#include "b2super_element.H"
Public Member Functions | |
int | get_number_of_dof () const |
size_t | set_global_dof_numbering (size_t index) |
std::pair< size_t, size_t > | get_global_dof_numbering () const |
void | set_nodes (std::pair< int, Node *const * > nodes_) |
std::pair< int, Node *const * > | get_nodes () const |
void | set_property (ElementProperty *property_) |
const ElementProperty * | get_property () const |
const std::vector< Element::VariableInfo > | get_value_info () const |
int | face_field_order (const int face_id, const std::string &field_name) |
bool | face_field_linear_on_dof (const int face_id, const std::string &field_name) |
int | face_field_polynomial_sub_face (const int face_id, const std::string &field_name, b2linalg::Matrix< double, b2linalg::Mrectangle > &sub_nodes, std::vector< Element::Triangle > &sub_faces) |
void | face_geom (const int face_id, const b2linalg::Vector< double, b2linalg::Vdense_constref > &internal_coor, b2linalg::Vector< double > &geom, b2linalg::Matrix< double, b2linalg::Mrectangle > &d_geom_d_icoor) |
![]() | |
virtual void | get_constraint (Model &model, const bool linear, const double time, const b2linalg::Matrix< T, b2linalg::Mrectangle_constref > &dof, b2linalg::Index &dof_numbering, b2linalg::Vector< T, b2linalg::Vdense > &constraint, b2linalg::Matrix< T, b2linalg::Mcompressed_col > &trans_d_constraint_d_dof, b2linalg::Vector< T, b2linalg::Vdense > &d_constraint_d_time) |
virtual void | edge_field_value (const int edge_id, const std::string &field_name, const double internal_coor, const b2linalg::Matrix< T, b2linalg::Mrectangle_constref > &dof, const double time, b2linalg::Vector< T, b2linalg::Vdense > &value, b2linalg::Vector< T, b2linalg::Vdense > &d_value_d_icoor, b2linalg::Index &dof_numbering, b2linalg::Matrix< T, b2linalg::Mrectangle > &d_value_d_dof, b2linalg::Index &d_value_d_dof_dep_col=b2linalg::Index::null) |
virtual void | face_field_value (const int face_id, const std::string &field_name, const b2linalg::Vector< double, b2linalg::Vdense_constref > &internal_coor, const b2linalg::Matrix< T, b2linalg::Mrectangle_constref > &dof, const double time, b2linalg::Vector< T, b2linalg::Vdense > &value, b2linalg::Matrix< T, b2linalg::Mrectangle > &d_value_d_icoor, b2linalg::Index &dof_numbering, b2linalg::Matrix< T, b2linalg::Mrectangle > &d_value_d_dof, b2linalg::Index &d_value_d_dof_dep_col) |
virtual void | body_field_value (const std::string &field_name, const b2linalg::Vector< double, b2linalg::Vdense_constref > &internal_coor, const b2linalg::Matrix< T, b2linalg::Mrectangle_constref > &dof, const double time, b2linalg::Vector< T, b2linalg::Vdense > &value, b2linalg::Matrix< T, b2linalg::Mrectangle > &d_value_d_icoor, b2linalg::Index &dof_numbering, b2linalg::Matrix< T, b2linalg::Mrectangle > &d_value_d_dof) |
virtual void | field_volume_integration (const b2linalg::Vector< T, b2linalg::Vdense_constref > &dof, const b2linalg::Vector< T, b2linalg::Vdense_constref > &dofdot, const b2linalg::Vector< T, b2linalg::Vdense_constref > &dofdotdot, const Field< T > &f, b2linalg::Index &dof_numbering, b2linalg::Vector< T, b2linalg::Vdense > &discretised_field, b2linalg::Matrix< T, b2linalg::Mpacked > &d_discretised_field_d_dof) |
virtual void | field_face_integration (const b2linalg::Vector< T, b2linalg::Vdense_constref > &dof, const b2linalg::Vector< T, b2linalg::Vdense_constref > &dofdot, const b2linalg::Vector< T, b2linalg::Vdense_constref > &dofdotdot, const Field< T > &field, int face, b2linalg::Index &dof_numbering, b2linalg::Vector< T, b2linalg::Vdense > &discretised_field, b2linalg::Matrix< T, b2linalg::Mpacked > &d_discretised_field_d_dof) |
virtual void | field_edge_integration (const b2linalg::Vector< T, b2linalg::Vdense_constref > &dof, const b2linalg::Vector< T, b2linalg::Vdense_constref > &dofdot, const b2linalg::Vector< T, b2linalg::Vdense_constref > &dofdotdot, const Field< T > &f, int edge, b2linalg::Index &dof_numbering, b2linalg::Vector< T, b2linalg::Vdense > &discretised_field, b2linalg::Matrix< T, b2linalg::Mpacked > &d_discretised_field_d_dof) |
template<typename MATRIX_FORMAT > | |
void | AssembleElementEffectiveSystem (const solver::TypedSolver< T, MATRIX_FORMAT > &solver, b2linalg::Matrix< T, typename MATRIX_FORMAT::dense > &k_eff, b2linalg::Vector< T, b2linalg::Vdense > &f_eff, b2linalg::Index &index) |
Assembles the element effective matrix and vector. | |
template<typename MATRIX_FORMAT > | |
b2linalg::Matrix< T, typename MATRIX_FORMAT::dense > | AssembleElementEffectiveMatrix (const solver::TypedSolver< T, MATRIX_FORMAT > &solver, b2linalg::Index &index) |
Assembles the element effective matrix. | |
template<typename MATRIX_FORMAT > | |
T | ComputeElementError (const solver::TypedSolver< T, MATRIX_FORMAT > &solver) |
TODO. | |
template<typename MATRIX_FORMAT > | |
b2linalg::Vector< T, b2linalg::Vdense > | AssembleElementEffectiveVector (const solver::TypedSolver< T, MATRIX_FORMAT > &solver, b2linalg::Index index) |
Assembles the element effective vector. | |
template<typename MATRIX_FORMAT > | |
void | ComputeElementGradient (const solver::TypedSolver< T, MATRIX_FORMAT > &solver, GradientContainer *const gradient_container) |
Compute the gradients for post-processing. | |
virtual void | get_value (Model &model, const bool linear, const EquilibriumSolution equilibrium_solution, const double time, const double delta_time, const b2linalg::Matrix< T, b2linalg::Mrectangle_constref > &dof, GradientContainer *gradient_container, SolverHints *solver_hints, b2linalg::Index &dof_numbering, b2linalg::Vector< T, b2linalg::Vdense > &value, const std::vector< bool > &d_value_d_dof_flags, std::vector< b2linalg::Matrix< T, b2linalg::Mrectangle > > &d_value_d_dof, b2linalg::Vector< T, b2linalg::Vdense > &d_value_d_time) |
Compute the internal forces and their derivatives of the element. | |
virtual void | get_value (Model &model, const bool linear, const EquilibriumSolution equilibrium_solution, const double time, const double delta_time, const b2linalg::Matrix< T, b2linalg::Mrectangle_constref > &dof, GradientContainer *gradient_container, SolverHints *solver_hints, b2linalg::Index &dof_numbering, b2linalg::Vector< T, b2linalg::Vdense > &value, const std::vector< bool > &d_value_d_dof_flags, std::vector< b2linalg::Matrix< T, b2linalg::Mpacked > > &d_value_d_dof, b2linalg::Vector< T, b2linalg::Vdense > &d_value_d_time) |
Compute the internal forces and their derivatives of the element. | |
![]() | |
virtual std::pair< int, VariableInfo > | get_constraint_info () |
Element (const bool symmetric) | |
Sets the symmetry property of this element. Call this constructor if your newly defined element is not symmetric. | |
size_t | get_id () const |
void | set_id (size_t id_) |
int | get_elem_type_number () const |
void | set_elem_type_number (int elno) |
const std::string & | get_object_name () const override |
void | set_object_name (const std::string &name) |
virtual void | get_dof_numbering (b2linalg::Index &dof_numbering) |
virtual void | set_additional_properties (const RTable &rtable) |
virtual void | init (Model &model) |
virtual std::pair< size_t, Element ** > | get_subelements () |
iterator | begin () |
iterator | end () const |
virtual int | edge_field_order (const int edge_id, const std::string &field_name) |
virtual bool | edge_field_linear_on_dof (const int edge_id, const std::string &field_name) |
virtual int | edge_field_polynomial_sub_edge (const int edge_id, const std::string &field_name, b2linalg::Vector< double, b2linalg::Vdense > &sub_nodes) |
virtual void | edge_geom (const int edge_id, const double internal_coor, b2linalg::Vector< double > &geom, b2linalg::Vector< double > &d_geom_d_icoor) |
virtual int | body_field_order (const std::string &field_name) |
virtual bool | body_field_linear_on_dof (const std::string &field_name) |
virtual int | body_field_polynomial_sub_volume (const std::string &field_name, b2linalg::Matrix< double, b2linalg::Mrectangle > &sub_nodes, std::vector< Tetrahedral > &sub_volumes) |
virtual void | body_geom (const b2linalg::Vector< double, b2linalg::Vdense_constref > &internal_coor, b2linalg::Vector< double > &geom, b2linalg::Matrix< double, b2linalg::Mrectangle > &d_geom_d_icoor) |
virtual std::pair< int, int > | get_state_buffer_size () const |
virtual std::pair< int *, double * > | get_state_buffer (std::pair< int *, double * > buffers) const |
virtual std::pair< const int *, const double * > | set_state_buffer (std::pair< const int *, const double * > buffers) |
![]() | |
virtual | ~Object () |
Additional Inherited Members | |
![]() | |
enum | VariableInfo { zero = 0 , constant = 1 , linear = 2 , nonlinear = 3 } |
![]() | |
static type_t | type |
Element interface typed by the dof. | |
![]() | |
static ObjectType | type |
Super element
|
inlinevirtual |
face_id | The face number, starting at 1, according to the B2000++ conventions. |
field_name | The field name. |
Solid elements depend linearly on the dof. Shell elements depend in nonlinear fashion on the side faces and on the lower and upper surfaces. Their reference surface depends linearly on the dof.
Reimplemented from b2000::Element.
|
inlinevirtual |
face_id | The face number, starting at 1, according to the B2000++ conventions. |
field_name | The field name. Usually, the integration order is independent of the field name. |
Reimplemented from b2000::Element.
|
inlinevirtual |
Sub-divide the face into triangular facets (sub-faces) and compute the face-internal coordinates of the sub-face vertices. For triangular faces, one triangle is sufficient. Quadrilateral faces can be sub-divided into two triangular facets.
face_id | The face number, starting at 1, according to the B2000++ conventions (input). |
field_name | The field name. |
sub_nodes | A matrix (2 x n) containing the face-internal coordinates of the n sub-nodes (output). |
sub_faces | A vector of Triangle structs with the connectivities of the sub-faces (output). |
Reimplemented from b2000::Element.
|
inlinevirtual |
Compute the face geometry at a specified integration point.
face_id | The face number, starting at 1, according to the B2000++ conventions (input). |
internal_coor | The face-internal coordinates of the integration point (input). |
geom | The branch-global coordinates of the integration point (vector, output). |
d_geom_d_icoor | The derivatives of the branch-global coordinates w.r.t. the internal coordinates; these are the tangential vectors (matrix, output). |
Reimplemented from b2000::Element.
|
inlinevirtual |
Elements without internal dof's do not need to override this function.
The variable internal_dof_index is a member of the element class and was set in set_global_dof_numbering().
Reimplemented from b2000::Element.
|
inlinevirtual |
Must be overridden.
Reimplemented from b2000::Element.
|
inlinevirtual |
Elements without internal dof's do not need to override this function.
A Finite Element with 4 internal dof's would return 4.
Reimplemented from b2000::Element.
|
inlinevirtual |
Must be overridden. This implementation returns a 0-pointer.
Reimplemented from b2000::Element.
|
inlinevirtual |
Must be overridden for elements with internal forces.
todo: not all elements override this function.
With "variable", the dof's and their time derivatives are meant. By default, the element specifies no response. Thus, this function must be overridden by elements producing a response in function of the dof's and/or their time derivatives. Elements formulating constraints only (such as rigid-body elements) do not need to override this function.
For stress elements, the variables are the displacements, the velocities, and the accelerations. A fully nonlinear element with support for stiffness, viscosity, and inertia effects specify:
whereas a stress element supporting only linear analysis without any viscosity and inertia effects specify:
Nonlinear heat elements specify:
Reimplemented from b2000::Element.
|
inlinevirtual |
Elements without internal dof's do not need to override this function.
Set the global degree-of-freedom numbers of the degrees-of-freedom belonging to the entity.
index | The index into the global solution vector of the first degree-of-freedom of the entity. |
A Finite Element with 4 internal dof's:
The variable internal_dof_index is a member of the element class.
Reimplemented from b2000::Element.
|
inlinevirtual |
Must be overridden. Set the nodes of the element. This function is called by the Domain object.
todo: not all elements override this function!
nodes | A pair (first, second) with first being the number of nodes, and second a pointer to a const array containing pointers to the Node instances. The length of the array corresponds to the number of nodes. |
Implementations should verify the number of nodes. The checking of the correct node type is optional.
It is recommended to perform a dynamic cast to the exact node type, this allows the use of the get_coor_s() and get_dof_numbering_s() functions without the overhead of the virtual functions Node::get_coor() and DegreesOfFreedom::get_global_dof_numbering().
A quadrilateral 2D stress element named Q4.S.2D.EXAMPLE. It has 3 coordinates (the Z-coordinate set to 0) and displacement dof's in X- and Y-direction. The following using statement gives the alias NodeType for the exact Node type:
The class has an array of NodeType pointers.
This array is initialized in the set_nodes() function:
The code to retrieve the coordinates of all 4 nodes of the element, without the overhead of virtual functions could then look like this:
The dof_numbering index vector (as given in the get_value(), get_constraint(), etc. functions) can then be initialized without the overhead of virtual functions:
Reimplemented from b2000::Element.
|
inlinevirtual |
Must be overridden. Set the element property. This function is called by the Domain object. This implementation performs no operation.
todo: not all elements override this function!
property | A pointer to the ElementProperty instance. |
ElementProperty instances contain material parameters and additional element properties (such as the thickness for shell elements).
Depending on the type of element (stress, heat-transfer, 2D, 3D, etc.), the appropriate property class should be used. The Domain object constructs a "generic" property object for each defined material, and passes a pointer to the Element instance by calling the set_property() function. The implementation of set_property() must perform a dynamic cast to the desired property class and check the result.
An element with 4 nodes and type Q4.S.2D.EXAMPLE. The class has a private variable for storing the pointer to the property. In this case, the element is a 2D stress element, and accordingly, the property class is SolidMechanics2D.
The function is implemented as follows:
Reimplemented from b2000::Element.