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