sampledPatch.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "sampledPatch.H"
27 #include "surfaceFields.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace sampledSurfaces
35 {
36  defineTypeNameAndDebug(patch, 0);
37  addToRunTimeSelectionTable(sampledSurface, patch, word);
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
45 (
46  const word& name,
47  const polyMesh& mesh,
48  const wordReList& patchNames,
49  const bool triangulate
50 )
51 :
52  sampledSurface(name, mesh),
53  patchNames_(patchNames),
54  triangulate_(triangulate),
55  needsUpdate_(true)
56 {}
57 
58 
60 (
61  const word& name,
62  const polyMesh& mesh,
63  const dictionary& dict
64 )
65 :
66  sampledSurface(name, mesh, dict),
67  patchNames_(dict.lookup("patches")),
68  triangulate_(dict.lookupOrDefault("triangulate", false)),
69  needsUpdate_(true)
70 {}
71 
72 
73 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
74 
76 {}
77 
78 
79 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
80 
82 {
83  if (patchIDs_.empty())
84  {
85  patchIDs_ = mesh().boundaryMesh().patchSet
86  (
87  patchNames_,
88  false
89  ).sortedToc();
90  }
91  return patchIDs_;
92 }
93 
94 
96 {
97  return needsUpdate_;
98 }
99 
100 
102 {
103  // already marked as expired
104  if (needsUpdate_)
105  {
106  return false;
107  }
108 
111  patchIDs_.clear();
112  patchIndex_.clear();
113  patchFaceLabels_.clear();
114  patchStart_.clear();
115 
116  needsUpdate_ = true;
117  return true;
118 }
119 
120 
122 {
123  if (!needsUpdate_)
124  {
125  return false;
126  }
127 
128  label sz = 0;
129  forAll(patchIDs(), i)
130  {
131  label patchi = patchIDs()[i];
132  const polyPatch& pp = mesh().boundaryMesh()[patchi];
133 
134  if (isA<emptyPolyPatch>(pp))
135  {
137  << "Cannot sample an empty patch. Patch " << pp.name()
138  << exit(FatalError);
139  }
140 
141  sz += pp.size();
142  }
143 
144  // For every face (or triangle) the originating patch and local face in the
145  // patch.
146  patchIndex_.setSize(sz);
147  patchFaceLabels_.setSize(sz);
148  patchStart_.setSize(patchIDs().size());
149  labelList meshFaceLabels(sz);
150 
151  sz = 0;
152 
153  forAll(patchIDs(), i)
154  {
155  label patchi = patchIDs()[i];
156 
157  patchStart_[i] = sz;
158 
159  const polyPatch& pp = mesh().boundaryMesh()[patchi];
160 
161  forAll(pp, j)
162  {
163  patchIndex_[sz] = i;
164  patchFaceLabels_[sz] = j;
165  meshFaceLabels[sz] = pp.start()+j;
166  sz++;
167  }
168  }
169 
170  indirectPrimitivePatch allPatches
171  (
172  IndirectList<face>(mesh().faces(), meshFaceLabels),
173  mesh().points()
174  );
175 
176  this->storedPoints() = allPatches.localPoints();
177  this->storedFaces() = allPatches.localFaces();
178 
179 
180  // triangulate uses remapFaces()
181  // - this is somewhat less efficient since it recopies the faces
182  // that we just created, but we probably don't want to do this
183  // too often anyhow.
184  if (triangulate_)
185  {
186  MeshStorage::triangulate();
187  }
188 
189  if (debug)
190  {
191  print(Pout);
192  Pout<< endl;
193  }
194 
195  needsUpdate_ = false;
196  return true;
197 }
198 
199 
200 // remap action on triangulation
201 void Foam::sampledSurfaces::patch::remapFaces(const labelUList& faceMap)
202 {
203  // recalculate the cells cut
204  if (notNull(faceMap) && faceMap.size())
205  {
206  MeshStorage::remapFaces(faceMap);
207  patchFaceLabels_ = labelList
208  (
209  UIndirectList<label>(patchFaceLabels_, faceMap)
210  );
211  patchIndex_ = labelList
212  (
213  UIndirectList<label>(patchIndex_, faceMap)
214  );
215 
216  // Redo patchStart.
217  if (patchIndex_.size() > 0)
218  {
219  patchStart_[patchIndex_[0]] = 0;
220  for (label i = 1; i < patchIndex_.size(); i++)
221  {
222  if (patchIndex_[i] != patchIndex_[i-1])
223  {
224  patchStart_[patchIndex_[i]] = i;
225  }
226  }
227  }
228  }
229 }
230 
231 
233 (
234  const volScalarField& vField
235 ) const
236 {
237  return sampleField(vField);
238 }
239 
240 
242 (
243  const volVectorField& vField
244 ) const
245 {
246  return sampleField(vField);
247 }
248 
249 
251 (
252  const volSphericalTensorField& vField
253 ) const
254 {
255  return sampleField(vField);
256 }
257 
258 
260 (
261  const volSymmTensorField& vField
262 ) const
263 {
264  return sampleField(vField);
265 }
266 
267 
269 (
270  const volTensorField& vField
271 ) const
272 {
273  return sampleField(vField);
274 }
275 
276 
278 (
279  const surfaceScalarField& sField
280 ) const
281 {
282  return sampleField(sField);
283 }
284 
285 
287 (
288  const surfaceVectorField& sField
289 ) const
290 {
291  return sampleField(sField);
292 }
293 
294 
296 (
297  const surfaceSphericalTensorField& sField
298 ) const
299 {
300  return sampleField(sField);
301 }
302 
303 
305 (
306  const surfaceSymmTensorField& sField
307 ) const
308 {
309  return sampleField(sField);
310 }
311 
312 
314 (
315  const surfaceTensorField& sField
316 ) const
317 {
318  return sampleField(sField);
319 }
320 
321 
323 (
324  const interpolation<scalar>& interpolator
325 ) const
326 {
327  return interpolateField(interpolator);
328 }
329 
330 
332 (
333  const interpolation<vector>& interpolator
334 ) const
335 {
336  return interpolateField(interpolator);
337 }
338 
339 
341 (
342  const interpolation<sphericalTensor>& interpolator
343 ) const
344 {
345  return interpolateField(interpolator);
346 }
347 
348 
350 (
351  const interpolation<symmTensor>& interpolator
352 ) const
353 {
354  return interpolateField(interpolator);
355 }
356 
357 
359 (
360  const interpolation<tensor>& interpolator
361 ) const
362 {
363  return interpolateField(interpolator);
364 }
365 
366 
368 {
369  os << "patch: " << name() << " :"
370  << " patches:" << patchNames()
371  << " faces:" << faces().size()
372  << " points:" << points().size();
373 }
374 
375 
376 // ************************************************************************* //
Foam::surfaceFields.
virtual void print(Ostream &) const
Write.
Definition: sampledPatch.C:367
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual bool expire()
Mark the surface as needing an update.
Definition: sampledPatch.C:101
tUEqn clear()
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const word & name() const
Return name.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const word & name() const
Return the name of this functionObject.
const labelList & patchIDs() const
Definition: sampledPatch.C:81
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
An abstract class for surfaces with sampling.
defineTypeNameAndDebug(distanceSurface, 0)
virtual ~patch()
Destructor.
Definition: sampledPatch.C:75
bool interpolate() const
Interpolation requested for surface.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual bool needsUpdate() const
Does the surface need an update?
Definition: sampledPatch.C:95
const Field< PointType > & localPoints() const
Return pointField of points in patch.
patch(const word &name, const polyMesh &mesh, const wordReList &patchNames, const bool triangulate=false)
Construct from components.
Definition: sampledPatch.C:45
Macros for easy insertion into run-time selection tables.
A list of faces which address into the list of points.
dynamicFvMesh & mesh
virtual void clearGeom() const
const pointField & points
A class for handling words, derived from string.
Definition: word.H:59
wordList patchNames(nPatches)
List< label > labelList
A List of labels.
Definition: labelList.H:56
virtual tmp< scalarField > sample(const volScalarField &) const
Sample field on surface.
Definition: sampledPatch.C:233
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
virtual bool update()
Update the surface as required.
Definition: sampledPatch.C:121
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
addToRunTimeSelectionTable(sampledSurface, distanceSurface, word)
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:52
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
label patchi
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
const List< FaceType > & localFaces() const
Return patch faces addressing into local point list.
A List with indirect addressing.
Definition: fvMatrix.H:106
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:303
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
A class for managing temporary objects.
Definition: PtrList.H:53
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
label size() const
Return the number of elements in the UList.
Definition: UListI.H:299
A List with indirect addressing.
Definition: IndirectList.H:101
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:812