From d10769f5a9cb606c67ac08a5b16b9709411cf783 Mon Sep 17 00:00:00 2001
From: mattijs <mattijs>
Date: Wed, 10 Mar 2021 17:55:15 +0000
Subject: [PATCH] ENH: checkMesh: write surface fields. Fixes #2023

---
 .../mesh/manipulation/checkMesh/checkMesh.C   | 22 ++++-
 .../mesh/manipulation/checkMesh/writeFields.C | 80 ++++++++++++++++++-
 .../mesh/manipulation/checkMesh/writeFields.H |  3 +-
 3 files changed, 100 insertions(+), 5 deletions(-)

diff --git a/applications/utilities/mesh/manipulation/checkMesh/checkMesh.C b/applications/utilities/mesh/manipulation/checkMesh/checkMesh.C
index ec1a2314b1f..c68bf78709b 100644
--- a/applications/utilities/mesh/manipulation/checkMesh/checkMesh.C
+++ b/applications/utilities/mesh/manipulation/checkMesh/checkMesh.C
@@ -6,7 +6,7 @@
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
     Copyright (C) 2011-2017 OpenFOAM Foundation
-    Copyright (C) 2015-2020 OpenCFD Ltd.
+    Copyright (C) 2015-2021 OpenCFD Ltd.
 -------------------------------------------------------------------------------
 License
     This file is part of OpenFOAM.
@@ -58,6 +58,9 @@ Usage
     \param -writeAllFields \n
     Writes all mesh quality measures as fields.
 
+    \param -writeAllSurfaceFields \n
+    Adds writing of surface fields when used in combination with writeAllFields.
+
     \param -writeFields '(\<fieldName\>)' \n
     Writes selected mesh quality measures as fields.
 
@@ -111,6 +114,11 @@ int main(int argc, char *argv[])
         "writeAllFields",
         "Write volFields with mesh quality parameters"
     );
+    argList::addBoolOption
+    (
+        "writeAllSurfaceFields",
+        "Write surfaceFields with mesh quality parameters"
+    );
     argList::addOption
     (
         "writeFields",
@@ -163,6 +171,7 @@ int main(int argc, char *argv[])
         "faceZone"
     });
 
+    const bool writeFaceFields = args.found("writeAllSurfaceFields");
     wordHashSet selectedFields;
     if (args.found("writeFields"))
     {
@@ -182,6 +191,13 @@ int main(int argc, char *argv[])
     {
         selectedFields = allFields;
     }
+    else if (writeFaceFields)
+    {
+        FatalErrorInFunction
+            << "Option 'writeAllSurfaceFields' only valid in combination"
+            << " with 'writeFields' or 'writeAllFields'"
+            << nl << exit(FatalError);
+    }
 
 
     if (noTopology)
@@ -306,7 +322,7 @@ int main(int argc, char *argv[])
 
 
             // Write selected fields
-            Foam::writeFields(mesh, selectedFields);
+            Foam::writeFields(mesh, selectedFields, writeFaceFields);
         }
         else if (state == polyMesh::POINTS_MOVED)
         {
@@ -338,7 +354,7 @@ int main(int argc, char *argv[])
 
 
             // Write selected fields
-            Foam::writeFields(mesh, selectedFields);
+            Foam::writeFields(mesh, selectedFields, writeFaceFields);
         }
     }
 
diff --git a/applications/utilities/mesh/manipulation/checkMesh/writeFields.C b/applications/utilities/mesh/manipulation/checkMesh/writeFields.C
index 92e6e73f475..e9c4956ad3c 100644
--- a/applications/utilities/mesh/manipulation/checkMesh/writeFields.C
+++ b/applications/utilities/mesh/manipulation/checkMesh/writeFields.C
@@ -117,10 +117,42 @@ void minFaceToCell
 }
 
 
+void writeSurfaceField
+(
+    const fvMesh& mesh,
+    const fileName& fName,
+    const scalarField& faceData
+)
+{
+    // Write single surfaceScalarField
+
+    surfaceScalarField fld
+    (
+        IOobject
+        (
+            fName,
+            mesh.time().timeName(),
+            mesh,
+            IOobject::NO_READ,
+            IOobject::AUTO_WRITE,
+            false
+        ),
+        mesh,
+        dimensionedScalar(dimless, Zero),
+        calculatedFvsPatchScalarField::typeName
+    );
+    fld.primitiveFieldRef() = faceData;
+    //fld.correctBoundaryConditions();
+    Info<< "    Writing face data to " << fName << endl;
+    fld.write();
+}
+
+
 void Foam::writeFields
 (
     const fvMesh& mesh,
-    const wordHashSet& selectedFields
+    const wordHashSet& selectedFields,
+    const bool writeFaceFields
 )
 {
     if (selectedFields.empty())
@@ -173,6 +205,16 @@ void Foam::writeFields
         Info<< "    Writing non-orthogonality (angle) to "
             << cellNonOrthoAngle.name() << endl;
         cellNonOrthoAngle.write();
+
+        if (writeFaceFields)
+        {
+            writeSurfaceField
+            (
+                mesh,
+                "face_nonOrthoAngle",
+                SubField<scalar>(nonOrthoAngle, mesh.nInternalFaces())
+            );
+        }
     }
 
     if (selectedFields.found("faceWeight"))
@@ -202,6 +244,11 @@ void Foam::writeFields
         Info<< "    Writing face interpolation weights (0..0.5) to "
             << cellWeights.name() << endl;
         cellWeights.write();
+
+        if (writeFaceFields)
+        {
+            writeSurfaceField(mesh, "face_faceWeight", mesh.weights());
+        }
     }
 
 
@@ -243,6 +290,16 @@ void Foam::writeFields
         maxFaceToCell(faceSkewness, cellSkewness);
         Info<< "    Writing face skewness to " << cellSkewness.name() << endl;
         cellSkewness.write();
+
+        if (writeFaceFields)
+        {
+            writeSurfaceField
+            (
+                mesh,
+                "face_skewness",
+                SubField<scalar>(faceSkewness, mesh.nInternalFaces())
+            );
+        }
     }
 
 
@@ -431,6 +488,16 @@ void Foam::writeFields
         Info<< "    Writing cell volume ratio to "
             << cellVolumeRatio.name() << endl;
         cellVolumeRatio.write();
+
+        if (writeFaceFields)
+        {
+            writeSurfaceField
+            (
+                mesh,
+                "face_cellVolumeRatio",
+                SubField<scalar>(faceVolumeRatio, mesh.nInternalFaces())
+            );
+        }
     }
 
     // minTetVolume
@@ -560,6 +627,17 @@ void Foam::writeFields
         minPyrVolume.correctBoundaryConditions();
         Info<< "    Writing minPyrVolume to " << minPyrVolume.name() << endl;
         minPyrVolume.write();
+
+        if (writeFaceFields)
+        {
+            scalarField minFacePyrVol(neiPyrVol);
+            minFacePyrVol = min
+            (
+                minFacePyrVol,
+                SubField<scalar>(ownPyrVol, mesh.nInternalFaces())
+            );
+            writeSurfaceField(mesh, "face_minPyrVolume", minFacePyrVol);
+        }
     }
 
     if (selectedFields.found("cellRegion"))
diff --git a/applications/utilities/mesh/manipulation/checkMesh/writeFields.H b/applications/utilities/mesh/manipulation/checkMesh/writeFields.H
index 675784c40d7..63e42eaffbb 100644
--- a/applications/utilities/mesh/manipulation/checkMesh/writeFields.H
+++ b/applications/utilities/mesh/manipulation/checkMesh/writeFields.H
@@ -5,6 +5,7 @@ namespace Foam
     void writeFields
     (
         const fvMesh& mesh,
-        const wordHashSet& selectedFields
+        const wordHashSet& selectedFields,
+        const bool writeFaceFields
     );
 }
-- 
GitLab