374 lines
10 KiB
C
374 lines
10 KiB
C
|
/*---------------------------------------------------------------------------*\
|
||
|
========= |
|
||
|
\\ / 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
|
||
|
|
||
|
\*---------------------------------------------------------------------------*/
|
||
|
|
||
|
#ifndef vtkPV3FoamAddVolumeMesh_H
|
||
|
#define vtkPV3FoamAddVolumeMesh_H
|
||
|
|
||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||
|
|
||
|
#include "vtkPV3Foam.H"
|
||
|
|
||
|
// Foam includes
|
||
|
#include "fvMesh.H"
|
||
|
#include "cellModeller.H"
|
||
|
#include "vtkPV3FoamInsertNextPoint.H"
|
||
|
|
||
|
// VTK includes
|
||
|
#include "vtkCellArray.h"
|
||
|
#include "vtkUnstructuredGrid.h"
|
||
|
|
||
|
|
||
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||
|
|
||
|
void Foam::vtkPV3Foam::addVolumeMesh
|
||
|
(
|
||
|
const fvMesh& mesh,
|
||
|
vtkUnstructuredGrid* vtkmesh,
|
||
|
labelList& superCells
|
||
|
)
|
||
|
{
|
||
|
if (debug)
|
||
|
{
|
||
|
Info<< "<beg> Foam::vtkPV3Foam::addVolumeMesh" << endl;
|
||
|
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;
|
||
|
|
||
|
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;
|
||
|
}
|
||
|
|
||
|
// Get access to faces for volume mesh analysis
|
||
|
const faceList& faces = mesh.faces();
|
||
|
|
||
|
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)
|
||
|
{
|
||
|
const face& f = faces[cFaces[cFaceI]];
|
||
|
|
||
|
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)
|
||
|
addPointCellLabels_.setSize(nAddPoints);
|
||
|
|
||
|
// Set size of additional cells mapping array
|
||
|
// (from added cell to original cell)
|
||
|
|
||
|
if (debug)
|
||
|
{
|
||
|
Info<<"mesh.nCells() = " << mesh.nCells() << nl
|
||
|
<<"mesh.nPoints() = " << mesh.nPoints() << nl
|
||
|
<<"nAddCells = " << nAddCells << nl
|
||
|
<<"nAddPoints = " << nAddPoints << endl;
|
||
|
}
|
||
|
|
||
|
superCells.setSize(mesh.nCells() + nAddCells);
|
||
|
|
||
|
if (debug)
|
||
|
{
|
||
|
Info<< "... converting points" << endl;
|
||
|
}
|
||
|
|
||
|
// Convert Foam mesh vertices to VTK
|
||
|
vtkPoints *vtkpoints = vtkPoints::New();
|
||
|
vtkpoints->Allocate(mesh.nPoints() + nAddPoints);
|
||
|
|
||
|
// Use all points: support for inactive points and faces.
|
||
|
// HJ, 28/Mar/2009
|
||
|
const pointField& points = mesh.allPoints();
|
||
|
|
||
|
forAll(points, i)
|
||
|
{
|
||
|
vtkPV3FoamInsertNextPoint(vtkpoints, points[i]);
|
||
|
}
|
||
|
|
||
|
|
||
|
if (debug)
|
||
|
{
|
||
|
Info<< "... converting cells" << endl;
|
||
|
}
|
||
|
|
||
|
vtkmesh->Allocate(mesh.nCells() + nAddCells);
|
||
|
|
||
|
// Set counters for additional points and additional cells
|
||
|
label api = 0, aci = 0;
|
||
|
|
||
|
// Create storage for points - needed for mapping from Foam to VTK
|
||
|
// data types - max 'order' = hex = 8 points
|
||
|
vtkIdType nodeIds[8];
|
||
|
|
||
|
forAll(cellShapes, celli)
|
||
|
{
|
||
|
const cellShape& cellShape = cellShapes[celli];
|
||
|
const cellModel& cellModel = cellShape.model();
|
||
|
|
||
|
superCells[aci++] = celli;
|
||
|
|
||
|
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
|
||
|
addPointCellLabels_[api] = celli;
|
||
|
|
||
|
// Insert the new vertex from the cell-centre
|
||
|
label newVertexLabel = mesh.nPoints() + api;
|
||
|
vtkPV3FoamInsertNextPoint(vtkpoints, mesh.C()[celli]);
|
||
|
|
||
|
// Whether to insert cell in place of original or not.
|
||
|
bool substituteCell = true;
|
||
|
|
||
|
const labelList& cFaces = mesh.cells()[celli];
|
||
|
|
||
|
forAll(cFaces, cFaceI)
|
||
|
{
|
||
|
const face& f = faces[cFaces[cFaceI]];
|
||
|
|
||
|
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)
|
||
|
{
|
||
|
thisCellI = celli;
|
||
|
substituteCell = false;
|
||
|
}
|
||
|
else
|
||
|
{
|
||
|
thisCellI = mesh.nCells() + aci;
|
||
|
superCells[aci++] = celli;
|
||
|
}
|
||
|
|
||
|
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)
|
||
|
{
|
||
|
thisCellI = celli;
|
||
|
substituteCell = false;
|
||
|
}
|
||
|
else
|
||
|
{
|
||
|
thisCellI = mesh.nCells() + aci;
|
||
|
superCells[aci++] = celli;
|
||
|
}
|
||
|
|
||
|
nodeIds[0] = f[0];
|
||
|
nodeIds[1] = f[qpi + 1];
|
||
|
nodeIds[2] = f[qpi + 2];
|
||
|
nodeIds[3] = newVertexLabel;
|
||
|
vtkmesh->InsertNextCell
|
||
|
(
|
||
|
VTK_TETRA,
|
||
|
4,
|
||
|
nodeIds
|
||
|
);
|
||
|
}
|
||
|
}
|
||
|
|
||
|
api++;
|
||
|
}
|
||
|
}
|
||
|
|
||
|
vtkmesh->SetPoints(vtkpoints);
|
||
|
vtkpoints->Delete();
|
||
|
|
||
|
if (debug)
|
||
|
{
|
||
|
Info<< "<end> Foam::vtkPV3Foam::addVolumeMesh" << endl;
|
||
|
printMemory();
|
||
|
}
|
||
|
}
|
||
|
|
||
|
|
||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||
|
|
||
|
#endif
|
||
|
|
||
|
// ************************************************************************* //
|