MechanicalFemPhysics
The GeMA Mechanical FEM Physics Plugin
gmpMaterialCreep.h
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 
25 #ifndef _GEMA_PLUGIN_MECHANICALMATERIAL_CREEP_H_
26 #define _GEMA_PLUGIN_MECHANICALMATERIAL_CREEP_H_
27 
28 #include "gmpMechanicalConfig.h"
29 #include <gmpFemPhysics.h>
30 #include <gmMathUtils.h>
31 
32 #include <gmpMaterialElastic.h>
33 #include "gmpMechanicPoint.h"
34 
35 class GMP_MECHANICAL_PHYSICS_API_EXPORT GmpMaterialCreep : public GmpMechanicalMaterialElastic
36 {
37 protected:
40  {
45 
46  // --- NO ADDING BELOW THIS LINE
48  };
49 
52  {
55 
56  // ------ NO ADDING BELOW THIS LINE
58  };
59 
62  {
64 
65  // ------ NO ADDING BELOW THIS LINE
67  };
68 
71  {
75 
76  // --- NO ADDING BELOW THIS LINE
78  };
79 
80 public:
81 
83  GmpMaterialCreep(int typeIndex, QString typeName, const GmLogCategory& logger)
84  : GmpMechanicalMaterialElastic(typeIndex, typeName, logger) {}
85 
87  virtual ~GmpMaterialCreep() {}
88 
89  // Returns a map with material associated properties.
90  virtual const QVariantMap* materialMetaDataMap();
91 
92  // Checks material loaded data.
93  virtual bool checkLoadedData(const GmElement* e) const;
94 
95  // Returns the stresses according to the material behavior adopted
96  virtual bool mechanicalConstitutiveModel(const GmElement* e, GmMatrix& Dep, const GmpMechanicPoint* mp, const GmVector* coord, const GmVector& Time, unsigned nc, bool ips) const;
97 
99  virtual double yieldStrengthRatio(const GmElement* e, const GmVector& S, const GmVector* coord, int ip, unsigned sc) const;
100 
102  virtual double fillTemperatureFromNodalAttr(const GmElement* e, const GmVector* coord) const;
103 
105  virtual double fillCreepStrainRate(const GmElement* e, const GmpMechanicPoint* mp, const GmVector* coord, const GmVector& Time, double Smises) const = 0;
106 
108  virtual void evaluatesCreepStrain(const GmElement* e, const GmpMechanicPoint* mp, const GmVector* coord, const GmVector& Stress, const GmVector& Time, GmVector& dStrainCreep, unsigned sc) const;
109 
111  virtual bool isIsotropic() const { return true; }
112 
114  virtual double elasticModulus(const GmElement* e, const GmVector* coord, int ip) const
115  {
116  return propertyAc(E_ID)->scalarValueAt(e, coord, ip);
117  }
119  virtual double poissonRatio(const GmElement* e, const GmVector* coord, int ip) const
120  {
121  return propertyAc(NU_ID)->scalarValueAt(e, coord, ip);
122  }
124  virtual double temperatureProperty(const GmElement* e, const GmVector* coord, int ip) const
125  {
126  return propertyAc(T_ID)->scalarValueAt(e, coord, ip);
127  }
129  virtual double stressTolerance(const GmElement* e, const GmVector* coord, int ip) const
130  {
131  S_TRACE();
132  assert(e);
133  GmCellAccessor* sTolAcc = propertyAc(STOL_ID);
134  // Reads accessor of Stress tolerance
135  if (sTolAcc == NULL)
136  {
137  return 1.0e-4;
138  }
139  return sTolAcc->scalarValueAt(e, coord, ip);
140  }
142  virtual int dilatancyCriterion(const GmElement* e, const GmVector* coord, int ip) const
143  {
144  S_TRACE(); assert(e);
145  GmCellAccessor* DilCrtAcc = propertyAc(DILATCRIT_ID);
146  if (DilCrtAcc == NULL)
147  {
148  return 0;
149  }
150  return DilCrtAcc->scalarValueAt(e, coord, ip);
151  }
153  virtual GmVector dilatancyParameters(const GmElement* e, const GmVector* coord, int ip) const
154  {
155  S_TRACE();
156  GmVector data;
157  GmCellAccessor* DPAcc = propertyAc(DILATPARAM_ID);
158  if (DPAcc != NULL) {
159  vectorPropertyValue(e, DPAcc, data, coord, ip);
160  }
161  return data;
162  }
163 };
164 
165 #endif
Damage criterion proposed by DeVries.
Definition: gmpMaterialCreep.h:73
Id for retrieving the dilatancy parameters.
Definition: gmpMaterialCreep.h:44
Id for retrieving the creep strain accessor at the previous state (old creep strain)
Definition: gmpMaterialCreep.h:54
virtual double elasticModulus(const GmElement *e, const GmVector *coord, int ip) const
Returns the material elastic modulus.
Definition: gmpMaterialCreep.h:114
The number of property ids above.
Definition: gmpMaterialElastic.h:46
Base Id for node attribute(s) used to store the calculated stress.
Definition: gmpMaterialCreep.h:63
ElementPropertyIds
IDs for material element properties.
Definition: gmpMaterialCreep.h:39
GaussAttributeIds
IDs for material Gauss attributes.
Definition: gmpMaterialCreep.h:51
The number of node attribute ids above.
Definition: gmpMaterialCreep.h:66
virtual int dilatancyCriterion(const GmElement *e, const GmVector *coord, int ip) const
Returns the dilatancy criterion.
Definition: gmpMaterialCreep.h:142
virtual double poissonRatio(const GmElement *e, const GmVector *coord, int ip) const
Returns the material poisson ratio.
Definition: gmpMaterialCreep.h:119
#define S_TRACE()
Definition: gmpMaterialElastic.h:34
Definition: gmpMaterialCreep.h:35
Definition: gmpMechanicPoint.h:32
The number of property ids above.
Definition: gmpMaterialCreep.h:47
virtual bool isIsotropic() const
Returns true if the material is isotropic, false otherwise.
Definition: gmpMaterialCreep.h:111
GmpMaterialCreep(int typeIndex, QString typeName, const GmLogCategory &logger)
Constructor. Gets as parameters the material index and its name.
Definition: gmpMaterialCreep.h:83
virtual ~GmpMaterialCreep()
Virtual destructor.
Definition: gmpMaterialCreep.h:87
virtual double temperatureProperty(const GmElement *e, const GmVector *coord, int ip) const
Returns the Temperature from property table.
Definition: gmpMaterialCreep.h:124
virtual double stressTolerance(const GmElement *e, const GmVector *coord, int ip) const
Returns the Stress tolerance.
Definition: gmpMaterialCreep.h:129
dilatancyCriterion
Dilatancy Criterion.
Definition: gmpMaterialCreep.h:70
The number of dilatancy criterion.
Definition: gmpMaterialCreep.h:77
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
Id for retrieving the stress tolerance.
Definition: gmpMaterialCreep.h:42
NodeAttributeIds
IDs for creep material node attributes.
Definition: gmpMaterialCreep.h:61
Declaration of the GmpMechanicPoint class.
Declaration of the GmpMechanicalMaterialElastic class.
Declaration of usefull configuration definitions for the plugin library.
The number of Gauss attribute ids above.
Definition: gmpMaterialCreep.h:57
Damage criterion proposed by VanSambeek.
Definition: gmpMaterialCreep.h:72
virtual GmVector dilatancyParameters(const GmElement *e, const GmVector *coord, int ip) const
Returns thedilatancy parameters.
Definition: gmpMaterialCreep.h:153
arma::vec GmVector
Damage criterion proposed by HunsCheMod.
Definition: gmpMaterialCreep.h:74
arma::mat GmMatrix
Base Id for Gauss attribute(s) used to store the calculated creep strain at the current state.
Definition: gmpMaterialCreep.h:53
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
Id for retrieving the dilatancy criterion.
Definition: gmpMaterialCreep.h:43