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.
plane_stress_d_plus_d_minus_damage_masonry_2d.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: Philip Kalkbrenner
12 // Massimo Petracca
13 // Alejandro Cornejo
14 //
15 //
16 #pragma once
17 
18 // System includes
19 
20 // External includes
21 
22 // Project includes
24 
25 namespace Kratos
26 {
31 class KRATOS_API(CONSTITUTIVE_LAWS_APPLICATION) DamageDPlusDMinusMasonry2DLaw
32  : public ConstitutiveLaw
33 {
34 public:
35 
37 
40 
42  static constexpr SizeType Dimension = 2;
43 
45  static constexpr SizeType VoigtSize = 3;
46 
48  static constexpr double tolerance = std::numeric_limits<double>::epsilon();
49 
51 
56 
61  {
62  }
63 
67  ConstitutiveLaw::Pointer Clone() const override;
68 
73  {
74  return Dimension;
75  };
76 
80  SizeType GetStrainSize() const override
81  {
82  return VoigtSize;
83  };
84 
86 
87  // Elastic Properties
88  double YoungModulus; double PoissonRatio;
90 
91  // Tension Damage Properties
92  double YieldStressTension; double FractureEnergyTension;
93 
94  // Compression Damage Properties
95  double DamageOnsetStressCompression; double YieldStressCompression;
96  double ResidualStressCompression; double YieldStrainCompression;
97  double BezierControllerC1; double BezierControllerC2;
98  double BezierControllerC3; double FractureEnergyCompression;
99  double BiaxialCompressionMultiplier; double ShearCompressionReductor;
100 
101  // Effective Stress Data
104  Matrix ProjectionTensorTension; Matrix ProjectionTensorCompression;
105 
106  // Misc
107  double CharacteristicLength; double DeltaTime;
109 
110  };
111 
117  bool Has(const Variable<double>& rThisVariable) override;
118 
124  bool Has(const Variable<Vector>& rThisVariable) override;
125 
131  bool Has(const Variable<Matrix>& rThisVariable) override;
132 
139  bool Has(const Variable<array_1d<double, 3 > >& rThisVariable) override;
140 
147  bool Has(const Variable<array_1d<double, 6 > >& rThisVariable) override;
148 
155  double& GetValue(
156  const Variable<double>& rThisVariable,
157  double& rValue) override;
158 
165  Vector& GetValue(
166  const Variable<Vector>& rThisVariable,
167  Vector& rValue) override;
168 
174  Matrix& GetValue(
175  const Variable<Matrix>& rThisVariable,
176  Matrix& rValue) override;
177 
185  const Variable<array_1d<double, 3 > >& rVariable,
186  array_1d<double, 3 > & rValue) override;
187 
195  const Variable<array_1d<double, 6 > >& rVariable,
196  array_1d<double, 6 > & rValue) override;
197 
204  void SetValue(
205  const Variable<double>& rVariable,
206  const double& rValue,
207  const ProcessInfo& rCurrentProcessInfo) override;
208 
215  void SetValue(
216  const Variable<Vector >& rVariable,
217  const Vector& rValue,
218  const ProcessInfo& rCurrentProcessInfo) override;
219 
226  void SetValue(
227  const Variable<Matrix >& rVariable,
228  const Matrix& rValue,
229  const ProcessInfo& rCurrentProcessInfo) override;
230 
237  void SetValue(
238  const Variable<array_1d<double, 3 > >& rVariable,
239  const array_1d<double, 3 > & rValue,
240  const ProcessInfo& rCurrentProcessInfo) override;
241 
248  void SetValue(
249  const Variable<array_1d<double, 6 > >& rVariable,
250  const array_1d<double, 6 > & rValue,
251  const ProcessInfo& rCurrentProcessInfo) override;
252 
261  bool ValidateInput(const Properties& rMaterialProperties) override;
262 
267  StrainMeasure GetStrainMeasure() override;
268 
273  StressMeasure GetStressMeasure() override;
274 
280  bool IsIncremental() override;
281 
290  void InitializeMaterial(
291  const Properties& rMaterialProperties,
292  const GeometryType& rElementGeometry,
293  const Vector& rShapeFunctionsValues) override;
294 
299  void InitializeMaterialResponsePK2 (
300  ConstitutiveLaw::Parameters& rValues) override;
301 
311  const Properties& rMaterialProperties,
312  const GeometryType& rElementGeometry,
313  const Vector& rShapeFunctionsValues,
314  const ProcessInfo& rCurrentProcessInfo) override;
315 
324  void FinalizeSolutionStep(
325  const Properties& rMaterialProperties,
326  const GeometryType& rElementGeometry,
327  const Vector& rShapeFunctionsValues,
328  const ProcessInfo& rCurrentProcessInfo) override;
329 
334  void CalculateMaterialResponsePK1 (ConstitutiveLaw::Parameters& rValues) override;
335 
340  void CalculateMaterialResponsePK2 (ConstitutiveLaw::Parameters& rValues) override;
341 
346  void CalculateMaterialResponseKirchhoff(ConstitutiveLaw::Parameters& rValues) override;
347 
352  void CalculateMaterialResponseCauchy (ConstitutiveLaw::Parameters& rValues) override;
353 
358  void FinalizeMaterialResponsePK1 (ConstitutiveLaw::Parameters& rValues) override;
359 
364  void FinalizeMaterialResponsePK2 (ConstitutiveLaw::Parameters& rValues) override;
365 
370  void FinalizeMaterialResponseKirchhoff (ConstitutiveLaw::Parameters& rValues) override;
371 
376  void FinalizeMaterialResponseCauchy (ConstitutiveLaw::Parameters& rValues) override;
377 
386  void ResetMaterial(
387  const Properties& rMaterialProperties,
388  const GeometryType& rElementGeometry,
389  const Vector& rShapeFunctionsValues) override;
390 
395  void GetLawFeatures(Features& rFeatures) override;
396 
406  int Check(
407  const Properties& rMaterialProperties,
408  const GeometryType& rElementGeometry,
409  const ProcessInfo& rCurrentProcessInfo) const override;
410 
411 
412  void CalculateMaterialResponse(const Vector& StrainVector,
413  const Matrix& DeformationGradient,
414  Vector& StressVector,
415  Matrix& AlgorithmicTangent,
416  const ProcessInfo& rCurrentProcessInfo,
417  const Properties& rMaterialProperties,
418  const GeometryType& rElementGeometry,
419  const Vector& rShapeFunctionsValues,
420  bool CalculateStresses = true,
421  int CalculateTangent = true,
422  bool SaveInternalVariables = true);
423 
424 protected:
425 
428 
429  // Initialization
430  bool InitializeDamageLaw = false;
431 
432  // Tension & Compression Thresholds
433 
434  // for IMPLEX_Integration:
435  double PreviousThresholdTension = 0.0; double PreviousThresholdCompression = 0.0; // at time step n - 1:
436  // end for IMPLEX_Integration
437 
438  double CurrentThresholdTension = 0.0; double CurrentThresholdCompression = 0.0; // at time step n:
439  double ThresholdTension = 0.0; double ThresholdCompression = 0.0; // at time step n + 1:
440 
441  // Damage Parameters & Uniaxial Stresses
442  double DamageParameterTension = 0.0; double DamageParameterCompression = 0.0;
443  double UniaxialStressTension = 0.0; double UniaxialStressCompression = 0.0;
444 
445  // Misc
446  double InitialCharacteristicLength = 0.0;
447 
448  // for IMPLEX Integration
449  double CurrentDeltaTime = 0.0; // at time step n + 1
450  double PreviousDeltaTime = 0.0; // at time step n
451 
452  // temporary internal variables saved for the implicit step in the finalize solution step of the IMPLEX
453  double TemporaryImplicitThresholdTension = 0.0;
454  double TemporaryImplicitThresholdTCompression = 0.0;
455  // end for IMPLEX Integration
456 
458 
459 
462 
470  void InitializeCalculationData(
471  const Properties& props,
472  const GeometryType& geom,
473  const ProcessInfo& pinfo,
475 
480  void CalculateElasticityMatrix(
482 
487  void TensionCompressionSplit(
489 
494  void ConstructProjectionTensors(
496 
502  void CalculateEquivalentStressTension(
504  double& UniaxialStressTension);
505 
511  void CalculateEquivalentStressCompression(
513  double& UniaxialStressCompression);
514 
522  void CalculateDamageTension(
524  double internal_variable,
525  double& rDamageTension);
526 
574  void CalculateDamageCompression(
576  double internal_variable,
577  double& rDamage);
578 
586  void ComputeBezierEnergy(
587  double& rBezierEnergy,
588  double& rBezierEnergy1,
589  double s_p, double s_k, double s_r,
590  double e_p, double e_j, double e_k, double e_r, double e_u);
591 
597  double EvaluateBezierArea(
598  double x1, double x2, double x3,
599  double y1, double y2, double y3);
600 
607  void ApplyBezierStretcherToStrains(
608  double stretcher, double e_p, double& e_j,
609  double& e_k, double& e_r, double& e_u);
610 
618  void EvaluateBezierCurve(
619  double& rDamageParameter, double xi,
620  double x1, double x2, double x3,
621  double y1, double y2, double y3);
622 
628  void ComputeCharacteristicLength(
629  const GeometryType& geom,
630  double& rCharacteristicLength);
631 
638  void CalculateMaterialResponseInternal(
639  const Vector& strain_vector,
640  Vector& stress_vector,
642  const Properties props);
643 
649  void CheckDamageLoadingUnloading(
650  bool& is_damaging_tension,
651  bool& is_damaging_compression);
652 
660  void CalculateTangentTensor(
661  Parameters& rValues,
662  Vector strain_vector,
663  Vector stress_vector,
665  const Properties& props);
666 
672  void CalculateSecantTensor(
673  Parameters& rValues,
675 
676 
678 
682 
683 
684 private:
685 
689 
693 
697 
701 
705 
708 
709  friend class Serializer;
710 
711  void save(Serializer& rSerializer) const override
712  {
714  }
715 
716  void load(Serializer& rSerializer) override
717  {
719  }
720 
722 
723 }; // Class DamageDPlusDMinusMasonry2DLaw
724 
725 } // namespace Kratos
Definition: constitutive_law.h:47
StrainMeasure
Definition: constitutive_law.h:52
StressMeasure
Definition: constitutive_law.h:69
Definition: plane_stress_d_plus_d_minus_damage_masonry_2d.h:33
~DamageDPlusDMinusMasonry2DLaw() override
Definition: plane_stress_d_plus_d_minus_damage_masonry_2d.h:60
SizeType WorkingSpaceDimension() override
returns the working space dimension of the current constitutive law
Definition: plane_stress_d_plus_d_minus_damage_masonry_2d.h:72
SizeType GetStrainSize() const override
returns the size of the strain vector of the current constitutive law
Definition: plane_stress_d_plus_d_minus_damage_masonry_2d.h:80
KRATOS_CLASS_POINTER_DEFINITION(DamageDPlusDMinusMasonry2DLaw)
Geometry base class.
Definition: geometry.h:71
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
#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
void InitializeSolutionStep(ConstructionUtility &rThisUtil, std::string ThermalSubModelPartName, std::string MechanicalSubModelPartName, std::string HeatFluxSubModelPartName, std::string HydraulicPressureSubModelPartName, bool thermal_conditions, bool mechanical_conditions, int phase)
Definition: add_custom_utilities_to_python.cpp:45
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
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
x2
Definition: generate_frictional_mortar_condition.py:122
x1
Definition: generate_frictional_mortar_condition.py:121
data
Definition: mesh_to_mdpa_converter.py:59
def load(f)
Definition: ode_solve.py:307
Definition: constitutive_law.h:137
Definition: constitutive_law.h:189
Definition: plane_stress_d_plus_d_minus_damage_masonry_2d.h:85
double BiaxialCompressionMultiplier
Definition: plane_stress_d_plus_d_minus_damage_masonry_2d.h:99
Matrix ProjectionTensorCompression
Definition: plane_stress_d_plus_d_minus_damage_masonry_2d.h:104
double ResidualStressCompression
Definition: plane_stress_d_plus_d_minus_damage_masonry_2d.h:96
double BezierControllerC3
Definition: plane_stress_d_plus_d_minus_damage_masonry_2d.h:98
double DamageOnsetStressCompression
Definition: plane_stress_d_plus_d_minus_damage_masonry_2d.h:95
array_1d< double, 3 > EffectiveStressVector
Definition: plane_stress_d_plus_d_minus_damage_masonry_2d.h:102
double FractureEnergyTension
Definition: plane_stress_d_plus_d_minus_damage_masonry_2d.h:92
array_1d< double, 3 > EffectiveCompressionStressVector
Definition: plane_stress_d_plus_d_minus_damage_masonry_2d.h:103
double PoissonRatio
Definition: plane_stress_d_plus_d_minus_damage_masonry_2d.h:88
Matrix ElasticityMatrix
Definition: plane_stress_d_plus_d_minus_damage_masonry_2d.h:89
double CharacteristicLength
Definition: plane_stress_d_plus_d_minus_damage_masonry_2d.h:107
int TensionYieldModel
Definition: plane_stress_d_plus_d_minus_damage_masonry_2d.h:108
double BezierControllerC1
Definition: plane_stress_d_plus_d_minus_damage_masonry_2d.h:97