32 template<
class BasePhaseSystem>
36 autoPtr<phaseSystem::dmdtfTable> totalDmdtfsPtr
38 new phaseSystem::dmdtfTable
40 phaseSystem::dmdtfTable& totalDmdtfs = totalDmdtfsPtr();
44 phaseTransferModelTable,
46 phaseTransferModelIter
49 const phaseInterface&
interface = phaseTransferModelIter()->interface();
51 totalDmdtfs.insert(interface, phaseSystem::dmdtf(interface).ptr());
53 if (phaseTransferModelIter()->mixture())
55 *totalDmdtfs[interface] += *dmdtfs_[interface];
60 HashPtrTable<volScalarField>,
65 *totalDmdtfs[interface] += *dmidtfIter();
69 return totalDmdtfsPtr;
75 template<
class BasePhaseSystem>
98 *eqns[Y1.
name()] += dmdtf21*Y2 +
fvm::Sp(dmdtf12, Y1);
99 *eqns[Y2.
name()] -= dmdtf12*Y1 +
fvm::Sp(dmdtf21, Y2);
105 template<
class BasePhaseSystem>
121 const word& member = dmidtfJter.key();
128 *eqns[Y1.
name()] += dmidtf;
134 *eqns[Y2.
name()] -= dmidtf;
143 template<
class BasePhaseSystem>
149 BasePhaseSystem(mesh)
151 this->generateInterfacialModels(phaseTransferModels_);
156 phaseTransferModels_,
157 phaseTransferModelIter
160 const phaseInterface&
interface = phaseTransferModelIter()->interface();
162 this->
template validateMassTransfer<phaseTransferModel>(interface);
164 if (phaseTransferModelIter()->
mixture())
175 "phaseTransfer:dmdtf",
178 this->mesh().time().name(),
195 "phaseTransfer:d2mdtdpf",
198 this->mesh().time().name(),
209 const hashedWordList species(phaseTransferModelIter()->species());
215 dmidtfs_[interface]->
insert
226 "phaseTransfer:dmidtf",
231 this->mesh().time().name(),
245 template<
class BasePhaseSystem>
252 template<
class BasePhaseSystem>
261 if (phaseTransferModels_.found(key))
263 if (phaseTransferModels_[key]->
mixture())
265 tDmdtf.
ref() += *dmdtfs_[key];
275 tDmdtf.
ref() += *dmidtfIter();
283 template<
class BasePhaseSystem>
304 template<
class BasePhaseSystem>
314 addField(interface.
phase1(),
"d2mdtdp", *d2mdtdpfIter(), d2mdtdps);
315 addField(interface.
phase2(),
"d2mdtdp", - *d2mdtdpfIter(), d2mdtdps);
322 template<
class BasePhaseSystem>
327 BasePhaseSystem::momentumTransfer();
331 this->addDmdtUfs(totalDmdtfs(), eqns);
337 template<
class BasePhaseSystem>
342 BasePhaseSystem::momentumTransferf();
346 this->addDmdtUfs(totalDmdtfs(), eqns);
352 template<
class BasePhaseSystem>
357 BasePhaseSystem::heatTransfer();
361 this->addDmdtHefs(dmdtfs_, eqns);
362 this->addDmidtHefs(dmidtfs_, eqns);
368 template<
class BasePhaseSystem>
380 forAll(this->phaseModels_, phasei)
382 const phaseModel& phase = this->phaseModels_[phasei];
396 this->addDmdtYfs(dmdtfs_, eqns);
397 this->addDmidtYf(dmidtfs_, eqns);
403 template<
class BasePhaseSystem>
412 phaseTransferModels_,
413 phaseTransferModelIter
416 const phaseInterface&
interface = phaseTransferModelIter()->interface();
418 if (phaseTransferModelIter()->
mixture())
420 *dmdtfs_[interface] =
Zero;
421 *d2mdtdpfs_[interface] =
Zero;
424 const hashedWordList species(phaseTransferModelIter()->species());
438 phaseTransferModels_,
439 phaseTransferModelIter
442 const phaseInterface&
interface = phaseTransferModelIter()->interface();
444 if (phaseTransferModelIter()->
mixture())
446 *dmdtfs_[interface] += phaseTransferModelIter()->dmdtf();
447 *d2mdtdpfs_[interface] += phaseTransferModelIter()->d2mdtdpf();
452 phaseTransferModelIter()->dmidtf()
457 *(*dmidtfs_[interface])[dmidtfIter.key()] += *dmidtfIter();
463 template<
class BasePhaseSystem>
#define forAll(list, i)
Loop across all elements in list.
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Generic GeometricField class.
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
static word member(const word &name)
Return member (name without the extension)
const word & name() const
Return name.
static word groupName(Name name, const word &group)
Class which models non-thermally-coupled or weakly thermally coupled mass transfers.
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransfer()
Return the momentum transfer matrices for the cell-based algorithm.
virtual autoPtr< phaseSystem::heatTransferTable > heatTransfer() const
Return the heat transfer matrices.
virtual void correct()
Correct the mass transfer rates.
void addDmidtYf(const phaseSystem::dmidtfTable &dmidtfs, phaseSystem::specieTransferTable &eqns) const
Add specie transfer terms which result from specie mass transfers.
virtual tmp< volScalarField > dmdtf(const phaseInterfaceKey &key) const
Return the mass transfer rate for an interface.
virtual autoPtr< phaseSystem::specieTransferTable > specieTransfer() const
Return the specie transfer matrices.
PhaseTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
virtual ~PhaseTransferPhaseSystem()
Destructor.
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransferf()
Return the momentum transfer matrices for the face-based algorithm.
void addDmdtYfs(const phaseSystem::dmdtfTable &dmdtfs, phaseSystem::specieTransferTable &eqns) const
Add specie transfer terms which result from bulk mass transfers.
virtual PtrList< volScalarField > d2mdtdps() const
Return the mass transfer pressure linearisation coefficients for.
virtual bool read()
Read base phaseProperties dictionary.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Mesh data needed to do the Finite Volume discretisation.
A wordList with hashed indices for faster lookup by name.
Word-pair based class used for keying interface models in hash tables.
Class to represent an interface between phases. Derivations can further specify the configuration of ...
const phaseModel & phase1() const
Return phase 1.
const phaseModel & phase2() const
Return phase 2.
virtual bool pure() const =0
Return whether the phase is pure (i.e., not multi-component)
virtual const PtrList< volScalarField > & Y() const =0
Return the species mass fractions.
Base class of the thermophysical property types.
A class for managing temporary objects.
T & ref() const
Return non-const reference or generate a fatal error.
A class for handling words, derived from string.
Calculate the matrix for implicit and explicit sources.
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su)
tmp< VolField< Type > > Sp(const volScalarField &sp, const VolField< Type > &vf)
const dimensionSet dimPressure
word name(const bool)
Return a word representation of a bool.
void addField(const label phasei, const word &name, tmp< GeoField > field, PtrList< GeoField > &fieldList)
fvMatrix< scalar > fvScalarMatrix
const dimensionSet dimTime
dimensionedScalar negPart(const dimensionedScalar &ds)
const dimensionSet dimDensity
const dimensionSet dimMass
dimensionedScalar posPart(const dimensionedScalar &ds)