pub struct ArrudaBoyce<'a> { /* private fields */ }
Expand description
The Arruda-Boyce hyperelastic constitutive model.1
Parameters
- The bulk modulus $
\kappa
$. - The shear modulus $
\mu
$. - The number of links $
N_b
$.
External variables
- The deformation gradient $
\mathbf{F}
$.
Internal variables
- None.
Notes
- The nondimensional end-to-end length per link of the chains is $
\gamma=\sqrt{\mathrm{tr}(\mathbf{B}^*)/3N_b}
$. - The nondimensional force is given by the inverse Langevin function as $
\eta=\mathcal{L}^{-1}(\gamma)
$. - The initial values are given by $
\gamma_0=\sqrt{1/3N_b}
$ and $\eta_0=\mathcal{L}^{-1}(\gamma_0)
$. - The Arruda-Boyce model reduces to the Neo-Hookean model when $
N_b\to\infty
$.
E.M. Arruda and M.C. Boyce, J. Mech. Phys. Solids 41, 389 (1993). ↩
Trait Implementations§
Source§impl<'a> Constitutive<'a> for ArrudaBoyce<'a>
impl<'a> Constitutive<'a> for ArrudaBoyce<'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 ArrudaBoyce<'a>
impl<'a> Debug for ArrudaBoyce<'a>
Source§impl<'a> Elastic<'a> for ArrudaBoyce<'a>
impl<'a> Elastic<'a> for ArrudaBoyce<'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\gamma_0\eta}{J\gamma\eta_0}\,{\mathbf{B}^*}' + \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\gamma_0\eta}{J^{5/3}\gamma\eta_0}\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{\mu\gamma_0\eta}{3J^{7/3}N_b\gamma^2\eta_0}\left(\frac{1}{\eta\mathcal{L}'(\eta)} - \frac{1}{\gamma}\right)B_{ij}'B_{km}'F_{mL}^{-T} + \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 ArrudaBoyce<'a>
impl<'a> Hyperelastic<'a> for ArrudaBoyce<'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{3\mu N_b\gamma_0}{\eta_0}\left[\gamma\eta - \gamma_0\eta_0 - \ln\left(\frac{\eta_0\sinh\eta}{\eta\sinh\eta_0}\right) \right] + \frac{\kappa}{2}\left[\frac{1}{2}\left(J^2 - 1\right) - \ln J\right]
Source§impl<'a> Solid<'a> for ArrudaBoyce<'a>
impl<'a> Solid<'a> for ArrudaBoyce<'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 ArrudaBoyce<'a>
impl<'a> RefUnwindSafe for ArrudaBoyce<'a>
impl<'a> Send for ArrudaBoyce<'a>
impl<'a> Sync for ArrudaBoyce<'a>
impl<'a> Unpin for ArrudaBoyce<'a>
impl<'a> UnwindSafe for ArrudaBoyce<'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