Builds finite element arrays for thermal diffusion and convection in 2-D domains using 3-Node triangles. More...
Public Member Functions | |
DC2DT3 () | |
Default Constructor. Constructs an empty equation. | |
DC2DT3 (Mesh &ms) | |
Constructor using Mesh data. More... | |
DC2DT3 (Mesh &ms, Vect< real_t > &u) | |
Constructor using Mesh and initial condition. More... | |
~DC2DT3 () | |
Destructor. | |
void | LCapacity (real_t coef=1) |
Add lumped capacity matrix to element matrix after multiplying it by coefficient coef More... | |
void | Capacity (real_t coef=1) |
Add Consistent capacity matrix to element matrix after multiplying it by coefficient coef More... | |
void | Diffusion (real_t coef=1) |
Add diffusion matrix to element matrix after multiplying it by coefficient coef More... | |
void | Diffusion (const LocalMatrix< real_t, 2, 2 > &diff, real_t coef=1) |
Add diffusion matrix to element matrix after multiplying it by coefficient coef More... | |
void | Convection (const Point< real_t > &v, real_t coef=1) |
Add convection matrix to element matrix after multiplying it by coefficient coef More... | |
void | Convection (const Vect< real_t > &v, real_t coef=1) |
Add convection matrix to element matrix after multiplying it by coefficient coef More... | |
void | Convection (real_t coef=1) |
Add convection matrix to element matrix after multiplying it by coefficient coef More... | |
void | LinearExchange (real_t coef, real_t T) |
Add an edge linear exchange term to left and right-hand sides. More... | |
void | BodyRHS (const Vect< real_t > &f) |
Add body right-hand side term to right hand side. More... | |
void | BodyRHS (real_t f) |
Add body right-hand side term to right hand side. More... | |
void | BoundaryRHS (real_t flux) |
Add boundary right-hand side flux to right hand side. More... | |
void | BoundaryRHS (const Vect< real_t > &f) |
Add boundary right-hand side term to right hand side after multiplying it by coefficient coef More... | |
void | Periodic (real_t coef=1.e20) |
Add contribution of periodic boundary condition (by a penalty technique). More... | |
Point< real_t > & | Flux () const |
Return (constant) heat flux in element. | |
void | Grad (Vect< Point< real_t > > &g) |
Compute gradient of solution. More... | |
Point< real_t > & | Grad (const Vect< real_t > &u) const |
Return gradient of a vector in element. More... | |
void | setInput (EType opt, Vect< real_t > &u) |
Set equation input data. More... | |
void | JouleHeating (const Vect< real_t > &sigma, const Vect< real_t > &psi) |
Set Joule heating term as source. More... | |
Public Member Functions inherited from Equa_Therm< 3, 3, 2, 2 > | |
Equa_Therm () | |
Default constructor. More... | |
virtual | ~Equa_Therm () |
Destructor. | |
virtual void | setStab () |
Set stabilized formulation. More... | |
void | build () |
Build the linear system of equations. More... | |
void | build (TimeStepping &s) |
Build the linear system of equations. More... | |
void | build (EigenProblemSolver &e) |
Build the linear system for an eigenvalue problem. More... | |
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. 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. | |
Mesh & | getMesh () const |
Return reference to Mesh instance. More... | |
LinearSolver & | getLinearSolver () |
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. | |
Additional Inherited Members | |
Protected Member Functions inherited from Equa_Therm< 3, 3, 2, 2 > | |
void | setMaterial () |
Set material properties. | |
Builds finite element arrays for thermal diffusion and convection in 2-D domains using 3-Node triangles.
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.
Add body right-hand side term to right hand side.
[in] | f | Vector containing source at nodes. |
Reimplemented from Equa_Therm< 3, 3, 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.
[in] | f | Value of thermal source (Constant in element). |
void BoundaryRHS | ( | real_t | flux | ) |
Add boundary right-hand side flux to right hand side.
[in] | flux | Vector containing source at side nodes. |
Add boundary right-hand side term to right hand side after multiplying it by coefficient coef
[in] | f | Vector containing source at nodes |
Reimplemented from Equa_Therm< 3, 3, 2, 2 >.
|
virtual |
Add Consistent capacity matrix to element matrix after multiplying it by coefficient coef
[in] | coef | Coefficient to multiply by added term [Default: 1 ] |
Reimplemented from Equa_Therm< 3, 3, 2, 2 >.
Add convection matrix to element matrix after multiplying it by coefficient coef
[in] | v | Constant velocity vector |
[in] | coef | Coefficient to multiply by added term [Default: 1 ] |
Add convection matrix to element matrix after multiplying it by coefficient coef
Case where velocity field is given by a vector v
[in] | v | Velocity vector |
[in] | coef | Coefficient to multiply by added term (Default: 1 ] |
|
virtual |
Add convection matrix to element matrix after multiplying it by coefficient coef
Case where velocity field has been previouly defined
[in] | coef | Coefficient to multiply by added term [Default: 1 ] |
Reimplemented from Equa_Therm< 3, 3, 2, 2 >.
|
virtual |
Add diffusion matrix to element matrix after multiplying it by coefficient coef
[in] | coef | Coefficient to multiply by added term [Default: 1 ] |
Reimplemented from Equa_Therm< 3, 3, 2, 2 >.
void Diffusion | ( | const LocalMatrix< real_t, 2, 2 > & | diff, |
real_t | coef = 1 |
||
) |
Add diffusion matrix to element matrix after multiplying it by coefficient coef
Case where the diffusivity matrix is given as an argument.
[in] | diff | Diffusion matrix (class LocalMatrix). |
[in] | coef | Coefficient to multiply by added term [Default: 1 ] |
Compute gradient of solution.
[in] | g | Elementwise vector containing gradient of solution. |
Return gradient of a vector in element.
[in] | u | Global vector for which gradient is computed. Vector u has as size the total number of nodes |
|
virtual |
Add lumped capacity matrix to element matrix after multiplying it by coefficient coef
[in] | coef | Coefficient to multiply by added term [Default: 1 ]. |
Reimplemented from Equa_Therm< 3, 3, 2, 2 >.
Add an edge linear exchange term to left and right-hand sides.
[in] | coef | Coefficient of exchange |
[in] | T | External value for exchange |
T
void Periodic | ( | real_t | coef = 1.e20 | ) |
Add contribution of periodic boundary condition (by a penalty technique).
Boundary nodes where periodic boundary conditions are to be imposed must have codes equal to PERIODIC_A
on one side and PERIODIC_B
on the opposite side.
[in] | coef | Value of penalty parameter [Default: 1.e20 ] |
Set equation input data.
[in] | opt | Parameter to select type of input (enumerated values)
|
[in] | u | Vector containing input data |