vtkPVFoamMeshLagrangian.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 "vtkPVFoamReader.h"
28 #include "vtkOpenFOAMPoints.H"
29 
30 // OpenFOAM includes
31 #include "domainDecomposition.H"
33 #include "LagrangianMesh.H"
35 #include "passiveParticleCloud.H"
36 
37 // VTK includes
38 #include "vtkCellArray.h"
39 #include "vtkPoints.h"
40 #include "vtkPolyData.h"
41 
42 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
43 
44 vtkPolyData* Foam::vtkPVFoam::lagrangianVTKMesh
45 (
46  const word& lagrangianName,
47  autoPtr<lagrangianFieldReconstructor>& lreconstructorPtr
48 )
49 {
50  vtkPolyData* vtkmesh = nullptr;
51 
52  if (reader_->GetDecomposedCase())
53  {
54  lreconstructorPtr.reset
55  (
56  new lagrangianFieldReconstructor
57  (
58  procMeshesPtr_->completeMesh(),
59  procMeshesPtr_->procMeshes(),
60  procMeshesPtr_->procFaceAddressing(),
61  procMeshesPtr_->procCellAddressing(),
62  lagrangianName
63  )
64  );
65  }
66 
67  autoPtr<passiveParticleCloud> tcloud
68  (
69  reader_->GetDecomposedCase()
70  ? lreconstructorPtr->completeCloud().ptr()
71  : new passiveParticleCloud
72  (
73  procMeshesPtr_->completeMesh(),
74  lagrangianName,
75  false
76  )
77  );
78  const passiveParticleCloud& cloud = tcloud();
79 
80  vtkmesh = vtkPolyData::New();
81  vtkPoints* vtkpoints = vtkPoints::New();
82  vtkCellArray* vtkcells = vtkCellArray::New();
83 
84  vtkpoints->Allocate(cloud.size());
85  vtkcells->Allocate(cloud.size());
86 
87  vtkIdType particleId = 0;
88  forAllConstIter(lagrangian::Cloud<passiveParticle>, cloud, iter)
89  {
90  vtkInsertNextOpenFOAMPoint(vtkpoints, iter().position(cloud.pMesh()));
91 
92  vtkcells->InsertNextCell(1, &particleId);
93  particleId++;
94  }
95 
96  vtkmesh->SetPoints(vtkpoints);
97  vtkpoints->Delete();
98 
99  vtkmesh->SetVerts(vtkcells);
100  vtkcells->Delete();
101 
102  return vtkmesh;
103 }
104 
105 
106 vtkPolyData* Foam::vtkPVFoam::LagrangianVTKMesh
107 (
108  const word& LagrangianName,
109  autoPtr<LagrangianMesh>& LmeshPtr,
110  autoPtr<LagrangianFieldReconstructor>& LreconstructorPtr
111 )
112 {
113  vtkPolyData* vtkmesh = nullptr;
114 
115  if (reader_->GetDecomposedCase())
116  {
117  LreconstructorPtr.reset
118  (
119  new LagrangianFieldReconstructor
120  (
121  procMeshesPtr_->completeMesh(),
122  procMeshesPtr_->procMeshes(),
123  procMeshesPtr_->procFaceAddressing(),
124  procMeshesPtr_->procCellAddressing(),
125  LagrangianName
126  )
127  );
128  }
129  else
130  {
131  LmeshPtr.reset
132  (
133  new LagrangianMesh(procMeshesPtr_->completeMesh(), LagrangianName)
134  );
135  }
136 
137  const LagrangianMesh& Lmesh =
138  reader_->GetDecomposedCase()
139  ? LreconstructorPtr().completeMesh()
140  : LmeshPtr();
141 
142  const LagrangianVectorInternalField Lpositions(Lmesh.position());
143 
144  vtkmesh = vtkPolyData::New();
145 
146  vtkPoints* vtkpoints = vtkPoints::New();
147  vtkCellArray* vtkcells = vtkCellArray::New();
148 
149  vtkpoints->Allocate(Lmesh.size());
150  vtkcells->Allocate(Lmesh.size());
151 
152  for (vtkIdType i = 0; i < Lmesh.size(); ++ i)
153  {
154  vtkInsertNextOpenFOAMPoint(vtkpoints, Lpositions[i]);
155 
156  vtkcells->InsertNextCell(1, &i);
157  }
158 
159  vtkmesh->SetPoints(vtkpoints);
160  vtkpoints->Delete();
161 
162  vtkmesh->SetVerts(vtkcells);
163  vtkcells->Delete();
164 
165  return vtkmesh;
166 }
167 
168 
169 // ************************************************************************* //
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
point position(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label facei, const label faceTrii, const scalar stepFraction)
Return the position given the coordinates and tet topology.
Definition: trackingI.H:224
LagrangianInternalField< vector > LagrangianVectorInternalField
void vtkInsertNextOpenFOAMPoint(vtkPoints *points, const Foam::point &p)