Builds finite element arrays for thermal diffusion and convection in 1-D using 2-Node elements. More...

#include <DC1DL2.h>

Inheritance diagram for DC1DL2:
Equa_Therm< 2, 2, 1, 1 > Equation< NEN_, NEE_, NSN_, NSE_ > Equa

Public Member Functions

 DC1DL2 ()
 Default Constructor.
 
 DC1DL2 (Mesh &ms)
 
 DC1DL2 (Mesh &ms, Vect< real_t > &u)
 
 ~DC1DL2 ()
 Destructor.
 
void LCapacity (real_t coef=1)
 Add lumped capacity matrix to element matrix after multiplying it by coefficient coef
 
void Capacity (real_t coef=1)
 Add Consistent capacity matrix to element matrix after multiplying it by coefficient coef.
 
void Diffusion (real_t coef=1)
 Add diffusion matrix to element matrix after multiplying it by coefficient coef
 
void Convection (const real_t &v, real_t coef=1)
 Add convection matrix to element matrix after multiplying it by coefficient coef
 
void Convection (const Vect< real_t > &v, real_t coef=1)
 Add convection matrix to element matrix after multiplying it by coefficient coef
 
void Convection (real_t coef=1)
 Add convection matrix to element matrix after multiplying it by coefficient coef
 
void BodyRHS (const Vect< real_t > &f)
 Add body right-hand side term to right hand side.
 
real_t Flux () const
 Return (constant) heat flux in element.
 
void setInput (EType opt, Vect< real_t > &u)
 Set equation input data.
 
- Public Member Functions inherited from Equa_Therm< 2, 2, 1, 1 >
 Equa_Therm ()
 Default constructor.
 
virtual ~Equa_Therm ()
 Destructor.
 
virtual void setStab ()
 Set stabilized formulation.
 
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.
 
void build (TimeStepping &s)
 Build the linear system of equations.
 
void build (EigenProblemSolver &e)
 Build the linear system for an eigenvalue problem.
 
void setRho (const real_t &rho)
 Set Density (constant)
 
void setCp (const real_t &cp)
 Set Specific heat (constant)
 
void setConductivity (const real_t &diff)
 Set (constant) thermal conductivity.
 
- Public Member Functions inherited from Equation< NEN_, NEE_, NSN_, NSE_ >
 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 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 ElementNodeVector (const Vect< real_t > &b, LocalVect< real_t, NEN_ > &be, int dof)
 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 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 ElementAssembly (Vect< real_t > &v)
 Assemble element vector 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.
 

Additional Inherited Members

- Protected Member Functions inherited from Equa_Therm< 2, 2, 1, 1 >
void setMaterial ()
 Set material properties.
 

Detailed Description

Builds finite element arrays for thermal diffusion and convection in 1-D using 2-Node elements.

Note that members calculating element arrays have as an argument a real coef that will be multiplied by the contribution of the current element. This makes possible testing different algorithms.

Constructor & Destructor Documentation

◆ DC1DL2() [1/3]

DC1DL2 ( )

Default Constructor.

Constructs an empty equation.

◆ DC1DL2() [2/3]

DC1DL2 ( Mesh ms)

Constructor using mesh instance

Parameters
[in]msMesh instance

◆ DC1DL2() [3/3]

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

Constructor using mesh instance and solution vector

Parameters
[in]msMesh instance
[in,out]uVect instance containing solution vector

Member Function Documentation

◆ BodyRHS()

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_Therm< 2, 2, 1, 1 >.

◆ Capacity()

void Capacity ( real_t  coef = 1)
virtual

Add Consistent capacity matrix to element matrix after multiplying it by coefficient coef.

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

Reimplemented from Equa_Therm< 2, 2, 1, 1 >.

◆ Convection() [1/3]

void Convection ( const real_t &  v,
real_t  coef = 1 
)

Add convection matrix to element matrix after multiplying it by coefficient coef

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

◆ Convection() [2/3]

void Convection ( const Vect< real_t > &  v,
real_t  coef = 1 
)

Add convection matrix to element matrix after multiplying it by coefficient coef

Case where velocity field is given by a vector v

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

◆ Convection() [3/3]

void Convection ( real_t  coef = 1)
virtual

Add convection matrix to element matrix after multiplying it by coefficient coef

Case where velocity field has been previouly defined

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

Reimplemented from Equa_Therm< 2, 2, 1, 1 >.

◆ Diffusion()

void Diffusion ( real_t  coef = 1)
virtual

Add diffusion matrix to element matrix after multiplying it by coefficient coef

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

Reimplemented from Equa_Therm< 2, 2, 1, 1 >.

◆ LCapacity()

void LCapacity ( real_t  coef = 1)
virtual

Add lumped capacity matrix to element matrix after multiplying it by coefficient coef

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

Reimplemented from Equa_Therm< 2, 2, 1, 1 >.

◆ setInput()

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

Set equation input data.

Parameters
[in]optParameter that selects data type for input. This parameter is to be chosen in the enumerated variable EqDataType
  • INITIAL_FIELD: Initial temperature
  • BOUNDARY_CONDITION_DATA: Boundary condition (Dirichlet)
  • SOURCE_DATA: Heat source
  • FLUX_DATA: Heat flux (Neumann boundary condition)
  • VELOCITY: Velocity vector (for the convection term)
[in]uVector containing input data