13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, 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 14*c9c2c079SJeremy 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 // ----------------------------------------------------------------------------- 295754ecacSJeremy L Thompson CEED_QFUNCTION(ElasLinearF)(void *ctx, CeedInt Q, const CeedScalar *const *in, 305754ecacSJeremy L Thompson CeedScalar *const *out) { 315754ecacSJeremy L Thompson // *INDENT-OFF* 325754ecacSJeremy L Thompson // Inputs 335754ecacSJeremy L Thompson const CeedScalar (*ug)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[0], 345754ecacSJeremy L Thompson (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 355754ecacSJeremy L Thompson 365754ecacSJeremy L Thompson // Outputs 375754ecacSJeremy L Thompson CeedScalar (*dvdX)[3][CEED_Q_VLA] = (CeedScalar(*)[3][CEED_Q_VLA])out[0]; 385754ecacSJeremy L Thompson // grad_u not used for linear elasticity 395754ecacSJeremy L Thompson // (*grad_u)[3][CEED_Q_VLA] = (CeedScalar(*)[3][CEED_Q_VLA])out[1]; 405754ecacSJeremy L Thompson // *INDENT-ON* 415754ecacSJeremy L Thompson 425754ecacSJeremy L Thompson // Context 435754ecacSJeremy L Thompson const Physics context = (Physics)ctx; 445754ecacSJeremy L Thompson const CeedScalar E = context->E; 455754ecacSJeremy L Thompson const CeedScalar nu = context->nu; 465754ecacSJeremy L Thompson 475754ecacSJeremy L Thompson // Quadrature Point Loop 485754ecacSJeremy L Thompson CeedPragmaSIMD 495754ecacSJeremy L Thompson for (CeedInt i=0; i<Q; i++) { 505754ecacSJeremy L Thompson // Read spatial derivatives of u 515754ecacSJeremy L Thompson // *INDENT-OFF* 525754ecacSJeremy L Thompson const CeedScalar du[3][3] = {{ug[0][0][i], 535754ecacSJeremy L Thompson ug[1][0][i], 545754ecacSJeremy L Thompson ug[2][0][i]}, 555754ecacSJeremy L Thompson {ug[0][1][i], 565754ecacSJeremy L Thompson ug[1][1][i], 575754ecacSJeremy L Thompson ug[2][1][i]}, 585754ecacSJeremy L Thompson {ug[0][2][i], 595754ecacSJeremy L Thompson ug[1][2][i], 605754ecacSJeremy L Thompson ug[2][2][i]} 615754ecacSJeremy L Thompson }; 625754ecacSJeremy L Thompson // -- Qdata 635754ecacSJeremy L Thompson const CeedScalar wdetJ = q_data[0][i]; 645754ecacSJeremy L Thompson const CeedScalar dXdx[3][3] = {{q_data[1][i], 655754ecacSJeremy L Thompson q_data[2][i], 665754ecacSJeremy L Thompson q_data[3][i]}, 675754ecacSJeremy L Thompson {q_data[4][i], 685754ecacSJeremy L Thompson q_data[5][i], 695754ecacSJeremy L Thompson q_data[6][i]}, 705754ecacSJeremy L Thompson {q_data[7][i], 715754ecacSJeremy L Thompson q_data[8][i], 725754ecacSJeremy L Thompson q_data[9][i]} 735754ecacSJeremy L Thompson }; 745754ecacSJeremy L Thompson // *INDENT-ON* 755754ecacSJeremy L Thompson 765754ecacSJeremy L Thompson // Compute grad_u 775754ecacSJeremy L Thompson // dXdx = (dx/dX)^(-1) 785754ecacSJeremy L Thompson // Apply dXdx to du = grad_u 795754ecacSJeremy L Thompson CeedScalar grad_u[3][3]; 805754ecacSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) // Component 815754ecacSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { // Derivative 825754ecacSJeremy L Thompson grad_u[j][k] = 0; 835754ecacSJeremy L Thompson for (CeedInt m = 0; m < 3; m++) 845754ecacSJeremy L Thompson grad_u[j][k] += dXdx[m][k] * du[j][m]; 855754ecacSJeremy L Thompson } 865754ecacSJeremy L Thompson 875754ecacSJeremy L Thompson // Compute Strain : e (epsilon) 885754ecacSJeremy L Thompson // e = 1/2 (grad u + (grad u)^T) 895754ecacSJeremy L Thompson 905754ecacSJeremy L Thompson // *INDENT-OFF* 915754ecacSJeremy L Thompson const CeedScalar e[3][3] = {{(grad_u[0][0] + grad_u[0][0])/2., 925754ecacSJeremy L Thompson (grad_u[0][1] + grad_u[1][0])/2., 935754ecacSJeremy L Thompson (grad_u[0][2] + grad_u[2][0])/2.}, 945754ecacSJeremy L Thompson {(grad_u[1][0] + grad_u[0][1])/2., 955754ecacSJeremy L Thompson (grad_u[1][1] + grad_u[1][1])/2., 965754ecacSJeremy L Thompson (grad_u[1][2] + grad_u[2][1])/2.}, 975754ecacSJeremy L Thompson {(grad_u[2][0] + grad_u[0][2])/2., 985754ecacSJeremy L Thompson (grad_u[2][1] + grad_u[1][2])/2., 995754ecacSJeremy L Thompson (grad_u[2][2] + grad_u[2][2])/2.} 1005754ecacSJeremy L Thompson }; 1015754ecacSJeremy L Thompson // *INDENT-ON* 1025754ecacSJeremy L Thompson 1035754ecacSJeremy L Thompson // 1045754ecacSJeremy L Thompson // Formulation uses Voigt notation: 1055754ecacSJeremy L Thompson // stress (sigma) strain (epsilon) 1065754ecacSJeremy L Thompson // 1075754ecacSJeremy L Thompson // [sigma00] [e00] 1085754ecacSJeremy L Thompson // [sigma11] [e11] 1095754ecacSJeremy L Thompson // [sigma22] = S * [e22] 1105754ecacSJeremy L Thompson // [sigma12] [e12] 1115754ecacSJeremy L Thompson // [sigma02] [e02] 1125754ecacSJeremy L Thompson // [sigma01] [e01] 1135754ecacSJeremy L Thompson // 1145754ecacSJeremy L Thompson // where 1155754ecacSJeremy L Thompson // [1-nu nu nu ] 1165754ecacSJeremy L Thompson // [ nu 1-nu nu ] 1175754ecacSJeremy L Thompson // [ nu nu 1-nu ] 1185754ecacSJeremy L Thompson // S = E/((1+nu)*(1-2*nu)) [ (1-2*nu)/2 ] 1195754ecacSJeremy L Thompson // [ (1-2*nu)/2 ] 1205754ecacSJeremy L Thompson // [ (1-2*nu)/2 ] 1215754ecacSJeremy L Thompson 1225754ecacSJeremy L Thompson // Above Voigt Notation is placed in a 3x3 matrix: 1235754ecacSJeremy L Thompson const CeedScalar ss = E / ((1 + nu)*(1 - 2*nu)); 1245754ecacSJeremy L Thompson const CeedScalar sigma00 = ss*((1 - nu)*e[0][0] + nu*e[1][1] + nu*e[2][2]), 1255754ecacSJeremy L Thompson sigma11 = ss*(nu*e[0][0] + (1 - nu)*e[1][1] + nu*e[2][2]), 1265754ecacSJeremy L Thompson sigma22 = ss*(nu*e[0][0] + nu*e[1][1] + (1 - nu)*e[2][2]), 1275754ecacSJeremy L Thompson sigma12 = ss*(1 - 2*nu)*e[1][2]*0.5, 1285754ecacSJeremy L Thompson sigma02 = ss*(1 - 2*nu)*e[0][2]*0.5, 1295754ecacSJeremy L Thompson sigma01 = ss*(1 - 2*nu)*e[0][1]*0.5; 1305754ecacSJeremy L Thompson // *INDENT-OFF* 1315754ecacSJeremy L Thompson const CeedScalar sigma[3][3] = {{sigma00, sigma01, sigma02}, 1325754ecacSJeremy L Thompson {sigma01, sigma11, sigma12}, 1335754ecacSJeremy L Thompson {sigma02, sigma12, sigma22} 1345754ecacSJeremy L Thompson }; 1355754ecacSJeremy L Thompson // *INDENT-ON* 1365754ecacSJeremy L Thompson 1375754ecacSJeremy L Thompson // Apply dXdx^T and weight to sigma 1385754ecacSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) // Component 1395754ecacSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { // Derivative 1405754ecacSJeremy L Thompson dvdX[k][j][i] = 0; 1415754ecacSJeremy L Thompson for (CeedInt m = 0; m < 3; m++) 1425754ecacSJeremy L Thompson dvdX[k][j][i] += dXdx[k][m] * sigma[j][m] * wdetJ; 1435754ecacSJeremy L Thompson } 1445754ecacSJeremy L Thompson 1455754ecacSJeremy L Thompson } // End of Quadrature Point Loop 1465754ecacSJeremy L Thompson 1475754ecacSJeremy L Thompson return 0; 1485754ecacSJeremy L Thompson } 1495754ecacSJeremy L Thompson 1505754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 1515754ecacSJeremy L Thompson // Jacobian evaluation for linear elasticity 1525754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 1535754ecacSJeremy L Thompson CEED_QFUNCTION(ElasLineardF)(void *ctx, CeedInt Q, const CeedScalar *const *in, 1545754ecacSJeremy L Thompson CeedScalar *const *out) { 1555754ecacSJeremy L Thompson // *INDENT-OFF* 1565754ecacSJeremy L Thompson // Inputs 1575754ecacSJeremy L Thompson const CeedScalar (*deltaug)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[0], 1585754ecacSJeremy L Thompson (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 1595754ecacSJeremy L Thompson // grad_u not used for linear elasticity 1605754ecacSJeremy L Thompson // (*grad_u)[3][Q] = (CeedScalar(*)[3][Q])in[2]; 1615754ecacSJeremy L Thompson 1625754ecacSJeremy L Thompson // Outputs 1635754ecacSJeremy L Thompson CeedScalar (*deltadvdX)[3][CEED_Q_VLA] = (CeedScalar(*)[3][CEED_Q_VLA])out[0]; 1645754ecacSJeremy L Thompson // *INDENT-ON* 1655754ecacSJeremy L Thompson 1665754ecacSJeremy L Thompson // Context 1675754ecacSJeremy L Thompson const Physics context = (Physics)ctx; 1685754ecacSJeremy L Thompson const CeedScalar E = context->E; 1695754ecacSJeremy L Thompson const CeedScalar nu = context->nu; 1705754ecacSJeremy L Thompson 1715754ecacSJeremy L Thompson // Quadrature Point Loop 1725754ecacSJeremy L Thompson CeedPragmaSIMD 1735754ecacSJeremy L Thompson for (CeedInt i=0; i<Q; i++) { 1745754ecacSJeremy L Thompson // Read spatial derivatives of u 1755754ecacSJeremy L Thompson // *INDENT-OFF* 1765754ecacSJeremy L Thompson const CeedScalar deltadu[3][3] = {{deltaug[0][0][i], 1775754ecacSJeremy L Thompson deltaug[1][0][i], 1785754ecacSJeremy L Thompson deltaug[2][0][i]}, 1795754ecacSJeremy L Thompson {deltaug[0][1][i], 1805754ecacSJeremy L Thompson deltaug[1][1][i], 1815754ecacSJeremy L Thompson deltaug[2][1][i]}, 1825754ecacSJeremy L Thompson {deltaug[0][2][i], 1835754ecacSJeremy L Thompson deltaug[1][2][i], 1845754ecacSJeremy L Thompson deltaug[2][2][i]} 1855754ecacSJeremy L Thompson }; 1865754ecacSJeremy L Thompson // -- Qdata 1875754ecacSJeremy L Thompson const CeedScalar wdetJ = q_data[0][i]; 1885754ecacSJeremy L Thompson const CeedScalar dXdx[3][3] = {{q_data[1][i], 1895754ecacSJeremy L Thompson q_data[2][i], 1905754ecacSJeremy L Thompson q_data[3][i]}, 1915754ecacSJeremy L Thompson {q_data[4][i], 1925754ecacSJeremy L Thompson q_data[5][i], 1935754ecacSJeremy L Thompson q_data[6][i]}, 1945754ecacSJeremy L Thompson {q_data[7][i], 1955754ecacSJeremy L Thompson q_data[8][i], 1965754ecacSJeremy L Thompson q_data[9][i]} 1975754ecacSJeremy L Thompson }; 1985754ecacSJeremy L Thompson // *INDENT-ON* 1995754ecacSJeremy L Thompson 2005754ecacSJeremy L Thompson // Compute graddeltau 2015754ecacSJeremy L Thompson // dXdx = (dx/dX)^(-1) 2025754ecacSJeremy L Thompson // Apply dXdx to deltadu = graddeltau 2035754ecacSJeremy L Thompson CeedScalar graddeltau[3][3]; 2045754ecacSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) // Component 2055754ecacSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { // Derivative 2065754ecacSJeremy L Thompson graddeltau[j][k] = 0; 2075754ecacSJeremy L Thompson for (CeedInt m = 0; m < 3; m++) 2085754ecacSJeremy L Thompson graddeltau[j][k] += dXdx[m][k] * deltadu[j][m]; 2095754ecacSJeremy L Thompson } 2105754ecacSJeremy L Thompson 2115754ecacSJeremy L Thompson // Compute Strain : e (epsilon) 2125754ecacSJeremy L Thompson // e = 1/2 (grad u + (grad u)^T) 2135754ecacSJeremy L Thompson // *INDENT-OFF* 2145754ecacSJeremy L Thompson const CeedScalar de[3][3] = {{(graddeltau[0][0] + graddeltau[0][0])/2., 2155754ecacSJeremy L Thompson (graddeltau[0][1] + graddeltau[1][0])/2., 2165754ecacSJeremy L Thompson (graddeltau[0][2] + graddeltau[2][0])/2.}, 2175754ecacSJeremy L Thompson {(graddeltau[1][0] + graddeltau[0][1])/2., 2185754ecacSJeremy L Thompson (graddeltau[1][1] + graddeltau[1][1])/2., 2195754ecacSJeremy L Thompson (graddeltau[1][2] + graddeltau[2][1])/2.}, 2205754ecacSJeremy L Thompson {(graddeltau[2][0] + graddeltau[0][2])/2., 2215754ecacSJeremy L Thompson (graddeltau[2][1] + graddeltau[1][2])/2., 2225754ecacSJeremy L Thompson (graddeltau[2][2] + graddeltau[2][2])/2.} 2235754ecacSJeremy L Thompson }; 2245754ecacSJeremy L Thompson 2255754ecacSJeremy L Thompson // *INDENT-ON* 2265754ecacSJeremy L Thompson // Formulation uses Voigt notation: 2275754ecacSJeremy L Thompson // stress (sigma) strain (epsilon) 2285754ecacSJeremy L Thompson // 2295754ecacSJeremy L Thompson // [dsigma00] [de00] 2305754ecacSJeremy L Thompson // [dsigma11] [de11] 2315754ecacSJeremy L Thompson // [dsigma22] = S * [de22] 2325754ecacSJeremy L Thompson // [dsigma12] [de12] 2335754ecacSJeremy L Thompson // [dsigma02] [de02] 2345754ecacSJeremy L Thompson // [dsigma01] [de01] 2355754ecacSJeremy L Thompson // 2365754ecacSJeremy L Thompson // where 2375754ecacSJeremy L Thompson // [1-nu nu nu ] 2385754ecacSJeremy L Thompson // [ nu 1-nu nu ] 2395754ecacSJeremy L Thompson // [ nu nu 1-nu ] 2405754ecacSJeremy L Thompson // S = E/((1+nu)*(1-2*nu)) [ (1-2*nu)/2 ] 2415754ecacSJeremy L Thompson // [ (1-2*nu)/2 ] 2425754ecacSJeremy L Thompson // [ (1-2*nu)/2 ] 2435754ecacSJeremy L Thompson 2445754ecacSJeremy L Thompson // Above Voigt Notation is placed in a 3x3 matrix: 2455754ecacSJeremy L Thompson const CeedScalar ss = E / ((1 + nu)*(1 - 2*nu)); 2465754ecacSJeremy L Thompson const CeedScalar dsigma00 = ss*((1 - nu)*de[0][0]+nu*de[1][1]+nu*de[2][2]), 2475754ecacSJeremy L Thompson dsigma11 = ss*(nu*de[0][0]+(1 - nu)*de[1][1]+nu*de[2][2]), 2485754ecacSJeremy L Thompson dsigma22 = ss*(nu*de[0][0]+nu*de[1][1]+(1 - nu)*de[2][2]), 2495754ecacSJeremy L Thompson dsigma12 = ss*(1 - 2*nu)*de[1][2] / 2, 2505754ecacSJeremy L Thompson dsigma02 = ss*(1 - 2*nu)*de[0][2] / 2, 2515754ecacSJeremy L Thompson dsigma01 = ss*(1 - 2*nu)*de[0][1] / 2; 2525754ecacSJeremy L Thompson // *INDENT-OFF* 2535754ecacSJeremy L Thompson const CeedScalar dsigma[3][3] = {{dsigma00, dsigma01, dsigma02}, 2545754ecacSJeremy L Thompson {dsigma01, dsigma11, dsigma12}, 2555754ecacSJeremy L Thompson {dsigma02, dsigma12, dsigma22} 2565754ecacSJeremy L Thompson }; 2575754ecacSJeremy L Thompson // *INDENT-ON* 2585754ecacSJeremy L Thompson 2595754ecacSJeremy L Thompson // Apply dXdx^T and weight 2605754ecacSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) // Component 2615754ecacSJeremy L Thompson for (CeedInt k = 0; k < 3 ; k++) { // Derivative 2625754ecacSJeremy L Thompson deltadvdX[k][j][i] = 0; 2635754ecacSJeremy L Thompson for (CeedInt m = 0; m < 3; m++) 2645754ecacSJeremy L Thompson deltadvdX[k][j][i] += dXdx[k][m] * dsigma[j][m] * wdetJ; 2655754ecacSJeremy L Thompson } 2665754ecacSJeremy L Thompson 2675754ecacSJeremy L Thompson } // End of Quadrature Point Loop 2685754ecacSJeremy L Thompson 2695754ecacSJeremy L Thompson return 0; 2705754ecacSJeremy L Thompson } 2715754ecacSJeremy L Thompson 2725754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 2735754ecacSJeremy L Thompson // Strain energy computation for linear elasticity 2745754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 2755754ecacSJeremy L Thompson CEED_QFUNCTION(ElasLinearEnergy)(void *ctx, CeedInt Q, 2765754ecacSJeremy L Thompson const CeedScalar *const *in, 2775754ecacSJeremy L Thompson CeedScalar *const *out) { 2785754ecacSJeremy L Thompson // *INDENT-OFF* 2795754ecacSJeremy L Thompson // Inputs 2805754ecacSJeremy L Thompson const CeedScalar (*ug)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[0], 2815754ecacSJeremy L Thompson (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 2825754ecacSJeremy L Thompson 2835754ecacSJeremy L Thompson // Outputs 2845754ecacSJeremy L Thompson CeedScalar (*energy) = (CeedScalar(*))out[0]; 2855754ecacSJeremy L Thompson // *INDENT-ON* 2865754ecacSJeremy L Thompson 2875754ecacSJeremy L Thompson // Context 2885754ecacSJeremy L Thompson const Physics context = (Physics)ctx; 2895754ecacSJeremy L Thompson const CeedScalar E = context->E; 2905754ecacSJeremy L Thompson const CeedScalar nu = context->nu; 2915754ecacSJeremy L Thompson 2925754ecacSJeremy L Thompson // Constants 2935754ecacSJeremy L Thompson const CeedScalar TwoMu = E / (1 + nu); 2945754ecacSJeremy L Thompson const CeedScalar mu = TwoMu / 2; 2955754ecacSJeremy L Thompson const CeedScalar Kbulk = E / (3*(1 - 2*nu)); // Bulk Modulus 2965754ecacSJeremy L Thompson const CeedScalar lambda = (3*Kbulk - TwoMu) / 3; 2975754ecacSJeremy L Thompson 2985754ecacSJeremy L Thompson // Quadrature Point Loop 2995754ecacSJeremy L Thompson CeedPragmaSIMD 3005754ecacSJeremy L Thompson for (CeedInt i=0; i<Q; i++) { 3015754ecacSJeremy L Thompson // Read spatial derivatives of u 3025754ecacSJeremy L Thompson // *INDENT-OFF* 3035754ecacSJeremy L Thompson const CeedScalar du[3][3] = {{ug[0][0][i], 3045754ecacSJeremy L Thompson ug[1][0][i], 3055754ecacSJeremy L Thompson ug[2][0][i]}, 3065754ecacSJeremy L Thompson {ug[0][1][i], 3075754ecacSJeremy L Thompson ug[1][1][i], 3085754ecacSJeremy L Thompson ug[2][1][i]}, 3095754ecacSJeremy L Thompson {ug[0][2][i], 3105754ecacSJeremy L Thompson ug[1][2][i], 3115754ecacSJeremy L Thompson ug[2][2][i]} 3125754ecacSJeremy L Thompson }; 3135754ecacSJeremy L Thompson // -- Qdata 3145754ecacSJeremy L Thompson const CeedScalar wdetJ = q_data[0][i]; 3155754ecacSJeremy L Thompson const CeedScalar dXdx[3][3] = {{q_data[1][i], 3165754ecacSJeremy L Thompson q_data[2][i], 3175754ecacSJeremy L Thompson q_data[3][i]}, 3185754ecacSJeremy L Thompson {q_data[4][i], 3195754ecacSJeremy L Thompson q_data[5][i], 3205754ecacSJeremy L Thompson q_data[6][i]}, 3215754ecacSJeremy L Thompson {q_data[7][i], 3225754ecacSJeremy L Thompson q_data[8][i], 3235754ecacSJeremy L Thompson q_data[9][i]} 3245754ecacSJeremy L Thompson }; 3255754ecacSJeremy L Thompson // *INDENT-ON* 3265754ecacSJeremy L Thompson 3275754ecacSJeremy L Thompson // Compute grad_u 3285754ecacSJeremy L Thompson // dXdx = (dx/dX)^(-1) 3295754ecacSJeremy L Thompson // Apply dXdx to du = grad_u 3305754ecacSJeremy L Thompson CeedScalar grad_u[3][3]; 3315754ecacSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) // Component 3325754ecacSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { // Derivative 3335754ecacSJeremy L Thompson grad_u[j][k] = 0; 3345754ecacSJeremy L Thompson for (CeedInt m = 0; m < 3; m++) 3355754ecacSJeremy L Thompson grad_u[j][k] += dXdx[m][k] * du[j][m]; 3365754ecacSJeremy L Thompson } 3375754ecacSJeremy L Thompson 3385754ecacSJeremy L Thompson // Compute Strain : e (epsilon) 3395754ecacSJeremy L Thompson // e = 1/2 (grad u + (grad u)^T) 3405754ecacSJeremy L Thompson 3415754ecacSJeremy L Thompson // *INDENT-OFF* 3425754ecacSJeremy L Thompson const CeedScalar e[3][3] = {{(grad_u[0][0] + grad_u[0][0])/2., 3435754ecacSJeremy L Thompson (grad_u[0][1] + grad_u[1][0])/2., 3445754ecacSJeremy L Thompson (grad_u[0][2] + grad_u[2][0])/2.}, 3455754ecacSJeremy L Thompson {(grad_u[1][0] + grad_u[0][1])/2., 3465754ecacSJeremy L Thompson (grad_u[1][1] + grad_u[1][1])/2., 3475754ecacSJeremy L Thompson (grad_u[1][2] + grad_u[2][1])/2.}, 3485754ecacSJeremy L Thompson {(grad_u[2][0] + grad_u[0][2])/2., 3495754ecacSJeremy L Thompson (grad_u[2][1] + grad_u[1][2])/2., 3505754ecacSJeremy L Thompson (grad_u[2][2] + grad_u[2][2])/2.} 3515754ecacSJeremy L Thompson }; 3525754ecacSJeremy L Thompson // *INDENT-ON* 3535754ecacSJeremy L Thompson 3545754ecacSJeremy L Thompson // Strain energy 3555754ecacSJeremy L Thompson const CeedScalar strain_vol = e[0][0] + e[1][1] + e[2][2]; 3565754ecacSJeremy L Thompson energy[i] = (lambda*strain_vol*strain_vol/2. + strain_vol*mu + 3575754ecacSJeremy L Thompson (e[0][1]*e[0][1]+e[0][2]*e[0][2]+e[1][2]*e[1][2])*2*mu)*wdetJ; 3585754ecacSJeremy L Thompson 3595754ecacSJeremy L Thompson } // End of Quadrature Point Loop 3605754ecacSJeremy L Thompson 3615754ecacSJeremy L Thompson return 0; 3625754ecacSJeremy L Thompson } 3635754ecacSJeremy L Thompson 3645754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 3655754ecacSJeremy L Thompson // Nodal diagnostic quantities for linear elasticity 3665754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 3675754ecacSJeremy L Thompson CEED_QFUNCTION(ElasLinearDiagnostic)(void *ctx, CeedInt Q, 3685754ecacSJeremy L Thompson const CeedScalar *const *in, 3695754ecacSJeremy L Thompson CeedScalar *const *out) { 3705754ecacSJeremy L Thompson // *INDENT-OFF* 3715754ecacSJeremy L Thompson // Inputs 3725754ecacSJeremy L Thompson const CeedScalar (*u)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 3735754ecacSJeremy L Thompson (*ug)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[1], 3745754ecacSJeremy L Thompson (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 3755754ecacSJeremy L Thompson 3765754ecacSJeremy L Thompson // Outputs 3775754ecacSJeremy L Thompson CeedScalar (*diagnostic)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 3785754ecacSJeremy L Thompson // *INDENT-ON* 3795754ecacSJeremy L Thompson 3805754ecacSJeremy L Thompson // Context 3815754ecacSJeremy L Thompson const Physics context = (Physics)ctx; 3825754ecacSJeremy L Thompson const CeedScalar E = context->E; 3835754ecacSJeremy L Thompson const CeedScalar nu = context->nu; 3845754ecacSJeremy L Thompson 3855754ecacSJeremy L Thompson // Constants 3865754ecacSJeremy L Thompson const CeedScalar TwoMu = E / (1 + nu); 3875754ecacSJeremy L Thompson const CeedScalar mu = TwoMu / 2; 3885754ecacSJeremy L Thompson const CeedScalar Kbulk = E / (3*(1 - 2*nu)); // Bulk Modulus 3895754ecacSJeremy L Thompson const CeedScalar lambda = (3*Kbulk - TwoMu) / 3; 3905754ecacSJeremy L Thompson 3915754ecacSJeremy L Thompson // Quadrature Point Loop 3925754ecacSJeremy L Thompson CeedPragmaSIMD 3935754ecacSJeremy L Thompson for (CeedInt i=0; i<Q; i++) { 3945754ecacSJeremy L Thompson // Read spatial derivatives of u 3955754ecacSJeremy L Thompson // *INDENT-OFF* 3965754ecacSJeremy L Thompson const CeedScalar du[3][3] = {{ug[0][0][i], 3975754ecacSJeremy L Thompson ug[1][0][i], 3985754ecacSJeremy L Thompson ug[2][0][i]}, 3995754ecacSJeremy L Thompson {ug[0][1][i], 4005754ecacSJeremy L Thompson ug[1][1][i], 4015754ecacSJeremy L Thompson ug[2][1][i]}, 4025754ecacSJeremy L Thompson {ug[0][2][i], 4035754ecacSJeremy L Thompson ug[1][2][i], 4045754ecacSJeremy L Thompson ug[2][2][i]} 4055754ecacSJeremy L Thompson }; 4065754ecacSJeremy L Thompson // -- Qdata 4075754ecacSJeremy L Thompson const CeedScalar dXdx[3][3] = {{q_data[1][i], 4085754ecacSJeremy L Thompson q_data[2][i], 4095754ecacSJeremy L Thompson q_data[3][i]}, 4105754ecacSJeremy L Thompson {q_data[4][i], 4115754ecacSJeremy L Thompson q_data[5][i], 4125754ecacSJeremy L Thompson q_data[6][i]}, 4135754ecacSJeremy L Thompson {q_data[7][i], 4145754ecacSJeremy L Thompson q_data[8][i], 4155754ecacSJeremy L Thompson q_data[9][i]} 4165754ecacSJeremy L Thompson }; 4175754ecacSJeremy L Thompson // *INDENT-ON* 4185754ecacSJeremy L Thompson 4195754ecacSJeremy L Thompson // Compute grad_u 4205754ecacSJeremy L Thompson // dXdx = (dx/dX)^(-1) 4215754ecacSJeremy L Thompson // Apply dXdx to du = grad_u 4225754ecacSJeremy L Thompson CeedScalar grad_u[3][3]; 4235754ecacSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) // Component 4245754ecacSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { // Derivative 4255754ecacSJeremy L Thompson grad_u[j][k] = 0; 4265754ecacSJeremy L Thompson for (CeedInt m = 0; m < 3; m++) 4275754ecacSJeremy L Thompson grad_u[j][k] += dXdx[m][k] * du[j][m]; 4285754ecacSJeremy L Thompson } 4295754ecacSJeremy L Thompson 4305754ecacSJeremy L Thompson // Compute Strain : e (epsilon) 4315754ecacSJeremy L Thompson // e = 1/2 (grad u + (grad u)^T) 4325754ecacSJeremy L Thompson 4335754ecacSJeremy L Thompson // *INDENT-OFF* 4345754ecacSJeremy L Thompson const CeedScalar e[3][3] = {{(grad_u[0][0] + grad_u[0][0])/2., 4355754ecacSJeremy L Thompson (grad_u[0][1] + grad_u[1][0])/2., 4365754ecacSJeremy L Thompson (grad_u[0][2] + grad_u[2][0])/2.}, 4375754ecacSJeremy L Thompson {(grad_u[1][0] + grad_u[0][1])/2., 4385754ecacSJeremy L Thompson (grad_u[1][1] + grad_u[1][1])/2., 4395754ecacSJeremy L Thompson (grad_u[1][2] + grad_u[2][1])/2.}, 4405754ecacSJeremy L Thompson {(grad_u[2][0] + grad_u[0][2])/2., 4415754ecacSJeremy L Thompson (grad_u[2][1] + grad_u[1][2])/2., 4425754ecacSJeremy L Thompson (grad_u[2][2] + grad_u[2][2])/2.} 4435754ecacSJeremy L Thompson }; 4445754ecacSJeremy L Thompson // *INDENT-ON* 4455754ecacSJeremy L Thompson 4465754ecacSJeremy L Thompson // Displacement 4475754ecacSJeremy L Thompson diagnostic[0][i] = u[0][i]; 4485754ecacSJeremy L Thompson diagnostic[1][i] = u[1][i]; 4495754ecacSJeremy L Thompson diagnostic[2][i] = u[2][i]; 4505754ecacSJeremy L Thompson 4515754ecacSJeremy L Thompson // Pressure 4525754ecacSJeremy L Thompson const CeedScalar strain_vol = e[0][0] + e[1][1] + e[2][2]; 4535754ecacSJeremy L Thompson diagnostic[3][i] = -lambda*strain_vol; 4545754ecacSJeremy L Thompson 4555754ecacSJeremy L Thompson // Stress tensor invariants 4565754ecacSJeremy L Thompson diagnostic[4][i] = strain_vol; 4575754ecacSJeremy L Thompson diagnostic[5][i] = 0.; 4585754ecacSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) 4595754ecacSJeremy L Thompson for (CeedInt m = 0; m < 3; m++) 4605754ecacSJeremy L Thompson diagnostic[5][i] += e[j][m] * e[m][j]; 4615754ecacSJeremy L Thompson diagnostic[6][i] = 1 + strain_vol; 4625754ecacSJeremy L Thompson 4635754ecacSJeremy L Thompson // Strain energy 4645754ecacSJeremy L Thompson diagnostic[7][i] = (lambda*strain_vol*strain_vol/2. + strain_vol*mu + 4655754ecacSJeremy L Thompson (e[0][1]*e[0][1]+e[0][2]*e[0][2]+e[1][2]*e[1][2])*2*mu); 4665754ecacSJeremy L Thompson 4675754ecacSJeremy L Thompson } // End of Quadrature Point Loop 4685754ecacSJeremy L Thompson 4695754ecacSJeremy L Thompson return 0; 4705754ecacSJeremy L Thompson } 4715754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 4725754ecacSJeremy L Thompson 4735754ecacSJeremy L Thompson #endif //End of ELAS_LINEAR_H 474