diff --git a/applications/utilities/mesh/manipulation/checkMesh/checkMesh.C b/applications/utilities/mesh/manipulation/checkMesh/checkMesh.C index d02cd5f6edccb342c3ec61391d1f7f49f58215ad..25355afd6fb1ef1f8987db30df0929ff21acc429 100644 --- a/applications/utilities/mesh/manipulation/checkMesh/checkMesh.C +++ b/applications/utilities/mesh/manipulation/checkMesh/checkMesh.C @@ -149,6 +149,7 @@ int main(int argc, char *argv[]) selectedFields.insert("cellShapes"); selectedFields.insert("cellVolume"); selectedFields.insert("cellVolumeRatio"); + selectedFields.insert("minTetVolume"); } diff --git a/applications/utilities/mesh/manipulation/checkMesh/writeFields.C b/applications/utilities/mesh/manipulation/checkMesh/writeFields.C index d8c23f33a3a01a6e8fd1189020c4a9416b244031..44e39e7a56416f13cf32e777e3a29889f85647f3 100644 --- a/applications/utilities/mesh/manipulation/checkMesh/writeFields.C +++ b/applications/utilities/mesh/manipulation/checkMesh/writeFields.C @@ -3,6 +3,7 @@ #include "polyMeshTools.H" #include "zeroGradientFvPatchFields.H" #include "syncTools.H" +#include "tetPointRef.H" using namespace Foam; @@ -273,6 +274,7 @@ void Foam::writeFields // cell type // ~~~~~~~~~ + if (selectedFields.found("cellShapes")) { volScalarField shape @@ -356,5 +358,72 @@ void Foam::writeFields cellVolumeRatio.write(); } + // minTetVolume + if (selectedFields.found("minTetVolume")) + { + volScalarField minTetVolume + ( + IOobject + ( + "minTetVolume", + mesh.time().timeName(), + mesh, + IOobject::NO_READ, + IOobject::AUTO_WRITE, + false + ), + mesh, + dimensionedScalar("minTetVolume", dimless, GREAT), + zeroGradientFvPatchScalarField::typeName + ); + + + const labelList& own = mesh.faceOwner(); + const labelList& nei = mesh.faceNeighbour(); + const pointField& p = mesh.points(); + forAll(own, facei) + { + const face& f = mesh.faces()[facei]; + const point& fc = mesh.faceCentres()[facei]; + + { + const point& ownCc = mesh.cellCentres()[own[facei]]; + scalar& ownVol = minTetVolume[own[facei]]; + forAll(f, fp) + { + scalar tetQual = tetPointRef + ( + p[f[fp]], + p[f.nextLabel(fp)], + ownCc, + fc + ).quality(); + ownVol = min(ownVol, tetQual); + } + } + if (mesh.isInternalFace(facei)) + { + const point& neiCc = mesh.cellCentres()[nei[facei]]; + scalar& neiVol = minTetVolume[nei[facei]]; + forAll(f, fp) + { + scalar tetQual = tetPointRef + ( + p[f[fp]], + p[f.nextLabel(fp)], + fc, + neiCc + ).quality(); + neiVol = min(neiVol, tetQual); + } + } + } + + minTetVolume.correctBoundaryConditions(); + Info<< " Writing minTetVolume to " << minTetVolume.name() << endl; + minTetVolume.write(); + } + + Info<< endl; }