MechanicalFemPhysics
The GeMA Mechanical FEM Physics Plugin
gmpMaterialCapModel.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 
29 #ifndef _GEMA_PLUGIN_MECHANICALMATERIAL_CAP_MODEL_H_
30 #define _GEMA_PLUGIN_MECHANICALMATERIAL_CAP_MODEL_H_
31 
32 #include "gmpMaterialElastoplasticMultiSurface.h"
33 #include "gmpMechanicPoint.h"
34 #include "gmpFemPhysics.h"
35 #include "gmMathUtils.h"
36 
37 class GMP_MECHANICAL_PHYSICS_API_EXPORT GmpMaterialCapModel : public GmpMaterialElastoplasticMultiSurface
38 {
39 protected:
42  {
52 
53  // --- NO ADDING BELOW THIS LINE
55  };
56 
59  {
61 
62  // ------ NO ADDING BELOW THIS LINE
64  };
65 
66 public:
67 
69  GmpMaterialCapModel(int typeIndex, QString typeName, const GmLogCategory& logger)
70  : GmpMaterialElastoplasticMultiSurface(typeIndex, typeName, logger) {}
71 
73  virtual ~GmpMaterialCapModel() {}
74 
76  static GmpFemPhysicsCommonMaterial* instance(GmSimulationData* simulation, int typeIndex, QString typeName, const GmLogCategory& logger)
77  {
78  Q_UNUSED(simulation);
79  return new GmpMaterialCapModel(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  // Returns the stresses according to the implicit return algorithm
89  virtual bool implicitReturnAlgorithm(const GmElement* c, GmMatrix& Dep, const GmpMechanicPoint* mp, const GmVector* coord, unsigned sc, bool ips) const;
90 
91  // Returns the stresses according to a Newton-Raphson implicit state-update algorithm based on complementary functions
92  virtual bool cFunctionsReturnAlgorithm(const GmElement* c, GmMatrix& Dep, const GmpMechanicPoint* mp, const GmVector* coord, unsigned sc, bool ips) const;
93 
94  // Returns the stresses according to a Newton-Krylov implicit state-update algorithm based on complementary functions
95  virtual bool newtonKrylovReturnAlgorithm(const GmElement* c, GmMatrix& Dep, const GmpMechanicPoint* mp, const GmVector* coord, unsigned sc, bool ips) const;
96 
97  // Sets the initial conditions required by Cap model material
98  virtual bool setInitialConditions(const GmElement* e, GmpMechanicPoint* mp, const GmVector* coord, unsigned sc) const;
99 
101  virtual double yieldStrengthRatio(const GmElement* e, const GmVector& S, const GmVector* coord, int ip, unsigned sc) const;
102 
104  virtual bool isIsotropic() const { return true; }
105 
107  virtual double alphaParameter(const GmElement* e, const GmVector* coord, int ip) const
108  {
109  return propertyAc(ALPHA_ID)->scalarValueAt(e, coord, ip);
110  }
112  virtual double betaParameter(const GmElement* e, const GmVector* coord, int ip) const
113  {
114  return propertyAc(BETA_ID)->scalarValueAt(e, coord, ip);
115  }
117  virtual double lambdaParameter(const GmElement* e, const GmVector* coord, int ip) const
118  {
119  return propertyAc(LAMBDA_ID)->scalarValueAt(e, coord, ip);
120  }
122  virtual double thetaParameter(const GmElement* e, const GmVector* coord, int ip) const
123  {
124  return propertyAc(THETA_ID)->scalarValueAt(e, coord, ip);
125  }
127  virtual double initialCenterEllipse(const GmElement* e, const GmVector* coord, int ip) const
128  {
129  return propertyAc(INITIALCENTER_ID)->scalarValueAt(e, coord, ip);
130  }
132  virtual double ratioSemiaxes(const GmElement* e, const GmVector* coord, int ip) const
133  {
134  return propertyAc(RATIO_ID)->scalarValueAt(e, coord, ip);
135  }
137  virtual double tensionCutOff(const GmElement* e, const GmVector* coord, int ip) const
138  {
139  return propertyAc(CUTOFF_ID)->scalarValueAt(e, coord, ip);
140  }
142  virtual double maxVolumetricStrain(const GmElement* e, const GmVector* coord, int ip) const
143  {
144  S_TRACE(); assert(e);
145  // Get accessor for the ultimate compressive volumetric strain
146  GmCellAccessor* wAccc = propertyAc(MAXVSTRAIN_ID);
147 
148  if (wAccc == NULL)
149  {
150  return 1.0e10;
151  }
152 
153  return wAccc->scalarValueAt(e, coord, ip);
154  }
155 
157  virtual double volumetricStrainRate(const GmElement* e, const GmVector* coord, int ip) const
158  {
159  S_TRACE(); assert(e);
160  // Get accessor for the volumetric strain rate with respect to the compressive hydrostatic strain
161  GmCellAccessor* dAccc = propertyAc(VSTRAINRATE_ID);
162 
163  if (dAccc == NULL)
164  {
165  return 1.0e-3;
166  }
167 
168  return dAccc->scalarValueAt(e, coord, ip);
169  }
170 
171  // Yield criterion - Conical Surface
172  virtual double yieldCriterion(const GmElement*, const GmVector&, const GmVector*, int, unsigned) const;
173  virtual bool yieldStressGradient(const GmElement*, GmVector&, const GmVector&, const GmVector*, int, unsigned, bool) const;
174  virtual bool yieldGradient_Hardening(const GmElement*, GmVector&, const GmVector&, const GmVector*, int, unsigned) const;
175  virtual bool yieldHessian(const GmElement*, GmMatrix&, const GmVector&, const GmVector*, int, unsigned, bool) const;
176 
177  // Yield criterion - Cap Surface
178  virtual double capSurfaceCriterion(const GmElement*, const GmVector&, const double, const GmVector*, int, unsigned) const;
179  virtual bool capSurfaceGradient(const GmElement*, GmVector&, const GmVector&, const double, const GmVector*, int, unsigned, bool) const;
180  virtual bool capSurfaceGradient_Hardening(const GmElement*, GmVector&, const GmVector&, const double, const GmVector*, int, unsigned) const;
181  virtual bool capSurfaceHessian(const GmElement*, GmMatrix&, const GmVector&, const double, const GmVector*, int, unsigned, bool) const;
182 
183  // Yield criterion - Tension Cut-off
184  virtual double tensionCutOffCriterion(const GmElement*, const GmVector&, const double, const GmVector*, int, unsigned) const;
185  virtual bool tensionCutOffGradient(const GmElement*, GmVector&, const GmVector&, const double, const GmVector*, int, unsigned, bool) const;
186  virtual bool tensionCutOffGradient_Hardening(const GmElement*, GmVector&, const GmVector&, const double, const GmVector*, int, unsigned) const;
187  virtual bool tensionCutOffHessian(const GmElement*, GmMatrix&, const GmVector&, const double, const GmVector*, int, unsigned, bool) const;
188 
189  // Plastic potential - Conical Surface
190  virtual double plasticFPotential(const GmElement*, const GmVector&, const GmVector*, int, unsigned) const;
191  virtual bool flowVector(const GmElement*, GmVector&, const GmVector&, const GmVector*, int, unsigned, bool) const;
192  virtual bool flowVectorStressGradient(const GmElement*, GmMatrix&, const GmVector&, const GmVector*, int, unsigned, bool) const;
193  virtual bool flowVectorDerivative_Hardening(const GmElement*, GmVector&, const GmVector&, const GmVector*, int, unsigned) const;
194 
195  // Plastic potential - Cap Surface
196  virtual double capSurfacePlasticFPotential(const GmElement*, const GmVector&, const double, const GmVector*, int, unsigned) const;
197  virtual bool capSurfacePlasticFGradient(const GmElement*, GmVector&, const GmVector&, const double, const GmVector*, int, unsigned, bool) const;
198  virtual bool capSurfacePlasticFHessian(const GmElement*, GmMatrix&, const GmVector&, const double, const GmVector*, int, unsigned, bool) const;
199  virtual bool capSurfaceFlowVectorDerivative_Hardening(const GmElement*, GmVector&, const GmVector&, const double, const GmVector*, int, unsigned, bool) const;
200 
201  // Plastic potential - Tension Cut-off
202  virtual double tensionCutOffPlasticFPotential(const GmElement*, const GmVector&, const double, const GmVector*, int, unsigned) const;
203  virtual bool tensionCutOffPlasticFGradient(const GmElement*, GmVector&, const GmVector&, const double, const GmVector*, int, unsigned, bool) const;
204  virtual bool tensionCutOffPlasticFHessian(const GmElement*, GmMatrix&, const GmVector&, const double, const GmVector*, int, unsigned, bool) const;
205  virtual bool tensionCutOffFlowVectorDerivative_Hardening(const GmElement*, GmVector&, const GmVector&, const double, const GmVector*, int, unsigned) const;
206 
207  // Hardening Evolution - Conical Surface
208  virtual bool hardeningEvolution(const GmElement*, GmVector&, const GmVector&, const double, const GmVector*, int, unsigned) const;
209  virtual bool hardeningEvolutionDerivative_Stress(const GmElement*, GmVector&, const GmVector&, const double, const GmVector*, int, unsigned) const;
210  virtual bool hardeningEvolutionDerivative_Hardening(const GmElement*, GmMatrix&, const GmVector&, const double, const GmVector*, int, unsigned) const;
211 
212  // Hardening Evolution - Cap Surface
213  virtual bool capSurfaceHardeningEvolution(const GmElement*, GmVector&, const GmVector&, const double, const GmVector*, int, unsigned) const;
214  virtual bool capSurfaceHardeningEvolutionDerivative_Stress(const GmElement*, GmVector&, const GmVector&, const double, const GmVector*, int, unsigned, bool) const;
215  virtual bool capSurfaceHardeningEvolutionDerivative_Hardening(const GmElement*, GmMatrix&, const GmVector&, const double, const GmVector*, int, unsigned) const;
216 
217  // Hardening Evolution - Tension Cut-off
218  virtual bool tensionCutOffHardeningEvolution(const GmElement*, GmVector&, const GmVector&, const double, const GmVector*, int, unsigned) const;
219  virtual bool tensionCutOffHardeningEvolutionDerivative_Stress(const GmElement*, GmVector&, const GmVector&, const double, const GmVector*, int, unsigned) const;
220  virtual bool tensionCutOffHardeningEvolutionDerivative_Hardening(const GmElement*, GmMatrix&, const GmVector&, const double, const GmVector*, int, unsigned) const;
221 
222  // Complementary Methods
223  virtual bool yieldSurfaces(const GmElement*, GmVector&, GmVector, const double, const GmVector*, int, unsigned) const;
224  virtual bool flowVectors(const GmElement*, GmVector&, int, GmVector, const double, const GmVector*, int, unsigned, bool) const;
225  virtual bool hardeningLaws(const GmElement*, GmVector&, int, GmVector, const double, const GmVector*, int, unsigned) const;
226  virtual bool gradientVectors_Stress(const GmElement*, GmVector&, int, GmVector, const double, const GmVector*, int, unsigned, bool) const;
227  virtual bool gradientVectors_Hardening(const GmElement*, GmVector&, int, GmVector, const double, const GmVector*, int, unsigned) const;
228  virtual bool flowVectorDerivatives_Stress(const GmElement*, GmMatrix&, int, GmVector, const double, const GmVector*, int, unsigned, bool) const;
229  virtual bool flowVectorDerivatives_Hardening(const GmElement*, GmVector&, int, GmVector, const double, const GmVector*, int, unsigned, bool) const;
230  virtual bool hardeningLawDerivatives_Stress(const GmElement*, GmVector&, int, GmVector, const double, const GmVector*, int, unsigned, bool) const;
231  virtual bool hardeningLawDerivatives_Hardening(const GmElement*, GmMatrix&, int, GmVector, const double q, const GmVector*, int, unsigned) const;
232  void updateDLamb(GmVector&, GmVector, GmVector) const;
233  virtual int numberActiveSurface(GmVector, const double) const;
234  void indexActiveSurfaces(GmVector&, GmVector, int, const double) const;
235  void deactivateSurfaces(GmVector&, GmVector, int) const;
236 
237  // Line Search
238  virtual double lineSearch(const GmElement*, GmMatrix, GmVector, GmVector, GmVector, GmVector, const double,
239  const double, const double, GmVector, GmVector, const GmpMechanicPoint*, const GmVector*, unsigned, bool) const;
240  virtual double goldenSectionMethod(const GmElement*, GmMatrix, GmVector, GmVector, GmVector, GmVector, const double,
241  const double, const double, GmVector, GmVector, const GmpMechanicPoint*, const GmVector*, unsigned, bool) const;
242  virtual double quadraticInterpolation(const GmElement*, GmMatrix, GmVector, GmVector, GmVector, GmVector, const double,
243  const double, const double, GmVector, GmVector, const GmpMechanicPoint*, const GmVector*, unsigned, bool) const;
244  virtual double cubicInterpolation(const GmElement*, GmMatrix, GmVector, GmVector, GmVector, GmVector, const double,
245  const double, const double, GmVector, GmVector, const GmpMechanicPoint*, const GmVector*, unsigned, bool) const;
246 
247  // Residual Function
248  virtual double residualFunction(const GmElement*, GmMatrix, const GmpMechanicPoint*, const GmVector*, GmVector,
249  GmVector, GmVector, const double, const double, GmVector, unsigned, bool) const;
250 
251  // Method for the Newton-Krylov algorithm
252  // Apply Givens Rotations
253  virtual bool applyGivensRotation(const GmVector c, const GmVector s, GmVector& vrot, const int k) const;
254  // Vector function to be solved by the Newton-Krylov method
255  virtual bool vectorFunction(const GmElement* e, const GmpMechanicPoint* mp, const GmVector* coord, GmVector& r, const GmMatrix C, const GmVector de,
256  const GmVector s, const GmVector s_old, const double q, const double q_old, const GmVector dg, const double cphi, const double cdg, const double delta, unsigned sc, bool ips) const;
257  // GMRES method
258  virtual bool generalizedMinimumResidual(const GmVector fi, const GmMatrix M_inv, double tol, const int max_iter, const GmElement* e, const GmpMechanicPoint* mp,
259  const GmVector* coord, const GmMatrix C, const GmVector de, GmVector& d, double& res, int& counter, const GmVector s_current, const GmVector s_old,
260  const double q_current, const double q_old, const GmVector dg_current, const double cphi, const double cdg, const double delta, unsigned sc, bool ips) const;
261  // Numerical Jacobian
262  virtual bool numericalJacobian(const GmElement* e, const GmpMechanicPoint* mp, const GmVector* coord, GmMatrix& J, const GmMatrix C, const GmVector de,
263  const GmVector s, const GmVector s_old, const double q, const double q_old, const GmVector dg, const double cphi, const double cdg, const double delta, unsigned sc, bool ips) const;
264 
265 };
266 
267 #endif
Definition: gmpMaterialCapModel.h:37
virtual double betaParameter(const GmElement *e, const GmVector *coord, int ip) const
Returns the material constant beta.
Definition: gmpMaterialCapModel.h:112
virtual double lambdaParameter(const GmElement *e, const GmVector *coord, int ip) const
Returns the material constant lambda.
Definition: gmpMaterialCapModel.h:117
The number of gauss attributes.
Definition: gmpMaterialElastoplastic.h:92
virtual double ratioSemiaxes(const GmElement *e, const GmVector *coord, int ip) const
Returns the ratio between the semiaxes.
Definition: gmpMaterialCapModel.h:132
virtual bool isIsotropic() const
Returns true if the material is isotropic, false otherwise.
Definition: gmpMaterialCapModel.h:104
The number of Gauss attribute ids above.
Definition: gmpMaterialCapModel.h:63
#define S_TRACE()
Definition: gmpMaterialElastoplasticMultiSurface.h:33
virtual double alphaParameter(const GmElement *e, const GmVector *coord, int ip) const
Returns the material constant alpha.
Definition: gmpMaterialCapModel.h:107
The number of property ids above.
Definition: gmpMaterialElastoplastic.h:79
GaussAttributeIds
IDs for material Gauss attributes.
Definition: gmpMaterialCapModel.h:58
Definition: gmpMechanicPoint.h:32
virtual double maxVolumetricStrain(const GmElement *e, const GmVector *coord, int ip) const
Returns the ultimate compressive volumetric strain.
Definition: gmpMaterialCapModel.h:142
Id for retrieving the ratio between the semiaxes.
Definition: gmpMaterialCapModel.h:48
Id for retrieving the material constant theta accessor.
Definition: gmpMaterialCapModel.h:46
Id for retrieving the material constant lambda accessor.
Definition: gmpMaterialCapModel.h:45
Id for retrieving the initial center of the ellipse.
Definition: gmpMaterialCapModel.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)
The number of property ids above.
Definition: gmpMaterialCapModel.h:54
virtual bool cFunctionsReturnAlgorithm(const GmElement *, GmMatrix &, const GmpMechanicPoint *, const GmVector *, unsigned, bool) const
Returns the updated stresses using a Newton-Raphson implicit state-update algorithm based on compleme...
Definition: gmpMaterialElastoplastic.cpp:527
virtual double tensionCutOff(const GmElement *e, const GmVector *coord, int ip) const
Returns the tension cut-off.
Definition: gmpMaterialCapModel.h:137
Id for retrieving the tension cut-off.
Definition: gmpMaterialCapModel.h:49
virtual double initialCenterEllipse(const GmElement *e, const GmVector *coord, int ip) const
Returns the center of the ellipse.
Definition: gmpMaterialCapModel.h:127
Declaration of the GmpMechanicPoint class.
ElementPropertyIds
IDs for material element properties.
Definition: gmpMaterialCapModel.h:41
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 implicitReturnAlgorithm(const GmElement *, GmMatrix &, const GmpMechanicPoint *, const GmVector *, unsigned, bool) const
Returns the updated stresses using closet point return algorithm.
Definition: gmpMaterialElastoplasticMultiSurface.cpp:60
Id for retrieving the material constant beta accessor.
Definition: gmpMaterialCapModel.h:44
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
virtual ~GmpMaterialCapModel()
Virtual destructor.
Definition: gmpMaterialCapModel.h:73
GmpMaterialCapModel(int typeIndex, QString typeName, const GmLogCategory &logger)
Constructor. Gets as parameters the material index and its name.
Definition: gmpMaterialCapModel.h:69
arma::vec GmVector
Id for retrieving the volumetric strain rate with respect to the compressive hydrostatic strain.
Definition: gmpMaterialCapModel.h:51
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: gmpMaterialCapModel.h:76
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: gmpMaterialCapModel.h:157
Id for retrieving the ultimate compressive volumetric strain.
Definition: gmpMaterialCapModel.h:50
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
virtual double thetaParameter(const GmElement *e, const GmVector *coord, int ip) const
Returns the material constant theta.
Definition: gmpMaterialCapModel.h:122
virtual bool newtonKrylovReturnAlgorithm(const GmElement *, GmMatrix &, const GmpMechanicPoint *, const GmVector *, unsigned, bool) const
Returns the updated stresses using a Newton-Krylov return mapping algorithm.
Definition: gmpMaterialElastoplastic.cpp:539