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-2019 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
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 Namespace
27  Foam::vtk::Tools
28 
29 Description
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 
37 Note
38  The class is implemented as headers-only.
39 
40 SourceFiles
41  foamVtkToolsI.H
42  foamVtkToolsTemplates.C
43 
44 \*---------------------------------------------------------------------------*/
45 
46 #ifndef foamVtkTools_H
47 #define foamVtkTools_H
48 
49 #include "faceList.H"
50 #include "pointField.H"
51 #include "symmTensor.H"
52 #include "MinMax.H"
53 
54 // VTK includes
55 #include "vtkCellArray.h"
56 #include "vtkFloatArray.h"
57 #include "vtkDoubleArray.h"
58 #include "vtkIdTypeArray.h"
59 #include "vtkSmartPointer.h"
60 #include "vtkUnsignedCharArray.h"
61 #include "vtkPoints.h"
62 #include "vtkPolyData.h"
63 
64 #include <utility>
65 
66 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67 
68 // Forward Declarations
69 class vtkDataSet;
70 class vtkCellData;
71 class vtkPointData;
72 
73 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
74 
75 namespace Foam
76 {
77 namespace vtk
78 {
79 
80 /*---------------------------------------------------------------------------*\
81  Class vtk::Caching Declaration
82 \*---------------------------------------------------------------------------*/
83 
84 //- Bookkeeping for internal caching.
85 // Retain an original copy of the geometry as well as a shallow copy
86 // with the output fields.
87 // The original copy is reused for different timesteps etc.
88 template<class DataType>
89 struct Caching
90 {
91  typedef DataType dataType;
92 
93  //- The geometry, without any cell/point data
95 
96  //- The shallow-copy of geometry, plus additional data
98 
99  //- Number of points associated with the geometry
100  inline uint64_t nPoints() const
101  {
102  return vtkgeom ? vtkgeom->GetNumberOfPoints() : 0;
103  }
104 
105  //- Clear geometry and dataset
106  void clearGeom()
107  {
108  vtkgeom = nullptr;
109  dataset = nullptr;
110  }
111 
112  //- Return a shallow copy of vtkgeom for manipulation
114  {
115  auto copy = vtkSmartPointer<dataType>::New();
116  if (vtkgeom)
117  {
118  copy->ShallowCopy(vtkgeom);
119  }
120  return copy;
121  }
122 
123  //- Make a shallow copy of vtkgeom into dataset
124  void reuse()
125  {
127  if (vtkgeom)
128  {
129  dataset->ShallowCopy(vtkgeom);
130  }
131  }
132 
133  //- Set the geometry and make a shallow copy to dataset
135  {
136  vtkgeom = geom;
137  reuse();
138  }
139 
140  //- Report basic information to output
141  void PrintSelf(std::ostream& os) const
142  {
143  os << "geom" << nl;
144  if (vtkgeom)
145  {
146  vtkgeom->PrintSelf(std::cout, vtkIndent(2));
147  }
148  else
149  {
150  os << "nullptr";
151  }
152  os << nl;
153 
154  os << "copy" << nl;
155  if (dataset)
156  {
157  dataset->PrintSelf(std::cout, vtkIndent(2));
158  }
159  else
160  {
161  os << "nullptr";
162  }
163  os << nl;
164  }
165 };
166 
167 
168 /*---------------------------------------------------------------------------*\
169  Namespace vtk::Tools
170 \*---------------------------------------------------------------------------*/
171 
172 namespace Tools
173 {
174  //- Wrap vtkUnsignedCharArray as a UList
175  inline UList<uint8_t> asUList
176  (
177  vtkUnsignedCharArray* array,
178  const label size
179  );
180 
181  //- Wrap vtkIdTypeArray as a UList
183  (
184  vtkIdTypeArray* array,
185  const label size
186  );
187 
188  //- Wrap vtkCellArray as a UList
190  (
191  vtkCellArray* cells,
192  const label nCells,
193  const label size
194  );
195 
196 
197  //- Return a list of points as vtkPoints
199  (
200  const UList<point>& pts
201  );
202 
203  //- Return an indirect list of points as vtkPoints
205  (
206  const UList<point>& pts,
207  const labelUList& addr
208  );
209 
210  //- Convert a list of faces (or triFaces) to vtk polygon cells
211  template<class Face>
213 
214  //- Return vtkPolyData of vertices for each point
216  (
217  const UList<point>& pts
218  );
219 
220  //- Return vtkPolyData of vertices for each point
222  (
223  const UList<point>& pts,
224  const labelUList& addr
225  );
226 
227  //- Min/Max of scalar, or mag() of non-scalars. Includes nullptr check.
228  inline scalarMinMax rangeOf(vtkDataArray* data);
229 
230 
231  //- Convert OpenFOAM patch to vtkPolyData
232  struct Patch
233  {
234  //- Return local patch points as vtkPoints
235  template<class PatchType>
236  static vtkSmartPointer<vtkPoints> points(const PatchType& p);
237 
238  //- Convert patch faces to vtk polygon cells
239  template<class PatchType>
240  static vtkSmartPointer<vtkCellArray> faces(const PatchType& p);
241 
242  //- Convert patch points/faces to vtkPolyData
243  template<class PatchType>
244  static vtkSmartPointer<vtkPolyData> mesh(const PatchType& p);
245 
246  //- Convert patch face normals to vtkFloatArray
247  template<class PatchType>
248  static vtkSmartPointer<vtkFloatArray> faceNormals(const PatchType& p);
249 
250  //- Return patch face centres as vtkPoints
251  template<class PatchType>
252  static vtkSmartPointer<vtkPoints> faceCentres(const PatchType& p);
253 
254  //- Convert points/faces component to vtkPolyData
255  template<class Face>
257  (
258  const UList<point>& pts,
259  const UList<Face>& fcs
260  );
261  };
262 
263 
264  //- Remapping for some OpenFOAM data types (eg, symmTensor)
265  // \param data[in,out] The data to be remapped in-place
266  template<class Type>
267  inline void remapTuple(float data[]) {}
268 
269  //- Template specialization for symmTensor ordering
270  template<>
271  inline void remapTuple<symmTensor>(float data[]);
272 
273  //- Remapping for some OpenFOAM data types (eg, symmTensor)
274  // \param data[in,out] The data to be remapped in-place
275  template<class Type>
276  inline void remapTuple(double data[]) {}
277 
278  //- Template specialization for symmTensor ordering
279  template<>
280  inline void remapTuple<symmTensor>(double data[]);
281 
282  //- Copy/transcribe OpenFOAM data types to VTK format
283  // This allows a change of data type (float vs double) as well as
284  // addressing any swapping issues (eg, symmTensor)
285  //
286  // \param output[out] The output scratch space. Must be long enough
287  // to hold the result.
288  // \param val[in] The input data to copy/transcribe
289  template<class Type>
290  inline void foamToVtkTuple(float output[], const Type& val);
291 
292  //- Copy/transcribe OpenFOAM data types to VTK format
293  // This allows a change of data type (float vs double) as well as
294  // addressing any swapping issues (eg, symmTensor)
295  //
296  // \param output[out] The output scratch space. Must be long enough
297  // to hold the result.
298  // \param val[in] The input data to copy/transcribe
299  template<class Type>
300  inline void foamToVtkTuple(double output[], const Type& val);
301 
302 
303  // Field Conversion Functions
304 
305  //- Copy list to pre-allocated vtk array.
306  // \return number of input items copied
307  template<class Type>
309  (
310  vtkFloatArray* array,
311  const UList<Type>& input,
312  vtkIdType start = 0
313  );
314 
315  //- Create named field initialized to zero
316  template<class Type>
318  (
319  const word& name,
320  const label size
321  );
322 
323  //- Convert field data to a vtkFloatArray
324  template<class Type>
326  (
327  const word& name,
328  const UList<Type>& fld
329  );
330 
331  //- An identity list of VTK_VERTEX
333  (
334  const label size
335  );
336 
337 } // End namespace Tools
338 
339 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
340 
341 } // End namespace vtk
342 } // End namespace Foam
343 
344 
345 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
346 
347 // Specializations
348 
349 //- Template specialization for symmTensor ordering
350 template<>
351 void Foam::vtk::Tools::remapTuple<Foam::symmTensor>(float data[])
352 {
353  std::swap(data[1], data[3]); // swap XY <-> YY
354  std::swap(data[2], data[5]); // swap XZ <-> ZZ
355 }
356 
357 
358 //- Template specialization for symmTensor ordering
359 template<>
360 void Foam::vtk::Tools::remapTuple<Foam::symmTensor>(double data[])
361 {
362  std::swap(data[1], data[3]); // swap XY <-> YY
363  std::swap(data[2], data[5]); // swap XZ <-> ZZ
364 }
365 
366 
367 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
368 
369 #include "foamVtkToolsI.H"
370 
371 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
372 
373 #ifdef NoRepository
374  #include "foamVtkToolsTemplates.C"
375 #endif
376 
377 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
378 
379 #endif
380 
381 // ************************************************************************* //
Foam::vtk::Caching::nPoints
uint64_t nPoints() const
Number of points associated with the geometry.
Definition: foamVtkTools.H:100
Foam::val
label ListType::const_reference val
Definition: ListOps.H:407
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::vtk::Tools::Vertices
vtkSmartPointer< vtkPolyData > Vertices(const UList< point > &pts)
Return vtkPolyData of vertices for each point.
Definition: foamVtkToolsI.H:104
Foam::vtk::Tools::Faces
vtkSmartPointer< vtkCellArray > Faces(const UList< Face > &faces)
Convert a list of faces (or triFaces) to vtk polygon cells.
Definition: foamVtkToolsTemplates.C:44
Foam::vtk::Tools::remapTuple< symmTensor >
void remapTuple< symmTensor >(float data[])
Template specialization for symmTensor ordering.
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::vtk::Tools::Patch::mesh
static vtkSmartPointer< vtkPolyData > mesh(const PatchType &p)
Convert patch points/faces to vtkPolyData.
Definition: foamVtkToolsTemplates.C:92
Foam::vtk::Tools::Patch::faceNormals
static vtkSmartPointer< vtkFloatArray > faceNormals(const PatchType &p)
Convert patch face normals to vtkFloatArray.
Definition: foamVtkToolsTemplates.C:122
Foam::vtk::Tools::Patch::faceCentres
static vtkSmartPointer< vtkPoints > faceCentres(const PatchType &p)
Return patch face centres as vtkPoints.
Definition: foamVtkToolsTemplates.C:156
Foam::vtk::Tools::zeroField
vtkSmartPointer< vtkFloatArray > zeroField(const word &name, const label size)
Create named field initialized to zero.
Definition: foamVtkToolsTemplates.C:258
Foam::vtk::Tools::asUList
UList< uint8_t > asUList(vtkUnsignedCharArray *array, const label size)
Wrap vtkUnsignedCharArray as a UList.
Definition: foamVtkToolsI.H:31
MinMax.H
foamVtkToolsTemplates.C
vtkSmartPointer< dataType >
Foam::vtk::Tools::Patch
Convert OpenFOAM patch to vtkPolyData.
Definition: foamVtkTools.H:232
faceList.H
Foam::vtk::Caching::PrintSelf
void PrintSelf(std::ostream &os) const
Report basic information to output.
Definition: foamVtkTools.H:141
Foam::vtk::Tools::transcribeFloatData
label transcribeFloatData(vtkFloatArray *array, const UList< Type > &input, vtkIdType start=0)
Copy list to pre-allocated vtk array.
Foam::vtk::Caching::clearGeom
void clearGeom()
Clear geometry and dataset.
Definition: foamVtkTools.H:106
Foam::vtk::Caching::vtkgeom
vtkSmartPointer< dataType > vtkgeom
The geometry, without any cell/point data.
Definition: foamVtkTools.H:94
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
symmTensor.H
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::vtk::Tools::identityVertices
vtkSmartPointer< vtkCellArray > identityVertices(const label size)
An identity list of VTK_VERTEX.
Definition: foamVtkToolsI.H:149
fld
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputLagrangian.H:23
Foam::vtk::Tools::convertFieldToVTK
vtkSmartPointer< vtkFloatArray > convertFieldToVTK(const word &name, const UList< Type > &fld)
Convert field data to a vtkFloatArray.
Definition: foamVtkToolsTemplates.C:286
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::vtk::Caching::getCopy
vtkSmartPointer< dataType > getCopy() const
Return a shallow copy of vtkgeom for manipulation.
Definition: foamVtkTools.H:113
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
pointField.H
foamVtkToolsI.H
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::start
label ListType::const_reference const label start
Definition: ListOps.H:408
Foam::vtk::Caching::dataType
DataType dataType
Definition: foamVtkTools.H:91
Foam::UList
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:103
Foam::vtk::Tools::Patch::points
static vtkSmartPointer< vtkPoints > points(const PatchType &p)
Return local patch points as vtkPoints.
Definition: foamVtkToolsTemplates.C:75
Foam::vtk::Caching::reuse
void reuse()
Make a shallow copy of vtkgeom into dataset.
Definition: foamVtkTools.H:124
Foam::vtk::Caching::dataset
vtkSmartPointer< dataType > dataset
The shallow-copy of geometry, plus additional data.
Definition: foamVtkTools.H:97
Foam::vtk::Tools::Points
vtkSmartPointer< vtkPoints > Points(const UList< point > &pts)
Return a list of points as vtkPoints.
Definition: foamVtkToolsI.H:70
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::vtk::Caching::set
void set(vtkSmartPointer< dataType > geom)
Set the geometry and make a shallow copy to dataset.
Definition: foamVtkTools.H:134
Foam::vtk::Tools::foamToVtkTuple
void foamToVtkTuple(float output[], const Type &val)
Copy/transcribe OpenFOAM data types to VTK format.
Definition: foamVtkToolsI.H:173
Foam::MinMax< scalar >
Foam::vtk::Tools::rangeOf
scalarMinMax rangeOf(vtkDataArray *data)
Min/Max of scalar, or mag() of non-scalars. Includes nullptr check.
Definition: foamVtkToolsI.H:127
Foam::vtk::Tools::remapTuple
void remapTuple(float data[])
Remapping for some OpenFOAM data types (eg, symmTensor)
Definition: foamVtkTools.H:267
Foam::data
Database for solution data, solver performance and other reduced data.
Definition: data.H:54
Foam::vtk::Caching
Bookkeeping for internal caching.
Definition: foamVtkTools.H:89
Foam::vtk::Tools::Patch::faces
static vtkSmartPointer< vtkCellArray > faces(const PatchType &p)
Convert patch faces to vtk polygon cells.
Definition: foamVtkToolsTemplates.C:84