LifeV
StabilizationSD< MeshType, DofType > Class Template Reference

StabilizationSD Class. More...

#include <StabilizationSD.hpp>

+ Collaboration diagram for StabilizationSD< MeshType, DofType >:
typedef MeshType mesh_Type
 
typedef DofType dof_Type
 

Constructor and Destructor

 StabilizationSD (const GetPot &dataFile, const mesh_Type &mesh, const dof_Type &dof, const ReferenceFE &refFE, const QuadratureRule &quadRule, Real viscosity)
 Constructor. More...
 
virtual ~StabilizationSD ()
 ~Destructor More...
 

Methods

template<typename MatrixType , typename VectorType >
void applySUPG (const Real dt, MatrixType &matrix, const VectorType &state)
 compute the SUPG stabilization terms and add them into the monolithic N.S. matrix More...
 
template<typename MatrixType , typename VectorType >
void applySD (const Real dt, MatrixType &matrix, const VectorType &state)
 compute the SD stabilization terms and add them into the Momentum matrix More...
 
template<typename VectorType , typename SourceType >
void applyRHS (const Real dt, VectorType &vector, const VectorType &state, const SourceType &source, const Real &time)
 compute the SUPG stabilization terms and add them into the right and side More...
 
void showMe (std::ostream &output=std::cout) const
 Display class informations. More...
 

Set Methods

void setGammaBeta (const Real &gammaBeta)
 Set the stabilization parameter $\gamma_\beta$ for $( \beta \nabla \mathbf{u} , \beta \nabla \mathbf{v})$. More...
 
void setGammaDiv (const Real &gammaDiv)
 Set the stabilization parameter $\gamma_d$ for $( div \mathbf{u} , div \nabla \mathbf{v})$. More...
 

Private Constructor

 StabilizationSD ()
 Default Constructor. More...
 
 StabilizationSD (const StabilizationSD< mesh_Type, dof_Type > &original)
 Copy Constructor. More...
 

Private Methods

template<typename VectorType >
void computeParameters (const Real dt, const UInt iVol, const VectorType &state, VectorElemental &beta, Real &coeffBeta, Real &coeffDiv) const
 Compute the stabilization coefficients for each element. More...
 
void bgradu_bgradv (const Real &coef, VectorElemental &vel, ElemMat &elmat, const CurrentFE &fe, UInt iblock, UInt jblock, UInt nb) const
 Evaluate the varf $(\beta \nabla \mathbf{u}, \beta \nabla \mathbf{v})$. More...
 
void lapu_bgradv (const Real &coef, VectorElemental &vel, MatrixElemental &elmat, const CurrentFE &fe, UInt iblock, UInt jblock, UInt nb) const
 Evaluate the varf $(\Delta \mathbf{u}, \beta \nabla \mathbf{v})$. More...
 
void gradp_bgradv (const Real &coef, VectorElemental &vel, MatrixElemental &elmat, const CurrentFE &fe) const
 Evaluate the varf $(\nabla p, \beta \nabla \mathbf{v})$. More...
 
void lapu_gradq (const Real &coef, MatrixElemental &elmat, const CurrentFE &fe) const
 Evaluate the varf $(\Delta \mathbf{u}, \nabla q)$. More...
 
template<typename SourceType >
void f_bgradv (const Real &coef, SourceType &source, VectorElemental &vel, VectorElemental &elvec, const CurrentFE &fe, UInt iblock, const Real &time) const
 Evaluate the varf $(f, \beta \nabla \mathbf{v})$. More...
 
template<typename SourceType >
void f_gradq (const Real &coef, SourceType &source, VectorElemental &elvec, const CurrentFE &fe, UInt iblock, const Real &time) const
 Evaluate the varf $(f, \nabla q)$. More...
 

Private Attributes

const mesh_TypeM_mesh
 the mesh object More...
 
const dof_TypeM_dof
 the dof object More...
 
CurrentFE M_fe
 current fe for the assembling of stabilization terms. More...
 
Real M_viscosity
 fluid viscosity $\nu$ More...
 
Real M_gammaBeta
 Stabilization coefficient of $(c(h,dt, |\beta|, nu) \beta \nabla \mathbf{u} , \beta \nabla \mathbf{v})$. More...
 
Real M_gammaDiv
 Stabilization coefficient of $(c(h,dt) div \mathbf{u} , div \nabla \mathbf{v})$. More...
 
MatrixElemental M_elMat
 Elementary Matrix for assembling the stabilization terms. More...
 
VectorElemental M_elVec
 Elementary Vector for assembling the stabilization terms. More...
 

Detailed Description

template<typename MeshType, typename DofType>
class LifeV::StabilizationSD< MeshType, DofType >

StabilizationSD Class.

Streamline diffusion and SUPG for Navier-Stokes

Author
M.A. Fernandez

Implementation of streamline diffusion (SD) and SUPG for the Navier-Stokes equations.
SUPG can be used only with equal order finite element for the discretization of velocity and pressure fields.

Definition at line 59 of file StabilizationSD.hpp.

Member Typedef Documentation

◆ mesh_Type

Definition at line 65 of file StabilizationSD.hpp.

◆ dof_Type

typedef DofType dof_Type

Definition at line 66 of file StabilizationSD.hpp.

Constructor & Destructor Documentation

◆ StabilizationSD() [1/3]

StabilizationSD ( const GetPot dataFile,
const mesh_Type mesh,
const dof_Type dof,
const ReferenceFE refFE,
const QuadratureRule quadRule,
Real  viscosity 
)

Constructor.

Parameters
dataFileGetPot dataFile where to read the stabilization parameter
meshmesh_Type mesh
dofdof_Type velocity field degree of freedom
refFErefFE velocity field reference finite element
quadRuleQuadratureRule quadrature rule for the integration of the stabilization variational forms
viscosityviscosity fluid viscosity $\nu$

Definition at line 236 of file StabilizationSD.hpp.

+ Here is the caller graph for this function:

◆ ~StabilizationSD()

virtual ~StabilizationSD ( )
inlinevirtual

~Destructor

Definition at line 87 of file StabilizationSD.hpp.

◆ StabilizationSD() [2/3]

StabilizationSD ( )
private

Default Constructor.

◆ StabilizationSD() [3/3]

StabilizationSD ( const StabilizationSD< mesh_Type, dof_Type > &  original)
private

Copy Constructor.

Member Function Documentation

◆ applySUPG()

void applySUPG ( const Real  dt,
MatrixType &  matrix,
const VectorType &  state 
)

compute the SUPG stabilization terms and add them into the monolithic N.S. matrix

This function adds the following stabilization terms to the Navier-Stokes monolithic matrix:

  1. Block(1,1): $c_\beta(\beta \nabla \mathbf{u} , \beta \nabla \mathbf{v}) + c_\beta( - \mu L \mathbf{u} , \beta \nabla \mathbf{v} )+ c_d( div \mathbf{u} , div \nabla \mathbf{v})$
  2. Block(1,2): $c_\beta(\nabla p , \beta \nabla \mathbf{v})$
  3. Block(2,1): $c_\beta(\beta \nabla \mathbf{u} , \nabla q) + c_\beta( - \mu L \mathbf{u} , \nabla q )$
  4. Block(2,2): $c_\beta(\nabla p , \nabla q)$

where $(\cdot, \cdot)$ represents the $L^2$ scalar product, and

  1. $\displaystyle c_\beta = \frac{\gamma_\beta}{\sqrt{ 4/( dt^2) + 4\|\beta\|_\infty/h^2 + 16*\nu/h^4 }} $, where $h$ is the diameter of the element;
  2. $\displaystyle c_d = \gamma_d h \|\beta\|_\infty $.

Both high Pechlet numbers and inf-sup incompatible FEM are stabilized.

PREREQUISITE: The velocity and the pressure field should belong to the same finite element space

Parameters are the followings:

Parameters
dtReal timestep (INPUT)
matrixMatrixType where the stabilization terms are added into. (OUTPUT)
stateVectorType velocity field for the linearization of the stabilization (INPUT)

Definition at line 257 of file StabilizationSD.hpp.

◆ applySD()

void applySD ( const Real  dt,
MatrixType &  matrix,
const VectorType &  state 
)

compute the SD stabilization terms and add them into the Momentum matrix

The following stabilization term is added: $c_\beta(\beta \nabla \mathbf{u} , \beta \nabla \mathbf{v})$

Parameters
dtReal timestep
matrixMatrixType where the stabilization terms are added into
stateVectorType velocity field for the linearization of the stabilization

Definition at line 344 of file StabilizationSD.hpp.

◆ applyRHS()

void applyRHS ( const Real  dt,
VectorType &  vector,
const VectorType &  state,
const SourceType &  source,
const Real time 
)

compute the SUPG stabilization terms and add them into the right and side

Add the following stabilization terms to rhs

  1. Block(1,1): $(c_\beta f , \beta \nabla \mathbf{v})$;
  2. Block(1,2): $(c_\beta f , \nabla q )$.
Parameters
dtReal timestep
vectorVectorType where the stabilization terms are added into
stateVectorType velocity field for the linearization of the stabilization
sourceSourceType a functor f(time, x, y, z, ic) which represents the forcing term
timeREAL the actual time in which source should be evaluated

Definition at line 417 of file StabilizationSD.hpp.

◆ showMe()

void showMe ( std::ostream &  output = std::cout) const

Display class informations.

Definition at line 487 of file StabilizationSD.hpp.

◆ setGammaBeta()

void setGammaBeta ( const Real gammaBeta)
inline

Set the stabilization parameter $\gamma_\beta$ for $( \beta \nabla \mathbf{u} , \beta \nabla \mathbf{v})$.

Definition at line 156 of file StabilizationSD.hpp.

◆ setGammaDiv()

void setGammaDiv ( const Real gammaDiv)
inline

Set the stabilization parameter $\gamma_d$ for $( div \mathbf{u} , div \nabla \mathbf{v})$.

Definition at line 161 of file StabilizationSD.hpp.

◆ computeParameters()

void computeParameters ( const Real  dt,
const UInt  iVol,
const VectorType &  state,
VectorElemental beta,
Real coeffBeta,
Real coeffDiv 
) const
private

Compute the stabilization coefficients for each element.

Definition at line 502 of file StabilizationSD.hpp.

◆ bgradu_bgradv()

void bgradu_bgradv ( const Real coef,
VectorElemental vel,
ElemMat &  elmat,
const CurrentFE fe,
UInt  iblock,
UInt  jblock,
UInt  nb 
) const
private

Evaluate the varf $(\beta \nabla \mathbf{u}, \beta \nabla \mathbf{v})$.

Definition at line 592 of file StabilizationSD.hpp.

◆ lapu_bgradv()

void lapu_bgradv ( const Real coef,
VectorElemental vel,
MatrixElemental elmat,
const CurrentFE fe,
UInt  iblock,
UInt  jblock,
UInt  nb 
) const
private

Evaluate the varf $(\Delta \mathbf{u}, \beta \nabla \mathbf{v})$.

Definition at line 651 of file StabilizationSD.hpp.

◆ gradp_bgradv()

void gradp_bgradv ( const Real coef,
VectorElemental vel,
MatrixElemental elmat,
const CurrentFE fe 
) const
private

Evaluate the varf $(\nabla p, \beta \nabla \mathbf{v})$.

Definition at line 545 of file StabilizationSD.hpp.

◆ lapu_gradq()

void lapu_gradq ( const Real coef,
MatrixElemental elmat,
const CurrentFE fe 
) const
private

Evaluate the varf $(\Delta \mathbf{u}, \nabla q)$.

Definition at line 714 of file StabilizationSD.hpp.

◆ f_bgradv()

void f_bgradv ( const Real coef,
SourceType &  source,
VectorElemental vel,
VectorElemental elvec,
const CurrentFE fe,
UInt  iblock,
const Real time 
) const
private

Evaluate the varf $(f, \beta \nabla \mathbf{v})$.

Definition at line 748 of file StabilizationSD.hpp.

◆ f_gradq()

void f_gradq ( const Real coef,
SourceType &  source,
VectorElemental elvec,
const CurrentFE fe,
UInt  iblock,
const Real time 
) const
private

Evaluate the varf $(f, \nabla q)$.

Definition at line 793 of file StabilizationSD.hpp.

Field Documentation

◆ M_mesh

const mesh_Type& M_mesh
private

the mesh object

Definition at line 213 of file StabilizationSD.hpp.

◆ M_dof

const dof_Type& M_dof
private

the dof object

Definition at line 215 of file StabilizationSD.hpp.

◆ M_fe

CurrentFE M_fe
private

current fe for the assembling of stabilization terms.

Definition at line 217 of file StabilizationSD.hpp.

◆ M_viscosity

Real M_viscosity
private

fluid viscosity $\nu$

Definition at line 219 of file StabilizationSD.hpp.

◆ M_gammaBeta

Real M_gammaBeta
private

Stabilization coefficient of $(c(h,dt, |\beta|, nu) \beta \nabla \mathbf{u} , \beta \nabla \mathbf{v})$.

Definition at line 221 of file StabilizationSD.hpp.

◆ M_gammaDiv

Real M_gammaDiv
private

Stabilization coefficient of $(c(h,dt) div \mathbf{u} , div \nabla \mathbf{v})$.

Definition at line 223 of file StabilizationSD.hpp.

◆ M_elMat

MatrixElemental M_elMat
private

Elementary Matrix for assembling the stabilization terms.

Definition at line 225 of file StabilizationSD.hpp.

◆ M_elVec

VectorElemental M_elVec
private

Elementary Vector for assembling the stabilization terms.

Definition at line 227 of file StabilizationSD.hpp.


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