oomph::GeneralisedHookean Class Reference

#include <constitutive_laws.h>

+ Inheritance diagram for oomph::GeneralisedHookean:

Public Member Functions

 GeneralisedHookean (double *nu_pt, double *e_pt)
 
 GeneralisedHookean (double *nu_pt)
 
virtual ~GeneralisedHookean ()
 Virtual destructor. More...
 
void calculate_second_piola_kirchhoff_stress (const DenseMatrix< double > &g, const DenseMatrix< double > &G, DenseMatrix< double > &sigma)
 
void calculate_second_piola_kirchhoff_stress (const DenseMatrix< double > &g, const DenseMatrix< double > &G, DenseMatrix< double > &sigma_dev, DenseMatrix< double > &G_contra, double &Gdet)
 
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)
 
bool requires_incompressibility_constraint ()
 
- Public Member Functions inherited from oomph::ConstitutiveLaw
 ConstitutiveLaw ()
 Empty constructor. More...
 
virtual ~ConstitutiveLaw ()
 Empty virtual destructor. More...
 
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)
 
virtual void calculate_d_second_piola_kirchhoff_stress_dG (const DenseMatrix< double > &g, const DenseMatrix< double > &G, const DenseMatrix< double > &sigma, const double &detG, const double &interpolated_solid_p, RankFourTensor< double > &d_sigma_dG, DenseMatrix< double > &d_detG_dG, const bool &symmetrize_tensor=true)
 
virtual void calculate_d_second_piola_kirchhoff_stress_dG (const DenseMatrix< double > &g, const DenseMatrix< double > &G, const DenseMatrix< double > &sigma, const double &gen_dil, const double &inv_kappa, const double &interpolated_solid_p, RankFourTensor< double > &d_sigma_dG, DenseMatrix< double > &d_gen_dil_dG, const bool &symmetrize_tensor=true)
 

Private Attributes

doubleNu_pt
 Poisson ratio. More...
 
doubleE_pt
 Young's modulus. More...
 
bool Must_delete_e
 

Additional Inherited Members

- Protected Member Functions inherited from oomph::ConstitutiveLaw
bool is_matrix_square (const DenseMatrix< double > &M)
 Test whether a matrix is square. More...
 
bool are_matrices_of_equal_dimensions (const DenseMatrix< double > &M1, const DenseMatrix< double > &M2)
 Test whether two matrices are of equal dimensions. More...
 
void error_checking_in_input (const DenseMatrix< double > &g, const DenseMatrix< double > &G, DenseMatrix< double > &sigma)
 
double calculate_contravariant (const DenseMatrix< double > &Gcov, DenseMatrix< double > &Gcontra)
 The function to calculate the contravariant tensor from a covariant one. More...
 
void calculate_d_contravariant_dG (const DenseMatrix< double > &Gcov, RankFourTensor< double > &dGcontra_dG, DenseMatrix< double > &d_detG_dG)
 

Detailed Description

////////////////////////////////////////////////////////////////////// ////////////////////////////////////////////////////////////////////// ////////////////////////////////////////////////////////////////////// Class for a "non-rational" extension of classical linear elasticity to large displacements:

\[ \sigma^{ij} = E^{ijkl} \gamma_{kl} \]

where

\[ E^{ijkl} = \frac{E}{(1+\nu)} \left( \frac{\nu}{(1-2\nu)} G^{ij} G^{kl} + \frac{1}{2} \left( G^{ik} G^{jl} + G^{il} G^{jk} \right) \right) \]

For small strains \( (| G_{ij} - g_{ij} | \ll 1)\) this approaches the version appropriate for linear elasticity, obtained by replacing \( G^{ij}\) with \( g^{ij}\).

We provide three versions of calculate_second_piola_kirchhoff_stress():

  1. If \( \nu \ne 1/2 \) (and not close to \( 1/2 \)), the constitutive law can be used directly in the above form, using the deformed and undeformed metric tensors as input.
  2. If the material is incompressible ( \( \nu = 1/2 \)), the first term in the above expression for \( E^{ijkl} \) is singular. We re-write the constitutive equation for this case as

    \[ \sigma^{ij} = -p G^{ij} + \frac{E}{3} \left( G^{ik} G^{jl} + G^{il} G^{jk} \right) \gamma_{kl} \]

    where the pressure \( p \) needs to be determined independently via the incompressibility constraint. In this case, the stress returned by calculate_second_piola_kirchhoff_stress() contains only the deviatoric part of the 2nd Piola Kirchhoff stress,

    \[ \overline{\sigma}^{ij} = \frac{E}{3} \left( G^{ik} G^{jl} + G^{il} G^{jk} \right) \gamma_{kl}. \]

    The function also returns the contravariant metric tensor \( G^{ij}\) (since it is needed to form the complete stress tensor), and the determinant of the deformed covariant metric tensor \( {\tt detG} = \det G_{ij} \) (since it is needed in the equation that enforces the incompressibility).

  3. If \( \nu \approx 1/2 \), the original form of the constitutive equation could be used, but the resulting equations tend to be ill-conditioned since they contain the product of the large "bulk modulus"

    \[ \kappa = \frac{E\nu}{(1+\nu)(1-2\nu)} \]

    and the small "generalised dilatation"

    \[ d = \frac{1}{2} G^{ij} (G_{ij}-g_{ij}). \]

    [ \( d \) represents the actual dilatation in the small strain limit; for large deformations it doesn't have any sensible interpretation (or does it?). It is simply the term that needs to go to zero as \( \kappa \to \infty\).] In this case, the stress returned by calculate_second_piola_kirchhoff_stress() contains only the deviatoric part of the 2nd Piola Kirchhoff stress,

    \[ \overline{\sigma}^{ij} = \frac{E}{3} \left( G^{ik} G^{jl} + G^{il} G^{jk} \right) \gamma_{kl}. \]

    The function also returns the contravariant metric tensor \( G^{ij}\) (since it is needed to form the complete stress tensor), the inverse of the bulk modulus, and the generalised dilatation (since they are needed in the equation that determines the pressure).

Constructor & Destructor Documentation

◆ GeneralisedHookean() [1/2]

oomph::GeneralisedHookean::GeneralisedHookean ( double nu_pt,
double e_pt 
)
inline

The constructor takes the pointers to values of material parameters: Poisson's ratio and Young's modulus.

704  : ConstitutiveLaw(), Nu_pt(nu_pt), E_pt(e_pt), Must_delete_e(false)
705  {
706  }
ConstitutiveLaw()
Empty constructor.
Definition: constitutive_laws.h:501
bool Must_delete_e
Definition: constitutive_laws.h:787
double * E_pt
Young's modulus.
Definition: constitutive_laws.h:783
double * Nu_pt
Poisson ratio.
Definition: constitutive_laws.h:780

◆ GeneralisedHookean() [2/2]

oomph::GeneralisedHookean::GeneralisedHookean ( double nu_pt)
inline

The constructor takes the pointers to value of Poisson's ratio . Young's modulus is set to E=1.0, implying that all stresses have been non-dimensionalised on on it.

713  : ConstitutiveLaw(),
714  Nu_pt(nu_pt),
715  E_pt(new double(1.0)),
716  Must_delete_e(true)
717  {
718  }

◆ ~GeneralisedHookean()

virtual oomph::GeneralisedHookean::~GeneralisedHookean ( )
inlinevirtual

Virtual destructor.

723  {
724  if (Must_delete_e) delete E_pt;
725  }

References E_pt, and Must_delete_e.

Member Function Documentation

◆ calculate_second_piola_kirchhoff_stress() [1/3]

void oomph::GeneralisedHookean::calculate_second_piola_kirchhoff_stress ( const DenseMatrix< double > &  g,
const DenseMatrix< double > &  G,
DenseMatrix< double > &  sigma 
)
virtual

Calculate the contravariant 2nd Piola Kirchhoff stress tensor. Arguments are the covariant undeformed and deformed metric tensor and the matrix in which to return the stress tensor

////////////////////////////////////////////////////////////////// ////////////////////////////////////////////////////////////////// ////////////////////////////////////////////////////////////////// Calculate the contravariant 2nd Piola Kirchhoff stress tensor. Arguments are the covariant undeformed (stress-free) and deformed metric tensors, g and G, and the matrix in which to return the stress tensor.

Implements oomph::ConstitutiveLaw.

686  {
687  // Error checking
688 #ifdef PARANOID
690 #endif
691 
692  // Find the dimension of the problem
693  unsigned dim = G.nrow();
694 
695  // Calculate the contravariant Deformed metric tensor
696  DenseMatrix<double> Gup(dim);
697  // We don't need the Jacobian so cast the function to void
698  (void)calculate_contravariant(G, Gup);
699 
700  // Premultiply some constants
701  double C1 = (*E_pt) / (2.0 * (1.0 + (*Nu_pt)));
702  double C2 = 2.0 * (*Nu_pt) / (1.0 - 2.0 * (*Nu_pt));
703 
704  // Strain tensor
705  DenseMatrix<double> strain(dim, dim);
706 
707  // Upper triangle
708  for (unsigned i = 0; i < dim; i++)
709  {
710  for (unsigned j = i; j < dim; j++)
711  {
712  strain(i, j) = 0.5 * (G(i, j) - g(i, j));
713  }
714  }
715 
716  // Copy across
717  for (unsigned i = 0; i < dim; i++)
718  {
719  for (unsigned j = 0; j < i; j++)
720  {
721  strain(i, j) = strain(j, i);
722  }
723  }
724 
725  // Compute upper triangle of stress
726  for (unsigned i = 0; i < dim; i++)
727  {
728  for (unsigned j = i; j < dim; j++)
729  {
730  // Initialise this component of sigma
731  sigma(i, j) = 0.0;
732  for (unsigned k = 0; k < dim; k++)
733  {
734  for (unsigned l = 0; l < dim; l++)
735  {
736  sigma(i, j) += C1 *
737  (Gup(i, k) * Gup(j, l) + Gup(i, l) * Gup(j, k) +
738  C2 * Gup(i, j) * Gup(k, l)) *
739  strain(k, l);
740  }
741  }
742  }
743  }
744 
745  // Copy across
746  for (unsigned i = 0; i < dim; i++)
747  {
748  for (unsigned j = 0; j < i; j++)
749  {
750  sigma(i, j) = sigma(j, i);
751  }
752  }
753  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
JacobiRotation< float > G
Definition: Jacobi_makeGivens.cpp:2
void error_checking_in_input(const DenseMatrix< double > &g, const DenseMatrix< double > &G, DenseMatrix< double > &sigma)
Definition: constitutive_laws.cc:71
double calculate_contravariant(const DenseMatrix< double > &Gcov, DenseMatrix< double > &Gcontra)
The function to calculate the contravariant tensor from a covariant one.
Definition: constitutive_laws.cc:113
char char char int int * k
Definition: level2_impl.h:374
double C1
"Mooney Rivlin" coefficient for generalised Mooney Rivlin law
Definition: TwenteMeshGluing.cpp:74
double C2
"Mooney Rivlin" coefficient for generalised Mooney Rivlin law
Definition: mpi/distribution/airy_cantilever/airy_cantilever2.cc:156
int sigma
Definition: calibrate.py:179
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References Global_Physical_Variables::C1, Global_Physical_Variables::C2, oomph::ConstitutiveLaw::calculate_contravariant(), oomph::ConstitutiveLaw::error_checking_in_input(), G, i, j, k, and calibrate::sigma.

Referenced by oomph::AnisotropicHookean::calculate_second_piola_kirchhoff_stress(), and calculate_second_piola_kirchhoff_stress().

◆ calculate_second_piola_kirchhoff_stress() [2/3]

void oomph::GeneralisedHookean::calculate_second_piola_kirchhoff_stress ( const DenseMatrix< double > &  g,
const DenseMatrix< double > &  G,
DenseMatrix< double > &  sigma_dev,
DenseMatrix< double > &  Gup,
double detG 
)
virtual

Calculate the deviatoric part \( \overline{ \sigma^{ij}}\) of the contravariant 2nd Piola Kirchhoff stress tensor \( \sigma^{ij}\). Also return the contravariant deformed metric tensor and the determinant of the deformed metric tensor. This form is appropriate for truly-incompressible materials for which \( \sigma^{ij} = - p G^{ij} +\overline{ \sigma^{ij}} \) where the "pressure" \( p \) is determined by \( \det G_{ij} - \det g_{ij} = 0 \).

Reimplemented from oomph::ConstitutiveLaw.

817  {
818  // Error checking
819 #ifdef PARANOID
820  error_checking_in_input(g, G, sigma_dev);
821 #endif
822 
823  // Find the dimension of the problem
824  unsigned dim = G.nrow();
825 
826  // Calculate the contravariant Deformed metric tensor
827  detG = calculate_contravariant(G, Gup);
828 
829  // Premultiply the appropriate physical constant
830  double C1 = (*E_pt) / (2.0 * (1.0 + (*Nu_pt)));
831 
832  // Strain tensor
833  DenseMatrix<double> strain(dim, dim);
834 
835  // Upper triangle
836  for (unsigned i = 0; i < dim; i++)
837  {
838  for (unsigned j = i; j < dim; j++)
839  {
840  strain(i, j) = 0.5 * (G(i, j) - g(i, j));
841  }
842  }
843 
844  // Copy across
845  for (unsigned i = 0; i < dim; i++)
846  {
847  for (unsigned j = 0; j < i; j++)
848  {
849  strain(i, j) = strain(j, i);
850  }
851  }
852 
853  // Compute upper triangle of stress
854  for (unsigned i = 0; i < dim; i++)
855  {
856  for (unsigned j = i; j < dim; j++)
857  {
858  // Initialise this component of sigma
859  sigma_dev(i, j) = 0.0;
860  for (unsigned k = 0; k < dim; k++)
861  {
862  for (unsigned l = 0; l < dim; l++)
863  {
864  sigma_dev(i, j) += C1 *
865  (Gup(i, k) * Gup(j, l) + Gup(i, l) * Gup(j, k)) *
866  strain(k, l);
867  }
868  }
869  }
870  }
871 
872  // Copy across
873  for (unsigned i = 0; i < dim; i++)
874  {
875  for (unsigned j = 0; j < i; j++)
876  {
877  sigma_dev(i, j) = sigma_dev(j, i);
878  }
879  }
880  }

References Global_Physical_Variables::C1, oomph::ConstitutiveLaw::calculate_contravariant(), oomph::ConstitutiveLaw::error_checking_in_input(), G, i, j, and k.

◆ calculate_second_piola_kirchhoff_stress() [3/3]

void oomph::GeneralisedHookean::calculate_second_piola_kirchhoff_stress ( const DenseMatrix< double > &  g,
const DenseMatrix< double > &  G,
DenseMatrix< double > &  sigma_dev,
DenseMatrix< double > &  Gup,
double gen_dil,
double inv_kappa 
)
virtual

Calculate the deviatoric part of the contravariant 2nd Piola Kirchoff stress tensor. Also return the contravariant deformed metric tensor, the generalised dilatation, \( d, \) and the inverse of the bulk modulus \( \kappa\). This form is appropriate for near-incompressible materials for which \( \sigma^{ij} = -p G^{ij} + \overline{ \sigma^{ij}} \) where the "pressure" \( p \) is determined from \( p / \kappa - d =0 \).

Reimplemented from oomph::ConstitutiveLaw.

772  {
773  // Find the dimension of the problem
774  unsigned dim = G.nrow();
775 
776  // Assign memory for the determinant of the deformed metric tensor
777  double detG = 0.0;
778 
779  // Compute deviatoric stress by calling the incompressible
780  // version of this function
781  calculate_second_piola_kirchhoff_stress(g, G, sigma_dev, Gup, detG);
782 
783  // Calculate the inverse of the "bulk" modulus
784  inv_kappa =
785  (1.0 - 2.0 * (*Nu_pt)) * (1.0 + (*Nu_pt)) / ((*E_pt) * (*Nu_pt));
786 
787  // Finally compute the generalised dilatation (i.e. the term that
788  // must be zero if \kappa \to \infty
789  gen_dil = 0.0;
790  for (unsigned i = 0; i < dim; i++)
791  {
792  for (unsigned j = 0; j < dim; j++)
793  {
794  gen_dil += Gup(i, j) * 0.5 * (G(i, j) - g(i, j));
795  }
796  }
797  }
void calculate_second_piola_kirchhoff_stress(const DenseMatrix< double > &g, const DenseMatrix< double > &G, DenseMatrix< double > &sigma)
Definition: constitutive_laws.cc:682

References calculate_second_piola_kirchhoff_stress(), E_pt, G, i, j, and Nu_pt.

◆ requires_incompressibility_constraint()

bool oomph::GeneralisedHookean::requires_incompressibility_constraint ( )
inlinevirtual

Pure virtual function in which the writer must declare if the constitutive equation requires an incompressible formulation in which the volume constraint is enforced explicitly. Used as a sanity check in PARANOID mode. False.

Implements oomph::ConstitutiveLaw.

774  {
775  return false;
776  }

Member Data Documentation

◆ E_pt

double* oomph::GeneralisedHookean::E_pt
private

◆ Must_delete_e

bool oomph::GeneralisedHookean::Must_delete_e
private

Boolean flag to indicate if storage for elastic modulus must be deleted in destructor

Referenced by ~GeneralisedHookean().

◆ Nu_pt

double* oomph::GeneralisedHookean::Nu_pt
private

Poisson ratio.

Referenced by calculate_second_piola_kirchhoff_stress().


The documentation for this class was generated from the following files: