3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
Public Member Functions | List of all members
Dumux::NavierStokesMomentumFluxVariables< TypeTag > Class Template Reference

The flux variables class for the Navier-Stokes model using the staggered grid discretization. More...

#include <dumux/freeflow/navierstokes/momentum/fluxvariables.hh>

Description

template<class TypeTag>
class Dumux::NavierStokesMomentumFluxVariables< TypeTag >

The flux variables class for the Navier-Stokes model using the staggered grid discretization.

Public Member Functions

 NavierStokesMomentumFluxVariables (const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvFace, const ElementVolumeVariables &elemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const ElementBoundaryTypes &elemBcTypes)
 
const Problem & problem () const
 
const Element & element () const
 
const SubControlVolumeFace & scvFace () const
 
const FVElementGeometry & fvGeometry () const
 
const ElementVolumeVariables & elemVolVars () const
 
const ElementFluxVariablesCache & elemFluxVarsCache () const
 
const ElementBoundaryTypes & elemBcTypes () const
 
NumEqVector advectiveMomentumFlux () const
 Returns the diffusive momentum flux due to viscous forces. More...
 
NumEqVector diffusiveMomentumFlux () const
 Returns the diffusive momentum flux due to viscous forces. More...
 
NumEqVector frontalDiffusiveMomentumFlux () const
 Returns the frontal part of the momentum flux. This treats the flux over the staggered face at the center of an element, parallel to the current scvf where the velocity dof of interest lives. More...
 
NumEqVector lateralDiffusiveMomentumFlux () const
 Returns the diffusive momentum flux over the staggered face perpendicular to the scvf where the velocity dof of interest lives (coinciding with the element's scvfs). More...
 
NumEqVector pressureContribution () const
 Returns the frontal pressure contribution. More...
 
NumEqVector frontalAdvectiveMomentumFlux () const
 Returns the frontal part of the momentum flux. This treats the flux over the staggered face at the center of an element, parallel to the current scvf where the velocity dof of interest lives. More...
 
NumEqVector lateralAdvectiveMomentumFlux () const
 Returns the advective momentum flux over the staggered face perpendicular to the scvf where the velocity dof of interest lives (coinciding with the element's scvfs). More...
 

Constructor & Destructor Documentation

◆ NavierStokesMomentumFluxVariables()

template<class TypeTag >
Dumux::NavierStokesMomentumFluxVariables< TypeTag >::NavierStokesMomentumFluxVariables ( const Problem &  problem,
const Element &  element,
const FVElementGeometry &  fvGeometry,
const SubControlVolumeFace &  scvFace,
const ElementVolumeVariables &  elemVolVars,
const ElementFluxVariablesCache &  elemFluxVarsCache,
const ElementBoundaryTypes &  elemBcTypes 
)
inline

Member Function Documentation

◆ advectiveMomentumFlux()

template<class TypeTag >
NumEqVector Dumux::NavierStokesMomentumFluxVariables< TypeTag >::advectiveMomentumFlux ( ) const
inline

Returns the diffusive momentum flux due to viscous forces.

◆ diffusiveMomentumFlux()

template<class TypeTag >
NumEqVector Dumux::NavierStokesMomentumFluxVariables< TypeTag >::diffusiveMomentumFlux ( ) const
inline

Returns the diffusive momentum flux due to viscous forces.

◆ elemBcTypes()

template<class TypeTag >
const ElementBoundaryTypes & Dumux::NavierStokesMomentumFluxVariables< TypeTag >::elemBcTypes ( ) const
inline

◆ element()

template<class TypeTag >
const Element & Dumux::NavierStokesMomentumFluxVariables< TypeTag >::element ( ) const
inline

◆ elemFluxVarsCache()

template<class TypeTag >
const ElementFluxVariablesCache & Dumux::NavierStokesMomentumFluxVariables< TypeTag >::elemFluxVarsCache ( ) const
inline

◆ elemVolVars()

template<class TypeTag >
const ElementVolumeVariables & Dumux::NavierStokesMomentumFluxVariables< TypeTag >::elemVolVars ( ) const
inline

◆ frontalAdvectiveMomentumFlux()

template<class TypeTag >
NumEqVector Dumux::NavierStokesMomentumFluxVariables< TypeTag >::frontalAdvectiveMomentumFlux ( ) const
inline

Returns the frontal part of the momentum flux. This treats the flux over the staggered face at the center of an element, parallel to the current scvf where the velocity dof of interest lives.

*                    scvf
*              ---------=======                 == and # staggered half-control-volume
*              |       #      | current scv
*              |       #      |                 # staggered face over which fluxes are calculated
*      vel.Opp |~~>    O~~~>  x~~~~> vel.Self
*              |       #      |                 x dof position
*              |       #      |
*              --------========                 -- element
*                   scvf
*                                               O integration point
* 

◆ frontalDiffusiveMomentumFlux()

template<class TypeTag >
NumEqVector Dumux::NavierStokesMomentumFluxVariables< TypeTag >::frontalDiffusiveMomentumFlux ( ) const
inline

Returns the frontal part of the momentum flux. This treats the flux over the staggered face at the center of an element, parallel to the current scvf where the velocity dof of interest lives.

*                    scvf
*              ---------=======                 == and # staggered half-control-volume
*              |       #      | current scv
*              |       #      |                 # staggered face over which fluxes are calculated
*      vel.Opp |~~>    O~~~>  x~~~~> vel.Self
*              |       #      |                 x dof position
*              |       #      |
*              --------========                 -- element
*                   scvf
*                                               O integration point
* 

◆ fvGeometry()

template<class TypeTag >
const FVElementGeometry & Dumux::NavierStokesMomentumFluxVariables< TypeTag >::fvGeometry ( ) const
inline

◆ lateralAdvectiveMomentumFlux()

template<class TypeTag >
NumEqVector Dumux::NavierStokesMomentumFluxVariables< TypeTag >::lateralAdvectiveMomentumFlux ( ) const
inline

Returns the advective momentum flux over the staggered face perpendicular to the scvf where the velocity dof of interest lives (coinciding with the element's scvfs).

*              ----------------
*              |     inner    |
*              |    transp.   |
*              |      vel.    |~~~~> outer vel.
*              |       ^      |
*              |       |      |
*              ---------######O                 || and # staggered half-control-volume
*              |      ||      | scv
*              |      ||      |                 # lateral staggered faces over which fluxes are calculated
*              |      ||      x~~~~> inner vel.
*              |      ||      |                 x dof position
*              |      ||      |
*              ---------#######                -- elements
*
*                                               O integration point
* 

◆ lateralDiffusiveMomentumFlux()

template<class TypeTag >
NumEqVector Dumux::NavierStokesMomentumFluxVariables< TypeTag >::lateralDiffusiveMomentumFlux ( ) const
inline

Returns the diffusive momentum flux over the staggered face perpendicular to the scvf where the velocity dof of interest lives (coinciding with the element's scvfs).

*              ----------------
*              |              |vel.
*              |    in.norm.  |Parallel
*              |       vel.   |~~~~>
*              |       ^      |        ^ out.norm.vel.
*              |       |      |        |
*       scvf   ---------#######:::::::::       || and # staggered half-control-volume (own element)
*              |      ||      | curr. ::
*              |      ||      | scvf  ::       :: staggered half-control-volume (neighbor element)
*              |      ||      x~~~~>  ::
*              |      ||      | vel.  ::       # lateral staggered faces over which fluxes are calculated
*        scvf  |      ||      | Self  ::
*              ---------#######:::::::::       x dof position
*                 scvf
*                                              -- elements
* 

◆ pressureContribution()

template<class TypeTag >
NumEqVector Dumux::NavierStokesMomentumFluxVariables< TypeTag >::pressureContribution ( ) const
inline

Returns the frontal pressure contribution.

*
*              ---------=======                 == and # staggered half-control-volume
*              |       #      | current scv
*              |       #      |                 # frontal staggered face for which pressure contribution is calculated
*              |       P      x
*              |       #      |                 x dof position
*              |       #      |
*              --------========                 -- element
*
*                                               P integration point, pressure DOF lives here
* 

◆ problem()

template<class TypeTag >
const Problem & Dumux::NavierStokesMomentumFluxVariables< TypeTag >::problem ( ) const
inline

◆ scvFace()

template<class TypeTag >
const SubControlVolumeFace & Dumux::NavierStokesMomentumFluxVariables< TypeTag >::scvFace ( ) const
inline

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