diff --git a/applications/solvers/heatTransfer/buoyantPimpleFoam/Make/options b/applications/solvers/heatTransfer/buoyantPimpleFoam/Make/options index d3224df411..a9b16d2f25 100644 --- a/applications/solvers/heatTransfer/buoyantPimpleFoam/Make/options +++ b/applications/solvers/heatTransfer/buoyantPimpleFoam/Make/options @@ -1,4 +1,5 @@ EXE_INC = \ + -I../../compressible/rhoPimpleFoam \ -I$(LIB_SRC)/finiteVolume/lnInclude \ -I$(LIB_SRC)/meshTools/lnInclude \ -I$(LIB_SRC)/sampling/lnInclude \ @@ -7,6 +8,8 @@ EXE_INC = \ -I$(LIB_SRC)/thermophysicalModels/radiation/lnInclude \ -I$(LIB_SRC)/TurbulenceModels/turbulenceModels/lnInclude \ -I$(LIB_SRC)/TurbulenceModels/compressible/lnInclude \ + -I$(LIB_SRC)/dynamicMesh/lnInclude \ + -I$(LIB_SRC)/dynamicFvMesh/lnInclude \ -I$(LIB_SRC)/regionFaModels\lnInclude EXE_LIBS = \ @@ -21,4 +24,7 @@ EXE_LIBS = \ -lturbulenceModels \ -lcompressibleTurbulenceModels \ -latmosphericModels \ + -ldynamicMesh \ + -ldynamicFvMesh \ + -ltopoChangerFvMesh \ -lregionFaModels diff --git a/applications/solvers/heatTransfer/buoyantPimpleFoam/buoyantPimpleFoam.C b/applications/solvers/heatTransfer/buoyantPimpleFoam/buoyantPimpleFoam.C index 691916c17e..70595ebba3 100644 --- a/applications/solvers/heatTransfer/buoyantPimpleFoam/buoyantPimpleFoam.C +++ b/applications/solvers/heatTransfer/buoyantPimpleFoam/buoyantPimpleFoam.C @@ -6,6 +6,7 @@ \\/ M anipulation | ------------------------------------------------------------------------------- Copyright (C) 2011-2017 OpenFOAM Foundation + Copyright (C) 2021 OpenCFD Ltd. ------------------------------------------------------------------------------- License This file is part of OpenFOAM. @@ -30,8 +31,9 @@ Group grpHeatTransferSolvers Description - Transient solver for buoyant, turbulent flow of compressible fluids for - ventilation and heat-transfer. + Transient solver for buoyant, turbulent flow of compressible fluids + for ventilation and heat-transfer, with optional mesh motion + and mesh topology changes. Turbulence is modelled using a run-time selectable compressible RAS or LES model. @@ -39,12 +41,16 @@ Description \*---------------------------------------------------------------------------*/ #include "fvCFD.H" +#include "dynamicFvMesh.H" #include "rhoThermo.H" #include "turbulentFluidThermoModel.H" #include "radiationModel.H" +#include "CorrectPhi.H" #include "fvOptions.H" #include "pimpleControl.H" #include "pressureControl.H" +#include "localEulerDdtScheme.H" +#include "fvcSmooth.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // @@ -53,7 +59,8 @@ int main(int argc, char *argv[]) argList::addNote ( "Transient solver for buoyant, turbulent fluid flow" - " of compressible fluids, including radiation." + " of compressible fluids, including radiation," + " with optional mesh motion and mesh topology changes." ); #include "postProcess.H" @@ -61,36 +68,105 @@ int main(int argc, char *argv[]) #include "addCheckCaseOptions.H" #include "setRootCaseLists.H" #include "createTime.H" - #include "createMesh.H" - #include "createControl.H" + #include "createDynamicFvMesh.H" + #include "createDyMControls.H" + #include "initContinuityErrs.H" #include "createFields.H" #include "createFieldRefs.H" - #include "initContinuityErrs.H" - #include "createTimeControls.H" - #include "compressibleCourantNo.H" - #include "setInitialDeltaT.H" + #include "createRhoUfIfPresent.H" turbulence->validate(); + if (!LTS) + { + #include "compressibleCourantNo.H" + #include "setInitialDeltaT.H" + } + // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // Info<< "\nStarting time loop\n" << endl; while (runTime.run()) { - #include "readTimeControls.H" - #include "compressibleCourantNo.H" - #include "setDeltaT.H" + #include "readDyMControls.H" + + // Store divrhoU from the previous mesh + // so that it can be mapped and used in correctPhi + // to ensure the corrected phi has the same divergence + autoPtr divrhoU; + if (correctPhi) + { + divrhoU.reset + ( + new volScalarField + ( + "divrhoU", + fvc::div(fvc::absolute(phi, rho, U)) + ) + ); + } + + if (LTS) + { + #include "setRDeltaT.H" + } + else + { + #include "compressibleCourantNo.H" + #include "setDeltaT.H" + } ++runTime; Info<< "Time = " << runTime.timeName() << nl << endl; - #include "rhoEqn.H" - // --- Pressure-velocity PIMPLE corrector loop while (pimple.loop()) { + if (pimple.firstIter() || moveMeshOuterCorrectors) + { + // Store momentum to set rhoUf for introduced faces. + autoPtr rhoU; + if (rhoUf.valid()) + { + rhoU.reset(new volVectorField("rhoU", rho*U)); + } + + // Do any mesh changes + mesh.update(); + + if (mesh.changing()) + { + gh = (g & mesh.C()) - ghRef; + ghf = (g & mesh.Cf()) - ghRef; + + MRF.update(); + + if (correctPhi) + { + // Calculate absolute flux + // from the mapped surface velocity + phi = mesh.Sf() & rhoUf(); + + #include "correctPhi.H" + + // Make the fluxes relative to the mesh-motion + fvc::makeRelative(phi, rho, U); + } + + if (checkMeshCourantNo) + { + #include "meshCourantNo.H" + } + } + } + + if (pimple.firstIter() && !pimple.SIMPLErho()) + { + #include "rhoEqn.H" + } + #include "UEqn.H" #include "EEqn.H" diff --git a/applications/solvers/heatTransfer/buoyantPimpleFoam/createFields.H b/applications/solvers/heatTransfer/buoyantPimpleFoam/createFields.H index 22dc129a87..813d1c3240 100644 --- a/applications/solvers/heatTransfer/buoyantPimpleFoam/createFields.H +++ b/applications/solvers/heatTransfer/buoyantPimpleFoam/createFields.H @@ -1,3 +1,5 @@ +#include "createRDeltaT.H" + Info<< "Reading thermophysical properties\n" << endl; autoPtr pThermo(rhoThermo::New(mesh)); diff --git a/applications/solvers/heatTransfer/buoyantPimpleFoam/pEqn.H b/applications/solvers/heatTransfer/buoyantPimpleFoam/pEqn.H index 886f9c6a5f..733e97028d 100644 --- a/applications/solvers/heatTransfer/buoyantPimpleFoam/pEqn.H +++ b/applications/solvers/heatTransfer/buoyantPimpleFoam/pEqn.H @@ -28,6 +28,8 @@ MRF.makeRelative(fvc::interpolate(rho), phiHbyA); // Update the pressure BCs to ensure flux consistency constrainPressure(p_rgh, rho, U, phiHbyA, rhorAUf, MRF); +fvc::makeRelative(phiHbyA, rho, U); + fvScalarMatrix p_rghDDtEqn ( fvc::ddt(rho) + psi*correction(fvm::ddt(p_rgh)) @@ -104,7 +106,15 @@ else rho = thermo.rho(); +// Correct rhoUf if the mesh is moving +fvc::correctRhoUf(rhoUf, rho, U, phi); + if (thermo.dpdt()) { dpdt = fvc::ddt(p); + + if (mesh.moving()) + { + dpdt -= fvc::div(fvc::meshPhi(rho, U), p); + } }