sampledSurfaces.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 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 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(sampledSurfaces, 0);
39 }
40 
41 bool Foam::sampledSurfaces::verbose_ = false;
42 Foam::scalar Foam::sampledSurfaces::mergeTol_ = 1e-10;
43 
44 
45 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46 
47 void Foam::sampledSurfaces::writeGeometry() const
48 {
49  // Write to time directory under outputPath_
50  // Skip surface without faces (eg, a failed cut-plane)
51 
52  const fileName outputDir = outputPath_/mesh_.time().timeName();
53 
54  forAll(*this, surfI)
55  {
56  const sampledSurface& s = operator[](surfI);
57 
58  if (Pstream::parRun())
59  {
60  if (Pstream::master() && mergeList_[surfI].faces.size())
61  {
62  formatter_->write
63  (
64  outputDir,
65  s.name(),
66  mergeList_[surfI].points,
67  mergeList_[surfI].faces
68  );
69  }
70  }
71  else if (s.faces().size())
72  {
73  formatter_->write
74  (
75  outputDir,
76  s.name(),
77  s.points(),
78  s.faces()
79  );
80  }
81  }
82 }
83 
84 
85 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
86 
87 Foam::sampledSurfaces::sampledSurfaces
88 (
89  const word& name,
90  const objectRegistry& obr,
91  const dictionary& dict,
92  const bool loadFromFiles
93 )
94 :
96  name_(name),
97  mesh_(refCast<const fvMesh>(obr)),
98  loadFromFiles_(loadFromFiles),
99  outputPath_(fileName::null),
100  fieldSelection_(),
101  interpolationScheme_(word::null),
102  mergeList_(),
103  formatter_(NULL)
104 {
105  if (Pstream::parRun())
106  {
107  outputPath_ = mesh_.time().path()/".."/"postProcessing"/name_;
108  }
109  else
110  {
111  outputPath_ = mesh_.time().path()/"postProcessing"/name_;
112  }
113 
114  read(dict);
115 }
116 
117 
118 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
119 
121 {}
122 
123 
124 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
125 
126 void Foam::sampledSurfaces::verbose(const bool verbosity)
127 {
128  verbose_ = verbosity;
129 }
130 
131 
133 {
134  // Do nothing - only valid on write
135 }
136 
137 
139 {
140  // Do nothing - only valid on write
141 }
142 
143 
145 {
146  // Do nothing - only valid on write
147 }
148 
149 
151 {
152  if (size())
153  {
154  // Finalize surfaces, merge points etc.
155  update();
156 
157  const label nFields = classifyFields();
158 
159  if (Pstream::master())
160  {
161  if (debug)
162  {
163  Pout<< "Creating directory "
164  << outputPath_/mesh_.time().timeName() << nl << endl;
165 
166  }
167 
168  mkDir(outputPath_/mesh_.time().timeName());
169  }
170 
171  // Write geometry first if required,
172  // or when no fields would otherwise be written
173  if (nFields == 0 || formatter_->separateGeometry())
174  {
175  writeGeometry();
176  }
177 
178  const IOobjectList objects(mesh_, mesh_.time().timeName());
179 
180  sampleAndWrite<volScalarField>(objects);
181  sampleAndWrite<volVectorField>(objects);
182  sampleAndWrite<volSphericalTensorField>(objects);
183  sampleAndWrite<volSymmTensorField>(objects);
184  sampleAndWrite<volTensorField>(objects);
185 
186  sampleAndWrite<surfaceScalarField>(objects);
187  sampleAndWrite<surfaceVectorField>(objects);
188  sampleAndWrite<surfaceSphericalTensorField>(objects);
189  sampleAndWrite<surfaceSymmTensorField>(objects);
190  sampleAndWrite<surfaceTensorField>(objects);
191  }
192 }
193 
194 
196 {
197  bool surfacesFound = dict.found("surfaces");
198 
199  if (surfacesFound)
200  {
201  dict.lookup("fields") >> fieldSelection_;
202 
203  dict.lookup("interpolationScheme") >> interpolationScheme_;
204  const word writeType(dict.lookup("surfaceFormat"));
205 
206  // Define the surface formatter
207  // Optionally defined extra controls for the output formats
208  formatter_ = surfaceWriter::New
209  (
210  writeType,
211  dict.subOrEmptyDict("formatOptions").subOrEmptyDict(writeType)
212  );
213 
215  (
216  dict.lookup("surfaces"),
217  sampledSurface::iNew(mesh_)
218  );
219  transfer(newList);
220 
221  if (Pstream::parRun())
222  {
223  mergeList_.setSize(size());
224  }
225 
226  // Ensure all surfaces and merge information are expired
227  expire();
228 
229  if (this->size())
230  {
231  Info<< "Reading surface description:" << nl;
232  forAll(*this, surfI)
233  {
234  Info<< " " << operator[](surfI).name() << nl;
235  }
236  Info<< endl;
237  }
238  }
239 
240  if (Pstream::master() && debug)
241  {
242  Pout<< "sample fields:" << fieldSelection_ << nl
243  << "sample surfaces:" << nl << "(" << nl;
244 
245  forAll(*this, surfI)
246  {
247  Pout<< " " << operator[](surfI) << endl;
248  }
249  Pout<< ")" << endl;
250  }
251 }
252 
253 
255 {
256  expire();
257 
258  // pointMesh and interpolation will have been reset in mesh.update
259 }
260 
261 
263 {
264  expire();
265 }
266 
267 
269 {
270  if (state != polyMesh::UNCHANGED)
271  {
272  expire();
273  }
274 }
275 
276 
278 {
279  forAll(*this, surfI)
280  {
281  if (operator[](surfI).needsUpdate())
282  {
283  return true;
284  }
285  }
286 
287  return false;
288 }
289 
290 
292 {
293  bool justExpired = false;
294 
295  forAll(*this, surfI)
296  {
297  if (operator[](surfI).expire())
298  {
299  justExpired = true;
300  }
301 
302  // Clear merge information
303  if (Pstream::parRun())
304  {
305  mergeList_[surfI].clear();
306  }
307  }
308 
309  // true if any surfaces just expired
310  return justExpired;
311 }
312 
313 
315 {
316  bool updated = false;
317 
318  if (!needsUpdate())
319  {
320  return updated;
321  }
322 
323  // Serial: quick and easy, no merging required
324  if (!Pstream::parRun())
325  {
326  forAll(*this, surfI)
327  {
328  if (operator[](surfI).update())
329  {
330  updated = true;
331  }
332  }
333 
334  return updated;
335  }
336 
337  // Dimension as fraction of mesh bounding box
338  scalar mergeDim = mergeTol_ * mesh_.bounds().mag();
339 
340  if (Pstream::master() && debug)
341  {
342  Pout<< nl << "Merging all points within "
343  << mergeDim << " metre" << endl;
344  }
345 
346  forAll(*this, surfI)
347  {
348  sampledSurface& s = operator[](surfI);
349 
350  if (s.update())
351  {
352  updated = true;
353  }
354  else
355  {
356  continue;
357  }
358 
360  (
361  mergeDim,
363  (
364  SubList<face>(s.faces(), s.faces().size()),
365  s.points()
366  ),
367  mergeList_[surfI].points,
368  mergeList_[surfI].faces,
369  mergeList_[surfI].pointsMap
370  );
371  }
372 
373  return updated;
374 }
375 
376 
377 // ************************************************************************* //
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:380
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:306
const pointField & points
virtual bool update()
Update the surfaces as required and merge surface points (parallel).
virtual void timeSet()
Called when time was set at the end of the Time::operator++.
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
An abstract class for surfaces with sampling.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
virtual void write()
Sample and write.
virtual void movePoints(const polyMesh &)
Update for mesh point-motion - expires the surfaces.
A class for handling words, derived from string.
Definition: word.H:59
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
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
virtual void execute()
Execute, currently does nothing.
virtual bool needsUpdate() const
Does any of the surfaces need an update?
messageStream Info
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
virtual void readUpdate(const polyMesh::readUpdateState state)
Update for changes of mesh due to readUpdate - expires the surfaces.
Namespace for OpenFOAM.
virtual bool expire()
Mark the surfaces as needing an update.
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.
static const char nl
Definition: Ostream.H:260
const double e
Elementary charge.
Definition: doubleFloat.H:78
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual void end()
Execute at the final time-loop, currently does nothing.
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:88
#define forAll(list, i)
Definition: UList.H:421
void verbose(const bool verbosity=true)
Set verbosity level.
static const fileName null
An empty fileName.
Definition: fileName.H:97
virtual bool update()=0
Update the surface as required.
Istream and Ostream manipulators taking arguments.
virtual void read(const dictionary &)
Read the sampledSurfaces dictionary.
A list of faces which address into the list of points.
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh - expires the surfaces.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:452
Registry of regIOobjects.
Class used for the PtrLists read-construction.
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:398
bool read(const char *, int32_t &)
Definition: int32IO.C:87
virtual const pointField & points() const =0
Points of surface.
A List obtained as a section of another List.
Definition: SubList.H:53
virtual const faceList & faces() const =0
Faces of surface.
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:420
virtual ~sampledSurfaces()
Destructor.
static const word null
An empty word.
Definition: word.H:77
static autoPtr< surfaceWriter > New(const word &writeType)
Return a reference to the selected surfaceWriter.
Definition: surfaceWriter.C:55
defineTypeNameAndDebug(combustionModel, 0)
dictionary subOrEmptyDict(const word &, const bool mustRead=false) const
Find and return a sub-dictionary as a copy, or.
Definition: dictionary.C:675
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53