1 // Copyright (c) 2017-2026, Lawrence Livermore National Security, LLC and other CEED contributors. 2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3 // 4 // SPDX-License-Identifier: BSD-2-Clause 5 // 6 // This file is part of CEED: http://github.com/ceed 7 8 #include <ceed/types.h> 9 10 /// A structure used to pass additional data to f_build_diff 11 struct BuildContext { 12 CeedInt dim, space_dim; 13 }; 14 15 /// libCEED Q-function for building quadrature data for a diffusion operator 16 CEED_QFUNCTION(build_diff)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 17 // in[0] is Jacobians with shape [dim, dim, Q] 18 // in[1] is quadrature weights, size (Q) 19 const CeedScalar *w = in[1]; 20 CeedScalar(*q_data)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 21 struct BuildContext *build_data = (struct BuildContext *)ctx; 22 23 // At every quadrature point, compute w/det(J).adj(J).adj(J)^T and store 24 // the symmetric part of the result. 25 switch (build_data->dim + 10 * build_data->space_dim) { 26 case 11: { 27 const CeedScalar(*J)[1][CEED_Q_VLA] = (const CeedScalar(*)[1][CEED_Q_VLA])in[0]; 28 29 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { q_data[0][i] = w[i] / J[0][0][i]; } // End of Quadrature Point Loop 30 } break; 31 case 22: { 32 const CeedScalar(*J)[2][CEED_Q_VLA] = (const CeedScalar(*)[2][CEED_Q_VLA])in[0]; 33 34 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 35 // J: 0 2 q_data: 0 2 adj(J): J11 -J01 36 // 1 3 2 1 -J10 J00 37 const CeedScalar J00 = J[0][0][i]; 38 const CeedScalar J10 = J[0][1][i]; 39 const CeedScalar J01 = J[1][0][i]; 40 const CeedScalar J11 = J[1][1][i]; 41 const CeedScalar qw = w[i] / (J00 * J11 - J10 * J01); 42 43 q_data[0][i] = qw * (J01 * J01 + J11 * J11); 44 q_data[1][i] = qw * (J00 * J00 + J10 * J10); 45 q_data[2][i] = -qw * (J00 * J01 + J10 * J11); 46 } // End of Quadrature Point Loop 47 } break; 48 case 33: { 49 const CeedScalar(*J)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[0]; 50 51 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 52 // Compute the adjoint 53 CeedScalar A[3][3]; 54 55 for (CeedInt j = 0; j < 3; j++) { 56 for (CeedInt k = 0; k < 3; k++) { 57 // Equivalent code with J as a VLA and no mod operations: 58 // A[k][j] = J[j+1][k+1]*J[j+2][k+2] - J[j+1][k+2]*J[j+2][k+1] 59 A[k][j] = 60 J[(k + 1) % 3][(j + 1) % 3][i] * J[(k + 2) % 3][(j + 2) % 3][i] - J[(k + 2) % 3][(j + 1) % 3][i] * J[(k + 1) % 3][(j + 2) % 3][i]; 61 } 62 } 63 64 // Compute quadrature weight / det(J) 65 const CeedScalar qw = w[i] / (J[0][0][i] * A[0][0] + J[0][1][i] * A[0][1] + J[0][2][i] * A[0][2]); 66 67 // Compute geometric factors 68 // Stored in Voigt convention 69 // 0 5 4 70 // 5 1 3 71 // 4 3 2 72 q_data[0][i] = qw * (A[0][0] * A[0][0] + A[0][1] * A[0][1] + A[0][2] * A[0][2]); 73 q_data[1][i] = qw * (A[1][0] * A[1][0] + A[1][1] * A[1][1] + A[1][2] * A[1][2]); 74 q_data[2][i] = qw * (A[2][0] * A[2][0] + A[2][1] * A[2][1] + A[2][2] * A[2][2]); 75 q_data[3][i] = qw * (A[1][0] * A[2][0] + A[1][1] * A[2][1] + A[1][2] * A[2][2]); 76 q_data[4][i] = qw * (A[0][0] * A[2][0] + A[0][1] * A[2][1] + A[0][2] * A[2][2]); 77 q_data[5][i] = qw * (A[0][0] * A[1][0] + A[0][1] * A[1][1] + A[0][2] * A[1][2]); 78 } // End of Quadrature Point Loop 79 } break; 80 } 81 return CEED_ERROR_SUCCESS; 82 } 83 84 /// libCEED Q-function for applying a diff operator 85 CEED_QFUNCTION(apply_diff)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 86 struct BuildContext *build_data = (struct BuildContext *)ctx; 87 // in[0], out[0] solution gradients with shape [dim, 1, Q] 88 // in[1] is quadrature data with shape [num_components, Q] 89 const CeedScalar(*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 90 91 switch (build_data->dim) { 92 case 1: { 93 const CeedScalar *ug = in[0]; 94 CeedScalar *vg = out[0]; 95 96 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { vg[i] = ug[i] * q_data[0][i]; } // End of Quadrature Point Loop 97 } break; 98 case 2: { 99 const CeedScalar(*ug)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 100 CeedScalar(*vg)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 101 102 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 103 // Read q_data (dXdxdXdx_T symmetric matrix) 104 // Stored in Voigt convention 105 // 0 2 106 // 2 1 107 const CeedScalar dXdxdXdx_T[2][2] = { 108 {q_data[0][i], q_data[2][i]}, 109 {q_data[2][i], q_data[1][i]} 110 }; 111 112 // j = direction of vg 113 for (int j = 0; j < 2; j++) vg[j][i] = (ug[0][i] * dXdxdXdx_T[0][j] + ug[1][i] * dXdxdXdx_T[1][j]); 114 } // End of Quadrature Point Loop 115 } break; 116 case 3: { 117 const CeedScalar(*ug)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 118 CeedScalar(*vg)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 119 120 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 121 // Read q_data (dXdxdXdx_T symmetric matrix) 122 // Stored in Voigt convention 123 // 0 5 4 124 // 5 1 3 125 // 4 3 2 126 const CeedScalar dXdxdXdx_T[3][3] = { 127 {q_data[0][i], q_data[5][i], q_data[4][i]}, 128 {q_data[5][i], q_data[1][i], q_data[3][i]}, 129 {q_data[4][i], q_data[3][i], q_data[2][i]} 130 }; 131 132 // j = direction of vg 133 for (int j = 0; j < 3; j++) vg[j][i] = (ug[0][i] * dXdxdXdx_T[0][j] + ug[1][i] * dXdxdXdx_T[1][j] + ug[2][i] * dXdxdXdx_T[2][j]); 134 } // End of Quadrature Point Loop 135 } break; 136 } 137 return CEED_ERROR_SUCCESS; 138 } 139