1 // Copyright (c) 2017-2025, 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 #pragma once 8 9 #include <ceed.h> 10 #include "ex-common.h" 11 12 /// libCEED Q-function for building quadrature data for a mass + diffusion operator 13 CEED_QFUNCTION(build_mass_diff)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 14 // in[0] is Jacobians with shape [dim, dim, Q] 15 // in[1] is quadrature weights, size (Q) 16 const CeedScalar *w = in[1]; 17 CeedScalar(*q_data)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 18 struct BuildContext *build_data = (struct BuildContext *)ctx; 19 20 // At every quadrature point, compute w/det(J).adj(J).adj(J)^T and store 21 // the symmetric part of the result. 22 switch (build_data->dim + 10 * build_data->space_dim) { 23 case 11: { // dim = 1, space_dim = 1 24 const CeedScalar(*J)[1][CEED_Q_VLA] = (const CeedScalar(*)[1][CEED_Q_VLA])in[0]; 25 26 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 27 // Mass 28 q_data[0][i] = w[i] * J[0][0][i]; 29 30 // Diffusion 31 q_data[1][i] = w[i] / J[0][0][i]; 32 } 33 } break; 34 case 22: { // dim = 2, space_dim = 2 35 const CeedScalar(*J)[2][CEED_Q_VLA] = (const CeedScalar(*)[2][CEED_Q_VLA])in[0]; 36 37 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 38 // J: 0 2 q_data: 0 2 adj(J): J22 -J12 39 // 1 3 2 1 -J10 J00 40 const CeedScalar J00 = J[0][0][i]; 41 const CeedScalar J10 = J[0][1][i]; 42 const CeedScalar J01 = J[1][0][i]; 43 const CeedScalar J11 = J[1][1][i]; 44 const CeedScalar qw = w[i] / (J00 * J11 - J10 * J01); 45 46 // Mass 47 q_data[0][i] = w[i] * (J00 * J11 - J10 * J01); 48 49 // Diffusion 50 q_data[1][i] = qw * (J01 * J01 + J11 * J11); 51 q_data[2][i] = qw * (J00 * J00 + J10 * J10); 52 q_data[3][i] = -qw * (J00 * J01 + J10 * J11); 53 } 54 } break; 55 case 33: { // dim = 3, space_dim = 3 56 const CeedScalar(*J)[3][CEED_Q_VLA] = (const CeedScalar(*)[3][CEED_Q_VLA])in[0]; 57 58 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 59 // Compute the adjoint 60 CeedScalar A[3][3]; 61 for (CeedInt j = 0; j < 3; j++) { 62 for (CeedInt k = 0; k < 3; k++) { 63 A[k][j] = 64 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]; 65 } 66 } 67 68 // Compute quadrature weight / det(J) 69 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]); 70 71 // Mass 72 q_data[0][i] = w[i] * (J[0][0][i] * A[0][0] + J[0][1][i] * A[0][1] + J[0][2][i] * A[0][2]); 73 74 // Diffusion 75 // Stored in Voigt convention 76 // 1 6 5 77 // 6 2 4 78 // 5 4 3 79 q_data[1][i] = qw * (A[0][0] * A[0][0] + A[0][1] * A[0][1] + A[0][2] * A[0][2]); 80 q_data[2][i] = qw * (A[1][0] * A[1][0] + A[1][1] * A[1][1] + A[1][2] * A[1][2]); 81 q_data[3][i] = qw * (A[2][0] * A[2][0] + A[2][1] * A[2][1] + A[2][2] * A[2][2]); 82 q_data[4][i] = qw * (A[1][0] * A[2][0] + A[1][1] * A[2][1] + A[1][2] * A[2][2]); 83 q_data[5][i] = qw * (A[0][0] * A[2][0] + A[0][1] * A[2][1] + A[0][2] * A[2][2]); 84 q_data[6][i] = qw * (A[0][0] * A[1][0] + A[0][1] * A[1][1] + A[0][2] * A[1][2]); 85 } 86 } break; 87 } 88 return CEED_ERROR_SUCCESS; 89 } 90 91 /// libCEED Q-function for applying a mass + diffusion operator 92 CEED_QFUNCTION(apply_mass_diff)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 93 struct BuildContext *build_data = (struct BuildContext *)ctx; 94 // in[0], out[0] solution values with shape [1, 1, Q] 95 // in[1], out[1] solution gradients with shape [dim, 1, Q] 96 // in[2] is quadrature data with shape [num_components, Q] 97 const CeedScalar(*q_data)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[2]; 98 99 switch (build_data->dim) { 100 case 1: { 101 const CeedScalar *u = in[0], *ug = in[1]; 102 CeedScalar *v = out[0], *vg = out[1]; 103 104 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 105 // Mass 106 v[i] = q_data[0][i] * u[i]; 107 108 // Diffusion 109 vg[i] = q_data[1][i] * ug[i]; 110 } 111 } break; 112 case 2: { 113 const CeedScalar *u = in[0]; 114 const CeedScalar(*ug)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 115 CeedScalar *v = out[0]; 116 CeedScalar(*vg)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[1]; 117 118 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 119 // Mass 120 v[i] = q_data[0][i] * u[i]; 121 122 // Diffusion 123 // Read q_data (dXdxdXdx_T symmetric matrix) 124 // Stored in Voigt convention 125 // 1 3 126 // 3 2 127 const CeedScalar dXdxdXdx_T[2][2] = { 128 {q_data[1][i], q_data[3][i]}, 129 {q_data[3][i], q_data[2][i]} 130 }; 131 132 // j = direction of vg 133 for (int j = 0; j < 2; j++) { 134 vg[j][i] = (ug[0][i] * dXdxdXdx_T[0][j] + ug[1][i] * dXdxdXdx_T[1][j]); 135 } 136 } 137 } break; 138 case 3: { 139 const CeedScalar *u = in[0]; 140 const CeedScalar(*ug)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[1]; 141 CeedScalar *v = out[0]; 142 CeedScalar(*vg)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[1]; 143 144 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 145 // Mass 146 v[i] = q_data[0][i] * u[i]; 147 148 // Diffusion 149 // Read q_data (dXdxdXdx_T symmetric matrix) 150 // Stored in Voigt convention 151 // 1 6 5 152 // 6 2 4 153 // 5 4 3 154 const CeedScalar dXdxdXdx_T[3][3] = { 155 {q_data[1][i], q_data[6][i], q_data[5][i]}, 156 {q_data[6][i], q_data[2][i], q_data[4][i]}, 157 {q_data[5][i], q_data[4][i], q_data[3][i]} 158 }; 159 160 // j = direction of vg 161 for (int j = 0; j < 3; j++) { 162 vg[j][i] = (ug[0][i] * dXdxdXdx_T[0][j] + ug[1][i] * dXdxdXdx_T[1][j] + ug[2][i] * dXdxdXdx_T[2][j]); 163 } 164 } 165 } break; 166 } 167 return CEED_ERROR_SUCCESS; 168 } 169