Builds finite element arrays for incompressible Navier-Stokes equations in 2-D domains using Q1/P0 element and a penaly formulation for the incompressibility condition. More...
#include <NSP2DQ41.h>

Public Member Functions | |
NSP2DQ41 (Mesh &ms) | |
Constructor using mesh data. More... | |
NSP2DQ41 (Mesh &ms, Vect< real_t > &u) | |
Constructor using mesh data and velocity vector. More... | |
~NSP2DQ41 () | |
Destructor. | |
void | setPenalty (real_t lambda) |
Define penalty parameter. More... | |
void | setInput (EqDataType opt, Vect< real_t > &u) |
Set equation input data. More... | |
void | Periodic (real_t coef=1.e20) |
Add contribution of periodic boundary condition (by a penalty technique). More... | |
void | build () |
Build the linear system of equations. More... | |
int | runOneTimeStep () |
Run one time step. More... | |
![]() | |
Equa_Fluid () | |
Default constructor. More... | |
virtual | ~Equa_Fluid () |
Destructor. | |
void | Reynolds (const real_t &Re) |
Set Reynolds number. | |
void | Viscosity (const real_t &visc) |
Set (constant) Viscosity. | |
void | Viscosity (const string &exp) |
Set viscosity given by an algebraic expression. | |
void | Density (const real_t &dens) |
Set (constant) Viscosity. | |
void | Density (const string &exp) |
Set Density given by an algebraic expression. | |
void | ThermalExpansion (const real_t *e) |
Set (constant) thermal expansion coefficient. | |
void | ThermalExpansion (const string &exp) |
Set thermal expansion coefficient given by an algebraic expression. | |
void | setMaterial () |
Set material properties. | |
![]() | |
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, std::valarray< 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, std::valarray< 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... | |
![]() | |
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. | |
Detailed Description
Builds finite element arrays for incompressible Navier-Stokes equations in 2-D domains using Q1/P0 element and a penaly formulation for the incompressibility condition.
- Copyright
- GNU Lesser Public License
Constructor & Destructor Documentation
◆ NSP2DQ41() [1/2]
◆ NSP2DQ41() [2/2]
Constructor using mesh data and velocity vector.
- Parameters
-
[in] ms Mesh instance [in,out] u Velocity vector
Member Function Documentation
◆ build()
void build | ( | ) |
Build the linear system of equations.
Before using this function, one must have properly selected appropriate options for:
- The choice of a steady state or transient analysis. By default, the analysis is stationary
- In the case of transient analysis, the choice of a time integration scheme and a lumped or consistent capacity matrix. If transient analysis is chosen, the lumped capacity matrix option is chosen by default, and the implicit Euler scheme is used by default for time integration.
◆ Periodic()
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
].
◆ runOneTimeStep()
int runOneTimeStep | ( | ) |
Run one time step.
This function performs one time step, once a time integration scheme has been selected.
◆ setInput()
void setInput | ( | EqDataType | opt, |
Vect< real_t > & | u | ||
) |
Set equation input data.
- Parameters
-
[in] opt Parameter that selects data type for input. This parameter is to be chosen in the enumerated variable EqDataType [in] u Vect instance that contains input vector data List of data types contains INITIAL_FIELD
,BOUNDARY_CONDITION_DATA
,SOURCE_DATA
orFLUX
with obvious meaning
◆ setPenalty()
void setPenalty | ( | real_t | lambda | ) |
Define penalty parameter.
Penalty parameter is used to enforce the incompressibility constraint
- Parameters
-
[in] lambda Penaly parameter: Large value [Default: 1.e07
]