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

BUG: compressibleInterFoam family: Corrected transonic option

Resolves bug-report https://bugs.openfoam.org/view.php?id=2785

ENH: compressibleInterFoam family: merged two-phase momentum stress modelling from compressibleInterPhaseTransportFoam

The new momentum stress model selector class
compressibleInterPhaseTransportModel is now used to select between the options:

Description
    Transport model selection class for the compressibleInterFoam family of
    solvers.

    By default the standard mixture transport modelling approach is used in
    which a single momentum stress model (laminar, non-Newtonian, LES or RAS) is
    constructed for the mixture.  However if the \c simulationType in
    constant/turbulenceProperties is set to \c twoPhaseTransport the alternative
    Euler-Euler two-phase transport modelling approach is used in which separate
    stress models (laminar, non-Newtonian, LES or RAS) are instantiated for each
    of the two phases allowing for different modeling for the phases.

Mixture and two-phase momentum stress modelling is now supported in
compressibleInterFoam, compressibleInterDyMFoam and compressibleInterFilmFoam.
The prototype compressibleInterPhaseTransportFoam solver is no longer needed and
has been removed.
parent e061de2c
......@@ -10,6 +10,11 @@ IOobject alphaPhi10Header
const bool alphaRestart =
alphaPhi10Header.typeHeaderOk<surfaceScalarField>(true);
if (alphaRestart)
{
Info << "Restarting alpha" << endl;
}
// MULES flux from previous time-step
surfaceScalarField alphaPhi10
(
......
......@@ -3,6 +3,8 @@ cd ${0%/*} || exit 1 # Run from this directory
wclean libso twoPhaseMixtureThermo
wclean libso surfaceTensionModels
wclean libso VoFphaseCompressibleTurbulenceModels
wclean
wclean compressibleInterDyMFoam
wclean compressibleInterFilmFoam
......
......@@ -6,10 +6,10 @@ cd ${0%/*} || exit 1 # Run from this directory
wmake $targetType twoPhaseMixtureThermo
wmake $targetType surfaceTensionModels
wmake $targetType VoFphaseCompressibleTurbulenceModels
wmake $targetType
wmake $targetType compressibleInterDyMFoam
wmake $targetType compressibleInterFilmFoam
compressibleInterPhaseTransportFoam/Allwmake $targetType $*
#------------------------------------------------------------------------------
......@@ -8,6 +8,8 @@ EXE_INC = \
-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
......@@ -22,6 +24,7 @@ EXE_LIBS = \
-linterfaceProperties \
-lturbulenceModels \
-lcompressibleTurbulenceModels \
-lVoFphaseCompressibleTurbulenceModels \
-lfiniteVolume \
-lfvOptions \
-lmeshTools
......@@ -3,7 +3,7 @@
(
fvm::ddt(rho, T) + fvm::div(rhoPhi, T)
- fvm::Sp(contErr, T)
- fvm::laplacian(mixture.alphaEff(turbulence->mut()), T)
- fvm::laplacian(turbulence.alphaEff(), T)
+ (
divU*p
+ fvc::ddt(rho, K) + fvc::div(rhoPhi, K)
......
......@@ -3,7 +3,7 @@
fvm::ddt(rho, U) + fvm::div(rhoPhi, U)
- fvm::Sp(contErr, U)
+ MRF.DDt(rho, U)
+ turbulence->divDevRhoReff(U)
+ turbulence.divDevRhoReff(U)
==
fvOptions(rho, U)
);
......
EXE_INC = \
-I../twoPhaseMixtureThermo \
-I$(LIB_SRC)/transportModels/compressible/lnInclude \
-I$(LIB_SRC)/thermophysicalModels/basic/lnInclude \
-I$(LIB_SRC)/transportModels/incompressible/transportModel \
-I$(LIB_SRC)/transportModels/twoPhaseMixture/lnInclude \
-I$(LIB_SRC)/transportModels/interfaceProperties/lnInclude \
-I$(LIB_SRC)/TurbulenceModels/compressible/lnInclude \
-I$(LIB_SRC)/TurbulenceModels/turbulenceModels/lnInclude \
-I$(LIB_SRC)/TurbulenceModels/phaseCompressible/lnInclude \
......@@ -9,12 +11,15 @@ EXE_INC = \
-I$(LIB_SRC)/meshTools/lnInclude
LIB_LIBS = \
-ltwoPhaseMixtureThermo \
-lcompressibleTransportModels \
-lfluidThermophysicalModels \
-lspecie \
-ltwoPhaseMixture \
-ltwoPhaseProperties \
-linterfaceProperties \
-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/>.
\*---------------------------------------------------------------------------*/
#include "compressibleInterPhaseTransportModel.H"
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::compressibleInterPhaseTransportModel::compressibleInterPhaseTransportModel
(
const volScalarField& rho,
const volVectorField& U,
const surfaceScalarField& phi,
const surfaceScalarField& rhoPhi,
const surfaceScalarField& alphaPhi10,
const twoPhaseMixtureThermo& mixture
)
:
twoPhaseTransport_(false),
mixture_(mixture),
phi_(phi),
alphaPhi10_(alphaPhi10)
{
{
IOdictionary turbulenceProperties
(
IOobject
(
turbulenceModel::propertiesName,
U.time().constant(),
U.db(),
IOobject::MUST_READ,
IOobject::NO_WRITE
)
);
word simulationType
(
turbulenceProperties.lookup("simulationType")
);
if (simulationType == "twoPhaseTransport")
{
twoPhaseTransport_ = true;
}
}
if (twoPhaseTransport_)
{
const volScalarField& alpha1(mixture_.alpha1());
const volScalarField& alpha2(mixture_.alpha2());
const volScalarField& rho1 = mixture_.thermo1().rho();
const volScalarField& rho2 = mixture_.thermo2().rho();
alphaRhoPhi1_ =
(
new surfaceScalarField
(
IOobject::groupName("alphaRhoPhi", alpha1.group()),
fvc::interpolate(rho1)*alphaPhi10_
)
);
alphaRhoPhi2_ =
(
new surfaceScalarField
(
IOobject::groupName("alphaRhoPhi", alpha2.group()),
fvc::interpolate(rho2)*(phi_ - alphaPhi10_)
)
);
turbulence1_ =
(
PhaseCompressibleTurbulenceModel<fluidThermo>::New
(
alpha1,
rho1,
U,
alphaRhoPhi1_(),
phi,
mixture.thermo1()
)
);
turbulence2_ =
(
PhaseCompressibleTurbulenceModel<fluidThermo>::New
(
alpha2,
rho2,
U,
alphaRhoPhi2_(),
phi,
mixture.thermo2()
)
);
}
else
{
turbulence_ = compressible::turbulenceModel::New
(
rho,
U,
rhoPhi,
mixture
);
turbulence_->validate();
}
}
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
Foam::tmp<Foam::volScalarField>
Foam::compressibleInterPhaseTransportModel::alphaEff() const
{
if (twoPhaseTransport_)
{
return
mixture_.alpha1()*mixture_.thermo1().alphaEff(turbulence1_->mut())
+ mixture_.alpha2()*mixture_.thermo2().alphaEff(turbulence2_->mut());
}
else
{
return turbulence_->mut();
}
}
Foam::tmp<Foam::fvVectorMatrix>
Foam::compressibleInterPhaseTransportModel::divDevRhoReff
(
volVectorField& U
) const
{
if (twoPhaseTransport_)
{
return
turbulence1_->divDevRhoReff(U)
+ turbulence2_->divDevRhoReff(U);
}
else
{
return turbulence_->divDevRhoReff(U);
}
}
void Foam::compressibleInterPhaseTransportModel::correctPhasePhi()
{
if (twoPhaseTransport_)
{
const volScalarField& rho1 = mixture_.thermo1().rho();
const volScalarField& rho2 = mixture_.thermo2().rho();
alphaRhoPhi1_.ref() = fvc::interpolate(rho1)*alphaPhi10_;
alphaRhoPhi2_.ref() = fvc::interpolate(rho2)*(phi_ - alphaPhi10_);
}
}
void Foam::compressibleInterPhaseTransportModel::correct()
{
if (twoPhaseTransport_)
{
turbulence1_->correct();
turbulence2_->correct();
}
else
{
turbulence_->correct();
}
}
// ************************************************************************* //
......@@ -21,160 +21,126 @@ License
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
Class
Foam::compressibleInterPhaseTransportModel
Description
Solver for 2 compressible, non-isothermal immiscible fluids using a VOF
(volume of fluid) phase-fraction based interface capturing approach.
Transport model selection class for the compressibleInterFoam family of
solvers.
The momentum and other fluid properties are of the "mixture" and a single
momentum equation is solved.
By default the standard mixture transport modelling approach is used in
which a single momentum stress model (laminar, non-Newtonian, LES or RAS) is
constructed for the mixture. However if the \c simulationType in
constant/turbulenceProperties is set to \c twoPhaseTransport the alternative
Euler-Euler two-phase transport modelling approach is used in which separate
stress models (laminar, non-Newtonian, LES or RAS) are instantiated for each
of the two phases allowing for different modeling for the phases.
The fluid stress modelling is generic Euler-Euler two-phase in which
separate laminar, RAS or LES are selected for each of the phases.
SourceFiles
compressibleInterPhaseTransportModel.C
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
#include "CMULES.H"
#include "EulerDdtScheme.H"
#include "localEulerDdtScheme.H"
#include "CrankNicolsonDdtScheme.H"
#include "subCycle.H"
#include "rhoThermo.H"
#ifndef compressibleInterPhaseTransportModel_H
#define compressibleInterPhaseTransportModel_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[])
namespace Foam
{
#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
/*---------------------------------------------------------------------------*\
Class compressibleInterPhaseTransportModel Declaration
\*---------------------------------------------------------------------------*/
class compressibleInterPhaseTransportModel
{
// Private data
//- Switch to select two-phase or mixture transport modelling
Switch twoPhaseTransport_;
//- Two-phase mixture
const twoPhaseMixtureThermo& mixture_;
//- Mixture volumetric flux
const surfaceScalarField& phi_;
//- Phase volumetric flux
const surfaceScalarField& alphaPhi10_;
//- Phase-1 mass-flux (constructed for two-phase transport)
tmp<surfaceScalarField> alphaRhoPhi1_;
//- Phase-2 mass-flux (constructed for two-phase transport)
tmp<surfaceScalarField> alphaRhoPhi2_;
//- Mixture transport model (constructed for mixture transport)
autoPtr<compressible::turbulenceModel> turbulence_;
//- Phase-1 transport model (constructed for two-phase transport)
autoPtr<PhaseCompressibleTurbulenceModel<fluidThermo>> turbulence1_;
//- Phase-2 transport model (constructed for two-phase transport)
autoPtr<PhaseCompressibleTurbulenceModel<fluidThermo>> turbulence2_;
// Private Member Functions
//- Disallow default bitwise copy construct
compressibleInterPhaseTransportModel
(
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
const compressibleInterPhaseTransportModel&
);
//- Disallow default bitwise assignment
void operator=(const compressibleInterPhaseTransportModel&);
public:
// Constructors
//- Construct from components
compressibleInterPhaseTransportModel
(
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;
}
const volScalarField& rho,
const volVectorField& U,
const surfaceScalarField& phi,
const surfaceScalarField& rhoPhi,
const surfaceScalarField& alphaPhi10,
const twoPhaseMixtureThermo& mixture
);