Go to the documentation of this file.
39 #ifndef MomentumTransferPhaseSystem_H
40 #define MomentumTransferPhaseSystem_H
50 template <
class modelType>
51 class BlendedInterfacialModel;
55 class virtualMassModel;
57 class wallLubricationModel;
58 class turbulentDispersionModel;
64 template<
class BasePhaseSystem>
67 public BasePhaseSystem
liftModelTable liftModels_
Lift models.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
HashTable< autoPtr< BlendedInterfacialModel< virtualMassModel > >, phasePairKey, phasePairKey::hash > virtualMassModelTable
virtual tmp< volScalarField > D(const phasePairKey &key) const
Return the turbulent diffusivity.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
A class for managing temporary objects.
wallLubricationModelTable wallLubricationModels_
Wall lubrication models.
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransfer() const
Return the momentum transfer matrices.
virtual tmp< volScalarField > Kd(const phasePairKey &key) const
Return the drag coefficient.
virtual ~MomentumTransferPhaseSystem()
Destructor.
HashTable< autoPtr< BlendedInterfacialModel< wallLubricationModel > >, phasePairKey, phasePairKey::hash > wallLubricationModelTable
dragModelTable dragModels_
Drag models.
virtual autoPtr< PtrList< volVectorField > > Fs() const
Return the combined force (lift + wall-lubrication)
virtual const phaseSystem::KdTable & Kds() const
Constant access to drag coefficients.
volVectorField & setF(PtrList< volVectorField > &Fs, const label phasei) const
Construct element phasei of Fs if not set and return.
virtual tmp< surfaceScalarField > Kdf(const phasePairKey &key) const
Return the face drag coefficient.
HashTable< autoPtr< BlendedInterfacialModel< dragModel > >, phasePairKey, phasePairKey::hash > dragModelTable
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
HashTable< autoPtr< BlendedInterfacialModel< turbulentDispersionModel > >, phasePairKey, phasePairKey::hash > turbulentDispersionModelTable
MomentumTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
Mesh data needed to do the Finite Volume discretisation.
phaseSystem::KdTable Kds_
Drag coefficients.
virtualMassModelTable virtualMassModels_
Virtual mass models.
phaseSystem::VmTable Vms_
Virtual mass coefficients.
virtual bool read()
Read base phaseProperties dictionary.
An STL-conforming hash table.
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)
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.
HashTable< autoPtr< BlendedInterfacialModel< liftModel > >, phasePairKey, phasePairKey::hash > liftModelTable
virtual tmp< volScalarField > Vm(const phasePairKey &key) const
Return the virtual mass coefficient.
PtrList< volScalarField > rAUs(fluid.phases().size())
Generic GeometricField class.
Class which models interfacial momenum transfer between a number of phases. Drag, virtual mass,...
virtual autoPtr< PtrList< surfaceScalarField > > phiDs(const PtrList< volScalarField > &rAUs) const
Return the turbulent dispersion force on faces for phase pair.
turbulentDispersionModelTable turbulentDispersionModels_
Turbulent dispersion models.
virtual tmp< surfaceScalarField > Ff(const phasePairKey &key) const
Return the combined face-force (lift + wall-lubrication)