1*5754ecacSJeremy L Thompson // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 2*5754ecacSJeremy L Thompson // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 3*5754ecacSJeremy L Thompson // reserved. See files LICENSE and NOTICE for details. 4*5754ecacSJeremy L Thompson // 5*5754ecacSJeremy L Thompson // This file is part of CEED, a collection of benchmarks, miniapps, software 6*5754ecacSJeremy L Thompson // libraries and APIs for efficient high-order finite element and spectral 7*5754ecacSJeremy L Thompson // element discretizations for exascale applications. For more information and 8*5754ecacSJeremy L Thompson // source code availability see http://github.com/ceed. 9*5754ecacSJeremy L Thompson // 10*5754ecacSJeremy L Thompson // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 11*5754ecacSJeremy L Thompson // a collaborative effort of two U.S. Department of Energy organizations (Office 12*5754ecacSJeremy L Thompson // of Science and the National Nuclear Security Administration) responsible for 13*5754ecacSJeremy L Thompson // the planning and preparation of a capable exascale ecosystem, including 14*5754ecacSJeremy L Thompson // software, applications, hardware, advanced system engineering and early 15*5754ecacSJeremy L Thompson // testbed platforms, in support of the nation's exascale computing imperative. 16*5754ecacSJeremy L Thompson 17*5754ecacSJeremy L Thompson /// @file 18*5754ecacSJeremy L Thompson /// Linear elasticity for solid mechanics example using PETSc 19*5754ecacSJeremy L Thompson 20*5754ecacSJeremy L Thompson #ifndef ELAS_LINEAR_H 21*5754ecacSJeremy L Thompson #define ELAS_LINEAR_H 22*5754ecacSJeremy L Thompson 23*5754ecacSJeremy L Thompson #ifndef __CUDACC__ 24*5754ecacSJeremy L Thompson # include <math.h> 25*5754ecacSJeremy L Thompson #endif 26*5754ecacSJeremy L Thompson 27*5754ecacSJeremy L Thompson #ifndef PHYSICS_STRUCT 28*5754ecacSJeremy L Thompson #define PHYSICS_STRUCT 29*5754ecacSJeremy L Thompson typedef struct Physics_private *Physics; 30*5754ecacSJeremy L Thompson struct Physics_private { 31*5754ecacSJeremy L Thompson CeedScalar nu; // Poisson's ratio 32*5754ecacSJeremy L Thompson CeedScalar E; // Young's Modulus 33*5754ecacSJeremy L Thompson }; 34*5754ecacSJeremy L Thompson #endif 35*5754ecacSJeremy L Thompson 36*5754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 37*5754ecacSJeremy L Thompson // Residual evaluation for linear elasticity 38*5754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 39*5754ecacSJeremy L Thompson CEED_QFUNCTION(ElasLinearF)(void *ctx, CeedInt Q, const CeedScalar *const *in, 40*5754ecacSJeremy L Thompson CeedScalar *const *out) { 41*5754ecacSJeremy L Thompson // *INDENT-OFF* 42*5754ecacSJeremy L Thompson // Inputs 43*5754ecacSJeremy L Thompson const CeedScalar (*ug)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[0], 44*5754ecacSJeremy L Thompson (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 45*5754ecacSJeremy L Thompson 46*5754ecacSJeremy L Thompson // Outputs 47*5754ecacSJeremy L Thompson CeedScalar (*dvdX)[3][CEED_Q_VLA] = (CeedScalar(*)[3][CEED_Q_VLA])out[0]; 48*5754ecacSJeremy L Thompson // grad_u not used for linear elasticity 49*5754ecacSJeremy L Thompson // (*grad_u)[3][CEED_Q_VLA] = (CeedScalar(*)[3][CEED_Q_VLA])out[1]; 50*5754ecacSJeremy L Thompson // *INDENT-ON* 51*5754ecacSJeremy L Thompson 52*5754ecacSJeremy L Thompson // Context 53*5754ecacSJeremy L Thompson const Physics context = (Physics)ctx; 54*5754ecacSJeremy L Thompson const CeedScalar E = context->E; 55*5754ecacSJeremy L Thompson const CeedScalar nu = context->nu; 56*5754ecacSJeremy L Thompson 57*5754ecacSJeremy L Thompson // Quadrature Point Loop 58*5754ecacSJeremy L Thompson CeedPragmaSIMD 59*5754ecacSJeremy L Thompson for (CeedInt i=0; i<Q; i++) { 60*5754ecacSJeremy L Thompson // Read spatial derivatives of u 61*5754ecacSJeremy L Thompson // *INDENT-OFF* 62*5754ecacSJeremy L Thompson const CeedScalar du[3][3] = {{ug[0][0][i], 63*5754ecacSJeremy L Thompson ug[1][0][i], 64*5754ecacSJeremy L Thompson ug[2][0][i]}, 65*5754ecacSJeremy L Thompson {ug[0][1][i], 66*5754ecacSJeremy L Thompson ug[1][1][i], 67*5754ecacSJeremy L Thompson ug[2][1][i]}, 68*5754ecacSJeremy L Thompson {ug[0][2][i], 69*5754ecacSJeremy L Thompson ug[1][2][i], 70*5754ecacSJeremy L Thompson ug[2][2][i]} 71*5754ecacSJeremy L Thompson }; 72*5754ecacSJeremy L Thompson // -- Qdata 73*5754ecacSJeremy L Thompson const CeedScalar wdetJ = q_data[0][i]; 74*5754ecacSJeremy L Thompson const CeedScalar dXdx[3][3] = {{q_data[1][i], 75*5754ecacSJeremy L Thompson q_data[2][i], 76*5754ecacSJeremy L Thompson q_data[3][i]}, 77*5754ecacSJeremy L Thompson {q_data[4][i], 78*5754ecacSJeremy L Thompson q_data[5][i], 79*5754ecacSJeremy L Thompson q_data[6][i]}, 80*5754ecacSJeremy L Thompson {q_data[7][i], 81*5754ecacSJeremy L Thompson q_data[8][i], 82*5754ecacSJeremy L Thompson q_data[9][i]} 83*5754ecacSJeremy L Thompson }; 84*5754ecacSJeremy L Thompson // *INDENT-ON* 85*5754ecacSJeremy L Thompson 86*5754ecacSJeremy L Thompson // Compute grad_u 87*5754ecacSJeremy L Thompson // dXdx = (dx/dX)^(-1) 88*5754ecacSJeremy L Thompson // Apply dXdx to du = grad_u 89*5754ecacSJeremy L Thompson CeedScalar grad_u[3][3]; 90*5754ecacSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) // Component 91*5754ecacSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { // Derivative 92*5754ecacSJeremy L Thompson grad_u[j][k] = 0; 93*5754ecacSJeremy L Thompson for (CeedInt m = 0; m < 3; m++) 94*5754ecacSJeremy L Thompson grad_u[j][k] += dXdx[m][k] * du[j][m]; 95*5754ecacSJeremy L Thompson } 96*5754ecacSJeremy L Thompson 97*5754ecacSJeremy L Thompson // Compute Strain : e (epsilon) 98*5754ecacSJeremy L Thompson // e = 1/2 (grad u + (grad u)^T) 99*5754ecacSJeremy L Thompson 100*5754ecacSJeremy L Thompson // *INDENT-OFF* 101*5754ecacSJeremy L Thompson const CeedScalar e[3][3] = {{(grad_u[0][0] + grad_u[0][0])/2., 102*5754ecacSJeremy L Thompson (grad_u[0][1] + grad_u[1][0])/2., 103*5754ecacSJeremy L Thompson (grad_u[0][2] + grad_u[2][0])/2.}, 104*5754ecacSJeremy L Thompson {(grad_u[1][0] + grad_u[0][1])/2., 105*5754ecacSJeremy L Thompson (grad_u[1][1] + grad_u[1][1])/2., 106*5754ecacSJeremy L Thompson (grad_u[1][2] + grad_u[2][1])/2.}, 107*5754ecacSJeremy L Thompson {(grad_u[2][0] + grad_u[0][2])/2., 108*5754ecacSJeremy L Thompson (grad_u[2][1] + grad_u[1][2])/2., 109*5754ecacSJeremy L Thompson (grad_u[2][2] + grad_u[2][2])/2.} 110*5754ecacSJeremy L Thompson }; 111*5754ecacSJeremy L Thompson // *INDENT-ON* 112*5754ecacSJeremy L Thompson 113*5754ecacSJeremy L Thompson // 114*5754ecacSJeremy L Thompson // Formulation uses Voigt notation: 115*5754ecacSJeremy L Thompson // stress (sigma) strain (epsilon) 116*5754ecacSJeremy L Thompson // 117*5754ecacSJeremy L Thompson // [sigma00] [e00] 118*5754ecacSJeremy L Thompson // [sigma11] [e11] 119*5754ecacSJeremy L Thompson // [sigma22] = S * [e22] 120*5754ecacSJeremy L Thompson // [sigma12] [e12] 121*5754ecacSJeremy L Thompson // [sigma02] [e02] 122*5754ecacSJeremy L Thompson // [sigma01] [e01] 123*5754ecacSJeremy L Thompson // 124*5754ecacSJeremy L Thompson // where 125*5754ecacSJeremy L Thompson // [1-nu nu nu ] 126*5754ecacSJeremy L Thompson // [ nu 1-nu nu ] 127*5754ecacSJeremy L Thompson // [ nu nu 1-nu ] 128*5754ecacSJeremy L Thompson // S = E/((1+nu)*(1-2*nu)) [ (1-2*nu)/2 ] 129*5754ecacSJeremy L Thompson // [ (1-2*nu)/2 ] 130*5754ecacSJeremy L Thompson // [ (1-2*nu)/2 ] 131*5754ecacSJeremy L Thompson 132*5754ecacSJeremy L Thompson // Above Voigt Notation is placed in a 3x3 matrix: 133*5754ecacSJeremy L Thompson const CeedScalar ss = E / ((1 + nu)*(1 - 2*nu)); 134*5754ecacSJeremy L Thompson const CeedScalar sigma00 = ss*((1 - nu)*e[0][0] + nu*e[1][1] + nu*e[2][2]), 135*5754ecacSJeremy L Thompson sigma11 = ss*(nu*e[0][0] + (1 - nu)*e[1][1] + nu*e[2][2]), 136*5754ecacSJeremy L Thompson sigma22 = ss*(nu*e[0][0] + nu*e[1][1] + (1 - nu)*e[2][2]), 137*5754ecacSJeremy L Thompson sigma12 = ss*(1 - 2*nu)*e[1][2]*0.5, 138*5754ecacSJeremy L Thompson sigma02 = ss*(1 - 2*nu)*e[0][2]*0.5, 139*5754ecacSJeremy L Thompson sigma01 = ss*(1 - 2*nu)*e[0][1]*0.5; 140*5754ecacSJeremy L Thompson // *INDENT-OFF* 141*5754ecacSJeremy L Thompson const CeedScalar sigma[3][3] = {{sigma00, sigma01, sigma02}, 142*5754ecacSJeremy L Thompson {sigma01, sigma11, sigma12}, 143*5754ecacSJeremy L Thompson {sigma02, sigma12, sigma22} 144*5754ecacSJeremy L Thompson }; 145*5754ecacSJeremy L Thompson // *INDENT-ON* 146*5754ecacSJeremy L Thompson 147*5754ecacSJeremy L Thompson // Apply dXdx^T and weight to sigma 148*5754ecacSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) // Component 149*5754ecacSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { // Derivative 150*5754ecacSJeremy L Thompson dvdX[k][j][i] = 0; 151*5754ecacSJeremy L Thompson for (CeedInt m = 0; m < 3; m++) 152*5754ecacSJeremy L Thompson dvdX[k][j][i] += dXdx[k][m] * sigma[j][m] * wdetJ; 153*5754ecacSJeremy L Thompson } 154*5754ecacSJeremy L Thompson 155*5754ecacSJeremy L Thompson } // End of Quadrature Point Loop 156*5754ecacSJeremy L Thompson 157*5754ecacSJeremy L Thompson return 0; 158*5754ecacSJeremy L Thompson } 159*5754ecacSJeremy L Thompson 160*5754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 161*5754ecacSJeremy L Thompson // Jacobian evaluation for linear elasticity 162*5754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 163*5754ecacSJeremy L Thompson CEED_QFUNCTION(ElasLineardF)(void *ctx, CeedInt Q, const CeedScalar *const *in, 164*5754ecacSJeremy L Thompson CeedScalar *const *out) { 165*5754ecacSJeremy L Thompson // *INDENT-OFF* 166*5754ecacSJeremy L Thompson // Inputs 167*5754ecacSJeremy L Thompson const CeedScalar (*deltaug)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[0], 168*5754ecacSJeremy L Thompson (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 169*5754ecacSJeremy L Thompson // grad_u not used for linear elasticity 170*5754ecacSJeremy L Thompson // (*grad_u)[3][Q] = (CeedScalar(*)[3][Q])in[2]; 171*5754ecacSJeremy L Thompson 172*5754ecacSJeremy L Thompson // Outputs 173*5754ecacSJeremy L Thompson CeedScalar (*deltadvdX)[3][CEED_Q_VLA] = (CeedScalar(*)[3][CEED_Q_VLA])out[0]; 174*5754ecacSJeremy L Thompson // *INDENT-ON* 175*5754ecacSJeremy L Thompson 176*5754ecacSJeremy L Thompson // Context 177*5754ecacSJeremy L Thompson const Physics context = (Physics)ctx; 178*5754ecacSJeremy L Thompson const CeedScalar E = context->E; 179*5754ecacSJeremy L Thompson const CeedScalar nu = context->nu; 180*5754ecacSJeremy L Thompson 181*5754ecacSJeremy L Thompson // Quadrature Point Loop 182*5754ecacSJeremy L Thompson CeedPragmaSIMD 183*5754ecacSJeremy L Thompson for (CeedInt i=0; i<Q; i++) { 184*5754ecacSJeremy L Thompson // Read spatial derivatives of u 185*5754ecacSJeremy L Thompson // *INDENT-OFF* 186*5754ecacSJeremy L Thompson const CeedScalar deltadu[3][3] = {{deltaug[0][0][i], 187*5754ecacSJeremy L Thompson deltaug[1][0][i], 188*5754ecacSJeremy L Thompson deltaug[2][0][i]}, 189*5754ecacSJeremy L Thompson {deltaug[0][1][i], 190*5754ecacSJeremy L Thompson deltaug[1][1][i], 191*5754ecacSJeremy L Thompson deltaug[2][1][i]}, 192*5754ecacSJeremy L Thompson {deltaug[0][2][i], 193*5754ecacSJeremy L Thompson deltaug[1][2][i], 194*5754ecacSJeremy L Thompson deltaug[2][2][i]} 195*5754ecacSJeremy L Thompson }; 196*5754ecacSJeremy L Thompson // -- Qdata 197*5754ecacSJeremy L Thompson const CeedScalar wdetJ = q_data[0][i]; 198*5754ecacSJeremy L Thompson const CeedScalar dXdx[3][3] = {{q_data[1][i], 199*5754ecacSJeremy L Thompson q_data[2][i], 200*5754ecacSJeremy L Thompson q_data[3][i]}, 201*5754ecacSJeremy L Thompson {q_data[4][i], 202*5754ecacSJeremy L Thompson q_data[5][i], 203*5754ecacSJeremy L Thompson q_data[6][i]}, 204*5754ecacSJeremy L Thompson {q_data[7][i], 205*5754ecacSJeremy L Thompson q_data[8][i], 206*5754ecacSJeremy L Thompson q_data[9][i]} 207*5754ecacSJeremy L Thompson }; 208*5754ecacSJeremy L Thompson // *INDENT-ON* 209*5754ecacSJeremy L Thompson 210*5754ecacSJeremy L Thompson // Compute graddeltau 211*5754ecacSJeremy L Thompson // dXdx = (dx/dX)^(-1) 212*5754ecacSJeremy L Thompson // Apply dXdx to deltadu = graddeltau 213*5754ecacSJeremy L Thompson CeedScalar graddeltau[3][3]; 214*5754ecacSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) // Component 215*5754ecacSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { // Derivative 216*5754ecacSJeremy L Thompson graddeltau[j][k] = 0; 217*5754ecacSJeremy L Thompson for (CeedInt m = 0; m < 3; m++) 218*5754ecacSJeremy L Thompson graddeltau[j][k] += dXdx[m][k] * deltadu[j][m]; 219*5754ecacSJeremy L Thompson } 220*5754ecacSJeremy L Thompson 221*5754ecacSJeremy L Thompson // Compute Strain : e (epsilon) 222*5754ecacSJeremy L Thompson // e = 1/2 (grad u + (grad u)^T) 223*5754ecacSJeremy L Thompson // *INDENT-OFF* 224*5754ecacSJeremy L Thompson const CeedScalar de[3][3] = {{(graddeltau[0][0] + graddeltau[0][0])/2., 225*5754ecacSJeremy L Thompson (graddeltau[0][1] + graddeltau[1][0])/2., 226*5754ecacSJeremy L Thompson (graddeltau[0][2] + graddeltau[2][0])/2.}, 227*5754ecacSJeremy L Thompson {(graddeltau[1][0] + graddeltau[0][1])/2., 228*5754ecacSJeremy L Thompson (graddeltau[1][1] + graddeltau[1][1])/2., 229*5754ecacSJeremy L Thompson (graddeltau[1][2] + graddeltau[2][1])/2.}, 230*5754ecacSJeremy L Thompson {(graddeltau[2][0] + graddeltau[0][2])/2., 231*5754ecacSJeremy L Thompson (graddeltau[2][1] + graddeltau[1][2])/2., 232*5754ecacSJeremy L Thompson (graddeltau[2][2] + graddeltau[2][2])/2.} 233*5754ecacSJeremy L Thompson }; 234*5754ecacSJeremy L Thompson 235*5754ecacSJeremy L Thompson // *INDENT-ON* 236*5754ecacSJeremy L Thompson // Formulation uses Voigt notation: 237*5754ecacSJeremy L Thompson // stress (sigma) strain (epsilon) 238*5754ecacSJeremy L Thompson // 239*5754ecacSJeremy L Thompson // [dsigma00] [de00] 240*5754ecacSJeremy L Thompson // [dsigma11] [de11] 241*5754ecacSJeremy L Thompson // [dsigma22] = S * [de22] 242*5754ecacSJeremy L Thompson // [dsigma12] [de12] 243*5754ecacSJeremy L Thompson // [dsigma02] [de02] 244*5754ecacSJeremy L Thompson // [dsigma01] [de01] 245*5754ecacSJeremy L Thompson // 246*5754ecacSJeremy L Thompson // where 247*5754ecacSJeremy L Thompson // [1-nu nu nu ] 248*5754ecacSJeremy L Thompson // [ nu 1-nu nu ] 249*5754ecacSJeremy L Thompson // [ nu nu 1-nu ] 250*5754ecacSJeremy L Thompson // S = E/((1+nu)*(1-2*nu)) [ (1-2*nu)/2 ] 251*5754ecacSJeremy L Thompson // [ (1-2*nu)/2 ] 252*5754ecacSJeremy L Thompson // [ (1-2*nu)/2 ] 253*5754ecacSJeremy L Thompson 254*5754ecacSJeremy L Thompson // Above Voigt Notation is placed in a 3x3 matrix: 255*5754ecacSJeremy L Thompson const CeedScalar ss = E / ((1 + nu)*(1 - 2*nu)); 256*5754ecacSJeremy L Thompson const CeedScalar dsigma00 = ss*((1 - nu)*de[0][0]+nu*de[1][1]+nu*de[2][2]), 257*5754ecacSJeremy L Thompson dsigma11 = ss*(nu*de[0][0]+(1 - nu)*de[1][1]+nu*de[2][2]), 258*5754ecacSJeremy L Thompson dsigma22 = ss*(nu*de[0][0]+nu*de[1][1]+(1 - nu)*de[2][2]), 259*5754ecacSJeremy L Thompson dsigma12 = ss*(1 - 2*nu)*de[1][2] / 2, 260*5754ecacSJeremy L Thompson dsigma02 = ss*(1 - 2*nu)*de[0][2] / 2, 261*5754ecacSJeremy L Thompson dsigma01 = ss*(1 - 2*nu)*de[0][1] / 2; 262*5754ecacSJeremy L Thompson // *INDENT-OFF* 263*5754ecacSJeremy L Thompson const CeedScalar dsigma[3][3] = {{dsigma00, dsigma01, dsigma02}, 264*5754ecacSJeremy L Thompson {dsigma01, dsigma11, dsigma12}, 265*5754ecacSJeremy L Thompson {dsigma02, dsigma12, dsigma22} 266*5754ecacSJeremy L Thompson }; 267*5754ecacSJeremy L Thompson // *INDENT-ON* 268*5754ecacSJeremy L Thompson 269*5754ecacSJeremy L Thompson // Apply dXdx^T and weight 270*5754ecacSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) // Component 271*5754ecacSJeremy L Thompson for (CeedInt k = 0; k < 3 ; k++) { // Derivative 272*5754ecacSJeremy L Thompson deltadvdX[k][j][i] = 0; 273*5754ecacSJeremy L Thompson for (CeedInt m = 0; m < 3; m++) 274*5754ecacSJeremy L Thompson deltadvdX[k][j][i] += dXdx[k][m] * dsigma[j][m] * wdetJ; 275*5754ecacSJeremy L Thompson } 276*5754ecacSJeremy L Thompson 277*5754ecacSJeremy L Thompson } // End of Quadrature Point Loop 278*5754ecacSJeremy L Thompson 279*5754ecacSJeremy L Thompson return 0; 280*5754ecacSJeremy L Thompson } 281*5754ecacSJeremy L Thompson 282*5754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 283*5754ecacSJeremy L Thompson // Strain energy computation for linear elasticity 284*5754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 285*5754ecacSJeremy L Thompson CEED_QFUNCTION(ElasLinearEnergy)(void *ctx, CeedInt Q, 286*5754ecacSJeremy L Thompson const CeedScalar *const *in, 287*5754ecacSJeremy L Thompson CeedScalar *const *out) { 288*5754ecacSJeremy L Thompson // *INDENT-OFF* 289*5754ecacSJeremy L Thompson // Inputs 290*5754ecacSJeremy L Thompson const CeedScalar (*ug)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[0], 291*5754ecacSJeremy L Thompson (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 292*5754ecacSJeremy L Thompson 293*5754ecacSJeremy L Thompson // Outputs 294*5754ecacSJeremy L Thompson CeedScalar (*energy) = (CeedScalar(*))out[0]; 295*5754ecacSJeremy L Thompson // *INDENT-ON* 296*5754ecacSJeremy L Thompson 297*5754ecacSJeremy L Thompson // Context 298*5754ecacSJeremy L Thompson const Physics context = (Physics)ctx; 299*5754ecacSJeremy L Thompson const CeedScalar E = context->E; 300*5754ecacSJeremy L Thompson const CeedScalar nu = context->nu; 301*5754ecacSJeremy L Thompson 302*5754ecacSJeremy L Thompson // Constants 303*5754ecacSJeremy L Thompson const CeedScalar TwoMu = E / (1 + nu); 304*5754ecacSJeremy L Thompson const CeedScalar mu = TwoMu / 2; 305*5754ecacSJeremy L Thompson const CeedScalar Kbulk = E / (3*(1 - 2*nu)); // Bulk Modulus 306*5754ecacSJeremy L Thompson const CeedScalar lambda = (3*Kbulk - TwoMu) / 3; 307*5754ecacSJeremy L Thompson 308*5754ecacSJeremy L Thompson // Quadrature Point Loop 309*5754ecacSJeremy L Thompson CeedPragmaSIMD 310*5754ecacSJeremy L Thompson for (CeedInt i=0; i<Q; i++) { 311*5754ecacSJeremy L Thompson // Read spatial derivatives of u 312*5754ecacSJeremy L Thompson // *INDENT-OFF* 313*5754ecacSJeremy L Thompson const CeedScalar du[3][3] = {{ug[0][0][i], 314*5754ecacSJeremy L Thompson ug[1][0][i], 315*5754ecacSJeremy L Thompson ug[2][0][i]}, 316*5754ecacSJeremy L Thompson {ug[0][1][i], 317*5754ecacSJeremy L Thompson ug[1][1][i], 318*5754ecacSJeremy L Thompson ug[2][1][i]}, 319*5754ecacSJeremy L Thompson {ug[0][2][i], 320*5754ecacSJeremy L Thompson ug[1][2][i], 321*5754ecacSJeremy L Thompson ug[2][2][i]} 322*5754ecacSJeremy L Thompson }; 323*5754ecacSJeremy L Thompson // -- Qdata 324*5754ecacSJeremy L Thompson const CeedScalar wdetJ = q_data[0][i]; 325*5754ecacSJeremy L Thompson const CeedScalar dXdx[3][3] = {{q_data[1][i], 326*5754ecacSJeremy L Thompson q_data[2][i], 327*5754ecacSJeremy L Thompson q_data[3][i]}, 328*5754ecacSJeremy L Thompson {q_data[4][i], 329*5754ecacSJeremy L Thompson q_data[5][i], 330*5754ecacSJeremy L Thompson q_data[6][i]}, 331*5754ecacSJeremy L Thompson {q_data[7][i], 332*5754ecacSJeremy L Thompson q_data[8][i], 333*5754ecacSJeremy L Thompson q_data[9][i]} 334*5754ecacSJeremy L Thompson }; 335*5754ecacSJeremy L Thompson // *INDENT-ON* 336*5754ecacSJeremy L Thompson 337*5754ecacSJeremy L Thompson // Compute grad_u 338*5754ecacSJeremy L Thompson // dXdx = (dx/dX)^(-1) 339*5754ecacSJeremy L Thompson // Apply dXdx to du = grad_u 340*5754ecacSJeremy L Thompson CeedScalar grad_u[3][3]; 341*5754ecacSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) // Component 342*5754ecacSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { // Derivative 343*5754ecacSJeremy L Thompson grad_u[j][k] = 0; 344*5754ecacSJeremy L Thompson for (CeedInt m = 0; m < 3; m++) 345*5754ecacSJeremy L Thompson grad_u[j][k] += dXdx[m][k] * du[j][m]; 346*5754ecacSJeremy L Thompson } 347*5754ecacSJeremy L Thompson 348*5754ecacSJeremy L Thompson // Compute Strain : e (epsilon) 349*5754ecacSJeremy L Thompson // e = 1/2 (grad u + (grad u)^T) 350*5754ecacSJeremy L Thompson 351*5754ecacSJeremy L Thompson // *INDENT-OFF* 352*5754ecacSJeremy L Thompson const CeedScalar e[3][3] = {{(grad_u[0][0] + grad_u[0][0])/2., 353*5754ecacSJeremy L Thompson (grad_u[0][1] + grad_u[1][0])/2., 354*5754ecacSJeremy L Thompson (grad_u[0][2] + grad_u[2][0])/2.}, 355*5754ecacSJeremy L Thompson {(grad_u[1][0] + grad_u[0][1])/2., 356*5754ecacSJeremy L Thompson (grad_u[1][1] + grad_u[1][1])/2., 357*5754ecacSJeremy L Thompson (grad_u[1][2] + grad_u[2][1])/2.}, 358*5754ecacSJeremy L Thompson {(grad_u[2][0] + grad_u[0][2])/2., 359*5754ecacSJeremy L Thompson (grad_u[2][1] + grad_u[1][2])/2., 360*5754ecacSJeremy L Thompson (grad_u[2][2] + grad_u[2][2])/2.} 361*5754ecacSJeremy L Thompson }; 362*5754ecacSJeremy L Thompson // *INDENT-ON* 363*5754ecacSJeremy L Thompson 364*5754ecacSJeremy L Thompson // Strain energy 365*5754ecacSJeremy L Thompson const CeedScalar strain_vol = e[0][0] + e[1][1] + e[2][2]; 366*5754ecacSJeremy L Thompson energy[i] = (lambda*strain_vol*strain_vol/2. + strain_vol*mu + 367*5754ecacSJeremy L Thompson (e[0][1]*e[0][1]+e[0][2]*e[0][2]+e[1][2]*e[1][2])*2*mu)*wdetJ; 368*5754ecacSJeremy L Thompson 369*5754ecacSJeremy L Thompson } // End of Quadrature Point Loop 370*5754ecacSJeremy L Thompson 371*5754ecacSJeremy L Thompson return 0; 372*5754ecacSJeremy L Thompson } 373*5754ecacSJeremy L Thompson 374*5754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 375*5754ecacSJeremy L Thompson // Nodal diagnostic quantities for linear elasticity 376*5754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 377*5754ecacSJeremy L Thompson CEED_QFUNCTION(ElasLinearDiagnostic)(void *ctx, CeedInt Q, 378*5754ecacSJeremy L Thompson const CeedScalar *const *in, 379*5754ecacSJeremy L Thompson CeedScalar *const *out) { 380*5754ecacSJeremy L Thompson // *INDENT-OFF* 381*5754ecacSJeremy L Thompson // Inputs 382*5754ecacSJeremy L Thompson const CeedScalar (*u)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0], 383*5754ecacSJeremy L Thompson (*ug)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[1], 384*5754ecacSJeremy L Thompson (*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 385*5754ecacSJeremy L Thompson 386*5754ecacSJeremy L Thompson // Outputs 387*5754ecacSJeremy L Thompson CeedScalar (*diagnostic)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 388*5754ecacSJeremy L Thompson // *INDENT-ON* 389*5754ecacSJeremy L Thompson 390*5754ecacSJeremy L Thompson // Context 391*5754ecacSJeremy L Thompson const Physics context = (Physics)ctx; 392*5754ecacSJeremy L Thompson const CeedScalar E = context->E; 393*5754ecacSJeremy L Thompson const CeedScalar nu = context->nu; 394*5754ecacSJeremy L Thompson 395*5754ecacSJeremy L Thompson // Constants 396*5754ecacSJeremy L Thompson const CeedScalar TwoMu = E / (1 + nu); 397*5754ecacSJeremy L Thompson const CeedScalar mu = TwoMu / 2; 398*5754ecacSJeremy L Thompson const CeedScalar Kbulk = E / (3*(1 - 2*nu)); // Bulk Modulus 399*5754ecacSJeremy L Thompson const CeedScalar lambda = (3*Kbulk - TwoMu) / 3; 400*5754ecacSJeremy L Thompson 401*5754ecacSJeremy L Thompson // Quadrature Point Loop 402*5754ecacSJeremy L Thompson CeedPragmaSIMD 403*5754ecacSJeremy L Thompson for (CeedInt i=0; i<Q; i++) { 404*5754ecacSJeremy L Thompson // Read spatial derivatives of u 405*5754ecacSJeremy L Thompson // *INDENT-OFF* 406*5754ecacSJeremy L Thompson const CeedScalar du[3][3] = {{ug[0][0][i], 407*5754ecacSJeremy L Thompson ug[1][0][i], 408*5754ecacSJeremy L Thompson ug[2][0][i]}, 409*5754ecacSJeremy L Thompson {ug[0][1][i], 410*5754ecacSJeremy L Thompson ug[1][1][i], 411*5754ecacSJeremy L Thompson ug[2][1][i]}, 412*5754ecacSJeremy L Thompson {ug[0][2][i], 413*5754ecacSJeremy L Thompson ug[1][2][i], 414*5754ecacSJeremy L Thompson ug[2][2][i]} 415*5754ecacSJeremy L Thompson }; 416*5754ecacSJeremy L Thompson // -- Qdata 417*5754ecacSJeremy L Thompson const CeedScalar dXdx[3][3] = {{q_data[1][i], 418*5754ecacSJeremy L Thompson q_data[2][i], 419*5754ecacSJeremy L Thompson q_data[3][i]}, 420*5754ecacSJeremy L Thompson {q_data[4][i], 421*5754ecacSJeremy L Thompson q_data[5][i], 422*5754ecacSJeremy L Thompson q_data[6][i]}, 423*5754ecacSJeremy L Thompson {q_data[7][i], 424*5754ecacSJeremy L Thompson q_data[8][i], 425*5754ecacSJeremy L Thompson q_data[9][i]} 426*5754ecacSJeremy L Thompson }; 427*5754ecacSJeremy L Thompson // *INDENT-ON* 428*5754ecacSJeremy L Thompson 429*5754ecacSJeremy L Thompson // Compute grad_u 430*5754ecacSJeremy L Thompson // dXdx = (dx/dX)^(-1) 431*5754ecacSJeremy L Thompson // Apply dXdx to du = grad_u 432*5754ecacSJeremy L Thompson CeedScalar grad_u[3][3]; 433*5754ecacSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) // Component 434*5754ecacSJeremy L Thompson for (CeedInt k = 0; k < 3; k++) { // Derivative 435*5754ecacSJeremy L Thompson grad_u[j][k] = 0; 436*5754ecacSJeremy L Thompson for (CeedInt m = 0; m < 3; m++) 437*5754ecacSJeremy L Thompson grad_u[j][k] += dXdx[m][k] * du[j][m]; 438*5754ecacSJeremy L Thompson } 439*5754ecacSJeremy L Thompson 440*5754ecacSJeremy L Thompson // Compute Strain : e (epsilon) 441*5754ecacSJeremy L Thompson // e = 1/2 (grad u + (grad u)^T) 442*5754ecacSJeremy L Thompson 443*5754ecacSJeremy L Thompson // *INDENT-OFF* 444*5754ecacSJeremy L Thompson const CeedScalar e[3][3] = {{(grad_u[0][0] + grad_u[0][0])/2., 445*5754ecacSJeremy L Thompson (grad_u[0][1] + grad_u[1][0])/2., 446*5754ecacSJeremy L Thompson (grad_u[0][2] + grad_u[2][0])/2.}, 447*5754ecacSJeremy L Thompson {(grad_u[1][0] + grad_u[0][1])/2., 448*5754ecacSJeremy L Thompson (grad_u[1][1] + grad_u[1][1])/2., 449*5754ecacSJeremy L Thompson (grad_u[1][2] + grad_u[2][1])/2.}, 450*5754ecacSJeremy L Thompson {(grad_u[2][0] + grad_u[0][2])/2., 451*5754ecacSJeremy L Thompson (grad_u[2][1] + grad_u[1][2])/2., 452*5754ecacSJeremy L Thompson (grad_u[2][2] + grad_u[2][2])/2.} 453*5754ecacSJeremy L Thompson }; 454*5754ecacSJeremy L Thompson // *INDENT-ON* 455*5754ecacSJeremy L Thompson 456*5754ecacSJeremy L Thompson // Displacement 457*5754ecacSJeremy L Thompson diagnostic[0][i] = u[0][i]; 458*5754ecacSJeremy L Thompson diagnostic[1][i] = u[1][i]; 459*5754ecacSJeremy L Thompson diagnostic[2][i] = u[2][i]; 460*5754ecacSJeremy L Thompson 461*5754ecacSJeremy L Thompson // Pressure 462*5754ecacSJeremy L Thompson const CeedScalar strain_vol = e[0][0] + e[1][1] + e[2][2]; 463*5754ecacSJeremy L Thompson diagnostic[3][i] = -lambda*strain_vol; 464*5754ecacSJeremy L Thompson 465*5754ecacSJeremy L Thompson // Stress tensor invariants 466*5754ecacSJeremy L Thompson diagnostic[4][i] = strain_vol; 467*5754ecacSJeremy L Thompson diagnostic[5][i] = 0.; 468*5754ecacSJeremy L Thompson for (CeedInt j = 0; j < 3; j++) 469*5754ecacSJeremy L Thompson for (CeedInt m = 0; m < 3; m++) 470*5754ecacSJeremy L Thompson diagnostic[5][i] += e[j][m] * e[m][j]; 471*5754ecacSJeremy L Thompson diagnostic[6][i] = 1 + strain_vol; 472*5754ecacSJeremy L Thompson 473*5754ecacSJeremy L Thompson // Strain energy 474*5754ecacSJeremy L Thompson diagnostic[7][i] = (lambda*strain_vol*strain_vol/2. + strain_vol*mu + 475*5754ecacSJeremy L Thompson (e[0][1]*e[0][1]+e[0][2]*e[0][2]+e[1][2]*e[1][2])*2*mu); 476*5754ecacSJeremy L Thompson 477*5754ecacSJeremy L Thompson } // End of Quadrature Point Loop 478*5754ecacSJeremy L Thompson 479*5754ecacSJeremy L Thompson return 0; 480*5754ecacSJeremy L Thompson } 481*5754ecacSJeremy L Thompson // ----------------------------------------------------------------------------- 482*5754ecacSJeremy L Thompson 483*5754ecacSJeremy L Thompson #endif //End of ELAS_LINEAR_H 484