29#ifndef OOMPH_CONSTITUTIVE_LAWS_HEADER
30#define OOMPH_CONSTITUTIVE_LAWS_HEADER
34#include <oomph-lib-config.h>
60 std::string error_message =
61 "The strain-energy function as a function of the strain-tensor,\n";
63 "gamma, is not implemented for this strain energy function.\n";
74 std::string error_message =
75 "The strain-energy function as a function of the strain\n ";
77 "invariants, I1, I2, I3, is not implemented for this strain\n ";
78 error_message +=
"energy function\n";
93 "Sorry, the FD setup of dW/dgamma hasn't been implemented yet",
110 for (
unsigned i = 0;
i < 3;
i++)
165 return (*
C1_pt) * (I[0] - 3.0) + (*
C2_pt) * (I[1] - 3.0);
257 double G = (*E_pt) / (2.0 * (1.0 + (*Nu_pt)));
258 return 0.5 * ((*C1_pt) * (I[0] - 3.0) + (
G - (*
C1_pt)) * (I[1] - 3.0) +
259 ((*C1_pt) - 2.0 *
G) * (I[2] - 1.0) +
260 (1.0 - (*Nu_pt)) *
G * (I[2] - 1.0) * (I[2] - 1.0) /
261 (2.0 * (1.0 - 2.0 * (*
Nu_pt))));
269 double G = (*E_pt) / (2.0 * (1.0 + (*Nu_pt)));
270 dWdI[0] = 0.5 * (*C1_pt);
271 dWdI[1] = 0.5 * (
G - (*C1_pt));
272 dWdI[2] = 0.5 * ((*C1_pt) - 2.0 *
G +
273 2.0 * (1.0 - (*Nu_pt)) *
G * (I[2] - 1.0) /
274 (2.0 * (1.0 - 2.0 * (*Nu_pt))));
555 "Incompressible formulation not implemented for this constitutive law",
578 const double& interpolated_solid_p,
601 "Near-incompressible formulation not implemented for constitutive law",
622 const double& interpolated_solid_p,
A class for constitutive laws for elements that solve the equations of solid mechanics based upon the...
void error_checking_in_input(const DenseMatrix< double > &g, const DenseMatrix< double > &G, DenseMatrix< double > &sigma)
Check for errors in the input, i.e. check that the dimensions of the arrays are all consistent.
virtual void calculate_second_piola_kirchhoff_stress(const DenseMatrix< double > &g, const DenseMatrix< double > &G, DenseMatrix< double > &sigma_dev, DenseMatrix< double > &Gcontra, double &gen_dil, double &inv_kappa)
Calculate the deviatoric part of the contravariant 2nd Piola Kirchoff stress tensor....
virtual bool requires_incompressibility_constraint()=0
Pure virtual function in which the user must declare if the constitutive equation requires an incompr...
double calculate_contravariant(const DenseMatrix< double > &Gcov, DenseMatrix< double > &Gcontra)
Calculate a contravariant tensor from a covariant tensor, and return the determinant of the covariant...
virtual void calculate_d_second_piola_kirchhoff_stress_dG(const DenseMatrix< double > &g, const DenseMatrix< double > &G, const DenseMatrix< double > &sigma, RankFourTensor< double > &d_sigma_dG, const bool &symmetrize_tensor=true)
Calculate the derivatives of the contravariant 2nd Piola Kirchhoff stress tensor with respect to the ...
void calculate_d_contravariant_dG(const DenseMatrix< double > &Gcov, RankFourTensor< double > &dGcontra_dG, DenseMatrix< double > &d_detG_dG)
Calculate the derivatives of the contravariant tensor and the derivatives of the determinant of the c...
bool is_matrix_square(const DenseMatrix< double > &M)
Test whether a matrix is square.
ConstitutiveLaw()
Empty constructor.
bool are_matrices_of_equal_dimensions(const DenseMatrix< double > &M1, const DenseMatrix< double > &M2)
Test whether two matrices are of equal dimensions.
virtual void calculate_second_piola_kirchhoff_stress(const DenseMatrix< double > &g, const DenseMatrix< double > &G, DenseMatrix< double > &sigma)=0
Calculate the contravariant 2nd Piola Kirchhoff stress tensor. Arguments are the covariant undeformed...
virtual ~ConstitutiveLaw()
Empty virtual destructor.
virtual void calculate_second_piola_kirchhoff_stress(const DenseMatrix< double > &g, const DenseMatrix< double > &G, DenseMatrix< double > &sigma_dev, DenseMatrix< double > &G_contra, double &Gdet)
Calculate the deviatoric part of the contravariant 2nd Piola Kirchhoff stress tensor ....
Class for a "non-rational" extension of classical linear elasticity to large displacements:
virtual ~GeneralisedHookean()
Virtual destructor.
bool Must_delete_e
Boolean flag to indicate if storage for elastic modulus must be deleted in destructor.
GeneralisedHookean(double *nu_pt)
The constructor takes the pointers to value of Poisson's ratio . Young's modulus is set to E=1....
void calculate_second_piola_kirchhoff_stress(const DenseMatrix< double > &g, const DenseMatrix< double > &G, DenseMatrix< double > &sigma)
Calculate the contravariant 2nd Piola Kirchhoff stress tensor. Arguments are the covariant undeformed...
double * E_pt
Young's modulus.
bool requires_incompressibility_constraint()
Pure virtual function in which the writer must declare if the constitutive equation requires an incom...
double * Nu_pt
Poisson ratio.
GeneralisedHookean(double *nu_pt, double *e_pt)
The constructor takes the pointers to values of material parameters: Poisson's ratio and Young's modu...
Generalisation of Mooney Rivlin constitutive law to compressible media as suggested on p....
GeneralisedMooneyRivlin(double *nu_pt, double *c1_pt)
Constructor takes the pointers to the constitutive parameters: Poisson's ratio, the Mooney-Rivlin par...
bool requires_incompressibility_constraint()
Pure virtual function in which the user must declare if the constitutive equation requires an incompr...
double W(const Vector< double > &I)
Return the strain energy in terms of the strain invariants.
bool Must_delete_e
Boolean flag to indicate if storage for elastic modulus must be deleted in destructor.
double * E_pt
Young's modulus.
GeneralisedMooneyRivlin(double *nu_pt, double *c1_pt, double *e_pt)
Constructor takes the pointers to the constitutive parameters: Poisson's ratio, the Mooney-Rivlin par...
virtual ~GeneralisedMooneyRivlin()
Virtual destructor.
double W(const DenseMatrix< double > &gamma)
Return the strain energy in terms of strain tensor.
void derivatives(Vector< double > &I, Vector< double > &dWdI)
Return the derivatives of the strain energy function with respect to the strain invariants.
double * Nu_pt
Poisson's ratio.
double * C1_pt
Mooney-Rivlin parameter.
A class for constitutive laws derived from strain-energy functions. Theory is in Green and Zerna.
StrainEnergyFunction * Strain_energy_function_pt
Pointer to the strain energy function.
bool requires_incompressibility_constraint()
State if the constitutive equation requires an incompressible formulation in which the volume constra...
IsotropicStrainEnergyFunctionConstitutiveLaw(StrainEnergyFunction *const &strain_energy_function_pt)
Constructor takes a pointer to the strain energy function.
void calculate_second_piola_kirchhoff_stress(const DenseMatrix< double > &g, const DenseMatrix< double > &G, DenseMatrix< double > &sigma)
Calculate the contravariant 2nd Piola Kirchhoff stress tensor. Arguments are the covariant undeformed...
MooneyRivlin strain-energy function. with constitutive parameters C1 and C2:
void derivatives(Vector< double > &I, Vector< double > &dWdI)
Return the derivatives of the strain energy function with respect to the strain invariants.
double * C2_pt
Pointer to second Mooney Rivlin constant.
virtual ~MooneyRivlin()
Empty Virtual destructor.
double W(const DenseMatrix< double > &gamma)
Return the strain energy in terms of strain tensor.
double W(const Vector< double > &I)
Return the strain energy in terms of the strain invariants.
double * C1_pt
Pointer to first Mooney Rivlin constant.
MooneyRivlin(double *c1_pt, double *c2_pt)
Constructor takes the pointer to the value of the constants.
bool requires_incompressibility_constraint()
Pure virtual function in which the user must declare if the constitutive equation requires an incompr...
An OomphLibError object which should be thrown when an run-time error is encountered....
Base class for strain energy functions to be used in solid mechanics computations.
StrainEnergyFunction()
Constructor takes no arguments.
virtual bool requires_incompressibility_constraint()=0
Pure virtual function in which the user must declare if the constitutive equation requires an incompr...
virtual void derivative(const DenseMatrix< double > &gamma, DenseMatrix< double > &dWdgamma)
Return the derivatives of the strain energy function with respect to the components of the strain ten...
virtual ~StrainEnergyFunction()
Empty virtual destructor.
virtual double W(const DenseMatrix< double > &gamma)
Return the strain energy in terms of the strain tensor.
virtual void derivatives(Vector< double > &I, Vector< double > &dWdI)
Return the derivatives of the strain energy function with respect to the strain invariants....
virtual double W(const Vector< double > &I)
Return the strain energy in terms of the strain invariants.
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).