vtkPVFoamTemplates.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 "polyPatch.H"
30 #include "primitivePatch.H"
31 #include "vtkOpenFOAMPoints.H"
32 
33 // VTK includes
34 #include "vtkCellArray.h"
35 #include "vtkPoints.h"
36 #include "vtkPolyData.h"
37 
38 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
39 
40 template<class PatchType>
41 vtkPolyData* Foam::vtkPVFoam::patchVTKMesh
42 (
43  const word& name,
44  const PatchType& p
45 )
46 {
47  vtkPolyData* vtkmesh = vtkPolyData::New();
48 
49  if (debug)
50  {
51  Info<< "<beg> Foam::vtkPVFoam::patchVTKMesh - " << name << endl;
52  printMemory();
53  }
54 
55  // Convert OpenFOAM mesh vertices to VTK
56  const Foam::pointField& points = p.localPoints();
57 
58  vtkPoints* vtkpoints = vtkPoints::New();
59  vtkpoints->Allocate(points.size());
60  forAll(points, i)
61  {
62  vtkInsertNextOpenFOAMPoint(vtkpoints, points[i]);
63  }
64 
65  vtkmesh->SetPoints(vtkpoints);
66  vtkpoints->Delete();
67 
68 
69  // Add faces as polygons
70  const faceList& faces = p.localFaces();
71 
72  vtkCellArray* vtkcells = vtkCellArray::New();
73  vtkcells->Allocate(faces.size());
74  forAll(faces, facei)
75  {
76  const face& f = faces[facei];
77  vtkIdType nodeIds[f.size()];
78 
79  forAll(f, fp)
80  {
81  nodeIds[fp] = f[fp];
82  }
83  vtkcells->InsertNextCell(f.size(), nodeIds);
84  }
85 
86  vtkmesh->SetPolys(vtkcells);
87  vtkcells->Delete();
88 
89  if (debug)
90  {
91  Info<< "<end> Foam::vtkPVFoam::patchVTKMesh - " << name << endl;
92  printMemory();
93  }
94 
95  return vtkmesh;
96 }
97 
98 
99 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
static void printMemory()
Simple memory used debugging information.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
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)
messageStream Info