b2api
B2000++ API Reference Manual, VERSION 4.6
 
Loading...
Searching...
No Matches
b2000::Element Class Reference

Defines the complete interface for Element instances (C++ representations of Finite Elements). More...

#include "b2element.H"

Inheritance diagram for b2000::Element:
Collaboration diagram for b2000::Element:

Classes

class  iterator
 
struct  Tetrahedral
 
struct  Triangle
 

Element response

Implementations for elements with internal forces need to override the get_value_info() and get_value() functions, implementation for elements with constraints (e.g. rigid-body and coupling elements) need to override the get_constraint_info() and get_constraint() functions.

Elements may have at the same time internal forces and constraints. Note that TypedEssentialBoundaryCondition provides an alternative method for specifying constraints.

enum  VariableInfo { zero = 0 , constant = 1 , linear = 2 , nonlinear = 3 }
 
virtual const std::vector< VariableInfoget_value_info () const
 
virtual std::pair< int, VariableInfoget_constraint_info ()
 

Constructor

Element objects are constructed by the Domain object, using the default constructor.

 Element (const bool symmetric)
 Sets the symmetry property of this element. Call this constructor if your newly defined element is not symmetric.
 

Element Identifiers and Names

Derived classes should not override these functions. The set-functions are called by the Domain object.

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)
 

Element-internal degrees-of-freedom

Elements without internal dof's do not need to override these functions.

int get_number_of_dof () const override
 
size_t set_global_dof_numbering (size_t index) override
 
std::pair< size_t, size_t > get_global_dof_numbering () const override
 
virtual void get_dof_numbering (b2linalg::Index &dof_numbering)
 

Nodes, ElementProperty, and initialization

Implementations must override these functions.

virtual std::pair< int, Node *const * > get_nodes () const
 
virtual void set_nodes (std::pair< int, Node *const * > nodes)
 
virtual const ElementPropertyget_property () const
 
virtual void set_property (ElementProperty *property)
 
virtual void set_additional_properties (const RTable &rtable)
 
virtual void init (Model &model)
 

Sub-elements

These functions do not need not be overwritten for most implementations.

Sub-elements are useful for elements with a large number of degrees-of-freedom and a sparse element stiffness matrix. In this case, it is computationally more effective to execute the assembly by means of sub-elements, where the number of dofs for each sub-element are small.

Sub-elements should be created in the init() function and destroyed in the destructor.

virtual std::pair< size_t, Element ** > get_subelements ()
 
iterator begin ()
 
iterator end () const
 

Functions to integrate along edges

Elements that support edge integrations should override these functions. Note that, while shell element edges are actually surfaces, integration by means of the edge_* functions is performed along the reference surface.

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)
 

Functions to integrate over faces

Elements that support face integration (e.g. traction forces) should override these functions. Solid elements have faces, and the sides of shell elements and the upper and lower and reference surfaces of shell elements are considered faces, too.

virtual int face_field_order (const int face_id, const std::string &field_name)
 
virtual bool face_field_linear_on_dof (const int face_id, const std::string &field_name)
 
virtual int face_field_polynomial_sub_face (const int face_id, const std::string &field_name, b2linalg::Matrix< double, b2linalg::Mrectangle > &sub_nodes, std::vector< Triangle > &sub_faces)
 
virtual 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)
 

Functions to integrate over the element volume

Elements that support volume integration (e.g. body forces) should override these functions.

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)
 

Save and restore element state

Elements without need for restart functionality need not override these functions.

These functions are called by the Domain::save_state() and Domain::load_state() functions, which themselves are called by the nonlinear Solver(s).

The element state comprises the path-dependent state of all ElementProperty objects owned by the Element plus the path-dependent state of the Element itself (e.g. configuration of the last converged solution for elements using the Update Lagrange formulation).

There are two buffers: One for integer numbers, and one for floating-point numbers. Elements write to and read from both buffers. The Domain object then writes/reads the buffers to/from the database as positional datasets "STATE.GLOB.*".

The builtin B2000++ elements do not implement these functions and therefore, do not support restart when used in conjunction with path-dependent ElementProperty instances (the elements themselves are Total-Lagrangian, thus path-independent).

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)
 

Additional Inherited Members

- Public Member Functions inherited from b2000::Object
virtual ~Object ()
 
- Static Public Attributes inherited from b2000::Object
static ObjectType type
 

Detailed Description

Defines the complete interface for Element instances (C++ representations of Finite Elements).

Note
Implementations should not directly derive from Element, but instead from one of:
  • b2000::TypedElement<double> for real-valued elements.
  • b2000::TypedElement<b2000::csda<double>> for csda-valued elements.
  • b2000::TypedElement<std::complex<double>> for complex-valued elements.

Member Enumeration Documentation

◆ VariableInfo

How the element response varies in function of a variable.

Enumerator
zero 

No response is computed for this variable.

constant 

The element response is constant for this variable.

linear 

The element response is linear in function of this variable.

nonlinear 

The element response is nonlinear in function of this variable.

Constructor & Destructor Documentation

◆ Element()

b2000::Element::Element ( const bool  symmetric)
inlineexplicit

Sets the symmetry property of this element. Call this constructor if your newly defined element is not symmetric.

Parameters
symmetriceither true or false

Member Function Documentation

◆ begin()

iterator b2000::Element::begin ( )
inline

Returns an iterator over the element and its subelements.

◆ body_field_linear_on_dof()

virtual bool b2000::Element::body_field_linear_on_dof ( 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
field_nameThe field name.

Solid elements depend linearly on the dof. Shell elements depend in nonlinear fashion on the dof.

◆ body_field_order()

virtual int b2000::Element::body_field_order ( const std::string &  field_name)
inlinevirtual
Returns
The integration order over the element volume. This depends on the shape functions. For instance, linear elements return 1, quadratic elements return 2, etc.
Parameters
field_nameThe field name. Usually, the integration order is independent of the field name.

◆ body_field_polynomial_sub_volume()

virtual int b2000::Element::body_field_polynomial_sub_volume ( const std::string &  field_name,
b2linalg::Matrix< double, b2linalg::Mrectangle > &  sub_nodes,
std::vector< Tetrahedral > &  sub_volumes 
)
inlinevirtual

Sub-divide the face into tetrahedron (sub-volumes) and compute the element-internal coordinates of the sub-volume vertices. For tetrahedral elements, specifying one tetrahedron is sufficient. Hexahedral elements can be subdivided in 5 or 6 tetrahedrons.

Parameters
field_nameThe field name.
sub_nodesA matrix (3 x n) containing the element-internal coordinates of the n sub-nodes (output).
sub_volumesA vector of Tetrahedral structs with the connectivities of the sub-volume (output).
Returns
The number of sub-volumes.

◆ body_geom()

virtual void b2000::Element::body_geom ( 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 body geometry at a specified integration point.

Parameters
internal_coorThe element-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; this is the transposed of the Jacobian matrix (matrix, output).

◆ edge_field_linear_on_dof()

virtual bool b2000::Element::edge_field_linear_on_dof ( const int  edge_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). Most elements should return true.
Parameters
edge_idThe edge number, starting at 1, according to the B2000++ conventions.
field_nameThe field name.

Reimplemented in b2000::ElementCFDHeat2DFieldL2Constant< FLUX, AXISYMMETRIC >.

◆ edge_field_order()

virtual int b2000::Element::edge_field_order ( const int  edge_id,
const std::string &  field_name 
)
inlinevirtual
Returns
The integration order along the edge. This depends on the shape functions. Linear elements return 1, quadratic elements return 2, etc.
Parameters
edge_idThe edge number, starting at 1, according to the B2000++ conventions.
field_nameThe field name. Usually, the integration order is independent of the field name.

Reimplemented in b2000::ElementCFDHeat2DFieldL2Constant< FLUX, AXISYMMETRIC >.

◆ edge_field_polynomial_sub_edge()

virtual int b2000::Element::edge_field_polynomial_sub_edge ( const int  edge_id,
const std::string &  field_name,
b2linalg::Vector< double, b2linalg::Vdense > &  sub_nodes 
)
inlinevirtual

Sub-divide the edge into polynomial sub-edges and compute the edge-internal coordinates of the nodes at the sub-edge boundaries. Most elements should put [-1, +1] (one sub-edge).

Parameters
edge_idThe edge number, starting at 1, according to the B2000++ conventions.
field_nameThe field name.
sub_nodesA vector containing the edge-internal coordinates of the sub-nodes. The coordinates must be strictly increasing, and the first sub-node must be at -1, and the last sub-node must be at +1.
Returns
The number of sub-edges, it must be equal to sub_nodes.size() - 1. Most elements should return 1.

Reimplemented in b2000::ElementCFDHeat2DFieldL2Constant< FLUX, AXISYMMETRIC >.

◆ edge_geom()

virtual void b2000::Element::edge_geom ( const int  edge_id,
const double  internal_coor,
b2linalg::Vector< double > &  geom,
b2linalg::Vector< double > &  d_geom_d_icoor 
)
inlinevirtual

Compute the edge geometry at a specified integration point.

Parameters
edge_idThe edge number, starting at 1, according to the B2000++ conventions (input).
internal_coorThe edge-internal coordinate (from -1 to +1) 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 coordinate; this is the edge tangential vector (vector, output).

Reimplemented in b2000::ElementCFDHeat2DFieldL2Constant< FLUX, AXISYMMETRIC >.

◆ end()

iterator b2000::Element::end ( ) const
inline

Returns the 'end' iterator.

◆ face_field_linear_on_dof()

virtual bool b2000::Element::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 in b2000::ElementCFD3DField< SHAPE, EULER >, and b2000::SuperElement< T, MATRIX_TYPE >.

◆ face_field_order()

virtual int b2000::Element::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 in b2000::SuperElement< T, MATRIX_TYPE >, and b2000::ElementCFD3DField< SHAPE, EULER >.

◆ face_field_polynomial_sub_face()

virtual int b2000::Element::face_field_polynomial_sub_face ( const int  face_id,
const std::string &  field_name,
b2linalg::Matrix< double, b2linalg::Mrectangle > &  sub_nodes,
std::vector< 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 in b2000::SuperElement< T, MATRIX_TYPE >.

◆ face_geom()

virtual void b2000::Element::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 in b2000::SuperElement< T, MATRIX_TYPE >, and b2000::ElementCFD3DField< SHAPE, EULER >.

◆ get_constraint_info()

virtual std::pair< int, VariableInfo > b2000::Element::get_constraint_info ( )
inlinevirtual

Must be overridden for elements with constraints. Return information about the constraints specified by the element.

todo: not all elements override this function.

In B2000++, elements are allowed, in addition to evaluate first and second variations (e.g. internal force vector and stiffness matrix), to specify constraints between dof's. This is useful for e.g. rigid-body elements or interface elements.

Returns
a pair (first, second), where 'first' is the number of constraints, and 'second' specifies the variation of the constraints w.r.t. the dof's.

◆ get_dof_numbering()

virtual void b2000::Element::get_dof_numbering ( b2linalg::Index &  dof_numbering)
inlinevirtual
Deprecated:
This function is deprecated.

◆ get_elem_type_number()

int b2000::Element::get_elem_type_number ( ) const
inline
Returns
element Nr. (elno) of an element in ELEMENT_PARAMETER table in MemCom database.
See also
set_elem_type_number()

◆ get_global_dof_numbering()

std::pair< size_t, size_t > b2000::Element::get_global_dof_numbering ( ) const
inlineoverridevirtual

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 override
Definition b2element.H:302

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

Reimplemented from b2000::DegreesOfFreedom.

Reimplemented in b2000::SuperElement< T, MATRIX_TYPE >, and b2000::ElementCFDHeat2DFieldL2Constant< FLUX, AXISYMMETRIC >.

◆ get_id()

size_t b2000::Element::get_id ( ) const
inline
Returns
The internal element identifier assigned by the Domain object.
See also
set_id()

◆ get_nodes()

virtual std::pair< int, Node *const * > b2000::Element::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 in b2000::SuperElement< T, MATRIX_TYPE >, b2000::ElementCFD3DField< SHAPE, EULER >, and b2000::ElementCFDHeat2DFieldL2Constant< FLUX, AXISYMMETRIC >.

◆ get_number_of_dof()

int b2000::Element::get_number_of_dof ( ) const
inlineoverridevirtual

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::DegreesOfFreedom.

Reimplemented in b2000::SuperElement< T, MATRIX_TYPE >, and b2000::ElementCFDHeat2DFieldL2Constant< FLUX, AXISYMMETRIC >.

◆ get_object_name()

const std::string & b2000::Element::get_object_name ( ) const
inlineoverridevirtual
Returns
The name of the element assigned by the Domain object.
See also
set_object_name()

Reimplemented from b2000::Object.

◆ get_property()

virtual const ElementProperty * b2000::Element::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 in b2000::SuperElement< T, MATRIX_TYPE >, and b2000::ElementCFD3DField< SHAPE, EULER >.

◆ get_state_buffer()

virtual std::pair< int *, double * > b2000::Element::get_state_buffer ( std::pair< int *, double * >  buffers) const
inlinevirtual

Write the element state into buffers.

Parameters
buffersA pair of pointers to an integer array and to a floating-point array. The buffer sizes of the arrays correspond to get_state_buffer_size().
Returns
A pair of pointers to an integer array and to a floating-point array, but incremented by get_state_buffer_size().

◆ get_state_buffer_size()

virtual std::pair< int, int > b2000::Element::get_state_buffer_size ( ) const
inlinevirtual

Return the state buffers size in integer and double to store the state data of the element.

◆ get_subelements()

virtual std::pair< size_t, Element ** > b2000::Element::get_subelements ( )
inlinevirtual

Returns an array of pointers to sub-elements. By default, no sub-elements are present.

◆ get_value_info()

virtual const std::vector< VariableInfo > b2000::Element::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
virtual const std::vector< VariableInfo > get_value_info() const
Definition b2element.H:670

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 in b2000::SuperElement< T, MATRIX_TYPE >.

◆ init()

virtual void b2000::Element::init ( Model model)
inlinevirtual

Initialisation of the element. This method is called by the Domain object, after set_nodes() and set_property(). This implementation performs no operation.

This allows for setting up private data structures that depend on the nodes and on the property.

Examples

  • Stress elements: Determine the element-internal positions of the integration points and compute the Jacobian matrices of the initial configuration at these integration points.
  • Elements with many dof's: Set-up sub-elements.

◆ set_additional_properties()

virtual void b2000::Element::set_additional_properties ( const RTable rtable)
inlinevirtual

Set additional element properties. Work-around until the new element properties are implemented. This implementation performs no operation.

Parameters
rtableA relational table containing the element definition as read from database.

◆ set_elem_type_number()

void b2000::Element::set_elem_type_number ( int  elno)
inline

This function is used to set the element type Nr. read from MemCom database. In the database file, the element type Nr. is stored as ITYP under ETAB table. ITYP is set in b2ip_branch.C and elno is stored as ITYP in b2ip_elements.C. It is used to create the element object.

◆ set_global_dof_numbering()

size_t b2000::Element::set_global_dof_numbering ( size_t  index)
inlineoverridevirtual

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) override
Definition b2element.H:279

The variable internal_dof_index is a member of the element class.

Reimplemented from b2000::DegreesOfFreedom.

Reimplemented in b2000::SuperElement< T, MATRIX_TYPE >, and b2000::ElementCFDHeat2DFieldL2Constant< FLUX, AXISYMMETRIC >.

◆ set_id()

void b2000::Element::set_id ( size_t  id_)
inline

This function is called by the Domain object. In B2000++, internal element identifiers start from 0 and are consecutive.

◆ set_nodes()

virtual void b2000::Element::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
virtual void set_nodes(std::pair< int, Node *const * > nodes)
Definition b2element.H:409
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 in b2000::SuperElement< T, MATRIX_TYPE >, b2000::ElementCFD3DField< SHAPE, EULER >, and b2000::ElementCFDHeat2DFieldL2Constant< FLUX, AXISYMMETRIC >.

◆ set_object_name()

void b2000::Element::set_object_name ( const std::string &  name)
inline

This function is called by the Domain object. In B2000++, the naming is "branch.eid", where 'branch' is the branch number and 'eid' is the element number as specified in the MDL input file and corresponding to the EID field in the ETAB.branch dataset of the database.

◆ set_property()

virtual void b2000::Element::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
virtual void set_property(ElementProperty *property)
Definition b2element.H:465

Reimplemented in b2000::SuperElement< T, MATRIX_TYPE >, and b2000::ElementCFD3DField< SHAPE, EULER >.

◆ set_state_buffer()

virtual std::pair< const int *, const double * > b2000::Element::set_state_buffer ( std::pair< const int *, const double * >  buffers)
inlinevirtual

Initialize the element state from buffers (for restart). Actual implementations may need to perform further initializations of the internal state (e.g. re-computation of Jacobian matrices).

Parameters
buffersA pair of pointers to an integer array and to a floating-point array. The buffer sizes of the arrays correspond to get_state_buffer_size().
Returns
A pair of pointers to an integer array and to a floating-point array, but incremented by get_state_buffer_size().

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