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