The cell-base momentum/pressure algorithm in the multiphaseEuler solver module
has been substantially updated to improve consistency, conservation and reduce
drag generated staggering patterns at sharp interfaces and the boundaries with
stationary phases. For most if not all cases this new algorithm can be used to
provide well resolved and reliable solutions where the faceMomentum algorithm
would have been chosen previously in order to obtain sufficiently smooth
solutions but at the expense of a noticeable loss in accuracy and resolution.
The first significant change in the momentum/pressure algorithm is in the
interpolation practice used to construct the flux predictor equation from the
cell momentum equation: rather than interpolating the H/A ratio to the faces
i.e. (H/A)_f the terms in the momentum equation are interpolated separately so
that H_f/A_f is used. The same approach is used for the drag i.e. (D_f/A_f) and
virtual mass contributions. The advantage of this change is that the phase
forces are now consistent in both the momentum and flux equations, i.e. sum to
zero for each pair of phases.
The second significant change is in the handling of ddtCorr which converts the
old-time time-derivative contributions in H from velocity to flux which is now
consistent due to the change to H/A interpolation and also generalised to use
the fvc::ddtCorr function which has been updated for multiphase. Additionally
ddtCorr may optionally be applied to the time-derivative in the virtual mass
term in a consistent manner so that the contributions to the flux equation sum
to zero for each pair of phases.
The third significant change is the addition of an optional drag correction term
to the momentum corrector to reduce the staggering patters generated in the
velocity field due to sudden changes in drag force between phase, e.g. at sharp
interfaces between phases or at the boundaries with stationary phases. This is
particularly beneficial for fluidised bed simulations. However this correction
is not and cannot be phase consistent, i.e. the correction does not sum to zero
for pairs of phases it is applied to so a small drag error is introduced, but
tests so far have shown that the error is small and outweighed by the benefit in
the reduction in numerical artefacts in the solution.
The final significant change is in the handling of residualAlpha for drag and
virtual mass to provide stable and physical phase velocities in the limit of the
phase-fraction -> 0. The new approach is phase asymmetric such that the
residual drag is applied only to the phase with a phase-fraction less than
residualAlpha and not to the carrier phase. This change ensures that the flow
of a pure phase is unaffected by the residualAlpha and residual drag of the
other phases that are stabilised in pure phase region.
There are now four options in the PIMPLE section of the fvSolutions dictionary
relating to the multiphase momentum/pressure algorithm:
PIMPLE
{
faceMomentum no;
VmDdtCorrection yes;
dragCorrection yes;
partialElimination no;
}
faceMomentum:
Switches between the cell and face momentum equation algorithms.
Provides much smoother and reliable solutions for even the most challenging
multiphase cases at the expense of a noticeable loss in accuracy and resolution.
Defaults to 'no'.
VmDdtCorrection:
Includes the ddtCorr correction term to the time-derivative part of the
virtual-mass term in the flux equation which ensures consistency between the
phase virtual mass force on the faces but generates solutions which are
slightly less smooth and more likely to contain numerical artefacts.
Defaults to 'no'.
Testing so far has shown that the loss in smoothness is small and there is
some noticeable improvement is some cases so in the future the default may
be changed to 'yes'.
dragCorrection:
Includes the momentum corrector drag correction term to reduce the
staggering patters generated in the velocity field due to sudden changes in
drag force at the expense of a small error in drag consistency.
Defaults to 'no'
partialElimination:
Switches the partial-elimination momentum corrector which inverts the drag
matrix for both the momentum equations and/or flux equations to provide a
drag implicit correction to the phase velocity and flux fields. The
algorithm is the same as previously but updated for the new consistent drag
interpolation.
All the tutorials/modules/multiphaseEuler tutorial cases have been updated and
tested with the above developments and the four options set appropriately for
each.
586 lines
13 KiB
C++
586 lines
13 KiB
C++
/*---------------------------------------------------------------------------*\
|
|
========= |
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
\\ / O peration | Website: https://openfoam.org
|
|
\\ / A nd | Copyright (C) 2015-2023 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 "MovingPhaseModel.H"
|
|
#include "phaseSystem.H"
|
|
#include "fixedValueFvPatchFields.H"
|
|
#include "slipFvPatchFields.H"
|
|
#include "partialSlipFvPatchFields.H"
|
|
|
|
#include "fvmDdt.H"
|
|
#include "fvmDiv.H"
|
|
#include "fvmSup.H"
|
|
#include "fvcDdt.H"
|
|
#include "fvcDiv.H"
|
|
#include "fvcFlux.H"
|
|
|
|
// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
|
|
|
|
template<class BasePhaseModel>
|
|
Foam::tmp<Foam::surfaceScalarField>
|
|
Foam::MovingPhaseModel<BasePhaseModel>::phi(const volVectorField& U) const
|
|
{
|
|
word phiName(IOobject::groupName("phi", this->name()));
|
|
|
|
typeIOobject<surfaceScalarField> phiHeader
|
|
(
|
|
phiName,
|
|
U.mesh().time().name(),
|
|
U.mesh(),
|
|
IOobject::NO_READ
|
|
);
|
|
|
|
if (phiHeader.headerOk())
|
|
{
|
|
Info<< "Reading face flux field " << phiName << endl;
|
|
|
|
return tmp<surfaceScalarField>
|
|
(
|
|
new surfaceScalarField
|
|
(
|
|
IOobject
|
|
(
|
|
phiName,
|
|
U.mesh().time().name(),
|
|
U.mesh(),
|
|
IOobject::MUST_READ,
|
|
IOobject::AUTO_WRITE
|
|
),
|
|
U.mesh()
|
|
)
|
|
);
|
|
}
|
|
else
|
|
{
|
|
Info<< "Calculating face flux field " << phiName << endl;
|
|
|
|
wordList phiTypes
|
|
(
|
|
U.boundaryField().size(),
|
|
calculatedFvPatchScalarField::typeName
|
|
);
|
|
|
|
forAll(U.boundaryField(), patchi)
|
|
{
|
|
if (!U.boundaryField()[patchi].assignable())
|
|
{
|
|
phiTypes[patchi] = fixedValueFvPatchScalarField::typeName;
|
|
}
|
|
}
|
|
|
|
return tmp<surfaceScalarField>
|
|
(
|
|
new surfaceScalarField
|
|
(
|
|
IOobject
|
|
(
|
|
phiName,
|
|
U.mesh().time().name(),
|
|
U.mesh(),
|
|
IOobject::NO_READ,
|
|
IOobject::AUTO_WRITE
|
|
),
|
|
fvc::flux(U),
|
|
phiTypes
|
|
)
|
|
);
|
|
}
|
|
}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
|
|
|
template<class BasePhaseModel>
|
|
Foam::MovingPhaseModel<BasePhaseModel>::MovingPhaseModel
|
|
(
|
|
const phaseSystem& fluid,
|
|
const word& phaseName,
|
|
const bool referencePhase,
|
|
const label index
|
|
)
|
|
:
|
|
BasePhaseModel(fluid, phaseName, referencePhase, index),
|
|
U_
|
|
(
|
|
IOobject
|
|
(
|
|
IOobject::groupName("U", this->name()),
|
|
fluid.mesh().time().name(),
|
|
fluid.mesh(),
|
|
IOobject::MUST_READ,
|
|
IOobject::AUTO_WRITE
|
|
),
|
|
fluid.mesh()
|
|
),
|
|
phi_(phi(U_)),
|
|
alphaPhi_
|
|
(
|
|
IOobject
|
|
(
|
|
IOobject::groupName("alphaPhi", this->name()),
|
|
fluid.mesh().time().name(),
|
|
fluid.mesh()
|
|
),
|
|
fluid.mesh(),
|
|
dimensionedScalar(dimensionSet(0, 3, -1, 0, 0), 0)
|
|
),
|
|
alphaRhoPhi_
|
|
(
|
|
IOobject
|
|
(
|
|
IOobject::groupName("alphaRhoPhi", this->name()),
|
|
fluid.mesh().time().name(),
|
|
fluid.mesh()
|
|
),
|
|
fluid.mesh(),
|
|
dimensionedScalar(dimensionSet(1, 0, -1, 0, 0), 0)
|
|
),
|
|
Uf_(nullptr),
|
|
DUDt_(nullptr),
|
|
DUDtf_(nullptr),
|
|
divU_(nullptr),
|
|
momentumTransport_
|
|
(
|
|
phaseCompressible::momentumTransportModel::New
|
|
(
|
|
*this,
|
|
this->thermo().rho(),
|
|
U_,
|
|
alphaRhoPhi_,
|
|
phi_,
|
|
*this
|
|
)
|
|
),
|
|
thermophysicalTransport_
|
|
(
|
|
PhaseThermophysicalTransportModel
|
|
<
|
|
phaseCompressible::momentumTransportModel,
|
|
transportThermoModel
|
|
>::New(momentumTransport_, this->thermo_)
|
|
),
|
|
continuityError_
|
|
(
|
|
IOobject
|
|
(
|
|
IOobject::groupName("continuityError", this->name()),
|
|
fluid.mesh().time().name(),
|
|
fluid.mesh()
|
|
),
|
|
fluid.mesh(),
|
|
dimensionedScalar(dimDensity/dimTime, 0)
|
|
),
|
|
K_(nullptr)
|
|
{
|
|
phi_.writeOpt() = IOobject::AUTO_WRITE;
|
|
|
|
if (fluid.mesh().dynamic() || this->fluid().MRF().size())
|
|
{
|
|
Uf_ = new surfaceVectorField
|
|
(
|
|
IOobject
|
|
(
|
|
IOobject::groupName("Uf", this->name()),
|
|
fluid.mesh().time().name(),
|
|
fluid.mesh(),
|
|
IOobject::READ_IF_PRESENT,
|
|
IOobject::AUTO_WRITE
|
|
),
|
|
fvc::interpolate(U_)
|
|
);
|
|
}
|
|
|
|
correctKinematics();
|
|
}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
|
|
|
|
template<class BasePhaseModel>
|
|
Foam::MovingPhaseModel<BasePhaseModel>::~MovingPhaseModel()
|
|
{}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
|
|
|
template<class BasePhaseModel>
|
|
void Foam::MovingPhaseModel<BasePhaseModel>::correctContinuityError
|
|
(
|
|
const volScalarField& source
|
|
)
|
|
{
|
|
volScalarField& rho = this->thermoRef().rho();
|
|
|
|
continuityError_ = fvc::ddt(*this, rho) + fvc::div(alphaRhoPhi_) - source;
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
void Foam::MovingPhaseModel<BasePhaseModel>::correct()
|
|
{
|
|
BasePhaseModel::correct();
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
void Foam::MovingPhaseModel<BasePhaseModel>::correctKinematics()
|
|
{
|
|
BasePhaseModel::correctKinematics();
|
|
|
|
if (DUDt_.valid())
|
|
{
|
|
DUDt_.clear();
|
|
DUDt();
|
|
}
|
|
|
|
if (DUDtf_.valid())
|
|
{
|
|
DUDtf_.clear();
|
|
DUDtf();
|
|
}
|
|
|
|
if (K_.valid())
|
|
{
|
|
K_.ref() = 0.5*magSqr(this->U());
|
|
}
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
void Foam::MovingPhaseModel<BasePhaseModel>::predictMomentumTransport()
|
|
{
|
|
BasePhaseModel::predictMomentumTransport();
|
|
momentumTransport_->predict();
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
void Foam::MovingPhaseModel<BasePhaseModel>::predictThermophysicalTransport()
|
|
{
|
|
BasePhaseModel::predictThermophysicalTransport();
|
|
thermophysicalTransport_->predict();
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
void Foam::MovingPhaseModel<BasePhaseModel>::correctMomentumTransport()
|
|
{
|
|
BasePhaseModel::correctMomentumTransport();
|
|
momentumTransport_->correct();
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
void Foam::MovingPhaseModel<BasePhaseModel>::correctThermophysicalTransport()
|
|
{
|
|
BasePhaseModel::correctThermophysicalTransport();
|
|
thermophysicalTransport_->correct();
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
void Foam::MovingPhaseModel<BasePhaseModel>::correctUf()
|
|
{
|
|
const fvMesh& mesh = this->fluid().mesh();
|
|
|
|
if (Uf_.valid())
|
|
{
|
|
Uf_() = fvc::interpolate(U_);
|
|
surfaceVectorField n(mesh.Sf()/mesh.magSf());
|
|
Uf_() +=
|
|
n*(
|
|
this->fluid().MRF().absolute(fvc::absolute(phi_, U_))
|
|
/mesh.magSf()
|
|
- (n & Uf_())
|
|
);
|
|
}
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
bool Foam::MovingPhaseModel<BasePhaseModel>::stationary() const
|
|
{
|
|
return false;
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
Foam::tmp<Foam::fvVectorMatrix>
|
|
Foam::MovingPhaseModel<BasePhaseModel>::UEqn()
|
|
{
|
|
const volScalarField& alpha = *this;
|
|
const volScalarField& rho = this->thermo().rho();
|
|
|
|
return
|
|
(
|
|
fvm::ddt(alpha, rho, U_)
|
|
+ fvm::div(alphaRhoPhi_, U_)
|
|
+ fvm::SuSp(-this->continuityError(), U_)
|
|
+ this->fluid().MRF().DDt(alpha*rho, U_)
|
|
+ momentumTransport_->divDevTau(U_)
|
|
);
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
Foam::tmp<Foam::fvVectorMatrix>
|
|
Foam::MovingPhaseModel<BasePhaseModel>::UfEqn()
|
|
{
|
|
// As the "normal" U-eqn but without the ddt terms
|
|
|
|
const volScalarField& alpha = *this;
|
|
const volScalarField& rho = this->thermo().rho();
|
|
|
|
return
|
|
(
|
|
fvm::div(alphaRhoPhi_, U_)
|
|
+ fvm::SuSp(fvc::ddt(*this, rho) - this->continuityError(), U_)
|
|
+ this->fluid().MRF().DDt(alpha*rho, U_)
|
|
+ momentumTransport_->divDevTau(U_)
|
|
);
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
Foam::tmp<Foam::volVectorField>
|
|
Foam::MovingPhaseModel<BasePhaseModel>::U() const
|
|
{
|
|
return U_;
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
Foam::volVectorField&
|
|
Foam::MovingPhaseModel<BasePhaseModel>::URef()
|
|
{
|
|
return U_;
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
Foam::tmp<Foam::surfaceScalarField>
|
|
Foam::MovingPhaseModel<BasePhaseModel>::phi() const
|
|
{
|
|
return phi_;
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
Foam::surfaceScalarField&
|
|
Foam::MovingPhaseModel<BasePhaseModel>::phiRef()
|
|
{
|
|
return phi_;
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
const Foam::autoPtr<Foam::surfaceVectorField>&
|
|
Foam::MovingPhaseModel<BasePhaseModel>::Uf() const
|
|
{
|
|
return Uf_;
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
Foam::surfaceVectorField&
|
|
Foam::MovingPhaseModel<BasePhaseModel>::UfRef()
|
|
{
|
|
if (Uf_.valid())
|
|
{
|
|
return Uf_();
|
|
}
|
|
else
|
|
{
|
|
FatalErrorInFunction
|
|
<< "Uf has not been allocated."
|
|
<< exit(FatalError);
|
|
|
|
return const_cast<surfaceVectorField&>(surfaceVectorField::null());
|
|
}
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
Foam::tmp<Foam::surfaceScalarField>
|
|
Foam::MovingPhaseModel<BasePhaseModel>::alphaPhi() const
|
|
{
|
|
return alphaPhi_;
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
Foam::surfaceScalarField&
|
|
Foam::MovingPhaseModel<BasePhaseModel>::alphaPhiRef()
|
|
{
|
|
return alphaPhi_;
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
Foam::tmp<Foam::surfaceScalarField>
|
|
Foam::MovingPhaseModel<BasePhaseModel>::alphaRhoPhi() const
|
|
{
|
|
return alphaRhoPhi_;
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
Foam::surfaceScalarField&
|
|
Foam::MovingPhaseModel<BasePhaseModel>::alphaRhoPhiRef()
|
|
{
|
|
return alphaRhoPhi_;
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
Foam::tmp<Foam::volVectorField>
|
|
Foam::MovingPhaseModel<BasePhaseModel>::DUDt() const
|
|
{
|
|
if (!DUDt_.valid())
|
|
{
|
|
const tmp<surfaceScalarField> taphi(fvc::absolute(phi_, U_));
|
|
const surfaceScalarField& aphi(taphi());
|
|
DUDt_ =
|
|
new volVectorField
|
|
(
|
|
IOobject::groupName("DUDt", this->name()),
|
|
fvc::ddt(U_) + fvc::div(aphi, U_) - fvc::div(aphi)*U_
|
|
);
|
|
}
|
|
|
|
return tmp<volVectorField>(DUDt_());
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
Foam::tmp<Foam::surfaceScalarField>
|
|
Foam::MovingPhaseModel<BasePhaseModel>::DUDtf() const
|
|
{
|
|
if (!DUDtf_.valid())
|
|
{
|
|
DUDtf_ =
|
|
new surfaceScalarField
|
|
(
|
|
IOobject::groupName("DUDtf", this->name()),
|
|
byDt(phi_ - phi_.oldTime())
|
|
);
|
|
}
|
|
|
|
return tmp<surfaceScalarField>(DUDtf_());
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
Foam::tmp<Foam::volScalarField>
|
|
Foam::MovingPhaseModel<BasePhaseModel>::continuityError() const
|
|
{
|
|
return continuityError_;
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
Foam::tmp<Foam::volScalarField>
|
|
Foam::MovingPhaseModel<BasePhaseModel>::K() const
|
|
{
|
|
if (!K_.valid())
|
|
{
|
|
K_ =
|
|
new volScalarField
|
|
(
|
|
IOobject::groupName("K", this->name()),
|
|
0.5*magSqr(this->U())
|
|
);
|
|
}
|
|
|
|
return tmp<volScalarField>(K_());
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
Foam::tmp<Foam::volScalarField>
|
|
Foam::MovingPhaseModel<BasePhaseModel>::divU() const
|
|
{
|
|
return divU_.valid() ? tmp<volScalarField>(divU_()) : tmp<volScalarField>();
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
void Foam::MovingPhaseModel<BasePhaseModel>::divU(tmp<volScalarField> divU)
|
|
{
|
|
if (!divU_.valid())
|
|
{
|
|
divU_ = divU;
|
|
divU_.ref().rename(IOobject::groupName("divU", this->name()));
|
|
divU_.ref().checkIn();
|
|
}
|
|
else
|
|
{
|
|
divU_.ref() = divU;
|
|
}
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
Foam::tmp<Foam::volScalarField>
|
|
Foam::MovingPhaseModel<BasePhaseModel>::k() const
|
|
{
|
|
return momentumTransport_->k();
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
Foam::tmp<Foam::volScalarField>
|
|
Foam::MovingPhaseModel<BasePhaseModel>::pPrime() const
|
|
{
|
|
return momentumTransport_->pPrime();
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
Foam::tmp<Foam::scalarField>
|
|
Foam::MovingPhaseModel<BasePhaseModel>::kappaEff(const label patchi) const
|
|
{
|
|
return thermophysicalTransport_->kappaEff(patchi);
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
Foam::tmp<Foam::fvScalarMatrix>
|
|
Foam::MovingPhaseModel<BasePhaseModel>::divq(volScalarField& he) const
|
|
{
|
|
return thermophysicalTransport_->divq(he);
|
|
}
|
|
|
|
|
|
template<class BasePhaseModel>
|
|
Foam::tmp<Foam::fvScalarMatrix>
|
|
Foam::MovingPhaseModel<BasePhaseModel>::divj(volScalarField& Yi) const
|
|
{
|
|
return thermophysicalTransport_->divj(Yi);
|
|
}
|
|
|
|
|
|
// ************************************************************************* //
|