vtkPVFoamMeshSet.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "vtkPVFoam.H"
27 
28 // OpenFOAM includes
29 #include "faceSet.H"
30 #include "pointSet.H"
31 #include "vtkOpenFOAMPoints.H"
32 
33 // VTK includes
34 #include "vtkPoints.h"
35 #include "vtkPolyData.h"
36 #include "vtkCellArray.h"
37 
38 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
39 
40 vtkPolyData* Foam::vtkPVFoam::faceSetVTKMesh
41 (
42  const fvMesh& mesh,
43  const faceSet& fSet
44 )
45 {
46  vtkPolyData* vtkmesh = vtkPolyData::New();
47 
48  if (debug)
49  {
50  Info<< "<beg> Foam::vtkPVFoam::faceSetVTKMesh" << endl;
51  printMemory();
52  }
53 
54  // Construct primitivePatch of faces in fSet.
55 
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());
64 
65 
66  // The balance of this routine should be identical to patchVTKMesh
67 
68  // Convert OpenFOAM mesh vertices to VTK
69  const pointField& points = p.localPoints();
70 
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();
79 
80  // Add faces as polygons
81  const faceList& faces = p.localFaces();
82 
83  vtkCellArray* vtkcells = vtkCellArray::New();
84  vtkcells->Allocate(faces.size());
85 
86  forAll(faces, facei)
87  {
88  const face& f = faces[facei];
89  vtkIdType nodeIds[f.size()];
90 
91  forAll(f, fp)
92  {
93  nodeIds[fp] = f[fp];
94  }
95  vtkcells->InsertNextCell(f.size(), nodeIds);
96  }
97 
98  vtkmesh->SetPolys(vtkcells);
99  vtkcells->Delete();
100 
101  if (debug)
102  {
103  Info<< "<end> Foam::vtkPVFoam::faceSetVTKMesh" << endl;
104  printMemory();
105  }
106 
107  return vtkmesh;
108 }
109 
110 
111 vtkPolyData* Foam::vtkPVFoam::pointSetVTKMesh
112 (
113  const fvMesh& mesh,
114  const pointSet& pSet
115 )
116 {
117  vtkPolyData* vtkmesh = vtkPolyData::New();
118 
119  if (debug)
120  {
121  Info<< "<beg> Foam::vtkPVFoam::pointSetVTKMesh" << endl;
122  printMemory();
123  }
124 
125  const pointField& meshPoints = mesh.points();
126 
127  vtkPoints* vtkpoints = vtkPoints::New();
128  vtkpoints->Allocate(pSet.size());
129 
130  forAllConstIter(pointSet, pSet, iter)
131  {
132  vtkInsertNextOpenFOAMPoint(vtkpoints, meshPoints[iter.key()]);
133  }
134 
135  vtkmesh->SetPoints(vtkpoints);
136  vtkpoints->Delete();
137 
138  if (debug)
139  {
140  Info<< "<end> Foam::vtkPVFoam::pointSetVTKMesh" << endl;
141  printMemory();
142  }
143 
144  return vtkmesh;
145 }
146 
147 
148 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
static void printMemory()
Simple memory used debugging information.
List< face > faceList
Definition: faceListFwd.H:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
void vtkInsertNextOpenFOAMPoint(vtkPoints *points, const Foam::point &p)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
PrimitivePatch< face, List, const pointField & > primitiveFacePatch
Foam::primitiveFacePatch.
messageStream Info
volScalarField & p