GranOO  3.0
A robust and versatile workbench to build 3D dynamic simulations based on the Discrete Element Method
XBeam.hpp
Go to the documentation of this file.
1 // This file is part of GranOO, a workbench for DEM simulation.
2 //
3 // Author(s) : - Damien Andre IRCER/UNILIM, Limoges France
4 // <damien.andre@unilim.fr>
5 // - Jean-luc Charles Arts et Metiers ParisTech, CNRS, I2M, Bordeaux France
6 // <jean-luc.charles@ensam.eu>
7 // - Jeremie Girardot Arts et Metiers ParisTech, CNRS, I2M, Bordeaux France
8 // <jeremie.girardot@ensam.eu>
9 // - Cedric Hubert LAMIH/UPHF, Valenciennes France
10 // <cedric.hubert@uphf.fr>
11 // - Ivan Iordanoff Arts et Metiers ParisTech, CNRS, I2M, Bordeaux France
12 // <ivan.iordanoff@ensam.eu>
13 //
14 // Copyright (C) 2008-2019 D. Andre, JL. Charles, J. Girardot, C. Hubert, I. Iordanoff
15 //
16 // This program is free software: you can redistribute it and/or modify
17 // it under the terms of the GNU General Public License as published by
18 // the Free Software Foundation, either version 3 of the License, or
19 // (at your option) any later version.
20 //
21 // This program is distributed in the hope that it will be useful,
22 // but WITHOUT ANY WARRANTY; without even the implied warranty of
23 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
24 // GNU General Public License for more details.
25 //
26 // You should have received a copy of the GNU General Public License
27 // along with this program. If not, see <http://www.gnu.org/licenses/>.
28 
29 #ifndef XBeam_hpp
30 #define XBeam_hpp
31 
32 #include "GranOO3/3rdParty/Eigen/Dense"
33 
34 #include "GranOO3/DEM/Bond.hpp"
35 #include "GranOO3/DEM/Element.hpp"
36 
37 namespace GranOO3
38 {
39  namespace DEM
40  {
41 
42  class XBeam : public Bond,
43  public Core::Register<XBeam>,
45 
46  GRANOO_CLASS(DEM::XBeam, Bond);
47 
48  public:
49  XBeam(Element& el1, Element& el2);
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);
52  virtual ~XBeam();
53 
54  public:
55  virtual void compute_load();
56  virtual double compute_critical_time_step() const;
57  virtual void stiffness_matrix(Eigen::MatrixXd& kMatrix);
58  virtual void update_dof();
59  virtual std::string info() const;
60  virtual void draw();
61 
62  //ACCESSORS
63  GRANOO_ACCESS(young_modulus , double, _young_modulus);
64  GRANOO_ACCESS(damping_factor, double, _damping_factor);
65  GRANOO_ACCESS(radius , double, _radius);
66  GRANOO_ACCESS(poisson_ratio , double, _poisson_ratio);
67  GRANOO_ACCESS(max_stress , double, _max_stress);
68 
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.);
71  GRANOO_ACCESS_GET_COPY(traction_stiffness, double, _young_modulus*get_section()/_relaxed_length);
72  GRANOO_ACCESS_GET_COPY(flexion_stiffness , double, 2.*_young_modulus*get_inertia()/_relaxed_length);
73 
74  GRANOO_ACCESS_GET(tensile_energy , double, _tensile_energy );
75  GRANOO_ACCESS_GET(shear_energy , double, _shear_energy );
76  GRANOO_ACCESS_GET(torsion_energy , double, _torsion_energy );
77  GRANOO_ACCESS_GET(bending_energy , double, _bending_energy );
78  GRANOO_ACCESS_GET_COPY(total_energy, double, _tensile_energy+_torsion_energy+_bending_energy+_shear_energy);
79 
80  GRANOO_ACCESS_GET(current_stress , double, _current_stress );
81  GRANOO_ACCESS_GET(current_tensile_stress, double, _current_tensile_stress );
82  GRANOO_ACCESS_GET(current_bending_stress, double, _current_bending_stress );
83  GRANOO_ACCESS_GET(current_torsion_stress, double, _current_torsion_stress );
84  GRANOO_ACCESS_GET(current_normal_stress , double, _current_normal_stress );
85 
86  private:
87  void update();
90  void on_matrix_from_xdir(const Geom::Vector& xDir, const Geom::Vector& yDir, Eigen::Matrix3d& rotMatrix) const;
91  void rot_matrix_to_quaternion(const Eigen::Matrix3d& rotMatrix, Geom::Quaternion& quaternion) const;
92  void quaternion_to_rot_matrix(const Geom::Quaternion& quaternion, Eigen::Matrix3d& rotMatrix) const;
93 
94 
95  private:
96  XBeam() = delete;
97  XBeam(const XBeam &frame) = delete;
98  XBeam & operator = (const XBeam &) = delete;
99 
100  private:
101  //SERIALIZATION
103  template<class Archive> void save(Archive &, const unsigned int ) const;
104  template<class Archive> void load(Archive &, const unsigned int);
106 
107  protected:
108  //Parameter
111  double _radius;
113  double _max_stress;
114 
120 
121  //Deformation energies
126 
127  // We store here coefficient values to speedup the computation
128  // Note that this value must be updated if the setting of the beam
129  // are changed (radius, relaxed length and so on...)
130  double _Ig, _Io, _G, _S, _K;
131 
132  private:
133  bool _nl_geom;
135  double _ldFactor;
136  double _adFactor;
137 
142  Eigen::MatrixXd _local_kmat;
143  Eigen::MatrixXd _tmat;
144  Eigen::MatrixXd _global_kmat;
145  Eigen::VectorXd _reaction;
146  Eigen::Matrix3d _mbeam_abs_rot;
147  Eigen::Matrix3d _mbeam_ref_rot;
148  };
149 
150  template<class Archive> void
151  XBeam::save(Archive &ar, const unsigned int v) const {
152  ar << boost::serialization::base_object<Bond>(*this);
153  ar << _young_modulus;
154  ar << _damping_factor;
155  ar << _radius;
156  ar << _poisson_ratio;
157  ar << _max_stress;
158  ar << _current_stress;
163  ar << _tensile_energy;
164  ar << _shear_energy;
165  ar << _torsion_energy;
166  ar << _bending_energy;
167  }
168 
169  template<class Archive> void
170  XBeam::load(Archive &ar, const unsigned int v) {
171  ar >> boost::serialization::base_object<Bond>(*this);
172  ar >> _young_modulus;
173  ar >> _damping_factor;
174  ar >> _radius;
175  ar >> _poisson_ratio;
176  ar >> _max_stress;
177  ar >> _current_stress;
182  ar >> _tensile_energy;
183  ar >> _shear_energy;
184  ar >> _torsion_energy;
185  ar >> _bending_energy;
186  }
187 
188  inline void
190 #ifndef SERVER
191  Bond::draw();
192 #endif
193  }
194 
195  }
196 }
197 
198 #include <boost/serialization/version.hpp>
200 
201 namespace boost
202 {
203 
204  namespace serialization
205  {
206 
207  template<class Archive> void
208  save_construct_data(Archive &ar, const GranOO3::DEM::XBeam *t, const unsigned int) {
209  const GranOO3::DEM::Element* el1 = &t->get_element1();
210  const GranOO3::DEM::Element* el2 = &t->get_element2();
211  ar << el1;
212  ar << el2;
213  }
214 
215  template<class Archive> void
216  load_construct_data(Archive &ar, GranOO3::DEM::XBeam *t, const unsigned int) {
217  GranOO3::DEM::Element *el1 = nullptr;
218  GranOO3::DEM::Element *el2 = nullptr;
219  ar >> el1;
220  ar >> el2;
221  ::new(t)GranOO3::DEM::XBeam(*el1, *el2);
222  }
223 
224  }
225 }
226 
227 namespace GranOO3 {
228  GRANOO_CLASS_DECLARE_TPL(DEM::XBeam);
229 }
230 
231 #endif
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
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
Definition: Pair.hpp:202
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