2010-05-12 13:27:55 +00:00
|
|
|
/*---------------------------------------------------------------------------*\
|
|
|
|
========= |
|
|
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
|
|
\\ / O peration |
|
|
|
|
\\ / A nd | Copyright held by original author
|
|
|
|
\\/ M anipulation |
|
|
|
|
-------------------------------------------------------------------------------
|
|
|
|
License
|
|
|
|
This file is part of OpenFOAM.
|
|
|
|
|
|
|
|
OpenFOAM is free software; you can redistribute it and/or modify it
|
|
|
|
under the terms of the GNU General Public License as published by the
|
|
|
|
Free Software Foundation; either version 2 of the License, or (at your
|
|
|
|
option) any later version.
|
|
|
|
|
|
|
|
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
|
|
|
|
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
|
|
|
|
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
|
|
|
|
for more details.
|
|
|
|
|
|
|
|
You should have received a copy of the GNU General Public License
|
|
|
|
along with OpenFOAM; if not, write to the Free Software Foundation,
|
|
|
|
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
|
|
|
|
|
|
|
|
Description
|
|
|
|
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
|
|
|
|
#include "vtkPV3Foam.H"
|
|
|
|
|
|
|
|
// Foam includes
|
|
|
|
#include "fvMesh.H"
|
|
|
|
#include "cellModeller.H"
|
2010-08-26 14:22:03 +00:00
|
|
|
#include "vtkPV3FoamPoints.H"
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
// VTK includes
|
|
|
|
#include "vtkCellArray.h"
|
|
|
|
#include "vtkUnstructuredGrid.h"
|
|
|
|
|
|
|
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
|
|
|
|
2010-08-26 14:22:03 +00:00
|
|
|
vtkUnstructuredGrid* Foam::vtkPV3Foam::volumeVTKMesh
|
2010-05-12 13:27:55 +00:00
|
|
|
(
|
|
|
|
const fvMesh& mesh,
|
2010-08-26 14:22:03 +00:00
|
|
|
polyDecomp& decompInfo
|
2010-05-12 13:27:55 +00:00
|
|
|
)
|
|
|
|
{
|
2010-08-26 14:22:03 +00:00
|
|
|
vtkUnstructuredGrid* vtkmesh = vtkUnstructuredGrid::New();
|
|
|
|
|
2010-05-12 13:27:55 +00:00
|
|
|
if (debug)
|
|
|
|
{
|
2010-08-26 14:22:03 +00:00
|
|
|
Info<< "<beg> Foam::vtkPV3Foam::volumeVTKMesh" << endl;
|
2010-05-12 13:27:55 +00:00
|
|
|
printMemory();
|
|
|
|
}
|
|
|
|
|
|
|
|
// Number of additional points needed by the decomposition of polyhedra
|
|
|
|
label nAddPoints = 0;
|
|
|
|
|
|
|
|
// Number of additional cells generated by the decomposition of polyhedra
|
|
|
|
label nAddCells = 0;
|
|
|
|
|
2010-08-26 14:22:03 +00:00
|
|
|
labelList& superCells = decompInfo.superCells();
|
|
|
|
labelList& addPointCellLabels = decompInfo.addPointCellLabels();
|
|
|
|
|
2010-05-12 13:27:55 +00:00
|
|
|
const cellModel& tet = *(cellModeller::lookup("tet"));
|
|
|
|
const cellModel& pyr = *(cellModeller::lookup("pyr"));
|
|
|
|
const cellModel& prism = *(cellModeller::lookup("prism"));
|
|
|
|
const cellModel& wedge = *(cellModeller::lookup("wedge"));
|
|
|
|
const cellModel& tetWedge = *(cellModeller::lookup("tetWedge"));
|
|
|
|
const cellModel& hex = *(cellModeller::lookup("hex"));
|
|
|
|
|
|
|
|
// Scan for cells which need to be decomposed and count additional points
|
|
|
|
// and cells
|
|
|
|
if (debug)
|
|
|
|
{
|
|
|
|
Info<< "... building cell-shapes" << endl;
|
|
|
|
}
|
|
|
|
const cellShapeList& cellShapes = mesh.cellShapes();
|
|
|
|
|
|
|
|
if (debug)
|
|
|
|
{
|
|
|
|
Info<< "... scanning" << endl;
|
|
|
|
}
|
|
|
|
forAll(cellShapes, cellI)
|
|
|
|
{
|
|
|
|
const cellModel& model = cellShapes[cellI].model();
|
|
|
|
|
|
|
|
if
|
|
|
|
(
|
|
|
|
model != hex
|
|
|
|
&& model != wedge
|
|
|
|
&& model != prism
|
|
|
|
&& model != pyr
|
|
|
|
&& model != tet
|
|
|
|
&& model != tetWedge
|
|
|
|
)
|
|
|
|
{
|
|
|
|
const cell& cFaces = mesh.cells()[cellI];
|
|
|
|
|
|
|
|
forAll(cFaces, cFaceI)
|
|
|
|
{
|
2010-08-26 14:22:03 +00:00
|
|
|
const face& f = mesh.faces()[cFaces[cFaceI]];
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
label nFacePoints = f.size();
|
|
|
|
|
|
|
|
label nQuads = (nFacePoints - 2)/2;
|
|
|
|
label nTris = (nFacePoints - 2)%2;
|
|
|
|
nAddCells += nQuads + nTris;
|
|
|
|
}
|
|
|
|
|
|
|
|
nAddCells--;
|
|
|
|
nAddPoints++;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
// Set size of additional point addressing array
|
|
|
|
// (from added point to original cell)
|
2010-08-26 14:22:03 +00:00
|
|
|
addPointCellLabels.setSize(nAddPoints);
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
// Set size of additional cells mapping array
|
|
|
|
// (from added cell to original cell)
|
|
|
|
|
|
|
|
if (debug)
|
|
|
|
{
|
2010-08-26 14:22:03 +00:00
|
|
|
Info<<" mesh nCells = " << mesh.nCells() << nl
|
|
|
|
<<" nPoints = " << mesh.nPoints() << nl
|
|
|
|
<<" nAddCells = " << nAddCells << nl
|
|
|
|
<<" nAddPoints = " << nAddPoints << endl;
|
2010-05-12 13:27:55 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
superCells.setSize(mesh.nCells() + nAddCells);
|
|
|
|
|
|
|
|
if (debug)
|
|
|
|
{
|
|
|
|
Info<< "... converting points" << endl;
|
|
|
|
}
|
|
|
|
|
|
|
|
// Convert Foam mesh vertices to VTK
|
|
|
|
vtkPoints *vtkpoints = vtkPoints::New();
|
2010-08-26 14:22:03 +00:00
|
|
|
vtkpoints->Allocate( mesh.nPoints() + nAddPoints );
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2010-08-26 14:22:03 +00:00
|
|
|
const Foam::pointField& points = mesh.points();
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
forAll(points, i)
|
|
|
|
{
|
|
|
|
vtkPV3FoamInsertNextPoint(vtkpoints, points[i]);
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
if (debug)
|
|
|
|
{
|
|
|
|
Info<< "... converting cells" << endl;
|
|
|
|
}
|
|
|
|
|
2010-08-26 14:22:03 +00:00
|
|
|
vtkmesh->Allocate( mesh.nCells() + nAddCells );
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
// Set counters for additional points and additional cells
|
2010-08-26 14:22:03 +00:00
|
|
|
label addPointI = 0, addCellI = 0;
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
// Create storage for points - needed for mapping from Foam to VTK
|
|
|
|
// data types - max 'order' = hex = 8 points
|
|
|
|
vtkIdType nodeIds[8];
|
|
|
|
|
2010-08-26 14:22:03 +00:00
|
|
|
forAll(cellShapes, cellI)
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
2010-08-26 14:22:03 +00:00
|
|
|
const cellShape& cellShape = cellShapes[cellI];
|
2010-05-12 13:27:55 +00:00
|
|
|
const cellModel& cellModel = cellShape.model();
|
|
|
|
|
2010-08-26 14:22:03 +00:00
|
|
|
superCells[addCellI++] = cellI;
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
if (cellModel == tet)
|
|
|
|
{
|
|
|
|
for (int j = 0; j < 4; j++)
|
|
|
|
{
|
|
|
|
nodeIds[j] = cellShape[j];
|
|
|
|
}
|
|
|
|
vtkmesh->InsertNextCell
|
|
|
|
(
|
|
|
|
VTK_TETRA,
|
|
|
|
4,
|
|
|
|
nodeIds
|
|
|
|
);
|
|
|
|
}
|
|
|
|
else if (cellModel == pyr)
|
|
|
|
{
|
|
|
|
for (int j = 0; j < 5; j++)
|
|
|
|
{
|
|
|
|
nodeIds[j] = cellShape[j];
|
|
|
|
}
|
|
|
|
vtkmesh->InsertNextCell
|
|
|
|
(
|
|
|
|
VTK_PYRAMID,
|
|
|
|
5,
|
|
|
|
nodeIds
|
|
|
|
);
|
|
|
|
}
|
|
|
|
else if (cellModel == prism)
|
|
|
|
{
|
|
|
|
for (int j = 0; j < 6; j++)
|
|
|
|
{
|
|
|
|
nodeIds[j] = cellShape[j];
|
|
|
|
}
|
|
|
|
vtkmesh->InsertNextCell
|
|
|
|
(
|
|
|
|
VTK_WEDGE,
|
|
|
|
6,
|
|
|
|
nodeIds
|
|
|
|
);
|
|
|
|
}
|
|
|
|
else if (cellModel == tetWedge)
|
|
|
|
{
|
|
|
|
// Treat as squeezed prism
|
|
|
|
|
|
|
|
nodeIds[0] = cellShape[0];
|
|
|
|
nodeIds[1] = cellShape[2];
|
|
|
|
nodeIds[2] = cellShape[1];
|
|
|
|
nodeIds[3] = cellShape[3];
|
|
|
|
nodeIds[4] = cellShape[4];
|
|
|
|
nodeIds[5] = cellShape[4];
|
|
|
|
|
|
|
|
vtkmesh->InsertNextCell
|
|
|
|
(
|
|
|
|
VTK_WEDGE,
|
|
|
|
6,
|
|
|
|
nodeIds
|
|
|
|
);
|
|
|
|
}
|
|
|
|
else if (cellModel == wedge)
|
|
|
|
{
|
|
|
|
// Treat as squeezed hex
|
|
|
|
|
|
|
|
nodeIds[0] = cellShape[0];
|
|
|
|
nodeIds[1] = cellShape[1];
|
|
|
|
nodeIds[2] = cellShape[2];
|
|
|
|
nodeIds[3] = cellShape[2];
|
|
|
|
nodeIds[4] = cellShape[3];
|
|
|
|
nodeIds[5] = cellShape[4];
|
|
|
|
nodeIds[6] = cellShape[5];
|
|
|
|
nodeIds[7] = cellShape[6];
|
|
|
|
|
|
|
|
vtkmesh->InsertNextCell
|
|
|
|
(
|
|
|
|
VTK_HEXAHEDRON,
|
|
|
|
8,
|
|
|
|
nodeIds
|
|
|
|
);
|
|
|
|
}
|
|
|
|
else if (cellModel == hex)
|
|
|
|
{
|
|
|
|
for (int j = 0; j < 8; j++)
|
|
|
|
{
|
|
|
|
nodeIds[j] = cellShape[j];
|
|
|
|
}
|
|
|
|
vtkmesh->InsertNextCell
|
|
|
|
(
|
|
|
|
VTK_HEXAHEDRON,
|
|
|
|
8,
|
|
|
|
nodeIds
|
|
|
|
);
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
// Polyhedral cell. Decompose into tets + prisms.
|
|
|
|
|
|
|
|
// Mapping from additional point to cell
|
2010-08-26 14:22:03 +00:00
|
|
|
addPointCellLabels[addPointI] = cellI;
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
// Insert the new vertex from the cell-centre
|
2010-08-26 14:22:03 +00:00
|
|
|
label newVertexLabel = mesh.nPoints() + addPointI;
|
|
|
|
vtkPV3FoamInsertNextPoint(vtkpoints, mesh.C()[cellI]);
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
// Whether to insert cell in place of original or not.
|
|
|
|
bool substituteCell = true;
|
|
|
|
|
2010-08-26 14:22:03 +00:00
|
|
|
const labelList& cFaces = mesh.cells()[cellI];
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
forAll(cFaces, cFaceI)
|
|
|
|
{
|
2010-08-26 14:22:03 +00:00
|
|
|
const face& f = mesh.faces()[cFaces[cFaceI]];
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
label nFacePoints = f.size();
|
|
|
|
|
|
|
|
label nQuads = (nFacePoints - 2)/2;
|
|
|
|
label nTris = (nFacePoints - 2)%2;
|
|
|
|
|
|
|
|
label qpi = 0;
|
|
|
|
|
|
|
|
for (label quadi=0; quadi<nQuads; quadi++)
|
|
|
|
{
|
|
|
|
label thisCellI = -1;
|
|
|
|
|
|
|
|
if (substituteCell)
|
|
|
|
{
|
2010-08-26 14:22:03 +00:00
|
|
|
thisCellI = cellI;
|
2010-05-12 13:27:55 +00:00
|
|
|
substituteCell = false;
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
2010-08-26 14:22:03 +00:00
|
|
|
thisCellI = mesh.nCells() + addCellI;
|
|
|
|
superCells[addCellI++] = cellI;
|
2010-05-12 13:27:55 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
nodeIds[0] = f[0];
|
|
|
|
nodeIds[1] = f[qpi + 1];
|
|
|
|
nodeIds[2] = f[qpi + 2];
|
|
|
|
nodeIds[3] = f[qpi + 3];
|
|
|
|
nodeIds[4] = newVertexLabel;
|
|
|
|
vtkmesh->InsertNextCell
|
|
|
|
(
|
|
|
|
VTK_PYRAMID,
|
|
|
|
5,
|
|
|
|
nodeIds
|
|
|
|
);
|
|
|
|
|
|
|
|
qpi += 2;
|
|
|
|
}
|
|
|
|
|
|
|
|
if (nTris)
|
|
|
|
{
|
|
|
|
label thisCellI = -1;
|
|
|
|
|
|
|
|
if (substituteCell)
|
|
|
|
{
|
2010-08-26 14:22:03 +00:00
|
|
|
thisCellI = cellI;
|
2010-05-12 13:27:55 +00:00
|
|
|
substituteCell = false;
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
2010-08-26 14:22:03 +00:00
|
|
|
thisCellI = mesh.nCells() + addCellI;
|
|
|
|
superCells[addCellI++] = cellI;
|
2010-05-12 13:27:55 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
nodeIds[0] = f[0];
|
|
|
|
nodeIds[1] = f[qpi + 1];
|
|
|
|
nodeIds[2] = f[qpi + 2];
|
|
|
|
nodeIds[3] = newVertexLabel;
|
|
|
|
vtkmesh->InsertNextCell
|
|
|
|
(
|
|
|
|
VTK_TETRA,
|
|
|
|
4,
|
|
|
|
nodeIds
|
|
|
|
);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2010-08-26 14:22:03 +00:00
|
|
|
addPointI++;
|
2010-05-12 13:27:55 +00:00
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
vtkmesh->SetPoints(vtkpoints);
|
|
|
|
vtkpoints->Delete();
|
|
|
|
|
|
|
|
if (debug)
|
|
|
|
{
|
2010-08-26 14:22:03 +00:00
|
|
|
Info<< "<end> Foam::vtkPV3Foam::volumeVTKMesh" << endl;
|
2010-05-12 13:27:55 +00:00
|
|
|
printMemory();
|
|
|
|
}
|
|
|
|
|
2010-08-26 14:22:03 +00:00
|
|
|
return vtkmesh;
|
|
|
|
}
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
|
|
|
|
// ************************************************************************* //
|