conspire::constitutive::solid::hyperelastic

Struct SaintVenantKirchhoff

Source
pub struct SaintVenantKirchhoff<'a> { /* private fields */ }
Expand description

The Saint Venant-Kirchhoff hyperelastic constitutive model.

Parameters

  • The bulk modulus $\kappa$.
  • The shear modulus $\mu$.

External variables

  • The deformation gradient $\mathbf{F}$.

Internal variables

  • None.

Notes

  • The Green-Saint Venant strain measure is given by $\mathbf{E}=\tfrac{1}{2}(\mathbf{C}-\mathbf{1})$.

Trait Implementations§

Source§

impl<'a> Constitutive<'a> for SaintVenantKirchhoff<'a>

Source§

fn new(parameters: Parameters<'a>) -> 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<'a> Debug for SaintVenantKirchhoff<'a>

Source§

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

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

impl<'a> Elastic<'a> for SaintVenantKirchhoff<'a>

Source§

fn second_piola_kirchhoff_stress( &self, deformation_gradient: &DeformationGradient, ) -> Result<SecondPiolaKirchhoffStress, ConstitutiveError>

\mathbf{S}(\mathbf{F}) = 2\mu\mathbf{E}' + \kappa\,\mathrm{tr}(\mathbf{E})\mathbf{1}
Source§

fn second_piola_kirchhoff_tangent_stiffness( &self, deformation_gradient: &DeformationGradient, ) -> Result<SecondPiolaKirchhoffTangentStiffness, ConstitutiveError>

\mathcal{G}_{IJkL}(\mathbf{F}) = \mu\,\delta_{JL}F_{kI} + \mu\,\delta_{IL}F_{kJ} + \left(\kappa - \frac{2}{3}\,\mu\right)\delta_{IJ}F_{kL}
Source§

fn cauchy_stress( &self, deformation_gradient: &DeformationGradient, ) -> Result<CauchyStress, ConstitutiveError>

Calculates and returns the Cauchy stress. Read more
Source§

fn cauchy_tangent_stiffness( &self, deformation_gradient: &DeformationGradient, ) -> Result<CauchyTangentStiffness, ConstitutiveError>

Calculates and returns the tangent stiffness associated with the Cauchy stress. Read more
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 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 SaintVenantKirchhoff<'a>

Source§

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

a(\mathbf{F}) = \mu\,\mathrm{tr}(\mathbf{E}^2) + \frac{1}{2}\left(\kappa - \frac{2}{3}\,\mu\right)\mathrm{tr}(\mathbf{E})^2
Source§

impl<'a> Solid<'a> for SaintVenantKirchhoff<'a>

Source§

fn bulk_modulus(&self) -> &Scalar

Returns the bulk modulus.
Source§

fn shear_modulus(&self) -> &Scalar

Returns the shear modulus.

Auto Trait Implementations§

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