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