Loading...
Searching...
No Matches
foamVtkTools.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-2025 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
26Namespace
27 Foam::vtk::Tools
28
29Description
30 A collection of static methods to assist converting OpenFOAM data
31 structures into VTK internal data structures.
32
33 Remapping of the symmTensor order is required in input or output
34 directions. OpenFOAM uses (XX, XY, XZ, YY, YZ, ZZ) order,
35 VTK uses (XX, YY, ZZ, XY, YZ, XZ) order.
36
37Note
38 This file is implemented as headers-only.
39
40SourceFiles
41 foamVtkToolsI.H
42 foamVtkTools.txx
43
44\*---------------------------------------------------------------------------*/
45
46#ifndef Foam_vtk_vtkTools_H
47#define Foam_vtk_vtkTools_H
48
49#include "faceList.H"
50#include "pointField.H"
51#include "symmTensor.H"
52#include "MinMax.H"
53#include "foamVtkCore.H"
54
55// VTK includes
56#include "vtkCellArray.h"
57#include "vtkFloatArray.h"
58#include "vtkDoubleArray.h"
59#include "vtkIdTypeArray.h"
60#include "vtkSmartPointer.h"
61#include "vtkUnsignedCharArray.h"
62#include "vtkPoints.h"
63#include "vtkPolyData.h"
64
65// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66
67// Forward Declarations
68class vtkDataSet;
69class vtkCellData;
70class vtkPointData;
71
72namespace Foam
73{
74namespace vtk
75{
76
77/*---------------------------------------------------------------------------*\
78 Class vtk::Caching Declaration
79\*---------------------------------------------------------------------------*/
80
81//- Bookkeeping for internal caching.
82// Retain an original copy of the geometry as well as a shallow copy
83// with the output fields.
84// The original copy is reused for different timesteps etc.
85template<class DataType>
86struct Caching
87{
88 //- The VTK geometry data type
89 typedef DataType dataType;
91 //- The geometry, without any cell/point data
92 vtkSmartPointer<dataType> vtkgeom;
93
94 //- The shallow-copy of geometry, plus additional data
95 vtkSmartPointer<dataType> dataset;
96
97 //- Number of points associated with the geometry
98 inline uint64_t nPoints() const
99 {
100 return vtkgeom ? vtkgeom->GetNumberOfPoints() : 0;
101 }
102
103 //- Clear geometry and dataset
104 void clearGeom()
106 vtkgeom = nullptr;
107 dataset = nullptr;
108 }
109
110 //- Return a shallow copy of vtkgeom for manipulation
111 vtkSmartPointer<dataType> getCopy() const
112 {
113 auto copy = vtkSmartPointer<dataType>::New();
114 if (vtkgeom)
115 {
116 copy->ShallowCopy(vtkgeom);
117 }
118 return copy;
119 }
120
121 //- Make a shallow copy of vtkgeom into dataset
122 void reuse()
123 {
124 dataset = vtkSmartPointer<dataType>::New();
125 if (vtkgeom)
126 {
127 dataset->ShallowCopy(vtkgeom);
128 }
129 }
130
131 //- Set the geometry and make a shallow copy to dataset
132 void set(vtkSmartPointer<dataType> geom)
133 {
134 vtkgeom = geom;
136 }
137
138 //- Report basic information to output
139 void PrintSelf(std::ostream& os) const
140 {
141 os << "geom" << nl;
142 if (vtkgeom)
143 {
144 vtkgeom->PrintSelf(std::cout, vtkIndent(2));
145 }
146 else
148 os << "nullptr";
149 }
150 os << nl;
151
152 os << "copy" << nl;
153 if (dataset)
154 {
155 dataset->PrintSelf(std::cout, vtkIndent(2));
157 else
158 {
159 os << "nullptr";
160 }
161 os << nl;
162 }
163};
164
165
166/*---------------------------------------------------------------------------*\
167 Namespace vtk::Tools
168\*---------------------------------------------------------------------------*/
169
170namespace Tools
171{
172 //- Wrap vtkUnsignedCharArray as a UList
173 inline UList<uint8_t> asUList
174 (
175 vtkUnsignedCharArray* array,
176 const label size
177 );
178
179 //- Wrap vtkIdTypeArray as a UList
180 inline UList<vtkIdType> asUList
181 (
182 vtkIdTypeArray* array,
183 const label size
184 );
185
186 //- Allocate vtkPoints
187 inline vtkSmartPointer<vtkPoints>
188 NewPoints(vtkIdType count);
189
190 //- Return a list of points as vtkPoints
191 inline vtkSmartPointer<vtkPoints>
192 Points(const UList<point>& pts);
193
194 //- Return an indirect list of points as vtkPoints
195 inline vtkSmartPointer<vtkPoints>
196 Points
197 (
198 const UList<point>& pts,
199 const labelUList& addr
200 );
201
202 //- Convert a list of faces (or triFaces) to vtk polygon cells
203 template<class Face>
204 vtkSmartPointer<vtkCellArray>
205 Faces(const UList<Face>& faces);
206
207 //- Return vtkPolyData of vertices for each point
208 inline vtkSmartPointer<vtkPolyData>
210 (
211 const UList<point>& pts
212 );
213
214 //- Return vtkPolyData of vertices for each point
215 inline vtkSmartPointer<vtkPolyData>
217 (
218 const UList<point>& pts,
219 const labelUList& addr
220 );
221
222 //- Min/Max of scalar, or mag() of non-scalars. Includes nullptr check.
223 inline scalarMinMax rangeOf(vtkDataArray* data);
224
225
226 //- Convert OpenFOAM patch to vtkPolyData
227 struct Patch
228 {
229 //- Return local patch points as vtkPoints
230 template<class PatchType>
231 static vtkSmartPointer<vtkPoints>
232 points(const PatchType& p);
233
234 //- Convert patch faces to vtk polygon cells
235 template<class PatchType>
236 static vtkSmartPointer<vtkCellArray>
237 faces(const PatchType& p);
238
239 //- Convert patch points/faces to vtkPolyData
240 template<class PatchType>
241 static vtkSmartPointer<vtkPolyData>
242 mesh(const PatchType& p);
243
244 //- Convert patch face normals to vtkFloatArray
245 template<class PatchType>
246 static vtkSmartPointer<vtkFloatArray>
247 faceNormals(const PatchType& p);
248
249 //- Return patch face centres as vtkPoints
250 template<class PatchType>
251 static vtkSmartPointer<vtkPoints>
252 faceCentres(const PatchType& p);
253
254 //- Convert points/faces component to vtkPolyData
255 template<class Face>
256 static vtkSmartPointer<vtkPolyData>
257 mesh
258 (
259 const UList<point>& pts,
260 const UList<Face>& fcs
261 );
262 };
263
265 //- Older name for copyTuple
266 template<class Type>
267 void foamToVtkTuple(float output[], const Type& value)
268 {
269 vtk::Tools::copyTuple(output, value);
270 }
272 //- Older name for copyTuple
273 template<class Type>
274 void foamToVtkTuple(double output[], const Type& value)
275 {
276 vtk::Tools::copyTuple(output, value);
277 }
279
280 // Field Conversion Functions
281
282 //- Copy list to pre-allocated vtk array.
283 // \return number of input items copied
284 template<class Type>
286 (
289 vtkDataArray* array,
290 const UList<Type>& input,
291 vtkIdType start = 0
292 );
293
294 //- Create named field initialized to zero
295 template<class Type, class DataArrayType = vtkFloatArray>
296 vtkSmartPointer<DataArrayType>
298 (
299 const word& name,
300 const label size
301 );
302
303 //- Convert field data to a vtkFloatArray
304 template<class Type, class DataArrayType = vtkFloatArray>
305 vtkSmartPointer<DataArrayType>
307 (
308 const word& name,
309 const UList<Type>& fld
310 );
311
312 //- Convert field data to a vtkFloatArray, using indirect addressing
313 template<class Type, class DataArrayType = vtkFloatArray>
314 vtkSmartPointer<DataArrayType>
316 (
317 const word& name,
319 const labelUList& addr
320 );
321
322 //- An identity list of VTK_VERTEX
323 inline vtkSmartPointer<vtkCellArray>
324 identityVertices(const vtkIdType size);
325
326} // End namespace Tools
328// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
329
330} // End namespace vtk
331} // End namespace Foam
332
333// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
334
335#include "foamVtkToolsI.H"
336
337// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
338
339#ifdef NoRepository
340 #include "foamVtkTools.txx"
341#endif
342
343// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
344
345#endif
346
347// ************************************************************************* //
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
A collection of static methods to assist converting OpenFOAM data structures into VTK internal data s...
label transcribeFloatData(vtkDataArray *array, const UList< Type > &input, vtkIdType start=0)
Copy list to pre-allocated vtk array.
vtkSmartPointer< DataArrayType > convertFieldToVTK(const word &name, const UList< Type > &fld)
Convert field data to a vtkFloatArray.
UList< uint8_t > asUList(vtkUnsignedCharArray *array, const label size)
Wrap vtkUnsignedCharArray as a UList.
scalarMinMax rangeOf(vtkDataArray *data)
Min/Max of scalar, or mag() of non-scalars. Includes nullptr check.
vtkSmartPointer< vtkPoints > Points(const UList< point > &pts)
Return a list of points as vtkPoints.
void foamToVtkTuple(float output[], const Type &value)
Older name for copyTuple.
vtkSmartPointer< vtkCellArray > identityVertices(const vtkIdType size)
An identity list of VTK_VERTEX.
vtkSmartPointer< vtkCellArray > Faces(const UList< Face > &faces)
Convert a list of faces (or triFaces) to vtk polygon cells.
vtkSmartPointer< vtkPoints > NewPoints(vtkIdType count)
Allocate vtkPoints.
vtkSmartPointer< DataArrayType > zeroField(const word &name, const label size)
Create named field initialized to zero.
vtkSmartPointer< vtkPolyData > Vertices(const UList< point > &pts)
Return vtkPolyData of vertices for each point.
Type copyTuple(const Type &value)
Return a copy of the value, with components changed from OpenFOAM order to VTK order (for symmTensor)...
Namespace for handling VTK output. Contains classes and functions for writing VTK file content.
Namespace for OpenFOAM.
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition MinMax.H:97
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
UList< label > labelUList
A UList of labels.
Definition UList.H:75
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
const pointField & pts
Bookkeeping for internal caching.
void set(vtkSmartPointer< dataType > geom)
Set the geometry and make a shallow copy to dataset.
void PrintSelf(std::ostream &os) const
Report basic information to output.
void clearGeom()
Clear geometry and dataset.
vtkSmartPointer< dataType > vtkgeom
The geometry, without any cell/point data.
DataType dataType
The VTK geometry data type.
void reuse()
Make a shallow copy of vtkgeom into dataset.
vtkSmartPointer< dataType > dataset
The shallow-copy of geometry, plus additional data.
uint64_t nPoints() const
Number of points associated with the geometry.
vtkSmartPointer< dataType > getCopy() const
Return a shallow copy of vtkgeom for manipulation.
Convert OpenFOAM patch to vtkPolyData.
static vtkSmartPointer< vtkCellArray > faces(const PatchType &p)
Convert patch faces to vtk polygon cells.
static vtkSmartPointer< vtkPoints > points(const PatchType &p)
Return local patch points as vtkPoints.
static vtkSmartPointer< vtkFloatArray > faceNormals(const PatchType &p)
Convert patch face normals to vtkFloatArray.
static vtkSmartPointer< vtkPolyData > mesh(const UList< point > &pts, const UList< Face > &fcs)
Convert points/faces component to vtkPolyData.
static vtkSmartPointer< vtkPolyData > mesh(const PatchType &p)
Convert patch points/faces to vtkPolyData.
static vtkSmartPointer< vtkPoints > faceCentres(const PatchType &p)
Return patch face centres as vtkPoints.