initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / applications / utilities / postProcessing / graphics / PV3FoamReader / vtkPV3Foam / vtkPV3FoamAddSetMesh.C
blob18835fe427f38c36d7a9aa1ac13f6a47af54fb9b
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 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 "faceSet.H"
33 #include "pointSet.H"
34 #include "vtkPV3FoamInsertNextPoint.H"
36 // VTK includes
37 #include "vtkPoints.h"
38 #include "vtkPolyData.h"
39 #include "vtkCellArray.h"
41 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
43 void Foam::vtkPV3Foam::addFaceSetMesh
45     const fvMesh& mesh,
46     const faceSet& fSet,
47     vtkPolyData* vtkmesh
50     if (debug)
51     {
52         Info<< "<beg> Foam::vtkPV3Foam::addFaceSetMesh" << endl;
53         printMemory();
54     }
56     // Construct primitivePatch of faces in fSet.
58     const faceList& meshFaces = mesh.faces();
59     faceList patchFaces(fSet.size());
60     label faceI = 0;
61     forAllConstIter(faceSet, fSet, iter)
62     {
63         patchFaces[faceI++] = meshFaces[iter.key()];
64     }
65     primitiveFacePatch p(patchFaces, mesh.points());
68     // The balance of this routine should be identical to addPatchMesh
70     // Convert Foam mesh vertices to VTK
71     const pointField& points = p.localPoints();
73     vtkPoints *vtkpoints = vtkPoints::New();
74     vtkpoints->Allocate(points.size());
75     forAll(points, i)
76     {
77         vtkPV3FoamInsertNextPoint(vtkpoints, points[i]);
78     }
79     vtkmesh->SetPoints(vtkpoints);
80     vtkpoints->Delete();
82     // Add faces as polygons
83     const faceList& faces = p.localFaces();
85     vtkCellArray* vtkcells = vtkCellArray::New();
86     vtkcells->Allocate(p.size());
87     forAll(faces, faceI)
88     {
89         const face& f = faces[faceI];
90         vtkIdType nodeIds[f.size()];
92         forAll (f, fp)
93         {
94             nodeIds[fp] = f[fp];
95         }
96         vtkcells->InsertNextCell(f.size(), nodeIds);
97     }
99     vtkmesh->SetPolys(vtkcells);
100     vtkcells->Delete();
102     if (debug)
103     {
104         Info<< "<end> Foam::vtkPV3Foam::addFaceSetMesh" << endl;
105         printMemory();
106     }
110 void Foam::vtkPV3Foam::addPointSetMesh
112     const fvMesh& mesh,
113     const pointSet& pSet,
114     vtkPolyData* vtkmesh
117     if (debug)
118     {
119         Info<< "<beg> Foam::vtkPV3Foam::addPointSetMesh" << endl;
120         printMemory();
121     }
123     const pointField& meshPoints = mesh.points();
125     vtkPoints *vtkpoints = vtkPoints::New();
126     vtkpoints->Allocate(pSet.size());
128     forAllConstIter(pointSet, pSet, iter)
129     {
130         vtkPV3FoamInsertNextPoint(vtkpoints, meshPoints[iter.key()]);
131     }
133     vtkmesh->SetPoints(vtkpoints);
134     vtkpoints->Delete();
136     if (debug)
137     {
138         Info<< "<end> Foam::vtkPV3Foam::addPointSetMesh" << endl;
139         printMemory();
140     }
143 // ************************************************************************* //