BUGFIX: Uninitialised member variables
[foam-extend-3.2.git] / applications / utilities / postProcessing / graphics / PV3FoamReader / vtkPV3Foam / vtkPV3FoamMeshZone.C
blob2edc0b97245ba70d4f7fb242c7b4bb3f67880cf6
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
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 #include "vtkPV3Foam.H"
31 // Foam includes
32 #include "vtkPV3FoamPoints.H"
34 // VTK includes
35 #include "vtkPoints.h"
36 #include "vtkPolyData.h"
37 #include "vtkCellArray.h"
39 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
41 vtkPolyData* Foam::vtkPV3Foam::faceZoneVTKMesh
43     const fvMesh& mesh,
44     const labelList& faceLabels
47     vtkPolyData* vtkmesh = vtkPolyData::New();
49     if (debug)
50     {
51         Info<< "<beg> Foam::vtkPV3Foam::faceZoneVTKMesh" << endl;
52         printMemory();
53     }
55     // Construct primitivePatch of faces in faceZone
57     const faceList& meshFaces = mesh.allFaces();
58     faceList patchFaces(faceLabels.size());
59     label npf = 0;
61     // Filter faces that are not in live mesh
62     // Bug fix.  HJ, 21/Mar/2011
63     forAll(faceLabels, faceI)
64     {
65         if (faceLabels[faceI] < mesh.nFaces())
66         {
67             patchFaces[npf] = meshFaces[faceLabels[faceI]];
68             npf++;
69         }
70     }
71     patchFaces.setSize(npf);
73     primitiveFacePatch p(patchFaces, mesh.points());
76     // The balance of this routine should be identical to patchVTKMesh
78     // Convert Foam mesh vertices to VTK
79     const pointField& points = p.localPoints();
81     vtkPoints* vtkpoints = vtkPoints::New();
82     vtkpoints->Allocate(points.size());
83     forAll(points, i)
84     {
85         vtkPV3FoamInsertNextPoint(vtkpoints, points[i]);
86     }
88     vtkmesh->SetPoints(vtkpoints);
89     vtkpoints->Delete();
92     // Add faces as polygons
93     const faceList& faces = p.localFaces();
95     vtkCellArray* vtkcells = vtkCellArray::New();
96     vtkcells->Allocate( faces.size() );
98     forAll(faces, faceI)
99     {
100         const face& f = faces[faceI];
101         vtkIdType nodeIds[f.size()];
103         forAll(f, fp)
104         {
105             nodeIds[fp] = f[fp];
106         }
107         vtkcells->InsertNextCell(f.size(), nodeIds);
108     }
110     vtkmesh->SetPolys(vtkcells);
111     vtkcells->Delete();
113     if (debug)
114     {
115         Info<< "<end> Foam::vtkPV3Foam::faceZoneVTKMesh" << endl;
116         printMemory();
117     }
119     return vtkmesh;
123 vtkPolyData* Foam::vtkPV3Foam::pointZoneVTKMesh
125     const fvMesh& mesh,
126     const labelList& pointLabels
129     vtkPolyData* vtkmesh = vtkPolyData::New();
131     if (debug)
132     {
133         Info<< "<beg> Foam::vtkPV3Foam::pointZoneVTKMesh" << endl;
134         printMemory();
135     }
137     const pointField& meshPoints = mesh.allPoints();
139     // Filter point labels to include only live points
140     labelList pl(pointLabels.size());
141     label npl = 0;
143     forAll (pointLabels, pointI)
144     {
145         if (pointLabels[pointI] < mesh.nPoints())
146         {
147             pl[npl] = pointLabels[pointI];
148             npl++;
149         }
150     }
151     pl.setSize(npl);
153     vtkPoints *vtkpoints = vtkPoints::New();
154     vtkpoints->Allocate( pl.size());
156     forAll(pointLabels, pointI)
157     {
158         vtkPV3FoamInsertNextPoint(vtkpoints, meshPoints[pl[pointI]]);
159     }
161     vtkmesh->SetPoints(vtkpoints);
162     vtkpoints->Delete();
164     if (debug)
165     {
166         Info<< "<beg> Foam::vtkPV3Foam::pointZoneVTKMesh" << endl;
167         printMemory();
168     }
170     return vtkmesh;
174 // ************************************************************************* //