3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
Static Public Member Functions | List of all members
Dumux::StaggeredVelocityGradients< Scalar, GridGeometry, BoundaryTypes, Indices > Class Template Reference

Helper class for calculating the velocity gradients for the Navier-Stokes model using the staggered grid discretization. More...

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

Description

template<class Scalar, class GridGeometry, class BoundaryTypes, class Indices>
class Dumux::StaggeredVelocityGradients< Scalar, GridGeometry, BoundaryTypes, Indices >

Helper class for calculating the velocity gradients for the Navier-Stokes model using the staggered grid discretization.

Static Public Member Functions

template<class FaceVariables >
static Scalar velocityGradII (const SubControlVolumeFace &scvf, const FaceVariables &faceVars)
 Returns the in-axis velocity gradient. More...
 
template<class Problem , class FaceVariables >
static Scalar velocityGradIJ (const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf, const FaceVariables &faceVars, const std::optional< BoundaryTypes > &currentScvfBoundaryTypes, const std::optional< BoundaryTypes > &lateralFaceBoundaryTypes, const std::size_t localSubFaceIdx)
 Returns the velocity gradient perpendicular to the orientation of our current scvf. More...
 
template<class Problem , class FaceVariables >
static Scalar velocityGradJI (const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf, const FaceVariables &faceVars, const std::optional< BoundaryTypes > &currentScvfBoundaryTypes, const std::optional< BoundaryTypes > &lateralFaceBoundaryTypes, const std::size_t localSubFaceIdx)
 Returns the velocity gradient in line with our current scvf. More...
 
template<class Problem , class FaceVariables >
static Scalar beaversJosephVelocityAtCurrentScvf (const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf, const FaceVariables &faceVars, const std::optional< BoundaryTypes > &currentScvfBoundaryTypes, const std::optional< BoundaryTypes > &lateralFaceBoundaryTypes, const std::size_t localSubFaceIdx)
 Returns the Beavers-Jospeh slip velocity for a scvf which lies on the boundary itself. More...
 
template<class Problem , class FaceVariables >
static Scalar beaversJosephVelocityAtLateralScvf (const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf, const FaceVariables &faceVars, const std::optional< BoundaryTypes > &currentScvfBoundaryTypes, const std::optional< BoundaryTypes > &lateralFaceBoundaryTypes, const std::size_t localSubFaceIdx)
 Returns the Beavers-Jospeh slip velocity for a lateral scvf which lies on the boundary. More...
 

Member Function Documentation

◆ beaversJosephVelocityAtCurrentScvf()

template<class Scalar , class GridGeometry , class BoundaryTypes , class Indices >
template<class Problem , class FaceVariables >
static Scalar Dumux::StaggeredVelocityGradients< Scalar, GridGeometry, BoundaryTypes, Indices >::beaversJosephVelocityAtCurrentScvf ( const Problem &  problem,
const Element &  element,
const FVElementGeometry &  fvGeometry,
const SubControlVolumeFace &  scvf,
const FaceVariables &  faceVars,
const std::optional< BoundaryTypes > &  currentScvfBoundaryTypes,
const std::optional< BoundaryTypes > &  lateralFaceBoundaryTypes,
const std::size_t  localSubFaceIdx 
)
inlinestatic

Returns the Beavers-Jospeh slip velocity for a scvf which lies on the boundary itself.

*                  in.norm.  B-J slip
*                     vel.   vel.
*                     ^       ^
*                     |       |
*       scvf   ---------######|*               * boundary
*              |      ||      |* curr.
*              |      ||      |* scvf          || and # staggered half-control-volume (own element)
*              |      ||      x~~~~>
*              |      ||      |* vel.          # lateral staggered faces
*        scvf  |      ||      |* Self
*              ---------#######*                x dof position
*                 scvf
*                                              -- element
* 

◆ beaversJosephVelocityAtLateralScvf()

template<class Scalar , class GridGeometry , class BoundaryTypes , class Indices >
template<class Problem , class FaceVariables >
static Scalar Dumux::StaggeredVelocityGradients< Scalar, GridGeometry, BoundaryTypes, Indices >::beaversJosephVelocityAtLateralScvf ( const Problem &  problem,
const Element &  element,
const FVElementGeometry &  fvGeometry,
const SubControlVolumeFace &  scvf,
const FaceVariables &  faceVars,
const std::optional< BoundaryTypes > &  currentScvfBoundaryTypes,
const std::optional< BoundaryTypes > &  lateralFaceBoundaryTypes,
const std::size_t  localSubFaceIdx 
)
inlinestatic

Returns the Beavers-Jospeh slip velocity for a lateral scvf which lies on the boundary.

*                             B-J slip                  * boundary
*              ************** vel. *****
*       scvf   ---------##### ~~~~> ::::                || and # staggered half-control-volume (own element)
*              |      ||      | curr. ::
*              |      ||      | scvf  ::                :: staggered half-control-volume (neighbor element)
*              |      ||      x~~~~>  ::
*              |      ||      | vel.  ::                # lateral staggered faces
*        scvf  |      ||      | Self  ::
*              ---------#######:::::::::                x dof position
*                 scvf
*                                                       -- elements
* 

◆ velocityGradII()

template<class Scalar , class GridGeometry , class BoundaryTypes , class Indices >
template<class FaceVariables >
static Scalar Dumux::StaggeredVelocityGradients< Scalar, GridGeometry, BoundaryTypes, Indices >::velocityGradII ( const SubControlVolumeFace &  scvf,
const FaceVariables &  faceVars 
)
inlinestatic

Returns the in-axis velocity gradient.

*              ---------=======                 == and # staggered half-control-volume
*              |       #      | current scvf
*              |       #      |                 # staggered face over which fluxes are calculated
*   vel.Opp <~~|       O~~>   x~~~~> vel.Self
*              |       #      |                 x dof position
*        scvf  |       #      |
*              --------========                 -- element
*
*                                               O position at which gradient is evaluated
* 

◆ velocityGradIJ()

template<class Scalar , class GridGeometry , class BoundaryTypes , class Indices >
template<class Problem , class FaceVariables >
static Scalar Dumux::StaggeredVelocityGradients< Scalar, GridGeometry, BoundaryTypes, Indices >::velocityGradIJ ( const Problem &  problem,
const Element &  element,
const FVElementGeometry &  fvGeometry,
const SubControlVolumeFace &  scvf,
const FaceVariables &  faceVars,
const std::optional< BoundaryTypes > &  currentScvfBoundaryTypes,
const std::optional< BoundaryTypes > &  lateralFaceBoundaryTypes,
const std::size_t  localSubFaceIdx 
)
inlinestatic

Returns the velocity gradient perpendicular to the orientation of our current scvf.

*              ----------------
*              |              |vel.
*              |              |Parallel
*              |              |~~~~>       ------->
*              |              |             ------>     * gradient
*              |              |              ----->
*       scvf   ---------######O:::::::::      ---->     || 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
*
*                                                       O position at which gradient is evaluated
* 

| xxxx s | xxxx a | xxxx s --------—O--------— | yyyy s zzzz | | yyyy b zzzz |

| yyyy s zzzz |

In a corner geometry (scvf is sas or sbs), we calculate the velocity gradient at O, by (velocity(a)-velocity(b))/distance(a,b) for the half-control volumes x and y, but by (velocity(O)-velocity(b))/distance(O,b) for z. This does not harm flux continuity (x and y use the same formulation). We do this different formulation for y (approximate gradient by central differncing) and z (approximate gradient by forward/backward differencing), because it is the natural way of implementing it and it is not clear which gradient is the better approximation in this case anyway. Particularly, for the frequent case of no-slip, no-flow boundaries, the velocity would be zero at O and a and thus, the gradient within the flow domain might be better approximated by velocity(b)/distanc(O,b) than by velocity(b)/distance(a,b).

◆ velocityGradJI()

template<class Scalar , class GridGeometry , class BoundaryTypes , class Indices >
template<class Problem , class FaceVariables >
static Scalar Dumux::StaggeredVelocityGradients< Scalar, GridGeometry, BoundaryTypes, Indices >::velocityGradJI ( const Problem &  problem,
const Element &  element,
const FVElementGeometry &  fvGeometry,
const SubControlVolumeFace &  scvf,
const FaceVariables &  faceVars,
const std::optional< BoundaryTypes > &  currentScvfBoundaryTypes,
const std::optional< BoundaryTypes > &  lateralFaceBoundaryTypes,
const std::size_t  localSubFaceIdx 
)
inlinestatic

Returns the velocity gradient in line with our current scvf.

*                      ^       gradient
*                      |  ^
*                      |  |  ^
*                      |  |  |  ^
*                      |  |  |  |  ^
*                      |  |  |  |  |  ^
*                      |  |  |  |  |  |
*
*              ----------------
*              |              |
*              |    in.norm.  |
*              |       vel.   |
*              |       ^      |        ^ out.norm.vel.
*              |       |      |        |
*       scvf   ---------######O:::::::::       || 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
*
*                                              O position at which gradient is evaluated
* 

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