GranOO  3.0
A robust and versatile workbench to build 3D dynamic simulations based on the Discrete Element Method
Element.hpp
Go to the documentation of this file.
1 // This file is part of GranOO, a workbench for FEM simulation.
2 //
3 // Author(s) : - Jean-luc CHARLES I2M-DuMAS/ENSAM Talence France
4 // <jean-luc.charles@ensam.eu>
5 // - Damien ANDRE SPCTS/ENS Ceramique industrielle, Limoges France
6 // <damien.andre@unilim.fr>
7 // - Jeremie GIRARDOT I2M-DuMAS/ENSAM Talence France
8 // <jeremie.girardot@ensam.eu>
9 // - Cedric Hubert LAMIH/UVHC Valenciennes France
10 // <cedric.hubert@univ-valenciennes.fr>
11 // - Ivan IORDANOFF I2M-DuMAS-MPI/ENSAM Talence France
12 // <ivan.iordanoff@ensam.eu>
13 //
14 // Copyright (C) 2008-2016 JL. CHARLES, D. ANDRE, I. IORDANOFF, J. GIRARDOT
15 //
16 //
17 //
18 //
19 //
20 // This program is free software: you can redistribute it and/or modify
21 // it under the terms of the GNU General Public License as published by
22 // the Free Software Foundation, either version 3 of the License, or
23 // (at your option) any later version.
24 //
25 // This program is distributed in the hope that it will be useful,
26 // but WITHOUT ANY WARRANTY; without even the implied warranty of
27 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
28 // GNU General Public License for more details.
29 //
30 // You should have received a copy of the GNU General Public License
31 // along with this program. If not, see <http://www.gnu.org/licenses/>.
32 
33 #ifndef _LibFEM_Element
34 #define _LibFEM_Element
35 
36 #include <string>
37 #include <iostream>
38 #include <fstream>
39 
40 #include "GranOO3/Physic/ExtendedSetOf/Node.hpp"
41 #include "GranOO3/Physic/Node.hpp"
42 #include "GranOO3/3rdParty/Eigen/Sparse"
43 #include "GranOO3/3rdParty/Eigen/Dense"
44 #include "GranOO3/Common.hpp"
45 
46 namespace GranOO3 {
47 
48  namespace FEM {
49 
50  class Element : public Core::Base, public Core::Register<Element> {
51 
52  public:
53  GRANOO_CLASS(FEM::Element, Core::Base);
54 
55  Element();
56 
57  virtual ~Element();
58 
59  virtual Eigen::VectorXd shape_function_at_coordinate(double xi, double eta, double zeta) = 0;
60  virtual Eigen::MatrixXd shape_function_derivative_at_coordinate(double xi, double eta, double zeta) = 0;
61 
62  virtual unsigned int node_count() = 0;
63  virtual unsigned int dof_per_node_count() = 0;
64  virtual unsigned int gauss_point_count() = 0;
65  virtual std::string name() = 0;
66  virtual std::string vtk_name() = 0;
67  virtual Eigen::MatrixXd d_ndx_at_gauss_point(unsigned int gaussPointIndex);
68  virtual Eigen::MatrixXd shape_function_at_gauss_point(unsigned int gaussPointIndex);
69  virtual Eigen::MatrixXd shape_function_derivative_at_gauss_point(unsigned int gaussPointIndex);
70  virtual Eigen::VectorXd strain_at_gauss_point(unsigned int gaussPointIndex);
71  virtual Eigen::VectorXd internal_force();
72  virtual Eigen::VectorXd internal_force_at_gauss_point(unsigned int gaussPointIndex);
73  virtual Eigen::MatrixXd mass_matrix();
74  virtual Eigen::MatrixXd mass_matrix_at_gauss_point(unsigned int gaussPointIndex);
75 
76  void setup_element();
77  void reset_variable();
78  Eigen::MatrixXd b_matrix_for_d_ndX(Eigen::MatrixXd dNdx);
79  void update_state_at_gauss_point(unsigned int gaussPointIndex, const Eigen::VectorXd &dStrain);
80  void update_variable();
81 
82  inline std::vector<GranOO3::Physic::Node *> &node_list() { return _node_list; };
83 
84  inline const Eigen::MatrixXd &gauss_point() { return _gauss_point; };
85  inline const Eigen::MatrixXd &node_coord() { return _node_coord; };
86  inline const std::vector<unsigned int> &connectivity() { return _connectivity; };
87  inline const Eigen::VectorXd & stress_at_gauss_point(unsigned int igaussPoint) { return _oldStress[igaussPoint]; };
88 
89  public:
90  double _youngModulus;
91  double _poissonRatio;
92  double _density;
93 
94  protected:
95  Eigen::MatrixXd _gauss_point;
96  Eigen::VectorXd _weights;
97  std::vector<GranOO3::Physic::Node *> _node_list;
98 
99  Eigen::MatrixXd _node_coord;
100  std::vector<Eigen::MatrixXd> _strainDisplacementMatrices;
101  std::vector<double> _jacobians;
102  std::vector<Eigen::VectorXd> _oldStrain;
103  std::vector<Eigen::VectorXd> _newStrain;
104  std::vector<Eigen::VectorXd> _oldStress;
105  std::vector<Eigen::VectorXd> _newStress;
106  std::vector<Eigen::VectorXd> _oldStateVariables;
107  std::vector<Eigen::VectorXd> _newStateVariables;
108  std::vector<Eigen::MatrixXd> _ctMatrix;
109  std::vector<unsigned int> _connectivity;
110  };
111  }
112 
113  GRANOO_CLASS_DECLARE_TPL(FEM::Element);
114 }
115 
116 #endif
Definition: Base.hpp:61
Definition: SetOf.hpp:346
Definition: Element.hpp:50
void setup_element()
Definition: Element.cpp:57
const std::vector< unsigned int > & connectivity()
Definition: Element.hpp:86
double _density
Definition: Element.hpp:92
const Eigen::MatrixXd & gauss_point()
Definition: Element.hpp:84
std::vector< Eigen::MatrixXd > _ctMatrix
Definition: Element.hpp:108
Eigen::MatrixXd b_matrix_for_d_ndX(Eigen::MatrixXd dNdx)
Definition: Element.cpp:111
double _poissonRatio
Definition: Element.hpp:91
std::vector< GranOO3::Physic::Node * > & node_list()
Definition: Element.hpp:82
Eigen::MatrixXd _gauss_point
Definition: Element.hpp:95
virtual std::string vtk_name()=0
void update_state_at_gauss_point(unsigned int gaussPointIndex, const Eigen::VectorXd &dStrain)
Definition: Element.cpp:238
virtual Eigen::MatrixXd mass_matrix_at_gauss_point(unsigned int gaussPointIndex)
Definition: Element.cpp:181
Element()
Definition: Element.cpp:45
virtual unsigned int gauss_point_count()=0
virtual std::string name()=0
virtual Eigen::VectorXd internal_force()
Definition: Element.cpp:217
std::vector< Eigen::VectorXd > _newStress
Definition: Element.hpp:105
virtual Eigen::VectorXd internal_force_at_gauss_point(unsigned int gaussPointIndex)
Definition: Element.cpp:229
std::vector< Eigen::VectorXd > _oldStress
Definition: Element.hpp:104
virtual Eigen::MatrixXd shape_function_derivative_at_coordinate(double xi, double eta, double zeta)=0
virtual Eigen::VectorXd strain_at_gauss_point(unsigned int gaussPointIndex)
Definition: Element.cpp:165
virtual Eigen::VectorXd shape_function_at_coordinate(double xi, double eta, double zeta)=0
void reset_variable()
Definition: Element.cpp:94
const Eigen::MatrixXd & node_coord()
Definition: Element.hpp:85
double _youngModulus
Definition: Element.hpp:87
std::vector< Eigen::VectorXd > _newStateVariables
Definition: Element.hpp:107
virtual unsigned int dof_per_node_count()=0
std::vector< double > _jacobians
Definition: Element.hpp:101
virtual Eigen::MatrixXd d_ndx_at_gauss_point(unsigned int gaussPointIndex)
Definition: Element.cpp:134
std::vector< Eigen::VectorXd > _newStrain
Definition: Element.hpp:103
Eigen::VectorXd _weights
Definition: Element.hpp:96
std::vector< unsigned int > _connectivity
Definition: Element.hpp:109
virtual ~Element()
Definition: Element.cpp:51
std::vector< Eigen::VectorXd > _oldStrain
Definition: Element.hpp:102
std::vector< Eigen::VectorXd > _oldStateVariables
Definition: Element.hpp:106
virtual Eigen::MatrixXd shape_function_at_gauss_point(unsigned int gaussPointIndex)
Definition: Element.cpp:149
Eigen::MatrixXd _node_coord
Definition: Element.hpp:99
const Eigen::VectorXd & stress_at_gauss_point(unsigned int igaussPoint)
Definition: Element.hpp:87
void update_variable()
Definition: Element.cpp:269
std::vector< GranOO3::Physic::Node * > _node_list
Definition: Element.hpp:97
std::vector< Eigen::MatrixXd > _strainDisplacementMatrices
Definition: Element.hpp:100
virtual Eigen::MatrixXd mass_matrix()
Definition: Element.cpp:204
virtual Eigen::MatrixXd shape_function_derivative_at_gauss_point(unsigned int gaussPointIndex)
Definition: Element.cpp:157
virtual unsigned int node_count()=0
Definition: Common.hpp:198