GranOO  3.0
A robust and versatile workbench to build 3D dynamic simulations based on the Discrete Element Method
Public Member Functions | Protected Attributes | Private Member Functions | Private Attributes | Friends | List of all members
GranOO3::DEM::XBeam Class Reference

A bond XBeam that works in tension, compression, bending and torsion very similar to the Beam bond. More...

#include <XBeam.hpp>

Inheritance diagram for GranOO3::DEM::XBeam:
Inheritance graph
[legend]

Public Member Functions

 XBeam (Element &el1, Element &el2)
 constructor More...
 
 XBeam (Element &el1, Element &el2, double young_mod, double radius, double poisson_ratio, double max_stress=0.)
 constructor More...
 
 XBeam (Element &el1, Element &el2, double young_mod, double damping, double radius, double poisson_ratio, double max_stress, bool nl_geom=false)
 constructor More...
 
virtual ~XBeam ()
 
virtual void compute_load ()
 Compute the reaction force and torque of the beam. More...
 
virtual double compute_critical_time_step () const
 compute the critical time of the bond More...
 
virtual void stiffness_matrix (Eigen::MatrixXd &kMatrix)
 Assigns the global stiffness matrix of the XBeam kMatrix. More...
 
virtual void update_dof ()
 Builds the list of global Degrees Of Freedom (DOFs), required by the solver to assemble the global (system-level) stiffness matrix and the displacement and forces vectors. More...
 
virtual std::string info () const
 Displays some useful info in the terminal. More...
 
virtual void draw ()
 OpenGL draw of the bond
More...
 
- Public Member Functions inherited from GranOO3::DEM::Bond
 Bond (Element &el1, Element &el2, bool contact=false)
 Constructor. More...
 
virtual ~Bond ()=0
 Destructor. More...
 
Elementget_item1 ()
 
const Elementget_item1 () const
 
Elementget_item2 ()
 
const Elementget_item2 () const
 
Elementitem1 ()
 
Elementitem2 ()
 
double get_normal_force_intensity () const
 
double get_tangential_force_intensity () const
 
Elementbonded_element (const Element &)
 get the bonded element More...
 
const Elementget_bonded_element (const Element &) const
 same as bonded_element(const Element&) (provided for convenience) More...
 
virtual double get_linear_stiffness () const
 compute the linear stiffness of the bond More...
 
virtual double get_angular_stiffness () const
 similar as get_linear_stiffness() for rotation More...
 
void apply_linear_damping_factor ()
 compute and apply the linear (translation) damping factor of the bond More...
 
void apply_angular_damping_factor ()
 similar as apply_linear_damping_factor() for rotation More...
 
void apply_load ()
 apply the computed reaction force and torque to the bonded Element (_element1 and _element2) More...
 
void disable (bool manage_slave=true)
 disable the bond More...
 
void disable_fast_mode (bool manage_slave=true)
 simply disable the bond by switching the _disabled attribute to true without doing the lot of things done by the disable() method More...
 
void enable (bool manage_slave=true)
 simply enable the bond by switching the _disabled attribute to false More...
 
bool is_disabled () const
 check if the bond is enable or not More...
 
void crack (bool manage_slave=true)
 
void uncrack (bool manage_slave=true)
 
bool is_cracked () const
 
void init_local_frame (bool update_initial_param=false)
 initialize the bond local frame More...
 
void update_local_frame ()
 update the bond local frame More...
 
void update_current_length ()
 update the current value of the bond length (see _current_length attribute) More...
 
void build_voronoi_bond ()
 build a voronoi bond from the current bond More...
 
const std::vector< unsigned int > & get_dof () const
 get the dof vector More...
 
const Geom::Pointinitial_pos (const Element &) const
 get the initial position of the el element More...
 
Geom::Vector get_disp (const Element &) const
 get the displacement of the el element More...
 
Geom::Vector get_disp1 () const
 get the displacement of the first bonded element (_element1) More...
 
Geom::Vector get_disp2 () const
 similar as get_disp1() More...
 
void swap (Element &from, Element &to)
 swap the bond from the Element from to to More...
 
void add_slave (Bond &b)
 add a slave bond (useful for periodic simulation) More...
 
virtual std::ostream & export_to_povray (std::ostream &out) const
 exporting to povray format More...
 
- Public Member Functions inherited from GranOO3::Physic::BodyInteraction
 BodyInteraction (Body &b1, Body &b2, bool contact_interaction)
 constructor More...
 
virtual ~BodyInteraction ()
 destructor More...
 
void apply_torque ()
 apply the current torque value More...
 
void clear_torque ()
 set the torques to zero More...
 
Geom::Vectortorque_on (const Body &b)
 get the torque value applied on b More...
 
const Geom::Vectortorque_on (const Body &b) const
 same as BodyInteraction::torque_on More...
 
void add_label_torque (const std::string &label, const Geom::Vector &t1, const Geom::Vector &t2)
 add a label for the torque More...
 
void swap (Body &from, Body &to)
 swap one body More...
 
- Public Member Functions inherited from GranOO3::Physic::NodeInteraction
 NodeInteraction (Node &de1, Node &de2, bool build_node_interaction)
 
virtual ~NodeInteraction ()
 
bool contact_interaction () const
 
void apply_force ()
 
void clear_force ()
 
Geom::Vectorforce_on (const Node &)
 
const Geom::Vectorforce_on (const Node &) const
 
void swap (Node &from, Node &to)
 
virtual std::ostream & write_ascii (std::ostream &out) const
 
void add_label_force (const std::string &label, const Geom::Vector &f1, const Geom::Vector &f2)
 
- Public Member Functions inherited from GranOO3::Core::Base
virtual ~Base ()
 
size_t numID () const
 
size_t uID () const
 
void set_numID (size_t val)
 
void clear_numID ()
 
Physic::Materialget_mat () const
 
void set_mat (Physic::Material *)
 
Baseitem ()
 
const Baseitem () const
 
bool is_same (const Base &) const
 
template<class T >
T & cast_to ()
 
template<class T >
const T & cast_to () const
 
template<class T >
bool is () const
 
virtual bool is (size_t) const
 
virtual const std::string & get_ID () const
 
virtual Baseclone ()
 
template<class T >
T & clone_to ()
 
virtual std::istream & read_ascii (std::istream &in)
 
Signal< Base & > & deleted_signal ()
 
- Public Member Functions inherited from GranOO3::Core::Null
 Null ()
 
virtual ~Null ()
 
- Public Member Functions inherited from GranOO3::Core::Drawable
 Drawable ()
 
virtual ~Drawable ()
 
virtual void draw_edge ()
 
virtual void init_default_color ()
 
virtual std::ostream & get_info (std::ostream &os) const
 
virtual const Colordefault_color () const
 
Colorget_color ()
 
const Colorget_color () const
 
virtual void set_color (const Color &)
 
virtual void set_alpha (float alpha)
 
void apply_color () const
 
void apply_edge_color () const
 
void apply_default_color ()
 
void apply_selected_color ()
 
void set_line_width (float)
 
float get_line_width () const
 
float & get_line_width ()
 
void apply_line_width () const
 
bool is_visible () const
 
void set_visible (bool)
 
void paint ()
 
void paint_edge ()
 
unsigned int get_item_glkey () const
 
- Public Member Functions inherited from GranOO3::Core::Register< Base >
 Register ()
 
virtual ~Register ()
 
void erase_from_all_setof ()
 
bool belong_to_setof (const std::string &setOfId) const
 
bool belong_to_setof (const SetOf< Base > &set) const
 
std::list< SetOf< Base > * > & get_setof_list ()
 
unsigned long long int get_numeric_ID () const
 
- Public Member Functions inherited from GranOO3::Core::Register< NodeInteraction >
 Register ()
 
virtual ~Register ()
 
void erase_from_all_setof ()
 
bool belong_to_setof (const std::string &setOfId) const
 
bool belong_to_setof (const SetOf< NodeInteraction > &set) const
 
std::list< SetOf< NodeInteraction > * > & get_setof_list ()
 
unsigned long long int get_numeric_ID () const
 
- Public Member Functions inherited from GranOO3::Core::Register< BodyInteraction >
 Register ()
 
virtual ~Register ()
 
void erase_from_all_setof ()
 
bool belong_to_setof (const std::string &setOfId) const
 
bool belong_to_setof (const SetOf< BodyInteraction > &set) const
 
std::list< SetOf< BodyInteraction > * > & get_setof_list ()
 
unsigned long long int get_numeric_ID () const
 
- Public Member Functions inherited from GranOO3::Core::PropClass< Bond >
 PropClass ()
 
virtual ~PropClass ()
 
T & new_object ()
 
T & get ()
 
const T & get () const
 
bool prop_exist () const
 
void add_prop (Core::Prop< Bond > *)
 
void remove_prop (Core::Prop< Bond > *)
 
std::string info () const
 
- Public Member Functions inherited from GranOO3::Core::Register< Bond >
 Register ()
 
virtual ~Register ()
 
void erase_from_all_setof ()
 
bool belong_to_setof (const std::string &setOfId) const
 
bool belong_to_setof (const SetOf< Bond > &set) const
 
std::list< SetOf< Bond > * > & get_setof_list ()
 
unsigned long long int get_numeric_ID () const
 
- Public Member Functions inherited from GranOO3::Core::Register< XBeam >
 Register ()
 
virtual ~Register ()
 
void erase_from_all_setof ()
 
bool belong_to_setof (const std::string &setOfId) const
 
bool belong_to_setof (const SetOf< XBeam > &set) const
 
std::list< SetOf< XBeam > * > & get_setof_list ()
 
unsigned long long int get_numeric_ID () const
 

Protected Attributes

double _young_modulus
 The Young's modulus value of the beam in [Pa]. More...
 
double _damping_factor
 The damping factor of the beam expressed as a ratio of the critical damping. More...
 
double _radius
 The radius of the beam in [m]. More...
 
double _poisson_ratio
 The Poisson's ratio of the beam in [-]. More...
 
double _max_stress
 The maximal tensile stress of the beam in [Pa]. More...
 
double _current_stress
 The current maximal stress (including normal, bending and torsion stresses) in [Pa]. More...
 
double _current_tensile_stress
 The current tensile stress into the beam. More...
 
double _current_bending_stress
 The current maximal bending stress into the beam. More...
 
double _current_torsion_stress
 The current maximal torsion stress into the beam. More...
 
double _current_normal_stress
 The current maximal normal stress into the beam (sum of tensile and bending stresses) More...
 
double _tensile_energy
 
double _shear_energy
 
double _torsion_energy
 
double _bending_energy
 
double _Ig
 For storing some intermediate values to speed-up the computation (quadratic inertia) More...
 
double _Io
 For storing some intermediate values to speed-up the computation (quadratic polar inertia) More...
 
double _G
 For storing some intermediate values to speed-up the computation (Coulomb's modulus) More...
 
double _S
 For storing some intermediate values to speed-up the computation (normal surface) More...
 
double _K
 For storing some intermediate values to speed-up the computation (normal stiffness) More...
 
- Protected Attributes inherited from GranOO3::DEM::Bond
Element_element1
 the first of the two bonded Element More...
 
Element_element2
 the second of the two bonded Element More...
 
double _linear_damping_factor
 the linear damping factor of the bond More...
 
double _angular_damping_factor
 the angular damping factor of the bond More...
 
double _relaxed_length
 the relaxed length of the bond. More...
 
double _current_length
 the current length of the bond that correspond to the norm of the vector between the center of the bonded elements More...
 
double _max_relative_elongation
 a threshold corresponding to a maximal admissible relative elongation (the longitudinal strain) More...
 
double _surface
 this attribute stores the corresponding normal surface of the bond More...
 
Core::SetOfBase< Bond_slave
 a container that list the slave bond More...
 
std::vector< unsigned int > _dofs
 a vector that stores the dof of the bond More...
 
- Protected Attributes inherited from GranOO3::Physic::BodyInteraction
Body_body1
 a pointer the first body More...
 
Body_body2
 a pointer the second body More...
 
Geom::Vector _torque_on1
 the current value of the applied torque on the body 1 More...
 
Geom::Vector _torque_on2
 the current value of the applied torque on the body 2 More...
 
- Protected Attributes inherited from GranOO3::Physic::NodeInteraction
const bool _contact_interaction
 
Node_node1
 
Node_node2
 
Geom::Vector _force_on1
 
Geom::Vector _force_on2
 
- Protected Attributes inherited from GranOO3::Core::PropClass< Bond >
std::vector< Core::Prop< Bond > * > _prop
 

Private Member Functions

void update ()
 Updates the beam's current configuration. This is required in case of large transformations, because the beam's current frame may have bee rotated since the previous configuration. More...
 
void build_transformation_matrix ()
 Builds the rotation matrix used to rotate the beam's local stiffness matrix (expressed along the x-axis) to its initial frame (the one usually defined by domain construction). More...
 
void build_local_stiffness_matrix ()
 Builds the beam's local stiffness matrix (expressed along the x-axis). More...
 
void on_matrix_from_xdir (const Geom::Vector &xDir, const Geom::Vector &yDir, Eigen::Matrix3d &rotMatrix) const
 Builds an orthonormal matrix given two vectors, x- and y-oriented. More...
 
void rot_matrix_to_quaternion (const Eigen::Matrix3d &rotMatrix, Geom::Quaternion &quaternion) const
 Builds a quaternion from a rotation matrix. More...
 
void quaternion_to_rot_matrix (const Geom::Quaternion &quaternion, Eigen::Matrix3d &rotMatrix) const
 Builds a rotation matrix from a quaternion. More...
 
 XBeam ()=delete
 
 XBeam (const XBeam &frame)=delete
 
XBeamoperator= (const XBeam &)=delete
 
template<class Archive >
void save (Archive &, const unsigned int) const
 complete serializing of the bond in the *.gdd format More...
 
template<class Archive >
void load (Archive &, const unsigned int)
 Complete serializing of the bond in the *.gdd format. More...
 
 BOOST_SERIALIZATION_SPLIT_MEMBER ()
 

Private Attributes

bool _nl_geom
 Defines whether non linear geometric effects should be accounted for, i.e. finite rotations. More...
 
Geom::Vector _xAxis
 The current x-axis of thje XBeam, define as the direction between both ends of the beam. More...
 
double _ldFactor
 The linear damping factor, precomputed at the initialization of the XBeam. More...
 
double _adFactor
 The angular damping factor, precomputed at the initialization of the XBeam. More...
 
Geom::Quaternion _qbeam_ref_rot
 The beam's reference (or initial) rotation quaternion. More...
 
Geom::Quaternion _qbeam_abs_rot
 The beam's absolute (or current) rotation quaternion. More...
 
Geom::Quaternion _qel1_ref_rot
 The node 1 reference (or initial) rotation quaternion. More...
 
Geom::Quaternion _qel2_ref_rot
 The node 2 reference (or initial) rotation quaternion. More...
 
Eigen::MatrixXd _local_kmat
 The local stiffness matrix, defined in the x-direction. More...
 
Eigen::MatrixXd _tmat
 The transformation (rotation) matrix, used to express the beam's stiffness matrix in its current frame. More...
 
Eigen::MatrixXd _global_kmat
 The beam's stiffness matrix in its current frame. More...
 
Eigen::VectorXd _reaction
 The reaction (or internal) forces vector, updated at each time increment, at the XBeam level. More...
 
Eigen::Matrix3d _mbeam_abs_rot
 The beam's absolute (or current) rotation matrix. More...
 
Eigen::Matrix3d _mbeam_ref_rot
 The beam's reference (or initial) rotation matrix. More...
 

Friends

class boost::serialization::access
 

Additional Inherited Members

- Public Types inherited from GranOO3::DEM::Bond
enum  FailureMode { NONE =0 , TENSION =1 , COMPRESSION =2 , SHEAR =3 }
 
enum  DrawMode { SOLID , LINE , NO }
 
- Static Public Member Functions inherited from GranOO3::DEM::Bond
static double get_disabled_cumulative_surface ()
 get the total surface disabled surface More...
 
- Static Public Member Functions inherited from GranOO3::Core::Base
static Baseget_by_numID (size_t)
 
static void clear_all_numID ()
 
static unsigned int get_sub_class_number ()
 
- Static Public Member Functions inherited from GranOO3::Core::Drawable
static Drawableget_drawable_item_by_glkey (int)
 
static void set_draw_precision (unsigned int p)
 
static void increase_draw_precision ()
 
static void decrease_draw_precision ()
 
static unsigned int get_draw_precision ()
 
- Static Public Member Functions inherited from GranOO3::Physic::CriticalTimeStep
static const std::set< CriticalTimeStep * > & get_all ()
 
- Static Public Attributes inherited from GranOO3::DEM::Bond
static double epsilon
 a very small value used as a threshold More...
 
static double draw_solid_factor
 
static DrawMode draw_mode = LINE
 
static Core::Signal< Bond & > disable_signal = Core::Signal<Bond&>()
 get the signal when a bond is disabled More...
 
- Static Public Attributes inherited from GranOO3::Core::Null
static Null null = Null()
 
- Protected Member Functions inherited from GranOO3::DEM::Bond
template<class Archive >
void save (Archive &, const unsigned int) const
 
template<class Archive >
void load (Archive &, const unsigned int)
 
 BOOST_SERIALIZATION_SPLIT_MEMBER ()
 
- Protected Member Functions inherited from GranOO3::Core::Base
 Base ()
 
- Protected Member Functions inherited from GranOO3::Physic::CriticalTimeStep
 CriticalTimeStep ()
 constructor More...
 
virtual ~CriticalTimeStep ()
 destructor More...
 
- Static Protected Member Functions inherited from GranOO3::Core::Base
static unsigned int affect_class_rank_ID ()
 

Detailed Description

A bond XBeam that works in tension, compression, bending and torsion very similar to the Beam bond.

The difference between Beam and XBeam is that XBeam uses a FEM formulation. So, XBeam can be used in FEM problem and it could advantageous speed-up a lattice calculation in comparison with the standard Beam.

Constructor & Destructor Documentation

◆ XBeam() [1/5]

GranOO3::DEM::XBeam::XBeam ( Element el1,
Element el2 
)

constructor

Parameters
[in]el1: the first bonded element
[in]el2: the second bonded element

Builds a new XBeam between the Element el1 and el2. Mechanical attributes are set to zero.

◆ XBeam() [2/5]

GranOO3::DEM::XBeam::XBeam ( Element el1,
Element el2,
double  young_mod,
double  radius,
double  poisson_ratio,
double  max_stress = 0. 
)

constructor

Parameters
[in]el1: the first bonded element
[in]el2: the second bonded element
[in]young_mod: the Young's modulus of the beam in [Pa]
[in]radius: the radius of the beam in [m]
[in]poisson_ratio: the Poisson's ratio value of the beam
[in]max_stress: the maximal stress threshold of the beam in [Pa] (note that zero means infinite strength)

Build sa new XBeam between the Element el1 and el2 with the required mechanical attributes. In case of continuum modeling, note that the max_stress is not the best way for managing fracture. You may prefer to use the stress computation (see the ApplyFracture PlugIn) which gives more accurate results.

◆ XBeam() [3/5]

GranOO3::DEM::XBeam::XBeam ( Element el1,
Element el2,
double  young_mod,
double  damping,
double  radius,
double  poisson_ratio,
double  max_stress,
bool  nl_geom = false 
)

constructor

Parameters
[in]el1: the first bonded element
[in]el2: the second bonded element
[in]young_mod: the Young's modulus of the beam in [Pa]
[in]damping: the damping factor (in the [0,1] range) which is expressed as the ratio of the critical damping
[in]radius: the radius of the beam in [m]
[in]poisson_ratio: the Poisson's ratio value of the beam
[in]max_stress: the maximal stress threshold of the beam in [Pa] (note that zero means infinite strength)
[in]nl_geom: a flag to speed-up the computation (for simulation with small strain, displacement and rotation)

Builds a new XBeam between the Element el1 and el2 with the required mechanical attributes. In case of continuum modeling, note that the max_stress is not the best way for managing fracture. You may prefer to use the stress computation (see the ApplyFracture PlugIn) which gives more accurate results. You can set the nl_geom flag to true if you know that your simulation will handle only small deformation, displacement and rotation. In such a case, passing true to the nl_geom flag can speed-up the computation.

◆ ~XBeam()

GranOO3::DEM::XBeam::~XBeam ( )
virtual

◆ XBeam() [4/5]

GranOO3::DEM::XBeam::XBeam ( )
privatedelete

◆ XBeam() [5/5]

GranOO3::DEM::XBeam::XBeam ( const XBeam frame)
privatedelete

Member Function Documentation

◆ BOOST_SERIALIZATION_SPLIT_MEMBER()

GranOO3::DEM::XBeam::BOOST_SERIALIZATION_SPLIT_MEMBER ( )
private

◆ build_local_stiffness_matrix()

void GranOO3::DEM::XBeam::build_local_stiffness_matrix ( )
private

Builds the beam's local stiffness matrix (expressed along the x-axis).

◆ build_transformation_matrix()

void GranOO3::DEM::XBeam::build_transformation_matrix ( )
private

Builds the rotation matrix used to rotate the beam's local stiffness matrix (expressed along the x-axis) to its initial frame (the one usually defined by domain construction).

◆ compute_critical_time_step()

double GranOO3::DEM::XBeam::compute_critical_time_step ( ) const
virtual

compute the critical time of the bond

Returns
the maximal time step

Note that this method computes on-the-fly the critical time step.

Implements GranOO3::Physic::CriticalTimeStep.

◆ compute_load()

void GranOO3::DEM::XBeam::compute_load ( )
virtual

Compute the reaction force and torque of the beam.

Reimplemented from GranOO3::DEM::Bond.

◆ draw()

void GranOO3::DEM::XBeam::draw ( )
inlinevirtual

OpenGL draw of the bond

This method is used for drawing the bond. It is used by the granoo-viewer vizualization tool.

Reimplemented from GranOO3::DEM::Bond.

◆ info()

std::string GranOO3::DEM::XBeam::info ( ) const
virtual

Displays some useful info in the terminal.

Reimplemented from GranOO3::DEM::Bond.

◆ load()

template<class Archive >
void GranOO3::DEM::XBeam::load ( Archive &  ar,
const unsigned int  v 
)
private

Complete serializing of the bond in the *.gdd format.

See also
the Domain classes to get additional info about I/O.

◆ on_matrix_from_xdir()

void GranOO3::DEM::XBeam::on_matrix_from_xdir ( const Geom::Vector xDir,
const Geom::Vector yDir,
Eigen::Matrix3d &  rotMatrix 
) const
private

Builds an orthonormal matrix given two vectors, x- and y-oriented.

Parameters
[in]xDir: A const reference to a Geom::Vector object, defining the x-direction.
[in]yDir: A const reference to a Geom::Vector object, defining the y-direction.
[out]rotMatrix: The orthonormal matrix, built from the x- and y-oriented vectors.

◆ operator=()

XBeam& GranOO3::DEM::XBeam::operator= ( const XBeam )
privatedelete

◆ quaternion_to_rot_matrix()

void GranOO3::DEM::XBeam::quaternion_to_rot_matrix ( const Geom::Quaternion quaternion,
Eigen::Matrix3d &  rotMatrix 
) const
private

Builds a rotation matrix from a quaternion.

Parameters
[in]quaternion: A const reference to a Geom::Quaterion.
[out]rotMatrix: The resulting rotation matrix.

◆ rot_matrix_to_quaternion()

void GranOO3::DEM::XBeam::rot_matrix_to_quaternion ( const Eigen::Matrix3d &  rotMatrix,
Geom::Quaternion quaternion 
) const
private

Builds a quaternion from a rotation matrix.

Parameters
[in]rotMatrix: A const reference to an Eigen::Matrix3d rotation matrix.
[out]quaternion: The resulting quaternion.

◆ save()

template<class Archive >
void GranOO3::DEM::XBeam::save ( Archive &  ar,
const unsigned int  v 
) const
private

complete serializing of the bond in the *.gdd format

See also
the Domain classes to get additional info about I/O.

◆ stiffness_matrix()

void GranOO3::DEM::XBeam::stiffness_matrix ( Eigen::MatrixXd &  kMatrix)
virtual

Assigns the global stiffness matrix of the XBeam kMatrix.

Parameters
[in,out]kMatrix: the stiffness to be updated.

Reimplemented from GranOO3::DEM::Bond.

◆ update()

void GranOO3::DEM::XBeam::update ( )
private

Updates the beam's current configuration. This is required in case of large transformations, because the beam's current frame may have bee rotated since the previous configuration.

◆ update_dof()

void GranOO3::DEM::XBeam::update_dof ( )
virtual

Builds the list of global Degrees Of Freedom (DOFs), required by the solver to assemble the global (system-level) stiffness matrix and the displacement and forces vectors.

Reimplemented from GranOO3::DEM::Bond.

Friends And Related Function Documentation

◆ boost::serialization::access

friend class boost::serialization::access
friend

Member Data Documentation

◆ _adFactor

double GranOO3::DEM::XBeam::_adFactor
private

The angular damping factor, precomputed at the initialization of the XBeam.

◆ _bending_energy

double GranOO3::DEM::XBeam::_bending_energy
protected

◆ _current_bending_stress

double GranOO3::DEM::XBeam::_current_bending_stress
protected

The current maximal bending stress into the beam.

◆ _current_normal_stress

double GranOO3::DEM::XBeam::_current_normal_stress
protected

The current maximal normal stress into the beam (sum of tensile and bending stresses)

◆ _current_stress

double GranOO3::DEM::XBeam::_current_stress
protected

The current maximal stress (including normal, bending and torsion stresses) in [Pa].

Note that the _max_stress threshold uses this value as computational reference for determining the maximal strength.

◆ _current_tensile_stress

double GranOO3::DEM::XBeam::_current_tensile_stress
protected

The current tensile stress into the beam.

◆ _current_torsion_stress

double GranOO3::DEM::XBeam::_current_torsion_stress
protected

The current maximal torsion stress into the beam.

◆ _damping_factor

double GranOO3::DEM::XBeam::_damping_factor
protected

The damping factor of the beam expressed as a ratio of the critical damping.

◆ _G

double GranOO3::DEM::XBeam::_G
protected

For storing some intermediate values to speed-up the computation (Coulomb's modulus)

◆ _global_kmat

Eigen::MatrixXd GranOO3::DEM::XBeam::_global_kmat
private

The beam's stiffness matrix in its current frame.

◆ _Ig

double GranOO3::DEM::XBeam::_Ig
protected

For storing some intermediate values to speed-up the computation (quadratic inertia)

◆ _Io

double GranOO3::DEM::XBeam::_Io
protected

For storing some intermediate values to speed-up the computation (quadratic polar inertia)

◆ _K

double GranOO3::DEM::XBeam::_K
protected

For storing some intermediate values to speed-up the computation (normal stiffness)

◆ _ldFactor

double GranOO3::DEM::XBeam::_ldFactor
private

The linear damping factor, precomputed at the initialization of the XBeam.

◆ _local_kmat

Eigen::MatrixXd GranOO3::DEM::XBeam::_local_kmat
private

The local stiffness matrix, defined in the x-direction.

◆ _max_stress

double GranOO3::DEM::XBeam::_max_stress
protected

The maximal tensile stress of the beam in [Pa].

◆ _mbeam_abs_rot

Eigen::Matrix3d GranOO3::DEM::XBeam::_mbeam_abs_rot
private

The beam's absolute (or current) rotation matrix.

◆ _mbeam_ref_rot

Eigen::Matrix3d GranOO3::DEM::XBeam::_mbeam_ref_rot
private

The beam's reference (or initial) rotation matrix.

◆ _nl_geom

bool GranOO3::DEM::XBeam::_nl_geom
private

Defines whether non linear geometric effects should be accounted for, i.e. finite rotations.

◆ _poisson_ratio

double GranOO3::DEM::XBeam::_poisson_ratio
protected

The Poisson's ratio of the beam in [-].

◆ _qbeam_abs_rot

Geom::Quaternion GranOO3::DEM::XBeam::_qbeam_abs_rot
private

The beam's absolute (or current) rotation quaternion.

◆ _qbeam_ref_rot

Geom::Quaternion GranOO3::DEM::XBeam::_qbeam_ref_rot
private

The beam's reference (or initial) rotation quaternion.

◆ _qel1_ref_rot

Geom::Quaternion GranOO3::DEM::XBeam::_qel1_ref_rot
private

The node 1 reference (or initial) rotation quaternion.

◆ _qel2_ref_rot

Geom::Quaternion GranOO3::DEM::XBeam::_qel2_ref_rot
private

The node 2 reference (or initial) rotation quaternion.

◆ _radius

double GranOO3::DEM::XBeam::_radius
protected

The radius of the beam in [m].

◆ _reaction

Eigen::VectorXd GranOO3::DEM::XBeam::_reaction
private

The reaction (or internal) forces vector, updated at each time increment, at the XBeam level.

◆ _S

double GranOO3::DEM::XBeam::_S
protected

For storing some intermediate values to speed-up the computation (normal surface)

◆ _shear_energy

double GranOO3::DEM::XBeam::_shear_energy
protected

◆ _tensile_energy

double GranOO3::DEM::XBeam::_tensile_energy
protected

◆ _tmat

Eigen::MatrixXd GranOO3::DEM::XBeam::_tmat
private

The transformation (rotation) matrix, used to express the beam's stiffness matrix in its current frame.

◆ _torsion_energy

double GranOO3::DEM::XBeam::_torsion_energy
protected

◆ _xAxis

Geom::Vector GranOO3::DEM::XBeam::_xAxis
private

The current x-axis of thje XBeam, define as the direction between both ends of the beam.

◆ _young_modulus

double GranOO3::DEM::XBeam::_young_modulus
protected

The Young's modulus value of the beam in [Pa].


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