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.
serial_parallel_rule_of_mixtures_law.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
12 // Vicente Mataix
13 // Fernando Rastellini
14 // Collaborator: Lucia Barbu
15 //
16 
17 #pragma once
18 
19 // System includes
20 
21 // External includes
22 
23 // Project includes
25 
26 
27 namespace Kratos
28 {
31 
35 
39 
43 
47 
54 class KRATOS_API(CONSTITUTIVE_LAWS_APPLICATION) SerialParallelRuleOfMixturesLaw
55  : public ConstitutiveLaw
56 {
57  public:
60 
62  typedef Node NodeType;
63 
66 
68  static constexpr double machine_tolerance = std::numeric_limits<double>::epsilon();
69 
72 
76 
81  {
82  }
83 
87  SerialParallelRuleOfMixturesLaw(double FiberVolParticipation, const Vector& rParallelDirections)
88  : mFiberVolumetricParticipation(FiberVolParticipation), mParallelDirections(rParallelDirections)
89  {
90  }
91 
95  ConstitutiveLaw::Pointer Clone() const override
96  {
97  return Kratos::make_shared<SerialParallelRuleOfMixturesLaw>(*this);
98  }
99 
100  // Copy constructor
102  : ConstitutiveLaw(rOther), mpMatrixConstitutiveLaw(rOther.mpMatrixConstitutiveLaw), mpFiberConstitutiveLaw(rOther.mpFiberConstitutiveLaw),
103  mFiberVolumetricParticipation(rOther.mFiberVolumetricParticipation), mParallelDirections(rOther.mParallelDirections) ,
104  mPreviousStrainVector(rOther.mPreviousStrainVector) , mPreviousSerialStrainMatrix(rOther.mPreviousSerialStrainMatrix) , mIsPrestressed(rOther.mIsPrestressed)
105  {
106  }
107 
112  {
113  }
114 
118 
122 
128  ConstitutiveLaw::Pointer Create(Kratos::Parameters NewParameters) const override;
129 
134  {
135  return 3;
136  };
137 
141  SizeType GetStrainSize() const override
142  {
143  return 6;
144  };
145 
150  void CalculateMaterialResponsePK1(ConstitutiveLaw::Parameters& rValues) override;
151 
156  void CalculateMaterialResponsePK2(ConstitutiveLaw::Parameters& rValues) override;
157 
162  void CalculateMaterialResponseKirchhoff(ConstitutiveLaw::Parameters& rValues) override;
163 
168  void CalculateMaterialResponseCauchy(ConstitutiveLaw::Parameters& rValues) override;
169 
174  void FinalizeMaterialResponsePK1(ConstitutiveLaw::Parameters& rValues) override;
175 
180  void FinalizeMaterialResponsePK2(ConstitutiveLaw::Parameters& rValues) override;
181 
186  void FinalizeMaterialResponseKirchhoff(ConstitutiveLaw::Parameters& rValues) override;
187 
192  void FinalizeMaterialResponseCauchy(ConstitutiveLaw::Parameters& rValues) override;
193 
200  Vector& GetValue(const Variable<Vector>& rThisVariable, Vector& rValue) override;
201 
208  bool& GetValue(const Variable<bool>& rThisVariable, bool& rValue) override;
209 
216  int& GetValue(const Variable<int>& rThisVariable, int& rValue) override;
217 
224  double& GetValue(const Variable<double>& rThisVariable, double& rValue) override;
225 
231  Matrix& GetValue(const Variable<Matrix>& rThisVariable, Matrix& rValue) override;
232 
239  void SetValue(
240  const Variable<bool>& rThisVariable,
241  const bool& rValue,
242  const ProcessInfo& rCurrentProcessInfo
243  ) override;
244 
251  void SetValue(
252  const Variable<int>& rThisVariable,
253  const int& rValue,
254  const ProcessInfo& rCurrentProcessInfo
255  ) override;
256 
263  void SetValue(
264  const Variable<double>& rThisVariable,
265  const double& rValue,
266  const ProcessInfo& rCurrentProcessInfo
267  ) override;
268 
274  bool Has(const Variable<bool>& rThisVariable) override;
275 
281  bool Has(const Variable<int>& rThisVariable) override;
282 
288  bool Has(const Variable<double>& rThisVariable) override;
289 
295  bool Has(const Variable<Vector>& rThisVariable) override;
296 
302  bool Has(const Variable<Matrix>& rThisVariable) override;
303 
311  bool& CalculateValue(
312  Parameters& rParameterValues,
313  const Variable<bool>& rThisVariable,
314  bool& rValue) override;
315 
323  int& CalculateValue(
324  Parameters& rParameterValues,
325  const Variable<int>& rThisVariable,
326  int& rValue) override;
327 
335  double& CalculateValue(
336  Parameters& rParameterValues,
337  const Variable<double>& rThisVariable,
338  double& rValue) override;
339 
347  Vector& CalculateValue(
348  Parameters& rParameterValues,
349  const Variable<Vector>& rThisVariable,
350  Vector& rValue) override;
351 
359  Matrix& CalculateValue(
360  Parameters& rParameterValues,
361  const Variable<Matrix>& rThisVariable,
362  Matrix& rValue) override;
363 
372  void InitializeMaterial(
373  const Properties& rMaterialProperties,
374  const GeometryType& rElementGeometry,
375  const Vector& rShapeFunctionsValues) override;
376 
377 
387  void IntegrateStrainSerialParallelBehaviour(
388  const Vector& rStrainVector,
389  Vector& FiberStressVector,
390  Vector& MatrixStressVector,
391  const Properties& rMaterialProperties,
393  Vector& rSerialStrainMatrix,
395 
401  void CalculateSerialParallelProjectionMatrices(
402  Matrix& rParallelProjector,
403  Matrix& rSerialProjector);
404 
409  void InitializeMaterialResponsePK2(ConstitutiveLaw::Parameters& rValues) override;
410 
415  void InitializeMaterialResponseCauchy(ConstitutiveLaw::Parameters& rValues) override;
416 
422  void InitializeMaterialResponsePK1(ConstitutiveLaw::Parameters& rValues) override;
423 
429  void InitializeMaterialResponseKirchhoff(ConstitutiveLaw::Parameters& rValues) override;
430 
435  void CalculateGreenLagrangeStrain(Parameters &rValues);
436 
441  void CalculateAlmansiStrain(Parameters &rValues);
442 
453  void CalculateStrainsOnEachComponent(
454  const Vector& rStrainVector,
455  const Matrix& rParallelProjector,
456  const Matrix& rSerialProjector,
457  const Vector& rSerialStrainMatrix,
458  Vector& rStrainVectorMatrix,
459  Vector& rStrainVectorFiber,
461  const int Iteration = 1);
462 
475  void CalculateInitialApproximationSerialStrainMatrix(
476  const Vector& rStrainVector,
477  const Vector& rPreviousStrainVector,
478  const Properties& rMaterialProperties,
479  const Matrix& rParallelProjector,
480  const Matrix& rSerialProjector,
481  Matrix& rConstitutiveTensorMatrixSS,
482  Matrix& rConstitutiveTensorFiberSS,
483  Vector& rInitialApproximationSerialStrainMatrix,
485  const ConstitutiveLaw::StressMeasure& rStressMeasure);
486 
494  void IntegrateStressesOfFiberAndMatrix(
496  Vector& rMatrixStrainVector,
497  Vector& rFiberStrainVector,
498  Vector& rMatrixStressVector,
499  Vector& rFiberStressVector,
500  const ConstitutiveLaw::StressMeasure& rStressMeasure);
501 
513  void CheckStressEquilibrium(
515  const Vector& rStrainVector,
516  const Matrix& rSerialProjector,
517  const Vector& rMatrixStressVector,
518  const Vector& rFiberStressVector,
519  Vector& rStressResidual,
520  bool& rIsConverged,
521  const Matrix& rConstitutiveTensorMatrixSS,
522  const Matrix& rConstitutiveTensorFiberSS);
523 
531  void CorrectSerialStrainMatrix(
533  const Vector& rResidualStresses,
534  Vector& rSerialStrainMatrix,
535  const Matrix& rSerialProjector,
536  const ConstitutiveLaw::StressMeasure& rStressMeasure);
537 
542  void CalculateTangentTensor(ConstitutiveLaw::Parameters& rValues,
544 
549  {
550  return true;
551  }
552 
557  {
558  return true;
559  }
560 
569  int Check(
570  const Properties& rMaterialProperties,
571  const GeometryType& rElementGeometry,
572  const ProcessInfo& rCurrentProcessInfo
573  ) const override;
574 
578 
582 
586 
590 
592 
593  protected:
596 
600 
604 
608  ConstitutiveLaw::Pointer GetMatrixConstitutiveLaw()
609  {
610  return mpMatrixConstitutiveLaw;
611  }
612 
616  void SetMatrixConstitutiveLaw(ConstitutiveLaw::Pointer pMatrixConstitutiveLaw)
617  {
618  mpMatrixConstitutiveLaw = pMatrixConstitutiveLaw;
619  }
620 
624  ConstitutiveLaw::Pointer GetFiberConstitutiveLaw()
625  {
626  return mpFiberConstitutiveLaw;
627  }
628 
632  void SetFiberConstitutiveLaw(ConstitutiveLaw::Pointer pFiberConstitutiveLaw)
633  {
634  mpFiberConstitutiveLaw = pFiberConstitutiveLaw;
635  }
636 
642  {
643  const int parallel_components = inner_prod(mParallelDirections, mParallelDirections);
644  return this->GetStrainSize() - parallel_components;
645  }
646 
650 
654 
658 
662 
664  private:
667 
671 
672  ConstitutiveLaw::Pointer mpMatrixConstitutiveLaw;
673  ConstitutiveLaw::Pointer mpFiberConstitutiveLaw;
674  double mFiberVolumetricParticipation;
675  array_1d<double, 6> mParallelDirections = ZeroVector(6);
676  array_1d<double, 6> mPreviousStrainVector = ZeroVector(6);
677  Vector mPreviousSerialStrainMatrix = ZeroVector(GetNumberOfSerialComponents());
678  bool mIsPrestressed = false;
679 
683 
687 
691 
695 
699 
700  // Serialization
701 
702  friend class Serializer;
703 
704  void save(Serializer& rSerializer) const override
705  {
707  rSerializer.save("MatrixConstitutiveLaw", mpMatrixConstitutiveLaw);
708  rSerializer.save("FiberConstitutiveLaw", mpFiberConstitutiveLaw);
709  rSerializer.save("FiberVolumetricParticipation", mFiberVolumetricParticipation);
710  rSerializer.save("ParallelDirections", mParallelDirections);
711  rSerializer.save("PreviousStrainVector", mPreviousStrainVector);
712  rSerializer.save("PreviousSerialStrainMatrix", mPreviousSerialStrainMatrix);
713  rSerializer.save("IsPrestressed", mIsPrestressed);
714  }
715 
716  void load(Serializer& rSerializer) override
717  {
719  rSerializer.load("MatrixConstitutiveLaw", mpMatrixConstitutiveLaw);
720  rSerializer.load("FiberConstitutiveLaw", mpFiberConstitutiveLaw);
721  rSerializer.load("FiberVolumetricParticipation", mFiberVolumetricParticipation);
722  rSerializer.load("ParallelDirections", mParallelDirections);
723  rSerializer.load("PreviousStrainVector", mPreviousStrainVector);
724  rSerializer.load("PreviousSerialStrainMatrix", mPreviousSerialStrainMatrix);
725  rSerializer.load("IsPrestressed", mIsPrestressed);
726  }
727 
729 
730 }; // Class SerialParallelRuleOfMixturesLaw
731 
732 } // namespace Kratos
Definition: constitutive_law.h:47
StressMeasure
Definition: constitutive_law.h:69
@ StressMeasure_Cauchy
Definition: constitutive_law.h:73
std::size_t SizeType
Definition: constitutive_law.h:82
Geometry base class.
Definition: geometry.h:71
This class defines the node.
Definition: node.h:65
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
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
This CL implements the serial-parallel rule of mixtures developed by F.Rastellini.
Definition: serial_parallel_rule_of_mixtures_law.h:56
void SetMatrixConstitutiveLaw(ConstitutiveLaw::Pointer pMatrixConstitutiveLaw)
This method sets the constitutive law of the matrix material.
Definition: serial_parallel_rule_of_mixtures_law.h:616
Node NodeType
The node definition.
Definition: serial_parallel_rule_of_mixtures_law.h:62
SizeType WorkingSpaceDimension() override
Dimension of the law:
Definition: serial_parallel_rule_of_mixtures_law.h:133
SerialParallelRuleOfMixturesLaw()
Definition: serial_parallel_rule_of_mixtures_law.h:80
bool RequiresInitializeMaterialResponse() override
If the CL requires to initialize the material response, called by the element in InitializeSolutionSt...
Definition: serial_parallel_rule_of_mixtures_law.h:548
Geometry< NodeType > GeometryType
The geometry definition.
Definition: serial_parallel_rule_of_mixtures_law.h:65
int GetNumberOfSerialComponents()
This method returns the number of directions with serial behaviour (iso-stress behaviour)
Definition: serial_parallel_rule_of_mixtures_law.h:641
bool RequiresFinalizeMaterialResponse() override
If the CL requires to initialize the material response, called by the element in InitializeSolutionSt...
Definition: serial_parallel_rule_of_mixtures_law.h:556
KRATOS_CLASS_POINTER_DEFINITION(SerialParallelRuleOfMixturesLaw)
Counted pointer of SerialParallelRuleOfMixturesLaw.
SizeType GetStrainSize() const override
Voigt tensor size:
Definition: serial_parallel_rule_of_mixtures_law.h:141
SerialParallelRuleOfMixturesLaw(double FiberVolParticipation, const Vector &rParallelDirections)
Definition: serial_parallel_rule_of_mixtures_law.h:87
~SerialParallelRuleOfMixturesLaw() override
Definition: serial_parallel_rule_of_mixtures_law.h:111
ConstitutiveLaw::Pointer Clone() const override
Definition: serial_parallel_rule_of_mixtures_law.h:95
ConstitutiveLaw::Pointer GetMatrixConstitutiveLaw()
This method the constitutive law of the matrix material.
Definition: serial_parallel_rule_of_mixtures_law.h:608
void SetFiberConstitutiveLaw(ConstitutiveLaw::Pointer pFiberConstitutiveLaw)
This method sets the constitutive law of the fiber material.
Definition: serial_parallel_rule_of_mixtures_law.h:632
ConstitutiveLaw::Pointer GetFiberConstitutiveLaw()
This method the constitutive law of the fiber material.
Definition: serial_parallel_rule_of_mixtures_law.h:624
SerialParallelRuleOfMixturesLaw(SerialParallelRuleOfMixturesLaw const &rOther)
Definition: serial_parallel_rule_of_mixtures_law.h:101
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
Variable class contains all information needed to store and retrive data from a data container.
Definition: variable.h:63
#define KRATOS_SERIALIZE_SAVE_BASE_CLASS(Serializer, BaseType)
Definition: define.h:812
#define KRATOS_SERIALIZE_LOAD_BASE_CLASS(Serializer, BaseType)
Definition: define.h:815
Modeler::Pointer Create(const std::string &ModelerName, Model &rModel, const Parameters ModelParameters)
Checks if the modeler is registered.
Definition: modeler_factory.cpp:30
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
TExpression1Type::data_type inner_prod(AMatrix::MatrixExpression< TExpression1Type, TCategory1 > const &First, AMatrix::MatrixExpression< TExpression2Type, TCategory2 > const &Second)
Definition: amatrix_interface.h:592
def SetValue(entity, variable, value)
Definition: coupling_interface_data.py:256
def load(f)
Definition: ode_solve.py:307
Definition: constitutive_law.h:189