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 11*c0b5abf0SJeremy L Thompson #include <ceed/types.h> 12*c0b5abf0SJeremy L Thompson #ifndef CEED_RUNNING_JIT_PASS 135754ecacSJeremy L Thompson #include <math.h> 14*c0b5abf0SJeremy L Thompson #endif 155754ecacSJeremy L Thompson 165754ecacSJeremy L Thompson #ifndef PHYSICS_STRUCT 175754ecacSJeremy L Thompson #define PHYSICS_STRUCT 185754ecacSJeremy L Thompson typedef struct Physics_private *Physics; 195754ecacSJeremy L Thompson struct Physics_private { 205754ecacSJeremy L Thompson CeedScalar nu; // Poisson's ratio 215754ecacSJeremy L Thompson CeedScalar E; // Young's Modulus 225754ecacSJeremy L Thompson }; 235754ecacSJeremy L Thompson #endif 245754ecacSJeremy L Thompson 255754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 265754ecacSJeremy L Thompson // Residual evaluation for linear elasticity 275754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 28c8565611SJeremy L Thompson CEED_QFUNCTION(ElasResidual_Linear)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 295754ecacSJeremy L Thompson // Inputs 302b730f8bSJeremy 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]; 315754ecacSJeremy L Thompson 325754ecacSJeremy L Thompson // Outputs 335754ecacSJeremy L Thompson CeedScalar(*dvdX)[3][CEED_Q_VLA] = (CeedScalar(*)[3][CEED_Q_VLA])out[0]; 345754ecacSJeremy L Thompson // grad_u not used for linear elasticity 355754ecacSJeremy L Thompson // (*grad_u)[3][CEED_Q_VLA] = (CeedScalar(*)[3][CEED_Q_VLA])out[1]; 365754ecacSJeremy L Thompson 375754ecacSJeremy L Thompson // Context 385754ecacSJeremy L Thompson const Physics context = (Physics)ctx; 395754ecacSJeremy L Thompson const CeedScalar E = context->E; 405754ecacSJeremy L Thompson const CeedScalar nu = context->nu; 415754ecacSJeremy L Thompson 425754ecacSJeremy L Thompson // Quadrature Point Loop 432b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 445754ecacSJeremy L Thompson // Read spatial derivatives of u 452b730f8bSJeremy L Thompson const CeedScalar du[3][3] = { 462b730f8bSJeremy L Thompson {ug[0][0][i], ug[1][0][i], ug[2][0][i]}, 472b730f8bSJeremy L Thompson {ug[0][1][i], ug[1][1][i], ug[2][1][i]}, 482b730f8bSJeremy L Thompson {ug[0][2][i], ug[1][2][i], ug[2][2][i]} 495754ecacSJeremy L Thompson }; 505754ecacSJeremy L Thompson // -- Qdata 515754ecacSJeremy L Thompson const CeedScalar wdetJ = q_data[0][i]; 522b730f8bSJeremy L Thompson const CeedScalar dXdx[3][3] = { 532b730f8bSJeremy L Thompson {q_data[1][i], q_data[2][i], q_data[3][i]}, 542b730f8bSJeremy L Thompson {q_data[4][i], q_data[5][i], q_data[6][i]}, 552b730f8bSJeremy L Thompson {q_data[7][i], q_data[8][i], q_data[9][i]} 565754ecacSJeremy L Thompson }; 575754ecacSJeremy L Thompson 585754ecacSJeremy L Thompson // Compute grad_u 595754ecacSJeremy L Thompson // dXdx = (dx/dX)^(-1) 605754ecacSJeremy L Thompson // Apply dXdx to du = grad_u 615754ecacSJeremy L Thompson CeedScalar grad_u[3][3]; 622b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { // Component 635754ecacSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { // Derivative 645754ecacSJeremy L Thompson grad_u[j][k] = 0; 652b730f8bSJeremy L Thompson for (CeedInt m = 0; m < 3; m++) grad_u[j][k] += dXdx[m][k] * du[j][m]; 662b730f8bSJeremy L Thompson } 675754ecacSJeremy L Thompson } 685754ecacSJeremy L Thompson 695754ecacSJeremy L Thompson // Compute Strain : e (epsilon) 705754ecacSJeremy L Thompson // e = 1/2 (grad u + (grad u)^T) 715754ecacSJeremy L Thompson 722b730f8bSJeremy L Thompson const CeedScalar e[3][3] = { 732b730f8bSJeremy 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.}, 742b730f8bSJeremy 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.}, 752b730f8bSJeremy 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.} 765754ecacSJeremy L Thompson }; 775754ecacSJeremy L Thompson 785754ecacSJeremy L Thompson // 795754ecacSJeremy L Thompson // Formulation uses Voigt notation: 805754ecacSJeremy L Thompson // stress (sigma) strain (epsilon) 815754ecacSJeremy L Thompson // 825754ecacSJeremy L Thompson // [sigma00] [e00] 835754ecacSJeremy L Thompson // [sigma11] [e11] 845754ecacSJeremy L Thompson // [sigma22] = S * [e22] 855754ecacSJeremy L Thompson // [sigma12] [e12] 865754ecacSJeremy L Thompson // [sigma02] [e02] 875754ecacSJeremy L Thompson // [sigma01] [e01] 885754ecacSJeremy L Thompson // 895754ecacSJeremy L Thompson // where 905754ecacSJeremy L Thompson // [1-nu nu nu ] 915754ecacSJeremy L Thompson // [ nu 1-nu nu ] 925754ecacSJeremy L Thompson // [ nu nu 1-nu ] 935754ecacSJeremy L Thompson // S = E/((1+nu)*(1-2*nu)) [ (1-2*nu)/2 ] 945754ecacSJeremy L Thompson // [ (1-2*nu)/2 ] 955754ecacSJeremy L Thompson // [ (1-2*nu)/2 ] 965754ecacSJeremy L Thompson 975754ecacSJeremy L Thompson // Above Voigt Notation is placed in a 3x3 matrix: 985754ecacSJeremy L Thompson const CeedScalar ss = E / ((1 + nu) * (1 - 2 * nu)); 995754ecacSJeremy L Thompson const CeedScalar sigma00 = ss * ((1 - nu) * e[0][0] + nu * e[1][1] + nu * e[2][2]), 1005754ecacSJeremy L Thompson sigma11 = ss * (nu * e[0][0] + (1 - nu) * e[1][1] + nu * e[2][2]), 1012b730f8bSJeremy 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, 1022b730f8bSJeremy L Thompson sigma02 = ss * (1 - 2 * nu) * e[0][2] * 0.5, sigma01 = ss * (1 - 2 * nu) * e[0][1] * 0.5; 1032b730f8bSJeremy L Thompson const CeedScalar sigma[3][3] = { 1042b730f8bSJeremy L Thompson {sigma00, sigma01, sigma02}, 1055754ecacSJeremy L Thompson {sigma01, sigma11, sigma12}, 1065754ecacSJeremy L Thompson {sigma02, sigma12, sigma22} 1075754ecacSJeremy L Thompson }; 1085754ecacSJeremy L Thompson 1095754ecacSJeremy L Thompson // Apply dXdx^T and weight to sigma 1102b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { // Component 1115754ecacSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { // Derivative 1125754ecacSJeremy L Thompson dvdX[k][j][i] = 0; 1132b730f8bSJeremy L Thompson for (CeedInt m = 0; m < 3; m++) dvdX[k][j][i] += dXdx[k][m] * sigma[j][m] * wdetJ; 1145754ecacSJeremy L Thompson } 1152b730f8bSJeremy L Thompson } 1165754ecacSJeremy L Thompson } // End of Quadrature Point Loop 1175754ecacSJeremy L Thompson 118c8565611SJeremy L Thompson return CEED_ERROR_SUCCESS; 1195754ecacSJeremy L Thompson } 1205754ecacSJeremy L Thompson 1215754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 1225754ecacSJeremy L Thompson // Jacobian evaluation for linear elasticity 1235754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 124c8565611SJeremy L Thompson CEED_QFUNCTION(ElasJacobian_Linear)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 1255754ecacSJeremy L Thompson // Inputs 1265754ecacSJeremy L Thompson const CeedScalar(*deltaug)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[0], 1275754ecacSJeremy L Thompson (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 1285754ecacSJeremy L Thompson // grad_u not used for linear elasticity 1295754ecacSJeremy L Thompson // (*grad_u)[3][Q] = (CeedScalar(*)[3][Q])in[2]; 1305754ecacSJeremy L Thompson 1315754ecacSJeremy L Thompson // Outputs 1325754ecacSJeremy L Thompson CeedScalar(*deltadvdX)[3][CEED_Q_VLA] = (CeedScalar(*)[3][CEED_Q_VLA])out[0]; 1335754ecacSJeremy L Thompson 1345754ecacSJeremy L Thompson // Context 1355754ecacSJeremy L Thompson const Physics context = (Physics)ctx; 1365754ecacSJeremy L Thompson const CeedScalar E = context->E; 1375754ecacSJeremy L Thompson const CeedScalar nu = context->nu; 1385754ecacSJeremy L Thompson 1395754ecacSJeremy L Thompson // Quadrature Point Loop 1402b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 1415754ecacSJeremy L Thompson // Read spatial derivatives of u 1422b730f8bSJeremy L Thompson const CeedScalar deltadu[3][3] = { 1432b730f8bSJeremy L Thompson {deltaug[0][0][i], deltaug[1][0][i], deltaug[2][0][i]}, 1442b730f8bSJeremy L Thompson {deltaug[0][1][i], deltaug[1][1][i], deltaug[2][1][i]}, 1452b730f8bSJeremy L Thompson {deltaug[0][2][i], deltaug[1][2][i], deltaug[2][2][i]} 1465754ecacSJeremy L Thompson }; 1475754ecacSJeremy L Thompson // -- Qdata 1485754ecacSJeremy L Thompson const CeedScalar wdetJ = q_data[0][i]; 1492b730f8bSJeremy L Thompson const CeedScalar dXdx[3][3] = { 1502b730f8bSJeremy L Thompson {q_data[1][i], q_data[2][i], q_data[3][i]}, 1512b730f8bSJeremy L Thompson {q_data[4][i], q_data[5][i], q_data[6][i]}, 1522b730f8bSJeremy L Thompson {q_data[7][i], q_data[8][i], q_data[9][i]} 1535754ecacSJeremy L Thompson }; 1545754ecacSJeremy L Thompson 1555754ecacSJeremy L Thompson // Compute graddeltau 1565754ecacSJeremy L Thompson // dXdx = (dx/dX)^(-1) 1575754ecacSJeremy L Thompson // Apply dXdx to deltadu = graddeltau 1585754ecacSJeremy L Thompson CeedScalar graddeltau[3][3]; 1592b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { // Component 1605754ecacSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { // Derivative 1615754ecacSJeremy L Thompson graddeltau[j][k] = 0; 1622b730f8bSJeremy L Thompson for (CeedInt m = 0; m < 3; m++) graddeltau[j][k] += dXdx[m][k] * deltadu[j][m]; 1632b730f8bSJeremy L Thompson } 1645754ecacSJeremy L Thompson } 1655754ecacSJeremy L Thompson 1665754ecacSJeremy L Thompson // Compute Strain : e (epsilon) 1675754ecacSJeremy L Thompson // e = 1/2 (grad u + (grad u)^T) 1682b730f8bSJeremy L Thompson const CeedScalar de[3][3] = { 1692b730f8bSJeremy L Thompson {(graddeltau[0][0] + graddeltau[0][0]) / 2., (graddeltau[0][1] + graddeltau[1][0]) / 2., (graddeltau[0][2] + graddeltau[2][0]) / 2.}, 1702b730f8bSJeremy L Thompson {(graddeltau[1][0] + graddeltau[0][1]) / 2., (graddeltau[1][1] + graddeltau[1][1]) / 2., (graddeltau[1][2] + graddeltau[2][1]) / 2.}, 1712b730f8bSJeremy L Thompson {(graddeltau[2][0] + graddeltau[0][2]) / 2., (graddeltau[2][1] + graddeltau[1][2]) / 2., (graddeltau[2][2] + graddeltau[2][2]) / 2.} 1725754ecacSJeremy L Thompson }; 1735754ecacSJeremy L Thompson 1745754ecacSJeremy L Thompson // Formulation uses Voigt notation: 1755754ecacSJeremy L Thompson // stress (sigma) strain (epsilon) 1765754ecacSJeremy L Thompson // 1775754ecacSJeremy L Thompson // [dsigma00] [de00] 1785754ecacSJeremy L Thompson // [dsigma11] [de11] 1795754ecacSJeremy L Thompson // [dsigma22] = S * [de22] 1805754ecacSJeremy L Thompson // [dsigma12] [de12] 1815754ecacSJeremy L Thompson // [dsigma02] [de02] 1825754ecacSJeremy L Thompson // [dsigma01] [de01] 1835754ecacSJeremy L Thompson // 1845754ecacSJeremy L Thompson // where 1855754ecacSJeremy L Thompson // [1-nu nu nu ] 1865754ecacSJeremy L Thompson // [ nu 1-nu nu ] 1875754ecacSJeremy L Thompson // [ nu nu 1-nu ] 1885754ecacSJeremy L Thompson // S = E/((1+nu)*(1-2*nu)) [ (1-2*nu)/2 ] 1895754ecacSJeremy L Thompson // [ (1-2*nu)/2 ] 1905754ecacSJeremy L Thompson // [ (1-2*nu)/2 ] 1915754ecacSJeremy L Thompson 1925754ecacSJeremy L Thompson // Above Voigt Notation is placed in a 3x3 matrix: 1935754ecacSJeremy L Thompson const CeedScalar ss = E / ((1 + nu) * (1 - 2 * nu)); 1945754ecacSJeremy L Thompson const CeedScalar dsigma00 = ss * ((1 - nu) * de[0][0] + nu * de[1][1] + nu * de[2][2]), 1955754ecacSJeremy L Thompson dsigma11 = ss * (nu * de[0][0] + (1 - nu) * de[1][1] + nu * de[2][2]), 1962b730f8bSJeremy 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, 1972b730f8bSJeremy L Thompson dsigma02 = ss * (1 - 2 * nu) * de[0][2] / 2, dsigma01 = ss * (1 - 2 * nu) * de[0][1] / 2; 1982b730f8bSJeremy L Thompson const CeedScalar dsigma[3][3] = { 1992b730f8bSJeremy L Thompson {dsigma00, dsigma01, dsigma02}, 2005754ecacSJeremy L Thompson {dsigma01, dsigma11, dsigma12}, 2015754ecacSJeremy L Thompson {dsigma02, dsigma12, dsigma22} 2025754ecacSJeremy L Thompson }; 2035754ecacSJeremy L Thompson 2045754ecacSJeremy L Thompson // Apply dXdx^T and weight 2052b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { // Component 2065754ecacSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { // Derivative 2075754ecacSJeremy L Thompson deltadvdX[k][j][i] = 0; 2082b730f8bSJeremy L Thompson for (CeedInt m = 0; m < 3; m++) deltadvdX[k][j][i] += dXdx[k][m] * dsigma[j][m] * wdetJ; 2095754ecacSJeremy L Thompson } 2102b730f8bSJeremy L Thompson } 2115754ecacSJeremy L Thompson } // End of Quadrature Point Loop 2125754ecacSJeremy L Thompson 213c8565611SJeremy L Thompson return CEED_ERROR_SUCCESS; 2145754ecacSJeremy L Thompson } 2155754ecacSJeremy L Thompson 2165754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 2175754ecacSJeremy L Thompson // Strain energy computation for linear elasticity 2185754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 219c8565611SJeremy L Thompson CEED_QFUNCTION(ElasEnergy_Linear)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 2205754ecacSJeremy L Thompson // Inputs 2212b730f8bSJeremy 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]; 2225754ecacSJeremy L Thompson 2235754ecacSJeremy L Thompson // Outputs 2245754ecacSJeremy L Thompson CeedScalar(*energy) = (CeedScalar(*))out[0]; 2255754ecacSJeremy L Thompson 2265754ecacSJeremy L Thompson // Context 2275754ecacSJeremy L Thompson const Physics context = (Physics)ctx; 2285754ecacSJeremy L Thompson const CeedScalar E = context->E; 2295754ecacSJeremy L Thompson const CeedScalar nu = context->nu; 2305754ecacSJeremy L Thompson 2315754ecacSJeremy L Thompson // Constants 2325754ecacSJeremy L Thompson const CeedScalar TwoMu = E / (1 + nu); 2335754ecacSJeremy L Thompson const CeedScalar mu = TwoMu / 2; 2345754ecacSJeremy L Thompson const CeedScalar Kbulk = E / (3 * (1 - 2 * nu)); // Bulk Modulus 2355754ecacSJeremy L Thompson const CeedScalar lambda = (3 * Kbulk - TwoMu) / 3; 2365754ecacSJeremy L Thompson 2375754ecacSJeremy L Thompson // Quadrature Point Loop 2382b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 2395754ecacSJeremy L Thompson // Read spatial derivatives of u 2402b730f8bSJeremy L Thompson const CeedScalar du[3][3] = { 2412b730f8bSJeremy L Thompson {ug[0][0][i], ug[1][0][i], ug[2][0][i]}, 2422b730f8bSJeremy L Thompson {ug[0][1][i], ug[1][1][i], ug[2][1][i]}, 2432b730f8bSJeremy L Thompson {ug[0][2][i], ug[1][2][i], ug[2][2][i]} 2445754ecacSJeremy L Thompson }; 2455754ecacSJeremy L Thompson // -- Qdata 2465754ecacSJeremy L Thompson const CeedScalar wdetJ = q_data[0][i]; 2472b730f8bSJeremy L Thompson const CeedScalar dXdx[3][3] = { 2482b730f8bSJeremy L Thompson {q_data[1][i], q_data[2][i], q_data[3][i]}, 2492b730f8bSJeremy L Thompson {q_data[4][i], q_data[5][i], q_data[6][i]}, 2502b730f8bSJeremy L Thompson {q_data[7][i], q_data[8][i], q_data[9][i]} 2515754ecacSJeremy L Thompson }; 2525754ecacSJeremy L Thompson 2535754ecacSJeremy L Thompson // Compute grad_u 2545754ecacSJeremy L Thompson // dXdx = (dx/dX)^(-1) 2555754ecacSJeremy L Thompson // Apply dXdx to du = grad_u 2565754ecacSJeremy L Thompson CeedScalar grad_u[3][3]; 2572b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { // Component 2585754ecacSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { // Derivative 2595754ecacSJeremy L Thompson grad_u[j][k] = 0; 2602b730f8bSJeremy L Thompson for (CeedInt m = 0; m < 3; m++) grad_u[j][k] += dXdx[m][k] * du[j][m]; 2612b730f8bSJeremy L Thompson } 2625754ecacSJeremy L Thompson } 2635754ecacSJeremy L Thompson 2645754ecacSJeremy L Thompson // Compute Strain : e (epsilon) 2655754ecacSJeremy L Thompson // e = 1/2 (grad u + (grad u)^T) 2665754ecacSJeremy L Thompson 2672b730f8bSJeremy L Thompson const CeedScalar e[3][3] = { 2682b730f8bSJeremy 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.}, 2692b730f8bSJeremy 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.}, 2702b730f8bSJeremy 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.} 2715754ecacSJeremy L Thompson }; 2725754ecacSJeremy L Thompson 2735754ecacSJeremy L Thompson // Strain energy 2745754ecacSJeremy L Thompson const CeedScalar strain_vol = e[0][0] + e[1][1] + e[2][2]; 2752b730f8bSJeremy L Thompson energy[i] = 2762b730f8bSJeremy 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; 2775754ecacSJeremy L Thompson 2785754ecacSJeremy L Thompson } // End of Quadrature Point Loop 2795754ecacSJeremy L Thompson 280c8565611SJeremy L Thompson return CEED_ERROR_SUCCESS; 2815754ecacSJeremy L Thompson } 2825754ecacSJeremy L Thompson 2835754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 2845754ecacSJeremy L Thompson // Nodal diagnostic quantities for linear elasticity 2855754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 286c8565611SJeremy L Thompson CEED_QFUNCTION(ElasDiagnostic_Linear)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 2875754ecacSJeremy L Thompson // Inputs 2882b730f8bSJeremy 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], 2895754ecacSJeremy L Thompson (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 2905754ecacSJeremy L Thompson 2915754ecacSJeremy L Thompson // Outputs 2925754ecacSJeremy L Thompson CeedScalar(*diagnostic)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 2935754ecacSJeremy L Thompson 2945754ecacSJeremy L Thompson // Context 2955754ecacSJeremy L Thompson const Physics context = (Physics)ctx; 2965754ecacSJeremy L Thompson const CeedScalar E = context->E; 2975754ecacSJeremy L Thompson const CeedScalar nu = context->nu; 2985754ecacSJeremy L Thompson 2995754ecacSJeremy L Thompson // Constants 3005754ecacSJeremy L Thompson const CeedScalar TwoMu = E / (1 + nu); 3015754ecacSJeremy L Thompson const CeedScalar mu = TwoMu / 2; 3025754ecacSJeremy L Thompson const CeedScalar Kbulk = E / (3 * (1 - 2 * nu)); // Bulk Modulus 3035754ecacSJeremy L Thompson const CeedScalar lambda = (3 * Kbulk - TwoMu) / 3; 3045754ecacSJeremy L Thompson 3055754ecacSJeremy L Thompson // Quadrature Point Loop 3062b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 3075754ecacSJeremy L Thompson // Read spatial derivatives of u 3082b730f8bSJeremy L Thompson const CeedScalar du[3][3] = { 3092b730f8bSJeremy L Thompson {ug[0][0][i], ug[1][0][i], ug[2][0][i]}, 3102b730f8bSJeremy L Thompson {ug[0][1][i], ug[1][1][i], ug[2][1][i]}, 3112b730f8bSJeremy L Thompson {ug[0][2][i], ug[1][2][i], ug[2][2][i]} 3125754ecacSJeremy L Thompson }; 3135754ecacSJeremy L Thompson // -- Qdata 3142b730f8bSJeremy L Thompson const CeedScalar dXdx[3][3] = { 3152b730f8bSJeremy L Thompson {q_data[1][i], q_data[2][i], q_data[3][i]}, 3162b730f8bSJeremy L Thompson {q_data[4][i], q_data[5][i], q_data[6][i]}, 3172b730f8bSJeremy L Thompson {q_data[7][i], q_data[8][i], q_data[9][i]} 3185754ecacSJeremy L Thompson }; 3195754ecacSJeremy L Thompson 3205754ecacSJeremy L Thompson // Compute grad_u 3215754ecacSJeremy L Thompson // dXdx = (dx/dX)^(-1) 3225754ecacSJeremy L Thompson // Apply dXdx to du = grad_u 3235754ecacSJeremy L Thompson CeedScalar grad_u[3][3]; 3242b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { // Component 3255754ecacSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { // Derivative 3265754ecacSJeremy L Thompson grad_u[j][k] = 0; 3272b730f8bSJeremy L Thompson for (CeedInt m = 0; m < 3; m++) grad_u[j][k] += dXdx[m][k] * du[j][m]; 3282b730f8bSJeremy L Thompson } 3295754ecacSJeremy L Thompson } 3305754ecacSJeremy L Thompson 3315754ecacSJeremy L Thompson // Compute Strain : e (epsilon) 3325754ecacSJeremy L Thompson // e = 1/2 (grad u + (grad u)^T) 3335754ecacSJeremy L Thompson 3342b730f8bSJeremy L Thompson const CeedScalar e[3][3] = { 3352b730f8bSJeremy 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.}, 3362b730f8bSJeremy 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.}, 3372b730f8bSJeremy 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.} 3385754ecacSJeremy L Thompson }; 3395754ecacSJeremy L Thompson 3405754ecacSJeremy L Thompson // Displacement 3415754ecacSJeremy L Thompson diagnostic[0][i] = u[0][i]; 3425754ecacSJeremy L Thompson diagnostic[1][i] = u[1][i]; 3435754ecacSJeremy L Thompson diagnostic[2][i] = u[2][i]; 3445754ecacSJeremy L Thompson 3455754ecacSJeremy L Thompson // Pressure 3465754ecacSJeremy L Thompson const CeedScalar strain_vol = e[0][0] + e[1][1] + e[2][2]; 3475754ecacSJeremy L Thompson diagnostic[3][i] = -lambda * strain_vol; 3485754ecacSJeremy L Thompson 3495754ecacSJeremy L Thompson // Stress tensor invariants 3505754ecacSJeremy L Thompson diagnostic[4][i] = strain_vol; 3515754ecacSJeremy L Thompson diagnostic[5][i] = 0.; 3522b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 3532b730f8bSJeremy L Thompson for (CeedInt m = 0; m < 3; m++) diagnostic[5][i] += e[j][m] * e[m][j]; 3542b730f8bSJeremy L Thompson } 3555754ecacSJeremy L Thompson diagnostic[6][i] = 1 + strain_vol; 3565754ecacSJeremy L Thompson 3575754ecacSJeremy L Thompson // Strain energy 3582b730f8bSJeremy L Thompson diagnostic[7][i] = 3592b730f8bSJeremy 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); 3605754ecacSJeremy L Thompson } // End of Quadrature Point Loop 3615754ecacSJeremy L Thompson 362c8565611SJeremy L Thompson return CEED_ERROR_SUCCESS; 3635754ecacSJeremy L Thompson } 3645754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 365