Loading...
Searching...
No Matches
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...

#include <MomentumTransferPhaseSystem.H>

Inheritance diagram for MomentumTransferPhaseSystem< BasePhaseSystem >:
Collaboration diagram for MomentumTransferPhaseSystem< BasePhaseSystem >:

Public Member Functions

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

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 62 of file MomentumTransferPhaseSystem.H.

Member Typedef Documentation

◆ KdTable

template<class BasePhaseSystem>
typedef HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > KdTable
protected

Definition at line 75 of file MomentumTransferPhaseSystem.H.

◆ KdfTable

template<class BasePhaseSystem>
typedef HashPtrTable< surfaceScalarField, phasePairKey, phasePairKey::hash > KdfTable
protected

Definition at line 82 of file MomentumTransferPhaseSystem.H.

◆ VmTable

template<class BasePhaseSystem>
typedef HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > VmTable
protected

Definition at line 89 of file MomentumTransferPhaseSystem.H.

◆ VmfTable

template<class BasePhaseSystem>
typedef HashPtrTable< surfaceScalarField, phasePairKey, phasePairKey::hash > VmfTable
protected

Definition at line 96 of file MomentumTransferPhaseSystem.H.

◆ dragModelTable

template<class BasePhaseSystem>
typedef HashTable< autoPtr<BlendedInterfacialModel<dragModel> >, phasePairKey, phasePairKey::hash > dragModelTable
protected

Definition at line 103 of file MomentumTransferPhaseSystem.H.

◆ virtualMassModelTable

Definition at line 110 of file MomentumTransferPhaseSystem.H.

◆ liftModelTable

template<class BasePhaseSystem>
typedef HashTable< autoPtr<BlendedInterfacialModel<liftModel> >, phasePairKey, phasePairKey::hash > liftModelTable
protected

Definition at line 117 of file MomentumTransferPhaseSystem.H.

◆ wallLubricationModelTable

◆ turbulentDispersionModelTable

Constructor & Destructor Documentation

◆ MomentumTransferPhaseSystem()

template<class BasePhaseSystem>
MomentumTransferPhaseSystem ( const fvMesh & mesh)

Construct from fvMesh.

Definition at line 158 of file MomentumTransferPhaseSystem.C.

References forAllConstIter, IOobject::groupName(), mesh, and phasePair::name().

Here is the call graph for this function:

◆ ~MomentumTransferPhaseSystem()

template<class BasePhaseSystem>
~MomentumTransferPhaseSystem ( )
virtual

Destructor.

Definition at line 261 of file MomentumTransferPhaseSystem.C.

Member Function Documentation

◆ momentumTransfer()

template<class BasePhaseSystem>
Foam::autoPtr< Foam::phaseSystem::momentumTransferTable > momentumTransfer ( )
virtual

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 270 of file MomentumTransferPhaseSystem.C.

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

Here is the call graph for this function:

◆ momentumTransferf()

template<class BasePhaseSystem>
Foam::autoPtr< Foam::phaseSystem::momentumTransferTable > momentumTransferf ( )
virtual

◆ AFfs()

template<class BasePhaseSystem>
Foam::PtrList< Foam::surfaceScalarField > AFfs ( ) const
virtual

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

algorithm.

Definition at line 450 of file MomentumTransferPhaseSystem.C.

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

Referenced by AFfs().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ phiFs()

template<class BasePhaseSystem>
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 489 of file MomentumTransferPhaseSystem.C.

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

Referenced by phiFs().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ phiFfs()

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

As phiFs, but for the face-based algorithm.

Definition at line 630 of file MomentumTransferPhaseSystem.C.

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

Referenced by phiFfs().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ phiKdPhis()

template<class BasePhaseSystem>
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 793 of file MomentumTransferPhaseSystem.C.

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

Referenced by phiKdPhis().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ phiKdPhifs()

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

As phiKdPhis, but for the face-based algorithm.

Definition at line 835 of file MomentumTransferPhaseSystem.C.

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

Referenced by phiKdPhifs().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ KdUByAs()

template<class BasePhaseSystem>
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 877 of file MomentumTransferPhaseSystem.C.

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

Referenced by KdUByAs().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ partialElimination()

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

◆ partialEliminationf()

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

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

for the fluxes.

Definition at line 1143 of file MomentumTransferPhaseSystem.C.

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

Here is the call graph for this function:

◆ ddtCorrByAs()

template<class BasePhaseSystem>
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 913 of file MomentumTransferPhaseSystem.C.

References alpha, Foam::fvc::average(), Foam::byDt(), ddtCorrByAs(), Foam::fvc::flux(), forAll, forAllConstIter, phaseModel::index(), Foam::fvc::interpolate(), Foam::isA(), MRF, GeometricField< Type, PatchField, GeoMesh >::oldTime(), phaseModel::otherPhase(), phasei, phaseModel::phi(), Foam::pos0(), rAUs, tmp< T >::ref(), phase::rho(), PtrList< T >::set(), and phaseModel::U().

Referenced by ddtCorrByAs().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ DByAfs()

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

Return the phase diffusivities divided by the momentum coefficients.

Definition at line 1243 of file MomentumTransferPhaseSystem.C.

◆ read()

template<class BasePhaseSystem>
bool read ( )
virtual

Read base phaseProperties dictionary.

Definition at line 1250 of file MomentumTransferPhaseSystem.C.


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