Struct MooneyRivlin

Source
pub struct MooneyRivlin<P> { /* 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

Implementations§

Source§

impl<P> MooneyRivlin<P>
where P: Parameters,

Source

pub fn extra_modulus(&self) -> &Scalar

Returns the extra modulus.

Trait Implementations§

Source§

impl<P> Constitutive<P> for MooneyRivlin<P>
where P: Parameters,

Source§

fn new(parameters: P) -> Self

Constructs and returns a new constitutive model.
Source§

fn jacobian( &self, deformation_gradient: &DeformationGradient, ) -> Result<Scalar, ConstitutiveError>

Calculates and returns the Jacobian.
Source§

impl<P: Debug> Debug for MooneyRivlin<P>

Source§

fn fmt(&self, f: &mut Formatter<'_>) -> Result

Formats the value using the given formatter. Read more
Source§

impl<P> Elastic for MooneyRivlin<P>
where P: Parameters,

Source§

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>

\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>

Calculates and returns the first Piola-Kirchhoff stress. Read more
Source§

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>

Calculates and returns the second Piola-Kirchhoff stress. Read more
Source§

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 root( &self, applied_load: AppliedLoad, ) -> Result<DeformationGradient, OptimizeError>

Solve for the unknown components of the deformation gradient under an applied load. Read more
Source§

impl<P> Hyperelastic for MooneyRivlin<P>
where P: Parameters,

Source§

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§

fn minimize( &self, applied_load: AppliedLoad, ) -> Result<DeformationGradient, OptimizeError>

Solve for the unknown components of the deformation gradient under an applied load. Read more
Source§

impl<P> Solid for MooneyRivlin<P>
where P: Parameters,

Source§

fn bulk_modulus(&self) -> &Scalar

Returns the bulk modulus.
Source§

fn shear_modulus(&self) -> &Scalar

Returns the shear modulus.

Auto Trait Implementations§

§

impl<P> Freeze for MooneyRivlin<P>
where P: Freeze,

§

impl<P> RefUnwindSafe for MooneyRivlin<P>
where P: RefUnwindSafe,

§

impl<P> Send for MooneyRivlin<P>
where P: Send,

§

impl<P> Sync for MooneyRivlin<P>
where P: Sync,

§

impl<P> Unpin for MooneyRivlin<P>
where P: Unpin,

§

impl<P> UnwindSafe for MooneyRivlin<P>
where P: UnwindSafe,

Blanket Implementations§

Source§

impl<T> Any for T
where T: 'static + ?Sized,

Source§

fn type_id(&self) -> TypeId

Gets the TypeId of self. Read more
Source§

impl<T> Borrow<T> for T
where T: ?Sized,

Source§

fn borrow(&self) -> &T

Immutably borrows from an owned value. Read more
Source§

impl<T> BorrowMut<T> for T
where T: ?Sized,

Source§

fn borrow_mut(&mut self) -> &mut T

Mutably borrows from an owned value. Read more
Source§

impl<T> From<T> for T

Source§

fn from(t: T) -> T

Returns the argument unchanged.

Source§

impl<T, U> Into<U> for T
where U: From<T>,

Source§

fn into(self) -> U

Calls U::from(self).

That is, this conversion is whatever the implementation of From<T> for U chooses to do.

Source§

impl<T, U> TryFrom<U> for T
where U: Into<T>,

Source§

type Error = Infallible

The type returned in the event of a conversion error.
Source§

fn try_from(value: U) -> Result<T, <T as TryFrom<U>>::Error>

Performs the conversion.
Source§

impl<T, U> TryInto<U> for T
where U: TryFrom<T>,

Source§

type Error = <U as TryFrom<T>>::Error

The type returned in the event of a conversion error.
Source§

fn try_into(self) -> Result<U, <U as TryFrom<T>>::Error>

Performs the conversion.
Source§

impl<A, Y, U> OdeSolver<Y, U> for A
where A: Debug, Y: Tensor, U: TensorVec<Item = Y>,