xref: /honee/qfunctions/utils.h (revision 8c85b83518484fc9eaa5bccfe38a7cce91b34691)
1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3 #pragma once
4 
5 #include <ceed.h>
6 #include <math.h>
7 
8 #ifndef M_PI
9 #define M_PI 3.14159265358979323846
10 #endif
11 
12 CEED_QFUNCTION_HELPER CeedScalar Max(CeedScalar a, CeedScalar b) { return a < b ? b : a; }
13 CEED_QFUNCTION_HELPER CeedScalar Min(CeedScalar a, CeedScalar b) { return a < b ? a : b; }
14 
15 CEED_QFUNCTION_HELPER void SwapScalar(CeedScalar *a, CeedScalar *b) {
16   CeedScalar temp = *a;
17   *a              = *b;
18   *b              = temp;
19 }
20 
21 CEED_QFUNCTION_HELPER CeedScalar Square(CeedScalar x) { return x * x; }
22 CEED_QFUNCTION_HELPER CeedScalar Cube(CeedScalar x) { return x * x * x; }
23 
24 // @brief Scale vector of length N by scalar alpha
25 CEED_QFUNCTION_HELPER void ScaleN(CeedScalar *u, const CeedScalar alpha, const CeedInt N) {
26   CeedPragmaSIMD for (CeedInt i = 0; i < N; i++) u[i] *= alpha;
27 }
28 
29 // @brief Set vector of length N to a value alpha
30 CEED_QFUNCTION_HELPER void SetValueN(CeedScalar *u, const CeedScalar alpha, const CeedInt N) {
31   CeedPragmaSIMD for (CeedInt i = 0; i < N; i++) u[i] = alpha;
32 }
33 
34 // @brief Copy N elements from x to y
35 CEED_QFUNCTION_HELPER void CopyN(const CeedScalar *x, CeedScalar *y, const CeedInt N) { CeedPragmaSIMD for (CeedInt i = 0; i < N; i++) y[i] = x[i]; }
36 
37 // @brief Copy 3x3 matrix from A to B
38 CEED_QFUNCTION_HELPER void CopyMat3(const CeedScalar A[3][3], CeedScalar B[3][3]) { CopyN((const CeedScalar *)A, (CeedScalar *)B, 9); }
39 
40 // @brief Dot product of vectors with N elements
41 CEED_QFUNCTION_HELPER CeedScalar DotN(const CeedScalar *u, const CeedScalar *v, const CeedInt N) {
42   CeedScalar output = 0;
43   CeedPragmaSIMD for (CeedInt i = 0; i < N; i++) output += u[i] * v[i];
44   return output;
45 }
46 
47 // @brief Dot product of 3 element vectors
48 CEED_QFUNCTION_HELPER CeedScalar Dot3(const CeedScalar *u, const CeedScalar *v) { return u[0] * v[0] + u[1] * v[1] + u[2] * v[2]; }
49 
50 // @brief \ell^2 norm of 3 element vectors
51 CEED_QFUNCTION_HELPER CeedScalar Norm3(const CeedScalar *u) { return sqrt(u[0] * u[0] + u[1] * u[1] + u[2] * u[2]); }
52 
53 // @brief Cross product of vectors with 3 elements
54 CEED_QFUNCTION_HELPER void Cross3(const CeedScalar u[3], const CeedScalar v[3], CeedScalar w[3]) {
55   w[0] = (u[1] * v[2]) - (u[2] * v[1]);
56   w[1] = (u[2] * v[0]) - (u[0] * v[2]);
57   w[2] = (u[0] * v[1]) - (u[1] * v[0]);
58 }
59 
60 // @brief Curl of vector given its gradient
61 CEED_QFUNCTION_HELPER void Curl3(const CeedScalar gradient[3][3], CeedScalar v[3]) {
62   v[0] = gradient[2][1] - gradient[1][2];
63   v[1] = gradient[0][2] - gradient[2][0];
64   v[2] = gradient[1][0] - gradient[0][1];
65 }
66 
67 // @brief Matrix vector product, b = Ax + b. A is NxM, x is M, b is N
68 CEED_QFUNCTION_HELPER void MatVecNM(const CeedScalar *A, const CeedScalar *x, const CeedInt N, const CeedInt M, const CeedTransposeMode transpose_A,
69                                     CeedScalar *b) {
70   switch (transpose_A) {
71     case CEED_NOTRANSPOSE:
72       CeedPragmaSIMD for (CeedInt i = 0; i < N; i++) b[i] += DotN(&A[i * M], x, M);
73       break;
74     case CEED_TRANSPOSE:
75       CeedPragmaSIMD for (CeedInt i = 0; i < M; i++) { CeedPragmaSIMD for (CeedInt j = 0; j < N; j++) b[i] += A[j * M + i] * x[j]; }
76       break;
77   }
78 }
79 
80 // @brief 3x3 Matrix vector product  b = Ax + b.
81 CEED_QFUNCTION_HELPER void MatVec3(const CeedScalar A[3][3], const CeedScalar x[3], const CeedTransposeMode transpose_A, CeedScalar b[3]) {
82   MatVecNM((const CeedScalar *)A, (const CeedScalar *)x, 3, 3, transpose_A, (CeedScalar *)b);
83 }
84 
85 // @brief Matrix-Matrix product, B = DA + B, where D is diagonal.
86 // @details A is NxM, D is diagonal NxN, represented by a vector of length N, and B is NxM. Optionally, A may be transposed.
87 CEED_QFUNCTION_HELPER void MatDiagNM(const CeedScalar *A, const CeedScalar *D, const CeedInt N, const CeedInt M, const CeedTransposeMode transpose_A,
88                                      CeedScalar *B) {
89   switch (transpose_A) {
90     case CEED_NOTRANSPOSE:
91       CeedPragmaSIMD for (CeedInt i = 0; i < N; i++) { CeedPragmaSIMD for (CeedInt j = 0; j < M; j++) B[i * M + j] += D[i] * A[i * M + j]; }
92       break;
93     case CEED_TRANSPOSE:
94       CeedPragmaSIMD for (CeedInt i = 0; i < M; i++) { CeedPragmaSIMD for (CeedInt j = 0; j < N; j++) B[i * N + j] += D[i] * A[j * M + i]; }
95       break;
96   }
97 }
98 
99 // @brief 3x3 Matrix-Matrix product, B = DA + B, where D is diagonal.
100 // @details Optionally, A may be transposed.
101 CEED_QFUNCTION_HELPER void MatDiag3(const CeedScalar A[3][3], const CeedScalar D[3], const CeedTransposeMode transpose_A, CeedScalar B[3][3]) {
102   MatDiagNM((const CeedScalar *)A, (const CeedScalar *)D, 3, 3, transpose_A, (CeedScalar *)B);
103 }
104 // @brief NxN Matrix-Matrix product, C = AB + C
105 CEED_QFUNCTION_HELPER void MatMatN(const CeedScalar *A, const CeedScalar *B, const CeedInt N, const CeedTransposeMode transpose_A,
106                                    const CeedTransposeMode transpose_B, CeedScalar *C) {
107   switch (transpose_A) {
108     case CEED_NOTRANSPOSE:
109       switch (transpose_B) {
110         case CEED_NOTRANSPOSE:
111           CeedPragmaSIMD for (CeedInt i = 0; i < N; i++) {
112             CeedPragmaSIMD for (CeedInt j = 0; j < N; j++) {
113               CeedPragmaSIMD for (CeedInt k = 0; k < N; k++) C[i * N + j] += A[i * N + k] * B[k * N + j];
114             }
115           }
116           break;
117         case CEED_TRANSPOSE:
118           CeedPragmaSIMD for (CeedInt i = 0; i < N; i++) {
119             CeedPragmaSIMD for (CeedInt j = 0; j < N; j++) {
120               CeedPragmaSIMD for (CeedInt k = 0; k < N; k++) C[i * N + j] += A[i * N + k] * B[j * N + k];
121             }
122           }
123           break;
124       }
125       break;
126     case CEED_TRANSPOSE:
127       switch (transpose_B) {
128         case CEED_NOTRANSPOSE:
129           CeedPragmaSIMD for (CeedInt i = 0; i < N; i++) {
130             CeedPragmaSIMD for (CeedInt j = 0; j < N; j++) {
131               CeedPragmaSIMD for (CeedInt k = 0; k < N; k++) C[i * N + j] += A[k * N + i] * B[k * N + j];
132             }
133           }
134           break;
135         case CEED_TRANSPOSE:
136           CeedPragmaSIMD for (CeedInt i = 0; i < N; i++) {
137             CeedPragmaSIMD for (CeedInt j = 0; j < N; j++) {
138               CeedPragmaSIMD for (CeedInt k = 0; k < N; k++) C[i * N + j] += A[k * N + i] * B[j * N + k];
139             }
140           }
141           break;
142       }
143       break;
144   }
145 }
146 
147 // @brief 3x3 Matrix-Matrix product, C = AB + C
148 CEED_QFUNCTION_HELPER void MatMat3(const CeedScalar A[3][3], const CeedScalar B[3][3], const CeedTransposeMode transpose_A,
149                                    const CeedTransposeMode transpose_B, CeedScalar C[3][3]) {
150   MatMatN((const CeedScalar *)A, (const CeedScalar *)B, 3, transpose_A, transpose_B, (CeedScalar *)C);
151 }
152 
153 /**
154   @brief MxN Matrix-Matrix product, C = AB + C
155 
156   C is NxM, A is NxP, B is PxM
157 
158   @param[in]  mat_A Row-major matrix `A`
159   @param[in]  mat_B Row-major matrix `B`
160   @param[out] mat_C Row-major output matrix `C`
161   @param[in]  N     Number of rows of `C`
162   @param[in]  M     Number of columns of `C`
163   @param[in]  P     Number of columns of `A`/rows of `B`
164 **/
165 CEED_QFUNCTION_HELPER void MatMatNM(const CeedScalar *mat_A, const CeedScalar *mat_B, CeedScalar *mat_C, CeedInt N, CeedInt M, CeedInt P) {
166   for (CeedInt i = 0; i < N; i++) {
167     for (CeedInt j = 0; j < M; j++) {
168       for (CeedInt k = 0; k < P; k++) mat_C[i * M + j] += mat_A[i * P + k] * mat_B[k * M + j];
169     }
170   }
171 }
172 
173 // @brief Unpack Kelvin-Mandel notation symmetric tensor into full tensor
174 CEED_QFUNCTION_HELPER void KMUnpack(const CeedScalar v[6], CeedScalar A[3][3]) {
175   const CeedScalar weight = 1 / sqrt(2.);
176   A[0][0]                 = v[0];
177   A[1][1]                 = v[1];
178   A[2][2]                 = v[2];
179   A[2][1] = A[1][2] = weight * v[3];
180   A[2][0] = A[0][2] = weight * v[4];
181   A[1][0] = A[0][1] = weight * v[5];
182 }
183 
184 // @brief Pack full tensor into Kelvin-Mandel notation symmetric tensor
185 CEED_QFUNCTION_HELPER void KMPack(const CeedScalar A[3][3], CeedScalar v[6]) {
186   const CeedScalar weight = sqrt(2.);
187   v[0]                    = A[0][0];
188   v[1]                    = A[1][1];
189   v[2]                    = A[2][2];
190   v[3]                    = A[2][1] * weight;
191   v[4]                    = A[2][0] * weight;
192   v[5]                    = A[1][0] * weight;
193 }
194 
195 // @brief Calculate metric tensor from mapping, g_{ij} = xi_{k,i} xi_{k,j} = dXdx^T dXdx
196 CEED_QFUNCTION_HELPER void KMMetricTensor(const CeedScalar dXdx[3][3], CeedScalar km_g_ij[6]) {
197   CeedScalar g_ij[3][3] = {{0.}};
198   MatMat3(dXdx, dXdx, CEED_TRANSPOSE, CEED_NOTRANSPOSE, g_ij);
199   KMPack(g_ij, km_g_ij);
200 }
201 
202 // @brief Linear ramp evaluation
203 CEED_QFUNCTION_HELPER CeedScalar LinearRampCoefficient(CeedScalar amplitude, CeedScalar length, CeedScalar start, CeedScalar x) {
204   if (x < start) {
205     return amplitude;
206   } else if (x < start + length) {
207     return amplitude * ((x - start) * (-1 / length) + 1);
208   } else {
209     return 0;
210   }
211 }
212 
213 /**
214   @brief Pack stored values at quadrature point
215 
216   @param[in]   Q              Number of quadrature points
217   @param[in]   i              Current quadrature point
218   @param[in]   start          Starting index to store components
219   @param[in]   num_comp       Number of components to store
220   @param[in]   values_at_qpnt Local values for quadrature point i
221   @param[out]  stored         Stored values
222 
223   @return An error code: 0 - success, otherwise - failure
224 **/
225 CEED_QFUNCTION_HELPER int StoredValuesPack(CeedInt Q, CeedInt i, CeedInt start, CeedInt num_comp, const CeedScalar *values_at_qpnt,
226                                            CeedScalar *stored) {
227   for (CeedInt j = 0; j < num_comp; j++) stored[(start + j) * Q + i] = values_at_qpnt[j];
228 
229   return CEED_ERROR_SUCCESS;
230 }
231 
232 /**
233   @brief Unpack stored values at quadrature point
234 
235   @param[in]   Q              Number of quadrature points
236   @param[in]   i              Current quadrature point
237   @param[in]   start          Starting index to store components
238   @param[in]   num_comp       Number of components to store
239   @param[in]   stored         Stored values
240   @param[out]  values_at_qpnt Local values for quadrature point i
241 
242   @return An error code: 0 - success, otherwise - failure
243 **/
244 CEED_QFUNCTION_HELPER int StoredValuesUnpack(CeedInt Q, CeedInt i, CeedInt start, CeedInt num_comp, const CeedScalar *stored,
245                                              CeedScalar *values_at_qpnt) {
246   for (CeedInt j = 0; j < num_comp; j++) values_at_qpnt[j] = stored[(start + j) * Q + i];
247 
248   return CEED_ERROR_SUCCESS;
249 }
250 
251 /**
252   @brief Unpack N-D element q_data at quadrature point
253 
254   @param[in]   dim       Dimension of the element
255   @param[in]   Q         Number of quadrature points
256   @param[in]   i         Current quadrature point
257   @param[in]   q_data    Pointer to q_data (generated by `setupgeo.h:Setup`)
258   @param[out]  wdetJ     Quadrature weight times determinant of the mapping Jacobian, or `NULL`
259   @param[out]  dXdx      Inverse of the mapping Jacobian (shape [dim][dim]), or `NULL`
260 **/
261 CEED_QFUNCTION_HELPER void QdataUnpack_ND(CeedInt dim, CeedInt Q, CeedInt i, const CeedScalar *q_data, CeedScalar *wdetJ, CeedScalar *dXdx) {
262   switch (dim) {
263     case 2:
264       if (wdetJ) StoredValuesUnpack(Q, i, 0, 1, q_data, wdetJ);
265       if (dXdx) StoredValuesUnpack(Q, i, 1, 4, q_data, dXdx);
266       break;
267     case 3:
268       if (wdetJ) StoredValuesUnpack(Q, i, 0, 1, q_data, wdetJ);
269       if (dXdx) StoredValuesUnpack(Q, i, 1, 9, q_data, dXdx);
270       break;
271   }
272 }
273 
274 /**
275   @brief Unpack boundary element q_data for N-D problem at quadrature point
276 
277   @param[in]   Q         Number of quadrature points
278   @param[in]   i         Current quadrature point
279   @param[in]   q_data    Pointer to q_data (generated by `setupgeo.h:SetupBoundary`)
280   @param[out]  wdetJ     Quadrature weight times determinant of the mapping Jacobian, or `NULL`
281   @param[out]  dXdx      Inverse of the mapping Jacobian (shape [dim - 1][dim]), or `NULL`
282   @param[out]  normal    Components of the normal vector (shape [dim]), or `NULL`
283 **/
284 CEED_QFUNCTION_HELPER void QdataBoundaryUnpack_ND(CeedInt dim, CeedInt Q, CeedInt i, const CeedScalar *q_data, CeedScalar *wdetJ, CeedScalar *dXdx,
285                                                   CeedScalar *normal) {
286   switch (dim) {
287     case 2:
288       if (wdetJ) StoredValuesUnpack(Q, i, 0, 1, q_data, wdetJ);
289       if (normal) StoredValuesUnpack(Q, i, 1, 2, q_data, normal);
290       break;
291     case 3:
292       if (wdetJ) StoredValuesUnpack(Q, i, 0, 1, q_data, wdetJ);
293       if (normal) StoredValuesUnpack(Q, i, 1, 3, q_data, normal);
294       if (dXdx) StoredValuesUnpack(Q, i, 4, 6, q_data, (CeedScalar *)dXdx);
295       break;
296   }
297 }
298 
299 /**
300   @brief Unpack 3D element q_data at quadrature point
301 
302   @param[in]   Q         Number of quadrature points
303   @param[in]   i         Current quadrature point
304   @param[in]   q_data    Pointer to q_data (generated by `setupgeo.h:Setup`)
305   @param[out]  wdetJ     Quadrature weight times determinant of the mapping Jacobian
306   @param[out]  dXdx      Inverse of the mapping Jacobian (shape [3][3])
307 
308   @return An error code: 0 - success, otherwise - failure
309 **/
310 CEED_QFUNCTION_HELPER int QdataUnpack_3D(CeedInt Q, CeedInt i, const CeedScalar *q_data, CeedScalar *wdetJ, CeedScalar dXdx[3][3]) {
311   QdataUnpack_ND(3, Q, i, q_data, wdetJ, (CeedScalar *)dXdx);
312   return CEED_ERROR_SUCCESS;
313 }
314 
315 /**
316   @brief Unpack boundary element q_data for 3D problem at quadrature point
317 
318   @param[in]   Q         Number of quadrature points
319   @param[in]   i         Current quadrature point
320   @param[in]   q_data    Pointer to q_data (generated by `setupgeo.h:SetupBoundary`)
321   @param[out]  wdetJ     Quadrature weight times determinant of the mapping Jacobian, or `NULL`
322   @param[out]  dXdx      Inverse of the mapping Jacobian (shape [2][3]), or `NULL`
323   @param[out]  normal    Components of the normal vector (shape [3]), or `NULL`
324 
325   @return An error code: 0 - success, otherwise - failure
326 **/
327 CEED_QFUNCTION_HELPER int QdataBoundaryUnpack_3D(CeedInt Q, CeedInt i, const CeedScalar *q_data, CeedScalar *wdetJ, CeedScalar dXdx[2][3],
328                                                  CeedScalar normal[3]) {
329   QdataBoundaryUnpack_ND(3, Q, i, q_data, wdetJ, (CeedScalar *)dXdx, normal);
330   return CEED_ERROR_SUCCESS;
331 }
332 
333 /**
334   @brief Unpack boundary element q_data for 3D problem at quadrature point
335 
336   @param[in]   Q         Number of quadrature points
337   @param[in]   i         Current quadrature point
338   @param[in]   q_data    Pointer to q_data (generated by `setupgeo.h:SetupBoundary`)
339   @param[out]  wdetJ     Quadrature weight times determinant of the mapping Jacobian, or `NULL`
340   @param[out]  dXdx      Inverse of the mapping Jacobian (shape [2][3]), or `NULL`
341   @param[out]  normal    Components of the normal vector (shape [3]), or `NULL`
342 
343   @return An error code: 0 - success, otherwise - failure
344 **/
345 CEED_QFUNCTION_HELPER int QdataBoundaryGradientUnpack_3D(CeedInt Q, CeedInt i, const CeedScalar *q_data, CeedScalar *wdetJ, CeedScalar dXdx[2][3],
346                                                          CeedScalar normal[3]) {
347   if (wdetJ) StoredValuesUnpack(Q, i, 0, 1, q_data, wdetJ);
348   if (dXdx) StoredValuesUnpack(Q, i, 1, 9, q_data, (CeedScalar *)dXdx);
349   if (normal) StoredValuesUnpack(Q, i, 10, 3, q_data, normal);
350   return CEED_ERROR_SUCCESS;
351 }
352 
353 /**
354   @brief Unpack 2D element q_data at quadrature point
355 
356   @param[in]   Q         Number of quadrature points
357   @param[in]   i         Current quadrature point
358   @param[in]   q_data    Pointer to q_data (generated by `setupgeo.h:Setup`)
359   @param[out]  wdetJ     Quadrature weight times determinant of the mapping Jacobian
360   @param[out]  dXdx      Inverse of the mapping Jacobian (shape [2][2])
361 
362   @return An error code: 0 - success, otherwise - failure
363 **/
364 CEED_QFUNCTION_HELPER int QdataUnpack_2D(CeedInt Q, CeedInt i, const CeedScalar *q_data, CeedScalar *wdetJ, CeedScalar dXdx[2][2]) {
365   QdataUnpack_ND(2, Q, i, q_data, wdetJ, (CeedScalar *)dXdx);
366   return CEED_ERROR_SUCCESS;
367 }
368 
369 /**
370   @brief Unpack boundary element q_data for 2D problem at quadrature point
371 
372   @param[in]   Q         Number of quadrature points
373   @param[in]   i         Current quadrature point
374   @param[in]   q_data    Pointer to q_data (generated by `setupgeo.h:SetupBoundary2d`)
375   @param[out]  wdetJ     Quadrature weight times determinant of the mapping Jacobian, or `NULL`
376   @param[out]  normal    Components of the normal vector (shape [2]), or `NULL`
377 
378   @return An error code: 0 - success, otherwise - failure
379 **/
380 CEED_QFUNCTION_HELPER int QdataBoundaryUnpack_2D(CeedInt Q, CeedInt i, const CeedScalar *q_data, CeedScalar *wdetJ, CeedScalar normal[2]) {
381   QdataBoundaryUnpack_ND(3, Q, i, q_data, wdetJ, NULL, normal);
382   return CEED_ERROR_SUCCESS;
383 }
384 
385 /**
386   @brief Unpack `CEED_EVAL_GRAD` QF input into quadrature-point local array
387 
388   @param[in]  Q        Number of quadrature points
389   @param[in]  i        Current quadrature point
390   @param[in]  num_comp Number of components of the input
391   @param[in]  dim      Topological dimension of the element (ie. number of derivative terms per component)
392   @param[in]  grad     QF gradient input, shape `[dim][num_comp][Q]`
393   @param[out] local    Gradient array at quadrature point Q, shape `[num_comp][dim]`
394 **/
395 CEED_QFUNCTION_HELPER void GradUnpackN(CeedInt Q, CeedInt i, CeedInt num_comp, CeedInt dim, const CeedScalar *grad, CeedScalar *local) {
396   for (CeedInt d = 0; d < dim; d++) {
397     for (CeedInt c = 0; c < num_comp; c++) {
398       local[dim * c + d] = grad[(Q * num_comp) * d + Q * c + i];
399     }
400   }
401 }
402 
403 /**
404   @brief Unpack `CEED_EVAL_GRAD` QF input into quadrature-point local array for 3D elements
405 
406   @param[in]  Q        Number of quadrature points
407   @param[in]  i        Current quadrature point
408   @param[in]  num_comp Number of components of the input
409   @param[in]  grad     QF gradient input, shape `[dim][num_comp][Q]`
410   @param[out] local    Gradient array at quadrature point Q, shape `[num_comp][dim]`
411 **/
412 CEED_QFUNCTION_HELPER void GradUnpack3(CeedInt Q, CeedInt i, CeedInt num_comp, const CeedScalar *grad, CeedScalar (*local)[3]) {
413   GradUnpackN(Q, i, num_comp, 3, grad, (CeedScalar *)local);
414 }
415 
416 /**
417   @brief Calculate divergence from reference gradient
418 
419   Given gradient array G_{ij} and inverse element mapping X_{ij}, then the divergence is
420 
421   G_{ij} X{ji}
422 
423   @param[in]  grad_qn    Gradient array, orientation [vector component][gradient direction]
424   @param[in]  dXdx       Inverse of the mapping Jacobian (shape [dim][dim])
425   @param[in]  dim        Dimension of the problem
426   @param[out] divergence The divergence
427 **/
428 CEED_QFUNCTION_HELPER void DivergenceND(const CeedScalar *grad_qn, const CeedScalar *dXdx, const CeedInt dim, CeedScalar *divergence) {
429   for (CeedInt i = 0; i < dim; i++) {
430     for (CeedInt j = 0; j < dim; j++) {
431       *divergence += grad_qn[i * dim + j] * dXdx[j * dim + i];
432     }
433   }
434 }
435 
436 /**
437   @brief Calculate divergence from reference gradient for 3D problem
438 
439   Given gradient array G_{ij} and inverse element mapping X_{ij}, then the divergence is
440 
441   G_{ij} X{ji}
442 
443   @param[in]  grad_qn    Gradient array, orientation [vector component][gradient direction]
444   @param[in]  dXdx       Inverse of the mapping Jacobian (shape [3][3])
445   @param[out] divergence The divergence
446 **/
447 CEED_QFUNCTION_HELPER void Divergence3D(const CeedScalar grad_qn[3][3], const CeedScalar dXdx[3][3], CeedScalar *divergence) {
448   DivergenceND((const CeedScalar *)grad_qn, (const CeedScalar *)dXdx, 3, divergence);
449 }
450