This base class provides the template for assembling the matrix and right-hand-side vector when solving any partial differential equations using finite element method.
More...
#include <assemble_system.h>
|
| | AssembleSystem (EquationSystems &es) |
| | Constructor. More...
|
| |
| | ~AssembleSystem () |
| | Destructor. More...
|
| |
| virtual void | assemble_global_K (const std::string &system_name, const std::string &option)=0 |
| | Assemble the Global Matrix K. More...
|
| |
| virtual void | assemble_global_F (const std::string &system_name, const std::string &option)=0 |
| | Assemble the Global force vector F. More...
|
| |
| virtual void | assemble_element_KIJ (const std::vector< Real > &JxW, const std::vector< std::vector< RealGradient > > &dphi, const unsigned int n_u_dofs, const unsigned int I, const unsigned int J, DenseMatrix< Number > &Kij)=0 |
| | Assemble the element matrix K_IJ. More...
|
| |
| virtual void | compute_element_rhs (const Elem *elem, const unsigned int n_u_dofs, FEBase &fe_v, const std::vector< std::size_t > n_list, const bool &pf_flag, const std::string &option, const Real &alpha, DenseVector< Number > &Fe)=0 |
| | Assemble function for calculating each element's contribution to the right-hand-side vector. It is used in assemble_global_F() More...
|
| |
| virtual void | apply_bc_by_penalty (const Elem *elem, const std::string &matrix_or_vector, DenseMatrix< Number > &Ke, DenseVector< Number > &Fe, const std::string &option)=0 |
| | Apply Boundary Conditions by penalty method. More...
|
| |
| void | assemble_int_force (const Elem *elem, const unsigned int n_u_dofs, FEBase &fe_v) |
| | Assemble int_force matrix for every element, this includes Gaussian quadrature weights multiplied by shape functions. The product is calculated once and is stored in _int_force. More...
|
| |
| virtual void | select_boundary_side (const Elem *elem)=0 |
| | select sides on the boundary for all elements More...
|
| |
| void | penalize_elem_matrix_vector (DenseMatrix< Number > &Ke, DenseVector< Number > &Fe, const std::string &matrix_or_vector, const unsigned int &var_number, const unsigned int &local_node_id, const unsigned int &n_nodes_elem, const Real &penalty, const Real &value) |
| | Universal function to penalize element matrix or vector with a large number. More...
|
| |
|
| EquationSystems & | _eqn_sys |
| |
| unsigned int | _dim |
| |
| MeshBase & | _mesh |
| |
| const std::vector< boundary_id_type > | _boundary_id_3D = {4,2,1,3,0,5} |
| |
| const std::vector< std::string > | _boundary_name_3D = {"left", "right", "bottom", "top", "back", "front"} |
| |
| std::vector< std::vector< Real > > | _int_force |
| |
| std::vector< std::vector< Point > > | _q_xyz |
| |
| std::vector< unsigned int > | _n_dofs |
| |
| std::vector< unsigned int > | _n_u_dofs |
| |
| std::vector< unsigned int > | _n_p_dofs |
| |
| std::vector< unsigned int > | _n_uvw_dofs |
| |
| std::vector< std::vector< dof_id_type > > | _dof_indices |
| |
| std::vector< std::vector< dof_id_type > > | _dof_indices_u |
| |
| std::vector< std::vector< dof_id_type > > | _dof_indices_p |
| |
| std::vector< std::vector< unsigned int > > | _boundary_sides |
| |
This base class provides the template for assembling the matrix and right-hand-side vector when solving any partial differential equations using finite element method.
| AssembleSystem::AssembleSystem |
( |
EquationSystems & |
es | ) |
|
| AssembleSystem::~AssembleSystem |
( |
| ) |
|
| virtual void AssembleSystem::apply_bc_by_penalty |
( |
const Elem * |
elem, |
|
|
const std::string & |
matrix_or_vector, |
|
|
DenseMatrix< Number > & |
Ke, |
|
|
DenseVector< Number > & |
Fe, |
|
|
const std::string & |
option |
|
) |
| |
|
pure virtual |
Apply Boundary Conditions by penalty method.
Implemented in AssembleStokes.
| virtual void AssembleSystem::assemble_element_KIJ |
( |
const std::vector< Real > & |
JxW, |
|
|
const std::vector< std::vector< RealGradient > > & |
dphi, |
|
|
const unsigned int |
n_u_dofs, |
|
|
const unsigned int |
I, |
|
|
const unsigned int |
J, |
|
|
DenseMatrix< Number > & |
Kij |
|
) |
| |
|
pure virtual |
Assemble the element matrix K_IJ.
Reinit and compute the element matrix K_ij, which will be added into K matrix after calling assemble_global_K(). For Stokes equation, size of this submatrix is n_u_dofs * n_u_dofs = n_v_dofs * n_v_dofs = n_w_dofs * n_w_dofs
Implemented in AssembleStokes.
| virtual void AssembleSystem::assemble_global_F |
( |
const std::string & |
system_name, |
|
|
const std::string & |
option |
|
) |
| |
|
pure virtual |
Assemble the Global force vector F.
- Parameters
-
| [in] | system_name | Name of the system (could be "Stokes", "Poisson", or "NP") |
| [in] | option | Options of assembling the system ("disturbed" or "undisturbed") for Stokes equation |
| [out] | Fe | Add rhs vector to system. |
Implemented in AssembleStokes.
| virtual void AssembleSystem::assemble_global_K |
( |
const std::string & |
system_name, |
|
|
const std::string & |
option |
|
) |
| |
|
pure virtual |
Assemble the Global Matrix K.
- Parameters
-
| [in] | system_name | Name of the system (could be "Stokes", "Poisson", or "NP") |
| [in] | Option | options of assembling the system ("disturbed" or "undisturbed") for Stokes equation |
| [out] | Ke | Add element matrix to system |
Implemented in AssembleStokes.
| void AssembleSystem::assemble_int_force |
( |
const Elem * |
elem, |
|
|
const unsigned int |
n_u_dofs, |
|
|
FEBase & |
fe_v |
|
) |
| |
Assemble int_force matrix for every element, this includes Gaussian quadrature weights multiplied by shape functions. The product is calculated once and is stored in _int_force.
| virtual void AssembleSystem::compute_element_rhs |
( |
const Elem * |
elem, |
|
|
const unsigned int |
n_u_dofs, |
|
|
FEBase & |
fe_v, |
|
|
const std::vector< std::size_t > |
n_list, |
|
|
const bool & |
pf_flag, |
|
|
const std::string & |
option, |
|
|
const Real & |
alpha, |
|
|
DenseVector< Number > & |
Fe |
|
) |
| |
|
pure virtual |
| void AssembleSystem::penalize_elem_matrix_vector |
( |
DenseMatrix< Number > & |
Ke, |
|
|
DenseVector< Number > & |
Fe, |
|
|
const std::string & |
matrix_or_vector, |
|
|
const unsigned int & |
var_number, |
|
|
const unsigned int & |
local_node_id, |
|
|
const unsigned int & |
n_nodes_elem, |
|
|
const Real & |
penalty, |
|
|
const Real & |
value |
|
) |
| |
Universal function to penalize element matrix or vector with a large number.
| virtual void AssembleSystem::select_boundary_side |
( |
const Elem * |
elem | ) |
|
|
pure virtual |
select sides on the boundary for all elements
Implemented in AssembleStokes.
| const std::vector<boundary_id_type> AssembleSystem::_boundary_id_3D = {4,2,1,3,0,5} |
|
protected |
| const std::vector<std::string> AssembleSystem::_boundary_name_3D = {"left", "right", "bottom", "top", "back", "front"} |
|
protected |
| std::vector<std::vector<unsigned int> > AssembleSystem::_boundary_sides |
|
protected |
| unsigned int AssembleSystem::_dim |
|
protected |
| std::vector<std::vector<dof_id_type> > AssembleSystem::_dof_indices |
|
protected |
| std::vector<std::vector<dof_id_type> > AssembleSystem::_dof_indices_p |
|
protected |
| std::vector<std::vector<dof_id_type> > AssembleSystem::_dof_indices_u |
|
protected |
| EquationSystems& AssembleSystem::_eqn_sys |
|
protected |
| std::vector<std::vector<Real> > AssembleSystem::_int_force |
|
protected |
| MeshBase& AssembleSystem::_mesh |
|
protected |
| std::vector<unsigned int> AssembleSystem::_n_dofs |
|
protected |
| std::vector<unsigned int> AssembleSystem::_n_p_dofs |
|
protected |
| std::vector<unsigned int> AssembleSystem::_n_u_dofs |
|
protected |
| std::vector<unsigned int> AssembleSystem::_n_uvw_dofs |
|
protected |
| std::vector<std::vector<Point> > AssembleSystem::_q_xyz |
|
protected |
The documentation for this class was generated from the following files: