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