MechanicalFemPhysics
The GeMA Mechanical FEM Physics Plugin
gmpMaterialIsotropicDamage.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 
27 #ifndef _GEMA_PLUGIN_MECHANICALMATERIAL_ISOTROPICDAMAGE_H_
28 #define _GEMA_PLUGIN_MECHANICALMATERIAL_ISOTROPICDAMAGE_H_
29 
30 #include "gmpMechanicalConfig.h"
31 #include "gmpFemPhysics.h"
32 #include "gmMathUtils.h"
33 
34 #include "gmpMaterialElastic.h"
35 #include "gmpMechanicPoint.h"
36 
37 
38 class GMP_MECHANICAL_PHYSICS_API_EXPORT GmpMaterialIsotropicDamage : public GmpMechanicalMaterialElastic
39 {
40 protected:
43  {
44  // Damage Law
49  A_ID,
50  B_ID,
52  // Damage Criterion
54 
55  // --- NO ADDING BELOW THIS LINE
57  };
58 
59 public:
60 
62  GmpMaterialIsotropicDamage(int typeIndex, QString typeName, const GmLogCategory& logger)
63  : GmpMechanicalMaterialElastic(typeIndex, typeName, logger) {}
64 
67 
69  static GmpFemPhysicsCommonMaterial* instance(GmSimulationData* simulation, int typeIndex, QString typeName, const GmLogCategory& logger)
70  {
71  Q_UNUSED(simulation);
72  return new GmpMaterialIsotropicDamage(typeIndex, typeName, logger);
73  }
74 
75  // Returns a map with material associated properties.
76  virtual const QVariantMap* materialMetaDataMap();
77 
78  // Checks material loaded data.
79  virtual bool checkLoadedData(const GmElement* e) const;
80 
83  {
89 
90  // ------ NO ADDING BELOW THIS LINE
92  };
93 
96  {
100 
101  // --- NO ADDING BELOW THIS LINE
103  };
104 
107  {
113 
114  // --- NO ADDING BELOW THIS LINE
116  };
117 
118  // Sets the initial conditions required by isotropic damage material
119  virtual bool setInitialConditions(const GmElement* e, GmpMechanicPoint* mp, const GmVector* coord, unsigned sc) const;
120 
121  // Returns the stresses according to the material behavior adopted
122  virtual bool mechanicalConstitutiveModel(const GmElement* e, GmMatrix& Ded, const GmpMechanicPoint* mp, const GmVector* coord, const GmVector& Time, unsigned nc, bool ips) const;
123 
125  virtual int damageCriterion(const GmElement* e, const GmVector* coord, int ip) const
126  {
127  S_TRACE(); assert(e);
128  GmCellAccessor* DCritAcc = propertyAc(DAMCRITERION_ID);
129  if (DCritAcc == NULL)
130  {
131  return 0;
132  }
133  return DCritAcc->scalarValueAt(e, coord, ip);
134  }
135 
137  virtual int damageLaw(const GmElement* e, const GmVector* coord, int ip) const
138  {
139  S_TRACE(); assert(e);
140  GmCellAccessor* DLawAcc = propertyAc(DAMLAW_ID);
141  if (DLawAcc == NULL)
142  {
143  return 0;
144  }
145  return DLawAcc->scalarValueAt(e, coord, ip);
146  }
147 
149  virtual double damageThreshold(const GmElement* e, const GmVector* coord, int ip) const
150  {
151  return propertyAc(R0_ID)->scalarValueAt(e, coord, ip);
152  }
153 
155  virtual double maximumInternalVariable(const GmElement* e, const GmVector* coord, int ip) const
156  {
157  return propertyAc(RMAX_ID)->scalarValueAt(e, coord, ip);
158  }
159 
161  virtual double parameterA(const GmElement* e, const GmVector* coord, int ip) const
162  {
163  return propertyAc(A_ID)->scalarValueAt(e, coord, ip);
164  }
165 
167  virtual double parameterB(const GmElement* e, const GmVector* coord, int ip) const
168  {
169  return propertyAc(B_ID)->scalarValueAt(e, coord, ip);
170  }
171 
173  virtual double fractureEnergy(const GmElement* e, const GmVector* coord, int ip) const
174  {
175  return propertyAc(GF_ID)->scalarValueAt(e, coord, ip);
176  }
177 
179  virtual double ratioCompressiveTensile(const GmElement* e, const GmVector* coord, int ip) const
180  {
181  return propertyAc(KAPPA_ID)->scalarValueAt(e, coord, ip);
182  }
183 
185  virtual bool damageCriterion(const GmElement*, const GmpMechanicPoint*, const GmVector*, unsigned, GmVector, double&, GmVector&, unsigned) const;
186 
188  virtual bool damageLaw(const GmElement*, const GmVector*, const GmpMechanicPoint*, unsigned, GmVector, double, GmVector, double&, double&) const;
189 
191  virtual bool engineeringToCauchyShearStrain(unsigned, GmVector&) const;
192 
194  virtual bool tensorialProduct(bool, unsigned, GmVector, GmVector&) const;
195 
197  virtual double characteristicLength(const GmElement*, GmValueAccessor*, bool type = false) const;
198 };
199 
200 #endif
Damage evolution proposed by Oliver.
Definition: gmpMaterialIsotropicDamage.h:112
damageLaw
Damage evolution law.
Definition: gmpMaterialIsotropicDamage.h:106
Id for retrieving the material point state accessor at the previous state.
Definition: gmpMaterialIsotropicDamage.h:87
GmpMaterialIsotropicDamage(int typeIndex, QString typeName, const GmLogCategory &logger)
Constructor. Gets as parameters the material index and its name.
Definition: gmpMaterialIsotropicDamage.h:62
The number of property ids above.
Definition: gmpMaterialElastic.h:46
virtual double parameterB(const GmElement *e, const GmVector *coord, int ip) const
Returns the exponential or polynomial parameter B.
Definition: gmpMaterialIsotropicDamage.h:167
The number of gauss attributes.
Definition: gmpMaterialIsotropicDamage.h:91
Id for retrieving the exponential or polynomial parameter B accessor.
Definition: gmpMaterialIsotropicDamage.h:50
Definition: gmpMaterialIsotropicDamage.h:38
virtual ~GmpMaterialIsotropicDamage()
Virtual destructor.
Definition: gmpMaterialIsotropicDamage.h:66
Id for retrieving the fracture energy accessor.
Definition: gmpMaterialIsotropicDamage.h:51
Damage criterion based on modified Von Mises.
Definition: gmpMaterialIsotropicDamage.h:98
Id for retrieving the damage law accessor.
Definition: gmpMaterialIsotropicDamage.h:46
The number of property ids above.
Definition: gmpMaterialIsotropicDamage.h:56
#define S_TRACE()
ElementPropertyIds
IDs for material element properties.
Definition: gmpMaterialElastic.h:38
Definition: gmpMaterialElastic.h:34
Damage evolution proposed by Kurumatani.
Definition: gmpMaterialIsotropicDamage.h:111
virtual int damageCriterion(const GmElement *e, const GmVector *coord, int ip) const
Returns the damage criterion.
Definition: gmpMaterialIsotropicDamage.h:125
Id for retrieving the damage threshold accessor.
Definition: gmpMaterialIsotropicDamage.h:47
virtual double parameterA(const GmElement *e, const GmVector *coord, int ip) const
Returns the exponential or polynomial parameter A.
Definition: gmpMaterialIsotropicDamage.h:161
Id for retrieving the exponential or polynomial parameter A accessor.
Definition: gmpMaterialIsotropicDamage.h:49
elasticDamageGaussAttrIds
IDs for Gauss attributes of elastic damage material.
Definition: gmpMaterialIsotropicDamage.h:82
Definition: gmpMechanicPoint.h:32
Id for retrieving the ratio between the uniaxial compressive strength and the uniaxial tensile streng...
Definition: gmpMaterialIsotropicDamage.h:53
Linear evolution.
Definition: gmpMaterialIsotropicDamage.h:108
Id for retrieving the geostatic strain accessor.
Definition: gmpMaterialIsotropicDamage.h:88
Exponential evolution.
Definition: gmpMaterialIsotropicDamage.h:109
virtual double fractureEnergy(const GmElement *e, const GmVector *coord, int ip) const
Returns the fracture energy.
Definition: gmpMaterialIsotropicDamage.h:173
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
static GmpFemPhysicsCommonMaterial * instance(GmSimulationData *simulation, int typeIndex, QString typeName, const GmLogCategory &logger)
A "factory" function used to register the material with the physics material factory.
Definition: gmpMaterialIsotropicDamage.h:69
Damage criterion proposed by Mazars.
Definition: gmpMaterialIsotropicDamage.h:99
Base Id for Gauss attribute(s) used to store the damage variable at the current state.
Definition: gmpMaterialIsotropicDamage.h:84
The number of damage initiation criterion.
Definition: gmpMaterialIsotropicDamage.h:115
virtual int damageLaw(const GmElement *e, const GmVector *coord, int ip) const
Returns the damage law.
Definition: gmpMaterialIsotropicDamage.h:137
virtual double ratioCompressiveTensile(const GmElement *e, const GmVector *coord, int ip) const
Returns the ratio between the uniaxial compressive strength and the uniaxial tensile strength.
Definition: gmpMaterialIsotropicDamage.h:179
Declaration of the GmpMechanicPoint class.
Declaration of the GmpMechanicalMaterialElastic class.
damageCriterion
Damage Criterion.
Definition: gmpMaterialIsotropicDamage.h:95
Id for retrieving the maximum admissible internal variable acessor.
Definition: gmpMaterialIsotropicDamage.h:48
Declaration of usefull configuration definitions for the plugin library.
The number of damage initiation criterion.
Definition: gmpMaterialIsotropicDamage.h:102
virtual double damageThreshold(const GmElement *e, const GmVector *coord, int ip) const
Returns the damage threshold.
Definition: gmpMaterialIsotropicDamage.h:149
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
Id for retrieving the material point state accessor at the previous state.
Definition: gmpMaterialIsotropicDamage.h:85
Base Id for Gauss attribute(s) used to store the internal state variable at the current state.
Definition: gmpMaterialIsotropicDamage.h:86
Polynomial evolution.
Definition: gmpMaterialIsotropicDamage.h:110
virtual double maximumInternalVariable(const GmElement *e, const GmVector *coord, int ip) const
Returns the maximum admissible internal variable.
Definition: gmpMaterialIsotropicDamage.h:155
arma::vec GmVector
Damage criterion based on normal strain.
Definition: gmpMaterialIsotropicDamage.h:97
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