Loading...
Searching...
No Matches
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-2024 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\*---------------------------------------------------------------------------*/
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>
43 (
44 sampler.psi().name(),
46 )
47 );
48
49 const labelList& elements = sampleElements_;
50
51
52 if (!onBoundary())
53 {
54 // Sample cells
55
57 (
58 sampler,
59 elements,
60 faces(),
61 points(),
62 deflt
63 );
64 }
65
66
67 //
68 // Sample boundary faces
69 //
70
71 auto tvalues = tmp<Field<Type>>::New(elements.size());
72 auto& values = tvalues.ref();
73
74 // Create flat boundary field
75
77
78 Field<Type> bVals(mesh().nBoundaryFaces(), deflt);
79
80 const auto& bField = sampler.psi().boundaryField();
81
82 forAll(bField, patchi)
83 {
84 // Note: restrict transcribing to actual size of the patch field
85 // - handles "empty" patch type etc.
86 const auto& pfld = bField[patchi];
87 SubList<Type>(bVals, pfld.size(), pbm[patchi].offset()) = pfld;
88 }
89
90 // Sample within the flat boundary field
91
92 forAll(elements, i)
93 {
94 const label bFacei = (elements[i] - mesh().nInternalFaces());
95
96 if (bFacei < 0)
97 {
98 values[i] = deflt;
99 }
100 else
101 {
102 values[i] = bVals[bFacei];
103 }
104 }
106 return tvalues;
107}
108
109
110template<class Type>
112Foam::sampledMeshedSurface::sampleOnPoints
113(
114 const interpolation<Type>& interpolator
115) const
116{
117 const Type deflt
118 (
119 defaultValues_.getOrDefault<Type>
120 (
121 interpolator.psi().name(),
122 Foam::zero{}
123 )
124 );
125
126 const labelList& elements = sampleElements_;
127
128
129 // One value per vertex.
130 // - sampleElements_ and samplePoints_ have the identical size
131 auto tvalues = tmp<Field<Type>>::New(elements.size());
132 auto& values = tvalues.ref();
133
134 if (!onBoundary())
135 {
136 // Sample cells
137
138 forAll(elements, i)
139 {
140 const label celli = elements[i];
141
142 if (celli < 0)
143 {
144 values[i] = deflt;
145 }
146 else
147 {
148 values[i] = interpolator.interpolate
149 (
150 samplePoints_[i],
151 celli
152 );
153 }
154 }
155
156 return tvalues;
157 }
158
159
160 //
161 // Sample boundary faces
162 //
163
164 forAll(elements, i)
165 {
166 const label facei = elements[i];
167
168 if (facei < 0)
169 {
170 values[i] = deflt;
171 }
172 else
173 {
174 values[i] = interpolator.interpolate
175 (
176 samplePoints_[i],
177 mesh().faceOwner()[facei],
178 facei
179 );
180 }
181 }
182
183 return tvalues;
184}
185
186
187// ************************************************************************* //
const polyBoundaryMesh & pbm
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
A non-owning sub-view of a List (allocated or unallocated storage).
Definition SubList.H:61
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Abstract base class for volume field interpolation.
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 registered IO, a reference to the associated polyMesh,...
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
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:75
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
dynamicFvMesh & mesh
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition HashOps.H:164
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
List< label > labelList
A List of labels.
Definition List.H:62
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299