LifeV
OseenSolver< MeshType, SolverType > Class Template Reference

This class contains an Oseen equation solver. More...

#include <OseenSolver.hpp>

+ Inheritance diagram for OseenSolver< MeshType, SolverType >:
+ Collaboration diagram for OseenSolver< MeshType, SolverType >:

Public Member Functions

preconditionerPtr_Typepreconditioner ()
 Return a shared pointer to the preconditioner (of type derived from EpetraPreconditioner) More...
 

Protected Attributes

dataPtr_Type M_oseenData
 data for Navier-Stokes solvers More...
 
FESpace< mesh_Type, MapEpetra > & M_velocityFESpace
 
FESpace< mesh_Type, MapEpetra > & M_pressureFESpace
 
Displayer M_Displayer
 MPI communicator. More...
 
MapEpetra M_localMap
 
matrixPtr_Type M_velocityMatrixMass
 mass matrix More...
 
matrixPtr_Type M_pressureMatrixMass
 mass matrix More...
 
matrixPtr_Type M_matrixStokes
 Stokes matrix: nu*stiff. More...
 
matrixPtr_Type M_matrixNoBC
 matrix to be solved More...
 
matrixPtr_Type M_matrixStabilization
 stabilization matrix More...
 
source_Type M_source
 source term for Navier-Stokes equations More...
 
vectorPtr_Type M_rightHandSideNoBC
 Right hand side for the velocity component. More...
 
vectorPtr_Type M_solution
 Global solution. More...
 
vectorPtr_Type M_residual
 residual More...
 
linearSolverPtr_Type M_linearSolver
 
bool M_steady
 
std::shared_ptr< PostProcessingBoundary< mesh_Type > > M_postProcessing
 Postprocessing class. More...
 
bool M_stabilization
 Stabilization. More...
 
bool M_reuseStabilization
 
bool M_resetStabilization
 
Int M_iterReuseStabilization
 
details::StabilizationIP< mesh_Type, DOFM_ipStabilization
 
Real M_gammaBeta
 
Real M_gammaDiv
 
Real M_gammaPress
 
const function_TypeM_betaFunction
 
bool M_divBetaUv
 
bool M_stiffStrain
 
Real M_diagonalize
 
UInt M_count
 
bool M_recomputeMatrix
 
bool M_isDiagonalBlockPreconditioner
 
MatrixElemental M_elementMatrixStiff
 Elementary matrices and vectors. More...
 
MatrixElemental M_elementMatrixMass
 
MatrixElemental M_elementMatrixPreconditioner
 
MatrixElemental M_elementMatrixDivergence
 
MatrixElemental M_elementMatrixGradient
 
VectorElemental M_elementRightHandSide
 
matrixPtr_Type M_blockPreconditioner
 
VectorElemental M_wLoc
 
VectorElemental M_uLoc
 
std::shared_ptr< vector_TypeM_un
 

Public Types

typedef MeshType mesh_Type
 
typedef SolverType linearSolver_Type
 
typedef std::shared_ptr< linearSolver_TypelinearSolverPtr_Type
 
typedef OseenData data_Type
 
typedef std::shared_ptr< data_TypedataPtr_Type
 
typedef std::function< Real(const Real &t, const Real &x, const Real &y, const Real &z, const ID &i) > function_Type
 
typedef std::function< Real(const Real &t, const Real &x, const Real &y, const Real &z, const ID &i) > source_Type
 
typedef BCHandler bcHandler_Type
 
typedef std::shared_ptr< bcHandler_TypebcHandlerPtr_Type
 
typedef linearSolver_Type::matrix_type matrix_Type
 
typedef std::shared_ptr< matrix_TypematrixPtr_Type
 
typedef linearSolver_Type::vector_type vector_Type
 
typedef std::shared_ptr< vector_TypevectorPtr_Type
 
typedef vector_Type solution_Type
 
typedef std::shared_ptr< solution_TypesolutionPtr_Type
 
typedef linearSolver_Type::prec_raw_type preconditioner_Type
 
typedef linearSolver_Type::prec_type preconditionerPtr_Type
 

Constructors & Destructor

 OseenSolver ()
 Empty constructor. More...
 
 OseenSolver (std::shared_ptr< data_Type > dataType, FESpace< mesh_Type, MapEpetra > &velocityFESpace, FESpace< mesh_Type, MapEpetra > &pressureFESpace, std::shared_ptr< Epetra_Comm > &communicator, const Int lagrangeMultiplier=0)
 Constructor. More...
 
 OseenSolver (std::shared_ptr< data_Type > dataType, FESpace< mesh_Type, MapEpetra > &velocityFESpace, FESpace< mesh_Type, MapEpetra > &pressureFESpace, std::shared_ptr< Epetra_Comm > &communicator, const MapEpetra monolithicMap, const UInt offset=0)
 Constructor. More...
 
 OseenSolver (std::shared_ptr< data_Type > dataType, FESpace< mesh_Type, MapEpetra > &velocityFESpace, FESpace< mesh_Type, MapEpetra > &pressureFESpace, const std::vector< Int > &lagrangeMultipliers, std::shared_ptr< Epetra_Comm > &communicator)
 Constructor. More...
 
virtual ~OseenSolver ()
 virtual destructor More...
 

Methods

virtual void setUp (const GetPot &dataFile)
 Set up data from GetPot. More...
 
void initialize (const function_Type &velocityFunction, const function_Type &pressureFunction)
 Initialize with velocityFunction and pressureFunction. More...
 
void initialize (const vector_Type &velocityInitialGuess, const vector_Type &pressureInitialGuess)
 Initialize with velocityInitialGuess and pressureInitialGuess. More...
 
void initialize (const vector_Type &velocityAndPressure)
 Initialize with velocityAndPressure. More...
 
virtual void buildSystem ()
 Build linear system. More...
 
virtual void updateSystem (const Real alpha, const vector_Type &betaVector, const vector_Type &sourceVector)
 Update system. More...
 
virtual void updateSystem (const Real alpha, const vector_Type &betaVector, const vector_Type &sourceVector, matrixPtr_Type matrix, const vector_Type &un)
 Update system. More...
 
void updateStabilization (matrix_Type &matrixFull)
 Update stabilization term. More...
 
virtual void updateRightHandSide (const vector_Type &rightHandSide)
 Update the right hand side. More...
 
void updateSourceTerm (const source_Type &source)
 Update the source term. More...
 
virtual void iterate (bcHandler_Type &bcHandler)
 Update convective term, boundary condition and solve the linearized ns system. More...
 
void reduceSolution (Vector &velocity, Vector &pressure)
 Reduce the local solution in global vectors. More...
 
void reduceResidual (Vector &residual)
 Reduce the residual. More...
 
void setBlockPreconditioner (matrixPtr_Type blockPreconditioner)
 Set a block preconditioner. More...
 
void getFluidMatrix (matrix_Type &matrixFull)
 Update and return the coefficient matrix. More...
 
void setupPostProc ()
 Set up post processing. More...
 
Real area (const markerID_Type &flag)
 Compute area on a boundary face with given flag. More...
 
Vector normal (const markerID_Type &flag)
 Compute the outgoing normal of a boundary face with given flag. More...
 
Vector geometricCenter (const markerID_Type &flag)
 Compute the geometric center of a boundary face with given flag. More...
 
Real flux (const markerID_Type &flag, const vector_Type &solution)
 Compute flux on a boundary face with given flag and a given solution. More...
 
Real flux (const markerID_Type &flag)
 Compute flux on a boundary face with given flag. More...
 
Real kineticNormalStress (const markerID_Type &flag, const vector_Type &solution)
 Compute the kinetic normal stress (i.e., the normal stress due to the kinetic energy) on a boundary face with a given flag and a given solution. More...
 
Real kineticNormalStress (const markerID_Type &flag)
 Compute the kinetic normal stress (i.e., the normal stress due to the kinetic energy) on a boundary face with a given flag. More...
 
Real pressure (const markerID_Type &flag, const vector_Type &solution)
 Compute average pressure on a boundary face with given flag and a given solution. More...
 
Real pressure (const markerID_Type &flag)
 Compute average pressure on a boundary face with given flag. More...
 
Real meanNormalStress (const markerID_Type &flag, bcHandler_Type &bcHandler, const vector_Type &solution)
 Compute the mean normal stress on a boundary face with a given flag and a given solution. More...
 
Real meanNormalStress (const markerID_Type &flag, bcHandler_Type &bcHandler)
 Compute the mean normal stress on a boundary face with a given flag. More...
 
Real meanTotalNormalStress (const markerID_Type &flag, bcHandler_Type &bcHandler, const vector_Type &solution)
 Compute the mean total normal stress on a boundary face with a given flag and a given solution. More...
 
Real meanTotalNormalStress (const markerID_Type &flag, bcHandler_Type &bcHandler)
 Compute the mean total normal stress on a boundary face with a given flag. More...
 
Real lagrangeMultiplier (const markerID_Type &flag, bcHandler_Type &bcHandler)
 Get the Lagrange multiplier related to a flux imposed on a given part of the boundary. More...
 
Real lagrangeMultiplier (const markerID_Type &flag, bcHandler_Type &bcHandler, const vector_Type &solution)
 Get the Lagrange multiplier related to a flux imposed on a given part of the boundary. More...
 
void resetPreconditioner (bool reset=true)
 Reset the preconditioner. More...
 
void resetStabilization ()
 Reset stabilization matrix at the same time as the preconditioner. More...
 
void updateUn ()
 Update. More...
 
void updateUn (const vector_Type &solution)
 Update for the monolithic. More...
 
void showMe (std::ostream &output=std::cout) const
 Display general information about the content of the class. More...
 

Set Methods

void setRecomputeMatrix (const bool &recomputeMatrix)
 Set. More...
 
void setSourceTerm (source_Type source)
 set the source term functor More...
 
void setTolMaxIteration (const Real &tolerance, const Int &maxIteration=-1)
 Set the tolerance and the maximum number of iterations of the linear solver. More...
 

Get Methods

const dataPtr_Typedata () const
 Return the data container of the fluid. More...
 
const Realdensity () const
 Return the density of the fluid. More...
 
const Realviscosity () const
 Return the viscosity of the fluid. More...
 
const vectorPtr_Typesolution () const
 Return the local solution vector. More...
 
const vector_Typeresidual () const
 Return the local residual vector. More...
 
FESpace< mesh_Type, MapEpetra > & velocityFESpace ()
 Return velocity FE space. More...
 
const FESpace< mesh_Type, MapEpetra > & velocityFESpace () const
 
FESpace< mesh_Type, MapEpetra > & pressureFESpace ()
 Return pressure FE space. More...
 
const FESpace< mesh_Type, MapEpetra > & pressureFESpace () const
 
const source_TypesourceTerm () const
 Get the source term. More...
 
PostProcessingBoundary< mesh_Type > & postProcessing ()
 Returns the post processing structure. More...
 
const PostProcessingBoundary< mesh_Type > & postProcessing () const
 
const MapEpetragetMap () const
 Return MapEpetra. More...
 
const std::shared_ptr< Epetra_Comm > & comm () const
 Return Epetra communicator. More...
 
const DisplayergetDisplayer () const
 Return displayer. More...
 
const bool & recomputeMatrix () const
 Return. More...
 
matrix_TypematrixNoBC ()
 Return matrix without boundary conditions. More...
 
const matrix_TypematrixNoBC () const
 
matrix_TypematrixMass ()
 Return mass matrix. More...
 
const matrix_TypematrixMass () const
 
const matrixPtr_Type matrixMassPtr () const
 
void postProcessingSetArea ()
 Set up post processing structures. More...
 
void postProcessingSetNormal ()
 Set up post processing. More...
 
void postProcessingSetPhi ()
 Set up post processing. More...
 
bool getIsDiagonalBlockPreconditioner ()
 Return a bool value if using diagonal block preconditioner. More...
 
const bool & getIsDiagonalBlockPreconditioner () const
 

Constructor

 OseenSolver (const OseenSolver &oseen)
 Empty copy constructor. More...
 

Private Methods

Real removeMean (vector_Type &x)
 Removes mean of component of vector x. More...
 
void applyBoundaryConditions (matrix_Type &matrix, vector_Type &rightHandSide, bcHandler_Type &bcHandler)
 Apply boundary conditions. More...
 
void echo (std::string message)
 Echo message. More...
 
const UIntdimVelocity () const
 Return the dim of velocity FE space. More...
 
const UIntdimPressure () const
 Return the dim of pressure FE space. More...
 

Detailed Description

template<typename MeshType, typename SolverType = LifeV::SolverAztecOO>
class LifeV::OseenSolver< MeshType, SolverType >

This class contains an Oseen equation solver.

Author
Gilles Fourestey gille.nosp@m.s.fo.nosp@m.urest.nosp@m.ey@i.nosp@m.mag.f.nosp@m.r Simone Deparis simon.nosp@m.e.de.nosp@m.paris.nosp@m.@epf.nosp@m.l.ch
Contributor:
Zhen Wang zhen..nosp@m.wang.nosp@m.@emor.nosp@m.y.ed.nosp@m.u

Definition at line 88 of file OseenSolver.hpp.

Member Typedef Documentation

◆ mesh_Type

Definition at line 96 of file OseenSolver.hpp.

◆ linearSolver_Type

typedef SolverType linearSolver_Type

Definition at line 97 of file OseenSolver.hpp.

◆ linearSolverPtr_Type

typedef std::shared_ptr<linearSolver_Type> linearSolverPtr_Type

Definition at line 98 of file OseenSolver.hpp.

◆ data_Type

Definition at line 99 of file OseenSolver.hpp.

◆ dataPtr_Type

typedef std::shared_ptr< data_Type > dataPtr_Type

Definition at line 100 of file OseenSolver.hpp.

◆ function_Type

typedef std::function< Real ( const Real& t, const Real& x, const Real& y, const Real& z, const ID& i ) > function_Type

Definition at line 103 of file OseenSolver.hpp.

◆ source_Type

typedef std::function< Real ( const Real& t, const Real& x, const Real& y, const Real& z, const ID& i ) > source_Type

Definition at line 106 of file OseenSolver.hpp.

◆ bcHandler_Type

Definition at line 108 of file OseenSolver.hpp.

◆ bcHandlerPtr_Type

typedef std::shared_ptr<bcHandler_Type> bcHandlerPtr_Type

Definition at line 109 of file OseenSolver.hpp.

◆ matrix_Type

Definition at line 114 of file OseenSolver.hpp.

◆ matrixPtr_Type

typedef std::shared_ptr<matrix_Type> matrixPtr_Type

Definition at line 116 of file OseenSolver.hpp.

◆ vector_Type

Definition at line 117 of file OseenSolver.hpp.

◆ vectorPtr_Type

typedef std::shared_ptr<vector_Type> vectorPtr_Type

Definition at line 118 of file OseenSolver.hpp.

◆ solution_Type

Definition at line 120 of file OseenSolver.hpp.

◆ solutionPtr_Type

typedef std::shared_ptr<solution_Type> solutionPtr_Type

Definition at line 121 of file OseenSolver.hpp.

◆ preconditioner_Type

typedef linearSolver_Type::prec_raw_type preconditioner_Type

Definition at line 123 of file OseenSolver.hpp.

◆ preconditionerPtr_Type

typedef linearSolver_Type::prec_type preconditionerPtr_Type

Definition at line 124 of file OseenSolver.hpp.

Constructor & Destructor Documentation

◆ OseenSolver() [1/5]

Empty constructor.

◆ OseenSolver() [2/5]

OseenSolver ( std::shared_ptr< data_Type dataType,
FESpace< mesh_Type, MapEpetra > &  velocityFESpace,
FESpace< mesh_Type, MapEpetra > &  pressureFESpace,
std::shared_ptr< Epetra_Comm > &  communicator,
const Int  lagrangeMultiplier = 0 
)

Constructor.

Parameters
dataTypeOseenData class
velocityFESpaceVelocity FE space
pressureFESpacePressure FE space
communicatorMPI communicator
lagrangeMultiplierLagrange multiplier

Definition at line 899 of file OseenSolver.hpp.

◆ OseenSolver() [3/5]

OseenSolver ( std::shared_ptr< data_Type dataType,
FESpace< mesh_Type, MapEpetra > &  velocityFESpace,
FESpace< mesh_Type, MapEpetra > &  pressureFESpace,
std::shared_ptr< Epetra_Comm > &  communicator,
const MapEpetra  monolithicMap,
const UInt  offset = 0 
)

Constructor.

Parameters
dataTypeOseenData class
velocityFESpaceVelocity FE space
pressureFESpacePressure FE space
communicatorMPI communicator
monolithicMapMapEpetra class
offset

Definition at line 964 of file OseenSolver.hpp.

◆ OseenSolver() [4/5]

OseenSolver ( std::shared_ptr< data_Type dataType,
FESpace< mesh_Type, MapEpetra > &  velocityFESpace,
FESpace< mesh_Type, MapEpetra > &  pressureFESpace,
const std::vector< Int > &  lagrangeMultipliers,
std::shared_ptr< Epetra_Comm > &  communicator 
)

Constructor.

Parameters
dataTypeOseenData class
velocityFESpaceVelocity FE space
pressureFESpacePressure FE space
lagrangeMultipliers(lagrange multipliers for the flux problem with rufaec flag)
communicatorMPI communicator

Definition at line 1021 of file OseenSolver.hpp.

◆ ~OseenSolver()

~OseenSolver ( )
virtual

virtual destructor

Definition at line 1077 of file OseenSolver.hpp.

◆ OseenSolver() [5/5]

OseenSolver ( const OseenSolver< MeshType, SolverType > &  oseen)
protected

Empty copy constructor.

Member Function Documentation

◆ setUp()

void setUp ( const GetPot dataFile)
virtual

Set up data from GetPot.

Parameters
dataFileGetPot object

Reimplemented in OseenSolverShapeDerivative< MeshType, SolverType >.

Definition at line 1089 of file OseenSolver.hpp.

◆ initialize() [1/3]

void initialize ( const function_Type velocityFunction,
const function_Type pressureFunction 
)

Initialize with velocityFunction and pressureFunction.

Parameters
velocityFunction
pressureFunction

Definition at line 1127 of file OseenSolver.hpp.

◆ initialize() [2/3]

void initialize ( const vector_Type velocityInitialGuess,
const vector_Type pressureInitialGuess 
)

Initialize with velocityInitialGuess and pressureInitialGuess.

Parameters
velocityInitialGuess
pressureInitialGuess

Definition at line 1146 of file OseenSolver.hpp.

◆ initialize() [3/3]

void initialize ( const vector_Type velocityAndPressure)

Initialize with velocityAndPressure.

Parameters
velocityAndPressure

Definition at line 1159 of file OseenSolver.hpp.

◆ buildSystem()

void buildSystem ( )
virtual

Build linear system.

Definition at line 1172 of file OseenSolver.hpp.

◆ updateSystem() [1/2]

void updateSystem ( const Real  alpha,
const vector_Type betaVector,
const vector_Type sourceVector 
)
virtual

Update system.

Parameters
alpha
betaVector
sourceVector

Definition at line 1386 of file OseenSolver.hpp.

◆ updateSystem() [2/2]

void updateSystem ( const Real  alpha,
const vector_Type betaVector,
const vector_Type sourceVector,
matrixPtr_Type  matrix,
const vector_Type un 
)
virtual

Update system.

Parameters
alpha
betaVector
sourceVector
matrix
un

managing the convective term

managing the convective term : semi-implicit approximation of the convective term

Definition at line 1410 of file OseenSolver.hpp.

◆ updateStabilization()

void updateStabilization ( matrix_Type matrixFull)

Update stabilization term.

Parameters
matrixFull

Definition at line 1623 of file OseenSolver.hpp.

◆ updateRightHandSide()

virtual void updateRightHandSide ( const vector_Type rightHandSide)
inlinevirtual

Update the right hand side.

Parameters
rightHandSideright hand side

Definition at line 250 of file OseenSolver.hpp.

◆ updateSourceTerm()

void updateSourceTerm ( const source_Type source)

Update the source term.

Parameters
sourceTerm

Definition at line 1634 of file OseenSolver.hpp.

◆ iterate()

void iterate ( bcHandler_Type bcHandler)
virtual

Update convective term, boundary condition and solve the linearized ns system.

Parameters
bcHandlerBC handler

Definition at line 1663 of file OseenSolver.hpp.

◆ reduceSolution()

void reduceSolution ( Vector velocity,
Vector pressure 
)

Reduce the local solution in global vectors.

Parameters
velocity
pressure

Definition at line 1722 of file OseenSolver.hpp.

◆ reduceResidual()

void reduceResidual ( Vector residual)

Reduce the residual.

Parameters
residual

Definition at line 1743 of file OseenSolver.hpp.

◆ setBlockPreconditioner()

void setBlockPreconditioner ( matrixPtr_Type  blockPreconditioner)

Set a block preconditioner.

Block preconditioner

Definition at line 1759 of file OseenSolver.hpp.

◆ getFluidMatrix()

void getFluidMatrix ( matrix_Type matrixFull)

Update and return the coefficient matrix.

Parameters
matrixFullThe coefficient matrix

Definition at line 1767 of file OseenSolver.hpp.

◆ setupPostProc()

void setupPostProc ( )

Set up post processing.

Definition at line 2051 of file OseenSolver.hpp.

◆ area()

Real area ( const markerID_Type flag)

Compute area on a boundary face with given flag.

Parameters
flag
Returns
area

Definition at line 1834 of file OseenSolver.hpp.

◆ normal()

Vector normal ( const markerID_Type flag)

Compute the outgoing normal of a boundary face with given flag.

Parameters
flag
Returns
boundary normal

Definition at line 1841 of file OseenSolver.hpp.

◆ geometricCenter()

Vector geometricCenter ( const markerID_Type flag)

Compute the geometric center of a boundary face with given flag.

Parameters
flag
Returns
geometric center

Definition at line 1848 of file OseenSolver.hpp.

◆ flux() [1/2]

Real flux ( const markerID_Type flag,
const vector_Type solution 
)

Compute flux on a boundary face with given flag and a given solution.

Parameters
flag
solution
Returns
flux

Definition at line 1803 of file OseenSolver.hpp.

◆ flux() [2/2]

Real flux ( const markerID_Type flag)

Compute flux on a boundary face with given flag.

Parameters
flag
Returns
flux

Definition at line 1796 of file OseenSolver.hpp.

◆ kineticNormalStress() [1/2]

Real kineticNormalStress ( const markerID_Type flag,
const vector_Type solution 
)

Compute the kinetic normal stress (i.e., the normal stress due to the kinetic energy) on a boundary face with a given flag and a given solution.

See also
[2] [15]

This method computes the following quantity:

\[ \mathcal{K} = \frac{1}{2}\rho_\textrm{F}\frac{1}{\left|\Gamma^t_{\textrm{F},j}\right|}\displaystyle\int_{\Gamma^t_{\textrm{F},j}}\left({\mathbf{u}}_\textrm{F} \mathbf{\cdot} {\mathbf{n}}_\textrm{F}\right)^2 \textrm{d} \Gamma \]

Parameters
flagboundary flag
solutionproblem solution
Returns
kinetic normal stress

Definition at line 1822 of file OseenSolver.hpp.

◆ kineticNormalStress() [2/2]

Real kineticNormalStress ( const markerID_Type flag)

Compute the kinetic normal stress (i.e., the normal stress due to the kinetic energy) on a boundary face with a given flag.

See also
[2] [15]

This method computes the following quantity:

\[ \mathcal{K} = \frac{1}{2}\rho_\textrm{F}\frac{1}{\left|\Gamma^t_{\textrm{F},j}\right|}\displaystyle\int_{\Gamma^t_{\textrm{F},j}}\left({\mathbf{u}}_\textrm{F} \mathbf{\cdot} {\mathbf{n}}_\textrm{F}\right)^2 \textrm{d} \Gamma \]

Parameters
flagboundary flag
Returns
kinetic normal stress

Definition at line 1815 of file OseenSolver.hpp.

◆ pressure() [1/2]

Real pressure ( const markerID_Type flag,
const vector_Type solution 
)

Compute average pressure on a boundary face with given flag and a given solution.

Parameters
flag
solution
Returns
average pressure

Definition at line 1862 of file OseenSolver.hpp.

◆ pressure() [2/2]

Real pressure ( const markerID_Type flag)

Compute average pressure on a boundary face with given flag.

Parameters
flag
Returns
average pressure

Definition at line 1855 of file OseenSolver.hpp.

◆ meanNormalStress() [1/2]

Real meanNormalStress ( const markerID_Type flag,
bcHandler_Type bcHandler,
const vector_Type solution 
)

Compute the mean normal stress on a boundary face with a given flag and a given solution.

See also
[2] [15]

The mean normal stress is defined as the average of the normal component of the traction vector.

\[ \mathcal{S}^{3-D}_j = \frac{1}{\left|\Gamma^t_{\textrm{F},j}\right|}\displaystyle\int_{\Gamma^t_{\textrm{F},j}} \left(\sigma_\textrm{F} \mathbf{\cdot} \mathbf{n}_\textrm{F}\right) \mathbf{\cdot} \mathbf{n}_\textrm{F} \: \textrm{d} \Gamma \]

Parameters
flagFlag of the boundary face
bcHandlerBChandler containing the boundary conditions of the problem.
solutionVector containing the solution of the problem (and also the Lagrange multipliers at the end).
Returns
mean normal stress

Definition at line 1883 of file OseenSolver.hpp.

◆ meanNormalStress() [2/2]

Real meanNormalStress ( const markerID_Type flag,
bcHandler_Type bcHandler 
)

Compute the mean normal stress on a boundary face with a given flag.

See also
[2] [15]

The mean normal stress is defined as the average of the normal component of the traction vector.

\[ \mathcal{S}^{3-D}_j = \frac{1}{\left|\Gamma^t_{\textrm{F},j}\right|}\displaystyle\int_{\Gamma^t_{\textrm{F},j}} \left(\sigma_\textrm{F} \mathbf{\cdot} \mathbf{n}_\textrm{F}\right) \mathbf{\cdot} \mathbf{n}_\textrm{F} \: \textrm{d} \Gamma \]

TODO The current version returns the exact mean normal stress if a flow rate boundary condition is imposed on the chosen boundary face. On the contrary, if other boundary conditions are applied, the mean normal stress is approximated with the mean pressure, which is a reasonable approximation for several applications.

Parameters
flagFlag of the boundary face
bcHandlerBChandler containing the boundary conditions of the problem.
Returns
mean normal stress

Definition at line 1876 of file OseenSolver.hpp.

◆ meanTotalNormalStress() [1/2]

Real meanTotalNormalStress ( const markerID_Type flag,
bcHandler_Type bcHandler,
const vector_Type solution 
)

Compute the mean total normal stress on a boundary face with a given flag and a given solution.

See also
[2] [15]

The mean total normal stress is defined as the average of the normal component of the traction vector minus the kinetic contribution.

\[ \mathcal{T}^{3-D}_j = \frac{1}{\left|\Gamma^t_{\textrm{F},j}\right|}\displaystyle\int_{\Gamma^t_{\textrm{F},j}} \left(\sigma_\textrm{F} \mathbf{\cdot} \mathbf{n}_\textrm{F}\right) \mathbf{\cdot} \mathbf{n}_\textrm{F} \: \textrm{d} \Gamma - \frac{1}{2}\rho_\textrm{F}\frac{1}{\left|\Gamma^t_{\textrm{F},j}\right|}\displaystyle\int_{\Gamma^t_{\textrm{F},j}}\left(\mathbf{u}_\textrm{F} \mathbf{\cdot} \mathbf{n}_\textrm{F}\right)^2 \: \textrm{d} \Gamma \]

TODO The current version returns the exact mean normal stress if a flow rate boundary condition is imposed on the chosen boundary face. On the contrary, if other boundary conditions are applied, the mean normal stress is approximated with the mean pressure, which is a reasonable approximation for several applications.

Parameters
flagFlag of the boundary face
bcHandlerBChandler containing the boundary conditions of the problem.
solutionVector containing the solution of the problem (and also the Lagrange multipliers at the end).
Returns
mean total normal stress

Definition at line 1909 of file OseenSolver.hpp.

◆ meanTotalNormalStress() [2/2]

Real meanTotalNormalStress ( const markerID_Type flag,
bcHandler_Type bcHandler 
)

Compute the mean total normal stress on a boundary face with a given flag.

See also
[2] [15]

The mean total normal stress is defined as the average of the normal component of the traction vector minus the kinetic contribution.

\[ \mathcal{T}^{3-D}_j = \frac{1}{\left|\Gamma^t_{\textrm{F},j}\right|}\displaystyle\int_{\Gamma^t_{\textrm{F},j}} \left(\sigma_\textrm{F} \mathbf{\cdot} \mathbf{n}_\textrm{F}\right) \mathbf{\cdot} \mathbf{n}_\textrm{F} \: \textrm{d} \Gamma - \frac{1}{2}\rho_\textrm{F}\frac{1}{\left|\Gamma^t_{\textrm{F},j}\right|}\displaystyle\int_{\Gamma^t_{\textrm{F},j}}\left(\mathbf{u}_\textrm{F} \mathbf{\cdot} \mathbf{n}_\textrm{F}\right)^2 \: \textrm{d} \Gamma \]

Parameters
flagFlag of the boundary face
bcHandlerBChandler containing the boundary conditions of the problem.
Returns
mean total normal stress

Definition at line 1902 of file OseenSolver.hpp.

◆ lagrangeMultiplier() [1/2]

Real lagrangeMultiplier ( const markerID_Type flag,
bcHandler_Type bcHandler 
)

Get the Lagrange multiplier related to a flux imposed on a given part of the boundary.

Parameters
flagFlag of the boundary face associated with the flux and the Lagrange multiplier we want.
bcHandlerBChandler containing the boundary conditions of the problem.
Returns
Lagrange multiplier

Definition at line 1916 of file OseenSolver.hpp.

◆ lagrangeMultiplier() [2/2]

Real lagrangeMultiplier ( const markerID_Type flag,
bcHandler_Type bcHandler,
const vector_Type solution 
)

Get the Lagrange multiplier related to a flux imposed on a given part of the boundary.

Parameters
flagFlag of the boundary face associated with the flux and the Lagrange multiplier we want.
bcHandlerBChandler containing the boundary conditions of the problem.
solutionVector containing the solution of the problem (and also the Lagrange multipliers at the end).
Returns
Lagrange multiplier

Definition at line 1923 of file OseenSolver.hpp.

◆ resetPreconditioner()

void resetPreconditioner ( bool  reset = true)
inline

Reset the preconditioner.

Parameters
resetReset preconditioner.

Definition at line 486 of file OseenSolver.hpp.

◆ resetStabilization()

void resetStabilization ( )
inline

Reset stabilization matrix at the same time as the preconditioner.

Definition at line 495 of file OseenSolver.hpp.

◆ updateUn() [1/2]

void updateUn ( )
inline

Update.

Definition at line 501 of file OseenSolver.hpp.

◆ updateUn() [2/2]

void updateUn ( const vector_Type solution)
inline

Update for the monolithic.

Definition at line 507 of file OseenSolver.hpp.

◆ showMe()

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

Display general information about the content of the class.

Parameters
outputspecify the output format (std::cout by default)

◆ setRecomputeMatrix()

void setRecomputeMatrix ( const bool &  recomputeMatrix)
inline

Set.

Parameters
recomputeMatrix

Definition at line 527 of file OseenSolver.hpp.

◆ setSourceTerm()

void setSourceTerm ( source_Type  source)
inline

set the source term functor

Parameters
source

Definition at line 536 of file OseenSolver.hpp.

◆ setTolMaxIteration()

void setTolMaxIteration ( const Real tolerance,
const Int maxIteration = -1 
)

Set the tolerance and the maximum number of iterations of the linear solver.

Parameters
toleranceTolerance
maxIterationmaximum number of iterations

Definition at line 2063 of file OseenSolver.hpp.

◆ data()

const dataPtr_Type& data ( ) const
inline

Return the data container of the fluid.

Returns
data container of the fluid

Definition at line 557 of file OseenSolver.hpp.

◆ density()

const Real& density ( ) const
inline

Return the density of the fluid.

Returns
Density of the fluid

Definition at line 566 of file OseenSolver.hpp.

◆ viscosity()

const Real& viscosity ( ) const
inline

Return the viscosity of the fluid.

Returns
Viscosity of the fluid

Definition at line 575 of file OseenSolver.hpp.

◆ solution()

const vectorPtr_Type& solution ( ) const
inline

Return the local solution vector.

Returns
vectorPtr_Type Solution vector

Definition at line 584 of file OseenSolver.hpp.

◆ residual()

const vector_Type& residual ( ) const
inline

Return the local residual vector.

Returns
Residual vector

Definition at line 593 of file OseenSolver.hpp.

◆ velocityFESpace() [1/2]

FESpace<mesh_Type, MapEpetra>& velocityFESpace ( )
inline

Return velocity FE space.

Returns
velocity FE space

Definition at line 602 of file OseenSolver.hpp.

◆ velocityFESpace() [2/2]

const FESpace<mesh_Type, MapEpetra>& velocityFESpace ( ) const
inline

Definition at line 607 of file OseenSolver.hpp.

◆ pressureFESpace() [1/2]

FESpace<mesh_Type, MapEpetra>& pressureFESpace ( )
inline

Return pressure FE space.

Returns
pressure FE space

Definition at line 616 of file OseenSolver.hpp.

◆ pressureFESpace() [2/2]

const FESpace<mesh_Type, MapEpetra>& pressureFESpace ( ) const
inline

Definition at line 621 of file OseenSolver.hpp.

◆ sourceTerm()

const source_Type& sourceTerm ( ) const
inline

Get the source term.

Returns
Source term

Definition at line 630 of file OseenSolver.hpp.

◆ postProcessing() [1/2]

PostProcessingBoundary<mesh_Type>& postProcessing ( )
inline

Returns the post processing structure.

Returns
Post processing

Definition at line 639 of file OseenSolver.hpp.

◆ postProcessing() [2/2]

const PostProcessingBoundary<mesh_Type>& postProcessing ( ) const
inline

Definition at line 644 of file OseenSolver.hpp.

◆ getMap()

const MapEpetra& getMap ( ) const
inline

Return MapEpetra.

Returns
MapEpetra

Definition at line 653 of file OseenSolver.hpp.

◆ comm()

const std::shared_ptr<Epetra_Comm>& comm ( ) const
inline

Return Epetra communicator.

Returns
Epetra communicator

Definition at line 662 of file OseenSolver.hpp.

◆ getDisplayer()

const Displayer& getDisplayer ( ) const
inline

Return displayer.

Returns

Definition at line 671 of file OseenSolver.hpp.

◆ recomputeMatrix()

const bool& recomputeMatrix ( ) const
inline

Return.

Returns
recomputeMatrix

Definition at line 680 of file OseenSolver.hpp.

◆ matrixNoBC() [1/2]

matrix_Type& matrixNoBC ( )
inline

Return matrix without boundary conditions.

Returns
Matrix without boundary conditions

Definition at line 689 of file OseenSolver.hpp.

◆ matrixNoBC() [2/2]

const matrix_Type& matrixNoBC ( ) const
inline

Definition at line 694 of file OseenSolver.hpp.

◆ matrixMass() [1/2]

matrix_Type& matrixMass ( )
inline

Return mass matrix.

Returns
Mass matrix

Definition at line 703 of file OseenSolver.hpp.

◆ matrixMass() [2/2]

const matrix_Type& matrixMass ( ) const
inline

Definition at line 708 of file OseenSolver.hpp.

◆ matrixMassPtr()

const matrixPtr_Type matrixMassPtr ( ) const
inline

Definition at line 713 of file OseenSolver.hpp.

◆ postProcessingSetArea()

void postProcessingSetArea ( )

Set up post processing structures.

Definition at line 1775 of file OseenSolver.hpp.

◆ postProcessingSetNormal()

void postProcessingSetNormal ( )

Set up post processing.

Definition at line 1782 of file OseenSolver.hpp.

◆ postProcessingSetPhi()

void postProcessingSetPhi ( )

Set up post processing.

Definition at line 1789 of file OseenSolver.hpp.

◆ getIsDiagonalBlockPreconditioner() [1/2]

bool getIsDiagonalBlockPreconditioner ( )
inline

Return a bool value if using diagonal block preconditioner.

Definition at line 732 of file OseenSolver.hpp.

◆ getIsDiagonalBlockPreconditioner() [2/2]

const bool& getIsDiagonalBlockPreconditioner ( ) const
inline

Definition at line 737 of file OseenSolver.hpp.

◆ preconditioner()

preconditionerPtr_Type& preconditioner ( )
inline

Return a shared pointer to the preconditioner (of type derived from EpetraPreconditioner)

Definition at line 745 of file OseenSolver.hpp.

◆ removeMean()

Real removeMean ( vector_Type x)
protected

Removes mean of component of vector x.

Parameters
x
Returns

Definition at line 1947 of file OseenSolver.hpp.

◆ applyBoundaryConditions()

void applyBoundaryConditions ( matrix_Type matrix,
vector_Type rightHandSide,
bcHandler_Type bcHandler 
)
protected

Apply boundary conditions.

Parameters
matrix
rightHandSide
bcHandler

Definition at line 2005 of file OseenSolver.hpp.

◆ echo()

void echo ( std::string  message)
protected

Echo message.

Parameters
message

◆ dimVelocity()

const UInt& dimVelocity ( ) const
inlineprotected

Return the dim of velocity FE space.

Definition at line 787 of file OseenSolver.hpp.

◆ dimPressure()

const UInt& dimPressure ( ) const
inlineprotected

Return the dim of pressure FE space.

Definition at line 793 of file OseenSolver.hpp.

Field Documentation

◆ M_oseenData

dataPtr_Type M_oseenData
protected

data for Navier-Stokes solvers

Definition at line 803 of file OseenSolver.hpp.

◆ M_velocityFESpace

FESpace<mesh_Type, MapEpetra>& M_velocityFESpace
protected

Definition at line 806 of file OseenSolver.hpp.

◆ M_pressureFESpace

FESpace<mesh_Type, MapEpetra>& M_pressureFESpace
protected

Definition at line 807 of file OseenSolver.hpp.

◆ M_Displayer

Displayer M_Displayer
protected

MPI communicator.

Definition at line 810 of file OseenSolver.hpp.

◆ M_localMap

MapEpetra M_localMap
protected

Definition at line 812 of file OseenSolver.hpp.

◆ M_velocityMatrixMass

matrixPtr_Type M_velocityMatrixMass
protected

mass matrix

Definition at line 815 of file OseenSolver.hpp.

◆ M_pressureMatrixMass

matrixPtr_Type M_pressureMatrixMass
protected

mass matrix

Definition at line 818 of file OseenSolver.hpp.

◆ M_matrixStokes

matrixPtr_Type M_matrixStokes
protected

Stokes matrix: nu*stiff.

Definition at line 821 of file OseenSolver.hpp.

◆ M_matrixNoBC

matrixPtr_Type M_matrixNoBC
protected

matrix to be solved

matrix without boundary conditions

Definition at line 827 of file OseenSolver.hpp.

◆ M_matrixStabilization

matrixPtr_Type M_matrixStabilization
protected

stabilization matrix

Definition at line 830 of file OseenSolver.hpp.

◆ M_source

source_Type M_source
protected

source term for Navier-Stokes equations

Definition at line 833 of file OseenSolver.hpp.

◆ M_rightHandSideNoBC

vectorPtr_Type M_rightHandSideNoBC
protected

Right hand side for the velocity component.

Definition at line 836 of file OseenSolver.hpp.

◆ M_solution

vectorPtr_Type M_solution
protected

Global solution.

Definition at line 839 of file OseenSolver.hpp.

◆ M_residual

vectorPtr_Type M_residual
protected

residual

Definition at line 842 of file OseenSolver.hpp.

◆ M_linearSolver

linearSolverPtr_Type M_linearSolver
protected

Definition at line 844 of file OseenSolver.hpp.

◆ M_steady

bool M_steady
protected

Definition at line 846 of file OseenSolver.hpp.

◆ M_postProcessing

std::shared_ptr<PostProcessingBoundary<mesh_Type> > M_postProcessing
protected

Postprocessing class.

Definition at line 849 of file OseenSolver.hpp.

◆ M_stabilization

bool M_stabilization
protected

Stabilization.

Definition at line 852 of file OseenSolver.hpp.

◆ M_reuseStabilization

bool M_reuseStabilization
protected

Definition at line 853 of file OseenSolver.hpp.

◆ M_resetStabilization

bool M_resetStabilization
protected

Definition at line 854 of file OseenSolver.hpp.

◆ M_iterReuseStabilization

Int M_iterReuseStabilization
protected

Definition at line 855 of file OseenSolver.hpp.

◆ M_ipStabilization

details::StabilizationIP<mesh_Type, DOF> M_ipStabilization
protected

Definition at line 857 of file OseenSolver.hpp.

◆ M_gammaBeta

Real M_gammaBeta
protected

Definition at line 858 of file OseenSolver.hpp.

◆ M_gammaDiv

Real M_gammaDiv
protected

Definition at line 859 of file OseenSolver.hpp.

◆ M_gammaPress

Real M_gammaPress
protected

Definition at line 860 of file OseenSolver.hpp.

◆ M_betaFunction

const function_Type* M_betaFunction
protected

Definition at line 862 of file OseenSolver.hpp.

◆ M_divBetaUv

bool M_divBetaUv
protected

Definition at line 864 of file OseenSolver.hpp.

◆ M_stiffStrain

bool M_stiffStrain
protected

Definition at line 866 of file OseenSolver.hpp.

◆ M_diagonalize

Real M_diagonalize
protected

Definition at line 869 of file OseenSolver.hpp.

◆ M_count

UInt M_count
protected

Definition at line 871 of file OseenSolver.hpp.

◆ M_recomputeMatrix

bool M_recomputeMatrix
protected

Definition at line 873 of file OseenSolver.hpp.

◆ M_isDiagonalBlockPreconditioner

bool M_isDiagonalBlockPreconditioner
protected

Definition at line 875 of file OseenSolver.hpp.

◆ M_elementMatrixStiff

MatrixElemental M_elementMatrixStiff
protected

Elementary matrices and vectors.

Definition at line 878 of file OseenSolver.hpp.

◆ M_elementMatrixMass

MatrixElemental M_elementMatrixMass
protected

Definition at line 879 of file OseenSolver.hpp.

◆ M_elementMatrixPreconditioner

MatrixElemental M_elementMatrixPreconditioner
protected

Definition at line 880 of file OseenSolver.hpp.

◆ M_elementMatrixDivergence

MatrixElemental M_elementMatrixDivergence
protected

Definition at line 881 of file OseenSolver.hpp.

◆ M_elementMatrixGradient

MatrixElemental M_elementMatrixGradient
protected

Definition at line 882 of file OseenSolver.hpp.

◆ M_elementRightHandSide

VectorElemental M_elementRightHandSide
protected

Definition at line 883 of file OseenSolver.hpp.

◆ M_blockPreconditioner

matrixPtr_Type M_blockPreconditioner
protected

Definition at line 884 of file OseenSolver.hpp.

◆ M_wLoc

VectorElemental M_wLoc
protected

Definition at line 885 of file OseenSolver.hpp.

◆ M_uLoc

VectorElemental M_uLoc
protected

Definition at line 886 of file OseenSolver.hpp.

◆ M_un

std::shared_ptr<vector_Type> M_un
protected

Definition at line 887 of file OseenSolver.hpp.


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