Go to the documentation of this file.
28 #include "BlendedInterfacialModel.H"
29 #include "dragModel.H"
30 #include "virtualMassModel.H"
31 #include "liftModel.H"
32 #include "wallLubricationModel.H"
33 #include "turbulentDispersionModel.H"
46 template<
class BasePhaseSystem>
55 this->generatePairsAndSubModels
61 this->generatePairsAndSubModels
67 this->generatePairsAndSubModels
73 this->generatePairsAndSubModels
76 wallLubricationModels_
79 this->generatePairsAndSubModels
81 "turbulentDispersion",
82 turbulentDispersionModels_
92 const phasePair& pair(this->phasePairs_[dragModelIter.key()]);
107 virtualMassModelTable,
112 const phasePair& pair(this->phasePairs_[virtualMassModelIter.key()]);
120 virtualMassModelIter()->
K()
129 template<
class BasePhaseSystem>
137 template<
class BasePhaseSystem>
141 const phasePairKey& key
144 return dragModels_[key]->K();
148 template<
class BasePhaseSystem>
152 const phasePairKey& key
155 return dragModels_[key]->Kf();
159 template<
class BasePhaseSystem>
166 tmp<volScalarField> tKd
173 this->mesh_.time().timeName(),
180 dimensionSet(1, -3, -1, 0, 0),
195 const phasePair& pair(this->phasePairs_[KdIter.key()]);
197 const phaseModel*
phase1 = &pair.phase1();
198 const phaseModel*
phase2 = &pair.phase2();
215 template<
class BasePhaseSystem>
219 const phasePairKey& key
222 if (virtualMassModels_.found(key))
224 return virtualMassModels_[key]->K();
228 return tmp<volScalarField>
234 virtualMassModel::typeName +
":K",
235 this->mesh_.time().timeName(),
249 template<
class BasePhaseSystem>
253 const phasePairKey& key
256 if (virtualMassModels_.found(key))
258 return virtualMassModels_[key]->Kf();
262 return tmp<surfaceScalarField>
268 virtualMassModel::typeName +
":Kf",
269 this->mesh_.time().timeName(),
283 template<
class BasePhaseSystem>
287 const phasePairKey& key
290 if (liftModels_.found(key) && wallLubricationModels_.found(key))
293 liftModels_[key]->template F<vector>()
294 + wallLubricationModels_[key]->template F<vector>();
296 else if (liftModels_.found(key))
298 return liftModels_[key]->template F<vector>();
300 else if (wallLubricationModels_.found(key))
302 return wallLubricationModels_[key]->template F<vector>();
306 return tmp<volVectorField>
312 liftModel::typeName +
":F",
313 this->mesh_.time().timeName(),
327 template<
class BasePhaseSystem>
331 const phasePairKey& key
334 if (liftModels_.found(key) && wallLubricationModels_.found(key))
337 liftModels_[key]->Ff()
338 + wallLubricationModels_[key]->Ff();
340 else if (liftModels_.found(key))
342 return liftModels_[key]->Ff();
344 else if (wallLubricationModels_.found(key))
346 return wallLubricationModels_[key]->Ff();
350 return tmp<surfaceScalarField>
356 liftModel::typeName +
":Ff",
357 this->mesh_.time().timeName(),
371 template<
class BasePhaseSystem>
375 const phasePairKey& key
378 if (turbulentDispersionModels_.found(key))
380 return turbulentDispersionModels_[key]->D();
384 return tmp<volScalarField>
390 turbulentDispersionModel::typeName +
":D",
391 this->mesh_.time().timeName(),
405 template<
class BasePhaseSystem>
410 autoPtr<phaseSystem::momentumTransferTable> eqnsPtr
419 const phaseModel& phase = this->phaseModels_[
phasei];
436 *Kds_[dragModelIter.key()] = dragModelIter()->K();
449 const phasePair& pair(this->phasePairs_[KdIter.key()]);
451 const phaseModel* phase = &pair.phase1();
452 const phaseModel* otherPhase = &pair.phase2();
460 Swap(phase, otherPhase);
467 virtualMassModelTable,
472 *Vms_[virtualMassModelIter.key()] = virtualMassModelIter()->K();
485 const phasePair& pair(this->phasePairs_[VmIter.key()]);
487 const phaseModel* phase = &pair.phase1();
488 const phaseModel* otherPhase = &pair.phase2();
495 *eqns[phase->name()] -=
503 + this->MRF_.DDt(Vm,
U - otherPhase->U());
505 Swap(phase, otherPhase);
513 template<
class BasePhaseSystem>
528 liftModel::typeName +
":F",
529 this->mesh_.time().timeName(),
545 template<
class BasePhaseSystem>
549 autoPtr<PtrList<volVectorField> > tFs
551 new PtrList<volVectorField>(this->
phases().size())
553 PtrList<volVectorField>& Fs = tFs();
565 const phasePair& pair(this->phasePairs_[liftModelIter.key()]);
567 setF(Fs, pair.phase1().index()) +=
F;
568 setF(Fs, pair.phase2().index()) -=
F;
574 wallLubricationModelTable,
575 wallLubricationModels_,
576 wallLubricationModelIter
582 pair(this->phasePairs_[wallLubricationModelIter.key()]);
584 setF(Fs, pair.phase1().index()) +=
F;
585 setF(Fs, pair.phase2().index()) -=
F;
592 template<
class BasePhaseSystem>
596 PtrList<surfaceScalarField>& phiDs,
const label phasei
608 turbulentDispersionModel::typeName +
":phiD",
609 this->mesh_.time().timeName(),
630 template<
class BasePhaseSystem>
634 const PtrList<volScalarField>&
rAUs
637 autoPtr<PtrList<surfaceScalarField> > tphiDs
639 new PtrList<surfaceScalarField>(this->
phases().size())
641 PtrList<surfaceScalarField>& phiDs = tphiDs();
646 turbulentDispersionModelTable,
647 turbulentDispersionModels_,
648 turbulentDispersionModelIter
652 pair(this->phasePairs_[turbulentDispersionModelIter.key()]);
660 setPhiD(phiDs, pair.phase1().index()) +=
662 setPhiD(phiDs, pair.phase2().index()) -=
670 template<
class BasePhaseSystem>
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
static word groupName(Name name, const word &group)
virtual tmp< volScalarField > D(const phasePairKey &key) const
Return the turbulent diffusivity.
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
const dimensionSet dimVelocity
const dimensionSet dimDensity
surfaceScalarField snGradAlpha1(fvc::snGrad(alpha1) *mesh.magSf())
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransfer() const
Return the momentum transfer matrices.
bool read(const char *, int32_t &)
Calculate the snGrad of the given volField.
virtual tmp< volScalarField > Kd(const phasePairKey &key) const
Return the drag coefficient.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
virtual ~MomentumTransferPhaseSystem()
Destructor.
static const dimensionSet dimF
Force dimensions.
Calculate the divergence of the given field.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
virtual autoPtr< PtrList< volVectorField > > Fs() const
Return the combined force (lift + wall-lubrication)
const dimensionedScalar F
Faraday constant: default SI units: [C/mol].
Calculate the matrix for the divergence of the given field and flux.
volVectorField & setF(PtrList< volVectorField > &Fs, const label phasei) const
Construct element phasei of Fs if not set and return.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
virtual tmp< surfaceScalarField > Kdf(const phasePairKey &key) const
Return the face drag coefficient.
CGAL::Exact_predicates_exact_constructions_kernel K
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
fvMatrix< vector > fvVectorMatrix
const dimensionSet dimArea(sqr(dimLength))
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const word & name() const
GeometricField< scalar, fvPatchField, volMesh > volScalarField
static const dimensionSet dimD
Diffusivity dimensions.
MomentumTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
static const dimensionSet dimF
Force dimensions.
HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > KdTable
Calculate the matrix for implicit and explicit sources.
virtual bool read()
Read base phaseProperties dictionary.
GeometricField< vector, fvPatchField, volMesh > volVectorField
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
virtual tmp< volVectorField > F(const phasePairKey &key) const
Return the combined force (lift + wall-lubrication)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
HashPtrTable< fvVectorMatrix, word, string::hash > momentumTransferTable
HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > VmTable
multiphaseSystem::phaseModelList & phases
virtual tmp< surfaceScalarField > Vmf(const phasePairKey &key) const
Return the face virtual mass coefficient.
surfaceScalarField & setPhiD(PtrList< surfaceScalarField > &phiDs, const label phasei) const
Construct element phasei of phiDs if not set and return.
virtual tmp< volScalarField > Vm(const phasePairKey &key) const
Return the virtual mass coefficient.
PtrList< volScalarField > rAUs(fluid.phases().size())
Generic GeometricField class.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
virtual autoPtr< PtrList< surfaceScalarField > > phiDs(const PtrList< volScalarField > &rAUs) const
Return the turbulent dispersion force on faces for phase pair.
Calulate the matrix for the first temporal derivative.
static const dimensionSet dimK
Coefficient dimensions.
virtual tmp< surfaceScalarField > Ff(const phasePairKey &key) const
Return the combined face-force (lift + wall-lubrication)