initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / applications / utilities / postProcessing / graphics / PV3FoamReader / vtkPV3Foam / vtkPV3FoamMeshZone.C
blob14c5da92c2ab27bdb462d2acb0b22aa4f5913582
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 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 #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.faces();
58     faceList patchFaces(faceLabels.size());
59     forAll(faceLabels, faceI)
60     {
61         patchFaces[faceI] = meshFaces[faceLabels[faceI]];
62     }
63     primitiveFacePatch p(patchFaces, mesh.points());
66     // The balance of this routine should be identical to patchVTKMesh
68     // Convert Foam 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         vtkPV3FoamInsertNextPoint(vtkpoints, points[i]);
76     }
78     vtkmesh->SetPoints(vtkpoints);
79     vtkpoints->Delete();
82     // Add faces as polygons
83     const faceList& faces = p.localFaces();
85     vtkCellArray* vtkcells = vtkCellArray::New();
86     vtkcells->Allocate( faces.size() );
88     forAll(faces, faceI)
89     {
90         const face& f = faces[faceI];
91         vtkIdType nodeIds[f.size()];
93         forAll(f, fp)
94         {
95             nodeIds[fp] = f[fp];
96         }
97         vtkcells->InsertNextCell(f.size(), nodeIds);
98     }
100     vtkmesh->SetPolys(vtkcells);
101     vtkcells->Delete();
103     if (debug)
104     {
105         Info<< "<end> Foam::vtkPV3Foam::faceZoneVTKMesh" << endl;
106         printMemory();
107     }
109     return vtkmesh;
113 vtkPolyData* Foam::vtkPV3Foam::pointZoneVTKMesh
115     const fvMesh& mesh,
116     const labelList& pointLabels
119     vtkPolyData* vtkmesh = vtkPolyData::New();
121     if (debug)
122     {
123         Info<< "<beg> Foam::vtkPV3Foam::pointZoneVTKMesh" << endl;
124         printMemory();
125     }
127     const pointField& meshPoints = mesh.points();
129     vtkPoints *vtkpoints = vtkPoints::New();
130     vtkpoints->Allocate( pointLabels.size() );
132     forAll(pointLabels, pointI)
133     {
134         vtkPV3FoamInsertNextPoint(vtkpoints, meshPoints[pointLabels[pointI]]);
135     }
137     vtkmesh->SetPoints(vtkpoints);
138     vtkpoints->Delete();
140     if (debug)
141     {
142         Info<< "<beg> Foam::vtkPV3Foam::pointZoneVTKMesh" << endl;
143         printMemory();
144     }
146     return vtkmesh;
150 // ************************************************************************* //