mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
ENH: Added tet volume check to checkMesh and snappyHexMesh
This commit is contained in:
@ -26,6 +26,7 @@ License
|
||||
|
||||
#include "polyMeshGeometry.H"
|
||||
#include "pyramidPointFaceRef.H"
|
||||
#include "tetPointRef.H"
|
||||
#include "syncTools.H"
|
||||
#include "unitConversion.H"
|
||||
|
||||
@ -308,6 +309,59 @@ Foam::scalar Foam::polyMeshGeometry::calcSkewness
|
||||
}
|
||||
|
||||
|
||||
// Create the neighbour pyramid - it will have positive volume
|
||||
bool Foam::polyMeshGeometry::checkFaceTet
|
||||
(
|
||||
const polyMesh& mesh,
|
||||
const bool report,
|
||||
const scalar minTetVol,
|
||||
const pointField& p,
|
||||
const label faceI,
|
||||
const point& fc, // face centre
|
||||
const point& cc, // cell centre
|
||||
|
||||
labelHashSet* setPtr
|
||||
)
|
||||
{
|
||||
const face& f = mesh.faces()[faceI];
|
||||
|
||||
forAll(f, fp)
|
||||
{
|
||||
scalar tetVol = tetPointRef
|
||||
(
|
||||
p[f[fp]],
|
||||
p[f.nextLabel(fp)],
|
||||
fc,
|
||||
cc
|
||||
).mag();
|
||||
|
||||
if (tetVol < minTetVol)
|
||||
{
|
||||
if (report)
|
||||
{
|
||||
Pout<< "bool polyMeshGeometry::checkFaceTets("
|
||||
<< "const bool, const scalar, const pointField&"
|
||||
<< ", const pointField&, const labelList&,"
|
||||
<< " labelHashSet*): "
|
||||
<< "face " << faceI
|
||||
<< " has a triangle that points the wrong way."
|
||||
<< endl
|
||||
<< "Tet volume: " << tetVol
|
||||
<< " Face " << faceI
|
||||
<< endl;
|
||||
}
|
||||
|
||||
if (setPtr)
|
||||
{
|
||||
setPtr->insert(faceI);
|
||||
}
|
||||
return true;
|
||||
}
|
||||
}
|
||||
return false;
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||
|
||||
// Construct from components
|
||||
@ -719,6 +773,139 @@ bool Foam::polyMeshGeometry::checkFacePyramids
|
||||
}
|
||||
|
||||
|
||||
bool Foam::polyMeshGeometry::checkFaceTets
|
||||
(
|
||||
const bool report,
|
||||
const scalar minTetVol,
|
||||
const polyMesh& mesh,
|
||||
const vectorField& cellCentres,
|
||||
const vectorField& faceCentres,
|
||||
const pointField& p,
|
||||
const labelList& checkFaces,
|
||||
const List<labelPair>& baffles,
|
||||
labelHashSet* setPtr
|
||||
)
|
||||
{
|
||||
// check whether face area vector points to the cell with higher label
|
||||
const labelList& own = mesh.faceOwner();
|
||||
const labelList& nei = mesh.faceNeighbour();
|
||||
|
||||
label nErrorPyrs = 0;
|
||||
|
||||
forAll(checkFaces, i)
|
||||
{
|
||||
label faceI = checkFaces[i];
|
||||
|
||||
// Create the owner pyramid - note: exchange cell and face centre
|
||||
// to get positive volume.
|
||||
bool tetError = checkFaceTet
|
||||
(
|
||||
mesh,
|
||||
report,
|
||||
minTetVol,
|
||||
p,
|
||||
faceI,
|
||||
cellCentres[own[faceI]], // face centre
|
||||
faceCentres[faceI], // cell centre
|
||||
setPtr
|
||||
);
|
||||
|
||||
if (tetError)
|
||||
{
|
||||
nErrorPyrs++;
|
||||
}
|
||||
|
||||
if (mesh.isInternalFace(faceI))
|
||||
{
|
||||
// Create the neighbour pyramid - it will have positive volume
|
||||
bool tetError = checkFaceTet
|
||||
(
|
||||
mesh,
|
||||
report,
|
||||
minTetVol,
|
||||
p,
|
||||
faceI,
|
||||
faceCentres[faceI], // face centre
|
||||
cellCentres[nei[faceI]], // cell centre
|
||||
setPtr
|
||||
);
|
||||
if (tetError)
|
||||
{
|
||||
nErrorPyrs++;
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
forAll(baffles, i)
|
||||
{
|
||||
label face0 = baffles[i].first();
|
||||
label face1 = baffles[i].second();
|
||||
|
||||
bool tetError = checkFaceTet
|
||||
(
|
||||
mesh,
|
||||
report,
|
||||
minTetVol,
|
||||
p,
|
||||
face0,
|
||||
cellCentres[own[face0]], // face centre
|
||||
faceCentres[face0], // cell centre
|
||||
setPtr
|
||||
);
|
||||
|
||||
if (tetError)
|
||||
{
|
||||
nErrorPyrs++;
|
||||
}
|
||||
|
||||
// Create the neighbour pyramid - it will have positive volume
|
||||
tetError = checkFaceTet
|
||||
(
|
||||
mesh,
|
||||
report,
|
||||
minTetVol,
|
||||
p,
|
||||
face0,
|
||||
faceCentres[face0], // face centre
|
||||
cellCentres[own[face1]], // cell centre
|
||||
setPtr
|
||||
);
|
||||
|
||||
if (tetError)
|
||||
{
|
||||
nErrorPyrs++;
|
||||
}
|
||||
}
|
||||
|
||||
reduce(nErrorPyrs, sumOp<label>());
|
||||
|
||||
if (nErrorPyrs > 0)
|
||||
{
|
||||
if (report)
|
||||
{
|
||||
SeriousErrorIn
|
||||
(
|
||||
"polyMeshGeometry::checkFaceTets("
|
||||
"const bool, const scalar, const pointField&, const pointField&"
|
||||
", const labelList&, labelHashSet*)"
|
||||
) << "Error in face pyramids: faces pointing the wrong way!"
|
||||
<< endl;
|
||||
}
|
||||
|
||||
return true;
|
||||
}
|
||||
else
|
||||
{
|
||||
if (report)
|
||||
{
|
||||
Info<< "Face tets OK.\n" << endl;
|
||||
}
|
||||
|
||||
return false;
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
bool Foam::polyMeshGeometry::checkFaceSkewness
|
||||
(
|
||||
const bool report,
|
||||
@ -1946,6 +2133,31 @@ bool Foam::polyMeshGeometry::checkFacePyramids
|
||||
}
|
||||
|
||||
|
||||
bool Foam::polyMeshGeometry::checkFaceTets
|
||||
(
|
||||
const bool report,
|
||||
const scalar minTetVol,
|
||||
const pointField& p,
|
||||
const labelList& checkFaces,
|
||||
const List<labelPair>& baffles,
|
||||
labelHashSet* setPtr
|
||||
) const
|
||||
{
|
||||
return checkFaceTets
|
||||
(
|
||||
report,
|
||||
minTetVol,
|
||||
mesh_,
|
||||
cellCentres_,
|
||||
faceCentres_,
|
||||
p,
|
||||
checkFaces,
|
||||
baffles,
|
||||
setPtr
|
||||
);
|
||||
}
|
||||
|
||||
|
||||
bool Foam::polyMeshGeometry::checkFaceSkewness
|
||||
(
|
||||
const bool report,
|
||||
|
||||
Reference in New Issue
Block a user