sampledSurfaces.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-2022 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 "sampledSurfaces.H"
27 #include "PatchTools.H"
28 #include "polyTopoChangeMap.H"
29 #include "polyMeshMap.H"
30 #include "polyDistributionMap.H"
31 #include "OSspecific.H"
32 #include "writeFile.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace functionObjects
40 {
42 
44  (
48  );
49 }
50 }
51 
52 Foam::scalar Foam::functionObjects::sampledSurfaces::mergeTol_ = 1e-10;
53 
54 
55 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
56 
57 bool Foam::functionObjects::sampledSurfaces::needsUpdate() const
58 {
59  forAll(*this, si)
60  {
61  if (operator[](si).needsUpdate())
62  {
63  return true;
64  }
65  }
66 
67  return false;
68 }
69 
70 
71 bool Foam::functionObjects::sampledSurfaces::expire()
72 {
73  bool justExpired = false;
74 
75  forAll(*this, si)
76  {
77  if (operator[](si).expire())
78  {
79  justExpired = true;
80  }
81 
82  // Clear merge information
83  if (Pstream::parRun())
84  {
85  mergeList_[si].clear();
86  }
87  }
88 
89  // true if any surfaces just expired
90  return justExpired;
91 }
92 
93 
94 bool Foam::functionObjects::sampledSurfaces::update()
95 {
96  bool updated = false;
97 
98  if (!needsUpdate())
99  {
100  return updated;
101  }
102 
103  // Serial: quick and easy, no merging required
104  if (!Pstream::parRun())
105  {
106  forAll(*this, si)
107  {
108  if (operator[](si).update())
109  {
110  updated = true;
111  }
112  }
113 
114  return updated;
115  }
116 
117  // Dimension as fraction of mesh bounding box
118  scalar mergeDim = mergeTol_ * mesh_.bounds().mag();
119 
120  if (Pstream::master() && debug)
121  {
122  Pout<< nl << "Merging all points within "
123  << mergeDim << " metre" << endl;
124  }
125 
126  forAll(*this, si)
127  {
128  sampledSurface& s = operator[](si);
129 
130  if (s.update())
131  {
132  updated = true;
133  }
134  else
135  {
136  continue;
137  }
138 
140  (
141  mergeDim,
143  (
144  SubList<face>(s.faces(), s.faces().size()),
145  s.points()
146  ),
147  mergeList_[si].points,
148  mergeList_[si].faces,
149  mergeList_[si].pointsMap
150  );
151  }
152 
153  return updated;
154 }
155 
156 
157 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
158 
160 (
161  const word& name,
162  const Time& t,
163  const dictionary& dict
164 )
165 :
168  outputPath_
169  (
170  mesh_.time().globalPath()
171  /writeFile::outputPrefix
172  /(mesh_.name() != polyMesh::defaultRegion ? mesh_.name() : word())
173  /name
174  ),
175  fields_(),
176  interpolationScheme_(word::null),
177  writeEmpty_(false),
178  mergeList_(),
179  formatter_(nullptr)
180 {
181  read(dict);
182 }
183 
184 
185 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
186 
188 {}
189 
190 
191 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
192 
194 {
195  bool surfacesFound = dict.found("surfaces");
196 
197  if (surfacesFound)
198  {
199  dict.lookup("fields") >> fields_;
200 
201  dict.lookup("interpolationScheme") >> interpolationScheme_;
202 
203  dict.readIfPresent("writeEmpty", writeEmpty_);
204 
205  const word writeType(dict.lookup("surfaceFormat"));
206 
207  // Define the surface formatter
208  formatter_ = surfaceWriter::New(writeType, dict);
209 
211  (
212  dict.lookup("surfaces"),
213  sampledSurface::iNew(mesh_)
214  );
215  transfer(newList);
216 
217  if (Pstream::parRun())
218  {
219  mergeList_.setSize(size());
220  }
221 
222  // Ensure all surfaces and merge information are expired
223  expire();
224 
225  if (this->size())
226  {
227  Info<< "Reading surface description:" << nl;
228  forAll(*this, si)
229  {
230  Info<< " " << operator[](si).name() << nl;
231  }
232  Info<< endl;
233  }
234  }
235 
236  if (Pstream::master() && debug)
237  {
238  Pout<< "sample fields:" << fields_ << nl
239  << "sample surfaces:" << nl << "(" << nl;
240 
241  forAll(*this, si)
242  {
243  Pout<< " " << operator[](si) << endl;
244  }
245  Pout<< ")" << endl;
246  }
247 
248  return true;
249 }
250 
251 
253 {
254  wordList fields(fields_);
255 
256  forAll(*this, si)
257  {
258  fields.append(operator[](si).fields());
259  }
260 
261  return fields;
262 }
263 
264 
266 {
267  return true;
268 }
269 
270 
272 {
273  if (size())
274  {
275  // Finalise surfaces, merge points etc.
276  update();
277 
278  // Create the output directory
279  if (Pstream::master())
280  {
281  if (debug)
282  {
283  Pout<< "Creating directory "
284  << outputPath_/mesh_.time().name() << nl << endl;
285 
286  }
287 
288  mkDir(outputPath_/mesh_.time().name());
289  }
290 
291  // Create a list of names of fields that are actually available
293  forAll(fields_, fieldi)
294  {
295  #define FoundFieldType(Type, nullArg) \
296  || foundObject<VolField<Type>>(fields_[fieldi]) \
297  || foundObject<SurfaceField<Type>>(fields_[fieldi])
299  {
300  fieldNames.append(fields_[fieldi]);
301  }
302  else
303  {
304  cannotFindObject(fields_[fieldi]);
305  }
306  #undef FoundFieldType
307  }
308 
309  // Create table of cached interpolations, to prevent unnecessary work
310  // when interpolating fields over multiple surfaces
311  #define DeclareInterpolations(Type, nullArg) \
312  HashPtrTable<interpolation<Type>> interpolation##Type##s;
314  #undef DeclareInterpolations
315 
316  // Sample and write the surfaces
317  forAll(*this, surfi)
318  {
319  const sampledSurface& s = operator[](surfi);
320 
321  #define GenerateFieldTypeValues(Type, nullArg) \
322  PtrList<Field<Type>> field##Type##Values = \
323  sampleType<Type>(surfi, fieldNames, interpolation##Type##s);
325  #undef GenerateFieldTypeValues
326 
327  if (Pstream::parRun())
328  {
329  if
330  (
332  && (mergeList_[surfi].faces.size() || writeEmpty_)
333  )
334  {
335  formatter_->write
336  (
337  outputPath_/mesh_.time().name(),
338  s.name(),
339  mergeList_[surfi].points,
340  mergeList_[surfi].faces,
341  fieldNames,
342  s.interpolate()
343  #define FieldTypeValuesParameter(Type, nullArg) \
344  , field##Type##Values
347  );
348  }
349  }
350  else
351  {
352  if (s.faces().size() || writeEmpty_)
353  {
354  formatter_->write
355  (
356  outputPath_/mesh_.time().name(),
357  s.name(),
358  s.points(),
359  s.faces(),
360  fieldNames,
361  s.interpolate()
362  #define FieldTypeValuesParameter(Type, nullArg) \
363  , field##Type##Values
366  );
367  }
368  }
369  }
370  }
371 
372  return true;
373 }
374 
375 
377 {
378  if (&mesh == &mesh_)
379  {
380  expire();
381  }
382 }
383 
384 
386 (
387  const polyTopoChangeMap& map
388 )
389 {
390  if (&map.mesh() == &mesh_)
391  {
392  expire();
393  }
394 }
395 
396 
398 (
399  const polyMeshMap& map
400 )
401 {
402  if (&map.mesh() == &mesh_)
403  {
404  expire();
405  }
406 }
407 
408 
410 (
411  const polyDistributionMap& map
412 )
413 {
414  if (&map.mesh() == &mesh_)
415  {
416  expire();
417  }
418 }
419 
420 
421 // ************************************************************************* //
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:85
static void gatherAndMerge(const scalar mergeDist, const PrimitivePatch< FaceList, PointField > &p, Field< typename PrimitivePatch< FaceList, PointField >::PointType > &mergedPoints, List< typename PrimitivePatch< FaceList, PointField >::FaceType > &mergedFaces, labelList &pointMergeMap)
Gather points and faces onto master and merge into single patch.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Abstract base-class for Time/database functionObjects.
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
virtual wordList fields() const
Return the list of fields required.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map - expires the surfaces.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
virtual void movePoints(const polyMesh &)
Update for mesh point-motion - expires the surfaces.
sampledSurfaces(const word &name, const Time &time, const dictionary &dict)
Construct from Time and dictionary.
virtual bool execute()
Execute, currently does nothing.
virtual bool write()
Sample and write.
virtual bool read(const dictionary &)
Read the sampledSurfaces dictionary.
Writes run time, CPU time and clock time and optionally the CPU and clock times per time step.
functionObject base class for writing single files
Definition: writeFile.H:56
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
const polyMesh & mesh() const
Return polyMesh.
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
const polyMesh & mesh() const
Return polyMesh.
Definition: polyMeshMap.H:75
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const polyMesh & mesh() const
Return polyMesh.
Class used for the PtrLists read-construction.
An abstract class for surfaces with sampling.
static autoPtr< surfaceWriter > New(const word &writeType, const IOstream::streamFormat writeFormat, const IOstream::compressionType writeCompression)
Select given write options.
Definition: surfaceWriter.C:75
A class for handling words, derived from string.
Definition: word.H:62
static List< word > fieldNames
Definition: globalFoam.H:46
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:230
const dimensionedScalar e
Elementary charge.
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
Namespace for OpenFOAM.
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:290
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
FOR_ALL_FIELD_TYPES(DefineContiguousFvWallLocationDataType)
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
Addressing for a faceList slice.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
static const char nl
Definition: Ostream.H:260
#define FoundFieldType(Type, nullArg)
#define FieldTypeValuesParameter(Type, nullArg)
#define DeclareInterpolations(Type, nullArg)
#define GenerateFieldTypeValues(Type, nullArg)
dictionary dict