MechanicalFemPhysics
The GeMA Mechanical FEM Physics Plugin
gmpMaterialElastoplastic.h
Go to the documentation of this file.
1 /************************************************************************
2 **
3 ** Copyright (C) 2014 by Nilthson NoreƱa
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 
27 #ifndef _GEMA_PLUGIN_MECHANICALMATERIAL_ELASTOPLASTIC_H_
28 #define _GEMA_PLUGIN_MECHANICALMATERIAL_ELASTOPLASTIC_H_
29 
30 #include "gmpMaterialElastic.h"
31 #include "gmpMechanicPoint.h"
32 #include "gmpFemPhysics.h"
33 #include "gmMathUtils.h"
34 
35 class GMP_MECHANICAL_PHYSICS_API_EXPORT GmpMaterialElastoplastic : public GmpMechanicalMaterialElastic
36 {
37 public:
38  //constructors
39  GmpMaterialElastoplastic(int typeIndex, QString typeName, const GmLogCategory& logger)
40  : GmpMechanicalMaterialElastic(typeIndex, typeName, logger) {}
41 
42  //destructor
43  virtual ~GmpMaterialElastoplastic() {}
44 
45  // Returns a map with material associated properties.
46  virtual const QVariantMap* materialMetaDataMap();
47 
50  {
57  };
58 
61  {
62  none,
66  };
67 
70  {
77 
78  // --- NO ADDING BELOW THIS LINE
80  };
81 
84  {
90 
91  // ------ NO ADDING BELOW THIS LINE
93  };
94 
95  //Returns the stresses according to the material behavior adopted
96  virtual bool mechanicalConstitutiveModel(const GmElement*, GmMatrix&, const GmpMechanicPoint*, const GmVector*, const GmVector&, unsigned, bool) const;
97 
98  // Returns the stresses according to the explicit return algorithm
99  virtual bool explicitReturnAlgorithm(const GmElement*, GmMatrix&, const GmpMechanicPoint*, const GmVector*, unsigned, bool) const ;
100 
101  // Returns the stresses according to the cutting plane return algorithm
102  virtual bool semiImplicitReturnAlgorithm(const GmElement*, GmMatrix&, const GmpMechanicPoint*, const GmVector*, unsigned, bool) const;
103 
104  // Returns the stresses according to the implicit return algorithm
105  virtual bool implicitReturnAlgorithm(const GmElement*, GmMatrix&, const GmpMechanicPoint*, const GmVector*, unsigned, bool) const ;
106 
107  // Returns the stresses according to a Newton-Krylov return mapping algorithm
108  virtual bool newtonKrylovReturnAlgorithm(const GmElement*, GmMatrix&, const GmpMechanicPoint*, const GmVector*, unsigned, bool) const;
109 
110  // Returns the stresses according to a Newton-Raphson implicit state-update algorithm based on complementary functions
111  virtual bool cFunctionsReturnAlgorithm(const GmElement*, GmMatrix&, const GmpMechanicPoint*, const GmVector*, unsigned, bool) const;
112 
113  // Returns the stresses according to an alternative return mapping algorithm
114  virtual bool alternativeReturnAlgorithm(const GmElement*, GmMatrix&, const GmpMechanicPoint*, const GmVector*, unsigned, bool) const;
115 
116  // Yield criterion
117  virtual double yieldCriterion(const GmElement* e, const GmVector& S, const GmVector* coord, int ip, unsigned sc) const = 0;
118  // Stress gradient of Yield function
119  virtual bool yieldStressGradient(const GmElement* e, GmVector& dfs, const GmVector& S, const GmVector* coord, int ip, unsigned sc, bool ips) const = 0;
120 
121  // Flow vector (Stress gradient of potencial function)
122  virtual bool flowVector(const GmElement* e, GmVector& dgs, const GmVector& s, const GmVector* coord, int ip, unsigned sc, bool ips) const = 0;
123  // Stress gradient of Flow vector
124  virtual bool flowVectorStressGradient(const GmElement* e, GmMatrix& dgss, const GmVector& s, const GmVector* coord, int ip, unsigned sc, bool ips) const = 0;
125 
126  // Hardening law
127  virtual double hardeningLaw(const GmElement* , const GmpMechanicPoint*, const GmVector*, int, unsigned) const { return 0.0; }
128  // Stress gradient of Hardening law
129  virtual double hardeningStressGradient(const GmElement* , const GmpMechanicPoint*, const GmVector*, int, unsigned) const { return 0.0; };
130  // hardening gradient of hardening law
131  virtual double hardeningHessian(const GmElement* , const GmpMechanicPoint*, const GmVector*, int, unsigned) const { return 0.0; };
132 
134  virtual double yieldStrengthRatio(const GmElement* e, const GmVector& S, const GmVector* coord, int ip, unsigned sc) const = 0;
135 
136  // Line Search for the optimized return algorithm
137  virtual double lineSearch(const GmElement*, GmMatrix&, GmVector&, GmVector&, const GmpMechanicPoint*, const GmVector*,
138  double, double, unsigned, bool) const;
139  // Golden section method
140  virtual double goldenSectionMethod(const GmElement*, GmMatrix&, GmVector&, GmVector&, const GmpMechanicPoint*, const GmVector*,
141  double, double, unsigned, bool) const;
142  // Quadratic interpolation
143  virtual double quadraticInterpolation(const GmElement*, GmMatrix&, GmVector&, GmVector&, const GmpMechanicPoint*, const GmVector*,
144  double, double, unsigned, bool) const;
145  // Quadratic interpolation
146  virtual double cubicInterpolation(const GmElement*, GmMatrix&, GmVector&, GmVector&, const GmpMechanicPoint*, const GmVector*,
147  double, double, unsigned, bool) const;
148 
149  // Residual Function for the optimized return algorithm
150  virtual double residualFunction(const GmElement*, GmMatrix&, GmVector&, const GmpMechanicPoint*, const GmVector*, double, unsigned, bool) const;
151 
153  virtual double returnAlgorithm(const GmElement* e, const GmVector* coord, int ip) const
154  {
155  S_TRACE();
156  assert(e);
157  GmCellAccessor* rmAcc = propertyAc(RMA_ID);
158 
159  // Read return mapping option
160  double value = 2;
161  if (rmAcc != NULL)
162  {
163  value = rmAcc->scalarValueAt(e, coord, ip);
164  }
165 
166  return value;
167  }
168 
170  virtual bool substepping(const GmElement* e, const GmVector* coord, int ip) const
171  {
172  S_TRACE();
173  assert(e);
174  GmCellAccessor* subStepAcc = propertyAc(SUBSTEP_ID);
175 
176  // Read substepping option
177  bool value = false;
178  if (subStepAcc != NULL)
179  {
180  value = (subStepAcc->scalarValueAt(e, coord, ip) > 0)? true : false;
181  }
182 
183  return value;
184  }
185 
187  virtual double lineSearchStrategy(const GmElement* e, const GmVector* coord, int ip) const
188  {
189  S_TRACE();
190  assert(e);
191  GmCellAccessor* lineSearchAcc = propertyAc(LINESEARCH_ID);
192 
193  // Read line search option
194  double value = 1;
195  if (lineSearchAcc != NULL)
196  {
197  value = lineSearchAcc->scalarValueAt(e, coord, ip);
198  }
199 
200  return value;
201  }
202 
204  virtual double yieldTolerance(const GmElement* e, const GmVector* coord, int ip) const
205  {
206  S_TRACE();
207  assert(e);
208  // Reads accessor of Yield function tolerance
209  GmCellAccessor* fTolAcc = propertyAc(FTOL_ID);
210  if (fTolAcc == NULL)
211  {
212  return 1.0e-6;
213  }
214  return fTolAcc->scalarValueAt(e, coord, ip);
215  }
216 
218  virtual double stressTolerance(const GmElement* e, const GmVector* coord, int ip) const
219  {
220  S_TRACE();
221  assert(e);
222  GmCellAccessor* sTolAcc = propertyAc(STOL_ID);
223  // Reads accessor of Stress tolerance
224  if (sTolAcc == NULL)
225  {
226  return 1.0e-6;
227  }
228  return sTolAcc->scalarValueAt(e, coord, ip);
229  }
230 
232  virtual double hardeningTolerance(const GmElement* e, const GmVector* coord, int ip) const
233  {
234  S_TRACE();
235  assert(e);
236  // Reads accessor of hardening tolerance
237  GmCellAccessor* hTolAcc = propertyAc(HTOL_ID);
238  if (hTolAcc == NULL)
239  {
240  return 1.0e-6;
241  }
242  return hTolAcc->scalarValueAt(e, coord, ip);
243  }
244 };
245 
246 #endif
Id for explicit state-update algorithm.
Definition: gmpMaterialElastoplastic.h:51
Id for retrieving the hardening variable tolerance.
Definition: gmpMaterialElastoplastic.h:76
The number of property ids above.
Definition: gmpMaterialElastic.h:46
The number of gauss attributes.
Definition: gmpMaterialElastoplastic.h:92
virtual double returnAlgorithm(const GmElement *e, const GmVector *coord, int ip) const
Returns the type of return mapping algorithm.
Definition: gmpMaterialElastoplastic.h:153
elastoplasticGaussAttrIds
IDs for Gauss attributes of elastoplastic material.
Definition: gmpMaterialElastoplastic.h:83
#define S_TRACE()
Id for Return mapping algorithm with golden section method.
Definition: gmpMaterialElastoplastic.h:63
ElementPropertyIds
IDs for material element properties.
Definition: gmpMaterialElastic.h:38
Definition: gmpMaterialElastic.h:34
virtual bool substepping(const GmElement *e, const GmVector *coord, int ip) const
Returns the substepping algorithm.
Definition: gmpMaterialElastoplastic.h:170
Id for retrieving the plastic strain accessor at the previous state (old plastic strain)
Definition: gmpMaterialElastoplastic.h:86
The number of property ids above.
Definition: gmpMaterialElastoplastic.h:79
Id for retrieving the material point state accessor at the previous state.
Definition: gmpMaterialElastoplastic.h:88
Id for a Newton-Raphson implicit state-update algorithm based on complementary functions.
Definition: gmpMaterialElastoplastic.h:55
Id for retrieving the line search method.
Definition: gmpMaterialElastoplastic.h:73
Id for retrieving the plastic multiplier at the plastic state.
Definition: gmpMaterialElastoplastic.h:89
Definition: gmpMechanicPoint.h:32
Id for Return mapping algorithm with quadratic line search.
Definition: gmpMaterialElastoplastic.h:64
virtual double stressTolerance(const GmElement *e, const GmVector *coord, int ip) const
Returns the Stress tolerance.
Definition: gmpMaterialElastoplastic.h:218
Id for retrieving the yield function tolerance.
Definition: gmpMaterialElastoplastic.h:74
Id for a Newton-Krylov implicit state-update algorithm.
Definition: gmpMaterialElastoplastic.h:54
GmpMechanicalMaterialElastic(int typeIndex, QString typeName, const GmLogCategory &logger)
Constructor. Gets as parameters the material index and its name.
Definition: gmpMaterialElastic.h:52
virtual bool mechanicalConstitutiveModel(const GmElement *e, GmMatrix &Dep, const GmpMechanicPoint *mp, const GmVector *coord, const GmVector &Time, unsigned nc, bool ips) const
Evaluates stress and tangent matrix according to the material behavior adopted.
Definition: gmpMaterialElastic.cpp:83
virtual double yieldTolerance(const GmElement *e, const GmVector *coord, int ip) const
Returns the yield function tolerance.
Definition: gmpMaterialElastoplastic.h:204
virtual double lineSearchStrategy(const GmElement *e, const GmVector *coord, int ip) const
Returns the line search strategy.
Definition: gmpMaterialElastoplastic.h:187
virtual double hardeningTolerance(const GmElement *e, const GmVector *coord, int ip) const
Returns the hardening variable tolerance.
Definition: gmpMaterialElastoplastic.h:232
Id for retrieving the substepping algorithm.
Definition: gmpMaterialElastoplastic.h:72
Id for implicit state-update algorithm.
Definition: gmpMaterialElastoplastic.h:53
Declaration of the GmpMechanicPoint class.
Declaration of the GmpMechanicalMaterialElastic class.
returnMappingStrategy
IDs for return mapping strategies.
Definition: gmpMaterialElastoplastic.h:49
Base Id for Gauss attribute(s) used to store the calculated hardening at the current state.
Definition: gmpMaterialElastoplastic.h:87
Id for Return mapping algorithm without line search.
Definition: gmpMaterialElastoplastic.h:62
ID for Return mapping algorithm with cubic line search.
Definition: gmpMaterialElastoplastic.h:65
Id for semi-implicit state-update algorithm.
Definition: gmpMaterialElastoplastic.h:52
arma::vec GmVector
Id for an alternative state-update algorithm.
Definition: gmpMaterialElastoplastic.h:56
lineSearchStrategy
IDs for line search strategies.
Definition: gmpMaterialElastoplastic.h:60
Id for retrieving the stress tolerance.
Definition: gmpMaterialElastoplastic.h:75
Definition: gmpMaterialElastoplastic.h:35
arma::mat GmMatrix
virtual const QVariantMap * materialMetaDataMap()
Returns a pointer to the material attribute map, built when the function is called for the first time...
Definition: gmpMaterialElastic.cpp:40
Base Id for Gauss attribute(s) used to store the calculated plastic strain at the current state.
Definition: gmpMaterialElastoplastic.h:85