sampledPatch.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 | Copyright (C) 2011-2016 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 "dictionary.H"
28 #include "polyMesh.H"
29 #include "polyPatch.H"
30 #include "volFields.H"
31 #include "surfaceFields.H"
32 
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(sampledPatch, 0);
40  addNamedToRunTimeSelectionTable(sampledSurface, sampledPatch, word, patch);
41 }
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
47 (
48  const word& name,
49  const polyMesh& mesh,
50  const wordReList& patchNames,
51  const bool triangulate
52 )
53 :
54  sampledSurface(name, mesh),
55  patchNames_(patchNames),
56  triangulate_(triangulate),
57  needsUpdate_(true)
58 {}
59 
60 
62 (
63  const word& name,
64  const polyMesh& mesh,
65  const dictionary& dict
66 )
67 :
68  sampledSurface(name, mesh, dict),
69  patchNames_(dict.lookup("patches")),
70  triangulate_(dict.lookupOrDefault("triangulate", false)),
71  needsUpdate_(true)
72 {}
73 
74 
75 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
76 
78 {}
79 
80 
81 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
82 
84 {
85  if (patchIDs_.empty())
86  {
87  patchIDs_ = mesh().boundaryMesh().patchSet
88  (
89  patchNames_,
90  false
91  ).sortedToc();
92  }
93  return patchIDs_;
94 }
95 
96 
98 {
99  return needsUpdate_;
100 }
101 
102 
104 {
105  // already marked as expired
106  if (needsUpdate_)
107  {
108  return false;
109  }
110 
113  patchIDs_.clear();
114  patchIndex_.clear();
115  patchFaceLabels_.clear();
116  patchStart_.clear();
117 
118  needsUpdate_ = true;
119  return true;
120 }
121 
122 
124 {
125  if (!needsUpdate_)
126  {
127  return false;
128  }
129 
130  label sz = 0;
131  forAll(patchIDs(), i)
132  {
133  label patchi = patchIDs()[i];
134  const polyPatch& pp = mesh().boundaryMesh()[patchi];
135 
136  if (isA<emptyPolyPatch>(pp))
137  {
139  << "Cannot sample an empty patch. Patch " << pp.name()
140  << exit(FatalError);
141  }
142 
143  sz += pp.size();
144  }
145 
146  // For every face (or triangle) the originating patch and local face in the
147  // patch.
148  patchIndex_.setSize(sz);
149  patchFaceLabels_.setSize(sz);
150  patchStart_.setSize(patchIDs().size());
151  labelList meshFaceLabels(sz);
152 
153  sz = 0;
154 
155  forAll(patchIDs(), i)
156  {
157  label patchi = patchIDs()[i];
158 
159  patchStart_[i] = sz;
160 
161  const polyPatch& pp = mesh().boundaryMesh()[patchi];
162 
163  forAll(pp, j)
164  {
165  patchIndex_[sz] = i;
166  patchFaceLabels_[sz] = j;
167  meshFaceLabels[sz] = pp.start()+j;
168  sz++;
169  }
170  }
171 
172  indirectPrimitivePatch allPatches
173  (
174  IndirectList<face>(mesh().faces(), meshFaceLabels),
175  mesh().points()
176  );
177 
178  this->storedPoints() = allPatches.localPoints();
179  this->storedFaces() = allPatches.localFaces();
180 
181 
182  // triangulate uses remapFaces()
183  // - this is somewhat less efficient since it recopies the faces
184  // that we just created, but we probably don't want to do this
185  // too often anyhow.
186  if (triangulate_)
187  {
188  MeshStorage::triangulate();
189  }
190 
191  if (debug)
192  {
193  print(Pout);
194  Pout<< endl;
195  }
196 
197  needsUpdate_ = false;
198  return true;
199 }
200 
201 
202 // remap action on triangulation
203 void Foam::sampledPatch::remapFaces(const labelUList& faceMap)
204 {
205  // recalculate the cells cut
206  if (notNull(faceMap) && faceMap.size())
207  {
208  MeshStorage::remapFaces(faceMap);
209  patchFaceLabels_ = labelList
210  (
211  UIndirectList<label>(patchFaceLabels_, faceMap)
212  );
213  patchIndex_ = labelList
214  (
215  UIndirectList<label>(patchIndex_, faceMap)
216  );
217 
218  // Redo patchStart.
219  if (patchIndex_.size() > 0)
220  {
221  patchStart_[patchIndex_[0]] = 0;
222  for (label i = 1; i < patchIndex_.size(); i++)
223  {
224  if (patchIndex_[i] != patchIndex_[i-1])
225  {
226  patchStart_[patchIndex_[i]] = i;
227  }
228  }
229  }
230  }
231 }
232 
233 
235 (
236  const volScalarField& vField
237 ) const
238 {
239  return sampleField(vField);
240 }
241 
242 
244 (
245  const volVectorField& vField
246 ) const
247 {
248  return sampleField(vField);
249 }
250 
251 
253 (
254  const volSphericalTensorField& vField
255 ) const
256 {
257  return sampleField(vField);
258 }
259 
260 
262 (
263  const volSymmTensorField& vField
264 ) const
265 {
266  return sampleField(vField);
267 }
268 
269 
271 (
272  const volTensorField& vField
273 ) const
274 {
275  return sampleField(vField);
276 }
277 
278 
280 (
281  const surfaceScalarField& sField
282 ) const
283 {
284  return sampleField(sField);
285 }
286 
287 
289 (
290  const surfaceVectorField& sField
291 ) const
292 {
293  return sampleField(sField);
294 }
295 
296 
298 (
299  const surfaceSphericalTensorField& sField
300 ) const
301 {
302  return sampleField(sField);
303 }
304 
305 
307 (
308  const surfaceSymmTensorField& sField
309 ) const
310 {
311  return sampleField(sField);
312 }
313 
314 
316 (
317  const surfaceTensorField& sField
318 ) const
319 {
320  return sampleField(sField);
321 }
322 
323 
325 (
326  const interpolation<scalar>& interpolator
327 ) const
328 {
329  return interpolateField(interpolator);
330 }
331 
332 
334 (
335  const interpolation<vector>& interpolator
336 ) const
337 {
338  return interpolateField(interpolator);
339 }
340 
341 
343 (
344  const interpolation<sphericalTensor>& interpolator
345 ) const
346 {
347  return interpolateField(interpolator);
348 }
349 
350 
352 (
353  const interpolation<symmTensor>& interpolator
354 ) const
355 {
356  return interpolateField(interpolator);
357 }
358 
359 
361 (
362  const interpolation<tensor>& interpolator
363 ) const
364 {
365  return interpolateField(interpolator);
366 }
367 
368 
370 {
371  os << "sampledPatch: " << name() << " :"
372  << " patches:" << patchNames()
373  << " faces:" << faces().size()
374  << " points:" << points().size();
375 }
376 
377 
378 // ************************************************************************* //
Foam::surfaceFields.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:424
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
virtual tmp< scalarField > sample(const volScalarField &) const
Sample field on surface.
Definition: sampledPatch.C:235
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
An abstract class for surfaces with sampling.
tUEqn clear()
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
bool interpolate() const
Interpolation requested for surface.
virtual ~sampledPatch()
Destructor.
Definition: sampledPatch.C:77
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool usePatchGroups=true) const
Return the set of patch IDs corresponding to the given names.
Macros for easy insertion into run-time selection tables.
const List< Face > & localFaces() const
Return patch faces addressing into local point list.
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
virtual bool update()
Update the surface as required.
Definition: sampledPatch.C:123
wordList patchNames(nPatches)
virtual bool expire()
Mark the surface as needing an update.
Definition: sampledPatch.C:103
List< label > labelList
A List of labels.
Definition: labelList.H:56
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:61
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
const Field< PointType > & localPoints() const
Return pointField of points in patch.
defineTypeNameAndDebug(combustionModel, 0)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
virtual void print(Ostream &) const
Write.
Definition: sampledPatch.C:369
addNamedToRunTimeSelectionTable(GAMGProcAgglomeration, noneGAMGProcAgglomeration, GAMGAgglomeration, none)
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:46
const labelList & patchIDs() const
Definition: sampledPatch.C:83
label patchi
virtual bool needsUpdate() const
Does the surface need an update?
Definition: sampledPatch.C:97
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
List< Key > sortedToc() const
Return the table of contents as a sorted list.
Definition: HashTable.C:217
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:300
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:102
sampledPatch(const word &name, const polyMesh &mesh, const wordReList &patchNames, const bool triangulate=false)
Construct from components.
Definition: sampledPatch.C:47
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576