MechanicalFemPhysics
The GeMA Mechanical FEM Physics Plugin
gmpMaterialDMCreep.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 
25 #ifndef _GEMA_PLUGIN_MECHANICALMATERIAL_DMCREEP_H_
26 #define _GEMA_PLUGIN_MECHANICALMATERIAL_DMCREEP_H_
27 
28 #include "gmpMechanicalConfig.h"
29 #include <gmpFemPhysics.h>
30 #include <gmMathUtils.h>
31 
32 #include <gmpMaterialCreep.h>
33 #include "gmpMechanicPoint.h"
34 
35 class GMP_MECHANICAL_PHYSICS_API_EXPORT GmpMaterialDMCreep : public GmpMaterialCreep
36 {
37 protected:
40  {
41  QT_ID = GmpMaterialCreep::ElementPropertyIds::NUM_PROPERTY_IDS,
42  R_ID,
48 
49  // --- NO ADDING BELOW THIS LINE
51  };
52 
53 public:
54 
56  GmpMaterialDMCreep(int typeIndex, QString typeName, const GmLogCategory& logger)
57  : GmpMaterialCreep(typeIndex, typeName, logger) {}
58 
60  virtual ~GmpMaterialDMCreep() {}
61 
63  static GmpFemPhysicsCommonMaterial* instance(GmSimulationData* simulation, int typeIndex, QString typeName, const GmLogCategory& logger)
64  {
65  Q_UNUSED(simulation);
66  return new GmpMaterialDMCreep(typeIndex, typeName, logger);
67  }
68 
69  // Returns a map with material associated properties.
70  virtual const QVariantMap* materialMetaDataMap()
71  {
72  S_TRACE();
73  static QVariantMap m;
74  if (m.isEmpty())
75  {
76  // read the material attribute map from elastic interface material
78 
79  // Adds material properties
80  GmpFemPhysicsCommon::ValueList matDMProp = m["properties"].value<GmpFemPhysicsCommon::ValueList>();
81  matDMProp.append(GmpFemPhysicsCommon::ScalarValue(QT_ID, "QT", QObject::tr("Thermal activation energy"), "J/mol", true));
82  matDMProp.append(GmpFemPhysicsCommon::ScalarValue(R_ID, "R", QObject::tr("Universal gas constant"), "J/(mol*K)", true));
83  matDMProp.append(GmpFemPhysicsCommon::ScalarValue(DEDT_ID, "DEDT", QObject::tr("Threshold creep rate"), "1/s", true));
84  matDMProp.append(GmpFemPhysicsCommon::ScalarValue(So_ID, "So", QObject::tr("Threshold deviatoric stress"), "kPa", true));
85  matDMProp.append(GmpFemPhysicsCommon::ScalarValue(To_ID, "To", QObject::tr("Threshold temperature"), "K", true));
86  matDMProp.append(GmpFemPhysicsCommon::ScalarValue(AN1_ID, "AN1", QObject::tr("Dislocation creep"), "", true));
87  matDMProp.append(GmpFemPhysicsCommon::ScalarValue(AN2_ID, "AN2", QObject::tr("Steady state cracking"), "", true));
88  m["properties"] = matDMProp;
89  }
90  return &m;
91  }
92 
93  // Returns the stresses according to the material behavior adopted
94  virtual double fillCreepStrainRate(const GmElement* e, const GmpMechanicPoint* mp, const GmVector* coord, const GmVector& Time, double MISES) const
95  {
96  S_TRACE();
97  assert(e); Q_UNUSED(Time);
98 
99  double eCreepRate, NP, So, Eo, QT, Rg, To, TEMP;
100 
101  // Reads DM creep Properties
102  So = thresholdDevStress(e, coord, mp->_index);
103  Eo = thresholdCreepRate(e, coord, mp->_index);
104  QT = thermalActEnergy(e, coord, mp->_index);
105  Rg = universalGasConstant(e, coord, mp->_index);
106  To = thresholdTemperature(e, coord, mp->_index);
107 
108  // Read the Temperature at integration point
109  if (nodeAttrAc(T_NA_ID)!= NULL)
110  {
111  TEMP = fillTemperatureFromNodalAttr(e, coord);
112  }
113  else
114  {
115  TEMP = temperatureProperty(e, coord, mp->_index);
116  }
117 
118  // Set the suitable stress power
119  if (MISES < So)
120  {
121  NP = disCreep(e, coord, mp->_index);
122  }
123  else
124  {
125  NP = steadyStateCracking(e, coord, mp->_index);
126  }
127 
128  // Computes the creep strain rate
129  eCreepRate = Eo*exp(QT / (Rg*To))*exp((-1.0)*QT / (Rg*TEMP))* pow(MISES / So, NP);
130 
131  return eCreepRate;
132  }
133 
135  double thermalActEnergy(const GmElement* e, const GmVector* coord, int ip) const
136  {
137  return propertyAc(QT_ID)->scalarValueAt(e, coord, ip);
138  }
139 
141  double universalGasConstant(const GmElement* e, const GmVector* coord, int ip) const
142  {
143  return propertyAc(R_ID)->scalarValueAt(e, coord, ip);
144  }
145 
147  double thresholdCreepRate(const GmElement* e, const GmVector* coord, int ip) const
148  {
149  return propertyAc(DEDT_ID)->scalarValueAt(e, coord, ip);
150  }
151 
153  double thresholdDevStress(const GmElement* e, const GmVector* coord, int ip) const
154  {
155  return propertyAc(So_ID)->scalarValueAt(e, coord, ip);
156  }
157 
159  double thresholdTemperature(const GmElement* e, const GmVector* coord, int ip) const
160  {
161  return propertyAc(To_ID)->scalarValueAt(e, coord, ip);
162  }
163 
165  double disCreep(const GmElement* e, const GmVector* coord, int ip) const
166  {
167  return propertyAc(AN1_ID)->scalarValueAt(e, coord, ip);
168  }
169 
171  double steadyStateCracking(const GmElement* e, const GmVector* coord, int ip) const
172  {
173  return propertyAc(AN2_ID)->scalarValueAt(e, coord, ip);
174  }
175 };
176 
177 #endif
The number of property ids above.
Definition: gmpMaterialDMCreep.h:50
ElementPropertyIds
IDs for material element properties.
Definition: gmpMaterialCreep.h:39
double disCreep(const GmElement *e, const GmVector *coord, int ip) const
Returns the Dislocation creep.
Definition: gmpMaterialDMCreep.h:165
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: gmpMaterialDMCreep.h:63
Id for retrieving the Dislocation creep accessor.
Definition: gmpMaterialDMCreep.h:46
#define S_TRACE()
Definition: gmpMaterialDMCreep.h:35
Definition: gmpMaterialCreep.h:35
QString tr(const char *sourceText, const char *disambiguation, int n)
double steadyStateCracking(const GmElement *e, const GmVector *coord, int ip) const
Returns the Steady state cracking.
Definition: gmpMaterialDMCreep.h:171
Definition: gmpMechanicPoint.h:32
Id for retrieving the Universal gas constant accessor.
Definition: gmpMaterialDMCreep.h:42
double universalGasConstant(const GmElement *e, const GmVector *coord, int ip) const
Returns the Universal gas constant R.
Definition: gmpMaterialDMCreep.h:141
virtual double fillTemperatureFromNodalAttr(const GmElement *e, const GmVector *coord) const
Return temperature at integration point using the nodal Temperature field.
Definition: gmpMaterialCreep.cpp:400
GmpMaterialDMCreep(int typeIndex, QString typeName, const GmLogCategory &logger)
Constructor. Gets as parameters the material index and its name.
Definition: gmpMaterialDMCreep.h:56
double thresholdTemperature(const GmElement *e, const GmVector *coord, int ip) const
Returns the test temperature To.
Definition: gmpMaterialDMCreep.h:159
virtual double temperatureProperty(const GmElement *e, const GmVector *coord, int ip) const
Returns the Temperature from property table.
Definition: gmpMaterialCreep.h:124
Id for retrieving the Threshold deviatoric stress accessor.
Definition: gmpMaterialDMCreep.h:44
Id for retrieving the Temperature accessor.
Definition: gmpMaterialDMCreep.h:45
double thermalActEnergy(const GmElement *e, const GmVector *coord, int ip) const
Returns the Thermal activation energy QT.
Definition: gmpMaterialDMCreep.h:135
double thresholdCreepRate(const GmElement *e, const GmVector *coord, int ip) const
Returns the Threshold creep rate.
Definition: gmpMaterialDMCreep.h:147
virtual ~GmpMaterialDMCreep()
Virtual destructor.
Definition: gmpMaterialDMCreep.h:60
virtual double fillCreepStrainRate(const GmElement *e, const GmpMechanicPoint *mp, const GmVector *coord, const GmVector &Time, double MISES) const
Computes the creep deformation rate.
Definition: gmpMaterialDMCreep.h:94
Declaration of the GmpMechanicPoint class.
Declaration of usefull configuration definitions for the plugin library.
Id for retrieving the Threshold creep rate accessor.
Definition: gmpMaterialDMCreep.h:43
double thresholdDevStress(const GmElement *e, const GmVector *coord, int ip) const
Returns the Threshold deviatoric stress.
Definition: gmpMaterialDMCreep.h:153
Id for retrieving the Steady state cracking accessor.
Definition: gmpMaterialDMCreep.h:47
arma::vec GmVector
virtual const QVariantMap * materialMetaDataMap()
Returns a pointer to the material attribute map, built when the function is called for the first time...
Definition: gmpMaterialDMCreep.h:70
virtual const QVariantMap * materialMetaDataMap()
Returns a pointer to the material attribute map, built when the function is called for the first time...
Definition: gmpMaterialCreep.cpp:38