39 #ifndef MomentumTransferPhaseSystem_H
40 #define MomentumTransferPhaseSystem_H
51 class blendedDragModel;
52 class blendedVirtualMassModel;
53 class blendedLiftModel;
54 class blendedWallLubricationModel;
55 class blendedTurbulentDispersionModel;
61 template<
class BasePhaseSystem>
64 public BasePhaseSystem
A HashTable specialisation for hashing pointers.
An STL-conforming hash table.
Class which models interfacial momentum transfer between a number of phases. Drag,...
virtual void partialElimination(const PtrList< volScalarField > &rAUs, const PtrList< volVectorField > &KdUs, const PtrList< surfaceScalarField > &alphafs, const PtrList< surfaceScalarField > &rAUfs, const PtrList< surfaceScalarField > &KdPhis)
Solve the drag system for the velocities and fluxes.
virtual void dragCorrs(PtrList< volVectorField > &dragCorrs, PtrList< surfaceScalarField > &dragCorrf) const
Set the cell and faces drag correction fields.
virtual PtrList< surfaceScalarField > ddtCorrs() const
Return the flux corrections for the cell-based algorithm. These.
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransfer()
Return the momentum transfer matrices for the cell-based algorithm.
virtual tmp< surfaceScalarField > alphaDByAf(const PtrList< volScalarField > &rAUs, const PtrList< surfaceScalarField > &rAUfs) const
Return the phase diffusivity.
virtual PtrList< volVectorField > KdUs() const
Return the explicit part of the drag force for the cell-based.
virtual void partialEliminationf(const PtrList< surfaceScalarField > &rAUfs, const PtrList< surfaceScalarField > &alphafs, const PtrList< surfaceScalarField > &KdPhifs)
As partialElimination, but for the face-based algorithm. Only solves.
virtual PtrList< volScalarField > Kds() const
Return the implicit part of the drag force.
void addDmdtUfs(const phaseSystem::dmdtfTable &dmdtfs, phaseSystem::momentumTransferTable &eqns)
Add momentum transfer terms which result from bulk mass transfers.
virtual ~MomentumTransferPhaseSystem()
Destructor.
virtual PtrList< surfaceScalarField > Ffs() const
As Fs, but for the face-based algorithm.
virtual bool implicitPhasePressure() const
Returns true if the phase pressure is treated implicitly.
virtual PtrList< surfaceScalarField > Fs() const
Return the explicit force fluxes for the cell-based algorithm, that.
virtual PtrList< surfaceScalarField > KdVmfs() const
Return implicit force coefficients on the faces, for the face-based.
MomentumTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransferf()
As momentumTransfer, but for the face-based algorithm.
virtual PtrList< surfaceScalarField > KdPhifs() const
As KdPhis, but for the face-based algorithm.
virtual PtrList< surfaceScalarField > KdPhis() const
Return the explicit drag force fluxes for the cell-based algorithm.
virtual bool read()
Read base phaseProperties dictionary.
void addTmpField(tmp< surfaceScalarField > &result, const tmp< surfaceScalarField > &field) const
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.
Word-pair based class used for keying interface models in hash tables.
A class for managing temporary objects.
SurfaceField< scalar > surfaceScalarField
VolField< scalar > volScalarField