Mka3D  1.0
Discrete Element method for solid mechanics
solide.hpp
Go to the documentation of this file.
1 //Copyright 2017 Laurent Monasse
2 
3 /*
4  This file is part of Mka3D.
5 
6  Mka3D is free software: you can redistribute it and/or modify
7  it under the terms of the GNU General Public License as published by
8  the Free Software Foundation, either version 3 of the License, or
9  (at your option) any later version.
10 
11  Mka3D is distributed in the hope that it will be useful,
12  but WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  GNU General Public License for more details.
15 
16  You should have received a copy of the GNU General Public License
17  along with Mka3D. If not, see <http://www.gnu.org/licenses/>.
18 */
19 
26 #ifndef SOLIDE_HPP
27 #define SOLIDE_HPP
28 
29 #include "geometry.hpp"
30 
31 
32 
34 class Vertex
35 {
36 public:
37  Vertex();
38  Vertex(const Point_3 &p, const std::vector<int> & parts);
39  Vertex & operator=(const Vertex &V); // op&eacute;rateur = surcharge pour l'affectation
40  Point_3 pos;
41  int num;
42 
43  int size(){
44  return particules.size();
45  }
46  std::vector<int> particules;
47 };
48 
50 class Face
51 {
52 public:
53  Face();//:vertex(std::vector<Vertex>(1)){}
54  Face(const std::vector<Vertex> & v, const int &part);
55  Face(const std::vector<Vertex> & v, const int &part, const double &dist);
56  Face & operator=(const Face &F); // op&eacute;rateur = surcharge pour l'affectation
57  int size(){
58  return vertex.size();
59  }
60  void compFaceIntegrals(double &Fa, double &Fb, double &Fc, double &Faa, double &Fbb, double &Fcc, double &Faaa, double &Fbbb, double &Fccc, double &Faab, double &Fbbc, double &Fcca, const double& na, const double& nb, const double& nc, const int& a, const int& b, const int& c);
61  void compProjectionIntegrals(double &P1, double &Pa, double &Pb, double &Paa, double &Pab, double &Pbb, double &Paaa, double &Paab, double &Pabb, double &Pbbb, const int &a, const int &b, const int &c);
62  void Inertie();
63  Point_3 centre;
64  Vector_3 normale;
65  double S; //Surface de la face
66  double Is;
67  double It;
68  Vector_3 s;
69  Vector_3 t;
70  std::vector<Vertex> vertex;
71  int voisin;
72  double D0;
73 };
74 
75 
77 class Particule
78 {
79 
80  public:
81 
82  Particule();//:faces(std::vector<Face>(1)){}
83 
84  Particule(const double &x_min, const double &y_min, const double &z_min,
85  const double &x_max, const double &y_max,const double &z_max);
86 
87  Particule(const Point_3 &c, const double &x_min, const double &y_min, const double &z_min,
88  const double &x_max, const double &y_max,const double &z_max,
89  const std::vector<Face> & F);
90  ~Particule();
91  Particule & operator=(const Particule &P); // opérateur = surcharge pour l'affectation
92  void Affiche(); //fonction auxilaire utile pour les tests
93  double volume();
94  void CompVolumeIntegrals(double &T1, double &Tx, double &Ty, double &Tz, double &Txx, double &Tyy, double &Tzz, double &Txy, double &Tyz, double &Tzx);
95  void Inertie(const double &rho);
96  void Volume_libre();
97 void solve_position(const double &dt, const bool &flag_2d);
98 void solve_vitesse(const double &dt, const bool &flag_2d);
99  Vector_3 vitesse_parois(const Point_3& X_f);
100  Vector_3 vitesse_parois_prev(const Point_3& X_f);
101  //double min_x; //!< la plus petite coordonn&eacute;e de la particule selon x
102  //double min_y; //!< la plus petite coordonn&eacute;e de la particule selon y
103  //double min_z; //!< la plus petite coordonn&eacute;e de la particule selon z
104  //double max_x; //!< la plus grande coordonn&eacute;e de la particule selon x
105  //double max_y; //!< la plus petite coordonn&eacute;e de la particule selon y
106  //double max_z; //!< la plus petite coordonn&eacute;e de la particule selon z
107 bool cube;
108 Bbox bbox;
109 
110  std::vector<Face> faces;
111 
112  std::vector<Point_3> vertices;
113 
117  Triangles triangles;
118 
122  Triangles triangles_prev;
123 
127  std::vector<Vector_3> normales;
128 
132  std::vector<Vector_3> normales_prev;
133 
137  std::vector<bool> vide;
138 
141  std::vector<bool> fluide;
142 
146  std::vector<bool> fluide_prev;
147 
150  std::vector< std::vector<Point_3> > Points_interface;
151 
155  std::vector< std::vector<Point_3> > Points_interface_prev;
156 
159  std::vector< std::vector<Triangle_3> > Triangles_interface;
160 
164  std::vector< std::vector< std::vector<int> > > Position_Triangles_interface;
165 
169  std::vector< std::vector<Triangle_3> > Triangles_interface_prev;
170 
174  std::vector< std::vector<std::vector<int> > > Position_Triangles_interface_prev;
175 
176  int fixe;
177  double m;
178  double V;
179  double Vl;
180  double epsilon;
181  double I[3];
182  double rotref[3][3];
183  Point_3 x0;
184  Vector_3 Dx;
185  Vector_3 Dxprev;
186  Vector_3 Fi;
187 
190  Vector_3 Ff;
191 
194  Vector_3 Ffprev;
195  Vector_3 Mi;
196 
199  Vector_3 Mf;
200 
203  Vector_3 Mfprev;
204  Vector_3 u;
205  Vector_3 u_half;
206  Vector_3 omega;
207  Vector_3 omega_half;
208  Vector_3 e;
209  Vector_3 eref;
210  Vector_3 eprev;
211  Aff_transformation_3 mvt_t;
212  Aff_transformation_3 mvt_tprev;
213 };
214 
216 class Solide
217 {
218 
219 public:
220 
221  Solide();//:solide(std::vector<Particule>(1)){}
222  Solide(const std::vector<Particule> & Part);
223  ~Solide();
224  Solide & operator=(const Solide &S); // opérateur = surcharge pour l'affectation
225  void Affiche(); //fonction auxilaire utile pour les test
226  int size(){
227  return solide.size();
228  }
229  void Impression(const int &n, const bool &reconstruction);
230  void Init(const char* s, const bool &rep, const int &numrep, const double &rho);
231  void Solve_position(const double &dt, const bool &flag_2d);
232  void Solve_vitesse(const double &dt, const bool &flag_2d);
233  void Forces(const int &N_dim, const double &nu, const double &E);
234  void Forces_internes(const int &N_dim, const double &nu, const double &E);
235  void update_triangles();
236  //void breaking_criterion();
237  double Energie(const int &N_dim, const double &nu, const double &E);
238  double Energie_potentielle(const int &N_dim, const double &nu, const double &E);
239  double Energie_cinetique();
240  double pas_temps(const double &t, const double &T, const double &cfls, const double &E, const double &nu, const double &rhos);
241  // private :
242  std::vector<Particule> solide;
243 };
244 
245 
246 #endif
std::vector< Vertex > vertex
Les sommets de la face.
Definition: solide.hpp:70
std::vector< std::vector< Point_3 > > Points_interface_prev
Liste de points d&#39;intersections de Particule.triangles_prev avec la grille fluide au temps t-dt...
Definition: solide.hpp:155
Vector_3 Dxprev
Déplacement du centre de la particule en t-dt.
Definition: solide.hpp:185
Vector_3 u_half
Vitesse de la particule au temps t-dt/2.
Definition: solide.hpp:205
Point_3 pos
Coordonnées du sommet.
Definition: solide.hpp:40
std::vector< std::vector< std::vector< int > > > Position_Triangles_interface_prev
index de la cellule où se trouve Triangles_interface au temps t-dt
Definition: solide.hpp:174
std::vector< bool > fluide
=true si Particule.triangles en contact avec le fluide
Definition: solide.hpp:141
int size()
Definition: solide.hpp:57
Aff_transformation_3 mvt_t
Transformation affine de la particule au temps t.
Definition: solide.hpp:211
std::vector< std::vector< Triangle_3 > > Triangles_interface_prev
Triangulation des Particule.triangles_prev au temps t.
Definition: solide.hpp:169
Vertex & operator=(const Vertex &V)
operator =: affectation overload.
Definition: solide.cpp:75
int size()
Definition: solide.hpp:43
Vector_3 Ff
Forces fluides exercées sur le solide entre t et t+dt/2.
Definition: solide.hpp:190
double S
Definition: solide.hpp:65
double V
Volume de la particule.
Definition: solide.hpp:178
Point_3 x0
Position du centre de la particule à t=0.
Definition: solide.hpp:183
Point_3 centre
Centre de la face.
Definition: solide.hpp:63
Vector_3 u
Vitesse de la particule au temps t.
Definition: solide.hpp:204
Vertex()
Default constructor.
Definition: solide.cpp:50
Vector_3 Mf
Moments fluides exercés sur le solide entre t et t+dt/2.
Definition: solide.hpp:199
Vector_3 Ffprev
Forces fluides exercées sur le solide entre t-dt/2 et t.
Definition: solide.hpp:194
Triangles triangles_prev
Triangulation des faces de la particule au temps t-dt.
Definition: solide.hpp:122
double It
Second moment d&#39;inertie de la face.
Definition: solide.hpp:67
std::vector< Vector_3 > normales
normales extérieures aux Particule.triangles
Definition: solide.hpp:127
Vector_3 Mi
Moments intérieurs du solide.
Definition: solide.hpp:195
std::vector< std::vector< std::vector< int > > > Position_Triangles_interface
index de la cellule où se trouve Triangles_interface au temps t
Definition: solide.hpp:164
int num
Numéro du point dans le maillage de construction.
Definition: solide.hpp:41
std::vector< int > particules
Vecteur de particules auxquelles pos appartient.
Definition: solide.hpp:46
Aff_transformation_3 mvt_tprev
Transformation affine de la particule au temps t-dt.
Definition: solide.hpp:212
Vector_3 eprev
Vecteur de rotation de la particule au temps t-dt.
Definition: solide.hpp:210
std::vector< bool > fluide_prev
=true si Particule.triangles_prev en contact avec le fluide
Definition: solide.hpp:146
int size()
Definition: solide.hpp:226
bool cube
= true si la particule est un cube, false sinon
Definition: solide.hpp:107
Vector_3 omega_half
Vecteur rotation au temps t-dt/2.
Definition: solide.hpp:207
double Is
Premier moment d&#39;inertie de la face.
Definition: solide.hpp:66
Vector_3 Dx
Déplacement du centre de la particule en t.
Definition: solide.hpp:184
double D0
Distance à l&#39;équilibre avec la particule voisine.
Definition: solide.hpp:72
double m
Masse de la particule.
Definition: solide.hpp:177
std::vector< bool > vide
=true si Particule.triangles en contact avec le fluide
Definition: solide.hpp:137
Vector_3 s
Vecteur selon le premier axe principal d&#39;inertie de la face.
Definition: solide.hpp:68
int voisin
Le numéro de la particule voisine. -1 si le voisin est le fluide.
Definition: solide.hpp:71
std::vector< Vector_3 > normales_prev
normales extérieures aux Particule.triangles_prev
Definition: solide.hpp:132
double epsilon
Déformation volumique globale de la particule.
Definition: solide.hpp:180
Définition de la classe Particule.
Definition: solide.hpp:77
std::vector< std::vector< Triangle_3 > > Triangles_interface
Triangulation des Particule.triangles au temps t.
Definition: solide.hpp:159
Définition de la classe Solide.
Definition: solide.hpp:216
Vector_3 e
Vecteur de rotation de la particule au temps t.
Definition: solide.hpp:208
double Vl
Volume libre de la particule (pour le calcul d&#39;epsilon)
Definition: solide.hpp:179
Bbox bbox
Definition: solide.hpp:108
Vector_3 eref
Definition: solide.hpp:209
std::vector< std::vector< Point_3 > > Points_interface
Liste de points d&#39;intersections de Particule.triangles avec la grille fluide au temps t...
Definition: solide.hpp:150
Vector_3 t
Vecteur selon le second axe principal d&#39;inertie de la face.
Definition: solide.hpp:69
std::vector< Point_3 > vertices
liste des sommets de la particule
Definition: solide.hpp:112
Vector_3 Fi
Forces intérieures du solide.
Definition: solide.hpp:186
int fixe
=true si la particule est fixée, false sinon
Definition: solide.hpp:176
Vector_3 omega
Vecteur rotation au temps t.
Definition: solide.hpp:206
Triangles triangles
Triangulation des faces de la particule au temps t.
Definition: solide.hpp:117
std::vector< Particule > solide
Maillage solide.
Definition: solide.hpp:242
Définition de la classe Face.
Definition: solide.hpp:50
Vector_3 Mfprev
Moments fluides exercés sur le solide entre t-dt/2 et t.
Definition: solide.hpp:203
Définition de la classe Vertex.
Definition: solide.hpp:34
Vector_3 normale
Normale sortante à la face.
Definition: solide.hpp:64
std::vector< Face > faces
liste de faces de la particule
Definition: solide.hpp:110