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