Loading...
Searching...
No Matches
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-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
26Class
27 Foam::vtk::vtuSizing
28
29Description
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, id0,id1,..., nFace1Pts, id0,...]
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, id0,id1,..., nFace1Pts, id0,...]
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, id0,id1,..., nFace1Pts, id0,...]
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, id0,id1,..., nFace1Pts, id0,...]
85 \c facelocation | begin location for each of \c faces, with -1 for primitive cells
86 \endtable
87
88 \table
89 hdf storage (similar to internal2)
90 \c Types | vtk cell type (1-255)
91 \c Connectivity | unique vertex labels used by the cell
92 \c Offsets | begin/end offsets for \c Connectivity
93 \c FaceConnectivity | faces for polyhedral cells (their vertex labels)
94 | [id0,id1,...]
95 \c FaceOffsets | begin/end offsets for FaceConnectivity.
96 \c PolyhedronToFaces | face ids (faces are defined by FaceOffsets/FaceConnectivity)
97 \c PolyhedronOffsets | begin/end offsets for \c PolyhedronToFaces
98 | Has size nCells+1
99 \endtable
100
101 The VTK storage concept for "connectivity" and "faces" somewhat resemble
102 a CompactListList.
103
104Note
105 It is possible to specify a global point offset (via the globalIndex)
106 so that the cell point labels will use global numbering.
107 There is no support for point renumbering with merged mesh points,
108 since it likely more efficient to use VTK point-blanking to mark duplicate
109 points instead of merging points ourselves.
110
111Note
112 The VTK_CELL_ARRAY_V2 define (from vtkCellArray.h) indicates if the new
113 (internal2) new format is being used.
114
115SourceFiles
116 foamVtuSizing.cxx
117 foamVtuSizingImpl.cxx
118 foamVtuSizingI.H
119
120\*---------------------------------------------------------------------------*/
121
122#ifndef Foam_vtk_vtuSizing_H
123#define Foam_vtk_vtuSizing_H
124
125#include "foamVtkMeshMaps.H"
126
127// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
128
129namespace Foam
130{
131
132// Forward Declarations
133class cellShape;
134class polyMesh;
135
136namespace vtk
137{
138
139/*---------------------------------------------------------------------------*\
140 Class Foam::vtk::vtuSizing Declaration
141\*---------------------------------------------------------------------------*/
142
143class vtuSizing
144{
145public:
146
147 // Public Data
148
149 //- Types of content that the storage may represent
150 enum class contentType : unsigned char
151 {
152 LEGACY,
153 XML,
154 INTERNAL1,
155 INTERNAL2,
156 HDF,
157 };
158
159 //- The possible storage 'slots' that can be used
160 enum class slotType : unsigned char
161 {
163 CELLS,
169 FACES,
179 };
180
181 //- How the mesh cells have been selected or defined
182 enum class selectionModeType : unsigned char
183 {
184 FULL_MESH,
187 };
188
189
190private:
191
192 // Private Member Data
193
194 //- Polyhedral decomposition requested
195 bool decompose_;
196
197 //- Manifold cells detected
198 bool manifold_;
199
200 //- How the mesh cells have been selected or defined
201 selectionModeType selectionMode_;
202
203 //- Number of cells in the mesh
204 label nCells_;
205
206 //- Number of points in the mesh
207 label nPoints_;
208
209 //- Number of vertex labels to represent the mesh
210 label nVertLabels_;
211
212
213 // Polyhedrals
214
215 //- Number of polyhedral face labels for the mesh
216 label nFaceLabels_;
217
218 //- Number of polyhedral cells
219 label nCellsPoly_;
220
221 //- Number of (non-unique) faces used by polyhedral cells
222 label nFacesPoly_;
223
224 //- Number of vertex labels used by polyhedrals
225 label nVertPoly_;
226
227
228 // Decomposed polyhedrals
229
230 //- Number of additional (decomposed) cells for the mesh
231 label nAddCells_;
232
233 //- Number of additional (decomposed) points for the mesh
234 label nAddPoints_;
235
236 //- Number of additional (decomposed) vertices for the mesh
237 label nAddVerts_;
239
240 // Private Member Functions
241
242 //- Allocate for cellMap and additionalIds
243 void presizeMaps(foamVtkMeshMaps& maps) const;
244
245 //- Verify list sizes
246 static void checkSizes
248 const vtk::vtuSizing& sizing,
250 const label cellTypes_size,
251 const label vertLabels_size, const label vertOffset_size,
252 const label faceLabels_size, const label faceOffset_size,
254 // HDF only
255 const label polyFaceIds_size, const label polyFaceOffset_size,
256
257 const enum contentType output,
258
259 const label cellMap_size,
260 const label addPointsIds_size
261 );
263 //- Populate lists.
264 //- For (legacy | xml | internal | hdf) VTK representations
265 template<class LabelType>
266 static void populateArrays
267 (
269 const vtk::vtuSizing& sizing,
271 UList<LabelType>& vertLabels,
272 UList<LabelType>& vertOffset,
274 UList<LabelType>& faceOffset,
275 UList<LabelType>& polyFaceIds, // HDF only
276 UList<LabelType>& polyFaceOffsets, // HDF only
277 const enum contentType output,
278 labelUList& cellMap,
279 labelUList& addPointsIds
280 );
281
282 //- Populate lists - Primitive mesh shapes only!!
283 template<class LabelType>
284 static void populateArrays
286 const UList<cellShape>& shapes,
287 const vtk::vtuSizing& sizing,
289 UList<LabelType>& vertLabels,
290 UList<LabelType>& vertOffset,
291 UList<LabelType>& faceLabels, // unused
292 UList<LabelType>& faceOffset, // unused
293 UList<LabelType>& polyFaceIds, // unused
294 UList<LabelType>& polyFaceOffsets, // unused
295 const enum contentType output,
296 labelUList& cellMap,
297 labelUList& addPointsIds // unused
298 );
299
300
301protected:
302
303 // Protected Member Functions
304
305 //- Reset list for primitive shapes only (ADVANCED USAGE)
307 (
308 const UList<cellShape>& shapes,
310 labelUList& connectivity,
311 foamVtkMeshMaps& maps
312 ) const;
313
314 //- Reset list for primitive shapes only (ADVANCED USAGE)
316 (
317 const UList<cellShape>& shapes,
319 labelUList& connectivity,
320 labelUList& offsets,
321 labelUList& faces,
322 labelUList& facesOffsets,
323 foamVtkMeshMaps& maps
324 ) const;
325
326 //- Return a list of dummy faceOffsets
328 (
329 const label numCells,
330 const enum contentType output,
331 label beginOffset = 0
332 );
333
334public:
335
336 // Constructors
337
338 //- Default construct
340
341 //- Construct sizing by analyzing the mesh.
342 // No polyhedral decomposition.
343 explicit vtuSizing(const polyMesh& mesh);
344
345 //- Construct sizing by analyzing the mesh.
346 // Optionally with polyhedral decomposition.
347 vtuSizing(const polyMesh& mesh, const bool decompose);
348
349
350 // Member Functions
351
352 // Edit
353
354 //- Reset sizing by analyzing the mesh.
355 // Optionally with polyhedral decomposition.
356 void reset(const polyMesh& mesh, const bool decompose=false);
357
358 //- Reset sizing by analyzing a subset of the mesh.
359 // A labelUList::null() selection corresponds to the entire mesh.
360 // Polyhedral decomposition is disabled for subsets.
361 void reset
362 (
363 const polyMesh& mesh,
364 const labelUList& subsetCellsIds,
365 const bool decompose = false
366 );
367
368 //- Reset sizing using primitive shapes only (ADVANCED USAGE)
369 // Effectively removes any polyhedrals!
370 void resetShapes(const UList<cellShape>& shapes);
371
372 //- Reset all sizes to zero.
373 void clear() noexcept;
374
375
376 // Access
377
378 //- Query the decompose flag (normally off)
379 inline bool decompose() const noexcept;
380
381 //- Manifold mesh cells detected? Globally consistent quantity.
382 inline bool manifold() const noexcept;
383
384 //- Query how the mesh cells have been selected or defined
386
387 //- Has polyhedral cells?
388 bool hasPolyCells() const noexcept { return nCellsPoly_; }
389
390 //- Number of cells for the mesh
391 inline label nCells() const noexcept;
392
393 //- Number of points for the mesh
394 inline label nPoints() const noexcept;
395
396 //- Number of vertex labels for the mesh
397 inline label nVertLabels() const noexcept;
398
399 //- Number of polyhedral face labels for the mesh,
400 //- without any prefixing information
401 inline label nFaceLabels() const noexcept;
402
403 //- Number of polyhedral face labels for the mesh
404 //- with format-dependent prefixing information
405 label nFaceLabels(contentType output) const;
406
407 //- Number of output polyhedral cells for the mesh
408 inline label nCellsPoly() const noexcept;
409
410 //- Number of (non-unique) faces used by polyhedral cells
411 inline label nFacesPoly() const noexcept;
412
413 //- Number of vertex labels for polyhedral cells of the mesh
414 inline label nVertPoly() const noexcept;
415
416 //- Number of additional (decomposed) cells for the mesh
417 inline label nAddCells() const noexcept;
418
419 //- Number of additional (decomposed) points for the mesh
420 inline label nAddPoints() const noexcept;
421
422 //- Number of additional (decomposed) vertices for the mesh
423 inline label nAddVerts() const noexcept;
424
425
426 //- Number of field cells = nCells + nAddCells
427 inline label nFieldCells() const noexcept;
428
429 //- Number of field points = nPoints + nAddPoints
430 inline label nFieldPoints() const noexcept;
431
432
433 // Edit
434
435 //- Alter number of mesh points (ADVANCED USAGE)
436 void setNumPoints(label n) noexcept { nPoints_ = n; }
437
438 //- Alter number of additional (cell-centre) points (ADVANCED USAGE)
439 void setNumAddPoints(label n) noexcept { nAddPoints_ = n; }
440
441
442 // Derived Sizes
444 //- Return the required size for the storage slot
445 label sizeOf
446 (
447 const enum contentType output,
448 const enum slotType slot
449 ) const noexcept;
450
451 //- Return the required size for the storage slot
452 template<vtuSizing::slotType slot>
453 label sizeOf(contentType output) const noexcept
455 return sizeOf(output, slot);
456 }
457
458
459 // Routines for populating the output lists
460
461 //- Populate lists for Legacy output
462 void populateLegacy
463 (
464 const polyMesh& mesh,
466 labelUList& connectivity,
467 foamVtkMeshMaps& maps
468 ) const;
469
470 //- Populate lists for XML output
471 void populateXml
472 (
473 const polyMesh& mesh,
475 labelUList& connectivity,
476 labelUList& offsets,
477 labelUList& faces,
478 labelUList& facesOffsets,
479 foamVtkMeshMaps& maps
480 ) const;
481
482 //- Populate lists for HDF output
483 void populateHdf
484 (
485 const polyMesh& mesh,
487 labelUList& connectivity,
488 labelUList& offsets,
490 labelUList& facesOffsets,
491 labelUList& polyFaceIds,
492 labelUList& polyFaceOffsets,
493 foamVtkMeshMaps& maps
494 ) const;
495
497 // Internal types. The size of vtkIdType is unknown here
498
499#undef declarePopulateInternalMethod
500#define declarePopulateInternalMethod(Type) \
501 \
502 \
503 void populateInternal \
504 ( \
505 const polyMesh& mesh, \
506 UList<uint8_t>& cellTypes, \
507 UList<Type>& connectivity, \
508 UList<Type>& offsets, \
509 UList<Type>& faces, \
510 UList<Type>& facesOffsets, \
511 foamVtkMeshMaps& maps, \
512 const enum contentType output \
513 ) const; \
514 \
515 \
516 void populateInternal \
517 ( \
518 const polyMesh& mesh, \
519 UList<uint8_t>& cellTypes, \
520 UList<Type>& connectivity, \
521 UList<Type>& offsets, \
522 UList<Type>& faces, \
523 UList<Type>& facesOffsets, \
524 labelUList& cellMap, \
525 labelUList& addPointsIds, \
526 const enum contentType output \
527 ) const
529
534 #undef declarePopulateInternalMethod
535
536
537 // Routines for renumber vertices with a global point offset
538 // Legacy and xml only, internal version less likely to be needed
539
540 //- Renumber vertex labels by (global) point offset
541 static void renumberVertLabels
542 (
543 labelUList& connectivity,
544 const label pointOffset,
545 const contentType output
546 );
547
548 //- Renumber faces stream labels by (global) point offset
549 static void renumberFaceLabels
550 (
552 const label pointOffset,
553 const contentType output
554 );
555
556 //- Renumber face offsets with a specified (global) begin offset
557 static void renumberFaceOffsets
558 (
559 labelUList& faceOffsets,
560 const label beginOffset,
561 const contentType output
562 );
563
564
565 // Write
566
567 //- Report some information
568 void info(Ostream& os) const;
569
570
571 // Member Operators
572
573 //- Test equality
574 bool operator==(const vtuSizing& rhs) const;
575
576 //- Test inequality
577 bool operator!=(const vtuSizing& rhs) const;
578
579
580 // Housekeeping
581
582 //- Copy vertex labels with a (global) point offset
584 (
585 const labelUList& connectivity,
586 const label pointOffset,
587 const contentType output
588 )
589 {
590 labelList result(connectivity);
591 renumberVertLabels(result, pointOffset, output);
592 return result;
593 }
594
595 //- Copy faces stream labels with a (global) point offset
597 (
598 const labelUList& faceLabels,
599 const label pointOffset,
600 const contentType output
601 )
602 {
603 labelList result(faceLabels);
604 renumberFaceLabels(result, pointOffset, output);
605 return result;
606 }
607
608 //- Copy face offsets with a specified (global) begin offset
610 (
611 const labelUList& faceOffsets,
612 const label beginOffset,
613 const contentType output
614 )
615 {
616 labelList result(faceOffsets);
617 renumberFaceOffsets(result, beginOffset, output);
618 return result;
619 }
620
621 //- Renumber vertex labels by (global) point offset - legacy format
622 FOAM_DEPRECATED_FOR(2025-11, "renumberVertLabels(...) method")
623 static void renumberVertLabelsLegacy
624 (
625 labelUList& connectivity,
626 const label pointOffset
627 )
628 {
629 renumberVertLabels(connectivity, pointOffset, contentType::LEGACY);
630 }
631
632 //- Renumber vertex labels by (global) point offset - XML format
633 FOAM_DEPRECATED_FOR(2025-11, "renumberVertLabels(...) method")
634 static void renumberVertLabelsXml
635 (
636 labelUList& connectivity,
637 const label pointOffset
638 )
639 {
640 renumberVertLabels(connectivity, pointOffset, contentType::XML);
641 }
642
643 //- Renumber faces stream labels by (global) point offset - XML format
644 FOAM_DEPRECATED_FOR(2025-11, "renumberFaceLabels(...) method")
646 (
648 const label pointOffset
649 )
650 {
652 }
653
654 //- Renumber face offsets with a specified begin offset - XML format
655 FOAM_DEPRECATED_FOR(2025-11, "renumberFaceOffsets(...) method")
656 static void renumberFaceOffsetsXml
657 (
658 labelUList& faceOffsets,
659 const label beginOffset
660 )
661 {
662 renumberFaceOffsets(faceOffsets, beginOffset, contentType::XML);
663 }
664
665 //- Copy vertex labels with a (global) point offset - legacy format
666 FOAM_DEPRECATED_FOR(2025-11, "copyVertLabels(...) method")
668 (
669 const labelUList& connectivity,
670 const label pointOffset
671 )
672 {
673 return
674 copyVertLabels(connectivity, pointOffset, contentType::LEGACY);
675 }
676
677 //- Copy vertex labels with a (global) point offset - XML format
678 FOAM_DEPRECATED_FOR(2025-11, "copyVertLabels(...) method")
680 (
681 const labelUList& connectivity,
682 const label pointOffset
683 )
684 {
685 return
686 copyVertLabels(connectivity, pointOffset, contentType::XML);
687 }
688
689 //- Copy faces stream labels with a (global) point offset - XML format
690 FOAM_DEPRECATED_FOR(2025-11, "copyFaceLabels(...) method")
692 (
693 const labelUList& faceLabels,
694 const label pointOffset
695 )
696 {
697 return
699 }
700
701 //- Copy face offsets with an offset from previous - XML format
702 FOAM_DEPRECATED_FOR(2025-11, "copyFaceOffsets(...) method")
704 (
705 const labelUList& faceOffsets,
706 const label beginOffset
707 )
709 return
710 copyFaceOffsets(faceOffsets, beginOffset, contentType::XML);
711 }
712
713 //- The calculated size for legacy storage
714 FOAM_DEPRECATED_FOR(2025-11, "sizeOf(LEGACY, CELLS)")
715 label sizeLegacy() const
716 {
718 }
719
720 //- The calculated size for legacy storage of the specified slot
721 FOAM_DEPRECATED_FOR(2025-11, "sizeOf(LEGACY, ...)")
722 label sizeLegacy(slotType slot) const
723 {
724 return sizeOf(contentType::LEGACY, slot);
725 }
726
727 //- The calculated size for xml storage of the specified slot
728 FOAM_DEPRECATED_FOR(2025-11, "sizeOf(XML, ...)")
729 label sizeXml(slotType slot) const
730 {
731 return sizeOf(contentType::XML, slot);
732 }
733
734 //- The calculated size for vtk-internal storage of the specified slot
735 FOAM_DEPRECATED_FOR(2025-11, "sizeOf(INTERNAL1, ...)")
736 label sizeInternal1(slotType slot) const
737 {
741 //- The calculated size for vtk-internal storage of the specified slot
742 FOAM_DEPRECATED_FOR(2025-11, "sizeOf(INTERNAL2, ...)")
743 label sizeInternal2(slotType slot) const
744 {
745 return sizeOf(contentType::INTERNAL2, slot);
746 }
747};
748
749
750// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
752} // End namespace vtk
753} // End namespace Foam
754
755// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
756
757#include "foamVtuSizingI.H"
758
759// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
760
761#endif
762
763// ************************************************************************* //
label n
label nFaceLabels
labelList faceLabels(nFaceLabels)
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
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
An analytical geometric cellShape.
Definition cellShape.H:71
Bookkeeping for mesh subsetting and/or polyhedral cell decomposition. Although the main use case is f...
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
Sizing descriptions and routines for transcribing an OpenFOAM volume mesh into a VTK unstructured gri...
static labelList dummyFaceOffsets(const label numCells, const enum contentType output, label beginOffset=0)
Return a list of dummy faceOffsets.
label nAddCells() const noexcept
Number of additional (decomposed) cells for the mesh.
void clear() noexcept
Reset all sizes to zero.
vtuSizing() noexcept
Default construct.
void resetShapes(const UList< cellShape > &shapes)
Reset sizing using primitive shapes only (ADVANCED USAGE).
void setNumPoints(label n) noexcept
Alter number of mesh points (ADVANCED USAGE).
static labelList copyFaceOffsetsXml(const labelUList &faceOffsets, const label beginOffset)
Copy face offsets with an offset from previous - XML format.
static void renumberFaceLabels(labelUList &faceLabels, const label pointOffset, const contentType output)
Renumber faces stream labels by (global) point offset.
label sizeInternal1(slotType slot) const
The calculated size for vtk-internal storage of the specified slot.
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).
label nCellsPoly() const noexcept
Number of output polyhedral cells for the mesh.
bool operator!=(const vtuSizing &rhs) const
Test inequality.
slotType
The possible storage 'slots' that can be used.
@ FACES
Polyhedral face-stream (XML, INTERNAL) or face vertices (HDF).
@ CELLS
Cell connectivity (ALL).
@ POLY_FACEIDS
Polyhedral face lookup (HDF only).
static labelList copyFaceLabelsXml(const labelUList &faceLabels, const label pointOffset)
Copy faces stream labels with a (global) point offset - XML format.
label nFacesPoly() const noexcept
Number of (non-unique) faces used by polyhedral cells.
void populateHdf(const polyMesh &mesh, UList< uint8_t > &cellTypes, labelUList &connectivity, labelUList &offsets, labelUList &faces, labelUList &facesOffsets, labelUList &polyFaceIds, labelUList &polyFaceOffsets, foamVtkMeshMaps &maps) const
Populate lists for HDF output.
selectionModeType selectionMode() const noexcept
Query how the mesh cells have been selected or defined.
static void renumberFaceOffsets(labelUList &faceOffsets, const label beginOffset, const contentType output)
Renumber face offsets with a specified (global) begin offset.
label nAddPoints() const noexcept
Number of additional (decomposed) points for the mesh.
bool decompose() const noexcept
Query the decompose flag (normally off).
bool manifold() const noexcept
Manifold mesh cells detected? Globally consistent quantity.
label sizeInternal2(slotType slot) const
The calculated size for vtk-internal storage of the specified slot.
label nVertPoly() const noexcept
Number of vertex labels for polyhedral cells of the mesh.
void populateShapesLegacy(const UList< cellShape > &shapes, UList< uint8_t > &cellTypes, labelUList &connectivity, foamVtkMeshMaps &maps) const
Reset list for primitive shapes only (ADVANCED USAGE).
contentType
Types of content that the storage may represent.
@ INTERNAL2
Internal vtkUnstructuredGrid content, VTK_CELL_ARRAY_V2.
@ HDF
Similar to INTERNAL2, but with VTKHDF specifics.
@ INTERNAL1
Internal vtkUnstructuredGrid content (old).
static labelList copyVertLabelsLegacy(const labelUList &connectivity, const label pointOffset)
Copy vertex labels with a (global) point offset - legacy format.
bool hasPolyCells() const noexcept
Has polyhedral cells?
static void renumberVertLabelsLegacy(labelUList &connectivity, const label pointOffset)
Renumber vertex labels by (global) point offset - legacy format.
label sizeOf(contentType output) const noexcept
Return the required size for the storage slot.
void reset(const polyMesh &mesh, const bool decompose=false)
Reset sizing by analyzing the mesh.
label nAddVerts() const noexcept
Number of additional (decomposed) vertices for the mesh.
void info(Ostream &os) const
Report some information.
label sizeXml(slotType slot) const
The calculated size for xml storage of the specified slot.
static void renumberFaceLabelsXml(labelUList &faceLabels, const label pointOffset)
Renumber faces stream labels by (global) point offset - XML format.
static labelList copyVertLabels(const labelUList &connectivity, const label pointOffset, const contentType output)
Copy vertex labels with a (global) point offset.
label nVertLabels() const noexcept
Number of vertex labels for the mesh.
static void renumberVertLabels(labelUList &connectivity, const label pointOffset, const contentType output)
Renumber vertex labels by (global) point offset.
label sizeOf(const enum contentType output, const enum slotType slot) const noexcept
Return the required size for the storage slot.
void setNumAddPoints(label n) noexcept
Alter number of additional (cell-centre) points (ADVANCED USAGE).
label nFieldPoints() const noexcept
Number of field points = nPoints + nAddPoints.
label nCells() const noexcept
Number of cells for the mesh.
bool operator==(const vtuSizing &rhs) const
Test equality.
label nFaceLabels() const noexcept
Number of polyhedral face labels for the mesh, without any prefixing information.
selectionModeType
How the mesh cells have been selected or defined.
@ SHAPE_MESH
A limited subset of primitive shapes.
@ SUBSET_MESH
A mesh subset, using all mesh points.
static labelList copyVertLabelsXml(const labelUList &connectivity, const label pointOffset)
Copy vertex labels with a (global) point offset - XML format.
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 populateLegacy(const polyMesh &mesh, UList< uint8_t > &cellTypes, labelUList &connectivity, foamVtkMeshMaps &maps) const
Populate lists for Legacy output.
static labelList copyFaceOffsets(const labelUList &faceOffsets, const label beginOffset, const contentType output)
Copy face offsets with a specified (global) begin offset.
static void renumberVertLabelsXml(labelUList &connectivity, const label pointOffset)
Renumber vertex labels by (global) point offset - XML format.
label nFieldCells() const noexcept
Number of field cells = nCells + nAddCells.
label sizeLegacy() const
The calculated size for legacy storage.
static labelList copyFaceLabels(const labelUList &faceLabels, const label pointOffset, const contentType output)
Copy faces stream labels with a (global) point offset.
static void renumberFaceOffsetsXml(labelUList &faceOffsets, const label beginOffset)
Renumber face offsets with a specified begin offset - XML format.
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
#define declarePopulateInternalMethod(Type)
label nPoints
Namespace for handling VTK output. Contains classes and functions for writing VTK file content.
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition List.H:62
void rhs(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
const direction noexcept
Definition scalarImpl.H:265
UList< label > labelUList
A UList of labels.
Definition UList.H:75
const labelUList & cellTypes
Definition setCellMask.H:27
#define FOAM_DEPRECATED_FOR(since, replacement)
Definition stdFoam.H:43