27 #include "vtkPVFoamReader.h" 32 #include "vtkOpenFOAMPoints.H" 37 #include "vtkCellArray.h" 38 #include "vtkIdTypeArray.h" 39 #include "vtkUnstructuredGrid.h" 43 vtkUnstructuredGrid* Foam::vtkPVFoam::volumeVTKMesh
46 polyDecomp& decompInfo
73 const labelList& owner = mesh.faceOwner();
75 labelList& superCells = decompInfo.superCells();
76 labelList& addPointCellLabels = decompInfo.addPointCellLabels();
80 if (!reader_->GetUseVTKPolyhedron())
85 <<
"... scanning for polyhedra" <<
endl;
90 const cellModel& model = cellShapes[celli].model();
102 const cell& cFaces = mesh.cells()[celli];
106 const face& f = mesh.faces()[cFaces[cFacei]];
110 nAddCells += f.nTriangles();
126 addPointCellLabels.setSize(nAddPoints);
133 Info<<
" mesh nCells = " << mesh.nCells() <<
nl 134 <<
" nPoints = " << mesh.nPoints() <<
nl 135 <<
" nAddCells = " << nAddCells <<
nl 136 <<
" nAddPoints = " << nAddPoints <<
endl;
139 superCells.setSize(mesh.nCells() + nAddCells);
143 Info<<
" ... converting points" <<
endl;
148 vtkpoints->Allocate(mesh.nPoints() + nAddPoints);
160 Info<<
" ... converting cells" <<
endl;
163 vtkmesh->Allocate(mesh.nCells() + nAddCells);
166 polygonTriangulate triEngine;
169 label addPointi = 0, addCelli = 0;
173 vtkIdType nodeIds[8];
177 DynamicList<vtkIdType> faceStream(256);
181 const cellShape& cellShape = cellShapes[celli];
182 const cellModel& cellModel = cellShape.model();
184 superCells[addCelli++] = celli;
186 if (cellModel == tet)
188 for (
int j = 0; j < 4; j++)
190 nodeIds[j] = cellShape[j];
192 vtkmesh->InsertNextCell
199 else if (cellModel == pyr)
201 for (
int j = 0; j < 5; j++)
203 nodeIds[j] = cellShape[j];
205 vtkmesh->InsertNextCell
212 else if (cellModel == prism)
216 nodeIds[0] = cellShape[0];
217 nodeIds[1] = cellShape[2];
218 nodeIds[2] = cellShape[1];
219 nodeIds[3] = cellShape[3];
220 nodeIds[4] = cellShape[5];
221 nodeIds[5] = cellShape[4];
223 vtkmesh->InsertNextCell
230 else if (cellModel == tetWedge && !reader_->GetUseVTKPolyhedron())
234 nodeIds[0] = cellShape[0];
235 nodeIds[1] = cellShape[2];
236 nodeIds[2] = cellShape[1];
237 nodeIds[3] = cellShape[3];
238 nodeIds[4] = cellShape[4];
239 nodeIds[5] = cellShape[3];
241 vtkmesh->InsertNextCell
248 else if (cellModel == wedge)
252 nodeIds[0] = cellShape[0];
253 nodeIds[1] = cellShape[1];
254 nodeIds[2] = cellShape[2];
255 nodeIds[3] = cellShape[2];
256 nodeIds[4] = cellShape[3];
257 nodeIds[5] = cellShape[4];
258 nodeIds[6] = cellShape[5];
259 nodeIds[7] = cellShape[6];
261 vtkmesh->InsertNextCell
268 else if (cellModel == hex)
270 for (
int j = 0; j < 8; j++)
272 nodeIds[j] = cellShape[j];
274 vtkmesh->InsertNextCell
281 else if (reader_->GetUseVTKPolyhedron())
284 const labelList& cFaces = mesh.cells()[celli];
286 #ifdef HAS_VTK_POLYHEDRON 287 vtkIdType nFaces = cFaces.
size();
288 vtkIdType nLabels = nFaces;
293 const face& f = mesh.faces()[cFaces[cFacei]];
301 faceStream.reserve(nLabels + nFaces);
305 const face& f = mesh.faces()[cFaces[cFacei]];
306 const bool isOwner = (owner[cFaces[cFacei]] == celli);
307 const label nFacePoints = f.size();
310 faceStream.append(nFacePoints);
316 faceStream.append(f[fp]);
325 faceStream.append(f[fp]);
330 vtkmesh->InsertNextCell(VTK_POLYHEDRON, nFaces, faceStream.data());
336 HashSet<vtkIdType, Hash<label>> hashUniqId(2*256);
340 const face& f = mesh.faces()[cFaces[cFacei]];
344 hashUniqId.insert(f[fp]);
349 faceStream = hashUniqId.sortedToc();
351 vtkmesh->InsertNextCell
353 VTK_CONVEX_POINT_SET,
354 vtkIdType(faceStream.size()),
364 addPointCellLabels[addPointi] = celli;
367 const label newVertexLabel = mesh.nPoints() + addPointi;
371 bool substituteCell =
true;
373 const labelList& cFaces = mesh.cells()[celli];
376 const face& f = mesh.faces()[cFaces[cFacei]];
377 const bool isOwner = (owner[cFaces[cFacei]] == celli);
380 faceList triFcs(f.size() == 4 ? 0 : f.nTriangles());
381 faceList quadFcs(f.size() == 4 ? 1 : 0);
384 triEngine.triangulate
386 UIndirectList<point>(mesh.points(),
f)
388 forAll(triEngine.triPoints(), trii)
390 triFcs[trii] = triEngine.triPoints(trii, f);
402 substituteCell =
false;
406 superCells[addCelli++] = celli;
409 const face& quad = quadFcs[quadI];
420 nodeIds[0] = quad[3];
421 nodeIds[1] = quad[2];
422 nodeIds[2] = quad[1];
423 nodeIds[3] = quad[0];
427 nodeIds[0] = quad[0];
428 nodeIds[1] = quad[1];
429 nodeIds[2] = quad[2];
430 nodeIds[3] = quad[3];
432 nodeIds[4] = newVertexLabel;
433 vtkmesh->InsertNextCell
445 substituteCell =
false;
449 superCells[addCelli++] = celli;
452 const face& tri = triFcs[triI];
467 nodeIds[3] = newVertexLabel;
469 vtkmesh->InsertNextCell
482 vtkmesh->SetPoints(vtkpoints);
#define forAll(list, i)
Loop across all elements in list.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
static const cellModel * lookup(const word &)
Look up a model by name and return a pointer to the model or nullptr.
static void printMemory()
Simple memory used debugging information.
void vtkInsertNextOpenFOAMPoint(vtkPoints *points, const Foam::point &p)
void size(const label)
Override size to be inconsistent with allocated storage.
Ostream & endl(Ostream &os)
Add newline and flush stream.
#define forAllReverse(list, i)
Reverse loop across all elements in list.
List< cellShape > cellShapeList
List of cellShapes and PtrList of List of cellShape.
Pre-declare SubField and related Field type.
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
List< label > labelList
A List of labels.
#define InfoInFunction
Report an information message using Foam::Info.