sampledMeshedSurfaceTemplates.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2016-2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
30
31// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32
33template<class Type>
35Foam::sampledMeshedSurface::sampleOnFaces
36(
37 const interpolation<Type>& sampler
38) const
39{
40 const Type deflt
41 (
42 defaultValues_.getOrDefault<Type>(sampler.psi().name(), Zero)
43 );
44
45 const labelList& elements = sampleElements_;
46
47
48 if (!onBoundary())
49 {
50 // Sample cells
51
53 (
54 sampler,
55 elements,
56 faces(),
57 points(),
58 deflt
59 );
60 }
61
62
63 //
64 // Sample boundary faces
65 //
66
67 auto tvalues = tmp<Field<Type>>::New(elements.size());
68 auto& values = tvalues.ref();
69
70 // Create flat boundary field
71
72 const polyBoundaryMesh& pbm = mesh().boundaryMesh();
73
74 Field<Type> bVals(mesh().nBoundaryFaces(), Zero);
75
76 const auto& bField = sampler.psi().boundaryField();
77
78 forAll(bField, patchi)
79 {
80 SubList<Type>(bVals, pbm[patchi].range()) = bField[patchi];
81 }
82
83 // Sample within the flat boundary field
84
85 forAll(elements, i)
86 {
87 const label bFacei = (elements[i] - mesh().nInternalFaces());
88
89 if (bFacei < 0)
90 {
91 values[i] = deflt;
92 }
93 else
94 {
95 values[i] = bVals[bFacei];
96 }
97 }
98
99 return tvalues;
100}
101
102
103template<class Type>
105Foam::sampledMeshedSurface::sampleOnPoints
106(
107 const interpolation<Type>& interpolator
108) const
109{
110 const Type deflt
111 (
112 defaultValues_.getOrDefault<Type>(interpolator.psi().name(), Zero)
113 );
114
115 const labelList& elements = sampleElements_;
116
117
118 // One value per vertex.
119 // - sampleElements_ and samplePoints_ have the identical size
120 auto tvalues = tmp<Field<Type>>::New(elements.size());
121 auto& values = tvalues.ref();
122
123 if (!onBoundary())
124 {
125 // Sample cells
126
127 forAll(elements, i)
128 {
129 const label celli = elements[i];
130
131 if (celli < 0)
132 {
133 values[i] = deflt;
134 }
135 else
136 {
137 values[i] = interpolator.interpolate
138 (
139 samplePoints_[i],
140 celli
141 );
142 }
143 }
144
145 return tvalues;
146 }
147
148
149 //
150 // Sample boundary faces
151 //
152
153 forAll(elements, i)
154 {
155 const label facei = elements[i];
156
157 if (facei < 0)
158 {
159 values[i] = deflt;
160 }
161 else
162 {
163 values[i] = interpolator.interpolate
164 (
165 samplePoints_[i],
166 mesh().faceOwner()[facei],
167 facei
168 );
169 }
170 }
171
172 return tvalues;
173}
174
175
176// ************************************************************************* //
scalar range
Generic templated field type.
Definition: Field.H:82
A List obtained as a section of another List.
Definition: SubList.H:70
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Abstract base class for volume field interpolation.
Definition: interpolation.H:60
const GeometricField< Type, fvPatchField, volMesh > & psi() const noexcept
Return the field to be interpolated.
virtual Type interpolate(const vector &position, const label celli, const label facei=-1) const =0
Interpolate field to the given point in the given cell.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
label nInternalFaces() const noexcept
Number of internal faces.
virtual const pointField & points() const
Points of surface.
bool onBoundary() const
Sampling boundary values instead of cell values.
virtual const faceList & faces() const
Faces of surface.
static tmp< Field< Type > > sampleOnFaces(const interpolation< Type > &sampler, const labelUList &elements, const faceList &fcs, const pointField &pts, const Type &defaultValue=Type(Zero))
Loop for sampling volume elements to faces.
static autoPtr< sampledSurface > New(const word &name, const polyMesh &mesh, const dictionary &dict)
Return a reference to the selected surface.
const polyMesh & mesh() const noexcept
Access to the underlying mesh.
A class for managing temporary objects.
Definition: tmp.H:65
dynamicFvMesh & mesh
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333