multiphaseEulerFoam::populationBalanceModel: Removed temporary dilatation correction

and updated tutorials to work with the current phase limit stabilisation.
This commit is contained in:
Henry Weller
2022-04-12 10:23:42 +01:00
parent 043b06f129
commit 5e99344348
4 changed files with 2 additions and 4 deletions

View File

@ -206,7 +206,6 @@ void Foam::diameterModels::shapeModels::fractal::correct()
fvScalarMatrix kappaEqn
(
fvm::ddt(alpha, fi, kappa_) + fvm::div(fAlphaPhi, kappa_)
- fvm::Sp(fvc::ddt(alpha, fi) + fvc::div(fAlphaPhi), kappa_)
==
- sinteringModel_->R()
+ Su_

View File

@ -1072,7 +1072,6 @@ void Foam::diameterModels::populationBalanceModel::solve()
fvScalarMatrix sizeGroupEqn
(
fvm::ddt(alpha, fi) + fvm::div(phase.alphaPhi(), fi)
- fvm::Sp(fvc::ddt(alpha) + fvc::div(phase.alphaPhi()), fi)
==
Su_[i]
- fvm::Sp(Sp_[i], fi)

View File

@ -28,7 +28,7 @@ solvers
tolerance 1e-4;
scale true;
solveOnFinalIterOnly true;
sourceUpdateInterval 20;
sourceUpdateInterval 1;
}
p_rgh

View File

@ -28,7 +28,7 @@ solvers
tolerance 1e-4;
scale true;
solveOnFinalIterOnly true;
sourceUpdateInterval 20;
sourceUpdateInterval 1;
}
p_rgh