Equa_LinearPDE< NEN_, NSN_ > Class Template Reference

Abstract class for Finite Element classes for lienar PDEs'. More...

#include <Equa_LinearPDE.h>

Inheritance diagram for Equa_LinearPDE< NEN_, NSN_ >:
Equation< NEN_, NEN_, NSN_, NSN_ > Equa

Public Member Functions

 Equa_LinearPDE ()
 Default constructor.
 
virtual ~Equa_LinearPDE ()
 Destructor.
 
void setNoLumping ()
 Set no lumping.
 
void setStab ()
 Set stabilization for convection term.
 
void set_00 (real_t a=1.0)
 Set coefficient for term (0,0): 0th order in time and space.
 
void set_00 (Fct &f)
 Set coefficient for term (0,0): 0th order in time and space.
 
void set_00 (const string &f)
 Set coefficient for term (0,0): 0th order in time and space.
 
void set_10 (real_t a=1.0)
 Set coefficient for term (1,0): 1st order in time, 0th order in space.
 
void set_10 (Fct &f)
 Set coefficient for term (1,0): 1st order in time, 0th order in space.
 
void set_20 (real_t a=1.0)
 Set coefficient for term (2,0): 2nd order in time, 0th order in space.
 
void set_20 (Fct &f)
 Set coefficient for term (2,0): 2nd order in time, 0th order in space.
 
void set_01 (real_t a=1.0)
 Set coefficient for term (0,1): 0th order in time, 1st order in space.
 
void set_01 (Point< real_t > &a)
 Set coefficient for term (0,1): 0th order in time, 1st order in space.
 
void set_01 (Fct &f)
 Set coefficient for term (0,1): 0th order in time, 1st order in time and space.
 
void set_01 (Vect< real_t > &a)
 Set coefficient for term (0,1): 0th order in time, 1st order in space.
 
void set_02 (real_t a=1.0)
 Set coefficient for term (0,2): 0th order in time, 2nd order in space.
 
void set_02 (Fct &f)
 Set coefficient for term (0,2): 0th order in time, 2nd order in time and space.
 
virtual void Mat_00 (real_t coef=1.0)
 Add 0th order term, in time and space, to left-hand side.
 
virtual void Mat_10 (real_t coef=1.0)
 Add 1st order term in time, 0th in space to left-hand side.
 
virtual void Mat_20 (real_t coef=1.0)
 Add 2nd order term in time, 0th in space to left-hand side.
 
virtual void Mat_01 (real_t coef=1.0)
 Add 0th order term in time, 1st in space to left-hand side.
 
virtual void Mat_02 (real_t coef=1.0)
 Add 0th order term in time, 2nd in space to left-hand side.
 
virtual void BodyRHS (const Vect< real_t > &f)
 Add body right-hand side term to right-hand side.
 
virtual void BoundaryRHS (const Vect< real_t > &f)
 Add boundary right-hand side term to right-hand side.
 
void build ()
 Build the linear system of equations for the steady state case.
 
void build (TimeStepping &s)
 Build the linear system of equations.
 
void build (EigenProblemSolver &e)
 Build the linear system for an eigenvalue problem.
 
- Public Member Functions inherited from Equation< NEN_, NEN_, NSN_, NSN_ >
 Equation ()
 
 Equation (Mesh &mesh)
 Constructor with mesh instance.
 
 Equation (Mesh &mesh, Vect< real_t > &u)
 Constructor with mesh instance and solution vector.
 
 Equation (Mesh &mesh, Vect< real_t > &u, real_t &init_time, real_t &final_time, real_t &time_step)
 Constructor with mesh instance, matrix and right-hand side.
 
 ~Equation ()
 Destructor.
 
void updateBC (const Element &el, const Vect< real_t > &bc)
 Update Right-Hand side by taking into account essential boundary conditions.
 
void DiagBC (DOFSupport dof_type=NODE_DOF, int dof=0)
 Update element matrix to impose bc by diagonalization technique.
 
void LocalNodeVector (Vect< real_t > &b)
 Localize element vector from a Vect instance.
 
void ElementNodeVector (const Vect< real_t > &b, LocalVect< real_t, NEE_ > &be)
 Localize element vector from a Vect instance.
 
void ElementNodeVector (const Vect< real_t > &b, LocalVect< real_t, NEN_ > &be, int dof)
 Localize Element Vector from a Vect instance.
 
void SideNodeVector (const Vect< real_t > &b, LocalVect< real_t, NSE_ > &bs)
 Localize side vector from a Vect instance.
 
void SideSideVector (const Vect< real_t > &b, vector< real_t > &bs)
 Localize side vector from a Vect instance.
 
void ElementNodeVectorSingleDOF (const Vect< real_t > &b, LocalVect< real_t, NEN_ > &be)
 Localize Element Vector from a Vect instance.
 
void ElementSideVector (const Vect< real_t > &b, LocalVect< real_t, NSE_ > &be)
 Localize Element Vector from a Vect instance.
 
void ElementVector (const Vect< real_t > &b, DOFSupport dof_type=NODE_DOF, int flag=0)
 Localize element vector.
 
void SideVector (const Vect< real_t > &b, vector< real_t > &sb)
 Localize side vector.
 
void ElementNodeCoordinates ()
 Localize coordinates of element nodes.
 
void SideNodeCoordinates ()
 Localize coordinates of side nodes.
 
void ElementAssembly (Matrix< real_t > *A)
 Assemble element matrix into global one.
 
void ElementAssembly (BMatrix< real_t > &A)
 Assemble element matrix into global one.
 
void ElementAssembly (SkSMatrix< real_t > &A)
 Assemble element matrix into global one.
 
void ElementAssembly (SkMatrix< real_t > &A)
 Assemble element matrix into global one.
 
void ElementAssembly (SpMatrix< real_t > &A)
 Assemble element matrix into global one.
 
void ElementAssembly (TrMatrix< real_t > &A)
 Assemble element matrix into global one.
 
void ElementAssembly (Vect< real_t > &v)
 Assemble element vector into global one.
 
void DGElementAssembly (Matrix< real_t > *A)
 Assemble element matrix into global one for the Discontinuous Galerkin approximation.
 
void DGElementAssembly (SkSMatrix< real_t > &A)
 Assemble element matrix into global one for the Discontinuous Galerkin approximation.
 
void DGElementAssembly (SkMatrix< real_t > &A)
 Assemble element matrix into global one for the Discontinuous Galerkin approximation.
 
void DGElementAssembly (SpMatrix< real_t > &A)
 Assemble element matrix into global one for the Discontinuous Galerkin approximation.
 
void DGElementAssembly (TrMatrix< real_t > &A)
 Assemble element matrix into global one for the Discontinuous Galerkin approximation.
 
void SideAssembly (Matrix< real_t > *A)
 Assemble side (edge or face) matrix into global one.
 
void SideAssembly (SkSMatrix< real_t > &A)
 Assemble side (edge or face) matrix into global one.
 
void SideAssembly (SkMatrix< real_t > &A)
 Assemble side (edge or face) matrix into global one.
 
void SideAssembly (SpMatrix< real_t > &A)
 Assemble side (edge or face) matrix into global one.
 
void SideAssembly (Vect< real_t > &v)
 Assemble side (edge or face) vector into global one.
 
void AxbAssembly (const Element &el, const Vect< real_t > &x, Vect< real_t > &b)
 Assemble product of element matrix by element vector into global vector.
 
void AxbAssembly (const Side &sd, const Vect< real_t > &x, Vect< real_t > &b)
 Assemble product of side matrix by side vector into global vector.
 
size_t getNbNodes () const
 Return number of element nodes.
 
size_t getNbEq () const
 Return number of element equations.
 
real_t setMaterialProperty (const string &exp, const string &prop)
 Define a material property by an algebraic expression.
 
- Public Member Functions inherited from Equa
 Equa ()
 Default constructor.
 
virtual ~Equa ()
 Destructor.
 
void setMesh (Mesh &m)
 Define mesh and renumber DOFs after removing imposed ones.
 
MeshgetMesh () const
 Return reference to Mesh instance.
 
LinearSolvergetLinearSolver ()
 Return reference to linear solver instance.
 
Matrix< real_t > * getMatrix () const
 Return pointer to matrix.
 
void setSolver (Iteration ls, Preconditioner pc=IDENT_PREC)
 Choose solver for the linear system.
 
void setMatrixType (int t)
 Choose type of matrix.
 
int solveLinearSystem (Matrix< real_t > *A, Vect< real_t > &b, Vect< real_t > &x)
 Solve the linear system with given matrix and right-hand side.
 
int solveLinearSystem (Vect< real_t > &b, Vect< real_t > &x)
 Solve the linear system with given right-hand side.
 
void LinearSystemInfo ()
 Print info on linear system solver.
 

Detailed Description

template<size_t NEN_, size_t NSN_>
class OFELI::Equa_LinearPDE< NEN_, NSN_ >

Abstract class for Finite Element classes for lienar PDEs'.

Template Parameters
<NEN_>Number of element nodes
<NSN_>Number of side nodes

Constructor & Destructor Documentation

◆ Equa_LinearPDE()

template<size_t NEN_, size_t NSN_>
Equa_LinearPDE ( )

Default constructor.

Constructs an empty equation.

Member Function Documentation

◆ BodyRHS()

template<size_t NEN_, size_t NSN_>
virtual void BodyRHS ( const Vect< real_t > &  f)
virtual

Add body right-hand side term to right-hand side.

Parameters
[in]fVector containing source at nodes.

Reimplemented in LinearPDE1D, and LinearPDE2D.

◆ BoundaryRHS()

template<size_t NEN_, size_t NSN_>
virtual void BoundaryRHS ( const Vect< real_t > &  f)
virtual

Add boundary right-hand side term to right-hand side.

Parameters
[in]fVector containing source at nodes.

Reimplemented in LinearPDE2D.

◆ build() [1/2]

template<size_t NEN_, size_t NSN_>
void build ( EigenProblemSolver e)

Build the linear system for an eigenvalue problem.

Parameters
[in]eReference to used EigenProblemSolver instance

◆ build() [2/2]

template<size_t NEN_, size_t NSN_>
void build ( TimeStepping s)

Build the linear system of equations.

Before using this function, one must have properly selected appropriate options for:

  • The choice of a steady state or transient analysis. By default, the analysis is stationary
  • In the case of transient analysis, the choice of a time integration scheme and a lumped or consistent capacity matrix. If transient analysis is chosen, the lumped capacity matrix option is chosen by default, and the implicit Euler scheme is used by default for time integration.
Parameters
[in]sReference to used TimeStepping instance

◆ Mat_00()

template<size_t NEN_, size_t NSN_>
virtual void Mat_00 ( real_t  coef = 1.0)
virtual

Add 0th order term, in time and space, to left-hand side.

Parameters
[in]coefcoefficient to multiply by the matrix before adding [Default: 1]

Reimplemented in LinearPDE1D, and LinearPDE2D.

◆ Mat_01()

template<size_t NEN_, size_t NSN_>
virtual void Mat_01 ( real_t  coef = 1.0)
virtual

Add 0th order term in time, 1st in space to left-hand side.

Parameters
[in]coefcoefficient to multiply by the matrix before adding [Default: 1]

Reimplemented in LinearPDE1D, and LinearPDE2D.

◆ Mat_02()

template<size_t NEN_, size_t NSN_>
virtual void Mat_02 ( real_t  coef = 1.0)
virtual

Add 0th order term in time, 2nd in space to left-hand side.

Parameters
[in]coefcoefficient to multiply by the matrix before adding [Default: 1]

Reimplemented in LinearPDE1D, and LinearPDE2D.

◆ Mat_10()

template<size_t NEN_, size_t NSN_>
virtual void Mat_10 ( real_t  coef = 1.0)
virtual

Add 1st order term in time, 0th in space to left-hand side.

Parameters
[in]coefcoefficient to multiply by the matrix before adding [Default: 1]

Reimplemented in LinearPDE1D, and LinearPDE2D.

◆ Mat_20()

template<size_t NEN_, size_t NSN_>
virtual void Mat_20 ( real_t  coef = 1.0)
virtual

Add 2nd order term in time, 0th in space to left-hand side.

Parameters
[in]coefcoefficient to multiply by the matrix before adding [Default: 1]

Reimplemented in LinearPDE1D, and LinearPDE2D.

◆ set_00() [1/3]

template<size_t NEN_, size_t NSN_>
void set_00 ( const string &  f)

Set coefficient for term (0,0): 0th order in time and space.

Parameters
[in]fFunction to multiply by 0-th order term (Function of x and t)

◆ set_00() [2/3]

template<size_t NEN_, size_t NSN_>
void set_00 ( Fct &  f)

Set coefficient for term (0,0): 0th order in time and space.

Parameters
[in]fFunction to multiply by 0-th order term (Function of x and t )

◆ set_00() [3/3]

template<size_t NEN_, size_t NSN_>
void set_00 ( real_t  a = 1.0)

Set coefficient for term (0,0): 0th order in time and space.

Parameters
[in]aConstant coefficient to multiply by 0-th order term [Default: 1.]

◆ set_01() [1/4]

template<size_t NEN_, size_t NSN_>
void set_01 ( Fct &  f)

Set coefficient for term (0,1): 0th order in time, 1st order in time and space.

Parameters
[in]fFunction to multiply by (0,1)-order term (Function of x and t)

◆ set_01() [2/4]

template<size_t NEN_, size_t NSN_>
void set_01 ( Point< real_t > &  a)

Set coefficient for term (0,1): 0th order in time, 1st order in space.

Parameters
[in]aConstant coefficient to multiply by (0,1)-order term [Default: 1.]

◆ set_01() [3/4]

template<size_t NEN_, size_t NSN_>
void set_01 ( real_t  a = 1.0)

Set coefficient for term (0,1): 0th order in time, 1st order in space.

Parameters
[in]aConstant coefficient to multiply by (0,1)-order term [Default: 1.]

◆ set_01() [4/4]

template<size_t NEN_, size_t NSN_>
void set_01 ( Vect< real_t > &  a)

Set coefficient for term (0,1): 0th order in time, 1st order in space.

Parameters
[in]aConstant coefficient to multiply by (0,1)-order term [Default: 1.]

◆ set_02() [1/2]

template<size_t NEN_, size_t NSN_>
void set_02 ( Fct &  f)

Set coefficient for term (0,2): 0th order in time, 2nd order in time and space.

Parameters
[in]fFunction to multiply by (0,2)-order term (Function of x and t)

◆ set_02() [2/2]

template<size_t NEN_, size_t NSN_>
void set_02 ( real_t  a = 1.0)

Set coefficient for term (0,2): 0th order in time, 2nd order in space.

Parameters
[in]aConstant coefficient to multiply by (0,2)-order term [Default: 1.]

◆ set_10() [1/2]

template<size_t NEN_, size_t NSN_>
void set_10 ( Fct &  f)

Set coefficient for term (1,0): 1st order in time, 0th order in space.

Parameters
[in]fFunction to multiply by (1,0)-order term (Function of x and t )

◆ set_10() [2/2]

template<size_t NEN_, size_t NSN_>
void set_10 ( real_t  a = 1.0)

Set coefficient for term (1,0): 1st order in time, 0th order in space.

Parameters
[in]aConstant coefficient to multiply by (1,0)-order term [Default: 1.]

◆ set_20() [1/2]

template<size_t NEN_, size_t NSN_>
void set_20 ( Fct &  f)

Set coefficient for term (2,0): 2nd order in time, 0th order in space.

Parameters
[in]fFunction to multiply by (2,0)-order term (Function of x and t)

◆ set_20() [2/2]

template<size_t NEN_, size_t NSN_>
void set_20 ( real_t  a = 1.0)

Set coefficient for term (2,0): 2nd order in time, 0th order in space.

Parameters
[in]aConstant coefficient to multiply by (2,0)-order term [Default: 1.]

◆ setNoLumping()

template<size_t NEN_, size_t NSN_>
void setNoLumping ( )

Set no lumping.

Remarks
Default is lumping

◆ setStab()

template<size_t NEN_, size_t NSN_>
void setStab ( )

Set stabilization for convection term.

Remarks
Default is no stabilization