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 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 // ----------------------------------------------------------------------------- 29*2b730f8bSJeremy L Thompson CEED_QFUNCTION(ElasLinearF)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 305754ecacSJeremy L Thompson // *INDENT-OFF* 315754ecacSJeremy L Thompson // Inputs 32*2b730f8bSJeremy 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]; 335754ecacSJeremy L Thompson 345754ecacSJeremy L Thompson // Outputs 355754ecacSJeremy L Thompson CeedScalar(*dvdX)[3][CEED_Q_VLA] = (CeedScalar(*)[3][CEED_Q_VLA])out[0]; 365754ecacSJeremy L Thompson // grad_u not used for linear elasticity 375754ecacSJeremy L Thompson // (*grad_u)[3][CEED_Q_VLA] = (CeedScalar(*)[3][CEED_Q_VLA])out[1]; 385754ecacSJeremy L Thompson // *INDENT-ON* 395754ecacSJeremy L Thompson 405754ecacSJeremy L Thompson // Context 415754ecacSJeremy L Thompson const Physics context = (Physics)ctx; 425754ecacSJeremy L Thompson const CeedScalar E = context->E; 435754ecacSJeremy L Thompson const CeedScalar nu = context->nu; 445754ecacSJeremy L Thompson 455754ecacSJeremy L Thompson // Quadrature Point Loop 46*2b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 475754ecacSJeremy L Thompson // Read spatial derivatives of u 485754ecacSJeremy L Thompson // *INDENT-OFF* 49*2b730f8bSJeremy L Thompson const CeedScalar du[3][3] = { 50*2b730f8bSJeremy L Thompson {ug[0][0][i], ug[1][0][i], ug[2][0][i]}, 51*2b730f8bSJeremy L Thompson {ug[0][1][i], ug[1][1][i], ug[2][1][i]}, 52*2b730f8bSJeremy L Thompson {ug[0][2][i], ug[1][2][i], ug[2][2][i]} 535754ecacSJeremy L Thompson }; 545754ecacSJeremy L Thompson // -- Qdata 555754ecacSJeremy L Thompson const CeedScalar wdetJ = q_data[0][i]; 56*2b730f8bSJeremy L Thompson const CeedScalar dXdx[3][3] = { 57*2b730f8bSJeremy L Thompson {q_data[1][i], q_data[2][i], q_data[3][i]}, 58*2b730f8bSJeremy L Thompson {q_data[4][i], q_data[5][i], q_data[6][i]}, 59*2b730f8bSJeremy L Thompson {q_data[7][i], q_data[8][i], q_data[9][i]} 605754ecacSJeremy L Thompson }; 615754ecacSJeremy L Thompson // *INDENT-ON* 625754ecacSJeremy L Thompson 635754ecacSJeremy L Thompson // Compute grad_u 645754ecacSJeremy L Thompson // dXdx = (dx/dX)^(-1) 655754ecacSJeremy L Thompson // Apply dXdx to du = grad_u 665754ecacSJeremy L Thompson CeedScalar grad_u[3][3]; 67*2b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { // Component 685754ecacSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { // Derivative 695754ecacSJeremy L Thompson grad_u[j][k] = 0; 70*2b730f8bSJeremy L Thompson for (CeedInt m = 0; m < 3; m++) grad_u[j][k] += dXdx[m][k] * du[j][m]; 71*2b730f8bSJeremy L Thompson } 725754ecacSJeremy L Thompson } 735754ecacSJeremy L Thompson 745754ecacSJeremy L Thompson // Compute Strain : e (epsilon) 755754ecacSJeremy L Thompson // e = 1/2 (grad u + (grad u)^T) 765754ecacSJeremy L Thompson 775754ecacSJeremy L Thompson // *INDENT-OFF* 78*2b730f8bSJeremy L Thompson const CeedScalar e[3][3] = { 79*2b730f8bSJeremy 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.}, 80*2b730f8bSJeremy 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.}, 81*2b730f8bSJeremy 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.} 825754ecacSJeremy L Thompson }; 835754ecacSJeremy L Thompson // *INDENT-ON* 845754ecacSJeremy L Thompson 855754ecacSJeremy L Thompson // 865754ecacSJeremy L Thompson // Formulation uses Voigt notation: 875754ecacSJeremy L Thompson // stress (sigma) strain (epsilon) 885754ecacSJeremy L Thompson // 895754ecacSJeremy L Thompson // [sigma00] [e00] 905754ecacSJeremy L Thompson // [sigma11] [e11] 915754ecacSJeremy L Thompson // [sigma22] = S * [e22] 925754ecacSJeremy L Thompson // [sigma12] [e12] 935754ecacSJeremy L Thompson // [sigma02] [e02] 945754ecacSJeremy L Thompson // [sigma01] [e01] 955754ecacSJeremy L Thompson // 965754ecacSJeremy L Thompson // where 975754ecacSJeremy L Thompson // [1-nu nu nu ] 985754ecacSJeremy L Thompson // [ nu 1-nu nu ] 995754ecacSJeremy L Thompson // [ nu nu 1-nu ] 1005754ecacSJeremy L Thompson // S = E/((1+nu)*(1-2*nu)) [ (1-2*nu)/2 ] 1015754ecacSJeremy L Thompson // [ (1-2*nu)/2 ] 1025754ecacSJeremy L Thompson // [ (1-2*nu)/2 ] 1035754ecacSJeremy L Thompson 1045754ecacSJeremy L Thompson // Above Voigt Notation is placed in a 3x3 matrix: 1055754ecacSJeremy L Thompson const CeedScalar ss = E / ((1 + nu) * (1 - 2 * nu)); 1065754ecacSJeremy L Thompson const CeedScalar sigma00 = ss * ((1 - nu) * e[0][0] + nu * e[1][1] + nu * e[2][2]), 1075754ecacSJeremy L Thompson sigma11 = ss * (nu * e[0][0] + (1 - nu) * e[1][1] + nu * e[2][2]), 108*2b730f8bSJeremy 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, 109*2b730f8bSJeremy L Thompson sigma02 = ss * (1 - 2 * nu) * e[0][2] * 0.5, sigma01 = ss * (1 - 2 * nu) * e[0][1] * 0.5; 1105754ecacSJeremy L Thompson // *INDENT-OFF* 111*2b730f8bSJeremy L Thompson const CeedScalar sigma[3][3] = { 112*2b730f8bSJeremy L Thompson {sigma00, sigma01, sigma02}, 1135754ecacSJeremy L Thompson {sigma01, sigma11, sigma12}, 1145754ecacSJeremy L Thompson {sigma02, sigma12, sigma22} 1155754ecacSJeremy L Thompson }; 1165754ecacSJeremy L Thompson // *INDENT-ON* 1175754ecacSJeremy L Thompson 1185754ecacSJeremy L Thompson // Apply dXdx^T and weight to sigma 119*2b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { // Component 1205754ecacSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { // Derivative 1215754ecacSJeremy L Thompson dvdX[k][j][i] = 0; 122*2b730f8bSJeremy L Thompson for (CeedInt m = 0; m < 3; m++) dvdX[k][j][i] += dXdx[k][m] * sigma[j][m] * wdetJ; 1235754ecacSJeremy L Thompson } 124*2b730f8bSJeremy L Thompson } 1255754ecacSJeremy L Thompson } // End of Quadrature Point Loop 1265754ecacSJeremy L Thompson 1275754ecacSJeremy L Thompson return 0; 1285754ecacSJeremy L Thompson } 1295754ecacSJeremy L Thompson 1305754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 1315754ecacSJeremy L Thompson // Jacobian evaluation for linear elasticity 1325754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 133*2b730f8bSJeremy L Thompson CEED_QFUNCTION(ElasLineardF)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 1345754ecacSJeremy L Thompson // *INDENT-OFF* 1355754ecacSJeremy L Thompson // Inputs 1365754ecacSJeremy L Thompson const CeedScalar(*deltaug)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[0], 1375754ecacSJeremy L Thompson (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 1385754ecacSJeremy L Thompson // grad_u not used for linear elasticity 1395754ecacSJeremy L Thompson // (*grad_u)[3][Q] = (CeedScalar(*)[3][Q])in[2]; 1405754ecacSJeremy L Thompson 1415754ecacSJeremy L Thompson // Outputs 1425754ecacSJeremy L Thompson CeedScalar(*deltadvdX)[3][CEED_Q_VLA] = (CeedScalar(*)[3][CEED_Q_VLA])out[0]; 1435754ecacSJeremy L Thompson // *INDENT-ON* 1445754ecacSJeremy L Thompson 1455754ecacSJeremy L Thompson // Context 1465754ecacSJeremy L Thompson const Physics context = (Physics)ctx; 1475754ecacSJeremy L Thompson const CeedScalar E = context->E; 1485754ecacSJeremy L Thompson const CeedScalar nu = context->nu; 1495754ecacSJeremy L Thompson 1505754ecacSJeremy L Thompson // Quadrature Point Loop 151*2b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 1525754ecacSJeremy L Thompson // Read spatial derivatives of u 1535754ecacSJeremy L Thompson // *INDENT-OFF* 154*2b730f8bSJeremy L Thompson const CeedScalar deltadu[3][3] = { 155*2b730f8bSJeremy L Thompson {deltaug[0][0][i], deltaug[1][0][i], deltaug[2][0][i]}, 156*2b730f8bSJeremy L Thompson {deltaug[0][1][i], deltaug[1][1][i], deltaug[2][1][i]}, 157*2b730f8bSJeremy L Thompson {deltaug[0][2][i], deltaug[1][2][i], deltaug[2][2][i]} 1585754ecacSJeremy L Thompson }; 1595754ecacSJeremy L Thompson // -- Qdata 1605754ecacSJeremy L Thompson const CeedScalar wdetJ = q_data[0][i]; 161*2b730f8bSJeremy L Thompson const CeedScalar dXdx[3][3] = { 162*2b730f8bSJeremy L Thompson {q_data[1][i], q_data[2][i], q_data[3][i]}, 163*2b730f8bSJeremy L Thompson {q_data[4][i], q_data[5][i], q_data[6][i]}, 164*2b730f8bSJeremy L Thompson {q_data[7][i], q_data[8][i], q_data[9][i]} 1655754ecacSJeremy L Thompson }; 1665754ecacSJeremy L Thompson // *INDENT-ON* 1675754ecacSJeremy L Thompson 1685754ecacSJeremy L Thompson // Compute graddeltau 1695754ecacSJeremy L Thompson // dXdx = (dx/dX)^(-1) 1705754ecacSJeremy L Thompson // Apply dXdx to deltadu = graddeltau 1715754ecacSJeremy L Thompson CeedScalar graddeltau[3][3]; 172*2b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { // Component 1735754ecacSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { // Derivative 1745754ecacSJeremy L Thompson graddeltau[j][k] = 0; 175*2b730f8bSJeremy L Thompson for (CeedInt m = 0; m < 3; m++) graddeltau[j][k] += dXdx[m][k] * deltadu[j][m]; 176*2b730f8bSJeremy L Thompson } 1775754ecacSJeremy L Thompson } 1785754ecacSJeremy L Thompson 1795754ecacSJeremy L Thompson // Compute Strain : e (epsilon) 1805754ecacSJeremy L Thompson // e = 1/2 (grad u + (grad u)^T) 1815754ecacSJeremy L Thompson // *INDENT-OFF* 182*2b730f8bSJeremy L Thompson const CeedScalar de[3][3] = { 183*2b730f8bSJeremy L Thompson {(graddeltau[0][0] + graddeltau[0][0]) / 2., (graddeltau[0][1] + graddeltau[1][0]) / 2., (graddeltau[0][2] + graddeltau[2][0]) / 2.}, 184*2b730f8bSJeremy L Thompson {(graddeltau[1][0] + graddeltau[0][1]) / 2., (graddeltau[1][1] + graddeltau[1][1]) / 2., (graddeltau[1][2] + graddeltau[2][1]) / 2.}, 185*2b730f8bSJeremy L Thompson {(graddeltau[2][0] + graddeltau[0][2]) / 2., (graddeltau[2][1] + graddeltau[1][2]) / 2., (graddeltau[2][2] + graddeltau[2][2]) / 2.} 1865754ecacSJeremy L Thompson }; 1875754ecacSJeremy L Thompson 1885754ecacSJeremy L Thompson // *INDENT-ON* 1895754ecacSJeremy L Thompson // Formulation uses Voigt notation: 1905754ecacSJeremy L Thompson // stress (sigma) strain (epsilon) 1915754ecacSJeremy L Thompson // 1925754ecacSJeremy L Thompson // [dsigma00] [de00] 1935754ecacSJeremy L Thompson // [dsigma11] [de11] 1945754ecacSJeremy L Thompson // [dsigma22] = S * [de22] 1955754ecacSJeremy L Thompson // [dsigma12] [de12] 1965754ecacSJeremy L Thompson // [dsigma02] [de02] 1975754ecacSJeremy L Thompson // [dsigma01] [de01] 1985754ecacSJeremy L Thompson // 1995754ecacSJeremy L Thompson // where 2005754ecacSJeremy L Thompson // [1-nu nu nu ] 2015754ecacSJeremy L Thompson // [ nu 1-nu nu ] 2025754ecacSJeremy L Thompson // [ nu nu 1-nu ] 2035754ecacSJeremy L Thompson // S = E/((1+nu)*(1-2*nu)) [ (1-2*nu)/2 ] 2045754ecacSJeremy L Thompson // [ (1-2*nu)/2 ] 2055754ecacSJeremy L Thompson // [ (1-2*nu)/2 ] 2065754ecacSJeremy L Thompson 2075754ecacSJeremy L Thompson // Above Voigt Notation is placed in a 3x3 matrix: 2085754ecacSJeremy L Thompson const CeedScalar ss = E / ((1 + nu) * (1 - 2 * nu)); 2095754ecacSJeremy L Thompson const CeedScalar dsigma00 = ss * ((1 - nu) * de[0][0] + nu * de[1][1] + nu * de[2][2]), 2105754ecacSJeremy L Thompson dsigma11 = ss * (nu * de[0][0] + (1 - nu) * de[1][1] + nu * de[2][2]), 211*2b730f8bSJeremy 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, 212*2b730f8bSJeremy L Thompson dsigma02 = ss * (1 - 2 * nu) * de[0][2] / 2, dsigma01 = ss * (1 - 2 * nu) * de[0][1] / 2; 2135754ecacSJeremy L Thompson // *INDENT-OFF* 214*2b730f8bSJeremy L Thompson const CeedScalar dsigma[3][3] = { 215*2b730f8bSJeremy L Thompson {dsigma00, dsigma01, dsigma02}, 2165754ecacSJeremy L Thompson {dsigma01, dsigma11, dsigma12}, 2175754ecacSJeremy L Thompson {dsigma02, dsigma12, dsigma22} 2185754ecacSJeremy L Thompson }; 2195754ecacSJeremy L Thompson // *INDENT-ON* 2205754ecacSJeremy L Thompson 2215754ecacSJeremy L Thompson // Apply dXdx^T and weight 222*2b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { // Component 2235754ecacSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { // Derivative 2245754ecacSJeremy L Thompson deltadvdX[k][j][i] = 0; 225*2b730f8bSJeremy L Thompson for (CeedInt m = 0; m < 3; m++) deltadvdX[k][j][i] += dXdx[k][m] * dsigma[j][m] * wdetJ; 2265754ecacSJeremy L Thompson } 227*2b730f8bSJeremy L Thompson } 2285754ecacSJeremy L Thompson } // End of Quadrature Point Loop 2295754ecacSJeremy L Thompson 2305754ecacSJeremy L Thompson return 0; 2315754ecacSJeremy L Thompson } 2325754ecacSJeremy L Thompson 2335754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 2345754ecacSJeremy L Thompson // Strain energy computation for linear elasticity 2355754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 236*2b730f8bSJeremy L Thompson CEED_QFUNCTION(ElasLinearEnergy)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 2375754ecacSJeremy L Thompson // *INDENT-OFF* 2385754ecacSJeremy L Thompson // Inputs 239*2b730f8bSJeremy 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]; 2405754ecacSJeremy L Thompson 2415754ecacSJeremy L Thompson // Outputs 2425754ecacSJeremy L Thompson CeedScalar(*energy) = (CeedScalar(*))out[0]; 2435754ecacSJeremy L Thompson // *INDENT-ON* 2445754ecacSJeremy L Thompson 2455754ecacSJeremy L Thompson // Context 2465754ecacSJeremy L Thompson const Physics context = (Physics)ctx; 2475754ecacSJeremy L Thompson const CeedScalar E = context->E; 2485754ecacSJeremy L Thompson const CeedScalar nu = context->nu; 2495754ecacSJeremy L Thompson 2505754ecacSJeremy L Thompson // Constants 2515754ecacSJeremy L Thompson const CeedScalar TwoMu = E / (1 + nu); 2525754ecacSJeremy L Thompson const CeedScalar mu = TwoMu / 2; 2535754ecacSJeremy L Thompson const CeedScalar Kbulk = E / (3 * (1 - 2 * nu)); // Bulk Modulus 2545754ecacSJeremy L Thompson const CeedScalar lambda = (3 * Kbulk - TwoMu) / 3; 2555754ecacSJeremy L Thompson 2565754ecacSJeremy L Thompson // Quadrature Point Loop 257*2b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 2585754ecacSJeremy L Thompson // Read spatial derivatives of u 2595754ecacSJeremy L Thompson // *INDENT-OFF* 260*2b730f8bSJeremy L Thompson const CeedScalar du[3][3] = { 261*2b730f8bSJeremy L Thompson {ug[0][0][i], ug[1][0][i], ug[2][0][i]}, 262*2b730f8bSJeremy L Thompson {ug[0][1][i], ug[1][1][i], ug[2][1][i]}, 263*2b730f8bSJeremy L Thompson {ug[0][2][i], ug[1][2][i], ug[2][2][i]} 2645754ecacSJeremy L Thompson }; 2655754ecacSJeremy L Thompson // -- Qdata 2665754ecacSJeremy L Thompson const CeedScalar wdetJ = q_data[0][i]; 267*2b730f8bSJeremy L Thompson const CeedScalar dXdx[3][3] = { 268*2b730f8bSJeremy L Thompson {q_data[1][i], q_data[2][i], q_data[3][i]}, 269*2b730f8bSJeremy L Thompson {q_data[4][i], q_data[5][i], q_data[6][i]}, 270*2b730f8bSJeremy L Thompson {q_data[7][i], q_data[8][i], q_data[9][i]} 2715754ecacSJeremy L Thompson }; 2725754ecacSJeremy L Thompson // *INDENT-ON* 2735754ecacSJeremy L Thompson 2745754ecacSJeremy L Thompson // Compute grad_u 2755754ecacSJeremy L Thompson // dXdx = (dx/dX)^(-1) 2765754ecacSJeremy L Thompson // Apply dXdx to du = grad_u 2775754ecacSJeremy L Thompson CeedScalar grad_u[3][3]; 278*2b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { // Component 2795754ecacSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { // Derivative 2805754ecacSJeremy L Thompson grad_u[j][k] = 0; 281*2b730f8bSJeremy L Thompson for (CeedInt m = 0; m < 3; m++) grad_u[j][k] += dXdx[m][k] * du[j][m]; 282*2b730f8bSJeremy L Thompson } 2835754ecacSJeremy L Thompson } 2845754ecacSJeremy L Thompson 2855754ecacSJeremy L Thompson // Compute Strain : e (epsilon) 2865754ecacSJeremy L Thompson // e = 1/2 (grad u + (grad u)^T) 2875754ecacSJeremy L Thompson 2885754ecacSJeremy L Thompson // *INDENT-OFF* 289*2b730f8bSJeremy L Thompson const CeedScalar e[3][3] = { 290*2b730f8bSJeremy 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.}, 291*2b730f8bSJeremy 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.}, 292*2b730f8bSJeremy 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.} 2935754ecacSJeremy L Thompson }; 2945754ecacSJeremy L Thompson // *INDENT-ON* 2955754ecacSJeremy L Thompson 2965754ecacSJeremy L Thompson // Strain energy 2975754ecacSJeremy L Thompson const CeedScalar strain_vol = e[0][0] + e[1][1] + e[2][2]; 298*2b730f8bSJeremy L Thompson energy[i] = 299*2b730f8bSJeremy 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; 3005754ecacSJeremy L Thompson 3015754ecacSJeremy L Thompson } // End of Quadrature Point Loop 3025754ecacSJeremy L Thompson 3035754ecacSJeremy L Thompson return 0; 3045754ecacSJeremy L Thompson } 3055754ecacSJeremy L Thompson 3065754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 3075754ecacSJeremy L Thompson // Nodal diagnostic quantities for linear elasticity 3085754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 309*2b730f8bSJeremy L Thompson CEED_QFUNCTION(ElasLinearDiagnostic)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 3105754ecacSJeremy L Thompson // *INDENT-OFF* 3115754ecacSJeremy L Thompson // Inputs 312*2b730f8bSJeremy 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], 3135754ecacSJeremy L Thompson (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 3145754ecacSJeremy L Thompson 3155754ecacSJeremy L Thompson // Outputs 3165754ecacSJeremy L Thompson CeedScalar(*diagnostic)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 3175754ecacSJeremy L Thompson // *INDENT-ON* 3185754ecacSJeremy L Thompson 3195754ecacSJeremy L Thompson // Context 3205754ecacSJeremy L Thompson const Physics context = (Physics)ctx; 3215754ecacSJeremy L Thompson const CeedScalar E = context->E; 3225754ecacSJeremy L Thompson const CeedScalar nu = context->nu; 3235754ecacSJeremy L Thompson 3245754ecacSJeremy L Thompson // Constants 3255754ecacSJeremy L Thompson const CeedScalar TwoMu = E / (1 + nu); 3265754ecacSJeremy L Thompson const CeedScalar mu = TwoMu / 2; 3275754ecacSJeremy L Thompson const CeedScalar Kbulk = E / (3 * (1 - 2 * nu)); // Bulk Modulus 3285754ecacSJeremy L Thompson const CeedScalar lambda = (3 * Kbulk - TwoMu) / 3; 3295754ecacSJeremy L Thompson 3305754ecacSJeremy L Thompson // Quadrature Point Loop 331*2b730f8bSJeremy L Thompson CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 3325754ecacSJeremy L Thompson // Read spatial derivatives of u 3335754ecacSJeremy L Thompson // *INDENT-OFF* 334*2b730f8bSJeremy L Thompson const CeedScalar du[3][3] = { 335*2b730f8bSJeremy L Thompson {ug[0][0][i], ug[1][0][i], ug[2][0][i]}, 336*2b730f8bSJeremy L Thompson {ug[0][1][i], ug[1][1][i], ug[2][1][i]}, 337*2b730f8bSJeremy L Thompson {ug[0][2][i], ug[1][2][i], ug[2][2][i]} 3385754ecacSJeremy L Thompson }; 3395754ecacSJeremy L Thompson // -- Qdata 340*2b730f8bSJeremy L Thompson const CeedScalar dXdx[3][3] = { 341*2b730f8bSJeremy L Thompson {q_data[1][i], q_data[2][i], q_data[3][i]}, 342*2b730f8bSJeremy L Thompson {q_data[4][i], q_data[5][i], q_data[6][i]}, 343*2b730f8bSJeremy L Thompson {q_data[7][i], q_data[8][i], q_data[9][i]} 3445754ecacSJeremy L Thompson }; 3455754ecacSJeremy L Thompson // *INDENT-ON* 3465754ecacSJeremy L Thompson 3475754ecacSJeremy L Thompson // Compute grad_u 3485754ecacSJeremy L Thompson // dXdx = (dx/dX)^(-1) 3495754ecacSJeremy L Thompson // Apply dXdx to du = grad_u 3505754ecacSJeremy L Thompson CeedScalar grad_u[3][3]; 351*2b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { // Component 3525754ecacSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { // Derivative 3535754ecacSJeremy L Thompson grad_u[j][k] = 0; 354*2b730f8bSJeremy L Thompson for (CeedInt m = 0; m < 3; m++) grad_u[j][k] += dXdx[m][k] * du[j][m]; 355*2b730f8bSJeremy L Thompson } 3565754ecacSJeremy L Thompson } 3575754ecacSJeremy L Thompson 3585754ecacSJeremy L Thompson // Compute Strain : e (epsilon) 3595754ecacSJeremy L Thompson // e = 1/2 (grad u + (grad u)^T) 3605754ecacSJeremy L Thompson 3615754ecacSJeremy L Thompson // *INDENT-OFF* 362*2b730f8bSJeremy L Thompson const CeedScalar e[3][3] = { 363*2b730f8bSJeremy 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.}, 364*2b730f8bSJeremy 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.}, 365*2b730f8bSJeremy 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.} 3665754ecacSJeremy L Thompson }; 3675754ecacSJeremy L Thompson // *INDENT-ON* 3685754ecacSJeremy L Thompson 3695754ecacSJeremy L Thompson // Displacement 3705754ecacSJeremy L Thompson diagnostic[0][i] = u[0][i]; 3715754ecacSJeremy L Thompson diagnostic[1][i] = u[1][i]; 3725754ecacSJeremy L Thompson diagnostic[2][i] = u[2][i]; 3735754ecacSJeremy L Thompson 3745754ecacSJeremy L Thompson // Pressure 3755754ecacSJeremy L Thompson const CeedScalar strain_vol = e[0][0] + e[1][1] + e[2][2]; 3765754ecacSJeremy L Thompson diagnostic[3][i] = -lambda * strain_vol; 3775754ecacSJeremy L Thompson 3785754ecacSJeremy L Thompson // Stress tensor invariants 3795754ecacSJeremy L Thompson diagnostic[4][i] = strain_vol; 3805754ecacSJeremy L Thompson diagnostic[5][i] = 0.; 381*2b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) { 382*2b730f8bSJeremy L Thompson for (CeedInt m = 0; m < 3; m++) diagnostic[5][i] += e[j][m] * e[m][j]; 383*2b730f8bSJeremy L Thompson } 3845754ecacSJeremy L Thompson diagnostic[6][i] = 1 + strain_vol; 3855754ecacSJeremy L Thompson 3865754ecacSJeremy L Thompson // Strain energy 387*2b730f8bSJeremy L Thompson diagnostic[7][i] = 388*2b730f8bSJeremy 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); 3895754ecacSJeremy L Thompson } // End of Quadrature Point Loop 3905754ecacSJeremy L Thompson 3915754ecacSJeremy L Thompson return 0; 3925754ecacSJeremy L Thompson } 3935754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 3945754ecacSJeremy L Thompson 3955754ecacSJeremy L Thompson #endif // End of ELAS_LINEAR_H 396