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 "OSspecific.H"
30 #include "writeFile.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace functionObjects
38 {
39  defineTypeNameAndDebug(sampledSurfaces, 0);
40 
42  (
43  functionObject,
44  sampledSurfaces,
45  dictionary
46  );
47 }
48 }
49 
50 bool Foam::functionObjects::sampledSurfaces::verbose_ = false;
51 Foam::scalar Foam::functionObjects::sampledSurfaces::mergeTol_ = 1e-10;
52 
53 
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55 
57 (
58  const word& name,
59  const Time& t,
60  const dictionary& dict
61 )
62 :
63  fvMeshFunctionObject(name, t, dict),
65  outputPath_(fileName::null),
66  fields_(),
67  interpolationScheme_(word::null),
68  writeEmpty_(false),
69  mergeList_(),
70  formatter_(nullptr)
71 {
72  outputPath_ =
73  mesh_.time().globalPath()/functionObjects::writeFile::outputPrefix/name;
74 
75  if (mesh_.name() != fvMesh::defaultRegion)
76  {
77  outputPath_ = outputPath_/mesh_.name();
78  }
79 
80  read(dict);
81 }
82 
83 
84 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
85 
87 {}
88 
89 
90 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
91 
93 {
94  verbose_ = verbosity;
95 }
96 
97 
99 {
100  return true;
101 }
102 
103 
105 {
106  if (size())
107  {
108  // Finalise surfaces, merge points etc.
109  update();
110 
111  // Create the output directory
112  if (Pstream::master())
113  {
114  if (debug)
115  {
116  Pout<< "Creating directory "
117  << outputPath_/mesh_.time().timeName() << nl << endl;
118 
119  }
120 
121  mkDir(outputPath_/mesh_.time().timeName());
122  }
123 
124  // Create a list of names of fields that are actually available
126  forAll(fields_, fieldi)
127  {
128  #define FoundFieldType(Type, nullArg) \
129  || foundObject<VolField<Type>>(fields_[fieldi]) \
130  || foundObject<SurfaceField<Type>>(fields_[fieldi])
132  {
133  fieldNames.append(fields_[fieldi]);
134  }
135  else
136  {
137  cannotFindObject(fields_[fieldi]);
138  }
139  #undef FoundFieldType
140  }
141 
142  // Create table of cached interpolations, to prevent unnecessary work
143  // when interpolating fields over multiple surfaces
144  #define DeclareInterpolations(Type, nullArg) \
145  HashPtrTable<interpolation<Type>> interpolation##Type##s;
147  #undef DeclareInterpolations
148 
149  // Sample and write the surfaces
150  forAll(*this, surfi)
151  {
152  const sampledSurface& s = operator[](surfi);
153 
154  #define GenerateFieldTypeValues(Type, nullArg) \
155  PtrList<Field<Type>> field##Type##Values = \
156  sampleType<Type>(surfi, fieldNames, interpolation##Type##s);
158  #undef GenerateFieldTypeValues
159 
160  if (Pstream::parRun())
161  {
162  if
163  (
165  && (mergeList_[surfi].faces.size() || writeEmpty_)
166  )
167  {
168  formatter_->write
169  (
170  outputPath_/mesh_.time().timeName(),
171  s.name(),
172  mergeList_[surfi].points,
173  mergeList_[surfi].faces,
174  fieldNames,
175  s.interpolate()
176  #define FieldTypeValuesParameter(Type, nullArg) \
177  , field##Type##Values
179  #undef FieldTypeValuesParameter
180  );
181  }
182  }
183  else
184  {
185  if (s.faces().size() || writeEmpty_)
186  {
187  formatter_->write
188  (
189  outputPath_/mesh_.time().timeName(),
190  s.name(),
191  s.points(),
192  s.faces(),
193  fieldNames,
194  s.interpolate()
195  #define FieldTypeValuesParameter(Type, nullArg) \
196  , field##Type##Values
198  #undef FieldTypeValuesParameter
199  );
200  }
201  }
202  }
203  }
204 
205  return true;
206 }
207 
208 
210 {
211  bool surfacesFound = dict.found("surfaces");
212 
213  if (surfacesFound)
214  {
215  dict.lookup("fields") >> fields_;
216 
217  dict.lookup("interpolationScheme") >> interpolationScheme_;
218 
219  dict.readIfPresent("writeEmpty", writeEmpty_);
220 
221  const word writeType(dict.lookup("surfaceFormat"));
222 
223  // Define the surface formatter
224  formatter_ = surfaceWriter::New(writeType, dict);
225 
227  (
228  dict.lookup("surfaces"),
229  sampledSurface::iNew(mesh_)
230  );
231  transfer(newList);
232 
233  if (Pstream::parRun())
234  {
235  mergeList_.setSize(size());
236  }
237 
238  // Ensure all surfaces and merge information are expired
239  expire();
240 
241  if (this->size())
242  {
243  Info<< "Reading surface description:" << nl;
244  forAll(*this, si)
245  {
246  Info<< " " << operator[](si).name() << nl;
247  }
248  Info<< endl;
249  }
250  }
251 
252  if (Pstream::master() && debug)
253  {
254  Pout<< "sample fields:" << fields_ << nl
255  << "sample surfaces:" << nl << "(" << nl;
256 
257  forAll(*this, si)
258  {
259  Pout<< " " << operator[](si) << endl;
260  }
261  Pout<< ")" << endl;
262  }
263 
264  return true;
265 }
266 
267 
269 {
270  wordList fields(fields_);
271 
272  forAll(*this, si)
273  {
274  fields.append(operator[](si).fields());
275  }
276 
277  return fields;
278 }
279 
280 
282 {
283  if (&mesh == &mesh_)
284  {
285  expire();
286  }
287 }
288 
289 
291 (
292  const polyTopoChangeMap& map
293 )
294 {
295  if (&map.mesh() == &mesh_)
296  {
297  expire();
298  }
299 
300  // pointMesh and interpolation will have been reset in mesh.update
301 }
302 
303 
305 (
306  const polyMeshMap& map
307 )
308 {
309  expire();
310 }
311 
312 
314 (
315  const polyMesh::readUpdateState state
316 )
317 {
318  if (state != polyMesh::UNCHANGED)
319  {
320  expire();
321  }
322 }
323 
324 
326 {
327  forAll(*this, si)
328  {
329  if (operator[](si).needsUpdate())
330  {
331  return true;
332  }
333  }
334 
335  return false;
336 }
337 
338 
340 {
341  bool justExpired = false;
342 
343  forAll(*this, si)
344  {
345  if (operator[](si).expire())
346  {
347  justExpired = true;
348  }
349 
350  // Clear merge information
351  if (Pstream::parRun())
352  {
353  mergeList_[si].clear();
354  }
355  }
356 
357  // true if any surfaces just expired
358  return justExpired;
359 }
360 
361 
363 {
364  bool updated = false;
365 
366  if (!needsUpdate())
367  {
368  return updated;
369  }
370 
371  // Serial: quick and easy, no merging required
372  if (!Pstream::parRun())
373  {
374  forAll(*this, si)
375  {
376  if (operator[](si).update())
377  {
378  updated = true;
379  }
380  }
381 
382  return updated;
383  }
384 
385  // Dimension as fraction of mesh bounding box
386  scalar mergeDim = mergeTol_ * mesh_.bounds().mag();
387 
388  if (Pstream::master() && debug)
389  {
390  Pout<< nl << "Merging all points within "
391  << mergeDim << " metre" << endl;
392  }
393 
394  forAll(*this, si)
395  {
396  sampledSurface& s = operator[](si);
397 
398  if (s.update())
399  {
400  updated = true;
401  }
402  else
403  {
404  continue;
405  }
406 
408  (
409  mergeDim,
411  (
412  SubList<face>(s.faces(), s.faces().size()),
413  s.points()
414  ),
415  mergeList_[si].points,
416  mergeList_[si].faces,
417  mergeList_[si].pointsMap
418  );
419  }
420 
421  return updated;
422 }
423 
424 
425 // ************************************************************************* //
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:663
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual bool execute()
Execute, currently does nothing.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
An abstract class for surfaces with sampling.
const word & name() const
Name of surface.
#define DeclareInterpolations(Type, nullArg)
sampledSurfaces(const word &name, const Time &time, const dictionary &dict)
Construct from Time and dictionary.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
addToRunTimeSelectionTable(functionObject, Qdot, dictionary)
bool interpolate() const
Interpolation requested for surface.
static const fileName null
An empty fileName.
Definition: fileName.H:97
virtual bool needsUpdate() const
Does any of the surfaces need an update?
virtual wordList fields() const
Return the list of fields required.
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:325
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual void movePoints(const polyMesh &)
Update for mesh point-motion - expires the surfaces.
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
fvMesh & mesh
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Macros for easy insertion into run-time selection tables.
static List< word > fieldNames
Definition: globalFoam.H:46
virtual const faceList & faces() const =0
Faces of surface.
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
A list of faces which address into the list of points.
A List obtained as a section of another List.
Definition: SubList.H:53
bool read(const char *, int32_t &)
Definition: int32IO.C:85
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.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
virtual bool write()
Sample and write.
const pointField & points
static const word outputPrefix
Directory prefix.
Definition: writeFile.H:72
A class for handling words, derived from string.
Definition: word.H:59
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
#define FieldTypeValuesParameter(Type, nullArg)
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
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.
#define GenerateFieldTypeValues(Type, nullArg)
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
static const word null
An empty word.
Definition: word.H:77
static autoPtr< surfaceWriter > New(const word &writeType, const IOstream::streamFormat writeFormat, const IOstream::compressionType writeCompression)
Select given write options.
Definition: surfaceWriter.C:75
const polyMesh & mesh() const
Return polyMesh.
virtual bool update()=0
Update the surface as required.
static const char nl
Definition: Ostream.H:260
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map - expires the surfaces.
virtual bool expire()
Mark the surfaces as needing an update.
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:290
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
defineTypeNameAndDebug(Qdot, 0)
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
Info<< "Reading field p_rgh\"<< endl;volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);pressureReference pressureReference(p, p_rgh, pimple.dict(), thermo.incompressible());mesh.schemes().setFluxRequired(p_rgh.name());hydrostaticInitialisation(p_rgh, p, rho, U, gh, ghf, pRef, thermo, pimple.dict());Info<< "Creating field dpdt\"<< endl;volScalarField dpdt(IOobject("dpdt", runTime.timeName(), mesh), mesh, dimensionedScalar(p.dimensions()/dimTime, 0));Info<< "Creating field kinetic energy K\"<< endl;volScalarField K("K", 0.5 *magSqr(U));dimensionedScalar initialMass=fvc::domainIntegrate(rho);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:131
void verbose(const bool verbosity=true)
Set verbosity level.
virtual bool update()
Update the surfaces as required and merge surface points (parallel).
virtual const pointField & points() const =0
Points of surface.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
virtual bool read(const dictionary &)
Read the sampledSurfaces dictionary.
messageStream Info
FOR_ALL_FIELD_TYPES(DefineFvWallInfoType)
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:90
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:50
#define FoundFieldType(Type, nullArg)
virtual void readUpdate(const polyMesh::readUpdateState state)
Update topology using the given map due to readUpdate.
Class used for the PtrLists read-construction.
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864