reactingFoam: Added frozenFlow option

with frozenFlow set true in the PIMPLE sub-dictionary of fvSolution the p-U
system is not solved and the energy-composition system including reactions is
solved with the fixed flow.
This commit is contained in:
Henry Weller
2020-12-16 11:24:00 +00:00
parent e27f719aed
commit a620e93fa7
3 changed files with 81 additions and 74 deletions

View File

@ -9,34 +9,32 @@ tmp<fv::convectionScheme<scalar>> mvConvection
) )
); );
reaction->correct();
forAll(Y, i)
{ {
reaction->correct(); if (composition.solve(i))
forAll(Y, i)
{ {
if (composition.solve(i)) volScalarField& Yi = Y[i];
{
volScalarField& Yi = Y[i];
fvScalarMatrix YiEqn fvScalarMatrix YiEqn
( (
fvm::ddt(rho, Yi) fvm::ddt(rho, Yi)
+ mvConvection->fvmDiv(phi, Yi) + mvConvection->fvmDiv(phi, Yi)
+ thermophysicalTransport->divj(Yi) + thermophysicalTransport->divj(Yi)
== ==
reaction->R(Yi) reaction->R(Yi)
+ fvOptions(rho, Yi) + fvOptions(rho, Yi)
); );
YiEqn.relax(); YiEqn.relax();
fvOptions.constrain(YiEqn); fvOptions.constrain(YiEqn);
YiEqn.solve("Yi"); YiEqn.solve("Yi");
fvOptions.correct(Yi); fvOptions.correct(Yi);
}
} }
composition.normalise();
} }
composition.normalise();

View File

@ -0,0 +1 @@
#include "../../compressible/rhoPimpleFoam/correctPhi.H"

View File

@ -101,63 +101,71 @@ int main(int argc, char *argv[])
Info<< "Time = " << runTime.timeName() << nl << endl; Info<< "Time = " << runTime.timeName() << nl << endl;
// --- Pressure-velocity PIMPLE corrector loop if (pimple.frozenFlow())
while (pimple.loop())
{ {
if (pimple.firstPimpleIter() || moveMeshOuterCorrectors)
{
// Store momentum to set rhoUf for introduced faces.
autoPtr<volVectorField> rhoU;
if (rhoUf.valid())
{
rhoU = new volVectorField("rhoU", rho*U);
}
// Do any mesh changes
mesh.update();
if (mesh.changing())
{
MRF.update();
if (correctPhi)
{
// Calculate absolute flux
// from the mapped surface velocity
phi = mesh.Sf() & rhoUf();
#include "../../compressible/rhoPimpleFoam/correctPhi.H"
// Make the fluxes relative to the mesh-motion
fvc::makeRelative(phi, rho, U);
}
if (checkMeshCourantNo)
{
#include "meshCourantNo.H"
}
}
}
if (pimple.firstPimpleIter() && !pimple.simpleRho())
{
#include "rhoEqn.H"
}
#include "UEqn.H"
#include "YEqn.H" #include "YEqn.H"
#include "EEqn.H" #include "EEqn.H"
}
// --- Pressure corrector loop else
while (pimple.correct()) {
// --- Pressure-velocity PIMPLE corrector loop
while (pimple.loop())
{ {
#include "../../compressible/rhoPimpleFoam/pEqn.H" if (pimple.firstPimpleIter() || moveMeshOuterCorrectors)
} {
// Store momentum to set rhoUf for introduced faces.
autoPtr<volVectorField> rhoU;
if (rhoUf.valid())
{
rhoU = new volVectorField("rhoU", rho*U);
}
if (pimple.turbCorr()) // Do any mesh changes
{ mesh.update();
turbulence->correct();
thermophysicalTransport->correct(); if (mesh.changing())
{
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.firstPimpleIter() && !pimple.simpleRho())
{
#include "rhoEqn.H"
}
#include "UEqn.H"
#include "YEqn.H"
#include "EEqn.H"
// --- Pressure corrector loop
while (pimple.correct())
{
#include "../../compressible/rhoPimpleFoam/pEqn.H"
}
if (pimple.turbCorr())
{
turbulence->correct();
thermophysicalTransport->correct();
}
} }
} }