MomentumTransferPhaseSystem< BasePhaseSystem > Class Template Reference

Class which models interfacial momenum transfer between a number of phases. Drag, virtual mass, lift, wall lubrication and turbulent dispersion are all modelled. The explicit contribution from the drag is omitted from the transfer matrices, as this forms part of the solution of the pressure equation. More...

Inheritance diagram for MomentumTransferPhaseSystem< BasePhaseSystem >:
[legend]
Collaboration diagram for MomentumTransferPhaseSystem< BasePhaseSystem >:
[legend]

Public Member Functions

 MomentumTransferPhaseSystem (const fvMesh &)
 Construct from fvMesh. More...
 
virtual ~MomentumTransferPhaseSystem ()
 Destructor. More...
 
virtual autoPtr< phaseSystem::momentumTransferTablemomentumTransfer ()
 Return the momentum transfer matrices for the cell-based algorithm. More...
 
virtual autoPtr< phaseSystem::momentumTransferTablemomentumTransferf ()
 As momentumTransfer, but for the face-based algorithm. More...
 
virtual PtrList< surfaceScalarFieldAFfs () const
 Return implicit force coefficients on the faces, for the face-based. More...
 
virtual PtrList< surfaceScalarFieldphiFs (const PtrList< volScalarField > &rAUs)
 Return the explicit force fluxes for the cell-based algorithm, that. More...
 
virtual PtrList< surfaceScalarFieldphiFfs (const PtrList< surfaceScalarField > &rAUfs)
 As phiFs, but for the face-based algorithm. More...
 
virtual PtrList< surfaceScalarFieldphiKdPhis (const PtrList< volScalarField > &rAUs) const
 Return the explicit drag force fluxes for the cell-based algorithm. More...
 
virtual PtrList< surfaceScalarFieldphiKdPhifs (const PtrList< surfaceScalarField > &rAUfs) const
 As phiKdPhis, but for the face-based algorithm. More...
 
virtual PtrList< volVectorFieldKdUByAs (const PtrList< volScalarField > &rAUs) const
 Return the explicit part of the drag force for the cell-based. More...
 
virtual void partialElimination (const PtrList< volScalarField > &rAUs)
 Solve the drag system for the velocities and fluxes. More...
 
virtual void partialEliminationf (const PtrList< surfaceScalarField > &rAUfs)
 As partialElimination, but for the face-based algorithm. Only solves. More...
 
virtual PtrList< surfaceScalarFieldddtCorrByAs (const PtrList< volScalarField > &rAUs, const bool includeVirtualMass=false) const
 Return the flux corrections for the cell-based algorithm. These. More...
 
virtual const HashPtrTable< surfaceScalarField > & DByAfs () const
 Return the phase diffusivities divided by the momentum coefficients. More...
 
virtual bool read ()
 Read base phaseProperties dictionary. More...
 

Protected Types

typedef HashPtrTable< volScalarField, phasePairKey, phasePairKey::hashKdTable
 
typedef HashPtrTable< surfaceScalarField, phasePairKey, phasePairKey::hashKdfTable
 
typedef HashPtrTable< volScalarField, phasePairKey, phasePairKey::hashVmTable
 
typedef HashPtrTable< surfaceScalarField, phasePairKey, phasePairKey::hashVmfTable
 
typedef HashTable< autoPtr< BlendedInterfacialModel< dragModel > >, phasePairKey, phasePairKey::hashdragModelTable
 
typedef HashTable< autoPtr< BlendedInterfacialModel< virtualMassModel > >, phasePairKey, phasePairKey::hashvirtualMassModelTable
 
typedef HashTable< autoPtr< BlendedInterfacialModel< liftModel > >, phasePairKey, phasePairKey::hashliftModelTable
 
typedef HashTable< autoPtr< BlendedInterfacialModel< wallLubricationModel > >, phasePairKey, phasePairKey::hashwallLubricationModelTable
 
typedef HashTable< autoPtr< BlendedInterfacialModel< turbulentDispersionModel > >, phasePairKey, phasePairKey::hashturbulentDispersionModelTable
 

Detailed Description

template<class BasePhaseSystem>
class Foam::MomentumTransferPhaseSystem< BasePhaseSystem >

Class which models interfacial momenum transfer between a number of phases. Drag, virtual mass, lift, wall lubrication and turbulent dispersion are all modelled. The explicit contribution from the drag is omitted from the transfer matrices, as this forms part of the solution of the pressure equation.

Source files

Definition at line 66 of file MomentumTransferPhaseSystem.H.

Member Typedef Documentation

◆ KdTable

◆ KdfTable

◆ VmTable

◆ VmfTable

◆ dragModelTable

◆ virtualMassModelTable

◆ liftModelTable

◆ wallLubricationModelTable

◆ turbulentDispersionModelTable

Constructor & Destructor Documentation

◆ MomentumTransferPhaseSystem()

Construct from fvMesh.

Definition at line 170 of file MomentumTransferPhaseSystem.C.

References forAllConstIter.

◆ ~MomentumTransferPhaseSystem()

Destructor.

Definition at line 272 of file MomentumTransferPhaseSystem.C.

Member Function Documentation

◆ momentumTransfer()

Return the momentum transfer matrices for the cell-based algorithm.

This includes implicit and explicit forces that add into the cell UEqn in the normal way.

Definition at line 280 of file MomentumTransferPhaseSystem.C.

References Foam::fac::ddt(), Foam::dimMass, Foam::dimTime, Foam::dimVelocity, Foam::fac::div(), phaseModel::DUDt(), forAll, forAllConstIter, phase::name(), phaseModel::otherPhase(), phi, fvMatrix< Type >::psi(), HashPtrTable< T, Key, Hash >::set(), Sp, U, and phaseModel::U().

Here is the call graph for this function:

◆ momentumTransferf()

Foam::autoPtr< Foam::phaseSystem::momentumTransferTable > momentumTransferf ( )
virtual

As momentumTransfer, but for the face-based algorithm.

Definition at line 382 of file MomentumTransferPhaseSystem.C.

References Foam::dimMass, Foam::dimTime, Foam::dimVelocity, Foam::fac::div(), forAll, forAllConstIter, phaseModel::index(), MRF, phase::name(), phaseModel::otherPhase(), phasei, HashPtrTable< T, Key, Hash >::set(), Sp, U, and phaseModel::U().

Here is the call graph for this function:

◆ AFfs()

Foam::PtrList< Foam::surfaceScalarField > AFfs ( ) const
virtual

Return implicit force coefficients on the faces, for the face-based.

algorithm.

Definition at line 458 of file MomentumTransferPhaseSystem.C.

References AFfs(), Foam::byDt(), Foam::dimDensity, Foam::dimTime, fillFields(), forAllConstIter, and Vmf().

Here is the call graph for this function:

◆ phiFs()

Foam::PtrList< Foam::surfaceScalarField > phiFs ( const PtrList< volScalarField > &  rAUs)
virtual

Return the explicit force fluxes for the cell-based algorithm, that.

do not depend on phase mass/volume fluxes, and can therefore be evaluated outside the corrector loop. This includes things like lift, turbulent dispersion, and wall lubrication.

Definition at line 498 of file MomentumTransferPhaseSystem.C.

References D, Foam::dimDensity, Foam::dimForce, Foam::dimVelocity, F(), fillFields(), Foam::fvc::flux(), forAll, forAllConstIter, Foam::fac::interpolate(), phasei, phiFs(), rAUs, Foam::fvc::snGrad(), and snGradAlpha1().

Here is the call graph for this function:

◆ phiFfs()

Foam::PtrList< Foam::surfaceScalarField > phiFfs ( const PtrList< surfaceScalarField > &  rAUfs)
virtual

As phiFs, but for the face-based algorithm.

Definition at line 639 of file MomentumTransferPhaseSystem.C.

References Foam::fvc::absolute(), Foam::byDt(), D, Foam::dimDensity, Foam::dimForce, Foam::dimVelocity, Ff(), fillFields(), forAll, forAllConstIter, Foam::fac::interpolate(), MRF, oldTime, phasei, phi, phiFfs(), rAUfs, Foam::fvc::snGrad(), snGradAlpha1(), and Vmf().

Here is the call graph for this function:

◆ phiKdPhis()

Foam::PtrList< Foam::surfaceScalarField > phiKdPhis ( const PtrList< volScalarField > &  rAUs) const
virtual

Return the explicit drag force fluxes for the cell-based algorithm.

These depend on phase mass/volume fluxes, and must therefore be evaluated inside the corrector loop.

Definition at line 802 of file MomentumTransferPhaseSystem.C.

References Foam::fvc::absolute(), Foam::dimDensity, Foam::dimForce, Foam::dimVelocity, fillFields(), forAllConstIter, Foam::fac::interpolate(), MRF, and rAUs.

Here is the call graph for this function:

◆ phiKdPhifs()

Foam::PtrList< Foam::surfaceScalarField > phiKdPhifs ( const PtrList< surfaceScalarField > &  rAUfs) const
virtual

As phiKdPhis, but for the face-based algorithm.

Definition at line 844 of file MomentumTransferPhaseSystem.C.

References Foam::fvc::absolute(), Foam::dimDensity, Foam::dimForce, Foam::dimVelocity, fillFields(), forAllConstIter, MRF, and rAUfs.

Here is the call graph for this function:

◆ KdUByAs()

Foam::PtrList< Foam::volVectorField > KdUByAs ( const PtrList< volScalarField > &  rAUs) const
virtual

Return the explicit part of the drag force for the cell-based.

algorithm. This is the cell-equivalent of phiKdPhis. These depend on phase velocities, and must therefore be evaluated inside the corrector loop.

Definition at line 886 of file MomentumTransferPhaseSystem.C.

References Foam::dimVelocity, fillFields(), forAllConstIter, and rAUs.

Here is the call graph for this function:

◆ partialElimination()

void partialElimination ( const PtrList< volScalarField > &  rAUs)
virtual

Solve the drag system for the velocities and fluxes.

Definition at line 1021 of file MomentumTransferPhaseSystem.C.

References Foam::dimless, Foam::endl(), fillFields(), forAll, forAllConstIter, Foam::gMin(), phaseModel::index(), Foam::Info, Foam::fac::interpolate(), k, phasePair::phase1(), phasei, phases, GeometricField< Type, PatchField, GeoMesh >::primitiveField(), and rAUs.

Here is the call graph for this function:

◆ partialEliminationf()

void partialEliminationf ( const PtrList< surfaceScalarField > &  rAUfs)
virtual

As partialElimination, but for the face-based algorithm. Only solves.

for the fluxes.

Definition at line 1152 of file MomentumTransferPhaseSystem.C.

References Foam::dimless, Foam::endl(), fillFields(), forAll, forAllConstIter, Foam::gMin(), phaseModel::index(), Foam::Info, k, phasePair::phase1(), phasei, phases, GeometricField< Type, PatchField, GeoMesh >::primitiveField(), and rAUfs.

Here is the call graph for this function:

◆ ddtCorrByAs()

Foam::PtrList< Foam::surfaceScalarField > ddtCorrByAs ( const PtrList< volScalarField > &  rAUs,
const bool  includeVirtualMass = false 
) const
virtual

Return the flux corrections for the cell-based algorithm. These.

depend on phase mass/volume fluxes, and must therefore be evaluated inside the corrector loop.

Definition at line 922 of file MomentumTransferPhaseSystem.C.

References Foam::fvc::absolute(), Foam::constant::atomic::alpha, Foam::fac::average(), Foam::byDt(), Foam::fvc::flux(), forAll, forAllConstIter, phaseModel::index(), Foam::fac::interpolate(), MRF, oldTime, phaseModel::otherPhase(), phasei, phaseModel::phi(), Foam::pos0(), rAUs, tmp< T >::ref(), phase::rho(), and phaseModel::U().

Here is the call graph for this function:

◆ DByAfs()

const Foam::HashPtrTable< Foam::surfaceScalarField > & DByAfs ( ) const
virtual

Return the phase diffusivities divided by the momentum coefficients.

Definition at line 1251 of file MomentumTransferPhaseSystem.C.

◆ read()

bool read ( )
virtual

Read base phaseProperties dictionary.

Definition at line 1258 of file MomentumTransferPhaseSystem.C.

References Foam::blockMeshTools::read().

Here is the call graph for this function:

The documentation for this class was generated from the following files: