conspire/constitutive/solid/hyperelastic_viscoplastic/hencky/
mod.rs1#[cfg(test)]
2mod test;
3
4use crate::{
5 constitutive::{
6 ConstitutiveError,
7 fluid::{plastic::Plastic, viscoplastic::Viscoplastic},
8 solid::{
9 Solid, TWO_THIRDS, elastic_viscoplastic::ElasticViscoplastic,
10 hyperelastic_viscoplastic::HyperelasticViscoplastic,
11 },
12 },
13 math::{ContractThirdFourthIndicesWithFirstSecondIndicesOf, IDENTITY, Rank2},
14 mechanics::{
15 CauchyStress, CauchyTangentStiffness, CauchyTangentStiffnessElastic, Deformation,
16 DeformationGradient, DeformationGradientPlastic, Scalar,
17 },
18};
19
20#[doc = include_str!("doc.md")]
21#[derive(Clone, Debug)]
22pub struct Hencky {
23 pub bulk_modulus: Scalar,
25 pub shear_modulus: Scalar,
27 pub initial_yield_stress: Scalar,
29 pub hardening_slope: Scalar,
31 pub rate_sensitivity: Scalar,
33 pub reference_flow_rate: Scalar,
35}
36
37impl Solid for Hencky {
38 fn bulk_modulus(&self) -> Scalar {
39 self.bulk_modulus
40 }
41 fn shear_modulus(&self) -> Scalar {
42 self.shear_modulus
43 }
44}
45
46impl Plastic for Hencky {
47 fn initial_yield_stress(&self) -> Scalar {
48 self.initial_yield_stress
49 }
50 fn hardening_slope(&self) -> Scalar {
51 self.hardening_slope
52 }
53}
54
55impl Viscoplastic for Hencky {
56 fn rate_sensitivity(&self) -> Scalar {
57 self.rate_sensitivity
58 }
59 fn reference_flow_rate(&self) -> Scalar {
60 self.reference_flow_rate
61 }
62}
63
64impl ElasticViscoplastic for Hencky {
65 #[doc = include_str!("cauchy_stress.md")]
66 fn cauchy_stress(
67 &self,
68 deformation_gradient: &DeformationGradient,
69 deformation_gradient_p: &DeformationGradientPlastic,
70 ) -> Result<CauchyStress, ConstitutiveError> {
71 let jacobian = self.jacobian(deformation_gradient)?;
72 let deformation_gradient_e = deformation_gradient * deformation_gradient_p.inverse();
73 let (deviatoric_strain_e, strain_trace_e) =
74 (deformation_gradient_e.left_cauchy_green().logm()? * 0.5).deviatoric_and_trace();
75 Ok(
76 deviatoric_strain_e * (2.0 * self.shear_modulus() / jacobian)
77 + IDENTITY * (self.bulk_modulus() * strain_trace_e / jacobian),
78 )
79 }
80 #[doc = include_str!("cauchy_tangent_stiffness.md")]
81 fn cauchy_tangent_stiffness(
82 &self,
83 deformation_gradient: &DeformationGradient,
84 deformation_gradient_p: &DeformationGradientPlastic,
85 ) -> Result<CauchyTangentStiffness, ConstitutiveError> {
86 let jacobian = self.jacobian(deformation_gradient)?;
87 let deformation_gradient_inverse_p = deformation_gradient_p.inverse();
88 let deformation_gradient_e = deformation_gradient * &deformation_gradient_inverse_p;
89 let left_cauchy_green_e = deformation_gradient_e.left_cauchy_green();
90 let (deviatoric_strain_e, strain_trace_e) =
91 (left_cauchy_green_e.logm()? * 0.5).deviatoric_and_trace();
92 let scaled_deformation_gradient_e =
93 &deformation_gradient_e * self.shear_modulus() / jacobian;
94 Ok((left_cauchy_green_e
95 .dlogm()?
96 .contract_third_fourth_indices_with_first_second_indices_of(
97 &(CauchyTangentStiffnessElastic::dyad_il_jk(
98 &scaled_deformation_gradient_e,
99 &IDENTITY,
100 ) + CauchyTangentStiffnessElastic::dyad_ik_jl(
101 &IDENTITY,
102 &scaled_deformation_gradient_e,
103 )),
104 ))
105 * deformation_gradient_inverse_p.transpose()
106 + (CauchyTangentStiffness::dyad_ij_kl(
107 &(IDENTITY
108 * ((self.bulk_modulus() - TWO_THIRDS * self.shear_modulus()) / jacobian)
109 - deviatoric_strain_e * (2.0 * self.shear_modulus() / jacobian)
110 - IDENTITY * (self.bulk_modulus() * strain_trace_e / jacobian)),
111 &deformation_gradient.inverse_transpose(),
112 )))
113 }
114}
115
116impl HyperelasticViscoplastic for Hencky {
117 #[doc = include_str!("helmholtz_free_energy_density.md")]
118 fn helmholtz_free_energy_density(
119 &self,
120 deformation_gradient: &DeformationGradient,
121 deformation_gradient_p: &DeformationGradientPlastic,
122 ) -> Result<Scalar, ConstitutiveError> {
123 let _jacobian = self.jacobian(deformation_gradient)?;
124 let deformation_gradient_e = deformation_gradient * deformation_gradient_p.inverse();
125 let strain_e = deformation_gradient_e.left_cauchy_green().logm()? * 0.5;
126 Ok(self.shear_modulus() * strain_e.squared_trace()
127 + 0.5
128 * (self.bulk_modulus() - TWO_THIRDS * self.shear_modulus())
129 * strain_e.trace().powi(2))
130 }
131}