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-2021 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 (vtk DataFile Version 2.0)
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 Foam_vtk_vtuSizing_H
109 #define Foam_vtk_vtuSizing_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 cellShape;
122 class polyMesh;
123 
124 namespace vtk
125 {
126 
127 /*---------------------------------------------------------------------------*\
128  Class Foam::vtk::vtuSizing Declaration
129 \*---------------------------------------------------------------------------*/
130 
131 class vtuSizing
132 {
133 public:
134 
135  // Public Data
136 
137  //- Types of content that the storage may represent
138  enum contentType : char
139  {
140  LEGACY,
141  XML,
142  INTERNAL1,
143  INTERNAL2
144  };
145 
146  //- The possible storage 'slots' that can be used
147  enum slotType : char
148  {
149  CELLS,
150  CELLS_OFFSETS,
151  FACES,
154  };
155 
156  //- How the mesh cells have been selected or defined
157  enum selectionModeType : char
158  {
159  FULL_MESH,
160  SUBSET_MESH,
161  SHAPE_MESH
162  };
163 
164 
165 private:
166 
167  // Private Member Data
168 
169  //- Polyhedral decomposition requested
170  bool decompose_;
171 
172  //- Manifold cells detected
173  bool manifold_;
174 
175  //- How the mesh cells have been selected or defined
176  selectionModeType selectionMode_;
177 
178  //- Number of cells in the mesh
179  label nCells_;
180 
181  //- Number of points in the mesh
182  label nPoints_;
183 
184  //- Number of vertex labels to represent the mesh
185  label nVertLabels_;
186 
187 
188  // Polyhedrals
189 
190  //- Number of polyhedral face labels for the mesh
191  label nFaceLabels_;
192 
193  //- Number of polyhedral cells (informational)
194  label nCellsPoly_;
195 
196  //- Number of vertex labels used by polyhedrals
197  label nVertPoly_;
198 
199 
200  // Decomposed polyhedrals
201 
202  //- Number of additional (decomposed) cells for the mesh
203  label nAddCells_;
204 
205  //- Number of additional (decomposed) points for the mesh
206  label nAddPoints_;
207 
208  //- Number of additional (decomposed) vertices for the mesh
209  label nAddVerts_;
212  // Private Member Functions
213 
214  //- set-size for cellMap and additionalIds
215  void presizeMaps(foamVtkMeshMaps& maps) const;
216 
217  //- Verify list sizes
218  static void checkSizes
219  (
220  const vtk::vtuSizing& sizing,
221 
222  const label cellTypes_size,
223  const label vertLabels_size,
224  const label vertOffset_size,
225  const label faceLabels_size,
226  const label faceOffset_size,
227 
228  const enum contentType output,
230  const label cellMap_size,
231  const label addPointsIds_size
232  );
234  //- Adjust cell/face offsets
235  template<class LabelType>
236  static void adjustOffsets
237  (
238  UList<LabelType>& vertOffset,
239  UList<LabelType>& faceOffset,
240  const enum contentType output,
241  const bool hasFaceStream
242  );
243 
244  //- Populate lists. For (legacy | xml | internal) VTK representations
245  template<class LabelType>
246  static void populateArrays
247  (
248  const polyMesh& mesh,
249  const vtk::vtuSizing& sizing,
251  UList<LabelType>& vertLabels,
252  UList<LabelType>& vertOffset,
254  UList<LabelType>& faceOffset,
255  const enum contentType output,
256  labelUList& cellMap,
257  labelUList& addPointsIds
258  );
259 
260  //- Populate lists - Primitive mesh shapes only!!
261  template<class LabelType>
262  static void populateArrays
263  (
264  const UList<cellShape>& shapes,
265  const vtk::vtuSizing& sizing,
267  UList<LabelType>& vertLabels,
268  UList<LabelType>& vertOffset,
269  UList<LabelType>& faceLabels, // unused
270  UList<LabelType>& faceOffset, // unused
271  const enum contentType output,
272  labelUList& cellMap,
273  labelUList& addPointsIds // unused
274  );
275 
276 
277 protected:
278 
279  // Protected Member Functions
280 
281  //- Reset list for primitive shapes only (ADVANCED USAGE)
283  (
284  const UList<cellShape>& shapes,
286  labelUList& connectivity,
287  foamVtkMeshMaps& maps
288  ) const;
289 
290  //- Reset list for primitive shapes only (ADVANCED USAGE)
291  void populateShapesXml
292  (
293  const UList<cellShape>& shapes,
295  labelUList& connectivity,
296  labelUList& offsets,
297  labelUList& faces,
298  labelUList& facesOffsets,
299  foamVtkMeshMaps& maps
300  ) const;
301 
302 
303 public:
304 
305  // Constructors
306 
307  //- Default construct
309 
310  //- Construct sizing by analyzing the mesh.
311  // No polyhedral decomposition.
312  explicit vtuSizing(const polyMesh& mesh);
313 
314  //- Construct sizing by analyzing the mesh.
315  // Optionally with polyhedral decomposition.
316  vtuSizing(const polyMesh& mesh, const bool decompose);
317 
318 
319  // Member Functions
320 
321  // Edit
322 
323  //- Reset sizing by analyzing the mesh.
324  // Optionally with polyhedral decomposition.
325  void reset(const polyMesh& mesh, const bool decompose=false);
326 
327  //- Reset sizing by analyzing a subset of the mesh.
328  // A labelUList::null() selection corresponds to the entire mesh.
329  // Polyhedral decomposition is disabled for subsets.
330  void reset
331  (
332  const polyMesh& mesh,
333  const labelUList& subsetCellsIds,
334  const bool decompose = false
335  );
336 
337  //- Reset sizing using primitive shapes only (ADVANCED USAGE)
338  // Effectively removes any polyhedrals!
339  void resetShapes(const UList<cellShape>& shapes);
340 
341  //- Reset all sizes to zero.
342  void clear() noexcept;
343 
344 
345  // Access
346 
347  //- Query the decompose flag (normally off)
348  inline bool decompose() const noexcept;
349 
350  //- Manifold mesh cells detected? Globally consistent quantity.
351  inline bool manifold() const noexcept;
352 
353  //- Query how the mesh cells have been selected or defined
355 
356  //- Number of cells for the mesh
357  inline label nCells() const noexcept;
358 
359  //- Number of points for the mesh
360  inline label nPoints() const noexcept;
361 
362  //- Number of vertex labels for the mesh
363  inline label nVertLabels() const noexcept;
364 
365  //- Number of polyhedral face labels for the mesh
366  inline label nFaceLabels() const noexcept;
367 
368  //- Number of polyhedral cells for the mesh
369  inline label nCellsPoly() const noexcept;
370 
371  //- Number of vertex labels for polyhedral cells of the mesh
372  inline label nVertPoly() const noexcept;
373 
374  //- Number of additional (decomposed) cells for the mesh
375  inline label nAddCells() const noexcept;
376 
377  //- Number of additional (decomposed) points for the mesh
378  inline label nAddPoints() const noexcept;
379 
380  //- Number of additional (decomposed) vertices for the mesh
381  inline label nAddVerts() const noexcept;
382 
383 
384  //- Number of field cells = nCells + nAddCells
385  inline label nFieldCells() const noexcept;
386 
387  //- Number of field points = nPoints + nAddPoints
388  inline label nFieldPoints() const noexcept;
389 
390 
391  // Edit
392 
393  //- Alter number of mesh points (ADVANCED USAGE)
394  void setNumPoints(label n) noexcept { nPoints_ = n; }
395 
396  //- Alter number of additional (cell-centre) points (ADVANCED USAGE)
397  void setNumAddPoints(label n) noexcept { nAddPoints_ = n; }
398 
399 
400  // Derived Sizes
401 
402  //- Return the required size for the storage slot
403  label sizeOf
404  (
405  const enum contentType output,
406  const enum slotType slot
407  ) const;
408 
409 
410  //- The calculated size for legacy storage
411  inline label sizeLegacy() const;
412 
413  //- The calculated size for legacy storage of the specified slot
414  inline label sizeLegacy(const enum slotType slot) const;
415 
416  //- The calculated size for xml storage of the specified slot
417  inline label sizeXml(const enum slotType slot) const;
418 
419  //- The calculated size for vtk-internal storage of the specified slot
420  inline label sizeInternal1(const enum slotType slot) const;
421 
422  //- The calculated size for vtk-internal storage of the specified slot
423  inline label sizeInternal2(const enum slotType slot) const;
424 
425 
426  // Routines for populating the output lists
427 
428  //- Populate lists for Legacy output
429  void populateLegacy
430  (
431  const polyMesh& mesh,
433  labelUList& connectivity,
434  foamVtkMeshMaps& maps
435  ) const;
436 
437  //- Populate lists for XML output
438  void populateXml
439  (
440  const polyMesh& mesh,
442  labelUList& connectivity,
443  labelUList& offsets,
444  labelUList& faces,
445  labelUList& facesOffsets,
446  foamVtkMeshMaps& maps
447  ) const;
448 
449 
450  // Internal types. The size of vtkIdType is unknown here
451 
452 #undef declarePopulateInternalMethod
453 #define declarePopulateInternalMethod(Type) \
454  \
455  \
456  void populateInternal \
457  ( \
458  const polyMesh& mesh, \
459  UList<uint8_t>& cellTypes, \
460  UList<Type>& connectivity, \
461  UList<Type>& offsets, \
462  UList<Type>& faces, \
463  UList<Type>& facesOffsets, \
464  foamVtkMeshMaps& maps, \
465  const enum contentType output \
466  ) const; \
467  \
468  \
469  void populateInternal \
470  ( \
471  const polyMesh& mesh, \
472  UList<uint8_t>& cellTypes, \
473  UList<Type>& connectivity, \
474  UList<Type>& offsets, \
475  UList<Type>& faces, \
476  UList<Type>& facesOffsets, \
477  labelUList& cellMap, \
478  labelUList& addPointsIds, \
479  const enum contentType output \
480  ) const
481 
482 
486 
487  #undef declarePopulateInternalMethod
488 
489 
490  // Routines for renumber vertices with a global point offset
491  // Legacy and xml only, internal version less likely to be needed
492 
493  //- Copy vertex labels with a global point offset - legacy format
495  (
496  const labelUList& connectivity,
497  const label globalPointOffset
498  );
499 
500  //- Copy vertex labels with a global point offset - XML format
502  (
503  const labelUList& connectivity,
504  const label globalPointOffset
505  );
506 
507  //- Copy faces stream labels with a global point offset - XML format
509  (
510  const labelUList& faceLabels,
511  const label globalPointOffset
512  );
513 
514  //- Copy face offsets with an offset from previous - XML format
516  (
517  const labelUList& faceOffsets,
518  const label prevOffset
519  );
520 
521  //- Renumber vertex labels by global point offset - legacy format
522  static void renumberVertLabelsLegacy
523  (
524  labelUList& connectivity,
525  const label globalPointOffset
526  );
527 
528  //- Renumber vertex labels by global point offset - XML format
529  static void renumberVertLabelsXml
530  (
531  labelUList& connectivity,
532  const label globalPointOffset
533  );
534 
535  //- Renumber faces stream labels by global point offset - XML format
536  static void renumberFaceLabelsXml
537  (
539  const label globalPointOffset
540  );
541 
542  //- Renumber face offsets with an offset from previous - XML format
543  static void renumberFaceOffsetsXml
544  (
545  labelUList& faceOffsets,
546  const label prevOffset
547  );
548 
549 
550  // Write
551 
552  //- Report some information
553  void info(Ostream& os) const;
554 
555 
556  // Member Operators
557 
558  //- Test equality
559  bool operator==(const vtuSizing& rhs) const;
560 
561  //- Test inequality
562  bool operator!=(const vtuSizing& rhs) const;
563 };
564 
565 
566 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
567 
568 } // End namespace vtk
569 } // End namespace Foam
570 
571 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
572 
573 #include "foamVtuSizingI.H"
574 
575 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
576 
577 #endif
578 
579 // ************************************************************************* //
vtuSizing() noexcept
Default construct.
void populateShapesLegacy(const UList< cellShape > &shapes, UList< uint8_t > &cellTypes, labelUList &connectivity, foamVtkMeshMaps &maps) const
Reset list for primitive shapes only (ADVANCED USAGE)
selectionModeType
How the mesh cells have been selected or defined.
label nVertLabels() const noexcept
Number of vertex labels for the mesh.
contentType
Types of content that the storage may represent.
static void renumberFaceOffsetsXml(labelUList &faceOffsets, const label prevOffset)
Renumber face offsets with an offset from previous - XML format.
bool operator==(const vtuSizing &rhs) const
Test equality.
label nAddVerts() const noexcept
Number of additional (decomposed) vertices for the mesh.
void populateLegacy(const polyMesh &mesh, UList< uint8_t > &cellTypes, labelUList &connectivity, foamVtkMeshMaps &maps) const
Populate lists for Legacy output.
An analytical geometric cellShape.
Definition: cellShape.H:68
void reset(const polyMesh &mesh, const bool decompose=false)
Reset sizing by analyzing the mesh.
bool decompose() const noexcept
Query the decompose flag (normally off)
label nVertPoly() const noexcept
Number of vertex labels for polyhedral cells of the mesh.
#define declarePopulateInternalMethod(Type)
Internal vtkUnstructuredGrid content.
label nFaceLabels() const noexcept
Number of polyhedral face labels for the mesh.
Face-stream (XML, INTERNAL)
label nFieldPoints() const noexcept
Number of field points = nPoints + nAddPoints.
labelList faceLabels(nFaceLabels)
Bookkeeping for mesh subsetting and/or polyhedral cell decomposition. Although the main use case is f...
slotType
The possible storage &#39;slots&#39; that can be used.
XML (VTU) content.
void populateShapesXml(const UList< cellShape > &shapes, UList< uint8_t > &cellTypes, labelUList &connectivity, labelUList &offsets, labelUList &faces, labelUList &facesOffsets, foamVtkMeshMaps &maps) const
Reset list for primitive shapes only (ADVANCED USAGE)
Faces end-offsets (XML) or locations (INTERNAL1)
dynamicFvMesh & mesh
static labelList copyFaceLabelsXml(const labelUList &faceLabels, const label globalPointOffset)
Copy faces stream labels with a global point offset - XML format.
static void renumberVertLabelsLegacy(labelUList &connectivity, const label globalPointOffset)
Renumber vertex labels by global point offset - legacy format.
bool operator!=(const vtuSizing &rhs) const
Test inequality.
void setNumPoints(label n) noexcept
Alter number of mesh points (ADVANCED USAGE)
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.
void info(Ostream &os) const
Report some information.
label nAddPoints() const noexcept
Number of additional (decomposed) points for the mesh.
static labelList copyVertLabelsLegacy(const labelUList &connectivity, const label globalPointOffset)
Copy vertex labels with a global point offset - legacy format.
label sizeInternal2(const enum slotType slot) const
The calculated size for vtk-internal storage of the specified slot.
label nCells() const noexcept
Number of cells for the mesh.
bool manifold() const noexcept
Manifold mesh cells detected? Globally consistent quantity.
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:105
const labelList & cellTypes
Definition: setCellMask.H:27
static labelList copyVertLabelsXml(const labelUList &connectivity, const label globalPointOffset)
Copy vertex labels with a global point offset - XML format.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
const direction noexcept
Definition: Scalar.H:258
OBJstream os(runTime.globalPath()/outputName)
Cell connectivity (ALL)
label nFieldCells() const noexcept
Number of field cells = nCells + nAddCells.
Legacy VTK content.
Sizing descriptions and routines for transcribing an OpenFOAM volume mesh into a VTK unstructured gri...
label sizeInternal1(const enum slotType slot) const
The calculated size for vtk-internal storage of the specified slot.
static void renumberFaceLabelsXml(labelUList &faceLabels, const label globalPointOffset)
Renumber faces stream labels by global point offset - XML format.
static labelList copyFaceOffsetsXml(const labelUList &faceOffsets, const label prevOffset)
Copy face offsets with an offset from previous - XML format.
label nPoints() const noexcept
Number of points for the mesh.
label sizeOf(const enum contentType output, const enum slotType slot) const
Return the required size for the storage slot.
label sizeLegacy() const
The calculated size for legacy storage.
static Ostream & output(Ostream &os, const IntRange< T > &range)
Definition: IntRanges.C:44
label n
label nAddCells() const noexcept
Number of additional (decomposed) cells for the mesh.
label nCellsPoly() const noexcept
Number of polyhedral cells for the mesh.
selectionModeType selectionMode() const noexcept
Query how the mesh cells have been selected or defined.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
void clear() noexcept
Reset all sizes to zero.
void resetShapes(const UList< cellShape > &shapes)
Reset sizing using primitive shapes only (ADVANCED USAGE)
void setNumAddPoints(label n) noexcept
Alter number of additional (cell-centre) points (ADVANCED USAGE)
Internal vtkUnstructuredGrid content, VTK_CELL_ARRAY_V2.
static void renumberVertLabelsXml(labelUList &connectivity, const label globalPointOffset)
Renumber vertex labels by global point offset - XML format.
label sizeXml(const enum slotType slot) const
The calculated size for xml storage of the specified slot.
Namespace for OpenFOAM.