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-2021 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 "mapPolyMesh.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 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
55 
56 void Foam::functionObjects::sampledSurfaces::writeGeometry() const
57 {
58  // Write to time directory under outputPath_
59  // Skip surface without faces (eg, a failed cut-plane)
60 
61  const fileName outputDir = outputPath_/mesh_.time().timeName();
62 
63  forAll(*this, surfI)
64  {
65  const sampledSurface& s = operator[](surfI);
66 
67  if (Pstream::parRun())
68  {
69  if (Pstream::master() && mergeList_[surfI].faces.size())
70  {
71  formatter_->write
72  (
73  outputDir,
74  s.name(),
75  mergeList_[surfI].points,
76  mergeList_[surfI].faces
77  );
78  }
79  }
80  else if (s.faces().size())
81  {
82  formatter_->write
83  (
84  outputDir,
85  s.name(),
86  s.points(),
87  s.faces()
88  );
89  }
90  }
91 }
92 
93 
94 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
95 
97 (
98  const word& name,
99  const Time& t,
100  const dictionary& dict
101 )
102 :
103  functionObject(name),
105  mesh_
106  (
107  refCast<const fvMesh>
108  (
110  (
112  )
113  )
114  ),
115  loadFromFiles_(false),
116  outputPath_(fileName::null),
117  fieldSelection_(),
118  interpolationScheme_(word::null),
119  mergeList_(),
120  formatter_(nullptr)
121 {
122  outputPath_ =
123  mesh_.time().globalPath()/functionObjects::writeFile::outputPrefix/name;
124 
125  if (mesh_.name() != fvMesh::defaultRegion)
126  {
127  outputPath_ = outputPath_/mesh_.name();
128  }
129 
130  read(dict);
131 }
132 
133 
135 (
136  const word& name,
137  const objectRegistry& obr,
138  const dictionary& dict,
139  const bool loadFromFiles
140 )
141 :
142  functionObject(name),
144  mesh_(refCast<const fvMesh>(obr)),
145  loadFromFiles_(loadFromFiles),
146  outputPath_(fileName::null),
147  fieldSelection_(),
148  interpolationScheme_(word::null),
149  mergeList_(),
150  formatter_(nullptr)
151 {
152  outputPath_ =
153  mesh_.time().globalPath()/functionObjects::writeFile::outputPrefix/name;
154 
155  if (mesh_.name() != fvMesh::defaultRegion)
156  {
157  outputPath_ = outputPath_/mesh_.name();
158  }
159 
160  read(dict);
161 }
162 
163 
164 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
165 
167 {}
168 
169 
170 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
171 
173 {
174  verbose_ = verbosity;
175 }
176 
177 
179 {
180  return true;
181 }
182 
183 
185 {
186  if (size())
187  {
188  // Finalise surfaces, merge points etc.
189  update();
190 
191  const label nFields = classifyFields();
192 
193  if (Pstream::master())
194  {
195  if (debug)
196  {
197  Pout<< "Creating directory "
198  << outputPath_/mesh_.time().timeName() << nl << endl;
199 
200  }
201 
202  mkDir(outputPath_/mesh_.time().timeName());
203  }
204 
205  // Write geometry first if required,
206  // or when no fields would otherwise be written
207  if (nFields == 0 || formatter_->separateGeometry())
208  {
209  writeGeometry();
210  }
211 
212  const IOobjectList objects(mesh_, mesh_.time().timeName());
213 
214  sampleAndWrite<volScalarField>(objects);
215  sampleAndWrite<volVectorField>(objects);
216  sampleAndWrite<volSphericalTensorField>(objects);
217  sampleAndWrite<volSymmTensorField>(objects);
218  sampleAndWrite<volTensorField>(objects);
219 
220  sampleAndWrite<surfaceScalarField>(objects);
221  sampleAndWrite<surfaceVectorField>(objects);
222  sampleAndWrite<surfaceSphericalTensorField>(objects);
223  sampleAndWrite<surfaceSymmTensorField>(objects);
224  sampleAndWrite<surfaceTensorField>(objects);
225  }
226 
227  return true;
228 }
229 
230 
232 {
233  bool surfacesFound = dict.found("surfaces");
234 
235  if (surfacesFound)
236  {
237  dict.lookup("fields") >> fieldSelection_;
238 
239  dict.lookup("interpolationScheme") >> interpolationScheme_;
240  const word writeType(dict.lookup("surfaceFormat"));
241 
242  // Define the surface formatter
243  formatter_ = surfaceWriter::New(writeType, dict);
244 
246  (
247  dict.lookup("surfaces"),
248  sampledSurface::iNew(mesh_)
249  );
250  transfer(newList);
251 
252  if (Pstream::parRun())
253  {
254  mergeList_.setSize(size());
255  }
256 
257  // Ensure all surfaces and merge information are expired
258  expire();
259 
260  if (this->size())
261  {
262  Info<< "Reading surface description:" << nl;
263  forAll(*this, surfI)
264  {
265  Info<< " " << operator[](surfI).name() << nl;
266  }
267  Info<< endl;
268  }
269  }
270 
271  if (Pstream::master() && debug)
272  {
273  Pout<< "sample fields:" << fieldSelection_ << nl
274  << "sample surfaces:" << nl << "(" << nl;
275 
276  forAll(*this, surfI)
277  {
278  Pout<< " " << operator[](surfI) << endl;
279  }
280  Pout<< ")" << endl;
281  }
282 
283  return true;
284 }
285 
286 
288 {
289  if (&mpm.mesh() == &mesh_)
290  {
291  expire();
292  }
293 
294  // pointMesh and interpolation will have been reset in mesh.update
295 }
296 
297 
299 {
300  if (&mesh == &mesh_)
301  {
302  expire();
303  }
304 }
305 
306 
308 (
309  const polyMesh::readUpdateState state
310 )
311 {
312  if (state != polyMesh::UNCHANGED)
313  {
314  expire();
315  }
316 }
317 
318 
320 {
321  forAll(*this, surfI)
322  {
323  if (operator[](surfI).needsUpdate())
324  {
325  return true;
326  }
327  }
328 
329  return false;
330 }
331 
332 
334 {
335  bool justExpired = false;
336 
337  forAll(*this, surfI)
338  {
339  if (operator[](surfI).expire())
340  {
341  justExpired = true;
342  }
343 
344  // Clear merge information
345  if (Pstream::parRun())
346  {
347  mergeList_[surfI].clear();
348  }
349  }
350 
351  // true if any surfaces just expired
352  return justExpired;
353 }
354 
355 
357 {
358  bool updated = false;
359 
360  if (!needsUpdate())
361  {
362  return updated;
363  }
364 
365  // Serial: quick and easy, no merging required
366  if (!Pstream::parRun())
367  {
368  forAll(*this, surfI)
369  {
370  if (operator[](surfI).update())
371  {
372  updated = true;
373  }
374  }
375 
376  return updated;
377  }
378 
379  // Dimension as fraction of mesh bounding box
380  scalar mergeDim = mergeTol_ * mesh_.bounds().mag();
381 
382  if (Pstream::master() && debug)
383  {
384  Pout<< nl << "Merging all points within "
385  << mergeDim << " metre" << endl;
386  }
387 
388  forAll(*this, surfI)
389  {
390  sampledSurface& s = operator[](surfI);
391 
392  if (s.update())
393  {
394  updated = true;
395  }
396  else
397  {
398  continue;
399  }
400 
402  (
403  mergeDim,
405  (
406  SubList<face>(s.faces(), s.faces().size()),
407  s.points()
408  ),
409  mergeList_[surfI].points,
410  mergeList_[surfI].faces,
411  mergeList_[surfI].pointsMap
412  );
413  }
414 
415  return updated;
416 }
417 
418 
419 // ************************************************************************* //
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:643
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
virtual bool execute()
Execute, currently does nothing.
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
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.
static autoPtr< surfaceWriter > New(const word &writeType, const IOstream::streamFormat writeFormat)
Return a reference to the selected surfaceWriter.
Definition: surfaceWriter.C:45
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)
static const fileName null
An empty fileName.
Definition: fileName.H:97
virtual bool needsUpdate() const
Does any of the surfaces need an update?
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:309
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
Abstract base-class for Time/database functionObjects.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh - expires the surfaces.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Macros for easy insertion into run-time selection tables.
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
dynamicFvMesh & mesh
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
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.
static const word null
An empty word.
Definition: word.H:77
virtual bool update()=0
Update the surface as required.
static const char nl
Definition: Ostream.H:260
objects
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
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
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
const polyMesh & mesh() const
Return polyMesh.
Definition: mapPolyMesh.H:352
virtual bool read(const dictionary &)
Read the sampledSurfaces dictionary.
messageStream Info
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:88
Registry of regIOobjects.
virtual void readUpdate(const polyMesh::readUpdateState state)
Update for changes of mesh due to readUpdate - expires the surfaces.
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:844