foamToVtk: Changed -poly option to -polyhedra
The new option takes a value indicating which cell types should be
written out as polyhedra. The values are as follows:
none: No polyhedral cells are written. Cells which match specific
shapes (hex, pyramid, prism, ...) are written as their
corresponding VTK cell types. Arbitrary polyhedral cells
that do not match a specific shape are decomposed into
tetrahedra.
polyhedra: Only arbitrary polyhedral cells are written as a VTK
polyhedron. Cells that match specific shapes are written as
their corresponding VTK cell types.
all: All cells are written as a VTK polyhedron.
The default is 'none', which retains the previous default behaviour.
This commit is contained in:
@ -2,7 +2,7 @@
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration | Website: https://openfoam.org
|
||||
\\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
|
||||
\\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
|
||||
\\/ M anipulation |
|
||||
-------------------------------------------------------------------------------
|
||||
License
|
||||
@ -78,8 +78,8 @@ Usage
|
||||
- \par -noLinks
|
||||
(in parallel) do not link processor files to master
|
||||
|
||||
- \par poly
|
||||
write polyhedral cells without tet/pyramid decomposition
|
||||
- \par -polyhedra
|
||||
Which cell types to write as polyhedra - 'none', 'polyhedra', or 'all'
|
||||
|
||||
- \par -allPatches
|
||||
Combine all patches into a single file
|
||||
@ -269,10 +269,11 @@ int main(int argc, char *argv[])
|
||||
"ascii",
|
||||
"write in ASCII format instead of binary"
|
||||
);
|
||||
argList::addBoolOption
|
||||
argList::addOption
|
||||
(
|
||||
"poly",
|
||||
"write polyhedral cells without tet/pyramid decomposition"
|
||||
"polyhedra",
|
||||
"types",
|
||||
"cell types to write as polyhedra - 'none', 'polyhedra', or 'all'"
|
||||
);
|
||||
argList::addBoolOption
|
||||
(
|
||||
@ -329,9 +330,15 @@ int main(int argc, char *argv[])
|
||||
const bool doLinks = !args.optionFound("noLinks");
|
||||
bool binary = !args.optionFound("ascii");
|
||||
const bool useTimeName = args.optionFound("useTimeName");
|
||||
|
||||
// Decomposition of polyhedral cells into tets/pyramids cells
|
||||
vtkTopo::decomposePoly = !args.optionFound("poly");
|
||||
const vtkTopo::vtkPolyhedra polyhedra =
|
||||
vtkTopo::vtkPolyhedraNames_
|
||||
[
|
||||
args.optionLookupOrDefault<word>
|
||||
(
|
||||
"polyhedra",
|
||||
vtkTopo::vtkPolyhedraNames_[vtkTopo::vtkPolyhedra::none]
|
||||
)
|
||||
];
|
||||
|
||||
if (binary && (sizeof(floatScalar) != 4 || sizeof(label) != 4))
|
||||
{
|
||||
@ -435,7 +442,7 @@ int main(int argc, char *argv[])
|
||||
|
||||
|
||||
// Mesh wrapper; does subsetting and decomposition
|
||||
vtkMesh vMesh(mesh, cellSetName);
|
||||
vtkMesh vMesh(mesh, polyhedra, cellSetName);
|
||||
|
||||
|
||||
// Scan for all possible lagrangian clouds
|
||||
|
||||
@ -2,7 +2,7 @@
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration | Website: https://openfoam.org
|
||||
\\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
|
||||
\\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
|
||||
\\/ M anipulation |
|
||||
-------------------------------------------------------------------------------
|
||||
License
|
||||
@ -33,14 +33,16 @@ License
|
||||
Foam::vtkMesh::vtkMesh
|
||||
(
|
||||
fvMesh& baseMesh,
|
||||
const vtkTopo::vtkPolyhedra polyhedra,
|
||||
const word& setName
|
||||
)
|
||||
:
|
||||
baseMesh_(baseMesh),
|
||||
subsetter_(baseMesh),
|
||||
polyhedra_(polyhedra),
|
||||
setName_(setName)
|
||||
{
|
||||
if (setName.size())
|
||||
if (setName != word::null)
|
||||
{
|
||||
// Read cellSet using whole mesh
|
||||
cellSet currentSet(baseMesh_, setName_);
|
||||
|
||||
@ -2,7 +2,7 @@
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration | Website: https://openfoam.org
|
||||
\\ / A nd | Copyright (C) 2011-2019 OpenFOAM Foundation
|
||||
\\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
|
||||
\\/ M anipulation |
|
||||
-------------------------------------------------------------------------------
|
||||
License
|
||||
@ -58,9 +58,12 @@ class vtkMesh
|
||||
//- Reference to mesh
|
||||
fvMesh& baseMesh_;
|
||||
|
||||
//- Subsetting engine + sub-fvMesh
|
||||
//- Subsetting engine
|
||||
fvMeshSubset subsetter_;
|
||||
|
||||
//- Cell types to retain as polyhedra
|
||||
const vtkTopo::vtkPolyhedra polyhedra_;
|
||||
|
||||
//- Current cellSet (or empty)
|
||||
const word setName_;
|
||||
|
||||
@ -73,7 +76,12 @@ public:
|
||||
// Constructors
|
||||
|
||||
//- Construct from components
|
||||
vtkMesh(fvMesh& baseMesh, const word& setName = "");
|
||||
vtkMesh
|
||||
(
|
||||
fvMesh& baseMesh,
|
||||
const vtkTopo::vtkPolyhedra polyhedra = vtkTopo::vtkPolyhedra::none,
|
||||
const word& setName = word::null
|
||||
);
|
||||
|
||||
//- Disallow default bitwise copy construction
|
||||
vtkMesh(const vtkMesh&) = delete;
|
||||
@ -89,6 +97,7 @@ public:
|
||||
return baseMesh_;
|
||||
}
|
||||
|
||||
//- Subsetting engine
|
||||
const fvMeshSubset& subsetter() const
|
||||
{
|
||||
return subsetter_;
|
||||
@ -100,12 +109,12 @@ public:
|
||||
return setName_.size();
|
||||
}
|
||||
|
||||
//- topology
|
||||
//- VTK topology
|
||||
const vtkTopo& topo() const
|
||||
{
|
||||
if (topoPtr_.empty())
|
||||
{
|
||||
topoPtr_.reset(new vtkTopo(mesh()));
|
||||
topoPtr_.reset(new vtkTopo(mesh(), polyhedra_));
|
||||
}
|
||||
return topoPtr_();
|
||||
}
|
||||
|
||||
@ -2,7 +2,7 @@
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration | Website: https://openfoam.org
|
||||
\\ / A nd | Copyright (C) 2011-2021 OpenFOAM Foundation
|
||||
\\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
|
||||
\\/ M anipulation |
|
||||
-------------------------------------------------------------------------------
|
||||
License
|
||||
@ -35,16 +35,27 @@ Description
|
||||
#include "Swap.H"
|
||||
#include "polygonTriangulate.H"
|
||||
|
||||
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
||||
// * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * * * //
|
||||
|
||||
bool Foam::vtkTopo::decomposePoly = true;
|
||||
namespace Foam
|
||||
{
|
||||
template<>
|
||||
const char*
|
||||
NamedEnum<vtkTopo::vtkPolyhedra, 3>::names[] =
|
||||
{"none", "polyhedra", "all"};
|
||||
}
|
||||
|
||||
const
|
||||
Foam::NamedEnum<Foam::vtkTopo::vtkPolyhedra, 3>
|
||||
Foam::vtkTopo::vtkPolyhedraNames_;
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||
|
||||
Foam::vtkTopo::vtkTopo(const polyMesh& mesh)
|
||||
Foam::vtkTopo::vtkTopo(const polyMesh& mesh, const vtkPolyhedra& polyhedra)
|
||||
:
|
||||
mesh_(mesh),
|
||||
polyhedra_(polyhedra),
|
||||
vertLabels_(),
|
||||
cellTypes_(),
|
||||
addPointCellLabels_(),
|
||||
@ -59,18 +70,13 @@ Foam::vtkTopo::vtkTopo(const polyMesh& mesh)
|
||||
|
||||
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;
|
||||
|
||||
// face owner is needed to determine the face orientation
|
||||
const labelList& owner = mesh.faceOwner();
|
||||
// Number of additional points and cells needed by the decomposition of
|
||||
// polyhedra
|
||||
label nAddPoints = 0, nAddCells = 0;
|
||||
|
||||
// Scan for cells which need to be decomposed and count additional points
|
||||
// and cells
|
||||
if (decomposePoly)
|
||||
if (polyhedra_ == vtkPolyhedra::none)
|
||||
{
|
||||
forAll(cellShapes, celli)
|
||||
{
|
||||
@ -108,7 +114,6 @@ Foam::vtkTopo::vtkTopo(const polyMesh& mesh)
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
// Set size of additional point addressing array
|
||||
// (from added point to original cell)
|
||||
addPointCellLabels_.setSize(nAddPoints);
|
||||
@ -126,9 +131,9 @@ Foam::vtkTopo::vtkTopo(const polyMesh& mesh)
|
||||
// Create a triangulation engine
|
||||
polygonTriangulate triEngine;
|
||||
|
||||
// Set counters for additional points and additional cells
|
||||
// Construct vtk cells and types
|
||||
label addPointi = 0, addCelli = 0;
|
||||
|
||||
const labelList& faceOwner = mesh.faceOwner();
|
||||
forAll(cellShapes, celli)
|
||||
{
|
||||
const cellShape& cellShape = cellShapes[celli];
|
||||
@ -136,19 +141,19 @@ Foam::vtkTopo::vtkTopo(const polyMesh& mesh)
|
||||
|
||||
labelList& vtkVerts = vertLabels_[celli];
|
||||
|
||||
if (cellModel == tet)
|
||||
if (polyhedra_ != vtkPolyhedra::all && cellModel == tet)
|
||||
{
|
||||
vtkVerts = cellShape;
|
||||
|
||||
cellTypes_[celli] = VTK_TETRA;
|
||||
}
|
||||
else if (cellModel == pyr)
|
||||
else if (polyhedra_ != vtkPolyhedra::all && cellModel == pyr)
|
||||
{
|
||||
vtkVerts = cellShape;
|
||||
|
||||
cellTypes_[celli] = VTK_PYRAMID;
|
||||
}
|
||||
else if (cellModel == prism)
|
||||
else if (polyhedra_ != vtkPolyhedra::all && cellModel == prism)
|
||||
{
|
||||
// VTK has a different node order for VTK_WEDGE
|
||||
// their triangles point outwards!
|
||||
@ -159,7 +164,7 @@ Foam::vtkTopo::vtkTopo(const polyMesh& mesh)
|
||||
|
||||
cellTypes_[celli] = VTK_WEDGE;
|
||||
}
|
||||
else if (cellModel == tetWedge && decomposePoly)
|
||||
else if (polyhedra_ == vtkPolyhedra::none && cellModel == tetWedge)
|
||||
{
|
||||
// Treat as squeezed prism (VTK_WEDGE)
|
||||
vtkVerts.setSize(6);
|
||||
@ -172,7 +177,7 @@ Foam::vtkTopo::vtkTopo(const polyMesh& mesh)
|
||||
|
||||
cellTypes_[celli] = VTK_WEDGE;
|
||||
}
|
||||
else if (cellModel == wedge)
|
||||
else if (polyhedra_ != vtkPolyhedra::all && cellModel == wedge)
|
||||
{
|
||||
// Treat as squeezed hex
|
||||
vtkVerts.setSize(8);
|
||||
@ -187,13 +192,13 @@ Foam::vtkTopo::vtkTopo(const polyMesh& mesh)
|
||||
|
||||
cellTypes_[celli] = VTK_HEXAHEDRON;
|
||||
}
|
||||
else if (cellModel == hex)
|
||||
else if (polyhedra_ != vtkPolyhedra::all && cellModel == hex)
|
||||
{
|
||||
vtkVerts = cellShape;
|
||||
|
||||
cellTypes_[celli] = VTK_HEXAHEDRON;
|
||||
}
|
||||
else if (decomposePoly)
|
||||
else if (polyhedra_ == vtkPolyhedra::none)
|
||||
{
|
||||
// Polyhedral cell. Decompose into tets + pyramids.
|
||||
|
||||
@ -210,7 +215,7 @@ Foam::vtkTopo::vtkTopo(const polyMesh& mesh)
|
||||
forAll(cFaces, cFacei)
|
||||
{
|
||||
const face& f = mesh_.faces()[cFaces[cFacei]];
|
||||
const bool isOwner = (owner[cFaces[cFacei]] == celli);
|
||||
const bool isOwner = faceOwner[cFaces[cFacei]] == celli;
|
||||
|
||||
// Do decomposition into triFcs and quadFcs.
|
||||
faceList triFcs(f.size() == 4 ? 0 : f.nTriangles());
|
||||
@ -347,7 +352,7 @@ Foam::vtkTopo::vtkTopo(const polyMesh& mesh)
|
||||
forAll(cFaces, cFacei)
|
||||
{
|
||||
const face& f = mesh.faces()[cFaces[cFacei]];
|
||||
const bool isOwner = (owner[cFaces[cFacei]] == celli);
|
||||
const bool isOwner = faceOwner[cFaces[cFacei]] == celli;
|
||||
|
||||
// number of labels for this face
|
||||
vtkVerts[nData++] = f.size();
|
||||
@ -372,7 +377,7 @@ Foam::vtkTopo::vtkTopo(const polyMesh& mesh)
|
||||
}
|
||||
}
|
||||
|
||||
if (decomposePoly)
|
||||
if (polyhedra_ == vtkPolyhedra::none)
|
||||
{
|
||||
Pout<< " Original cells:" << mesh_.nCells()
|
||||
<< " points:" << mesh_.nPoints()
|
||||
@ -380,7 +385,6 @@ Foam::vtkTopo::vtkTopo(const polyMesh& mesh)
|
||||
<< " additional points:" << addPointCellLabels_.size()
|
||||
<< nl << endl;
|
||||
}
|
||||
|
||||
}
|
||||
|
||||
|
||||
|
||||
@ -2,7 +2,7 @@
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration | Website: https://openfoam.org
|
||||
\\ / A nd | Copyright (C) 2011-2019 OpenFOAM Foundation
|
||||
\\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
|
||||
\\/ M anipulation |
|
||||
-------------------------------------------------------------------------------
|
||||
License
|
||||
@ -51,24 +51,20 @@ class polyMesh;
|
||||
|
||||
class vtkTopo
|
||||
{
|
||||
// Private Data
|
||||
|
||||
const polyMesh& mesh_;
|
||||
|
||||
//- Vertices per cell (including added cells) in vtk ordering
|
||||
labelListList vertLabels_;
|
||||
|
||||
//- Cell types (including added cells) in vtk numbering
|
||||
labelList cellTypes_;
|
||||
|
||||
labelList addPointCellLabels_;
|
||||
|
||||
labelList superCells_;
|
||||
|
||||
|
||||
public:
|
||||
|
||||
// Public static data
|
||||
// Public Enumerations
|
||||
|
||||
//- Groups of cell types retain as polyhedra
|
||||
enum class vtkPolyhedra
|
||||
{
|
||||
none,
|
||||
polyhedra,
|
||||
all
|
||||
};
|
||||
|
||||
//- Names of groups of cell types retain as polyhedra
|
||||
static const NamedEnum<vtkPolyhedra, 3> vtkPolyhedraNames_;
|
||||
|
||||
//- Equivalent to enumeration in "vtkCellType.h"
|
||||
enum vtkTypes
|
||||
@ -76,7 +72,6 @@ public:
|
||||
VTK_TRIANGLE = 5,
|
||||
VTK_POLYGON = 7,
|
||||
VTK_QUAD = 9,
|
||||
|
||||
VTK_TETRA = 10,
|
||||
VTK_HEXAHEDRON = 12,
|
||||
VTK_WEDGE = 13,
|
||||
@ -84,14 +79,36 @@ public:
|
||||
VTK_POLYHEDRON = 42
|
||||
};
|
||||
|
||||
//- Enable/disable polyhedron decomposition. Default = true
|
||||
static bool decomposePoly;
|
||||
|
||||
private:
|
||||
|
||||
// Private Data
|
||||
|
||||
//- Reference to the mesh
|
||||
const polyMesh& mesh_;
|
||||
|
||||
//- Cell types to retain as polyhedra
|
||||
const vtkPolyhedra polyhedra_;
|
||||
|
||||
//- Vertices per cell (including added cells) in vtk ordering
|
||||
labelListList vertLabels_;
|
||||
|
||||
//- Cell types (including added cells) in vtk numbering
|
||||
labelList cellTypes_;
|
||||
|
||||
//- ...
|
||||
labelList addPointCellLabels_;
|
||||
|
||||
//- ...
|
||||
labelList superCells_;
|
||||
|
||||
|
||||
public:
|
||||
|
||||
// Constructors
|
||||
|
||||
//- Construct from components
|
||||
vtkTopo(const polyMesh&);
|
||||
vtkTopo(const polyMesh& mesh, const vtkPolyhedra& polyhedra);
|
||||
|
||||
//- Disallow default bitwise copy construction
|
||||
vtkTopo(const vtkTopo&) = delete;
|
||||
|
||||
Reference in New Issue
Block a user