Files
OpenFOAM-6/applications/solvers/multiphase/reactingEulerFoam/phaseSystems/PhaseSystems/PhaseTransferPhaseSystem/PhaseTransferPhaseSystem.C
Will Bainbridge c259eac3e7 reactingEulerFoam: Removed continuity error from species equations
This fixes a consistency issue in the interface-composition method, and
also seems to improve stability/convergence of the pimple iteration in
the presence of significant mass transfer.
2018-05-29 15:15:57 +01:00

248 lines
6.2 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();
// Note that the phase YiEqn does not contain a continuity error term,
// so these additions represent the entire mass transfer
const volScalarField dmdt(this->rDmdt(pair));
const volScalarField dmdt12(negPart(dmdt));
const volScalarField dmdt21(posPart(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(dmdt12, eqns[name]->psi());
*eqns[otherName] -=
dmdt12*eqns[name]->psi()
+ fvm::Sp(dmdt21, 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;
}
}
// ************************************************************************* //