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-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
33 template<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 
77 template<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 
127 template<class Type>
130 (
132 )
133 {
134  const fvMesh& mesh = dynamic_cast<const fvMesh&>(pfld.mesh()());
135 
137  (
138  IOobject
139  (
140  "cellAvg",
141  mesh.time().timeName(),
142  pfld.db(),
143  IOobject::NO_READ,
144  IOobject::NO_WRITE,
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 
180 template<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 
204 template<class Type, class GeoMeshType>
206 (
207  const objectRegistry& obr,
208  const word& fieldName,
209  const dimensionSet& dims,
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 template<class Type, class GeoMeshType>
230 (
231  const word& fieldName,
232  const dimensionSet& dims,
233  const Field<Type>& values,
234  word lookupName
235 ) const
236 {
237  surfMesh* surfptr = this->getSurfMesh(lookupName);
238 
239  if (surfptr)
240  {
241  surfptr->storeField<Type, GeoMeshType>
242  (
243  fieldName, dims, values
244  );
245  }
246 
247  return surfptr;
248 }
249 
250 
251 template<class Type, class GeoMeshType>
253 (
254  const word& fieldName,
255  const dimensionSet& dims,
257  word lookupName
258 ) const
259 {
260  surfMesh* surfptr = this->getSurfMesh(lookupName);
261 
262  if (surfptr)
263  {
264  surfptr->storeField<Type, GeoMeshType>
265  (
266  fieldName, dims, std::move(values)
267  );
268  }
269 
270  return surfptr;
271 }
272 
273 
274 // ************************************************************************* //
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:63
Foam::sampledSurface::sampleOnFaces
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.
Foam::surfMesh::storeField
void storeField(const word &fieldName, const dimensionSet &dims, const Field< Type > &values)
Copy/store named field as face or point data (template parameter).
Definition: surfMeshTemplates.C:35
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::HashTableOps::values
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:149
Foam::surfMesh
A surface mesh consisting of general polygon faces that has IO capabilities and a registry for storin...
Definition: surfMesh.H:63
Foam::polySurface::storeField
void storeField(const word &fieldName, const dimensionSet &dims, const Field< Type > &values)
Copy/store named field as face or point data (template parameter).
Definition: polySurfaceTemplates.C:96
Foam::dimensionSet
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:108
Foam::sampledSurface::pointAverage
static tmp< GeometricField< Type, fvPatchField, volMesh > > pointAverage(const GeometricField< Type, pointPatchField, pointMesh > &pfld)
Create cell values by averaging the point values.
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
Foam::polySurface
A surface mesh consisting of general polygon faces and capable of holding fields.
Definition: polySurface.H:67
Foam::Field< vector >
sampledSurface.H
Foam::sampledSurface::storeRegistryField
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)
Definition: sampledSurfaceTemplates.C:182
Foam::interpolation
Abstract base class for interpolation.
Definition: mappedPatchFieldBase.H:96
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:42
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::sampledSurface::sampleOnPoints
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.
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::List< face >
Foam::UList< label >
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::UList::size
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:53
Foam::interpolation::interpolate
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.
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
Foam::sampledSurface::storeSurfMeshField
bool storeSurfMeshField(const word &fieldName, const dimensionSet &dims, const Field< Type > &values, word lookupName="") const
Copy/store sampled Face field onto surfMesh (if it exists)
Definition: sampledSurfaceTemplates.C:230