This repository has been archived on 2023-11-20. You can view files and clone it, but cannot push or open issues or pull requests.
foam-extend4.1-coherent-io/applications/utilities/postProcessing/graphics/PVFoamReader/vtkFoam/vtkFoamAddInternalMesh.C
2013-07-18 10:15:54 +02:00

299 lines
8.6 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
\*---------------------------------------------------------------------------*/
#include "vtkFoam.H"
#include "fvMesh.H"
#include "cellModeller.H"
#include "vtkUnstructuredGrid.h"
#include "vtkCellArray.h"
#include "vtkFoamInsertNextPoint.H"
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
void Foam::vtkFoam::addInternalMesh
(
const fvMesh& mesh,
vtkUnstructuredGrid* vtkMesh
)
{
SetName(vtkMesh, "Internal Mesh");
// 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;
}
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 = mesh.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)
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);
const Foam::pointField& points = mesh.points();
forAll(points, i)
{
vtkFoamInsertNextPoint(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;
forAll(cellShapes, celli)
{
const cellShape& cellShape = cellShapes[celli];
const cellModel& cellModel = cellShape.model();
superCells_[aci++] = celli;
if (cellModel == tet)
{
vtkMesh->InsertNextCell
(
VTK_TETRA,
4,
const_cast<int*>(cellShape.begin())
);
}
else if (cellModel == pyr)
{
vtkMesh->InsertNextCell
(
VTK_PYRAMID,
5,
const_cast<int*>(cellShape.begin())
);
}
else if (cellModel == prism)
{
vtkMesh->InsertNextCell
(
VTK_WEDGE,
6,
const_cast<int*>(cellShape.begin())
);
}
else if (cellModel == tetWedge)
{
// Treat as squeezed prism
int vtkVerts[6];
vtkVerts[0] = cellShape[0];
vtkVerts[1] = cellShape[2];
vtkVerts[2] = cellShape[1];
vtkVerts[3] = cellShape[3];
vtkVerts[4] = cellShape[4];
vtkVerts[5] = cellShape[4];
vtkMesh->InsertNextCell(VTK_WEDGE, 6, vtkVerts);
}
else if (cellModel == wedge)
{
// Treat as squeezed hex
int vtkVerts[8];
vtkVerts[0] = cellShape[0];
vtkVerts[1] = cellShape[1];
vtkVerts[2] = cellShape[2];
vtkVerts[3] = cellShape[2];
vtkVerts[4] = cellShape[3];
vtkVerts[5] = cellShape[4];
vtkVerts[6] = cellShape[5];
vtkVerts[7] = cellShape[6];
vtkMesh->InsertNextCell(VTK_HEXAHEDRON, 8, vtkVerts);
}
else if (cellModel == hex)
{
vtkMesh->InsertNextCell
(
VTK_HEXAHEDRON,
8,
const_cast<int*>(cellShape.begin())
);
}
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;
vtkFoamInsertNextPoint(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 = mesh.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;
}
int addVtkVerts[5];
addVtkVerts[0] = f[0];
addVtkVerts[1] = f[qpi + 1];
addVtkVerts[2] = f[qpi + 2];
addVtkVerts[3] = f[qpi + 3];
addVtkVerts[4] = newVertexLabel;
vtkMesh->InsertNextCell(VTK_PYRAMID, 5, addVtkVerts);
qpi += 2;
}
if (nTris)
{
label thisCellI = -1;
if (substituteCell)
{
thisCellI = celli;
substituteCell = false;
}
else
{
thisCellI = mesh.nCells() + aci;
superCells_[aci++] = celli;
}
int addVtkVerts[4];
addVtkVerts[0] = f[0];
addVtkVerts[1] = f[qpi + 1];
addVtkVerts[2] = f[qpi + 2];
addVtkVerts[3] = newVertexLabel;
vtkMesh->InsertNextCell(VTK_TETRA, 4, addVtkVerts);
}
}
api++;
}
}
vtkMesh->SetPoints(vtkpoints);
vtkpoints->Delete();
}
// ************************************************************************* //