BUGFIX: Uninitialised member variables
[foam-extend-3.2.git] / applications / utilities / postProcessing / graphics / PV3FoamReader / vtkPV3Foam / vtkPV3FoamMeshSet.C
blob4c934f05dd2150b92ab77df43618b2a2a7089d4e
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 "faceSet.H"
33 #include "pointSet.H"
34 #include "vtkPV3FoamPoints.H"
36 // VTK includes
37 #include "vtkPoints.h"
38 #include "vtkPolyData.h"
39 #include "vtkCellArray.h"
41 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
43 vtkPolyData* Foam::vtkPV3Foam::faceSetVTKMesh
45     const fvMesh& mesh,
46     const faceSet& fSet
49     vtkPolyData* vtkmesh = vtkPolyData::New();
51     if (debug)
52     {
53         Info<< "<beg> Foam::vtkPV3Foam::faceSetVTKMesh" << endl;
54         printMemory();
55     }
57     // Construct primitivePatch of faces in fSet.
59     const faceList& meshFaces = mesh.faces();
60     faceList patchFaces(fSet.size());
61     label faceI = 0;
62     forAllConstIter(faceSet, fSet, iter)
63     {
64         patchFaces[faceI++] = meshFaces[iter.key()];
65     }
66     primitiveFacePatch p(patchFaces, mesh.points());
69     // The balance of this routine should be identical to patchVTKMesh
71     // Convert Foam mesh vertices to VTK
72     const pointField& points = p.localPoints();
74     vtkPoints *vtkpoints = vtkPoints::New();
75     vtkpoints->Allocate( points.size() );
76     forAll(points, i)
77     {
78         vtkPV3FoamInsertNextPoint(vtkpoints, points[i]);
79     }
80     vtkmesh->SetPoints(vtkpoints);
81     vtkpoints->Delete();
83     // Add faces as polygons
84     const faceList& faces = p.localFaces();
86     vtkCellArray* vtkcells = vtkCellArray::New();
87     vtkcells->Allocate( faces.size() );
89     forAll(faces, faceI)
90     {
91         const face& f = faces[faceI];
92         vtkIdType nodeIds[f.size()];
94         forAll(f, fp)
95         {
96             nodeIds[fp] = f[fp];
97         }
98         vtkcells->InsertNextCell(f.size(), nodeIds);
99     }
101     vtkmesh->SetPolys(vtkcells);
102     vtkcells->Delete();
104     if (debug)
105     {
106         Info<< "<end> Foam::vtkPV3Foam::faceSetVTKMesh" << endl;
107         printMemory();
108     }
110     return vtkmesh;
114 vtkPolyData* Foam::vtkPV3Foam::pointSetVTKMesh
116     const fvMesh& mesh,
117     const pointSet& pSet
120     vtkPolyData* vtkmesh = vtkPolyData::New();
122     if (debug)
123     {
124         Info<< "<beg> Foam::vtkPV3Foam::pointSetVTKMesh" << endl;
125         printMemory();
126     }
128     const pointField& meshPoints = mesh.points();
130     vtkPoints *vtkpoints = vtkPoints::New();
131     vtkpoints->Allocate( pSet.size() );
133     forAllConstIter(pointSet, pSet, iter)
134     {
135         vtkPV3FoamInsertNextPoint(vtkpoints, meshPoints[iter.key()]);
136     }
138     vtkmesh->SetPoints(vtkpoints);
139     vtkpoints->Delete();
141     if (debug)
142     {
143         Info<< "<end> Foam::vtkPV3Foam::pointSetVTKMesh" << endl;
144         printMemory();
145     }
147     return vtkmesh;
151 // ************************************************************************* //