#include <dumux/freeflow/navierstokes/staggered/fluxvariables.hh>
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... | |
using Dumux::NavierStokesFluxVariablesImpl< TypeTag, DiscretizationMethod::staggered >::HeatConductionType = GetPropType<TypeTag, Properties::HeatConductionType> |
|
inlinestatic |
Returns the advective flux over a sub control volume face.
problem | The object specifying the problem which ought to be simulated |
elemVolVars | All volume variables for the element |
elemFaceVars | The face variables |
scvf | The sub control volume face |
upwindTerm | The uwind term (i.e. the advectively transported quantity) |
|
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 *
|
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 *
|
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 *
|
inline |
Returns the momentum flux over all staggered faces.
|
inlineinherited |
|
inlineinherited |
|
inlineinherited |
|
inlineinherited |
|
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 *
|
inlineinherited |
Initialize the flux variables storing some temporary pointers.
|
inlineinherited |
|
inlineinherited |