mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
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:
@ -30,7 +30,7 @@ if (ign.ignited())
|
||||
n /= mgb;
|
||||
|
||||
|
||||
# include "StCorr.H"
|
||||
#include "StCorr.H"
|
||||
|
||||
// Calculate turbulent flame speed flux
|
||||
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
|
||||
@ -58,7 +58,7 @@ if (ign.ignited())
|
||||
|
||||
// Add ignition cell contribution to b-equation
|
||||
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
|
||||
# include "ignite.H"
|
||||
#include "ignite.H"
|
||||
|
||||
|
||||
// Solve for b
|
||||
@ -134,7 +134,7 @@ if (ign.ignited())
|
||||
(sigmas*SuInf*(Su0 - SuInf) + sqr(SuMin)*sigmaExt)
|
||||
/(sqr(Su0 - SuInf) + sqr(SuMin));
|
||||
|
||||
solve
|
||||
fvScalarMatrix SuEqn
|
||||
(
|
||||
fvm::ddt(rho, Su)
|
||||
+ fvm::div(phi + phiXi, Su, "div(phiXi,Su)")
|
||||
@ -144,6 +144,9 @@ if (ign.ignited())
|
||||
- fvm::SuSp(rho*(sigmas + Rc), Su)
|
||||
);
|
||||
|
||||
SuEqn.relax();
|
||||
SuEqn.solve();
|
||||
|
||||
// Limit the maximum Su
|
||||
// ~~~~~~~~~~~~~~~~~~~~
|
||||
Su.min(SuMax);
|
||||
@ -196,7 +199,7 @@ if (ign.ignited())
|
||||
|
||||
// Solve for the flame wrinkling
|
||||
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
|
||||
solve
|
||||
fvScalarMatrix XiEqn
|
||||
(
|
||||
fvm::ddt(rho, Xi)
|
||||
+ fvm::div(phi + phiXi, Xi, "div(phiXi,Xi)")
|
||||
@ -215,6 +218,8 @@ if (ign.ignited())
|
||||
)
|
||||
);
|
||||
|
||||
XiEqn.relax();
|
||||
XiEqn.solve();
|
||||
|
||||
// Correct boundedness of Xi
|
||||
// ~~~~~~~~~~~~~~~~~~~~~~~~~
|
||||
|
||||
Reference in New Issue
Block a user