Loading...
Searching...
No Matches
foamVtkVtuAdaptorI.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
26\*---------------------------------------------------------------------------*/
27
28// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
29
30inline vtkSmartPointer<vtkPoints>
32(
33 vtkIdType count
34) const
35{
36 auto vtkpoints = vtkSmartPointer<vtkPoints>::New();
37 // if (narrowing || (sizeof(float) == sizeof(Foam::scalar)))
38 // {
39 // vtkpoints->SetDataTypeToFloat();
40 // }
41 // else
42 // {
43 // vtkpoints->SetDataTypeToDouble();
44 // }
45 vtkpoints->SetNumberOfPoints(count);
46
47 return vtkpoints;
48}
49
50
51inline vtkSmartPointer<vtkPoints>
53(
54 const polyMesh& mesh
55) const
56{
57 // Convert OpenFOAM mesh vertices to VTK
58
59 // Combine mesh points and any additional cellCentre points
60 // into a single field
61
62 const auto& pts = mesh.points();
63 const auto& cc = mesh.cellCentres();
64
65 // The additional cellCentre points
66 const labelUList& addPoints = this->additionalIds();
67
68 auto vtkpoints = NewPoints(pts.size() + addPoints.size());
69
70 // Normal points
71 vtkIdType pointId = 0;
72 for (const point& p : pts)
73 {
74 vtkpoints->SetPoint(pointId++, p.cdata());
75 }
76
77 // Cell centres
78 for (const label celli : addPoints)
79 {
80 const auto& p = cc[celli];
81 vtkpoints->SetPoint(pointId++, p.cdata());
82 }
83
84 return vtkpoints;
85}
86
87
88inline vtkSmartPointer<vtkPoints>
90(
91 const polyMesh& mesh,
92 const labelUList& pointMap
93) const
94{
95 // Convert OpenFOAM mesh vertices to VTK
96
97 // Combine mapped points and additional cellCentre points
98 // into a single field
99
100 const auto& pts = mesh.points();
101 const auto& cc = mesh.cellCentres();
102
103 // The additional cellCentre points
104 const labelUList& addPoints = this->additionalIds();
105
106 auto vtkpoints = NewPoints(pointMap.size() + addPoints.size());
107
108 // Normal points
109 vtkIdType pointId = 0;
110 for (const label pointi : pointMap)
111 {
112 const auto& p = pts[pointi];
113 vtkpoints->SetPoint(pointId++, p.cdata());
114 }
115
116 // Cell centres
117 for (const label celli : addPoints)
118 {
119 const auto& p = cc[celli];
120 vtkpoints->SetPoint(pointId++, p.cdata());
122
123 return vtkpoints;
124}
125
126
127inline vtkSmartPointer<vtkUnstructuredGrid>
129(
130 const fvMesh& mesh,
131 const bool decompPoly
132)
133{
134 const vtk::vtuSizing::contentType output
135 (
136 #ifdef VTK_CELL_ARRAY_V2
138 #else
140 #endif
141 );
142
143 vtk::vtuSizing sizing(mesh, decompPoly);
144
145 auto vtkmesh = vtkSmartPointer<vtkUnstructuredGrid>::New();
146
147 auto cellTypes = vtkSmartPointer<vtkUnsignedCharArray>::New();
148
149 UList<uint8_t> cellTypesUL
150 (
152 );
153
154 auto cells = vtkSmartPointer<vtkCellArray>::New();
155 auto faces = vtkSmartPointer<vtkIdTypeArray>::New();
156
157 auto cellOffsets = vtkSmartPointer<vtkIdTypeArray>::New();
158 auto faceLocations = vtkSmartPointer<vtkIdTypeArray>::New();
159
160 const auto nConnect
161 (
163 );
164
165 UList<vtkIdType> cellOffsetsUL
166 (
168 (
169 cellOffsets,
171 )
172 );
173
174 #ifdef VTK_CELL_ARRAY_V2
175
176 auto cellConnect = vtkSmartPointer<vtkIdTypeArray>::New();
177
178 UList<vtkIdType> cellsUL
179 (
180 vtk::Tools::asUList(cellConnect, nConnect)
181 );
182
183 #else
184
185 cells->GetData()->SetNumberOfTuples(sizing.nFieldCells());
186
187 UList<vtkIdType> cellsUL
188 (
189 cells->WritePointer(sizing.nFieldCells(), nConnect),
190 nConnect
191 );
192
193 #endif
194
195 UList<vtkIdType> facesUL
196 (
198 (
199 faces,
201 )
202 );
203
204 UList<vtkIdType> faceLocationsUL
205 (
207 (
208 faceLocations,
210 )
211 );
212
213
214 sizing.populateInternal
215 (
216 mesh,
217 cellTypesUL,
218 cellsUL, cellOffsetsUL,
219 facesUL, faceLocationsUL,
220 static_cast<foamVtkMeshMaps&>(*this),
221 output
222 );
223
224
225 // Convert OpenFOAM mesh vertices to VTK
226 // - can only do this *after* populating the decompInfo with cell-ids
227 // for any additional points (ie, mesh cell-centres)
228 vtkmesh->SetPoints(this->points(mesh));
229
230 #ifdef VTK_CELL_ARRAY_V2
231
232 // Move into a vtkCellArray
233 cells->SetData(cellOffsets, cellConnect);
234
235 // FUTURE: after VTK-9.4.0
236 // (VTK_MAJOR_VERSION >= 9) && (VTK_NINOR_VERSION >= 4)
237 // SetPolyhedralCells()
238
239 if (facesUL.size())
240 {
241 vtkmesh->SetCells(cellTypes, cells, faceLocations, faces);
242 }
243 else
244 {
245 vtkmesh->SetCells(cellTypes, cells, nullptr, nullptr);
246 }
247 #else
248
249 if (facesUL.size())
250 {
251 vtkmesh->SetCells(cellTypes, cellOffsets, cells, faceLocations, faces);
252 }
253 else
254 {
255 vtkmesh->SetCells(cellTypes, cellOffsets, cells, nullptr, nullptr);
256 }
257
258 #endif
259
260 return vtkmesh;
261}
262
263
264// ************************************************************************* //
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
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Bookkeeping for mesh subsetting and/or polyhedral cell decomposition. Although the main use case is f...
const labelList & additionalIds() const noexcept
Any additional (user) labels.
foamVtkMeshMaps()=default
Default construct: zero-sized, no reserved sizes.
const labelList & pointMap() const noexcept
Point labels for subsetted meshes.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
vtkSmartPointer< vtkPoints > NewPoints(vtkIdType count) const
Allocate vtkPoints.
vtkSmartPointer< vtkUnstructuredGrid > internal(const fvMesh &mesh, const bool decompPoly=false)
Internal mesh as vtkUnstructuredGrid.
vtkSmartPointer< vtkPoints > points(const polyMesh &mesh) const
The vtk points for the mesh (and decomposition).
Sizing descriptions and routines for transcribing an OpenFOAM volume mesh into a VTK unstructured gri...
@ FACES
Polyhedral face-stream (XML, INTERNAL) or face vertices (HDF).
@ CELLS
Cell connectivity (ALL).
contentType
Types of content that the storage may represent.
@ INTERNAL2
Internal vtkUnstructuredGrid content, VTK_CELL_ARRAY_V2.
@ INTERNAL1
Internal vtkUnstructuredGrid content (old).
label sizeOf(const enum contentType output, const enum slotType slot) const noexcept
Return the required size for the storage slot.
void populateInternal(const polyMesh &mesh, UList< uint8_t > &cellTypes, UList< int > &connectivity, UList< int > &offsets, UList< int > &faces, UList< int > &facesOffsets, foamVtkMeshMaps &maps, const enum contentType output) const
Populate lists for Internal VTK format.
label nFieldCells() const noexcept
Number of field cells = nCells + nAddCells.
volScalarField & p
dynamicFvMesh & mesh
const pointField & points
const cellShapeList & cells
UList< uint8_t > asUList(vtkUnsignedCharArray *array, const label size)
Wrap vtkUnsignedCharArray as a UList.
vector point
Point is a vector.
Definition point.H:37
UList< label > labelUList
A UList of labels.
Definition UList.H:75
const pointField & pts
const labelUList & cellTypes
Definition setCellMask.H:27