Files
OpenFOAM-6/applications/solvers/multiphase/reactingEulerFoam/phaseSystems/phaseModel/MultiComponentPhaseModel/MultiComponentPhaseModel.C
Will Bainbridge e352828514 reactingMultiphaseEulerFoam: Stationary phase
Two new phase models have been added as selectable options for
reactingMultiphaseEulerFoam; pureStationaryPhaseModel and
pureStationaryIsothermalPhaseModel. These phases do not store a
velocity and their phase fractions remain constant throughout the
simulation. They are intended for use in modelling static particle beds
and other forms of porous media by means of the existing Euler-Euler
transfer models (drag, heat transfer, etc...).

Note that this functionality has not been extended to
reactingTwoPhaseEulerFoam, or the non-reacting *EulerFoam solvers.

Additional maintenance work has been carried out on the phase model
and phase system structure. The system can now loop over subsets of
phases with specific functionality (moving, multi-component, etc...) in
order to avoid testing for the existence of equations or variables in
the top level solver. The mass transfer handling and it's effect on
per-phase source terms has been refactored to reduce duplication. Const
and non-const access to phase properties has been formalised by renaming
non-const accessors with a "Ref" suffix, which is consistent with other
recent developments to classes including tmp and GeometricField, among
others. More sub-modelling details have been made private in order to
reduce the size of interfaces and improve abstraction.

This work was supported by Zhen Li, at Evonik
2018-03-23 09:08:52 +00:00

208 lines
5.0 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2015-2018 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 "MultiComponentPhaseModel.H"
#include "phaseSystem.H"
#include "fvmDdt.H"
#include "fvmDiv.H"
#include "fvmSup.H"
#include "fvmLaplacian.H"
#include "fvcDdt.H"
#include "fvcDiv.H"
#include "fvMatrix.H"
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
template<class BasePhaseModel>
Foam::MultiComponentPhaseModel<BasePhaseModel>::MultiComponentPhaseModel
(
const phaseSystem& fluid,
const word& phaseName,
const label index
)
:
BasePhaseModel(fluid, phaseName, index),
Sc_
(
"Sc",
dimless,
fluid.subDict(phaseName)
),
residualAlpha_
(
"residualAlpha",
dimless,
fluid.mesh().solverDict("Yi")
),
inertIndex_(-1)
{
const word inertSpecie
(
this->thermo_->lookupOrDefault("inertSpecie", word::null)
);
if (inertSpecie != word::null)
{
inertIndex_ = this->thermo_->composition().species()[inertSpecie];
}
PtrList<volScalarField>& Y = this->thermo_->composition().Y();
forAll(Y, i)
{
if (i != inertIndex_ && this->thermo_->composition().active(i))
{
const label j = YActive_.size();
YActive_.resize(j + 1);
YActive_.set(j, &Y[i]);
}
}
}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
template<class BasePhaseModel>
Foam::MultiComponentPhaseModel<BasePhaseModel>::~MultiComponentPhaseModel()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
template<class BasePhaseModel>
void Foam::MultiComponentPhaseModel<BasePhaseModel>::correctThermo()
{
volScalarField Yt
(
IOobject
(
IOobject::groupName("Yt", this->name()),
this->fluid().mesh().time().timeName(),
this->fluid().mesh()
),
this->fluid().mesh(),
dimensionedScalar("zero", dimless, 0)
);
PtrList<volScalarField>& Yi = YRef();
forAll(Yi, i)
{
if (i != inertIndex_)
{
Yt += Yi[i];
}
}
if (inertIndex_ != -1)
{
Yi[inertIndex_] = scalar(1) - Yt;
Yi[inertIndex_].max(0);
}
else
{
forAll(Yi, i)
{
Yi[i] /= Yt;
Yi[i].max(0);
}
}
BasePhaseModel::correctThermo();
}
template<class BasePhaseModel>
bool Foam::MultiComponentPhaseModel<BasePhaseModel>::pure() const
{
return false;
}
template<class BasePhaseModel>
Foam::tmp<Foam::fvScalarMatrix>
Foam::MultiComponentPhaseModel<BasePhaseModel>::YiEqn(volScalarField& Yi)
{
const volScalarField& alpha = *this;
const surfaceScalarField alphaRhoPhi(this->alphaRhoPhi());
const volScalarField rho(this->rho());
return
(
fvm::ddt(alpha, rho, Yi)
+ fvm::div(alphaRhoPhi, Yi, "div(" + alphaRhoPhi.name() + ",Yi)")
- fvm::Sp(this->continuityError(), Yi)
- fvm::laplacian
(
fvc::interpolate(alpha)
*fvc::interpolate(this->muEff()/Sc_),
Yi
)
==
alpha*this->R(Yi)
+ fvc::ddt(residualAlpha_*rho, Yi)
- fvm::ddt(residualAlpha_*rho, Yi)
);
}
template<class BasePhaseModel>
const Foam::PtrList<Foam::volScalarField>&
Foam::MultiComponentPhaseModel<BasePhaseModel>::Y() const
{
return this->thermo_->composition().Y();
}
template<class BasePhaseModel>
Foam::PtrList<Foam::volScalarField>&
Foam::MultiComponentPhaseModel<BasePhaseModel>::YRef()
{
return this->thermo_->composition().Y();
}
template<class BasePhaseModel>
const Foam::UPtrList<Foam::volScalarField>&
Foam::MultiComponentPhaseModel<BasePhaseModel>::YActive() const
{
return YActive_;
}
template<class BasePhaseModel>
Foam::UPtrList<Foam::volScalarField>&
Foam::MultiComponentPhaseModel<BasePhaseModel>::YActiveRef()
{
return YActive_;
}
// ************************************************************************* //