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

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

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

Inheritance diagram for Dumux::NavierStokesFluxVariablesImpl< TypeTag, DiscretizationMethod::staggered >:
Inheritance graph

Description

template<class TypeTag>
class Dumux::NavierStokesFluxVariablesImpl< TypeTag, DiscretizationMethod::staggered >

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

Public Types

using HeatConductionType = GetPropType< TypeTag, Properties::HeatConductionType >
 

Public Member Functions

CellCenterPrimaryVariables computeMassFlux (const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolumeFace &scvf, const FluxVariablesCache &fluxVarsCache)
 Computes the flux for the cell center residual (mass balance). More...
 
FacePrimaryVariables computeMomentumFlux (const Problem &problem, const Element &element, const SubControlVolumeFace &scvf, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const GridFluxVariablesCache &gridFluxVarsCache)
 Returns the momentum flux over all staggered faces. More...
 
FacePrimaryVariables computeFrontalMomentumFlux (const Problem &problem, const Element &element, const SubControlVolumeFace &scvf, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const GridFluxVariablesCache &gridFluxVarsCache)
 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...
 
FacePrimaryVariables computeLateralMomentumFlux (const Problem &problem, const Element &element, const SubControlVolumeFace &scvf, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const GridFluxVariablesCache &gridFluxVarsCache)
 Returns the momentum flux over the staggered faces perpendicular to the scvf where the velocity dof of interest lives (coinciding with the element's scvfs). More...
 
FacePrimaryVariables inflowOutflowBoundaryFlux (const Problem &problem, const Element &element, const SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars) const
 Returns the momentum flux over an inflow or outflow boundary face. More...
 
void init (const GetPropType< TypeTag, Properties::Problem > &problem, const Element &element, const GetPropType< TypeTag, Properties::GridGeometry >::LocalView &fvGeometry, const GetPropType< TypeTag, Properties::GridVolumeVariables >::LocalView &elemVolVars, const SubControlVolumeFace &scvFace, const GetPropType< TypeTag, Properties::GridFluxVariablesCache >::LocalView &elemFluxVarsCache)
 Initialize the flux variables storing some temporary pointers. More...
 
const GetPropType< TypeTag, Properties::Problem > & problem () const
 
const Element & element () const
 
const SubControlVolumeFace & scvFace () const
 
const GetPropType< TypeTag, Properties::GridGeometry >::LocalView & fvGeometry () const
 
const GetPropType< TypeTag, Properties::GridVolumeVariables >::LocalView & elemVolVars () const
 
const GetPropType< TypeTag, Properties::GridFluxVariablesCache >::LocalView & elemFluxVarsCache () const
 

Static Public Member Functions

template<class UpwindTerm >
static Scalar advectiveFluxForCellCenter (const Problem &problem, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolumeFace &scvf, UpwindTerm upwindTerm)
 Returns the advective flux over a sub control volume face. More...
 

Member Typedef Documentation

◆ HeatConductionType

template<class TypeTag >
using Dumux::NavierStokesFluxVariablesImpl< TypeTag, DiscretizationMethod::staggered >::HeatConductionType = GetPropType<TypeTag, Properties::HeatConductionType>

Member Function Documentation

◆ advectiveFluxForCellCenter()

template<class TypeTag >
template<class UpwindTerm >
static Scalar Dumux::NavierStokesFluxVariablesImpl< TypeTag, DiscretizationMethod::staggered >::advectiveFluxForCellCenter ( const Problem &  problem,
const ElementVolumeVariables &  elemVolVars,
const ElementFaceVariables &  elemFaceVars,
const SubControlVolumeFace &  scvf,
UpwindTerm  upwindTerm 
)
inlinestatic

Returns the advective flux over a sub control volume face.

Parameters
problemThe object specifying the problem which ought to be simulated
elemVolVarsAll volume variables for the element
elemFaceVarsThe face variables
scvfThe sub control volume face
upwindTermThe uwind term (i.e. the advectively transported quantity)

◆ computeFrontalMomentumFlux()

template<class TypeTag >
FacePrimaryVariables Dumux::NavierStokesFluxVariablesImpl< TypeTag, DiscretizationMethod::staggered >::computeFrontalMomentumFlux ( const Problem &  problem,
const Element &  element,
const SubControlVolumeFace &  scvf,
const FVElementGeometry &  fvGeometry,
const ElementVolumeVariables &  elemVolVars,
const ElementFaceVariables &  elemFaceVars,
const GridFluxVariablesCache &  gridFluxVarsCache 
)
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 scvf
*              |       #      |                 # staggered face over which fluxes are calculated
*   vel.Opp <~~|       #~~>   x~~~~> vel.Self
*              |       #      |                 x dof position
*        scvf  |       #      |
*              --------========                 -- element
*                   scvf
* 

◆ computeLateralMomentumFlux()

template<class TypeTag >
FacePrimaryVariables Dumux::NavierStokesFluxVariablesImpl< TypeTag, DiscretizationMethod::staggered >::computeLateralMomentumFlux ( const Problem &  problem,
const Element &  element,
const SubControlVolumeFace &  scvf,
const FVElementGeometry &  fvGeometry,
const ElementVolumeVariables &  elemVolVars,
const ElementFaceVariables &  elemFaceVars,
const GridFluxVariablesCache &  gridFluxVarsCache 
)
inline

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

*                scvf
*              ---------#######                 || and # staggered half-control-volume
*              |      ||      | current scvf
*              |      ||      |                 # normal staggered sub faces over which fluxes are calculated
*              |      ||      x~~~~> vel.Self
*              |      ||      |                 x dof position
*        scvf  |      ||      |
*              --------########                -- element
*                 scvf
* 

◆ computeMassFlux()

template<class TypeTag >
CellCenterPrimaryVariables Dumux::NavierStokesFluxVariablesImpl< TypeTag, DiscretizationMethod::staggered >::computeMassFlux ( const Problem &  problem,
const Element &  element,
const FVElementGeometry &  fvGeometry,
const ElementVolumeVariables &  elemVolVars,
const ElementFaceVariables &  elemFaceVars,
const SubControlVolumeFace &  scvf,
const FluxVariablesCache &  fluxVarsCache 
)
inline

Computes the flux for the cell center residual (mass balance).

*                    scvf
*              ----------------
*              |              # current scvf    # scvf over which fluxes are calculated
*              |              #
*              |      x       #~~~~> vel.Self   x dof position
*              |              #
*        scvf  |              #                 -- element
*              ----------------
*                   scvf
* 

◆ computeMomentumFlux()

template<class TypeTag >
FacePrimaryVariables Dumux::NavierStokesFluxVariablesImpl< TypeTag, DiscretizationMethod::staggered >::computeMomentumFlux ( const Problem &  problem,
const Element &  element,
const SubControlVolumeFace &  scvf,
const FVElementGeometry &  fvGeometry,
const ElementVolumeVariables &  elemVolVars,
const ElementFaceVariables &  elemFaceVars,
const GridFluxVariablesCache &  gridFluxVarsCache 
)
inline

Returns the momentum flux over all staggered faces.

◆ element()

const Element & Dumux::FluxVariablesBase< GetPropType< TypeTag, Properties::Problem > , GetPropType< TypeTag, Properties::GridGeometry >::LocalView , GetPropType< TypeTag, Properties::GridVolumeVariables >::LocalView , GetPropType< TypeTag, Properties::GridFluxVariablesCache >::LocalView >::element ( ) const
inlineinherited

◆ elemFluxVarsCache()

const GetPropType< TypeTag, Properties::GridFluxVariablesCache >::LocalView & Dumux::FluxVariablesBase< GetPropType< TypeTag, Properties::Problem > , GetPropType< TypeTag, Properties::GridGeometry >::LocalView , GetPropType< TypeTag, Properties::GridVolumeVariables >::LocalView , GetPropType< TypeTag, Properties::GridFluxVariablesCache >::LocalView >::elemFluxVarsCache ( ) const
inlineinherited

◆ elemVolVars()

const GetPropType< TypeTag, Properties::GridVolumeVariables >::LocalView & Dumux::FluxVariablesBase< GetPropType< TypeTag, Properties::Problem > , GetPropType< TypeTag, Properties::GridGeometry >::LocalView , GetPropType< TypeTag, Properties::GridVolumeVariables >::LocalView , GetPropType< TypeTag, Properties::GridFluxVariablesCache >::LocalView >::elemVolVars ( ) const
inlineinherited

◆ fvGeometry()

const GetPropType< TypeTag, Properties::GridGeometry >::LocalView & Dumux::FluxVariablesBase< GetPropType< TypeTag, Properties::Problem > , GetPropType< TypeTag, Properties::GridGeometry >::LocalView , GetPropType< TypeTag, Properties::GridVolumeVariables >::LocalView , GetPropType< TypeTag, Properties::GridFluxVariablesCache >::LocalView >::fvGeometry ( ) const
inlineinherited

◆ inflowOutflowBoundaryFlux()

template<class TypeTag >
FacePrimaryVariables Dumux::NavierStokesFluxVariablesImpl< TypeTag, DiscretizationMethod::staggered >::inflowOutflowBoundaryFlux ( const Problem &  problem,
const Element &  element,
const SubControlVolumeFace &  scvf,
const ElementVolumeVariables &  elemVolVars,
const ElementFaceVariables &  elemFaceVars 
) const
inline

Returns the momentum flux over an inflow or outflow boundary face.

*                    scvf      //
*              ---------=======//               == and # staggered half-control-volume
*              |      ||      #// current scvf
*              |      ||      #//               # staggered boundary face over which fluxes are calculated
*              |      ||      x~~~~> vel.Self
*              |      ||      #//               x dof position
*        scvf  |      ||      #//
*              --------========//               -- element
*                   scvf       //
*                                              // boundary
* 

◆ init()

void Dumux::FluxVariablesBase< GetPropType< TypeTag, Properties::Problem > , GetPropType< TypeTag, Properties::GridGeometry >::LocalView , GetPropType< TypeTag, Properties::GridVolumeVariables >::LocalView , GetPropType< TypeTag, Properties::GridFluxVariablesCache >::LocalView >::init ( const GetPropType< TypeTag, Properties::Problem > &  problem,
const Element &  element,
const GetPropType< TypeTag, Properties::GridGeometry >::LocalView &  fvGeometry,
const GetPropType< TypeTag, Properties::GridVolumeVariables >::LocalView &  elemVolVars,
const SubControlVolumeFace &  scvFace,
const GetPropType< TypeTag, Properties::GridFluxVariablesCache >::LocalView &  elemFluxVarsCache 
)
inlineinherited

Initialize the flux variables storing some temporary pointers.

◆ problem()

const GetPropType< TypeTag, Properties::Problem > & Dumux::FluxVariablesBase< GetPropType< TypeTag, Properties::Problem > , GetPropType< TypeTag, Properties::GridGeometry >::LocalView , GetPropType< TypeTag, Properties::GridVolumeVariables >::LocalView , GetPropType< TypeTag, Properties::GridFluxVariablesCache >::LocalView >::problem ( ) const
inlineinherited

◆ scvFace()

const SubControlVolumeFace & Dumux::FluxVariablesBase< GetPropType< TypeTag, Properties::Problem > , GetPropType< TypeTag, Properties::GridGeometry >::LocalView , GetPropType< TypeTag, Properties::GridVolumeVariables >::LocalView , GetPropType< TypeTag, Properties::GridFluxVariablesCache >::LocalView >::scvFace ( ) const
inlineinherited

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