xref: /libCEED/rust/libceed-sys/c-src/interface/ceed-basis.c (revision ca62d558e488f6cc806d30b1949b4c18ff816b3a)
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 
8 #include <ceed-impl.h>
9 #include <ceed.h>
10 #include <ceed/backend.h>
11 #include <math.h>
12 #include <stdbool.h>
13 #include <stdio.h>
14 #include <string.h>
15 
16 /// @file
17 /// Implementation of CeedBasis interfaces
18 
19 /// @cond DOXYGEN_SKIP
20 static struct CeedBasis_private ceed_basis_none;
21 /// @endcond
22 
23 /// @addtogroup CeedBasisUser
24 /// @{
25 
26 /// Argument for @ref CeedOperatorSetField() indicating that the field does not require a `CeedBasis`
27 const CeedBasis CEED_BASIS_NONE = &ceed_basis_none;
28 
29 /// @}
30 
31 /// ----------------------------------------------------------------------------
32 /// CeedBasis Library Internal Functions
33 /// ----------------------------------------------------------------------------
34 /// @addtogroup CeedBasisDeveloper
35 /// @{
36 
37 /**
38   @brief Compute Chebyshev polynomial values at a point
39 
40   @param[in]  x           Coordinate to evaluate Chebyshev polynomials at
41   @param[in]  n           Number of Chebyshev polynomials to evaluate, `n >= 2`
42   @param[out] chebyshev_x Array of Chebyshev polynomial values
43 
44   @return An error code: 0 - success, otherwise - failure
45 
46   @ref Developer
47 **/
48 static int CeedChebyshevPolynomialsAtPoint(CeedScalar x, CeedInt n, CeedScalar *chebyshev_x) {
49   chebyshev_x[0] = 1.0;
50   chebyshev_x[1] = 2 * x;
51   for (CeedInt i = 2; i < n; i++) chebyshev_x[i] = 2 * x * chebyshev_x[i - 1] - chebyshev_x[i - 2];
52   return CEED_ERROR_SUCCESS;
53 }
54 
55 /**
56   @brief Compute values of the derivative of Chebyshev polynomials at a point
57 
58   @param[in]  x            Coordinate to evaluate derivative of Chebyshev polynomials at
59   @param[in]  n            Number of Chebyshev polynomials to evaluate, `n >= 2`
60   @param[out] chebyshev_dx Array of Chebyshev polynomial derivative values
61 
62   @return An error code: 0 - success, otherwise - failure
63 
64   @ref Developer
65 **/
66 static int CeedChebyshevDerivativeAtPoint(CeedScalar x, CeedInt n, CeedScalar *chebyshev_dx) {
67   CeedScalar chebyshev_x[3];
68 
69   chebyshev_x[1]  = 1.0;
70   chebyshev_x[2]  = 2 * x;
71   chebyshev_dx[0] = 0.0;
72   chebyshev_dx[1] = 2.0;
73   for (CeedInt i = 2; i < n; i++) {
74     chebyshev_x[0]  = chebyshev_x[1];
75     chebyshev_x[1]  = chebyshev_x[2];
76     chebyshev_x[2]  = 2 * x * chebyshev_x[1] - chebyshev_x[0];
77     chebyshev_dx[i] = 2 * x * chebyshev_dx[i - 1] + 2 * chebyshev_x[1] - chebyshev_dx[i - 2];
78   }
79   return CEED_ERROR_SUCCESS;
80 }
81 
82 /**
83   @brief Compute Householder reflection.
84 
85   Computes \f$A = (I - b v v^T) A\f$, where \f$A\f$ is an \f$m \times n\f$ matrix indexed as `A[i*row + j*col]`.
86 
87   @param[in,out] A   Matrix to apply Householder reflection to, in place
88   @param[in]     v   Householder vector
89   @param[in]     b   Scaling factor
90   @param[in]     m   Number of rows in `A`
91   @param[in]     n   Number of columns in `A`
92   @param[in]     row Row stride
93   @param[in]     col Col stride
94 
95   @return An error code: 0 - success, otherwise - failure
96 
97   @ref Developer
98 **/
99 static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v, CeedScalar b, CeedInt m, CeedInt n, CeedInt row, CeedInt col) {
100   for (CeedInt j = 0; j < n; j++) {
101     CeedScalar w = A[0 * row + j * col];
102 
103     for (CeedInt i = 1; i < m; i++) w += v[i] * A[i * row + j * col];
104     A[0 * row + j * col] -= b * w;
105     for (CeedInt i = 1; i < m; i++) A[i * row + j * col] -= b * w * v[i];
106   }
107   return CEED_ERROR_SUCCESS;
108 }
109 
110 /**
111   @brief Compute Givens rotation
112 
113   Computes \f$A = G A\f$ (or \f$G^T A\f$ in transpose mode), where \f$A\f$ is an \f$m \times n\f$ matrix indexed as `A[i*n + j*m]`.
114 
115   @param[in,out] A      Row major matrix to apply Givens rotation to, in place
116   @param[in]     c      Cosine factor
117   @param[in]     s      Sine factor
118   @param[in]     t_mode @ref CEED_NOTRANSPOSE to rotate the basis counter-clockwise, which has the effect of rotating columns of `A` clockwise;
119                           @ref CEED_TRANSPOSE for the opposite rotation
120   @param[in]     i      First row/column to apply rotation
121   @param[in]     k      Second row/column to apply rotation
122   @param[in]     m      Number of rows in `A`
123   @param[in]     n      Number of columns in `A`
124 
125   @return An error code: 0 - success, otherwise - failure
126 
127   @ref Developer
128 **/
129 static int CeedGivensRotation(CeedScalar *A, CeedScalar c, CeedScalar s, CeedTransposeMode t_mode, CeedInt i, CeedInt k, CeedInt m, CeedInt n) {
130   CeedInt stride_j = 1, stride_ik = m, num_its = n;
131 
132   if (t_mode == CEED_NOTRANSPOSE) {
133     stride_j  = n;
134     stride_ik = 1;
135     num_its   = m;
136   }
137 
138   // Apply rotation
139   for (CeedInt j = 0; j < num_its; j++) {
140     CeedScalar tau1 = A[i * stride_ik + j * stride_j], tau2 = A[k * stride_ik + j * stride_j];
141 
142     A[i * stride_ik + j * stride_j] = c * tau1 - s * tau2;
143     A[k * stride_ik + j * stride_j] = s * tau1 + c * tau2;
144   }
145   return CEED_ERROR_SUCCESS;
146 }
147 
148 /**
149   @brief View an array stored in a `CeedBasis`
150 
151   @param[in] name   Name of array
152   @param[in] fp_fmt Printing format
153   @param[in] m      Number of rows in array
154   @param[in] n      Number of columns in array
155   @param[in] a      Array to be viewed
156   @param[in] stream Stream to view to, e.g., `stdout`
157 
158   @return An error code: 0 - success, otherwise - failure
159 
160   @ref Developer
161 **/
162 static int CeedScalarView(const char *name, const char *fp_fmt, CeedInt m, CeedInt n, const CeedScalar *a, FILE *stream) {
163   if (m > 1) {
164     fprintf(stream, "  %s:\n", name);
165   } else {
166     char padded_name[12];
167 
168     snprintf(padded_name, 11, "%s:", name);
169     fprintf(stream, "  %-10s", padded_name);
170   }
171   for (CeedInt i = 0; i < m; i++) {
172     if (m > 1) fprintf(stream, "    [%" CeedInt_FMT "]", i);
173     for (CeedInt j = 0; j < n; j++) fprintf(stream, fp_fmt, fabs(a[i * n + j]) > 1E-14 ? a[i * n + j] : 0);
174     fputs("\n", stream);
175   }
176   return CEED_ERROR_SUCCESS;
177 }
178 
179 /**
180   @brief Create the interpolation and gradient matrices for projection from the nodes of `basis_from` to the nodes of `basis_to`.
181 
182   The interpolation is given by `interp_project = interp_to^+ * interp_from`, where the pseudoinverse `interp_to^+` is given by QR factorization.
183   The gradient is given by `grad_project = interp_to^+ * grad_from`, and is only computed for \f$H^1\f$ spaces otherwise it should not be used.
184 
185   Note: `basis_from` and `basis_to` must have compatible quadrature spaces.
186 
187   @param[in]  basis_from     `CeedBasis` to project from
188   @param[in]  basis_to       `CeedBasis` to project to
189   @param[out] interp_project Address of the variable where the newly created interpolation matrix will be stored
190   @param[out] grad_project   Address of the variable where the newly created gradient matrix will be stored
191 
192   @return An error code: 0 - success, otherwise - failure
193 
194   @ref Developer
195 **/
196 static int CeedBasisCreateProjectionMatrices(CeedBasis basis_from, CeedBasis basis_to, CeedScalar **interp_project, CeedScalar **grad_project) {
197   bool    are_both_tensor;
198   CeedInt Q, Q_to, Q_from, P_to, P_from;
199 
200   // Check for compatible quadrature spaces
201   CeedCall(CeedBasisGetNumQuadraturePoints(basis_to, &Q_to));
202   CeedCall(CeedBasisGetNumQuadraturePoints(basis_from, &Q_from));
203   CeedCheck(Q_to == Q_from, CeedBasisReturnCeed(basis_to), CEED_ERROR_DIMENSION,
204             "Bases must have compatible quadrature spaces."
205             " 'basis_from' has %" CeedInt_FMT " points and 'basis_to' has %" CeedInt_FMT,
206             Q_from, Q_to);
207   Q = Q_to;
208 
209   // Check for matching tensor or non-tensor
210   {
211     bool is_tensor_to, is_tensor_from;
212 
213     CeedCall(CeedBasisIsTensor(basis_to, &is_tensor_to));
214     CeedCall(CeedBasisIsTensor(basis_from, &is_tensor_from));
215     are_both_tensor = is_tensor_to && is_tensor_from;
216   }
217   if (are_both_tensor) {
218     CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_to));
219     CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_from));
220     CeedCall(CeedBasisGetNumQuadraturePoints1D(basis_from, &Q));
221   } else {
222     CeedCall(CeedBasisGetNumNodes(basis_to, &P_to));
223     CeedCall(CeedBasisGetNumNodes(basis_from, &P_from));
224   }
225 
226   // Check for matching FE space
227   CeedFESpace fe_space_to, fe_space_from;
228 
229   CeedCall(CeedBasisGetFESpace(basis_to, &fe_space_to));
230   CeedCall(CeedBasisGetFESpace(basis_from, &fe_space_from));
231   CeedCheck(fe_space_to == fe_space_from, CeedBasisReturnCeed(basis_to), CEED_ERROR_MINOR,
232             "Bases must both be the same FE space type."
233             " 'basis_from' is a %s and 'basis_to' is a %s",
234             CeedFESpaces[fe_space_from], CeedFESpaces[fe_space_to]);
235 
236   // Get source matrices
237   CeedInt           dim, q_comp = 1;
238   CeedScalar       *interp_to_inv, *interp_from;
239   const CeedScalar *interp_to_source = NULL, *interp_from_source = NULL, *grad_from_source = NULL;
240 
241   CeedCall(CeedBasisGetDimension(basis_from, &dim));
242   if (are_both_tensor) {
243     CeedCall(CeedBasisGetInterp1D(basis_to, &interp_to_source));
244     CeedCall(CeedBasisGetInterp1D(basis_from, &interp_from_source));
245   } else {
246     CeedCall(CeedBasisGetNumQuadratureComponents(basis_from, CEED_EVAL_INTERP, &q_comp));
247     CeedCall(CeedBasisGetInterp(basis_to, &interp_to_source));
248     CeedCall(CeedBasisGetInterp(basis_from, &interp_from_source));
249   }
250   CeedCall(CeedMalloc(Q * P_from * q_comp, &interp_from));
251   CeedCall(CeedCalloc(P_to * P_from, interp_project));
252 
253   // `grad_project = interp_to^+ * grad_from` is computed for the H^1 space case so the
254   // projection basis will have a gradient operation (allocated even if not H^1 for the
255   // basis construction later on)
256   if (fe_space_to == CEED_FE_SPACE_H1) {
257     if (are_both_tensor) {
258       CeedCall(CeedBasisGetGrad1D(basis_from, &grad_from_source));
259     } else {
260       CeedCall(CeedBasisGetGrad(basis_from, &grad_from_source));
261     }
262   }
263   CeedCall(CeedCalloc(P_to * P_from * (are_both_tensor ? 1 : dim), grad_project));
264 
265   // Compute interp_to^+, pseudoinverse of interp_to
266   CeedCall(CeedCalloc(Q * q_comp * P_to, &interp_to_inv));
267   CeedCall(CeedMatrixPseudoinverse(CeedBasisReturnCeed(basis_to), interp_to_source, Q * q_comp, P_to, interp_to_inv));
268   // Build matrices
269   CeedInt     num_matrices = 1 + (fe_space_to == CEED_FE_SPACE_H1) * (are_both_tensor ? 1 : dim);
270   CeedScalar *input_from[num_matrices], *output_project[num_matrices];
271 
272   input_from[0]     = (CeedScalar *)interp_from_source;
273   output_project[0] = *interp_project;
274   for (CeedInt m = 1; m < num_matrices; m++) {
275     input_from[m]     = (CeedScalar *)&grad_from_source[(m - 1) * Q * P_from];
276     output_project[m] = &((*grad_project)[(m - 1) * P_to * P_from]);
277   }
278   for (CeedInt m = 0; m < num_matrices; m++) {
279     // output_project = interp_to^+ * interp_from
280     memcpy(interp_from, input_from[m], Q * P_from * q_comp * sizeof(input_from[m][0]));
281     CeedCall(CeedMatrixMatrixMultiply(CeedBasisReturnCeed(basis_to), interp_to_inv, input_from[m], output_project[m], P_to, P_from, Q * q_comp));
282     // Round zero to machine precision
283     for (CeedInt i = 0; i < P_to * P_from; i++) {
284       if (fabs(output_project[m][i]) < 10 * CEED_EPSILON) output_project[m][i] = 0.0;
285     }
286   }
287 
288   // Cleanup
289   CeedCall(CeedFree(&interp_to_inv));
290   CeedCall(CeedFree(&interp_from));
291   return CEED_ERROR_SUCCESS;
292 }
293 
294 /**
295   @brief Check input vector dimensions for CeedBasisApply[Add]
296 
297   @param[in]  basis     `CeedBasis` to evaluate
298   @param[in]  num_elem  The number of elements to apply the basis evaluation to;
299                           the backend will specify the ordering in @ref CeedElemRestrictionCreate()
300   @param[in]  t_mode    @ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature points;
301                           @ref CEED_TRANSPOSE to apply the transpose, mapping from quadrature points to nodes
302   @param[in]  eval_mode @ref CEED_EVAL_NONE to use values directly,
303                           @ref CEED_EVAL_INTERP to use interpolated values,
304                           @ref CEED_EVAL_GRAD to use gradients,
305                           @ref CEED_EVAL_DIV to use divergence,
306                           @ref CEED_EVAL_CURL to use curl,
307                           @ref CEED_EVAL_WEIGHT to use quadrature weights
308   @param[in]  u         Input `CeedVector`
309   @param[out] v         Output `CeedVector`
310 
311   @return An error code: 0 - success, otherwise - failure
312 
313   @ref Developer
314 **/
315 static int CeedBasisApplyCheckDims(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
316   CeedInt  dim, num_comp, q_comp, num_nodes, num_qpts;
317   CeedSize u_length = 0, v_length;
318 
319   CeedCall(CeedBasisGetDimension(basis, &dim));
320   CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
321   CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp));
322   CeedCall(CeedBasisGetNumNodes(basis, &num_nodes));
323   CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts));
324   CeedCall(CeedVectorGetLength(v, &v_length));
325   if (u) CeedCall(CeedVectorGetLength(u, &u_length));
326 
327   // Check vector lengths to prevent out of bounds issues
328   bool has_good_dims = true;
329   switch (eval_mode) {
330     case CEED_EVAL_NONE:
331     case CEED_EVAL_INTERP:
332     case CEED_EVAL_GRAD:
333     case CEED_EVAL_DIV:
334     case CEED_EVAL_CURL:
335       has_good_dims = ((t_mode == CEED_TRANSPOSE && u_length >= (CeedSize)num_elem * (CeedSize)num_comp * (CeedSize)num_qpts * (CeedSize)q_comp &&
336                         v_length >= (CeedSize)num_elem * (CeedSize)num_comp * (CeedSize)num_nodes) ||
337                        (t_mode == CEED_NOTRANSPOSE && v_length >= (CeedSize)num_elem * (CeedSize)num_qpts * (CeedSize)num_comp * (CeedSize)q_comp &&
338                         u_length >= (CeedSize)num_elem * (CeedSize)num_comp * (CeedSize)num_nodes));
339       break;
340     case CEED_EVAL_WEIGHT:
341       has_good_dims = v_length >= (CeedSize)num_elem * (CeedSize)num_qpts;
342       break;
343   }
344   CeedCheck(has_good_dims, CeedBasisReturnCeed(basis), CEED_ERROR_DIMENSION, "Input/output vectors too short for basis and evaluation mode");
345   return CEED_ERROR_SUCCESS;
346 }
347 
348 /**
349   @brief Check input vector dimensions for CeedBasisApply[Add]AtPoints
350 
351   @param[in]  basis      `CeedBasis` to evaluate
352   @param[in]  num_elem   The number of elements to apply the basis evaluation to;
353                           the backend will specify the ordering in @ref CeedElemRestrictionCreate()
354   @param[in]  num_points Array of the number of points to apply the basis evaluation to in each element, size `num_elem`
355   @param[in]  t_mode     @ref CEED_NOTRANSPOSE to evaluate from nodes to points;
356                            @ref CEED_TRANSPOSE to apply the transpose, mapping from points to nodes
357   @param[in]  eval_mode  @ref CEED_EVAL_INTERP to use interpolated values,
358                            @ref CEED_EVAL_GRAD to use gradients,
359                            @ref CEED_EVAL_WEIGHT to use quadrature weights
360   @param[in]  x_ref      `CeedVector` holding reference coordinates of each point
361   @param[in]  u          Input `CeedVector`, of length `num_nodes * num_comp` for @ref CEED_NOTRANSPOSE
362   @param[out] v          Output `CeedVector`, of length `num_points * num_q_comp` for @ref CEED_NOTRANSPOSE with @ref CEED_EVAL_INTERP
363 
364   @return An error code: 0 - success, otherwise - failure
365 
366   @ref Developer
367 **/
368 static int CeedBasisApplyAtPointsCheckDims(CeedBasis basis, CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode,
369                                            CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) {
370   CeedInt  dim, num_comp, num_q_comp, num_nodes, P_1d = 1, Q_1d = 1, total_num_points = 0;
371   CeedSize x_length = 0, u_length = 0, v_length;
372 
373   CeedCall(CeedBasisGetDimension(basis, &dim));
374   CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d));
375   CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
376   CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
377   CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &num_q_comp));
378   CeedCall(CeedBasisGetNumNodes(basis, &num_nodes));
379   CeedCall(CeedVectorGetLength(v, &v_length));
380   if (x_ref != CEED_VECTOR_NONE) CeedCall(CeedVectorGetLength(x_ref, &x_length));
381   if (u != CEED_VECTOR_NONE) CeedCall(CeedVectorGetLength(u, &u_length));
382 
383   // Check compatibility coordinates vector
384   for (CeedInt i = 0; i < num_elem; i++) total_num_points += num_points[i];
385   CeedCheck((x_length >= (CeedSize)total_num_points * (CeedSize)dim) || (eval_mode == CEED_EVAL_WEIGHT), CeedBasisReturnCeed(basis),
386             CEED_ERROR_DIMENSION,
387             "Length of reference coordinate vector incompatible with basis dimension and number of points."
388             " Found reference coordinate vector of length %" CeedSize_FMT ", not of length %" CeedSize_FMT ".",
389             x_length, (CeedSize)total_num_points * (CeedSize)dim);
390 
391   // Check CEED_EVAL_WEIGHT only on CEED_NOTRANSPOSE
392   CeedCheck(eval_mode != CEED_EVAL_WEIGHT || t_mode == CEED_NOTRANSPOSE, CeedBasisReturnCeed(basis), CEED_ERROR_UNSUPPORTED,
393             "CEED_EVAL_WEIGHT only supported with CEED_NOTRANSPOSE");
394 
395   // Check vector lengths to prevent out of bounds issues
396   bool has_good_dims = true;
397   switch (eval_mode) {
398     case CEED_EVAL_INTERP:
399       has_good_dims = ((t_mode == CEED_TRANSPOSE && (u_length >= (CeedSize)total_num_points * (CeedSize)num_q_comp ||
400                                                      v_length >= (CeedSize)num_elem * (CeedSize)num_nodes * (CeedSize)num_comp)) ||
401                        (t_mode == CEED_NOTRANSPOSE && (v_length >= (CeedSize)total_num_points * (CeedSize)num_q_comp ||
402                                                        u_length >= (CeedSize)num_elem * (CeedSize)num_nodes * (CeedSize)num_comp)));
403       break;
404     case CEED_EVAL_GRAD:
405       has_good_dims = ((t_mode == CEED_TRANSPOSE && (u_length >= (CeedSize)total_num_points * (CeedSize)num_q_comp * (CeedSize)dim ||
406                                                      v_length >= (CeedSize)num_elem * (CeedSize)num_nodes * (CeedSize)num_comp)) ||
407                        (t_mode == CEED_NOTRANSPOSE && (v_length >= (CeedSize)total_num_points * (CeedSize)num_q_comp * (CeedSize)dim ||
408                                                        u_length >= (CeedSize)num_elem * (CeedSize)num_nodes * (CeedSize)num_comp)));
409       break;
410     case CEED_EVAL_WEIGHT:
411       has_good_dims = t_mode == CEED_NOTRANSPOSE && (v_length >= total_num_points);
412       break;
413       // LCOV_EXCL_START
414     case CEED_EVAL_NONE:
415     case CEED_EVAL_DIV:
416     case CEED_EVAL_CURL:
417       return CeedError(CeedBasisReturnCeed(basis), CEED_ERROR_UNSUPPORTED, "Evaluation at arbitrary points not supported for %s",
418                        CeedEvalModes[eval_mode]);
419       // LCOV_EXCL_STOP
420   }
421   CeedCheck(has_good_dims, CeedBasisReturnCeed(basis), CEED_ERROR_DIMENSION, "Input/output vectors too short for basis and evaluation mode");
422   return CEED_ERROR_SUCCESS;
423 }
424 
425 /**
426   @brief Default implimentation to apply basis evaluation from nodes to arbitrary points
427 
428   @param[in]  basis      `CeedBasis` to evaluate
429   @param[in]  apply_add  Sum result into target vector or overwrite
430   @param[in]  num_elem   The number of elements to apply the basis evaluation to;
431                           the backend will specify the ordering in @ref CeedElemRestrictionCreate()
432   @param[in]  num_points Array of the number of points to apply the basis evaluation to in each element, size `num_elem`
433   @param[in]  t_mode     @ref CEED_NOTRANSPOSE to evaluate from nodes to points;
434                            @ref CEED_TRANSPOSE to apply the transpose, mapping from points to nodes
435   @param[in]  eval_mode  @ref CEED_EVAL_INTERP to use interpolated values,
436                            @ref CEED_EVAL_GRAD to use gradients,
437                            @ref CEED_EVAL_WEIGHT to use quadrature weights
438   @param[in]  x_ref      `CeedVector` holding reference coordinates of each point
439   @param[in]  u          Input `CeedVector`, of length `num_nodes * num_comp` for @ref CEED_NOTRANSPOSE
440   @param[out] v          Output `CeedVector`, of length `num_points * num_q_comp` for @ref CEED_NOTRANSPOSE with @ref CEED_EVAL_INTERP
441 
442   @return An error code: 0 - success, otherwise - failure
443 
444   @ref Developer
445 **/
446 static int CeedBasisApplyAtPoints_Core(CeedBasis basis, bool apply_add, CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode,
447                                        CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) {
448   CeedInt dim, num_comp, P_1d = 1, Q_1d = 1, total_num_points = num_points[0];
449 
450   CeedCall(CeedBasisGetDimension(basis, &dim));
451   // Inserting check because clang-tidy doesn't understand this cannot occur
452   CeedCheck(dim > 0, CeedBasisReturnCeed(basis), CEED_ERROR_UNSUPPORTED, "Malformed CeedBasis, dim > 0 is required");
453   CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d));
454   CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
455   CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
456 
457   // Default implementation
458   {
459     bool is_tensor_basis;
460 
461     CeedCall(CeedBasisIsTensor(basis, &is_tensor_basis));
462     CeedCheck(is_tensor_basis, CeedBasisReturnCeed(basis), CEED_ERROR_UNSUPPORTED,
463               "Evaluation at arbitrary points only supported for tensor product bases");
464   }
465   CeedCheck(num_elem == 1, CeedBasisReturnCeed(basis), CEED_ERROR_UNSUPPORTED,
466             "Evaluation at arbitrary  points only supported for a single element at a time");
467   if (eval_mode == CEED_EVAL_WEIGHT) {
468     CeedCall(CeedVectorSetValue(v, 1.0));
469     return CEED_ERROR_SUCCESS;
470   }
471   if (!basis->basis_chebyshev) {
472     // Build basis mapping from nodes to Chebyshev coefficients
473     CeedScalar       *chebyshev_interp_1d, *chebyshev_grad_1d, *chebyshev_q_weight_1d;
474     const CeedScalar *q_ref_1d;
475     Ceed              ceed;
476 
477     CeedCall(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
478     CeedCall(CeedCalloc(P_1d * Q_1d, &chebyshev_grad_1d));
479     CeedCall(CeedCalloc(Q_1d, &chebyshev_q_weight_1d));
480     CeedCall(CeedBasisGetQRef(basis, &q_ref_1d));
481     CeedCall(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d));
482 
483     CeedCall(CeedBasisGetCeed(basis, &ceed));
484     CeedCall(CeedVectorCreate(ceed, num_comp * CeedIntPow(Q_1d, dim), &basis->vec_chebyshev));
485     CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P_1d, Q_1d, chebyshev_interp_1d, chebyshev_grad_1d, q_ref_1d, chebyshev_q_weight_1d,
486                                      &basis->basis_chebyshev));
487 
488     // Cleanup
489     CeedCall(CeedFree(&chebyshev_interp_1d));
490     CeedCall(CeedFree(&chebyshev_grad_1d));
491     CeedCall(CeedFree(&chebyshev_q_weight_1d));
492     CeedCall(CeedDestroy(&ceed));
493   }
494 
495   // Create TensorContract object if needed, such as a basis from the GPU backends
496   if (!basis->contract) {
497     Ceed      ceed_ref;
498     CeedBasis basis_ref = NULL;
499 
500     CeedCall(CeedInit("/cpu/self", &ceed_ref));
501     // Only need matching tensor contraction dimensions, any type of basis will work
502     CeedCall(CeedBasisCreateTensorH1Lagrange(ceed_ref, dim, num_comp, P_1d, Q_1d, CEED_GAUSS, &basis_ref));
503     // Note - clang-tidy doesn't know basis_ref->contract must be valid here
504     CeedCheck(basis_ref && basis_ref->contract, CeedBasisReturnCeed(basis), CEED_ERROR_UNSUPPORTED,
505               "Reference CPU ceed failed to create a tensor contraction object");
506     CeedCall(CeedTensorContractReferenceCopy(basis_ref->contract, &basis->contract));
507     CeedCall(CeedBasisDestroy(&basis_ref));
508     CeedCall(CeedDestroy(&ceed_ref));
509   }
510 
511   // Basis evaluation
512   switch (t_mode) {
513     case CEED_NOTRANSPOSE: {
514       // Nodes to arbitrary points
515       CeedScalar       *v_array;
516       const CeedScalar *chebyshev_coeffs, *x_array_read;
517 
518       // -- Interpolate to Chebyshev coefficients
519       CeedCall(CeedBasisApply(basis->basis_chebyshev, 1, CEED_NOTRANSPOSE, CEED_EVAL_INTERP, u, basis->vec_chebyshev));
520 
521       // -- Evaluate Chebyshev polynomials at arbitrary points
522       CeedCall(CeedVectorGetArrayRead(basis->vec_chebyshev, CEED_MEM_HOST, &chebyshev_coeffs));
523       CeedCall(CeedVectorGetArrayRead(x_ref, CEED_MEM_HOST, &x_array_read));
524       CeedCall(CeedVectorGetArrayWrite(v, CEED_MEM_HOST, &v_array));
525       switch (eval_mode) {
526         case CEED_EVAL_INTERP: {
527           CeedScalar tmp[2][num_comp * CeedIntPow(Q_1d, dim)], chebyshev_x[Q_1d];
528 
529           // ---- Values at point
530           for (CeedInt p = 0; p < total_num_points; p++) {
531             CeedInt pre = num_comp * CeedIntPow(Q_1d, dim - 1), post = 1;
532 
533             for (CeedInt d = 0; d < dim; d++) {
534               // ------ Tensor contract with current Chebyshev polynomial values
535               CeedCall(CeedChebyshevPolynomialsAtPoint(x_array_read[d * total_num_points + p], Q_1d, chebyshev_x));
536               CeedCall(CeedTensorContractApply(basis->contract, pre, Q_1d, post, 1, chebyshev_x, t_mode, false,
537                                                d == 0 ? chebyshev_coeffs : tmp[d % 2], tmp[(d + 1) % 2]));
538               pre /= Q_1d;
539               post *= 1;
540             }
541             for (CeedInt c = 0; c < num_comp; c++) v_array[c * total_num_points + p] = tmp[dim % 2][c];
542           }
543           break;
544         }
545         case CEED_EVAL_GRAD: {
546           CeedScalar tmp[2][num_comp * CeedIntPow(Q_1d, dim)], chebyshev_x[Q_1d];
547 
548           // ---- Values at point
549           for (CeedInt p = 0; p < total_num_points; p++) {
550             // Dim**2 contractions, apply grad when pass == dim
551             for (CeedInt pass = 0; pass < dim; pass++) {
552               CeedInt pre = num_comp * CeedIntPow(Q_1d, dim - 1), post = 1;
553 
554               for (CeedInt d = 0; d < dim; d++) {
555                 // ------ Tensor contract with current Chebyshev polynomial values
556                 if (pass == d) CeedCall(CeedChebyshevDerivativeAtPoint(x_array_read[d * total_num_points + p], Q_1d, chebyshev_x));
557                 else CeedCall(CeedChebyshevPolynomialsAtPoint(x_array_read[d * total_num_points + p], Q_1d, chebyshev_x));
558                 CeedCall(CeedTensorContractApply(basis->contract, pre, Q_1d, post, 1, chebyshev_x, t_mode, false,
559                                                  d == 0 ? chebyshev_coeffs : tmp[d % 2], tmp[(d + 1) % 2]));
560                 pre /= Q_1d;
561                 post *= 1;
562               }
563               for (CeedInt c = 0; c < num_comp; c++) v_array[(pass * num_comp + c) * total_num_points + p] = tmp[dim % 2][c];
564             }
565           }
566           break;
567         }
568         default:
569           // Nothing to do, excluded above
570           break;
571       }
572       CeedCall(CeedVectorRestoreArrayRead(basis->vec_chebyshev, &chebyshev_coeffs));
573       CeedCall(CeedVectorRestoreArrayRead(x_ref, &x_array_read));
574       CeedCall(CeedVectorRestoreArray(v, &v_array));
575       break;
576     }
577     case CEED_TRANSPOSE: {
578       // Note: No switch on e_mode here because only CEED_EVAL_INTERP is supported at this time
579       // Arbitrary points to nodes
580       CeedScalar       *chebyshev_coeffs;
581       const CeedScalar *u_array, *x_array_read;
582 
583       // -- Transpose of evaluation of Chebyshev polynomials at arbitrary points
584       CeedCall(CeedVectorGetArrayWrite(basis->vec_chebyshev, CEED_MEM_HOST, &chebyshev_coeffs));
585       CeedCall(CeedVectorGetArrayRead(x_ref, CEED_MEM_HOST, &x_array_read));
586       CeedCall(CeedVectorGetArrayRead(u, CEED_MEM_HOST, &u_array));
587 
588       switch (eval_mode) {
589         case CEED_EVAL_INTERP: {
590           CeedScalar tmp[2][num_comp * CeedIntPow(Q_1d, dim)], chebyshev_x[Q_1d];
591 
592           // ---- Values at point
593           for (CeedInt p = 0; p < total_num_points; p++) {
594             CeedInt pre = num_comp * 1, post = 1;
595 
596             for (CeedInt c = 0; c < num_comp; c++) tmp[0][c] = u_array[c * total_num_points + p];
597             for (CeedInt d = 0; d < dim; d++) {
598               // ------ Tensor contract with current Chebyshev polynomial values
599               CeedCall(CeedChebyshevPolynomialsAtPoint(x_array_read[d * total_num_points + p], Q_1d, chebyshev_x));
600               CeedCall(CeedTensorContractApply(basis->contract, pre, 1, post, Q_1d, chebyshev_x, t_mode, p > 0 && d == (dim - 1), tmp[d % 2],
601                                                d == (dim - 1) ? chebyshev_coeffs : tmp[(d + 1) % 2]));
602               pre /= 1;
603               post *= Q_1d;
604             }
605           }
606           break;
607         }
608         case CEED_EVAL_GRAD: {
609           CeedScalar tmp[2][num_comp * CeedIntPow(Q_1d, dim)], chebyshev_x[Q_1d];
610 
611           // ---- Values at point
612           for (CeedInt p = 0; p < total_num_points; p++) {
613             // Dim**2 contractions, apply grad when pass == dim
614             for (CeedInt pass = 0; pass < dim; pass++) {
615               CeedInt pre = num_comp * 1, post = 1;
616 
617               for (CeedInt c = 0; c < num_comp; c++) tmp[0][c] = u_array[(pass * num_comp + c) * total_num_points + p];
618               for (CeedInt d = 0; d < dim; d++) {
619                 // ------ Tensor contract with current Chebyshev polynomial values
620                 if (pass == d) CeedCall(CeedChebyshevDerivativeAtPoint(x_array_read[d * total_num_points + p], Q_1d, chebyshev_x));
621                 else CeedCall(CeedChebyshevPolynomialsAtPoint(x_array_read[d * total_num_points + p], Q_1d, chebyshev_x));
622                 CeedCall(CeedTensorContractApply(basis->contract, pre, 1, post, Q_1d, chebyshev_x, t_mode,
623                                                  (p > 0 || (p == 0 && pass > 0)) && d == (dim - 1), tmp[d % 2],
624                                                  d == (dim - 1) ? chebyshev_coeffs : tmp[(d + 1) % 2]));
625                 pre /= 1;
626                 post *= Q_1d;
627               }
628             }
629           }
630           break;
631         }
632         default:
633           // Nothing to do, excluded above
634           break;
635       }
636       CeedCall(CeedVectorRestoreArray(basis->vec_chebyshev, &chebyshev_coeffs));
637       CeedCall(CeedVectorRestoreArrayRead(x_ref, &x_array_read));
638       CeedCall(CeedVectorRestoreArrayRead(u, &u_array));
639 
640       // -- Interpolate transpose from Chebyshev coefficients
641       if (apply_add) CeedCall(CeedBasisApplyAdd(basis->basis_chebyshev, 1, CEED_TRANSPOSE, CEED_EVAL_INTERP, basis->vec_chebyshev, v));
642       else CeedCall(CeedBasisApply(basis->basis_chebyshev, 1, CEED_TRANSPOSE, CEED_EVAL_INTERP, basis->vec_chebyshev, v));
643       break;
644     }
645   }
646   return CEED_ERROR_SUCCESS;
647 }
648 
649 /// @}
650 
651 /// ----------------------------------------------------------------------------
652 /// Ceed Backend API
653 /// ----------------------------------------------------------------------------
654 /// @addtogroup CeedBasisBackend
655 /// @{
656 
657 /**
658   @brief Fallback to a reference implementation for a non tensor-product basis for \f$H^1\f$ discretizations.
659     This function may only be called inside of a backend `BasisCreateH1` function.
660     This is used by a backend when the specific parameters for a `CeedBasis` exceed the backend's support, such as
661     when a `interp` and `grad` matrices require too many bytes to fit into shared memory on a GPU.
662 
663   @param[in]  ceed      `Ceed` object used to create the `CeedBasis`
664   @param[in]  topo      Topology of element, e.g. hypercube, simplex, etc
665   @param[in]  num_comp  Number of field components (1 for scalar fields)
666   @param[in]  num_nodes Total number of nodes
667   @param[in]  num_qpts  Total number of quadrature points
668   @param[in]  interp    Row-major (`num_qpts * num_nodes`) matrix expressing the values of nodal basis functions at quadrature points
669   @param[in]  grad      Row-major (`dim * num_qpts * num_nodes`) matrix expressing derivatives of nodal basis functions at quadrature points
670   @param[in]  q_ref     Array of length `num_qpts * dim` holding the locations of quadrature points on the reference element
671   @param[in]  q_weight  Array of length `num_qpts` holding the quadrature weights on the reference element
672   @param[out] basis     Newly created `CeedBasis`
673 
674   @return An error code: 0 - success, otherwise - failure
675 
676   @ref User
677 **/
678 int CeedBasisCreateH1Fallback(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
679                               const CeedScalar *grad, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis basis) {
680   CeedInt P = num_nodes, Q = num_qpts, dim = 0;
681   Ceed    delegate;
682 
683   CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
684   CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateH1");
685 
686   CeedCall(CeedReferenceCopy(delegate, &(basis)->ceed));
687   CeedCall(CeedBasisGetTopologyDimension(topo, &dim));
688   CeedCall(delegate->BasisCreateH1(topo, dim, P, Q, interp, grad, q_ref, q_weight, basis));
689   CeedCall(CeedDestroy(&delegate));
690   return CEED_ERROR_SUCCESS;
691 }
692 
693 /**
694   @brief Return collocated gradient matrix
695 
696   @param[in]  basis         `CeedBasis`
697   @param[out] collo_grad_1d Row-major (`Q_1d * Q_1d`) matrix expressing derivatives of basis functions at quadrature points
698 
699   @return An error code: 0 - success, otherwise - failure
700 
701   @ref Backend
702 **/
703 int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *collo_grad_1d) {
704   Ceed              ceed;
705   CeedInt           P_1d, Q_1d;
706   CeedScalar       *interp_1d_pinv;
707   const CeedScalar *grad_1d, *interp_1d;
708 
709   // Note: This function is for backend use, so all errors are terminal and we do not need to clean up memory on failure.
710   CeedCall(CeedBasisGetCeed(basis, &ceed));
711   CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d));
712   CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
713 
714   // Compute interp_1d^+, pseudoinverse of interp_1d
715   CeedCall(CeedCalloc(P_1d * Q_1d, &interp_1d_pinv));
716   CeedCall(CeedBasisGetInterp1D(basis, &interp_1d));
717   CeedCall(CeedMatrixPseudoinverse(ceed, interp_1d, Q_1d, P_1d, interp_1d_pinv));
718   CeedCall(CeedBasisGetGrad1D(basis, &grad_1d));
719   CeedCall(CeedMatrixMatrixMultiply(ceed, grad_1d, (const CeedScalar *)interp_1d_pinv, collo_grad_1d, Q_1d, Q_1d, P_1d));
720 
721   CeedCall(CeedFree(&interp_1d_pinv));
722   CeedCall(CeedDestroy(&ceed));
723   return CEED_ERROR_SUCCESS;
724 }
725 
726 /**
727   @brief Return 1D interpolation matrix to Chebyshev polynomial coefficients on quadrature space
728 
729   @param[in]  basis               `CeedBasis`
730   @param[out] chebyshev_interp_1d Row-major (`P_1d * Q_1d`) matrix interpolating from basis nodes to Chebyshev polynomial coefficients
731 
732   @return An error code: 0 - success, otherwise - failure
733 
734   @ref Backend
735 **/
736 int CeedBasisGetChebyshevInterp1D(CeedBasis basis, CeedScalar *chebyshev_interp_1d) {
737   CeedInt           P_1d, Q_1d;
738   CeedScalar       *C, *chebyshev_coeffs_1d_inv;
739   const CeedScalar *interp_1d, *q_ref_1d;
740   Ceed              ceed;
741 
742   CeedCall(CeedBasisGetCeed(basis, &ceed));
743   CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d));
744   CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
745 
746   // Build coefficient matrix
747   // -- Note: Clang-tidy needs this check
748   CeedCheck((P_1d > 0) && (Q_1d > 0), ceed, CEED_ERROR_INCOMPATIBLE, "CeedBasis dimensions are malformed");
749   CeedCall(CeedCalloc(Q_1d * Q_1d, &C));
750   CeedCall(CeedBasisGetQRef(basis, &q_ref_1d));
751   for (CeedInt i = 0; i < Q_1d; i++) CeedCall(CeedChebyshevPolynomialsAtPoint(q_ref_1d[i], Q_1d, &C[i * Q_1d]));
752 
753   // Compute C^+, pseudoinverse of coefficient matrix
754   CeedCall(CeedCalloc(Q_1d * Q_1d, &chebyshev_coeffs_1d_inv));
755   CeedCall(CeedMatrixPseudoinverse(ceed, C, Q_1d, Q_1d, chebyshev_coeffs_1d_inv));
756 
757   // Build mapping from nodes to Chebyshev coefficients
758   CeedCall(CeedBasisGetInterp1D(basis, &interp_1d));
759   CeedCall(CeedMatrixMatrixMultiply(ceed, chebyshev_coeffs_1d_inv, interp_1d, chebyshev_interp_1d, Q_1d, P_1d, Q_1d));
760 
761   // Cleanup
762   CeedCall(CeedFree(&C));
763   CeedCall(CeedFree(&chebyshev_coeffs_1d_inv));
764   CeedCall(CeedDestroy(&ceed));
765   return CEED_ERROR_SUCCESS;
766 }
767 
768 /**
769   @brief Get tensor status for given `CeedBasis`
770 
771   @param[in]  basis     `CeedBasis`
772   @param[out] is_tensor Variable to store tensor status
773 
774   @return An error code: 0 - success, otherwise - failure
775 
776   @ref Backend
777 **/
778 int CeedBasisIsTensor(CeedBasis basis, bool *is_tensor) {
779   *is_tensor = basis->is_tensor_basis;
780   return CEED_ERROR_SUCCESS;
781 }
782 
783 /**
784   @brief Determine if given `CeedBasis` has nodes collocated with quadrature points
785 
786   @param[in]  basis     `CeedBasis`
787   @param[out] is_tensor Variable to store collocated status
788 
789   @return An error code: 0 - success, otherwise - failure
790 
791   @ref Backend
792 **/
793 int CeedBasisIsCollocated(CeedBasis basis, bool *is_collocated) {
794   if (basis->is_tensor_basis && (basis->Q_1d == basis->P_1d)) {
795     *is_collocated = true;
796 
797     for (CeedInt i = 0; i < basis->P_1d; i++) {
798       *is_collocated = *is_collocated && (fabs(basis->interp_1d[i + basis->P_1d * i] - 1.0) < 10 * CEED_EPSILON);
799       for (CeedInt j = 0; j < basis->Q_1d; j++) {
800         if (j != i) *is_collocated = *is_collocated && (fabs(basis->interp_1d[j + basis->P_1d * i]) < 10 * CEED_EPSILON);
801       }
802     }
803   } else {
804     *is_collocated = false;
805   }
806   return CEED_ERROR_SUCCESS;
807 }
808 
809 /**
810   @brief Get backend data of a `CeedBasis`
811 
812   @param[in]  basis `CeedBasis`
813   @param[out] data  Variable to store data
814 
815   @return An error code: 0 - success, otherwise - failure
816 
817   @ref Backend
818 **/
819 int CeedBasisGetData(CeedBasis basis, void *data) {
820   *(void **)data = basis->data;
821   return CEED_ERROR_SUCCESS;
822 }
823 
824 /**
825   @brief Set backend data of a `CeedBasis`
826 
827   @param[in,out] basis  `CeedBasis`
828   @param[in]     data   Data to set
829 
830   @return An error code: 0 - success, otherwise - failure
831 
832   @ref Backend
833 **/
834 int CeedBasisSetData(CeedBasis basis, void *data) {
835   basis->data = data;
836   return CEED_ERROR_SUCCESS;
837 }
838 
839 /**
840   @brief Increment the reference counter for a `CeedBasis`
841 
842   @param[in,out] basis `CeedBasis` to increment the reference counter
843 
844   @return An error code: 0 - success, otherwise - failure
845 
846   @ref Backend
847 **/
848 int CeedBasisReference(CeedBasis basis) {
849   basis->ref_count++;
850   return CEED_ERROR_SUCCESS;
851 }
852 
853 /**
854   @brief Get number of Q-vector components for given `CeedBasis`
855 
856   @param[in]  basis     `CeedBasis`
857   @param[in]  eval_mode @ref CEED_EVAL_INTERP to use interpolated values,
858                           @ref CEED_EVAL_GRAD to use gradients,
859                           @ref CEED_EVAL_DIV to use divergence,
860                           @ref CEED_EVAL_CURL to use curl
861   @param[out] q_comp    Variable to store number of Q-vector components of basis
862 
863   @return An error code: 0 - success, otherwise - failure
864 
865   @ref Backend
866 **/
867 int CeedBasisGetNumQuadratureComponents(CeedBasis basis, CeedEvalMode eval_mode, CeedInt *q_comp) {
868   CeedInt dim;
869 
870   CeedCall(CeedBasisGetDimension(basis, &dim));
871   switch (eval_mode) {
872     case CEED_EVAL_INTERP: {
873       CeedFESpace fe_space;
874 
875       CeedCall(CeedBasisGetFESpace(basis, &fe_space));
876       *q_comp = (fe_space == CEED_FE_SPACE_H1) ? 1 : dim;
877     } break;
878     case CEED_EVAL_GRAD:
879       *q_comp = dim;
880       break;
881     case CEED_EVAL_DIV:
882       *q_comp = 1;
883       break;
884     case CEED_EVAL_CURL:
885       *q_comp = (dim < 3) ? 1 : dim;
886       break;
887     case CEED_EVAL_NONE:
888     case CEED_EVAL_WEIGHT:
889       *q_comp = 1;
890       break;
891   }
892   return CEED_ERROR_SUCCESS;
893 }
894 
895 /**
896   @brief Estimate number of FLOPs required to apply `CeedBasis` in `t_mode` and `eval_mode`
897 
898   @param[in]  basis        `CeedBasis` to estimate FLOPs for
899   @param[in]  t_mode       Apply basis or transpose
900   @param[in]  eval_mode    @ref CeedEvalMode
901   @param[in]  is_at_points Evaluate the basis at points or quadrature points
902   @param[in]  num_points   Number of points basis is evaluated at
903   @param[out] flops        Address of variable to hold FLOPs estimate
904 
905   @ref Backend
906 **/
907 int CeedBasisGetFlopsEstimate(CeedBasis basis, CeedTransposeMode t_mode, CeedEvalMode eval_mode, bool is_at_points, CeedInt num_points,
908                               CeedSize *flops) {
909   bool is_tensor;
910 
911   CeedCall(CeedBasisIsTensor(basis, &is_tensor));
912   CeedCheck(!is_at_points || is_tensor, CeedBasisReturnCeed(basis), CEED_ERROR_INCOMPATIBLE, "Can only evaluate tensor-product bases at points");
913   if (is_tensor) {
914     CeedInt dim, num_comp, P_1d, Q_1d;
915 
916     CeedCall(CeedBasisGetDimension(basis, &dim));
917     CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
918     CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d));
919     CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
920     if (t_mode == CEED_TRANSPOSE) {
921       P_1d = Q_1d;
922       Q_1d = P_1d;
923     }
924     CeedInt tensor_flops = 0, pre = num_comp * CeedIntPow(P_1d, dim - 1), post = 1;
925 
926     for (CeedInt d = 0; d < dim; d++) {
927       tensor_flops += 2 * pre * P_1d * post * Q_1d;
928       pre /= P_1d;
929       post *= Q_1d;
930     }
931     if (is_at_points) {
932       bool is_gpu = false;
933 
934       {
935         CeedMemType mem_type;
936 
937         CeedCall(CeedGetPreferredMemType(CeedBasisReturnCeed(basis), &mem_type));
938         is_gpu = mem_type == CEED_MEM_DEVICE;
939       }
940 
941       CeedInt chebyshev_flops = (Q_1d - 2) * 3 + 1, d_chebyshev_flops = (Q_1d - 2) * 8 + 1;
942       CeedInt point_tensor_flops = 0, pre = CeedIntPow(Q_1d, dim - 1), post = 1;
943 
944       for (CeedInt d = 0; d < dim; d++) {
945         point_tensor_flops += 2 * pre * Q_1d * post * 1;
946         pre /= P_1d;
947         post *= Q_1d;
948       }
949 
950       switch (eval_mode) {
951         case CEED_EVAL_NONE:
952           *flops = 0;
953           break;
954         case CEED_EVAL_INTERP: {
955           *flops = tensor_flops + num_points * num_comp * (point_tensor_flops + (t_mode == CEED_TRANSPOSE ? CeedIntPow(Q_1d, dim) : 0));
956           if (dim == 3 && is_gpu) {
957             *flops += num_points * Q_1d *
958                       (chebyshev_flops + num_comp * (2 * chebyshev_flops + 2 * Q_1d * Q_1d + (t_mode == CEED_TRANSPOSE ? 2 * Q_1d + 1 : 3 * Q_1d)));
959           } else {
960             *flops += num_points * (is_gpu ? num_comp : 1) * dim * chebyshev_flops;
961           }
962           break;
963         }
964         case CEED_EVAL_GRAD: {
965           *flops = tensor_flops + num_points * num_comp * (point_tensor_flops + (t_mode == CEED_TRANSPOSE ? CeedIntPow(Q_1d, dim) : 0));
966           if (dim == 3 && is_gpu) {
967             CeedInt inner_flops =
968                 dim * (2 * Q_1d * Q_1d + (t_mode == CEED_TRANSPOSE ? 2 : 3) * Q_1d) + (dim - 1) * (2 * chebyshev_flops + d_chebyshev_flops);
969 
970             *flops += num_points * Q_1d * (chebyshev_flops + d_chebyshev_flops + num_comp * (inner_flops + (t_mode == CEED_TRANSPOSE ? 1 : 0)));
971           } else {
972             *flops += num_points * (is_gpu ? num_comp : 1) * dim * (d_chebyshev_flops + (dim - 1) * chebyshev_flops);
973           }
974           break;
975         }
976         case CEED_EVAL_DIV:
977         case CEED_EVAL_CURL: {
978           // LCOV_EXCL_START
979           return CeedError(CeedBasisReturnCeed(basis), CEED_ERROR_INCOMPATIBLE, "Tensor basis evaluation for %s not supported at points",
980                            CeedEvalModes[eval_mode]);
981           break;
982           // LCOV_EXCL_STOP
983         }
984         case CEED_EVAL_WEIGHT:
985           *flops = num_points;
986           break;
987       }
988     } else {
989       switch (eval_mode) {
990         case CEED_EVAL_NONE:
991           *flops = 0;
992           break;
993         case CEED_EVAL_INTERP:
994           *flops = tensor_flops;
995           break;
996         case CEED_EVAL_GRAD:
997           *flops = tensor_flops * 2;
998           break;
999         case CEED_EVAL_DIV:
1000         case CEED_EVAL_CURL: {
1001           // LCOV_EXCL_START
1002           return CeedError(CeedBasisReturnCeed(basis), CEED_ERROR_INCOMPATIBLE, "Tensor basis evaluation for %s not supported",
1003                            CeedEvalModes[eval_mode]);
1004           break;
1005           // LCOV_EXCL_STOP
1006         }
1007         case CEED_EVAL_WEIGHT:
1008           *flops = dim * CeedIntPow(Q_1d, dim);
1009           break;
1010       }
1011     }
1012   } else {
1013     CeedInt dim, num_comp, q_comp, num_nodes, num_qpts;
1014 
1015     CeedCall(CeedBasisGetDimension(basis, &dim));
1016     CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
1017     CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp));
1018     CeedCall(CeedBasisGetNumNodes(basis, &num_nodes));
1019     CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts));
1020     switch (eval_mode) {
1021       case CEED_EVAL_NONE:
1022         *flops = 0;
1023         break;
1024       case CEED_EVAL_INTERP:
1025       case CEED_EVAL_GRAD:
1026       case CEED_EVAL_DIV:
1027       case CEED_EVAL_CURL:
1028         *flops = num_nodes * num_qpts * num_comp * q_comp;
1029         break;
1030       case CEED_EVAL_WEIGHT:
1031         *flops = 0;
1032         break;
1033     }
1034   }
1035   return CEED_ERROR_SUCCESS;
1036 }
1037 
1038 /**
1039   @brief Get `CeedFESpace` for a `CeedBasis`
1040 
1041   @param[in]  basis    `CeedBasis`
1042   @param[out] fe_space Variable to store `CeedFESpace`
1043 
1044   @return An error code: 0 - success, otherwise - failure
1045 
1046   @ref Backend
1047 **/
1048 int CeedBasisGetFESpace(CeedBasis basis, CeedFESpace *fe_space) {
1049   *fe_space = basis->fe_space;
1050   return CEED_ERROR_SUCCESS;
1051 }
1052 
1053 /**
1054   @brief Get dimension for given `CeedElemTopology`
1055 
1056   @param[in]  topo `CeedElemTopology`
1057   @param[out] dim  Variable to store dimension of topology
1058 
1059   @return An error code: 0 - success, otherwise - failure
1060 
1061   @ref Backend
1062 **/
1063 int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) {
1064   *dim = (CeedInt)topo >> 16;
1065   return CEED_ERROR_SUCCESS;
1066 }
1067 
1068 /**
1069   @brief Get `CeedTensorContract` of a `CeedBasis`
1070 
1071   @param[in]  basis     `CeedBasis`
1072   @param[out] contract  Variable to store `CeedTensorContract`
1073 
1074   @return An error code: 0 - success, otherwise - failure
1075 
1076   @ref Backend
1077 **/
1078 int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) {
1079   *contract = basis->contract;
1080   return CEED_ERROR_SUCCESS;
1081 }
1082 
1083 /**
1084   @brief Set `CeedTensorContract` of a `CeedBasis`
1085 
1086   @param[in,out] basis    `CeedBasis`
1087   @param[in]     contract `CeedTensorContract` to set
1088 
1089   @return An error code: 0 - success, otherwise - failure
1090 
1091   @ref Backend
1092 **/
1093 int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract contract) {
1094   basis->contract = contract;
1095   CeedCall(CeedTensorContractReference(contract));
1096   return CEED_ERROR_SUCCESS;
1097 }
1098 
1099 /**
1100   @brief Return a reference implementation of matrix multiplication \f$C = A B\f$.
1101 
1102   Note: This is a reference implementation for CPU `CeedScalar` pointers that is not intended for high performance.
1103 
1104   @param[in]  ceed  `Ceed` context for error handling
1105   @param[in]  mat_A Row-major matrix `A`
1106   @param[in]  mat_B Row-major matrix `B`
1107   @param[out] mat_C Row-major output matrix `C`
1108   @param[in]  m     Number of rows of `C`
1109   @param[in]  n     Number of columns of `C`
1110   @param[in]  kk    Number of columns of `A`/rows of `B`
1111 
1112   @return An error code: 0 - success, otherwise - failure
1113 
1114   @ref Utility
1115 **/
1116 int CeedMatrixMatrixMultiply(Ceed ceed, const CeedScalar *mat_A, const CeedScalar *mat_B, CeedScalar *mat_C, CeedInt m, CeedInt n, CeedInt kk) {
1117   for (CeedInt i = 0; i < m; i++) {
1118     for (CeedInt j = 0; j < n; j++) {
1119       CeedScalar sum = 0;
1120 
1121       for (CeedInt k = 0; k < kk; k++) sum += mat_A[k + i * kk] * mat_B[j + k * n];
1122       mat_C[j + i * n] = sum;
1123     }
1124   }
1125   return CEED_ERROR_SUCCESS;
1126 }
1127 
1128 /**
1129   @brief Return QR Factorization of a matrix
1130 
1131   @param[in]     ceed `Ceed` context for error handling
1132   @param[in,out] mat  Row-major matrix to be factorized in place
1133   @param[in,out] tau  Vector of length `m` of scaling factors
1134   @param[in]     m    Number of rows
1135   @param[in]     n    Number of columns
1136 
1137   @return An error code: 0 - success, otherwise - failure
1138 
1139   @ref Utility
1140 **/
1141 int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau, CeedInt m, CeedInt n) {
1142   CeedScalar v[m];
1143 
1144   // Check matrix shape
1145   CeedCheck(n <= m, ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute QR factorization with n > m");
1146 
1147   for (CeedInt i = 0; i < n; i++) {
1148     CeedScalar sigma = 0.0;
1149 
1150     if (i >= m - 1) {  // last row of matrix, no reflection needed
1151       tau[i] = 0.;
1152       break;
1153     }
1154     // Calculate Householder vector, magnitude
1155     v[i] = mat[i + n * i];
1156     for (CeedInt j = i + 1; j < m; j++) {
1157       v[j] = mat[i + n * j];
1158       sigma += v[j] * v[j];
1159     }
1160     const CeedScalar norm = sqrt(v[i] * v[i] + sigma);  // norm of v[i:m]
1161     const CeedScalar R_ii = -copysign(norm, v[i]);
1162 
1163     v[i] -= R_ii;
1164     // norm of v[i:m] after modification above and scaling below
1165     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
1166     //   tau = 2 / (norm*norm)
1167     tau[i] = 2 * v[i] * v[i] / (v[i] * v[i] + sigma);
1168     for (CeedInt j = i + 1; j < m; j++) v[j] /= v[i];
1169 
1170     // Apply Householder reflector to lower right panel
1171     CeedHouseholderReflect(&mat[i * n + i + 1], &v[i], tau[i], m - i, n - i - 1, n, 1);
1172     // Save v
1173     mat[i + n * i] = R_ii;
1174     for (CeedInt j = i + 1; j < m; j++) mat[i + n * j] = v[j];
1175   }
1176   return CEED_ERROR_SUCCESS;
1177 }
1178 
1179 /**
1180   @brief Apply Householder Q matrix
1181 
1182   Compute `mat_A = mat_Q mat_A`, where `mat_Q` is \f$m \times m\f$ and `mat_A` is \f$m \times n\f$.
1183 
1184   @param[in,out] mat_A  Matrix to apply Householder Q to, in place
1185   @param[in]     mat_Q  Householder Q matrix
1186   @param[in]     tau    Householder scaling factors
1187   @param[in]     t_mode Transpose mode for application
1188   @param[in]     m      Number of rows in `A`
1189   @param[in]     n      Number of columns in `A`
1190   @param[in]     k      Number of elementary reflectors in Q, `k < m`
1191   @param[in]     row    Row stride in `A`
1192   @param[in]     col    Col stride in `A`
1193 
1194   @return An error code: 0 - success, otherwise - failure
1195 
1196   @ref Utility
1197 **/
1198 int CeedHouseholderApplyQ(CeedScalar *mat_A, const CeedScalar *mat_Q, const CeedScalar *tau, CeedTransposeMode t_mode, CeedInt m, CeedInt n,
1199                           CeedInt k, CeedInt row, CeedInt col) {
1200   CeedScalar *v;
1201 
1202   CeedCall(CeedMalloc(m, &v));
1203   for (CeedInt ii = 0; ii < k; ii++) {
1204     CeedInt i = t_mode == CEED_TRANSPOSE ? ii : k - 1 - ii;
1205     for (CeedInt j = i + 1; j < m; j++) v[j] = mat_Q[j * k + i];
1206     // Apply Householder reflector (I - tau v v^T) collo_grad_1d^T
1207     CeedCall(CeedHouseholderReflect(&mat_A[i * row], &v[i], tau[i], m - i, n, row, col));
1208   }
1209   CeedCall(CeedFree(&v));
1210   return CEED_ERROR_SUCCESS;
1211 }
1212 
1213 /**
1214   @brief Return pseudoinverse of a matrix
1215 
1216   @param[in]     ceed      Ceed context for error handling
1217   @param[in]     mat       Row-major matrix to compute pseudoinverse of
1218   @param[in]     m         Number of rows
1219   @param[in]     n         Number of columns
1220   @param[out]    mat_pinv  Row-major pseudoinverse matrix
1221 
1222   @return An error code: 0 - success, otherwise - failure
1223 
1224   @ref Utility
1225 **/
1226 int CeedMatrixPseudoinverse(Ceed ceed, const CeedScalar *mat, CeedInt m, CeedInt n, CeedScalar *mat_pinv) {
1227   CeedScalar *tau, *I, *mat_copy;
1228 
1229   CeedCall(CeedCalloc(m, &tau));
1230   CeedCall(CeedCalloc(m * m, &I));
1231   CeedCall(CeedCalloc(m * n, &mat_copy));
1232   memcpy(mat_copy, mat, m * n * sizeof mat[0]);
1233 
1234   // QR Factorization, mat = Q R
1235   CeedCall(CeedQRFactorization(ceed, mat_copy, tau, m, n));
1236 
1237   // -- Apply Q^T, I = Q^T * I
1238   for (CeedInt i = 0; i < m; i++) I[i * m + i] = 1.0;
1239   CeedCall(CeedHouseholderApplyQ(I, mat_copy, tau, CEED_TRANSPOSE, m, m, n, m, 1));
1240   // -- Apply R_inv, mat_pinv = R_inv * Q^T
1241   for (CeedInt j = 0; j < m; j++) {  // Column j
1242     mat_pinv[j + m * (n - 1)] = I[j + m * (n - 1)] / mat_copy[n * n - 1];
1243     for (CeedInt i = n - 2; i >= 0; i--) {  // Row i
1244       mat_pinv[j + m * i] = I[j + m * i];
1245       for (CeedInt k = i + 1; k < n; k++) mat_pinv[j + m * i] -= mat_copy[k + n * i] * mat_pinv[j + m * k];
1246       mat_pinv[j + m * i] /= mat_copy[i + n * i];
1247     }
1248   }
1249 
1250   // Cleanup
1251   CeedCall(CeedFree(&I));
1252   CeedCall(CeedFree(&tau));
1253   CeedCall(CeedFree(&mat_copy));
1254   return CEED_ERROR_SUCCESS;
1255 }
1256 
1257 /**
1258   @brief Return symmetric Schur decomposition of the symmetric matrix mat via symmetric QR factorization
1259 
1260   @param[in]     ceed   `Ceed` context for error handling
1261   @param[in,out] mat    Row-major matrix to be factorized in place
1262   @param[out]    lambda Vector of length n of eigenvalues
1263   @param[in]     n      Number of rows/columns
1264 
1265   @return An error code: 0 - success, otherwise - failure
1266 
1267   @ref Utility
1268 **/
1269 CeedPragmaOptimizeOff
1270 int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat, CeedScalar *lambda, CeedInt n) {
1271   // Check bounds for clang-tidy
1272   CeedCheck(n > 1, ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute symmetric Schur decomposition of scalars");
1273 
1274   CeedScalar v[n - 1], tau[n - 1], mat_T[n * n];
1275 
1276   // Copy mat to mat_T and set mat to I
1277   memcpy(mat_T, mat, n * n * sizeof(mat[0]));
1278   for (CeedInt i = 0; i < n; i++) {
1279     for (CeedInt j = 0; j < n; j++) mat[j + n * i] = (i == j) ? 1 : 0;
1280   }
1281 
1282   // Reduce to tridiagonal
1283   for (CeedInt i = 0; i < n - 1; i++) {
1284     // Calculate Householder vector, magnitude
1285     CeedScalar sigma = 0.0;
1286 
1287     v[i] = mat_T[i + n * (i + 1)];
1288     for (CeedInt j = i + 1; j < n - 1; j++) {
1289       v[j] = mat_T[i + n * (j + 1)];
1290       sigma += v[j] * v[j];
1291     }
1292     const CeedScalar norm = sqrt(v[i] * v[i] + sigma);  // norm of v[i:n-1]
1293     const CeedScalar R_ii = -copysign(norm, v[i]);
1294 
1295     v[i] -= R_ii;
1296     // norm of v[i:m] after modification above and scaling below
1297     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
1298     //   tau = 2 / (norm*norm)
1299     tau[i] = i == n - 2 ? 2 : 2 * v[i] * v[i] / (v[i] * v[i] + sigma);
1300     for (CeedInt j = i + 1; j < n - 1; j++) v[j] /= v[i];
1301 
1302     // Update sub and super diagonal
1303     for (CeedInt j = i + 2; j < n; j++) {
1304       mat_T[i + n * j] = 0;
1305       mat_T[j + n * i] = 0;
1306     }
1307     // Apply symmetric Householder reflector to lower right panel
1308     CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1);
1309     CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), 1, n);
1310 
1311     // Save v
1312     mat_T[i + n * (i + 1)] = R_ii;
1313     mat_T[(i + 1) + n * i] = R_ii;
1314     for (CeedInt j = i + 1; j < n - 1; j++) {
1315       mat_T[i + n * (j + 1)] = v[j];
1316     }
1317   }
1318   // Backwards accumulation of Q
1319   for (CeedInt i = n - 2; i >= 0; i--) {
1320     if (tau[i] > 0.0) {
1321       v[i] = 1;
1322       for (CeedInt j = i + 1; j < n - 1; j++) {
1323         v[j]                   = mat_T[i + n * (j + 1)];
1324         mat_T[i + n * (j + 1)] = 0;
1325       }
1326       CeedHouseholderReflect(&mat[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1);
1327     }
1328   }
1329 
1330   // Reduce sub and super diagonal
1331   CeedInt    p = 0, q = 0, itr = 0, max_itr = n * n * n * n;
1332   CeedScalar tol = CEED_EPSILON;
1333 
1334   while (itr < max_itr) {
1335     // Update p, q, size of reduced portions of diagonal
1336     p = 0;
1337     q = 0;
1338     for (CeedInt i = n - 2; i >= 0; i--) {
1339       if (fabs(mat_T[i + n * (i + 1)]) < tol) q += 1;
1340       else break;
1341     }
1342     for (CeedInt i = 0; i < n - q - 1; i++) {
1343       if (fabs(mat_T[i + n * (i + 1)]) < tol) p += 1;
1344       else break;
1345     }
1346     if (q == n - 1) break;  // Finished reducing
1347 
1348     // Reduce tridiagonal portion
1349     CeedScalar t_nn = mat_T[(n - 1 - q) + n * (n - 1 - q)], t_nnm1 = mat_T[(n - 2 - q) + n * (n - 1 - q)];
1350     CeedScalar d  = (mat_T[(n - 2 - q) + n * (n - 2 - q)] - t_nn) / 2;
1351     CeedScalar mu = t_nn - t_nnm1 * t_nnm1 / (d + copysign(sqrt(d * d + t_nnm1 * t_nnm1), d));
1352     CeedScalar x  = mat_T[p + n * p] - mu;
1353     CeedScalar z  = mat_T[p + n * (p + 1)];
1354 
1355     for (CeedInt k = p; k < n - q - 1; k++) {
1356       // Compute Givens rotation
1357       CeedScalar c = 1, s = 0;
1358 
1359       if (fabs(z) > tol) {
1360         if (fabs(z) > fabs(x)) {
1361           const CeedScalar tau = -x / z;
1362 
1363           s = 1 / sqrt(1 + tau * tau);
1364           c = s * tau;
1365         } else {
1366           const CeedScalar tau = -z / x;
1367 
1368           c = 1 / sqrt(1 + tau * tau);
1369           s = c * tau;
1370         }
1371       }
1372 
1373       // Apply Givens rotation to T
1374       CeedGivensRotation(mat_T, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n);
1375       CeedGivensRotation(mat_T, c, s, CEED_TRANSPOSE, k, k + 1, n, n);
1376 
1377       // Apply Givens rotation to Q
1378       CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n);
1379 
1380       // Update x, z
1381       if (k < n - q - 2) {
1382         x = mat_T[k + n * (k + 1)];
1383         z = mat_T[k + n * (k + 2)];
1384       }
1385     }
1386     itr++;
1387   }
1388 
1389   // Save eigenvalues
1390   for (CeedInt i = 0; i < n; i++) lambda[i] = mat_T[i + n * i];
1391 
1392   // Check convergence
1393   CeedCheck(itr < max_itr || q > n, ceed, CEED_ERROR_MINOR, "Symmetric QR failed to converge");
1394   return CEED_ERROR_SUCCESS;
1395 }
1396 CeedPragmaOptimizeOn
1397 
1398 /**
1399   @brief Return Simultaneous Diagonalization of two matrices.
1400 
1401   This solves the generalized eigenvalue problem `A x = lambda B x`, where `A` and `B` are symmetric and `B` is positive definite.
1402   We generate the matrix `X` and vector `Lambda` such that `X^T A X = Lambda` and `X^T B X = I`.
1403   This is equivalent to the LAPACK routine 'sygv' with `TYPE = 1`.
1404 
1405   @param[in]  ceed   `Ceed` context for error handling
1406   @param[in]  mat_A  Row-major matrix to be factorized with eigenvalues
1407   @param[in]  mat_B  Row-major matrix to be factorized to identity
1408   @param[out] mat_X  Row-major orthogonal matrix
1409   @param[out] lambda Vector of length `n` of generalized eigenvalues
1410   @param[in]  n      Number of rows/columns
1411 
1412   @return An error code: 0 - success, otherwise - failure
1413 
1414   @ref Utility
1415 **/
1416 CeedPragmaOptimizeOff
1417 int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *mat_A, CeedScalar *mat_B, CeedScalar *mat_X, CeedScalar *lambda, CeedInt n) {
1418   CeedScalar *mat_C, *mat_G, *vec_D;
1419 
1420   CeedCall(CeedCalloc(n * n, &mat_C));
1421   CeedCall(CeedCalloc(n * n, &mat_G));
1422   CeedCall(CeedCalloc(n, &vec_D));
1423 
1424   // Compute B = G D G^T
1425   memcpy(mat_G, mat_B, n * n * sizeof(mat_B[0]));
1426   CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_G, vec_D, n));
1427 
1428   // Sort eigenvalues
1429   for (CeedInt i = n - 1; i >= 0; i--) {
1430     for (CeedInt j = 0; j < i; j++) {
1431       if (fabs(vec_D[j]) > fabs(vec_D[j + 1])) {
1432         CeedScalarSwap(vec_D[j], vec_D[j + 1]);
1433         for (CeedInt k = 0; k < n; k++) CeedScalarSwap(mat_G[k * n + j], mat_G[k * n + j + 1]);
1434       }
1435     }
1436   }
1437 
1438   // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T
1439   //           = D^-1/2 G^T A G D^-1/2
1440   // -- D = D^-1/2
1441   for (CeedInt i = 0; i < n; i++) vec_D[i] = 1. / sqrt(vec_D[i]);
1442   // -- G = G D^-1/2
1443   // -- C = D^-1/2 G^T
1444   for (CeedInt i = 0; i < n; i++) {
1445     for (CeedInt j = 0; j < n; j++) {
1446       mat_G[i * n + j] *= vec_D[j];
1447       mat_C[j * n + i] = mat_G[i * n + j];
1448     }
1449   }
1450   // -- X = (D^-1/2 G^T) A
1451   CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_C, (const CeedScalar *)mat_A, mat_X, n, n, n));
1452   // -- C = (D^-1/2 G^T A) (G D^-1/2)
1453   CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_X, (const CeedScalar *)mat_G, mat_C, n, n, n));
1454 
1455   // Compute Q^T C Q = lambda
1456   CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_C, lambda, n));
1457 
1458   // Sort eigenvalues
1459   for (CeedInt i = n - 1; i >= 0; i--) {
1460     for (CeedInt j = 0; j < i; j++) {
1461       if (fabs(lambda[j]) > fabs(lambda[j + 1])) {
1462         CeedScalarSwap(lambda[j], lambda[j + 1]);
1463         for (CeedInt k = 0; k < n; k++) CeedScalarSwap(mat_C[k * n + j], mat_C[k * n + j + 1]);
1464       }
1465     }
1466   }
1467 
1468   // Set X = (G D^1/2)^-T Q
1469   //       = G D^-1/2 Q
1470   CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_G, (const CeedScalar *)mat_C, mat_X, n, n, n));
1471 
1472   // Cleanup
1473   CeedCall(CeedFree(&mat_C));
1474   CeedCall(CeedFree(&mat_G));
1475   CeedCall(CeedFree(&vec_D));
1476   return CEED_ERROR_SUCCESS;
1477 }
1478 CeedPragmaOptimizeOn
1479 
1480 /// @}
1481 
1482 /// ----------------------------------------------------------------------------
1483 /// CeedBasis Public API
1484 /// ----------------------------------------------------------------------------
1485 /// @addtogroup CeedBasisUser
1486 /// @{
1487 
1488 /**
1489   @brief Create a tensor-product basis for \f$H^1\f$ discretizations
1490 
1491   @param[in]  ceed        `Ceed` object used to create the `CeedBasis`
1492   @param[in]  dim         Topological dimension
1493   @param[in]  num_comp    Number of field components (1 for scalar fields)
1494   @param[in]  P_1d        Number of nodes in one dimension
1495   @param[in]  Q_1d        Number of quadrature points in one dimension
1496   @param[in]  interp_1d   Row-major (`Q_1d * P_1d`) matrix expressing the values of nodal basis functions at quadrature points
1497   @param[in]  grad_1d     Row-major (`Q_1d * P_1d`) matrix expressing derivatives of nodal basis functions at quadrature points
1498   @param[in]  q_ref_1d    Array of length `Q_1d` holding the locations of quadrature points on the 1D reference element `[-1, 1]`
1499   @param[in]  q_weight_1d Array of length `Q_1d` holding the quadrature weights on the reference element
1500   @param[out] basis       Address of the variable where the newly created `CeedBasis` will be stored
1501 
1502   @return An error code: 0 - success, otherwise - failure
1503 
1504   @ref User
1505 **/
1506 int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d,
1507                             const CeedScalar *grad_1d, const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis *basis) {
1508   if (!ceed->BasisCreateTensorH1) {
1509     Ceed delegate;
1510 
1511     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
1512     CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateTensorH1");
1513     CeedCall(CeedBasisCreateTensorH1(delegate, dim, num_comp, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis));
1514     CeedCall(CeedDestroy(&delegate));
1515     return CEED_ERROR_SUCCESS;
1516   }
1517 
1518   CeedCheck(dim > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis dimension must be a positive value");
1519   CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 component");
1520   CeedCheck(P_1d > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 node");
1521   CeedCheck(Q_1d > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 quadrature point");
1522 
1523   CeedElemTopology topo = dim == 1 ? CEED_TOPOLOGY_LINE : dim == 2 ? CEED_TOPOLOGY_QUAD : CEED_TOPOLOGY_HEX;
1524 
1525   CeedCall(CeedCalloc(1, basis));
1526   CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed));
1527   (*basis)->ref_count       = 1;
1528   (*basis)->is_tensor_basis = true;
1529   (*basis)->dim             = dim;
1530   (*basis)->topo            = topo;
1531   (*basis)->num_comp        = num_comp;
1532   (*basis)->P_1d            = P_1d;
1533   (*basis)->Q_1d            = Q_1d;
1534   (*basis)->P               = CeedIntPow(P_1d, dim);
1535   (*basis)->Q               = CeedIntPow(Q_1d, dim);
1536   (*basis)->fe_space        = CEED_FE_SPACE_H1;
1537   CeedCall(CeedCalloc(Q_1d, &(*basis)->q_ref_1d));
1538   CeedCall(CeedCalloc(Q_1d, &(*basis)->q_weight_1d));
1539   if (q_ref_1d) memcpy((*basis)->q_ref_1d, q_ref_1d, Q_1d * sizeof(q_ref_1d[0]));
1540   if (q_weight_1d) memcpy((*basis)->q_weight_1d, q_weight_1d, Q_1d * sizeof(q_weight_1d[0]));
1541   CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->interp_1d));
1542   CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->grad_1d));
1543   if (interp_1d) memcpy((*basis)->interp_1d, interp_1d, Q_1d * P_1d * sizeof(interp_1d[0]));
1544   if (grad_1d) memcpy((*basis)->grad_1d, grad_1d, Q_1d * P_1d * sizeof(grad_1d[0]));
1545   CeedCall(ceed->BasisCreateTensorH1(dim, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, *basis));
1546   return CEED_ERROR_SUCCESS;
1547 }
1548 
1549 /**
1550   @brief Create a tensor-product \f$H^1\f$ Lagrange basis
1551 
1552   @param[in]  ceed      `Ceed` object used to create the `CeedBasis`
1553   @param[in]  dim       Topological dimension of element
1554   @param[in]  num_comp  Number of field components (1 for scalar fields)
1555   @param[in]  P         Number of Gauss-Lobatto nodes in one dimension.
1556                           The polynomial degree of the resulting `Q_k` element is `k = P - 1`.
1557   @param[in]  Q         Number of quadrature points in one dimension.
1558   @param[in]  quad_mode Distribution of the `Q` quadrature points (affects order of accuracy for the quadrature)
1559   @param[out] basis     Address of the variable where the newly created `CeedBasis` will be stored
1560 
1561   @return An error code: 0 - success, otherwise - failure
1562 
1563   @ref User
1564 **/
1565 int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P, CeedInt Q, CeedQuadMode quad_mode, CeedBasis *basis) {
1566   // Allocate
1567   int        ierr = CEED_ERROR_SUCCESS;
1568   CeedScalar c1, c2, c3, c4, dx, *nodes, *interp_1d, *grad_1d, *q_ref_1d, *q_weight_1d;
1569 
1570   CeedCheck(dim > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis dimension must be a positive value");
1571   CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 component");
1572   CeedCheck(P > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 node");
1573   CeedCheck(Q > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 quadrature point");
1574 
1575   // Get Nodes and Weights
1576   CeedCall(CeedCalloc(P * Q, &interp_1d));
1577   CeedCall(CeedCalloc(P * Q, &grad_1d));
1578   CeedCall(CeedCalloc(P, &nodes));
1579   CeedCall(CeedCalloc(Q, &q_ref_1d));
1580   CeedCall(CeedCalloc(Q, &q_weight_1d));
1581   if (CeedLobattoQuadrature(P, nodes, NULL) != CEED_ERROR_SUCCESS) goto cleanup;
1582   switch (quad_mode) {
1583     case CEED_GAUSS:
1584       ierr = CeedGaussQuadrature(Q, q_ref_1d, q_weight_1d);
1585       break;
1586     case CEED_GAUSS_LOBATTO:
1587       ierr = CeedLobattoQuadrature(Q, q_ref_1d, q_weight_1d);
1588       break;
1589   }
1590   if (ierr != CEED_ERROR_SUCCESS) goto cleanup;
1591 
1592   // Build B, D matrix
1593   // Fornberg, 1998
1594   for (CeedInt i = 0; i < Q; i++) {
1595     c1                   = 1.0;
1596     c3                   = nodes[0] - q_ref_1d[i];
1597     interp_1d[i * P + 0] = 1.0;
1598     for (CeedInt j = 1; j < P; j++) {
1599       c2 = 1.0;
1600       c4 = c3;
1601       c3 = nodes[j] - q_ref_1d[i];
1602       for (CeedInt k = 0; k < j; k++) {
1603         dx = nodes[j] - nodes[k];
1604         c2 *= dx;
1605         if (k == j - 1) {
1606           grad_1d[i * P + j]   = c1 * (interp_1d[i * P + k] - c4 * grad_1d[i * P + k]) / c2;
1607           interp_1d[i * P + j] = -c1 * c4 * interp_1d[i * P + k] / c2;
1608         }
1609         grad_1d[i * P + k]   = (c3 * grad_1d[i * P + k] - interp_1d[i * P + k]) / dx;
1610         interp_1d[i * P + k] = c3 * interp_1d[i * P + k] / dx;
1611       }
1612       c1 = c2;
1613     }
1614   }
1615   // Pass to CeedBasisCreateTensorH1
1616   CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P, Q, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis));
1617 cleanup:
1618   CeedCall(CeedFree(&interp_1d));
1619   CeedCall(CeedFree(&grad_1d));
1620   CeedCall(CeedFree(&nodes));
1621   CeedCall(CeedFree(&q_ref_1d));
1622   CeedCall(CeedFree(&q_weight_1d));
1623   return CEED_ERROR_SUCCESS;
1624 }
1625 
1626 /**
1627   @brief Create a non tensor-product basis for \f$H^1\f$ discretizations
1628 
1629   @param[in]  ceed      `Ceed` object used to create the `CeedBasis`
1630   @param[in]  topo      Topology of element, e.g. hypercube, simplex, etc
1631   @param[in]  num_comp  Number of field components (1 for scalar fields)
1632   @param[in]  num_nodes Total number of nodes
1633   @param[in]  num_qpts  Total number of quadrature points
1634   @param[in]  interp    Row-major (`num_qpts * num_nodes`) matrix expressing the values of nodal basis functions at quadrature points
1635   @param[in]  grad      Row-major (`dim * num_qpts * num_nodes`) matrix expressing derivatives of nodal basis functions at quadrature points
1636   @param[in]  q_ref     Array of length `num_qpts * dim` holding the locations of quadrature points on the reference element
1637   @param[in]  q_weight  Array of length `num_qpts` holding the quadrature weights on the reference element
1638   @param[out] basis     Address of the variable where the newly created `CeedBasis` will be stored
1639 
1640   @return An error code: 0 - success, otherwise - failure
1641 
1642   @ref User
1643 **/
1644 int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
1645                       const CeedScalar *grad, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) {
1646   CeedInt P = num_nodes, Q = num_qpts, dim = 0;
1647 
1648   if (!ceed->BasisCreateH1) {
1649     Ceed delegate;
1650 
1651     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
1652     CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateH1");
1653     CeedCall(CeedBasisCreateH1(delegate, topo, num_comp, num_nodes, num_qpts, interp, grad, q_ref, q_weight, basis));
1654     CeedCall(CeedDestroy(&delegate));
1655     return CEED_ERROR_SUCCESS;
1656   }
1657 
1658   CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 component");
1659   CeedCheck(num_nodes > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 node");
1660   CeedCheck(num_qpts > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 quadrature point");
1661 
1662   CeedCall(CeedBasisGetTopologyDimension(topo, &dim));
1663 
1664   CeedCall(CeedCalloc(1, basis));
1665   CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed));
1666   (*basis)->ref_count       = 1;
1667   (*basis)->is_tensor_basis = false;
1668   (*basis)->dim             = dim;
1669   (*basis)->topo            = topo;
1670   (*basis)->num_comp        = num_comp;
1671   (*basis)->P               = P;
1672   (*basis)->Q               = Q;
1673   (*basis)->fe_space        = CEED_FE_SPACE_H1;
1674   CeedCall(CeedCalloc(Q * dim, &(*basis)->q_ref_1d));
1675   CeedCall(CeedCalloc(Q, &(*basis)->q_weight_1d));
1676   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0]));
1677   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0]));
1678   CeedCall(CeedCalloc(Q * P, &(*basis)->interp));
1679   CeedCall(CeedCalloc(dim * Q * P, &(*basis)->grad));
1680   if (interp) memcpy((*basis)->interp, interp, Q * P * sizeof(interp[0]));
1681   if (grad) memcpy((*basis)->grad, grad, dim * Q * P * sizeof(grad[0]));
1682   CeedCall(ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, q_ref, q_weight, *basis));
1683   return CEED_ERROR_SUCCESS;
1684 }
1685 
1686 /**
1687   @brief Create a non tensor-product basis for \f$H(\mathrm{div})\f$ discretizations
1688 
1689   @param[in]  ceed      `Ceed` object used to create the `CeedBasis`
1690   @param[in]  topo      Topology of element (`CEED_TOPOLOGY_QUAD`, `CEED_TOPOLOGY_PRISM`, etc.), dimension of which is used in some array sizes below
1691   @param[in]  num_comp  Number of components (usually 1 for vectors in H(div) bases)
1692   @param[in]  num_nodes Total number of nodes (DoFs per element)
1693   @param[in]  num_qpts  Total number of quadrature points
1694   @param[in]  interp    Row-major (`dim * num_qpts * num_nodes`) matrix expressing the values of basis functions at quadrature points
1695   @param[in]  div       Row-major (`num_qpts * num_nodes`) matrix expressing divergence of basis functions at quadrature points
1696   @param[in]  q_ref     Array of length `num_qpts` * dim holding the locations of quadrature points on the reference element
1697   @param[in]  q_weight  Array of length `num_qpts` holding the quadrature weights on the reference element
1698   @param[out] basis     Address of the variable where the newly created `CeedBasis` will be stored
1699 
1700   @return An error code: 0 - success, otherwise - failure
1701 
1702   @ref User
1703 **/
1704 int CeedBasisCreateHdiv(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
1705                         const CeedScalar *div, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) {
1706   CeedInt Q = num_qpts, P = num_nodes, dim = 0;
1707 
1708   if (!ceed->BasisCreateHdiv) {
1709     Ceed delegate;
1710 
1711     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
1712     CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateHdiv");
1713     CeedCall(CeedBasisCreateHdiv(delegate, topo, num_comp, num_nodes, num_qpts, interp, div, q_ref, q_weight, basis));
1714     CeedCall(CeedDestroy(&delegate));
1715     return CEED_ERROR_SUCCESS;
1716   }
1717 
1718   CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 component");
1719   CeedCheck(num_nodes > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 node");
1720   CeedCheck(num_qpts > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 quadrature point");
1721 
1722   CeedCall(CeedBasisGetTopologyDimension(topo, &dim));
1723 
1724   CeedCall(CeedCalloc(1, basis));
1725   CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed));
1726   (*basis)->ref_count       = 1;
1727   (*basis)->is_tensor_basis = false;
1728   (*basis)->dim             = dim;
1729   (*basis)->topo            = topo;
1730   (*basis)->num_comp        = num_comp;
1731   (*basis)->P               = P;
1732   (*basis)->Q               = Q;
1733   (*basis)->fe_space        = CEED_FE_SPACE_HDIV;
1734   CeedCall(CeedMalloc(Q * dim, &(*basis)->q_ref_1d));
1735   CeedCall(CeedMalloc(Q, &(*basis)->q_weight_1d));
1736   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0]));
1737   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0]));
1738   CeedCall(CeedMalloc(dim * Q * P, &(*basis)->interp));
1739   CeedCall(CeedMalloc(Q * P, &(*basis)->div));
1740   if (interp) memcpy((*basis)->interp, interp, dim * Q * P * sizeof(interp[0]));
1741   if (div) memcpy((*basis)->div, div, Q * P * sizeof(div[0]));
1742   CeedCall(ceed->BasisCreateHdiv(topo, dim, P, Q, interp, div, q_ref, q_weight, *basis));
1743   return CEED_ERROR_SUCCESS;
1744 }
1745 
1746 /**
1747   @brief Create a non tensor-product basis for \f$H(\mathrm{curl})\f$ discretizations
1748 
1749   @param[in]  ceed      `Ceed` object used to create the `CeedBasis`
1750   @param[in]  topo      Topology of element (`CEED_TOPOLOGY_QUAD`, `CEED_TOPOLOGY_PRISM`, etc.), dimension of which is used in some array sizes below
1751   @param[in]  num_comp  Number of components (usually 1 for vectors in \f$H(\mathrm{curl})\f$ bases)
1752   @param[in]  num_nodes Total number of nodes (DoFs per element)
1753   @param[in]  num_qpts  Total number of quadrature points
1754   @param[in]  interp    Row-major (`dim * num_qpts * num_nodes`) matrix expressing the values of basis functions at quadrature points
1755   @param[in]  curl      Row-major (`curl_comp * num_qpts * num_nodes`, `curl_comp = 1` if `dim < 3` otherwise `curl_comp = dim`) matrix expressing curl of basis functions at quadrature points
1756   @param[in]  q_ref     Array of length `num_qpts * dim` holding the locations of quadrature points on the reference element
1757   @param[in]  q_weight  Array of length `num_qpts` holding the quadrature weights on the reference element
1758   @param[out] basis     Address of the variable where the newly created `CeedBasis` will be stored
1759 
1760   @return An error code: 0 - success, otherwise - failure
1761 
1762   @ref User
1763 **/
1764 int CeedBasisCreateHcurl(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
1765                          const CeedScalar *curl, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) {
1766   CeedInt Q = num_qpts, P = num_nodes, dim = 0, curl_comp = 0;
1767 
1768   if (!ceed->BasisCreateHcurl) {
1769     Ceed delegate;
1770 
1771     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
1772     CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateHcurl");
1773     CeedCall(CeedBasisCreateHcurl(delegate, topo, num_comp, num_nodes, num_qpts, interp, curl, q_ref, q_weight, basis));
1774     CeedCall(CeedDestroy(&delegate));
1775     return CEED_ERROR_SUCCESS;
1776   }
1777 
1778   CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 component");
1779   CeedCheck(num_nodes > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 node");
1780   CeedCheck(num_qpts > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 quadrature point");
1781 
1782   CeedCall(CeedBasisGetTopologyDimension(topo, &dim));
1783   curl_comp = (dim < 3) ? 1 : dim;
1784 
1785   CeedCall(CeedCalloc(1, basis));
1786   CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed));
1787   (*basis)->ref_count       = 1;
1788   (*basis)->is_tensor_basis = false;
1789   (*basis)->dim             = dim;
1790   (*basis)->topo            = topo;
1791   (*basis)->num_comp        = num_comp;
1792   (*basis)->P               = P;
1793   (*basis)->Q               = Q;
1794   (*basis)->fe_space        = CEED_FE_SPACE_HCURL;
1795   CeedCall(CeedMalloc(Q * dim, &(*basis)->q_ref_1d));
1796   CeedCall(CeedMalloc(Q, &(*basis)->q_weight_1d));
1797   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0]));
1798   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0]));
1799   CeedCall(CeedMalloc(dim * Q * P, &(*basis)->interp));
1800   CeedCall(CeedMalloc(curl_comp * Q * P, &(*basis)->curl));
1801   if (interp) memcpy((*basis)->interp, interp, dim * Q * P * sizeof(interp[0]));
1802   if (curl) memcpy((*basis)->curl, curl, curl_comp * Q * P * sizeof(curl[0]));
1803   CeedCall(ceed->BasisCreateHcurl(topo, dim, P, Q, interp, curl, q_ref, q_weight, *basis));
1804   return CEED_ERROR_SUCCESS;
1805 }
1806 
1807 /**
1808   @brief Create a `CeedBasis` for projection from the nodes of `basis_from` to the nodes of `basis_to`.
1809 
1810   Only @ref CEED_EVAL_INTERP will be valid for the new basis, `basis_project`.
1811   For \f$H^1\f$ spaces, @ref CEED_EVAL_GRAD will also be valid.
1812   The interpolation is given by `interp_project = interp_to^+ * interp_from`, where the pseudoinverse `interp_to^+` is given by QR factorization.
1813   The gradient (for the \f$H^1\f$ case) is given by `grad_project = interp_to^+ * grad_from`.
1814 
1815   Note: `basis_from` and `basis_to` must have compatible quadrature spaces.
1816 
1817   Note: `basis_project` will have the same number of components as `basis_from`, regardless of the number of components that `basis_to` has.
1818         If `basis_from` has 3 components and `basis_to` has 5 components, then `basis_project` will have 3 components.
1819 
1820   Note: If either `basis_from` or `basis_to` are non-tensor, then `basis_project` will also be non-tensor
1821 
1822   @param[in]  basis_from    `CeedBasis` to prolong from
1823   @param[in]  basis_to      `CeedBasis` to prolong to
1824   @param[out] basis_project Address of the variable where the newly created `CeedBasis` will be stored
1825 
1826   @return An error code: 0 - success, otherwise - failure
1827 
1828   @ref User
1829 **/
1830 int CeedBasisCreateProjection(CeedBasis basis_from, CeedBasis basis_to, CeedBasis *basis_project) {
1831   Ceed        ceed;
1832   bool        create_tensor;
1833   CeedInt     dim, num_comp;
1834   CeedScalar *interp_project, *grad_project;
1835 
1836   CeedCall(CeedBasisGetCeed(basis_to, &ceed));
1837 
1838   // Create projection matrix
1839   CeedCall(CeedBasisCreateProjectionMatrices(basis_from, basis_to, &interp_project, &grad_project));
1840 
1841   // Build basis
1842   {
1843     bool is_tensor_to, is_tensor_from;
1844 
1845     CeedCall(CeedBasisIsTensor(basis_to, &is_tensor_to));
1846     CeedCall(CeedBasisIsTensor(basis_from, &is_tensor_from));
1847     create_tensor = is_tensor_from && is_tensor_to;
1848   }
1849   CeedCall(CeedBasisGetDimension(basis_to, &dim));
1850   CeedCall(CeedBasisGetNumComponents(basis_from, &num_comp));
1851   if (create_tensor) {
1852     CeedInt P_1d_to, P_1d_from;
1853 
1854     CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_1d_from));
1855     CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_1d_to));
1856     CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P_1d_from, P_1d_to, interp_project, grad_project, NULL, NULL, basis_project));
1857   } else {
1858     // Even if basis_to and basis_from are not H1, the resulting basis is H1 for interpolation to work
1859     CeedInt          num_nodes_to, num_nodes_from;
1860     CeedElemTopology topo;
1861 
1862     CeedCall(CeedBasisGetTopology(basis_from, &topo));
1863     CeedCall(CeedBasisGetNumNodes(basis_from, &num_nodes_from));
1864     CeedCall(CeedBasisGetNumNodes(basis_to, &num_nodes_to));
1865     CeedCall(CeedBasisCreateH1(ceed, topo, num_comp, num_nodes_from, num_nodes_to, interp_project, grad_project, NULL, NULL, basis_project));
1866   }
1867 
1868   // Cleanup
1869   CeedCall(CeedFree(&interp_project));
1870   CeedCall(CeedFree(&grad_project));
1871   CeedCall(CeedDestroy(&ceed));
1872   return CEED_ERROR_SUCCESS;
1873 }
1874 
1875 /**
1876   @brief Copy the pointer to a `CeedBasis`.
1877 
1878   Note: If the value of `*basis_copy` passed into this function is non-`NULL`, then it is assumed that `*basis_copy` is a pointer to a `CeedBasis`.
1879         This `CeedBasis` will be destroyed if `*basis_copy` is the only reference to this `CeedBasis`.
1880 
1881   @param[in]     basis      `CeedBasis` to copy reference to
1882   @param[in,out] basis_copy Variable to store copied reference
1883 
1884   @return An error code: 0 - success, otherwise - failure
1885 
1886   @ref User
1887 **/
1888 int CeedBasisReferenceCopy(CeedBasis basis, CeedBasis *basis_copy) {
1889   if (basis != CEED_BASIS_NONE) CeedCall(CeedBasisReference(basis));
1890   CeedCall(CeedBasisDestroy(basis_copy));
1891   *basis_copy = basis;
1892   return CEED_ERROR_SUCCESS;
1893 }
1894 
1895 /**
1896   @brief View a `CeedBasis`
1897 
1898   @param[in] basis  `CeedBasis` to view
1899   @param[in] stream Stream to view to, e.g., `stdout`
1900 
1901   @return An error code: 0 - success, otherwise - failure
1902 
1903   @ref User
1904 **/
1905 int CeedBasisView(CeedBasis basis, FILE *stream) {
1906   bool             is_tensor_basis;
1907   CeedElemTopology topo;
1908   CeedFESpace      fe_space;
1909 
1910   // Basis data
1911   CeedCall(CeedBasisIsTensor(basis, &is_tensor_basis));
1912   CeedCall(CeedBasisGetTopology(basis, &topo));
1913   CeedCall(CeedBasisGetFESpace(basis, &fe_space));
1914 
1915   // Print FE space and element topology of the basis
1916   fprintf(stream, "CeedBasis in a %s on a %s element\n", CeedFESpaces[fe_space], CeedElemTopologies[topo]);
1917   if (is_tensor_basis) {
1918     fprintf(stream, "  P: %" CeedInt_FMT "\n  Q: %" CeedInt_FMT "\n", basis->P_1d, basis->Q_1d);
1919   } else {
1920     fprintf(stream, "  P: %" CeedInt_FMT "\n  Q: %" CeedInt_FMT "\n", basis->P, basis->Q);
1921   }
1922   fprintf(stream, "  dimension: %" CeedInt_FMT "\n  field components: %" CeedInt_FMT "\n", basis->dim, basis->num_comp);
1923   // Print quadrature data, interpolation/gradient/divergence/curl of the basis
1924   if (is_tensor_basis) {  // tensor basis
1925     CeedInt           P_1d, Q_1d;
1926     const CeedScalar *q_ref_1d, *q_weight_1d, *interp_1d, *grad_1d;
1927 
1928     CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d));
1929     CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
1930     CeedCall(CeedBasisGetQRef(basis, &q_ref_1d));
1931     CeedCall(CeedBasisGetQWeights(basis, &q_weight_1d));
1932     CeedCall(CeedBasisGetInterp1D(basis, &interp_1d));
1933     CeedCall(CeedBasisGetGrad1D(basis, &grad_1d));
1934 
1935     CeedCall(CeedScalarView("qref1d", "\t% 12.8f", 1, Q_1d, q_ref_1d, stream));
1936     CeedCall(CeedScalarView("qweight1d", "\t% 12.8f", 1, Q_1d, q_weight_1d, stream));
1937     CeedCall(CeedScalarView("interp1d", "\t% 12.8f", Q_1d, P_1d, interp_1d, stream));
1938     CeedCall(CeedScalarView("grad1d", "\t% 12.8f", Q_1d, P_1d, grad_1d, stream));
1939   } else {  // non-tensor basis
1940     CeedInt           P, Q, dim, q_comp;
1941     const CeedScalar *q_ref, *q_weight, *interp, *grad, *div, *curl;
1942 
1943     CeedCall(CeedBasisGetNumNodes(basis, &P));
1944     CeedCall(CeedBasisGetNumQuadraturePoints(basis, &Q));
1945     CeedCall(CeedBasisGetDimension(basis, &dim));
1946     CeedCall(CeedBasisGetQRef(basis, &q_ref));
1947     CeedCall(CeedBasisGetQWeights(basis, &q_weight));
1948     CeedCall(CeedBasisGetInterp(basis, &interp));
1949     CeedCall(CeedBasisGetGrad(basis, &grad));
1950     CeedCall(CeedBasisGetDiv(basis, &div));
1951     CeedCall(CeedBasisGetCurl(basis, &curl));
1952 
1953     CeedCall(CeedScalarView("qref", "\t% 12.8f", 1, Q * dim, q_ref, stream));
1954     CeedCall(CeedScalarView("qweight", "\t% 12.8f", 1, Q, q_weight, stream));
1955     CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp));
1956     CeedCall(CeedScalarView("interp", "\t% 12.8f", q_comp * Q, P, interp, stream));
1957     if (grad) {
1958       CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp));
1959       CeedCall(CeedScalarView("grad", "\t% 12.8f", q_comp * Q, P, grad, stream));
1960     }
1961     if (div) {
1962       CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_DIV, &q_comp));
1963       CeedCall(CeedScalarView("div", "\t% 12.8f", q_comp * Q, P, div, stream));
1964     }
1965     if (curl) {
1966       CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_CURL, &q_comp));
1967       CeedCall(CeedScalarView("curl", "\t% 12.8f", q_comp * Q, P, curl, stream));
1968     }
1969   }
1970   return CEED_ERROR_SUCCESS;
1971 }
1972 
1973 /**
1974   @brief Apply basis evaluation from nodes to quadrature points or vice versa
1975 
1976   @param[in]  basis     `CeedBasis` to evaluate
1977   @param[in]  num_elem  The number of elements to apply the basis evaluation to;
1978                           the backend will specify the ordering in @ref CeedElemRestrictionCreate()
1979   @param[in]  t_mode    @ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature points;
1980                           @ref CEED_TRANSPOSE to apply the transpose, mapping from quadrature points to nodes
1981   @param[in]  eval_mode @ref CEED_EVAL_NONE to use values directly,
1982                           @ref CEED_EVAL_INTERP to use interpolated values,
1983                           @ref CEED_EVAL_GRAD to use gradients,
1984                           @ref CEED_EVAL_DIV to use divergence,
1985                           @ref CEED_EVAL_CURL to use curl,
1986                           @ref CEED_EVAL_WEIGHT to use quadrature weights
1987   @param[in]  u         Input `CeedVector`
1988   @param[out] v         Output `CeedVector`
1989 
1990   @return An error code: 0 - success, otherwise - failure
1991 
1992   @ref User
1993 **/
1994 int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
1995   CeedCall(CeedBasisApplyCheckDims(basis, num_elem, t_mode, eval_mode, u, v));
1996   CeedCheck(basis->Apply, CeedBasisReturnCeed(basis), CEED_ERROR_UNSUPPORTED, "Backend does not support CeedBasisApply");
1997   CeedCall(basis->Apply(basis, num_elem, t_mode, eval_mode, u, v));
1998   return CEED_ERROR_SUCCESS;
1999 }
2000 
2001 /**
2002   @brief Apply basis evaluation from quadrature points to nodes and sum into target vector
2003 
2004   @param[in]  basis     `CeedBasis` to evaluate
2005   @param[in]  num_elem  The number of elements to apply the basis evaluation to;
2006                           the backend will specify the ordering in @ref CeedElemRestrictionCreate()
2007   @param[in]  t_mode    @ref CEED_TRANSPOSE to apply the transpose, mapping from quadrature points to nodes;
2008                            @ref CEED_NOTRANSPOSE is not valid for `CeedBasisApplyAdd()`
2009   @param[in]  eval_mode @ref CEED_EVAL_NONE to use values directly,
2010                           @ref CEED_EVAL_INTERP to use interpolated values,
2011                           @ref CEED_EVAL_GRAD to use gradients,
2012                           @ref CEED_EVAL_DIV to use divergence,
2013                           @ref CEED_EVAL_CURL to use curl,
2014                           @ref CEED_EVAL_WEIGHT to use quadrature weights
2015   @param[in]  u         Input `CeedVector`
2016   @param[out] v         Output `CeedVector` to sum into
2017 
2018   @return An error code: 0 - success, otherwise - failure
2019 
2020   @ref User
2021 **/
2022 int CeedBasisApplyAdd(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
2023   CeedCheck(t_mode == CEED_TRANSPOSE, CeedBasisReturnCeed(basis), CEED_ERROR_UNSUPPORTED, "CeedBasisApplyAdd only supports CEED_TRANSPOSE");
2024   CeedCall(CeedBasisApplyCheckDims(basis, num_elem, t_mode, eval_mode, u, v));
2025   CeedCheck(basis->ApplyAdd, CeedBasisReturnCeed(basis), CEED_ERROR_UNSUPPORTED, "Backend does not implement CeedBasisApplyAdd");
2026   CeedCall(basis->ApplyAdd(basis, num_elem, t_mode, eval_mode, u, v));
2027   return CEED_ERROR_SUCCESS;
2028 }
2029 
2030 /**
2031   @brief Apply basis evaluation from nodes to arbitrary points
2032 
2033   @param[in]  basis      `CeedBasis` to evaluate
2034   @param[in]  num_elem   The number of elements to apply the basis evaluation to;
2035                           the backend will specify the ordering in @ref CeedElemRestrictionCreate()
2036   @param[in]  num_points Array of the number of points to apply the basis evaluation to in each element, size `num_elem`
2037   @param[in]  t_mode     @ref CEED_NOTRANSPOSE to evaluate from nodes to points;
2038                            @ref CEED_TRANSPOSE to apply the transpose, mapping from points to nodes
2039   @param[in]  eval_mode  @ref CEED_EVAL_INTERP to use interpolated values,
2040                            @ref CEED_EVAL_GRAD to use gradients,
2041                            @ref CEED_EVAL_WEIGHT to use quadrature weights
2042   @param[in]  x_ref      `CeedVector` holding reference coordinates of each point
2043   @param[in]  u          Input `CeedVector`, of length `num_nodes * num_comp` for @ref CEED_NOTRANSPOSE
2044   @param[out] v          Output `CeedVector`, of length `num_points * num_q_comp` for @ref CEED_NOTRANSPOSE with @ref CEED_EVAL_INTERP
2045 
2046   @return An error code: 0 - success, otherwise - failure
2047 
2048   @ref User
2049 **/
2050 int CeedBasisApplyAtPoints(CeedBasis basis, CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode, CeedEvalMode eval_mode,
2051                            CeedVector x_ref, CeedVector u, CeedVector v) {
2052   CeedCall(CeedBasisApplyAtPointsCheckDims(basis, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
2053   if (basis->ApplyAtPoints) {
2054     CeedCall(basis->ApplyAtPoints(basis, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
2055   } else {
2056     CeedCall(CeedBasisApplyAtPoints_Core(basis, false, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
2057   }
2058   return CEED_ERROR_SUCCESS;
2059 }
2060 
2061 /**
2062   @brief Apply basis evaluation from nodes to arbitrary points and sum into target vector
2063 
2064   @param[in]  basis      `CeedBasis` to evaluate
2065   @param[in]  num_elem   The number of elements to apply the basis evaluation to;
2066                           the backend will specify the ordering in @ref CeedElemRestrictionCreate()
2067   @param[in]  num_points Array of the number of points to apply the basis evaluation to in each element, size `num_elem`
2068   @param[in]  t_mode     @ref CEED_NOTRANSPOSE to evaluate from nodes to points;
2069                            @ref CEED_NOTRANSPOSE is not valid for `CeedBasisApplyAddAtPoints()`
2070   @param[in]  eval_mode  @ref CEED_EVAL_INTERP to use interpolated values,
2071                            @ref CEED_EVAL_GRAD to use gradients,
2072                            @ref CEED_EVAL_WEIGHT to use quadrature weights
2073   @param[in]  x_ref      `CeedVector` holding reference coordinates of each point
2074   @param[in]  u          Input `CeedVector`, of length `num_nodes * num_comp` for @ref CEED_NOTRANSPOSE
2075   @param[out] v          Output `CeedVector`, of length `num_points * num_q_comp` for @ref CEED_NOTRANSPOSE with @ref CEED_EVAL_INTERP
2076 
2077   @return An error code: 0 - success, otherwise - failure
2078 
2079   @ref User
2080 **/
2081 int CeedBasisApplyAddAtPoints(CeedBasis basis, CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode, CeedEvalMode eval_mode,
2082                               CeedVector x_ref, CeedVector u, CeedVector v) {
2083   CeedCheck(t_mode == CEED_TRANSPOSE, CeedBasisReturnCeed(basis), CEED_ERROR_UNSUPPORTED, "CeedBasisApplyAddAtPoints only supports CEED_TRANSPOSE");
2084   CeedCall(CeedBasisApplyAtPointsCheckDims(basis, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
2085   if (basis->ApplyAddAtPoints) {
2086     CeedCall(basis->ApplyAddAtPoints(basis, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
2087   } else {
2088     CeedCall(CeedBasisApplyAtPoints_Core(basis, true, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
2089   }
2090   return CEED_ERROR_SUCCESS;
2091 }
2092 
2093 /**
2094   @brief Get the `Ceed` associated with a `CeedBasis`
2095 
2096   @param[in]  basis `CeedBasis`
2097   @param[out] ceed  Variable to store `Ceed`
2098 
2099   @return An error code: 0 - success, otherwise - failure
2100 
2101   @ref Advanced
2102 **/
2103 int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) {
2104   *ceed = NULL;
2105   CeedCall(CeedReferenceCopy(CeedBasisReturnCeed(basis), ceed));
2106   return CEED_ERROR_SUCCESS;
2107 }
2108 
2109 /**
2110   @brief Return the `Ceed` associated with a `CeedBasis`
2111 
2112   @param[in]  basis `CeedBasis`
2113 
2114   @return `Ceed` associated with the `basis`
2115 
2116   @ref Advanced
2117 **/
2118 Ceed CeedBasisReturnCeed(CeedBasis basis) { return basis->ceed; }
2119 
2120 /**
2121   @brief Get dimension for given `CeedBasis`
2122 
2123   @param[in]  basis `CeedBasis`
2124   @param[out] dim   Variable to store dimension of basis
2125 
2126   @return An error code: 0 - success, otherwise - failure
2127 
2128   @ref Advanced
2129 **/
2130 int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) {
2131   *dim = basis->dim;
2132   return CEED_ERROR_SUCCESS;
2133 }
2134 
2135 /**
2136   @brief Get topology for given `CeedBasis`
2137 
2138   @param[in]  basis `CeedBasis`
2139   @param[out] topo  Variable to store topology of basis
2140 
2141   @return An error code: 0 - success, otherwise - failure
2142 
2143   @ref Advanced
2144 **/
2145 int CeedBasisGetTopology(CeedBasis basis, CeedElemTopology *topo) {
2146   *topo = basis->topo;
2147   return CEED_ERROR_SUCCESS;
2148 }
2149 
2150 /**
2151   @brief Get number of components for given `CeedBasis`
2152 
2153   @param[in]  basis    `CeedBasis`
2154   @param[out] num_comp Variable to store number of components
2155 
2156   @return An error code: 0 - success, otherwise - failure
2157 
2158   @ref Advanced
2159 **/
2160 int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *num_comp) {
2161   *num_comp = basis->num_comp;
2162   return CEED_ERROR_SUCCESS;
2163 }
2164 
2165 /**
2166   @brief Get total number of nodes (in `dim` dimensions) of a `CeedBasis`
2167 
2168   @param[in]  basis `CeedBasis`
2169   @param[out] P     Variable to store number of nodes
2170 
2171   @return An error code: 0 - success, otherwise - failure
2172 
2173   @ref Utility
2174 **/
2175 int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) {
2176   *P = basis->P;
2177   return CEED_ERROR_SUCCESS;
2178 }
2179 
2180 /**
2181   @brief Get total number of nodes (in 1 dimension) of a `CeedBasis`
2182 
2183   @param[in]  basis `CeedBasis`
2184   @param[out] P_1d  Variable to store number of nodes
2185 
2186   @return An error code: 0 - success, otherwise - failure
2187 
2188   @ref Advanced
2189 **/
2190 int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P_1d) {
2191   CeedCheck(basis->is_tensor_basis, CeedBasisReturnCeed(basis), CEED_ERROR_MINOR, "Cannot supply P_1d for non-tensor CeedBasis");
2192   *P_1d = basis->P_1d;
2193   return CEED_ERROR_SUCCESS;
2194 }
2195 
2196 /**
2197   @brief Get total number of quadrature points (in `dim` dimensions) of a `CeedBasis`
2198 
2199   @param[in]  basis `CeedBasis`
2200   @param[out] Q     Variable to store number of quadrature points
2201 
2202   @return An error code: 0 - success, otherwise - failure
2203 
2204   @ref Utility
2205 **/
2206 int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) {
2207   *Q = basis->Q;
2208   return CEED_ERROR_SUCCESS;
2209 }
2210 
2211 /**
2212   @brief Get total number of quadrature points (in 1 dimension) of a `CeedBasis`
2213 
2214   @param[in]  basis `CeedBasis`
2215   @param[out] Q_1d  Variable to store number of quadrature points
2216 
2217   @return An error code: 0 - success, otherwise - failure
2218 
2219   @ref Advanced
2220 **/
2221 int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q_1d) {
2222   CeedCheck(basis->is_tensor_basis, CeedBasisReturnCeed(basis), CEED_ERROR_MINOR, "Cannot supply Q_1d for non-tensor CeedBasis");
2223   *Q_1d = basis->Q_1d;
2224   return CEED_ERROR_SUCCESS;
2225 }
2226 
2227 /**
2228   @brief Get reference coordinates of quadrature points (in `dim` dimensions) of a `CeedBasis`
2229 
2230   @param[in]  basis `CeedBasis`
2231   @param[out] q_ref Variable to store reference coordinates of quadrature points
2232 
2233   @return An error code: 0 - success, otherwise - failure
2234 
2235   @ref Advanced
2236 **/
2237 int CeedBasisGetQRef(CeedBasis basis, const CeedScalar **q_ref) {
2238   *q_ref = basis->q_ref_1d;
2239   return CEED_ERROR_SUCCESS;
2240 }
2241 
2242 /**
2243   @brief Get quadrature weights of quadrature points (in `dim` dimensions) of a `CeedBasis`
2244 
2245   @param[in]  basis    `CeedBasis`
2246   @param[out] q_weight Variable to store quadrature weights
2247 
2248   @return An error code: 0 - success, otherwise - failure
2249 
2250   @ref Advanced
2251 **/
2252 int CeedBasisGetQWeights(CeedBasis basis, const CeedScalar **q_weight) {
2253   *q_weight = basis->q_weight_1d;
2254   return CEED_ERROR_SUCCESS;
2255 }
2256 
2257 /**
2258   @brief Get interpolation matrix of a `CeedBasis`
2259 
2260   @param[in]  basis  `CeedBasis`
2261   @param[out] interp Variable to store interpolation matrix
2262 
2263   @return An error code: 0 - success, otherwise - failure
2264 
2265   @ref Advanced
2266 **/
2267 int CeedBasisGetInterp(CeedBasis basis, const CeedScalar **interp) {
2268   if (!basis->interp && basis->is_tensor_basis) {
2269     // Allocate
2270     CeedCall(CeedMalloc(basis->Q * basis->P, &basis->interp));
2271 
2272     // Initialize
2273     for (CeedInt i = 0; i < basis->Q * basis->P; i++) basis->interp[i] = 1.0;
2274 
2275     // Calculate
2276     for (CeedInt d = 0; d < basis->dim; d++) {
2277       for (CeedInt qpt = 0; qpt < basis->Q; qpt++) {
2278         for (CeedInt node = 0; node < basis->P; node++) {
2279           CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d;
2280           CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d;
2281 
2282           basis->interp[qpt * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p];
2283         }
2284       }
2285     }
2286   }
2287   *interp = basis->interp;
2288   return CEED_ERROR_SUCCESS;
2289 }
2290 
2291 /**
2292   @brief Get 1D interpolation matrix of a tensor product `CeedBasis`
2293 
2294   @param[in]  basis     `CeedBasis`
2295   @param[out] interp_1d Variable to store interpolation matrix
2296 
2297   @return An error code: 0 - success, otherwise - failure
2298 
2299   @ref Backend
2300 **/
2301 int CeedBasisGetInterp1D(CeedBasis basis, const CeedScalar **interp_1d) {
2302   bool is_tensor_basis;
2303 
2304   CeedCall(CeedBasisIsTensor(basis, &is_tensor_basis));
2305   CeedCheck(is_tensor_basis, CeedBasisReturnCeed(basis), CEED_ERROR_MINOR, "CeedBasis is not a tensor product CeedBasis");
2306   *interp_1d = basis->interp_1d;
2307   return CEED_ERROR_SUCCESS;
2308 }
2309 
2310 /**
2311   @brief Get gradient matrix of a `CeedBasis`
2312 
2313   @param[in]  basis `CeedBasis`
2314   @param[out] grad  Variable to store gradient matrix
2315 
2316   @return An error code: 0 - success, otherwise - failure
2317 
2318   @ref Advanced
2319 **/
2320 int CeedBasisGetGrad(CeedBasis basis, const CeedScalar **grad) {
2321   if (!basis->grad && basis->is_tensor_basis) {
2322     // Allocate
2323     CeedCall(CeedMalloc(basis->dim * basis->Q * basis->P, &basis->grad));
2324 
2325     // Initialize
2326     for (CeedInt i = 0; i < basis->dim * basis->Q * basis->P; i++) basis->grad[i] = 1.0;
2327 
2328     // Calculate
2329     for (CeedInt d = 0; d < basis->dim; d++) {
2330       for (CeedInt i = 0; i < basis->dim; i++) {
2331         for (CeedInt qpt = 0; qpt < basis->Q; qpt++) {
2332           for (CeedInt node = 0; node < basis->P; node++) {
2333             CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d;
2334             CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d;
2335 
2336             if (i == d) basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->grad_1d[q * basis->P_1d + p];
2337             else basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p];
2338           }
2339         }
2340       }
2341     }
2342   }
2343   *grad = basis->grad;
2344   return CEED_ERROR_SUCCESS;
2345 }
2346 
2347 /**
2348   @brief Get 1D gradient matrix of a tensor product `CeedBasis`
2349 
2350   @param[in]  basis   `CeedBasis`
2351   @param[out] grad_1d Variable to store gradient matrix
2352 
2353   @return An error code: 0 - success, otherwise - failure
2354 
2355   @ref Advanced
2356 **/
2357 int CeedBasisGetGrad1D(CeedBasis basis, const CeedScalar **grad_1d) {
2358   bool is_tensor_basis;
2359 
2360   CeedCall(CeedBasisIsTensor(basis, &is_tensor_basis));
2361   CeedCheck(is_tensor_basis, CeedBasisReturnCeed(basis), CEED_ERROR_MINOR, "CeedBasis is not a tensor product CeedBasis");
2362   *grad_1d = basis->grad_1d;
2363   return CEED_ERROR_SUCCESS;
2364 }
2365 
2366 /**
2367   @brief Get divergence matrix of a `CeedBasis`
2368 
2369   @param[in]  basis `CeedBasis`
2370   @param[out] div   Variable to store divergence matrix
2371 
2372   @return An error code: 0 - success, otherwise - failure
2373 
2374   @ref Advanced
2375 **/
2376 int CeedBasisGetDiv(CeedBasis basis, const CeedScalar **div) {
2377   *div = basis->div;
2378   return CEED_ERROR_SUCCESS;
2379 }
2380 
2381 /**
2382   @brief Get curl matrix of a `CeedBasis`
2383 
2384   @param[in]  basis `CeedBasis`
2385   @param[out] curl  Variable to store curl matrix
2386 
2387   @return An error code: 0 - success, otherwise - failure
2388 
2389   @ref Advanced
2390 **/
2391 int CeedBasisGetCurl(CeedBasis basis, const CeedScalar **curl) {
2392   *curl = basis->curl;
2393   return CEED_ERROR_SUCCESS;
2394 }
2395 
2396 /**
2397   @brief Destroy a @ref  CeedBasis
2398 
2399   @param[in,out] basis `CeedBasis` to destroy
2400 
2401   @return An error code: 0 - success, otherwise - failure
2402 
2403   @ref User
2404 **/
2405 int CeedBasisDestroy(CeedBasis *basis) {
2406   if (!*basis || *basis == CEED_BASIS_NONE || --(*basis)->ref_count > 0) {
2407     *basis = NULL;
2408     return CEED_ERROR_SUCCESS;
2409   }
2410   if ((*basis)->Destroy) CeedCall((*basis)->Destroy(*basis));
2411   CeedCall(CeedTensorContractDestroy(&(*basis)->contract));
2412   CeedCall(CeedFree(&(*basis)->q_ref_1d));
2413   CeedCall(CeedFree(&(*basis)->q_weight_1d));
2414   CeedCall(CeedFree(&(*basis)->interp));
2415   CeedCall(CeedFree(&(*basis)->interp_1d));
2416   CeedCall(CeedFree(&(*basis)->grad));
2417   CeedCall(CeedFree(&(*basis)->grad_1d));
2418   CeedCall(CeedFree(&(*basis)->div));
2419   CeedCall(CeedFree(&(*basis)->curl));
2420   CeedCall(CeedVectorDestroy(&(*basis)->vec_chebyshev));
2421   CeedCall(CeedBasisDestroy(&(*basis)->basis_chebyshev));
2422   CeedCall(CeedDestroy(&(*basis)->ceed));
2423   CeedCall(CeedFree(basis));
2424   return CEED_ERROR_SUCCESS;
2425 }
2426 
2427 /**
2428   @brief Construct a Gauss-Legendre quadrature
2429 
2430   @param[in]  Q           Number of quadrature points (integrates polynomials of degree `2*Q-1` exactly)
2431   @param[out] q_ref_1d    Array of length `Q` to hold the abscissa on `[-1, 1]`
2432   @param[out] q_weight_1d Array of length `Q` to hold the weights
2433 
2434   @return An error code: 0 - success, otherwise - failure
2435 
2436   @ref Utility
2437 **/
2438 int CeedGaussQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) {
2439   CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0 * atan(1.0);
2440 
2441   // Build q_ref_1d, q_weight_1d
2442   for (CeedInt i = 0; i <= Q / 2; i++) {
2443     // Guess
2444     xi = cos(PI * (CeedScalar)(2 * i + 1) / ((CeedScalar)(2 * Q)));
2445     // Pn(xi)
2446     P0 = 1.0;
2447     P1 = xi;
2448     P2 = 0.0;
2449     for (CeedInt j = 2; j <= Q; j++) {
2450       P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
2451       P0 = P1;
2452       P1 = P2;
2453     }
2454     // First Newton Step
2455     dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
2456     xi  = xi - P2 / dP2;
2457     // Newton to convergence
2458     for (CeedInt k = 0; k < 100 && fabs(P2) > 10 * CEED_EPSILON; k++) {
2459       P0 = 1.0;
2460       P1 = xi;
2461       for (CeedInt j = 2; j <= Q; j++) {
2462         P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
2463         P0 = P1;
2464         P1 = P2;
2465       }
2466       dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
2467       xi  = xi - P2 / dP2;
2468     }
2469     // Save xi, wi
2470     wi                     = 2.0 / ((1.0 - xi * xi) * dP2 * dP2);
2471     q_weight_1d[i]         = wi;
2472     q_weight_1d[Q - 1 - i] = wi;
2473     q_ref_1d[i]            = -xi;
2474     q_ref_1d[Q - 1 - i]    = xi;
2475   }
2476   return CEED_ERROR_SUCCESS;
2477 }
2478 
2479 /**
2480   @brief Construct a Gauss-Legendre-Lobatto quadrature
2481 
2482   @param[in]  Q           Number of quadrature points (integrates polynomials of degree `2*Q-3` exactly)
2483   @param[out] q_ref_1d    Array of length `Q` to hold the abscissa on `[-1, 1]`
2484   @param[out] q_weight_1d Array of length `Q` to hold the weights
2485 
2486   @return An error code: 0 - success, otherwise - failure
2487 
2488   @ref Utility
2489 **/
2490 int CeedLobattoQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) {
2491   CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0 * atan(1.0);
2492 
2493   // Build q_ref_1d, q_weight_1d
2494   // Set endpoints
2495   CeedCheck(Q > 1, NULL, CEED_ERROR_DIMENSION, "Cannot create Lobatto quadrature with Q=%" CeedInt_FMT " < 2 points", Q);
2496   wi = 2.0 / ((CeedScalar)(Q * (Q - 1)));
2497   if (q_weight_1d) {
2498     q_weight_1d[0]     = wi;
2499     q_weight_1d[Q - 1] = wi;
2500   }
2501   q_ref_1d[0]     = -1.0;
2502   q_ref_1d[Q - 1] = 1.0;
2503   // Interior
2504   for (CeedInt i = 1; i <= (Q - 1) / 2; i++) {
2505     // Guess
2506     xi = cos(PI * (CeedScalar)(i) / (CeedScalar)(Q - 1));
2507     // Pn(xi)
2508     P0 = 1.0;
2509     P1 = xi;
2510     P2 = 0.0;
2511     for (CeedInt j = 2; j < Q; j++) {
2512       P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
2513       P0 = P1;
2514       P1 = P2;
2515     }
2516     // First Newton step
2517     dP2  = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
2518     d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi);
2519     xi   = xi - dP2 / d2P2;
2520     // Newton to convergence
2521     for (CeedInt k = 0; k < 100 && fabs(dP2) > 10 * CEED_EPSILON; k++) {
2522       P0 = 1.0;
2523       P1 = xi;
2524       for (CeedInt j = 2; j < Q; j++) {
2525         P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
2526         P0 = P1;
2527         P1 = P2;
2528       }
2529       dP2  = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
2530       d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi);
2531       xi   = xi - dP2 / d2P2;
2532     }
2533     // Save xi, wi
2534     wi = 2.0 / (((CeedScalar)(Q * (Q - 1))) * P2 * P2);
2535     if (q_weight_1d) {
2536       q_weight_1d[i]         = wi;
2537       q_weight_1d[Q - 1 - i] = wi;
2538     }
2539     q_ref_1d[i]         = -xi;
2540     q_ref_1d[Q - 1 - i] = xi;
2541   }
2542   return CEED_ERROR_SUCCESS;
2543 }
2544 
2545 /// @}
2546