foamVtuSizing.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) 2016-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 Class
27  Foam::vtk::vtuSizing
28 
29 Description
30  Sizing descriptions and routines for transcribing an OpenFOAM volume mesh
31  into a VTK unstructured grid, with possible decomposition of polyhedral
32  cells into primitive cell types.
33 
34  This class is intended to populate externally allocated arrays with content
35  that is compatible with what VTK expects. This approach allows an improved
36  separation of the OpenFOAM mesh description and the storage, and allows
37  support of alternative storage containers (eg, std::vector, vtkDataArray).
38  The ideal goal would be a zero-copy mechanism, but this does not work for
39  several reasons:
40  \par
41  - OpenFOAM and VTK have different point ordering for prism
42  - polyhedral decomposition
43  - face-stream are required for VTK
44  - VTK internal storage includes list size as part of the data
45  - VTK includes storage may be a different base size (eg, long long)
46  compared to the OpenFOAM label.
47 
48  \par Data Entries (slots)
49 
50  These are the storage entries normally associate with each output-type:
51  \table
52  legacy output
53  \c types | vtk cell type (1-255)
54  \c cells | nLabels and unique vertex labels used by the cell, or
55  | [nLabels nFaces, nFace0Pts, id1, id2, ..., nFace1Pts, id1, id2, ...]
56  \endtable
57 
58  \table
59  xml output
60  \c types | vtk cell type (1-255)
61  \c connectivity | unique vertex labels used by the cell
62  \c offsets | end offset for each of \c connectivity
63  \c faces | face stream for polyhedral cells
64  | [nFaces, nFace0Pts, id1, id2, ..., nFace1Pts, id1, id2, ...]
65  \c faceoffsets | end offset for each of \c faces, with -1 for primitive cells
66  \endtable
67 
68  \table
69  internal1 storage (VTK-8 and earlier)
70  \c types | vtk cell type (1-255)
71  \c connectivity | nLabels and unique vertex labels used by the cell
72  \c location | begin location for each of \c connectivity
73  \c faces | face stream for polyhedral cells
74  | [nFaces, nFace0Pts, id1, id2, ..., nFace1Pts, id1, id2, ...]
75  \c facelocation | begin location for each of \c faces, with -1 for primitive cells
76  \endtable
77 
78  \table
79  internal2 storage (with VTK_CELL_ARRAY_V2)
80  \c types | vtk cell type (1-255)
81  \c connectivity | unique vertex labels used by the cell
82  \c offsets | begin/end offsets for \c connectivity
83  \c faces | face stream for polyhedral cells
84  | [nFaces, nFace0Pts, id1, id2, ..., nFace1Pts, id1, id2, ...]
85  \c facelocation | begin location for each of \c faces, with -1 for primitive cells
86  \endtable
87 
88  The VTK storage concept for "connectivity" and "faces" somewhat resemble
89  a CompactListList.
90 
91 Note
92  It is possible to specify a global point offset (via the globalIndex)
93  so that the cell point labels will use global numbering.
94  There is no support for point renumbering with merged mesh points,
95  since it likely more efficient to use VTK point-blanking to mark duplicate
96  points instead of merging points ourselves.
97 
98 Note
99  The VTK_CELL_ARRAY_V2 define (from vtkCellArray.h) indicates if the new
100  (internal2) new format is being used.
101 
102 SourceFiles
103  foamVtuSizing.C
104  foamVtuSizingI.H
105 
106 \*---------------------------------------------------------------------------*/
107 
108 #ifndef foamVtuSizing_H
109 #define foamVtuSizing_H
110 
111 #include "label.H"
112 #include "labelList.H"
113 #include "foamVtkMeshMaps.H"
114 
115 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
116 
117 namespace Foam
118 {
119 
120 // Forward Declarations
121 class polyMesh;
122 
123 namespace vtk
124 {
125 
126 /*---------------------------------------------------------------------------*\
127  Class Foam::vtk::vtuSizing Declaration
128 \*---------------------------------------------------------------------------*/
129 
130 class vtuSizing
131 {
132 public:
133 
134  // Public Data
135 
136  //- Types of content that the storage may represent
137  enum contentType
138  {
139  LEGACY,
140  XML,
141  INTERNAL1,
142  INTERNAL2
143  };
144 
145  //- The possible storage 'slots' that can be used
146  enum slotType
147  {
148  CELLS,
149  CELLS_OFFSETS,
150  FACES,
153  };
154 
155 
156 private:
157 
158  // Private Member Data
159 
160  //- Polyhedral decomposition requested
161  bool decompose_;
162 
163  //- Number of cells in the mesh
164  label nCells_;
165 
166  //- Number of points in the mesh
167  label nPoints_;
168 
169  //- Number of vertex labels to represent the mesh
170  label nVertLabels_;
171 
172 
173  // Polyhedrals
174 
175  //- Number of polyhedral face labels for the mesh
176  label nFaceLabels_;
177 
178  //- Number of polyhedral cells (informational)
179  label nCellsPoly_;
180 
181  //- Number of vertex labels used by polyhedrals
182  label nVertPoly_;
183 
184 
185  // Decomposed polyhedrals
186 
187  //- Number of additional (decomposed) cells for the mesh
188  label nAddCells_;
189 
190  //- Number of additional (decomposed) points for the mesh
191  label nAddPoints_;
192 
193  //- Number of additional (decomposed) vertices for the mesh
194  label nAddVerts_;
195 
196 
197  // Private Member Functions
198 
199  //- set-size for cellMap and additionalIds
200  void presizeMaps(foamVtkMeshMaps& maps) const;
201 
202  //- Populate lists. For (legacy | xml | internal) VTK representations
203  template<class LabelType, class LabelType2>
204  static void populateArrays
205  (
206  const polyMesh& mesh,
207  const vtk::vtuSizing& sizing,
209  UList<LabelType>& vertLabels,
210  UList<LabelType>& vertOffset,
211  UList<LabelType>& faceLabels,
212  UList<LabelType>& faceOffset,
213  const enum contentType output,
214  UList<LabelType2>& cellMap,
215  UList<LabelType2>& addPointsIds
216  );
217 
218 
219 public:
220 
221  // Constructors
222 
223  //- Default construct
224  vtuSizing() noexcept;
225 
226  //- Construct sizing by analyzing the mesh.
227  // No polyhedral decomposition.
228  explicit vtuSizing(const polyMesh& mesh);
229 
230  //- Construct sizing by analyzing the mesh.
231  // Optionally with polyhedral decomposition.
232  vtuSizing(const polyMesh& mesh, const bool decompose);
233 
234 
235  // Member Functions
236 
237  // Edit
238 
239  //- Reset sizing by analyzing the mesh.
240  // Optionally with polyhedral decomposition.
241  void reset(const polyMesh& mesh, const bool decompose=false);
242 
243  //- Reset all sizes to zero.
244  void clear() noexcept;
245 
246 
247  // Access
248 
249  //- Query the decompose flag (normally off)
250  inline bool decompose() const;
251 
252  //- Number of cells for the mesh
253  inline label nCells() const;
254 
255  //- Number of points for the mesh
256  inline label nPoints() const;
257 
258  //- Number of vertex labels for the mesh
259  inline label nVertLabels() const;
260 
261  //- Number of polyhedral face labels for the mesh
262  inline label nFaceLabels() const;
263 
264  //- Number of polyhedral cells for the mesh
265  inline label nCellsPoly() const;
266 
267  //- Number of vertex labels for polyhedral cells of the mesh
268  inline label nVertPoly() const;
269 
270  //- Number of additional (decomposed) cells for the mesh
271  inline label nAddCells() const;
272 
273  //- Number of additional (decomposed) points for the mesh
274  inline label nAddPoints() const;
275 
276  //- Number of additional (decomposed) vertices for the mesh
277  inline label nAddVerts() const;
278 
279 
280  //- Number of field cells = nCells + nAddCells
281  inline label nFieldCells() const;
282 
283  //- Number of field points = nPoints + nAddPoints
284  inline label nFieldPoints() const;
285 
286 
287  // Derived sizes
288 
289  //- Return the required size for the storage slot
290  label sizeOf
291  (
292  const enum contentType output,
293  const enum slotType slot
294  ) const;
295 
296 
297  //- The calculated size for legacy storage
298  inline label sizeLegacy() const;
299 
300  //- The calculated size for legacy storage of the specified slot
301  inline label sizeLegacy(const enum slotType slot) const;
302 
303  //- The calculated size for xml storage of the specified slot
304  inline label sizeXml(const enum slotType slot) const;
305 
306  //- The calculated size for vtk-internal storage of the specified slot
307  inline label sizeInternal1(const enum slotType slot) const;
308 
309  //- The calculated size for vtk-internal storage of the specified slot
310  inline label sizeInternal2(const enum slotType slot) const;
311 
312 
313  // Routines for populating the output lists
314 
315  //- Populate lists for Legacy output
316  void populateLegacy
317  (
318  const polyMesh& mesh,
319  UList<uint8_t>& cellTypes,
320  labelUList& connectivity,
321  foamVtkMeshMaps& maps
322  ) const;
323 
324  //- Populate lists for XML output
325  void populateXml
326  (
327  const polyMesh& mesh,
328  UList<uint8_t>& cellTypes,
329  labelUList& connectivity,
330  labelUList& offsets,
331  labelUList& faces,
332  labelUList& facesOffsets,
333  foamVtkMeshMaps& maps
334  ) const;
335 
336 
337  // Internal types. The size of vtkIdType is unknown here
338 
339 #undef declarePopulateInternalMethod
340 #define declarePopulateInternalMethod(Type) \
341  \
342  \
343  void populateInternal \
344  ( \
345  const polyMesh& mesh, \
346  UList<uint8_t>& cellTypes, \
347  UList<Type>& connectivity, \
348  UList<Type>& offsets, \
349  UList<Type>& faces, \
350  UList<Type>& facesOffsets, \
351  foamVtkMeshMaps& maps, \
352  const enum contentType output \
353  ) const; \
354  \
355  \
356  void populateInternal \
357  ( \
358  const polyMesh& mesh, \
359  UList<uint8_t>& cellTypes, \
360  UList<Type>& connectivity, \
361  UList<Type>& offsets, \
362  UList<Type>& faces, \
363  UList<Type>& facesOffsets, \
364  labelUList& cellMap, \
365  labelUList& addPointsIds, \
366  const enum contentType output \
367  ) const
368 
369 
373 
374  #undef declarePopulateInternalMethod
375 
376 
377  // Routines for renumber vertices with a global point offset
378  // Legacy and xml only, internal version less likely to be needed
379 
380  //- Copy vertex labels with a global point offset - legacy format
382  (
383  const labelUList& connectivity,
384  const label globalPointOffset
385  );
386 
387  //- Copy vertex labels with a global point offset - XML format
389  (
390  const labelUList& connectivity,
391  const label globalPointOffset
392  );
393 
394  //- Copy faces stream labels with a global point offset - XML format
396  (
397  const labelUList& faceLabels,
398  const label globalPointOffset
399  );
400 
401  //- Copy face offsets with an offset from previous - XML format
403  (
404  const labelUList& faceOffsets,
405  const label prevOffset
406  );
407 
408  //- Renumber vertex labels by global point offset - legacy format
409  static void renumberVertLabelsLegacy
410  (
411  labelUList& connectivity,
412  const label globalPointOffset
413  );
414 
415  //- Renumber vertex labels by global point offset - XML format
416  static void renumberVertLabelsXml
417  (
418  labelUList& connectivity,
419  const label globalPointOffset
420  );
421 
422  //- Renumber faces stream labels by global point offset - XML format
423  static void renumberFaceLabelsXml
424  (
425  labelUList& faceLabels,
426  const label globalPointOffset
427  );
428 
429  //- Renumber face offsets with an offset from previous - XML format
430  static void renumberFaceOffsetsXml
431  (
432  labelUList& faceOffsets,
433  const label prevOffset
434  );
435 
436 
437  // Write
438 
439  //- Report some information
440  void info(Ostream& os) const;
441 
442 
443  // Member Operators
444 
445  //- Test equality
446  bool operator==(const vtuSizing& rhs) const;
447 
448  //- Test inequality
449  bool operator!=(const vtuSizing& rhs) const;
450 };
451 
452 
453 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
454 
455 } // End namespace vtk
456 } // End namespace Foam
457 
458 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
459 
460 #include "foamVtuSizingI.H"
461 
462 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
463 
464 #endif
465 
466 // ************************************************************************* //
Foam::vtk::vtuSizing::nCells
label nCells() const
Number of cells for the mesh.
Definition: foamVtuSizingI.H:38
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:71
Foam::vtk::vtuSizing::renumberFaceLabelsXml
static void renumberFaceLabelsXml(labelUList &faceLabels, const label globalPointOffset)
Renumber faces stream labels by global point offset - XML format.
Definition: foamVtuSizing.C:583
Foam::vtk::vtuSizing::CELLS_OFFSETS
Definition: foamVtuSizing.H:219
Foam::output
static Ostream & output(Ostream &os, const IntRange< T > &range)
Definition: IntRanges.C:66
Foam::vtk::vtuSizing::decompose
bool decompose() const
Query the decompose flag (normally off)
Definition: foamVtuSizingI.H:32
declarePopulateInternalMethod
#define declarePopulateInternalMethod(Type)
Definition: foamVtuSizing.H:410
Foam::vtk::vtuSizing::nCellsPoly
label nCellsPoly() const
Number of polyhedral cells for the mesh.
Definition: foamVtuSizingI.H:62
Foam::vtk::vtuSizing::populateLegacy
void populateLegacy(const polyMesh &mesh, UList< uint8_t > &cellTypes, labelUList &connectivity, foamVtkMeshMaps &maps) const
Populate lists for Legacy output.
Definition: foamVtuSizing.C:314
Foam::vtk::vtuSizing::nAddVerts
label nAddVerts() const
Number of additional (decomposed) vertices for the mesh.
Definition: foamVtuSizingI.H:86
Foam::vtk::vtuSizing::info
void info(Ostream &os) const
Report some information.
Definition: foamVtuSizing.C:663
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::vtk::vtuSizing::INTERNAL2
Internal vtkUnstructuredGrid content, VTK_CELL_ARRAY_V2.
Definition: foamVtuSizing.H:212
Foam::vtk::vtuSizing::copyVertLabelsLegacy
static labelList copyVertLabelsLegacy(const labelUList &connectivity, const label globalPointOffset)
Copy vertex labels with a global point offset - legacy format.
Definition: foamVtuSizing.C:443
Foam::vtk::vtuSizing::operator!=
bool operator!=(const vtuSizing &rhs) const
Test inequality.
Definition: foamVtuSizing.C:717
Foam::vtk::vtuSizing::sizeLegacy
label sizeLegacy() const
The calculated size for legacy storage.
Definition: foamVtuSizingI.H:104
Foam::vtk::vtuSizing::nFieldCells
label nFieldCells() const
Number of field cells = nCells + nAddCells.
Definition: foamVtuSizingI.H:92
Foam::vtk::vtuSizing::nVertPoly
label nVertPoly() const
Number of vertex labels for polyhedral cells of the mesh.
Definition: foamVtuSizingI.H:68
Foam::vtk::vtuSizing::populateXml
void populateXml(const polyMesh &mesh, UList< uint8_t > &cellTypes, labelUList &connectivity, labelUList &offsets, labelUList &faces, labelUList &facesOffsets, foamVtkMeshMaps &maps) const
Populate lists for XML output.
Definition: foamVtuSizing.C:343
Foam::vtk::vtuSizing::XML
XML (VTU) content.
Definition: foamVtuSizing.H:210
labelList.H
foamVtuSizingI.H
Foam::vtk::vtuSizing::nFaceLabels
label nFaceLabels() const
Number of polyhedral face labels for the mesh.
Definition: foamVtuSizingI.H:56
Foam::vtk::vtuSizing::FACES_OFFSETS
Faces end-offsets (XML) or locations (INTERNAL1)
Definition: foamVtuSizing.H:222
Foam::vtk::vtuSizing::renumberVertLabelsLegacy
static void renumberVertLabelsLegacy(labelUList &connectivity, const label globalPointOffset)
Renumber vertex labels by global point offset - legacy format.
Definition: foamVtuSizing.C:461
Foam::vtk::vtuSizing::nFieldPoints
label nFieldPoints() const
Number of field points = nPoints + nAddPoints.
Definition: foamVtuSizingI.H:98
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::vtk::vtuSizing::nPoints
label nPoints() const
Number of points for the mesh.
Definition: foamVtuSizingI.H:44
Foam::vtk::vtuSizing::CELLS
Cell connectivity (ALL)
Definition: foamVtuSizing.H:218
Foam::vtk::vtuSizing::nVertLabels
label nVertLabels() const
Number of vertex labels for the mesh.
Definition: foamVtuSizingI.H:50
Foam::vtk::vtuSizing::contentType
contentType
Types of content that the storage may represent.
Definition: foamVtuSizing.H:207
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::vtk::vtuSizing
Sizing descriptions and routines for transcribing an OpenFOAM volume mesh into a VTK unstructured gri...
Definition: foamVtuSizing.H:200
Foam::vtk::vtuSizing::renumberFaceOffsetsXml
static void renumberFaceOffsetsXml(labelUList &faceOffsets, const label prevOffset)
Renumber face offsets with an offset from previous - XML format.
Definition: foamVtuSizing.C:638
Foam::vtk::vtuSizing::sizeInternal2
label sizeInternal2(const enum slotType slot) const
The calculated size for vtk-internal storage of the specified slot.
Definition: foamVtuSizingI.H:138
Foam::vtk::vtuSizing::FACES
Face-stream (XML, INTERNAL)
Definition: foamVtuSizing.H:221
Foam::vtk::vtuSizing::clear
void clear() noexcept
Reset all sizes to zero.
Definition: foamVtuSizing.C:66
Foam::vtk::vtuSizing::copyVertLabelsXml
static labelList copyVertLabelsXml(const labelUList &connectivity, const label globalPointOffset)
Copy vertex labels with a global point offset - XML format.
Definition: foamVtuSizing.C:526
cellTypes
const labelList & cellTypes
Definition: setCellMask.H:33
Foam::vtk::vtuSizing::copyFaceOffsetsXml
static labelList copyFaceOffsetsXml(const labelUList &faceOffsets, const label prevOffset)
Copy face offsets with an offset from previous - XML format.
Definition: foamVtuSizing.C:620
label.H
Foam::List< label >
Foam::vtk::vtuSizing::renumberVertLabelsXml
static void renumberVertLabelsXml(labelUList &connectivity, const label globalPointOffset)
Renumber vertex labels by global point offset - XML format.
Definition: foamVtuSizing.C:544
Foam::vtk::vtuSizing::copyFaceLabelsXml
static labelList copyFaceLabelsXml(const labelUList &faceLabels, const label globalPointOffset)
Copy faces stream labels with a global point offset - XML format.
Definition: foamVtuSizing.C:565
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::vtuSizing::LEGACY
Legacy VTK content.
Definition: foamVtuSizing.H:209
Foam::vtk::vtuSizing::vtuSizing
vtuSizing() noexcept
Default construct.
Definition: foamVtuSizing.C:47
Foam::vtk::vtuSizing::slotType
slotType
The possible storage 'slots' that can be used.
Definition: foamVtuSizing.H:216
Foam::vtk::vtuSizing::nAddPoints
label nAddPoints() const
Number of additional (decomposed) points for the mesh.
Definition: foamVtuSizingI.H:80
Foam::vtk::vtuSizing::sizeInternal1
label sizeInternal1(const enum slotType slot) const
The calculated size for vtk-internal storage of the specified slot.
Definition: foamVtuSizingI.H:129
Foam::vtk::vtuSizing::sizeXml
label sizeXml(const enum slotType slot) const
The calculated size for xml storage of the specified slot.
Definition: foamVtuSizingI.H:120
Foam::labelUList
UList< label > labelUList
A UList of labels.
Definition: UList.H:80
Foam::vtk::vtuSizing::nAddCells
label nAddCells() const
Number of additional (decomposed) cells for the mesh.
Definition: foamVtuSizingI.H:74
Foam::vtk::vtuSizing::operator==
bool operator==(const vtuSizing &rhs) const
Test equality.
Definition: foamVtuSizing.C:698
Foam::vtk::vtuSizing::reset
void reset(const polyMesh &mesh, const bool decompose=false)
Reset sizing by analyzing the mesh.
Definition: foamVtuSizing.C:84
foamVtkMeshMaps.H
Foam::foamVtkMeshMaps
Bookkeeping for mesh subsetting and/or polyhedral cell decomposition. Although the main use case is f...
Definition: foamVtkMeshMaps.H:58
Foam::vtk::vtuSizing::INTERNAL1
Internal vtkUnstructuredGrid content.
Definition: foamVtuSizing.H:211
Foam::vtk::vtuSizing::sizeOf
label sizeOf(const enum contentType output, const enum slotType slot) const
Return the required size for the storage slot.
Definition: foamVtuSizing.C:207