vtkPVFoamMeshSet.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2025 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 #include "vtkOpenFOAMPoints.H"
28 
29 // OpenFOAM includes
30 #include "fvMesh.H"
31 #include "faceSet.H"
32 #include "pointSet.H"
33 
34 // VTK includes
35 #include "vtkPoints.h"
36 #include "vtkPolyData.h"
37 #include "vtkCellArray.h"
38 
39 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
40 
41 vtkPolyData* Foam::vtkPVFoam::faceSetVTKMesh
42 (
43  const fvMesh& mesh,
44  const faceSet& fSet
45 )
46 {
47  vtkPolyData* vtkmesh = vtkPolyData::New();
48 
49  // Construct primitivePatch of faces in fSet.
50  const faceList& meshFaces = mesh.faces();
51  faceList patchFaces(fSet.size());
52  label facei = 0;
53  forAllConstIter(faceSet, fSet, iter)
54  {
55  patchFaces[facei++] = meshFaces[iter.key()];
56  }
57  primitiveFacePatch p(patchFaces, mesh.points());
58 
59  // The balance of this routine should be identical to patchVTKMesh
60 
61  // Convert OpenFOAM mesh vertices to VTK
62  const pointField& points = p.localPoints();
63 
64  vtkPoints* vtkpoints = vtkPoints::New();
65  vtkpoints->Allocate(points.size());
66  forAll(points, i)
67  {
68  vtkInsertNextOpenFOAMPoint(vtkpoints, points[i]);
69  }
70  vtkmesh->SetPoints(vtkpoints);
71  vtkpoints->Delete();
72 
73  // Add faces as polygons
74  const faceList& faces = p.localFaces();
75 
76  vtkCellArray* vtkcells = vtkCellArray::New();
77  vtkcells->Allocate(faces.size());
78 
79  DynamicList<vtkIdType> nodeIds(4);
80 
81  forAll(faces, facei)
82  {
83  const face& f = faces[facei];
84 
85  nodeIds.resize(f.size());
86  forAll(f, fp)
87  {
88  nodeIds[fp] = f[fp];
89  }
90 
91  vtkcells->InsertNextCell(f.size(), nodeIds.data());
92  }
93 
94  vtkmesh->SetPolys(vtkcells);
95  vtkcells->Delete();
96 
97  return vtkmesh;
98 }
99 
100 
101 vtkPolyData* Foam::vtkPVFoam::pointSetVTKMesh
102 (
103  const fvMesh& mesh,
104  const pointSet& pSet
105 )
106 {
107  vtkPolyData* vtkmesh = vtkPolyData::New();
108 
109  const pointField& meshPoints = mesh.points();
110 
111  vtkPoints* vtkpoints = vtkPoints::New();
112  vtkpoints->Allocate(pSet.size());
113 
114  forAllConstIter(pointSet, pSet, iter)
115  {
116  vtkInsertNextOpenFOAMPoint(vtkpoints, meshPoints[iter.key()]);
117  }
118 
119  vtkmesh->SetPoints(vtkpoints);
120  vtkpoints->Delete();
121 
122  return vtkmesh;
123 }
124 
125 
126 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1344
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1331
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
const pointField & points
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
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
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
PrimitivePatch< List< face >, const pointField & > primitiveFacePatch
Foam::primitiveFacePatch.
List< face > faceList
Definition: faceListFwd.H:41
labelList f(nPoints)
volScalarField & p
void vtkInsertNextOpenFOAMPoint(vtkPoints *points, const Foam::point &p)