vtkTopo.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 Description
25  Note: bug in vtk displaying wedges? Seems to display ok if we decompose
26  them. Should be thoroughly tested!
27  (they appear rarely in polyhedral meshes, do appear in some cut meshes)
28 
29 \*---------------------------------------------------------------------------*/
30 
31 #include "vtkTopo.H"
32 #include "polyMesh.H"
33 #include "cellShape.H"
34 #include "cellModeller.H"
35 #include "Swap.H"
36 #include "polygonTriangulate.H"
37 
38 // * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * * * //
39 
42 {"none", "polyhedra", "all"};
43 
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
47 Foam::vtkTopo::vtkTopo(const polyMesh& mesh, const vtkPolyhedra& polyhedra)
48 :
49  mesh_(mesh),
50  polyhedra_(polyhedra),
51  vertLabels_(),
52  cellTypes_(),
53  addPointCellLabels_(),
54  superCells_()
55 {
56  const cellModel& tet = *(cellModeller::lookup("tet"));
57  const cellModel& pyr = *(cellModeller::lookup("pyr"));
58  const cellModel& prism = *(cellModeller::lookup("prism"));
59  const cellModel& wedge = *(cellModeller::lookup("wedge"));
60  const cellModel& tetWedge = *(cellModeller::lookup("tetWedge"));
61  const cellModel& hex = *(cellModeller::lookup("hex"));
62 
63  const cellShapeList& cellShapes = mesh_.cellShapes();
64 
65  // Number of additional points and cells needed by the decomposition of
66  // polyhedra
67  label nAddPoints = 0, nAddCells = 0;
68 
69  // Scan for cells which need to be decomposed and count additional points
70  // and cells
71  if (polyhedra_ == vtkPolyhedra::none)
72  {
73  forAll(cellShapes, celli)
74  {
75  const cellModel& model = cellShapes[celli].model();
76 
77  if
78  (
79  model != hex
80  && model != wedge // See above.
81  && model != prism
82  && model != pyr
83  && model != tet
84  && model != tetWedge
85  )
86  {
87  const cell& cFaces = mesh_.cells()[celli];
88 
89  forAll(cFaces, cFacei)
90  {
91  const face& f = mesh_.faces()[cFaces[cFacei]];
92 
93  if (f.size() != 4)
94  {
95  nAddCells += f.nTriangles();
96  }
97  else
98  {
99  nAddCells ++;
100  }
101  }
102 
103  nAddCells--;
104  nAddPoints++;
105  }
106  }
107  }
108 
109  // Set size of additional point addressing array
110  // (from added point to original cell)
111  addPointCellLabels_.setSize(nAddPoints);
112 
113  // Set size of additional cells mapping array
114  // (from added cell to original cell)
115  superCells_.setSize(nAddCells);
116 
117  // List of vertex labels in VTK ordering
118  vertLabels_.setSize(cellShapes.size() + nAddCells);
119 
120  // Label of vtk type
121  cellTypes_.setSize(cellShapes.size() + nAddCells);
122 
123  // Create a triangulation engine
124  polygonTriangulate triEngine;
125 
126  // Construct vtk cells and types
127  label addPointi = 0, addCelli = 0;
128  const labelList& faceOwner = mesh.faceOwner();
129  forAll(cellShapes, celli)
130  {
131  const cellShape& cellShape = cellShapes[celli];
132  const cellModel& cellModel = cellShape.model();
133 
134  labelList& vtkVerts = vertLabels_[celli];
135 
136  if (polyhedra_ != vtkPolyhedra::all && cellModel == tet)
137  {
138  vtkVerts = cellShape;
139 
140  cellTypes_[celli] = VTK_TETRA;
141  }
142  else if (polyhedra_ != vtkPolyhedra::all && cellModel == pyr)
143  {
144  vtkVerts = cellShape;
145 
146  cellTypes_[celli] = VTK_PYRAMID;
147  }
148  else if (polyhedra_ != vtkPolyhedra::all && cellModel == prism)
149  {
150  // VTK has a different node order for VTK_WEDGE
151  // their triangles point outwards!
152  vtkVerts = cellShape;
153 
154  Foam::Swap(vtkVerts[1], vtkVerts[2]);
155  Foam::Swap(vtkVerts[4], vtkVerts[5]);
156 
157  cellTypes_[celli] = VTK_WEDGE;
158  }
159  else if (polyhedra_ == vtkPolyhedra::none && cellModel == tetWedge)
160  {
161  // Treat as squeezed prism (VTK_WEDGE)
162  vtkVerts.setSize(6);
163  vtkVerts[0] = cellShape[0];
164  vtkVerts[1] = cellShape[2];
165  vtkVerts[2] = cellShape[1];
166  vtkVerts[3] = cellShape[3];
167  vtkVerts[4] = cellShape[4];
168  vtkVerts[5] = cellShape[3];
169 
170  cellTypes_[celli] = VTK_WEDGE;
171  }
172  else if (polyhedra_ != vtkPolyhedra::all && cellModel == wedge)
173  {
174  // Treat as squeezed hex
175  vtkVerts.setSize(8);
176  vtkVerts[0] = cellShape[0];
177  vtkVerts[1] = cellShape[1];
178  vtkVerts[2] = cellShape[2];
179  vtkVerts[3] = cellShape[2];
180  vtkVerts[4] = cellShape[3];
181  vtkVerts[5] = cellShape[4];
182  vtkVerts[6] = cellShape[5];
183  vtkVerts[7] = cellShape[6];
184 
185  cellTypes_[celli] = VTK_HEXAHEDRON;
186  }
187  else if (polyhedra_ != vtkPolyhedra::all && cellModel == hex)
188  {
189  vtkVerts = cellShape;
190 
191  cellTypes_[celli] = VTK_HEXAHEDRON;
192  }
193  else if (polyhedra_ == vtkPolyhedra::none)
194  {
195  // Polyhedral cell. Decompose into tets + pyramids.
196 
197  // Mapping from additional point to cell
198  addPointCellLabels_[addPointi] = celli;
199 
200  // The new vertex from the cell-centre
201  const label newVertexLabel = mesh_.nPoints() + addPointi;
202 
203  // Whether to insert cell in place of original or not.
204  bool substituteCell = true;
205 
206  const labelList& cFaces = mesh_.cells()[celli];
207  forAll(cFaces, cFacei)
208  {
209  const face& f = mesh_.faces()[cFaces[cFacei]];
210  const bool isOwner = faceOwner[cFaces[cFacei]] == celli;
211 
212  // Do decomposition into triFcs and quadFcs.
213  faceList triFcs(f.size() == 4 ? 0 : f.nTriangles());
214  faceList quadFcs(f.size() == 4 ? 1 : 0);
215  if (f.size() != 4)
216  {
217  triEngine.triangulate
218  (
219  UIndirectList<point>(mesh.points(), f)
220  );
221  forAll(triEngine.triPoints(), trii)
222  {
223  triFcs[trii] = triEngine.triPoints(trii, f);
224  }
225  }
226  else
227  {
228  quadFcs[0] = f;
229  }
230 
231  forAll(quadFcs, quadI)
232  {
233  label thisCelli;
234 
235  if (substituteCell)
236  {
237  thisCelli = celli;
238  substituteCell = false;
239  }
240  else
241  {
242  thisCelli = mesh_.nCells() + addCelli;
243  superCells_[addCelli++] = celli;
244  }
245 
246  labelList& addVtkVerts = vertLabels_[thisCelli];
247 
248  addVtkVerts.setSize(5);
249 
250  const face& quad = quadFcs[quadI];
251 
252  // Ensure we have the correct orientation for the
253  // base of the primitive cell shape.
254  // If the cell is face owner, the orientation needs to be
255  // flipped.
256  // At the moment, VTK doesn't actually seem to care if
257  // negative cells are defined, but we'll do it anyhow
258  // (for safety).
259  if (isOwner)
260  {
261  addVtkVerts[0] = quad[3];
262  addVtkVerts[1] = quad[2];
263  addVtkVerts[2] = quad[1];
264  addVtkVerts[3] = quad[0];
265  }
266  else
267  {
268  addVtkVerts[0] = quad[0];
269  addVtkVerts[1] = quad[1];
270  addVtkVerts[2] = quad[2];
271  addVtkVerts[3] = quad[3];
272  }
273  addVtkVerts[4] = newVertexLabel;
274 
275  cellTypes_[thisCelli] = VTK_PYRAMID;
276  }
277 
278  forAll(triFcs, triI)
279  {
280  label thisCelli;
281 
282  if (substituteCell)
283  {
284  thisCelli = celli;
285  substituteCell = false;
286  }
287  else
288  {
289  thisCelli = mesh_.nCells() + addCelli;
290  superCells_[addCelli++] = celli;
291  }
292 
293 
294  labelList& addVtkVerts = vertLabels_[thisCelli];
295 
296  const face& tri = triFcs[triI];
297 
298  addVtkVerts.setSize(4);
299 
300  // See note above about the orientation.
301  if (isOwner)
302  {
303  addVtkVerts[0] = tri[2];
304  addVtkVerts[1] = tri[1];
305  addVtkVerts[2] = tri[0];
306  }
307  else
308  {
309  addVtkVerts[0] = tri[0];
310  addVtkVerts[1] = tri[1];
311  addVtkVerts[2] = tri[2];
312  }
313  addVtkVerts[3] = newVertexLabel;
314 
315  cellTypes_[thisCelli] = VTK_TETRA;
316  }
317  }
318 
319  addPointi++;
320  }
321  else
322  {
323  // Polyhedral cell - not decomposed
324  cellTypes_[celli] = VTK_POLYHEDRON;
325 
326  const labelList& cFaces = mesh_.cells()[celli];
327 
328  // space for the number of faces and size of each face
329  label nData = 1 + cFaces.size();
330 
331  // count total number of face points
332  forAll(cFaces, cFacei)
333  {
334  const face& f = mesh.faces()[cFaces[cFacei]];
335  nData += f.size(); // space for the face labels
336  }
337 
338  vtkVerts.setSize(nData);
339 
340  nData = 0;
341  vtkVerts[nData++] = cFaces.size();
342 
343  // build face stream
344  forAll(cFaces, cFacei)
345  {
346  const face& f = mesh.faces()[cFaces[cFacei]];
347  const bool isOwner = faceOwner[cFaces[cFacei]] == celli;
348 
349  // number of labels for this face
350  vtkVerts[nData++] = f.size();
351 
352  if (isOwner)
353  {
354  forAll(f, fp)
355  {
356  vtkVerts[nData++] = f[fp];
357  }
358  }
359  else
360  {
361  // fairly immaterial if we reverse the list
362  // or use face::reverseFace()
363  forAllReverse(f, fp)
364  {
365  vtkVerts[nData++] = f[fp];
366  }
367  }
368  }
369  }
370  }
371 
372  if (polyhedra_ == vtkPolyhedra::none)
373  {
374  Pout<< " Original cells:" << mesh_.nCells()
375  << " points:" << mesh_.nPoints()
376  << " Additional cells:" << superCells_.size()
377  << " additional points:" << addPointCellLabels_.size()
378  << nl << endl;
379  }
380 }
381 
382 
383 // ************************************************************************* //
Swap its arguments.
#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 setSize(const label)
Reset size of List.
Definition: List.C:281
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:55
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
vtkTopo(const polyMesh &mesh, const vtkPolyhedra &polyhedra)
Construct from components.
static const NamedEnum< vtkPolyhedra, 3 > vtkPolyhedraNames_
Names of groups of cell types retain as polyhedra.
Definition: vtkTopo.H:66
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
const cellShapeList & cellShapes
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
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
void Swap(T &a, T &b)
Definition: Swap.H:43
List< face > faceList
Definition: faceListFwd.H:41
static const char nl
Definition: Ostream.H:267
labelList f(nPoints)