From e00a56cdbb78f3a5bef4e67711f6a23f76e49206 Mon Sep 17 00:00:00 2001 From: Henry Weller Date: Tue, 6 Dec 2016 09:30:50 +0000 Subject: [PATCH] reactingTwoPhaseEulerFoam: Corrected LTS support Resolves bug-report http://bugs.openfoam.org/view.php?id=2374 --- .../reactingTwoPhaseEulerFoam/pUf/pEqn.H | 10 ++++------ .../reactingTwoPhaseEulerFoam.C | 16 ++++++++++++++++ 2 files changed, 20 insertions(+), 6 deletions(-) diff --git a/applications/solvers/multiphase/reactingEulerFoam/reactingTwoPhaseEulerFoam/pUf/pEqn.H b/applications/solvers/multiphase/reactingEulerFoam/reactingTwoPhaseEulerFoam/pUf/pEqn.H index b76932216..7950d76b9 100644 --- a/applications/solvers/multiphase/reactingEulerFoam/reactingTwoPhaseEulerFoam/pUf/pEqn.H +++ b/applications/solvers/multiphase/reactingEulerFoam/reactingTwoPhaseEulerFoam/pUf/pEqn.H @@ -27,14 +27,12 @@ surfaceScalarField Kdf("Kdf", fluid.Kdf()); // Virtual-mass coefficient surfaceScalarField Vmf("Vmf", fluid.Vmf()); -const surfaceScalarField& rDeltaTf = fv::localEulerDdt::localRDeltaTf(mesh); - surfaceScalarField rAUf1 ( IOobject::groupName("rAUf", phase1.name()), 1.0 /( - rDeltaTf*(alphaRhof10 + Vmf) + byDt(alphaRhof10 + Vmf) + fvc::interpolate(U1Eqn.A()) + Kdf ) @@ -45,7 +43,7 @@ surfaceScalarField rAUf2 IOobject::groupName("rAUf", phase2.name()), 1.0 /( - rDeltaTf*(alphaRhof20 + Vmf) + byDt(alphaRhof20 + Vmf) + fvc::interpolate(U2Eqn.A()) + Kdf ) @@ -159,7 +157,7 @@ while (pimple.correct()) rAUf1 *( (alphaRhof10 + Vmf) - *rDeltaTf*MRF.absolute(phi1.oldTime()) + *byDt(MRF.absolute(phi1.oldTime())) + fvc::flux(U1Eqn.H()) + Vmf*ddtPhi2 + Kdf*MRF.absolute(phi2) @@ -177,7 +175,7 @@ while (pimple.correct()) rAUf2 *( (alphaRhof20 + Vmf) - *rDeltaTf*MRF.absolute(phi2.oldTime()) + *byDt(MRF.absolute(phi2.oldTime())) + fvc::flux(U2Eqn.H()) + Vmf*ddtPhi1 + Kdf*MRF.absolute(phi1) diff --git a/applications/solvers/multiphase/reactingEulerFoam/reactingTwoPhaseEulerFoam/reactingTwoPhaseEulerFoam.C b/applications/solvers/multiphase/reactingEulerFoam/reactingTwoPhaseEulerFoam/reactingTwoPhaseEulerFoam.C index 6e9d4e4be..bebe56200 100644 --- a/applications/solvers/multiphase/reactingEulerFoam/reactingTwoPhaseEulerFoam/reactingTwoPhaseEulerFoam.C +++ b/applications/solvers/multiphase/reactingEulerFoam/reactingTwoPhaseEulerFoam/reactingTwoPhaseEulerFoam.C @@ -40,6 +40,22 @@ Description #include "localEulerDdtScheme.H" #include "fvcSmooth.H" +namespace Foam +{ + tmp byDt(const surfaceScalarField& sf) + { + if (fv::localEulerDdt::enabled(sf.mesh())) + { + return fv::localEulerDdt::localRDeltaTf(sf.mesh())*sf; + } + else + { + return sf/sf.mesh().time().deltaT(); + } + } +} + + // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // int main(int argc, char *argv[])