1*bcb2dfaeSJed Brown(example-petsc-elasticity)= 2*bcb2dfaeSJed Brown 3*bcb2dfaeSJed Brown# Solid mechanics mini-app 4*bcb2dfaeSJed Brown 5*bcb2dfaeSJed BrownThis example is located in the subdirectory {file}`examples/solids`. 6*bcb2dfaeSJed BrownIt solves the steady-state static momentum balance equations using unstructured high-order finite/spectral element spatial discretizations. 7*bcb2dfaeSJed BrownAs for the {ref}`example-petsc-navier-stokes` case, the solid mechanics elasticity example has been developed using PETSc, so that the pointwise physics (defined at quadrature points) is separated from the parallelization and meshing concerns. 8*bcb2dfaeSJed Brown 9*bcb2dfaeSJed BrownIn this mini-app, we consider three formulations used in solid mechanics applications: linear elasticity, Neo-Hookean hyperelasticity at small strain, and Neo-Hookean hyperelasticity at finite strain. 10*bcb2dfaeSJed BrownWe provide the strong and weak forms of static balance of linear momentum in the small strain and finite strain regimes. 11*bcb2dfaeSJed BrownThe stress-strain relationship (constitutive law) for each of the material models is provided. 12*bcb2dfaeSJed BrownDue to the nonlinearity of material models in Neo-Hookean hyperelasticity, the Newton linearization of the material models is provided. 13*bcb2dfaeSJed Brown 14*bcb2dfaeSJed Brown:::{note} 15*bcb2dfaeSJed BrownLinear elasticity and small-strain hyperelasticity can both by obtained from the finite-strain hyperelastic formulation by linearization of geometric and constitutive nonlinearities. 16*bcb2dfaeSJed BrownThe effect of these linearizations is sketched in the diagram below, where $\bm \sigma$ and $\bm \epsilon$ are stress and strain, respectively, in the small strain regime, while $\bm S$ and $\bm E$ are their finite-strain generalizations (second Piola-Kirchoff tensor and Green-Lagrange strain tensor, respectively) defined in the initial configuration, and $\mathsf C$ is a linearized constitutive model. 17*bcb2dfaeSJed Brown 18*bcb2dfaeSJed Brown$$ 19*bcb2dfaeSJed Brown\begin{CD} 20*bcb2dfaeSJed Brown {\overbrace{\bm S(\bm E)}^{\text{Finite Strain Hyperelastic}}} 21*bcb2dfaeSJed Brown @>{\text{constitutive}}>{\text{linearization}}> 22*bcb2dfaeSJed Brown {\overbrace{\bm S = \mathsf C \bm E}^{\text{St. Venant-Kirchoff}}} \\ 23*bcb2dfaeSJed Brown @V{\text{geometric}}V{\begin{smallmatrix}\bm E \to \bm \epsilon \\ \bm S \to \bm \sigma \end{smallmatrix}}V 24*bcb2dfaeSJed Brown @V{\begin{smallmatrix}\bm E \to \bm \epsilon \\ \bm S \to \bm \sigma \end{smallmatrix}}V{\text{geometric}}V \\ 25*bcb2dfaeSJed Brown {\underbrace{\bm \sigma(\bm \epsilon)}_\text{Small Strain Hyperelastic}} 26*bcb2dfaeSJed Brown @>{\text{constitutive}}>\text{linearization}> 27*bcb2dfaeSJed Brown {\underbrace{\bm \sigma = \mathsf C \bm \epsilon}_\text{Linear Elastic}} 28*bcb2dfaeSJed Brown\end{CD} 29*bcb2dfaeSJed Brown$$ (hyperelastic-cd) 30*bcb2dfaeSJed Brown::: 31*bcb2dfaeSJed Brown 32*bcb2dfaeSJed Brown(running-elasticity)= 33*bcb2dfaeSJed Brown 34*bcb2dfaeSJed Brown## Running the mini-app 35*bcb2dfaeSJed Brown 36*bcb2dfaeSJed Brown```{include} README.md 37*bcb2dfaeSJed Brown:start-after: inclusion-solids-marker 38*bcb2dfaeSJed Brown``` 39*bcb2dfaeSJed Brown 40*bcb2dfaeSJed Brown(problem-linear-elasticity)= 41*bcb2dfaeSJed Brown 42*bcb2dfaeSJed Brown## Linear Elasticity 43*bcb2dfaeSJed Brown 44*bcb2dfaeSJed BrownThe strong form of the static balance of linear momentum at small strain for the three-dimensional linear elasticity problem is given by {cite}`hughes2012finite`: 45*bcb2dfaeSJed Brown 46*bcb2dfaeSJed Brown$$ 47*bcb2dfaeSJed Brown\nabla \cdot \bm{\sigma} + \bm{g} = \bm{0} 48*bcb2dfaeSJed Brown$$ (lin-elas) 49*bcb2dfaeSJed Brown 50*bcb2dfaeSJed Brownwhere $\bm{\sigma}$ and $\bm{g}$ are stress and forcing functions, respectively. 51*bcb2dfaeSJed BrownWe multiply {math:numref}`lin-elas` by a test function $\bm v$ and integrate the divergence term by parts to arrive at the weak form: find $\bm u \in \mathcal V \subset H^1(\Omega)$ such that 52*bcb2dfaeSJed Brown 53*bcb2dfaeSJed Brown$$ 54*bcb2dfaeSJed Brown\int_{\Omega}{ \nabla \bm{v} \tcolon \bm{\sigma}} \, dV 55*bcb2dfaeSJed Brown- \int_{\partial \Omega}{\bm{v} \cdot \left(\bm{\sigma} \cdot \hat{\bm{n}}\right)} \, dS 56*bcb2dfaeSJed Brown- \int_{\Omega}{\bm{v} \cdot \bm{g}} \, dV 57*bcb2dfaeSJed Brown= 0, \quad \forall \bm v \in \mathcal V, 58*bcb2dfaeSJed Brown$$ (lin-elas-weak) 59*bcb2dfaeSJed Brown 60*bcb2dfaeSJed Brownwhere $\bm{\sigma} \cdot \hat{\bm{n}}|_{\partial \Omega}$ is replaced by an applied force/traction boundary condition written in terms of the initial configuration. 61*bcb2dfaeSJed BrownWhen inhomogeneous Dirichlet boundary conditions are present, $\mathcal V$ is an affine space that satisfies those boundary conditions. 62*bcb2dfaeSJed Brown 63*bcb2dfaeSJed Brown### Constitutive modeling 64*bcb2dfaeSJed Brown 65*bcb2dfaeSJed BrownIn their most general form, constitutive models define $\bm \sigma$ in terms of state variables. 66*bcb2dfaeSJed BrownIn the model taken into consideration in the present mini-app, the state variables are constituted by the vector displacement field $\bm u$, and its gradient $\nabla \bm u$. 67*bcb2dfaeSJed BrownWe begin by defining the symmetric (small/infintesimal) strain tensor as 68*bcb2dfaeSJed Brown 69*bcb2dfaeSJed Brown$$ 70*bcb2dfaeSJed Brown\bm{\epsilon} = \dfrac{1}{2}\left(\nabla \bm{u} + \nabla \bm{u}^T \right). 71*bcb2dfaeSJed Brown$$ (small-strain) 72*bcb2dfaeSJed Brown 73*bcb2dfaeSJed BrownThis constitutive model $\bm \sigma(\bm \epsilon)$ is a linear tensor-valued function of a tensor-valued input, but we will consider the more general nonlinear case in other models below. 74*bcb2dfaeSJed BrownIn these cases, an arbitrary choice of such a function will generally not be invariant under orthogonal transformations and thus will not admissible as a physical model must not depend on the coordinate system chosen to express it. 75*bcb2dfaeSJed BrownIn particular, given an orthogonal transformation $Q$, we desire 76*bcb2dfaeSJed Brown 77*bcb2dfaeSJed Brown$$ 78*bcb2dfaeSJed BrownQ \bm \sigma(\bm \epsilon) Q^T = \bm \sigma(Q \bm \epsilon Q^T), 79*bcb2dfaeSJed Brown$$ (elastic-invariance) 80*bcb2dfaeSJed Brown 81*bcb2dfaeSJed Brownwhich means that we can change our reference frame before or after computing $\bm \sigma$, and get the same result either way. 82*bcb2dfaeSJed BrownConstitutive relations in which $\bm \sigma$ is uniquely determined by $\bm \epsilon$ while satisfying the invariance property {math:numref}`elastic-invariance` are known as Cauchy elastic materials. 83*bcb2dfaeSJed BrownHere, we define a strain energy density functional $\Phi(\bm \epsilon) \in \mathbb R$ and obtain the strain energy from its gradient, 84*bcb2dfaeSJed Brown 85*bcb2dfaeSJed Brown$$ 86*bcb2dfaeSJed Brown\bm \sigma(\bm \epsilon) = \frac{\partial \Phi}{\partial \bm \epsilon}. 87*bcb2dfaeSJed Brown$$ (strain-energy-grad) 88*bcb2dfaeSJed Brown 89*bcb2dfaeSJed Brown:::{note} 90*bcb2dfaeSJed BrownThe strain energy density functional cannot be an arbitrary function $\Phi(\bm \epsilon)$; it can only depend on *invariants*, scalar-valued functions $\gamma$ satisfying 91*bcb2dfaeSJed Brown 92*bcb2dfaeSJed Brown$$ 93*bcb2dfaeSJed Brown\gamma(\bm \epsilon) = \gamma(Q \bm \epsilon Q^T) 94*bcb2dfaeSJed Brown$$ 95*bcb2dfaeSJed Brown 96*bcb2dfaeSJed Brownfor all orthogonal matrices $Q$. 97*bcb2dfaeSJed Brown::: 98*bcb2dfaeSJed Brown 99*bcb2dfaeSJed BrownFor the linear elasticity model, the strain energy density is given by 100*bcb2dfaeSJed Brown 101*bcb2dfaeSJed Brown$$ 102*bcb2dfaeSJed Brown\bm{\Phi} = \frac{\lambda}{2} (\operatorname{trace} \bm{\epsilon})^2 + \mu \bm{\epsilon} : \bm{\epsilon} . 103*bcb2dfaeSJed Brown$$ 104*bcb2dfaeSJed Brown 105*bcb2dfaeSJed BrownThe constitutive law (stress-strain relationship) is therefore given by its gradient, 106*bcb2dfaeSJed Brown 107*bcb2dfaeSJed Brown$$ 108*bcb2dfaeSJed Brown\bm\sigma = \lambda (\operatorname{trace} \bm\epsilon) \bm I_3 + 2 \mu \bm\epsilon, 109*bcb2dfaeSJed Brown$$ 110*bcb2dfaeSJed Brown 111*bcb2dfaeSJed Brownwhere $\bm I_3$ is the $3 \times 3$ identity matrix, the colon represents a double contraction (over both indices of $\bm \epsilon$), and the Lamé parameters are given by 112*bcb2dfaeSJed Brown 113*bcb2dfaeSJed Brown$$ 114*bcb2dfaeSJed Brown\begin{aligned} \lambda &= \frac{E \nu}{(1 + \nu)(1 - 2 \nu)} \\ \mu &= \frac{E}{2(1 + \nu)} \end{aligned}. 115*bcb2dfaeSJed Brown$$ 116*bcb2dfaeSJed Brown 117*bcb2dfaeSJed BrownThe constitutive law (stress-strain relationship) can also be written as 118*bcb2dfaeSJed Brown 119*bcb2dfaeSJed Brown$$ 120*bcb2dfaeSJed Brown\bm{\sigma} = \mathsf{C} \!:\! \bm{\epsilon}. 121*bcb2dfaeSJed Brown$$ (linear-stress-strain) 122*bcb2dfaeSJed Brown 123*bcb2dfaeSJed BrownFor notational convenience, we express the symmetric second order tensors $\bm \sigma$ and $\bm \epsilon$ as vectors of length 6 using the [Voigt notation](https://en.wikipedia.org/wiki/Voigt_notation). 124*bcb2dfaeSJed BrownHence, the fourth order elasticity tensor $\mathsf C$ (also known as elastic moduli tensor or material stiffness tensor) can be represented as 125*bcb2dfaeSJed Brown 126*bcb2dfaeSJed Brown$$ 127*bcb2dfaeSJed Brown\mathsf C = \begin{pmatrix} 128*bcb2dfaeSJed Brown\lambda + 2\mu & \lambda & \lambda & & & \\ 129*bcb2dfaeSJed Brown\lambda & \lambda + 2\mu & \lambda & & & \\ 130*bcb2dfaeSJed Brown\lambda & \lambda & \lambda + 2\mu & & & \\ 131*bcb2dfaeSJed Brown& & & \mu & & \\ 132*bcb2dfaeSJed Brown& & & & \mu & \\ 133*bcb2dfaeSJed Brown& & & & & \mu 134*bcb2dfaeSJed Brown\end{pmatrix}. 135*bcb2dfaeSJed Brown$$ (linear-elasticity-tensor) 136*bcb2dfaeSJed Brown 137*bcb2dfaeSJed BrownNote that the incompressible limit $\nu \to \frac 1 2$ causes $\lambda \to \infty$, and thus $\mathsf C$ becomes singular. 138*bcb2dfaeSJed Brown 139*bcb2dfaeSJed Brown(problem-hyper-small-strain)= 140*bcb2dfaeSJed Brown 141*bcb2dfaeSJed Brown## Hyperelasticity at Small Strain 142*bcb2dfaeSJed Brown 143*bcb2dfaeSJed BrownThe strong and weak forms given above, in {math:numref}`lin-elas` and {math:numref}`lin-elas-weak`, are valid for Neo-Hookean hyperelasticity at small strain. 144*bcb2dfaeSJed BrownHowever, the strain energy density differs and is given by 145*bcb2dfaeSJed Brown 146*bcb2dfaeSJed Brown$$ 147*bcb2dfaeSJed Brown\bm{\Phi} = \lambda (1 + \operatorname{trace} \bm{\epsilon}) (\log(1 + \operatorname{trace} \bm\epsilon) - 1) + \mu \bm{\epsilon} : \bm{\epsilon} . 148*bcb2dfaeSJed Brown$$ 149*bcb2dfaeSJed Brown 150*bcb2dfaeSJed BrownAs above, we have the corresponding constitutive law given by 151*bcb2dfaeSJed Brown 152*bcb2dfaeSJed Brown$$ 153*bcb2dfaeSJed Brown\bm{\sigma} = \lambda \log(1 + \operatorname{trace} \bm\epsilon) \bm{I}_3 + 2\mu \bm{\epsilon} 154*bcb2dfaeSJed Brown$$ (eq-neo-hookean-small-strain) 155*bcb2dfaeSJed Brown 156*bcb2dfaeSJed Brownwhere $\bm{\epsilon}$ is defined as in {math:numref}`small-strain`. 157*bcb2dfaeSJed Brown 158*bcb2dfaeSJed Brown### Newton linearization 159*bcb2dfaeSJed Brown 160*bcb2dfaeSJed BrownDue to nonlinearity in the constitutive law, we require a Newton linearization of {math:numref}`eq-neo-hookean-small-strain`. 161*bcb2dfaeSJed BrownTo derive the Newton linearization, we begin by expressing the derivative, 162*bcb2dfaeSJed Brown 163*bcb2dfaeSJed Brown$$ 164*bcb2dfaeSJed Brown\diff \bm{\sigma} = \dfrac{\partial \bm{\sigma}}{\partial \bm{\epsilon}} \tcolon \diff \bm{\epsilon} 165*bcb2dfaeSJed Brown$$ 166*bcb2dfaeSJed Brown 167*bcb2dfaeSJed Brownwhere 168*bcb2dfaeSJed Brown 169*bcb2dfaeSJed Brown$$ 170*bcb2dfaeSJed Brown\diff \bm{\epsilon} = \dfrac{1}{2}\left( \nabla \diff \bm{u} + \nabla \diff \bm{u}^T \right) 171*bcb2dfaeSJed Brown$$ 172*bcb2dfaeSJed Brown 173*bcb2dfaeSJed Brownand 174*bcb2dfaeSJed Brown 175*bcb2dfaeSJed Brown$$ 176*bcb2dfaeSJed Brown\diff \nabla \bm{u} = \nabla \diff \bm{u} . 177*bcb2dfaeSJed Brown$$ 178*bcb2dfaeSJed Brown 179*bcb2dfaeSJed BrownTherefore, 180*bcb2dfaeSJed Brown 181*bcb2dfaeSJed Brown$$ 182*bcb2dfaeSJed Brown\diff \bm{\sigma} = \bar{\lambda} \cdot \operatorname{trace} \diff \bm{\epsilon} \cdot \bm{I}_3 + 2\mu \diff \bm{\epsilon} 183*bcb2dfaeSJed Brown$$ (derss) 184*bcb2dfaeSJed Brown 185*bcb2dfaeSJed Brownwhere we have introduced the symbol 186*bcb2dfaeSJed Brown 187*bcb2dfaeSJed Brown$$ 188*bcb2dfaeSJed Brown\bar{\lambda} = \dfrac{\lambda}{1 + \epsilon_v } 189*bcb2dfaeSJed Brown$$ 190*bcb2dfaeSJed Brown 191*bcb2dfaeSJed Brownwhere volumetric strain is given by $\epsilon_v = \sum_i \epsilon_{ii}$. 192*bcb2dfaeSJed Brown 193*bcb2dfaeSJed BrownEquation {math:numref}`derss` can be written in Voigt matrix notation as follows: 194*bcb2dfaeSJed Brown 195*bcb2dfaeSJed Brown$$ 196*bcb2dfaeSJed Brown\begin{pmatrix} 197*bcb2dfaeSJed Brown \diff \sigma_{11} \\ 198*bcb2dfaeSJed Brown \diff \sigma_{22} \\ 199*bcb2dfaeSJed Brown \diff \sigma_{33} \\ 200*bcb2dfaeSJed Brown \diff \sigma_{23} \\ 201*bcb2dfaeSJed Brown \diff \sigma_{13} \\ 202*bcb2dfaeSJed Brown \diff \sigma_{12} 203*bcb2dfaeSJed Brown\end{pmatrix} = 204*bcb2dfaeSJed Brown\begin{pmatrix} 205*bcb2dfaeSJed Brown 2 \mu +\bar{\lambda} & \bar{\lambda} & \bar{\lambda} & & & \\ 206*bcb2dfaeSJed Brown \bar{\lambda} & 2 \mu +\bar{\lambda} & \bar{\lambda} & & & \\ 207*bcb2dfaeSJed Brown \bar{\lambda} & \bar{\lambda} & 2 \mu +\bar{\lambda} & & & \\ 208*bcb2dfaeSJed Brown & & & \mu & & \\ 209*bcb2dfaeSJed Brown & & & & \mu & \\ 210*bcb2dfaeSJed Brown & & & & & \mu \\ 211*bcb2dfaeSJed Brown\end{pmatrix} 212*bcb2dfaeSJed Brown\begin{pmatrix} 213*bcb2dfaeSJed Brown \diff \epsilon_{11} \\ 214*bcb2dfaeSJed Brown \diff \epsilon_{22} \\ 215*bcb2dfaeSJed Brown \diff \epsilon_{33} \\ 216*bcb2dfaeSJed Brown 2 \diff \epsilon_{23} \\ 217*bcb2dfaeSJed Brown 2 \diff \epsilon_{13} \\ 218*bcb2dfaeSJed Brown 2 \diff \epsilon_{12} 219*bcb2dfaeSJed Brown\end{pmatrix}. 220*bcb2dfaeSJed Brown$$ (mdss) 221*bcb2dfaeSJed Brown 222*bcb2dfaeSJed Brown(problem-hyperelasticity-finite-strain)= 223*bcb2dfaeSJed Brown 224*bcb2dfaeSJed Brown## Hyperelasticity at Finite Strain 225*bcb2dfaeSJed Brown 226*bcb2dfaeSJed BrownIn the *total Lagrangian* approach for the Neo-Hookean hyperelasticity problem, the discrete equations are formulated with respect to the initial configuration. 227*bcb2dfaeSJed BrownIn this formulation, we solve for displacement $\bm u(\bm X)$ in the reference frame $\bm X$. 228*bcb2dfaeSJed BrownThe notation for elasticity at finite strain is inspired by {cite}`holzapfel2000nonlinear` to distinguish between the current and initial configurations. 229*bcb2dfaeSJed BrownAs explained in the {ref}`common-notation` section, we denote by capital letters the reference frame and by small letters the current one. 230*bcb2dfaeSJed Brown 231*bcb2dfaeSJed BrownThe strong form of the static balance of linear-momentum at *finite strain* (total Lagrangian) is given by: 232*bcb2dfaeSJed Brown 233*bcb2dfaeSJed Brown$$ 234*bcb2dfaeSJed Brown- \nabla_X \cdot \bm{P} - \rho_0 \bm{g} = \bm{0} 235*bcb2dfaeSJed Brown$$ (sblFinS) 236*bcb2dfaeSJed Brown 237*bcb2dfaeSJed Brownwhere the $_X$ in $\nabla_X$ indicates that the gradient is calculated with respect to the initial configuration in the finite strain regime. 238*bcb2dfaeSJed Brown$\bm{P}$ and $\bm{g}$ are the *first Piola-Kirchhoff stress* tensor and the prescribed forcing function, respectively. 239*bcb2dfaeSJed Brown$\rho_0$ is known as the *initial* mass density. 240*bcb2dfaeSJed BrownThe tensor $\bm P$ is not symmetric, living in the current configuration on the left and the initial configuration on the right. 241*bcb2dfaeSJed Brown 242*bcb2dfaeSJed Brown$\bm{P}$ can be decomposed as 243*bcb2dfaeSJed Brown 244*bcb2dfaeSJed Brown$$ 245*bcb2dfaeSJed Brown\bm{P} = \bm{F} \, \bm{S}, 246*bcb2dfaeSJed Brown$$ (1st2nd) 247*bcb2dfaeSJed Brown 248*bcb2dfaeSJed Brownwhere $\bm S$ is the *second Piola-Kirchhoff stress* tensor, a symmetric tensor defined entirely in the initial configuration, and $\bm{F} = \bm I_3 + \nabla_X \bm u$ is the deformation gradient. 249*bcb2dfaeSJed BrownDifferent constitutive models can define $\bm S$. 250*bcb2dfaeSJed Brown 251*bcb2dfaeSJed Brown### Constitutive modeling 252*bcb2dfaeSJed Brown 253*bcb2dfaeSJed BrownFor the constitutive modeling of hyperelasticity at finite strain, we begin by defining two symmetric tensors in the initial configuration, the right Cauchy-Green tensor 254*bcb2dfaeSJed Brown 255*bcb2dfaeSJed Brown$$ 256*bcb2dfaeSJed Brown\bm C = \bm F^T \bm F 257*bcb2dfaeSJed Brown$$ 258*bcb2dfaeSJed Brown 259*bcb2dfaeSJed Brownand the Green-Lagrange strain tensor 260*bcb2dfaeSJed Brown 261*bcb2dfaeSJed Brown$$ 262*bcb2dfaeSJed Brown\bm E = \frac 1 2 (\bm C - \bm I_3) = \frac 1 2 \Big( \nabla_X \bm u + (\nabla_X \bm u)^T + (\nabla_X \bm u)^T \nabla_X \bm u \Big), 263*bcb2dfaeSJed Brown$$ (eq-green-lagrange-strain) 264*bcb2dfaeSJed Brown 265*bcb2dfaeSJed Brownthe latter of which converges to the linear strain tensor $\bm \epsilon$ in the small-deformation limit. 266*bcb2dfaeSJed BrownThe constitutive models considered, appropriate for large deformations, express $\bm S$ as a function of $\bm E$, similar to the linear case, shown in equation {math:numref}`linear-stress-strain`, which expresses the relationship between $\bm\sigma$ and $\bm\epsilon$. 267*bcb2dfaeSJed Brown 268*bcb2dfaeSJed BrownRecall that the strain energy density functional can only depend upon invariants. 269*bcb2dfaeSJed BrownWe will assume without loss of generality that $\bm E$ is diagonal and take its set of eigenvalues as the invariants. 270*bcb2dfaeSJed BrownIt is clear that there can be only three invariants, and there are many alternate choices, such as $\operatorname{trace}(\bm E), \operatorname{trace}(\bm E^2), \lvert \bm E \rvert$, and combinations thereof. 271*bcb2dfaeSJed BrownIt is common in the literature for invariants to be taken from $\bm C = \bm I_3 + 2 \bm E$ instead of $\bm E$. 272*bcb2dfaeSJed Brown 273*bcb2dfaeSJed BrownFor example, if we take the compressible Neo-Hookean model, 274*bcb2dfaeSJed Brown 275*bcb2dfaeSJed Brown$$ 276*bcb2dfaeSJed Brown\begin{aligned} 277*bcb2dfaeSJed Brown\Phi(\bm E) &= \frac{\lambda}{2}(\log J)^2 - \mu \log J + \frac \mu 2 (\operatorname{trace} \bm C - 3) \\ 278*bcb2dfaeSJed Brown &= \frac{\lambda}{2}(\log J)^2 - \mu \log J + \mu \operatorname{trace} \bm E, 279*bcb2dfaeSJed Brown\end{aligned} 280*bcb2dfaeSJed Brown$$ (neo-hookean-energy) 281*bcb2dfaeSJed Brown 282*bcb2dfaeSJed Brownwhere $J = \lvert \bm F \rvert = \sqrt{\lvert \bm C \rvert}$ is the determinant of deformation (i.e., volume change) and $\lambda$ and $\mu$ are the Lamé parameters in the infinitesimal strain limit. 283*bcb2dfaeSJed Brown 284*bcb2dfaeSJed BrownTo evaluate {math:numref}`strain-energy-grad`, we make use of 285*bcb2dfaeSJed Brown 286*bcb2dfaeSJed Brown$$ 287*bcb2dfaeSJed Brown\frac{\partial J}{\partial \bm E} = \frac{\partial \sqrt{\lvert \bm C \rvert}}{\partial \bm E} = \lvert \bm C \rvert^{-1/2} \lvert \bm C \rvert \bm C^{-1} = J \bm C^{-1}, 288*bcb2dfaeSJed Brown$$ 289*bcb2dfaeSJed Brown 290*bcb2dfaeSJed Brownwhere the factor of $\frac 1 2$ has been absorbed due to $\bm C = \bm I_3 + 2 \bm E.$ 291*bcb2dfaeSJed BrownCarrying through the differentiation {math:numref}`strain-energy-grad` for the model {math:numref}`neo-hookean-energy`, we arrive at 292*bcb2dfaeSJed Brown 293*bcb2dfaeSJed Brown$$ 294*bcb2dfaeSJed Brown\bm S = \lambda \log J \bm C^{-1} + \mu (\bm I_3 - \bm C^{-1}). 295*bcb2dfaeSJed Brown$$ (neo-hookean-stress) 296*bcb2dfaeSJed Brown 297*bcb2dfaeSJed Brown:::{tip} 298*bcb2dfaeSJed BrownAn equivalent form of {math:numref}`neo-hookean-stress` is 299*bcb2dfaeSJed Brown 300*bcb2dfaeSJed Brown$$ 301*bcb2dfaeSJed Brown\bm S = \lambda \log J \bm C^{-1} + 2 \mu \bm C^{-1} \bm E, 302*bcb2dfaeSJed Brown$$ (neo-hookean-stress-stable) 303*bcb2dfaeSJed Brown 304*bcb2dfaeSJed Brownwhich is more numerically stable for small $\bm E$, and thus preferred for computation. 305*bcb2dfaeSJed BrownNote that the product $\bm C^{-1} \bm E$ is also symmetric, and that $\bm E$ should be computed using {math:numref}`eq-green-lagrange-strain`. 306*bcb2dfaeSJed Brown 307*bcb2dfaeSJed BrownSimilarly, it is preferable to compute $\log J$ using `log1p`, especially in case of nearly incompressible materials. 308*bcb2dfaeSJed BrownTo sketch this idea, suppose we have the $2\times 2$ non-symmetric matrix $\bm{F} = \left( \begin{smallmatrix} 1 + u_{0,0} & u_{0,1} \\ u_{1,0} & 1 + u_{1,1} \end{smallmatrix} \right)$. 309*bcb2dfaeSJed BrownThen we compute 310*bcb2dfaeSJed Brown 311*bcb2dfaeSJed Brown$$ 312*bcb2dfaeSJed Brown\log J = \mathtt{log1p}(u_{0,0} + u_{1,1} + u_{0,0} u_{1,1} - u_{0,1} u_{1,0}), 313*bcb2dfaeSJed Brown$$ (log1p) 314*bcb2dfaeSJed Brown 315*bcb2dfaeSJed Brownwhich gives accurate results even in the limit when the entries $u_{i,j}$ are very small. 316*bcb2dfaeSJed BrownFor example, if $u_{i,j} \sim 10^{-8}$, then naive computation of $\bm I_3 - \bm C^{-1}$ and $\log J$ will have a relative accuracy of order $10^{-8}$ in double precision and no correct digits in single precision. 317*bcb2dfaeSJed BrownWhen using the stable choices above, these quantities retain full $\varepsilon_{\text{machine}}$ relative accuracy. 318*bcb2dfaeSJed Brown::: 319*bcb2dfaeSJed Brown 320*bcb2dfaeSJed Brown:::{dropdown} Mooney-Rivlin model 321*bcb2dfaeSJed BrownWhile the Neo-Hookean model depends on just two scalar invariants, $\mathbb I_1 = \trace \bm C = 3 + 2\trace \bm E$ and $J$, Mooney-Rivlin models depend on the additional invariant, $\mathbb I_2 = \frac 1 2 (\mathbb I_1^2 - \bm C \tcolon \bm C)$. 322*bcb2dfaeSJed BrownA coupled Mooney-Rivlin strain energy density (cf. Neo-Hookean {math:numref}`neo-hookean-energy`) is {cite}`holzapfel2000nonlinear` 323*bcb2dfaeSJed Brown 324*bcb2dfaeSJed Brown$$ 325*bcb2dfaeSJed Brown\Phi(\mathbb{I_1}, \mathbb{I_2}, J) = \frac{\lambda}{2}(\log J)^2 - (\mu_1 + 2\mu_2) \log J + \frac{\mu_1}{2}(\mathbb{I_1} - 3) + \frac{\mu_2}{2}(\mathbb{I_2} - 3). 326*bcb2dfaeSJed Brown$$ (mooney-rivlin-energy_coupled) 327*bcb2dfaeSJed Brown 328*bcb2dfaeSJed BrownWe differentiate $\Phi$ as in the Neo-Hookean case {math:numref}`neo-hookean-stress` to yield the second Piola-Kirchoff tensor, 329*bcb2dfaeSJed Brown 330*bcb2dfaeSJed Brown$$ 331*bcb2dfaeSJed Brown\begin{aligned} 332*bcb2dfaeSJed Brown\bm S &= \lambda \log J \bm{C}^{-1} - (\mu_1 + 2\mu_2) \bm{C}^{-1} + \mu_1\bm I_3 + \mu_2(\mathbb{I_1} \bm I_3 - \bm C) \\ 333*bcb2dfaeSJed Brown&= (\lambda \log J - \mu_1 - 2\mu_2) \bm C^{-1} + (\mu_1 + \mu_2 \mathbb I_1) \bm I_3 - \mu_2 \bm C, 334*bcb2dfaeSJed Brown\end{aligned} 335*bcb2dfaeSJed Brown$$ (mooney-rivlin-stress_coupled) 336*bcb2dfaeSJed Brown 337*bcb2dfaeSJed Brownwhere we have used 338*bcb2dfaeSJed Brown 339*bcb2dfaeSJed Brown$$ 340*bcb2dfaeSJed Brown\begin{aligned} 341*bcb2dfaeSJed Brown\frac{\partial \mathbb{I_1}}{\partial \bm E} &= 2 \bm I_3, & \frac{\partial \mathbb{I_2}}{\partial \bm E} &= 2 \mathbb I_1 \bm I_3 - 2 \bm C, & \frac{\partial \log J}{\partial \bm E} &= \bm{C}^{-1}. 342*bcb2dfaeSJed Brown\end{aligned} 343*bcb2dfaeSJed Brown$$ (None) 344*bcb2dfaeSJed Brown 345*bcb2dfaeSJed BrownThis is a common model for vulcanized rubber, with a shear modulus (defined for the small-strain limit) of $\mu_1 + \mu_2$ that should be significantly smaller than the first Lamé parameter $\lambda$. 346*bcb2dfaeSJed Brown::: 347*bcb2dfaeSJed Brown 348*bcb2dfaeSJed Brown:::{dropdown} Mooney-Rivlin strain energy comparison 349*bcb2dfaeSJed BrownWe apply traction to a block and plot integrated strain energy $\Phi$ as a function of the loading paramater. 350*bcb2dfaeSJed Brown 351*bcb2dfaeSJed Brown```{eval-rst} 352*bcb2dfaeSJed Brown.. altair-plot:: 353*bcb2dfaeSJed Brown :hide-code: 354*bcb2dfaeSJed Brown 355*bcb2dfaeSJed Brown import altair as alt 356*bcb2dfaeSJed Brown import pandas as pd 357*bcb2dfaeSJed Brown nh = pd.read_csv("source/examples/solids/output/NH-strain.csv") 358*bcb2dfaeSJed Brown nh["model"] = "Neo-Hookean" 359*bcb2dfaeSJed Brown nh["parameters"] = "E=2.8, nu=0.4" 360*bcb2dfaeSJed Brown 361*bcb2dfaeSJed Brown mr = pd.read_csv("source/examples/solids/output/MR-strain.csv") 362*bcb2dfaeSJed Brown mr["model"] = "Mooney-Rivlin; Neo-Hookean equivalent" 363*bcb2dfaeSJed Brown mr["parameters"] = "mu_1=1, mu_2=0, nu=.4" 364*bcb2dfaeSJed Brown 365*bcb2dfaeSJed Brown mr1 = pd.read_csv("source/examples/solids/output/MR-strain1.csv") 366*bcb2dfaeSJed Brown mr1["model"] = "Mooney-Rivlin" 367*bcb2dfaeSJed Brown mr1["parameters"] = "mu_1=0.5, mu_2=0.5, nu=.4" 368*bcb2dfaeSJed Brown 369*bcb2dfaeSJed Brown df = pd.concat([nh, mr, mr1]) 370*bcb2dfaeSJed Brown highlight = alt.selection_single( 371*bcb2dfaeSJed Brown on = "mouseover", 372*bcb2dfaeSJed Brown nearest = True, 373*bcb2dfaeSJed Brown fields=["model", "parameters"], 374*bcb2dfaeSJed Brown ) 375*bcb2dfaeSJed Brown base = alt.Chart(df).encode( 376*bcb2dfaeSJed Brown alt.X("increment"), 377*bcb2dfaeSJed Brown alt.Y("energy", scale=alt.Scale(type="sqrt")), 378*bcb2dfaeSJed Brown alt.Color("model"), 379*bcb2dfaeSJed Brown alt.Tooltip(("model", "parameters")), 380*bcb2dfaeSJed Brown opacity=alt.condition(highlight, alt.value(1), alt.value(.5)), 381*bcb2dfaeSJed Brown size=alt.condition(highlight, alt.value(2), alt.value(1)), 382*bcb2dfaeSJed Brown ) 383*bcb2dfaeSJed Brown base.mark_point().add_selection(highlight) + base.mark_line() 384*bcb2dfaeSJed Brown``` 385*bcb2dfaeSJed Brown::: 386*bcb2dfaeSJed Brown 387*bcb2dfaeSJed Brown:::{note} 388*bcb2dfaeSJed BrownOne can linearize {math:numref}`neo-hookean-stress` around $\bm E = 0$, for which $\bm C = \bm I_3 + 2 \bm E \to \bm I_3$ and $J \to 1 + \operatorname{trace} \bm E$, therefore {math:numref}`neo-hookean-stress` reduces to 389*bcb2dfaeSJed Brown 390*bcb2dfaeSJed Brown$$ 391*bcb2dfaeSJed Brown\bm S = \lambda (\trace \bm E) \bm I_3 + 2 \mu \bm E, 392*bcb2dfaeSJed Brown$$ (eq-st-venant-kirchoff) 393*bcb2dfaeSJed Brown 394*bcb2dfaeSJed Brownwhich is the St. Venant-Kirchoff model (constitutive linearization without geometric linearization; see {math:numref}`hyperelastic-cd`). 395*bcb2dfaeSJed Brown 396*bcb2dfaeSJed BrownThis model can be used for geometrically nonlinear mechanics (e.g., snap-through of thin structures), but is inappropriate for large strain. 397*bcb2dfaeSJed Brown 398*bcb2dfaeSJed BrownAlternatively, one can drop geometric nonlinearities, $\bm E \to \bm \epsilon$ and $\bm C \to \bm I_3$, while retaining the nonlinear dependence on $J \to 1 + \operatorname{trace} \bm \epsilon$, thereby yielding {math:numref}`eq-neo-hookean-small-strain` (see {math:numref}`hyperelastic-cd`). 399*bcb2dfaeSJed Brown::: 400*bcb2dfaeSJed Brown 401*bcb2dfaeSJed Brown### Weak form 402*bcb2dfaeSJed Brown 403*bcb2dfaeSJed BrownWe multiply {math:numref}`sblFinS` by a test function $\bm v$ and integrate by parts to obtain the weak form for finite-strain hyperelasticity: 404*bcb2dfaeSJed Brownfind $\bm u \in \mathcal V \subset H^1(\Omega_0)$ such that 405*bcb2dfaeSJed Brown 406*bcb2dfaeSJed Brown$$ 407*bcb2dfaeSJed Brown\int_{\Omega_0}{\nabla_X \bm{v} \tcolon \bm{P}} \, dV 408*bcb2dfaeSJed Brown - \int_{\Omega_0}{\bm{v} \cdot \rho_0 \bm{g}} \, dV 409*bcb2dfaeSJed Brown - \int_{\partial \Omega_0}{\bm{v} \cdot (\bm{P} \cdot \hat{\bm{N}})} \, dS 410*bcb2dfaeSJed Brown = 0, \quad \forall \bm v \in \mathcal V, 411*bcb2dfaeSJed Brown$$ (hyperelastic-weak-form-initial) 412*bcb2dfaeSJed Brown 413*bcb2dfaeSJed Brownwhere $\bm{P} \cdot \hat{\bm{N}}|_{\partial\Omega}$ is replaced by any prescribed force/traction boundary condition written in terms of the initial configuration. 414*bcb2dfaeSJed BrownThis equation contains material/constitutive nonlinearities in defining $\bm S(\bm E)$, as well as geometric nonlinearities through $\bm P = \bm F\, \bm S$, $\bm E(\bm F)$, and the body force $\bm g$, which must be pulled back from the current configuration to the initial configuration. 415*bcb2dfaeSJed BrownDiscretization of {math:numref}`hyperelastic-weak-form-initial` produces a finite-dimensional system of nonlinear algebraic equations, which we solve using Newton-Raphson methods. 416*bcb2dfaeSJed BrownOne attractive feature of Galerkin discretization is that we can arrive at the same linear system by discretizing the Newton linearization of the continuous form; that is, discretization and differentiation (Newton linearization) commute. 417*bcb2dfaeSJed Brown 418*bcb2dfaeSJed Brown### Newton linearization 419*bcb2dfaeSJed Brown 420*bcb2dfaeSJed BrownTo derive a Newton linearization of {math:numref}`hyperelastic-weak-form-initial`, we begin by expressing the derivative of {math:numref}`1st2nd` in incremental form, 421*bcb2dfaeSJed Brown 422*bcb2dfaeSJed Brown$$ 423*bcb2dfaeSJed Brown\diff \bm P = \frac{\partial \bm P}{\partial \bm F} \!:\! \diff \bm F = \diff \bm F\, \bm S + \bm F \underbrace{\frac{\partial \bm S}{\partial \bm E} \!:\! \diff \bm E}_{\diff \bm S} 424*bcb2dfaeSJed Brown$$ (eq-diff-P) 425*bcb2dfaeSJed Brown 426*bcb2dfaeSJed Brownwhere 427*bcb2dfaeSJed Brown 428*bcb2dfaeSJed Brown$$ 429*bcb2dfaeSJed Brown\diff \bm E = \frac{\partial \bm E}{\partial \bm F} \!:\! \diff \bm F = \frac 1 2 \Big( \diff \bm F^T \bm F + \bm F^T \diff \bm F \Big) 430*bcb2dfaeSJed Brown$$ 431*bcb2dfaeSJed Brown 432*bcb2dfaeSJed Brownand $\diff\bm F = \nabla_X\diff\bm u$. 433*bcb2dfaeSJed BrownThe quantity ${\partial \bm S} / {\partial \bm E}$ is known as the incremental elasticity tensor, and is analogous to the linear elasticity tensor $\mathsf C$ of {math:numref}`linear-elasticity-tensor`. 434*bcb2dfaeSJed BrownWe now evaluate $\diff \bm S$ for the Neo-Hookean model {math:numref}`neo-hookean-stress`, 435*bcb2dfaeSJed Brown 436*bcb2dfaeSJed Brown$$ 437*bcb2dfaeSJed Brown\diff\bm S = \frac{\partial \bm S}{\partial \bm E} \!:\! \diff \bm E 438*bcb2dfaeSJed Brown= \lambda (\bm C^{-1} \!:\! \diff\bm E) \bm C^{-1} 439*bcb2dfaeSJed Brown + 2 (\mu - \lambda \log J) \bm C^{-1} \diff\bm E \, \bm C^{-1}, 440*bcb2dfaeSJed Brown$$ (eq-neo-hookean-incremental-stress) 441*bcb2dfaeSJed Brown 442*bcb2dfaeSJed Brownwhere we have used 443*bcb2dfaeSJed Brown 444*bcb2dfaeSJed Brown$$ 445*bcb2dfaeSJed Brown\diff \bm C^{-1} = \frac{\partial \bm C^{-1}}{\partial \bm E} \!:\! \diff\bm E = -2 \bm C^{-1} \diff \bm E \, \bm C^{-1} . 446*bcb2dfaeSJed Brown$$ 447*bcb2dfaeSJed Brown 448*bcb2dfaeSJed Brown:::{note} 449*bcb2dfaeSJed BrownIn the small-strain limit, $\bm C \to \bm I_3$ and $\log J \to 0$, thereby reducing {math:numref}`eq-neo-hookean-incremental-stress` to the St. Venant-Kirchoff model {math:numref}`eq-st-venant-kirchoff`. 450*bcb2dfaeSJed Brown::: 451*bcb2dfaeSJed Brown 452*bcb2dfaeSJed Brown:::{dropdown} Newton linearization of Mooney-Rivlin 453*bcb2dfaeSJed BrownSimilar to {math:numref}`eq-neo-hookean-incremental-stress`, we differentiate {math:numref}`mooney-rivlin-stress_coupled` using variational notation, 454*bcb2dfaeSJed Brown 455*bcb2dfaeSJed Brown$$ 456*bcb2dfaeSJed Brown\begin{aligned} 457*bcb2dfaeSJed Brown\diff\bm S &= \lambda (\bm C^{-1} \tcolon \diff\bm E) \bm C^{-1} \\ 458*bcb2dfaeSJed Brown&\quad + 2(\mu_1 + 2\mu_2 - \lambda \log J) \bm C^{-1} \diff\bm E \bm C^{-1} \\ 459*bcb2dfaeSJed Brown&\quad + 2 \mu_2 \Big[ \trace (\diff\bm E) \bm I_3 - \diff\bm E\Big] . 460*bcb2dfaeSJed Brown\end{aligned} 461*bcb2dfaeSJed Brown$$ (mooney-rivlin-dS-coupled) 462*bcb2dfaeSJed Brown 463*bcb2dfaeSJed BrownNote that this agrees with {math:numref}`eq-neo-hookean-incremental-stress` if $\mu_1 = \mu, \mu_2 = 0$. 464*bcb2dfaeSJed BrownMoving from Neo-Hookean to Mooney-Rivlin modifies the second term and adds the third. 465*bcb2dfaeSJed Brown::: 466*bcb2dfaeSJed Brown 467*bcb2dfaeSJed Brown:::{dropdown} Cancellation vs symmetry 468*bcb2dfaeSJed BrownSome cancellation is possible (at the expense of symmetry) if we substitute {math:numref}`eq-neo-hookean-incremental-stress` into {math:numref}`eq-diff-P`, 469*bcb2dfaeSJed Brown 470*bcb2dfaeSJed Brown$$ 471*bcb2dfaeSJed Brown\begin{aligned} 472*bcb2dfaeSJed Brown\diff \bm P &= \diff \bm F\, \bm S 473*bcb2dfaeSJed Brown + \lambda (\bm C^{-1} : \diff \bm E) \bm F^{-T} + 2(\mu - \lambda \log J) \bm F^{-T} \diff\bm E \, \bm C^{-1} \\ 474*bcb2dfaeSJed Brown&= \diff \bm F\, \bm S 475*bcb2dfaeSJed Brown + \lambda (\bm F^{-T} : \diff \bm F) \bm F^{-T} + (\mu - \lambda \log J) \bm F^{-T} (\bm F^T \diff \bm F + \diff \bm F^T \bm F) \bm C^{-1} \\ 476*bcb2dfaeSJed Brown&= \diff \bm F\, \bm S 477*bcb2dfaeSJed Brown + \lambda (\bm F^{-T} : \diff \bm F) \bm F^{-T} + (\mu - \lambda \log J) \Big( \diff \bm F\, \bm C^{-1} + \bm F^{-T} \diff \bm F^T \bm F^{-T} \Big), 478*bcb2dfaeSJed Brown\end{aligned} 479*bcb2dfaeSJed Brown$$ (eq-diff-P-dF) 480*bcb2dfaeSJed Brown 481*bcb2dfaeSJed Brownwhere we have exploited $\bm F \bm C^{-1} = \bm F^{-T}$ and 482*bcb2dfaeSJed Brown 483*bcb2dfaeSJed Brown$$ 484*bcb2dfaeSJed Brown\begin{aligned} \bm C^{-1} \!:\! \diff \bm E = \bm C_{IJ}^{-1} \diff \bm E_{IJ} &= \frac 1 2 \bm F_{Ik}^{-1} \bm F_{Jk}^{-1} (\bm F_{\ell I} \diff \bm F_{\ell J} + \diff \bm F_{\ell I} \bm F_{\ell J}) \\ &= \frac 1 2 \Big( \delta_{\ell k} \bm F_{Jk}^{-1} \diff \bm F_{\ell J} + \delta_{\ell k} \bm F_{Ik}^{-1} \diff \bm F_{\ell I} \Big) \\ &= \bm F_{Ik}^{-1} \diff \bm F_{kI} = \bm F^{-T} \!:\! \diff \bm F. \end{aligned} 485*bcb2dfaeSJed Brown$$ 486*bcb2dfaeSJed Brown 487*bcb2dfaeSJed BrownWe prefer to compute with {math:numref}`eq-neo-hookean-incremental-stress` because {math:numref}`eq-diff-P-dF` is more expensive, requiring access to (non-symmetric) $\bm F^{-1}$ in addition to (symmetric) $\bm C^{-1} = \bm F^{-1} \bm F^{-T}$, having fewer symmetries to exploit in contractions, and being less numerically stable. 488*bcb2dfaeSJed Brown::: 489*bcb2dfaeSJed Brown 490*bcb2dfaeSJed Brown:::{dropdown} $\diff\bm S$ in index notation 491*bcb2dfaeSJed BrownIt is sometimes useful to express {math:numref}`eq-neo-hookean-incremental-stress` in index notation, 492*bcb2dfaeSJed Brown 493*bcb2dfaeSJed Brown$$ 494*bcb2dfaeSJed Brown\begin{aligned} 495*bcb2dfaeSJed Brown\diff\bm S_{IJ} &= \frac{\partial \bm S_{IJ}}{\partial \bm E_{KL}} \diff \bm E_{KL} \\ 496*bcb2dfaeSJed Brown &= \lambda (\bm C^{-1}_{KL} \diff\bm E_{KL}) \bm C^{-1}_{IJ} + 2 (\mu - \lambda \log J) \bm C^{-1}_{IK} \diff\bm E_{KL} \bm C^{-1}_{LJ} \\ 497*bcb2dfaeSJed Brown &= \underbrace{\Big( \lambda \bm C^{-1}_{IJ} \bm C^{-1}_{KL} + 2 (\mu - \lambda \log J) \bm C^{-1}_{IK} \bm C^{-1}_{JL} \Big)}_{\mathsf C_{IJKL}} \diff \bm E_{KL} \,, 498*bcb2dfaeSJed Brown\end{aligned} 499*bcb2dfaeSJed Brown$$ (eq-neo-hookean-incremental-stress-index) 500*bcb2dfaeSJed Brown 501*bcb2dfaeSJed Brownwhere we have identified the effective elasticity tensor $\mathsf C = \mathsf C_{IJKL}$. 502*bcb2dfaeSJed BrownIt is generally not desirable to store $\mathsf C$, but rather to use the earlier expressions so that only $3\times 3$ tensors (most of which are symmetric) must be manipulated. 503*bcb2dfaeSJed BrownThat is, given the linearization point $\bm F$ and solution increment $\diff \bm F = \nabla_X (\diff \bm u)$ (which we are solving for in the Newton step), we compute $\diff \bm P$ via 504*bcb2dfaeSJed Brown 505*bcb2dfaeSJed Brown1. recover $\bm C^{-1}$ and $\log J$ (either stored at quadrature points or recomputed), 506*bcb2dfaeSJed Brown2. proceed with $3\times 3$ matrix products as in {math:numref}`eq-neo-hookean-incremental-stress` or the second line of {math:numref}`eq-neo-hookean-incremental-stress-index` to compute $\diff \bm S$ while avoiding computation or storage of higher order tensors, and 507*bcb2dfaeSJed Brown3. conclude by {math:numref}`eq-diff-P`, where $\bm S$ is either stored or recomputed from its definition exactly as in the nonlinear residual evaluation. 508*bcb2dfaeSJed Brown::: 509*bcb2dfaeSJed Brown 510*bcb2dfaeSJed BrownNote that the Newton linearization of {math:numref}`hyperelastic-weak-form-initial` may be written as a weak form for linear operators: find $\diff\bm u \in \mathcal V_0$ such that 511*bcb2dfaeSJed Brown 512*bcb2dfaeSJed Brown$$ 513*bcb2dfaeSJed Brown\int_{\Omega_0} \nabla_X \bm v \!:\! \diff\bm P dV = \text{rhs}, \quad \forall \bm v \in \mathcal V_0, 514*bcb2dfaeSJed Brown$$ 515*bcb2dfaeSJed Brown 516*bcb2dfaeSJed Brownwhere $\diff \bm P$ is defined by {math:numref}`eq-diff-P` and {math:numref}`eq-neo-hookean-incremental-stress`, and $\mathcal V_0$ is the homogeneous space corresponding to $\mathcal V$. 517*bcb2dfaeSJed Brown 518*bcb2dfaeSJed Brown:::{note} 519*bcb2dfaeSJed BrownThe decision of whether to recompute or store functions of the current state $\bm F$ depends on a roofline analysis {cite}`williams2009roofline,brown2010` of the computation and the cost of the constitutive model. 520*bcb2dfaeSJed BrownFor low-order elements where flops tend to be in surplus relative to memory bandwidth, recomputation is likely to be preferable, where as the opposite may be true for high-order elements. 521*bcb2dfaeSJed BrownSimilarly, analysis with a simple constitutive model may see better performance while storing little or nothing while an expensive model such as Arruda-Boyce {cite}`arruda1993largestretch`, which contains many special functions, may be faster when using more storage to avoid recomputation. 522*bcb2dfaeSJed BrownIn the case where complete linearization is preferred, note the symmetry $\mathsf C_{IJKL} = \mathsf C_{KLIJ}$ evident in {math:numref}`eq-neo-hookean-incremental-stress-index`, thus $\mathsf C$ can be stored as a symmetric $6\times 6$ matrix, which has 21 unique entries. 523*bcb2dfaeSJed BrownAlong with 6 entries for $\bm S$, this totals 27 entries of overhead compared to computing everything from $\bm F$. 524*bcb2dfaeSJed BrownThis compares with 13 entries of overhead for direct storage of $\{ \bm S, \bm C^{-1}, \log J \}$, which is sufficient for the Neo-Hookean model to avoid all but matrix products. 525*bcb2dfaeSJed Brown::: 526*bcb2dfaeSJed Brown 527*bcb2dfaeSJed Brown(problem-hyperelasticity-finite-strain-current-configuration)= 528*bcb2dfaeSJed Brown 529*bcb2dfaeSJed Brown## Hyperelasticity in current configuration 530*bcb2dfaeSJed Brown 531*bcb2dfaeSJed BrownIn the preceeding discussion, all equations have been formulated in the initial configuration. 532*bcb2dfaeSJed BrownThis may feel convenient in that the computational domain is clearly independent of the solution, but there are some advantages to defining the equations in the current configuration. 533*bcb2dfaeSJed Brown 534*bcb2dfaeSJed Brown1. Body forces (like gravity), traction, and contact are more easily defined in the current configuration. 535*bcb2dfaeSJed Brown2. Mesh quality in the initial configuration can be very bad for large deformation. 536*bcb2dfaeSJed Brown3. The required storage and numerical representation can be smaller in the current configuration. 537*bcb2dfaeSJed Brown 538*bcb2dfaeSJed BrownMost of the benefit in case 3 can be attained solely by moving the Jacobian representation to the current configuration {cite}`davydov2020matrix`, though residual evaluation may also be slightly faster in current configuration. 539*bcb2dfaeSJed BrownThere are multiple commuting paths from the nonlinear weak form in initial configuration {math:numref}`hyperelastic-weak-form-initial` to the Jacobian weak form in current configuration {math:numref}`jacobian-weak-form-current`. 540*bcb2dfaeSJed BrownOne may push forward to the current configuration and then linearize or linearize in initial configuration and then push forward, as summarized below. 541*bcb2dfaeSJed Brown 542*bcb2dfaeSJed Brown$$ 543*bcb2dfaeSJed Brown\begin{CD} 544*bcb2dfaeSJed Brown {\overbrace{\nabla_X \bm{v} \tcolon \bm{FS}}^{\text{Initial Residual}}} 545*bcb2dfaeSJed Brown @>{\text{push forward}}>{}> 546*bcb2dfaeSJed Brown {\overbrace{\nabla_x \bm{v} \tcolon \bm{\tau}}^{\text{Current Residual}}} \\ 547*bcb2dfaeSJed Brown @V{\text{linearize}}V{\begin{smallmatrix} \diff\bm F = \nabla_X\diff\bm u \\ \diff\bm S(\diff\bm E) \end{smallmatrix}}V 548*bcb2dfaeSJed Brown @V{\begin{smallmatrix} \diff\nabla_x\bm v = -\nabla_x\bm v \nabla_x \diff\bm u \\ \diff\bm\tau(\diff\bm\epsilon) \end{smallmatrix}}V{\text{linearize}}V \\ 549*bcb2dfaeSJed Brown {\underbrace{\nabla_X\bm{v}\tcolon \Big(\diff\bm{F}\bm{S} + \bm{F}\diff\bm{S}\Big)}_\text{Initial Jacobian}} 550*bcb2dfaeSJed Brown @>{\text{push forward}}>{}> 551*bcb2dfaeSJed Brown {\underbrace{\nabla_x\bm{v}\tcolon \Big(\diff\bm{\tau} -\bm{\tau}(\nabla_x \diff\bm{u})^T \Big)}_\text{Current Jacobian}} 552*bcb2dfaeSJed Brown\end{CD} 553*bcb2dfaeSJed Brown$$ (initial-current-linearize) 554*bcb2dfaeSJed Brown 555*bcb2dfaeSJed BrownWe will follow both paths for consistency and because both intermediate representations may be useful for implementation. 556*bcb2dfaeSJed Brown 557*bcb2dfaeSJed Brown### Push forward, then linearize 558*bcb2dfaeSJed Brown 559*bcb2dfaeSJed BrownThe first term of {math:numref}`hyperelastic-weak-form-initial` can be rewritten in terms of the symmetric Kirchhoff stress tensor 560*bcb2dfaeSJed Brown$\bm{\tau}=J\bm{\sigma}=\bm{P}\bm{F}^T = \bm F \bm S \bm F^T$ as 561*bcb2dfaeSJed Brown 562*bcb2dfaeSJed Brown$$ 563*bcb2dfaeSJed Brown\nabla_X \bm{v} \tcolon \bm{P} = \nabla_X \bm{v} \tcolon \bm{\tau}\bm{F}^{-T} = \nabla_X \bm{v}\bm{F}^{-1} \tcolon \bm{\tau} = \nabla_x \bm{v} \tcolon \bm{\tau} 564*bcb2dfaeSJed Brown$$ 565*bcb2dfaeSJed Brown 566*bcb2dfaeSJed Browntherefore, the weak form in terms of $\bm{\tau}$ and $\nabla_x$ with integral over $\Omega_0$ is 567*bcb2dfaeSJed Brown 568*bcb2dfaeSJed Brown$$ 569*bcb2dfaeSJed Brown\int_{\Omega_0}{\nabla_x \bm{v} \tcolon \bm{\tau}} \, dV 570*bcb2dfaeSJed Brown - \int_{\Omega_0}{\bm{v} \cdot \rho_0 \bm{g}} \, dV 571*bcb2dfaeSJed Brown - \int_{\partial \Omega_0}{\bm{v}\cdot(\bm{P}\cdot\hat{\bm{N}})} \, dS 572*bcb2dfaeSJed Brown = 0, \quad \forall \bm v \in \mathcal V. 573*bcb2dfaeSJed Brown$$ (hyperelastic-weak-form-current) 574*bcb2dfaeSJed Brown 575*bcb2dfaeSJed Brown#### Linearize in current configuration 576*bcb2dfaeSJed Brown 577*bcb2dfaeSJed BrownTo derive a Newton linearization of {math:numref}`hyperelastic-weak-form-current`, first we define 578*bcb2dfaeSJed Brown 579*bcb2dfaeSJed Brown$$ 580*bcb2dfaeSJed Brown\nabla_x \diff \bm{u} = \nabla_X \diff \bm{u} \ \bm{F}^{-1} = \diff \bm{F} \bm{F}^{-1} 581*bcb2dfaeSJed Brown$$ (nabla_xdu) 582*bcb2dfaeSJed Brown 583*bcb2dfaeSJed Brownand $\bm{\tau}$ for Neo-Hookean materials as the push forward of {math:numref}`neo-hookean-stress` 584*bcb2dfaeSJed Brown 585*bcb2dfaeSJed Brown$$ 586*bcb2dfaeSJed Brown\bm{\tau} = \bm{F}\bm{S}\bm{F}^T = \mu (\bm{b} - \bm I_3) + \lambda \log J \bm{I}_3, 587*bcb2dfaeSJed Brown$$ (tau-neo-hookean) 588*bcb2dfaeSJed Brown 589*bcb2dfaeSJed Brownwhere $\bm{b} = \bm{F} \bm{F}^T$, is the left Cauchy-Green tensor. 590*bcb2dfaeSJed BrownThen by expanding the directional derivative of $\nabla_x \bm{v} \tcolon \bm{\tau}$, we arrive at 591*bcb2dfaeSJed Brown 592*bcb2dfaeSJed Brown$$ 593*bcb2dfaeSJed Brown\diff \ (\nabla_x \bm{v} \tcolon \bm{\tau}) = \diff \ (\nabla_x \bm{v})\tcolon \bm{\tau} + \nabla_x \bm{v} \tcolon \diff \bm{\tau} . 594*bcb2dfaeSJed Brown$$ (hyperelastic-linearization-current1) 595*bcb2dfaeSJed Brown 596*bcb2dfaeSJed BrownThe first term of {math:numref}`hyperelastic-linearization-current1` can be written as 597*bcb2dfaeSJed Brown 598*bcb2dfaeSJed Brown$$ 599*bcb2dfaeSJed Brown\begin{aligned} \diff \ (\nabla_x \bm{v})\tcolon \bm{\tau} &= \diff \ (\nabla_X \bm{v} \bm{F}^{-1})\tcolon \bm{\tau} = \Big(\underbrace{\nabla_X (\diff \bm{v})}_{0}\bm{F}^{-1} + \nabla_X \bm{v}\diff \bm{F}^{-1}\Big)\tcolon \bm{\tau}\\ &= \Big(-\nabla_X \bm{v} \bm{F}^{-1}\diff\bm{F}\bm{F}^{-1}\Big)\tcolon \bm{\tau}=\Big(-\nabla_x \bm{v} \diff\bm{F}\bm{F}^{-1}\Big)\tcolon \bm{\tau}\\ &= \Big(-\nabla_x \bm{v} \nabla_x \diff\bm{u} \Big)\tcolon \bm{\tau}= -\nabla_x \bm{v}\tcolon\bm{\tau}(\nabla_x \diff\bm{u})^T \,, \end{aligned} 600*bcb2dfaeSJed Brown$$ 601*bcb2dfaeSJed Brown 602*bcb2dfaeSJed Brownwhere we have used $\diff \bm{F}^{-1}=-\bm{F}^{-1} \diff \bm{F} \bm{F}^{-1}$ and {math:numref}`nabla_xdu`. 603*bcb2dfaeSJed BrownUsing this and {math:numref}`hyperelastic-linearization-current1` in {math:numref}`hyperelastic-weak-form-current` yields the weak form in the current configuration 604*bcb2dfaeSJed Brown 605*bcb2dfaeSJed Brown$$ 606*bcb2dfaeSJed Brown\int_{\Omega_0} \nabla_x \bm v \tcolon \Big(\diff\bm\tau - \bm\tau (\nabla_x \diff\bm u)^T \Big) = \text{rhs}. 607*bcb2dfaeSJed Brown$$ (jacobian-weak-form-current) 608*bcb2dfaeSJed Brown 609*bcb2dfaeSJed BrownIn the following, we will sometimes make use of the incremental strain tensor in the current configuration, 610*bcb2dfaeSJed Brown 611*bcb2dfaeSJed Brown$$ 612*bcb2dfaeSJed Brown\diff\bm\epsilon \equiv \frac{1}{2}\Big(\nabla_x \diff\bm{u} + (\nabla_x \diff\bm{u})^T \Big) . 613*bcb2dfaeSJed Brown$$ 614*bcb2dfaeSJed Brown 615*bcb2dfaeSJed Brown:::{dropdown} Deriving $\diff\bm\tau$ for Neo-Hookean material 616*bcb2dfaeSJed BrownTo derive a useful expression of $\diff\bm\tau$ for Neo-Hookean materials, we will use the representations 617*bcb2dfaeSJed Brown 618*bcb2dfaeSJed Brown$$ 619*bcb2dfaeSJed Brown\begin{aligned} 620*bcb2dfaeSJed Brown\diff \bm{b} &= \diff \bm{F} \bm{F}^T + \bm{F} \diff \bm{F}^T \\ 621*bcb2dfaeSJed Brown&= \nabla_x \diff \bm{u} \ \bm{b} + \bm{b} \ (\nabla_x \diff \bm{u})^T \\ 622*bcb2dfaeSJed Brown&= (\nabla_x \diff\bm u)(\bm b - \bm I_3) + (\bm b - \bm I_3) (\nabla_x \diff\bm u)^T + 2 \diff\bm\epsilon 623*bcb2dfaeSJed Brown\end{aligned} 624*bcb2dfaeSJed Brown$$ 625*bcb2dfaeSJed Brown 626*bcb2dfaeSJed Brownand 627*bcb2dfaeSJed Brown 628*bcb2dfaeSJed Brown$$ 629*bcb2dfaeSJed Brown\begin{aligned} \diff\ (\log J) &= \frac{\partial \log J}{\partial \bm{b}}\tcolon \diff \bm{b} = \frac{\partial J}{J\partial \bm{b}}\tcolon \diff \bm{b}=\frac{1}{2}\bm{b}^{-1}\tcolon \diff \bm{b} \\ &= \frac 1 2 \bm b^{-1} \tcolon \Big(\nabla_x \diff\bm u \ \bm b + \bm b (\nabla_x \diff\bm u)^T \Big) \\ &= \trace (\nabla_x \diff\bm u) \\ &= \trace \diff\bm\epsilon . \end{aligned} 630*bcb2dfaeSJed Brown$$ 631*bcb2dfaeSJed Brown 632*bcb2dfaeSJed BrownSubstituting into {math:numref}`tau-neo-hookean` gives 633*bcb2dfaeSJed Brown 634*bcb2dfaeSJed Brown$$ 635*bcb2dfaeSJed Brown\begin{aligned} 636*bcb2dfaeSJed Brown\diff \bm{\tau} &= \mu \diff \bm{b} + \lambda \trace (\diff\bm\epsilon) \bm I_3 \\ 637*bcb2dfaeSJed Brown&= \underbrace{2 \mu \diff\bm\epsilon + \lambda \trace (\diff\bm\epsilon) \bm I_3 - 2\lambda \log J \diff\bm\epsilon}_{\bm F \diff\bm S \bm F^T} \\ 638*bcb2dfaeSJed Brown&\quad + (\nabla_x \diff\bm u)\underbrace{\Big( \mu (\bm b - \bm I_3) + \lambda \log J \bm I_3 \Big)}_{\bm\tau} \\ 639*bcb2dfaeSJed Brown&\quad + \underbrace{\Big( \mu (\bm b - \bm I_3) + \lambda \log J \bm I_3 \Big)}_{\bm\tau} (\nabla_x \diff\bm u)^T , 640*bcb2dfaeSJed Brown\end{aligned} 641*bcb2dfaeSJed Brown$$ (dtau-neo-hookean) 642*bcb2dfaeSJed Brown 643*bcb2dfaeSJed Brownwhere the final expression has been identified according to 644*bcb2dfaeSJed Brown 645*bcb2dfaeSJed Brown$$ 646*bcb2dfaeSJed Brown\diff\bm\tau = \diff\ (\bm F \bm S \bm F^T) = (\nabla_x \diff\bm u) \bm\tau + \bm F \diff\bm S \bm F^T + \bm\tau(\nabla_x \diff\bm u)^T. 647*bcb2dfaeSJed Brown$$ 648*bcb2dfaeSJed Brown::: 649*bcb2dfaeSJed Brown 650*bcb2dfaeSJed BrownCollecting terms, we may thus opt to use either of the two forms 651*bcb2dfaeSJed Brown 652*bcb2dfaeSJed Brown$$ 653*bcb2dfaeSJed Brown\begin{aligned} 654*bcb2dfaeSJed Brown\diff \bm{\tau} -\bm{\tau}(\nabla_x \diff\bm{u})^T &= (\nabla_x \diff\bm u)\bm\tau + \bm F \diff\bm S \bm F^T \\ 655*bcb2dfaeSJed Brown&= (\nabla_x \diff\bm u)\bm\tau + \lambda \trace(\diff\bm\epsilon) \bm I_3 + 2(\mu - \lambda \log J) \diff\bm\epsilon, 656*bcb2dfaeSJed Brown\end{aligned} 657*bcb2dfaeSJed Brown$$ (cur_simp_Jac) 658*bcb2dfaeSJed Brown 659*bcb2dfaeSJed Brownwith the last line showing the especially compact representation available for Neo-Hookean materials. 660*bcb2dfaeSJed Brown 661*bcb2dfaeSJed Brown### Linearize, then push forward 662*bcb2dfaeSJed Brown 663*bcb2dfaeSJed BrownWe can move the derivatives to the current configuration via 664*bcb2dfaeSJed Brown 665*bcb2dfaeSJed Brown$$ 666*bcb2dfaeSJed Brown\nabla_X \bm v \!:\! \diff\bm P = (\nabla_X \bm v) \bm F^{-1} \!:\! \diff \bm P \bm F^T = \nabla_x \bm v \!:\! \diff\bm P \bm F^T 667*bcb2dfaeSJed Brown$$ 668*bcb2dfaeSJed Brown 669*bcb2dfaeSJed Brownand expand 670*bcb2dfaeSJed Brown 671*bcb2dfaeSJed Brown$$ 672*bcb2dfaeSJed Brown\begin{aligned} 673*bcb2dfaeSJed Brown\diff\bm P \bm F^T &= \diff\bm F \bm S \bm F^T + \bm F \diff\bm S \bm F^T \\ 674*bcb2dfaeSJed Brown&= \underbrace{\diff\bm F \bm F^{-1}}_{\nabla_x \diff\bm u} \underbrace{\bm F \bm S \bm F^T}_{\bm\tau} + \bm F \diff\bm S \bm F^T . 675*bcb2dfaeSJed Brown\end{aligned} 676*bcb2dfaeSJed Brown$$ 677*bcb2dfaeSJed Brown 678*bcb2dfaeSJed Brown:::{dropdown} Representation of $\bm F \diff\bm S \bm F^T$ for Neo-Hookean materials 679*bcb2dfaeSJed BrownNow we push {math:numref}`eq-neo-hookean-incremental-stress` forward via 680*bcb2dfaeSJed Brown 681*bcb2dfaeSJed Brown$$ 682*bcb2dfaeSJed Brown\begin{aligned} 683*bcb2dfaeSJed Brown\bm F \diff\bm S \bm F^T &= \lambda (\bm C^{-1} \!:\! \diff\bm E) \bm F \bm C^{-1} \bm F^T 684*bcb2dfaeSJed Brown + 2 (\mu - \lambda \log J) \bm F \bm C^{-1} \diff\bm E \, \bm C^{-1} \bm F^T \\ 685*bcb2dfaeSJed Brown &= \lambda (\bm C^{-1} \!:\! \diff\bm E) \bm I_3 + 2 (\mu - \lambda \log J) \bm F^{-T} \diff\bm E \, \bm F^{-1} \\ 686*bcb2dfaeSJed Brown &= \lambda \operatorname{trace}(\nabla_x \diff\bm u) \bm I_3 + 2 (\mu - \lambda \log J) \diff\bm \epsilon 687*bcb2dfaeSJed Brown\end{aligned} 688*bcb2dfaeSJed Brown$$ 689*bcb2dfaeSJed Brown 690*bcb2dfaeSJed Brownwhere we have used 691*bcb2dfaeSJed Brown 692*bcb2dfaeSJed Brown$$ 693*bcb2dfaeSJed Brown\begin{aligned} 694*bcb2dfaeSJed Brown\bm C^{-1} \!:\! \diff\bm E &= \bm F^{-1} \bm F^{-T} \!:\! \bm F^T \diff\bm F \\ 695*bcb2dfaeSJed Brown&= \operatorname{trace}(\bm F^{-1} \bm F^{-T} \bm F^T \diff \bm F) \\ 696*bcb2dfaeSJed Brown&= \operatorname{trace}(\bm F^{-1} \diff\bm F) \\ 697*bcb2dfaeSJed Brown&= \operatorname{trace}(\diff \bm F \bm F^{-1}) \\ 698*bcb2dfaeSJed Brown&= \operatorname{trace}(\nabla_x \diff\bm u) 699*bcb2dfaeSJed Brown\end{aligned} 700*bcb2dfaeSJed Brown$$ 701*bcb2dfaeSJed Brown 702*bcb2dfaeSJed Brownand 703*bcb2dfaeSJed Brown 704*bcb2dfaeSJed Brown$$ 705*bcb2dfaeSJed Brown\begin{aligned} 706*bcb2dfaeSJed Brown\bm F^{-T} \diff\bm E \, \bm F^{-1} &= \frac 1 2 \bm F^{-T} (\bm F^T \diff\bm F + \diff\bm F^T \bm F) \bm F^{-1} \\ 707*bcb2dfaeSJed Brown&= \frac 1 2 (\diff \bm F \bm F^{-1} + \bm F^{-T} \diff\bm F^T) \\ 708*bcb2dfaeSJed Brown&= \frac 1 2 \Big(\nabla_x \diff\bm u + (\nabla_x\diff\bm u)^T \Big) \equiv \diff\bm\epsilon. 709*bcb2dfaeSJed Brown\end{aligned} 710*bcb2dfaeSJed Brown$$ 711*bcb2dfaeSJed Brown::: 712*bcb2dfaeSJed Brown 713*bcb2dfaeSJed BrownCollecting terms, the weak form of the Newton linearization for Neo-Hookean materials in the current configuration is 714*bcb2dfaeSJed Brown 715*bcb2dfaeSJed Brown$$ 716*bcb2dfaeSJed Brown\int_{\Omega_0} \nabla_x \bm v \!:\! \Big( (\nabla_x \diff\bm u) \bm\tau + \lambda \operatorname{trace}(\diff\bm\epsilon)\bm I_3 + 2(\mu - \lambda\log J)\diff \bm\epsilon \Big) dV = \text{rhs}, 717*bcb2dfaeSJed Brown$$ (jacobian-weak-form-current2) 718*bcb2dfaeSJed Brown 719*bcb2dfaeSJed Brownwhich equivalent to Algorithm 2 of {cite}`davydov2020matrix` and requires only derivatives with respect to the current configuration. Note that {math:numref}`cur_simp_Jac` and {math:numref}`jacobian-weak-form-current2` have recovered the same representation 720*bcb2dfaeSJed Brownusing different algebraic manipulations. 721*bcb2dfaeSJed Brown 722*bcb2dfaeSJed Brown:::{tip} 723*bcb2dfaeSJed BrownWe define a second order *Green-Euler* strain tensor (cf. Green-Lagrange strain {math:numref}`eq-green-lagrange-strain`) as 724*bcb2dfaeSJed Brown 725*bcb2dfaeSJed Brown$$ 726*bcb2dfaeSJed Brown\bm e = \frac 1 2 \Big(\bm{b} - \bm{I}_3 \Big) = \frac 1 2 \Big( \nabla_X \bm{u} + (\nabla_X \bm{u})^T + \nabla_X \bm{u} \, (\nabla_X \bm{u})^T \Big). 727*bcb2dfaeSJed Brown$$ (green-euler-strain) 728*bcb2dfaeSJed Brown 729*bcb2dfaeSJed BrownThen, the Kirchhoff stress tensor {math:numref}`tau-neo-hookean` can be written as 730*bcb2dfaeSJed Brown 731*bcb2dfaeSJed Brown$$ 732*bcb2dfaeSJed Brown\bm \tau = \lambda \log J \bm I_{3} + 2\mu \bm e, 733*bcb2dfaeSJed Brown$$ (tau-neo-hookean-stable) 734*bcb2dfaeSJed Brown 735*bcb2dfaeSJed Brownwhich is more numerically stable for small strain, and thus preferred for computation. Note that the $\log J$ is computed via `log1p` {math:numref}`log1p`, as we discussed in the previous tip. 736*bcb2dfaeSJed Brown::: 737*bcb2dfaeSJed Brown 738*bcb2dfaeSJed Brown### Jacobian representation 739*bcb2dfaeSJed Brown 740*bcb2dfaeSJed BrownWe have implemented four storage variants for the Jacobian in our finite strain hyperelasticity. In each case, some variables are computed during residual evaluation and used during Jacobian application. 741*bcb2dfaeSJed Brown 742*bcb2dfaeSJed Brown```{eval-rst} 743*bcb2dfaeSJed Brown.. list-table:: Four algorithms for Jacobian action in finite strain hyperelasticity problem 744*bcb2dfaeSJed Brown :header-rows: 1 745*bcb2dfaeSJed Brown 746*bcb2dfaeSJed Brown * - Option :code:`-problem` 747*bcb2dfaeSJed Brown - Static storage 748*bcb2dfaeSJed Brown - Computed storage 749*bcb2dfaeSJed Brown - # scalars 750*bcb2dfaeSJed Brown - Equations 751*bcb2dfaeSJed Brown 752*bcb2dfaeSJed Brown 753*bcb2dfaeSJed Brown * - :code:`FSInitial-NH1` 754*bcb2dfaeSJed Brown - :math:`\nabla_{X} \hat X, \operatorname{det}\nabla_{\hat X} X` 755*bcb2dfaeSJed Brown - :math:`\nabla_X \bm u` 756*bcb2dfaeSJed Brown - 19 757*bcb2dfaeSJed Brown - :math:numref:`eq-diff-P` :math:numref:`eq-neo-hookean-incremental-stress` 758*bcb2dfaeSJed Brown 759*bcb2dfaeSJed Brown * - :code:`FSInitial-NH2` 760*bcb2dfaeSJed Brown - :math:`\nabla_{X} \hat X, \operatorname{det}\nabla_{\hat X} X` 761*bcb2dfaeSJed Brown - :math:`\nabla_X \bm u, \bm C^{-1}, \lambda \log J` 762*bcb2dfaeSJed Brown - 26 763*bcb2dfaeSJed Brown - :math:numref:`eq-diff-P` :math:numref:`eq-neo-hookean-incremental-stress` 764*bcb2dfaeSJed Brown 765*bcb2dfaeSJed Brown * - :code:`FSCurrent-NH1` 766*bcb2dfaeSJed Brown - :math:`\nabla_{X} \hat X, \operatorname{det}\nabla_{\hat X} X` 767*bcb2dfaeSJed Brown - :math:`\nabla_X \bm u` 768*bcb2dfaeSJed Brown - 19 769*bcb2dfaeSJed Brown - :math:numref:`jacobian-weak-form-current` :math:numref:`nabla_xdu` 770*bcb2dfaeSJed Brown 771*bcb2dfaeSJed Brown * - :code:`FSCurrent-NH2` 772*bcb2dfaeSJed Brown - :math:`\operatorname{det}\nabla_{\hat X} X` 773*bcb2dfaeSJed Brown - :math:`\nabla_x \hat X, \bm \tau, \lambda \log J` 774*bcb2dfaeSJed Brown - 17 775*bcb2dfaeSJed Brown - :math:numref:`jacobian-weak-form-current` :math:numref:`jacobian-weak-form-current2` 776*bcb2dfaeSJed Brown``` 777