27 #include "vtkPVFoamReader.h"
28 #include "vtkOpenFOAMPoints.H"
36 #include "vtkCellArray.h"
37 #include "vtkIdTypeArray.h"
38 #include "vtkUnstructuredGrid.h"
42 vtkUnstructuredGrid* Foam::vtkPVFoam::volumeVTKMesh
45 polyDecomp& decompInfo
68 labelList& superCells = decompInfo.superCells();
69 labelList& addPointCellLabels = decompInfo.addPointCellLabels();
73 if (!reader_->GetUseVTKPolyhedron())
77 const cellModel& model =
cellShapes[celli].model();
97 nAddCells +=
f.nTriangles();
113 addPointCellLabels.setSize(nAddPoints);
117 superCells.setSize(
mesh.
nCells() + nAddCells);
122 <<
", nAddCells=" << nAddCells
123 <<
", nAddPoints=" << nAddPoints <<
endl;
139 polygonTriangulate triEngine;
142 label addPointi = 0, addCelli = 0;
146 vtkIdType nodeIds[8];
150 DynamicList<vtkIdType> faceStream(256);
154 const cellShape& cellShape =
cellShapes[celli];
155 const cellModel& cellModel = cellShape.model();
157 superCells[addCelli++] = celli;
159 if (cellModel == tet)
161 for (
int j = 0; j < 4; j++)
163 nodeIds[j] = cellShape[j];
165 vtkmesh->InsertNextCell
172 else if (cellModel == pyr)
174 for (
int j = 0; j < 5; j++)
176 nodeIds[j] = cellShape[j];
178 vtkmesh->InsertNextCell
185 else if (cellModel == prism)
189 nodeIds[0] = cellShape[0];
190 nodeIds[1] = cellShape[2];
191 nodeIds[2] = cellShape[1];
192 nodeIds[3] = cellShape[3];
193 nodeIds[4] = cellShape[5];
194 nodeIds[5] = cellShape[4];
196 vtkmesh->InsertNextCell
203 else if (cellModel == tetWedge && !reader_->GetUseVTKPolyhedron())
207 nodeIds[0] = cellShape[0];
208 nodeIds[1] = cellShape[2];
209 nodeIds[2] = cellShape[1];
210 nodeIds[3] = cellShape[3];
211 nodeIds[4] = cellShape[4];
212 nodeIds[5] = cellShape[3];
214 vtkmesh->InsertNextCell
221 else if (cellModel == wedge)
225 nodeIds[0] = cellShape[0];
226 nodeIds[1] = cellShape[1];
227 nodeIds[2] = cellShape[2];
228 nodeIds[3] = cellShape[2];
229 nodeIds[4] = cellShape[3];
230 nodeIds[5] = cellShape[4];
231 nodeIds[6] = cellShape[5];
232 nodeIds[7] = cellShape[6];
234 vtkmesh->InsertNextCell
241 else if (cellModel ==
hex)
243 for (
int j = 0; j < 8; j++)
245 nodeIds[j] = cellShape[j];
247 vtkmesh->InsertNextCell
254 else if (reader_->GetUseVTKPolyhedron())
259 #ifdef HAS_VTK_POLYHEDRON
260 vtkIdType nFaces = cFaces.
size();
261 vtkIdType nLabels = nFaces;
274 faceStream.reserve(nLabels + nFaces);
279 const bool isOwner = (owner[cFaces[cFacei]] == celli);
280 const label nFacePoints =
f.size();
283 faceStream.append(nFacePoints);
289 faceStream.append(
f[fp]);
298 faceStream.append(
f[fp]);
303 vtkmesh->InsertNextCell(VTK_POLYHEDRON, nFaces, faceStream.data());
309 HashSet<vtkIdType, Hash<label>> hashUniqId(2*256);
317 hashUniqId.insert(
f[fp]);
322 faceStream = hashUniqId.sortedToc();
324 vtkmesh->InsertNextCell
326 VTK_CONVEX_POINT_SET,
327 vtkIdType(faceStream.size()),
337 addPointCellLabels[addPointi] = celli;
355 bool substituteCell =
true;
361 const bool isOwner = (owner[cFaces[cFacei]] == celli);
364 faceList triFcs(
f.size() == 4 ? 0 :
f.nTriangles());
368 triEngine.triangulate
372 forAll(triEngine.triPoints(), trii)
374 triFcs[trii] = triEngine.triPoints(trii,
f);
386 substituteCell =
false;
390 superCells[addCelli++] = celli;
393 const face& quad = quadFcs[quadI];
404 nodeIds[0] = quad[3];
405 nodeIds[1] = quad[2];
406 nodeIds[2] = quad[1];
407 nodeIds[3] = quad[0];
411 nodeIds[0] = quad[0];
412 nodeIds[1] = quad[1];
413 nodeIds[2] = quad[2];
414 nodeIds[3] = quad[3];
416 nodeIds[4] = newVertexLabel;
417 vtkmesh->InsertNextCell
429 substituteCell =
false;
433 superCells[addCelli++] = celli;
436 const face& tri = triFcs[triI];
451 nodeIds[3] = newVertexLabel;
453 vtkmesh->InsertNextCell
466 vtkmesh->SetPoints(vtkpoints);
#define forAll(list, i)
Loop across all elements in list.
#define forAllReverse(list, i)
Reverse loop across all elements in list.
void size(const label)
Override size to be inconsistent with allocated storage.
static const cellModel * lookup(const word &)
Look up a model by name and return a pointer to the model or nullptr.
virtual const faceList & faces() const
Return raw faces.
virtual const labelList & faceOwner() const
Return face owner.
virtual const pointField & points() const
Return raw points.
const vectorField & cellCentres() const
const cellShapeList & cellShapes() const
Return cell shapes.
bool has(const Type &(primitiveMesh::*method)() const) const
Return whether the result of the given method has been allocated.
const cellList & cells() const
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
const cellShapeList & cellShapes
#define DebugInFunction
Report an information message using Foam::Info.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
List< cellShape > cellShapeList
List of cellShapes and PtrList of List of cellShape.
List< label > labelList
A List of labels.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Ostream & endl(Ostream &os)
Add newline and flush stream.
IOstream & hex(IOstream &io)
void vtkInsertNextOpenFOAMPoint(vtkPoints *points, const Foam::point &p)