32 #include "GranOO3/3rdParty/Eigen/Dense"
50 XBeam(
Element& el1,
Element& el2,
double young_mod,
double radius,
double poisson_ratio,
double max_stress = 0.);
51 XBeam(
Element& el1,
Element& el2,
double young_mod,
double damping,
double radius,
double poisson_ratio,
double max_stress,
bool nl_geom =
false);
59 virtual std::string
info()
const;
65 GRANOO_ACCESS(radius ,
double,
_radius);
69 GRANOO_ACCESS_GET_COPY(section,
double, M_PI*
pow(
_radius,2));
70 GRANOO_ACCESS_GET_COPY(inertia,
double, M_PI*
pow(
_radius*2.,4)/64.);
103 template<
class Archive>
void save(Archive &,
const unsigned int )
const;
104 template<
class Archive>
void load(Archive &,
const unsigned int);
150 template<
class Archive>
void
152 ar << boost::serialization::base_object<Bond>(*
this);
169 template<
class Archive>
void
171 ar >> boost::serialization::base_object<Bond>(*
this);
198 #include <boost/serialization/version.hpp>
204 namespace serialization
207 template<
class Archive>
void
215 template<
class Archive>
void
228 GRANOO_CLASS_DECLARE_TPL(DEM::XBeam);
BOOST_CLASS_VERSION(GranOO3::DEM::XBeam, 0) namespace boost
Definition: XBeam.hpp:199
Definition: SetOf.hpp:346
the base class for all bonds between discrete elements.
Definition: Bond.hpp:49
virtual void draw()
OpenGL draw of the bond
Definition: Bond.cpp:428
double _relaxed_length
the relaxed length of the bond.
Definition: Bond.hpp:209
a base class that represents an element
Definition: Element.hpp:55
A bond XBeam that works in tension, compression, bending and torsion very similar to the Beam bond.
Definition: XBeam.hpp:44
virtual void compute_load()
Compute the reaction force and torque of the beam.
Definition: XBeam.cpp:234
double _G
For storing some intermediate values to speed-up the computation (Coulomb's modulus)
Definition: XBeam.hpp:130
Geom::Quaternion _qbeam_ref_rot
The beam's reference (or initial) rotation quaternion.
Definition: XBeam.hpp:138
Eigen::MatrixXd _tmat
The transformation (rotation) matrix, used to express the beam's stiffness matrix in its current fram...
Definition: XBeam.hpp:143
double _poisson_ratio
The Poisson's ratio of the beam in [-].
Definition: XBeam.hpp:112
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.
Definition: XBeam.cpp:410
Eigen::VectorXd _reaction
The reaction (or internal) forces vector, updated at each time increment, at the XBeam level.
Definition: XBeam.hpp:145
XBeam & operator=(const XBeam &)=delete
double _young_modulus
The Young's modulus value of the beam in [Pa].
Definition: XBeam.hpp:109
void save(Archive &, const unsigned int) const
complete serializing of the bond in the *.gdd format
Definition: XBeam.hpp:151
double _Ig
For storing some intermediate values to speed-up the computation (quadratic inertia)
Definition: XBeam.hpp:130
XBeam(const XBeam &frame)=delete
double _max_stress
The maximal tensile stress of the beam in [Pa].
Definition: XBeam.hpp:113
Eigen::MatrixXd _local_kmat
The local stiffness matrix, defined in the x-direction.
Definition: XBeam.hpp:142
Geom::Quaternion _qbeam_abs_rot
The beam's absolute (or current) rotation quaternion.
Definition: XBeam.hpp:139
void quaternion_to_rot_matrix(const Geom::Quaternion &quaternion, Eigen::Matrix3d &rotMatrix) const
Builds a rotation matrix from a quaternion.
Definition: XBeam.cpp:386
void build_local_stiffness_matrix()
Builds the beam's local stiffness matrix (expressed along the x-axis).
Definition: XBeam.cpp:137
double _tensile_energy
Definition: XBeam.hpp:122
double _damping_factor
The damping factor of the beam expressed as a ratio of the critical damping.
Definition: XBeam.hpp:110
double _current_tensile_stress
The current tensile stress into the beam.
Definition: XBeam.hpp:116
virtual void stiffness_matrix(Eigen::MatrixXd &kMatrix)
Assigns the global stiffness matrix of the XBeam kMatrix.
Definition: XBeam.cpp:370
bool _nl_geom
Defines whether non linear geometric effects should be accounted for, i.e. finite rotations.
Definition: XBeam.hpp:133
double _current_bending_stress
The current maximal bending stress into the beam.
Definition: XBeam.hpp:117
virtual std::string info() const
Displays some useful info in the terminal.
Definition: XBeam.cpp:556
Geom::Quaternion _qel1_ref_rot
The node 1 reference (or initial) rotation quaternion.
Definition: XBeam.hpp:140
virtual void draw()
OpenGL draw of the bond
Definition: XBeam.hpp:189
double _adFactor
The angular damping factor, precomputed at the initialization of the XBeam.
Definition: XBeam.hpp:136
Eigen::MatrixXd _global_kmat
The beam's stiffness matrix in its current frame.
Definition: XBeam.hpp:144
virtual void update_dof()
Builds the list of global Degrees Of Freedom (DOFs), required by the solver to assemble the global (s...
Definition: XBeam.cpp:375
double _K
For storing some intermediate values to speed-up the computation (normal stiffness)
Definition: XBeam.hpp:130
Eigen::Matrix3d _mbeam_abs_rot
The beam's absolute (or current) rotation matrix.
Definition: XBeam.hpp:146
double _bending_energy
Definition: XBeam.hpp:125
void update()
Updates the beam's current configuration. This is required in case of large transformations,...
Definition: XBeam.cpp:220
void load(Archive &, const unsigned int)
Complete serializing of the bond in the *.gdd format.
Definition: XBeam.hpp:170
double _Io
For storing some intermediate values to speed-up the computation (quadratic polar inertia)
Definition: XBeam.hpp:130
double _current_normal_stress
The current maximal normal stress into the beam (sum of tensile and bending stresses)
Definition: XBeam.hpp:119
double _shear_energy
Definition: XBeam.hpp:123
double _radius
The radius of the beam in [m].
Definition: XBeam.hpp:111
virtual double compute_critical_time_step() const
compute the critical time of the bond
Definition: XBeam.cpp:529
double _current_stress
The current maximal stress (including normal, bending and torsion stresses) in [Pa].
Definition: XBeam.hpp:115
double _current_torsion_stress
The current maximal torsion stress into the beam.
Definition: XBeam.hpp:118
Eigen::Matrix3d _mbeam_ref_rot
The beam's reference (or initial) rotation matrix.
Definition: XBeam.hpp:147
Geom::Vector _xAxis
The current x-axis of thje XBeam, define as the direction between both ends of the beam.
Definition: XBeam.hpp:134
void rot_matrix_to_quaternion(const Eigen::Matrix3d &rotMatrix, Geom::Quaternion &quaternion) const
Builds a quaternion from a rotation matrix.
Definition: XBeam.cpp:461
friend class boost::serialization::access
Definition: XBeam.hpp:102
double _S
For storing some intermediate values to speed-up the computation (normal surface)
Definition: XBeam.hpp:130
void build_transformation_matrix()
Builds the rotation matrix used to rotate the beam's local stiffness matrix (expressed along the x-ax...
Definition: XBeam.cpp:193
BOOST_SERIALIZATION_SPLIT_MEMBER()
virtual ~XBeam()
Definition: XBeam.cpp:133
Geom::Quaternion _qel2_ref_rot
The node 2 reference (or initial) rotation quaternion.
Definition: XBeam.hpp:141
double _ldFactor
The linear damping factor, precomputed at the initialization of the XBeam.
Definition: XBeam.hpp:135
double _torsion_energy
Definition: XBeam.hpp:124
Definition: Quaternion.hpp:54
Definition: Vector.hpp:75
pure virtual class for modeling classes able to compute a critical time step
Definition: CriticalTimeStep.hpp:50
Definition: Common.hpp:198
void save_construct_data(Archive &ar, const GranOO3::Core::Pair< type > *t, const unsigned int)
Definition: Pair.hpp:207
void load_construct_data(Archive &ar, GranOO3::Core::Pair< type > *t, const unsigned int)
Definition: Pair.hpp:217
T pow(const T v0, const T v1)
Definition: Exprtk.hpp:1491
x y t t *t x y t t t x y t t t x *y t *t t x *y t *t t x y t t t x y t t t t(t+t)") define_sfop3(16