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

#include "b2element.H"

Inheritance diagram for b2000::TypedElement< T >:
Collaboration diagram for b2000::TypedElement< T >:

Public Member Functions

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)
 
Old integration functions

The field_face_integration() and field_volume_integration() functions are currently in use by the implementation of the NaturalBoundaryCondition object implementing traction and body forces.

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)
 
- Public Member Functions inherited from b2000::Element
virtual const std::vector< VariableInfoget_value_info () const
 
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)
 
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)
 
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)
 
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 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)
 
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 ()
 

Bounding box

This function does not have to be overwritten for most implementations. Their purpose is to compute the bounding box of a deformed configuration, aligned with the axis of the branch-global coordinate system. This functionality can be used e.g. for contact overlay elements.

The bounding box may be larger than the element, but not smaller.

For elements using Lagrange shape functions, it is sufficient to compute the minimum and maximum deformed coordinates of all nodes in the given dimensions.

static type_t type
 Element interface typed by the 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.
 

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::Object
static ObjectType type
 

Detailed Description

template<typename T>
class b2000::TypedElement< T >

Element interface typed by the dof type.

Member Function Documentation

◆ AssembleElementEffectiveMatrix()

template<typename T >
template<typename MATRIX_FORMAT >
b2linalg::Matrix< T, typename MATRIX_FORMAT::dense > b2000::TypedElement< T >::AssembleElementEffectiveMatrix ( const solver::TypedSolver< T, MATRIX_FORMAT > &  solver,
b2linalg::Index &  index 
)
inline

Assembles the element effective matrix.

Template Parameters
MATRIX_FORMATEither b2000::b2linalg::Msymmetric or b2000::b2linalg::Munsymmetric.
Parameters
[in]solverSolver instance chosen by the user.
[out]indexLocal to global node mapping.
Returns
k_eff The element effective matrix.
Note
This is an adapter interface for the new solver design.

◆ AssembleElementEffectiveSystem()

template<typename T >
template<typename MATRIX_FORMAT >
void b2000::TypedElement< T >::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 
)
inline

Assembles the element effective matrix and vector.

Template Parameters
MATRIX_FORMATEither b2000::b2linalg::Msymmetric or b2000::b2linalg::Munsymmetric.
Parameters
[in]solverSolver instance chosen by the user.
[out]indexLocal to global node mapping.
[out]k_effThe element effective matrix.
[out]f_effThe element effective vector.
Note
This is an adapter interface for the new solver design.

◆ AssembleElementEffectiveVector()

template<typename T >
template<typename MATRIX_FORMAT >
b2linalg::Vector< T, b2linalg::Vdense > b2000::TypedElement< T >::AssembleElementEffectiveVector ( const solver::TypedSolver< T, MATRIX_FORMAT > &  solver,
b2linalg::Index  index 
)
inline

Assembles the element effective vector.

Template Parameters
MATRIX_FORMATEither b2000::b2linalg::Msymmetric or b2000::b2linalg::Munsymmetric.
Parameters
[in]solverinstance chosen by the user.
[out]indexLocal to global node mapping.
Returns
f_eff The element effective vector.
Note
This is an adapter interface for the new solver design.

◆ body_field_value()

template<typename T >
virtual void b2000::TypedElement< T >::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 
)
inlinevirtual

Interpolate the value of the field and its derivatives at a specified face integration point for real-valued elements.

Parameters
field_nameThe field name.
internal_coorThe element-internal coordinates (input).
dofThe solution vector (input).
timeThe current time (input).
valueThe value of the field interpolated at the given element-internal coordinates for the given dof and time.
d_value_d_icoorThe derivatives of the interpolated field value w.r.t. the element-internal coordinates (matrix, output).
dof_numberingArray of the global dof indices on which the interpolated field value depends (vector, output).
d_value_d_dofThe derivatives of the interpolated field values w.r.t. the dof's (matrix, output).

◆ ComputeElementError()

template<typename T >
template<typename MATRIX_FORMAT >
T b2000::TypedElement< T >::ComputeElementError ( const solver::TypedSolver< T, MATRIX_FORMAT > &  solver)
inline

TODO.

Template Parameters
MATRIX_FORMATEither b2000::b2linalg::Msymmetric or b2000::b2linalg::Munsymmetric.
Parameters
[in]solverSolver instance chosen by the user.
Returns
TODO
Note
This is an adapter interface for the new solver design.

◆ ComputeElementGradient()

template<typename T >
template<typename MATRIX_FORMAT >
void b2000::TypedElement< T >::ComputeElementGradient ( const solver::TypedSolver< T, MATRIX_FORMAT > &  solver,
GradientContainer *const  gradient_container 
)
inline

Compute the gradients for post-processing.

Template Parameters
MATRIX_FORMATEither b2000::b2linalg::Msymmetric or b2000::b2linalg::Munsymmetric.
Parameters
[in]solverinstance chosen by the user.
[in]gradient_containerconst pointer to non-const gradient solution.
Note
This is an adapter interface for the new solver design.

◆ edge_field_value()

template<typename T >
virtual void b2000::TypedElement< T >::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 
)
inlinevirtual

Interpolate the value of the field and its derivatives at a specified edge integration point for real-valued elements.

Parameters
edge_idThe edge number, starting at 1, according to the B2000++ conventions (input).
field_nameThe field name.
internal_coorThe edge-internal coordinate (from -1 to +1) of the integration point (input).
dofThe solution vector (input).
timeThe current time (input).
valueThe value of the field interpolated at the given edge-internal coordinate for the given dof and time.
d_value_d_icoorThe derivatives of the interpolated field value w.r.t. the internal coordinate (vector, output).
dof_numberingArray of the global dof indices on which the interpolated field value depends (vector, output).
d_value_d_dofThe derivatives of the interpolated field values w.r.t. the dof's (matrix, output).
d_value_d_dof_dep_colArray of the global dof indices on which the interpolated field value depends for the current deformation (vector, output).

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

◆ face_field_value()

template<typename T >
virtual void b2000::TypedElement< T >::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 
)
inlinevirtual

Interpolate the value of the field and its derivatives at a specified face integration point for real-valued elements.

Parameters
face_idThe face number, starting at 1, according to the B2000++ conventions (input).
field_nameThe field name.
internal_coorThe face-internal coordinates (input).
dofThe solution vector (input).
timeThe current time (input).
valueThe value of the field interpolated at the given face-internal coordinates for the given dof and time.
d_value_d_icoorThe derivatives of the interpolated field value w.r.t. the face-internal coordinates (matrix, output).
dof_numberingArray of the global dof indices on which the interpolated field value depends (vector, output).
d_value_d_dofThe derivatives of the interpolated field values w.r.t. the dof's (matrix, output).
d_value_d_dof_dep_colArray of the global dof indices on which the interpolated field value depends for the current deformation (vector, output).

◆ field_edge_integration()

template<typename T >
virtual void b2000::TypedElement< T >::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 
)
inlinevirtual

Field integration on the (potentially deformed) edge of the element.

Parameters
dofThe current dof used by the element (input).
dofdotThe current velocity dof used by the element (input).
dofdotdotThe current acceleration dof used by the element (input).
fThe field to integrate.
edgeThe edge over which to integrate the field.
dof_numberingThe dof numbering of the value returned by this function.
discretised_fieldThe value at dofs of the integration of the field on the edge of the element.
d_discretised_field_d_dofThe derivative w.r.t. the dof's of the value at dofs of the integration of the field on the edge element (matrix, output).

If a null index/vector/matrix is passed for one of the output parameters, this output shall not be computed.

◆ field_face_integration()

template<typename T >
virtual void b2000::TypedElement< T >::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 
)
inlinevirtual

Field integration on the (potentially deformed) surface of the element.

Parameters
dofThe current dof's used by the element (input).
dofdotThe current velocity dof's used by the element (input).
dofdotdotThe current acceleration dof's used by the element (input).
fieldThe field to integrate.
faceThe face on witch to integrate the field.
dof_numberingThe dof numbering of the value returned by this function.
discretised_fieldThe value at dofs of the integration of the field on the surface of the element.
d_discretised_field_d_dofThe derivative w.r.t. the dof of the value at dofs of the integration of the field on the surface element (matrix, output).

If a null index/vector/matrix is passed for one of the output parameters, this output shall not be computed.

◆ field_volume_integration()

template<typename T >
virtual void b2000::TypedElement< T >::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 
)
inlinevirtual

Field integration on the (potentially deformed) volume element.

Parameters
dofThe current dof used by the element (input).
dofdotThe current velocity dof used by the element (input).
dofdotdotThe current acceleration dof used by the element (input).
fThe field to integrate.
dof_numberingThe dof numbering of the value returned by this function.
discretised_fieldThe value at dofs of the integration of the field on the volume element.
d_discretised_field_d_dofThe derivatives w.r.t. the dof's of the value at dofs of the integration of the field on the volume element (matrix, output).

If a null index/vector/matrix is passed for one of the output parameters, this output shall not be computed.

◆ get_constraint()

template<typename T >
virtual void b2000::TypedElement< T >::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 
)
inlinevirtual

Must be overridden for real-valued elements with constraints: Compute the constraints of the element and their derivatives.

todo: not overridden by all elements

Parameters
modelThe model object to which the Element belongs (input).
linearA flag to indicate if the value returned must be linear in function of the dof and the time derivative of the dof (input).
timeThe current time (input).
dofThe values and the time derivatives of all dof's present in the FE model. The first column of the matrix dof[0] are the dofs. The i.th column of the matrix (dof[i]) is the i.th time derivative of dof[0] (matrix, input).
dof_numberingArray of the global dof indices on which the Element depends (vector, output).
constraintThe evaluated constraints, in function of dof and time. The size of the constraint vector must be equal to the number of constraints (vector, output).
trans_d_constraint_d_dofThe derivatives of constraint w.r.t. dof[0]. If trans_d_constraint_d_dof.is_null(), then the derivatives shall not be computed (transposed sparse matrix, output).
d_constraint_d_timeThe derivatives of constraints w.r.t. time (vector, output).

◆ get_value() [1/2]

template<typename T >
virtual void b2000::TypedElement< T >::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 
)
inlinevirtual

Compute the internal forces and their derivatives of the element.

This is an overload where "d_value_d_dof" is supposed to be symmetric. For more information on individual input parameters see get_value() above.

◆ get_value() [2/2]

template<typename T >
virtual void b2000::TypedElement< T >::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 
)
inlinevirtual

Compute the internal forces and their derivatives of the element.

This is an overlaod where d_value_d_dof is supposed to be unsymmetric.

Note
This is supposed to be overwritten by actual element classes.

todo Not all deriving classes implement this function and the implementation here is actually used for those elements.

Parameters
[in]modelThe model object to which the Element belongs (input).
[in]linearA flag to indicate if the value returned must be linear in function of the dof and the time derivative of the dof (input).
[in]equilibrium_solutionUsed by the elements with path dependent materials only. If true, the given dof and time derivatives represent an equilibrium solution (input).
[in]timeThe current time (input).
[in]step_sizeThe current time step size of the incremental solution solver (input).
[in]dofThe values and the time derivatives of all dofs present in the FE model. The first column of the matrix dof[0] are the dofs. The i.th column of the matrix (dof[i]) is the i.th time derivative of dof[0] (matrix, input).
[out]gradient_containerTo store the gradient data computed by the Element and/or ElementProperty object (output).
[in]solver_hintsCan be used to add hints for the Solver (e.g. for damage material models, fracture elements, etc.) (output).
[in]dof_numberingArray of the global dof indices on which the internal forces depend (vector, output).
[in,out]valueThe first variation (internal forces), in function of dof and time (vector, output).
[in]d_value_d_dof_flagsA vector of bool indicating for which derivatives of dof the second variation should be computed.
[in,out]d_value_d_dofThe second variation, that is, d_value_d_dof[i] contains the derivatives of value w.r.t. dof[i]. If d_value_d_dof[i].is_null(), then the matrix shall not be computed (vector of packed matrices, output).
[out]d_value_d_timeThe derivatives of value w.r.t. time (vector, output).

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