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>
78 const phaseSystem::dmdtfTable& dmdtfs,
79 phaseSystem::specieTransferTable& eqns
84 const phaseInterface interface(*
this, dmdtfIter.key());
90 const phaseModel& phase1 = interface.phase1();
91 const phaseModel& phase2 = interface.phase2();
98 *eqns[Y1.name()] += dmdtf21*Y2 +
fvm::Sp(dmdtf12, Y1);
99 *eqns[Y2.name()] -= dmdtf12*Y1 +
fvm::Sp(dmdtf21, Y2);
105 template<
class BasePhaseSystem>
108 const phaseSystem::dmidtfTable& dmidtfs,
109 phaseSystem::specieTransferTable& eqns
114 const phaseInterface interface(*
this, dmidtfIter.key());
116 const phaseModel& phase1 = interface.phase1();
117 const phaseModel& phase2 = interface.phase2();
119 forAllConstIter(HashPtrTable<volScalarField>, *dmidtfIter(), dmidtfJter)
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_);
155 phaseTransferModelTable,
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().timeName(),
187 dmidtfs_.insert(interface,
new HashPtrTable<volScalarField>());
189 const hashedWordList species(phaseTransferModelIter()->species());
193 const word& specie = *specieIter;
195 dmidtfs_[interface]->insert
206 "phaseTransfer:dmidtf",
211 this->
mesh().time().timeName(),
225 template<
class BasePhaseSystem>
232 template<
class BasePhaseSystem>
236 const phaseInterfaceKey& key
239 tmp<volScalarField> tDmdtf = BasePhaseSystem::dmdtf(key);
241 if (phaseTransferModels_.found(key))
243 if (phaseTransferModels_[key]->mixture())
245 tDmdtf.ref() += *dmdtfs_[key];
250 HashPtrTable<volScalarField>,
255 tDmdtf.ref() += *dmidtfIter();
263 template<
class BasePhaseSystem>
267 PtrList<volScalarField> dmdts(BasePhaseSystem::dmdts());
269 autoPtr<phaseSystem::dmdtfTable> totalDmdtfsPtr = this->totalDmdtfs();
274 const phaseInterface interface(*
this, totalDmdtfIter.key());
276 addField(interface.phase1(),
"dmdt", *totalDmdtfIter(), dmdts);
277 addField(interface.phase2(),
"dmdt", - *totalDmdtfIter(), dmdts);
284 template<
class BasePhaseSystem>
288 autoPtr<phaseSystem::momentumTransferTable> eqnsPtr =
293 this->addDmdtUfs(totalDmdtfs(), eqns);
299 template<
class BasePhaseSystem>
303 autoPtr<phaseSystem::momentumTransferTable> eqnsPtr =
304 BasePhaseSystem::momentumTransferf();
308 this->addDmdtUfs(totalDmdtfs(), eqns);
314 template<
class BasePhaseSystem>
318 autoPtr<phaseSystem::heatTransferTable> eqnsPtr =
319 BasePhaseSystem::heatTransfer();
330 template<
class BasePhaseSystem>
334 autoPtr<phaseSystem::specieTransferTable> eqnsPtr
344 const phaseModel& phase = this->phaseModels_[
phasei];
346 const PtrList<volScalarField>& Yi = phase.Y();
358 this->addDmdtYfs(dmdtfs_, eqns);
359 this->addDmidtYf(dmidtfs_, eqns);
365 template<
class BasePhaseSystem>
373 phaseTransferModelTable,
374 phaseTransferModels_,
375 phaseTransferModelIter
378 const phaseInterface&
interface = phaseTransferModelIter()->interface();
380 if (phaseTransferModelIter()->mixture())
382 *dmdtfs_[interface] =
Zero;
385 const hashedWordList species(phaseTransferModelIter()->species());
389 const word& specie = *specieIter;
391 *(*dmidtfs_[interface])[specie] =
Zero;
398 phaseTransferModelTable,
399 phaseTransferModels_,
400 phaseTransferModelIter
403 const phaseInterface&
interface = phaseTransferModelIter()->interface();
405 if (phaseTransferModelIter()->mixture())
407 *dmdtfs_[interface] += phaseTransferModelIter()->dmdtf();
410 const HashPtrTable<volScalarField> dmidtf
412 phaseTransferModelIter()->dmidtf()
417 *(*dmidtfs_[interface])[dmidtfIter.key()] += *dmidtfIter();
423 template<
class BasePhaseSystem>
fvMatrix< scalar > fvScalarMatrix
#define forAll(list, i)
Loop across all elements in list.
PhaseTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-10/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
HashPtrTable< fvVectorMatrix > momentumTransferTable
virtual autoPtr< phaseSystem::heatTransferTable > heatTransfer() const
Return the heat transfer matrices.
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransferf()
Return the momentum transfer matrices for the face-based algorithm.
Class which models non-thermally-coupled or weakly thermally coupled mass transfers.
HashPtrTable< fvScalarMatrix > heatTransferTable
HashPtrTable< fvScalarMatrix > specieTransferTable
dimensionedScalar posPart(const dimensionedScalar &ds)
virtual void correct()
Correct the mass transfer rates.
void addDmdtHefs(const phaseSystem::dmdtfTable &dmdtfs, phaseSystem::heatTransferTable &eqns) const
Add energy transfer terms which result from bulk mass transfers.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
void addDmidtHefs(const phaseSystem::dmidtfTable &dmidtfs, phaseSystem::heatTransferTable &eqns) const
Add energy transfer terms which result from specie mass transfers.
const dimensionSet dimTime
virtual bool read()
Read base phaseProperties dictionary.
static word groupName(Name name, const word &group)
const dimensionSet dimDensity
void addDmdtYfs(const phaseSystem::dmdtfTable &dmdtfs, phaseSystem::specieTransferTable &eqns) const
Add specie transfer terms which result from bulk mass transfers.
HashPtrTable< volScalarField, phaseInterfaceKey, phaseInterfaceKey::hash > dmdtfTable
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
const dimensionSet dimMass
word name(const complex &)
Return a string representation of a complex.
phaseSystem::momentumTransferTable & momentumTransfer(momentumTransferPtr())
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransfer()
Return the momentum transfer matrices for the cell-based algorithm.
const tmp< volScalarField::Internal > & Sp
virtual autoPtr< phaseSystem::specieTransferTable > specieTransfer() const
Return the specie transfer matrices.
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.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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...
A class for managing temporary objects.
Calculate the matrix for implicit and explicit sources.
dimensionedScalar negPart(const dimensionedScalar &ds)
virtual ~PhaseTransferPhaseSystem()
Destructor.