From 29e83f3958209cab89c8bb67d98ffb069eb64e3e Mon Sep 17 00:00:00 2001 From: Henry Weller <http://cfd.direct> Date: Tue, 7 Feb 2017 18:59:40 +0000 Subject: [PATCH] compressibleInterFoam: Added support for fvOptions in both the U and T equations --- .../solvers/multiphase/compressibleInterFoam/TEqn.H | 9 +++++++-- .../solvers/multiphase/compressibleInterFoam/UEqn.H | 6 ++++++ .../compressibleInterDyMFoam/compressibleInterDyMFoam.C | 4 +++- .../compressibleInterDyMFoam/pEqn.H | 4 +--- .../compressibleInterFoam/compressibleInterFoam.C | 4 +++- .../solvers/multiphase/compressibleInterFoam/pEqn.H | 4 +--- 6 files changed, 21 insertions(+), 10 deletions(-) diff --git a/applications/solvers/multiphase/compressibleInterFoam/TEqn.H b/applications/solvers/multiphase/compressibleInterFoam/TEqn.H index 13e1feb546e..7bd4ca35688 100644 --- a/applications/solvers/multiphase/compressibleInterFoam/TEqn.H +++ b/applications/solvers/multiphase/compressibleInterFoam/TEqn.H @@ -12,12 +12,17 @@ alpha1/mixture.thermo1().Cv() + alpha2/mixture.thermo2().Cv() ) + == + fvOptions(rho, T) ); TEqn.relax(); + + fvOptions.constrain(TEqn); + TEqn.solve(); - mixture.correct(); + fvOptions.correct(T); - Info<< "min(T) " << min(T).value() << endl; + mixture.correct(); } diff --git a/applications/solvers/multiphase/compressibleInterFoam/UEqn.H b/applications/solvers/multiphase/compressibleInterFoam/UEqn.H index 44be32aac22..18e947c8a7a 100644 --- a/applications/solvers/multiphase/compressibleInterFoam/UEqn.H +++ b/applications/solvers/multiphase/compressibleInterFoam/UEqn.H @@ -3,10 +3,14 @@ fvm::ddt(rho, U) + fvm::div(rhoPhi, U) + turbulence->divDevRhoReff(U) + == + fvOptions(rho, U) ); UEqn.relax(); + fvOptions.constrain(UEqn); + if (pimple.momentumPredictor()) { solve @@ -23,5 +27,7 @@ ) ); + fvOptions.correct(U); + K = 0.5*magSqr(U); } diff --git a/applications/solvers/multiphase/compressibleInterFoam/compressibleInterDyMFoam/compressibleInterDyMFoam.C b/applications/solvers/multiphase/compressibleInterFoam/compressibleInterDyMFoam/compressibleInterDyMFoam.C index 6f7ee7ffcea..157830ae6bb 100644 --- a/applications/solvers/multiphase/compressibleInterFoam/compressibleInterDyMFoam/compressibleInterDyMFoam.C +++ b/applications/solvers/multiphase/compressibleInterFoam/compressibleInterDyMFoam/compressibleInterDyMFoam.C @@ -2,7 +2,7 @@ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | - \\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation + \\ / A nd | Copyright (C) 2011-2017 OpenFOAM Foundation \\/ M anipulation | ------------------------------------------------------------------------------- License @@ -46,6 +46,7 @@ Description #include "twoPhaseMixtureThermo.H" #include "turbulentFluidThermoModel.H" #include "pimpleControl.H" +#include "fvOptions.H" #include "CorrectPhi.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // @@ -60,6 +61,7 @@ int main(int argc, char *argv[]) #include "initContinuityErrs.H" #include "createControl.H" #include "createFields.H" + #include "createFvOptions.H" #include "createUf.H" #include "createControls.H" #include "CourantNo.H" diff --git a/applications/solvers/multiphase/compressibleInterFoam/compressibleInterDyMFoam/pEqn.H b/applications/solvers/multiphase/compressibleInterFoam/compressibleInterDyMFoam/pEqn.H index fef73aa3659..661cb3a78eb 100644 --- a/applications/solvers/multiphase/compressibleInterFoam/compressibleInterDyMFoam/pEqn.H +++ b/applications/solvers/multiphase/compressibleInterFoam/compressibleInterDyMFoam/pEqn.H @@ -101,6 +101,7 @@ U = HbyA + rAU*fvc::reconstruct((phig + p_rghEqnIncomp.flux())/rAUf); U.correctBoundaryConditions(); + fvOptions.correct(U); } } @@ -121,7 +122,4 @@ p_rgh.correctBoundaryConditions(); K = 0.5*magSqr(U); - - Info<< "max(U) " << max(mag(U)).value() << endl; - Info<< "min(p_rgh) " << min(p_rgh).value() << endl; } diff --git a/applications/solvers/multiphase/compressibleInterFoam/compressibleInterFoam.C b/applications/solvers/multiphase/compressibleInterFoam/compressibleInterFoam.C index 150cc8bdfdb..8775d0fbd33 100644 --- a/applications/solvers/multiphase/compressibleInterFoam/compressibleInterFoam.C +++ b/applications/solvers/multiphase/compressibleInterFoam/compressibleInterFoam.C @@ -2,7 +2,7 @@ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | - \\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation + \\ / A nd | Copyright (C) 2011-2017 OpenFOAM Foundation \\/ M anipulation | ------------------------------------------------------------------------------- License @@ -44,6 +44,7 @@ Description #include "twoPhaseMixtureThermo.H" #include "turbulentFluidThermoModel.H" #include "pimpleControl.H" +#include "fvOptions.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // @@ -57,6 +58,7 @@ int main(int argc, char *argv[]) #include "createControl.H" #include "createTimeControls.H" #include "createFields.H" + #include "createFvOptions.H" #include "CourantNo.H" #include "setInitialDeltaT.H" diff --git a/applications/solvers/multiphase/compressibleInterFoam/pEqn.H b/applications/solvers/multiphase/compressibleInterFoam/pEqn.H index 78e46f46642..1328b0abd9d 100644 --- a/applications/solvers/multiphase/compressibleInterFoam/pEqn.H +++ b/applications/solvers/multiphase/compressibleInterFoam/pEqn.H @@ -98,6 +98,7 @@ U = HbyA + rAU*fvc::reconstruct((phig + p_rghEqnIncomp.flux())/rAUf); U.correctBoundaryConditions(); + fvOptions.correct(U); } } @@ -112,7 +113,4 @@ p_rgh.correctBoundaryConditions(); K = 0.5*magSqr(U); - - Info<< "max(U) " << max(mag(U)).value() << endl; - Info<< "min(p_rgh) " << min(p_rgh).value() << endl; } -- GitLab