pub struct SaintVenantKirchhoff<'a> { /* private fields */ }
Expand description
The Saint Venant-Kirchhoff thermohyperelastic constitutive model.
Parameters
- The bulk modulus $
\kappa
$. - The shear modulus $
\mu
$. - The coefficient of thermal expansion $
\alpha
$. - The reference temperature $
T_\mathrm{ref}
$.
External variables
- The deformation gradient $
\mathbf{F}
$. - The temperature $
T
$.
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>
impl<'a> Constitutive<'a> for SaintVenantKirchhoff<'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 SaintVenantKirchhoff<'a>
impl<'a> Debug for SaintVenantKirchhoff<'a>
Source§impl<'a> Solid<'a> for SaintVenantKirchhoff<'a>
impl<'a> Solid<'a> for SaintVenantKirchhoff<'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.
Source§impl<'a> Thermoelastic<'a> for SaintVenantKirchhoff<'a>
impl<'a> Thermoelastic<'a> for SaintVenantKirchhoff<'a>
Source§fn second_piola_kirchhoff_stress(
&self,
deformation_gradient: &DeformationGradient,
temperature: &Scalar,
) -> Result<SecondPiolaKirchhoffStress, ConstitutiveError>
fn second_piola_kirchhoff_stress( &self, deformation_gradient: &DeformationGradient, temperature: &Scalar, ) -> Result<SecondPiolaKirchhoffStress, ConstitutiveError>
Calculates and returns the second Piola-Kirchhoff stress.
\mathbf{S}(\mathbf{F}, T) = 2\mu\mathbf{E}' + \kappa\,\mathrm{tr}(\mathbf{E})\mathbf{1} - 3\alpha\kappa(T - T_\mathrm{ref})\mathbf{1}
Source§fn second_piola_kirchhoff_tangent_stiffness(
&self,
deformation_gradient: &DeformationGradient,
_: &Scalar,
) -> Result<SecondPiolaKirchhoffTangentStiffness, ConstitutiveError>
fn second_piola_kirchhoff_tangent_stiffness( &self, deformation_gradient: &DeformationGradient, _: &Scalar, ) -> Result<SecondPiolaKirchhoffTangentStiffness, ConstitutiveError>
Calculates and returns the tangent stiffness associated with the second Piola-Kirchhoff stress.
\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 coefficient_of_thermal_expansion(&self) -> &Scalar
fn coefficient_of_thermal_expansion(&self) -> &Scalar
Returns the coefficient of thermal expansion.
Source§fn reference_temperature(&self) -> &Scalar
fn reference_temperature(&self) -> &Scalar
Returns the reference temperature.
Source§fn cauchy_stress(
&self,
deformation_gradient: &DeformationGradient,
temperature: &Scalar,
) -> Result<CauchyStress, ConstitutiveError>
fn cauchy_stress( &self, deformation_gradient: &DeformationGradient, temperature: &Scalar, ) -> Result<CauchyStress, ConstitutiveError>
Calculates and returns the Cauchy stress. Read more
Source§fn cauchy_tangent_stiffness(
&self,
deformation_gradient: &DeformationGradient,
temperature: &Scalar,
) -> Result<CauchyTangentStiffness, ConstitutiveError>
fn cauchy_tangent_stiffness( &self, deformation_gradient: &DeformationGradient, temperature: &Scalar, ) -> 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,
temperature: &Scalar,
) -> Result<FirstPiolaKirchhoffStress, ConstitutiveError>
fn first_piola_kirchhoff_stress( &self, deformation_gradient: &DeformationGradient, temperature: &Scalar, ) -> Result<FirstPiolaKirchhoffStress, ConstitutiveError>
Calculates and returns the first Piola-Kirchhoff stress. Read more
Source§fn first_piola_kirchhoff_tangent_stiffness(
&self,
deformation_gradient: &DeformationGradient,
temperature: &Scalar,
) -> Result<FirstPiolaKirchhoffTangentStiffness, ConstitutiveError>
fn first_piola_kirchhoff_tangent_stiffness( &self, deformation_gradient: &DeformationGradient, temperature: &Scalar, ) -> Result<FirstPiolaKirchhoffTangentStiffness, ConstitutiveError>
Calculates and returns the tangent stiffness associated with the first Piola-Kirchhoff stress. Read more
Source§impl<'a> Thermohyperelastic<'a> for SaintVenantKirchhoff<'a>
impl<'a> Thermohyperelastic<'a> for SaintVenantKirchhoff<'a>
Source§fn helmholtz_free_energy_density(
&self,
deformation_gradient: &DeformationGradient,
temperature: &Scalar,
) -> Result<Scalar, ConstitutiveError>
fn helmholtz_free_energy_density( &self, deformation_gradient: &DeformationGradient, temperature: &Scalar, ) -> Result<Scalar, ConstitutiveError>
Calculates and returns the Helmholtz free energy density.
a(\mathbf{F}, T) = \mu\,\mathrm{tr}(\mathbf{E}^2) + \frac{1}{2}\left(\kappa - \frac{2}{3}\,\mu\right)\mathrm{tr}(\mathbf{E})^2 - 3\alpha\kappa\,\mathrm{tr}(\mathbf{E})(T - T_\mathrm{ref})
Auto Trait Implementations§
impl<'a> Freeze for SaintVenantKirchhoff<'a>
impl<'a> RefUnwindSafe for SaintVenantKirchhoff<'a>
impl<'a> Send for SaintVenantKirchhoff<'a>
impl<'a> Sync for SaintVenantKirchhoff<'a>
impl<'a> Unpin for SaintVenantKirchhoff<'a>
impl<'a> UnwindSafe for SaintVenantKirchhoff<'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