pub struct Gent<'a> { /* private fields */ }
Expand description
The Gent hyperelastic constitutive model.1
Parameters
- The bulk modulus $
\kappa
$. - The shear modulus $
\mu
$. - The extensibility $
J_m
$.
External variables
- The deformation gradient $
\mathbf{F}
$.
Internal variables
- None.
Notes
- The Gent model reduces to the Neo-Hookean model when $
J_m\to\infty
$.
A.N. Gent, Rubber Chem. Technol. 69, 59 (1996). ↩
Trait Implementations§
Source§impl<'a> Constitutive<'a> for Gent<'a>
impl<'a> Constitutive<'a> for Gent<'a>
Source§fn new(parameters: Parameters<'a>) -> Self
fn new(parameters: Parameters<'a>) -> Self
Constructs and returns a new constitutive model.
Source§fn jacobian(
&self,
deformation_gradient: &DeformationGradient,
) -> Result<Scalar, ConstitutiveError>
fn jacobian( &self, deformation_gradient: &DeformationGradient, ) -> Result<Scalar, ConstitutiveError>
Calculates and returns the Jacobian.
Source§impl<'a> Elastic<'a> for Gent<'a>
impl<'a> Elastic<'a> for Gent<'a>
Source§fn cauchy_stress(
&self,
deformation_gradient: &DeformationGradient,
) -> Result<CauchyStress, ConstitutiveError>
fn cauchy_stress( &self, deformation_gradient: &DeformationGradient, ) -> Result<CauchyStress, ConstitutiveError>
\boldsymbol{\sigma}(\mathbf{F}) = \frac{J^{-1}\mu J_m {\mathbf{B}^* }'}{J_m - \mathrm{tr}(\mathbf{B}^* ) + 3} + \frac{\kappa}{2}\left(J - \frac{1}{J}\right)\mathbf{1}
Source§fn cauchy_tangent_stiffness(
&self,
deformation_gradient: &DeformationGradient,
) -> Result<CauchyTangentStiffness, ConstitutiveError>
fn cauchy_tangent_stiffness( &self, deformation_gradient: &DeformationGradient, ) -> Result<CauchyTangentStiffness, ConstitutiveError>
\begin{aligned}
\mathcal{T}_{ijkL}(\mathbf{F}) =
\,& \frac{J^{-5/3}\mu J_m}{J_m - \mathrm{tr}(\mathbf{B}^* ) + 3}\Bigg[\delta_{ik}F_{jL} + \delta_{jk}F_{iL} - \frac{2}{3}\,\delta_{ij}F_{kL} + \frac{2{B_{ij}^* }' F_{kL}}{J_m - \mathrm{tr}(\mathbf{B}^* ) + 3}
\\
&- \left(\frac{5}{3} + \frac{2}{3}\frac{\mathrm{tr}(\mathbf{B}^* )}{J_m - \mathrm{tr}(\mathbf{B}^* ) + 3}\right) J^{2/3} {B_{ij}^* }' F_{kL}^{-T} \Bigg] + \frac{\kappa}{2} \left(J + \frac{1}{J}\right)\delta_{ij}F_{kL}^{-T}
\end{aligned}
Source§fn first_piola_kirchhoff_stress(
&self,
deformation_gradient: &DeformationGradient,
) -> Result<FirstPiolaKirchhoffStress, ConstitutiveError>
fn first_piola_kirchhoff_stress( &self, deformation_gradient: &DeformationGradient, ) -> Result<FirstPiolaKirchhoffStress, ConstitutiveError>
Calculates and returns the first Piola-Kirchhoff stress. Read more
Source§fn first_piola_kirchhoff_tangent_stiffness(
&self,
deformation_gradient: &DeformationGradient,
) -> Result<FirstPiolaKirchhoffTangentStiffness, ConstitutiveError>
fn first_piola_kirchhoff_tangent_stiffness( &self, deformation_gradient: &DeformationGradient, ) -> Result<FirstPiolaKirchhoffTangentStiffness, ConstitutiveError>
Calculates and returns the tangent stiffness associated with the first Piola-Kirchhoff stress. Read more
Source§fn second_piola_kirchhoff_stress(
&self,
deformation_gradient: &DeformationGradient,
) -> Result<SecondPiolaKirchhoffStress, ConstitutiveError>
fn second_piola_kirchhoff_stress( &self, deformation_gradient: &DeformationGradient, ) -> Result<SecondPiolaKirchhoffStress, ConstitutiveError>
Calculates and returns the second Piola-Kirchhoff stress. Read more
Source§fn second_piola_kirchhoff_tangent_stiffness(
&self,
deformation_gradient: &DeformationGradient,
) -> Result<SecondPiolaKirchhoffTangentStiffness, ConstitutiveError>
fn second_piola_kirchhoff_tangent_stiffness( &self, deformation_gradient: &DeformationGradient, ) -> Result<SecondPiolaKirchhoffTangentStiffness, ConstitutiveError>
Calculates and returns the tangent stiffness associated with the second Piola-Kirchhoff stress. Read more
Source§fn solve(
&self,
applied_load: AppliedLoad,
) -> Result<(DeformationGradient, CauchyStress), ConstitutiveError>
fn solve( &self, applied_load: AppliedLoad, ) -> Result<(DeformationGradient, CauchyStress), ConstitutiveError>
Solve for the unknown components of the Cauchy stress and deformation gradient under an applied load.
Source§impl<'a> Hyperelastic<'a> for Gent<'a>
impl<'a> Hyperelastic<'a> for Gent<'a>
Source§fn helmholtz_free_energy_density(
&self,
deformation_gradient: &DeformationGradient,
) -> Result<Scalar, ConstitutiveError>
fn helmholtz_free_energy_density( &self, deformation_gradient: &DeformationGradient, ) -> Result<Scalar, ConstitutiveError>
a(\mathbf{F}) = -\frac{\mu J_m}{2}\,\ln\left[1 - \frac{\mathrm{tr}(\mathbf{B}^* ) - 3}{J_m}\right] + \frac{\kappa}{2}\left[\frac{1}{2}\left(J^2 - 1\right) - \ln J\right]
Auto Trait Implementations§
impl<'a> Freeze for Gent<'a>
impl<'a> RefUnwindSafe for Gent<'a>
impl<'a> Send for Gent<'a>
impl<'a> Sync for Gent<'a>
impl<'a> Unpin for Gent<'a>
impl<'a> UnwindSafe for Gent<'a>
Blanket Implementations§
Source§impl<T> BorrowMut<T> for Twhere
T: ?Sized,
impl<T> BorrowMut<T> for Twhere
T: ?Sized,
Source§fn borrow_mut(&mut self) -> &mut T
fn borrow_mut(&mut self) -> &mut T
Mutably borrows from an owned value. Read more