areaWrite.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2019-2023 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "areaWrite.H"
29 #include "polySurface.H"
30 
31 #include "fvMesh.H"
32 #include "mapPolyMesh.H"
33 #include "areaFields.H"
34 #include "HashOps.H"
35 #include "ListOps.H"
36 #include "Time.H"
37 #include "IndirectList.H"
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44  defineTypeNameAndDebug(areaWrite, 0);
45 
47  (
48  functionObject,
49  areaWrite,
50  dictionary
51  );
52 }
53 
54 Foam::scalar Foam::areaWrite::mergeTol_ = 1e-10;
55 
56 
57 // Implementation
58 #include "areaWriteImpl.C"
59 
60 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
61 
62 Foam::areaWrite::areaWrite
63 (
64  const word& name,
65  const Time& runTime,
66  const dictionary& dict
67 )
68 :
69  functionObjects::fvMeshFunctionObject(name, runTime, dict),
70  loadFromFiles_(false),
71  verbose_(false),
72  outputPath_
73  (
74  time_.globalPath()/functionObject::outputPrefix/name
75  ),
76  selectAreas_(),
77  fieldSelection_(),
78  meshes_(),
79  surfaces_(nullptr),
80  writers_()
81 {
82  outputPath_.clean(); // Remove unneeded ".."
83 
84  read(dict);
85 }
86 
87 
88 Foam::areaWrite::areaWrite
89 (
90  const word& name,
91  const objectRegistry& obr,
92  const dictionary& dict,
93  const bool loadFromFiles
94 )
95 :
96  functionObjects::fvMeshFunctionObject(name, obr, dict),
97  loadFromFiles_(loadFromFiles),
98  verbose_(false),
99  outputPath_
100  (
101  time_.globalPath()/functionObject::outputPrefix/name
102  ),
103  selectAreas_(),
104  fieldSelection_(),
105  meshes_(),
106  surfaces_(nullptr),
107  writers_()
108 {
109  outputPath_.clean(); // Remove unneeded ".."
111  read(dict);
112 }
113 
114 
115 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
116 
117 bool Foam::areaWrite::verbose(const bool on) noexcept
118 {
119  bool old(verbose_);
120  verbose_ = on;
121  return old;
122 }
123 
124 
126 {
128 
129  writers_.clear();
130  selectAreas_.clear();
131  fieldSelection_.clear();
132 
133  surfaces_.reset
134  (
135  new objectRegistry
136  (
137  IOobject
138  (
139  "::areaWrite::",
140  obr_.time().constant(),
141  obr_,
145  )
146  )
147  );
148 
149  verbose_ = dict.getOrDefault("verbose", false);
150 
151  // Registry containing all finite-area meshes on the polyMesh
152  const auto* faRegistry = faMesh::registry(mesh_);
153 
154  dict.readIfPresent("areas", selectAreas_);
155 
156  if (selectAreas_.empty())
157  {
158  word areaName;
159 
160  if (!dict.readIfPresent("area", areaName))
161  {
162  if (faRegistry)
163  {
164  wordList available = faRegistry->sortedNames<faMesh>();
165  if (!available.empty())
166  {
167  areaName = available.front();
168  }
169  }
170  }
171 
172  if (!areaName.empty())
173  {
174  selectAreas_.resize(1);
175  selectAreas_.front() = areaName;
176  }
177  }
178 
179  // Restrict to specified meshes
180  meshes_.clear();
181 
182  if (faRegistry)
183  {
184  meshes_ = faRegistry->csorted<faMesh>(selectAreas_);
185  }
186 
187  dict.readEntry("fields", fieldSelection_);
188  fieldSelection_.uniq();
189 
190 
191  // Surface writer type and format options
192  const word writerType = dict.get<word>("surfaceFormat");
193 
194  const dictionary writerOptions
195  (
197  );
198 
199  for (const faMesh& areaMesh : meshes_)
200  {
201  const word& areaName = areaMesh.name();
202 
203  // Define surface writer, but do NOT yet attach a surface
204 
205  auto surfWriter = surfaceWriter::New(writerType, writerOptions);
206 
207  // Use outputDir/TIME/surface-name
208  surfWriter->useTimeDir(true);
209  surfWriter->verbose(verbose_);
210 
211  writers_.set(areaName, surfWriter);
212  }
213 
214  // Ensure all surfaces and merge information are expired
215  expire();
216 
217  return true;
218 }
219 
222 {
223  return true;
224 }
225 
226 
228 {
229  // Just needed for warnings
230  wordList allFields;
231  HashTable<wordHashSet> selected;
232  DynamicList<label> missed(fieldSelection_.size());
233 
234 
235  for (const faMesh& areaMesh : meshes_)
236  {
237  const word& areaName = areaMesh.name();
238 
239  polySurface* surfptr = surfaces_->getObjectPtr<polySurface>(areaName);
240 
241  if (!surfptr)
242  {
243  // Construct null and add to registry (owned by registry)
244  surfptr = new polySurface(areaName, *surfaces_, true);
245  }
246 
247  pointField pts(areaMesh.patch().localPoints());
248  faceList fcs(areaMesh.patch().localFaces());
249 
250  // Copy in geometry
251  surfptr->transfer(std::move(pts), std::move(fcs));
252 
253  surfaceWriter& outWriter = *writers_[areaName];
254 
255  if (outWriter.needsUpdate())
256  {
257  outWriter.setSurface(*surfptr);
258  }
259 
260 
261  // Determine the per-surface number of fields
262  // Only seems to be needed for VTK legacy
263 
264  selected.clear();
265 
266  IOobjectList objects;
267 
268  if (loadFromFiles_)
269  {
270  // Check files for a particular time
271  objects = IOobjectList(areaMesh.thisDb(), obr_.time().timeName());
272 
273  allFields = objects.names();
274  selected = objects.classes(fieldSelection_);
275  }
276  else
277  {
278  // Check currently available fields
279  allFields = areaMesh.thisDb().names();
280  selected = areaMesh.thisDb().classes(fieldSelection_);
281  }
282 
283  // Parallel consistency (no-op in serial)
284  Pstream::mapCombineReduce(selected, HashSetOps::plusEqOp<word>());
285 
286  missed.clear();
287 
288  // Detect missing fields
289  forAll(fieldSelection_, i)
290  {
291  if (!ListOps::found(allFields, fieldSelection_[i]))
292  {
293  missed.push_back(i);
294  }
295  }
296 
297  if (missed.size())
298  {
300  << nl
301  << "Cannot find "
302  << (loadFromFiles_ ? "field file" : "registered field")
303  << " matching "
304  << UIndirectList<wordRe>(fieldSelection_, missed) << endl;
305  }
306 
307 
308  // Currently only support area field types
309  label nAreaFields = 0;
310 
311  forAllConstIters(selected, iter)
312  {
313  const word& clsName = iter.key();
314  const label n = iter.val().size();
315 
316  if
317  (
318  fieldTypes::area.contains(clsName)
319  || fieldTypes::area_internal.contains(clsName)
320  )
321  {
322  nAreaFields += n;
323  }
324  }
325 
326 
327  // Propagate field counts (per surface)
328  outWriter.nFields(nAreaFields);
329 
330 
331  // Begin writing
332 
333  outWriter.open(outputPath_/areaName);
334 
335  outWriter.beginTime(obr_.time());
336 
337  // Write fields
338 
339  {
340  // Area fields
341  #undef doLocalCode
342  #define doLocalCode(Type) \
343  performAction \
344  < \
345  GeometricField<Type, Foam::faPatchField, Foam::areaMesh> \
346  > \
347  ( \
348  outWriter, areaMesh, objects \
349  ); \
350 
351  doLocalCode(scalar);
356 
357  // Area internal fields
358  #undef doLocalCode
359  #define doLocalCode(Type) \
360  performAction \
361  < \
362  DimensionedField<Type, Foam::areaMesh> \
363  > \
364  ( \
365  outWriter, areaMesh, objects \
366  );
367 
368  doLocalCode(scalar);
373 
374  #undef doLocalCode
375  }
376 
377  // Finish this time step
378 
379  // Write geometry if no fields were written so that we still
380  // can have something to look at
381 
382  if (!outWriter.wroteData())
383  {
384  outWriter.write();
385  }
386 
387  outWriter.endTime();
388  }
389 
390  return true;
391 }
392 
393 
394 void Foam::areaWrite::expire()
395 {
396  surfaces_->clear();
397 
398  // Dimension as fraction of mesh bounding box
399  const scalar mergeDim = mergeTol_ * mesh_.bounds().mag();
400 
401  forAllIters(writers_, iter)
402  {
403  surfaceWriter& writer = *(iter.val());
404  writer.expire();
405  writer.mergeDim(mergeDim);
406  }
407 }
408 
409 
410 void Foam::areaWrite::updateMesh(const mapPolyMesh& mpm)
411 {
412  if (&mpm.mesh() == &mesh_)
413  {
414  expire();
415  }
416 }
417 
418 
419 void Foam::areaWrite::movePoints(const polyMesh& mesh)
420 {
421  if (&mesh == &mesh_)
422  {
423  expire();
424  }
425 }
426 
427 
429 {
430  if (state != polyMesh::UNCHANGED)
431  {
432  expire();
433  }
434 }
435 
437 Foam::scalar Foam::areaWrite::mergeTol() noexcept
438 {
439  return mergeTol_;
440 }
441 
442 
443 Foam::scalar Foam::areaWrite::mergeTol(const scalar tol) noexcept
444 {
445  scalar old(mergeTol_);
446  mergeTol_ = tol;
447  return old;
448 }
449 
450 
451 // ************************************************************************* //
A surface mesh consisting of general polygon faces and capable of holding fields. ...
Definition: polySurface.H:62
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition: faMesh.H:87
dictionary dict
virtual void readUpdate(const polyMesh::readUpdateState state)
Update for changes of mesh due to readUpdate - expires the surfaces.
Definition: areaWrite.C:421
vtk::lineWriter writer(edgeCentres, edgeList::null(), fileName(aMesh.time().globalPath()/"finiteArea-edgesCentres"))
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
bool found(const ListType &input, const UnaryPredicate &pred, const label start=0)
Same as found_if.
Definition: ListOps.H:826
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
bool empty() const noexcept
True if List is empty (ie, size() is zero)
Definition: UList.H:666
T & front()
Access first element of the list, position [0].
Definition: UListI.H:237
engineTime & runTime
Tensor< scalar > tensor
Definition: symmTensor.H:57
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Ignore writing from objectRegistry::writeObject()
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
Abstract base-class for Time/database function objects.
static dictionary formatOptions(const dictionary &dict, const word &formatName, const word &entryName="formatOptions")
Same as fileFormats::getFormatOptions.
Definition: surfaceWriter.C:57
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.
virtual void updateMesh(const mapPolyMesh &mpm)
Update for changes of mesh - expires the surfaces.
Definition: areaWrite.C:403
Various functions to operate on Lists.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect...
List< face > faceList
List of faces.
Definition: faceListFwd.H:39
static scalar mergeTol() noexcept
Get merge tolerance.
Definition: areaWrite.C:430
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const wordList area
Standard area field types (scalar, vector, tensor, etc)
SymmTensor< scalar > symmTensor
SymmTensor of scalars, i.e. SymmTensor<scalar>.
Definition: symmTensor.H:55
A class for handling words, derived from Foam::string.
Definition: word.H:63
Type * getObjectPtr(const word &name, const bool recursive=false) const
Return non-const pointer to the object of the given Type, using a const-cast to have it behave like a...
virtual bool write()
Sample and write.
Definition: areaWrite.C:220
void clear()
Remove all entries from table.
Definition: HashTable.C:725
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:336
Vector< scalar > vector
Definition: vector.H:57
A HashTable similar to std::unordered_map.
Definition: HashTable.H:108
const direction noexcept
Definition: Scalar.H:258
#define doLocalCode(Type)
void read(Istream &, label &val, const dictionary &)
In-place read with dictionary lookup.
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
defineTypeNameAndDebug(combustionModel, 0)
static bool clean(std::string &str)
Cleanup filename string, possibly applies other transformations such as changing the path separator e...
Definition: fileName.C:192
bool verbose(const bool on) noexcept
Enable/disable verbose output.
Definition: areaWrite.C:110
void transfer(pointField &&points, faceList &&faces, labelList &&zoneIds=labelList())
Transfer the contents of the argument and annul the argument.
Definition: polySurface.C:364
#define WarningInFunction
Report a warning using Foam::Warning.
const wordList area_internal
Standard dimensioned field types (scalar, vector, tensor, etc)
virtual bool read(const dictionary &dict)
Read the areaWrite dictionary.
Definition: areaWrite.C:118
Mesh data needed to do the Finite Area discretisation.
Definition: areaFaMesh.H:47
Nothing to be read.
label n
virtual void movePoints(const polyMesh &mesh)
Update for mesh point-motion - expires the surfaces.
Definition: areaWrite.C:412
static const objectRegistry * registry(const polyMesh &pMesh)
The parent registry containing all finite-area meshes on the polyMesh.
Definition: faMesh.C:66
Base class for surface writers.
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
SphericalTensor< scalar > sphericalTensor
SphericalTensor of scalars, i.e. SphericalTensor<scalar>.
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:90
Registry of regIOobjects.
static void mapCombineReduce(Container &values, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) applying cop to inplace combine map values from different processo...
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
virtual bool execute()
Execute, currently does nothing.
Definition: areaWrite.C:214
Do not request registration (bool: false)
static autoPtr< surfaceWriter > New(const word &writeType)
Return a reference to the selected surfaceWriter.
Definition: surfaceWriter.C:80
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
const pointField & pts