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/dataConversion/foamToVTK/vtkTopo.C
2018-06-01 18:11:37 +02:00

298 lines
9.1 KiB
C++

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | foam-extend: Open Source CFD
\\ / O peration | Version: 4.1
\\ / A nd | Web: http://www.foam-extend.org
\\/ M anipulation | For copyright notice see file Copyright
-------------------------------------------------------------------------------
License
This file is part of foam-extend.
foam-extend 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 3 of the License, or (at your
option) any later version.
foam-extend 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 foam-extend. If not, see <http://www.gnu.org/licenses/>.
Description
Note: bug in vtk displaying wedges? Seems to display ok if we decompose
them. Should be thoroughly tested!
(they appear rarely in polyhedral meshes, do appear in some cut meshes)
\*---------------------------------------------------------------------------*/
#include "vtkTopo.H"
#include "polyMesh.H"
#include "cellShape.H"
#include "cellModeller.H"
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
// Construct from components
Foam::vtkTopo::vtkTopo(const polyMesh& mesh)
:
mesh_(mesh),
vertLabels_(),
cellTypes_(),
addPointCellLabels_(),
superCells_()
{
const cellModel& tet = *(cellModeller::lookup("tet"));
const cellModel& pyr = *(cellModeller::lookup("pyr"));
const cellModel& prism = *(cellModeller::lookup("prism"));
const cellModel& tetWedge = *(cellModeller::lookup("tetWedge"));
const cellModel& hex = *(cellModeller::lookup("hex"));
const cellShapeList& cellShapes = mesh_.cellShapes();
// 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;
// Scan for cells which need to be decomposed and count additional points
// and cells
forAll(cellShapes, cellI)
{
const cellModel& model = cellShapes[cellI].model();
if
(
model != hex
// && model != wedge // See above.
&& 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 nQuads = 0;
label nTris = 0;
f.nTrianglesQuads(mesh_.points(), nTris, nQuads);
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(nAddCells);
// List of vertex labels in VTK ordering
vertLabels_.setSize(cellShapes.size() + nAddCells);
// Label of vtk type
cellTypes_.setSize(cellShapes.size() + 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();
labelList& vtkVerts = vertLabels_[cellI];
if (cellModel == tet)
{
vtkVerts = cellShape;
cellTypes_[cellI] = VTK_TETRA;
}
else if (cellModel == pyr)
{
vtkVerts = cellShape;
cellTypes_[cellI] = VTK_PYRAMID;
}
else if (cellModel == prism)
{
vtkVerts.setSize(6);
vtkVerts[0] = cellShape[0];
vtkVerts[1] = cellShape[2];
vtkVerts[2] = cellShape[1];
vtkVerts[3] = cellShape[3];
vtkVerts[4] = cellShape[5];
vtkVerts[5] = cellShape[4];
// VTK calls this a wedge.
cellTypes_[cellI] = VTK_WEDGE;
}
else if (cellModel == tetWedge)
{
// Treat as squeezed prism
vtkVerts.setSize(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];
cellTypes_[cellI] = VTK_WEDGE;
}
// else if (cellModel == wedge)
// {
// // Treat as squeezed hex
// vtkVerts.setSize(8);
// vtkVerts[0] = cellShape[0];
// vtkVerts[1] = cellShape[1];
// vtkVerts[2] = cellShape[2];
// vtkVerts[3] = cellShape[0];
// vtkVerts[4] = cellShape[3];
// vtkVerts[5] = cellShape[4];
// vtkVerts[6] = cellShape[5];
// vtkVerts[7] = cellShape[6];
//
// cellTypes_[cellI] = VTK_HEXAHEDRON;
// }
else if (cellModel == hex)
{
vtkVerts.setSize(8);
vtkVerts[0] = cellShape[0];
vtkVerts[1] = cellShape[1];
vtkVerts[2] = cellShape[2];
vtkVerts[3] = cellShape[3];
vtkVerts[4] = cellShape[4];
vtkVerts[5] = cellShape[5];
vtkVerts[6] = cellShape[6];
vtkVerts[7] = cellShape[7];
cellTypes_[cellI] = VTK_HEXAHEDRON;
}
else
{
// Polyhedral cell. Decompose into tets + prisms.
// (see dxFoamExec/createDxConnections.C)
// Mapping from additional point to cell
addPointCellLabels_[api] = 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]];
// Number of triangles and quads in decomposition
label nTris = 0;
label nQuads = 0;
f.nTrianglesQuads(mesh_.points(), nTris, nQuads);
// Do actual decomposition into triFcs and quadFcs.
faceList triFcs(nTris);
faceList quadFcs(nQuads);
label trii = 0;
label quadi = 0;
f.trianglesQuads(mesh_.points(), trii, quadi, triFcs, quadFcs);
forAll(quadFcs, quadi)
{
label thisCellI = -1;
if (substituteCell)
{
thisCellI = cellI;
substituteCell = false;
}
else
{
thisCellI = mesh_.nCells() + aci;
superCells_[aci] = cellI;
aci++;
}
labelList& addVtkVerts = vertLabels_[thisCellI];
addVtkVerts.setSize(5);
const face& quad = quadFcs[quadi];
addVtkVerts[0] = quad[0];
addVtkVerts[1] = quad[1];
addVtkVerts[2] = quad[2];
addVtkVerts[3] = quad[3];
addVtkVerts[4] = mesh_.nPoints() + api;
cellTypes_[thisCellI] = VTK_PYRAMID;
}
forAll(triFcs, trii)
{
label thisCellI = -1;
if (substituteCell)
{
thisCellI = cellI;
substituteCell = false;
}
else
{
thisCellI = mesh_.nCells() + aci;
superCells_[aci] = cellI;
aci++;
}
labelList& addVtkVerts = vertLabels_[thisCellI];
const face& tri = triFcs[trii];
addVtkVerts.setSize(4);
addVtkVerts[0] = tri[0];
addVtkVerts[1] = tri[1];
addVtkVerts[2] = tri[2];
addVtkVerts[3] = mesh_.nPoints() + api;
cellTypes_[thisCellI] = VTK_TETRA;
}
}
api++;
}
}
Pout<< " Original cells:" << mesh_.nCells()
<< " points:" << mesh_.nPoints()
<< " Additional cells:" << superCells_.size()
<< " additional points:" << addPointCellLabels_.size()
<< nl << endl;
}
// ************************************************************************* //