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>
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 > ¤tScvfBoundaryTypes, 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 > ¤tScvfBoundaryTypes, 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 > ¤tScvfBoundaryTypes, 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 > ¤tScvfBoundaryTypes, 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... | |
|
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 *
|
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 *
|
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 *
|
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 *
|
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 *