Go to the documentation of this file.
41 #ifndef MomentumTransferPhaseSystem_H
42 #define MomentumTransferPhaseSystem_H
44 #include "phaseSystem.H"
52 template<
class modelType>
53 class BlendedInterfacialModel;
57 class virtualMassModel;
59 class wallLubricationModel;
60 class turbulentDispersionModel;
66 template<
class BasePhaseSystem>
69 public BasePhaseSystem
188 void addMassTransferMomentumTransfer
276 const bool includeVirtualMass =
false
virtual PtrList< surfaceScalarField > phiFfs(const PtrList< surfaceScalarField > &rAUfs)
PtrList< volScalarField > rAUs
virtual ~MomentumTransferPhaseSystem()
HashTable< autoPtr< BlendedInterfacialModel< virtualMassModel > >, phasePairKey, phasePairKey::hash > virtualMassModelTable
A class for managing temporary objects.
HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > VmTable
virtual PtrList< volVectorField > KdUByAs(const PtrList< volScalarField > &rAUs) const
auto key(const Type &t) -> typename std::enable_if< std::is_enum< Type >::value, typename std::underlying_type< Type >::type >::type
HashTable< autoPtr< BlendedInterfacialModel< wallLubricationModel > >, phasePairKey, phasePairKey::hash > wallLubricationModelTable
HashTable< autoPtr< BlendedInterfacialModel< dragModel > >, phasePairKey, phasePairKey::hash > dragModelTable
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
GeometricField< scalar, fvPatchField, volMesh > volScalarField
An ordered or unorder pair of phase names. Typically specified as follows.
virtual void partialElimination(const PtrList< volScalarField > &rAUs)
virtual PtrList< surfaceScalarField > phiKdPhis(const PtrList< volScalarField > &rAUs) const
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransferf()
HashTable< autoPtr< BlendedInterfacialModel< turbulentDispersionModel > >, phasePairKey, phasePairKey::hash > turbulentDispersionModelTable
MomentumTransferPhaseSystem(const fvMesh &)
Mesh data needed to do the Finite Volume discretisation.
virtual PtrList< surfaceScalarField > phiKdPhifs(const PtrList< surfaceScalarField > &rAUfs) const
A HashTable similar to std::unordered_map.
virtual PtrList< surfaceScalarField > phiFs(const PtrList< volScalarField > &rAUs)
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > KdTable
HashPtrTable< surfaceScalarField, phasePairKey, phasePairKey::hash > KdfTable
A HashTable of pointers to objects of type <T>, with deallocation management of the pointers.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
PtrList< surfaceScalarField > rAUfs
virtual void partialEliminationf(const PtrList< surfaceScalarField > &rAUfs)
HashTable< autoPtr< BlendedInterfacialModel< liftModel > >, phasePairKey, phasePairKey::hash > liftModelTable
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransfer()
virtual PtrList< surfaceScalarField > AFfs() const
HashPtrTable< surfaceScalarField, phasePairKey, phasePairKey::hash > VmfTable
Class which models interfacial momenum transfer between a number of phases. Drag, virtual mass,...
virtual PtrList< surfaceScalarField > ddtCorrByAs(const PtrList< volScalarField > &rAUs, const bool includeVirtualMass=false) const
virtual const HashPtrTable< surfaceScalarField > & DByAfs() const