From 673e0d1704cdc1660238dfba6121b03e28ed0fff Mon Sep 17 00:00:00 2001 From: Henry Weller Date: Sat, 23 Apr 2016 10:04:39 +0100 Subject: [PATCH] fireFoam: Added optional hydrostatic initialization of the pressure and density Also added the new prghTotalHydrostaticPressure p_rgh BC which uses the hydrostatic pressure field as the reference state for the far-field which provides much more accurate entrainment is large open domains typical of many fire simulations. The hydrostatic field solution is controlled by the optional entries in the fvSolution.PIMPLE dictionary, e.g. hydrostaticInitialization yes; nHydrostaticCorrectors 5; and the solver must also be specified for the hydrostatic p_rgh field ph_rgh e.g. ph_rgh { $p_rgh; } Suitable boundary conditions for ph_rgh cannot always be derived from those for p_rgh and so the ph_rgh is read to provide them. To avoid accuracy issues with IO, restart and post-processing the p_rgh and ph_rgh the option to specify a suitable reference pressure is provided via the optional pRef file in the constant directory, e.g. dimensions [1 -1 -2 0 0 0 0]; value 101325; which is used in the relationship between p_rgh and p: p = p_rgh + rho*gh + pRef; Note that if pRef is specified all pressure BC specifications in the p_rgh and ph_rgh files are relative to the reference to avoid round-off errors. For examples of suitable BCs for p_rgh and ph_rgh for a range of fireFoam cases please study the tutorials in tutorials/combustion/fireFoam/les which have all been updated. Henry G. Weller CFD Direct Ltd. --- .../combustion/fireFoam/createFields.H | 101 ++++----- .../solvers/combustion/fireFoam/fireFoam.C | 1 - .../solvers/combustion/fireFoam/pEqn.H | 9 +- .../solvers/combustion/fireFoam/phrghEqn.H | 62 ++++++ src/finiteVolume/Make/files | 1 + .../cfdTools/general/include/readpRef.H | 13 ++ .../basic/fixedValue/fixedValueFvPatchField.H | 2 +- ...talHydrostaticPressureFvPatchScalarField.C | 168 +++++++++++++++ ...talHydrostaticPressureFvPatchScalarField.H | 202 ++++++++++++++++++ .../prghTotalPressureFvPatchScalarField.H | 20 +- .../flameSpreadWaterSuppressionPanel/0/p_rgh | 20 +- .../flameSpreadWaterSuppressionPanel/Allclean | 2 +- .../flameSpreadWaterSuppressionPanel/Allrun | 4 +- .../constant/pRef | 21 ++ .../system/fvSolution | 8 + .../les/oppositeBurningPanels/0/p_rgh | 9 +- .../les/oppositeBurningPanels/Allclean | 1 + .../fireFoam/les/oppositeBurningPanels/Allrun | 2 + .../les/oppositeBurningPanels/constant/pRef | 21 ++ .../oppositeBurningPanels/system/fvSolution | 9 + .../fireFoam/les/smallPoolFire2D/0/p_rgh | 9 +- .../fireFoam/les/smallPoolFire2D/Allclean | 11 + .../fireFoam/les/smallPoolFire2D/Allrun | 2 + .../les/smallPoolFire2D/constant/pRef | 21 ++ .../les/smallPoolFire2D/system/fvSolution | 8 + .../fireFoam/les/smallPoolFire3D/0/p_rgh | 8 +- .../fireFoam/les/smallPoolFire3D/Allclean | 11 + .../fireFoam/les/smallPoolFire3D/Allrun | 3 + .../les/smallPoolFire3D/constant/pRef | 21 ++ .../les/smallPoolFire3D/system/fvSolution | 8 + 30 files changed, 673 insertions(+), 105 deletions(-) create mode 100644 applications/solvers/combustion/fireFoam/phrghEqn.H create mode 100644 src/finiteVolume/cfdTools/general/include/readpRef.H create mode 100644 src/finiteVolume/fields/fvPatchFields/derived/prghTotalHydrostaticPressure/prghTotalHydrostaticPressureFvPatchScalarField.C create mode 100644 src/finiteVolume/fields/fvPatchFields/derived/prghTotalHydrostaticPressure/prghTotalHydrostaticPressureFvPatchScalarField.H create mode 100644 tutorials/combustion/fireFoam/les/flameSpreadWaterSuppressionPanel/constant/pRef create mode 100644 tutorials/combustion/fireFoam/les/oppositeBurningPanels/constant/pRef create mode 100755 tutorials/combustion/fireFoam/les/smallPoolFire2D/Allclean create mode 100644 tutorials/combustion/fireFoam/les/smallPoolFire2D/constant/pRef create mode 100755 tutorials/combustion/fireFoam/les/smallPoolFire3D/Allclean create mode 100644 tutorials/combustion/fireFoam/les/smallPoolFire3D/constant/pRef diff --git a/applications/solvers/combustion/fireFoam/createFields.H b/applications/solvers/combustion/fireFoam/createFields.H index c5fcc21c90..ce621258cc 100644 --- a/applications/solvers/combustion/fireFoam/createFields.H +++ b/applications/solvers/combustion/fireFoam/createFields.H @@ -54,6 +54,9 @@ volVectorField U #include "compressibleCreatePhi.H" +#include "createMRF.H" + + Info<< "Creating turbulence model\n" << endl; autoPtr turbulence ( @@ -69,6 +72,55 @@ autoPtr turbulence // Set the turbulence into the combustion model combustion->setTurbulence(turbulence()); + +#include "readGravitationalAcceleration.H" +#include "readhRef.H" +#include "gh.H" +#include "readpRef.H" + +volScalarField p_rgh +( + IOobject + ( + "p_rgh", + runTime.timeName(), + mesh, + IOobject::MUST_READ, + IOobject::AUTO_WRITE + ), + mesh +); + +mesh.setFluxRequired(p_rgh.name()); + +#include "phrghEqn.H" + + +multivariateSurfaceInterpolationScheme::fieldTable fields; + +forAll(Y, i) +{ + fields.add(Y[i]); +} +fields.add(thermo.he()); + +IOdictionary additionalControlsDict +( + IOobject + ( + "additionalControls", + runTime.constant(), + mesh, + IOobject::MUST_READ_IF_MODIFIED, + IOobject::NO_WRITE + ) +); + +Switch solvePrimaryRegion +( + additionalControlsDict.lookup("solvePrimaryRegion") +); + volScalarField dQ ( IOobject @@ -99,52 +151,3 @@ volScalarField dpdt Info<< "Creating field kinetic energy K\n" << endl; volScalarField K("K", 0.5*magSqr(U)); - - -#include "readGravitationalAcceleration.H" -#include "readhRef.H" -#include "gh.H" - - -volScalarField p_rgh -( - IOobject - ( - "p_rgh", - runTime.timeName(), - mesh, - IOobject::MUST_READ, - IOobject::AUTO_WRITE - ), - mesh -); - -// Force p_rgh to be consistent with p -p_rgh = p - rho*gh; - -mesh.setFluxRequired(p_rgh.name()); - -multivariateSurfaceInterpolationScheme::fieldTable fields; - -forAll(Y, i) -{ - fields.add(Y[i]); -} -fields.add(thermo.he()); - -IOdictionary additionalControlsDict -( - IOobject - ( - "additionalControls", - runTime.constant(), - mesh, - IOobject::MUST_READ_IF_MODIFIED, - IOobject::NO_WRITE - ) -); - -Switch solvePrimaryRegion -( - additionalControlsDict.lookup("solvePrimaryRegion") -); diff --git a/applications/solvers/combustion/fireFoam/fireFoam.C b/applications/solvers/combustion/fireFoam/fireFoam.C index fd8e2a4357..3e1689f656 100644 --- a/applications/solvers/combustion/fireFoam/fireFoam.C +++ b/applications/solvers/combustion/fireFoam/fireFoam.C @@ -54,7 +54,6 @@ int main(int argc, char *argv[]) pimpleControl pimple(mesh); #include "createFields.H" - #include "createMRF.H" #include "createFvOptions.H" #include "createClouds.H" #include "createSurfaceFilmModel.H" diff --git a/applications/solvers/combustion/fireFoam/pEqn.H b/applications/solvers/combustion/fireFoam/pEqn.H index 7b5249d57e..0d42787fff 100644 --- a/applications/solvers/combustion/fireFoam/pEqn.H +++ b/applications/solvers/combustion/fireFoam/pEqn.H @@ -4,7 +4,7 @@ volScalarField rAU(1.0/UEqn.A()); surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho*rAU)); volVectorField HbyA(constrainHbyA(rAU*UEqn.H(), U, p)); -surfaceScalarField phig(-rhorAUf*ghf*fvc::snGrad(rho)*mesh.magSf()); +surfaceScalarField phig("phig", -rhorAUf*ghf*fvc::snGrad(rho)*mesh.magSf()); surfaceScalarField phiHbyA ( @@ -25,9 +25,10 @@ while (pimple.correctNonOrthogonal()) { fvScalarMatrix p_rghEqn ( - fvc::ddt(psi, rho)*gh + fvm::ddt(psi, p_rgh) + + fvc::ddt(psi, rho)*gh + + fvc::ddt(psi)*pRef + fvc::div(phiHbyA) - + fvm::ddt(psi, p_rgh) - fvm::laplacian(rhorAUf, p_rgh) == parcels.Srho() @@ -46,7 +47,7 @@ while (pimple.correctNonOrthogonal()) } } -p = p_rgh + rho*gh; +p = p_rgh + rho*gh + pRef; #include "rhoEqn.H" #include "compressibleContinuityErrs.H" diff --git a/applications/solvers/combustion/fireFoam/phrghEqn.H b/applications/solvers/combustion/fireFoam/phrghEqn.H new file mode 100644 index 0000000000..3c5175f41f --- /dev/null +++ b/applications/solvers/combustion/fireFoam/phrghEqn.H @@ -0,0 +1,62 @@ +if (pimple.dict().lookupOrDefault("hydrostaticInitialization", false)) +{ + volScalarField& ph_rgh = regIOobject::store + ( + new volScalarField + ( + IOobject + ( + "ph_rgh", + "0", + mesh, + IOobject::MUST_READ, + IOobject::NO_WRITE + ), + mesh + ) + ); + + if (equal(runTime.value(), 0)) + { + p = ph_rgh + rho*gh + pRef; + thermo.correct(); + rho = thermo.rho(); + + label nCorr + ( + pimple.dict().lookupOrDefault