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_isotropic_plasticity.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 & Lucia Barbu
12 // Collaborator: Vicente Mataix Ferrandiz
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 TConstLawIntegratorType>
58 class KRATOS_API(CONSTITUTIVE_LAWS_APPLICATION) GenericSmallStrainIsotropicPlasticity
59  : public std::conditional<TConstLawIntegratorType::VoigtSize == 6, ElasticIsotropic3D, LinearPlaneStrain >::type
60 {
61 public:
64 
66  static constexpr SizeType Dimension = TConstLawIntegratorType::Dimension;
67 
69  static constexpr SizeType VoigtSize = TConstLawIntegratorType::VoigtSize;
70 
73 
76 
79 
82 
84  typedef Node NodeType;
85 
88 
92 
97  {
98  }
99 
103  ConstitutiveLaw::Pointer Clone() const override
104  {
105  return Kratos::make_shared<GenericSmallStrainIsotropicPlasticity<TConstLawIntegratorType>>(*this);
106  }
107 
112  : BaseType(rOther),
113  mPlasticDissipation(rOther.mPlasticDissipation),
114  mThreshold(rOther.mThreshold),
115  mPlasticStrain(rOther.mPlasticStrain)
116  {
117  }
118 
123  {
124  }
125 
129 
133 
138  void CalculateMaterialResponsePK1(ConstitutiveLaw::Parameters &rValues) override;
139 
144  void CalculateMaterialResponsePK2(ConstitutiveLaw::Parameters &rValues) override;
145 
150  void CalculateMaterialResponseKirchhoff(ConstitutiveLaw::Parameters &rValues) override;
151 
156  void CalculateMaterialResponseCauchy(ConstitutiveLaw::Parameters &rValues) override;
157 
164  void InitializeMaterial(
165  const Properties& rMaterialProperties,
166  const GeometryType& rElementGeometry,
167  const Vector& rShapeFunctionsValues
168  ) override;
169 
174  void InitializeMaterialResponsePK1 (ConstitutiveLaw::Parameters& rValues) override;
175 
180  void InitializeMaterialResponsePK2 (ConstitutiveLaw::Parameters& rValues) override;
181 
186  void InitializeMaterialResponseKirchhoff (ConstitutiveLaw::Parameters& rValues) override;
187 
192  void InitializeMaterialResponseCauchy (ConstitutiveLaw::Parameters& rValues) override;
193 
198  void FinalizeMaterialResponsePK1(ConstitutiveLaw::Parameters &rValues) override;
199 
204  void FinalizeMaterialResponsePK2(ConstitutiveLaw::Parameters &rValues) override;
205 
210  void FinalizeMaterialResponseKirchhoff(ConstitutiveLaw::Parameters &rValues) override;
215  void FinalizeMaterialResponseCauchy(ConstitutiveLaw::Parameters &rValues) override;
216 
222  bool Has(const Variable<double> &rThisVariable) override;
223 
229  bool Has(const Variable<Vector> &rThisVariable) override;
230 
236  bool Has(const Variable<Matrix> &rThisVariable) override;
237 
244  void SetValue(
245  const Variable<double> &rThisVariable,
246  const double& rValue,
247  const ProcessInfo& rCurrentProcessInfo
248  ) override;
249 
256  void SetValue(
257  const Variable<Vector> &rThisVariable,
258  const Vector& rValue,
259  const ProcessInfo& rCurrentProcessInfo
260  ) override;
261 
268  double& GetValue(
269  const Variable<double> &rThisVariable,
270  double& rValue
271  ) override;
272 
279  Vector& GetValue(
280  const Variable<Vector> &rThisVariable,
281  Vector& rValue
282  ) override;
283 
290  Matrix& GetValue(
291  const Variable<Matrix>& rThisVariable,
292  Matrix& rValue
293  ) override;
294 
299  {
300  return false;
301  }
302 
307  {
308  return true;
309  }
310 
318  double& CalculateValue(
319  ConstitutiveLaw::Parameters& rParameterValues,
320  const Variable<double>& rThisVariable,
321  double& rValue) override;
322 
330  Vector& CalculateValue(
331  ConstitutiveLaw::Parameters& rParameterValues,
332  const Variable<Vector>& rThisVariable,
333  Vector& rValue
334  ) override;
335 
343  Matrix& CalculateValue(
344  ConstitutiveLaw::Parameters& rParameterValues,
345  const Variable<Matrix>& rThisVariable,
346  Matrix& rValue
347  ) override;
348 
358  int Check(
359  const Properties& rMaterialProperties,
360  const GeometryType& rElementGeometry,
361  const ProcessInfo& rCurrentProcessInfo
362  ) const override;
363 
367 
371 
375 
379 
381 
382 protected:
385 
389 
393 
394  double& GetThreshold() { return mThreshold; }
395  double& GetPlasticDissipation() { return mPlasticDissipation; }
396  Vector& GetPlasticStrain() { return mPlasticStrain; }
397 
398  void SetThreshold(const double Threshold) { mThreshold = Threshold; }
399  void SetPlasticDissipation(const double PlasticDissipation) { mPlasticDissipation = PlasticDissipation; }
400  void SetPlasticStrain(const array_1d<double, VoigtSize>& rPlasticStrain) { mPlasticStrain = rPlasticStrain; }
401 
405 
409 
413 
417 
419  private:
422 
426 
427  // Converged values
428  double mPlasticDissipation = 0.0;
429  double mThreshold = 0.0;
430  Vector mPlasticStrain = ZeroVector(VoigtSize);
431 
435 
439 
444  void CalculateTangentTensor(ConstitutiveLaw::Parameters &rValues);
445 
449 
453 
457 
458  // Serialization
459 
460  friend class Serializer;
461 
462  void save(Serializer &rSerializer) const override
463  {
465  rSerializer.save("PlasticDissipation", mPlasticDissipation);
466  rSerializer.save("Threshold", mThreshold);
467  rSerializer.save("PlasticStrain", mPlasticStrain);
468  }
469 
470  void load(Serializer &rSerializer) override
471  {
473  rSerializer.load("PlasticDissipation", mPlasticDissipation);
474  rSerializer.load("Threshold", mThreshold);
475  rSerializer.load("PlasticStrain", mPlasticStrain);
476  }
477 
479 
480 }; // Class GenericYieldSurface
481 
482 } // namespace Kratos
Definition: constitutive_law.h:47
This class is the base class which define all the constitutive laws for plasticity in small deformati...
Definition: generic_small_strain_isotropic_plasticity.h:60
Geometry< NodeType > GeometryType
The geometry definition.
Definition: generic_small_strain_isotropic_plasticity.h:87
void SetPlasticStrain(const array_1d< double, VoigtSize > &rPlasticStrain)
Definition: generic_small_strain_isotropic_plasticity.h:400
GenericSmallStrainIsotropicPlasticity(const GenericSmallStrainIsotropicPlasticity &rOther)
Definition: generic_small_strain_isotropic_plasticity.h:111
Node NodeType
The node definition.
Definition: generic_small_strain_isotropic_plasticity.h:84
std::conditional< VoigtSize==6, ElasticIsotropic3D, LinearPlaneStrain >::type BaseType
Definition of the base class.
Definition: generic_small_strain_isotropic_plasticity.h:72
BoundedMatrix< double, Dimension, Dimension > BoundedMatrixType
The definition of the bounded matrix type.
Definition: generic_small_strain_isotropic_plasticity.h:78
ConstitutiveLaw::Pointer Clone() const override
Definition: generic_small_strain_isotropic_plasticity.h:103
bool RequiresInitializeMaterialResponse() override
If the CL requires to initialize the material response, called by the element in InitializeSolutionSt...
Definition: generic_small_strain_isotropic_plasticity.h:298
GenericSmallStrainIsotropicPlasticity()
Definition: generic_small_strain_isotropic_plasticity.h:96
KRATOS_CLASS_POINTER_DEFINITION(GenericSmallStrainIsotropicPlasticity)
Counted pointer of GenericSmallStrainIsotropicPlasticity.
void SetPlasticDissipation(const double PlasticDissipation)
Definition: generic_small_strain_isotropic_plasticity.h:399
array_1d< double, VoigtSize > BoundedArrayType
The definition of the Voigt array type.
Definition: generic_small_strain_isotropic_plasticity.h:75
double & GetThreshold()
Definition: generic_small_strain_isotropic_plasticity.h:394
bool RequiresFinalizeMaterialResponse() override
If the CL requires to initialize the material response, called by the element in InitializeSolutionSt...
Definition: generic_small_strain_isotropic_plasticity.h:306
~GenericSmallStrainIsotropicPlasticity() override
Definition: generic_small_strain_isotropic_plasticity.h:122
Vector & GetPlasticStrain()
Definition: generic_small_strain_isotropic_plasticity.h:396
void SetThreshold(const double Threshold)
Definition: generic_small_strain_isotropic_plasticity.h:398
double & GetPlasticDissipation()
Definition: generic_small_strain_isotropic_plasticity.h:395
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