10 #if !defined (KRATOS_HYPERELASTIC_3D_LAW_H_INCLUDED)
11 #define KRATOS_HYPERELASTIC_3D_LAW_H_INCLUDED
31 class KRATOS_API(SOLID_MECHANICS_APPLICATION) HyperElastic3DLaw :
public ConstitutiveLaw
36 struct MaterialResponseVariables
43 double ThermalExpansionCoefficient;
44 double ReferenceTemperature;
51 Matrix DeformationGradientF;
55 const Vector* mpShapeFunctionsValues;
95 ConstitutiveLaw::Pointer
Clone()
const override;
152 const double& rValue,
165 const Vector& rShapeFunctionsValues )
override;
275 Matrix mInverseDeformationGradientF0;
277 double mDeterminantF0;
279 double mStrainEnergy;
327 Vector& rIsoStressVector);
335 Vector& rVolStressVector );
344 Matrix& rConstitutiveMatrix);
353 const unsigned int&
a,
const unsigned int&
b,
354 const unsigned int&
c,
const unsigned int&
d);
365 const Matrix & rIsoStressMatrix,
366 Matrix& rConstitutiveMatrix);
374 const Matrix & rIsoStressMatrix,
375 const unsigned int&
a,
const unsigned int&
b,
376 const unsigned int&
c,
const unsigned int&
d);
385 Matrix& rConstitutiveMatrix);
395 const unsigned int&
a,
const unsigned int&
b,
396 const unsigned int&
c,
const unsigned int&
d);
433 double & rTemperature);
494 void save(
Serializer& rSerializer)
const override
497 rSerializer.
save(
"mInverseDeformationGradientF0",mInverseDeformationGradientF0);
498 rSerializer.
save(
"mDeterminantF0",mDeterminantF0);
499 rSerializer.
save(
"mStrainEnergy",mStrainEnergy);
505 rSerializer.
load(
"mInverseDeformationGradientF0",mInverseDeformationGradientF0);
506 rSerializer.
load(
"mDeterminantF0",mDeterminantF0);
507 rSerializer.
load(
"mStrainEnergy",mStrainEnergy);
Definition: constitutive_law.h:47
StressMeasure
Definition: constitutive_law.h:69
std::size_t SizeType
Definition: constitutive_law.h:82
Geometry base class.
Definition: geometry.h:71
Definition: hyperelastic_3D_law.hpp:38
virtual void CalculateIsochoricStress(const MaterialResponseVariables &rElasticVariables, StressMeasure rStressMeasure, Vector &rIsoStressVector)
bool Has(const Variable< Matrix > &rThisVariable) override
Returns whether this constitutive Law has specified variable (Matrix)
int Check(const Properties &rMaterialProperties, const GeometryType &rElementGeometry, const ProcessInfo &rCurrentProcessInfo) const override
std::size_t SizeType
Definition: hyperelastic_3D_law.hpp:74
virtual void CalculateVolumetricStress(const MaterialResponseVariables &rElasticVariables, Vector &rVolStressVector)
void GetLawFeatures(Features &rFeatures) override
double & GetValue(const Variable< double > &rThisVariable, double &rValue) override
Returns the value of a specified variable (double)
virtual void CalculateConstitutiveMatrix(const MaterialResponseVariables &rElasticVariables, Matrix &rConstitutiveMatrix)
void FinalizeMaterialResponseKirchhoff(Parameters &rValues) override
virtual Vector & CalculateVolumetricPressureFactors(const MaterialResponseVariables &rElasticVariables, Vector &rFactors)
void FinalizeMaterialResponseCauchy(Parameters &rValues) override
Matrix & Transform2DTo3D(Matrix &rMatrix)
virtual double & CalculateVolumetricPressure(const MaterialResponseVariables &rElasticVariables, double &rPressure)
bool Has(const Variable< Vector > &rThisVariable) override
Returns whether this constitutive Law has specified variable (Vector)
~HyperElastic3DLaw() override
void CalculateMaterialResponseCauchy(Parameters &rValues) override
void CalculateStress(const MaterialResponseVariables &rElasticVariables, StressMeasure rStressMeasure, Vector &rStressVector)
void FinalizeMaterialResponsePK1(Parameters &rValues) override
virtual double & CalculateVolumetricFactor(const MaterialResponseVariables &rElasticVariables, double &rFactor)
void CalculateMaterialResponsePK2(Parameters &rValues) override
Matrix & GetValue(const Variable< Matrix > &rThisVariable, Matrix &rValue) override
Returns the value of a specified variable (Matrix)
void SetValue(const Variable< double > &rVariable, const double &rValue, const ProcessInfo &rCurrentProcessInfo) override
Sets the value of a specified variable (double)
HyperElastic3DLaw(const HyperElastic3DLaw &rOther)
void CalculateMaterialResponseKirchhoff(Parameters &rValues) override
KRATOS_CLASS_POINTER_DEFINITION(HyperElastic3DLaw)
virtual void CalculateAlmansiStrain(const Matrix &rLeftCauchyGreen, Vector &rStrainVector)
virtual void CalculateIsochoricConstitutiveMatrix(const MaterialResponseVariables &rElasticVariables, const Matrix &rIsoStressMatrix, Matrix &rConstitutiveMatrix)
void SetValue(const Variable< Vector > &rThisVariable, const Vector &rValue, const ProcessInfo &rCurrentProcessInfo) override
Sets the value of a specified variable (Vector)
void InitializeMaterial(const Properties &rMaterialProperties, const GeometryType &rElementGeometry, const Vector &rShapeFunctionsValues) override
void FinalizeMaterialResponsePK2(Parameters &rValues) override
virtual void CalculateGreenLagrangeStrain(const Matrix &rRightCauchyGreen, Vector &rStrainVector)
double & ConstitutiveComponent(double &rCabcd, const MaterialResponseVariables &rElasticVariables, const unsigned int &a, const unsigned int &b, const unsigned int &c, const unsigned int &d)
virtual void UpdateInternalVariables(Parameters &rValues)
double & CalculateValue(Parameters &rParameterValues, const Variable< double > &rThisVariable, double &rValue) override
Calculates the value of a specified variable (double)
ProcessInfo ProcessInfoType
Definition: hyperelastic_3D_law.hpp:72
double & VolumetricConstitutiveComponent(double &rCabcd, const MaterialResponseVariables &rElasticVariables, const Vector &rFactors, const unsigned int &a, const unsigned int &b, const unsigned int &c, const unsigned int &d)
bool Has(const Variable< double > &rThisVariable) override
Returns whether this constitutive Law has specified variable (double)
ConstitutiveLaw BaseType
Definition: hyperelastic_3D_law.hpp:73
virtual bool CheckParameters(Parameters &rValues)
SizeType GetStrainSize() const override
Definition: hyperelastic_3D_law.hpp:134
Vector & GetValue(const Variable< Vector > &rThisVariable, Vector &rValue) override
Returns the value of a specified variable (Vector)
SizeType WorkingSpaceDimension() override
Definition: hyperelastic_3D_law.hpp:126
void CalculateMaterialResponsePK1(Parameters &rValues) override
virtual double & CalculateDomainTemperature(const MaterialResponseVariables &rElasticVariables, double &rTemperature)
double & IsochoricConstitutiveComponent(double &rCabcd, const MaterialResponseVariables &rElasticVariables, const Matrix &rIsoStressMatrix, const unsigned int &a, const unsigned int &b, const unsigned int &c, const unsigned int &d)
void SetValue(const Variable< Matrix > &rThisVariable, const Matrix &rValue, const ProcessInfo &rCurrentProcessInfo) override
Sets the value of a specified variable (Matrix)
ConstitutiveLaw::Pointer Clone() const override
virtual void CalculateVolumetricConstitutiveMatrix(const MaterialResponseVariables &rElasticVariables, Matrix &rConstitutiveMatrix)
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
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
Internals::Matrix< double, AMatrix::dynamic, 1 > Vector
Definition: amatrix_interface.h:472
Internals::Matrix< double, AMatrix::dynamic, AMatrix::dynamic > Matrix
Definition: amatrix_interface.h:470
Geometry< Node > GeometryType
The definition of the geometry.
Definition: mortar_classes.h:37
a
Definition: generate_stokes_twofluid_element.py:77
b
Definition: generate_total_lagrangian_mixed_volumetric_strain_element.py:31
c
Definition: generate_weakly_compressible_navier_stokes_element.py:108
def load(f)
Definition: ode_solve.py:307
int d
Definition: ode_solve.py:397
Definition: constitutive_law.h:137
Definition: constitutive_law.h:189
Definition: hyperelastic_3D_law.hpp:43
const GeometryType & GetElementGeometry() const
Definition: hyperelastic_3D_law.hpp:62
void SetShapeFunctionsValues(const Vector &rShapeFunctionsValues)
Definition: hyperelastic_3D_law.hpp:59
void SetElementGeometry(const GeometryType &rElementGeometry)
Definition: hyperelastic_3D_law.hpp:60
const Vector & GetShapeFunctionsValues() const
Definition: hyperelastic_3D_law.hpp:61