Commit f966d205 authored by Henry Weller's avatar Henry Weller Committed by Andrew Heather
Browse files

ENH: compressibleInterPhaseTransportFoam: New variant of compressibleInterFoam...

ENH: compressibleInterPhaseTransportFoam: New variant of compressibleInterFoam supporting separate phase stress models

In this version of compressibleInterFoam separate stress models (laminar,
non-Newtonian, LES or RAS) are instantiated for each of the two phases allowing
for completely different modeling for the phases.

e.g. in the climbingRod tutorial case provided a Newtonian laminar model is
instantiated for the air and a Maxwell non-Newtonian model is instantiated for
the viscoelastic liquid.  To stabilize the Maxwell model in regions where the
liquid phase-fraction is 0 the new symmTensorPhaseLimitStabilization fvOption is
applied.

Other phase stress modeling combinations are also possible, e.g. the air may be
turbulent but the liquid laminar and an RAS or LES model applied to the air
only.  However, to stabilize this combination a suitable fvOption would need to
be applied to the turbulence properties where the air phase-fraction is 0.

Henry G. Weller, Chris Greenshields
CFD Direct Ltd.
parent 4f922181
......@@ -10,5 +10,6 @@ wmake $targetType surfaceTensionModels
wmake $targetType
wmake $targetType compressibleInterDyMFoam
wmake $targetType compressibleInterFilmFoam
compressibleInterPhaseTransportFoam/Allwmake $targetType $*
#------------------------------------------------------------------------------
#!/bin/sh
cd ${0%/*} || exit 1 # Run from this directory
wclean libso VoFphaseCompressibleTurbulenceModels
wclean
#------------------------------------------------------------------------------
#!/bin/sh
cd ${0%/*} || exit 1 # Run from this directory
# Parse arguments for library compilation
. $WM_PROJECT_DIR/wmake/scripts/AllwmakeParseArguments
wmake $targetType VoFphaseCompressibleTurbulenceModels
wmake $targetType
#------------------------------------------------------------------------------
compressibleInterPhaseTransportFoam.C
EXE = $(FOAM_APPBIN)/compressibleInterPhaseTransportFoam
EXE_INC = \
-I. \
-I$(FOAM_SOLVERS)/multiphase/VoF \
-I$(FOAM_SOLVERS)/multiphase/compressibleInterFoam/twoPhaseMixtureThermo \
-I$(FOAM_SOLVERS)/multiphase/compressibleInterFoam \
-I$(LIB_SRC)/transportModels/compressible/lnInclude \
-I$(LIB_SRC)/thermophysicalModels/basic/lnInclude \
-I$(LIB_SRC)/transportModels/twoPhaseMixture/lnInclude \
-I$(LIB_SRC)/transportModels/interfaceProperties/lnInclude \
-I$(LIB_SRC)/TurbulenceModels/turbulenceModels/lnInclude \
-I$(LIB_SRC)/TurbulenceModels/compressible/lnInclude \
-I$(LIB_SRC)/TurbulenceModels/phaseCompressible/lnInclude \
-IVoFphaseCompressibleTurbulenceModels/lnInclude \
-I$(LIB_SRC)/finiteVolume/lnInclude \
-I$(LIB_SRC)/meshTools/lnInclude
EXE_LIBS = \
-ltwoPhaseMixtureThermo \
-ltwoPhaseSurfaceTension \
-lcompressibleTransportModels \
-lfluidThermophysicalModels \
-lspecie \
-ltwoPhaseMixture \
-ltwoPhaseProperties \
-linterfaceProperties \
-lturbulenceModels \
-lcompressibleTurbulenceModels \
-lVoFphaseCompressibleTurbulenceModels \
-lfiniteVolume \
-lfvOptions \
-lmeshTools
{
fvScalarMatrix TEqn
(
fvm::ddt(rho, T) + fvm::div(rhoPhi, T)
- fvm::Sp(contErr, T)
- fvm::laplacian
(
mixture.alphaEff
(
alpha1*turbulence1->mut()
+ alpha2*turbulence2->mut()
),
T
)
+ (
fvc::div(fvc::absolute(phi, U), p)
+ fvc::ddt(rho, K) + fvc::div(rhoPhi, K)
)
*(
alpha1/mixture.thermo1().Cv()
+ alpha2/mixture.thermo2().Cv()
)
==
fvOptions(rho, T)
);
TEqn.relax();
fvOptions.constrain(TEqn);
TEqn.solve();
fvOptions.correct(T);
mixture.correctThermo();
mixture.correct();
}
fvVectorMatrix UEqn
(
fvm::ddt(rho, U) + fvm::div(rhoPhi, U)
- fvm::Sp(contErr, U)
+ MRF.DDt(rho, U)
+ turbulence1->divDevRhoReff(U)
+ turbulence2->divDevRhoReff(U)
==
fvOptions(rho, U)
);
UEqn.relax();
fvOptions.constrain(UEqn);
if (pimple.momentumPredictor())
{
solve
(
UEqn
==
fvc::reconstruct
(
(
mixture.surfaceTensionForce()
- ghf*fvc::snGrad(rho)
- fvc::snGrad(p_rgh)
) * mesh.magSf()
)
);
fvOptions.correct(U);
K = 0.5*magSqr(U);
}
VoFphaseCompressibleTurbulenceModels.C
LIB = $(FOAM_LIBBIN)/libVoFphaseCompressibleTurbulenceModels
EXE_INC = \
-I$(LIB_SRC)/transportModels/compressible/lnInclude \
-I$(LIB_SRC)/thermophysicalModels/basic/lnInclude \
-I$(LIB_SRC)/transportModels/incompressible/transportModel \
-I$(LIB_SRC)/TurbulenceModels/compressible/lnInclude \
-I$(LIB_SRC)/TurbulenceModels/turbulenceModels/lnInclude \
-I$(LIB_SRC)/TurbulenceModels/phaseCompressible/lnInclude \
-I$(LIB_SRC)/finiteVolume/lnInclude \
-I$(LIB_SRC)/meshTools/lnInclude
LIB_LIBS = \
-lcompressibleTransportModels \
-lfluidThermophysicalModels \
-lspecie \
-lturbulenceModels \
-lcompressibleTurbulenceModels \
-lincompressibleTransportModels \
-lfiniteVolume \
-lfvOptions \
-lmeshTools
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2017 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Typedef
Foam::VoFphaseCompressibleTurbulenceModel
Description
Typedef for VoFphaseCompressibleTurbulenceModel
\*---------------------------------------------------------------------------*/
#ifndef VoFphaseCompressibleTurbulenceModel_H
#define VoFphaseCompressibleTurbulenceModel_H
#include "VoFphaseCompressibleTurbulenceModelFwd.H"
#include "PhaseCompressibleTurbulenceModel.H"
#include "ThermalDiffusivity.H"
#include "fluidThermo.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
typedef ThermalDiffusivity<PhaseCompressibleTurbulenceModel<fluidThermo>>
VoFphaseCompressibleTurbulenceModel;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2017 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Typedef
Foam::VoFphaseCompressibleTurbulenceModel
Description
Forward declaration of typedef for VoFphaseCompressibleTurbulenceModel
\*---------------------------------------------------------------------------*/
#ifndef VoFphaseCompressibleTurbulenceModelFwd_H
#define VoFphaseCompressibleTurbulenceModelFwd_H
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
class fluidThermo;
template<class TransportModel>
class PhaseCompressibleTurbulenceModel;
template<class BasicTurbulenceModel>
class ThermalDiffusivity;
typedef ThermalDiffusivity<PhaseCompressibleTurbulenceModel<fluidThermo>>
VoFphaseCompressibleTurbulenceModel;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2017 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "PhaseCompressibleTurbulenceModel.H"
#include "fluidThermo.H"
#include "addToRunTimeSelectionTable.H"
#include "makeTurbulenceModel.H"
#include "ThermalDiffusivity.H"
#include "EddyDiffusivity.H"
#include "laminarModel.H"
#include "RASModel.H"
#include "LESModel.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
makeTurbulenceModelTypes
(
volScalarField,
volScalarField,
compressibleTurbulenceModel,
PhaseCompressibleTurbulenceModel,
ThermalDiffusivity,
fluidThermo
);
makeBaseTurbulenceModel
(
volScalarField,
volScalarField,
compressibleTurbulenceModel,
PhaseCompressibleTurbulenceModel,
ThermalDiffusivity,
fluidThermo
);
#define makeLaminarModel(Type) \
makeTemplatedLaminarModel \
(fluidThermoPhaseCompressibleTurbulenceModel, laminar, Type)
#define makeRASModel(Type) \
makeTemplatedTurbulenceModel \
(fluidThermoPhaseCompressibleTurbulenceModel, RAS, Type)
#define makeLESModel(Type) \
makeTemplatedTurbulenceModel \
(fluidThermoPhaseCompressibleTurbulenceModel, LES, Type)
#include "Stokes.H"
makeLaminarModel(Stokes);
#include "Maxwell.H"
makeLaminarModel(Maxwell);
#include "kEpsilon.H"
makeRASModel(kEpsilon);
#include "kOmegaSST.H"
makeRASModel(kOmegaSST);
#include "Smagorinsky.H"
makeLESModel(Smagorinsky);
#include "kEqn.H"
makeLESModel(kEqn);
// ************************************************************************* //
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2017 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Application
compressibleInterPhaseTransportFoam
Description
Solver for 2 compressible, non-isothermal immiscible fluids using a VOF
(volume of fluid) phase-fraction based interface capturing approach.
The momentum and other fluid properties are of the "mixture" and a single
momentum equation is solved.
The fluid stress modelling is generic Euler-Euler two-phase in which
separate laminar, RAS or LES are selected for each of the phases.
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
#include "CMULES.H"
#include "EulerDdtScheme.H"
#include "localEulerDdtScheme.H"
#include "CrankNicolsonDdtScheme.H"
#include "subCycle.H"
#include "rhoThermo.H"
#include "twoPhaseMixture.H"
#include "twoPhaseMixtureThermo.H"
#include "turbulentFluidThermoModel.H"
#include "pimpleControl.H"
#include "fvOptions.H"
#include "fvcSmooth.H"
#include "VoFphaseCompressibleTurbulenceModel.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
int main(int argc, char *argv[])
{
#include "postProcess.H"
#include "setRootCase.H"
#include "createTime.H"
#include "createMesh.H"
#include "createControl.H"
#include "createTimeControls.H"
#include "createFields.H"
#include "createAlphaFluxes.H"
volScalarField& p = mixture.p();
volScalarField& T = mixture.T();
const volScalarField& psi1 = mixture.thermo1().psi();
const volScalarField& psi2 = mixture.thermo2().psi();
surfaceScalarField alphaRhoPhi1
(
IOobject::groupName("alphaRhoPhi", alpha1.group()),
fvc::interpolate(rho1)*alphaPhi10
);
autoPtr<PhaseCompressibleTurbulenceModel<fluidThermo>> turbulence1
(
PhaseCompressibleTurbulenceModel<fluidThermo>::New
(
alpha1,
rho1,
U,
alphaRhoPhi1,
phi,
mixture.thermo1()
)
);
surfaceScalarField alphaRhoPhi2
(
IOobject::groupName("alphaRhoPhi", alpha2.group()),
fvc::interpolate(rho2)*(phi - alphaPhi10)
);
autoPtr<PhaseCompressibleTurbulenceModel<fluidThermo>> turbulence2
(
PhaseCompressibleTurbulenceModel<fluidThermo>::New
(
alpha2,
rho2,
U,
alphaRhoPhi2,
phi,
mixture.thermo2()
)
);
if (!LTS)
{
#include "readTimeControls.H"
#include "CourantNo.H"
#include "setInitialDeltaT.H"
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Info<< "\nStarting time loop\n" << endl;
while (runTime.run())
{
#include "readTimeControls.H"
if (LTS)
{
#include "setRDeltaT.H"
}
else
{
#include "CourantNo.H"
#include "alphaCourantNo.H"
#include "setDeltaT.H"
}
runTime++;
Info<< "Time = " << runTime.timeName() << nl << endl;
// --- Pressure-velocity PIMPLE corrector loop
while (pimple.loop())
{
#include "alphaControls.H"
#include "compressibleAlphaEqnSubCycle.H"
alphaRhoPhi1 = fvc::interpolate(rho1)*alphaPhi10;
alphaRhoPhi2 = fvc::interpolate(rho2)*(phi - alphaPhi10);
#include "UEqn.H"
#include "TEqn.H"
// --- Pressure corrector loop
while (pimple.correct())
{
#include "pEqn.H"
}
if (pimple.turbCorr())
{
turbulence1->correct();
turbulence2->correct();
}
}
runTime.write();
Info<< "ExecutionTime = "
<< runTime.elapsedCpuTime()
<< " s\n\n" << endl;
}
Info<< "End\n" << endl;
return 0;
}
// ************************************************************************* //