///////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////// More...
#include <constitutive_laws.h>
Public Member Functions | |
ConstitutiveLaw () | |
Empty constructor. More... | |
virtual | ~ConstitutiveLaw () |
Empty virtual destructor. More... | |
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 and deformed metric tensor and the matrix in which to return the stress tensor. 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) |
Calculate the derivatives of the contravariant 2nd Piola Kirchhoff stress tensor with respect to the deformed metric tensor. Arguments are the covariant undeformed and deformed metric tensor, the current value of the stress tensor and the rank four tensor in which to return the derivatives of the stress tensor The default implementation uses finite differences, but can be overloaded for constitutive laws in which an analytic formulation is possible. If the boolean flag symmetrize_tensor is false, only the "upper triangular" entries of the tensor will be filled in. This is a useful efficiency when using the derivatives in Jacobian calculations. More... | |
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 . 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 where the "pressure" is determined by . More... | |
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) |
Calculate the derivatives of the contravariant 2nd Piola Kirchhoff stress tensor . with respect to the deformed metric tensor. Also return the derivatives of the determinant of the deformed metric tensor with respect to the deformed metric tensor. This form is appropriate for truly-incompressible materials. The default implementation uses finite differences for the derivatives that depend on the constitutive law, but not for the derivatives of the determinant, which are generic. / If the boolean flag symmetrize_tensor is false, only the "upper triangular" entries of the tensor will be filled in. This is a useful efficiency when using the derivatives in Jacobian calculations. More... | |
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. Also return the contravariant deformed metric tensor, the generalised dilatation, and the inverse of the bulk modulus . This form is appropriate for near-incompressible materials for which where the "pressure" is determined from . More... | |
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) |
Calculate the derivatives of the contravariant 2nd Piola Kirchoff stress tensor with respect to the deformed metric tensor. Also return the derivatives of the generalised dilatation, with respect to the deformed metric tensor. This form is appropriate for near-incompressible materials. The default implementation uses finite differences. If the boolean flag symmetrize_tensor is false, only the "upper triangular" entries of the tensor will be filled in. This is a useful efficiency when using the derivatives in Jacobian calculations. More... | |
virtual bool | requires_incompressibility_constraint ()=0 |
Pure virtual function in which the user 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. More... | |
Protected Member Functions | |
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) |
Check for errors in the input, i.e. check that the dimensions of the arrays are all consistent. More... | |
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 tensor. More... | |
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 covariant tensor with respect to the components of the covariant tensor. More... | |
///////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////////
A class for constitutive laws for elements that solve the equations of solid mechanics based upon the principle of virtual displacements. In that formulation, the information required from a constitutive law is the (2nd Piola-Kirchhoff) stress tensor as a function of the (Green) strain :
The Green strain is defined as
where and are the metric tensors in the deformed and undeformed (stress-free) configurations, respectively. A specific ConstitutiveLaw needs to be implement the pure virtual function
Equation (1) shows that the strain may be calculated from the undeformed and deformed metric tensors. Frequently, these tensors are also required in the constitutive law itself. To avoid unnecessary re-computation of these quantities, we pass the deformed and undeformed metric tensor to calculate_second_piola_kirchhoff_stress(...)
rather than the strain tensor itself.
The functional form of the constitutive equation is different for compressible/incompressible/near-incompressible behaviour and we provide interfaces that are appropriate for all of these cases.
Compressible Behaviour:
If the material is compressible, the stress can be computed from the deformed and undeformed metric tensors,
using the interface
Incompressible Behaviour:
If the material is incompressible, its deformation is constrained by the condition that
which ensures that the volume of infinitesimal material elements remains constant during the deformation. This condition is typically enforced by a Lagrange multiplier which plays the role of a pressure. In such cases, the stress tensor has form
where only the deviatoric part of the stress tensor, depends directly on the strain. The pressure needs to be determined independently from (2). Given the deformed and undeformed metric tensors, the computation of the stress tensor for an incompressible material therefore requires the computation of the following quantities:
These quantities can be obtained from the following interface
Nearly Incompressible Behaviour:
If the material is nearly incompressible, it is advantageous to split the stress into its deviatoric and hydrostatic parts by writing the constitutive law in the form
where the deviatoric part of the stress tensor, depends on the strain. This form of the constitutive law is identical to that of the incompressible case and it involves a pressure which needs to be determined from an additional equation. In the incompressible case, this equation was given by the incompressibility constraint (2). Here, we need to augment the constitutive law (3) by a separate equation for the pressure. Generally this takes the form
where is the "bulk modulus", a material property that needs to be specified by the constitutive law. is the (generalised) dilatation, i.e. the relative change in the volume of an infinitesimal material element (or some suitable generalised quantitiy that is related to it). As the material approaches incompressibility, , so that infinitely large pressures would be required to achieve any change in volume. To facilitate the implementation of (4) as the equation for the pressure, we re-write it in the form
which only involves quantities that remain finite as we approach true incompressibility.
Given the deformed and undeformed metric tensors, the computation of the stress tensor for a nearly incompressible material therefore requires the computation of the following quantities:
These quantities can be obtained from the following interface
Definition at line 470 of file constitutive_laws.h.
|
inline |
Empty constructor.
Definition at line 501 of file constitutive_laws.h.
|
inlinevirtual |
Empty virtual destructor.
Definition at line 505 of file constitutive_laws.h.
|
protected |
Test whether two matrices are of equal dimensions.
This function is used to check whether matrices are of equal dimension.
Definition at line 52 of file constitutive_laws.cc.
References oomph::DenseMatrix< T >::ncol(), and oomph::DenseMatrix< T >::nrow().
Referenced by calculate_contravariant(), calculate_d_second_piola_kirchhoff_stress_dG(), and error_checking_in_input().
|
protected |
Calculate a contravariant tensor from a covariant tensor, and return the determinant of the covariant tensor.
The function to calculate the contravariant tensor from a covariant one.
Three dimensions
Definition at line 113 of file constitutive_laws.cc.
References are_matrices_of_equal_dimensions(), is_matrix_square(), oomph::DenseMatrix< T >::ncol(), and oomph::Global_string_for_annotation::string().
Referenced by oomph::GeneralisedHookean::calculate_second_piola_kirchhoff_stress(), and oomph::IsotropicStrainEnergyFunctionConstitutiveLaw::calculate_second_piola_kirchhoff_stress().
|
protected |
Calculate the derivatives of the contravariant tensor and the derivatives of the determinant of the covariant tensor with respect to the components of the covariant tensor.
The function to calculate the derivatives of the contravariant tensor and determinant of covariant tensor with respect to the components of the covariant tensor.
Three dimensions
Definition at line 225 of file constitutive_laws.cc.
References is_matrix_square(), oomph::DenseMatrix< T >::ncol(), and oomph::Global_string_for_annotation::string().
|
virtual |
Calculate the derivatives of the contravariant 2nd Piola Kirchhoff stress tensor . with respect to the deformed metric tensor. Also return the derivatives of the determinant of the deformed metric tensor with respect to the deformed metric tensor. This form is appropriate for truly-incompressible materials. The default implementation uses finite differences for the derivatives that depend on the constitutive law, but not for the derivatives of the determinant, which are generic. / If the boolean flag symmetrize_tensor is false, only the "upper triangular" entries of the tensor will be filled in. This is a useful efficiency when using the derivatives in Jacobian calculations.
Calculate the derivatives of the contravariant 2nd Piola Kirchhoff stress tensor . with respect to the deformed metric tensor. Also return the derivatives of the determinant of the deformed metric tensor with respect to the deformed metric tensor. This form is appropriate for truly-incompressible materials. The default implementation uses finite differences for the derivatives that depend on the constitutive law, but not for the derivatives of the determinant, which are generic.
Definition at line 453 of file constitutive_laws.cc.
References are_matrices_of_equal_dimensions(), calculate_second_piola_kirchhoff_stress(), oomph::GeneralisedElement::Default_fd_jacobian_step, i, and oomph::DenseMatrix< T >::ncol().
|
virtual |
Calculate the derivatives of the contravariant 2nd Piola Kirchoff stress tensor with respect to the deformed metric tensor. Also return the derivatives of the generalised dilatation, with respect to the deformed metric tensor. This form is appropriate for near-incompressible materials. The default implementation uses finite differences. If the boolean flag symmetrize_tensor is false, only the "upper triangular" entries of the tensor will be filled in. This is a useful efficiency when using the derivatives in Jacobian calculations.
Calculate the derivatives of the contravariant 2nd Piola Kirchoff stress tensor with respect to the deformed metric tensor. Also return the derivatives of the generalised dilatation, with respect to the deformed metric tensor. This form is appropriate for near-incompressible materials. The default implementation uses finite differences.
Definition at line 566 of file constitutive_laws.cc.
References are_matrices_of_equal_dimensions(), calculate_second_piola_kirchhoff_stress(), oomph::GeneralisedElement::Default_fd_jacobian_step, i, and oomph::DenseMatrix< T >::ncol().
|
virtual |
Calculate the derivatives of the contravariant 2nd Piola Kirchhoff stress tensor with respect to the deformed metric tensor. Arguments are the covariant undeformed and deformed metric tensor, the current value of the stress tensor and the rank four tensor in which to return the derivatives of the stress tensor The default implementation uses finite differences, but can be overloaded for constitutive laws in which an analytic formulation is possible. If the boolean flag symmetrize_tensor is false, only the "upper triangular" entries of the tensor will be filled in. This is a useful efficiency when using the derivatives in Jacobian calculations.
Calculate the derivatives of the contravariant 2nd Piola Kirchhoff stress tensor with respect to the deformed metric tensor. Arguments are the covariant undeformed and deformed metric tensor and the matrix in which to return the derivatives of the stress tensor The default implementation uses finite differences, but can be overloaded for constitutive laws in which an analytic formulation is possible.
Definition at line 352 of file constitutive_laws.cc.
References are_matrices_of_equal_dimensions(), calculate_second_piola_kirchhoff_stress(), oomph::GeneralisedElement::Default_fd_jacobian_step, i, and oomph::DenseMatrix< T >::ncol().
Referenced by oomph::PVDEquationsWithPressure< DIM >::get_d_stress_dG_upper(), and oomph::PVDEquations< DIM >::get_d_stress_dG_upper().
|
pure 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.
Implemented in oomph::IsotropicStrainEnergyFunctionConstitutiveLaw, and oomph::GeneralisedHookean.
Referenced by calculate_d_second_piola_kirchhoff_stress_dG(), oomph::AxisymmetricPVDEquations::get_stress(), oomph::PVDEquations< DIM >::get_stress(), oomph::AxisymmetricPVDEquationsWithPressure::get_stress(), and oomph::PVDEquationsWithPressure< DIM >::get_stress().
|
inlinevirtual |
Calculate the deviatoric part of the contravariant 2nd Piola Kirchhoff stress tensor . 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 where the "pressure" is determined by .
Reimplemented in oomph::IsotropicStrainEnergyFunctionConstitutiveLaw, and oomph::GeneralisedHookean.
Definition at line 547 of file constitutive_laws.h.
|
inlinevirtual |
Calculate the deviatoric part of the contravariant 2nd Piola Kirchoff stress tensor. Also return the contravariant deformed metric tensor, the generalised dilatation, and the inverse of the bulk modulus . This form is appropriate for near-incompressible materials for which where the "pressure" is determined from .
Reimplemented in oomph::IsotropicStrainEnergyFunctionConstitutiveLaw, and oomph::GeneralisedHookean.
Definition at line 592 of file constitutive_laws.h.
|
protected |
Check for errors in the input, i.e. check that the dimensions of the arrays are all consistent.
This function is used to provide simple error (bounce) checks on the input to any calculate_second_piola_kirchhoff_stress.
Definition at line 71 of file constitutive_laws.cc.
References are_matrices_of_equal_dimensions(), is_matrix_square(), and oomph::Global_string_for_annotation::string().
Referenced by oomph::GeneralisedHookean::calculate_second_piola_kirchhoff_stress(), and oomph::IsotropicStrainEnergyFunctionConstitutiveLaw::calculate_second_piola_kirchhoff_stress().
|
protected |
Test whether a matrix is square.
This function is used to check whether a matrix is square.
Definition at line 36 of file constitutive_laws.cc.
References oomph::DenseMatrix< T >::ncol(), and oomph::DenseMatrix< T >::nrow().
Referenced by calculate_contravariant(), calculate_d_contravariant_dG(), and error_checking_in_input().
|
pure virtual |
Pure virtual function in which the user 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.
Implemented in oomph::IsotropicStrainEnergyFunctionConstitutiveLaw, and oomph::GeneralisedHookean.
Referenced by oomph::RefineablePVDEquations< DIM >::fill_in_generic_contribution_to_residuals_pvd(), oomph::PVDEquations< DIM >::fill_in_generic_contribution_to_residuals_pvd(), oomph::RefineablePVDEquationsWithPressure< DIM >::fill_in_generic_residual_contribution_pvd_with_pressure(), and oomph::PVDEquationsWithPressure< DIM >::fill_in_generic_residual_contribution_pvd_with_pressure().