Builds finite element arrays for thermal diffusion and convection in 3-D domains using 4-Node tetrahedra. More...


Public Member Functions | |
DC3DT4 () | |
Default Constructor. More... | |
DC3DT4 (const Element *el) | |
Constructor for an element. More... | |
DC3DT4 (const Side *sd) | |
Constructor for a boundary side. More... | |
DC3DT4 (const Element *el, const Vect< real_t > &u, real_t time=0.) | |
Constructor for an element (transient case). More... | |
DC3DT4 (const Side *sd, const Vect< real_t > &u, real_t time=0.) | |
Constructor for a boundary side (transient case). More... | |
DC3DT4 (const Element *el, const Vect< real_t > &u, real_t time, real_t deltat, int scheme) | |
Constructor for an element (transient case) with specification of time integration scheme. More... | |
DC3DT4 (const Side *sd, const Vect< real_t > &u, real_t time, real_t deltat, int scheme) | |
Constructor for a side (transient case) with specification of time integration scheme. More... | |
~DC3DT4 () | |
Destructor. | |
void | build () |
Build the linear system without solving. | |
void | LCapacity (real_t coef=1.) |
Add lumped capacity contribution to left and right-hand sides after multiplying it by coefficient coef . More... | |
void | LCapacityToLHS (real_t coef=1) |
Add lumped capacity matrix to left-hand side after multiplying it by coefficient coef More... | |
void | LCapacityToRHS (real_t coef=1) |
Add lumped capacity contribution to right-hand side after multiplying it by coefficient coef More... | |
void | Capacity (real_t coef=1) |
Add Consistent capacity contribution to left and right-hand sides after multiplying it by coefficient coef . More... | |
void | CapacityToLHS (real_t coef=1) |
Add consistent capacity matrix to left-hand side after multiplying it by coefficient coef More... | |
void | CapacityToRHS (real_t coef=1) |
Add consistent capacity contribution to right-hand side after multiplying it by coefficient coef More... | |
void | Diffusion (real_t coef=1) |
Add diffusion matrix to left hand side after multiplying it by coefficient coef . More... | |
void | Diffusion (const DMatrix< real_t > &diff, real_t coef=1) |
Add diffusion matrix to left hand side after multiplying it by coefficient coef More... | |
void | DiffusionToRHS (real_t coef=1) |
Add diffusion contribution to right hand side after multiplying it by coefficient coef More... | |
void | Convection (real_t coef=1) |
Add convection matrix to left-hand side after multiplying it by coefficient coef More... | |
void | Convection (const Point< real_t > &v, real_t coef=1) |
Add convection matrix to left-hand side after multiplying it by coefficient coef More... | |
void | Convection (const Vect< Point< real_t > > &v, real_t coef=1) |
Add convection matrix to left-hand side after multiplying it by coefficient coef More... | |
void | RHS_Convection (const Point< real_t > &v, real_t coef=1.) |
Add convection contribution to right-hand side after multiplying it by coefficient coef More... | |
void | BodyRHS (UserData< real_t > &ud, real_t coef=1) |
Add body right-hand side term to right hand side after multiplying it by coefficient coef More... | |
void | BodyRHS (const Vect< real_t > &b, int opt=GLOBAL_ARRAY) |
Add body right-hand side term to right hand side. More... | |
void | BoundaryRHS (UserData< real_t > &ud, real_t coef=1) |
Add boundary right-hand side term to right hand side after multiplying it by coefficient coef More... | |
void | BoundaryRHS (const Vect< real_t > &b, int opt=GLOBAL_ARRAY) |
Add boundary right-hand side term to right hand side after multiplying it by coefficient coef More... | |
void | BoundaryRHS (real_t flux) |
Add boundary right-hand side flux to right hand side. More... | |
Point< real_t > | Flux () const |
Return (constant) heat flux in element. | |
Point< real_t > | Grad (const Vect< real_t > &u) const |
Return gradient of vector u in element. u is a local vector. | |
void | Periodic (real_t coef=1.e20) |
Add contribution of periodic boundary condition (by a penalty technique). More... | |
virtual void | setStab () |
Set stabilized formulation. More... | |
void | setLumpedCapacity () |
Add lumped capacity contribution to left and right-hand sides taking into account time integration scheme. | |
void | setCapacity () |
Add consistent capacity contribution to left and right-hand sides taking into account time integration scheme. | |
void | setDiffusion () |
Add diffusion contribution to left and/or right-hand side taking into account time integration scheme. | |
virtual void | ConvectionToRHS (real_t coef=1.) |
Add convection term to right-hand side. | |
void | setConvection () |
Add convection contribution to left and/or right-hand side taking into account time integration scheme. | |
int | runTransient () |
Run one time step. More... | |
int | run () |
Run the equation. More... | |
void | setRhoCp (const real_t &rhocp) |
Set product of Density by Specific heat (constants) | |
void | setConductivity (const real_t &diff) |
Set (constant) thermal conductivity. | |
void | RhoCp (const string &exp) |
Set product of Density by Specific heat given by an algebraic expression. | |
void | Conduc (const string &exp) |
Set thermal conductivity given by an algebraic expression. | |
void | buildEigen (SkSMatrix< real_t > &K, SkSMatrix< real_t > &M) |
Build global stiffness and mass matrices for the eigen system. More... | |
void | buildEigen (SkSMatrix< real_t > &K, Vect< real_t > &M) |
Build global diffusion and capacity matrices for the eigen system. More... | |
void | updateBC (const Element &el, const Vect< real_t > &bc) |
Update Right-Hand side by taking into account essential boundary conditions. More... | |
void | updateBC (const Vect< real_t > &bc) |
Update Right-Hand side by taking into account essential boundary conditions. More... | |
void | DiagBC (int 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 | ElementNodeVector (const Vect< real_t > &b, LocalVect< real_t, NEN_ > &be, int dof) |
Localize Element 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 | 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, int dof_type=NODE_FIELD, int flag=0) |
Localize Element Vector. More... | |
void | SideVector (const Vect< real_t > &b) |
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 (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 | ElementAssembly (Vect< real_t > &v) |
Assemble element vector 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 | 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 * | A () |
Return element matrix as a C-array. | |
real_t * | sA () |
Return side matrix as a C-array. | |
real_t * | b () |
Return element right-hand side as a C-array. | |
real_t * | sb () |
Return side right-hand side as a C-array. | |
real_t * | Prev () |
Return element matrix as a C-array. | |
LocalMatrix< real_t, NEE_, NEE_ > & | EA () |
Return element matrix as a LocalMatrix instance. | |
LocalMatrix< real_t, NSE_, NSE_ > & | SA () |
Return side matrix as a LocalMatrix instance. | |
LocalVect< real_t, NEE_ > & | Eb () |
Return element right-hand side as a LocalVect instance. | |
LocalVect< real_t, NEE_ > & | Ep () |
Return element matrix as a C-array. | |
void | setInitialSolution (const Vect< real_t > &u) |
Set initial solution (previous time step) | |
real_t | setMaterialProperty (const string &exp, const string &prop) |
Define a material property by an algebraic expression. More... | |
void | setMesh (class Mesh &m) |
Define mesh and renumber DOFs after removing imposed ones. | |
Mesh & | getMesh () const |
Return reference to Mesh instance. More... | |
LinearSolver< real_t > & | getLinearSolver () |
Return reference to linear solver instance. | |
void | setSolver (int ls, int pc=IDENT_PREC) |
Choose solver for the linear system. More... | |
int | solveEigenProblem (int nb_eigv, bool g=false) |
Compute eigenvalues and eigenvectors. More... | |
real_t | getEigenValue (int n) const |
Return the n-th eigenvalue. More... | |
void | getEigenVector (int n, Vect< real_t > &v) const |
Store the eigenvector corresponding to a given eigenvalue. More... | |
class Eigen & | getEigenSolver () |
Return reference to eigenproblem solver. | |
Protected Member Functions | |
void | setMaterial () |
Set material properties. | |
void | Init (const Element *el) |
Set element arrays to zero. | |
void | Init (const Side *sd) |
Set side arrays to zero. | |
Detailed Description
Builds finite element arrays for thermal diffusion and convection in 3-D domains using 4-Node tetrahedra.
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
DC3DT4 | ( | ) |
Default Constructor.
Constructs an empty equation.
Constructor for an element (transient case).
- Parameters
-
[in] el Pointer to element. [in] u Vect instance that contains solution at previous time step. [in] time Current time value [Default: 0
].
Constructor for a boundary side (transient case).
- Parameters
-
[in] sd Pointer to side. [in] u Vect instance that contains solution at previous time step. [in] time Current time value [Default: 0
].
Constructor for an element (transient case) with specification of time integration scheme.
- Parameters
-
[in] el Pointer to element. [in] u Vect instance that contains solution at previous time step. [in] time Current time value. [in] deltat Value of time step [in] scheme Time Integration Scheme: -
FORWARD_EULER
: Forward Euler scheme -
BACKWARD_EULER
: Backward Euler scheme -
CRANK_NICOLSON
: Crank-Nicolson Euler scheme
-
Constructor for a side (transient case) with specification of time integration scheme.
- Parameters
-
[in] sd Pointer to side. [in] u Vect instance that contains solution at previous time step. [in] time Current time value. [in] deltat Value of time step [in] scheme Time Integration Scheme (): - FORWARD_EULER: for Forward Euler scheme
- BACKWARD_EULER: for Backward Euler scheme
- CRANK_NICOLSON: for Crank-Nicolson Euler scheme
Member Function Documentation
void LCapacity | ( | real_t | coef = 1. | ) |
Add lumped capacity contribution to left and right-hand sides after multiplying it by coefficient coef
.
- Parameters
-
[in] coef Coefficient to multiply by added term [Default: 1
].
|
virtual |
Add lumped capacity matrix to left-hand side after multiplying it by coefficient coef
- Parameters
-
[in] coef Coefficient to multiply by added term [Default: 1
].
Reimplemented from Equa_Therm< real_t, 4, 4, 3, 3 >.
|
virtual |
Add lumped capacity contribution to right-hand side after multiplying it by coefficient coef
- Parameters
-
[in] coef Coefficient to multiply by added term [Default: 1
].
Reimplemented from Equa_Therm< real_t, 4, 4, 3, 3 >.
void Capacity | ( | real_t | coef = 1 | ) |
Add Consistent capacity contribution to left and right-hand sides after multiplying it by coefficient coef
.
- Parameters
-
[in] coef Coefficient to multiply by added term [Default: 1
].
|
virtual |
Add consistent capacity matrix to left-hand side after multiplying it by coefficient coef
- Parameters
-
[in] coef Coefficient to multiply by added term [Default: 1
].
Reimplemented from Equa_Therm< real_t, 4, 4, 3, 3 >.
|
virtual |
Add consistent capacity contribution to right-hand side after multiplying it by coefficient coef
- Parameters
-
[in] coef Coefficient to multiply by added term [Default: 1
].
Reimplemented from Equa_Therm< real_t, 4, 4, 3, 3 >.
|
virtual |
Add diffusion matrix to left hand side after multiplying it by coefficient coef
.
- Parameters
-
[in] coef Coefficient to multiply by added term (default value = 1).
Reimplemented from Equa_Therm< real_t, 4, 4, 3, 3 >.
Add diffusion matrix to left hand side after multiplying it by coefficient coef
Case where the diffusivity matrix is given as an argument.
- Parameters
-
[in] diff Diffusion matrix (class DMatrix). [in] coef Coefficient to multiply by added term [Default: 1
].
|
virtual |
Add diffusion contribution to right hand side after multiplying it by coefficient coef
To be used for explicit diffusion term
- Parameters
-
[in] coef Coefficient to multiply by added term [Default: 1
].
Reimplemented from Equa_Therm< real_t, 4, 4, 3, 3 >.
|
virtual |
Add convection matrix to left-hand side after multiplying it by coefficient coef
Case where velocity field has been previouly defined
- Parameters
-
[in] coef Coefficient to multiply by added term [Default: 1
].
Reimplemented from Equa_Therm< real_t, 4, 4, 3, 3 >.
Add convection matrix to left-hand side after multiplying it by coefficient coef
- Parameters
-
[in] v Constant velocity vector [in] coef Coefficient to multiply by added term [Default: 1
].
Add convection matrix to left-hand side after multiplying it by coefficient coef
Case where velocity field is given by a vector v
.
- Parameters
-
[in] v Velocity vector. [in] coef Coefficient to multiply by added term [Default: 1
].
Add convection contribution to right-hand side after multiplying it by coefficient coef
To be used for explicit convection term.
- Parameters
-
[in] v Velocity vector. [in] coef Coefficient to multiply by added term [Default: 1
].
Add body right-hand side term to right hand side after multiplying it by coefficient coef
- Parameters
-
[in] ud Instance of UserData or of an inherited class. Contains a member function that provides body source. [in] coef Coefficient to multiply by added term [Default: 1
].
void BodyRHS | ( | const Vect< real_t > & | b, |
int | opt = GLOBAL_ARRAY |
||
) |
Add body right-hand side term to right hand side.
- Parameters
-
[in] b Local vector containing source at element nodes. [in] opt Vector is local ( LOCAL_ARRAY
) with size4
or global (GLOBAL_ARRAY
) with size = Number of nodes [Default:GLOBAL_ARRAY
].
Add boundary right-hand side term to right hand side after multiplying it by coefficient coef
- Parameters
-
[in] ud Instance of UserData or of an inherited class. Contains a member function that provides body source. [in] coef Value by which the added term is multiplied [Default: 1
].
void BoundaryRHS | ( | const Vect< real_t > & | b, |
int | opt = GLOBAL_ARRAY |
||
) |
Add boundary right-hand side term to right hand side after multiplying it by coefficient coef
Case where body source is given by a vector
- Parameters
-
[in] b Vector containing source at side nodes. [in] opt Vector is local ( LOCAL_ARRAY
) with size3
or global (GLOBAL_ARRAY
) with size = Number of nodes [Default:GLOBAL_ARRAY
].
void BoundaryRHS | ( | real_t | flux | ) |
Add boundary right-hand side flux to right hand side.
- Parameters
-
[in] flux Vector containing source at side nodes.
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.
- Parameters
-
[in] coef Value of penalty parameter [Default: 1.e20
].
|
virtualinherited |
Set stabilized formulation.
Stabilized variational formulations are to be used when the Péclet number is large.
By default, no stabilization is used.
|
inherited |
Run one time step.
This function performs one time step in equation solving. It is to be used only if a TRANSIENT analysis is required.
- Returns
- Return error from the linear system solver
|
inherited |
Run the equation.
If the analysis (see function setAnalysis) is STEADY_STATE
, then the function solves the stationary equation.
If the analysis is TRANSIENT
, then the function performs time stepping until the final time is reached.
Build global stiffness and mass matrices for the eigen system.
Case where the mass matrix is consistent
- Parameters
-
[in] K Stiffness matrix [in] M Consistent mass matrix
Build global diffusion and capacity matrices for the eigen system.
Case where the capacity matrix is lumped
- Parameters
-
[in] K Diffusion matrix [in] M Vector containing diagonal capacity matrix
Update Right-Hand side by taking into account essential boundary conditions.
- Parameters
-
[in] el Reference to current element instance [in] bc Vector that contains imposed values at all DOFs
Update Right-Hand side by taking into account essential boundary conditions.
- Parameters
-
[in] bc Vector that contains imposed values at all DOFs
- Remarks
- The current element is pointed by
_theElement
|
inherited |
Update element matrix to impose bc by diagonalization technique.
- Parameters
-
[in] dof_type DOF type option. To choose among the enumerated values: -
NODE_FIELD
, DOFs are supported by nodes [Default] -
ELEMENT_FIELD
, DOFs are supported by elements -
SIDE_FIELD
, DOFs are supported by sides
[in] dof DOF setting: -
= 0
, All DOFs are taken into account [Default] -
!= 0
, Only DOF No.dof
is handled in the system
-
Localize Element Vector from a Vect instance.
- Parameters
-
[in] b Reference 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.
Localize Element Vector from a Vect instance.
- Parameters
-
[in] b Global vector to be localized. [out] be Local vector, the length of which is the total number of element equations.
- Remarks
- All degrees of freedom are transferred to the local vector
|
inherited |
Localize Element Vector from a Vect instance.
- Parameters
-
[in] b Global vector to be localized. [out] be Local vector, the length of which is the total number of element equations. [in] dof Degree of freedom to transfer to the local vector
- Remarks
- Only yhe dega dof is transferred to the local vector
|
inherited |
Localize Element Vector from a Vect instance.
- Parameters
-
[in] b Global vector to be localized. [out] be Local 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.
Localize Element Vector from a Vect instance.
- Parameters
-
[in] b Global vector to be localized. [out] be Local vector, the length of which is
Localize Element Vector.
- Parameters
-
[in] b Global vector to be localized [in] dof_type DOF type option. To choose among the enumerated values: -
NODE_FIELD
, DOFs are supported by nodes [Default] -
ELEMENT_FIELD
, DOFs are supported by elements -
SIDE_FIELD
, DOFs are supported by sides
[in] flag Option to set: -
= 0
, All DOFs are taken into account [Default] -
!= 0
, Only DOF numberdof
is handled in the system
ePrev
. -
- Remarks
- This member function is to be used if a constructor with Element was invoked. It uses the Element pointer
_theElement
Localize Side Vector.
- Parameters
-
[in] b Global vector to be localized -
NODE_FIELD
, DOFs are supported by nodes [ default ] -
ELEMENT_FIELD
, DOFs are supported by elements -
SIDE_FIELD
, DOFs are supported by sides
ePrev
. -
- Remarks
- This member function is to be used if a constructor with Side was invoked. It uses the Side pointer
_theSide
|
inherited |
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
|
inherited |
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
Assemble element matrix into global one.
- Parameters
-
A Pointer to global matrix (abstract class: can be any of classes SkSMatrix, SkMatrix, SpMatrix)
- Warning
- The element pointer is given by the global variable
theElement
Assemble element matrix into global one.
- Parameters
-
A Global matrix stored as an SkSMatrix instance
- Warning
- The element pointer is given by the global variable
theElement
Assemble element matrix into global one.
- Parameters
-
[in] A Global matrix stored as an SkMatrix instance
- Warning
- The element pointer is given by the global variable
theElement
Assemble element matrix into global one.
- Parameters
-
[in] A Global matrix stored as an SpMatrix instance
- Warning
- The element pointer is given by the global variable
theElement
Assemble element matrix into global one.
- Parameters
-
[in] A Global matrix stored as an TrMatrix instance
- Warning
- The element pointer is given by the global variable
theElement
Assemble element vector into global one.
- Parameters
-
[in] v Global vector (Vect instance)
- Warning
- The element pointer is given by the global variable
theElement
Assemble element matrix into global one for the Discontinuous Galerkin approximation.
- Parameters
-
A Pointer to global matrix (abstract class: can be any of classes SkSMatrix, SkMatrix, SpMatrix)
- Warning
- The element pointer is given by the global variable
theElement
Assemble element matrix into global one for the Discontinuous Galerkin approximation.
- Parameters
-
A Global matrix stored as an SkSMatrix instance
- Warning
- The element pointer is given by the global variable
theElement
Assemble element matrix into global one for the Discontinuous Galerkin approximation.
- Parameters
-
[in] A Global matrix stored as an SkMatrix instance
- Warning
- The element pointer is given by the global variable
theElement
Assemble element matrix into global one for the Discontinuous Galerkin approximation.
- Parameters
-
[in] A Global matrix stored as an SpMatrix instance
- Warning
- The element pointer is given by the global variable
theElement
Assemble element matrix into global one for the Discontinuous Galerkin approximation.
- Parameters
-
[in] A Global matrix stored as an TrMatrix instance
- Warning
- The element pointer is given by the global variable
theElement
Assemble side (edge or face) matrix into global one.
- Parameters
-
A Pointer to global matrix (abstract class: can be any of classes SkSMatrix, SkMatrix, SpMatrix)
- Warning
- The side pointer is given by the global variable
theSide
Assemble side (edge or face) matrix into global one.
- Parameters
-
[in] A Global matrix stored as an SkSMatrix instance
- Warning
- The side pointer is given by the global variable
theSide
Assemble side (edge or face) matrix into global one.
- Parameters
-
[in] A Global matrix stored as an SkMatrix instance
- Warning
- The side pointer is given by the global variable
theSide
Assemble side (edge or face) matrix into global one.
- Parameters
-
[in] A Global matrix stored as an SpMatrix instance
- Warning
- The side pointer is given by the global variable
theSide
Assemble side (edge or face) vector into global one.
- Parameters
-
[in] v Global vector (Vect instance)
- Warning
- The side pointer is given by the global variable
theSide
Assemble product of element matrix by element vector into global vector.
- Parameters
-
[in] el Reference to Element instance [in] x Global vector to multiply by (Vect instance) [out] b Global vector to add (Vect instance)
Assemble product of side matrix by side vector into global vector.
- Parameters
-
[in] sd Reference to Side instance [in] x Global vector to multiply by (Vect instance) [out] b Global vector (Vect instance)
|
inherited |
Define a material property by an algebraic expression.
- Parameters
-
[in] exp Algebraic expression [in] prop Property name
- Returns
- Return value in expression evaluation:
-
=0
, Normal evaluation -
!=0
, An error message is displayed
-
|
inherited |
Return reference to Mesh instance.
- Returns
- Reference to Mesh instance
|
inherited |
Choose solver for the linear system.
- Parameters
-
[in] ls Solver of the linear system. To choose among the enumerated values: DIRECT_SOLVER
,CG_SOLVER
,GMRES_SOLVER
-
DIRECT_SOLVER
, Use a facorization solver [default] -
CG_SOLVER
, Conjugate Gradient iterative solver -
CGS_SOLVER
, Squared Conjugate Gradient iterative solver -
BICG_SOLVER
, BiConjugate Gradient iterative solver -
BICG_STAB_SOLVER
, BiConjugate Gradient Stabilized iterative solver -
GMRES_SOLVER
, GMRES iterative solver -
QMR_SOLVER
, QMR iterative solver
[in] pc Preconditioner to associate to the iterative solver. If the direct solver was chosen for the first argument this argument is not used. Otherwise choose among the enumerated values: -
IDENT_PREC
, Identity preconditioner (no preconditioning [default]) -
DIAG_PREC
, Diagonal preconditioner -
ILU_PREC
, Incomplete LU factorization preconditioner
-
|
inherited |
Compute eigenvalues and eigenvectors.
Eigenvalues and vectors are computed using the Bathe's subspace iteration method.
- Parameters
-
[in] nb_eigv Number of eigenvalues to compute [in] g Option to choose whether to solve a generalized eigenvalue problem (true) or a standard one (false). The generalized eigenvalue problem corresponds to the case where a consistent mass matrix (rather than a lumped one) is computed. Default value is false.
|
inherited |
Return the n-th eigenvalue.
This functions works only if the member function getEigen was called with an argument nb_eigv
greater or equal to n
. Otherwise it returns 0
.