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