Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / applications / utilities / postProcessing / graphics / PV3Readers / PV3FoamReader / vtkPV3Foam / vtkPV3FoamMeshSet.C
blob38b3e87a47e8a9d20479acbf19a349177bfbdf88
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2004-2010 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
13     the Free Software Foundation, either version 3 of the License, or
14     (at your 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, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "vtkPV3Foam.H"
28 // OpenFOAM includes
29 #include "faceSet.H"
30 #include "pointSet.H"
31 #include "vtkOpenFOAMPoints.H"
33 // VTK includes
34 #include "vtkPoints.h"
35 #include "vtkPolyData.h"
36 #include "vtkCellArray.h"
38 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
40 vtkPolyData* Foam::vtkPV3Foam::faceSetVTKMesh
42     const fvMesh& mesh,
43     const faceSet& fSet
46     vtkPolyData* vtkmesh = vtkPolyData::New();
48     if (debug)
49     {
50         Info<< "<beg> Foam::vtkPV3Foam::faceSetVTKMesh" << endl;
51         printMemory();
52     }
54     // Construct primitivePatch of faces in fSet.
56     const faceList& meshFaces = mesh.faces();
57     faceList patchFaces(fSet.size());
58     label faceI = 0;
59     forAllConstIter(faceSet, fSet, iter)
60     {
61         patchFaces[faceI++] = meshFaces[iter.key()];
62     }
63     primitiveFacePatch p(patchFaces, mesh.points());
66     // The balance of this routine should be identical to patchVTKMesh
68     // Convert OpenFOAM mesh vertices to VTK
69     const pointField& points = p.localPoints();
71     vtkPoints* vtkpoints = vtkPoints::New();
72     vtkpoints->Allocate(points.size());
73     forAll(points, i)
74     {
75         vtkInsertNextOpenFOAMPoint(vtkpoints, points[i]);
76     }
77     vtkmesh->SetPoints(vtkpoints);
78     vtkpoints->Delete();
80     // Add faces as polygons
81     const faceList& faces = p.localFaces();
83     vtkCellArray* vtkcells = vtkCellArray::New();
84     vtkcells->Allocate(faces.size());
86     forAll(faces, faceI)
87     {
88         const face& f = faces[faceI];
89         vtkIdType nodeIds[f.size()];
91         forAll(f, fp)
92         {
93             nodeIds[fp] = f[fp];
94         }
95         vtkcells->InsertNextCell(f.size(), nodeIds);
96     }
98     vtkmesh->SetPolys(vtkcells);
99     vtkcells->Delete();
101     if (debug)
102     {
103         Info<< "<end> Foam::vtkPV3Foam::faceSetVTKMesh" << endl;
104         printMemory();
105     }
107     return vtkmesh;
111 vtkPolyData* Foam::vtkPV3Foam::pointSetVTKMesh
113     const fvMesh& mesh,
114     const pointSet& pSet
117     vtkPolyData* vtkmesh = vtkPolyData::New();
119     if (debug)
120     {
121         Info<< "<beg> Foam::vtkPV3Foam::pointSetVTKMesh" << endl;
122         printMemory();
123     }
125     const pointField& meshPoints = mesh.points();
127     vtkPoints* vtkpoints = vtkPoints::New();
128     vtkpoints->Allocate(pSet.size());
130     forAllConstIter(pointSet, pSet, iter)
131     {
132         vtkInsertNextOpenFOAMPoint(vtkpoints, meshPoints[iter.key()]);
133     }
135     vtkmesh->SetPoints(vtkpoints);
136     vtkpoints->Delete();
138     if (debug)
139     {
140         Info<< "<end> Foam::vtkPV3Foam::pointSetVTKMesh" << endl;
141         printMemory();
142     }
144     return vtkmesh;
148 // ************************************************************************* //