vtkWrite.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) 2017-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 "vtkWrite.H"
29#include "dictionary.H"
30#include "Time.H"
31#include "areaFields.H"
32#include "stringListOps.H"
34#include "foamVtkPatchWriter.H"
35#include "foamVtkSeriesWriter.H"
36#include "foamVtmWriter.H"
38
39// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40
41namespace Foam
42{
43namespace functionObjects
44{
47}
48}
49
50
51// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52
53Foam::label Foam::functionObjects::vtkWrite::writeAllVolFields
54(
55 autoPtr<vtk::internalWriter>& internalWriter,
56 UPtrList<vtk::patchWriter>& patchWriters,
57 const fvMeshSubset& proxy,
58 const wordHashSet& acceptField
59) const
60{
61 #undef vtkWrite_WRITE_FIELD
62 #define vtkWrite_WRITE_FIELD(FieldType) \
63 writeVolFields<FieldType> \
64 ( \
65 internalWriter, \
66 patchWriters, \
67 proxy, \
68 acceptField \
69 )
70
71
72 label count = 0;
78
79 #undef vtkWrite_WRITE_FIELD
80 return count;
81}
82
83
84Foam::label Foam::functionObjects::vtkWrite::writeAllVolFields
85(
86 autoPtr<vtk::internalWriter>& internalWriter,
87 const autoPtr<volPointInterpolation>& pInterp,
88
89 UPtrList<vtk::patchWriter>& patchWriters,
90 const UPtrList<PrimitivePatchInterpolation<primitivePatch>>& patchInterps,
91 const fvMeshSubset& proxy,
92 const wordHashSet& acceptField
93) const
94{
95 #undef vtkWrite_WRITE_FIELD
96 #define vtkWrite_WRITE_FIELD(FieldType) \
97 writeVolFields<FieldType> \
98 ( \
99 internalWriter, pInterp, \
100 patchWriters, patchInterps, \
101 proxy, \
102 acceptField \
103 )
104
105
106 label count = 0;
112
113 #undef vtkWrite_WRITE_FIELD
114 return count;
115}
116
117
118// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
119
121(
122 const word& name,
123 const Time& runTime,
124 const dictionary& dict
125)
126:
128 outputDir_(),
129 printf_(),
130 writeOpts_(vtk::formatType::INLINE_BASE64),
131 verbose_(true),
132 doInternal_(true),
133 doBoundary_(true),
134 oneBoundary_(false),
135 interpolate_(false),
136 decompose_(false),
137 writeIds_(false),
138 meshState_(polyMesh::TOPO_CHANGE),
139 selectRegions_(),
140 selectPatches_(),
141 selectFields_(),
142 selection_(),
143 meshes_(),
144 meshSubsets_(),
145 vtuMappings_(),
146 series_()
147{
148 // May still want this? (OCT-2018)
149 // if (postProcess)
150 // {
151 // // Disable for post-process mode.
152 // // Emit as FatalError for the try/catch in the caller.
153 // FatalError
154 // << type() << " disabled in post-process mode"
155 // << exit(FatalError);
156 // }
157
158 read(dict);
159}
160
161
162// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
163
165{
167
168 readSelection(dict);
169
170 // We probably cannot trust old information after a reread
171 series_.clear();
172
173 // verbose_ = dict.getOrDefault("verbose", true);
174 doInternal_ = dict.getOrDefault("internal", true);
175 doBoundary_ = dict.getOrDefault("boundary", true);
176 oneBoundary_ = dict.getOrDefault("single", false);
177 interpolate_ = dict.getOrDefault("interpolate", false);
178
179 //
180 // Writer options - default is xml base64
181 //
183
184 writeOpts_.ascii
185 (
188 );
189
190 writeOpts_.legacy(dict.getOrDefault("legacy", false));
191
192 writeOpts_.precision
193 (
194 dict.getOrDefault("precision", IOstream::defaultPrecision())
195 );
196
197 // Info<< type() << " " << name() << " output-format: "
198 // << writeOpts_.description() << nl;
199
200 const int padWidth = dict.getOrDefault<int>("width", 8);
201
202 // Appropriate printf format - Enforce min/max sanity limits
203 if (padWidth < 1 || padWidth > 31)
204 {
205 printf_.clear();
206 }
207 else
208 {
209 printf_ = "%0" + std::to_string(padWidth) + "d";
210 }
211
212 //
213 // Other options
214 //
215
216 decompose_ = dict.getOrDefault("decompose", false);
217 writeIds_ = dict.getOrDefault("writeIds", false);
218
219
220 // Output directory
221
222 outputDir_.clear();
223 dict.readIfPresent("directory", outputDir_);
224
225 if (outputDir_.size())
226 {
227 // User-defined output directory
228 outputDir_.expand();
229 if (!outputDir_.isAbsolute())
230 {
231 outputDir_ = time_.globalPath()/outputDir_;
232 }
233 }
234 else
235 {
236 // Standard postProcessing/ naming
237 outputDir_ = time_.globalPath()/functionObject::outputPrefix/name();
238 }
239 outputDir_.clean(); // Remove unneeded ".."
240
241 return true;
242}
243
244
246{
247 return true;
248}
249
250
252{
253 // const word timeDesc =
254 // useTimeName ? time_.timeName() : Foam::name(time_.timeIndex());
255
256 const word timeDesc = "_" +
257 (
258 printf_.empty()
259 ? Foam::name(time_.timeIndex())
260 : word::printf(printf_, time_.timeIndex())
261 );
262
263 const scalar timeValue = time_.value();
264
265 update();
266
267 if (meshes_.empty() || (!doInternal_ && !doBoundary_))
268 {
269 // Skip
270 return true;
271 }
272
273
274 fileName vtkName = time_.globalCaseName();
275
276 vtk::vtmWriter vtmMultiRegion;
277
278 Info<< name() << " output Time: " << time_.timeName() << nl;
279
280 label regioni = 0;
281 for (const word& regionName : meshes_.sortedToc())
282 {
284
285 auto& meshProxy = meshSubsets_[regioni];
286 auto& vtuMeshCells = vtuMappings_[regioni];
287 ++regioni;
288
289 const fvMesh& baseMesh = meshProxy.baseMesh();
290
291 wordHashSet acceptField(baseMesh.names<void>(selectFields_));
292
293 // Prune restart fields
294 acceptField.filterKeys
295 (
296 [](const word& k){ return k.ends_with("_0"); },
297 true // prune
298 );
299
300 const label nVolFields =
301 (
302 (doInternal_ || doBoundary_)
303 ? baseMesh.count
304 (
306 acceptField
307 )
308 : 0
309 );
310
311 // Undecided if we want to automatically support DimensionedFields
312 // or only on demand:
313 const label nDimFields = 0;
314 // (
315 // (doInternal_ || doBoundary_)
316 // ? baseMesh.count
317 // (
318 // stringListOps::foundOp<word>(fieldTypes::internal),
319 // acceptField
320 // )
321 // : 0
322 // );
323
324
325 // Setup for the vtm writer.
326 // For legacy format, the information added is simply ignored.
327
329 (
330 outputDir_/regionDir/vtkName + timeDesc
331 );
332
333 // Combined internal + boundary in a vtm file
335
336 // Collect individual boundaries into a vtm file
338
339 // Setup the internal writer
341
342 // Interpolator for volume and dimensioned fields
344
345 if (doInternal_)
346 {
347 if (vtuMeshCells.empty())
348 {
349 // Use the appropriate mesh (baseMesh or subMesh)
350 vtuMeshCells.reset(meshProxy.mesh());
351
352 if (interpolate_ && vtuMeshCells.manifold())
353 {
354 interpolate_ = false;
356 << "Manifold cells detected - disabling PointData"
357 << endl;
358 }
359 }
360
362 (
363 meshProxy.mesh(),
364 vtuMeshCells,
365 writeOpts_,
366 // Output name for internal
367 (
368 writeOpts_.legacy()
370 : (vtmOutputBase / "internal")
371 ),
373 );
374
375 Info<< " Internal : "
376 << time_.relativePath(internalWriter->output())
377 << endl;
378
379 // No sub-block for internal
380 vtmWriter.append_vtu
381 (
382 "internal",
383 vtmOutputBase.name()/"internal"
384 );
385
386 internalWriter->writeTimeValue(timeValue);
387 internalWriter->writeGeometry();
388
389 if (interpolate_)
390 {
391 pInterp.reset(new volPointInterpolation(meshProxy.mesh()));
392 }
393 }
394
395
396 // Setup the patch writers
397
398 const polyBoundaryMesh& patches = meshProxy.mesh().boundaryMesh();
399
402
404 if (doBoundary_)
405 {
406 patchIds = getSelectedPatches(patches);
407 }
408
409 if (oneBoundary_ && patchIds.size())
410 {
412 (
413 meshProxy.mesh(),
414 patchIds,
415 writeOpts_,
416 // Output name for one patch: "boundary"
417 (
418 writeOpts_.legacy()
419 ? (outputDir_/regionDir/"boundary"/"boundary" + timeDesc)
420 : (vtmOutputBase / "boundary")
421 ),
423 );
424
425 // No sub-block for one-patch
426 vtmWriter.append_vtp
427 (
428 "boundary",
429 vtmOutputBase.name()/"boundary"
430 );
431
432 Info<< " Boundaries: "
433 << time_.relativePath(writer->output()) << nl;
434
435
436 writer->writeTimeValue(timeValue);
437 writer->writeGeometry();
438
439 // Transfer writer to list for later use
440 patchWriters.resize(1);
441 patchWriters.set(0, writer);
442
443 // Avoid patchInterpolation for each sub-patch
444 patchInterps.resize(1); // == nullptr
445 }
446 else if (patchIds.size())
447 {
448 patchWriters.resize(patchIds.size());
449 if (interpolate_)
450 {
451 patchInterps.resize(patchIds.size());
452 }
453
454 label nPatchWriters = 0;
455 label nPatchInterps = 0;
456
457 for (const label patchId : patchIds)
458 {
459 const polyPatch& pp = patches[patchId];
460
462 (
463 meshProxy.mesh(),
464 labelList(one{}, pp.index()),
465 writeOpts_,
466 // Output name for patch: "boundary"/name
467 (
468 writeOpts_.legacy()
469 ?
470 (
471 outputDir_/regionDir/pp.name()
472 / (pp.name()) + timeDesc
473 )
474 : (vtmOutputBase / "boundary" / pp.name())
475 ),
477 );
478
479 if (!nPatchWriters)
480 {
481 vtmWriter.beginBlock("boundary");
482 vtmBoundaries.beginBlock("boundary");
483 }
484
485 vtmWriter.append_vtp
486 (
487 pp.name(),
488 vtmOutputBase.name()/"boundary"/pp.name()
489 );
490
491 vtmBoundaries.append_vtp
492 (
493 pp.name(),
494 "boundary"/pp.name()
495 );
496
497 Info<< " Boundary : "
498 << time_.relativePath(writer->output()) << nl;
499
500 writer->writeTimeValue(timeValue);
501 writer->writeGeometry();
502
503 // Transfer writer to list for later use
505
506 if (patchInterps.size())
507 {
508 patchInterps.set
509 (
510 nPatchInterps++,
512 );
513 }
514 }
515
516 if (nPatchWriters)
517 {
518 vtmWriter.endBlock("boundary");
519 vtmBoundaries.endBlock("boundary");
520 }
521
523 patchInterps.resize(nPatchInterps);
524 }
525
526 // CellData
527 {
528 if (internalWriter)
529 {
530 // Optionally with cellID and procID fields
531 internalWriter->beginCellData
532 (
533 (writeIds_ ? 1 + (internalWriter->parallel() ? 1 : 0) : 0)
535 );
536
537 if (writeIds_)
538 {
539 internalWriter->writeCellIDs();
540 internalWriter->writeProcIDs(); // parallel only
541 }
542 }
543
544 if (nVolFields)
545 {
547 {
548 // Optionally with patchID field
549 writer.beginCellData
550 (
551 (writeIds_ ? 1 : 0)
552 + nVolFields
553 );
554
555 if (writeIds_)
556 {
557 writer.writePatchIDs();
558 }
559 }
560 }
561
563 (
566 meshProxy,
567 acceptField
568 );
569
570 // writeAllDimFields
571 // (
572 // internalWriter,
573 // meshProxy,
574 // acceptField
575 // );
576
577 // End CellData is implicit
578 }
579
580
581 // PointData
582 // - only construct pointMesh on request since it constructs
583 // edge addressing
584 if (interpolate_)
585 {
586 // Begin PointData
587 if (internalWriter)
588 {
589 internalWriter->beginPointData
590 (
592 );
593 }
594
595 forAll(patchWriters, writeri)
596 {
597 const label nPatchFields =
598 (
599 writeri < patchInterps.size() && patchInterps.set(writeri)
600 ? nVolFields
601 : 0
602 );
603
604 if (nPatchFields)
605 {
606 patchWriters[writeri].beginPointData(nPatchFields);
607 }
608 }
609
611 (
614 meshProxy,
615 acceptField
616 );
617
618 // writeAllDimFields
619 // (
620 // internalWriter, pInterp,
621 // meshProxy,
622 // acceptField
623 // );
624
625 // writeAllPointFields
626 // (
627 // internalWriter,
628 // patchWriters,
629 // meshProxy,
630 // acceptField
631 // );
632
633 // End PointData is implicit
634 }
635
636
637 // Finish writers
638 if (internalWriter)
639 {
640 internalWriter->close();
641 }
642
644 {
645 writer.close();
646 }
647
648 pInterp.clear();
649 patchWriters.clear();
650 patchInterps.clear();
651
652
653 // Collective output
654
655 if (Pstream::master())
656 {
657 // Naming for vtm, file series etc.
659
660 if (writeOpts_.legacy())
661 {
662 if (doInternal_)
663 {
664 // Add to file-series and emit as JSON
665
667
669
670 vtk::seriesWriter& series = series_(seriesName);
671
672 // First time?
673 // Load from file, verify against filesystem,
674 // prune time >= currentTime
675 if (series.empty())
676 {
677 series.load(seriesName, true, timeValue);
678 }
679
680 series.append(timeValue, timeDesc);
681 series.write(seriesName);
682 }
683 }
684 else
685 {
686 if (vtmWriter.size())
687 {
688 // Emit ".vtm"
689
690 outputName.ext(vtmWriter.ext());
691
692 vtmWriter.setTime(timeValue);
693 vtmWriter.write(outputName);
694
696
697 vtk::seriesWriter& series = series_(seriesName);
698
699 // First time?
700 // Load from file, verify against filesystem,
701 // prune time >= currentTime
702 if (series.empty())
703 {
704 series.load(seriesName, true, timeValue);
705 }
706
707 series.append(timeValue, outputName);
708 series.write(seriesName);
709
710 // Add to multi-region vtm
711 vtmMultiRegion.add(regionName, regionDir, vtmWriter);
712 }
713
714 if (vtmBoundaries.size())
715 {
716 // Emit "boundary.vtm" with collection of boundaries
717
718 // Naming for vtm
719 fileName outputName(vtmOutputBase / "boundary");
720 outputName.ext(vtmBoundaries.ext());
721
722 vtmBoundaries.setTime(timeValue);
724 }
725 }
726 }
727 }
728
729
730 // Emit multi-region vtm
731 if (Pstream::master() && meshes_.size() > 1)
732 {
734 (
735 outputDir_/vtkName + "-regions" + timeDesc + ".vtm"
736 );
737
738 vtmMultiRegion.setTime(timeValue);
739 vtmMultiRegion.write(outputName);
740
742
743 vtk::seriesWriter& series = series_(seriesName);
744
745 // First time?
746 // Load from file, verify against filesystem,
747 // prune time >= currentTime
748 if (series.empty())
749 {
750 series.load(seriesName, true, timeValue);
751 }
752
753 series.append(timeValue, outputName);
754 series.write(seriesName);
755 }
756
757 return true;
758}
759
760
762{
763 meshSubsets_.clear();
764 vtuMappings_.clear();
765 meshes_.clear();
766
767 return true;
768}
769
770
771// ************************************************************************* //
label k
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")
label filterKeys(const UnaryPredicate &pred, const bool pruning=false)
Generalized means to filter table entries based on their keys.
@ ASCII
"ascii" (normal default)
static streamFormat formatEnum(const word &formatName, const streamFormat deflt=streamFormat::ASCII)
static unsigned int defaultPrecision() noexcept
Return the default precision.
Definition: IOstream.H:342
Interpolation class within a primitive patch. Allows interpolation from points to faces and vice vers...
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
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 bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Base functionality common to reader and writer classes.
Definition: ccmBase.H:63
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A class for handling file names.
Definition: fileName.H:76
Abstract base-class for Time/database function objects.
static word outputPrefix
Directory prefix.
Virtual base class for function objects with a reference to Time.
Writes fields in VTK (xml or legacy) format. Writes cell-values or point-interpolated values for volF...
Definition: vtkWrite.H:258
virtual bool read(const dictionary &dict)
Read the vtkWrite specification.
Definition: vtkWrite.C:164
virtual bool execute()
Execute - does nothing.
Definition: vtkWrite.C:245
virtual bool write()
Write fields.
Definition: vtkWrite.C:251
virtual bool end()
On end - cleanup internal allocations.
Definition: vtkWrite.C:761
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
wordList names() const
The unsorted names of all objects.
label count(const char *clsName) const
The number of objects of the given class name.
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition: one.H:62
label index() const noexcept
The index of this patch in the boundaryMesh.
const word & name() const noexcept
The patch name.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
const word & regionName() const
The mesh region name or word::null if polyMesh::defaultRegion.
Definition: polyMesh.C:854
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
splitCell * master() const
Definition: splitCell.H:113
Interpolate from cell centres to points (vertices) using inverse distance weighting.
Write OpenFOAM patches and patch fields in VTP or legacy vtk format.
Provides a means of accumulating and generating VTK file series.
bool append(const fileNameInstant &inst)
Append the specified file instant.
bool empty() const noexcept
True if there are no data sets.
label load(const fileName &seriesName, const bool checkFiles=false, const scalar restartTime=ROOTVGREAT)
Clear contents and reload by parsing the specified file.
static void write(const fileName &base, const UList< instant > &series, const char sep='_')
Write file series (JSON format) to disk, for specified instances.
Provides a means of accumulating file entries for generating a vtkMultiBlockDataSet (....
Definition: foamVtmWriter.H:92
label write(const fileName &file)
void setTime(scalar timeValue)
Define "TimeValue" for FieldData (name as per Catalyst output)
void add(const word &blockName, const vtmWriter &other)
A class for handling words, derived from Foam::string.
Definition: word.H:68
static word printf(const char *fmt, const PrimitiveType &val)
Use a printf-style formatter for a primitive.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
mesh update()
labelList patchIds
const label nVolFields
label nPatchWriters
fileName vtmOutputBase(outputDir/regionDir/vtkName+timeDesc)
autoPtr< vtk::internalWriter > internalWriter
PtrList< vtk::patchWriter > patchWriters
const polyBoundaryMesh & patches
engineTime & runTime
Foam::word regionName(Foam::polyMesh::defaultRegion)
word outputName("finiteArea-edges.obj")
const word & regionDir
PtrList< PrimitivePatchInterpolation< primitivePatch > > patchInterps
vtk::vtmWriter vtmWriter
const label nDimFields
autoPtr< volPointInterpolation > pInterp
vtk::vtmWriter vtmBoundaries
label patchId(-1)
#define WarningInFunction
Report a warning using Foam::Warning.
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition: BitOps.H:78
const wordList volume
Standard volume field types (scalar, vector, tensor, etc)
const word fileExtension
Legacy file extension ("vtk")
@ INLINE_BASE64
XML inline base64, base64Formatter.
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:83
label writeAllVolFields(ensightCase &ensCase, const ensightMesh &ensMesh, const IOobjectList &objects, const bool nearCellValue=false)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:87
GeometricField< sphericalTensor, fvPatchField, volMesh > volSphericalTensorField
Definition: volFieldsFwd.H:85
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:86
HashSet< word, Hash< word > > wordHashSet
A HashSet of words, uses string hasher.
Definition: HashSet.H:82
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
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Operations on lists of strings.
Functor to determine if a string is exists in a list of strings.
#define vtkWrite_WRITE_FIELD(FieldType)