b2api
B2000++ API Reference Manual, VERSION 4.6
 
Loading...
Searching...
No Matches
b2000::SuperElement< T, MATRIX_TYPE > Class Template Reference

#include "b2super_element.H"

Inheritance diagram for b2000::SuperElement< T, MATRIX_TYPE >:
Collaboration diagram for b2000::SuperElement< T, MATRIX_TYPE >:

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 ElementPropertyget_property () const
 
const std::vector< Element::VariableInfoget_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)
 
- Public Member Functions inherited from b2000::TypedElement< T >
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 >
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.
 
- Public Member Functions inherited from b2000::Element
virtual std::pair< int, VariableInfoget_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)
 
- Public Member Functions inherited from b2000::Object
virtual ~Object ()
 

Additional Inherited Members

- Public Types inherited from b2000::Element
enum  VariableInfo { zero = 0 , constant = 1 , linear = 2 , nonlinear = 3 }
 
- Static Public Attributes inherited from b2000::TypedElement< T >
static type_t type
 Element interface typed by the dof.
 
- Static Public Attributes inherited from b2000::Object
static ObjectType type
 

Detailed Description

template<typename T, typename MATRIX_TYPE>
class b2000::SuperElement< T, MATRIX_TYPE >

Super element

Member Function Documentation

◆ face_field_linear_on_dof()

template<typename T , typename MATRIX_TYPE >
bool b2000::SuperElement< T, MATRIX_TYPE >::face_field_linear_on_dof ( const int  face_id,
const std::string &  field_name 
)
inlinevirtual
Returns
Whether the integrated field depends linearly on the dof (true), or whether the field depends in nonlinear fashion (false).
Parameters
face_idThe face number, starting at 1, according to the B2000++ conventions.
field_nameThe 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.

◆ face_field_order()

template<typename T , typename MATRIX_TYPE >
int b2000::SuperElement< T, MATRIX_TYPE >::face_field_order ( const int  face_id,
const std::string &  field_name 
)
inlinevirtual
Returns
The integration order over the element face. This depends on the shape functions. For instance, linear elements return 1, quadratic elements return 2, etc.
Parameters
face_idThe face number, starting at 1, according to the B2000++ conventions.
field_nameThe field name. Usually, the integration order is independent of the field name.

Reimplemented from b2000::Element.

◆ face_field_polynomial_sub_face()

template<typename T , typename MATRIX_TYPE >
int b2000::SuperElement< T, MATRIX_TYPE >::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 
)
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.

Parameters
face_idThe face number, starting at 1, according to the B2000++ conventions (input).
field_nameThe field name.
sub_nodesA matrix (2 x n) containing the face-internal coordinates of the n sub-nodes (output).
sub_facesA vector of Triangle structs with the connectivities of the sub-faces (output).
Returns
The number of sub-faces.

Reimplemented from b2000::Element.

◆ face_geom()

template<typename T , typename MATRIX_TYPE >
void b2000::SuperElement< T, MATRIX_TYPE >::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 
)
inlinevirtual

Compute the face geometry at a specified integration point.

Parameters
face_idThe face number, starting at 1, according to the B2000++ conventions (input).
internal_coorThe face-internal coordinates of the integration point (input).
geomThe branch-global coordinates of the integration point (vector, output).
d_geom_d_icoorThe derivatives of the branch-global coordinates w.r.t. the internal coordinates; these are the tangential vectors (matrix, output).

Reimplemented from b2000::Element.

◆ get_global_dof_numbering()

template<typename T , typename MATRIX_TYPE >
std::pair< size_t, size_t > b2000::SuperElement< T, MATRIX_TYPE >::get_global_dof_numbering ( ) const
inlinevirtual

Elements without internal dof's do not need to override this function.

Returns
A pair (i, i+n), with i being the global number of the first dof of the element, and n equal to get_number_of_dof().

Example

std::pair<size_t, size_t>
return std::pair<size_t, size_t>(
internal_dof_index, internal_dof_index + 4);
}
std::pair< size_t, size_t > get_global_dof_numbering() const
Definition b2super_element.H:48

The variable internal_dof_index is a member of the element class and was set in set_global_dof_numbering().

Reimplemented from b2000::Element.

◆ get_nodes()

template<typename T , typename MATRIX_TYPE >
std::pair< int, Node *const * > b2000::SuperElement< T, MATRIX_TYPE >::get_nodes ( ) const
inlinevirtual

Must be overridden.

Returns
A pair containing the number of nodes and a const pointer to an array of pointers to Node instances.
See also
set_nodes()

Reimplemented from b2000::Element.

◆ get_number_of_dof()

template<typename T , typename MATRIX_TYPE >
int b2000::SuperElement< T, MATRIX_TYPE >::get_number_of_dof ( ) const
inlinevirtual

Elements without internal dof's do not need to override this function.

Returns
The number of degrees-of-freedom of the entity.

Example

A Finite Element with 4 internal dof's would return 4.

Reimplemented from b2000::Element.

◆ get_property()

template<typename T , typename MATRIX_TYPE >
const ElementProperty * b2000::SuperElement< T, MATRIX_TYPE >::get_property ( ) const
inlinevirtual

Must be overridden. This implementation returns a 0-pointer.

Returns
A pointer to the element property.
See also
set_property()

Reimplemented from b2000::Element.

◆ get_value_info()

template<typename T , typename MATRIX_TYPE >
const std::vector< Element::VariableInfo > b2000::SuperElement< T, MATRIX_TYPE >::get_value_info ( ) const
inlinevirtual

Must be overridden for elements with internal forces.

todo: not all elements override this function.

Returns
A vector containing the type of response for each variable.

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.

Examples

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:

const std::vector<VariableInfo>
get_value_info() const {
return std::vector<VariableInfo>(3, Element::nonlinear);
}
@ nonlinear
Definition b2element.H:619
const std::vector< Element::VariableInfo > get_value_info() const
Definition b2super_element.H:78

whereas a stress element supporting only linear analysis without any viscosity and inertia effects specify:

const std::vector<VariableInfo>
get_value_info() const {
return std::vector<VariableInfo>(1, Element::linear);
}
@ linear
Definition b2element.H:615

Nonlinear heat elements specify:

const std::vector<VariableInfo>
get_value_info() const {
return std::vector<VariableInfo>(2, Element::nonlinear);
}

Reimplemented from b2000::Element.

◆ set_global_dof_numbering()

template<typename T , typename MATRIX_TYPE >
size_t b2000::SuperElement< T, MATRIX_TYPE >::set_global_dof_numbering ( size_t  index)
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.

Parameters
indexThe index into the global solution vector of the first degree-of-freedom of the entity.
Returns
The index into the global solution vector past the last degree-of-freedom of the entity.

Example

A Finite Element with 4 internal dof's:

size_t
set_global_dof_numbering(size_t index) {
internal_dof_index = index;
return index + 4;
}
size_t set_global_dof_numbering(size_t index)
Definition b2super_element.H:43

The variable internal_dof_index is a member of the element class.

Reimplemented from b2000::Element.

◆ set_nodes()

template<typename T , typename MATRIX_TYPE >
void b2000::SuperElement< T, MATRIX_TYPE >::set_nodes ( std::pair< int, Node *const * >  nodes)
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!

Parameters
nodesA 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().

Example

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:

using NodeType = HNode<Coor<Trans<X,Y,Z>>,Dof<DTrans<DX,DY>>>;

The class has an array of NodeType pointers.

private:
NodeType* nodes[4];

This array is initialized in the set_nodes() function:

void
set_nodes(std::pair<int, Node* const *> nodes_) {
if (nodes_.first != 4)
SizeError() << "Element " << get_object_name() << " of type "
"Q4.S.2D.EXAMPLE needs 4 nodes." << THROW;
for (int i = 0; i != 4; ++i) {
nodes[i] = dynamic_cast<NodeType*>(nodes_.second[i]);
if (nodes[i] == 0)
TypeError() << "Element " << get_object_name() <<
" of type Q4.S.2D.EXAMPLE requires nodes of type " <<
typeid(NodeType) << "." << THROW;
}
}
#define THROW
Definition b2exception.H:198
const std::string & get_object_name() const override
Definition b2element.H:220
void set_nodes(std::pair< int, Node *const * > nodes_)
Definition b2super_element.H:60
GenericException< SizeError_name > SizeError
Definition b2exception.H:346
GenericException< TypeError_name > TypeError
Definition b2exception.H:325

The code to retrieve the coordinates of all 4 nodes of the element, without the overhead of virtual functions could then look like this:

double coor[4][3];
for (int i = 0; i != 4; ++i)
NodeType::get_coor_s(coor[i], nodes[i]);

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:

{
dof_numbering.resize(8);
size_t* ptr = &dof_numbering[0];
for (int i = 0; i != 4; ++i)
ptr = NodeType::get_global_dof_numbering_s(ptr, nodes[i]);
}

Reimplemented from b2000::Element.

◆ set_property()

template<typename T , typename MATRIX_TYPE >
void b2000::SuperElement< T, MATRIX_TYPE >::set_property ( ElementProperty property)
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!

Parameters
propertyA 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.

Example

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.

SolidMechanics2D* property;
Definition b2solid_mechanics.H:562

The function is implemented as follows:

void
property = dynamic_cast<SolidMechanics2D*>(property_);
if (property == 0)
TypeError() <<
"Element " << get_object_name() << " of type Q4.S.2D.EXAMPLE "
"does not accept materials of type " << typeid(*property_) <<
" but requires instead materials derived from the class " <<
typeid(SolidMechanics2D) << "." << THROW;
}
Definition b2element.H:71
void set_property(ElementProperty *property_)
Definition b2super_element.H:68

Reimplemented from b2000::Element.


The documentation for this class was generated from the following file: