MechanicalFemPhysics
The GeMA Mechanical FEM Physics Plugin
gmpMechanicalMaterial.h
Go to the documentation of this file.
1 /************************************************************************
2 **
3 ** Copyright (C) 2014 by Carlos Augusto Teixera Mendes
4 ** All rights reserved.
5 **
6 ** This file is part of the "GeMA" software. It's use should respect
7 ** the terms in the license agreement that can be found together
8 ** with this source code.
9 ** It is provided AS IS, with NO WARRANTY OF ANY KIND,
10 ** INCLUDING THE WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR
11 ** A PARTICULAR PURPOSE.
12 **
13 ************************************************************************/
14 
26 #ifndef _GEMA_PLUGIN_MECHANICAL_MATERIAL_H_
27 #define _GEMA_PLUGIN_MECHANICAL_MATERIAL_H_
28 
29 #include "gmpMechanicalConfig.h"
30 
31 #include <gmpFemPhysicsCommonMaterial.h>
32 #include "gmpFemPhysicsCommon.h"
33 #include "gmpMechanicPoint.h"
34 #include "gmpFemPhysics.h"
35 #include "gmVector.h"
36 #include "gmMatrix.h"
37 #include <gmTrace.h>
38 
41 class GMP_MECHANICAL_PHYSICS_API_EXPORT GmpMechanicalMaterial : public GmpFemPhysicsCommonMaterial
42 {
43 public:
44 
46  GmpMechanicalMaterial(int typeIndex, QString typeName, const GmLogCategory& logger)
47  : GmpFemPhysicsCommonMaterial(typeIndex, typeName, logger)
48  {
49  }
50 
53 
55  virtual bool mechanicalConstitutiveModel(const GmElement* e, GmMatrix& Dep, const GmpMechanicPoint* mp, const GmVector* coord, const GmVector& Time, unsigned nc, bool ips = false) const = 0;
56 
58  virtual void tangentModulus(const GmElement* e, GmMatrix& Dep, const GmpMechanicPoint* mp, const GmVector* coord, unsigned nc, unsigned ips = false) const
59  {
60  // TODO: Revisar para retirar da barra que ainda usa
61  Q_UNUSED(e); Q_UNUSED(Dep); Q_UNUSED(mp); Q_UNUSED(coord); Q_UNUSED(nc); Q_UNUSED(ips);
62  }
63 
65  virtual bool isIsotropic() const = 0;
66 
68  virtual bool setInitialConditions(const GmElement* e, GmpMechanicPoint* mp, const GmVector* coord, unsigned sc) const
69  {
70  Q_UNUSED(e); Q_UNUSED(mp); Q_UNUSED(coord); Q_UNUSED(sc);
71  return true;
72  }
73 
74  // Updates deformation gradient F according to material requirement
75  virtual bool updateDeformationGradient(GmMatrix& F, const GmElement* e, const GmVector* coord, int ip, bool ips = false) const
76  {
77  S_TRACE();
78  Q_UNUSED(F); Q_UNUSED(ips); Q_UNUSED(e); Q_UNUSED(coord); Q_UNUSED(ip);
79  // To check and implement...
80  // Plane strain condition
81  //F.at(2, 2) = 1;
82  return true;
83  }
84 
86  virtual bool calcDerivedResults(const GmElement* e, GmpMechanicPoint* mp, const GmVector* coord, unsigned sc) const
87  {
88  Q_UNUSED(e); Q_UNUSED(mp); Q_UNUSED(coord); Q_UNUSED(sc);
89  return true;
90  }
91 };
92 
94 namespace GmpMechanicUtils
95 {
96  // Returns the Shear modulus
97  double shearModulus(double E, double nu);
98  // Returns the Bulk modulus
99  double bulkModulus(double E, double nu);
100  // Returns the Lame modulus
101  double lameModulus(double E, double nu);
102 
103  //stress state
104  void setStressStateVector(const double[6], GmVector&, unsigned, double = 1.0);
105  void setStressStateMatrix(const double[6][6], GmMatrix&, unsigned, double = 1.0);
106 
107  //stress invariants
108  double stressInvariantI1(const GmVector&, unsigned);
109  double stressInvariantI2(const GmVector&, unsigned);
110  double stressInvariantI3(const GmVector&, unsigned);
111 
112  GmVector stressInvariants(const GmVector&, unsigned);
113 
114  void stressGradientI1(const GmVector&, unsigned, GmVector&);
115  void stressGradientI2(const GmVector&, unsigned, GmVector&);
116  void stressGradientI3(const GmVector&, unsigned, GmVector&);
117 
118  void stressHessianI2(const GmVector&, unsigned, GmMatrix&);
119  void stressHessianI3(const GmVector&, unsigned, GmMatrix&);
120 
121  //principal stress
122  GmVector stressPrincipal(const GmVector&, unsigned);
123 
124  //hydrostatic stress
125  double hydrostaticStress(const GmVector&, unsigned);
126 
127  void hydrostaticGradient(const GmVector&, unsigned);
128 
129  //deviatoric stress
130  GmVector deviatoricTensor(const GmVector& s, unsigned st);
131 
132  // Deviatoric invariant J2
133  double deviatoricInvariantJ2(const GmVector&, unsigned);
134  // Deviatoric invariant J3
135  double deviatoricInvariantJ3(const GmVector&, unsigned);
136  // Deviatoric invariants J = [J1, J2, J3]'
137  GmVector deviatoricInvariants(const GmVector&, unsigned);
138  // Deviatoric stress gradient dJ2/ds
139  void deviatoricGradientJ2(const GmVector&, unsigned, GmVector&);
140  // Deviatoric stress gradient dJ3/ds
141  void deviatoricGradientJ3(const GmVector&, unsigned, GmVector&);
142  // Deviatoric stress hessian ddJ2/dds
143  void deviatoricHessianJ2(const GmVector&, unsigned, GmMatrix&);
144  // Deviatoric stress hessian ddJ3/dds
145  void deviatoricHessianJ3(const GmVector&, unsigned, GmMatrix&);
146 
147  //deviatoric principal
148  GmVector deviatoricPrincipal(const GmVector&, unsigned);
149 
150  //von mises stress
151  double vonMisesStress(const GmVector&, unsigned);
152  // Von-Mises stress gradient
153  void vonMisesGradient(const GmVector&, unsigned, GmVector&);
154  // Von-Mises hessian
155  void vonMisesHessian(const GmVector&, unsigned, GmMatrix&);
156 
157  // Lode angle
158  double lodeAngle(const GmVector&, unsigned);
159  // Lode coordinates
160  GmVector lodeCoordinates(const GmVector&, unsigned);
161  // Lode MC angle
162  double lodeAngleMC(const GmVector&, unsigned);
163 
164  // Isotropic damage matrix
166  // Anisotropic damage matrix
168 
169  // Reorganize from vector to tensorial representation and vice versa
170  void vectorialToTensorial(GmVector& Av, GmMatrix& Am, unsigned sc, QString type = "strain", bool mode = true);
171 
172  // Fills the vectorial rotation matrix
173  void fillTransformationMatrix(GmMatrix R, GmMatrix& Te, GmMatrix& Tei, int ns);
174 
175  // Fills rotation matrix using the angle of dip and strike
176  void fillRotationMatrixFromDip(GmMatrix& R, double dip, double strike);
177 
178  // Fills rotation matrix using The Euler's angles
179  void fillRotationMatrixFromEulerAngle(GmMatrix& R, double phi, double theta, double psi);
180 
181  // Computes the The Macaulay brackets
182  GmVector vectorialMacaulay(GmVector x);
183 
184  // Returns the scalar Macaulay function
185  double macaulayFunction(double x);
186 
187  // Evaluates pseudo-inverse
188  void pseudoInverse(GmMatrix& M, double tol = 0.0);
189 
190  // Computes the elastic stiffness matrix using E and Nu
191  void elasticStiffness(GmMatrix& De, double E, double nu, unsigned sc, unsigned ips);
192 
193  // Computes the elastic stiffness matrix using lame and G
194  void elasticStiffnessLameG(GmMatrix& De, double lame, double G, unsigned sc, unsigned ips);
195 
196  // Computes the elastic stiffness matrix using K and G
197  void elasticStiffnessKG(GmMatrix& De, double K, double G, unsigned sc, unsigned ips);
198 
199  // Computes the elastic flexibility matrix using E and Nu
200  void elasticFlexibility(GmMatrix& Ce, double E, double nu, unsigned sc, unsigned ips);
201 };
202 #endif
void elasticStiffnessKG(GmMatrix &De, double K, double G, unsigned sc, unsigned ips)
Computes the elastic stiffness matrix using K and G.
Definition: gmpMechanicalMaterial.cpp:983
void hydrostaticGradient(const GmVector &, unsigned st, GmVector &dsm)
Returns the hydrostatic pressure gradient dP/ds.
Definition: gmpMechanicalMaterial.cpp:289
GmMatrix anisotropicDamageMatrix(GmVector Sdv)
Fills the anisotropic damage matrix in the local references which is aligned.
Definition: gmpMechanicalMaterial.cpp:611
#define S_TRACE()
Basic class providing the interface for a mechanical material.
Definition: gmpMechanicalMaterial.h:41
double hydrostaticStress(const GmVector &s, unsigned st)
Returns the hydrostatic stress.
Definition: gmpMechanicalMaterial.cpp:282
GmMatrix isotropicDamageMatrix(double d)
Fills the isotropic damage matrix in the local references which is aligned.
Definition: gmpMechanicalMaterial.cpp:601
Groups utilitary routines for working with stress.
Definition: gmpMechanicalMaterial.cpp:40
virtual ~GmpMechanicalMaterial()
Virtual destructor.
Definition: gmpMechanicalMaterial.h:52
Definition: gmpMechanicPoint.h:32
void elasticStiffnessLameG(GmMatrix &De, double lame, double G, unsigned sc, unsigned ips)
Computes the elastic stiffness matrix using lame and G.
Definition: gmpMechanicalMaterial.cpp:946
double shearModulus(double E, double nu)
Returns the shear modulus (G) from Young's modulus (E) and Poisson's Coefficient (nu)
Definition: gmpMechanicalMaterial.cpp:44
virtual bool calcDerivedResults(const GmElement *e, GmpMechanicPoint *mp, const GmVector *coord, unsigned sc) const
Calculates the derived results required by materials.
Definition: gmpMechanicalMaterial.h:86
GmpMechanicalMaterial(int typeIndex, QString typeName, const GmLogCategory &logger)
Constructor. Gets as parameters the material index and its name.
Definition: gmpMechanicalMaterial.h:46
virtual void tangentModulus(const GmElement *e, GmMatrix &Dep, const GmpMechanicPoint *mp, const GmVector *coord, unsigned nc, unsigned ips=false) const
Returns the constitutive tangent matrix.
Definition: gmpMechanicalMaterial.h:58
Declaration of the GmpMechanicPoint class.
Declaration of usefull configuration definitions for the plugin library.
void elasticFlexibility(GmMatrix &Ce, double E, double nu, unsigned sc, unsigned ips)
Computes the elastic flexibility matrix using E and Nu.
Definition: gmpMechanicalMaterial.cpp:1020
virtual bool setInitialConditions(const GmElement *e, GmpMechanicPoint *mp, const GmVector *coord, unsigned sc) const
Sets the initial conditions required by Solid materials.
Definition: gmpMechanicalMaterial.h:68
double lameModulus(double E, double nu)
Returns the lame modulus (lame) from Young's modulus (E) and Poisson's Coefficient (nu)
Definition: gmpMechanicalMaterial.cpp:58
arma::vec GmVector
void elasticStiffness(GmMatrix &De, double E, double nu, unsigned sc, unsigned ips)
Computes the elastic stiffness matrix using E and Nu.
Definition: gmpMechanicalMaterial.cpp:874
arma::mat GmMatrix
double bulkModulus(double E, double nu)
Returns the bulk modulus (K) from Young's modulus (E) and Poisson's Coefficient (nu)
Definition: gmpMechanicalMaterial.cpp:51