From 15e553b88109e6646cf66775e13f4c2ab786654d Mon Sep 17 00:00:00 2001 From: Henry <Henry> Date: Mon, 6 May 2013 22:25:57 +0100 Subject: [PATCH] compressibleTwoPhaseEulerFoam: Corrected pressure work due to interface motion in EEqns --- .../compressibleTwoPhaseEulerFoam/EEqns.H | 8 ++++---- .../createFields.H | 18 +++--------------- .../compressibleTwoPhaseEulerFoam/pEqn.H | 9 ++------- 3 files changed, 9 insertions(+), 26 deletions(-) diff --git a/applications/solvers/multiphase/compressibleTwoPhaseEulerFoam/EEqns.H b/applications/solvers/multiphase/compressibleTwoPhaseEulerFoam/EEqns.H index 1f12fc9a121..1810bdd4d9f 100644 --- a/applications/solvers/multiphase/compressibleTwoPhaseEulerFoam/EEqns.H +++ b/applications/solvers/multiphase/compressibleTwoPhaseEulerFoam/EEqns.H @@ -28,8 +28,8 @@ + ( he1.name() == thermo1.phasePropertyName("e") - ? fvc::div(alphaPhi1, p) - : -dalpha1pdt + ? fvc::ddt(alpha1)*p + fvc::div(alphaPhi1, p) + : -alpha1*dpdt )/rho1 - fvm::laplacian(k1, he1) @@ -50,8 +50,8 @@ + ( he2.name() == thermo2.phasePropertyName("e") - ? fvc::div(alphaPhi2, p) - : -dalpha2pdt + ? fvc::ddt(alpha2)*p + fvc::div(alphaPhi2, p) + : -alpha2*dpdt )/rho2 - fvm::laplacian(k2, he2) diff --git a/applications/solvers/multiphase/compressibleTwoPhaseEulerFoam/createFields.H b/applications/solvers/multiphase/compressibleTwoPhaseEulerFoam/createFields.H index cc5b3ebe219..5b030b95ea6 100644 --- a/applications/solvers/multiphase/compressibleTwoPhaseEulerFoam/createFields.H +++ b/applications/solvers/multiphase/compressibleTwoPhaseEulerFoam/createFields.H @@ -275,8 +275,8 @@ ); - Info<< "Creating field dalpha1pdt\n" << endl; - volScalarField dalpha1pdt + Info<< "Creating field dpdt\n" << endl; + volScalarField dpdt ( IOobject ( @@ -285,21 +285,9 @@ mesh ), mesh, - dimensionedScalar("dalpha1pdt", p.dimensions()/dimTime, 0) + dimensionedScalar("dpdt", p.dimensions()/dimTime, 0) ); - Info<< "Creating field dalpha2pdt\n" << endl; - volScalarField dalpha2pdt - ( - IOobject - ( - "dpdt", - runTime.timeName(), - mesh - ), - mesh, - dimensionedScalar("dalpha2pdt", p.dimensions()/dimTime, 0) - ); Info<< "Creating field kinetic energy K\n" << endl; volScalarField K1("K" + phase1Name, 0.5*magSqr(U1)); diff --git a/applications/solvers/multiphase/compressibleTwoPhaseEulerFoam/pEqn.H b/applications/solvers/multiphase/compressibleTwoPhaseEulerFoam/pEqn.H index 02e49b56618..756a231fc2f 100644 --- a/applications/solvers/multiphase/compressibleTwoPhaseEulerFoam/pEqn.H +++ b/applications/solvers/multiphase/compressibleTwoPhaseEulerFoam/pEqn.H @@ -208,13 +208,8 @@ K1 = 0.5*magSqr(U1); K2 = 0.5*magSqr(U2); - if (thermo1.dpdt()) + if (thermo1.dpdt() || thermo2.dpdt()) { - dalpha1pdt = fvc::ddt(alpha1, p); - } - - if (thermo2.dpdt()) - { - dalpha2pdt = fvc::ddt(alpha2, p); + dpdt = fvc::ddt(p); } } -- GitLab