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-2018 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 "volFields.H"
28 #include "dictionary.H"
29 #include "Time.H"
30 #include "IOmanip.H"
31 #include "volPointInterpolation.H"
32 #include "PatchTools.H"
33 #include "mapPolyMesh.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  defineTypeNameAndDebug(sampledSurfaces, 0);
41 
43  (
44  functionObject,
45  sampledSurfaces,
46  dictionary
47  );
48 }
49 
50 bool Foam::sampledSurfaces::verbose_ = false;
51 Foam::scalar Foam::sampledSurfaces::mergeTol_ = 1e-10;
52 
53 
54 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
55 
56 void Foam::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 
96 Foam::sampledSurfaces::sampledSurfaces
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  if (Pstream::parRun())
123  {
124  outputPath_ = mesh_.time().path()/".."/"postProcessing"/name;
125  }
126  else
127  {
128  outputPath_ = mesh_.time().path()/"postProcessing"/name;
129  }
130  // Remove ".."
131  outputPath_.clean();
132 
133  read(dict);
134 }
135 
136 
137 Foam::sampledSurfaces::sampledSurfaces
138 (
139  const word& name,
140  const objectRegistry& obr,
141  const dictionary& dict,
142  const bool loadFromFiles
143 )
144 :
145  functionObject(name),
147  mesh_(refCast<const fvMesh>(obr)),
148  loadFromFiles_(loadFromFiles),
149  outputPath_(fileName::null),
150  fieldSelection_(),
151  interpolationScheme_(word::null),
152  mergeList_(),
153  formatter_(nullptr)
154 {
155  if (Pstream::parRun())
156  {
157  outputPath_ = mesh_.time().path()/".."/"postProcessing"/name;
158  }
159  else
160  {
161  outputPath_ = mesh_.time().path()/"postProcessing"/name;
162  }
163  // Remove ".."
164  outputPath_.clean();
165 
166  read(dict);
167 }
168 
169 
170 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
171 
173 {}
174 
175 
176 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
177 
178 void Foam::sampledSurfaces::verbose(const bool verbosity)
179 {
180  verbose_ = verbosity;
181 }
182 
183 
185 {
186  return true;
187 }
188 
189 
191 {
192  if (size())
193  {
194  // Finalize surfaces, merge points etc.
195  update();
196 
197  const label nFields = classifyFields();
198 
199  if (Pstream::master())
200  {
201  if (debug)
202  {
203  Pout<< "Creating directory "
204  << outputPath_/mesh_.time().timeName() << nl << endl;
205 
206  }
207 
208  mkDir(outputPath_/mesh_.time().timeName());
209  }
210 
211  // Write geometry first if required,
212  // or when no fields would otherwise be written
213  if (nFields == 0 || formatter_->separateGeometry())
214  {
215  writeGeometry();
216  }
217 
218  const IOobjectList objects(mesh_, mesh_.time().timeName());
219 
220  sampleAndWrite<volScalarField>(objects);
221  sampleAndWrite<volVectorField>(objects);
222  sampleAndWrite<volSphericalTensorField>(objects);
223  sampleAndWrite<volSymmTensorField>(objects);
224  sampleAndWrite<volTensorField>(objects);
225 
226  sampleAndWrite<surfaceScalarField>(objects);
227  sampleAndWrite<surfaceVectorField>(objects);
228  sampleAndWrite<surfaceSphericalTensorField>(objects);
229  sampleAndWrite<surfaceSymmTensorField>(objects);
230  sampleAndWrite<surfaceTensorField>(objects);
231  }
232 
233  return true;
234 }
235 
236 
238 {
239  bool surfacesFound = dict.found("surfaces");
240 
241  if (surfacesFound)
242  {
243  dict.lookup("fields") >> fieldSelection_;
244 
245  dict.lookup("interpolationScheme") >> interpolationScheme_;
246  const word writeType(dict.lookup("surfaceFormat"));
247 
248  // Define the surface formatter
249  // Optionally defined extra controls for the output formats
250  formatter_ = surfaceWriter::New
251  (
252  writeType,
253  dict.subOrEmptyDict("formatOptions").subOrEmptyDict(writeType)
254  );
255 
257  (
258  dict.lookup("surfaces"),
259  sampledSurface::iNew(mesh_)
260  );
261  transfer(newList);
262 
263  if (Pstream::parRun())
264  {
265  mergeList_.setSize(size());
266  }
267 
268  // Ensure all surfaces and merge information are expired
269  expire();
270 
271  if (this->size())
272  {
273  Info<< "Reading surface description:" << nl;
274  forAll(*this, surfI)
275  {
276  Info<< " " << operator[](surfI).name() << nl;
277  }
278  Info<< endl;
279  }
280  }
281 
282  if (Pstream::master() && debug)
283  {
284  Pout<< "sample fields:" << fieldSelection_ << nl
285  << "sample surfaces:" << nl << "(" << nl;
286 
287  forAll(*this, surfI)
288  {
289  Pout<< " " << operator[](surfI) << endl;
290  }
291  Pout<< ")" << endl;
292  }
293 
294  return true;
295 }
296 
297 
299 {
300  if (&mpm.mesh() == &mesh_)
301  {
302  expire();
303  }
304 
305  // pointMesh and interpolation will have been reset in mesh.update
306 }
307 
308 
310 {
311  if (&mesh == &mesh_)
312  {
313  expire();
314  }
315 }
316 
317 
319 {
320  if (state != polyMesh::UNCHANGED)
321  {
322  expire();
323  }
324 }
325 
326 
328 {
329  forAll(*this, surfI)
330  {
331  if (operator[](surfI).needsUpdate())
332  {
333  return true;
334  }
335  }
336 
337  return false;
338 }
339 
340 
342 {
343  bool justExpired = false;
344 
345  forAll(*this, surfI)
346  {
347  if (operator[](surfI).expire())
348  {
349  justExpired = true;
350  }
351 
352  // Clear merge information
353  if (Pstream::parRun())
354  {
355  mergeList_[surfI].clear();
356  }
357  }
358 
359  // true if any surfaces just expired
360  return justExpired;
361 }
362 
363 
365 {
366  bool updated = false;
367 
368  if (!needsUpdate())
369  {
370  return updated;
371  }
372 
373  // Serial: quick and easy, no merging required
374  if (!Pstream::parRun())
375  {
376  forAll(*this, surfI)
377  {
378  if (operator[](surfI).update())
379  {
380  updated = true;
381  }
382  }
383 
384  return updated;
385  }
386 
387  // Dimension as fraction of mesh bounding box
388  scalar mergeDim = mergeTol_ * mesh_.bounds().mag();
389 
390  if (Pstream::master() && debug)
391  {
392  Pout<< nl << "Merging all points within "
393  << mergeDim << " metre" << endl;
394  }
395 
396  forAll(*this, surfI)
397  {
398  sampledSurface& s = operator[](surfI);
399 
400  if (s.update())
401  {
402  updated = true;
403  }
404  else
405  {
406  continue;
407  }
408 
410  (
411  mergeDim,
413  (
414  SubList<face>(s.faces(), s.faces().size()),
415  s.points()
416  ),
417  mergeList_[surfI].points,
418  mergeList_[surfI].faces,
419  mergeList_[surfI].pointsMap
420  );
421  }
422 
423  return updated;
424 }
425 
426 
427 // ************************************************************************* //
static void gatherAndMerge(const scalar mergeDist, const PrimitivePatch< Face, FaceList, PointField, PointType > &p, Field< PointType > &mergedPoints, List< Face > &mergedFaces, labelList &pointMergeMap)
Gather points and faces onto master and merge into single patch.
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:431
#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
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:137
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:421
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 const word null
An empty word.
Definition: word.H:77
virtual bool update()=0
Update the surface as required.
Istream and Ostream manipulators taking arguments.
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
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:289
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:397
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:357
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:55
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:727
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