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