MechanicalFemPhysics
The GeMA Mechanical FEM Physics Plugin
gmpMaterialVonMises.h
Go to the documentation of this file.
1 /************************************************************************
2 **
3 ** Copyright (C) 2014 by Nilthson NoreƱa
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_VONMISES_H_
28 #define _GEMA_PLUGIN_MECHANICALMATERIAL_VONMISES_H_
29 
32 #include "gmpMechanicPoint.h"
33 #include "gmpFemPhysics.h"
34 #include "gmMathUtils.h"
35 
36 class GMP_MECHANICAL_PHYSICS_API_EXPORT GmpMaterialVonMises : public GmpMaterialElastoplastic
37 {
38 public:
39 
41  GmpMaterialVonMises(int typeIndex, QString typeName, const GmLogCategory& logger)
42  : GmpMaterialElastoplastic(typeIndex, typeName, logger) {}
43 
45  virtual ~GmpMaterialVonMises() {}
46 
48  static GmpFemPhysicsCommonMaterial* instance(GmSimulationData* simulation, int typeIndex, QString typeName, const GmLogCategory& logger)
49  {
50  Q_UNUSED(simulation);
51  return new GmpMaterialVonMises(typeIndex, typeName, logger);
52  }
53 
54  // Returns a map with material associated properties.
55  virtual const QVariantMap* materialMetaDataMap();
56 
58  virtual bool isIsotropic() const { return true; }
60  virtual double initialYieldStress(const GmElement* e, const GmVector* coord, int ip) const
61  {
62  return propertyAc(IYST_ID)->scalarValueAt(e, coord, ip);
63  }
65  virtual double plasticModulus(const GmElement* e, const GmVector* coord, int ip) const
66  {
67  GmCellAccessor* pmAc = propertyAc(PMD_ID);
68  if (pmAc == NULL)
69  {
70  return 0.0;
71  }
72  return propertyAc(PMD_ID)->scalarValueAt(e, coord, ip);
73  }
74  // Elastoplastic Model Equations
75  virtual double vonMisesCriterion(const GmElement*, const GmVector, const double, const GmVector*, int, unsigned, bool) const;
76  virtual bool vonMisesStressGradient(const GmElement*, GmVector&, const GmVector, const double, const GmVector*, int, unsigned, bool) const;
77  virtual bool vonMisesHardeningGradient(const GmElement*, GmVector&, const GmVector, const double, const GmVector*, int, unsigned, bool) const;
78  virtual bool vonMisesflowVectorStressGradient(const GmElement*, GmMatrix&, const GmVector, const double, const GmVector*, int, unsigned, bool) const;
79  virtual bool vonMisesflowVectorHardeningGradient(const GmElement*, GmVector&, const GmVector, const double, const GmVector*, int, unsigned, bool) const;
80  virtual double hardeningLaw(const GmElement*, const GmVector, const double, const GmVector*, int, unsigned, bool) const;
81  virtual bool hardeningLawStressGradient(const GmElement*, GmVector&, const GmVector, const double, const GmVector*, int, unsigned, bool) const;
82  virtual bool hardeningLawHardeningGradient(const GmElement*, GmVector&, const GmVector, const double, const GmVector*, int, unsigned, bool) const;
83  virtual double yieldStrengthRatio(const GmElement*, const GmVector, const double, const GmVector*, int, unsigned, bool) const;
84 
85  //
86  virtual double yieldCriterion(const GmElement*, const GmVector&, const GmVector*, int, unsigned) const;
87  virtual bool yieldStressGradient(const GmElement*, GmVector&, const GmVector&, const GmVector*, int, unsigned, bool) const;
88  virtual bool yieldHessian(const GmElement*, GmMatrix&, const GmVector&, const GmVector*, int, unsigned, bool) const;
89  virtual double plasticFPotential(const GmElement*, const GmVector&, const GmVector*, int, unsigned) const;
90  virtual bool flowVector(const GmElement*, GmVector&, const GmVector&, const GmVector*, int, unsigned, bool) const;
91  virtual bool flowVectorStressGradient(const GmElement*, GmMatrix&, const GmVector&, const GmVector*, int, unsigned, bool) const;
92  virtual double plasticPotential(const GmElement*, const GmpMechanicPoint*, const GmVector*, int, unsigned) const;
93  virtual double plasticGradient(const GmElement*, const GmpMechanicPoint*, const GmVector*, int, unsigned) const;
94  virtual double plasticHessian(const GmElement*, const GmpMechanicPoint*, const GmVector*, int, unsigned) const;
95  virtual double yieldStrengthRatio(const GmElement*, const GmVector&, const GmVector*, int, unsigned) const { return 0.0; };
96 
97 
98  // Returns the stresses according to the material behavior adopted
99  virtual bool mechanicalConstitutiveModel(const GmElement*, GmMatrix&, const GmpMechanicPoint*, const GmVector*, const GmVector&, unsigned, bool) const;
100 
101 protected:
104  {
107 
108  // --- NO ADDING BELOW THIS LINE
110  };
111 
112  };
113 
114 #endif
virtual double initialYieldStress(const GmElement *e, const GmVector *coord, int ip) const
Returns the material yield stress.
Definition: gmpMaterialVonMises.h:60
virtual bool isIsotropic() const
Returns true if the material is isotropic, false otherwise.
Definition: gmpMaterialVonMises.h:58
virtual ~GmpMaterialVonMises()
Virtual destructor.
Definition: gmpMaterialVonMises.h:45
ElementPropertyIds
IDs for material element properties.
Definition: gmpMaterialElastic.h:38
The number of property ids above.
Definition: gmpMaterialElastoplastic.h:79
Declaration of the gmpMaterial1DPlasticity classes.
Definition: gmpMechanicPoint.h:32
Declaration of the gmpMaterialElastoplastic classes.
virtual double yieldStrengthRatio(const GmElement *e, const GmVector &S, const GmVector *coord, int ip, unsigned sc) const =0
Returns the Yield Strength Ratio (Ysr)
Declaration of the GmpMechanicPoint class.
virtual double yieldStrengthRatio(const GmElement *, const GmVector &, const GmVector *, int, unsigned) const
Returns the Yield Strength Ratio (Ysr)
Definition: gmpMaterialVonMises.h:95
virtual double plasticModulus(const GmElement *e, const GmVector *coord, int ip) const
Returns the material plastic modulus.
Definition: gmpMaterialVonMises.h:65
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: gmpMaterialVonMises.h:48
Definition: gmpMaterialVonMises.h:36
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
arma::vec GmVector
Definition: gmpMaterialElastoplastic.h:35
GmpMaterialVonMises(int typeIndex, QString typeName, const GmLogCategory &logger)
Constructor. Gets as parameters the material index and its name.
Definition: gmpMaterialVonMises.h:41
Id for retrieving the plastic modulus accessor.
Definition: gmpMaterialVonMises.h:106
arma::mat GmMatrix
The number of property ids above.
Definition: gmpMaterialVonMises.h:109
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