xref: /libCEED/examples/solids/qfunctions/linear.h (revision c8565611f4f88586c9ab8f49f4be6e8b5d8096a7)
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
35754ecacSJeremy L Thompson //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
55754ecacSJeremy L Thompson //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
75754ecacSJeremy L Thompson 
85754ecacSJeremy L Thompson /// @file
95754ecacSJeremy L Thompson /// Linear elasticity for solid mechanics example using PETSc
105754ecacSJeremy L Thompson 
11c9c2c079SJeremy L Thompson #include <ceed.h>
125754ecacSJeremy L Thompson #include <math.h>
135754ecacSJeremy L Thompson 
145754ecacSJeremy L Thompson #ifndef PHYSICS_STRUCT
155754ecacSJeremy L Thompson #define PHYSICS_STRUCT
165754ecacSJeremy L Thompson typedef struct Physics_private *Physics;
175754ecacSJeremy L Thompson struct Physics_private {
185754ecacSJeremy L Thompson   CeedScalar nu;  // Poisson's ratio
195754ecacSJeremy L Thompson   CeedScalar E;   // Young's Modulus
205754ecacSJeremy L Thompson };
215754ecacSJeremy L Thompson #endif
225754ecacSJeremy L Thompson 
235754ecacSJeremy L Thompson // -----------------------------------------------------------------------------
245754ecacSJeremy L Thompson // Residual evaluation for linear elasticity
255754ecacSJeremy L Thompson // -----------------------------------------------------------------------------
26*c8565611SJeremy L Thompson CEED_QFUNCTION(ElasResidual_Linear)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
275754ecacSJeremy L Thompson   // Inputs
282b730f8bSJeremy L Thompson   const CeedScalar(*ug)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[0], (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1];
295754ecacSJeremy L Thompson 
305754ecacSJeremy L Thompson   // Outputs
315754ecacSJeremy L Thompson   CeedScalar(*dvdX)[3][CEED_Q_VLA] = (CeedScalar(*)[3][CEED_Q_VLA])out[0];
325754ecacSJeremy L Thompson   // grad_u not used for linear elasticity
335754ecacSJeremy L Thompson   // (*grad_u)[3][CEED_Q_VLA] = (CeedScalar(*)[3][CEED_Q_VLA])out[1];
345754ecacSJeremy L Thompson 
355754ecacSJeremy L Thompson   // Context
365754ecacSJeremy L Thompson   const Physics    context = (Physics)ctx;
375754ecacSJeremy L Thompson   const CeedScalar E       = context->E;
385754ecacSJeremy L Thompson   const CeedScalar nu      = context->nu;
395754ecacSJeremy L Thompson 
405754ecacSJeremy L Thompson   // Quadrature Point Loop
412b730f8bSJeremy L Thompson   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
425754ecacSJeremy L Thompson     // Read spatial derivatives of u
432b730f8bSJeremy L Thompson     const CeedScalar du[3][3] = {
442b730f8bSJeremy L Thompson         {ug[0][0][i], ug[1][0][i], ug[2][0][i]},
452b730f8bSJeremy L Thompson         {ug[0][1][i], ug[1][1][i], ug[2][1][i]},
462b730f8bSJeremy L Thompson         {ug[0][2][i], ug[1][2][i], ug[2][2][i]}
475754ecacSJeremy L Thompson     };
485754ecacSJeremy L Thompson     // -- Qdata
495754ecacSJeremy L Thompson     const CeedScalar wdetJ      = q_data[0][i];
502b730f8bSJeremy L Thompson     const CeedScalar dXdx[3][3] = {
512b730f8bSJeremy L Thompson         {q_data[1][i], q_data[2][i], q_data[3][i]},
522b730f8bSJeremy L Thompson         {q_data[4][i], q_data[5][i], q_data[6][i]},
532b730f8bSJeremy L Thompson         {q_data[7][i], q_data[8][i], q_data[9][i]}
545754ecacSJeremy L Thompson     };
555754ecacSJeremy L Thompson 
565754ecacSJeremy L Thompson     // Compute grad_u
575754ecacSJeremy L Thompson     //   dXdx = (dx/dX)^(-1)
585754ecacSJeremy L Thompson     // Apply dXdx to du = grad_u
595754ecacSJeremy L Thompson     CeedScalar grad_u[3][3];
602b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {    // Component
615754ecacSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) {  // Derivative
625754ecacSJeremy L Thompson         grad_u[j][k] = 0;
632b730f8bSJeremy L Thompson         for (CeedInt m = 0; m < 3; m++) grad_u[j][k] += dXdx[m][k] * du[j][m];
642b730f8bSJeremy L Thompson       }
655754ecacSJeremy L Thompson     }
665754ecacSJeremy L Thompson 
675754ecacSJeremy L Thompson     // Compute Strain : e (epsilon)
685754ecacSJeremy L Thompson     // e = 1/2 (grad u + (grad u)^T)
695754ecacSJeremy L Thompson 
702b730f8bSJeremy L Thompson     const CeedScalar e[3][3] = {
712b730f8bSJeremy L Thompson         {(grad_u[0][0] + grad_u[0][0]) / 2., (grad_u[0][1] + grad_u[1][0]) / 2., (grad_u[0][2] + grad_u[2][0]) / 2.},
722b730f8bSJeremy L Thompson         {(grad_u[1][0] + grad_u[0][1]) / 2., (grad_u[1][1] + grad_u[1][1]) / 2., (grad_u[1][2] + grad_u[2][1]) / 2.},
732b730f8bSJeremy L Thompson         {(grad_u[2][0] + grad_u[0][2]) / 2., (grad_u[2][1] + grad_u[1][2]) / 2., (grad_u[2][2] + grad_u[2][2]) / 2.}
745754ecacSJeremy L Thompson     };
755754ecacSJeremy L Thompson 
765754ecacSJeremy L Thompson     //
775754ecacSJeremy L Thompson     // Formulation uses Voigt notation:
785754ecacSJeremy L Thompson     //  stress (sigma)      strain (epsilon)
795754ecacSJeremy L Thompson     //
805754ecacSJeremy L Thompson     //    [sigma00]             [e00]
815754ecacSJeremy L Thompson     //    [sigma11]             [e11]
825754ecacSJeremy L Thompson     //    [sigma22]   =  S   *  [e22]
835754ecacSJeremy L Thompson     //    [sigma12]             [e12]
845754ecacSJeremy L Thompson     //    [sigma02]             [e02]
855754ecacSJeremy L Thompson     //    [sigma01]             [e01]
865754ecacSJeremy L Thompson     //
875754ecacSJeremy L Thompson     //        where
885754ecacSJeremy L Thompson     //                         [1-nu   nu    nu                                    ]
895754ecacSJeremy L Thompson     //                         [ nu   1-nu   nu                                    ]
905754ecacSJeremy L Thompson     //                         [ nu    nu   1-nu                                   ]
915754ecacSJeremy L Thompson     // S = E/((1+nu)*(1-2*nu)) [                  (1-2*nu)/2                       ]
925754ecacSJeremy L Thompson     //                         [                             (1-2*nu)/2            ]
935754ecacSJeremy L Thompson     //                         [                                        (1-2*nu)/2 ]
945754ecacSJeremy L Thompson 
955754ecacSJeremy L Thompson     // Above Voigt Notation is placed in a 3x3 matrix:
965754ecacSJeremy L Thompson     const CeedScalar ss      = E / ((1 + nu) * (1 - 2 * nu));
975754ecacSJeremy L Thompson     const CeedScalar sigma00 = ss * ((1 - nu) * e[0][0] + nu * e[1][1] + nu * e[2][2]),
985754ecacSJeremy L Thompson                      sigma11 = ss * (nu * e[0][0] + (1 - nu) * e[1][1] + nu * e[2][2]),
992b730f8bSJeremy L Thompson                      sigma22 = ss * (nu * e[0][0] + nu * e[1][1] + (1 - nu) * e[2][2]), sigma12 = ss * (1 - 2 * nu) * e[1][2] * 0.5,
1002b730f8bSJeremy L Thompson                      sigma02 = ss * (1 - 2 * nu) * e[0][2] * 0.5, sigma01 = ss * (1 - 2 * nu) * e[0][1] * 0.5;
1012b730f8bSJeremy L Thompson     const CeedScalar sigma[3][3] = {
1022b730f8bSJeremy L Thompson         {sigma00, sigma01, sigma02},
1035754ecacSJeremy L Thompson         {sigma01, sigma11, sigma12},
1045754ecacSJeremy L Thompson         {sigma02, sigma12, sigma22}
1055754ecacSJeremy L Thompson     };
1065754ecacSJeremy L Thompson 
1075754ecacSJeremy L Thompson     // Apply dXdx^T and weight to sigma
1082b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {    // Component
1095754ecacSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) {  // Derivative
1105754ecacSJeremy L Thompson         dvdX[k][j][i] = 0;
1112b730f8bSJeremy L Thompson         for (CeedInt m = 0; m < 3; m++) dvdX[k][j][i] += dXdx[k][m] * sigma[j][m] * wdetJ;
1125754ecacSJeremy L Thompson       }
1132b730f8bSJeremy L Thompson     }
1145754ecacSJeremy L Thompson   }  // End of Quadrature Point Loop
1155754ecacSJeremy L Thompson 
116*c8565611SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1175754ecacSJeremy L Thompson }
1185754ecacSJeremy L Thompson 
1195754ecacSJeremy L Thompson // -----------------------------------------------------------------------------
1205754ecacSJeremy L Thompson // Jacobian evaluation for linear elasticity
1215754ecacSJeremy L Thompson // -----------------------------------------------------------------------------
122*c8565611SJeremy L Thompson CEED_QFUNCTION(ElasJacobian_Linear)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
1235754ecacSJeremy L Thompson   // Inputs
1245754ecacSJeremy L Thompson   const CeedScalar(*deltaug)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[0],
1255754ecacSJeremy L Thompson         (*q_data)[CEED_Q_VLA]               = (const CeedScalar(*)[CEED_Q_VLA])in[1];
1265754ecacSJeremy L Thompson   // grad_u not used for linear elasticity
1275754ecacSJeremy L Thompson   // (*grad_u)[3][Q] = (CeedScalar(*)[3][Q])in[2];
1285754ecacSJeremy L Thompson 
1295754ecacSJeremy L Thompson   // Outputs
1305754ecacSJeremy L Thompson   CeedScalar(*deltadvdX)[3][CEED_Q_VLA] = (CeedScalar(*)[3][CEED_Q_VLA])out[0];
1315754ecacSJeremy L Thompson 
1325754ecacSJeremy L Thompson   // Context
1335754ecacSJeremy L Thompson   const Physics    context = (Physics)ctx;
1345754ecacSJeremy L Thompson   const CeedScalar E       = context->E;
1355754ecacSJeremy L Thompson   const CeedScalar nu      = context->nu;
1365754ecacSJeremy L Thompson 
1375754ecacSJeremy L Thompson   // Quadrature Point Loop
1382b730f8bSJeremy L Thompson   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
1395754ecacSJeremy L Thompson     // Read spatial derivatives of u
1402b730f8bSJeremy L Thompson     const CeedScalar deltadu[3][3] = {
1412b730f8bSJeremy L Thompson         {deltaug[0][0][i], deltaug[1][0][i], deltaug[2][0][i]},
1422b730f8bSJeremy L Thompson         {deltaug[0][1][i], deltaug[1][1][i], deltaug[2][1][i]},
1432b730f8bSJeremy L Thompson         {deltaug[0][2][i], deltaug[1][2][i], deltaug[2][2][i]}
1445754ecacSJeremy L Thompson     };
1455754ecacSJeremy L Thompson     // -- Qdata
1465754ecacSJeremy L Thompson     const CeedScalar wdetJ      = q_data[0][i];
1472b730f8bSJeremy L Thompson     const CeedScalar dXdx[3][3] = {
1482b730f8bSJeremy L Thompson         {q_data[1][i], q_data[2][i], q_data[3][i]},
1492b730f8bSJeremy L Thompson         {q_data[4][i], q_data[5][i], q_data[6][i]},
1502b730f8bSJeremy L Thompson         {q_data[7][i], q_data[8][i], q_data[9][i]}
1515754ecacSJeremy L Thompson     };
1525754ecacSJeremy L Thompson 
1535754ecacSJeremy L Thompson     // Compute graddeltau
1545754ecacSJeremy L Thompson     //   dXdx = (dx/dX)^(-1)
1555754ecacSJeremy L Thompson     // Apply dXdx to deltadu = graddeltau
1565754ecacSJeremy L Thompson     CeedScalar graddeltau[3][3];
1572b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {    // Component
1585754ecacSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) {  // Derivative
1595754ecacSJeremy L Thompson         graddeltau[j][k] = 0;
1602b730f8bSJeremy L Thompson         for (CeedInt m = 0; m < 3; m++) graddeltau[j][k] += dXdx[m][k] * deltadu[j][m];
1612b730f8bSJeremy L Thompson       }
1625754ecacSJeremy L Thompson     }
1635754ecacSJeremy L Thompson 
1645754ecacSJeremy L Thompson     // Compute Strain : e (epsilon)
1655754ecacSJeremy L Thompson     // e = 1/2 (grad u + (grad u)^T)
1662b730f8bSJeremy L Thompson     const CeedScalar de[3][3] = {
1672b730f8bSJeremy L Thompson         {(graddeltau[0][0] + graddeltau[0][0]) / 2., (graddeltau[0][1] + graddeltau[1][0]) / 2., (graddeltau[0][2] + graddeltau[2][0]) / 2.},
1682b730f8bSJeremy L Thompson         {(graddeltau[1][0] + graddeltau[0][1]) / 2., (graddeltau[1][1] + graddeltau[1][1]) / 2., (graddeltau[1][2] + graddeltau[2][1]) / 2.},
1692b730f8bSJeremy L Thompson         {(graddeltau[2][0] + graddeltau[0][2]) / 2., (graddeltau[2][1] + graddeltau[1][2]) / 2., (graddeltau[2][2] + graddeltau[2][2]) / 2.}
1705754ecacSJeremy L Thompson     };
1715754ecacSJeremy L Thompson 
1725754ecacSJeremy L Thompson     // Formulation uses Voigt notation:
1735754ecacSJeremy L Thompson     //  stress (sigma)      strain (epsilon)
1745754ecacSJeremy L Thompson     //
1755754ecacSJeremy L Thompson     //    [dsigma00]             [de00]
1765754ecacSJeremy L Thompson     //    [dsigma11]             [de11]
1775754ecacSJeremy L Thompson     //    [dsigma22]   =  S   *  [de22]
1785754ecacSJeremy L Thompson     //    [dsigma12]             [de12]
1795754ecacSJeremy L Thompson     //    [dsigma02]             [de02]
1805754ecacSJeremy L Thompson     //    [dsigma01]             [de01]
1815754ecacSJeremy L Thompson     //
1825754ecacSJeremy L Thompson     //        where
1835754ecacSJeremy L Thompson     //                         [1-nu   nu    nu                                    ]
1845754ecacSJeremy L Thompson     //                         [ nu   1-nu   nu                                    ]
1855754ecacSJeremy L Thompson     //                         [ nu    nu   1-nu                                   ]
1865754ecacSJeremy L Thompson     // S = E/((1+nu)*(1-2*nu)) [                  (1-2*nu)/2                       ]
1875754ecacSJeremy L Thompson     //                         [                             (1-2*nu)/2            ]
1885754ecacSJeremy L Thompson     //                         [                                        (1-2*nu)/2 ]
1895754ecacSJeremy L Thompson 
1905754ecacSJeremy L Thompson     // Above Voigt Notation is placed in a 3x3 matrix:
1915754ecacSJeremy L Thompson     const CeedScalar ss       = E / ((1 + nu) * (1 - 2 * nu));
1925754ecacSJeremy L Thompson     const CeedScalar dsigma00 = ss * ((1 - nu) * de[0][0] + nu * de[1][1] + nu * de[2][2]),
1935754ecacSJeremy L Thompson                      dsigma11 = ss * (nu * de[0][0] + (1 - nu) * de[1][1] + nu * de[2][2]),
1942b730f8bSJeremy L Thompson                      dsigma22 = ss * (nu * de[0][0] + nu * de[1][1] + (1 - nu) * de[2][2]), dsigma12 = ss * (1 - 2 * nu) * de[1][2] / 2,
1952b730f8bSJeremy L Thompson                      dsigma02 = ss * (1 - 2 * nu) * de[0][2] / 2, dsigma01 = ss * (1 - 2 * nu) * de[0][1] / 2;
1962b730f8bSJeremy L Thompson     const CeedScalar dsigma[3][3] = {
1972b730f8bSJeremy L Thompson         {dsigma00, dsigma01, dsigma02},
1985754ecacSJeremy L Thompson         {dsigma01, dsigma11, dsigma12},
1995754ecacSJeremy L Thompson         {dsigma02, dsigma12, dsigma22}
2005754ecacSJeremy L Thompson     };
2015754ecacSJeremy L Thompson 
2025754ecacSJeremy L Thompson     // Apply dXdx^T and weight
2032b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {    // Component
2045754ecacSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) {  // Derivative
2055754ecacSJeremy L Thompson         deltadvdX[k][j][i] = 0;
2062b730f8bSJeremy L Thompson         for (CeedInt m = 0; m < 3; m++) deltadvdX[k][j][i] += dXdx[k][m] * dsigma[j][m] * wdetJ;
2075754ecacSJeremy L Thompson       }
2082b730f8bSJeremy L Thompson     }
2095754ecacSJeremy L Thompson   }  // End of Quadrature Point Loop
2105754ecacSJeremy L Thompson 
211*c8565611SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2125754ecacSJeremy L Thompson }
2135754ecacSJeremy L Thompson 
2145754ecacSJeremy L Thompson // -----------------------------------------------------------------------------
2155754ecacSJeremy L Thompson // Strain energy computation for linear elasticity
2165754ecacSJeremy L Thompson // -----------------------------------------------------------------------------
217*c8565611SJeremy L Thompson CEED_QFUNCTION(ElasEnergy_Linear)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
2185754ecacSJeremy L Thompson   // Inputs
2192b730f8bSJeremy L Thompson   const CeedScalar(*ug)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[0], (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1];
2205754ecacSJeremy L Thompson 
2215754ecacSJeremy L Thompson   // Outputs
2225754ecacSJeremy L Thompson   CeedScalar(*energy) = (CeedScalar(*))out[0];
2235754ecacSJeremy L Thompson 
2245754ecacSJeremy L Thompson   // Context
2255754ecacSJeremy L Thompson   const Physics    context = (Physics)ctx;
2265754ecacSJeremy L Thompson   const CeedScalar E       = context->E;
2275754ecacSJeremy L Thompson   const CeedScalar nu      = context->nu;
2285754ecacSJeremy L Thompson 
2295754ecacSJeremy L Thompson   // Constants
2305754ecacSJeremy L Thompson   const CeedScalar TwoMu  = E / (1 + nu);
2315754ecacSJeremy L Thompson   const CeedScalar mu     = TwoMu / 2;
2325754ecacSJeremy L Thompson   const CeedScalar Kbulk  = E / (3 * (1 - 2 * nu));  // Bulk Modulus
2335754ecacSJeremy L Thompson   const CeedScalar lambda = (3 * Kbulk - TwoMu) / 3;
2345754ecacSJeremy L Thompson 
2355754ecacSJeremy L Thompson   // Quadrature Point Loop
2362b730f8bSJeremy L Thompson   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
2375754ecacSJeremy L Thompson     // Read spatial derivatives of u
2382b730f8bSJeremy L Thompson     const CeedScalar du[3][3] = {
2392b730f8bSJeremy L Thompson         {ug[0][0][i], ug[1][0][i], ug[2][0][i]},
2402b730f8bSJeremy L Thompson         {ug[0][1][i], ug[1][1][i], ug[2][1][i]},
2412b730f8bSJeremy L Thompson         {ug[0][2][i], ug[1][2][i], ug[2][2][i]}
2425754ecacSJeremy L Thompson     };
2435754ecacSJeremy L Thompson     // -- Qdata
2445754ecacSJeremy L Thompson     const CeedScalar wdetJ      = q_data[0][i];
2452b730f8bSJeremy L Thompson     const CeedScalar dXdx[3][3] = {
2462b730f8bSJeremy L Thompson         {q_data[1][i], q_data[2][i], q_data[3][i]},
2472b730f8bSJeremy L Thompson         {q_data[4][i], q_data[5][i], q_data[6][i]},
2482b730f8bSJeremy L Thompson         {q_data[7][i], q_data[8][i], q_data[9][i]}
2495754ecacSJeremy L Thompson     };
2505754ecacSJeremy L Thompson 
2515754ecacSJeremy L Thompson     // Compute grad_u
2525754ecacSJeremy L Thompson     //   dXdx = (dx/dX)^(-1)
2535754ecacSJeremy L Thompson     // Apply dXdx to du = grad_u
2545754ecacSJeremy L Thompson     CeedScalar grad_u[3][3];
2552b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {    // Component
2565754ecacSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) {  // Derivative
2575754ecacSJeremy L Thompson         grad_u[j][k] = 0;
2582b730f8bSJeremy L Thompson         for (CeedInt m = 0; m < 3; m++) grad_u[j][k] += dXdx[m][k] * du[j][m];
2592b730f8bSJeremy L Thompson       }
2605754ecacSJeremy L Thompson     }
2615754ecacSJeremy L Thompson 
2625754ecacSJeremy L Thompson     // Compute Strain : e (epsilon)
2635754ecacSJeremy L Thompson     // e = 1/2 (grad u + (grad u)^T)
2645754ecacSJeremy L Thompson 
2652b730f8bSJeremy L Thompson     const CeedScalar e[3][3] = {
2662b730f8bSJeremy L Thompson         {(grad_u[0][0] + grad_u[0][0]) / 2., (grad_u[0][1] + grad_u[1][0]) / 2., (grad_u[0][2] + grad_u[2][0]) / 2.},
2672b730f8bSJeremy L Thompson         {(grad_u[1][0] + grad_u[0][1]) / 2., (grad_u[1][1] + grad_u[1][1]) / 2., (grad_u[1][2] + grad_u[2][1]) / 2.},
2682b730f8bSJeremy L Thompson         {(grad_u[2][0] + grad_u[0][2]) / 2., (grad_u[2][1] + grad_u[1][2]) / 2., (grad_u[2][2] + grad_u[2][2]) / 2.}
2695754ecacSJeremy L Thompson     };
2705754ecacSJeremy L Thompson 
2715754ecacSJeremy L Thompson     // Strain energy
2725754ecacSJeremy L Thompson     const CeedScalar strain_vol = e[0][0] + e[1][1] + e[2][2];
2732b730f8bSJeremy L Thompson     energy[i] =
2742b730f8bSJeremy L Thompson         (lambda * strain_vol * strain_vol / 2. + strain_vol * mu + (e[0][1] * e[0][1] + e[0][2] * e[0][2] + e[1][2] * e[1][2]) * 2 * mu) * wdetJ;
2755754ecacSJeremy L Thompson 
2765754ecacSJeremy L Thompson   }  // End of Quadrature Point Loop
2775754ecacSJeremy L Thompson 
278*c8565611SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2795754ecacSJeremy L Thompson }
2805754ecacSJeremy L Thompson 
2815754ecacSJeremy L Thompson // -----------------------------------------------------------------------------
2825754ecacSJeremy L Thompson // Nodal diagnostic quantities for linear elasticity
2835754ecacSJeremy L Thompson // -----------------------------------------------------------------------------
284*c8565611SJeremy L Thompson CEED_QFUNCTION(ElasDiagnostic_Linear)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
2855754ecacSJeremy L Thompson   // Inputs
2862b730f8bSJeremy L Thompson   const CeedScalar(*u)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], (*ug)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[1],
2875754ecacSJeremy L Thompson         (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
2885754ecacSJeremy L Thompson 
2895754ecacSJeremy L Thompson   // Outputs
2905754ecacSJeremy L Thompson   CeedScalar(*diagnostic)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
2915754ecacSJeremy L Thompson 
2925754ecacSJeremy L Thompson   // Context
2935754ecacSJeremy L Thompson   const Physics    context = (Physics)ctx;
2945754ecacSJeremy L Thompson   const CeedScalar E       = context->E;
2955754ecacSJeremy L Thompson   const CeedScalar nu      = context->nu;
2965754ecacSJeremy L Thompson 
2975754ecacSJeremy L Thompson   // Constants
2985754ecacSJeremy L Thompson   const CeedScalar TwoMu  = E / (1 + nu);
2995754ecacSJeremy L Thompson   const CeedScalar mu     = TwoMu / 2;
3005754ecacSJeremy L Thompson   const CeedScalar Kbulk  = E / (3 * (1 - 2 * nu));  // Bulk Modulus
3015754ecacSJeremy L Thompson   const CeedScalar lambda = (3 * Kbulk - TwoMu) / 3;
3025754ecacSJeremy L Thompson 
3035754ecacSJeremy L Thompson   // Quadrature Point Loop
3042b730f8bSJeremy L Thompson   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
3055754ecacSJeremy L Thompson     // Read spatial derivatives of u
3062b730f8bSJeremy L Thompson     const CeedScalar du[3][3] = {
3072b730f8bSJeremy L Thompson         {ug[0][0][i], ug[1][0][i], ug[2][0][i]},
3082b730f8bSJeremy L Thompson         {ug[0][1][i], ug[1][1][i], ug[2][1][i]},
3092b730f8bSJeremy L Thompson         {ug[0][2][i], ug[1][2][i], ug[2][2][i]}
3105754ecacSJeremy L Thompson     };
3115754ecacSJeremy L Thompson     // -- Qdata
3122b730f8bSJeremy L Thompson     const CeedScalar dXdx[3][3] = {
3132b730f8bSJeremy L Thompson         {q_data[1][i], q_data[2][i], q_data[3][i]},
3142b730f8bSJeremy L Thompson         {q_data[4][i], q_data[5][i], q_data[6][i]},
3152b730f8bSJeremy L Thompson         {q_data[7][i], q_data[8][i], q_data[9][i]}
3165754ecacSJeremy L Thompson     };
3175754ecacSJeremy L Thompson 
3185754ecacSJeremy L Thompson     // Compute grad_u
3195754ecacSJeremy L Thompson     //   dXdx = (dx/dX)^(-1)
3205754ecacSJeremy L Thompson     // Apply dXdx to du = grad_u
3215754ecacSJeremy L Thompson     CeedScalar grad_u[3][3];
3222b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {    // Component
3235754ecacSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) {  // Derivative
3245754ecacSJeremy L Thompson         grad_u[j][k] = 0;
3252b730f8bSJeremy L Thompson         for (CeedInt m = 0; m < 3; m++) grad_u[j][k] += dXdx[m][k] * du[j][m];
3262b730f8bSJeremy L Thompson       }
3275754ecacSJeremy L Thompson     }
3285754ecacSJeremy L Thompson 
3295754ecacSJeremy L Thompson     // Compute Strain : e (epsilon)
3305754ecacSJeremy L Thompson     // e = 1/2 (grad u + (grad u)^T)
3315754ecacSJeremy L Thompson 
3322b730f8bSJeremy L Thompson     const CeedScalar e[3][3] = {
3332b730f8bSJeremy L Thompson         {(grad_u[0][0] + grad_u[0][0]) / 2., (grad_u[0][1] + grad_u[1][0]) / 2., (grad_u[0][2] + grad_u[2][0]) / 2.},
3342b730f8bSJeremy L Thompson         {(grad_u[1][0] + grad_u[0][1]) / 2., (grad_u[1][1] + grad_u[1][1]) / 2., (grad_u[1][2] + grad_u[2][1]) / 2.},
3352b730f8bSJeremy L Thompson         {(grad_u[2][0] + grad_u[0][2]) / 2., (grad_u[2][1] + grad_u[1][2]) / 2., (grad_u[2][2] + grad_u[2][2]) / 2.}
3365754ecacSJeremy L Thompson     };
3375754ecacSJeremy L Thompson 
3385754ecacSJeremy L Thompson     // Displacement
3395754ecacSJeremy L Thompson     diagnostic[0][i] = u[0][i];
3405754ecacSJeremy L Thompson     diagnostic[1][i] = u[1][i];
3415754ecacSJeremy L Thompson     diagnostic[2][i] = u[2][i];
3425754ecacSJeremy L Thompson 
3435754ecacSJeremy L Thompson     // Pressure
3445754ecacSJeremy L Thompson     const CeedScalar strain_vol = e[0][0] + e[1][1] + e[2][2];
3455754ecacSJeremy L Thompson     diagnostic[3][i]            = -lambda * strain_vol;
3465754ecacSJeremy L Thompson 
3475754ecacSJeremy L Thompson     // Stress tensor invariants
3485754ecacSJeremy L Thompson     diagnostic[4][i] = strain_vol;
3495754ecacSJeremy L Thompson     diagnostic[5][i] = 0.;
3502b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {
3512b730f8bSJeremy L Thompson       for (CeedInt m = 0; m < 3; m++) diagnostic[5][i] += e[j][m] * e[m][j];
3522b730f8bSJeremy L Thompson     }
3535754ecacSJeremy L Thompson     diagnostic[6][i] = 1 + strain_vol;
3545754ecacSJeremy L Thompson 
3555754ecacSJeremy L Thompson     // Strain energy
3562b730f8bSJeremy L Thompson     diagnostic[7][i] =
3572b730f8bSJeremy L Thompson         (lambda * strain_vol * strain_vol / 2. + strain_vol * mu + (e[0][1] * e[0][1] + e[0][2] * e[0][2] + e[1][2] * e[1][2]) * 2 * mu);
3585754ecacSJeremy L Thompson   }  // End of Quadrature Point Loop
3595754ecacSJeremy L Thompson 
360*c8565611SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3615754ecacSJeremy L Thompson }
3625754ecacSJeremy L Thompson // -----------------------------------------------------------------------------
363