multiphaseEulerFoam::fvModels::phaseTurbulenceStabilisation: Phase turbulence stabilisation

Description
    Phase turbulence stabilisation

    In the limit of a phase-fraction->0 the turbulence properties cannot be
    obtained from the phase turbulence model, coupling to the other phase/phases
    is required.  The phaseTurbulenceStabilisation fvModel stabilises the phase
    turbulence properties by adding transfer terms from the corresponding
    properties of the other phases when the phase-fraction is less than the
    specified \c alphaInversion.  This implementation is a generalisation of
    the approach used in the Foam::RASModels::LaheyKEpsilon and
    Foam::RASModels::continuousGasKEpsilon models to handle phase-inversion and
    free-surface flow and can be used with any combination of RAS turbulence
    models.

    To stabilise the solution of the phase turbulence equations \c
    alphaInversion can be set to a small value e.g. 1e-2, but unless the phase
    turbulence model is specifically designed to handle phase-inversion and both
    continuous and dispersed regimes it may be useful to set \c alphaInversion
    to a higher value, corresponding to the phase-fraction at which transision
    from continuous to dispersed happens and effectively use the turbulence
    properties of the other phase when the phase is dispersed.  This is of
    course an approximation to the real system and if accurate handling of both
    the continuous and dispersed phase regimes is required specially developed
    models should be used.

Usage
    Example usage:
    \verbatim
    phaseTurbulenceStabilisation
    {
        type    phaseTurbulenceStabilisation;

        libs    ("libmultiphaseEulerFoamFvModels.so");

        phase   air;

        alphaInversion  0.1;
    }
    \endverbatim
This commit is contained in:
Henry Weller
2022-02-16 11:17:24 +00:00
parent 898924aa48
commit 9925df5407
6 changed files with 470 additions and 0 deletions

View File

@ -6,6 +6,7 @@ wclean libso interfacialModels
wclean libso interfacialCompositionModels
wclean libso multiphaseCompressibleMomentumTransportModels
multiphaseEulerFoam/Allwclean
wclean libso fvModels
wclean libso functionObjects
#------------------------------------------------------------------------------

View File

@ -10,6 +10,7 @@ wmake $targetType interfacialCompositionModels
wmake $targetType multiphaseCompressibleMomentumTransportModels
wmake $targetType multiphaseReactions
multiphaseEulerFoam/Allwmake $targetType $*
wmake $targetType fvModels
wmake $targetType functionObjects
#------------------------------------------------------------------------------

View File

@ -0,0 +1,4 @@
phaseTurbulenceStabilisation/phaseTurbulenceStabilisation.C
interfaceTurbulenceDamping/interfaceTurbulenceDamping.C
LIB = $(FOAM_LIBBIN)/libmultiphaseEulerFoamFvModels

View File

@ -0,0 +1,16 @@
EXE_INC = \
-I../phaseSystems/lnInclude \
-I$(LIB_SRC)/physicalProperties/lnInclude \
-I$(LIB_SRC)/thermophysicalModels/basic/lnInclude \
-I$(LIB_SRC)/MomentumTransportModels/momentumTransportModels/lnInclude \
-I$(LIB_SRC)/MomentumTransportModels/compressible/lnInclude \
-I$(LIB_SRC)/MomentumTransportModels/phaseCompressible/lnInclude \
-I$(LIB_SRC)/finiteVolume/lnInclude \
-I$(LIB_SRC)/meshTools/lnInclude
LIB_LIBS = \
-lphysicalProperties \
-lphaseSystem \
-lmultiphaseMomentumTransportModels \
-lfiniteVolume \
-lmeshTools

View File

@ -0,0 +1,255 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2022 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 "phaseTurbulenceStabilisation.H"
#include "phaseSystem.H"
#include "surfaceInterpolate.H"
#include "fvcGrad.H"
#include "fvmSup.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
namespace Foam
{
namespace fv
{
defineTypeNameAndDebug(phaseTurbulenceStabilisation, 0);
addToRunTimeSelectionTable
(
fvModel,
phaseTurbulenceStabilisation,
dictionary
);
}
}
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
void Foam::fv::phaseTurbulenceStabilisation::addSup
(
const volScalarField& alpha,
const volScalarField& rho,
fvMatrix<scalar>& eqn,
tmp<volScalarField>
(phaseCompressible::momentumTransportModel::*psi)() const
) const
{
if (debug)
{
Info<< type() << ": applying source to " << eqn.psi().name() << endl;
}
const fvMesh& mesh = this->mesh();
const phaseSystem::phaseModelPartialList& movingPhases =
phase_.fluid().movingPhases();
volScalarField::Internal transferRate
(
volScalarField::Internal::New
(
"transferRate",
mesh,
dimensionedScalar(dimless/dimTime, 0)
)
);
volScalarField::Internal psiTransferRate
(
volScalarField::Internal::New
(
"psiTransferRate",
mesh,
dimensionedScalar((turbulence_.*psi)()().dimensions()/dimTime, 0)
)
);
forAll(movingPhases, phasei)
{
if (movingPhases[phasei] != phase_)
{
const phaseCompressible::momentumTransportModel& turbulence =
mesh.lookupObject<phaseCompressible::momentumTransportModel>
(
IOobject::groupName
(
momentumTransportModel::typeName,
phaseName_
)
);
if (notNull(turbulence))
{
const volScalarField::Internal phaseTransferRate
(
movingPhases[phasei]
*min
(
turbulence.epsilon()/turbulence.k(),
1.0/phase_.time().deltaT()
)
);
transferRate += phaseTransferRate;
psiTransferRate += phaseTransferRate*(turbulence.*psi)()();
}
}
}
const volScalarField::Internal transferCoeff
(
max(alphaInversion_ - alpha(), scalar(0))*rho()
);
eqn += transferCoeff*psiTransferRate;
eqn -= fvm::Sp(transferCoeff*transferRate, eqn.psi());
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::fv::phaseTurbulenceStabilisation::phaseTurbulenceStabilisation
(
const word& sourceName,
const word& modelType,
const dictionary& dict,
const fvMesh& mesh
)
:
fvModel(sourceName, modelType, dict, mesh),
phaseName_(dict.lookup("phase")),
alphaInversion_("alphaInversion", dimless, dict),
phase_
(
mesh.lookupObject<phaseModel>(IOobject::groupName("alpha", phaseName_))
),
turbulence_
(
mesh.lookupObject<phaseCompressible::momentumTransportModel>
(
IOobject::groupName
(
momentumTransportModel::typeName,
phaseName_
)
)
)
{
const word kName(IOobject::groupName("k", phaseName_));
const word epsilonName(IOobject::groupName("epsilon", phaseName_));
const word omegaName(IOobject::groupName("omega", phaseName_));
if (mesh.foundObject<volScalarField>(kName))
{
fieldNames_.append(kName);
}
if (mesh.foundObject<volScalarField>(epsilonName))
{
fieldNames_.append(epsilonName);
}
if (mesh.foundObject<volScalarField>(omegaName))
{
fieldNames_.append(omegaName);
}
}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
Foam::wordList Foam::fv::phaseTurbulenceStabilisation::addSupFields() const
{
return fieldNames_;
}
void Foam::fv::phaseTurbulenceStabilisation::addSup
(
const volScalarField& alpha,
const volScalarField& rho,
fvMatrix<scalar>& eqn,
const word& fieldName
) const
{
if (fieldName == IOobject::groupName("k", phaseName_))
{
addSup
(
alpha,
rho,
eqn,
&phaseCompressible::momentumTransportModel::k
);
}
else if (fieldName == IOobject::groupName("epsilon", phaseName_))
{
addSup
(
alpha,
rho,
eqn,
&phaseCompressible::momentumTransportModel::epsilon
);
}
else if (fieldName == IOobject::groupName("omega", phaseName_))
{
addSup
(
alpha,
rho,
eqn,
&phaseCompressible::momentumTransportModel::omega
);
}
else
{
FatalErrorInFunction
<< "Support for field " << fieldName << " is not implemented"
<< exit(FatalError);
}
}
void Foam::fv::phaseTurbulenceStabilisation::updateMesh(const mapPolyMesh&)
{}
void Foam::fv::phaseTurbulenceStabilisation::distribute
(
const mapDistributePolyMesh&
)
{}
bool Foam::fv::phaseTurbulenceStabilisation::movePoints()
{
return true;
}
// ************************************************************************* //

View File

@ -0,0 +1,193 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2022 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/>.
Class
Foam::fv::phaseTurbulenceStabilisation
Description
Phase turbulence stabilisation
In the limit of a phase-fraction->0 the turbulence properties cannot be
obtained from the phase turbulence model, coupling to the other phase/phases
is required. The phaseTurbulenceStabilisation fvModel stabilises the phase
turbulence properties by adding transfer terms from the corresponding
properties of the other phases when the phase-fraction is less than the
specified \c alphaInversion. This implementation is a generalisation of
the approach used in the Foam::RASModels::LaheyKEpsilon and
Foam::RASModels::continuousGasKEpsilon models to handle phase-inversion and
free-surface flow and can be used with any combination of RAS turbulence
models.
To stabilise the solution of the phase turbulence equations \c
alphaInversion can be set to a small value e.g. 1e-2, but unless the phase
turbulence model is specifically designed to handle phase-inversion and both
continuous and dispersed regimes it may be useful to set \c alphaInversion
to a higher value, corresponding to the phase-fraction at which transision
from continuous to dispersed happens and effectively use the turbulence
properties of the other phase when the phase is dispersed. This is of
course an approximation to the real system and if accurate handling of both
the continuous and dispersed phase regimes is required specially developed
models should be used.
Usage
Example usage:
\verbatim
phaseTurbulenceStabilisation
{
type phaseTurbulenceStabilisation;
libs ("libmultiphaseEulerFoamFvModels.so");
phase air;
alphaInversion 0.1;
}
\endverbatim
SourceFiles
phaseTurbulenceStabilisation.C
\*---------------------------------------------------------------------------*/
#ifndef phaseTurbulenceStabilisation_H
#define phaseTurbulenceStabilisation_H
#include "fvModel.H"
#include "phaseModel.H"
#include "phaseCompressibleMomentumTransportModel.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace fv
{
/*---------------------------------------------------------------------------*\
Class phaseTurbulenceStabilisation Declaration
\*---------------------------------------------------------------------------*/
class phaseTurbulenceStabilisation
:
public fvModel
{
// Private Data
//- The name of the Lagrangian phase
word phaseName_;
//- Field names
wordList fieldNames_;
//- Phase-fraction below which turbulence property blending is applied
dimensionedScalar alphaInversion_;
//- Reference to the phase
const phaseModel& phase_;
//- Reference to the mixture turbulence model
const phaseCompressible::momentumTransportModel& turbulence_;
// Private Member Functions
//- Add contribution to phase psi equation
void addSup
(
const volScalarField& alpha,
const volScalarField& rho,
fvMatrix<scalar>& eqn,
tmp<volScalarField>
(phaseCompressible::momentumTransportModel::*psi)() const
) const;
public:
//- Runtime type information
TypeName("phaseTurbulenceStabilisation");
// Constructors
//- Construct from explicit source name and mesh
phaseTurbulenceStabilisation
(
const word& sourceName,
const word& modelType,
const dictionary& dict,
const fvMesh& mesh
);
//- Disallow default bitwise copy construction
phaseTurbulenceStabilisation
(
const phaseTurbulenceStabilisation&
) = delete;
// Member Functions
//- Return the list of fields for which the option adds source term
// to the transport equation
virtual wordList addSupFields() const;
//- Add contribution to phase k, epsilon or omega equation
virtual void addSup
(
const volScalarField& alpha,
const volScalarField& rho,
fvMatrix<scalar>& eqn,
const word& fieldName
) const;
// Mesh changes
//- Update for mesh changes
virtual void updateMesh(const mapPolyMesh&);
//- Update mesh corresponding to the given distribution map
virtual void distribute(const mapDistributePolyMesh&);
//- Update for mesh motion
virtual bool movePoints();
// Member Operators
//- Disallow default bitwise assignment
void operator=(const phaseTurbulenceStabilisation&) = delete;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace fv
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //