pub struct MooneyRivlin<'a> { /* private fields */ }
Expand description
The Mooney-Rivlin hyperelastic constitutive model.1,2
Parameters
- The bulk modulus $
\kappa
$. - The shear modulus $
\mu
$. - The extra modulus $
\mu_m
$.
External variables
- The deformation gradient $
\mathbf{F}
$.
Internal variables
- None.
Notes
- The Mooney-Rivlin model reduces to the Neo-Hookean model when $
\mu_m\to 0
$.
M. Mooney, J. Appl. Phys. 11, 582 (1940). ↩
R.S. Rivlin, Philos. Trans. R. Soc. London, Ser. A 241, 379 (1948). ↩
Trait Implementations§
Source§impl<'a> Constitutive<'a> for MooneyRivlin<'a>
impl<'a> Constitutive<'a> for MooneyRivlin<'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> Debug for MooneyRivlin<'a>
impl<'a> Debug for MooneyRivlin<'a>
Source§impl<'a> Elastic<'a> for MooneyRivlin<'a>
impl<'a> Elastic<'a> for MooneyRivlin<'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{\mu - \mu_m}{J}\, {\mathbf{B}^* }' - \frac{\mu_m}{J}\left(\mathbf{B}^{* -1}\right)' + \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{\mu-\mu_m}{J^{5/3}}\left(\delta_{ik}F_{jL} + \delta_{jk}F_{iL} - \frac{2}{3}\,\delta_{ij}F_{kL}- \frac{5}{3} \, B_{ij}'F_{kL}^{-T} \right)
+ \frac{\kappa}{2} \left(J + \frac{1}{J}\right)\delta_{ij}F_{kL}^{-T}
\\
&- \frac{\mu_m}{J}\left[ \frac{2}{3}\,B_{ij}^{* -1}F_{kL}^{-T} - B_{ik}^{* -1}F_{jL}^{-T} - B_{ik}^{* -1}F_{iL}^{-T} + \frac{2}{3}\,\delta_{ij}\left(B_{km}^{* -1}\right)'F_{mL}^{-T} - \left(B_{ij}^{* -1}\right)'F_{kL}^{-T} \right]
\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 MooneyRivlin<'a>
impl<'a> Hyperelastic<'a> for MooneyRivlin<'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 - \mu_m}{2}\left[\mathrm{tr}(\mathbf{B}^* ) - 3\right] + \frac{\mu_m}{2}\left[I_2(\mathbf{B}^*) - 3\right] + \frac{\kappa}{2}\left[\frac{1}{2}\left(J^2 - 1\right) - \ln J\right]
Source§impl<'a> Solid<'a> for MooneyRivlin<'a>
impl<'a> Solid<'a> for MooneyRivlin<'a>
Source§fn bulk_modulus(&self) -> &Scalar
fn bulk_modulus(&self) -> &Scalar
Returns the bulk modulus.
Source§fn shear_modulus(&self) -> &Scalar
fn shear_modulus(&self) -> &Scalar
Returns the shear modulus.
Auto Trait Implementations§
impl<'a> Freeze for MooneyRivlin<'a>
impl<'a> RefUnwindSafe for MooneyRivlin<'a>
impl<'a> Send for MooneyRivlin<'a>
impl<'a> Sync for MooneyRivlin<'a>
impl<'a> Unpin for MooneyRivlin<'a>
impl<'a> UnwindSafe for MooneyRivlin<'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