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-2021 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 
29 // OpenFOAM includes
30 #include "fvMesh.H"
31 #include "cellModeller.H"
32 #include "vtkOpenFOAMPoints.H"
33 #include "Swap.H"
34 #include "polygonTriangulate.H"
35 
36 // VTK includes
37 #include "vtkCellArray.h"
38 #include "vtkIdTypeArray.h"
39 #include "vtkUnstructuredGrid.h"
40 
41 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
42 
43 vtkUnstructuredGrid* Foam::vtkPVFoam::volumeVTKMesh
44 (
45  const fvMesh& mesh,
46  polyDecomp& decompInfo
47 )
48 {
49  const cellModel& tet = *(cellModeller::lookup("tet"));
50  const cellModel& pyr = *(cellModeller::lookup("pyr"));
51  const cellModel& prism = *(cellModeller::lookup("prism"));
52  const cellModel& wedge = *(cellModeller::lookup("wedge"));
53  const cellModel& tetWedge = *(cellModeller::lookup("tetWedge"));
54  const cellModel& hex = *(cellModeller::lookup("hex"));
55 
56  vtkUnstructuredGrid* vtkmesh = vtkUnstructuredGrid::New();
57 
58  if (debug)
59  {
61  printMemory();
62  }
63 
64  const cellShapeList& cellShapes = mesh.cellShapes();
65 
66  // Number of additional points needed by the decomposition of polyhedra
67  label nAddPoints = 0;
68 
69  // Number of additional cells generated by the decomposition of polyhedra
70  label nAddCells = 0;
71 
72  // face owner is needed to determine the face orientation
73  const labelList& owner = mesh.faceOwner();
74 
75  labelList& superCells = decompInfo.superCells();
76  labelList& addPointCellLabels = decompInfo.addPointCellLabels();
77 
78  // Scan for cells which need to be decomposed and count additional points
79  // and cells
80  if (!reader_->GetUseVTKPolyhedron())
81  {
82  if (debug)
83  {
85  << "... scanning for polyhedra" << endl;
86  }
87 
88  forAll(cellShapes, celli)
89  {
90  const cellModel& model = cellShapes[celli].model();
91 
92  if
93  (
94  model != hex
95  && model != wedge
96  && model != prism
97  && model != pyr
98  && model != tet
99  && model != tetWedge
100  )
101  {
102  const cell& cFaces = mesh.cells()[celli];
103 
104  forAll(cFaces, cFacei)
105  {
106  const face& f = mesh.faces()[cFaces[cFacei]];
107 
108  if (f.size() != 4)
109  {
110  nAddCells += f.nTriangles();
111  }
112  else
113  {
114  nAddCells ++;
115  }
116  }
117 
118  nAddCells--;
119  nAddPoints++;
120  }
121  }
122  }
123 
124  // Set size of additional point addressing array
125  // (from added point to original cell)
126  addPointCellLabels.setSize(nAddPoints);
127 
128  // Set size of additional cells mapping array
129  // (from added cell to original cell)
130 
131  if (debug)
132  {
133  Info<< " mesh nCells = " << mesh.nCells() << nl
134  << " nPoints = " << mesh.nPoints() << nl
135  << " nAddCells = " << nAddCells << nl
136  << " nAddPoints = " << nAddPoints << endl;
137  }
138 
139  superCells.setSize(mesh.nCells() + nAddCells);
140 
141  if (debug)
142  {
143  Info<< " ... converting points" << endl;
144  }
145 
146  // Convert OpenFOAM mesh vertices to VTK
147  vtkPoints* vtkpoints = vtkPoints::New();
148  vtkpoints->Allocate(mesh.nPoints() + nAddPoints);
149 
150  const Foam::pointField& points = mesh.points();
151 
152  forAll(points, i)
153  {
154  vtkInsertNextOpenFOAMPoint(vtkpoints, points[i]);
155  }
156 
157 
158  if (debug)
159  {
160  Info<< " ... converting cells" << endl;
161  }
162 
163  vtkmesh->Allocate(mesh.nCells() + nAddCells);
164 
165  // Create a triangulation engine
166  polygonTriangulate triEngine;
167 
168  // Set counters for additional points and additional cells
169  label addPointi = 0, addCelli = 0;
170 
171  // Create storage for points - needed for mapping from OpenFOAM to VTK
172  // data types - max 'order' = hex = 8 points
173  vtkIdType nodeIds[8];
174 
175  // face-stream for a polyhedral cell
176  // [numFace0Pts, id1, id2, id3, numFace1Pts, id1, id2, id3, ...]
177  DynamicList<vtkIdType> faceStream(256);
178 
179  forAll(cellShapes, celli)
180  {
181  const cellShape& cellShape = cellShapes[celli];
182  const cellModel& cellModel = cellShape.model();
183 
184  superCells[addCelli++] = celli;
185 
186  if (cellModel == tet)
187  {
188  for (int j = 0; j < 4; j++)
189  {
190  nodeIds[j] = cellShape[j];
191  }
192  vtkmesh->InsertNextCell
193  (
194  VTK_TETRA,
195  4,
196  nodeIds
197  );
198  }
199  else if (cellModel == pyr)
200  {
201  for (int j = 0; j < 5; j++)
202  {
203  nodeIds[j] = cellShape[j];
204  }
205  vtkmesh->InsertNextCell
206  (
207  VTK_PYRAMID,
208  5,
209  nodeIds
210  );
211  }
212  else if (cellModel == prism)
213  {
214  // VTK has a different node order for VTK_WEDGE
215  // their triangles point outwards!
216  nodeIds[0] = cellShape[0];
217  nodeIds[1] = cellShape[2];
218  nodeIds[2] = cellShape[1];
219  nodeIds[3] = cellShape[3];
220  nodeIds[4] = cellShape[5];
221  nodeIds[5] = cellShape[4];
222 
223  vtkmesh->InsertNextCell
224  (
225  VTK_WEDGE,
226  6,
227  nodeIds
228  );
229  }
230  else if (cellModel == tetWedge && !reader_->GetUseVTKPolyhedron())
231  {
232  // Treat as squeezed prism (VTK_WEDGE)
233 
234  nodeIds[0] = cellShape[0];
235  nodeIds[1] = cellShape[2];
236  nodeIds[2] = cellShape[1];
237  nodeIds[3] = cellShape[3];
238  nodeIds[4] = cellShape[4];
239  nodeIds[5] = cellShape[3];
240 
241  vtkmesh->InsertNextCell
242  (
243  VTK_WEDGE,
244  6,
245  nodeIds
246  );
247  }
248  else if (cellModel == wedge)
249  {
250  // Treat as squeezed hex
251 
252  nodeIds[0] = cellShape[0];
253  nodeIds[1] = cellShape[1];
254  nodeIds[2] = cellShape[2];
255  nodeIds[3] = cellShape[2];
256  nodeIds[4] = cellShape[3];
257  nodeIds[5] = cellShape[4];
258  nodeIds[6] = cellShape[5];
259  nodeIds[7] = cellShape[6];
260 
261  vtkmesh->InsertNextCell
262  (
263  VTK_HEXAHEDRON,
264  8,
265  nodeIds
266  );
267  }
268  else if (cellModel == hex)
269  {
270  for (int j = 0; j < 8; j++)
271  {
272  nodeIds[j] = cellShape[j];
273  }
274  vtkmesh->InsertNextCell
275  (
276  VTK_HEXAHEDRON,
277  8,
278  nodeIds
279  );
280  }
281  else if (reader_->GetUseVTKPolyhedron())
282  {
283  // Polyhedral cell - use VTK_POLYHEDRON
284  const labelList& cFaces = mesh.cells()[celli];
285 
286 #ifdef HAS_VTK_POLYHEDRON
287  vtkIdType nFaces = cFaces.size();
288  vtkIdType nLabels = nFaces;
289 
290  // count size for face stream
291  forAll(cFaces, cFacei)
292  {
293  const face& f = mesh.faces()[cFaces[cFacei]];
294  nLabels += f.size();
295  }
296 
297  // build face-stream
298  // [numFace0Pts, id1, id2, id3, numFace1Pts, id1, id2, id3, ...]
299  // point Ids are global
300  faceStream.clear();
301  faceStream.reserve(nLabels + nFaces);
302 
303  forAll(cFaces, cFacei)
304  {
305  const face& f = mesh.faces()[cFaces[cFacei]];
306  const bool isOwner = (owner[cFaces[cFacei]] == celli);
307  const label nFacePoints = f.size();
308 
309  // number of labels for this face
310  faceStream.append(nFacePoints);
311 
312  if (isOwner)
313  {
314  forAll(f, fp)
315  {
316  faceStream.append(f[fp]);
317  }
318  }
319  else
320  {
321  // Fairly immaterial if we reverse the list
322  // or use face::reverseFace()
323  forAllReverse(f, fp)
324  {
325  faceStream.append(f[fp]);
326  }
327  }
328  }
329 
330  vtkmesh->InsertNextCell(VTK_POLYHEDRON, nFaces, faceStream.data());
331 #else
332  // This is a horrible substitute
333  // but avoids crashes when there is no vtkPolyhedron support
334 
335  // Establish unique node ids used
336  HashSet<vtkIdType, Hash<label>> hashUniqId(2*256);
337 
338  forAll(cFaces, cFacei)
339  {
340  const face& f = mesh.faces()[cFaces[cFacei]];
341 
342  forAll(f, fp)
343  {
344  hashUniqId.insert(f[fp]);
345  }
346  }
347 
348  // Use face stream to store unique node ids:
349  faceStream = hashUniqId.sortedToc();
350 
351  vtkmesh->InsertNextCell
352  (
353  VTK_CONVEX_POINT_SET,
354  vtkIdType(faceStream.size()),
355  faceStream.data()
356  );
357 #endif
358  }
359  else
360  {
361  // Polyhedral cell. Decompose into tets + prisms.
362 
363  // Mapping from additional point to cell
364  addPointCellLabels[addPointi] = celli;
365 
366  // The new vertex from the cell-centre
367  const label newVertexLabel = mesh.nPoints() + addPointi;
368  vtkInsertNextOpenFOAMPoint(vtkpoints, mesh.C()[celli]);
369 
370  // Whether to insert cell in place of original or not.
371  bool substituteCell = true;
372 
373  const labelList& cFaces = mesh.cells()[celli];
374  forAll(cFaces, cFacei)
375  {
376  const face& f = mesh.faces()[cFaces[cFacei]];
377  const bool isOwner = (owner[cFaces[cFacei]] == celli);
378 
379  // Do decomposition into triFcs and quadFcs.
380  faceList triFcs(f.size() == 4 ? 0 : f.nTriangles());
381  faceList quadFcs(f.size() == 4 ? 1 : 0);
382  if (f.size() != 4)
383  {
384  triEngine.triangulate
385  (
386  UIndirectList<point>(mesh.points(), f)
387  );
388  forAll(triEngine.triPoints(), trii)
389  {
390  triFcs[trii] = triEngine.triPoints(trii, f);
391  }
392  }
393  else
394  {
395  quadFcs[0] = f;
396  }
397 
398  forAll(quadFcs, quadI)
399  {
400  if (substituteCell)
401  {
402  substituteCell = false;
403  }
404  else
405  {
406  superCells[addCelli++] = celli;
407  }
408 
409  const face& quad = quadFcs[quadI];
410 
411  // Ensure we have the correct orientation for the
412  // base of the primitive cell shape.
413  // If the cell is face owner, the orientation needs to be
414  // flipped.
415  // At the moment, VTK doesn't actually seem to care if
416  // negative cells are defined, but we'll do it anyhow
417  // (for safety).
418  if (isOwner)
419  {
420  nodeIds[0] = quad[3];
421  nodeIds[1] = quad[2];
422  nodeIds[2] = quad[1];
423  nodeIds[3] = quad[0];
424  }
425  else
426  {
427  nodeIds[0] = quad[0];
428  nodeIds[1] = quad[1];
429  nodeIds[2] = quad[2];
430  nodeIds[3] = quad[3];
431  }
432  nodeIds[4] = newVertexLabel;
433  vtkmesh->InsertNextCell
434  (
435  VTK_PYRAMID,
436  5,
437  nodeIds
438  );
439  }
440 
441  forAll(triFcs, triI)
442  {
443  if (substituteCell)
444  {
445  substituteCell = false;
446  }
447  else
448  {
449  superCells[addCelli++] = celli;
450  }
451 
452  const face& tri = triFcs[triI];
453 
454  // See note above about the orientation.
455  if (isOwner)
456  {
457  nodeIds[0] = tri[2];
458  nodeIds[1] = tri[1];
459  nodeIds[2] = tri[0];
460  }
461  else
462  {
463  nodeIds[0] = tri[0];
464  nodeIds[1] = tri[1];
465  nodeIds[2] = tri[2];
466  }
467  nodeIds[3] = newVertexLabel;
468 
469  vtkmesh->InsertNextCell
470  (
471  VTK_TETRA,
472  4,
473  nodeIds
474  );
475  }
476  }
477 
478  addPointi++;
479  }
480  }
481 
482  vtkmesh->SetPoints(vtkpoints);
483  vtkpoints->Delete();
484 
485  if (debug)
486  {
487  printMemory();
488  }
489 
490  return vtkmesh;
491 }
492 
493 
494 // ************************************************************************* //
Swap its arguments.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: UList.H:446
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
static void printMemory()
Simple memory used debugging information.
const cellShapeList & cellShapes
const pointField & points
#define InfoInFunction
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:251
IOstream & hex(IOstream &io)
Definition: IOstream.H:561
messageStream Info
List< face > faceList
Definition: faceListFwd.H:43
static const char nl
Definition: Ostream.H:260
labelList f(nPoints)
void vtkInsertNextOpenFOAMPoint(vtkPoints *points, const Foam::point &p)