sampledSurfaceTemplates.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) 2018-2022 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
29#include "sampledSurface.H"
30
31// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
32
33template<class Type>
36(
37 const interpolation<Type>& sampler,
38 const labelUList& elements,
39 const faceList& fcs,
40 const pointField& pts,
41 const Type& defaultValue
42)
43{
44 const label len = elements.size();
45
46 if (len != fcs.size())
47 {
49 << "size mismatch: "
50 << "sampled elements (" << len
51 << ") != faces (" << fcs.size() << ')'
52 << exit(FatalError);
53 }
54
55 auto tvalues = tmp<Field<Type>>::New(len);
56 auto& values = tvalues.ref();
57
58 for (label i=0; i < len; ++i)
59 {
60 const label celli = elements[i];
61 if (celli < 0)
62 {
63 values[i] = defaultValue;
64 }
65 else
66 {
67 const point pt = fcs[i].centre(pts);
68
69 values[i] = sampler.interpolate(pt, celli);
70 }
71 }
72
73 return tvalues;
74}
75
76
77template<class Type>
80(
81 const interpolation<Type>& interpolator,
82 const labelUList& elements,
83 const faceList& fcs,
84 const pointField& pts
85)
86{
87 const label len = elements.size();
88
89 if (len != fcs.size())
90 {
92 << "size mismatch: "
93 << "sampled elements (" << len
94 << ") != faces (" << fcs.size() << ')'
95 << exit(FatalError);
96 }
97
98 // One value per point
99 // Initialize with Zero to handle missed/degenerate faces
100 auto tvalues = tmp<Field<Type>>::New(pts.size(), Zero);
101 auto& values = tvalues.ref();
102
103 bitSet pointDone(pts.size());
104
105 forAll(fcs, facei)
106 {
107 const face& f = fcs[facei];
108 const label celli = elements[facei];
109
110 for (const label pointi : f)
111 {
112 if (pointDone.set(pointi))
113 {
114 values[pointi] = interpolator.interpolate
115 (
116 pts[pointi],
117 celli
118 );
119 }
120 }
121 }
122
123 return tvalues;
124}
125
126
127template<class Type>
130(
131 const PointField<Type>& pfld
132)
133{
134 const fvMesh& mesh = dynamic_cast<const fvMesh&>(pfld.mesh()());
135
136 auto tcellAvg = tmp<VolumeField<Type>>::New
137 (
139 (
140 "cellAvg",
141 mesh.time().timeName(),
142 pfld.db(),
145 false
146 ),
147 mesh,
149 );
150 auto& cellAvg = tcellAvg.ref();
151
152 labelField nPointCells(mesh.nCells(), Zero);
153
154 for (label pointi = 0; pointi < mesh.nPoints(); ++pointi)
155 {
156 const Type& val = pfld[pointi];
157 const labelList& pCells = mesh.pointCells(pointi);
158
159 for (const label celli : pCells)
160 {
161 cellAvg[celli] += val;
162 ++nPointCells[celli];
163 }
164 }
165
166 forAll(cellAvg, celli)
167 {
168 cellAvg[celli] /= nPointCells[celli];
169 }
170
171 // Give value to calculatedFvPatchFields
172 cellAvg.correctBoundaryConditions();
173
174 return tcellAvg;
175}
176
177
178// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
179
180template<class Type, class GeoMeshType>
182(
183 const objectRegistry& obr,
184 const word& fieldName,
185 const dimensionSet& dims,
186 const Field<Type>& values,
187 word lookupName
188) const
189{
190 polySurface* surfptr = this->getRegistrySurface(obr, lookupName);
191
192 if (surfptr)
193 {
194 surfptr->storeField<Type, GeoMeshType>
195 (
196 fieldName, dims, values
197 );
198 }
199
200 return surfptr;
201}
202
203
204template<class Type, class GeoMeshType>
206(
207 const objectRegistry& obr,
208 const word& fieldName,
209 const dimensionSet& dims,
210 Field<Type>&& values,
211 word lookupName
212) const
213{
214 polySurface* surfptr = this->getRegistrySurface(obr, lookupName);
215
216 if (surfptr)
217 {
218 surfptr->storeField<Type, GeoMeshType>
219 (
220 fieldName, dims, std::move(values)
221 );
222 }
223
224 return surfptr;
225}
226
227
228// ************************************************************************* //
const Mesh & mesh() const
Return mesh.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:500
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:66
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition: bitSetI.H:590
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:109
Generic dimensioned Type class.
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
Abstract base class for volume field interpolation.
Definition: interpolation.H:60
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.
Registry of regIOobjects.
A surface mesh consisting of general polygon faces and capable of holding fields.
Definition: polySurface.H:71
void storeField(const word &fieldName, const dimensionSet &dims, const Field< Type > &values)
Copy/store named field as face or point data (template parameter).
const labelListList & pointCells() const
label nPoints() const noexcept
Number of mesh points.
label nCells() const noexcept
Number of mesh cells.
static tmp< Field< Type > > sampleOnPoints(const interpolation< Type > &interpolator, const labelUList &elements, const faceList &fcs, const pointField &pts)
Loop for interpolating volume elements to face points.
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 tmp< VolumeField< Type > > pointAverage(const PointField< Type > &pfld)
Create cell values by averaging the point values.
bool storeRegistryField(const objectRegistry &obr, const word &fieldName, const dimensionSet &dims, const Field< Type > &values, word lookupName="") const
Copy/store sampled field onto registered surface (if it exists)
static autoPtr< sampledSurface > New(const word &name, const polyMesh &mesh, const dictionary &dict)
Return a reference to the selected surface.
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const dimensionSet dimless
Dimensionless.
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.
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
labelList f(nPoints)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333