From bc493180d98e9cf069e582befe9c4dc4ce611df5 Mon Sep 17 00:00:00 2001
From: Henry Weller <http://cfd.direct>
Date: Wed, 17 Jun 2015 08:46:46 +0100
Subject: [PATCH] rhoCentralFoam: Ensure fixed value boundary conditions are
 preserved Resolves bug-report
 http://www.openfoam.org/mantisbt/view.php?id=1748

---
 .../rhoCentralFoam/createFields.H             | 22 +++++++++++----
 .../rhoCentralFoam/directionInterpolate.H     | 27 ++++++++++++++++++-
 .../rhoCentralFoam/rhoBoundaryTypes.H         | 14 ----------
 .../rhoCentralDyMFoam/rhoCentralDyMFoam.C     |  6 ++---
 .../rhoCentralFoam/rhoCentralFoam.C           |  6 ++---
 5 files changed, 49 insertions(+), 26 deletions(-)
 delete mode 100644 applications/solvers/compressible/rhoCentralFoam/rhoBoundaryTypes.H

diff --git a/applications/solvers/compressible/rhoCentralFoam/createFields.H b/applications/solvers/compressible/rhoCentralFoam/createFields.H
index c1618d6346a..425ba78b877 100644
--- a/applications/solvers/compressible/rhoCentralFoam/createFields.H
+++ b/applications/solvers/compressible/rhoCentralFoam/createFields.H
@@ -32,7 +32,6 @@ volVectorField U
     mesh
 );
 
-#include "rhoBoundaryTypes.H"
 volScalarField rho
 (
     IOobject
@@ -44,7 +43,7 @@ volScalarField rho
         IOobject::AUTO_WRITE
     ),
     thermo.rho(),
-    rhoBoundaryTypes
+    derivedPatchFieldTypes(p)
 );
 
 volVectorField rhoU
@@ -57,7 +56,8 @@ volVectorField rhoU
         IOobject::NO_READ,
         IOobject::NO_WRITE
     ),
-    rho*U
+    rho*U,
+    derivedPatchFieldTypes(U)
 );
 
 volScalarField rhoE
@@ -70,7 +70,8 @@ volScalarField rhoE
         IOobject::NO_READ,
         IOobject::NO_WRITE
     ),
-    rho*(e + 0.5*magSqr(U))
+    rho*(e + 0.5*magSqr(U)),
+    derivedPatchFieldTypes(T)
 );
 
 surfaceScalarField pos
@@ -98,7 +99,18 @@ surfaceScalarField neg
 );
 
 
-surfaceScalarField phi("phi", mesh.Sf() & fvc::interpolate(rhoU));
+surfaceScalarField phi
+(
+    IOobject
+    (
+        "phi",
+        runTime.timeName(),
+        mesh,
+        IOobject::NO_READ,
+        IOobject::AUTO_WRITE
+    ),
+    mesh.Sf() & fvc::interpolate(rhoU)
+);
 
 Info<< "Creating turbulence model\n" << endl;
 autoPtr<compressible::turbulenceModel> turbulence
diff --git a/applications/solvers/compressible/rhoCentralFoam/directionInterpolate.H b/applications/solvers/compressible/rhoCentralFoam/directionInterpolate.H
index d818f2f0e5f..a223ecb5032 100644
--- a/applications/solvers/compressible/rhoCentralFoam/directionInterpolate.H
+++ b/applications/solvers/compressible/rhoCentralFoam/directionInterpolate.H
@@ -31,8 +31,9 @@ tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > interpolate
     {
         if
         (
-           !sf.boundaryField()[patchi].coupled()
+            !sf.boundaryField()[patchi].coupled()
          && sf.boundaryField()[patchi].size()
+         && !vf.boundaryField()[patchi].fixesValue()
          && dir.boundaryField()[patchi][0] > 0
         )
         {
@@ -44,4 +45,28 @@ tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > interpolate
     return tsf;
 }
 
+
+template<class Type>
+wordList derivedPatchFieldTypes
+(
+    const GeometricField<Type, fvPatchField, volMesh>& vf
+)
+{
+    wordList phiTypes
+    (
+        vf.boundaryField().size(),
+        calculatedFvPatchField<Type>::typeName
+    );
+
+    forAll(vf.boundaryField(), patchi)
+    {
+        if (vf.boundaryField()[patchi].fixesValue())
+        {
+            phiTypes[patchi] = fixedValueFvPatchField<Type>::typeName;
+        }
+    }
+
+    return phiTypes;
+}
+
 }
diff --git a/applications/solvers/compressible/rhoCentralFoam/rhoBoundaryTypes.H b/applications/solvers/compressible/rhoCentralFoam/rhoBoundaryTypes.H
deleted file mode 100644
index 7118d24fc92..00000000000
--- a/applications/solvers/compressible/rhoCentralFoam/rhoBoundaryTypes.H
+++ /dev/null
@@ -1,14 +0,0 @@
-const volScalarField::GeometricBoundaryField& pbf = p.boundaryField();
-wordList rhoBoundaryTypes = pbf.types();
-
-forAll(rhoBoundaryTypes, patchi)
-{
-    if (rhoBoundaryTypes[patchi] == "waveTransmissive")
-    {
-        rhoBoundaryTypes[patchi] = zeroGradientFvPatchScalarField::typeName;
-    }
-    else if (pbf[patchi].fixesValue())
-    {
-        rhoBoundaryTypes[patchi] = fixedRhoFvPatchScalarField::typeName;
-    }
-}
diff --git a/applications/solvers/compressible/rhoCentralFoam/rhoCentralDyMFoam/rhoCentralDyMFoam.C b/applications/solvers/compressible/rhoCentralFoam/rhoCentralDyMFoam/rhoCentralDyMFoam.C
index c2633b382fa..114d97fbf10 100644
--- a/applications/solvers/compressible/rhoCentralFoam/rhoCentralDyMFoam/rhoCentralDyMFoam.C
+++ b/applications/solvers/compressible/rhoCentralFoam/rhoCentralDyMFoam/rhoCentralDyMFoam.C
@@ -187,7 +187,7 @@ int main(int argc, char *argv[])
             rhoU.dimensionedInternalField()
            /rho.dimensionedInternalField();
         U.correctBoundaryConditions();
-        rhoU.boundaryField() = rho.boundaryField()*U.boundaryField();
+        rhoU.boundaryField() == rho.boundaryField()*U.boundaryField();
 
         if (!inviscid)
         {
@@ -221,7 +221,7 @@ int main(int argc, char *argv[])
         e = rhoE/rho - 0.5*magSqr(U);
         e.correctBoundaryConditions();
         thermo.correct();
-        rhoE.boundaryField() =
+        rhoE.boundaryField() ==
             rho.boundaryField()*
             (
                 e.boundaryField() + 0.5*magSqr(U.boundaryField())
@@ -242,7 +242,7 @@ int main(int argc, char *argv[])
             rho.dimensionedInternalField()
            /psi.dimensionedInternalField();
         p.correctBoundaryConditions();
-        rho.boundaryField() = psi.boundaryField()*p.boundaryField();
+        rho.boundaryField() == psi.boundaryField()*p.boundaryField();
 
         turbulence->correct();
 
diff --git a/applications/solvers/compressible/rhoCentralFoam/rhoCentralFoam.C b/applications/solvers/compressible/rhoCentralFoam/rhoCentralFoam.C
index 05a4bdfb72d..0b5cba1c78a 100644
--- a/applications/solvers/compressible/rhoCentralFoam/rhoCentralFoam.C
+++ b/applications/solvers/compressible/rhoCentralFoam/rhoCentralFoam.C
@@ -169,7 +169,7 @@ int main(int argc, char *argv[])
             rhoU.dimensionedInternalField()
            /rho.dimensionedInternalField();
         U.correctBoundaryConditions();
-        rhoU.boundaryField() = rho.boundaryField()*U.boundaryField();
+        rhoU.boundaryField() == rho.boundaryField()*U.boundaryField();
 
         if (!inviscid)
         {
@@ -203,7 +203,7 @@ int main(int argc, char *argv[])
         e = rhoE/rho - 0.5*magSqr(U);
         e.correctBoundaryConditions();
         thermo.correct();
-        rhoE.boundaryField() =
+        rhoE.boundaryField() ==
             rho.boundaryField()*
             (
                 e.boundaryField() + 0.5*magSqr(U.boundaryField())
@@ -224,7 +224,7 @@ int main(int argc, char *argv[])
             rho.dimensionedInternalField()
            /psi.dimensionedInternalField();
         p.correctBoundaryConditions();
-        rho.boundaryField() = psi.boundaryField()*p.boundaryField();
+        rho.boundaryField() == psi.boundaryField()*p.boundaryField();
 
         turbulence->correct();
 
-- 
GitLab