MechanicalFemPhysics
The GeMA Mechanical FEM Physics Plugin
gmpMaterialDPCreep.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_DPCREEP_H_
26 #define _GEMA_PLUGIN_MECHANICALMATERIAL_DPCREEP_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 // Double power creep material
36 class GMP_MECHANICAL_PHYSICS_API_EXPORT GmpMaterialDPCreep : public GmpMaterialCreep
37 {
38 protected:
41  {
42  A1_ID = GmpMaterialCreep::ElementPropertyIds::NUM_PROPERTY_IDS,
49  R_ID,
50 
51  // --- NO ADDING BELOW THIS LINE
53  };
54 
55 public:
56 
58  GmpMaterialDPCreep(int typeIndex, QString typeName, const GmLogCategory& logger)
59  : GmpMaterialCreep(typeIndex, typeName, logger) {}
60 
62  virtual ~GmpMaterialDPCreep() {}
63 
65  static GmpFemPhysicsCommonMaterial* instance(GmSimulationData* simulation, int typeIndex, QString typeName, const GmLogCategory& logger)
66  {
67  Q_UNUSED(simulation);
68  return new GmpMaterialDPCreep(typeIndex, typeName, logger);
69  }
70 
71  // Returns a map with material associated properties.
72  virtual const QVariantMap* materialMetaDataMap()
73  {
74  S_TRACE();
75  static QVariantMap m;
76  if (m.isEmpty())
77  {
78  // read the material attribute map from elastic interface material
80 
81  // Adds material properties
82  m["properties"] = m["properties"].value<GmpFemPhysicsCommon::ValueList>()
83  << GmpFemPhysicsCommon::ScalarValue(A1_ID, "A1", QObject::tr("Structure factor for DCL"), "1/s", true)
84  << GmpFemPhysicsCommon::ScalarValue(AN1_ID, "n1", QObject::tr("Stress power for DCL"), "", true)
85  << GmpFemPhysicsCommon::ScalarValue(Q1_ID, "Q1", QObject::tr("Thermal activation energy for DCL"), "J/mol", true)
86  << GmpFemPhysicsCommon::ScalarValue(A2_ID, "A2", QObject::tr("Structure factor for UMC"), "1/s", true)
87  << GmpFemPhysicsCommon::ScalarValue(AN2_ID, "n2", QObject::tr("Stress power for UMC"), "", true)
88  << GmpFemPhysicsCommon::ScalarValue(Q2_ID, "Q2", QObject::tr("Thermal activation energy for UMC"), "J/mol", true)
89  << GmpFemPhysicsCommon::ScalarValue(So_ID, "So", QObject::tr("Threshold deviatoric stress"), "kPa", true)
90  << GmpFemPhysicsCommon::ScalarValue(R_ID, "R", QObject::tr("Universal gas constant"), "J/(mol*K)", true);
91  }
92  return &m;
93  }
94 
95  // Returns the stresses according to the material behavior adopted
96  virtual double fillCreepStrainRate(const GmElement* e, const GmpMechanicPoint* mp, const GmVector* coord, const GmVector& Time, double MISES) const
97  {
98  S_TRACE();
99  assert(e); Q_UNUSED(Time);
100 
101  double eCreepRate, A1, Q1, N1, A2, Q2, N2, So, R, TEMP;
102 
103  // Reads Double Power creep Properties
104 
105  A1 = structureFactorDCL(e, coord, mp->_index);
106  Q1 = thermalActivEnergyDCL(e, coord, mp->_index);
107  N1 = stressPowerDCL(e, coord, mp->_index);
108  A2 = structureFactorUMC(e, coord, mp->_index);
109  Q2 = thermalActivEnergyUMC(e, coord, mp->_index);
110  N2 = stressPowerUMC(e, coord, mp->_index);
111  So = thresholdDevStress(e, coord, mp->_index);
112  R = universalGasConstant(e, coord, mp->_index);
113 
114  // Read the Temperature at integration point
115  if (nodeAttrAc(T_NA_ID)!= NULL)
116  {
117  TEMP = fillTemperatureFromNodalAttr(e, coord);
118  }
119  else
120  {
121  TEMP = temperatureProperty(e, coord, mp->_index);
122  }
123  // Computes the creep strain rate
124  eCreepRate = A1 * exp(-Q1 / (R*TEMP))*pow(MISES / So, N1) + A2 * exp(-Q2 / (R*TEMP))*pow(MISES / So, N2);
125 
126  return eCreepRate;
127  }
128 
130  double structureFactorDCL(const GmElement* e, const GmVector* coord, int ip) const
131  {
132  return propertyAc(A1_ID)->scalarValueAt(e, coord, ip);
133  }
134 
136  double thermalActivEnergyDCL(const GmElement* e, const GmVector* coord, int ip) const
137  {
138  return propertyAc(Q1_ID)->scalarValueAt(e, coord, ip);
139  }
140 
142  double stressPowerDCL(const GmElement* e, const GmVector* coord, int ip) const
143  {
144  return propertyAc(AN1_ID)->scalarValueAt(e, coord, ip);
145  }
146 
148  double structureFactorUMC(const GmElement* e, const GmVector* coord, int ip) const
149  {
150  return propertyAc(A2_ID)->scalarValueAt(e, coord, ip);
151  }
152 
154  double thermalActivEnergyUMC(const GmElement* e, const GmVector* coord, int ip) const
155  {
156  return propertyAc(Q2_ID)->scalarValueAt(e, coord, ip);
157  }
158 
160  double stressPowerUMC(const GmElement* e, const GmVector* coord, int ip) const
161  {
162  return propertyAc(AN2_ID)->scalarValueAt(e, coord, ip);
163  }
164 
166  double thresholdDevStress(const GmElement* e, const GmVector* coord, int ip) const
167  {
168  return propertyAc(So_ID)->scalarValueAt(e, coord, ip);
169  }
170 
172  double universalGasConstant(const GmElement* e, const GmVector* coord, int ip) const
173  {
174  return propertyAc(R_ID)->scalarValueAt(e, coord, ip);
175  }
176 
177 };
178 
179 #endif
double thermalActivEnergyDCL(const GmElement *e, const GmVector *coord, int ip) const
Returns the Thermal activation energy for DCL.
Definition: gmpMaterialDPCreep.h:136
ElementPropertyIds
IDs for material element properties.
Definition: gmpMaterialCreep.h:39
Id for retrieving the universal gas constant accessor.
Definition: gmpMaterialDPCreep.h:49
#define S_TRACE()
virtual const QVariantMap * materialMetaDataMap()
Returns a pointer to the material attribute map, built when the function is called for the first time...
Definition: gmpMaterialDPCreep.h:72
virtual double fillCreepStrainRate(const GmElement *e, const GmpMechanicPoint *mp, const GmVector *coord, const GmVector &Time, double MISES) const
Computes the creep deformation rate.
Definition: gmpMaterialDPCreep.h:96
virtual ~GmpMaterialDPCreep()
Virtual destructor.
Definition: gmpMaterialDPCreep.h:62
Definition: gmpMaterialCreep.h:35
QString tr(const char *sourceText, const char *disambiguation, int n)
The number of property ids above.
Definition: gmpMaterialDPCreep.h:52
double universalGasConstant(const GmElement *e, const GmVector *coord, int ip) const
Returns the Universal gas constant R.
Definition: gmpMaterialDPCreep.h:172
Definition: gmpMechanicPoint.h:32
Id for retrieving the threshold deviatoric stress accessor.
Definition: gmpMaterialDPCreep.h:48
Id for retrieving the structure factor for UMC accessor.
Definition: gmpMaterialDPCreep.h:45
virtual double fillTemperatureFromNodalAttr(const GmElement *e, const GmVector *coord) const
Return temperature at integration point using the nodal Temperature field.
Definition: gmpMaterialCreep.cpp:400
double stressPowerDCL(const GmElement *e, const GmVector *coord, int ip) const
Returns the stress power for DCL.
Definition: gmpMaterialDPCreep.h:142
virtual double temperatureProperty(const GmElement *e, const GmVector *coord, int ip) const
Returns the Temperature from property table.
Definition: gmpMaterialCreep.h:124
Definition: gmpMaterialDPCreep.h:36
double stressPowerUMC(const GmElement *e, const GmVector *coord, int ip) const
Returns the stress power for UMC.
Definition: gmpMaterialDPCreep.h:160
double structureFactorUMC(const GmElement *e, const GmVector *coord, int ip) const
Returns the structure factor for UMC.
Definition: gmpMaterialDPCreep.h:148
Id for retrieving the thermal activation energy for DCL accessor.
Definition: gmpMaterialDPCreep.h:44
Id for retrieving the thermal activation energy for UMC accessor.
Definition: gmpMaterialDPCreep.h:47
Declaration of the GmpMechanicPoint class.
Declaration of usefull configuration definitions for the plugin library.
GmpMaterialDPCreep(int typeIndex, QString typeName, const GmLogCategory &logger)
Constructor. Gets as parameters the material index and its name.
Definition: gmpMaterialDPCreep.h:58
Id for retrieving the stress power for DCL accessor.
Definition: gmpMaterialDPCreep.h:43
double thresholdDevStress(const GmElement *e, const GmVector *coord, int ip) const
Returns the Threshold deviatoric stress.
Definition: gmpMaterialDPCreep.h:166
double thermalActivEnergyUMC(const GmElement *e, const GmVector *coord, int ip) const
Returns the Thermal activation energy for UMC.
Definition: gmpMaterialDPCreep.h:154
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: gmpMaterialDPCreep.h:65
arma::vec GmVector
double structureFactorDCL(const GmElement *e, const GmVector *coord, int ip) const
Returns the structure factor for DCL.
Definition: gmpMaterialDPCreep.h:130
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
Id for retrieving the stress power for UMC accessor.
Definition: gmpMaterialDPCreep.h:46