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-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 "sampledSurfaces.H"
27 #include "PatchTools.H"
28 #include "mapPolyMesh.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  defineTypeNameAndDebug(sampledSurfaces, 0);
36 
38  (
39  functionObject,
40  sampledSurfaces,
41  dictionary
42  );
43 }
44 
45 bool Foam::sampledSurfaces::verbose_ = false;
46 Foam::scalar Foam::sampledSurfaces::mergeTol_ = 1e-10;
47 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
51 void Foam::sampledSurfaces::writeGeometry() const
52 {
53  // Write to time directory under outputPath_
54  // Skip surface without faces (eg, a failed cut-plane)
55 
56  const fileName outputDir = outputPath_/mesh_.time().timeName();
57 
58  forAll(*this, surfI)
59  {
60  const sampledSurface& s = operator[](surfI);
61 
62  if (Pstream::parRun())
63  {
64  if (Pstream::master() && mergeList_[surfI].faces.size())
65  {
66  formatter_->write
67  (
68  outputDir,
69  s.name(),
70  mergeList_[surfI].points,
71  mergeList_[surfI].faces
72  );
73  }
74  }
75  else if (s.faces().size())
76  {
77  formatter_->write
78  (
79  outputDir,
80  s.name(),
81  s.points(),
82  s.faces()
83  );
84  }
85  }
86 }
87 
88 
89 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
90 
92 (
93  const word& name,
94  const Time& t,
95  const dictionary& dict
96 )
97 :
98  functionObject(name),
100  mesh_
101  (
102  refCast<const fvMesh>
103  (
105  (
107  )
108  )
109  ),
110  loadFromFiles_(false),
111  outputPath_(fileName::null),
112  fieldSelection_(),
113  interpolationScheme_(word::null),
114  mergeList_(),
115  formatter_(nullptr)
116 {
117  if (Pstream::parRun())
118  {
119  outputPath_ = mesh_.time().path()/".."/"postProcessing"/name;
120  }
121  else
122  {
123  outputPath_ = mesh_.time().path()/"postProcessing"/name;
124  }
125  // Remove ".."
126  outputPath_.clean();
127 
128  read(dict);
129 }
130 
131 
133 (
134  const word& name,
135  const objectRegistry& obr,
136  const dictionary& dict,
137  const bool loadFromFiles
138 )
139 :
140  functionObject(name),
142  mesh_(refCast<const fvMesh>(obr)),
143  loadFromFiles_(loadFromFiles),
144  outputPath_(fileName::null),
145  fieldSelection_(),
146  interpolationScheme_(word::null),
147  mergeList_(),
148  formatter_(nullptr)
149 {
150  if (Pstream::parRun())
151  {
152  outputPath_ = mesh_.time().path()/".."/"postProcessing"/name;
153  }
154  else
155  {
156  outputPath_ = mesh_.time().path()/"postProcessing"/name;
157  }
158  // Remove ".."
159  outputPath_.clean();
160 
161  read(dict);
162 }
163 
164 
165 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
166 
168 {}
169 
170 
171 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
172 
173 void Foam::sampledSurfaces::verbose(const bool verbosity)
174 {
175  verbose_ = verbosity;
176 }
177 
178 
180 {
181  return true;
182 }
183 
184 
186 {
187  if (size())
188  {
189  // Finalize surfaces, merge points etc.
190  update();
191 
192  const label nFields = classifyFields();
193 
194  if (Pstream::master())
195  {
196  if (debug)
197  {
198  Pout<< "Creating directory "
199  << outputPath_/mesh_.time().timeName() << nl << endl;
200 
201  }
202 
203  mkDir(outputPath_/mesh_.time().timeName());
204  }
205 
206  // Write geometry first if required,
207  // or when no fields would otherwise be written
208  if (nFields == 0 || formatter_->separateGeometry())
209  {
210  writeGeometry();
211  }
212 
213  const IOobjectList objects(mesh_, mesh_.time().timeName());
214 
215  sampleAndWrite<volScalarField>(objects);
216  sampleAndWrite<volVectorField>(objects);
217  sampleAndWrite<volSphericalTensorField>(objects);
218  sampleAndWrite<volSymmTensorField>(objects);
219  sampleAndWrite<volTensorField>(objects);
220 
221  sampleAndWrite<surfaceScalarField>(objects);
222  sampleAndWrite<surfaceVectorField>(objects);
223  sampleAndWrite<surfaceSphericalTensorField>(objects);
224  sampleAndWrite<surfaceSymmTensorField>(objects);
225  sampleAndWrite<surfaceTensorField>(objects);
226  }
227 
228  return true;
229 }
230 
231 
233 {
234  bool surfacesFound = dict.found("surfaces");
235 
236  if (surfacesFound)
237  {
238  dict.lookup("fields") >> fieldSelection_;
239 
240  dict.lookup("interpolationScheme") >> interpolationScheme_;
241  const word writeType(dict.lookup("surfaceFormat"));
242 
243  // Define the surface formatter
244  // Optionally defined extra controls for the output formats
245  formatter_ = surfaceWriter::New
246  (
247  writeType,
248  dict.subOrEmptyDict("formatOptions").subOrEmptyDict(writeType)
249  );
250 
252  (
253  dict.lookup("surfaces"),
254  sampledSurface::iNew(mesh_)
255  );
256  transfer(newList);
257 
258  if (Pstream::parRun())
259  {
260  mergeList_.setSize(size());
261  }
262 
263  // Ensure all surfaces and merge information are expired
264  expire();
265 
266  if (this->size())
267  {
268  Info<< "Reading surface description:" << nl;
269  forAll(*this, surfI)
270  {
271  Info<< " " << operator[](surfI).name() << nl;
272  }
273  Info<< endl;
274  }
275  }
276 
277  if (Pstream::master() && debug)
278  {
279  Pout<< "sample fields:" << fieldSelection_ << nl
280  << "sample surfaces:" << nl << "(" << nl;
281 
282  forAll(*this, surfI)
283  {
284  Pout<< " " << operator[](surfI) << endl;
285  }
286  Pout<< ")" << endl;
287  }
288 
289  return true;
290 }
291 
292 
294 {
295  if (&mpm.mesh() == &mesh_)
296  {
297  expire();
298  }
299 
300  // pointMesh and interpolation will have been reset in mesh.update
301 }
302 
303 
305 {
306  if (&mesh == &mesh_)
307  {
308  expire();
309  }
310 }
311 
312 
314 {
315  if (state != polyMesh::UNCHANGED)
316  {
317  expire();
318  }
319 }
320 
321 
323 {
324  forAll(*this, surfI)
325  {
326  if (operator[](surfI).needsUpdate())
327  {
328  return true;
329  }
330  }
331 
332  return false;
333 }
334 
335 
337 {
338  bool justExpired = false;
339 
340  forAll(*this, surfI)
341  {
342  if (operator[](surfI).expire())
343  {
344  justExpired = true;
345  }
346 
347  // Clear merge information
348  if (Pstream::parRun())
349  {
350  mergeList_[surfI].clear();
351  }
352  }
353 
354  // true if any surfaces just expired
355  return justExpired;
356 }
357 
358 
360 {
361  bool updated = false;
362 
363  if (!needsUpdate())
364  {
365  return updated;
366  }
367 
368  // Serial: quick and easy, no merging required
369  if (!Pstream::parRun())
370  {
371  forAll(*this, surfI)
372  {
373  if (operator[](surfI).update())
374  {
375  updated = true;
376  }
377  }
378 
379  return updated;
380  }
381 
382  // Dimension as fraction of mesh bounding box
383  scalar mergeDim = mergeTol_ * mesh_.bounds().mag();
384 
385  if (Pstream::master() && debug)
386  {
387  Pout<< nl << "Merging all points within "
388  << mergeDim << " metre" << endl;
389  }
390 
391  forAll(*this, surfI)
392  {
393  sampledSurface& s = operator[](surfI);
394 
395  if (s.update())
396  {
397  updated = true;
398  }
399  else
400  {
401  continue;
402  }
403 
405  (
406  mergeDim,
408  (
409  SubList<face>(s.faces(), s.faces().size()),
410  s.points()
411  ),
412  mergeList_[surfI].points,
413  mergeList_[surfI].faces,
414  mergeList_[surfI].pointsMap
415  );
416  }
417 
418  return updated;
419 }
420 
421 
422 // ************************************************************************* //
virtual bool write()
Sample and write.
virtual bool read(const dictionary &)
Read the sampledSurfaces dictionary.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:438
#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 update()
Update the surfaces as required and merge surface points (parallel).
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:158
An abstract class for surfaces with sampling.
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh - expires the surfaces.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
static const fileName null
An empty fileName.
Definition: fileName.H:97
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:309
void verbose(const bool verbosity=true)
Set verbosity level.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
virtual void readUpdate(const polyMesh::readUpdateState state)
Update for changes of mesh due to readUpdate - expires the surfaces.
virtual ~sampledSurfaces()
Destructor.
Abstract base-class for Time/database function objects.
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
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.
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
const pointField & points
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
sampledSurfaces(const word &name, const Time &time, const dictionary &dict)
Construct from Time and dictionary.
virtual bool update()=0
Update the surface as required.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
virtual void movePoints(const polyMesh &)
Update for mesh point-motion - expires the surfaces.
static const char nl
Definition: Ostream.H:265
objects
defineTypeNameAndDebug(combustionModel, 0)
virtual bool execute()
Execute, currently does nothing.
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
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,.
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 expire()
Mark the surfaces as needing an update.
messageStream Info
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:98
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 bool needsUpdate() const
Does any of the surfaces need an update?
static autoPtr< surfaceWriter > New(const word &writeType)
Return a reference to the selected surfaceWriter.
Definition: surfaceWriter.C:44
Class used for the PtrLists read-construction.
dictionary subOrEmptyDict(const word &, const bool mustRead=false) const
Find and return a sub-dictionary as a copy, or.
Definition: dictionary.C:734
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:583