xref: /libCEED/examples/solids/qfunctions/linear.h (revision 5aed82e4fa97acf4ba24a7f10a35f5303a6798e0)
1*5aed82e4SJeremy 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 
115754ecacSJeremy L Thompson #ifndef ELAS_LINEAR_H
125754ecacSJeremy L Thompson #define ELAS_LINEAR_H
135754ecacSJeremy L Thompson 
14c9c2c079SJeremy L Thompson #include <ceed.h>
155754ecacSJeremy L Thompson #include <math.h>
165754ecacSJeremy L Thompson 
175754ecacSJeremy L Thompson #ifndef PHYSICS_STRUCT
185754ecacSJeremy L Thompson #define PHYSICS_STRUCT
195754ecacSJeremy L Thompson typedef struct Physics_private *Physics;
205754ecacSJeremy L Thompson struct Physics_private {
215754ecacSJeremy L Thompson   CeedScalar nu;  // Poisson's ratio
225754ecacSJeremy L Thompson   CeedScalar E;   // Young's Modulus
235754ecacSJeremy L Thompson };
245754ecacSJeremy L Thompson #endif
255754ecacSJeremy L Thompson 
265754ecacSJeremy L Thompson // -----------------------------------------------------------------------------
275754ecacSJeremy L Thompson // Residual evaluation for linear elasticity
285754ecacSJeremy L Thompson // -----------------------------------------------------------------------------
292b730f8bSJeremy L Thompson CEED_QFUNCTION(ElasLinearF)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
305754ecacSJeremy L Thompson   // Inputs
312b730f8bSJeremy 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];
325754ecacSJeremy L Thompson 
335754ecacSJeremy L Thompson   // Outputs
345754ecacSJeremy L Thompson   CeedScalar(*dvdX)[3][CEED_Q_VLA] = (CeedScalar(*)[3][CEED_Q_VLA])out[0];
355754ecacSJeremy L Thompson   // grad_u not used for linear elasticity
365754ecacSJeremy L Thompson   // (*grad_u)[3][CEED_Q_VLA] = (CeedScalar(*)[3][CEED_Q_VLA])out[1];
375754ecacSJeremy L Thompson 
385754ecacSJeremy L Thompson   // Context
395754ecacSJeremy L Thompson   const Physics    context = (Physics)ctx;
405754ecacSJeremy L Thompson   const CeedScalar E       = context->E;
415754ecacSJeremy L Thompson   const CeedScalar nu      = context->nu;
425754ecacSJeremy L Thompson 
435754ecacSJeremy L Thompson   // Quadrature Point Loop
442b730f8bSJeremy L Thompson   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
455754ecacSJeremy L Thompson     // Read spatial derivatives of u
462b730f8bSJeremy L Thompson     const CeedScalar du[3][3] = {
472b730f8bSJeremy L Thompson         {ug[0][0][i], ug[1][0][i], ug[2][0][i]},
482b730f8bSJeremy L Thompson         {ug[0][1][i], ug[1][1][i], ug[2][1][i]},
492b730f8bSJeremy L Thompson         {ug[0][2][i], ug[1][2][i], ug[2][2][i]}
505754ecacSJeremy L Thompson     };
515754ecacSJeremy L Thompson     // -- Qdata
525754ecacSJeremy L Thompson     const CeedScalar wdetJ      = q_data[0][i];
532b730f8bSJeremy L Thompson     const CeedScalar dXdx[3][3] = {
542b730f8bSJeremy L Thompson         {q_data[1][i], q_data[2][i], q_data[3][i]},
552b730f8bSJeremy L Thompson         {q_data[4][i], q_data[5][i], q_data[6][i]},
562b730f8bSJeremy L Thompson         {q_data[7][i], q_data[8][i], q_data[9][i]}
575754ecacSJeremy L Thompson     };
585754ecacSJeremy L Thompson 
595754ecacSJeremy L Thompson     // Compute grad_u
605754ecacSJeremy L Thompson     //   dXdx = (dx/dX)^(-1)
615754ecacSJeremy L Thompson     // Apply dXdx to du = grad_u
625754ecacSJeremy L Thompson     CeedScalar grad_u[3][3];
632b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {    // Component
645754ecacSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) {  // Derivative
655754ecacSJeremy L Thompson         grad_u[j][k] = 0;
662b730f8bSJeremy L Thompson         for (CeedInt m = 0; m < 3; m++) grad_u[j][k] += dXdx[m][k] * du[j][m];
672b730f8bSJeremy L Thompson       }
685754ecacSJeremy L Thompson     }
695754ecacSJeremy L Thompson 
705754ecacSJeremy L Thompson     // Compute Strain : e (epsilon)
715754ecacSJeremy L Thompson     // e = 1/2 (grad u + (grad u)^T)
725754ecacSJeremy L Thompson 
732b730f8bSJeremy L Thompson     const CeedScalar e[3][3] = {
742b730f8bSJeremy 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.},
752b730f8bSJeremy 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.},
762b730f8bSJeremy 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.}
775754ecacSJeremy L Thompson     };
785754ecacSJeremy L Thompson 
795754ecacSJeremy L Thompson     //
805754ecacSJeremy L Thompson     // Formulation uses Voigt notation:
815754ecacSJeremy L Thompson     //  stress (sigma)      strain (epsilon)
825754ecacSJeremy L Thompson     //
835754ecacSJeremy L Thompson     //    [sigma00]             [e00]
845754ecacSJeremy L Thompson     //    [sigma11]             [e11]
855754ecacSJeremy L Thompson     //    [sigma22]   =  S   *  [e22]
865754ecacSJeremy L Thompson     //    [sigma12]             [e12]
875754ecacSJeremy L Thompson     //    [sigma02]             [e02]
885754ecacSJeremy L Thompson     //    [sigma01]             [e01]
895754ecacSJeremy L Thompson     //
905754ecacSJeremy L Thompson     //        where
915754ecacSJeremy L Thompson     //                         [1-nu   nu    nu                                    ]
925754ecacSJeremy L Thompson     //                         [ nu   1-nu   nu                                    ]
935754ecacSJeremy L Thompson     //                         [ nu    nu   1-nu                                   ]
945754ecacSJeremy L Thompson     // S = E/((1+nu)*(1-2*nu)) [                  (1-2*nu)/2                       ]
955754ecacSJeremy L Thompson     //                         [                             (1-2*nu)/2            ]
965754ecacSJeremy L Thompson     //                         [                                        (1-2*nu)/2 ]
975754ecacSJeremy L Thompson 
985754ecacSJeremy L Thompson     // Above Voigt Notation is placed in a 3x3 matrix:
995754ecacSJeremy L Thompson     const CeedScalar ss      = E / ((1 + nu) * (1 - 2 * nu));
1005754ecacSJeremy L Thompson     const CeedScalar sigma00 = ss * ((1 - nu) * e[0][0] + nu * e[1][1] + nu * e[2][2]),
1015754ecacSJeremy L Thompson                      sigma11 = ss * (nu * e[0][0] + (1 - nu) * e[1][1] + nu * e[2][2]),
1022b730f8bSJeremy 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,
1032b730f8bSJeremy L Thompson                      sigma02 = ss * (1 - 2 * nu) * e[0][2] * 0.5, sigma01 = ss * (1 - 2 * nu) * e[0][1] * 0.5;
1042b730f8bSJeremy L Thompson     const CeedScalar sigma[3][3] = {
1052b730f8bSJeremy L Thompson         {sigma00, sigma01, sigma02},
1065754ecacSJeremy L Thompson         {sigma01, sigma11, sigma12},
1075754ecacSJeremy L Thompson         {sigma02, sigma12, sigma22}
1085754ecacSJeremy L Thompson     };
1095754ecacSJeremy L Thompson 
1105754ecacSJeremy L Thompson     // Apply dXdx^T and weight to sigma
1112b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {    // Component
1125754ecacSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) {  // Derivative
1135754ecacSJeremy L Thompson         dvdX[k][j][i] = 0;
1142b730f8bSJeremy L Thompson         for (CeedInt m = 0; m < 3; m++) dvdX[k][j][i] += dXdx[k][m] * sigma[j][m] * wdetJ;
1155754ecacSJeremy L Thompson       }
1162b730f8bSJeremy L Thompson     }
1175754ecacSJeremy L Thompson   }  // End of Quadrature Point Loop
1185754ecacSJeremy L Thompson 
1195754ecacSJeremy L Thompson   return 0;
1205754ecacSJeremy L Thompson }
1215754ecacSJeremy L Thompson 
1225754ecacSJeremy L Thompson // -----------------------------------------------------------------------------
1235754ecacSJeremy L Thompson // Jacobian evaluation for linear elasticity
1245754ecacSJeremy L Thompson // -----------------------------------------------------------------------------
1252b730f8bSJeremy L Thompson CEED_QFUNCTION(ElasLineardF)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
1265754ecacSJeremy L Thompson   // Inputs
1275754ecacSJeremy L Thompson   const CeedScalar(*deltaug)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[0],
1285754ecacSJeremy L Thompson         (*q_data)[CEED_Q_VLA]               = (const CeedScalar(*)[CEED_Q_VLA])in[1];
1295754ecacSJeremy L Thompson   // grad_u not used for linear elasticity
1305754ecacSJeremy L Thompson   // (*grad_u)[3][Q] = (CeedScalar(*)[3][Q])in[2];
1315754ecacSJeremy L Thompson 
1325754ecacSJeremy L Thompson   // Outputs
1335754ecacSJeremy L Thompson   CeedScalar(*deltadvdX)[3][CEED_Q_VLA] = (CeedScalar(*)[3][CEED_Q_VLA])out[0];
1345754ecacSJeremy L Thompson 
1355754ecacSJeremy L Thompson   // Context
1365754ecacSJeremy L Thompson   const Physics    context = (Physics)ctx;
1375754ecacSJeremy L Thompson   const CeedScalar E       = context->E;
1385754ecacSJeremy L Thompson   const CeedScalar nu      = context->nu;
1395754ecacSJeremy L Thompson 
1405754ecacSJeremy L Thompson   // Quadrature Point Loop
1412b730f8bSJeremy L Thompson   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
1425754ecacSJeremy L Thompson     // Read spatial derivatives of u
1432b730f8bSJeremy L Thompson     const CeedScalar deltadu[3][3] = {
1442b730f8bSJeremy L Thompson         {deltaug[0][0][i], deltaug[1][0][i], deltaug[2][0][i]},
1452b730f8bSJeremy L Thompson         {deltaug[0][1][i], deltaug[1][1][i], deltaug[2][1][i]},
1462b730f8bSJeremy L Thompson         {deltaug[0][2][i], deltaug[1][2][i], deltaug[2][2][i]}
1475754ecacSJeremy L Thompson     };
1485754ecacSJeremy L Thompson     // -- Qdata
1495754ecacSJeremy L Thompson     const CeedScalar wdetJ      = q_data[0][i];
1502b730f8bSJeremy L Thompson     const CeedScalar dXdx[3][3] = {
1512b730f8bSJeremy L Thompson         {q_data[1][i], q_data[2][i], q_data[3][i]},
1522b730f8bSJeremy L Thompson         {q_data[4][i], q_data[5][i], q_data[6][i]},
1532b730f8bSJeremy L Thompson         {q_data[7][i], q_data[8][i], q_data[9][i]}
1545754ecacSJeremy L Thompson     };
1555754ecacSJeremy L Thompson 
1565754ecacSJeremy L Thompson     // Compute graddeltau
1575754ecacSJeremy L Thompson     //   dXdx = (dx/dX)^(-1)
1585754ecacSJeremy L Thompson     // Apply dXdx to deltadu = graddeltau
1595754ecacSJeremy L Thompson     CeedScalar graddeltau[3][3];
1602b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {    // Component
1615754ecacSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) {  // Derivative
1625754ecacSJeremy L Thompson         graddeltau[j][k] = 0;
1632b730f8bSJeremy L Thompson         for (CeedInt m = 0; m < 3; m++) graddeltau[j][k] += dXdx[m][k] * deltadu[j][m];
1642b730f8bSJeremy L Thompson       }
1655754ecacSJeremy L Thompson     }
1665754ecacSJeremy L Thompson 
1675754ecacSJeremy L Thompson     // Compute Strain : e (epsilon)
1685754ecacSJeremy L Thompson     // e = 1/2 (grad u + (grad u)^T)
1692b730f8bSJeremy L Thompson     const CeedScalar de[3][3] = {
1702b730f8bSJeremy L Thompson         {(graddeltau[0][0] + graddeltau[0][0]) / 2., (graddeltau[0][1] + graddeltau[1][0]) / 2., (graddeltau[0][2] + graddeltau[2][0]) / 2.},
1712b730f8bSJeremy L Thompson         {(graddeltau[1][0] + graddeltau[0][1]) / 2., (graddeltau[1][1] + graddeltau[1][1]) / 2., (graddeltau[1][2] + graddeltau[2][1]) / 2.},
1722b730f8bSJeremy L Thompson         {(graddeltau[2][0] + graddeltau[0][2]) / 2., (graddeltau[2][1] + graddeltau[1][2]) / 2., (graddeltau[2][2] + graddeltau[2][2]) / 2.}
1735754ecacSJeremy L Thompson     };
1745754ecacSJeremy L Thompson 
1755754ecacSJeremy L Thompson     // Formulation uses Voigt notation:
1765754ecacSJeremy L Thompson     //  stress (sigma)      strain (epsilon)
1775754ecacSJeremy L Thompson     //
1785754ecacSJeremy L Thompson     //    [dsigma00]             [de00]
1795754ecacSJeremy L Thompson     //    [dsigma11]             [de11]
1805754ecacSJeremy L Thompson     //    [dsigma22]   =  S   *  [de22]
1815754ecacSJeremy L Thompson     //    [dsigma12]             [de12]
1825754ecacSJeremy L Thompson     //    [dsigma02]             [de02]
1835754ecacSJeremy L Thompson     //    [dsigma01]             [de01]
1845754ecacSJeremy L Thompson     //
1855754ecacSJeremy L Thompson     //        where
1865754ecacSJeremy L Thompson     //                         [1-nu   nu    nu                                    ]
1875754ecacSJeremy L Thompson     //                         [ nu   1-nu   nu                                    ]
1885754ecacSJeremy L Thompson     //                         [ nu    nu   1-nu                                   ]
1895754ecacSJeremy L Thompson     // S = E/((1+nu)*(1-2*nu)) [                  (1-2*nu)/2                       ]
1905754ecacSJeremy L Thompson     //                         [                             (1-2*nu)/2            ]
1915754ecacSJeremy L Thompson     //                         [                                        (1-2*nu)/2 ]
1925754ecacSJeremy L Thompson 
1935754ecacSJeremy L Thompson     // Above Voigt Notation is placed in a 3x3 matrix:
1945754ecacSJeremy L Thompson     const CeedScalar ss       = E / ((1 + nu) * (1 - 2 * nu));
1955754ecacSJeremy L Thompson     const CeedScalar dsigma00 = ss * ((1 - nu) * de[0][0] + nu * de[1][1] + nu * de[2][2]),
1965754ecacSJeremy L Thompson                      dsigma11 = ss * (nu * de[0][0] + (1 - nu) * de[1][1] + nu * de[2][2]),
1972b730f8bSJeremy 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,
1982b730f8bSJeremy L Thompson                      dsigma02 = ss * (1 - 2 * nu) * de[0][2] / 2, dsigma01 = ss * (1 - 2 * nu) * de[0][1] / 2;
1992b730f8bSJeremy L Thompson     const CeedScalar dsigma[3][3] = {
2002b730f8bSJeremy L Thompson         {dsigma00, dsigma01, dsigma02},
2015754ecacSJeremy L Thompson         {dsigma01, dsigma11, dsigma12},
2025754ecacSJeremy L Thompson         {dsigma02, dsigma12, dsigma22}
2035754ecacSJeremy L Thompson     };
2045754ecacSJeremy L Thompson 
2055754ecacSJeremy L Thompson     // Apply dXdx^T and weight
2062b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {    // Component
2075754ecacSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) {  // Derivative
2085754ecacSJeremy L Thompson         deltadvdX[k][j][i] = 0;
2092b730f8bSJeremy L Thompson         for (CeedInt m = 0; m < 3; m++) deltadvdX[k][j][i] += dXdx[k][m] * dsigma[j][m] * wdetJ;
2105754ecacSJeremy L Thompson       }
2112b730f8bSJeremy L Thompson     }
2125754ecacSJeremy L Thompson   }  // End of Quadrature Point Loop
2135754ecacSJeremy L Thompson 
2145754ecacSJeremy L Thompson   return 0;
2155754ecacSJeremy L Thompson }
2165754ecacSJeremy L Thompson 
2175754ecacSJeremy L Thompson // -----------------------------------------------------------------------------
2185754ecacSJeremy L Thompson // Strain energy computation for linear elasticity
2195754ecacSJeremy L Thompson // -----------------------------------------------------------------------------
2202b730f8bSJeremy L Thompson CEED_QFUNCTION(ElasLinearEnergy)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
2215754ecacSJeremy L Thompson   // Inputs
2222b730f8bSJeremy 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];
2235754ecacSJeremy L Thompson 
2245754ecacSJeremy L Thompson   // Outputs
2255754ecacSJeremy L Thompson   CeedScalar(*energy) = (CeedScalar(*))out[0];
2265754ecacSJeremy L Thompson 
2275754ecacSJeremy L Thompson   // Context
2285754ecacSJeremy L Thompson   const Physics    context = (Physics)ctx;
2295754ecacSJeremy L Thompson   const CeedScalar E       = context->E;
2305754ecacSJeremy L Thompson   const CeedScalar nu      = context->nu;
2315754ecacSJeremy L Thompson 
2325754ecacSJeremy L Thompson   // Constants
2335754ecacSJeremy L Thompson   const CeedScalar TwoMu  = E / (1 + nu);
2345754ecacSJeremy L Thompson   const CeedScalar mu     = TwoMu / 2;
2355754ecacSJeremy L Thompson   const CeedScalar Kbulk  = E / (3 * (1 - 2 * nu));  // Bulk Modulus
2365754ecacSJeremy L Thompson   const CeedScalar lambda = (3 * Kbulk - TwoMu) / 3;
2375754ecacSJeremy L Thompson 
2385754ecacSJeremy L Thompson   // Quadrature Point Loop
2392b730f8bSJeremy L Thompson   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
2405754ecacSJeremy L Thompson     // Read spatial derivatives of u
2412b730f8bSJeremy L Thompson     const CeedScalar du[3][3] = {
2422b730f8bSJeremy L Thompson         {ug[0][0][i], ug[1][0][i], ug[2][0][i]},
2432b730f8bSJeremy L Thompson         {ug[0][1][i], ug[1][1][i], ug[2][1][i]},
2442b730f8bSJeremy L Thompson         {ug[0][2][i], ug[1][2][i], ug[2][2][i]}
2455754ecacSJeremy L Thompson     };
2465754ecacSJeremy L Thompson     // -- Qdata
2475754ecacSJeremy L Thompson     const CeedScalar wdetJ      = q_data[0][i];
2482b730f8bSJeremy L Thompson     const CeedScalar dXdx[3][3] = {
2492b730f8bSJeremy L Thompson         {q_data[1][i], q_data[2][i], q_data[3][i]},
2502b730f8bSJeremy L Thompson         {q_data[4][i], q_data[5][i], q_data[6][i]},
2512b730f8bSJeremy L Thompson         {q_data[7][i], q_data[8][i], q_data[9][i]}
2525754ecacSJeremy L Thompson     };
2535754ecacSJeremy L Thompson 
2545754ecacSJeremy L Thompson     // Compute grad_u
2555754ecacSJeremy L Thompson     //   dXdx = (dx/dX)^(-1)
2565754ecacSJeremy L Thompson     // Apply dXdx to du = grad_u
2575754ecacSJeremy L Thompson     CeedScalar grad_u[3][3];
2582b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {    // Component
2595754ecacSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) {  // Derivative
2605754ecacSJeremy L Thompson         grad_u[j][k] = 0;
2612b730f8bSJeremy L Thompson         for (CeedInt m = 0; m < 3; m++) grad_u[j][k] += dXdx[m][k] * du[j][m];
2622b730f8bSJeremy L Thompson       }
2635754ecacSJeremy L Thompson     }
2645754ecacSJeremy L Thompson 
2655754ecacSJeremy L Thompson     // Compute Strain : e (epsilon)
2665754ecacSJeremy L Thompson     // e = 1/2 (grad u + (grad u)^T)
2675754ecacSJeremy L Thompson 
2682b730f8bSJeremy L Thompson     const CeedScalar e[3][3] = {
2692b730f8bSJeremy 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.},
2702b730f8bSJeremy 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.},
2712b730f8bSJeremy 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.}
2725754ecacSJeremy L Thompson     };
2735754ecacSJeremy L Thompson 
2745754ecacSJeremy L Thompson     // Strain energy
2755754ecacSJeremy L Thompson     const CeedScalar strain_vol = e[0][0] + e[1][1] + e[2][2];
2762b730f8bSJeremy L Thompson     energy[i] =
2772b730f8bSJeremy 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;
2785754ecacSJeremy L Thompson 
2795754ecacSJeremy L Thompson   }  // End of Quadrature Point Loop
2805754ecacSJeremy L Thompson 
2815754ecacSJeremy L Thompson   return 0;
2825754ecacSJeremy L Thompson }
2835754ecacSJeremy L Thompson 
2845754ecacSJeremy L Thompson // -----------------------------------------------------------------------------
2855754ecacSJeremy L Thompson // Nodal diagnostic quantities for linear elasticity
2865754ecacSJeremy L Thompson // -----------------------------------------------------------------------------
2872b730f8bSJeremy L Thompson CEED_QFUNCTION(ElasLinearDiagnostic)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
2885754ecacSJeremy L Thompson   // Inputs
2892b730f8bSJeremy 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],
2905754ecacSJeremy L Thompson         (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2];
2915754ecacSJeremy L Thompson 
2925754ecacSJeremy L Thompson   // Outputs
2935754ecacSJeremy L Thompson   CeedScalar(*diagnostic)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
2945754ecacSJeremy L Thompson 
2955754ecacSJeremy L Thompson   // Context
2965754ecacSJeremy L Thompson   const Physics    context = (Physics)ctx;
2975754ecacSJeremy L Thompson   const CeedScalar E       = context->E;
2985754ecacSJeremy L Thompson   const CeedScalar nu      = context->nu;
2995754ecacSJeremy L Thompson 
3005754ecacSJeremy L Thompson   // Constants
3015754ecacSJeremy L Thompson   const CeedScalar TwoMu  = E / (1 + nu);
3025754ecacSJeremy L Thompson   const CeedScalar mu     = TwoMu / 2;
3035754ecacSJeremy L Thompson   const CeedScalar Kbulk  = E / (3 * (1 - 2 * nu));  // Bulk Modulus
3045754ecacSJeremy L Thompson   const CeedScalar lambda = (3 * Kbulk - TwoMu) / 3;
3055754ecacSJeremy L Thompson 
3065754ecacSJeremy L Thompson   // Quadrature Point Loop
3072b730f8bSJeremy L Thompson   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
3085754ecacSJeremy L Thompson     // Read spatial derivatives of u
3092b730f8bSJeremy L Thompson     const CeedScalar du[3][3] = {
3102b730f8bSJeremy L Thompson         {ug[0][0][i], ug[1][0][i], ug[2][0][i]},
3112b730f8bSJeremy L Thompson         {ug[0][1][i], ug[1][1][i], ug[2][1][i]},
3122b730f8bSJeremy L Thompson         {ug[0][2][i], ug[1][2][i], ug[2][2][i]}
3135754ecacSJeremy L Thompson     };
3145754ecacSJeremy L Thompson     // -- Qdata
3152b730f8bSJeremy L Thompson     const CeedScalar dXdx[3][3] = {
3162b730f8bSJeremy L Thompson         {q_data[1][i], q_data[2][i], q_data[3][i]},
3172b730f8bSJeremy L Thompson         {q_data[4][i], q_data[5][i], q_data[6][i]},
3182b730f8bSJeremy L Thompson         {q_data[7][i], q_data[8][i], q_data[9][i]}
3195754ecacSJeremy L Thompson     };
3205754ecacSJeremy L Thompson 
3215754ecacSJeremy L Thompson     // Compute grad_u
3225754ecacSJeremy L Thompson     //   dXdx = (dx/dX)^(-1)
3235754ecacSJeremy L Thompson     // Apply dXdx to du = grad_u
3245754ecacSJeremy L Thompson     CeedScalar grad_u[3][3];
3252b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {    // Component
3265754ecacSJeremy L Thompson       for (CeedInt k = 0; k < 3; k++) {  // Derivative
3275754ecacSJeremy L Thompson         grad_u[j][k] = 0;
3282b730f8bSJeremy L Thompson         for (CeedInt m = 0; m < 3; m++) grad_u[j][k] += dXdx[m][k] * du[j][m];
3292b730f8bSJeremy L Thompson       }
3305754ecacSJeremy L Thompson     }
3315754ecacSJeremy L Thompson 
3325754ecacSJeremy L Thompson     // Compute Strain : e (epsilon)
3335754ecacSJeremy L Thompson     // e = 1/2 (grad u + (grad u)^T)
3345754ecacSJeremy L Thompson 
3352b730f8bSJeremy L Thompson     const CeedScalar e[3][3] = {
3362b730f8bSJeremy 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.},
3372b730f8bSJeremy 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.},
3382b730f8bSJeremy 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.}
3395754ecacSJeremy L Thompson     };
3405754ecacSJeremy L Thompson 
3415754ecacSJeremy L Thompson     // Displacement
3425754ecacSJeremy L Thompson     diagnostic[0][i] = u[0][i];
3435754ecacSJeremy L Thompson     diagnostic[1][i] = u[1][i];
3445754ecacSJeremy L Thompson     diagnostic[2][i] = u[2][i];
3455754ecacSJeremy L Thompson 
3465754ecacSJeremy L Thompson     // Pressure
3475754ecacSJeremy L Thompson     const CeedScalar strain_vol = e[0][0] + e[1][1] + e[2][2];
3485754ecacSJeremy L Thompson     diagnostic[3][i]            = -lambda * strain_vol;
3495754ecacSJeremy L Thompson 
3505754ecacSJeremy L Thompson     // Stress tensor invariants
3515754ecacSJeremy L Thompson     diagnostic[4][i] = strain_vol;
3525754ecacSJeremy L Thompson     diagnostic[5][i] = 0.;
3532b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < 3; j++) {
3542b730f8bSJeremy L Thompson       for (CeedInt m = 0; m < 3; m++) diagnostic[5][i] += e[j][m] * e[m][j];
3552b730f8bSJeremy L Thompson     }
3565754ecacSJeremy L Thompson     diagnostic[6][i] = 1 + strain_vol;
3575754ecacSJeremy L Thompson 
3585754ecacSJeremy L Thompson     // Strain energy
3592b730f8bSJeremy L Thompson     diagnostic[7][i] =
3602b730f8bSJeremy 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);
3615754ecacSJeremy L Thompson   }  // End of Quadrature Point Loop
3625754ecacSJeremy L Thompson 
3635754ecacSJeremy L Thompson   return 0;
3645754ecacSJeremy L Thompson }
3655754ecacSJeremy L Thompson // -----------------------------------------------------------------------------
3665754ecacSJeremy L Thompson 
3675754ecacSJeremy L Thompson #endif  // End of ELAS_LINEAR_H
368