vtkPVFoamMeshVolume.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "vtkPVFoam.H"
27 #include "vtkPVFoamReader.h"
28 #include "vtkOpenFOAMPoints.H"
29 
30 // OpenFOAM includes
31 #include "volFields.H"
32 #include "cellModeller.H"
33 #include "polygonTriangulate.H"
34 
35 // VTK includes
36 #include "vtkCellArray.h"
37 #include "vtkIdTypeArray.h"
38 #include "vtkUnstructuredGrid.h"
39 
40 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
41 
42 vtkUnstructuredGrid* Foam::vtkPVFoam::volumeVTKMesh
43 (
44  const fvMesh& mesh,
45  polyDecomp& decompInfo
46 )
47 {
48  const cellModel& tet = *(cellModeller::lookup("tet"));
49  const cellModel& pyr = *(cellModeller::lookup("pyr"));
50  const cellModel& prism = *(cellModeller::lookup("prism"));
51  const cellModel& wedge = *(cellModeller::lookup("wedge"));
52  const cellModel& tetWedge = *(cellModeller::lookup("tetWedge"));
53  const cellModel& hex = *(cellModeller::lookup("hex"));
54 
55  vtkUnstructuredGrid* vtkmesh = vtkUnstructuredGrid::New();
56 
58 
59  // Number of additional points needed by the decomposition of polyhedra
60  label nAddPoints = 0;
61 
62  // Number of additional cells generated by the decomposition of polyhedra
63  label nAddCells = 0;
64 
65  // face owner is needed to determine the face orientation
66  const labelList& owner = mesh.faceOwner();
67 
68  labelList& superCells = decompInfo.superCells();
69  labelList& addPointCellLabels = decompInfo.addPointCellLabels();
70 
71  // Scan for cells which need to be decomposed and count additional points
72  // and cells
73  if (!reader_->GetUseVTKPolyhedron())
74  {
75  forAll(cellShapes, celli)
76  {
77  const cellModel& model = cellShapes[celli].model();
78 
79  if
80  (
81  model != hex
82  && model != wedge
83  && model != prism
84  && model != pyr
85  && model != tet
86  && model != tetWedge
87  )
88  {
89  const cell& cFaces = mesh.cells()[celli];
90 
91  forAll(cFaces, cFacei)
92  {
93  const face& f = mesh.faces()[cFaces[cFacei]];
94 
95  if (f.size() != 4)
96  {
97  nAddCells += f.nTriangles();
98  }
99  else
100  {
101  nAddCells ++;
102  }
103  }
104 
105  nAddCells--;
106  nAddPoints++;
107  }
108  }
109  }
110 
111  // Set size of additional point addressing array
112  // (from added point to original cell)
113  addPointCellLabels.setSize(nAddPoints);
114 
115  // Set size of additional cells mapping array
116  // (from added cell to original cell)
117  superCells.setSize(mesh.nCells() + nAddCells);
118 
120  << "nCells=" << mesh.nCells()
121  << ", nPoints=" << mesh.nPoints()
122  << ", nAddCells=" << nAddCells
123  << ", nAddPoints=" << nAddPoints << endl;
124 
125  // Convert OpenFOAM mesh vertices to VTK
126  vtkPoints* vtkpoints = vtkPoints::New();
127  vtkpoints->Allocate(mesh.nPoints() + nAddPoints);
128 
129  const Foam::pointField& points = mesh.points();
130 
131  forAll(points, i)
132  {
133  vtkInsertNextOpenFOAMPoint(vtkpoints, points[i]);
134  }
135 
136  vtkmesh->Allocate(mesh.nCells() + nAddCells);
137 
138  // Create a triangulation engine
139  polygonTriangulate triEngine;
140 
141  // Set counters for additional points and additional cells
142  label addPointi = 0, addCelli = 0;
143 
144  // Create storage for points - needed for mapping from OpenFOAM to VTK
145  // data types - max 'order' = hex = 8 points
146  vtkIdType nodeIds[8];
147 
148  // face-stream for a polyhedral cell
149  // [numFace0Pts, id1, id2, id3, numFace1Pts, id1, id2, id3, ...]
150  DynamicList<vtkIdType> faceStream(256);
151 
152  forAll(cellShapes, celli)
153  {
154  const cellShape& cellShape = cellShapes[celli];
155  const cellModel& cellModel = cellShape.model();
156 
157  superCells[addCelli++] = celli;
158 
159  if (cellModel == tet)
160  {
161  for (int j = 0; j < 4; j++)
162  {
163  nodeIds[j] = cellShape[j];
164  }
165  vtkmesh->InsertNextCell
166  (
167  VTK_TETRA,
168  4,
169  nodeIds
170  );
171  }
172  else if (cellModel == pyr)
173  {
174  for (int j = 0; j < 5; j++)
175  {
176  nodeIds[j] = cellShape[j];
177  }
178  vtkmesh->InsertNextCell
179  (
180  VTK_PYRAMID,
181  5,
182  nodeIds
183  );
184  }
185  else if (cellModel == prism)
186  {
187  // VTK has a different node order for VTK_WEDGE
188  // their triangles point outwards!
189  nodeIds[0] = cellShape[0];
190  nodeIds[1] = cellShape[2];
191  nodeIds[2] = cellShape[1];
192  nodeIds[3] = cellShape[3];
193  nodeIds[4] = cellShape[5];
194  nodeIds[5] = cellShape[4];
195 
196  vtkmesh->InsertNextCell
197  (
198  VTK_WEDGE,
199  6,
200  nodeIds
201  );
202  }
203  else if (cellModel == tetWedge && !reader_->GetUseVTKPolyhedron())
204  {
205  // Treat as squeezed prism (VTK_WEDGE)
206 
207  nodeIds[0] = cellShape[0];
208  nodeIds[1] = cellShape[2];
209  nodeIds[2] = cellShape[1];
210  nodeIds[3] = cellShape[3];
211  nodeIds[4] = cellShape[4];
212  nodeIds[5] = cellShape[3];
213 
214  vtkmesh->InsertNextCell
215  (
216  VTK_WEDGE,
217  6,
218  nodeIds
219  );
220  }
221  else if (cellModel == wedge)
222  {
223  // Treat as squeezed hex
224 
225  nodeIds[0] = cellShape[0];
226  nodeIds[1] = cellShape[1];
227  nodeIds[2] = cellShape[2];
228  nodeIds[3] = cellShape[2];
229  nodeIds[4] = cellShape[3];
230  nodeIds[5] = cellShape[4];
231  nodeIds[6] = cellShape[5];
232  nodeIds[7] = cellShape[6];
233 
234  vtkmesh->InsertNextCell
235  (
236  VTK_HEXAHEDRON,
237  8,
238  nodeIds
239  );
240  }
241  else if (cellModel == hex)
242  {
243  for (int j = 0; j < 8; j++)
244  {
245  nodeIds[j] = cellShape[j];
246  }
247  vtkmesh->InsertNextCell
248  (
249  VTK_HEXAHEDRON,
250  8,
251  nodeIds
252  );
253  }
254  else if (reader_->GetUseVTKPolyhedron())
255  {
256  // Polyhedral cell - use VTK_POLYHEDRON
257  const labelList& cFaces = mesh.cells()[celli];
258 
259 #ifdef HAS_VTK_POLYHEDRON
260  vtkIdType nFaces = cFaces.size();
261  vtkIdType nLabels = nFaces;
262 
263  // count size for face stream
264  forAll(cFaces, cFacei)
265  {
266  const face& f = mesh.faces()[cFaces[cFacei]];
267  nLabels += f.size();
268  }
269 
270  // build face-stream
271  // [numFace0Pts, id1, id2, id3, numFace1Pts, id1, id2, id3, ...]
272  // point Ids are global
273  faceStream.clear();
274  faceStream.reserve(nLabels + nFaces);
275 
276  forAll(cFaces, cFacei)
277  {
278  const face& f = mesh.faces()[cFaces[cFacei]];
279  const bool isOwner = (owner[cFaces[cFacei]] == celli);
280  const label nFacePoints = f.size();
281 
282  // number of labels for this face
283  faceStream.append(nFacePoints);
284 
285  if (isOwner)
286  {
287  forAll(f, fp)
288  {
289  faceStream.append(f[fp]);
290  }
291  }
292  else
293  {
294  // Fairly immaterial if we reverse the list
295  // or use face::reverseFace()
296  forAllReverse(f, fp)
297  {
298  faceStream.append(f[fp]);
299  }
300  }
301  }
302 
303  vtkmesh->InsertNextCell(VTK_POLYHEDRON, nFaces, faceStream.data());
304 #else
305  // This is a horrible substitute
306  // but avoids crashes when there is no vtkPolyhedron support
307 
308  // Establish unique node ids used
309  HashSet<vtkIdType, Hash<label>> hashUniqId(2*256);
310 
311  forAll(cFaces, cFacei)
312  {
313  const face& f = mesh.faces()[cFaces[cFacei]];
314 
315  forAll(f, fp)
316  {
317  hashUniqId.insert(f[fp]);
318  }
319  }
320 
321  // Use face stream to store unique node ids:
322  faceStream = hashUniqId.sortedToc();
323 
324  vtkmesh->InsertNextCell
325  (
326  VTK_CONVEX_POINT_SET,
327  vtkIdType(faceStream.size()),
328  faceStream.data()
329  );
330 #endif
331  }
332  else
333  {
334  // Polyhedral cell. Decompose into tets + prisms.
335 
336  // Mapping from additional point to cell
337  addPointCellLabels[addPointi] = celli;
338 
339  // Create a new vertex from the cell-centre. If the mesh doesn't
340  // have cell-centres stored then compute the centre on the fly.
341  // These polyhedra should (hopefully) only be a small minority of
342  // the cells, and computing cell-centres throughout the mesh is
343  // quite expensive and involves a lot of memory (it also creates
344  // face areas and centres and cell volumes).
345  const label newVertexLabel = mesh.nPoints() + addPointi;
347  (
348  vtkpoints,
350  ? mesh.cellCentres()[celli]
351  : mesh.cells()[celli].centre(mesh.points(), mesh.faces())
352  );
353 
354  // Whether to insert cell in place of original or not.
355  bool substituteCell = true;
356 
357  const labelList& cFaces = mesh.cells()[celli];
358  forAll(cFaces, cFacei)
359  {
360  const face& f = mesh.faces()[cFaces[cFacei]];
361  const bool isOwner = (owner[cFaces[cFacei]] == celli);
362 
363  // Do decomposition into triFcs and quadFcs.
364  faceList triFcs(f.size() == 4 ? 0 : f.nTriangles());
365  faceList quadFcs(f.size() == 4 ? 1 : 0);
366  if (f.size() != 4)
367  {
368  triEngine.triangulate
369  (
370  UIndirectList<point>(mesh.points(), f)
371  );
372  forAll(triEngine.triPoints(), trii)
373  {
374  triFcs[trii] = triEngine.triPoints(trii, f);
375  }
376  }
377  else
378  {
379  quadFcs[0] = f;
380  }
381 
382  forAll(quadFcs, quadI)
383  {
384  if (substituteCell)
385  {
386  substituteCell = false;
387  }
388  else
389  {
390  superCells[addCelli++] = celli;
391  }
392 
393  const face& quad = quadFcs[quadI];
394 
395  // Ensure we have the correct orientation for the
396  // base of the primitive cell shape.
397  // If the cell is face owner, the orientation needs to be
398  // flipped.
399  // At the moment, VTK doesn't actually seem to care if
400  // negative cells are defined, but we'll do it anyhow
401  // (for safety).
402  if (isOwner)
403  {
404  nodeIds[0] = quad[3];
405  nodeIds[1] = quad[2];
406  nodeIds[2] = quad[1];
407  nodeIds[3] = quad[0];
408  }
409  else
410  {
411  nodeIds[0] = quad[0];
412  nodeIds[1] = quad[1];
413  nodeIds[2] = quad[2];
414  nodeIds[3] = quad[3];
415  }
416  nodeIds[4] = newVertexLabel;
417  vtkmesh->InsertNextCell
418  (
419  VTK_PYRAMID,
420  5,
421  nodeIds
422  );
423  }
424 
425  forAll(triFcs, triI)
426  {
427  if (substituteCell)
428  {
429  substituteCell = false;
430  }
431  else
432  {
433  superCells[addCelli++] = celli;
434  }
435 
436  const face& tri = triFcs[triI];
437 
438  // See note above about the orientation.
439  if (isOwner)
440  {
441  nodeIds[0] = tri[2];
442  nodeIds[1] = tri[1];
443  nodeIds[2] = tri[0];
444  }
445  else
446  {
447  nodeIds[0] = tri[0];
448  nodeIds[1] = tri[1];
449  nodeIds[2] = tri[2];
450  }
451  nodeIds[3] = newVertexLabel;
452 
453  vtkmesh->InsertNextCell
454  (
455  VTK_TETRA,
456  4,
457  nodeIds
458  );
459  }
460  }
461 
462  addPointi++;
463  }
464  }
465 
466  vtkmesh->SetPoints(vtkpoints);
467  vtkpoints->Delete();
468 
469  return vtkmesh;
470 }
471 
472 
473 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: UList.H:445
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
static const cellModel * lookup(const word &)
Look up a model by name and return a pointer to the model or nullptr.
Definition: cellModeller.C:100
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1344
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1357
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1331
const vectorField & cellCentres() const
const cellShapeList & cellShapes() const
Return cell shapes.
label nCells() const
bool has(const Type &(primitiveMesh::*method)() const) const
Return whether the result of the given method has been allocated.
label nPoints() const
const cellList & cells() const
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
const cellShapeList & cellShapes
const pointField & points
#define DebugInFunction
Report an information message using Foam::Info.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
List< cellShape > cellShapeList
List of cellShapes and PtrList of List of cellShape.
Definition: cellShapeList.H:43
List< label > labelList
A List of labels.
Definition: labelList.H:56
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
IOstream & hex(IOstream &io)
Definition: IOstream.H:576
List< face > faceList
Definition: faceListFwd.H:41
labelList f(nPoints)
void vtkInsertNextOpenFOAMPoint(vtkPoints *points, const Foam::point &p)