Merge branch 'master' of ssh://git.code.sf.net/p/foam-extend/foam-extend-3.2
[foam-extend-3.2.git] / applications / utilities / postProcessing / graphics / PV3FoamReader / vtkPV3Foam / vtkPV3FoamMeshVolume.C
blob13224596414270a71bd9c41f3456bedf0bdaf5c5
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend is free software: you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by the
13     Free Software Foundation, either version 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  If not, see <http://www.gnu.org/licenses/>.
24 Description
26 \*---------------------------------------------------------------------------*/
28 #include "vtkPV3Foam.H"
30 // Foam includes
31 #include "fvMesh.H"
32 #include "cellModeller.H"
33 #include "vtkPV3FoamPoints.H"
35 // VTK includes
36 #include "vtkCellArray.h"
37 #include "vtkUnstructuredGrid.h"
39 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
41 vtkUnstructuredGrid* Foam::vtkPV3Foam::volumeVTKMesh
43     const fvMesh& mesh,
44     polyDecomp& decompInfo
47     vtkUnstructuredGrid* vtkmesh = vtkUnstructuredGrid::New();
49     if (debug)
50     {
51         Info<< "<beg> Foam::vtkPV3Foam::volumeVTKMesh" << endl;
52         printMemory();
53     }
55     // Number of additional points needed by the decomposition of polyhedra
56     label nAddPoints = 0;
58     // Number of additional cells generated by the decomposition of polyhedra
59     label nAddCells = 0;
61     labelList& superCells = decompInfo.superCells();
62     labelList& addPointCellLabels = decompInfo.addPointCellLabels();
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"));
71     // Scan for cells which need to be decomposed and count additional points
72     // and cells
73     if (debug)
74     {
75         Info<< "... building cell-shapes" << endl;
76     }
77     const cellShapeList& cellShapes = mesh.cellShapes();
79     if (debug)
80     {
81         Info<< "... scanning" << endl;
82     }
83     forAll(cellShapes, cellI)
84     {
85         const cellModel& model = cellShapes[cellI].model();
87         if
88         (
89             model != hex
90          && model != wedge
91          && model != prism
92          && model != pyr
93          && model != tet
94          && model != tetWedge
95         )
96         {
97             const cell& cFaces = mesh.cells()[cellI];
99             forAll(cFaces, cFaceI)
100             {
101                 const face& f = mesh.faces()[cFaces[cFaceI]];
103                 label nFacePoints = f.size();
105                 label nQuads = (nFacePoints - 2)/2;
106                 label nTris = (nFacePoints - 2)%2;
107                 nAddCells += nQuads + nTris;
108             }
110             nAddCells--;
111             nAddPoints++;
112         }
113     }
115     // Set size of additional point addressing array
116     // (from added point to original cell)
117     addPointCellLabels.setSize(nAddPoints);
119     // Set size of additional cells mapping array
120     // (from added cell to original cell)
122     if (debug)
123     {
124         Info<<" mesh nCells     = " << mesh.nCells() << nl
125             <<"      nPoints    = " << mesh.nPoints() << nl
126             <<"      nAddCells  = " << nAddCells << nl
127             <<"      nAddPoints = " << nAddPoints << endl;
128     }
130     superCells.setSize(mesh.nCells() + nAddCells);
132     if (debug)
133     {
134         Info<< "... converting points" << endl;
135     }
137     // Convert Foam mesh vertices to VTK
138     vtkPoints *vtkpoints = vtkPoints::New();
139     vtkpoints->Allocate( mesh.nPoints() + nAddPoints );
141     const Foam::pointField& points = mesh.points();
143     forAll(points, i)
144     {
145         vtkPV3FoamInsertNextPoint(vtkpoints, points[i]);
146     }
149     if (debug)
150     {
151         Info<< "... converting cells" << endl;
152     }
154     vtkmesh->Allocate( mesh.nCells() + nAddCells );
156     // Set counters for additional points and additional cells
157     label addPointI = 0, addCellI = 0;
159     // Create storage for points - needed for mapping from Foam to VTK
160     // data types - max 'order' = hex = 8 points
161     vtkIdType nodeIds[8];
163     forAll(cellShapes, cellI)
164     {
165         const cellShape& cellShape = cellShapes[cellI];
166         const cellModel& cellModel = cellShape.model();
168         superCells[addCellI++] = cellI;
170         if (cellModel == tet)
171         {
172             for (int j = 0; j < 4; j++)
173             {
174                 nodeIds[j] = cellShape[j];
175             }
176             vtkmesh->InsertNextCell
177             (
178                 VTK_TETRA,
179                 4,
180                 nodeIds
181             );
182         }
183         else if (cellModel == pyr)
184         {
185             for (int j = 0; j < 5; j++)
186             {
187                 nodeIds[j] = cellShape[j];
188             }
189             vtkmesh->InsertNextCell
190             (
191                 VTK_PYRAMID,
192                 5,
193                 nodeIds
194             );
195         }
196         else if (cellModel == prism)
197         {
198             for (int j = 0; j < 6; j++)
199             {
200                 nodeIds[j] = cellShape[j];
201             }
202             vtkmesh->InsertNextCell
203             (
204                 VTK_WEDGE,
205                 6,
206                 nodeIds
207             );
208         }
209         else if (cellModel == tetWedge)
210         {
211             // Treat as squeezed prism
213             nodeIds[0] = cellShape[0];
214             nodeIds[1] = cellShape[2];
215             nodeIds[2] = cellShape[1];
216             nodeIds[3] = cellShape[3];
217             nodeIds[4] = cellShape[4];
218             nodeIds[5] = cellShape[4];
220             vtkmesh->InsertNextCell
221             (
222                 VTK_WEDGE,
223                 6,
224                 nodeIds
225             );
226         }
227         else if (cellModel == wedge)
228         {
229             // Treat as squeezed hex
231             nodeIds[0] = cellShape[0];
232             nodeIds[1] = cellShape[1];
233             nodeIds[2] = cellShape[2];
234             nodeIds[3] = cellShape[2];
235             nodeIds[4] = cellShape[3];
236             nodeIds[5] = cellShape[4];
237             nodeIds[6] = cellShape[5];
238             nodeIds[7] = cellShape[6];
240             vtkmesh->InsertNextCell
241             (
242                 VTK_HEXAHEDRON,
243                 8,
244                 nodeIds
245             );
246         }
247         else if (cellModel == hex)
248         {
249             for (int j = 0; j < 8; j++)
250             {
251                 nodeIds[j] = cellShape[j];
252             }
253             vtkmesh->InsertNextCell
254             (
255                 VTK_HEXAHEDRON,
256                 8,
257                 nodeIds
258             );
259         }
260         else
261         {
262             // Polyhedral cell. Decompose into tets + prisms.
264             // Mapping from additional point to cell
265             addPointCellLabels[addPointI] = cellI;
267             // Insert the new vertex from the cell-centre
268             label newVertexLabel = mesh.nPoints() + addPointI;
269             vtkPV3FoamInsertNextPoint(vtkpoints, mesh.C()[cellI]);
271             // Whether to insert cell in place of original or not.
272             bool substituteCell = true;
274             const labelList& cFaces = mesh.cells()[cellI];
276             forAll(cFaces, cFaceI)
277             {
278                 const face& f = mesh.faces()[cFaces[cFaceI]];
280                 label nFacePoints = f.size();
282                 label nQuads = (nFacePoints - 2)/2;
283                 label nTris = (nFacePoints - 2)%2;
285                 label qpi = 0;
287                 for (label quadi=0; quadi<nQuads; quadi++)
288                 {
289                     if (substituteCell)
290                     {
291                         substituteCell = false;
292                     }
293                     else
294                     {
295                         superCells[addCellI++] = cellI;
296                     }
298                     nodeIds[0] = f[0];
299                     nodeIds[1] = f[qpi + 1];
300                     nodeIds[2] = f[qpi + 2];
301                     nodeIds[3] = f[qpi + 3];
302                     nodeIds[4] = newVertexLabel;
303                     vtkmesh->InsertNextCell
304                     (
305                         VTK_PYRAMID,
306                         5,
307                         nodeIds
308                     );
310                     qpi += 2;
311                 }
313                 if (nTris)
314                 {
315                     if (substituteCell)
316                     {
317                         substituteCell = false;
318                     }
319                     else
320                     {
321                         superCells[addCellI++] = cellI;
322                     }
324                     nodeIds[0] = f[0];
325                     nodeIds[1] = f[qpi + 1];
326                     nodeIds[2] = f[qpi + 2];
327                     nodeIds[3] = newVertexLabel;
328                     vtkmesh->InsertNextCell
329                     (
330                         VTK_TETRA,
331                         4,
332                         nodeIds
333                     );
334                 }
335             }
337             addPointI++;
338         }
339     }
341     vtkmesh->SetPoints(vtkpoints);
342     vtkpoints->Delete();
344     if (debug)
345     {
346         Info<< "<end> Foam::vtkPV3Foam::volumeVTKMesh" << endl;
347         printMemory();
348     }
350     return vtkmesh;
354 // ************************************************************************* //