xref: /libCEED/rust/libceed-sys/c-src/interface/ceed-basis.c (revision 2ab8f9abfe7d30e4e4d04b8d85a29fe1ff83aed3)
1 // Copyright (c) 2017-2022, 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_collocated;
21 /// @endcond
22 
23 /// @addtogroup CeedBasisUser
24 /// @{
25 
26 /// Indicate that the quadrature points are collocated with the nodes
27 const CeedBasis CEED_BASIS_COLLOCATED = &ceed_basis_collocated;
28 
29 /// @}
30 
31 /// ----------------------------------------------------------------------------
32 /// CeedBasis Library Internal Functions
33 /// ----------------------------------------------------------------------------
34 /// @addtogroup CeedBasisDeveloper
35 /// @{
36 
37 /**
38   @brief Compute Householder reflection
39 
40   Computes A = (I - b v v^T) A, where A is an mxn matrix indexed as A[i*row + j*col]
41 
42   @param[in,out] A   Matrix to apply Householder reflection to, in place
43   @param[in]     v   Householder vector
44   @param[in]     b   Scaling factor
45   @param[in]     m   Number of rows in A
46   @param[in]     n   Number of columns in A
47   @param[in]     row Row stride
48   @param[in]     col Col stride
49 
50   @return An error code: 0 - success, otherwise - failure
51 
52   @ref Developer
53 **/
54 static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v, CeedScalar b, CeedInt m, CeedInt n, CeedInt row, CeedInt col) {
55   for (CeedInt j = 0; j < n; j++) {
56     CeedScalar w = A[0 * row + j * col];
57     for (CeedInt i = 1; i < m; i++) w += v[i] * A[i * row + j * col];
58     A[0 * row + j * col] -= b * w;
59     for (CeedInt i = 1; i < m; i++) A[i * row + j * col] -= b * w * v[i];
60   }
61   return CEED_ERROR_SUCCESS;
62 }
63 
64 /**
65   @brief Compute Givens rotation
66 
67   Computes A = G A (or G^T A in transpose mode), where A is an mxn matrix indexed as A[i*n + j*m]
68 
69   @param[in,out] A      Row major matrix to apply Givens rotation to, in place
70   @param[in]     c      Cosine factor
71   @param[in]     s      Sine factor
72   @param[in]     t_mode @ref CEED_NOTRANSPOSE to rotate the basis counter-clockwise, which has the effect of rotating columns of A clockwise;
73                           @ref CEED_TRANSPOSE for the opposite rotation
74   @param[in]     i      First row/column to apply rotation
75   @param[in]     k      Second row/column to apply rotation
76   @param[in]     m      Number of rows in A
77   @param[in]     n      Number of columns in A
78 
79   @return An error code: 0 - success, otherwise - failure
80 
81   @ref Developer
82 **/
83 static int CeedGivensRotation(CeedScalar *A, CeedScalar c, CeedScalar s, CeedTransposeMode t_mode, CeedInt i, CeedInt k, CeedInt m, CeedInt n) {
84   CeedInt stride_j = 1, stride_ik = m, num_its = n;
85   if (t_mode == CEED_NOTRANSPOSE) {
86     stride_j  = n;
87     stride_ik = 1;
88     num_its   = m;
89   }
90 
91   // Apply rotation
92   for (CeedInt j = 0; j < num_its; j++) {
93     CeedScalar tau1 = A[i * stride_ik + j * stride_j], tau2 = A[k * stride_ik + j * stride_j];
94     A[i * stride_ik + j * stride_j] = c * tau1 - s * tau2;
95     A[k * stride_ik + j * stride_j] = s * tau1 + c * tau2;
96   }
97   return CEED_ERROR_SUCCESS;
98 }
99 
100 /**
101   @brief View an array stored in a CeedBasis
102 
103   @param[in] name   Name of array
104   @param[in] fp_fmt Printing format
105   @param[in] m      Number of rows in array
106   @param[in] n      Number of columns in array
107   @param[in] a      Array to be viewed
108   @param[in] stream Stream to view to, e.g., stdout
109 
110   @return An error code: 0 - success, otherwise - failure
111 
112   @ref Developer
113 **/
114 static int CeedScalarView(const char *name, const char *fp_fmt, CeedInt m, CeedInt n, const CeedScalar *a, FILE *stream) {
115   for (CeedInt i = 0; i < m; i++) {
116     if (m > 1) fprintf(stream, "%12s[%" CeedInt_FMT "]:", name, i);
117     else fprintf(stream, "%12s:", name);
118     for (CeedInt j = 0; j < n; j++) fprintf(stream, fp_fmt, fabs(a[i * n + j]) > 1E-14 ? a[i * n + j] : 0);
119     fputs("\n", stream);
120   }
121   return CEED_ERROR_SUCCESS;
122 }
123 
124 /**
125   @brief Create the interpolation and gradient matrices for projection from the nodes of `basis_from` to the nodes of `basis_to`.
126 
127   The interpolation is given by `interp_project = interp_to^+ * interp_from`, where the pseudoinverse `interp_to^+` is given by QR factorization.
128   The gradient is given by `grad_project = interp_to^+ * grad_from`, and is only computed for H^1 spaces otherwise it should not be used.
129 
130   Note: `basis_from` and `basis_to` must have compatible quadrature spaces.
131 
132   @param[in]  basis_from     CeedBasis to project from
133   @param[in]  basis_to       CeedBasis to project to
134   @param[out] interp_project Address of the variable where the newly created interpolation matrix will be stored.
135   @param[out] grad_project   Address of the variable where the newly created gradient matrix will be stored.
136 
137   @return An error code: 0 - success, otherwise - failure
138 
139   @ref Developer
140 **/
141 static int CeedBasisCreateProjectionMatrices(CeedBasis basis_from, CeedBasis basis_to, CeedScalar **interp_project, CeedScalar **grad_project) {
142   Ceed ceed;
143   CeedCall(CeedBasisGetCeed(basis_to, &ceed));
144 
145   // Check for compatible quadrature spaces
146   CeedInt Q_to, Q_from;
147   CeedCall(CeedBasisGetNumQuadraturePoints(basis_to, &Q_to));
148   CeedCall(CeedBasisGetNumQuadraturePoints(basis_from, &Q_from));
149   CeedCheck(Q_to == Q_from, ceed, CEED_ERROR_DIMENSION, "Bases must have compatible quadrature spaces");
150 
151   // Check for matching tensor or non-tensor
152   CeedInt P_to, P_from, Q = Q_to;
153   bool    is_tensor_to, is_tensor_from;
154   CeedCall(CeedBasisIsTensor(basis_to, &is_tensor_to));
155   CeedCall(CeedBasisIsTensor(basis_from, &is_tensor_from));
156   CeedCheck(is_tensor_to == is_tensor_from, ceed, CEED_ERROR_MINOR, "Bases must both be tensor or non-tensor");
157   if (is_tensor_to) {
158     CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_to));
159     CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_from));
160     CeedCall(CeedBasisGetNumQuadraturePoints1D(basis_from, &Q));
161   } else {
162     CeedCall(CeedBasisGetNumNodes(basis_to, &P_to));
163     CeedCall(CeedBasisGetNumNodes(basis_from, &P_from));
164   }
165 
166   // Check for matching FE space
167   CeedFESpace fe_space_to, fe_space_from;
168   CeedCall(CeedBasisGetFESpace(basis_to, &fe_space_to));
169   CeedCall(CeedBasisGetFESpace(basis_from, &fe_space_from));
170   CeedCheck(fe_space_to == fe_space_from, ceed, CEED_ERROR_MINOR, "Bases must both be the same FE space type");
171 
172   // Get source matrices
173   CeedInt           dim, q_comp = 1;
174   const CeedScalar *interp_to_source = NULL, *interp_from_source = NULL;
175   CeedScalar       *interp_to, *interp_from, *tau;
176   CeedCall(CeedBasisGetDimension(basis_to, &dim));
177   if (is_tensor_to) {
178     CeedCall(CeedBasisGetInterp1D(basis_to, &interp_to_source));
179     CeedCall(CeedBasisGetInterp1D(basis_from, &interp_from_source));
180   } else {
181     CeedCall(CeedBasisGetNumQuadratureComponents(basis_from, CEED_EVAL_INTERP, &q_comp));
182     CeedCall(CeedBasisGetInterp(basis_to, &interp_to_source));
183     CeedCall(CeedBasisGetInterp(basis_from, &interp_from_source));
184   }
185   CeedCall(CeedMalloc(Q * P_from * q_comp, &interp_from));
186   CeedCall(CeedMalloc(Q * P_to * q_comp, &interp_to));
187   CeedCall(CeedCalloc(P_to * P_from, interp_project));
188   CeedCall(CeedMalloc(Q * q_comp, &tau));
189 
190   // `grad_project = interp_to^+ * grad_from` is computed for the H^1 space case so the
191   // projection basis will have a gradient operation (allocated even if not H^1 for the
192   // basis construction later on)
193   const CeedScalar *grad_from_source = NULL;
194   if (fe_space_to == CEED_FE_SPACE_H1) {
195     if (is_tensor_to) {
196       CeedCall(CeedBasisGetGrad1D(basis_from, &grad_from_source));
197     } else {
198       CeedCall(CeedBasisGetGrad(basis_from, &grad_from_source));
199     }
200   }
201   CeedCall(CeedCalloc(P_to * P_from * (is_tensor_to ? 1 : dim), grad_project));
202 
203   // QR Factorization, interp_to = Q R
204   memcpy(interp_to, interp_to_source, Q * P_to * q_comp * sizeof(interp_to_source[0]));
205   CeedCall(CeedQRFactorization(ceed, interp_to, tau, Q * q_comp, P_to));
206 
207   // Build matrices
208   CeedInt     num_matrices = 1 + (fe_space_to == CEED_FE_SPACE_H1) * (is_tensor_to ? 1 : dim);
209   CeedScalar *input_from[num_matrices], *output_project[num_matrices];
210   input_from[0]     = (CeedScalar *)interp_from_source;
211   output_project[0] = *interp_project;
212   for (CeedInt m = 1; m < num_matrices; m++) {
213     input_from[m]     = (CeedScalar *)&grad_from_source[(m - 1) * Q * P_from];
214     output_project[m] = &((*grad_project)[(m - 1) * P_to * P_from]);
215   }
216   for (CeedInt m = 0; m < num_matrices; m++) {
217     // Apply Q^T, interp_from = Q^T interp_from
218     memcpy(interp_from, input_from[m], Q * P_from * q_comp * sizeof(input_from[m][0]));
219     CeedCall(CeedHouseholderApplyQ(interp_from, interp_to, tau, CEED_TRANSPOSE, Q * q_comp, P_from, P_to, P_from, 1));
220 
221     // Apply Rinv, output_project = Rinv interp_from
222     for (CeedInt j = 0; j < P_from; j++) {  // Column j
223       output_project[m][j + P_from * (P_to - 1)] = interp_from[j + P_from * (P_to - 1)] / interp_to[P_to * P_to - 1];
224       for (CeedInt i = P_to - 2; i >= 0; i--) {  // Row i
225         output_project[m][j + P_from * i] = interp_from[j + P_from * i];
226         for (CeedInt k = i + 1; k < P_to; k++) {
227           output_project[m][j + P_from * i] -= interp_to[k + P_to * i] * output_project[m][j + P_from * k];
228         }
229         output_project[m][j + P_from * i] /= interp_to[i + P_to * i];
230       }
231     }
232   }
233 
234   // Cleanup
235   CeedCall(CeedFree(&tau));
236   CeedCall(CeedFree(&interp_to));
237   CeedCall(CeedFree(&interp_from));
238 
239   return CEED_ERROR_SUCCESS;
240 }
241 
242 /// @}
243 
244 /// ----------------------------------------------------------------------------
245 /// Ceed Backend API
246 /// ----------------------------------------------------------------------------
247 /// @addtogroup CeedBasisBackend
248 /// @{
249 
250 /**
251   @brief Return collocated grad matrix
252 
253   @param[in]  basis         CeedBasis
254   @param[out] collo_grad_1d Row-major (Q_1d * Q_1d) matrix expressing derivatives of basis functions at quadrature points
255 
256   @return An error code: 0 - success, otherwise - failure
257 
258   @ref Backend
259 **/
260 int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *collo_grad_1d) {
261   int         i, j, k;
262   Ceed        ceed;
263   CeedInt     P_1d = (basis)->P_1d, Q_1d = (basis)->Q_1d;
264   CeedScalar *interp_1d, *grad_1d, *tau;
265 
266   CeedCall(CeedMalloc(Q_1d * P_1d, &interp_1d));
267   CeedCall(CeedMalloc(Q_1d * P_1d, &grad_1d));
268   CeedCall(CeedMalloc(Q_1d, &tau));
269   memcpy(interp_1d, (basis)->interp_1d, Q_1d * P_1d * sizeof(basis)->interp_1d[0]);
270   memcpy(grad_1d, (basis)->grad_1d, Q_1d * P_1d * sizeof(basis)->interp_1d[0]);
271 
272   // QR Factorization, interp_1d = Q R
273   CeedCall(CeedBasisGetCeed(basis, &ceed));
274   CeedCall(CeedQRFactorization(ceed, interp_1d, tau, Q_1d, P_1d));
275   // Note: This function is for backend use, so all errors are terminal and we do not need to clean up memory on failure.
276 
277   // Apply Rinv, collo_grad_1d = grad_1d Rinv
278   for (i = 0; i < Q_1d; i++) {  // Row i
279     collo_grad_1d[Q_1d * i] = grad_1d[P_1d * i] / interp_1d[0];
280     for (j = 1; j < P_1d; j++) {  // Column j
281       collo_grad_1d[j + Q_1d * i] = grad_1d[j + P_1d * i];
282       for (k = 0; k < j; k++) collo_grad_1d[j + Q_1d * i] -= interp_1d[j + P_1d * k] * collo_grad_1d[k + Q_1d * i];
283       collo_grad_1d[j + Q_1d * i] /= interp_1d[j + P_1d * j];
284     }
285     for (j = P_1d; j < Q_1d; j++) collo_grad_1d[j + Q_1d * i] = 0;
286   }
287 
288   // Apply Q^T, collo_grad_1d = collo_grad_1d Q^T
289   CeedCall(CeedHouseholderApplyQ(collo_grad_1d, interp_1d, tau, CEED_NOTRANSPOSE, Q_1d, Q_1d, P_1d, 1, Q_1d));
290 
291   CeedCall(CeedFree(&interp_1d));
292   CeedCall(CeedFree(&grad_1d));
293   CeedCall(CeedFree(&tau));
294   return CEED_ERROR_SUCCESS;
295 }
296 
297 /**
298   @brief Get tensor status for given CeedBasis
299 
300   @param[in]  basis     CeedBasis
301   @param[out] is_tensor Variable to store tensor status
302 
303   @return An error code: 0 - success, otherwise - failure
304 
305   @ref Backend
306 **/
307 int CeedBasisIsTensor(CeedBasis basis, bool *is_tensor) {
308   *is_tensor = basis->is_tensor_basis;
309   return CEED_ERROR_SUCCESS;
310 }
311 
312 /**
313   @brief Get backend data of a CeedBasis
314 
315   @param[in]  basis CeedBasis
316   @param[out] data  Variable to store data
317 
318   @return An error code: 0 - success, otherwise - failure
319 
320   @ref Backend
321 **/
322 int CeedBasisGetData(CeedBasis basis, void *data) {
323   *(void **)data = basis->data;
324   return CEED_ERROR_SUCCESS;
325 }
326 
327 /**
328   @brief Set backend data of a CeedBasis
329 
330   @param[in,out] basis  CeedBasis
331   @param[in]     data   Data to set
332 
333   @return An error code: 0 - success, otherwise - failure
334 
335   @ref Backend
336 **/
337 int CeedBasisSetData(CeedBasis basis, void *data) {
338   basis->data = data;
339   return CEED_ERROR_SUCCESS;
340 }
341 
342 /**
343   @brief Increment the reference counter for a CeedBasis
344 
345   @param[in,out] basis Basis to increment the reference counter
346 
347   @return An error code: 0 - success, otherwise - failure
348 
349   @ref Backend
350 **/
351 int CeedBasisReference(CeedBasis basis) {
352   basis->ref_count++;
353   return CEED_ERROR_SUCCESS;
354 }
355 
356 /**
357   @brief Get number of Q-vector components for given CeedBasis
358 
359   @param[in]  basis  CeedBasis
360   @param[in]  eval_mode \ref CEED_EVAL_INTERP to use interpolated values,
361                           \ref CEED_EVAL_GRAD to use gradients,
362                           \ref CEED_EVAL_DIV to use divergence,
363                           \ref CEED_EVAL_CURL to use curl.
364   @param[out] q_comp Variable to store number of Q-vector components of basis
365 
366   @return An error code: 0 - success, otherwise - failure
367 
368   @ref Backend
369 **/
370 int CeedBasisGetNumQuadratureComponents(CeedBasis basis, CeedEvalMode eval_mode, CeedInt *q_comp) {
371   switch (eval_mode) {
372     case CEED_EVAL_INTERP:
373       *q_comp = (basis->fe_space == CEED_FE_SPACE_H1) ? 1 : basis->dim;
374       break;
375     case CEED_EVAL_GRAD:
376       *q_comp = basis->dim;
377       break;
378     case CEED_EVAL_DIV:
379       *q_comp = 1;
380       break;
381     case CEED_EVAL_CURL:
382       *q_comp = (basis->dim < 3) ? 1 : basis->dim;
383       break;
384     case CEED_EVAL_NONE:
385     case CEED_EVAL_WEIGHT:
386       *q_comp = 1;
387       break;
388   }
389   return CEED_ERROR_SUCCESS;
390 }
391 
392 /**
393   @brief Estimate number of FLOPs required to apply CeedBasis in t_mode and eval_mode
394 
395   @param[in]  basis     Basis to estimate FLOPs for
396   @param[in]  t_mode    Apply basis or transpose
397   @param[in]  eval_mode Basis evaluation mode
398   @param[out] flops     Address of variable to hold FLOPs estimate
399 
400   @ref Backend
401 **/
402 int CeedBasisGetFlopsEstimate(CeedBasis basis, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedSize *flops) {
403   bool is_tensor;
404 
405   CeedCall(CeedBasisIsTensor(basis, &is_tensor));
406   if (is_tensor) {
407     CeedInt dim, num_comp, P_1d, Q_1d;
408     CeedCall(CeedBasisGetDimension(basis, &dim));
409     CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
410     CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d));
411     CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
412     if (t_mode == CEED_TRANSPOSE) {
413       P_1d = Q_1d;
414       Q_1d = P_1d;
415     }
416     CeedInt tensor_flops = 0, pre = num_comp * CeedIntPow(P_1d, dim - 1), post = 1;
417     for (CeedInt d = 0; d < dim; d++) {
418       tensor_flops += 2 * pre * P_1d * post * Q_1d;
419       pre /= P_1d;
420       post *= Q_1d;
421     }
422     switch (eval_mode) {
423       case CEED_EVAL_NONE:
424         *flops = 0;
425         break;
426       case CEED_EVAL_INTERP:
427         *flops = tensor_flops;
428         break;
429       case CEED_EVAL_GRAD:
430         *flops = tensor_flops * 2;
431         break;
432       case CEED_EVAL_DIV:
433       case CEED_EVAL_CURL:
434         // LCOV_EXCL_START
435         return CeedError(basis->ceed, CEED_ERROR_INCOMPATIBLE, "Tensor basis evaluation for %s not supported", CeedEvalModes[eval_mode]);
436         break;
437       // LCOV_EXCL_STOP
438       case CEED_EVAL_WEIGHT:
439         *flops = dim * CeedIntPow(Q_1d, dim);
440         break;
441     }
442   } else {
443     CeedInt dim, num_comp, q_comp, num_nodes, num_qpts;
444     CeedCall(CeedBasisGetDimension(basis, &dim));
445     CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
446     CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp));
447     CeedCall(CeedBasisGetNumNodes(basis, &num_nodes));
448     CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts));
449     switch (eval_mode) {
450       case CEED_EVAL_NONE:
451         *flops = 0;
452         break;
453       case CEED_EVAL_INTERP:
454       case CEED_EVAL_GRAD:
455       case CEED_EVAL_DIV:
456       case CEED_EVAL_CURL:
457         *flops = num_nodes * num_qpts * num_comp * q_comp;
458         break;
459       case CEED_EVAL_WEIGHT:
460         *flops = 0;
461         break;
462     }
463   }
464 
465   return CEED_ERROR_SUCCESS;
466 }
467 
468 /**
469   @brief Get CeedFESpace for a CeedBasis
470 
471   @param[in]  basis     CeedBasis
472   @param[out] fe_space  Variable to store CeedFESpace
473 
474   @return An error code: 0 - success, otherwise - failure
475 
476   @ref Backend
477 **/
478 int CeedBasisGetFESpace(CeedBasis basis, CeedFESpace *fe_space) {
479   *fe_space = basis->fe_space;
480   return CEED_ERROR_SUCCESS;
481 }
482 
483 /**
484   @brief Get dimension for given CeedElemTopology
485 
486   @param[in]  topo CeedElemTopology
487   @param[out] dim  Variable to store dimension of topology
488 
489   @return An error code: 0 - success, otherwise - failure
490 
491   @ref Backend
492 **/
493 int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) {
494   *dim = (CeedInt)topo >> 16;
495   return CEED_ERROR_SUCCESS;
496 }
497 
498 /**
499   @brief Get CeedTensorContract of a CeedBasis
500 
501   @param[in]  basis     CeedBasis
502   @param[out] contract  Variable to store CeedTensorContract
503 
504   @return An error code: 0 - success, otherwise - failure
505 
506   @ref Backend
507 **/
508 int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) {
509   *contract = basis->contract;
510   return CEED_ERROR_SUCCESS;
511 }
512 
513 /**
514   @brief Set CeedTensorContract of a CeedBasis
515 
516   @param[in,out] basis    CeedBasis
517   @param[in]     contract CeedTensorContract to set
518 
519   @return An error code: 0 - success, otherwise - failure
520 
521   @ref Backend
522 **/
523 int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract contract) {
524   basis->contract = contract;
525   CeedCall(CeedTensorContractReference(contract));
526   return CEED_ERROR_SUCCESS;
527 }
528 
529 /**
530   @brief Return a reference implementation of matrix multiplication C = A B.
531 
532   Note: This is a reference implementation for CPU CeedScalar pointers that is not intended for high performance.
533 
534   @param[in]  ceed  Ceed context for error handling
535   @param[in]  mat_A Row-major matrix A
536   @param[in]  mat_B Row-major matrix B
537   @param[out] mat_C Row-major output matrix C
538   @param[in]  m     Number of rows of C
539   @param[in]  n     Number of columns of C
540   @param[in]  kk    Number of columns of A/rows of B
541 
542   @return An error code: 0 - success, otherwise - failure
543 
544   @ref Utility
545 **/
546 int CeedMatrixMatrixMultiply(Ceed ceed, const CeedScalar *mat_A, const CeedScalar *mat_B, CeedScalar *mat_C, CeedInt m, CeedInt n, CeedInt kk) {
547   for (CeedInt i = 0; i < m; i++) {
548     for (CeedInt j = 0; j < n; j++) {
549       CeedScalar sum = 0;
550       for (CeedInt k = 0; k < kk; k++) sum += mat_A[k + i * kk] * mat_B[j + k * n];
551       mat_C[j + i * n] = sum;
552     }
553   }
554   return CEED_ERROR_SUCCESS;
555 }
556 
557 /**
558   @brief Return QR Factorization of a matrix
559 
560   @param[in]     ceed Ceed context for error handling
561   @param[in,out] mat  Row-major matrix to be factorized in place
562   @param[in,out] tau  Vector of length m of scaling factors
563   @param[in]     m    Number of rows
564   @param[in]     n    Number of columns
565 
566   @return An error code: 0 - success, otherwise - failure
567 
568   @ref Utility
569 **/
570 int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau, CeedInt m, CeedInt n) {
571   CeedScalar v[m];
572 
573   // Check matrix shape
574   CeedCheck(n <= m, ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute QR factorization with n > m");
575 
576   for (CeedInt i = 0; i < n; i++) {
577     if (i >= m - 1) {  // last row of matrix, no reflection needed
578       tau[i] = 0.;
579       break;
580     }
581     // Calculate Householder vector, magnitude
582     CeedScalar sigma = 0.0;
583     v[i]             = mat[i + n * i];
584     for (CeedInt j = i + 1; j < m; j++) {
585       v[j] = mat[i + n * j];
586       sigma += v[j] * v[j];
587     }
588     CeedScalar norm = sqrt(v[i] * v[i] + sigma);  // norm of v[i:m]
589     CeedScalar R_ii = -copysign(norm, v[i]);
590     v[i] -= R_ii;
591     // norm of v[i:m] after modification above and scaling below
592     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
593     //   tau = 2 / (norm*norm)
594     tau[i] = 2 * v[i] * v[i] / (v[i] * v[i] + sigma);
595     for (CeedInt j = i + 1; j < m; j++) v[j] /= v[i];
596 
597     // Apply Householder reflector to lower right panel
598     CeedHouseholderReflect(&mat[i * n + i + 1], &v[i], tau[i], m - i, n - i - 1, n, 1);
599     // Save v
600     mat[i + n * i] = R_ii;
601     for (CeedInt j = i + 1; j < m; j++) mat[i + n * j] = v[j];
602   }
603   return CEED_ERROR_SUCCESS;
604 }
605 
606 /**
607   @brief Apply Householder Q matrix
608 
609   Compute mat_A = mat_Q mat_A, where mat_Q is mxm and mat_A is mxn.
610 
611   @param[in,out] mat_A  Matrix to apply Householder Q to, in place
612   @param[in]     mat_Q  Householder Q matrix
613   @param[in]     tau    Householder scaling factors
614   @param[in]     t_mode Transpose mode for application
615   @param[in]     m      Number of rows in A
616   @param[in]     n      Number of columns in A
617   @param[in]     k      Number of elementary reflectors in Q, k<m
618   @param[in]     row    Row stride in A
619   @param[in]     col    Col stride in A
620 
621   @return An error code: 0 - success, otherwise - failure
622 
623   @ref Utility
624 **/
625 int CeedHouseholderApplyQ(CeedScalar *mat_A, const CeedScalar *mat_Q, const CeedScalar *tau, CeedTransposeMode t_mode, CeedInt m, CeedInt n,
626                           CeedInt k, CeedInt row, CeedInt col) {
627   CeedScalar *v;
628   CeedCall(CeedMalloc(m, &v));
629   for (CeedInt ii = 0; ii < k; ii++) {
630     CeedInt i = t_mode == CEED_TRANSPOSE ? ii : k - 1 - ii;
631     for (CeedInt j = i + 1; j < m; j++) v[j] = mat_Q[j * k + i];
632     // Apply Householder reflector (I - tau v v^T) collo_grad_1d^T
633     CeedCall(CeedHouseholderReflect(&mat_A[i * row], &v[i], tau[i], m - i, n, row, col));
634   }
635   CeedCall(CeedFree(&v));
636   return CEED_ERROR_SUCCESS;
637 }
638 
639 /**
640   @brief Return symmetric Schur decomposition of the symmetric matrix mat via symmetric QR factorization
641 
642   @param[in]     ceed   Ceed context for error handling
643   @param[in,out] mat    Row-major matrix to be factorized in place
644   @param[out]    lambda Vector of length n of eigenvalues
645   @param[in]     n      Number of rows/columns
646 
647   @return An error code: 0 - success, otherwise - failure
648 
649   @ref Utility
650 **/
651 CeedPragmaOptimizeOff
652 int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat, CeedScalar *lambda, CeedInt n) {
653   // Check bounds for clang-tidy
654   CeedCheck(n > 1, ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute symmetric Schur decomposition of scalars");
655 
656   CeedScalar v[n - 1], tau[n - 1], mat_T[n * n];
657 
658   // Copy mat to mat_T and set mat to I
659   memcpy(mat_T, mat, n * n * sizeof(mat[0]));
660   for (CeedInt i = 0; i < n; i++) {
661     for (CeedInt j = 0; j < n; j++) mat[j + n * i] = (i == j) ? 1 : 0;
662   }
663 
664   // Reduce to tridiagonal
665   for (CeedInt i = 0; i < n - 1; i++) {
666     // Calculate Householder vector, magnitude
667     CeedScalar sigma = 0.0;
668     v[i]             = mat_T[i + n * (i + 1)];
669     for (CeedInt j = i + 1; j < n - 1; j++) {
670       v[j] = mat_T[i + n * (j + 1)];
671       sigma += v[j] * v[j];
672     }
673     CeedScalar norm = sqrt(v[i] * v[i] + sigma);  // norm of v[i:n-1]
674     CeedScalar R_ii = -copysign(norm, v[i]);
675     v[i] -= R_ii;
676     // norm of v[i:m] after modification above and scaling below
677     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
678     //   tau = 2 / (norm*norm)
679     tau[i] = i == n - 2 ? 2 : 2 * v[i] * v[i] / (v[i] * v[i] + sigma);
680     for (CeedInt j = i + 1; j < n - 1; j++) v[j] /= v[i];
681 
682     // Update sub and super diagonal
683     for (CeedInt j = i + 2; j < n; j++) {
684       mat_T[i + n * j] = 0;
685       mat_T[j + n * i] = 0;
686     }
687     // Apply symmetric Householder reflector to lower right panel
688     CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1);
689     CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), 1, n);
690 
691     // Save v
692     mat_T[i + n * (i + 1)] = R_ii;
693     mat_T[(i + 1) + n * i] = R_ii;
694     for (CeedInt j = i + 1; j < n - 1; j++) {
695       mat_T[i + n * (j + 1)] = v[j];
696     }
697   }
698   // Backwards accumulation of Q
699   for (CeedInt i = n - 2; i >= 0; i--) {
700     if (tau[i] > 0.0) {
701       v[i] = 1;
702       for (CeedInt j = i + 1; j < n - 1; j++) {
703         v[j]                   = mat_T[i + n * (j + 1)];
704         mat_T[i + n * (j + 1)] = 0;
705       }
706       CeedHouseholderReflect(&mat[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1);
707     }
708   }
709 
710   // Reduce sub and super diagonal
711   CeedInt    p = 0, q = 0, itr = 0, max_itr = n * n * n * n;
712   CeedScalar tol = CEED_EPSILON;
713 
714   while (itr < max_itr) {
715     // Update p, q, size of reduced portions of diagonal
716     p = 0;
717     q = 0;
718     for (CeedInt i = n - 2; i >= 0; i--) {
719       if (fabs(mat_T[i + n * (i + 1)]) < tol) q += 1;
720       else break;
721     }
722     for (CeedInt i = 0; i < n - q - 1; i++) {
723       if (fabs(mat_T[i + n * (i + 1)]) < tol) p += 1;
724       else break;
725     }
726     if (q == n - 1) break;  // Finished reducing
727 
728     // Reduce tridiagonal portion
729     CeedScalar t_nn = mat_T[(n - 1 - q) + n * (n - 1 - q)], t_nnm1 = mat_T[(n - 2 - q) + n * (n - 1 - q)];
730     CeedScalar d  = (mat_T[(n - 2 - q) + n * (n - 2 - q)] - t_nn) / 2;
731     CeedScalar mu = t_nn - t_nnm1 * t_nnm1 / (d + copysign(sqrt(d * d + t_nnm1 * t_nnm1), d));
732     CeedScalar x  = mat_T[p + n * p] - mu;
733     CeedScalar z  = mat_T[p + n * (p + 1)];
734     for (CeedInt k = p; k < n - q - 1; k++) {
735       // Compute Givens rotation
736       CeedScalar c = 1, s = 0;
737       if (fabs(z) > tol) {
738         if (fabs(z) > fabs(x)) {
739           CeedScalar tau = -x / z;
740           s = 1 / sqrt(1 + tau * tau), c = s * tau;
741         } else {
742           CeedScalar tau = -z / x;
743           c = 1 / sqrt(1 + tau * tau), s = c * tau;
744         }
745       }
746 
747       // Apply Givens rotation to T
748       CeedGivensRotation(mat_T, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n);
749       CeedGivensRotation(mat_T, c, s, CEED_TRANSPOSE, k, k + 1, n, n);
750 
751       // Apply Givens rotation to Q
752       CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n);
753 
754       // Update x, z
755       if (k < n - q - 2) {
756         x = mat_T[k + n * (k + 1)];
757         z = mat_T[k + n * (k + 2)];
758       }
759     }
760     itr++;
761   }
762 
763   // Save eigenvalues
764   for (CeedInt i = 0; i < n; i++) lambda[i] = mat_T[i + n * i];
765 
766   // Check convergence
767   CeedCheck(itr < max_itr || q > n, ceed, CEED_ERROR_MINOR, "Symmetric QR failed to converge");
768   return CEED_ERROR_SUCCESS;
769 }
770 CeedPragmaOptimizeOn
771 
772 /**
773   @brief Return Simultaneous Diagonalization of two matrices.
774 
775   This solves the generalized eigenvalue problem A x = lambda B x, where A and B are symmetric and B is positive definite.
776   We generate the matrix X and vector Lambda such that X^T A X = Lambda and X^T B X = I.
777   This is equivalent to the LAPACK routine 'sygv' with TYPE = 1.
778 
779   @param[in]  ceed   Ceed context for error handling
780   @param[in]  mat_A  Row-major matrix to be factorized with eigenvalues
781   @param[in]  mat_B  Row-major matrix to be factorized to identity
782   @param[out] mat_X  Row-major orthogonal matrix
783   @param[out] lambda Vector of length n of generalized eigenvalues
784   @param[in]  n      Number of rows/columns
785 
786   @return An error code: 0 - success, otherwise - failure
787 
788   @ref Utility
789 **/
790 CeedPragmaOptimizeOff
791 int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *mat_A, CeedScalar *mat_B, CeedScalar *mat_X, CeedScalar *lambda, CeedInt n) {
792   CeedScalar *mat_C, *mat_G, *vec_D;
793   CeedCall(CeedCalloc(n * n, &mat_C));
794   CeedCall(CeedCalloc(n * n, &mat_G));
795   CeedCall(CeedCalloc(n, &vec_D));
796 
797   // Compute B = G D G^T
798   memcpy(mat_G, mat_B, n * n * sizeof(mat_B[0]));
799   CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_G, vec_D, n));
800 
801   // Sort eigenvalues
802   for (CeedInt i = n - 1; i >= 0; i--) {
803     for (CeedInt j = 0; j < i; j++) {
804       if (fabs(vec_D[j]) > fabs(vec_D[j + 1])) {
805         CeedScalar temp;
806         temp         = vec_D[j];
807         vec_D[j]     = vec_D[j + 1];
808         vec_D[j + 1] = temp;
809         for (CeedInt k = 0; k < n; k++) {
810           temp                 = mat_G[k * n + j];
811           mat_G[k * n + j]     = mat_G[k * n + j + 1];
812           mat_G[k * n + j + 1] = temp;
813         }
814       }
815     }
816   }
817 
818   // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T
819   //           = D^-1/2 G^T A G D^-1/2
820   // -- D = D^-1/2
821   for (CeedInt i = 0; i < n; i++) vec_D[i] = 1. / sqrt(vec_D[i]);
822   // -- G = G D^-1/2
823   // -- C = D^-1/2 G^T
824   for (CeedInt i = 0; i < n; i++) {
825     for (CeedInt j = 0; j < n; j++) {
826       mat_G[i * n + j] *= vec_D[j];
827       mat_C[j * n + i] = mat_G[i * n + j];
828     }
829   }
830   // -- X = (D^-1/2 G^T) A
831   CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_C, (const CeedScalar *)mat_A, mat_X, n, n, n));
832   // -- C = (D^-1/2 G^T A) (G D^-1/2)
833   CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_X, (const CeedScalar *)mat_G, mat_C, n, n, n));
834 
835   // Compute Q^T C Q = lambda
836   CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_C, lambda, n));
837 
838   // Sort eigenvalues
839   for (CeedInt i = n - 1; i >= 0; i--) {
840     for (CeedInt j = 0; j < i; j++) {
841       if (fabs(lambda[j]) > fabs(lambda[j + 1])) {
842         CeedScalar temp;
843         temp          = lambda[j];
844         lambda[j]     = lambda[j + 1];
845         lambda[j + 1] = temp;
846         for (CeedInt k = 0; k < n; k++) {
847           temp                 = mat_C[k * n + j];
848           mat_C[k * n + j]     = mat_C[k * n + j + 1];
849           mat_C[k * n + j + 1] = temp;
850         }
851       }
852     }
853   }
854 
855   // Set X = (G D^1/2)^-T Q
856   //       = G D^-1/2 Q
857   CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_G, (const CeedScalar *)mat_C, mat_X, n, n, n));
858 
859   // Cleanup
860   CeedCall(CeedFree(&mat_C));
861   CeedCall(CeedFree(&mat_G));
862   CeedCall(CeedFree(&vec_D));
863   return CEED_ERROR_SUCCESS;
864 }
865 CeedPragmaOptimizeOn
866 
867 /// @}
868 
869 /// ----------------------------------------------------------------------------
870 /// CeedBasis Public API
871 /// ----------------------------------------------------------------------------
872 /// @addtogroup CeedBasisUser
873 /// @{
874 
875 /**
876   @brief Create a tensor-product basis for H^1 discretizations
877 
878   @param[in]  ceed        Ceed object where the CeedBasis will be created
879   @param[in]  dim         Topological dimension
880   @param[in]  num_comp    Number of field components (1 for scalar fields)
881   @param[in]  P_1d        Number of nodes in one dimension
882   @param[in]  Q_1d        Number of quadrature points in one dimension
883   @param[in]  interp_1d   Row-major (Q_1d * P_1d) matrix expressing the values of nodal basis functions at quadrature points
884   @param[in]  grad_1d     Row-major (Q_1d * P_1d) matrix expressing derivatives of nodal basis functions at quadrature points
885   @param[in]  q_ref_1d    Array of length Q_1d holding the locations of quadrature points on the 1D reference element [-1, 1]
886   @param[in]  q_weight_1d Array of length Q_1d holding the quadrature weights on the reference element
887   @param[out] basis       Address of the variable where the newly created CeedBasis will be stored.
888 
889   @return An error code: 0 - success, otherwise - failure
890 
891   @ref User
892 **/
893 int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d,
894                             const CeedScalar *grad_1d, const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis *basis) {
895   if (!ceed->BasisCreateTensorH1) {
896     Ceed delegate;
897 
898     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
899     CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisCreateTensorH1");
900     CeedCall(CeedBasisCreateTensorH1(delegate, dim, num_comp, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis));
901     return CEED_ERROR_SUCCESS;
902   }
903 
904   CeedCheck(dim > 0, ceed, CEED_ERROR_DIMENSION, "Basis dimension must be a positive value");
905   CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component");
906   CeedCheck(P_1d > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node");
907   CeedCheck(Q_1d > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point");
908 
909   CeedElemTopology topo = dim == 1 ? CEED_TOPOLOGY_LINE : dim == 2 ? CEED_TOPOLOGY_QUAD : CEED_TOPOLOGY_HEX;
910 
911   CeedCall(CeedCalloc(1, basis));
912   CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed));
913   (*basis)->ref_count       = 1;
914   (*basis)->is_tensor_basis = true;
915   (*basis)->dim             = dim;
916   (*basis)->topo            = topo;
917   (*basis)->num_comp        = num_comp;
918   (*basis)->P_1d            = P_1d;
919   (*basis)->Q_1d            = Q_1d;
920   (*basis)->P               = CeedIntPow(P_1d, dim);
921   (*basis)->Q               = CeedIntPow(Q_1d, dim);
922   (*basis)->fe_space        = CEED_FE_SPACE_H1;
923   CeedCall(CeedCalloc(Q_1d, &(*basis)->q_ref_1d));
924   CeedCall(CeedCalloc(Q_1d, &(*basis)->q_weight_1d));
925   if (q_ref_1d) memcpy((*basis)->q_ref_1d, q_ref_1d, Q_1d * sizeof(q_ref_1d[0]));
926   if (q_weight_1d) memcpy((*basis)->q_weight_1d, q_weight_1d, Q_1d * sizeof(q_weight_1d[0]));
927   CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->interp_1d));
928   CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->grad_1d));
929   if (interp_1d) memcpy((*basis)->interp_1d, interp_1d, Q_1d * P_1d * sizeof(interp_1d[0]));
930   if (grad_1d) memcpy((*basis)->grad_1d, grad_1d, Q_1d * P_1d * sizeof(grad_1d[0]));
931   CeedCall(ceed->BasisCreateTensorH1(dim, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, *basis));
932   return CEED_ERROR_SUCCESS;
933 }
934 
935 /**
936   @brief Create a tensor-product Lagrange basis
937 
938   @param[in]  ceed      Ceed object where the CeedBasis will be created
939   @param[in]  dim       Topological dimension of element
940   @param[in]  num_comp  Number of field components (1 for scalar fields)
941   @param[in]  P         Number of Gauss-Lobatto nodes in one dimension.
942                           The polynomial degree of the resulting Q_k element is k=P-1.
943   @param[in]  Q         Number of quadrature points in one dimension.
944   @param[in]  quad_mode Distribution of the Q quadrature points (affects order of accuracy for the quadrature)
945   @param[out] basis     Address of the variable where the newly created CeedBasis will be stored.
946 
947   @return An error code: 0 - success, otherwise - failure
948 
949   @ref User
950 **/
951 int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P, CeedInt Q, CeedQuadMode quad_mode, CeedBasis *basis) {
952   // Allocate
953   int        ierr = CEED_ERROR_SUCCESS, i, j, k;
954   CeedScalar c1, c2, c3, c4, dx, *nodes, *interp_1d, *grad_1d, *q_ref_1d, *q_weight_1d;
955 
956   CeedCheck(dim > 0, ceed, CEED_ERROR_DIMENSION, "Basis dimension must be a positive value");
957   CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component");
958   CeedCheck(P > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node");
959   CeedCheck(Q > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point");
960 
961   // Get Nodes and Weights
962   CeedCall(CeedCalloc(P * Q, &interp_1d));
963   CeedCall(CeedCalloc(P * Q, &grad_1d));
964   CeedCall(CeedCalloc(P, &nodes));
965   CeedCall(CeedCalloc(Q, &q_ref_1d));
966   CeedCall(CeedCalloc(Q, &q_weight_1d));
967   if (CeedLobattoQuadrature(P, nodes, NULL) != CEED_ERROR_SUCCESS) goto cleanup;
968   switch (quad_mode) {
969     case CEED_GAUSS:
970       ierr = CeedGaussQuadrature(Q, q_ref_1d, q_weight_1d);
971       break;
972     case CEED_GAUSS_LOBATTO:
973       ierr = CeedLobattoQuadrature(Q, q_ref_1d, q_weight_1d);
974       break;
975   }
976   if (ierr != CEED_ERROR_SUCCESS) goto cleanup;
977 
978   // Build B, D matrix
979   // Fornberg, 1998
980   for (i = 0; i < Q; i++) {
981     c1                   = 1.0;
982     c3                   = nodes[0] - q_ref_1d[i];
983     interp_1d[i * P + 0] = 1.0;
984     for (j = 1; j < P; j++) {
985       c2 = 1.0;
986       c4 = c3;
987       c3 = nodes[j] - q_ref_1d[i];
988       for (k = 0; k < j; k++) {
989         dx = nodes[j] - nodes[k];
990         c2 *= dx;
991         if (k == j - 1) {
992           grad_1d[i * P + j]   = c1 * (interp_1d[i * P + k] - c4 * grad_1d[i * P + k]) / c2;
993           interp_1d[i * P + j] = -c1 * c4 * interp_1d[i * P + k] / c2;
994         }
995         grad_1d[i * P + k]   = (c3 * grad_1d[i * P + k] - interp_1d[i * P + k]) / dx;
996         interp_1d[i * P + k] = c3 * interp_1d[i * P + k] / dx;
997       }
998       c1 = c2;
999     }
1000   }
1001   // Pass to CeedBasisCreateTensorH1
1002   CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P, Q, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis));
1003 cleanup:
1004   CeedCall(CeedFree(&interp_1d));
1005   CeedCall(CeedFree(&grad_1d));
1006   CeedCall(CeedFree(&nodes));
1007   CeedCall(CeedFree(&q_ref_1d));
1008   CeedCall(CeedFree(&q_weight_1d));
1009   return CEED_ERROR_SUCCESS;
1010 }
1011 
1012 /**
1013   @brief Create a non tensor-product basis for H^1 discretizations
1014 
1015   @param[in]  ceed      Ceed object where the CeedBasis will be created
1016   @param[in]  topo      Topology of element, e.g. hypercube, simplex, ect
1017   @param[in]  num_comp  Number of field components (1 for scalar fields)
1018   @param[in]  num_nodes Total number of nodes
1019   @param[in]  num_qpts  Total number of quadrature points
1020   @param[in]  interp    Row-major (num_qpts * num_nodes) matrix expressing the values of nodal basis functions at quadrature points
1021   @param[in]  grad      Row-major (dim * num_qpts * num_nodes) matrix expressing derivatives of nodal basis functions at quadrature points
1022   @param[in]  q_ref     Array of length num_qpts * dim holding the locations of quadrature points on the reference element
1023   @param[in]  q_weight  Array of length num_qpts holding the quadrature weights on the reference element
1024   @param[out] basis     Address of the variable where the newly created CeedBasis will be stored.
1025 
1026   @return An error code: 0 - success, otherwise - failure
1027 
1028   @ref User
1029 **/
1030 int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
1031                       const CeedScalar *grad, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) {
1032   CeedInt P = num_nodes, Q = num_qpts, dim = 0;
1033 
1034   if (!ceed->BasisCreateH1) {
1035     Ceed delegate;
1036 
1037     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
1038     CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisCreateH1");
1039     CeedCall(CeedBasisCreateH1(delegate, topo, num_comp, num_nodes, num_qpts, interp, grad, q_ref, q_weight, basis));
1040     return CEED_ERROR_SUCCESS;
1041   }
1042 
1043   CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component");
1044   CeedCheck(num_nodes > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node");
1045   CeedCheck(num_qpts > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point");
1046 
1047   CeedCall(CeedBasisGetTopologyDimension(topo, &dim));
1048 
1049   CeedCall(CeedCalloc(1, basis));
1050   CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed));
1051   (*basis)->ref_count       = 1;
1052   (*basis)->is_tensor_basis = false;
1053   (*basis)->dim             = dim;
1054   (*basis)->topo            = topo;
1055   (*basis)->num_comp        = num_comp;
1056   (*basis)->P               = P;
1057   (*basis)->Q               = Q;
1058   (*basis)->fe_space        = CEED_FE_SPACE_H1;
1059   CeedCall(CeedCalloc(Q * dim, &(*basis)->q_ref_1d));
1060   CeedCall(CeedCalloc(Q, &(*basis)->q_weight_1d));
1061   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0]));
1062   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0]));
1063   CeedCall(CeedCalloc(Q * P, &(*basis)->interp));
1064   CeedCall(CeedCalloc(dim * Q * P, &(*basis)->grad));
1065   if (interp) memcpy((*basis)->interp, interp, Q * P * sizeof(interp[0]));
1066   if (grad) memcpy((*basis)->grad, grad, dim * Q * P * sizeof(grad[0]));
1067   CeedCall(ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, q_ref, q_weight, *basis));
1068   return CEED_ERROR_SUCCESS;
1069 }
1070 
1071 /**
1072   @brief Create a non tensor-product basis for \f$H(\mathrm{div})\f$ discretizations
1073 
1074   @param[in]  ceed      Ceed object where the CeedBasis will be created
1075   @param[in]  topo      Topology of element (`CEED_TOPOLOGY_QUAD`, `CEED_TOPOLOGY_PRISM`, etc.), dimension of which is used in some array sizes below
1076   @param[in]  num_comp  Number of components (usually 1 for vectors in H(div) bases)
1077   @param[in]  num_nodes Total number of nodes (dofs per element)
1078   @param[in]  num_qpts  Total number of quadrature points
1079   @param[in]  interp    Row-major (dim * num_qpts * num_nodes) matrix expressing the values of basis functions at quadrature points
1080   @param[in]  div       Row-major (num_qpts * num_nodes) matrix expressing divergence of basis functions at quadrature points
1081   @param[in]  q_ref     Array of length num_qpts * dim holding the locations of quadrature points on the reference element
1082   @param[in]  q_weight  Array of length num_qpts holding the quadrature weights on the reference element
1083   @param[out] basis     Address of the variable where the newly created CeedBasis will be stored.
1084 
1085   @return An error code: 0 - success, otherwise - failure
1086 
1087   @ref User
1088 **/
1089 int CeedBasisCreateHdiv(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
1090                         const CeedScalar *div, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) {
1091   CeedInt Q = num_qpts, P = num_nodes, dim = 0;
1092 
1093   if (!ceed->BasisCreateHdiv) {
1094     Ceed delegate;
1095 
1096     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
1097     CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateHdiv");
1098     CeedCall(CeedBasisCreateHdiv(delegate, topo, num_comp, num_nodes, num_qpts, interp, div, q_ref, q_weight, basis));
1099     return CEED_ERROR_SUCCESS;
1100   }
1101 
1102   CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component");
1103   CeedCheck(num_nodes > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node");
1104   CeedCheck(num_qpts > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point");
1105 
1106   CeedCall(CeedBasisGetTopologyDimension(topo, &dim));
1107 
1108   CeedCall(CeedCalloc(1, basis));
1109   CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed));
1110   (*basis)->ref_count       = 1;
1111   (*basis)->is_tensor_basis = false;
1112   (*basis)->dim             = dim;
1113   (*basis)->topo            = topo;
1114   (*basis)->num_comp        = num_comp;
1115   (*basis)->P               = P;
1116   (*basis)->Q               = Q;
1117   (*basis)->fe_space        = CEED_FE_SPACE_HDIV;
1118   CeedCall(CeedMalloc(Q * dim, &(*basis)->q_ref_1d));
1119   CeedCall(CeedMalloc(Q, &(*basis)->q_weight_1d));
1120   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0]));
1121   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0]));
1122   CeedCall(CeedMalloc(dim * Q * P, &(*basis)->interp));
1123   CeedCall(CeedMalloc(Q * P, &(*basis)->div));
1124   if (interp) memcpy((*basis)->interp, interp, dim * Q * P * sizeof(interp[0]));
1125   if (div) memcpy((*basis)->div, div, Q * P * sizeof(div[0]));
1126   CeedCall(ceed->BasisCreateHdiv(topo, dim, P, Q, interp, div, q_ref, q_weight, *basis));
1127   return CEED_ERROR_SUCCESS;
1128 }
1129 
1130 /**
1131   @brief Create a non tensor-product basis for \f$H(\mathrm{curl})\f$ discretizations
1132 
1133   @param[in]  ceed      Ceed object where the CeedBasis will be created
1134   @param[in]  topo      Topology of element (`CEED_TOPOLOGY_QUAD`, `CEED_TOPOLOGY_PRISM`, etc.), dimension of which is used in some array sizes below
1135   @param[in]  num_comp  Number of components (usually 1 for vectors in H(curl) bases)
1136   @param[in]  num_nodes Total number of nodes (dofs per element)
1137   @param[in]  num_qpts  Total number of quadrature points
1138   @param[in]  interp    Row-major (dim * num_qpts * num_nodes) matrix expressing the values of basis functions at quadrature points
1139   @param[in]  curl      Row-major (curl_comp * num_qpts * num_nodes, curl_comp = 1 if dim < 3 else dim) matrix expressing curl of basis functions at
1140 quadrature points
1141   @param[in]  q_ref     Array of length num_qpts * dim holding the locations of quadrature points on the reference element
1142   @param[in]  q_weight  Array of length num_qpts holding the quadrature weights on the reference element
1143   @param[out] basis     Address of the variable where the newly created CeedBasis will be stored.
1144 
1145   @return An error code: 0 - success, otherwise - failure
1146 
1147   @ref User
1148 **/
1149 int CeedBasisCreateHcurl(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
1150                          const CeedScalar *curl, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) {
1151   CeedInt Q = num_qpts, P = num_nodes, dim = 0, curl_comp = 0;
1152 
1153   if (!ceed->BasisCreateHdiv) {
1154     Ceed delegate;
1155 
1156     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
1157     CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateHcurl");
1158     CeedCall(CeedBasisCreateHcurl(delegate, topo, num_comp, num_nodes, num_qpts, interp, curl, q_ref, q_weight, basis));
1159     return CEED_ERROR_SUCCESS;
1160   }
1161 
1162   CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component");
1163   CeedCheck(num_nodes > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node");
1164   CeedCheck(num_qpts > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point");
1165 
1166   CeedCall(CeedBasisGetTopologyDimension(topo, &dim));
1167   curl_comp = (dim < 3) ? 1 : dim;
1168 
1169   CeedCall(CeedCalloc(1, basis));
1170   CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed));
1171   (*basis)->ref_count       = 1;
1172   (*basis)->is_tensor_basis = false;
1173   (*basis)->dim             = dim;
1174   (*basis)->topo            = topo;
1175   (*basis)->num_comp        = num_comp;
1176   (*basis)->P               = P;
1177   (*basis)->Q               = Q;
1178   (*basis)->fe_space        = CEED_FE_SPACE_HCURL;
1179   CeedCall(CeedMalloc(Q * dim, &(*basis)->q_ref_1d));
1180   CeedCall(CeedMalloc(Q, &(*basis)->q_weight_1d));
1181   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0]));
1182   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0]));
1183   CeedCall(CeedMalloc(dim * Q * P, &(*basis)->interp));
1184   CeedCall(CeedMalloc(curl_comp * Q * P, &(*basis)->curl));
1185   if (interp) memcpy((*basis)->interp, interp, dim * Q * P * sizeof(interp[0]));
1186   if (curl) memcpy((*basis)->curl, curl, curl_comp * Q * P * sizeof(curl[0]));
1187   CeedCall(ceed->BasisCreateHcurl(topo, dim, P, Q, interp, curl, q_ref, q_weight, *basis));
1188   return CEED_ERROR_SUCCESS;
1189 }
1190 
1191 /**
1192   @brief Create a CeedBasis for projection from the nodes of `basis_from` to the nodes of `basis_to`.
1193 
1194   Only `CEED_EVAL_INTERP` will be valid for the new basis, `basis_project`.
1195   For H^1 spaces, `CEED_EVAL_GRAD` will also be valid.
1196   The interpolation is given by `interp_project = interp_to^+ * interp_from`, where the pseudoinverse `interp_to^+` is given by QR
1197 factorization.
1198   The gradient (for the H^1 case) is given by `grad_project = interp_to^+ * grad_from`.
1199 
1200   Note: `basis_from` and `basis_to` must have compatible quadrature spaces.
1201 
1202   Note: `basis_project` will have the same number of components as `basis_from`, regardless of the number of components that `basis_to` has.
1203         If `basis_from` has 3 components and `basis_to` has 5 components, then `basis_project` will have 3 components.
1204 
1205   @param[in]  basis_from    CeedBasis to prolong from
1206   @param[in]  basis_to      CeedBasis to prolong to
1207   @param[out] basis_project Address of the variable where the newly created CeedBasis will be stored.
1208 
1209   @return An error code: 0 - success, otherwise - failure
1210 
1211   @ref User
1212 **/
1213 int CeedBasisCreateProjection(CeedBasis basis_from, CeedBasis basis_to, CeedBasis *basis_project) {
1214   Ceed ceed;
1215   CeedCall(CeedBasisGetCeed(basis_to, &ceed));
1216 
1217   // Create projection matrix
1218   CeedScalar *interp_project, *grad_project;
1219   CeedCall(CeedBasisCreateProjectionMatrices(basis_from, basis_to, &interp_project, &grad_project));
1220 
1221   // Build basis
1222   bool        is_tensor;
1223   CeedInt     dim, num_comp;
1224   CeedScalar *q_ref, *q_weight;
1225   CeedCall(CeedBasisIsTensor(basis_to, &is_tensor));
1226   CeedCall(CeedBasisGetDimension(basis_to, &dim));
1227   CeedCall(CeedBasisGetNumComponents(basis_from, &num_comp));
1228   if (is_tensor) {
1229     CeedInt P_1d_to, P_1d_from;
1230     CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_1d_from));
1231     CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_1d_to));
1232     CeedCall(CeedCalloc(P_1d_to, &q_ref));
1233     CeedCall(CeedCalloc(P_1d_to, &q_weight));
1234     CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P_1d_from, P_1d_to, interp_project, grad_project, q_ref, q_weight, basis_project));
1235   } else {
1236     // Even if basis_to and basis_from are not H1, the resulting basis is H1 for interpolation to work
1237     CeedElemTopology topo;
1238     CeedCall(CeedBasisGetTopology(basis_to, &topo));
1239     CeedInt num_nodes_to, num_nodes_from;
1240     CeedCall(CeedBasisGetNumNodes(basis_from, &num_nodes_from));
1241     CeedCall(CeedBasisGetNumNodes(basis_to, &num_nodes_to));
1242     CeedCall(CeedCalloc(num_nodes_to * dim, &q_ref));
1243     CeedCall(CeedCalloc(num_nodes_to, &q_weight));
1244     CeedCall(CeedBasisCreateH1(ceed, topo, num_comp, num_nodes_from, num_nodes_to, interp_project, grad_project, q_ref, q_weight, basis_project));
1245   }
1246 
1247   // Cleanup
1248   CeedCall(CeedFree(&interp_project));
1249   CeedCall(CeedFree(&grad_project));
1250   CeedCall(CeedFree(&q_ref));
1251   CeedCall(CeedFree(&q_weight));
1252 
1253   return CEED_ERROR_SUCCESS;
1254 }
1255 
1256 /**
1257   @brief Copy the pointer to a CeedBasis.
1258 
1259   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.
1260         This CeedBasis will be destroyed if `basis_copy` is the only reference to this CeedBasis.
1261 
1262   @param[in]     basis      CeedBasis to copy reference to
1263   @param[in,out] basis_copy Variable to store copied reference
1264 
1265   @return An error code: 0 - success, otherwise - failure
1266 
1267   @ref User
1268 **/
1269 int CeedBasisReferenceCopy(CeedBasis basis, CeedBasis *basis_copy) {
1270   if (basis != CEED_BASIS_COLLOCATED) CeedCall(CeedBasisReference(basis));
1271   CeedCall(CeedBasisDestroy(basis_copy));
1272   *basis_copy = basis;
1273   return CEED_ERROR_SUCCESS;
1274 }
1275 
1276 /**
1277   @brief View a CeedBasis
1278 
1279   @param[in] basis  CeedBasis to view
1280   @param[in] stream Stream to view to, e.g., stdout
1281 
1282   @return An error code: 0 - success, otherwise - failure
1283 
1284   @ref User
1285 **/
1286 int CeedBasisView(CeedBasis basis, FILE *stream) {
1287   CeedElemTopology topo     = basis->topo;
1288   CeedFESpace      fe_space = basis->fe_space;
1289   CeedInt          q_comp   = 0;
1290 
1291   // Print FE space and element topology of the basis
1292   if (basis->is_tensor_basis) {
1293     fprintf(stream, "CeedBasis (%s on a %s element): dim=%" CeedInt_FMT " P=%" CeedInt_FMT " Q=%" CeedInt_FMT "\n", CeedFESpaces[fe_space],
1294             CeedElemTopologies[topo], basis->dim, basis->P_1d, basis->Q_1d);
1295   } else {
1296     fprintf(stream, "CeedBasis (%s on a %s element): dim=%" CeedInt_FMT " P=%" CeedInt_FMT " Q=%" CeedInt_FMT "\n", CeedFESpaces[fe_space],
1297             CeedElemTopologies[topo], basis->dim, basis->P, basis->Q);
1298   }
1299   // Print quadrature data, interpolation/gradient/divergence/curl of the basis
1300   if (basis->is_tensor_basis) {  // tensor basis
1301     CeedCall(CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_ref_1d, stream));
1302     CeedCall(CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_weight_1d, stream));
1303     CeedCall(CeedScalarView("interp1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, basis->interp_1d, stream));
1304     CeedCall(CeedScalarView("grad1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, basis->grad_1d, stream));
1305   } else {  // non-tensor basis
1306     CeedCall(CeedScalarView("qref", "\t% 12.8f", 1, basis->Q * basis->dim, basis->q_ref_1d, stream));
1307     CeedCall(CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->q_weight_1d, stream));
1308     CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp));
1309     CeedCall(CeedScalarView("interp", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->interp, stream));
1310     if (basis->grad) {
1311       CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp));
1312       CeedCall(CeedScalarView("grad", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->grad, stream));
1313     }
1314     if (basis->div) {
1315       CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_DIV, &q_comp));
1316       CeedCall(CeedScalarView("div", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->div, stream));
1317     }
1318     if (basis->curl) {
1319       CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_CURL, &q_comp));
1320       CeedCall(CeedScalarView("curl", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->curl, stream));
1321     }
1322   }
1323   return CEED_ERROR_SUCCESS;
1324 }
1325 
1326 /**
1327   @brief Apply basis evaluation from nodes to quadrature points or vice versa
1328 
1329   @param[in]  basis      CeedBasis to evaluate
1330   @param[in]  num_elem   The number of elements to apply the basis evaluation to;
1331                            the backend will specify the ordering in CeedElemRestrictionCreateBlocked()
1332   @param[in]  t_mode    \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature points;
1333                           \ref CEED_TRANSPOSE to apply the transpose, mapping from quadrature points to nodes
1334   @param[in]  eval_mode \ref CEED_EVAL_NONE to use values directly,
1335                           \ref CEED_EVAL_INTERP to use interpolated values,
1336                           \ref CEED_EVAL_GRAD to use gradients,
1337                           \ref CEED_EVAL_DIV to use divergence,
1338                           \ref CEED_EVAL_CURL to use curl,
1339                           \ref CEED_EVAL_WEIGHT to use quadrature weights.
1340   @param[in]  u        Input CeedVector
1341   @param[out] v        Output CeedVector
1342 
1343   @return An error code: 0 - success, otherwise - failure
1344 
1345   @ref User
1346 **/
1347 int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
1348   CeedSize u_length = 0, v_length;
1349   CeedInt  dim, num_comp, q_comp, num_nodes, num_qpts;
1350   CeedCall(CeedBasisGetDimension(basis, &dim));
1351   CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
1352   CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp));
1353   CeedCall(CeedBasisGetNumNodes(basis, &num_nodes));
1354   CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts));
1355   CeedCall(CeedVectorGetLength(v, &v_length));
1356   if (u) {
1357     CeedCall(CeedVectorGetLength(u, &u_length));
1358   }
1359 
1360   CeedCheck(basis->Apply, basis->ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisApply");
1361 
1362   // Check compatibility of topological and geometrical dimensions
1363   CeedCheck((t_mode == CEED_TRANSPOSE && v_length % num_nodes == 0 && u_length % num_qpts == 0) ||
1364                 (t_mode == CEED_NOTRANSPOSE && u_length % num_nodes == 0 && v_length % num_qpts == 0),
1365             basis->ceed, CEED_ERROR_DIMENSION, "Length of input/output vectors incompatible with basis dimensions");
1366 
1367   // Check vector lengths to prevent out of bounds issues
1368   bool good_dims = true;
1369   switch (eval_mode) {
1370     case CEED_EVAL_NONE:
1371     case CEED_EVAL_INTERP:
1372     case CEED_EVAL_GRAD:
1373     case CEED_EVAL_DIV:
1374     case CEED_EVAL_CURL:
1375       good_dims =
1376           ((t_mode == CEED_TRANSPOSE && u_length >= num_elem * num_comp * num_qpts * q_comp && v_length >= num_elem * num_comp * num_nodes) ||
1377            (t_mode == CEED_NOTRANSPOSE && v_length >= num_elem * num_qpts * num_comp * q_comp && u_length >= num_elem * num_comp * num_nodes));
1378       break;
1379     case CEED_EVAL_WEIGHT:
1380       good_dims = v_length >= num_elem * num_qpts;
1381       break;
1382   }
1383   CeedCheck(good_dims, basis->ceed, CEED_ERROR_DIMENSION, "Input/output vectors too short for basis and evaluation mode");
1384 
1385   CeedCall(basis->Apply(basis, num_elem, t_mode, eval_mode, u, v));
1386   return CEED_ERROR_SUCCESS;
1387 }
1388 
1389 /**
1390   @brief Get Ceed associated with a CeedBasis
1391 
1392   @param[in]  basis CeedBasis
1393   @param[out] ceed  Variable to store Ceed
1394 
1395   @return An error code: 0 - success, otherwise - failure
1396 
1397   @ref Advanced
1398 **/
1399 int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) {
1400   *ceed = basis->ceed;
1401   return CEED_ERROR_SUCCESS;
1402 }
1403 
1404 /**
1405   @brief Get dimension for given CeedBasis
1406 
1407   @param[in]  basis CeedBasis
1408   @param[out] dim   Variable to store dimension of basis
1409 
1410   @return An error code: 0 - success, otherwise - failure
1411 
1412   @ref Advanced
1413 **/
1414 int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) {
1415   *dim = basis->dim;
1416   return CEED_ERROR_SUCCESS;
1417 }
1418 
1419 /**
1420   @brief Get topology for given CeedBasis
1421 
1422   @param[in]  basis CeedBasis
1423   @param[out] topo  Variable to store topology of basis
1424 
1425   @return An error code: 0 - success, otherwise - failure
1426 
1427   @ref Advanced
1428 **/
1429 int CeedBasisGetTopology(CeedBasis basis, CeedElemTopology *topo) {
1430   *topo = basis->topo;
1431   return CEED_ERROR_SUCCESS;
1432 }
1433 
1434 /**
1435   @brief Get number of components for given CeedBasis
1436 
1437   @param[in]  basis    CeedBasis
1438   @param[out] num_comp Variable to store number of components of basis
1439 
1440   @return An error code: 0 - success, otherwise - failure
1441 
1442   @ref Advanced
1443 **/
1444 int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *num_comp) {
1445   *num_comp = basis->num_comp;
1446   return CEED_ERROR_SUCCESS;
1447 }
1448 
1449 /**
1450   @brief Get total number of nodes (in dim dimensions) of a CeedBasis
1451 
1452   @param[in]  basis CeedBasis
1453   @param[out] P     Variable to store number of nodes
1454 
1455   @return An error code: 0 - success, otherwise - failure
1456 
1457   @ref Utility
1458 **/
1459 int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) {
1460   *P = basis->P;
1461   return CEED_ERROR_SUCCESS;
1462 }
1463 
1464 /**
1465   @brief Get total number of nodes (in 1 dimension) of a CeedBasis
1466 
1467   @param[in]  basis CeedBasis
1468   @param[out] P_1d  Variable to store number of nodes
1469 
1470   @return An error code: 0 - success, otherwise - failure
1471 
1472   @ref Advanced
1473 **/
1474 int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P_1d) {
1475   CeedCheck(basis->is_tensor_basis, basis->ceed, CEED_ERROR_MINOR, "Cannot supply P_1d for non-tensor basis");
1476   *P_1d = basis->P_1d;
1477   return CEED_ERROR_SUCCESS;
1478 }
1479 
1480 /**
1481   @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis
1482 
1483   @param[in]  basis CeedBasis
1484   @param[out] Q     Variable to store number of quadrature points
1485 
1486   @return An error code: 0 - success, otherwise - failure
1487 
1488   @ref Utility
1489 **/
1490 int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) {
1491   *Q = basis->Q;
1492   return CEED_ERROR_SUCCESS;
1493 }
1494 
1495 /**
1496   @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis
1497 
1498   @param[in]  basis CeedBasis
1499   @param[out] Q_1d  Variable to store number of quadrature points
1500 
1501   @return An error code: 0 - success, otherwise - failure
1502 
1503   @ref Advanced
1504 **/
1505 int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q_1d) {
1506   CeedCheck(basis->is_tensor_basis, basis->ceed, CEED_ERROR_MINOR, "Cannot supply Q_1d for non-tensor basis");
1507   *Q_1d = basis->Q_1d;
1508   return CEED_ERROR_SUCCESS;
1509 }
1510 
1511 /**
1512   @brief Get reference coordinates of quadrature points (in dim dimensions) of a CeedBasis
1513 
1514   @param[in]  basis CeedBasis
1515   @param[out] q_ref Variable to store reference coordinates of quadrature points
1516 
1517   @return An error code: 0 - success, otherwise - failure
1518 
1519   @ref Advanced
1520 **/
1521 int CeedBasisGetQRef(CeedBasis basis, const CeedScalar **q_ref) {
1522   *q_ref = basis->q_ref_1d;
1523   return CEED_ERROR_SUCCESS;
1524 }
1525 
1526 /**
1527   @brief Get quadrature weights of quadrature points (in dim dimensions) of a CeedBasis
1528 
1529   @param[in]  basis    CeedBasis
1530   @param[out] q_weight Variable to store quadrature weights
1531 
1532   @return An error code: 0 - success, otherwise - failure
1533 
1534   @ref Advanced
1535 **/
1536 int CeedBasisGetQWeights(CeedBasis basis, const CeedScalar **q_weight) {
1537   *q_weight = basis->q_weight_1d;
1538   return CEED_ERROR_SUCCESS;
1539 }
1540 
1541 /**
1542   @brief Get interpolation matrix of a CeedBasis
1543 
1544   @param[in]  basis  CeedBasis
1545   @param[out] interp Variable to store interpolation matrix
1546 
1547   @return An error code: 0 - success, otherwise - failure
1548 
1549   @ref Advanced
1550 **/
1551 int CeedBasisGetInterp(CeedBasis basis, const CeedScalar **interp) {
1552   if (!basis->interp && basis->is_tensor_basis) {
1553     // Allocate
1554     CeedCall(CeedMalloc(basis->Q * basis->P, &basis->interp));
1555 
1556     // Initialize
1557     for (CeedInt i = 0; i < basis->Q * basis->P; i++) basis->interp[i] = 1.0;
1558 
1559     // Calculate
1560     for (CeedInt d = 0; d < basis->dim; d++) {
1561       for (CeedInt qpt = 0; qpt < basis->Q; qpt++) {
1562         for (CeedInt node = 0; node < basis->P; node++) {
1563           CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d;
1564           CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d;
1565           basis->interp[qpt * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p];
1566         }
1567       }
1568     }
1569   }
1570   *interp = basis->interp;
1571   return CEED_ERROR_SUCCESS;
1572 }
1573 
1574 /**
1575   @brief Get 1D interpolation matrix of a tensor product CeedBasis
1576 
1577   @param[in]  basis     CeedBasis
1578   @param[out] interp_1d Variable to store interpolation matrix
1579 
1580   @return An error code: 0 - success, otherwise - failure
1581 
1582   @ref Backend
1583 **/
1584 int CeedBasisGetInterp1D(CeedBasis basis, const CeedScalar **interp_1d) {
1585   CeedCheck(basis->is_tensor_basis, basis->ceed, CEED_ERROR_MINOR, "CeedBasis is not a tensor product basis.");
1586   *interp_1d = basis->interp_1d;
1587   return CEED_ERROR_SUCCESS;
1588 }
1589 
1590 /**
1591   @brief Get gradient matrix of a CeedBasis
1592 
1593   @param[in]  basis CeedBasis
1594   @param[out] grad  Variable to store gradient matrix
1595 
1596   @return An error code: 0 - success, otherwise - failure
1597 
1598   @ref Advanced
1599 **/
1600 int CeedBasisGetGrad(CeedBasis basis, const CeedScalar **grad) {
1601   if (!basis->grad && basis->is_tensor_basis) {
1602     // Allocate
1603     CeedCall(CeedMalloc(basis->dim * basis->Q * basis->P, &basis->grad));
1604 
1605     // Initialize
1606     for (CeedInt i = 0; i < basis->dim * basis->Q * basis->P; i++) basis->grad[i] = 1.0;
1607 
1608     // Calculate
1609     for (CeedInt d = 0; d < basis->dim; d++) {
1610       for (CeedInt i = 0; i < basis->dim; i++) {
1611         for (CeedInt qpt = 0; qpt < basis->Q; qpt++) {
1612           for (CeedInt node = 0; node < basis->P; node++) {
1613             CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d;
1614             CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d;
1615             if (i == d) basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->grad_1d[q * basis->P_1d + p];
1616             else basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p];
1617           }
1618         }
1619       }
1620     }
1621   }
1622   *grad = basis->grad;
1623   return CEED_ERROR_SUCCESS;
1624 }
1625 
1626 /**
1627   @brief Get 1D gradient matrix of a tensor product CeedBasis
1628 
1629   @param[in]  basis   CeedBasis
1630   @param[out] grad_1d Variable to store gradient matrix
1631 
1632   @return An error code: 0 - success, otherwise - failure
1633 
1634   @ref Advanced
1635 **/
1636 int CeedBasisGetGrad1D(CeedBasis basis, const CeedScalar **grad_1d) {
1637   CeedCheck(basis->is_tensor_basis, basis->ceed, CEED_ERROR_MINOR, "CeedBasis is not a tensor product basis.");
1638   *grad_1d = basis->grad_1d;
1639   return CEED_ERROR_SUCCESS;
1640 }
1641 
1642 /**
1643   @brief Get divergence matrix of a CeedBasis
1644 
1645   @param[in]  basis CeedBasis
1646   @param[out] div   Variable to store divergence matrix
1647 
1648   @return An error code: 0 - success, otherwise - failure
1649 
1650   @ref Advanced
1651 **/
1652 int CeedBasisGetDiv(CeedBasis basis, const CeedScalar **div) {
1653   CeedCheck(basis->div, basis->ceed, CEED_ERROR_MINOR, "CeedBasis does not have divergence matrix.");
1654   *div = basis->div;
1655   return CEED_ERROR_SUCCESS;
1656 }
1657 
1658 /**
1659   @brief Get curl matrix of a CeedBasis
1660 
1661   @param[in]  basis CeedBasis
1662   @param[out] curl  Variable to store curl matrix
1663 
1664   @return An error code: 0 - success, otherwise - failure
1665 
1666   @ref Advanced
1667 **/
1668 int CeedBasisGetCurl(CeedBasis basis, const CeedScalar **curl) {
1669   CeedCheck(basis->curl, basis->ceed, CEED_ERROR_MINOR, "CeedBasis does not have curl matrix.");
1670   *curl = basis->curl;
1671   return CEED_ERROR_SUCCESS;
1672 }
1673 
1674 /**
1675   @brief Destroy a CeedBasis
1676 
1677   @param[in,out] basis CeedBasis to destroy
1678 
1679   @return An error code: 0 - success, otherwise - failure
1680 
1681   @ref User
1682 **/
1683 int CeedBasisDestroy(CeedBasis *basis) {
1684   if (!*basis || *basis == CEED_BASIS_COLLOCATED || --(*basis)->ref_count > 0) {
1685     *basis = NULL;
1686     return CEED_ERROR_SUCCESS;
1687   }
1688   if ((*basis)->Destroy) CeedCall((*basis)->Destroy(*basis));
1689   if ((*basis)->contract) CeedCall(CeedTensorContractDestroy(&(*basis)->contract));
1690   CeedCall(CeedFree(&(*basis)->q_ref_1d));
1691   CeedCall(CeedFree(&(*basis)->q_weight_1d));
1692   CeedCall(CeedFree(&(*basis)->interp));
1693   CeedCall(CeedFree(&(*basis)->interp_1d));
1694   CeedCall(CeedFree(&(*basis)->grad));
1695   CeedCall(CeedFree(&(*basis)->grad_1d));
1696   CeedCall(CeedFree(&(*basis)->div));
1697   CeedCall(CeedFree(&(*basis)->curl));
1698   CeedCall(CeedDestroy(&(*basis)->ceed));
1699   CeedCall(CeedFree(basis));
1700   return CEED_ERROR_SUCCESS;
1701 }
1702 
1703 /**
1704   @brief Construct a Gauss-Legendre quadrature
1705 
1706   @param[in]  Q           Number of quadrature points (integrates polynomials of degree 2*Q-1 exactly)
1707   @param[out] q_ref_1d    Array of length Q to hold the abscissa on [-1, 1]
1708   @param[out] q_weight_1d Array of length Q to hold the weights
1709 
1710   @return An error code: 0 - success, otherwise - failure
1711 
1712   @ref Utility
1713 **/
1714 int CeedGaussQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) {
1715   // Allocate
1716   CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0 * atan(1.0);
1717   // Build q_ref_1d, q_weight_1d
1718   for (CeedInt i = 0; i <= Q / 2; i++) {
1719     // Guess
1720     xi = cos(PI * (CeedScalar)(2 * i + 1) / ((CeedScalar)(2 * Q)));
1721     // Pn(xi)
1722     P0 = 1.0;
1723     P1 = xi;
1724     P2 = 0.0;
1725     for (CeedInt j = 2; j <= Q; j++) {
1726       P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
1727       P0 = P1;
1728       P1 = P2;
1729     }
1730     // First Newton Step
1731     dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
1732     xi  = xi - P2 / dP2;
1733     // Newton to convergence
1734     for (CeedInt k = 0; k < 100 && fabs(P2) > 10 * CEED_EPSILON; k++) {
1735       P0 = 1.0;
1736       P1 = xi;
1737       for (CeedInt j = 2; j <= Q; j++) {
1738         P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
1739         P0 = P1;
1740         P1 = P2;
1741       }
1742       dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
1743       xi  = xi - P2 / dP2;
1744     }
1745     // Save xi, wi
1746     wi                     = 2.0 / ((1.0 - xi * xi) * dP2 * dP2);
1747     q_weight_1d[i]         = wi;
1748     q_weight_1d[Q - 1 - i] = wi;
1749     q_ref_1d[i]            = -xi;
1750     q_ref_1d[Q - 1 - i]    = xi;
1751   }
1752   return CEED_ERROR_SUCCESS;
1753 }
1754 
1755 /**
1756   @brief Construct a Gauss-Legendre-Lobatto quadrature
1757 
1758   @param[in]  Q           Number of quadrature points (integrates polynomials of degree 2*Q-3 exactly)
1759   @param[out] q_ref_1d    Array of length Q to hold the abscissa on [-1, 1]
1760   @param[out] q_weight_1d Array of length Q to hold the weights
1761 
1762   @return An error code: 0 - success, otherwise - failure
1763 
1764   @ref Utility
1765 **/
1766 int CeedLobattoQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) {
1767   // Allocate
1768   CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0 * atan(1.0);
1769   // Build q_ref_1d, q_weight_1d
1770   // Set endpoints
1771   CeedCheck(Q > 1, NULL, CEED_ERROR_DIMENSION, "Cannot create Lobatto quadrature with Q=%" CeedInt_FMT " < 2 points", Q);
1772   wi = 2.0 / ((CeedScalar)(Q * (Q - 1)));
1773   if (q_weight_1d) {
1774     q_weight_1d[0]     = wi;
1775     q_weight_1d[Q - 1] = wi;
1776   }
1777   q_ref_1d[0]     = -1.0;
1778   q_ref_1d[Q - 1] = 1.0;
1779   // Interior
1780   for (CeedInt i = 1; i <= (Q - 1) / 2; i++) {
1781     // Guess
1782     xi = cos(PI * (CeedScalar)(i) / (CeedScalar)(Q - 1));
1783     // Pn(xi)
1784     P0 = 1.0;
1785     P1 = xi;
1786     P2 = 0.0;
1787     for (CeedInt j = 2; j < Q; j++) {
1788       P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
1789       P0 = P1;
1790       P1 = P2;
1791     }
1792     // First Newton step
1793     dP2  = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
1794     d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi);
1795     xi   = xi - dP2 / d2P2;
1796     // Newton to convergence
1797     for (CeedInt k = 0; k < 100 && fabs(dP2) > 10 * CEED_EPSILON; k++) {
1798       P0 = 1.0;
1799       P1 = xi;
1800       for (CeedInt j = 2; j < Q; j++) {
1801         P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
1802         P0 = P1;
1803         P1 = P2;
1804       }
1805       dP2  = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
1806       d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi);
1807       xi   = xi - dP2 / d2P2;
1808     }
1809     // Save xi, wi
1810     wi = 2.0 / (((CeedScalar)(Q * (Q - 1))) * P2 * P2);
1811     if (q_weight_1d) {
1812       q_weight_1d[i]         = wi;
1813       q_weight_1d[Q - 1 - i] = wi;
1814     }
1815     q_ref_1d[i]         = -xi;
1816     q_ref_1d[Q - 1 - i] = xi;
1817   }
1818   return CEED_ERROR_SUCCESS;
1819 }
1820 
1821 /// @}
1822