KratosMultiphysics
KRATOS Multiphysics (Kratos) is a framework for building parallel, multi-disciplinary simulation software, aiming at modularity, extensibility, and high performance. Kratos is written in C++, and counts with an extensive Python interface.
generic_small_strain_plastic_damage_model.h
Go to the documentation of this file.
1 // KRATOS ___ _ _ _ _ _ __ _
2 // / __\___ _ __ ___| |_(_) |_ _ _| |_(_)_ _____ / / __ ___ _____ /_\ _ __ _ __
3 // / / / _ \| '_ \/ __| __| | __| | | | __| \ \ / / _ \/ / / _` \ \ /\ / / __| //_\\| '_ \| '_ |
4 // / /__| (_) | | | \__ \ |_| | |_| |_| | |_| |\ V / __/ /__| (_| |\ V V /\__ \/ _ \ |_) | |_) |
5 // \____/\___/|_| |_|___/\__|_|\__|\__,_|\__|_| \_/ \___\____/\__,_| \_/\_/ |___/\_/ \_/ .__/| .__/
6 // |_| |_|
7 //
8 // License: BSD License
9 // license: structural_mechanics_application/license.txt
10 //
11 // Main authors: Alejandro Cornejo & Sergio Jimenez
12 //
13 //
14 
15 #pragma once
16 
17 // System includes
18 
19 // External includes
20 
21 // Project includes
22 #include "custom_constitutive/elastic_isotropic_3d.h"
23 #include "custom_constitutive/linear_plane_strain.h"
24 
25 namespace Kratos
26 {
29 
33 
34 // The size type definition
35 typedef std::size_t SizeType;
36 
40 
44 
48 
57 template <class TPlasticityIntegratorType, class TDamageIntegratorType>
58 class KRATOS_API(CONSTITUTIVE_LAWS_APPLICATION) GenericSmallStrainPlasticDamageModel
59  : public std::conditional<TPlasticityIntegratorType::VoigtSize == 6, ElasticIsotropic3D, LinearPlaneStrain >::type
60 {
61 public:
64 
65  // The index type definition
66  typedef std::size_t IndexType;
67 
69  static constexpr SizeType Dimension = TPlasticityIntegratorType::Dimension;
70 
72  static constexpr SizeType VoigtSize = TPlasticityIntegratorType::VoigtSize;
73 
76 
79 
82 
85 
87  typedef Node NodeType;
88 
91 
93  static constexpr double tolerance = std::numeric_limits<double>::epsilon();
94 
99  double DamageIndicator = 0.0;
100  double PlasticityIndicator = 0.0;
104  double DamageIncrement = 0.0;
105  double PlasticConsistencyIncrement = 0.0;
106  double UniaxialStressPlasticity = 0.0;
107  double UniaxialStressDamage = 0.0;
108  double HardeningParameterDamage = 0.0;
109  double DamageDissipationIncrement = 0.0;
111  double CharacteristicLength = 0.0;
112  double Damage = 0.0;
113  double PlasticDissipation = 0.0;
114  double DamageDissipation = 0.0;
115  double DamageThreshold = 0.0;
116  double PlasticityThreshold = 0.0;
117  double PlasticDenominator = 0.0;
118  double UndamagedFreeEnergy = 0.0;
119  };
123 
128  {
129  }
130 
134  ConstitutiveLaw::Pointer Clone() const override
135  {
136  return Kratos::make_shared<GenericSmallStrainPlasticDamageModel<TPlasticityIntegratorType, TDamageIntegratorType>>(*this);
137  }
138 
143  : BaseType(rOther),
144  mPlasticDissipation(rOther.mPlasticDissipation),
145  mThresholdPlasticity(rOther.mThresholdPlasticity),
146  mPlasticStrain(rOther.mPlasticStrain),
147  mThresholdDamage(rOther.mThresholdDamage),
148  mDamage(rOther.mDamage)
149  {
150  }
151 
156  {
157  }
158 
162 
166 
171  void CalculateMaterialResponsePK1(ConstitutiveLaw::Parameters &rValues) override;
172 
177  void CalculateMaterialResponsePK2(ConstitutiveLaw::Parameters &rValues) override;
178 
183  void CalculateMaterialResponseKirchhoff(ConstitutiveLaw::Parameters &rValues) override;
184 
189  void CalculateMaterialResponseCauchy(ConstitutiveLaw::Parameters &rValues) override;
190 
197  void InitializeMaterial(
198  const Properties& rMaterialProperties,
199  const GeometryType& rElementGeometry,
200  const Vector& rShapeFunctionsValues
201  ) override;
202 
207  void InitializeMaterialResponsePK1 (ConstitutiveLaw::Parameters& rValues) override;
208 
213  void InitializeMaterialResponsePK2 (ConstitutiveLaw::Parameters& rValues) override;
214 
219  void InitializeMaterialResponseKirchhoff (ConstitutiveLaw::Parameters& rValues) override;
220 
225  void InitializeMaterialResponseCauchy (ConstitutiveLaw::Parameters& rValues) override;
226 
231  void FinalizeMaterialResponsePK1(ConstitutiveLaw::Parameters &rValues) override;
232 
237  void FinalizeMaterialResponsePK2(ConstitutiveLaw::Parameters &rValues) override;
238 
243  void FinalizeMaterialResponseKirchhoff(ConstitutiveLaw::Parameters &rValues) override;
248  void FinalizeMaterialResponseCauchy(ConstitutiveLaw::Parameters &rValues) override;
249 
255  bool Has(const Variable<double> &rThisVariable) override;
256 
262  bool Has(const Variable<Vector> &rThisVariable) override;
263 
269  bool Has(const Variable<Matrix> &rThisVariable) override;
270 
277  void SetValue(
278  const Variable<double> &rThisVariable,
279  const double& rValue,
280  const ProcessInfo& rCurrentProcessInfo
281  ) override;
282 
289  void SetValue(
290  const Variable<Vector> &rThisVariable,
291  const Vector& rValue,
292  const ProcessInfo& rCurrentProcessInfo
293  ) override;
294 
301  double& GetValue(
302  const Variable<double> &rThisVariable,
303  double& rValue
304  ) override;
305 
312  Vector& GetValue(
313  const Variable<Vector> &rThisVariable,
314  Vector& rValue
315  ) override;
316 
323  Matrix& GetValue(
324  const Variable<Matrix>& rThisVariable,
325  Matrix& rValue
326  ) override;
327 
332  {
333  return false;
334  }
335 
340  {
341  return true;
342  }
343 
351  double& CalculateValue(
352  ConstitutiveLaw::Parameters& rParameterValues,
353  const Variable<double>& rThisVariable,
354  double& rValue) override;
355 
363  Vector& CalculateValue(
364  ConstitutiveLaw::Parameters& rParameterValues,
365  const Variable<Vector>& rThisVariable,
366  Vector& rValue
367  ) override;
368 
376  Matrix& CalculateValue(
377  ConstitutiveLaw::Parameters& rParameterValues,
378  const Variable<Matrix>& rThisVariable,
379  Matrix& rValue
380  ) override;
381 
391  int Check(
392  const Properties& rMaterialProperties,
393  const GeometryType& rElementGeometry,
394  const ProcessInfo& rCurrentProcessInfo
395  ) const override;
396 
416  double CalculateDamageParameters(
417  PlasticDamageParameters& rParameters,
418  const Matrix& rElasticMatrix,
419  ConstitutiveLaw::Parameters& rValues);
420 
429  void CalculateIndicatorsFactors(
430  const array_1d<double, 6>& rPredictiveStressVector,
431  double& rTensileIndicatorFactor,
432  double& rCompressionIndicatorFactor,
433  double& rSumPrincipalStresses,
434  array_1d<double, 3>& rPrincipalStresses);
435 
440  void CheckInternalVariable(
441  double& rInternalVariable);
442 
461  void CalculateIncrementsPlasticDamageCase(
462  PlasticDamageParameters& rParameters,
463  const Matrix& rElasticMatrix);
464 
482  double CalculatePlasticParameters(
483  PlasticDamageParameters& rParameters,
484  const Matrix& rConstitutiveMatrix,
485  ConstitutiveLaw::Parameters& rValues);
486 
496  void CalculatePlasticDenominator(
497  const array_1d<double, VoigtSize>& rFFlux,
498  const array_1d<double, VoigtSize>& rGFlux,
499  const Matrix& rConstitutiveMatrix,
500  double& rHardeningParameter,
501  const double Damage,
502  double& rPlasticDenominator);
503 
507 
511 
515 
519 
521 
522 protected:
525 
529 
533 
534  double& GetThresholdPlasticity() { return mThresholdPlasticity; }
535  double& GetPlasticDissipation() { return mPlasticDissipation; }
536  Vector& GetPlasticStrain() { return mPlasticStrain; }
537 
538  void SetThresholdPlasticity(const double ThresholdPlasticity) { mThresholdPlasticity = ThresholdPlasticity; }
539  void SetPlasticDissipation(const double PlasticDissipation) { mPlasticDissipation = PlasticDissipation; }
540  void SetPlasticStrain(const array_1d<double, VoigtSize>& rPlasticStrain) { mPlasticStrain = rPlasticStrain; }
541 
545 
549 
553 
557 
559  private:
562 
566 
567  // Converged values
568  double mPlasticDissipation = 0.0;
569  double mThresholdPlasticity = 0.0;
570  Vector mPlasticStrain = ZeroVector(VoigtSize);
571  double mThresholdDamage = 0.0;
572  double mDamage = 0.0;
573  double mDamageDissipation = 0.0;
574 
575  // only for printing
576  double mUniaxialStress = 0.0;
577 
581 
585 
590  void CalculateTangentTensor(ConstitutiveLaw::Parameters &rValues);
591 
595 
599 
603 
604  // Serialization
605 
606  friend class Serializer;
607 
608  void save(Serializer &rSerializer) const override
609  {
611  rSerializer.save("PlasticDissipation", mPlasticDissipation);
612  rSerializer.save("ThresholdPlasticity", mThresholdPlasticity);
613  rSerializer.save("PlasticStrain", mPlasticStrain);
614  rSerializer.save("ThresholdDamage", mThresholdDamage);
615  rSerializer.save("Damage", mDamage);
616  rSerializer.save("DamageDissipation", mDamageDissipation);
617  }
618 
619  void load(Serializer &rSerializer) override
620  {
622  rSerializer.load("PlasticDissipation", mPlasticDissipation);
623  rSerializer.load("ThresholdPlasticity", mThresholdPlasticity);
624  rSerializer.load("PlasticStrain", mPlasticStrain);
625  rSerializer.load("ThresholdDamage", mThresholdDamage);
626  rSerializer.load("Damage", mDamage);
627  rSerializer.load("DamageDissipation", mDamageDissipation);
628  }
629 
631 
632 }; // Class GenericYieldSurface
633 
634 } // namespace Kratos
Definition: constitutive_law.h:47
This class is the base class which define the Plastic Damage model developed by Luccioni B....
Definition: generic_small_strain_plastic_damage_model.h:60
void SetPlasticDissipation(const double PlasticDissipation)
Definition: generic_small_strain_plastic_damage_model.h:539
void SetPlasticStrain(const array_1d< double, VoigtSize > &rPlasticStrain)
Definition: generic_small_strain_plastic_damage_model.h:540
KRATOS_CLASS_POINTER_DEFINITION(GenericSmallStrainPlasticDamageModel)
Counted pointer of GenericSmallStrainPlasticDamageModel.
double & GetPlasticDissipation()
Definition: generic_small_strain_plastic_damage_model.h:535
Node NodeType
The node definition.
Definition: generic_small_strain_plastic_damage_model.h:87
array_1d< double, VoigtSize > BoundedArrayType
The definition of the Voigt array type.
Definition: generic_small_strain_plastic_damage_model.h:78
~GenericSmallStrainPlasticDamageModel() override
Definition: generic_small_strain_plastic_damage_model.h:155
void SetThresholdPlasticity(const double ThresholdPlasticity)
Definition: generic_small_strain_plastic_damage_model.h:538
Vector & GetPlasticStrain()
Definition: generic_small_strain_plastic_damage_model.h:536
GenericSmallStrainPlasticDamageModel()
Definition: generic_small_strain_plastic_damage_model.h:127
bool RequiresFinalizeMaterialResponse() override
If the CL requires to initialize the material response, called by the element in InitializeSolutionSt...
Definition: generic_small_strain_plastic_damage_model.h:339
std::size_t IndexType
Definition: generic_small_strain_plastic_damage_model.h:66
Geometry< NodeType > GeometryType
The geometry definition.
Definition: generic_small_strain_plastic_damage_model.h:90
ConstitutiveLaw::Pointer Clone() const override
Definition: generic_small_strain_plastic_damage_model.h:134
std::conditional< VoigtSize==6, ElasticIsotropic3D, LinearPlaneStrain >::type BaseType
Definition of the base class.
Definition: generic_small_strain_plastic_damage_model.h:75
double & GetThresholdPlasticity()
Definition: generic_small_strain_plastic_damage_model.h:534
bool RequiresInitializeMaterialResponse() override
If the CL requires to initialize the material response, called by the element in InitializeSolutionSt...
Definition: generic_small_strain_plastic_damage_model.h:331
BoundedMatrix< double, Dimension, Dimension > BoundedMatrixType
The definition of the bounded matrix type.
Definition: generic_small_strain_plastic_damage_model.h:81
GenericSmallStrainPlasticDamageModel(const GenericSmallStrainPlasticDamageModel &rOther)
Definition: generic_small_strain_plastic_damage_model.h:142
Geometry base class.
Definition: geometry.h:71
Definition: amatrix_interface.h:41
This class defines the node.
Definition: node.h:65
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
Properties encapsulates data shared by different Elements or Conditions. It can store any type of dat...
Definition: properties.h:69
The serialization consists in storing the state of an object into a storage format like data file or ...
Definition: serializer.h:123
void load(std::string const &rTag, TDataType &rObject)
Definition: serializer.h:207
void save(std::string const &rTag, std::array< TDataType, TDataSize > const &rObject)
Definition: serializer.h:545
#define KRATOS_SERIALIZE_SAVE_BASE_CLASS(Serializer, BaseType)
Definition: define.h:812
#define KRATOS_SERIALIZE_LOAD_BASE_CLASS(Serializer, BaseType)
Definition: define.h:815
bool Has(const std::string &ModelerName)
Checks if the modeler is registered.
Definition: modeler_factory.cpp:24
Parameters GetValue(Parameters &rParameters, const std::string &rEntry)
Definition: add_kratos_parameters_to_python.cpp:53
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
KratosZeroVector< double > ZeroVector
Definition: amatrix_interface.h:561
std::size_t SizeType
The definition of the size type.
Definition: mortar_classes.h:43
def SetValue(entity, variable, value)
Definition: coupling_interface_data.py:256
type
Definition: generate_gid_list_file.py:35
def load(f)
Definition: ode_solve.py:307
Definition: constitutive_law.h:189
Definition: generic_small_strain_plastic_damage_model.h:95
array_1d< double, VoigtSize > PlasticityGFLux
Definition: generic_small_strain_plastic_damage_model.h:97
array_1d< double, VoigtSize > PlasticStrain
Definition: generic_small_strain_plastic_damage_model.h:101
array_1d< double, VoigtSize > StrainVector
Definition: generic_small_strain_plastic_damage_model.h:102
array_1d< double, VoigtSize > DamageYieldFLux
Definition: generic_small_strain_plastic_damage_model.h:98
array_1d< double, VoigtSize > PlasticStrainIncrement
Definition: generic_small_strain_plastic_damage_model.h:110
array_1d< double, VoigtSize > StressVector
Definition: generic_small_strain_plastic_damage_model.h:103
array_1d< double, VoigtSize > PlasticityFFLux
Definition: generic_small_strain_plastic_damage_model.h:96