Files
OpenFOAM-12/applications/solvers/multiphase/reactingEulerFoam/phaseSystems/PhaseSystems/PhaseTransferPhaseSystem/PhaseTransferPhaseSystem.C
Will Bainbridge 85a9e17dd5 reactingEulerFoam: Added phase transfer structure
An additional layer has been added into the phase system hierarchy which
facilitates the application of phase transfer modelling. These are
models which exchange mass between phases without the thermal coupling
that would be required to represent phase change. They can be thought of
as representation changes; e.g., between two phases representing
different droplet sizes of the same physical fluid.

To facilitate this, the heat transfer phase systems have been modified
and renamed and now both support mass transfer. The two sided version
is only required for derivations which support phase change.

The following changes to case settings have been made:

- The simplest instantiated phase systems have been renamed to
basicTwoPhaseSystem and basicMultiphaseSystem. The
heatAndMomentumTransfer*System entries in constant/phaseProperties files
will need updating accordingly.

- A phaseTransfer sub-model entry will be required in the
constant/phaseProperties file. This can be an empty list.

- The massTransfer switch in thermal phase change cases has been renamed
phaseTransfer, so as not to be confused with the mass transfer models
used by interface composition cases.

This work was supported by Georg Skillas and Zhen Li, at Evonik
2018-04-05 15:11:39 +01:00

245 lines
6.1 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 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 "PhaseTransferPhaseSystem.H"
#include "phaseTransferModel.H"
#include "fvmSup.H"
// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
template<class BasePhaseSystem>
Foam::tmp<Foam::volScalarField>
Foam::PhaseTransferPhaseSystem<BasePhaseSystem>::rDmdt
(
const phasePairKey& key
) const
{
if (!rDmdt_.found(key))
{
return phaseSystem::dmdt(key);
}
const scalar rDmdtSign(Pair<word>::compare(rDmdt_.find(key).key(), key));
return rDmdtSign**rDmdt_[key];
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
template<class BasePhaseSystem>
Foam::PhaseTransferPhaseSystem<BasePhaseSystem>::PhaseTransferPhaseSystem
(
const fvMesh& mesh
)
:
BasePhaseSystem(mesh)
{
this->generatePairsAndSubModels
(
"phaseTransfer",
phaseTransferModels_,
false
);
forAllConstIter
(
phaseTransferModelTable,
phaseTransferModels_,
phaseTransferModelIter
)
{
this->rDmdt_.insert
(
phaseTransferModelIter.key(),
phaseSystem::dmdt(phaseTransferModelIter.key()).ptr()
);
}
}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
template<class BasePhaseSystem>
Foam::PhaseTransferPhaseSystem<BasePhaseSystem>::
~PhaseTransferPhaseSystem()
{}
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
template<class BasePhaseSystem>
Foam::tmp<Foam::volScalarField>
Foam::PhaseTransferPhaseSystem<BasePhaseSystem>::dmdt
(
const phasePairKey& key
) const
{
return BasePhaseSystem::dmdt(key) + this->rDmdt(key);
}
template<class BasePhaseSystem>
Foam::Xfer<Foam::PtrList<Foam::volScalarField>>
Foam::PhaseTransferPhaseSystem<BasePhaseSystem>::dmdts() const
{
PtrList<volScalarField> dmdts(BasePhaseSystem::dmdts());
forAllConstIter(rDmdtTable, rDmdt_, rDmdtIter)
{
const phasePair& pair = this->phasePairs_[rDmdtIter.key()];
const volScalarField& rDmdt = *rDmdtIter();
this->addField(pair.phase1(), "dmdt", rDmdt, dmdts);
this->addField(pair.phase2(), "dmdt", - rDmdt, dmdts);
}
return dmdts.xfer();
}
template<class BasePhaseSystem>
Foam::autoPtr<Foam::phaseSystem::massTransferTable>
Foam::PhaseTransferPhaseSystem<BasePhaseSystem>::massTransfer() const
{
// Create a mass transfer matrix for each species of each phase
autoPtr<phaseSystem::massTransferTable> eqnsPtr
(
new phaseSystem::massTransferTable()
);
phaseSystem::massTransferTable& eqns = eqnsPtr();
forAll(this->phaseModels_, phasei)
{
const phaseModel& phase = this->phaseModels_[phasei];
const PtrList<volScalarField>& Yi = phase.Y();
forAll(Yi, i)
{
eqns.insert
(
Yi[i].name(),
new fvScalarMatrix(Yi[i], dimMass/dimTime)
);
}
}
// Mass transfer across the interface
forAllConstIter
(
phaseTransferModelTable,
phaseTransferModels_,
phaseTransferModelIter
)
{
const phasePair& pair(this->phasePairs_[phaseTransferModelIter.key()]);
const phaseModel& phase = pair.phase1();
const phaseModel& otherPhase = pair.phase2();
const volScalarField dmdt(this->rDmdt(pair));
const volScalarField dmdt12(posPart(dmdt));
const volScalarField dmdt21(negPart(dmdt));
const PtrList<volScalarField>& Yi = phase.Y();
forAll(Yi, i)
{
const word name
(
IOobject::groupName(Yi[i].member(), phase.name())
);
const word otherName
(
IOobject::groupName(Yi[i].member(), otherPhase.name())
);
*eqns[name] +=
dmdt21*eqns[otherName]->psi()
- fvm::Sp(dmdt21, eqns[name]->psi());
*eqns[otherName] -=
dmdt12*eqns[name]->psi()
- fvm::Sp(dmdt12, eqns[otherName]->psi());
}
}
return eqnsPtr;
}
template<class BasePhaseSystem>
void Foam::PhaseTransferPhaseSystem<BasePhaseSystem>::correct()
{
BasePhaseSystem::correct();
forAllConstIter
(
phaseTransferModelTable,
phaseTransferModels_,
phaseTransferModelIter
)
{
*rDmdt_[phaseTransferModelIter.key()] =
dimensionedScalar("zero", dimDensity/dimTime, 0);
}
forAllConstIter
(
phaseTransferModelTable,
phaseTransferModels_,
phaseTransferModelIter
)
{
*rDmdt_[phaseTransferModelIter.key()] +=
phaseTransferModelIter()->dmdt();
}
}
template<class BasePhaseSystem>
bool Foam::PhaseTransferPhaseSystem<BasePhaseSystem>::read()
{
if (BasePhaseSystem::read())
{
bool readOK = true;
// Models ...
return readOK;
}
else
{
return false;
}
}
// ************************************************************************* //