The flux variables class for the Navier-Stokes model using the staggered grid discretization.
#include <dumux/freeflow/navierstokes/momentum/fluxvariables.hh>
|
| 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.
|
|
NumEqVector | diffusiveMomentumFlux () const |
| Returns the diffusive momentum flux due to viscous forces.
|
|
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.
|
|
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).
|
|
NumEqVector | pressureContribution () const |
| Returns the frontal pressure contribution.
|
|
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.
|
|
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).
|
|
◆ 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 |
◆ advectiveMomentumFlux()
◆ diffusiveMomentumFlux()
◆ elemBcTypes()
◆ element()
◆ elemFluxVarsCache()
◆ elemVolVars()
◆ frontalAdvectiveMomentumFlux()
* 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()
* 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()
◆ lateralAdvectiveMomentumFlux()
* ----------------
* | 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()
* ----------------
* | |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()
*
* ---------======= == 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()
◆ scvFace()
The documentation for this class was generated from the following file: