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.
small_strain_orthotropic_3D_law.hpp
Go to the documentation of this file.
1 //
2 // Project Name: KratosConstitutiveModelsApplication $
3 // Created by: $Author: JMCarbonell $
4 // Last modified by: $Co-Author: $
5 // Date: $Date: April 2017 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 #if !defined(KRATOS_SMALL_STRAIN_ORTHOTROPIC_3D_LAW_H_INCLUDED )
11 #define KRATOS_SMALL_STRAIN_ORTHOTROPIC_3D_LAW_H_INCLUDED
12 
13 // System includes
14 
15 // External includes
16 
17 // Project includes
19 
20 namespace Kratos
21 {
24 
27 
31 
35 
39 
43 
45 
47  class KRATOS_API(CONSTITUTIVE_MODELS_APPLICATION) SmallStrainOrthotropic3DLaw : public SmallStrain3DLaw
48  {
49  public:
52 
55 
59 
62 
65 
68 
71  {
73  return *this;
74  }
75 
77  ConstitutiveLaw::Pointer Clone() const override
78  {
79  return Kratos::make_shared<SmallStrainOrthotropic3DLaw>(*this);
80  }
81 
84 
85 
89 
91  SizeType WorkingSpaceDimension() override { return 3; }
92 
94  SizeType GetStrainSize() const override { return 6; }
95 
97  void GetLawFeatures(Features& rFeatures) override
98  {
100 
101  //Set the type of law
102  rFeatures.mOptions.Set( THREE_DIMENSIONAL_LAW );
103  rFeatures.mOptions.Set( INFINITESIMAL_STRAINS );
104  rFeatures.mOptions.Set( ANISOTROPIC );
105 
106  //Get model features
107  GetModelFeatures(rFeatures);
108 
109  //Set strain measure required by the consitutive law
110  rFeatures.mStrainMeasures.push_back(StrainMeasure_Infinitesimal);
111  rFeatures.mStrainMeasures.push_back(StrainMeasure_Deformation_Gradient);
112 
113  //Set the strain size
114  rFeatures.mStrainSize = GetStrainSize();
115 
116  //Set the spacedimension
117  rFeatures.mSpaceDimension = WorkingSpaceDimension();
118 
119  KRATOS_CATCH(" ")
120  }
121 
125 
127  {
128  KRATOS_TRY
129 
130  //0.- Check if the constitutive parameters are passed correctly to the law calculation
131  //CheckParameters(rValues);
132 
133  const Flags& rOptions = rValues.GetOptions();
134 
135  const Properties& rProperties = rValues.GetMaterialProperties();
136 
137  Vector& rStrainVector = rValues.GetStrainVector();
138  Vector& rStressVector = rValues.GetStressVector();
139 
140 
141  // Calculate total Kirchhoff stress
142 
143  if( rOptions.Is( ConstitutiveLaw::COMPUTE_STRESS ) ){
144 
145  Matrix& rConstitutiveMatrix = rValues.GetConstitutiveMatrix();
146 
147  this->CalculateConstitutiveMatrix( rConstitutiveMatrix, rProperties);
148 
149  this->CalculateStress( rStrainVector, rConstitutiveMatrix, rStressVector );
150 
151  }
152  else if( rOptions.Is( ConstitutiveLaw::COMPUTE_CONSTITUTIVE_TENSOR ) ){
153 
154  Matrix& rConstitutiveMatrix = rValues.GetConstitutiveMatrix();
155  this->CalculateConstitutiveMatrix(rConstitutiveMatrix, rProperties);
156 
157  }
158 
159  // std::cout<<" StrainVector "<<rValues.GetStrainVector()<<std::endl;
160  // std::cout<<" StressVector "<<rValues.GetStressVector()<<std::endl;
161  // std::cout<<" ConstitutiveMatrix "<<rValues.GetConstitutiveMatrix()<<std::endl;
162 
163  KRATOS_CATCH(" ")
164 
165  }
166 
170 
171 
175 
176 
186  int Check(const Properties& rProperties,
187  const GeometryType& rElementGeometry,
188  const ProcessInfo& rCurrentProcessInfo) const override
189  {
190 
191  if(YOUNG_MODULUS_X.Key() == 0 || !rProperties.Has(YOUNG_MODULUS_X))
192  KRATOS_ERROR << "YOUNG_MODULUS_X has Key zero or invalid value" << std::endl;
193 
194  if(YOUNG_MODULUS_Y.Key() == 0 || !rProperties.Has(YOUNG_MODULUS_Y))
195  KRATOS_ERROR << "YOUNG_MODULUS_Y has Key zero or invalid value" << std::endl;
196 
197  if(YOUNG_MODULUS_Z.Key() == 0 || !rProperties.Has(YOUNG_MODULUS_Z))
198  KRATOS_ERROR << "YOUNG_MODULUS_Z has Key zero or invalid value" << std::endl;
199 
200  if(POISSON_RATIO_XY.Key() == 0 || !rProperties.Has(POISSON_RATIO_XY))
201  KRATOS_ERROR << "POISSON_RATIO_XY has Key zero invalid value" << std::endl;
202 
203  if(POISSON_RATIO_YZ.Key() == 0 || !rProperties.Has(POISSON_RATIO_YZ))
204  KRATOS_ERROR << "POISSON_RATIO_YZ has Key zero invalid value" << std::endl;
205 
206  if(POISSON_RATIO_XZ.Key() == 0 || !rProperties.Has(POISSON_RATIO_XZ))
207  KRATOS_ERROR << "POISSON_RATIO_XZ has Key zero invalid value" << std::endl;
208 
209  if(DENSITY.Key() == 0 || !rProperties.Has(DENSITY))
210  KRATOS_ERROR << "DENSITY has Key zero or invalid value" << std::endl;
211 
212  return 0;
213 
214  }
215 
216 
217 
221 
223  std::string Info() const override
224  {
225  std::stringstream buffer;
226  buffer << "SmallStrainOrthotropic3DLaw" ;
227  return buffer.str();
228  }
229 
231  void PrintInfo(std::ostream& rOStream) const override {rOStream << "SmallStrainOrthotropic3DLaw";}
232 
234  void PrintData(std::ostream& rOStream) const override {}
235 
236 
240 
241 
243 
244  protected:
247 
248 
252 
253 
257 
258 
266  void CalculateConstitutiveMatrix( Matrix& rConstitutiveMatrix,
267  const Properties& rProperties) override
268  {
269  KRATOS_TRY
270 
271  // Orthotropic constitutive matrix
272  double E1 = rProperties[YOUNG_MODULUS_X];
273  double E2 = rProperties[YOUNG_MODULUS_Y];
274  double E3 = rProperties[YOUNG_MODULUS_Z];
275 
276  double v12 = rProperties[POISSON_RATIO_XY];
277  double v23 = rProperties[POISSON_RATIO_YZ];
278  double v13 = rProperties[POISSON_RATIO_XZ];
279 
280  double P1 = 1.0/(E2*E2*v12*v12 + 2.0*E3*E2*v12*v13*v23 + E3*E2*v13*v13 - E1*E2 + E1*E3*v23*v23);
281  double P2 = E1*E1;
282  double P3 = E2*E2;
283  double P4 = E1*v23 + E2*v12*v13;
284  double P5 = E2*v12 + E3*v13*v23;
285  double P6 = E3*E3;
286 
287  rConstitutiveMatrix(0, 0) = -P1*P2*(- E3*v23*v23 + E2);
288  rConstitutiveMatrix(0, 1) = -E1*E2*P1*P5;
289  rConstitutiveMatrix(0, 2) = -E2*E3*P1*(E1*v13 + E1*v12*v23);
290  rConstitutiveMatrix(1, 0) = -E1*E2*P1*P5;
291  rConstitutiveMatrix(1, 1) = -P1*P3*(- E3*v13*v13 + E1);
292  rConstitutiveMatrix(1, 2) = -E2*E3*P1*P4;
293  rConstitutiveMatrix(2, 0) = -E1*E2*E3*P1*(v13 + v12*v23);
294  rConstitutiveMatrix(2, 1) = -E2*E3*P1*P4;
295  rConstitutiveMatrix(2, 2) = -E2*E3*P1*(- E2*v12*v12 + E1);
296  rConstitutiveMatrix(3, 3) = (E2*P2)/(P2 + v12*(P2 + P3) + E1*E2)/2.0;
297  rConstitutiveMatrix(4, 4) = (E3*P3)/(P3 + v23*(P3 + P6) + E2*E3)/2.0;
298  rConstitutiveMatrix(5, 5) = (E3*P2)/(P2 + v13*(P2 + P6) + E1*E3)/2.0;
299 
300  KRATOS_CATCH(" ")
301  }
302 
303 
307 
308 
312 
313 
317 
318 
322 
323 
325 
326  private:
329 
330 
334 
335 
339 
340 
344 
345 
349 
350 
354 
355 
359  friend class Serializer;
360 
361  void save(Serializer& rSerializer) const override
362  {
364  }
365 
366  void load(Serializer& rSerializer) override
367  {
369  }
370 
371 
375 
377 
378  }; // Class SmallStrainOrthotropic3DLaw
379 
381 
384 
385 
389 
390 
392 
394 
395 } // namespace Kratos.
396 
397 #endif // KRATOS_SMALL_STRAIN_ORTHOTROPIC_3D_LAW_H_INCLUDED defined
std::size_t SizeType
Definition: constitutive_law.h:82
Definition: flags.h:58
void Set(const Flags ThisFlag)
Definition: flags.cpp:33
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
bool Has(TVariableType const &rThisVariable) const
Definition: properties.h:578
The serialization consists in storing the state of an object into a storage format like data file or ...
Definition: serializer.h:123
Definition: small_strain_3D_law.hpp:33
SmallStrain3DLaw & operator=(const SmallStrain3DLaw &rOther)
Assignment operator.
Definition: small_strain_3D_law.cpp:59
ModelType::Pointer ModelTypePointer
Definition: small_strain_3D_law.hpp:40
Short class definition.
Definition: small_strain_orthotropic_3D_law.hpp:48
SmallStrainOrthotropic3DLaw(const SmallStrainOrthotropic3DLaw &rOther)
Copy constructor.
Definition: small_strain_orthotropic_3D_law.hpp:67
~SmallStrainOrthotropic3DLaw() override
Destructor.
Definition: small_strain_orthotropic_3D_law.hpp:83
KRATOS_CLASS_POINTER_DEFINITION(SmallStrainOrthotropic3DLaw)
Pointer definition of SmallStrainOrthotropic3DLaw.
SmallStrainOrthotropic3DLaw(ModelTypePointer pModel)
Constructor.
Definition: small_strain_orthotropic_3D_law.hpp:64
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: small_strain_orthotropic_3D_law.hpp:231
ConstitutiveLaw::Pointer Clone() const override
Clone.
Definition: small_strain_orthotropic_3D_law.hpp:77
SizeType WorkingSpaceDimension() override
Law Dimension.
Definition: small_strain_orthotropic_3D_law.hpp:91
void GetLawFeatures(Features &rFeatures) override
Law Features.
Definition: small_strain_orthotropic_3D_law.hpp:97
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: small_strain_orthotropic_3D_law.hpp:234
SmallStrainOrthotropic3DLaw()
Default constructor.
Definition: small_strain_orthotropic_3D_law.hpp:61
std::string Info() const override
Turn back information as a string.
Definition: small_strain_orthotropic_3D_law.hpp:223
void CalculateConstitutiveMatrix(Matrix &rConstitutiveMatrix, const Properties &rProperties) override
Definition: small_strain_orthotropic_3D_law.hpp:266
int Check(const Properties &rProperties, const GeometryType &rElementGeometry, const ProcessInfo &rCurrentProcessInfo) const override
Definition: small_strain_orthotropic_3D_law.hpp:186
void CalculateMaterialResponseKirchhoff(Parameters &rValues) override
Definition: small_strain_orthotropic_3D_law.hpp:126
SizeType GetStrainSize() const override
Law Voigt Strain Size.
Definition: small_strain_orthotropic_3D_law.hpp:94
SmallStrainOrthotropic3DLaw & operator=(SmallStrainOrthotropic3DLaw const &rOther)
Assignment operator.
Definition: small_strain_orthotropic_3D_law.hpp:70
#define KRATOS_SERIALIZE_SAVE_BASE_CLASS(Serializer, BaseType)
Definition: define.h:812
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
#define KRATOS_SERIALIZE_LOAD_BASE_CLASS(Serializer, BaseType)
Definition: define.h:815
#define KRATOS_ERROR
Definition: exception.h:161
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
def load(f)
Definition: ode_solve.py:307
Definition: constitutive_law.h:137
SizeType mStrainSize
Definition: constitutive_law.h:152
std::vector< StrainMeasure > mStrainMeasures
Definition: constitutive_law.h:154
SizeType mSpaceDimension
Definition: constitutive_law.h:153
Flags mOptions
Definition: constitutive_law.h:151
Definition: constitutive_law.h:189
Flags & GetOptions()
Definition: constitutive_law.h:412
StrainVectorType & GetStrainVector()
Definition: constitutive_law.h:435
const Properties & GetMaterialProperties()
Definition: constitutive_law.h:457
VoigtSizeMatrixType & GetConstitutiveMatrix()
Definition: constitutive_law.h:446
StressVectorType & GetStressVector()
Definition: constitutive_law.h:440