Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / applications / utilities / postProcessing / graphics / PV3Readers / PV3FoamReader / vtkPV3Foam / vtkPV3FoamMeshPatch.C
blobd66b690ff1a280196b043357100f8e2ace75f2c6
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 "polyPatch.H"
30 #include "primitivePatch.H"
31 #include "vtkOpenFOAMPoints.H"
33 // VTK includes
34 #include "vtkCellArray.h"
35 #include "vtkPoints.h"
36 #include "vtkPolyData.h"
38 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
40 vtkPolyData* Foam::vtkPV3Foam::patchVTKMesh(const polyPatch& p)
42     vtkPolyData* vtkmesh = vtkPolyData::New();
44     if (debug)
45     {
46         Info<< "<beg> Foam::vtkPV3Foam::patchVTKMesh - " << p.name() << endl;
47         printMemory();
48     }
50     // Convert OpenFOAM mesh vertices to VTK
51     const Foam::pointField& points = p.localPoints();
53     vtkPoints* vtkpoints = vtkPoints::New();
54     vtkpoints->Allocate(points.size());
55     forAll(points, i)
56     {
57         vtkInsertNextOpenFOAMPoint(vtkpoints, points[i]);
58     }
60     vtkmesh->SetPoints(vtkpoints);
61     vtkpoints->Delete();
64     // Add faces as polygons
65     const faceList& faces = p.localFaces();
67     vtkCellArray* vtkcells = vtkCellArray::New();
68     vtkcells->Allocate(faces.size());
69     forAll(faces, faceI)
70     {
71         const face& f = faces[faceI];
72         vtkIdType nodeIds[f.size()];
74         forAll(f, fp)
75         {
76             nodeIds[fp] = f[fp];
77         }
78         vtkcells->InsertNextCell(f.size(), nodeIds);
79     }
81     vtkmesh->SetPolys(vtkcells);
82     vtkcells->Delete();
84     if (debug)
85     {
86         Info<< "<end> Foam::vtkPV3Foam::patchVTKMesh - " << p.name() << endl;
87         printMemory();
88     }
90     return vtkmesh;
94 // ************************************************************************* //