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-2020 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 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 
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 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65 
66 // Forward Declarations
67 class vtkDataSet;
68 class vtkCellData;
69 class vtkPointData;
70 
71 namespace Foam
72 {
73 namespace vtk
74 {
75 
76 /*---------------------------------------------------------------------------*\
77  Class vtk::Caching Declaration
78 \*---------------------------------------------------------------------------*/
79 
80 //- Bookkeeping for internal caching.
81 // Retain an original copy of the geometry as well as a shallow copy
82 // with the output fields.
83 // The original copy is reused for different timesteps etc.
84 template<class DataType>
85 struct Caching
86 {
87  typedef DataType dataType;
88 
89  //- The geometry, without any cell/point data
90  vtkSmartPointer<dataType> vtkgeom;
91 
92  //- The shallow-copy of geometry, plus additional data
93  vtkSmartPointer<dataType> dataset;
94 
95  //- Number of points associated with the geometry
96  inline uint64_t nPoints() const
97  {
98  return vtkgeom ? vtkgeom->GetNumberOfPoints() : 0;
99  }
100 
101  //- Clear geometry and dataset
102  void clearGeom()
103  {
104  vtkgeom = nullptr;
105  dataset = nullptr;
106  }
107 
108  //- Return a shallow copy of vtkgeom for manipulation
109  vtkSmartPointer<dataType> getCopy() const
110  {
111  auto copy = vtkSmartPointer<dataType>::New();
112  if (vtkgeom)
113  {
114  copy->ShallowCopy(vtkgeom);
115  }
116  return copy;
117  }
119  //- Make a shallow copy of vtkgeom into dataset
120  void reuse()
121  {
123  if (vtkgeom)
124  {
125  dataset->ShallowCopy(vtkgeom);
126  }
127  }
128 
129  //- Set the geometry and make a shallow copy to dataset
130  void set(vtkSmartPointer<dataType> geom)
131  {
132  vtkgeom = geom;
133  reuse();
134  }
135 
136  //- Report basic information to output
137  void PrintSelf(std::ostream& os) const
138  {
139  os << "geom" << nl;
140  if (vtkgeom)
141  {
142  vtkgeom->PrintSelf(std::cout, vtkIndent(2));
143  }
144  else
145  {
146  os << "nullptr";
147  }
148  os << nl;
149 
150  os << "copy" << nl;
151  if (dataset)
152  {
153  dataset->PrintSelf(std::cout, vtkIndent(2));
154  }
155  else
156  {
157  os << "nullptr";
158  }
159  os << nl;
160  }
161 };
162 
163 
164 /*---------------------------------------------------------------------------*\
165  Namespace vtk::Tools
166 \*---------------------------------------------------------------------------*/
167 
168 namespace Tools
169 {
170  //- Wrap vtkUnsignedCharArray as a UList
171  inline UList<uint8_t> asUList
172  (
173  vtkUnsignedCharArray* array,
174  const label size
175  );
176 
177  //- Wrap vtkIdTypeArray as a UList
178  inline UList<vtkIdType> asUList
179  (
180  vtkIdTypeArray* array,
181  const label size
182  );
183 
184  //- Return a list of points as vtkPoints
185  inline vtkSmartPointer<vtkPoints> Points
186  (
187  const UList<point>& pts
188  );
189 
190  //- Return an indirect list of points as vtkPoints
191  inline vtkSmartPointer<vtkPoints> Points
192  (
193  const UList<point>& pts,
194  const labelUList& addr
195  );
196 
197  //- Convert a list of faces (or triFaces) to vtk polygon cells
198  template<class Face>
199  vtkSmartPointer<vtkCellArray> Faces(const UList<Face>& faces);
200 
201  //- Return vtkPolyData of vertices for each point
202  inline vtkSmartPointer<vtkPolyData> Vertices
203  (
204  const UList<point>& pts
205  );
206 
207  //- Return vtkPolyData of vertices for each point
208  inline vtkSmartPointer<vtkPolyData> Vertices
209  (
210  const UList<point>& pts,
211  const labelUList& addr
212  );
213 
214  //- Min/Max of scalar, or mag() of non-scalars. Includes nullptr check.
215  inline scalarMinMax rangeOf(vtkDataArray* data);
216 
217 
218  //- Convert OpenFOAM patch to vtkPolyData
219  struct Patch
220  {
221  //- Return local patch points as vtkPoints
222  template<class PatchType>
223  static vtkSmartPointer<vtkPoints> points(const PatchType& p);
224 
225  //- Convert patch faces to vtk polygon cells
226  template<class PatchType>
227  static vtkSmartPointer<vtkCellArray> faces(const PatchType& p);
228 
229  //- Convert patch points/faces to vtkPolyData
230  template<class PatchType>
231  static vtkSmartPointer<vtkPolyData> mesh(const PatchType& p);
232 
233  //- Convert patch face normals to vtkFloatArray
234  template<class PatchType>
235  static vtkSmartPointer<vtkFloatArray> faceNormals(const PatchType& p);
236 
237  //- Return patch face centres as vtkPoints
238  template<class PatchType>
239  static vtkSmartPointer<vtkPoints> faceCentres(const PatchType& p);
240 
241  //- Convert points/faces component to vtkPolyData
242  template<class Face>
243  static vtkSmartPointer<vtkPolyData> mesh
244  (
245  const UList<point>& pts,
246  const UList<Face>& fcs
247  );
248  };
249 
250 
251  //- Remapping for some OpenFOAM data types (eg, symmTensor)
252  // \param data[in,out] The data to be remapped in-place
253  template<class Type>
254  inline void remapTuple(float data[]) {}
255 
256  //- Template specialization for symmTensor ordering
257  template<>
258  inline void remapTuple<symmTensor>(float data[]);
259 
260  //- Remapping for some OpenFOAM data types (eg, symmTensor)
261  // \param data[in,out] The data to be remapped in-place
262  template<class Type>
263  inline void remapTuple(double data[]) {}
264 
265  //- Template specialization for symmTensor ordering
266  template<>
267  inline void remapTuple<symmTensor>(double data[]);
268 
269  //- Copy/transcribe OpenFOAM data types to VTK format
270  // This allows a change of data type (float vs double) as well as
271  // addressing any swapping issues (eg, symmTensor)
272  //
273  // \param output[out] The output scratch space. Must be long enough
274  // to hold the result.
275  // \param val[in] The input data to copy/transcribe
276  template<class Type>
277  inline void foamToVtkTuple(float output[], const Type& val);
278 
279  //- Copy/transcribe OpenFOAM data types to VTK format
280  // This allows a change of data type (float vs double) as well as
281  // addressing any swapping issues (eg, symmTensor)
282  //
283  // \param output[out] The output scratch space. Must be long enough
284  // to hold the result.
285  // \param val[in] The input data to copy/transcribe
286  template<class Type>
287  inline void foamToVtkTuple(double output[], const Type& val);
288 
289 
290  // Field Conversion Functions
291 
292  //- Copy list to pre-allocated vtk array.
293  // \return number of input items copied
294  template<class Type>
295  label transcribeFloatData
296  (
297  vtkFloatArray* array,
298  const UList<Type>& input,
299  vtkIdType start = 0
300  );
301 
302  //- Create named field initialized to zero
303  template<class Type>
304  vtkSmartPointer<vtkFloatArray> zeroField
305  (
306  const word& name,
307  const label size
308  );
309 
310  //- Convert field data to a vtkFloatArray
311  template<class Type>
312  vtkSmartPointer<vtkFloatArray> convertFieldToVTK
313  (
314  const word& name,
315  const UList<Type>& fld
316  );
317 
318  //- An identity list of VTK_VERTEX
319  inline vtkSmartPointer<vtkCellArray> identityVertices
320  (
321  const label size
322  );
323 
324 } // End namespace Tools
325 
326 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
327 
328 } // End namespace vtk
329 } // End namespace Foam
330 
331 
332 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
333 
334 // Specializations
335 
336 //- Template specialization for symmTensor ordering
337 template<>
338 void Foam::vtk::Tools::remapTuple<Foam::symmTensor>(float data[])
339 {
340  std::swap(data[1], data[3]); // swap XY <-> YY
341  std::swap(data[2], data[5]); // swap XZ <-> ZZ
342 }
343 
344 
345 //- Template specialization for symmTensor ordering
346 template<>
347 void Foam::vtk::Tools::remapTuple<Foam::symmTensor>(double data[])
348 {
349  std::swap(data[1], data[3]); // swap XY <-> YY
350  std::swap(data[2], data[5]); // swap XZ <-> ZZ
351 }
352 
353 
354 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
355 
356 #include "foamVtkToolsI.H"
357 
358 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
359 
360 #ifdef NoRepository
361  #include "foamVtkToolsTemplates.C"
362 #endif
363 
364 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
365 
366 #endif
367 
368 // ************************************************************************* //
vtkSmartPointer< vtkCellArray > identityVertices(const label size)
An identity list of VTK_VERTEX.
static vtkSmartPointer< vtkPoints > points(const PatchType &p)
Return local patch points as vtkPoints.
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition: MinMax.H:97
static vtkSmartPointer< vtkCellArray > faces(const PatchType &p)
Convert patch faces to vtk polygon cells.
vtkSmartPointer< vtkPoints > Points(const UList< point > &pts)
Return a list of points as vtkPoints.
Definition: foamVtkToolsI.H:52
vtkSmartPointer< vtkCellArray > Faces(const UList< Face > &faces)
Convert a list of faces (or triFaces) to vtk polygon cells.
uint64_t nPoints() const
Number of points associated with the geometry.
Definition: foamVtkTools.H:101
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
vtkSmartPointer< vtkFloatArray > zeroField(const word &name, const label size)
Create named field initialized to zero.
UList< uint8_t > asUList(vtkUnsignedCharArray *array, const label size)
Wrap vtkUnsignedCharArray as a UList.
Definition: foamVtkToolsI.H:26
static vtkSmartPointer< vtkFloatArray > faceNormals(const PatchType &p)
Convert patch face normals to vtkFloatArray.
void clearGeom()
Clear geometry and dataset.
Definition: foamVtkTools.H:109
void PrintSelf(std::ostream &os) const
Report basic information to output.
Definition: foamVtkTools.H:152
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
vtkSmartPointer< dataType > vtkgeom
The geometry, without any cell/point data.
Definition: foamVtkTools.H:91
vtkSmartPointer< vtkPolyData > Vertices(const UList< point > &pts)
Return vtkPolyData of vertices for each point.
Definition: foamVtkToolsI.H:86
label transcribeFloatData(vtkFloatArray *array, const UList< Type > &input, vtkIdType start=0)
Copy list to pre-allocated vtk array.
Bookkeeping for internal caching.
Definition: foamVtkTools.H:84
void remapTuple< symmTensor >(float data[])
Template specialization for symmTensor ordering.
static vtkSmartPointer< vtkPolyData > mesh(const PatchType &p)
Convert patch points/faces to vtkPolyData.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
A class for handling words, derived from Foam::string.
Definition: word.H:63
static Istream & input(Istream &is, IntRange< T > &range)
Definition: IntRanges.C:33
static vtkSmartPointer< vtkPoints > faceCentres(const PatchType &p)
Return patch face centres as vtkPoints.
void remapTuple(float data[])
Remapping for some OpenFOAM data types (eg, symmTensor)
Definition: foamVtkTools.H:302
OBJstream os(runTime.globalPath()/outputName)
scalarMinMax rangeOf(vtkDataArray *data)
Min/Max of scalar, or mag() of non-scalars. Includes nullptr check.
vtkSmartPointer< vtkFloatArray > convertFieldToVTK(const word &name, const UList< Type > &fld)
Convert field data to a vtkFloatArray.
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))
vtkSmartPointer< dataType > dataset
The shallow-copy of geometry, plus additional data.
Definition: foamVtkTools.H:96
vtkSmartPointer< dataType > getCopy() const
Return a shallow copy of vtkgeom for manipulation.
Definition: foamVtkTools.H:118
void foamToVtkTuple(float output[], const Type &val)
Copy/transcribe OpenFOAM data types to VTK format.
static Ostream & output(Ostream &os, const IntRange< T > &range)
Definition: IntRanges.C:44
volScalarField & p
void reuse()
Make a shallow copy of vtkgeom into dataset.
Definition: foamVtkTools.H:131
Namespace for OpenFOAM.
const pointField & pts