LinearPDE2D Class Reference

Solves a generic linear PDE in 2-D using 3-Node triangular finite elements. More...

#include <LinearPDE2D.h>

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

Public Member Functions

 LinearPDE2D ()
 Default Constructor. Constructs an empty equation.
 
 LinearPDE2D (Mesh &ms)
 Constructor using Mesh data.
 
 LinearPDE2D (Mesh &ms, Vect< real_t > &u)
 Constructor using Mesh and initial condition.
 
 ~LinearPDE2D ()
 Destructor.
 
void Mat_00 (real_t coef=1.0)
 Add 0th order term, in time and space, to left-hand side coef
 
void Mat_10 (real_t coef=1.0)
 Add 1st order term in time, 0th in space to left-hand side.
 
void Mat_01 (real_t coef=1.0)
 Add Oth order term in time, 1st in space to left-hand side.
 
void Mat_20 (real_t coef=1.0)
 Add 2nd order term in time, 0th in space to left-hand side.
 
void Mat_02 (real_t coef=1.0)
 Add 0th order term in time, 2nd in space to left-hand side.
 
void BodyRHS (const Vect< real_t > &f)
 Add body right-hand side term to right hand side.
 
void BodyRHS (real_t f)
 Add body right-hand side term to right hand side.
 
void BoundaryRHS (real_t flux)
 Add boundary right-hand side flux to right hand side.
 
void BoundaryRHS (const Vect< real_t > &f)
 Add boundary right-hand side term to right hand side after multiplying it by coefficient coef
 
Point< real_t > & Flux () const
 Return (constant) heat flux in element.
 
void Grad (Vect< Point< real_t > > &g)
 Compute gradient of solution.
 
Point< real_t > & Grad (const Vect< real_t > &u) const
 Return gradient of a vector in element.
 
void setInput (EType opt, Vect< real_t > &u)
 Set equation input data.
 
- Public Member Functions inherited from Equa_LinearPDE< 3, 2 >
 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.
 
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

Solves a generic linear PDE in 2-D using 3-Node triangular finite elements.

This class handles the scalar linear equation: a d^2u/dt^2 + b du/dt - div(c Grad u) + d.Grad u + eu = f in a given interval, with Dirichlet and/or Neumann boundary conditions. The coefficients a, b, c, d, e can be constants or given functions of space and time

Numerical solution uses the P1 finite element method for space discretization. For transient (time-dependent) problems, the class TimeStepping must be used. For eigenproblems, the class EigenProblemSolver must be used.

Constructor & Destructor Documentation

◆ LinearPDE2D() [1/2]

LinearPDE2D ( Mesh ms)

Constructor using Mesh data.

Parameters
[in]msMesh instance

◆ LinearPDE2D() [2/2]

LinearPDE2D ( Mesh ms,
Vect< real_t > &  u 
)

Constructor using Mesh and initial condition.

Parameters
[in]msMesh instance
[in]uVect instance containing initial solution

Member Function Documentation

◆ BodyRHS() [1/2]

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 from Equa_LinearPDE< 3, 2 >.

◆ BodyRHS() [2/2]

void BodyRHS ( real_t  f)

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

Case where the body right-hand side is piecewise constant.

Parameters
[in]fValue of thermal source (Constant in element).

◆ BoundaryRHS() [1/2]

void BoundaryRHS ( const Vect< real_t > &  f)
virtual

Add boundary right-hand side term to right hand side after multiplying it by coefficient coef

Parameters
[in]fVector containing source at nodes

Reimplemented from Equa_LinearPDE< 3, 2 >.

◆ BoundaryRHS() [2/2]

void BoundaryRHS ( real_t  flux)

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

Parameters
[in]fluxVector containing source at side nodes.

◆ Grad() [1/2]

Point< real_t > & Grad ( const Vect< real_t > &  u) const

Return gradient of a vector in element.

Parameters
[in]uGlobal vector for which gradient is computed. Vector u has as size the total number of nodes

◆ Grad() [2/2]

void Grad ( Vect< Point< real_t > > &  g)

Compute gradient of solution.

Parameters
[in]gElementwise vector containing gradient of solution.

◆ Mat_00()

void Mat_00 ( real_t  coef = 1.0)
virtual

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

Parameters
[in]coefCoefficient to multiply by added term [default: 1]

Reimplemented from Equa_LinearPDE< 3, 2 >.

◆ Mat_01()

void Mat_01 ( real_t  coef = 1.0)
virtual

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

Parameters
[in]coefCoefficient to multiply by added term [default: 1]

Reimplemented from Equa_LinearPDE< 3, 2 >.

◆ Mat_02()

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 added term [default: 1]

Reimplemented from Equa_LinearPDE< 3, 2 >.

◆ Mat_10()

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 added term [default: 1]

Reimplemented from Equa_LinearPDE< 3, 2 >.

◆ Mat_20()

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 added term [default: 1]

Reimplemented from Equa_LinearPDE< 3, 2 >.

◆ setInput()

void setInput ( EType  opt,
Vect< real_t > &  u 
)

Set equation input data.

Parameters
[in]optParameter to select type of input (enumerated values)
  • INITIAL_FIELD: Initial temperature
  • BOUNDARY_CONDITION_DATA: Boundary condition (Dirichlet)
  • SOURCE_DATA: Heat source
  • FLUX_DATA: Heat flux (Neumann boundary condition)
[in]uVector containing input data