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