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-2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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
42namespace Foam
43{
45
47 (
51 );
52}
53
54Foam::scalar Foam::areaWrite::mergeTol_ = 1e-10;
55
56// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
57
59(
60 const word& name,
61 const Time& runTime,
62 const dictionary& dict
63)
64:
65 functionObjects::fvMeshFunctionObject(name, runTime, dict),
66 loadFromFiles_(false),
67 verbose_(false),
68 outputPath_
69 (
70 time_.globalPath()/functionObject::outputPrefix/name
71 ),
72 selectAreas_(),
73 fieldSelection_(),
74 meshes_(),
75 surfaces_(nullptr),
76 writers_()
77{
78 outputPath_.clean(); // Remove unneeded ".."
79
80 read(dict);
81}
82
83
85(
86 const word& name,
87 const objectRegistry& obr,
88 const dictionary& dict,
89 const bool loadFromFiles
90)
91:
92 functionObjects::fvMeshFunctionObject(name, obr, dict),
93 loadFromFiles_(loadFromFiles),
94 verbose_(false),
95 outputPath_
96 (
97 time_.globalPath()/functionObject::outputPrefix/name
98 ),
99 selectAreas_(),
100 fieldSelection_(),
101 meshes_(),
102 surfaces_(nullptr),
103 writers_()
104{
105 outputPath_.clean(); // Remove unneeded ".."
106
107 read(dict);
108}
109
110
111// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
112
113bool Foam::areaWrite::verbose(const bool on)
114{
115 bool old(verbose_);
116 verbose_ = on;
117 return old;
118}
119
120
122{
123 fvMeshFunctionObject::read(dict);
124
125 writers_.clear();
126 selectAreas_.clear();
127 fieldSelection_.clear();
128
129 surfaces_.reset
130 (
132 (
134 (
135 "::areaWrite::",
136 obr_.time().constant(),
137 obr_,
140 false
141 )
142 )
143 );
144
145 verbose_ = dict.getOrDefault("verbose", false);
146
147 // All possible area meshes for the given fvMesh region
148 meshes_ = obr().lookupClass<faMesh>();
149
150 dict.readIfPresent("areas", selectAreas_);
151
152 if (selectAreas_.empty())
153 {
154 word areaName;
155 if (!dict.readIfPresent("area", areaName))
156 {
157 wordList available = obr().sortedNames<faMesh>();
158
159 if (available.size())
160 {
161 areaName = available.first();
162 }
163 }
164
165 if (!areaName.empty())
166 {
167 selectAreas_.resize(1);
168 selectAreas_.first() = areaName;
169 }
170 }
171
172 // Restrict to specified meshes
173 meshes_.filterKeys(selectAreas_);
174
175 dict.readEntry("fields", fieldSelection_);
176 fieldSelection_.uniq();
177
178
179 // Surface writer type and format options
180 const word writerType = dict.get<word>("surfaceFormat");
181
182 const dictionary writerOptions
183 (
184 dict.subOrEmptyDict("formatOptions").subOrEmptyDict(writerType)
185 );
186
187 for (const word& areaName : meshes_.keys())
188 {
189 // Define surface writer, but do NOT yet attach a surface
190
191 auto surfWriter = surfaceWriter::New(writerType, writerOptions);
192
193 // Use outputDir/TIME/surface-name
194 surfWriter->useTimeDir(true);
195 surfWriter->verbose(verbose_);
196
197 writers_.set(areaName, surfWriter);
198 }
199
200 // Ensure all surfaces and merge information are expired
201 expire();
202
203 return true;
204}
205
206
208{
209 return true;
210}
211
212
214{
215 // Just needed for warnings
216 wordList allFields;
217 HashTable<wordHashSet> selected;
218 DynamicList<label> missed(fieldSelection_.size());
219
220
221 for (const word& areaName : meshes_.sortedToc())
222 {
223 const faMesh& areaMesh = *meshes_[areaName];
224
225 polySurface* surfptr = surfaces_->getObjectPtr<polySurface>(areaName);
226
227 if (!surfptr)
228 {
229 // Construct null and add to registry (owned by registry)
230 surfptr = new polySurface(areaName, *surfaces_, true);
231 }
232
233 pointField pts(areaMesh.patch().localPoints());
234 faceList fcs(areaMesh.patch().localFaces());
235
236 // Copy in geometry
237 surfptr->transfer(std::move(pts), std::move(fcs));
238
239 surfaceWriter& outWriter = *writers_[areaName];
240
241 if (outWriter.needsUpdate())
242 {
243 outWriter.setSurface(*surfptr);
244 }
245
246
247 // Determine the per-surface number of fields
248 // Only seems to be needed for VTK legacy
249
250 selected.clear();
251
252 IOobjectList objects(0);
253
254 if (loadFromFiles_)
255 {
256 // Check files for a particular time
257 objects = IOobjectList(areaMesh.thisDb(), obr_.time().timeName());
258
259 allFields = objects.names();
260 selected = objects.classes(fieldSelection_);
261 }
262 else
263 {
264 // Check currently available fields
265 allFields = areaMesh.thisDb().names();
266 selected = areaMesh.thisDb().classes(fieldSelection_);
267 }
268
269 // Parallel consistency (no-op in serial)
271
272 missed.clear();
273
274 // Detect missing fields
275 forAll(fieldSelection_, i)
276 {
277 if (!ListOps::found(allFields, fieldSelection_[i]))
278 {
279 missed.append(i);
280 }
281 }
282
283 if (missed.size())
284 {
286 << nl
287 << "Cannot find "
288 << (loadFromFiles_ ? "field file" : "registered field")
289 << " matching "
290 << UIndirectList<wordRe>(fieldSelection_, missed) << endl;
291 }
292
293
294 // Currently only support area field types
295 label nAreaFields = 0;
296
297 forAllConstIters(selected, iter)
298 {
299 const word& clsName = iter.key();
300 const label n = iter.val().size();
301
302 if (fieldTypes::area.found(clsName))
303 {
304 nAreaFields += n;
305 }
306 }
307
308
309 // Propagate field counts (per surface)
310 outWriter.nFields(nAreaFields);
311
312
313 // Begin writing
314
315 outWriter.open(outputPath_/areaName);
316
317 outWriter.beginTime(obr_.time());
318
319 // Write fields
320
321 performAction<areaScalarField>(outWriter, areaMesh, objects);
322 performAction<areaVectorField>(outWriter, areaMesh, objects);
323 performAction<areaSphericalTensorField>(outWriter, areaMesh, objects);
324 performAction<areaSymmTensorField>(outWriter, areaMesh, objects);
325 performAction<areaTensorField>(outWriter, areaMesh, objects);
326
327 // Finish this time step
328
329 // Write geometry if no fields were written so that we still
330 // can have something to look at
331
332 if (!outWriter.wroteData())
333 {
334 outWriter.write();
335 }
336
337 outWriter.endTime();
338 }
339
340 return true;
341}
342
343
344void Foam::areaWrite::expire()
345{
346 surfaces_->clear();
347
348 // Dimension as fraction of mesh bounding box
349 const scalar mergeDim = mergeTol_ * mesh_.bounds().mag();
350
351 forAllIters(writers_, iter)
352 {
353 surfaceWriter& writer = *(iter.val());
354 writer.expire();
355 writer.mergeDim(mergeDim);
356 }
357}
358
359
361{
362 if (&mpm.mesh() == &mesh_)
363 {
364 expire();
365 }
366}
367
368
370{
371 if (&mesh == &mesh_)
372 {
373 expire();
374 }
375}
376
377
379{
380 if (state != polyMesh::UNCHANGED)
381 {
382 expire();
383 }
384}
385
386
388{
389 return mergeTol_;
390}
391
392
393Foam::scalar Foam::areaWrite::mergeTol(const scalar tol)
394{
395 const scalar prev(mergeTol_);
396 mergeTol_ = tol;
397 return prev;
398}
399
400
401// ************************************************************************* //
Various functions to operate on Lists.
bool found
label n
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
vtk::internalMeshWriter writer(topoMesh, topoCells, vtk::formatType::INLINE_ASCII, runTime.path()/"blockTopology")
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:391
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:503
const objectRegistry & thisDb() const
Return the object registry.
Definition: GeoMesh.H:84
A HashTable similar to std::unordered_map.
Definition: HashTable.H:123
void clear()
Clear all entries from table.
Definition: HashTable.C:678
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:59
HashTable< wordHashSet > classes() const
A summary hash of classes used and their associated object names.
Definition: IOobjectList.C:320
wordList names() const
The unsorted names of the IOobjects.
Definition: IOobjectList.C:351
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
static void mapCombineAllGather(const List< commsStruct > &comms, Container &values, const CombineOp &cop, const int tag, const label comm)
After completion all processors have the same data.
virtual bool read()
Re-read model coefficients if they have changed.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
A List with indirect addressing. Like IndirectList but does not store addressing.
Definition: IndirectList.H:79
T & first()
Return the first element of the list.
Definition: UListI.H:202
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Mesh data needed to do the Finite Area discretisation.
Definition: areaFaMesh.H:56
Write finite area mesh/fields to standard output formats.
Definition: areaWrite.H:143
virtual bool read(const dictionary &dict)
Read the areaWrite dictionary.
Definition: areaWrite.C:121
static scalar mergeTol()
Get merge tolerance.
Definition: areaWrite.C:387
virtual bool execute()
Execute, currently does nothing.
Definition: areaWrite.C:207
virtual bool write()
Sample and write.
Definition: areaWrite.C:213
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
Definition: dictionary.C:540
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
int verbose() const noexcept
Output verbosity level.
Definition: ensightMesh.C:125
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition: faMesh.H:100
static bool clean(std::string &str)
Definition: fileName.C:199
Abstract base-class for Time/database function objects.
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in time.
Definition: fvMesh.C:675
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:162
const polyMesh & mesh() const
Return polyMesh.
Definition: mapPolyMesh.H:363
void movePoints()
Update for new mesh geometry.
void updateMesh()
Update for new mesh topology.
Registry of regIOobjects.
HashTable< wordHashSet > classes() const
A summary hash of classes used and their associated object names.
const Time & time() const noexcept
Return time registry.
wordList names() const
The unsorted names of all objects.
Type * getObjectPtr(const word &name, const bool recursive=false) const
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:91
A surface mesh consisting of general polygon faces and capable of holding fields.
Definition: polySurface.H:71
void transfer(pointField &&points, faceList &&faces, labelList &&zoneIds=labelList())
Transfer the contents of the argument and annul the argument.
Definition: polySurface.C:371
Base class for surface writers.
virtual void open(const fileName &outputPath)
Open for output on specified path, using existing surface.
virtual void endTime()
End a time-step.
virtual fileName write()=0
Write separate surface geometry to file.
virtual void beginTime(const Time &t)
Begin a time-step.
virtual void setSurface(const meshedSurf &surf, bool parallel)
virtual bool needsUpdate() const
Does the writer need an update (eg, lagging behind surface changes)
virtual bool wroteData() const
Geometry or fields written since the last open?
label nFields() const noexcept
The number of expected output fields.
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
dynamicFvMesh & mesh
engineTime & runTime
#define WarningInFunction
Report a warning using Foam::Warning.
bool found(const ListType &input, const UnaryPredicate &pred, const label start=0)
Same as found_if.
Definition: ListOps.H:679
const wordList area
Standard area field types (scalar, vector, tensor, etc)
Namespace for OpenFOAM.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:260
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278
Combine HashSet operation. Equivalent to 'a |= b'.
Definition: HashOps.H:67