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

A in-core implementation of the b2dbv3::Domain. All the Element, Node and ElementProperty objects are loaded in memory at theInitialisation of the object with the database. More...

#include "b2domain_database.H"

Inheritance diagram for b2000::b2dbv3::Domain:
Collaboration diagram for b2000::b2dbv3::Domain:

Public Member Functions

void set_model (b2000::Model &model) override
 
void set_case (b2000::Case &case_) override
 
void set_subcase_id (const int subcase_id_) override
 
bool have_temperature () const override
 
size_t get_number_of_dof () const override
 
size_t get_number_of_elements () const override
 
size_t get_number_of_elements_and_subelements () const override
 
size_t get_number_of_nodes () const override
 
const std::vector< Element::VariableInfoget_value_info () const override
 
node_iterator get_node_iterator () override
 
element_iterator get_element_iterator () override
 
void save_field (const std::string &name, const b2linalg::Vector< double, b2linalg::Vdense_constref > &dof) override
 
void save_field (const std::string &name, const b2linalg::Vector< b2000::csda< double >, b2linalg::Vdense_constref > &dof) override
 
void save_field (const std::string &name, const b2linalg::Vector< std::complex< double >, b2linalg::Vdense_constref > &dof) override
 
void save_state (const std::string &state_id) override
 
void load_state (const std::string &state_id) override
 
b2linalg::SparseMatrixConnectivityType get_dof_connectivity_type () const override
 
Nodeget_node (const std::string &node_name) override
 
Elementget_element (const std::string &element_name) override
 
void get_adjacent_elements (std::pair< int, const Node *const * > node_list, std::vector< AdjacentElement > &adjacent_elements) override
 
void get_elements_of_same_type_near (const Element *e, double dist, std::vector< Element * > &element_list) override
 
void get_slave_node (std::vector< std::pair< const Node *, const Node * > > &slave_node) override
 
void set_dof (const b2linalg::Matrix< double, b2linalg::Mrectangle_constref > &dof) override
 
void set_dof (const b2linalg::Matrix< b2000::csda< double >, b2linalg::Mrectangle_constref > &dof) override
 
void set_dof (const b2linalg::Matrix< std::complex< double >, b2linalg::Mrectangle_constref > &dof) override
 
const b2linalg::Matrix< double, b2linalg::Mrectangle_constref > & get_dof (double) override
 
Elementget_parent_element_and_nodes_mapping (b2000::Domain &parent_domain, const b2000::Element &element, b2linalg::Matrix< double > &parent_nodes_internal_coor) override
 
Elementget_parent_element_mapping (b2000::Domain &parent_domain, const b2000::Node &node, b2linalg::Vector< double > &parent_internal_coor) override
 
const double * get_node_local_referential_double (const Node &node) const override
 
const b2000::csda< double > * get_node_local_referential_csda (const Node &node) const override
 
- Public Member Functions inherited from b2000::Domain
template<typename T >
const T * get_node_local_referential (const Node &node) const
 
virtual void get_adjacent_elements (const Node *node, std::vector< AdjacentElement > &adjacent_elements)
 
std::vector< std::reference_wrapper< Element > > & GetElementContainer ()
 Return private member domain_element_container_ holding all elements of this domain.
 
void FillElementContainer ()
 Fill domain_element_container_ with all elements of this domain.
 
- Public Member Functions inherited from b2000::Object
virtual const std::string & get_object_name () const
 
virtual ~Object ()
 

Additional Inherited Members

- Public Types inherited from b2000::Domain
- Static Public Attributes inherited from b2000::Object
static ObjectType type
 

Detailed Description

A in-core implementation of the b2dbv3::Domain. All the Element, Node and ElementProperty objects are loaded in memory at theInitialisation of the object with the database.

Note
that the processing of the element types defined inside enum element_type_t takes place using a binary interpretation of the different values. Do not assume a simple sequential ordering when adding further element types!

Member Function Documentation

◆ get_adjacent_elements()

void b2000::b2dbv3::Domain::get_adjacent_elements ( std::pair< int, const Node *const * >  node_list,
std::vector< AdjacentElement > &  adjacent_elements 
)
inlineoverridevirtual

Get the list of elements that are adjacent to all the Node instances given as argument.

Parameters
node_listA pair (first, second) with first being the number of nodes and second a pointer to an array of Node pointers (input).
adjacent_elementsA vector of containing information on adjacent elements. Each entry corresponds to one Element.

The AdjacentElement::internal_node_id_list contains the node-internal ID's of the nodes given in node_list that

Implements b2000::Domain.

◆ get_dof()

const b2linalg::Matrix< double, b2linalg::Mrectangle_constref > & b2000::b2dbv3::Domain::get_dof ( double  )
inlineoverridevirtual

Get a subset of the current degrees-of-freedom.

Implements b2000::Domain.

◆ get_dof_connectivity_type()

b2linalg::SparseMatrixConnectivityType b2000::b2dbv3::Domain::get_dof_connectivity_type ( ) const
inlineoverridevirtual
Returns
The type of connectivity of FE elements present in the domain (1D, 2D, 3D, etc.). This function is called by the Solver instance (for the sparse matrix solver). Depending on this value, the Solver decides what graph partitioning algorithm to use for the symbolic factorization.

Implements b2000::Domain.

◆ get_element()

Element * b2000::b2dbv3::Domain::get_element ( const std::string &  element_name)
inlineoverridevirtual
Returns
A pointer to the Element instance with Element::get_object_name() == element_name, or 0 if such a Element instance does not exist.
Note
The computational cost of this operation can be O(n) where n is the number of elements in the FE domain.

Reimplemented from b2000::Domain.

◆ get_element_iterator()

element_iterator b2000::b2dbv3::Domain::get_element_iterator ( )
inlineoverridevirtual
Returns
A std::unique_ptr to ElementIterator to iterate over all elements of the domain.

Implements b2000::Domain.

◆ get_elements_of_same_type_near()

void b2000::b2dbv3::Domain::get_elements_of_same_type_near ( const Element e,
double  dist,
std::vector< Element * > &  element_list 
)
overridevirtual

Get the list of elements near a given element and of the same type. This function is useful e.g. for contact overlay elements.

Parameters
eA pointer to the Element instance (input).
distThe maximum search distance (input).
element_listA vector of pointers to Element instances (output).

The list of returned elements contains all elements of the same type as Element e, and for with the distance with e is less or equal dist. The list may also contain some elements whose distance to e is slightly larger than dist.

Implements b2000::Domain.

◆ get_node()

Node * b2000::b2dbv3::Domain::get_node ( const std::string &  node_name)
inlineoverridevirtual
Returns
A pointer to the Node instance with Node::get_object_name() == node_name, or 0 if such a Node instance does not exist.
Note
The computational cost of this operation can be O(n) where n is the number of nodes in the FE domain.

Reimplemented from b2000::Domain.

◆ get_node_iterator()

node_iterator b2000::b2dbv3::Domain::get_node_iterator ( )
inlineoverridevirtual

Return a std::unique_ptr to NodeIterator to iterate over all nodes of the FE domain.

Implements b2000::Domain.

◆ get_node_local_referential_csda()

const b2000::csda< double > * b2000::b2dbv3::Domain::get_node_local_referential_csda ( const Node node) const
inlineoverridevirtual

Must be overwritten by implementations of Domain.

Implements b2000::Domain.

◆ get_node_local_referential_double()

const double * b2000::b2dbv3::Domain::get_node_local_referential_double ( const Node node) const
inlineoverridevirtual

Must be overwritten by implementations of Domain.

Implements b2000::Domain.

◆ get_number_of_dof()

size_t b2000::b2dbv3::Domain::get_number_of_dof ( ) const
inlineoverridevirtual
Returns
The number of degrees-of-freedom in the domain.

Implements b2000::Domain.

◆ get_number_of_elements()

size_t b2000::b2dbv3::Domain::get_number_of_elements ( ) const
inlineoverridevirtual
Returns
The number of elements in the domain.

Implements b2000::Domain.

◆ get_number_of_elements_and_subelements()

size_t b2000::b2dbv3::Domain::get_number_of_elements_and_subelements ( ) const
inlineoverridevirtual
Returns
The number of elements in the domain plus the numbers of their subelements (if any).

Implements b2000::Domain.

◆ get_number_of_nodes()

size_t b2000::b2dbv3::Domain::get_number_of_nodes ( ) const
inlineoverridevirtual
Returns
The number of nodes in the domain.

Implements b2000::Domain.

◆ get_parent_element_and_nodes_mapping()

Element * b2000::b2dbv3::Domain::get_parent_element_and_nodes_mapping ( b2000::Domain parent_domain,
const b2000::Element element,
b2linalg::Matrix< double > &  parent_nodes_internal_coor 
)
inlineoverridevirtual

Get the Element of the parent domain containing an element belonging to this domain (the local domain).

Parameters
parent_domainThe parent domain to search.
elementThe Element instance of the this domain (= the local element).
parent_nodes_internal_coorThe element-internal coordinates of the local element inside the parent element.
Returns
A pointer to the parent Element of the parent Domain or 0, if no parent element exists.

Reimplemented from b2000::Domain.

◆ get_parent_element_mapping()

Element * b2000::b2dbv3::Domain::get_parent_element_mapping ( b2000::Domain parent_domain,
const b2000::Node node,
b2linalg::Vector< double > &  parent_internal_coor 
)
inlineoverridevirtual

Get the first Element of the parent domain containing a node belonging to this domain (the local domain.

Parameters
parent_domainThe parent domain to search.
parent_domainThe parent domain to search.
nodeThe Node instance of the this domain (= the local node).
parent_internal_coorThe element-internal coordinates of the local node inside the parent element.
Returns
A pointer to the parent Element of the parent Domain or 0, if no parent element exists.

Reimplemented from b2000::Domain.

◆ get_slave_node()

void b2000::b2dbv3::Domain::get_slave_node ( std::vector< std::pair< const Node *, const Node * > > &  slave_node)
overridevirtual

In a domain, coinciding nodes of different branches may be merged together. In this case, one of the nodes becomes the master node and the other node becomes the slave node.

Parameters
slave_nodeA vector of (slave,master) pairs of pointers to Node instances (output).

Implements b2000::Domain.

◆ get_value_info()

const std::vector< Element::VariableInfo > b2000::b2dbv3::Domain::get_value_info ( ) const
inlineoverridevirtual
Returns
A vector containing the type of response for each variable of the whole domain (all Element instances). See also Element::get_value_info().

Implements b2000::Domain.

◆ have_temperature()

bool b2000::b2dbv3::Domain::have_temperature ( ) const
inlineoverridevirtual

Whether temperature conditions affecting the right-hand side through the element response (thermal expansion) are present for the current case and subcase.

Implements b2000::Domain.

◆ load_state()

void b2000::b2dbv3::Domain::load_state ( const std::string &  state_id)
inlineoverridevirtual

Load the state of the FE domain with the identifier state_id. This function is called by the Solver instance upon restart.

Implements b2000::Domain.

◆ save_field() [1/3]

void b2000::b2dbv3::Domain::save_field ( const std::string &  name_id,
const b2linalg::Vector< b2000::csda< double >, b2linalg::Vdense_constref > &  dof 
)
inlineoverridevirtual

Save a csda-valued field defined on all degrees-of-freedom. This function is called by the Solver instance.

Parameters
name_idThe name of the field.
dofA vector containing values for all degrees-of-freedom present in the FE domain.

Implements b2000::Domain.

◆ save_field() [2/3]

void b2000::b2dbv3::Domain::save_field ( const std::string &  name_id,
const b2linalg::Vector< double, b2linalg::Vdense_constref > &  dof 
)
inlineoverridevirtual

Save a real-valued field defined on all degrees-of-freedom. This function is called by the Solver instance. Such fields may be, e.g. for stress analysis, the displacement field ("DISP"), or the field containing the summed-up internal forces ("RCFO").

Parameters
name_idThe name of the field.
dofA vector containing values for all degrees-of-freedom present in the FE domain.

Implements b2000::Domain.

◆ save_field() [3/3]

void b2000::b2dbv3::Domain::save_field ( const std::string &  name_id,
const b2linalg::Vector< std::complex< double >, b2linalg::Vdense_constref > &  dof 
)
inlineoverridevirtual

Save a complex-valued field defined on all degrees-of-freedom. This function is called by the Solver instance.

Parameters
name_idThe name of the field.
dofA vector containing values for all degrees-of-freedom present in the FE domain.

Implements b2000::Domain.

◆ save_state()

void b2000::b2dbv3::Domain::save_state ( const std::string &  state_id)
inlineoverridevirtual

Save the state of the FE domain with the identifier state_id. This function is called by the Solver instance.

Implements b2000::Domain.

◆ set_case()

void b2000::b2dbv3::Domain::set_case ( b2000::Case case_)
overridevirtual

Set the current Case instance. This function is called by the Model instance.

Implements b2000::Domain.

◆ set_dof() [1/3]

void b2000::b2dbv3::Domain::set_dof ( const b2linalg::Matrix< b2000::csda< double >, b2linalg::Mrectangle_constref > &  dof)
inlineoverridevirtual

Set the current degrees-of- freedom.

Implements b2000::Domain.

◆ set_dof() [2/3]

void b2000::b2dbv3::Domain::set_dof ( const b2linalg::Matrix< double, b2linalg::Mrectangle_constref > &  dof)
inlineoverridevirtual

Set the current degrees-of-freedom.

Implements b2000::Domain.

◆ set_dof() [3/3]

void b2000::b2dbv3::Domain::set_dof ( const b2linalg::Matrix< std::complex< double >, b2linalg::Mrectangle_constref > &  dof)
inlineoverridevirtual

Set the current degrees-of- freedom.

Implements b2000::Domain.

◆ set_model()

void b2000::b2dbv3::Domain::set_model ( b2000::Model model)
overridevirtual

Set a reference to the Model instance. This function is called by the Model instance.

Implements b2000::Domain.

◆ set_subcase_id()

void b2000::b2dbv3::Domain::set_subcase_id ( const int  subcase_id_)
inlineoverridevirtual

Set the subcase ID.

Implements b2000::Domain.


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