(9)

Equation< NEN_, NEE_, NSN_, NSE_ > Class Template Reference

Abstract class for all equation classes. More...

Inheritance diagram for Equation< NEN_, NEE_, NSN_, NSE_ >:
Equa Equa_Acoustics< 2, 2, 1, 1 > Equa_Acoustics< 3, 3, 2, 2 > Equa_Electromagnetics< 3, 6, 2, 4 > Equa_Fluid< 3, 6, 2, 4 > Equa_Fluid< 4, 12, 3, 9 > Equa_Fluid< 4, 8, 2, 4 > Equa_Laplace< 2, 2, 1, 1 > Equa_Laplace< 3, 3, 1, 1 > Equa_Laplace< 3, 3, 2, 2 > Equa_Laplace< 4, 4, 3, 3 > Equa_Laplace< 6, 6, 3, 3 > Equa_Porous< 2, 2, 1, 1 > Equa_Porous< 3, 3, 2, 2 > Equa_Solid< 2, 12, 1, 6 > Equa_Solid< 2, 4, 1, 2 > Equa_Solid< 3, 6, 2, 4 > Equa_Solid< 4, 12, 3, 9 > Equa_Solid< 4, 8, 2, 4 > Equa_Solid< 8, 24, 4, 12 > Equa_Therm< 2, 2, 1, 1 > Equa_Therm< 3, 3, 2, 2 > Equa_Therm< 4, 4, 3, 3 > Equa_Therm< 6, 6, 3, 3 > Equa_Acoustics< NEN_, NEE_, NSN_, NSE_ > Equa_Electromagnetics< NEN_, NEE_, NSN_, NSE_ > Equa_Fluid< NEN_, NEE_, NSN_, NSE_ > Equa_Laplace< NEN_, NEE_, NSN_, NSE_ > Equa_Porous< NEN_, NEE_, NSN_, NSE_ > Equa_Solid< NEN_, NEE_, NSN_, NSE_ > Equa_Therm< NEN_, NEE_, NSN_, NSE_ >

Public Member Functions

 Equation ()
 
 Equation (Mesh &mesh)
 Constructor with mesh instance. More...
 
 Equation (Mesh &mesh, Vect< real_t > &u)
 Constructor with mesh instance and solution vector. More...
 
 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. More...
 
 ~Equation ()
 Destructor.
 
void updateBC (const Element &el, const Vect< real_t > &bc)
 Update Right-Hand side by taking into account essential boundary conditions. More...
 
void DiagBC (DOFSupport dof_type=NODE_DOF, int dof=0)
 Update element matrix to impose bc by diagonalization technique. More...
 
void LocalNodeVector (Vect< real_t > &b)
 Localize element vector from a Vect instance. More...
 
void ElementNodeVector (const Vect< real_t > &b, LocalVect< real_t, NEE_ > &be)
 Localize element vector from a Vect instance. More...
 
void SideNodeVector (const Vect< real_t > &b, LocalVect< real_t, NSE_ > &bs)
 Localize side vector from a Vect instance. More...
 
void SideSideVector (const Vect< real_t > &b, vector< real_t > &bs)
 Localize side vector from a Vect instance. More...
 
void ElementNodeVectorSingleDOF (const Vect< real_t > &b, LocalVect< real_t, NEN_ > &be)
 Localize Element Vector from a Vect instance. More...
 
void ElementNodeVector (const Vect< real_t > &b, LocalVect< real_t, NEN_ > &be, int dof)
 Localize Element Vector from a Vect instance. More...
 
void ElementSideVector (const Vect< real_t > &b, LocalVect< real_t, NSE_ > &be)
 Localize Element Vector from a Vect instance. More...
 
void ElementVector (const Vect< real_t > &b, DOFSupport dof_type=NODE_DOF, int flag=0)
 Localize element vector. More...
 
void SideVector (const Vect< real_t > &b, vector< real_t > &sb)
 Localize side vector. More...
 
void ElementNodeCoordinates ()
 Localize coordinates of element nodes. More...
 
void SideNodeCoordinates ()
 Localize coordinates of side nodes. More...
 
void ElementAssembly (Matrix< real_t > *A)
 Assemble element matrix into global one. More...
 
void ElementAssembly (BMatrix< real_t > &A)
 Assemble element matrix into global one. More...
 
void ElementAssembly (SkSMatrix< real_t > &A)
 Assemble element matrix into global one. More...
 
void ElementAssembly (SkMatrix< real_t > &A)
 Assemble element matrix into global one. More...
 
void ElementAssembly (SpMatrix< real_t > &A)
 Assemble element matrix into global one. More...
 
void ElementAssembly (TrMatrix< real_t > &A)
 Assemble element matrix into global one. More...
 
void DGElementAssembly (Matrix< real_t > *A)
 Assemble element matrix into global one for the Discontinuous Galerkin approximation. More...
 
void DGElementAssembly (SkSMatrix< real_t > &A)
 Assemble element matrix into global one for the Discontinuous Galerkin approximation. More...
 
void DGElementAssembly (SkMatrix< real_t > &A)
 Assemble element matrix into global one for the Discontinuous Galerkin approximation. More...
 
void DGElementAssembly (SpMatrix< real_t > &A)
 Assemble element matrix into global one for the Discontinuous Galerkin approximation. More...
 
void DGElementAssembly (TrMatrix< real_t > &A)
 Assemble element matrix into global one for the Discontinuous Galerkin approximation. More...
 
void SideAssembly (Matrix< real_t > *A)
 Assemble side (edge or face) matrix into global one. More...
 
void SideAssembly (SkSMatrix< real_t > &A)
 Assemble side (edge or face) matrix into global one. More...
 
void SideAssembly (SkMatrix< real_t > &A)
 Assemble side (edge or face) matrix into global one. More...
 
void SideAssembly (SpMatrix< real_t > &A)
 Assemble side (edge or face) matrix into global one. More...
 
void ElementAssembly (Vect< real_t > &v)
 Assemble element vector into global one. More...
 
void SideAssembly (Vect< real_t > &v)
 Assemble side (edge or face) vector into global one. More...
 
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. More...
 
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. More...
 
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. More...
 
- 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. More...
 
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. More...
 
void setMatrixType (int t)
 Choose type of matrix. More...
 
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. More...
 
int solveLinearSystem (Vect< real_t > &b, Vect< real_t > &x)
 Solve the linear system with given right-hand side. More...
 
void LinearSystemInfo ()
 Print info on linear system solver.
 

Detailed Description

template<size_t NEN_, size_t NEE_, size_t NSN_, size_t NSE_>
class OFELI::Equation< NEN_, NEE_, NSN_, NSE_ >

Abstract class for all equation classes.

Template Arguments:

  • NEN_ : Number of element nodes
  • NEE_ : Number of element equations
  • NSN_ : Number of side nodes
  • NSN_ : Number of side equations
Author
Rachid Touzani

Constructor & Destructor Documentation

◆ Equation() [1/4]

Equation ( )

Default constructor. Constructs an "empty" equation

◆ Equation() [2/4]

Equation ( Mesh mesh)

Constructor with mesh instance.

Parameters
[in]meshMesh instance

◆ Equation() [3/4]

Equation ( Mesh mesh,
Vect< real_t > &  u 
)

Constructor with mesh instance and solution vector.

Parameters
[in]meshMesh instance
[in]uVect instance containing solution.

◆ Equation() [4/4]

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.

Parameters
[in]meshMesh instance
[in]uVect instance containing Right-hand side.
[in]init_timeInitial Time value
[in]final_timeFinal Time value
[in]time_stepTime step value

Member Function Documentation

◆ AxbAssembly() [1/2]

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.

Parameters
[in]elReference to Element instance
[in]xGlobal vector to multiply by (Vect instance)
[out]bGlobal vector to add (Vect instance)

◆ AxbAssembly() [2/2]

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.

Parameters
[in]sdReference to Side instance
[in]xGlobal vector to multiply by (Vect instance)
[out]bGlobal vector (Vect instance)

◆ DGElementAssembly() [1/5]

void DGElementAssembly ( Matrix< real_t > *  A)

Assemble element matrix into global one for the Discontinuous Galerkin approximation.

Parameters
APointer to global matrix (abstract class: can be any of classes SkSMatrix, SkMatrix, SpMatrix)
Warning
The element pointer is given by the global variable theElement

◆ DGElementAssembly() [2/5]

void DGElementAssembly ( SkSMatrix< real_t > &  A)

Assemble element matrix into global one for the Discontinuous Galerkin approximation.

Parameters
AGlobal matrix stored as an SkSMatrix instance
Warning
The element pointer is given by the global variable theElement

◆ DGElementAssembly() [3/5]

void DGElementAssembly ( SkMatrix< real_t > &  A)

Assemble element matrix into global one for the Discontinuous Galerkin approximation.

Parameters
[in]AGlobal matrix stored as an SkMatrix instance
Warning
The element pointer is given by the global variable theElement

◆ DGElementAssembly() [4/5]

void DGElementAssembly ( SpMatrix< real_t > &  A)

Assemble element matrix into global one for the Discontinuous Galerkin approximation.

Parameters
[in]AGlobal matrix stored as an SpMatrix instance
Warning
The element pointer is given by the global variable theElement

◆ DGElementAssembly() [5/5]

void DGElementAssembly ( TrMatrix< real_t > &  A)

Assemble element matrix into global one for the Discontinuous Galerkin approximation.

Parameters
[in]AGlobal matrix stored as an TrMatrix instance
Warning
The element pointer is given by the global variable theElement

◆ DiagBC()

void DiagBC ( DOFSupport  dof_type = NODE_DOF,
int  dof = 0 
)

Update element matrix to impose bc by diagonalization technique.

Parameters
[in]dof_typeDOF type option. To choose among the enumerated values:
  • NODE_DOF, DOFs are supported by nodes [Default]
  • ELEMENT_DOF, DOFs are supported by elements
  • SIDE_DOF, DOFs are supported by sides
[in]dofDOF setting:
  • = 0, All DOFs are taken into account [Default]
  • != 0, Only DOF No. dof is handled in the system

◆ ElementAssembly() [1/7]

void ElementAssembly ( Matrix< real_t > *  A)

Assemble element matrix into global one.

Parameters
APointer to global matrix (abstract class: can be any of classes SkSMatrix, SkMatrix, SpMatrix)
Warning
The element pointer is given by the global variable theElement

◆ ElementAssembly() [2/7]

void ElementAssembly ( BMatrix< real_t > &  A)

Assemble element matrix into global one.

Parameters
AGlobal matrix stored as a BMatrix instance
Warning
The element pointer is given by the global variable theElement

◆ ElementAssembly() [3/7]

void ElementAssembly ( SkSMatrix< real_t > &  A)

Assemble element matrix into global one.

Parameters
AGlobal matrix stored as an SkSMatrix instance
Warning
The element pointer is given by the global variable theElement

◆ ElementAssembly() [4/7]

void ElementAssembly ( SkMatrix< real_t > &  A)

Assemble element matrix into global one.

Parameters
[in]AGlobal matrix stored as an SkMatrix instance
Warning
The element pointer is given by the global variable theElement

◆ ElementAssembly() [5/7]

void ElementAssembly ( SpMatrix< real_t > &  A)

Assemble element matrix into global one.

Parameters
[in]AGlobal matrix stored as an SpMatrix instance
Warning
The element pointer is given by the global variable theElement

◆ ElementAssembly() [6/7]

void ElementAssembly ( TrMatrix< real_t > &  A)

Assemble element matrix into global one.

Parameters
[in]AGlobal matrix stored as an TrMatrix instance
Warning
The element pointer is given by the global variable theElement

◆ ElementAssembly() [7/7]

void ElementAssembly ( Vect< real_t > &  v)

Assemble element vector into global one.

Parameters
[in]vGlobal vector (Vect instance)
Warning
The element pointer is given by the global variable theElement

◆ ElementNodeCoordinates()

void ElementNodeCoordinates ( )

Localize coordinates of element nodes.

Coordinates are stored in array _x[0], _x[1], ... which are instances of class Point<real_t>

Remarks
This member function uses the Side pointer _theSide

◆ ElementNodeVector() [1/2]

void ElementNodeVector ( const Vect< real_t > &  b,
LocalVect< real_t, NEE_ > &  be 
)

Localize element vector from a Vect instance.

Parameters
[in]bGlobal vector to be localized.
[out]beLocal vector, the length of which is the total number of element equations.
Remarks
All degrees of freedom are transferred to the local vector

◆ ElementNodeVector() [2/2]

void ElementNodeVector ( const Vect< real_t > &  b,
LocalVect< real_t, NEN_ > &  be,
int  dof 
)

Localize Element Vector from a Vect instance.

Parameters
[in]bGlobal vector to be localized.
[out]beLocal vector, the length of which is the total number of element equations.
[in]dofDegree of freedom to transfer to the local vector
Remarks
Only yhe dega dof is transferred to the local vector

◆ ElementNodeVectorSingleDOF()

void ElementNodeVectorSingleDOF ( const Vect< real_t > &  b,
LocalVect< real_t, NEN_ > &  be 
)

Localize Element Vector from a Vect instance.

Parameters
[in]bGlobal vector to be localized.
[out]beLocal vector, the length of which is the total number of element equations.
Remarks
Vector b is assumed to contain only one degree of freedom by node.

◆ ElementSideVector()

void ElementSideVector ( const Vect< real_t > &  b,
LocalVect< real_t, NSE_ > &  be 
)

Localize Element Vector from a Vect instance.

Parameters
[in]bGlobal vector to be localized.
[out]beLocal vector, the length of which is

◆ ElementVector()

void ElementVector ( const Vect< real_t > &  b,
DOFSupport  dof_type = NODE_DOF,
int  flag = 0 
)

Localize element vector.

Parameters
[in]bGlobal vector to be localized
[in]dof_typeDOF type option. To choose among the enumerated values:
  • NODE_DOF, DOFs are supported by nodes [Default]
  • ELEMENT_DOF, DOFs are supported by elements
  • SIDE_DOF, DOFs are supported by sides
[in]flagOption to set:
  • = 0, All DOFs are taken into account [Default]
  • != 0, Only DOF number dof is handled in the system
The resulting local vector can be accessed by attribute ePrev.
Remarks
This member function is to be used if a constructor with Element was invoked. It uses the Element pointer _theElement

◆ LocalNodeVector()

void LocalNodeVector ( Vect< real_t > &  b)

Localize element vector from a Vect instance.

Parameters
[in]bReference to global vector to be localized. The resulting local vector can be accessed by attribute ePrev. This member function is to be used if a constructor with Element was invoked.

◆ setMaterialProperty()

real_t setMaterialProperty ( const string &  exp,
const string &  prop 
)

Define a material property by an algebraic expression.

Parameters
[in]expAlgebraic expression
[in]propProperty name
Returns
Return value in expression evaluation:
  • =0, Normal evaluation
  • !=0, An error message is displayed

◆ SideAssembly() [1/5]

void SideAssembly ( Matrix< real_t > *  A)

Assemble side (edge or face) matrix into global one.

Parameters
APointer to global matrix (abstract class: can be any of classes SkSMatrix, SkMatrix, SpMatrix)
Warning
The side pointer is given by the global variable theSide

◆ SideAssembly() [2/5]

void SideAssembly ( SkSMatrix< real_t > &  A)

Assemble side (edge or face) matrix into global one.

Parameters
[in]AGlobal matrix stored as an SkSMatrix instance
Warning
The side pointer is given by the global variable theSide

◆ SideAssembly() [3/5]

void SideAssembly ( SkMatrix< real_t > &  A)

Assemble side (edge or face) matrix into global one.

Parameters
[in]AGlobal matrix stored as an SkMatrix instance
Warning
The side pointer is given by the global variable theSide

◆ SideAssembly() [4/5]

void SideAssembly ( SpMatrix< real_t > &  A)

Assemble side (edge or face) matrix into global one.

Parameters
[in]AGlobal matrix stored as an SpMatrix instance
Warning
The side pointer is given by the global variable theSide

◆ SideAssembly() [5/5]

void SideAssembly ( Vect< real_t > &  v)

Assemble side (edge or face) vector into global one.

Parameters
[in]vGlobal vector (Vect instance)
Warning
The side pointer is given by the global variable theSide

◆ SideNodeCoordinates()

void SideNodeCoordinates ( )

Localize coordinates of side nodes.

Coordinates are stored in array _x[0], _x[1], ... which are instances of class Point<real_t>

Remarks
This member function uses the Element pointer _theElement

◆ SideNodeVector()

void SideNodeVector ( const Vect< real_t > &  b,
LocalVect< real_t, NSE_ > &  bs 
)

Localize side vector from a Vect instance.

Parameters
[in]bGlobal vector to be localized.
[out]bsLocal vector, the length of which is the total number of side equations.
Remarks
All degrees of freedom are transferred to the local vector

◆ SideSideVector()

void SideSideVector ( const Vect< real_t > &  b,
vector< real_t > &  bs 
)

Localize side vector from a Vect instance.

Parameters
[in]bGlobal vector to be localized.
[out]bsLocal constant value of vector at given side.
Remarks
All degrees of freedom are transferred to the local vector

◆ SideVector()

void SideVector ( const Vect< real_t > &  b,
vector< real_t > &  sb 
)

Localize side vector.

Parameters
[in]bGlobal vector to be localized
  • NODE_DOF, DOFs are supported by nodes [ default ]
  • ELEMENT_DOF, DOFs are supported by elements
  • SIDE_DOF, DOFs are supported by sides
  • BOUNDARY_SIDE_DOF, DOFs are supported by boundary sides
[out]sbArray in which local vector is stored The resulting local vector can be accessed by attribute ePrev.
Remarks
This member function is to be used if a constructor with Side was invoked. It uses the Side pointer _theSide

◆ updateBC()

void updateBC ( const Element el,
const Vect< real_t > &  bc 
)

Update Right-Hand side by taking into account essential boundary conditions.

Parameters
[in]elReference to current element instance
[in]bcVector that contains imposed values at all DOFs