MechanicalFemPhysics
The GeMA Mechanical FEM Physics Plugin
gmpMaterialCapGeneralYieldSurface.h
Go to the documentation of this file.
1 /************************************************************************
2 **
3 ** Copyright (C) 2016 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_CAP_GENERAL_YIELD_SURFACE_H_
28 #define _GEMA_PLUGIN_MECHANICALMATERIAL_CAP_GENERAL_YIELD_SURFACE_H_
29 
30 #include "gmpMaterialElastoplasticMultiSurface.h"
31 #include "gmpMechanicPoint.h"
32 #include "gmpFemPhysics.h"
33 #include "gmMathUtils.h"
34 
35 class GMP_MECHANICAL_PHYSICS_API_EXPORT GmpMaterialCapGeneralYieldSurface : public GmpMaterialElastoplasticMultiSurface
36 {
37 protected:
40  {
47  A_ID,
52 
53  // --- NO ADDING BELOW THIS LINE
55  };
56 
59  {
61 
62  // ------ NO ADDING BELOW THIS LINE
64  };
65 
66 public:
67 
69  GmpMaterialCapGeneralYieldSurface(int typeIndex, QString typeName, const GmLogCategory& logger)
70  : GmpMaterialElastoplasticMultiSurface(typeIndex, typeName, logger) {}
71 
74 
76  static GmpFemPhysicsCommonMaterial* instance(GmSimulationData* simulation, int typeIndex, QString typeName, const GmLogCategory& logger)
77  {
78  Q_UNUSED(simulation);
79  return new GmpMaterialCapGeneralYieldSurface(typeIndex, typeName, logger);
80  }
81 
82  // Returns a map with material associated properties.
83  virtual const QVariantMap* materialMetaDataMap();
84 
85  // Returns the stresses according to the material behavior adopted
86  virtual bool mechanicalConstitutiveModel(const GmElement* c, GmMatrix& Dep, const GmpMechanicPoint* mp, const GmVector* coord, const GmVector& Time, unsigned sc, bool ips) const;
87 
88  // Sets the initial conditions required by Cap model material
89  virtual bool setInitialConditions(const GmElement* e, GmpMechanicPoint* mp, const GmVector* coord, unsigned sc) const;
90 
92  virtual bool isIsotropic() const { return true; }
93 
95  virtual double yieldStrengthRatio(const GmElement* e, const GmVector& S, const GmVector* coord, int ip, unsigned sc) const;
96 
98  virtual double cohesion(const GmElement* e, const GmVector* coord, int ip) const
99  {
100  return propertyAc(COH_ID)->scalarValueAt(e, coord, ip);
101  }
103  virtual double frictionAngle(const GmElement* e, const GmVector* coord, int ip) const
104  {
105  return propertyAc(PHI_ID)->scalarValueAt(e, coord, ip);
106  }
108  virtual double dilationAngle(const GmElement* e, const GmVector* coord, int ip) const
109  {
110  return propertyAc(PSI_ID)->scalarValueAt(e, coord, ip);
111  }
113  virtual double firstShapeFactor(const GmElement* e, const GmVector* coord, int ip) const
114  {
115  return propertyAc(ALPHA_ID)->scalarValueAt(e, coord, ip);
116  }
118  virtual double secondShapeFactor(const GmElement* e, const GmVector* coord, int ip) const
119  {
120  return propertyAc(BETA_ID)->scalarValueAt(e, coord, ip);
121  }
123  virtual double thirdShapeFactor(const GmElement* e, const GmVector* coord, int ip) const
124  {
125  return propertyAc(GAMMA_ID)->scalarValueAt(e, coord, ip);
126  }
128  virtual double fourthShapeFactor(const GmElement* e, const GmVector* coord, int ip) const
129  {
130  return propertyAc(A_ID)->scalarValueAt(e, coord, ip);
131  }
133  virtual double initialCenterEllipse(const GmElement* e, const GmVector* coord, int ip) const
134  {
135  return propertyAc(INITIALCENTER_ID)->scalarValueAt(e, coord, ip);
136  }
138  virtual double ratioSemiaxes(const GmElement* e, const GmVector* coord, int ip) const
139  {
140  return propertyAc(RATIO_ID)->scalarValueAt(e, coord, ip);
141  }
142 
144  virtual double maxVolumetricStrain(const GmElement* e, const GmVector* coord, int ip) const
145  {
146  S_TRACE(); assert(e);
147  // Get accessor for the ultimate compressive volumetric strain
148  GmCellAccessor* wAccc = propertyAc(MAXVSTRAIN_ID);
149 
150  if (wAccc == NULL)
151  {
152  return 1.0e10;
153  }
154 
155  return wAccc->scalarValueAt(e, coord, ip);
156  }
157 
159  virtual double volumetricStrainRate(const GmElement* e, const GmVector* coord, int ip) const
160  {
161  S_TRACE(); assert(e);
162  // Get accessor for the volumetric strain rate with respect to the compressive hydrostatic strain
163  GmCellAccessor* dAccc = propertyAc(VSTRAINRATE_ID);
164 
165  if (dAccc == NULL)
166  {
167  return 1.0e-3;
168  }
169 
170  return dAccc->scalarValueAt(e, coord, ip);
171  }
172 
173  // Yield criterion - Conical Surface
174  virtual double yieldCriterion(const GmElement*, const GmVector&, const GmVector*, int, unsigned) const;
175  virtual bool yieldStressGradient(const GmElement*, GmVector&, const GmVector&, const GmVector*, int, unsigned, bool) const;
176  virtual bool yieldGradient_Hardening(const GmElement*, GmVector&, const GmVector&, const GmVector*, int, unsigned) const;
177  virtual bool yieldHessian(const GmElement*, GmMatrix&, const GmVector&, const GmVector*, int, unsigned, bool) const;
178 
179  // Yield criterion - Cap Surface
180  virtual double capSurfaceCriterion(const GmElement*, const GmVector&, const double, const GmVector*, int, unsigned) const;
181  virtual bool capSurfaceGradient(const GmElement*, GmVector&, const GmVector&, const double, const GmVector*, int, unsigned, bool) const;
182  virtual bool capSurfaceGradient_Hardening(const GmElement*, GmVector&, const GmVector&, const double, const GmVector*, int, unsigned) const;
183  virtual bool capSurfaceHessian(const GmElement*, GmMatrix&, const GmVector&, const double, const GmVector*, int, unsigned, bool) const;
184 
185  // Plastic potential - Conical Surface
186  virtual double plasticFPotential(const GmElement*, const GmVector&, const GmVector*, int, unsigned) const;
187  virtual bool flowVector(const GmElement*, GmVector&, const GmVector&, const GmVector*, int, unsigned, bool) const;
188  virtual bool flowVectorStressGradient(const GmElement*, GmMatrix&, const GmVector&, const GmVector*, int, unsigned, bool) const;
189  virtual bool flowVectorDerivative_Hardening(const GmElement*, GmVector&, const GmVector&, const GmVector*, int, unsigned) const;
190 
191  // Plastic potential - Cap Surface
192  virtual double capSurfacePlasticFPotential(const GmElement*, const GmVector&, const double, const GmVector*, int, unsigned) const;
193  virtual bool capSurfacePlasticFGradient(const GmElement*, GmVector&, const GmVector&, const double, const GmVector*, int, unsigned, bool) const;
194  virtual bool capSurfacePlasticFHessian(const GmElement*, GmMatrix&, const GmVector&, const double, const GmVector*, int, unsigned, bool) const;
195  virtual bool capSurfaceFlowVectorDerivative_Hardening(const GmElement*, GmVector&, const GmVector&, const double, const GmVector*, int, unsigned, bool) const;
196 
197  // Hardening Evolution - Conical Surface
198  virtual bool hardeningEvolution(const GmElement*, GmVector&, const GmVector&, const double, const GmVector*, int, unsigned) const;
199  virtual bool hardeningEvolutionDerivative_Stress(const GmElement*, GmVector&, const GmVector&, const double, const GmVector*, int, unsigned) const;
200  virtual bool hardeningEvolutionDerivative_Hardening(const GmElement*, GmMatrix&, const GmVector&, const double, const GmVector*, int, unsigned) const;
201 
202  // Hardening Evolution - Cap Surface
203  virtual bool capSurfaceHardeningEvolution(const GmElement*, GmVector&, const GmVector&, const double, const GmVector*, int, unsigned) const;
204  virtual bool capSurfaceHardeningEvolutionDerivative_Stress(const GmElement*, GmVector&, const GmVector&, const double, const GmVector*, int, unsigned, bool) const;
205  virtual bool capSurfaceHardeningEvolutionDerivative_Hardening(const GmElement*, GmMatrix&, const GmVector&, const double, const GmVector*, int, unsigned) const;
206 
207  // Complementary Methods
208  virtual bool yieldSurfaces(const GmElement*, GmVector&, GmVector, const double, const GmVector*, int, unsigned) const;
209  virtual bool flowVectors(const GmElement*, GmVector&, int, GmVector, const double, const GmVector*, int, unsigned, bool) const;
210  virtual bool hardeningLaws(const GmElement*, GmVector&, int, GmVector, const double, const GmVector*, int, unsigned) const;
211  virtual bool gradientVectors_Stress(const GmElement*, GmVector&, int, GmVector, const double, const GmVector*, int, unsigned, bool) const;
212  virtual bool gradientVectors_Hardening(const GmElement*, GmVector&, int, GmVector, const double, const GmVector*, int, unsigned) const;
213  virtual bool flowVectorDerivatives_Stress(const GmElement*, GmMatrix&, int, GmVector, const double, const GmVector*, int, unsigned, bool) const;
214  virtual bool flowVectorDerivatives_Hardening(const GmElement*, GmVector&, int, GmVector, const double, const GmVector*, int, unsigned, bool) const;
215  virtual bool hardeningLawDerivatives_Stress(const GmElement*, GmVector&, int, GmVector, const double, const GmVector*, int, unsigned, bool) const;
216  virtual bool hardeningLawDerivatives_Hardening(const GmElement*, GmMatrix&, int, GmVector, const double q, const GmVector*, int, unsigned) const;
217  void updateDLamb(GmVector&, GmVector, GmVector) const;
218  virtual int numberActiveSurface(GmVector, const double) const;
219  void indexActiveSurfaces(GmVector&, GmVector, int, const double) const;
220  void deactivateSurfaces(GmVector&, GmVector, int) const;
221 
222  // Line Search
223  virtual double lineSearch(const GmElement*, GmMatrix, GmVector, GmVector, GmVector, GmVector, const double,
224  const double, const double, GmVector, GmVector, const GmpMechanicPoint*, const GmVector*, unsigned, bool) const;
225  virtual double goldenSectionMethod(const GmElement*, GmMatrix, GmVector, GmVector, GmVector, GmVector, const double,
226  const double, const double, GmVector, GmVector, const GmpMechanicPoint*, const GmVector*, unsigned, bool) const;
227  virtual double quadraticInterpolation(const GmElement*, GmMatrix, GmVector, GmVector, GmVector, GmVector, const double,
228  const double, const double, GmVector, GmVector, const GmpMechanicPoint*, const GmVector*, unsigned, bool) const;
229  virtual double cubicInterpolation(const GmElement*, GmMatrix, GmVector, GmVector, GmVector, GmVector, const double,
230  const double, const double, GmVector, GmVector, const GmpMechanicPoint*, const GmVector*, unsigned, bool) const;
231 
232  // Residual Function
233  virtual double residualFunction(const GmElement*, GmMatrix, const GmpMechanicPoint*, const GmVector*, GmVector,
234  GmVector, GmVector, const double, const double, GmVector, unsigned, bool) const;
235 
236 };
237 
238 #endif
The number of Gauss attribute ids above.
Definition: gmpMaterialCapGeneralYieldSurface.h:63
The number of property ids above.
Definition: gmpMaterialCapGeneralYieldSurface.h:54
virtual double dilationAngle(const GmElement *e, const GmVector *coord, int ip) const
Returns the material dilation angle.
Definition: gmpMaterialCapGeneralYieldSurface.h:108
The number of gauss attributes.
Definition: gmpMaterialElastoplastic.h:92
virtual double thirdShapeFactor(const GmElement *e, const GmVector *coord, int ip) const
Returns the third shape factor.
Definition: gmpMaterialCapGeneralYieldSurface.h:123
virtual double ratioSemiaxes(const GmElement *e, const GmVector *coord, int ip) const
Returns the ratio between the semiaxes.
Definition: gmpMaterialCapGeneralYieldSurface.h:138
Id for retrieving the ultimate compressive volumetric strain.
Definition: gmpMaterialCapGeneralYieldSurface.h:50
Id for retrieving the ratio between the semiaxes.
Definition: gmpMaterialCapGeneralYieldSurface.h:49
Id for retrieving the third shape factor.
Definition: gmpMaterialCapGeneralYieldSurface.h:46
#define S_TRACE()
Definition: gmpMaterialElastoplasticMultiSurface.h:33
The number of property ids above.
Definition: gmpMaterialElastoplastic.h:79
Id for retrieving the friction angle accessor.
Definition: gmpMaterialCapGeneralYieldSurface.h:42
virtual double secondShapeFactor(const GmElement *e, const GmVector *coord, int ip) const
Returns the second shape factor.
Definition: gmpMaterialCapGeneralYieldSurface.h:118
Definition: gmpMechanicPoint.h:32
virtual double cohesion(const GmElement *e, const GmVector *coord, int ip) const
Returns the material cohesion.
Definition: gmpMaterialCapGeneralYieldSurface.h:98
virtual double fourthShapeFactor(const GmElement *e, const GmVector *coord, int ip) const
Returns the fourth shape factor.
Definition: gmpMaterialCapGeneralYieldSurface.h:128
Id for retrieving the initial center of the ellipse.
Definition: gmpMaterialCapGeneralYieldSurface.h:48
Id for retrieving the fourth shape factor.
Definition: gmpMaterialCapGeneralYieldSurface.h:47
virtual double yieldStrengthRatio(const GmElement *e, const GmVector &S, const GmVector *coord, int ip, unsigned sc) const =0
Returns the Yield Strength Ratio (Ysr)
GmpMaterialCapGeneralYieldSurface(int typeIndex, QString typeName, const GmLogCategory &logger)
Constructor. Gets as parameters the material index and its name.
Definition: gmpMaterialCapGeneralYieldSurface.h:69
ElementPropertyIds
IDs for material element properties.
Definition: gmpMaterialCapGeneralYieldSurface.h:39
Id for retrieving the second shape factor.
Definition: gmpMaterialCapGeneralYieldSurface.h:45
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: gmpMaterialCapGeneralYieldSurface.h:76
Declaration of the GmpMechanicPoint class.
virtual double maxVolumetricStrain(const GmElement *e, const GmVector *coord, int ip) const
Returns the ultimate compressive volumetric strain.
Definition: gmpMaterialCapGeneralYieldSurface.h:144
virtual ~GmpMaterialCapGeneralYieldSurface()
Virtual destructor.
Definition: gmpMaterialCapGeneralYieldSurface.h:73
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
virtual bool isIsotropic() const
Returns true if the material is isotropic, false otherwise.
Definition: gmpMaterialCapGeneralYieldSurface.h:92
Id for retrieving the first shape factor.
Definition: gmpMaterialCapGeneralYieldSurface.h:44
virtual double initialCenterEllipse(const GmElement *e, const GmVector *coord, int ip) const
Returns the center of the ellipse.
Definition: gmpMaterialCapGeneralYieldSurface.h:133
virtual bool mechanicalConstitutiveModel(const GmElement *, GmMatrix &, const GmpMechanicPoint *, const GmVector *, const GmVector &, unsigned, bool) const
Returns the updated stresses after the return mapping process.
Definition: gmpMaterialElastoplastic.cpp:82
Id for retrieving the dilatance angle accessor.
Definition: gmpMaterialCapGeneralYieldSurface.h:43
GaussAttributeIds
IDs for material Gauss attributes.
Definition: gmpMaterialCapGeneralYieldSurface.h:58
arma::vec GmVector
virtual double firstShapeFactor(const GmElement *e, const GmVector *coord, int ip) const
Returns the first shape factor.
Definition: gmpMaterialCapGeneralYieldSurface.h:113
Id for retrieving the volumetric strain rate with respect to the compressive hydrostatic strain.
Definition: gmpMaterialCapGeneralYieldSurface.h:51
arma::mat GmMatrix
virtual double volumetricStrainRate(const GmElement *e, const GmVector *coord, int ip) const
Returns the volumetric strain rate with respect to the compressive hydrostatic strain.
Definition: gmpMaterialCapGeneralYieldSurface.h:159
virtual const QVariantMap * materialMetaDataMap()
Returns a pointer to the material/Gauus attribute map, built when the function is called for the firs...
Definition: gmpMaterialElastoplastic.cpp:45
Definition: gmpMaterialCapGeneralYieldSurface.h:35