Loading...
Searching...
No Matches
vtkWrite.H
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-2023 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
26Class
27 Foam::functionObjects::vtkWrite
28
29Group
30 grpUtilitiesFunctionObjects
31
32Description
33 Writes fields in VTK (xml or legacy) format.
34 Writes cell-values or point-interpolated values for volFields.
35
36 Example of function object specification:
37 \verbatim
38 vtkWrite1
39 {
40 type vtkWrite;
41 libs (utilityFunctionObjects);
42 writeControl writeTime;
43 writeInterval 1;
44 format binary;
45 legacy false;
46 decompose false;
47 ...
48 fields (U p);
49 // excludeFields ("force.*");
50
51 selection
52 {
53 box
54 {
55 action use;
56 source box;
57 box (-0.1 -0.01 -0.1) (0.1 0.30 0.1);
58 }
59 dome
60 {
61 action add;
62 source sphere;
63 origin (-0.1 -0.01 -0.1);
64 radius 0.25;
65 }
66 centre
67 {
68 action subtract;
69 source sphere;
70 origin (-0.1 -0.01 -0.1);
71 radius 0.1;
72 }
73 blob
74 {
75 action add;
76 source surface;
77 surface triSurfaceMesh;
78 name blob.stl;
79 }
80 }
81 }
82 \endverbatim
83
84 \heading Basic Usage
85 \table
86 Property | Description | Required | Default
87 type | Type name: vtkWrite | yes |
88 fields | Select fields to output (wordRe list) | yes |
89 excludeFields | Exclude fields from output (wordRe list) | no |
90 boundary | Convert boundary fields | no | true
91 internal | Convert internal fields | no | true
92 single | Combine patches into a single boundary | no | false
93 interpolate | Interpolate for point values | no | false
94 \endtable
95
96 \heading Output Options
97 \table
98 Property | Description | Required | Default
99 format | ascii or binary format | no | binary
100 legacy | Legacy VTK output | no | false
101 precision | Write precision in ascii | no | same as IOstream
102 directory | The output directory name | no | postProcessing/NAME
103 width | Padding width for file name | no | 8
104 decompose | Decompose polyhedral cells | no | false
105 writeIds | Write cell,patch,proc id fields | no | false
106 \endtable
107
108 \heading Output Selection
109 \table
110 Property | Description | Required | Default
111 region | Name for a single region | no | region0
112 regions | List of regions (wordRe list) | no |
113 patches | Limit to listed patches (wordRe list) | no |
114 excludePatches | Exclude specified patches | no |
115 selection | Cell selection (topoSet actions) | no | empty dict
116 \endtable
117
118Note
119 The region of interest is defined by the selection dictionary
120 as a set of actions (use,add,subtract,subset,invert).
121 Omitting the selection dictionary is the same as specifying the
122 conversion of all cells (in the selected regions).
123 Omitting the patches entry is the same as specifying the conversion of all
124 patches.
125
126See also
127 Foam::functionObjects::ensightWrite
128 Foam::functionObjects::fvMeshFunctionObject
129 Foam::functionObjects::timeControl
130 Foam::cellBitSet::select
131
132SourceFiles
133 vtkWrite.cxx
134 vtkWriteImpl.cxx
135 vtkWriteUpdate.cxx
136
137\*---------------------------------------------------------------------------*/
138
139#ifndef Foam_functionObjects_vtkWrite_H
140#define Foam_functionObjects_vtkWrite_H
141
142#include "timeFunctionObject.H"
144#include "foamVtkPatchWriter.H"
145#include "foamVtkSeriesWriter.H"
146#include "fvMeshSubsetProxy.H"
147#include "searchableSurfaces.H"
148
149// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
150
151namespace Foam
152{
153namespace functionObjects
154{
155
156/*---------------------------------------------------------------------------*\
157 Class vtkWrite Declaration
158\*---------------------------------------------------------------------------*/
159
160class vtkWrite
161:
162 public functionObjects::timeFunctionObject
163{
164 // Private Data
165
166 //- The output directory
167 fileName outputDir_;
168
169 //- The printf format for zero-padding names
170 string printf_;
171
172 //- VTK output options
173 vtk::outputOptions writeOpts_;
174
175 //- Verbose output
176 bool verbose_;
177
178 //- Convert internal mesh
179 bool doInternal_;
180
181 //- Convert boundary mesh
182 bool doBoundary_;
183
184 //- Combine patches into a single file
185 bool oneBoundary_;
186
187 //- Interpolate cell to point values
188 bool interpolate_;
189
190 //- Decompose polyhedra
191 bool decompose_;
192
193 //- Write cell ids field
194 bool writeIds_;
195
196 //- Track changes in mesh geometry
197 polyMesh::readUpdateState meshState_;
198
199 //- Requested names of regions to process
200 wordRes selectRegions_;
201
202 //- Requested names of patches to process
203 wordRes selectPatches_;
204
205 //- Selection of patches to block
206 wordRes blockPatches_;
207
208 //- Requested selection of fields to process
209 wordRes selectFields_;
210
211 //- Selection of fields to block
212 wordRes blockFields_;
213
214 //- Dictionary of volume selections
215 dictionary selection_;
216
217 //- Pointers to the requested mesh regions (sorted)
218 UPtrList<const fvMesh> meshes_;
219
220 //- Subsetting for meshes.
221 // Access index according to sorted mesh names.
222 PtrList<fvMeshSubset> meshSubsets_;
223
224 //- Storage for VTU cells, sizing.
225 // Access index according to sorted mesh names.
226 PtrList<vtk::vtuCells> vtuMappings_;
227
228 //- VTK file series
229 HashTable<vtk::seriesWriter, fileName> series_;
230
231
232 // Private Member Functions
233
234 //- Update mesh subset according to zones, geometry, bounds
235 bool updateSubset(fvMeshSubset& subsetter) const;
236
237 //- Get patchIds selected in list
238 labelList getSelectedPatches(const polyBoundaryMesh& patches) const;
239
240 //- Read information for selections
241 bool readSelection(const dictionary& dict);
242
243 //- Update meshes, subMeshes etc.
244 bool update();
245
246
247 // Write
248
249 //- Write all volume fields
250 label writeAllVolFields
251 (
252 autoPtr<vtk::internalWriter>& internalWriter,
253 UPtrList<vtk::patchWriter>& patchWriters,
254 const fvMeshSubset& proxy,
255 const wordHashSet& candidateNames
256 ) const;
257
258 //- Write all volume fields with point interpolation
259 label writeAllVolFields
260 (
261 autoPtr<vtk::internalWriter>& internalWriter,
262 const autoPtr<volPointInterpolation>& pInterp,
263 UPtrList<vtk::patchWriter>& patchWriters,
264 const UPtrList
267 >& patchInterps,
268 const fvMeshSubset& proxy,
269 const wordHashSet& candidateNames
270 ) const;
271
272 //- Write selected GeoField fields.
273 template<class GeoField>
274 label writeVolFieldsImpl
275 (
278 const fvMeshSubset& proxy,
279 const wordHashSet& candidateNames
280 ) const;
281
282 //- Write selected GeoField fields with point interpolation
283 template<class GeoField>
284 label writeVolFieldsImpl
285 (
289 const UPtrList
290 <
292 >& patchInterps,
293 const fvMeshSubset& proxy,
294 const wordHashSet& candidateNames
295 ) const;
296
297
298public:
299
300 //- Runtime type information
301 TypeName("vtkWrite");
302
303
304 // Constructors
305
306 //- Construct from name, Time and dictionary
308 (
309 const word& name,
310 const Time& runTime,
311 const dictionary& dict
312 );
313
314 //- No copy construct
315 vtkWrite(const vtkWrite&) = delete;
316
317 //- No copy assignment
318 void operator=(const vtkWrite&) = delete;
319
320
321 //- Destructor
322 virtual ~vtkWrite() = default;
323
324
325 // Member Functions
326
327 //- Read the vtkWrite specification
328 virtual bool read(const dictionary& dict);
329
330 //- Execute - does nothing
331 virtual bool execute();
332
333 //- Write fields
334 virtual bool write();
335
336 //- On end - cleanup internal allocations
337 virtual bool end();
338
339 //- Update for changes of mesh
340 virtual void updateMesh(const mapPolyMesh& mpm);
341
342 //- Update for mesh point-motion
343 virtual void movePoints(const polyMesh& mesh);
344};
345
346
347// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
348
349} // End namespace functionObjects
350} // End namespace Foam
351
352// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
353
354#endif
355
356// ************************************************************************* //
A HashTable similar to std::unordered_map.
Definition HashTable.H:124
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:67
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition UPtrList.H:101
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A class for handling file names.
Definition fileName.H:75
const word & name() const noexcept
Return the name of this functionObject.
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:268
virtual bool read(const dictionary &dict)
Read the vtkWrite specification.
virtual bool end()
On end - cleanup internal allocations.
virtual bool execute()
Execute - does nothing.
void operator=(const vtkWrite &)=delete
No copy assignment.
virtual bool write()
Write fields.
vtkWrite(const vtkWrite &)=delete
No copy construct.
virtual ~vtkWrite()=default
Destructor.
TypeName("vtkWrite")
Runtime type information.
vtkWrite(const word &name, const Time &runTime, const dictionary &dict)
Construct from name, Time and dictionary.
virtual void movePoints(const polyMesh &mesh)
Update for mesh point-motion.
virtual void updateMesh(const mapPolyMesh &mpm)
Update for changes of mesh.
Holds a reference to the original mesh (the baseMesh) and optionally to a subset of that mesh (the su...
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition polyMesh.H:92
Encapsulated combinations of output format options. This is primarily useful when defining the output...
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
mesh update()
autoPtr< vtk::internalWriter > internalWriter
PtrList< vtk::patchWriter > patchWriters
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
engineTime & runTime
PtrList< PrimitivePatchInterpolation< primitivePatch > > patchInterps
autoPtr< volPointInterpolation > pInterp
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
Namespace for OpenFOAM.
HashSet< word, Hash< word > > wordHashSet
A HashSet of words, uses string hasher.
Definition HashSet.H:80
List< label > labelList
A List of labels.
Definition List.H:62
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68