XiFoam: Added relaxation statements to Su and Xi equations

This makes the solver more robust when the "Courant number" of the flame
propagation fluctuates.
This commit is contained in:
henry
2010-05-06 14:45:53 +01:00
parent 08916e83a2
commit 911f50cebe
2 changed files with 13 additions and 4 deletions

View File

@ -30,7 +30,7 @@ if (ign.ignited())
n /= mgb; n /= mgb;
# include "StCorr.H" #include "StCorr.H"
// Calculate turbulent flame speed flux // Calculate turbulent flame speed flux
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
@ -58,7 +58,7 @@ if (ign.ignited())
// Add ignition cell contribution to b-equation // Add ignition cell contribution to b-equation
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
# include "ignite.H" #include "ignite.H"
// Solve for b // Solve for b
@ -134,7 +134,7 @@ if (ign.ignited())
(sigmas*SuInf*(Su0 - SuInf) + sqr(SuMin)*sigmaExt) (sigmas*SuInf*(Su0 - SuInf) + sqr(SuMin)*sigmaExt)
/(sqr(Su0 - SuInf) + sqr(SuMin)); /(sqr(Su0 - SuInf) + sqr(SuMin));
solve fvScalarMatrix SuEqn
( (
fvm::ddt(rho, Su) fvm::ddt(rho, Su)
+ fvm::div(phi + phiXi, Su, "div(phiXi,Su)") + fvm::div(phi + phiXi, Su, "div(phiXi,Su)")
@ -144,6 +144,9 @@ if (ign.ignited())
- fvm::SuSp(rho*(sigmas + Rc), Su) - fvm::SuSp(rho*(sigmas + Rc), Su)
); );
SuEqn.relax();
SuEqn.solve();
// Limit the maximum Su // Limit the maximum Su
// ~~~~~~~~~~~~~~~~~~~~ // ~~~~~~~~~~~~~~~~~~~~
Su.min(SuMax); Su.min(SuMax);
@ -196,7 +199,7 @@ if (ign.ignited())
// Solve for the flame wrinkling // Solve for the flame wrinkling
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
solve fvScalarMatrix XiEqn
( (
fvm::ddt(rho, Xi) fvm::ddt(rho, Xi)
+ fvm::div(phi + phiXi, Xi, "div(phiXi,Xi)") + fvm::div(phi + phiXi, Xi, "div(phiXi,Xi)")
@ -215,6 +218,8 @@ if (ign.ignited())
) )
); );
XiEqn.relax();
XiEqn.solve();
// Correct boundedness of Xi // Correct boundedness of Xi
// ~~~~~~~~~~~~~~~~~~~~~~~~~ // ~~~~~~~~~~~~~~~~~~~~~~~~~

View File

@ -48,5 +48,9 @@ PISO
nNonOrthogonalCorrectors 0; nNonOrthogonalCorrectors 0;
} }
relaxationFactors
{
"(Xi|Su)" 1;
}
// ************************************************************************* // // ************************************************************************* //