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