initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / applications / utilities / postProcessing / graphics / PV3FoamReader / vtkPV3Foam / vtkPV3FoamAddVolumeMesh.C
blob292fb8d2a0964cd6f75582ea356d67b59e00dbae
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Description
27 \*---------------------------------------------------------------------------*/
29 #ifndef vtkPV3FoamAddVolumeMesh_H
30 #define vtkPV3FoamAddVolumeMesh_H
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 #include "vtkPV3Foam.H"
36 // Foam includes
37 #include "fvMesh.H"
38 #include "cellModeller.H"
39 #include "vtkPV3FoamInsertNextPoint.H"
41 // VTK includes
42 #include "vtkCellArray.h"
43 #include "vtkUnstructuredGrid.h"
46 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
48 void Foam::vtkPV3Foam::addVolumeMesh
50     const fvMesh& mesh,
51     vtkUnstructuredGrid* vtkmesh,
52     labelList& superCells
55     if (debug)
56     {
57         Info<< "<beg> Foam::vtkPV3Foam::addVolumeMesh" << endl;
58         printMemory();
59     }
61     // Number of additional points needed by the decomposition of polyhedra
62     label nAddPoints = 0;
64     // Number of additional cells generated by the decomposition of polyhedra
65     label nAddCells = 0;
67     const cellModel& tet = *(cellModeller::lookup("tet"));
68     const cellModel& pyr = *(cellModeller::lookup("pyr"));
69     const cellModel& prism = *(cellModeller::lookup("prism"));
70     const cellModel& wedge = *(cellModeller::lookup("wedge"));
71     const cellModel& tetWedge = *(cellModeller::lookup("tetWedge"));
72     const cellModel& hex = *(cellModeller::lookup("hex"));
74     // Scan for cells which need to be decomposed and count additional points
75     // and cells
76     if (debug)
77     {
78         Info<< "... building cell-shapes" << endl;
79     }
80     const cellShapeList& cellShapes = mesh.cellShapes();
82     if (debug)
83     {
84         Info<< "... scanning" << endl;
85     }
86     forAll(cellShapes, cellI)
87     {
88         const cellModel& model = cellShapes[cellI].model();
90         if
91         (
92             model != hex
93          && model != wedge
94          && model != prism
95          && model != pyr
96          && model != tet
97          && model != tetWedge
98         )
99         {
100             const cell& cFaces = mesh.cells()[cellI];
102             forAll(cFaces, cFaceI)
103             {
104                 const face& f = mesh.faces()[cFaces[cFaceI]];
106                 label nFacePoints = f.size();
108                 label nQuads = (nFacePoints - 2)/2;
109                 label nTris = (nFacePoints - 2)%2;
110                 nAddCells += nQuads + nTris;
111             }
113             nAddCells--;
114             nAddPoints++;
115         }
116     }
118     // Set size of additional point addressing array
119     // (from added point to original cell)
120     addPointCellLabels_.setSize(nAddPoints);
122     // Set size of additional cells mapping array
123     // (from added cell to original cell)
125     if (debug)
126     {
127         Info<<"mesh.nCells()  = " << mesh.nCells() << nl
128             <<"mesh.nPoints() = " << mesh.nPoints() << nl
129             <<"nAddCells      = " << nAddCells << nl
130             <<"nAddPoints     = " << nAddPoints << endl;
131     }
133     superCells.setSize(mesh.nCells() + nAddCells);
135     if (debug)
136     {
137         Info<< "... converting points" << endl;
138     }
140     // Convert Foam mesh vertices to VTK
141     vtkPoints *vtkpoints = vtkPoints::New();
142     vtkpoints->Allocate(mesh.nPoints() + nAddPoints);
144     const Foam::pointField& points = mesh.points();
146     forAll(points, i)
147     {
148         vtkPV3FoamInsertNextPoint(vtkpoints, points[i]);
149     }
152     if (debug)
153     {
154         Info<< "... converting cells" << endl;
155     }
157     vtkmesh->Allocate(mesh.nCells() + nAddCells);
159     // Set counters for additional points and additional cells
160     label api = 0, aci = 0;
162     // Create storage for points - needed for mapping from Foam to VTK
163     // data types - max 'order' = hex = 8 points
164     vtkIdType nodeIds[8];
166     forAll(cellShapes, celli)
167     {
168         const cellShape& cellShape = cellShapes[celli];
169         const cellModel& cellModel = cellShape.model();
171         superCells[aci++] = celli;
173         if (cellModel == tet)
174         {
175             for (int j = 0; j < 4; j++)
176             {
177                 nodeIds[j] = cellShape[j];
178             }
179             vtkmesh->InsertNextCell
180             (
181                 VTK_TETRA,
182                 4,
183                 nodeIds
184             );
185         }
186         else if (cellModel == pyr)
187         {
188             for (int j = 0; j < 5; j++)
189             {
190                 nodeIds[j] = cellShape[j];
191             }
192             vtkmesh->InsertNextCell
193             (
194                 VTK_PYRAMID,
195                 5,
196                 nodeIds
197             );
198         }
199         else if (cellModel == prism)
200         {
201             for (int j = 0; j < 6; j++)
202             {
203                 nodeIds[j] = cellShape[j];
204             }
205             vtkmesh->InsertNextCell
206             (
207                 VTK_WEDGE,
208                 6,
209                 nodeIds
210             );
211         }
212         else if (cellModel == tetWedge)
213         {
214             // Treat as squeezed prism
216             nodeIds[0] = cellShape[0];
217             nodeIds[1] = cellShape[2];
218             nodeIds[2] = cellShape[1];
219             nodeIds[3] = cellShape[3];
220             nodeIds[4] = cellShape[4];
221             nodeIds[5] = cellShape[4];
223             vtkmesh->InsertNextCell
224             (
225                 VTK_WEDGE,
226                 6,
227                 nodeIds
228             );
229         }
230         else if (cellModel == wedge)
231         {
232             // Treat as squeezed hex
234             nodeIds[0] = cellShape[0];
235             nodeIds[1] = cellShape[1];
236             nodeIds[2] = cellShape[2];
237             nodeIds[3] = cellShape[2];
238             nodeIds[4] = cellShape[3];
239             nodeIds[5] = cellShape[4];
240             nodeIds[6] = cellShape[5];
241             nodeIds[7] = cellShape[6];
243             vtkmesh->InsertNextCell
244             (
245                 VTK_HEXAHEDRON,
246                 8,
247                 nodeIds
248             );
249         }
250         else if (cellModel == hex)
251         {
252             for (int j = 0; j < 8; j++)
253             {
254                 nodeIds[j] = cellShape[j];
255             }
256             vtkmesh->InsertNextCell
257             (
258                 VTK_HEXAHEDRON,
259                 8,
260                 nodeIds
261             );
262         }
263         else
264         {
265             // Polyhedral cell. Decompose into tets + prisms.
267             // Mapping from additional point to cell
268             addPointCellLabels_[api] = celli;
270             // Insert the new vertex from the cell-centre
271             label newVertexLabel = mesh.nPoints() + api;
272             vtkPV3FoamInsertNextPoint(vtkpoints, mesh.C()[celli]);
274             // Whether to insert cell in place of original or not.
275             bool substituteCell = true;
277             const labelList& cFaces = mesh.cells()[celli];
279             forAll(cFaces, cFaceI)
280             {
281                 const face& f = mesh.faces()[cFaces[cFaceI]];
283                 label nFacePoints = f.size();
285                 label nQuads = (nFacePoints - 2)/2;
286                 label nTris = (nFacePoints - 2)%2;
288                 label qpi = 0;
290                 for (label quadi=0; quadi<nQuads; quadi++)
291                 {
292                     label thisCellI = -1;
294                     if (substituteCell)
295                     {
296                         thisCellI = celli;
297                         substituteCell = false;
298                     }
299                     else
300                     {
301                         thisCellI = mesh.nCells() + aci;
302                         superCells[aci++] = celli;
303                     }
305                     nodeIds[0] = f[0];
306                     nodeIds[1] = f[qpi + 1];
307                     nodeIds[2] = f[qpi + 2];
308                     nodeIds[3] = f[qpi + 3];
309                     nodeIds[4] = newVertexLabel;
310                     vtkmesh->InsertNextCell
311                     (
312                         VTK_PYRAMID,
313                         5,
314                         nodeIds
315                     );
317                     qpi += 2;
318                 }
320                 if (nTris)
321                 {
322                     label thisCellI = -1;
324                     if (substituteCell)
325                     {
326                         thisCellI = celli;
327                         substituteCell = false;
328                     }
329                     else
330                     {
331                         thisCellI = mesh.nCells() + aci;
332                         superCells[aci++] = celli;
333                     }
335                     nodeIds[0] = f[0];
336                     nodeIds[1] = f[qpi + 1];
337                     nodeIds[2] = f[qpi + 2];
338                     nodeIds[3] = newVertexLabel;
339                     vtkmesh->InsertNextCell
340                     (
341                         VTK_TETRA,
342                         4,
343                         nodeIds
344                     );
345                 }
346             }
348             api++;
349         }
350     }
352     vtkmesh->SetPoints(vtkpoints);
353     vtkpoints->Delete();
355     if (debug)
356     {
357         Info<< "<end> Foam::vtkPV3Foam::addVolumeMesh" << endl;
358         printMemory();
359     }
363 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
365 #endif
367 // ************************************************************************* //