xref: /libCEED/rust/libceed-sys/c-src/interface/ceed-basis.c (revision 7a982d89c751e293e39d23a7c44a161cef1fcd06)
1 // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3 // reserved. See files LICENSE and NOTICE for details.
4 //
5 // This file is part of CEED, a collection of benchmarks, miniapps, software
6 // libraries and APIs for efficient high-order finite element and spectral
7 // element discretizations for exascale applications. For more information and
8 // source code availability see http://github.com/ceed.
9 //
10 // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11 // a collaborative effort of two U.S. Department of Energy organizations (Office
12 // of Science and the National Nuclear Security Administration) responsible for
13 // the planning and preparation of a capable exascale ecosystem, including
14 // software, applications, hardware, advanced system engineering and early
15 // testbed platforms, in support of the nation's exascale computing imperative.
16 
17 #include <ceed-impl.h>
18 #include <ceed-backend.h>
19 #include <math.h>
20 #include <stdio.h>
21 #include <stdlib.h>
22 #include <string.h>
23 
24 /// @file
25 /// Implementation of CeedBasis interfaces
26 
27 /// @cond DOXYGEN_SKIP
28 static struct CeedBasis_private ceed_basis_collocated;
29 /// @endcond
30 
31 /// @addtogroup CeedBasisUser
32 /// @{
33 
34 /// Indicate that the quadrature points are collocated with the nodes
35 const CeedBasis CEED_BASIS_COLLOCATED = &ceed_basis_collocated;
36 
37 /// @}
38 
39 /// ----------------------------------------------------------------------------
40 /// CeedBasis Library Internal Functions
41 /// ----------------------------------------------------------------------------
42 /// @addtogroup CeedBasisDeveloper
43 /// @{
44 
45 /**
46   @brief Compute Householder reflection
47 
48     Computes A = (I - b v v^T) A
49     where A is an mxn matrix indexed as A[i*row + j*col]
50 
51   @param[in,out] A  Matrix to apply Householder reflection to, in place
52   @param v          Householder vector
53   @param b          Scaling factor
54   @param m          Number of rows in A
55   @param n          Number of columns in A
56   @param row        Row stride
57   @param col        Col stride
58 
59   @return An error code: 0 - success, otherwise - failure
60 
61   @ref Developer
62 **/
63 static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v,
64                                   CeedScalar b, CeedInt m, CeedInt n,
65                                   CeedInt row, CeedInt col) {
66   for (CeedInt j=0; j<n; j++) {
67     CeedScalar w = A[0*row + j*col];
68     for (CeedInt i=1; i<m; i++)
69       w += v[i] * A[i*row + j*col];
70     A[0*row + j*col] -= b * w;
71     for (CeedInt i=1; i<m; i++)
72       A[i*row + j*col] -= b * w * v[i];
73   }
74   return 0;
75 }
76 
77 /**
78   @brief Apply Householder Q matrix
79 
80     Compute A = Q A where Q is mxm and A is mxn.
81 
82   @param[in,out] A  Matrix to apply Householder Q to, in place
83   @param Q          Householder Q matrix
84   @param tau        Householder scaling factors
85   @param tmode      Transpose mode for application
86   @param m          Number of rows in A
87   @param n          Number of columns in A
88   @param k          Number of elementary reflectors in Q, k<m
89   @param row        Row stride in A
90   @param col        Col stride in A
91 
92   @return An error code: 0 - success, otherwise - failure
93 
94   @ref Developer
95 **/
96 static int CeedHouseholderApplyQ(CeedScalar *A, const CeedScalar *Q,
97                                  const CeedScalar *tau, CeedTransposeMode tmode,
98                                  CeedInt m, CeedInt n, CeedInt k,
99                                  CeedInt row, CeedInt col) {
100   CeedScalar v[m];
101   for (CeedInt ii=0; ii<k; ii++) {
102     CeedInt i = tmode == CEED_TRANSPOSE ? ii : k-1-ii;
103     for (CeedInt j=i+1; j<m; j++)
104       v[j] = Q[j*k+i];
105     // Apply Householder reflector (I - tau v v^T) collograd1d^T
106     CeedHouseholderReflect(&A[i*row], &v[i], tau[i], m-i, n, row, col);
107   }
108   return 0;
109 }
110 
111 /**
112   @brief Compute Givens rotation
113 
114     Computes A = G A (or G^T A in transpose mode)
115     where A is an mxn matrix indexed as A[i*n + j*m]
116 
117   @param[in,out] A  Row major matrix to apply Givens rotation to, in place
118   @param c          Cosine factor
119   @param s          Sine factor
120   @param i          First row/column to apply rotation
121   @param k          Second row/column to apply rotation
122   @param m          Number of rows in A
123   @param n          Number of columns in A
124 
125   @return An error code: 0 - success, otherwise - failure
126 
127   @ref Developer
128 **/
129 static int CeedGivensRotation(CeedScalar *A, CeedScalar c, CeedScalar s,
130                               CeedTransposeMode tmode, CeedInt i, CeedInt k,
131                               CeedInt m, CeedInt n) {
132   CeedInt stridej = 1, strideik = m, numits = n;
133   if (tmode == CEED_NOTRANSPOSE) {
134     stridej = n; strideik = 1; numits = m;
135   }
136 
137   // Apply rotation
138   for (CeedInt j=0; j<numits; j++) {
139     CeedScalar tau1 = A[i*strideik+j*stridej], tau2 = A[k*strideik+j*stridej];
140     A[i*strideik+j*stridej] = c*tau1 - s*tau2;
141     A[k*strideik+j*stridej] = s*tau1 + c*tau2;
142   }
143 
144   return 0;
145 }
146 
147 /**
148   @brief View an array stored in a CeedBasis
149 
150   @param name      Name of array
151   @param fpformat  Printing format
152   @param m         Number of rows in array
153   @param n         Number of columns in array
154   @param a         Array to be viewed
155   @param stream    Stream to view to, e.g., stdout
156 
157   @return An error code: 0 - success, otherwise - failure
158 
159   @ref Developer
160 **/
161 static int CeedScalarView(const char *name, const char *fpformat, CeedInt m,
162                           CeedInt n, const CeedScalar *a, FILE *stream) {
163   for (int i=0; i<m; i++) {
164     if (m > 1)
165       fprintf(stream, "%12s[%d]:", name, i);
166     else
167       fprintf(stream, "%12s:", name);
168     for (int j=0; j<n; j++)
169       fprintf(stream, fpformat, fabs(a[i*n+j]) > 1E-14 ? a[i*n+j] : 0);
170     fputs("\n", stream);
171   }
172   return 0;
173 }
174 
175 /// @}
176 
177 /// ----------------------------------------------------------------------------
178 /// Ceed Backend API
179 /// ----------------------------------------------------------------------------
180 /// @addtogroup CeedBasisBackend
181 /// @{
182 
183 /**
184   @brief Return collocated grad matrix
185 
186   @param basis             CeedBasis
187   @param[out] collograd1d  Row-major (Q1d * Q1d) matrix expressing derivatives of
188                             basis functions at quadrature points
189 
190   @return An error code: 0 - success, otherwise - failure
191 
192   @ref Backend
193 **/
194 int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *collograd1d) {
195   int i, j, k;
196   Ceed ceed;
197   CeedInt ierr, P1d=(basis)->P1d, Q1d=(basis)->Q1d;
198   CeedScalar *interp1d, *grad1d, tau[Q1d];
199 
200   ierr = CeedMalloc(Q1d*P1d, &interp1d); CeedChk(ierr);
201   ierr = CeedMalloc(Q1d*P1d, &grad1d); CeedChk(ierr);
202   memcpy(interp1d, (basis)->interp1d, Q1d*P1d*sizeof(basis)->interp1d[0]);
203   memcpy(grad1d, (basis)->grad1d, Q1d*P1d*sizeof(basis)->interp1d[0]);
204 
205   // QR Factorization, interp1d = Q R
206   ierr = CeedBasisGetCeed(basis, &ceed); CeedChk(ierr);
207   ierr = CeedQRFactorization(ceed, interp1d, tau, Q1d, P1d); CeedChk(ierr);
208 
209   // Apply Rinv, collograd1d = grad1d Rinv
210   for (i=0; i<Q1d; i++) { // Row i
211     collograd1d[Q1d*i] = grad1d[P1d*i]/interp1d[0];
212     for (j=1; j<P1d; j++) { // Column j
213       collograd1d[j+Q1d*i] = grad1d[j+P1d*i];
214       for (k=0; k<j; k++)
215         collograd1d[j+Q1d*i] -= interp1d[j+P1d*k]*collograd1d[k+Q1d*i];
216       collograd1d[j+Q1d*i] /= interp1d[j+P1d*j];
217     }
218     for (j=P1d; j<Q1d; j++)
219       collograd1d[j+Q1d*i] = 0;
220   }
221 
222   // Apply Qtranspose, collograd = collograd Qtranspose
223   CeedHouseholderApplyQ(collograd1d, interp1d, tau, CEED_NOTRANSPOSE,
224                         Q1d, Q1d, P1d, 1, Q1d);
225 
226   ierr = CeedFree(&interp1d); CeedChk(ierr);
227   ierr = CeedFree(&grad1d); CeedChk(ierr);
228 
229   return 0;
230 }
231 
232 /**
233   @brief Get Ceed associated with a CeedBasis
234 
235   @param basis      CeedBasis
236   @param[out] ceed  Variable to store Ceed
237 
238   @return An error code: 0 - success, otherwise - failure
239 
240   @ref Backend
241 **/
242 int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) {
243   *ceed = basis->ceed;
244   return 0;
245 }
246 
247 /**
248   @brief Get tensor status for given CeedBasis
249 
250   @param basis        CeedBasis
251   @param[out] tensor  Variable to store tensor status
252 
253   @return An error code: 0 - success, otherwise - failure
254 
255   @ref Backend
256 **/
257 int CeedBasisGetTensorStatus(CeedBasis basis, bool *tensor) {
258   *tensor = basis->tensorbasis;
259   return 0;
260 }
261 
262 /**
263   @brief Get dimension for given CeedBasis
264 
265   @param basis     CeedBasis
266   @param[out] dim  Variable to store dimension of basis
267 
268   @return An error code: 0 - success, otherwise - failure
269 
270   @ref Backend
271 **/
272 int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) {
273   *dim = basis->dim;
274   return 0;
275 }
276 
277 /**
278   @brief Get number of components for given CeedBasis
279 
280   @param basis         CeedBasis
281   @param[out] numcomp  Variable to store number of components of basis
282 
283   @return An error code: 0 - success, otherwise - failure
284 
285   @ref Backend
286 **/
287 int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *numcomp) {
288   *numcomp = basis->ncomp;
289   return 0;
290 }
291 
292 /**
293   @brief Get total number of nodes (in 1 dimension) of a CeedBasis
294 
295   @param basis     CeedBasis
296   @param[out] P1d  Variable to store number of nodes
297 
298   @return An error code: 0 - success, otherwise - failure
299 
300   @ref Backend
301 **/
302 int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P1d) {
303   if (!basis->tensorbasis)
304     // LCOV_EXCL_START
305     return CeedError(basis->ceed, 1, "Cannot supply P1d for non-tensor basis");
306   // LCOV_EXCL_STOP
307 
308   *P1d = basis->P1d;
309   return 0;
310 }
311 
312 /**
313   @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis
314 
315   @param basis     CeedBasis
316   @param[out] Q1d  Variable to store number of quadrature points
317 
318   @return An error code: 0 - success, otherwise - failure
319 
320   @ref Backend
321 **/
322 int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q1d) {
323   if (!basis->tensorbasis)
324     // LCOV_EXCL_START
325     return CeedError(basis->ceed, 1, "Cannot supply Q1d for non-tensor basis");
326   // LCOV_EXCL_STOP
327 
328   *Q1d = basis->Q1d;
329   return 0;
330 }
331 
332 /**
333   @brief Get reference coordinates of quadrature points (in dim dimensions)
334          of a CeedBasis
335 
336   @param basis      CeedBasis
337   @param[out] qref  Variable to store reference coordinates of quadrature points
338 
339   @return An error code: 0 - success, otherwise - failure
340 
341   @ref Backend
342 **/
343 int CeedBasisGetQRef(CeedBasis basis, CeedScalar **qref) {
344   *qref = basis->qref1d;
345   return 0;
346 }
347 
348 /**
349   @brief Get quadrature weights of quadrature points (in dim dimensions)
350          of a CeedBasis
351 
352   @param basis         CeedBasis
353   @param[out] qweight  Variable to store quadrature weights
354 
355   @return An error code: 0 - success, otherwise - failure
356 
357   @ref Backend
358 **/
359 int CeedBasisGetQWeights(CeedBasis basis, CeedScalar **qweight) {
360   *qweight = basis->qweight1d;
361   return 0;
362 }
363 
364 /**
365   @brief Get interpolation matrix of a CeedBasis
366 
367   @param basis        CeedBasis
368   @param[out] interp  Variable to store interpolation matrix
369 
370   @return An error code: 0 - success, otherwise - failure
371 
372   @ref Backend
373 **/
374 int CeedBasisGetInterp(CeedBasis basis, CeedScalar **interp) {
375   if (!basis->interp && basis->tensorbasis) {
376     // Allocate
377     int ierr;
378     ierr = CeedMalloc(basis->Q*basis->P, &basis->interp); CeedChk(ierr);
379 
380     // Initialize
381     for (CeedInt i=0; i<basis->Q*basis->P; i++)
382       basis->interp[i] = 1.0;
383 
384     // Calculate
385     for (CeedInt d=0; d<basis->dim; d++)
386       for (CeedInt qpt=0; qpt<basis->Q; qpt++)
387         for (CeedInt node=0; node<basis->P; node++) {
388           CeedInt p = (node / CeedIntPow(basis->P1d, d)) % basis->P1d;
389           CeedInt q = (qpt / CeedIntPow(basis->Q1d, d)) % basis->Q1d;
390           basis->interp[qpt*(basis->P)+node] *= basis->interp1d[q*basis->P1d+p];
391         }
392   }
393 
394   *interp = basis->interp;
395 
396   return 0;
397 }
398 
399 /**
400   @brief Get 1D interpolation matrix of a tensor product CeedBasis
401 
402   @param basis          CeedBasis
403   @param[out] interp1d  Variable to store interpolation matrix
404 
405   @return An error code: 0 - success, otherwise - failure
406 
407   @ref Backend
408 **/
409 int CeedBasisGetInterp1D(CeedBasis basis, CeedScalar **interp1d) {
410   if (!basis->tensorbasis)
411     // LCOV_EXCL_START
412     return CeedError(basis->ceed, 1, "CeedBasis is not a tensor product basis.");
413   // LCOV_EXCL_STOP
414 
415   *interp1d = basis->interp1d;
416 
417   return 0;
418 }
419 
420 /**
421   @brief Get gradient matrix of a CeedBasis
422 
423   @param basis      CeedBasis
424   @param[out] grad  Variable to store gradient matrix
425 
426   @return An error code: 0 - success, otherwise - failure
427 
428   @ref Backend
429 **/
430 int CeedBasisGetGrad(CeedBasis basis, CeedScalar **grad) {
431   if (!basis->grad && basis->tensorbasis) {
432     // Allocate
433     int ierr;
434     ierr = CeedMalloc(basis->dim*basis->Q*basis->P, &basis->grad);
435     CeedChk(ierr);
436 
437     // Initialize
438     for (CeedInt i=0; i<basis->dim*basis->Q*basis->P; i++)
439       basis->grad[i] = 1.0;
440 
441     // Calculate
442     for (CeedInt d=0; d<basis->dim; d++)
443       for (CeedInt i=0; i<basis->dim; i++)
444         for (CeedInt qpt=0; qpt<basis->Q; qpt++)
445           for (CeedInt node=0; node<basis->P; node++) {
446             CeedInt p = (node / CeedIntPow(basis->P1d, d)) % basis->P1d;
447             CeedInt q = (qpt / CeedIntPow(basis->Q1d, d)) % basis->Q1d;
448             if (i == d)
449               basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *=
450                 basis->grad1d[q*basis->P1d+p];
451             else
452               basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *=
453                 basis->interp1d[q*basis->P1d+p];
454           }
455   }
456 
457   *grad = basis->grad;
458 
459   return 0;
460 }
461 
462 /**
463   @brief Get 1D gradient matrix of a tensor product CeedBasis
464 
465   @param basis        CeedBasis
466   @param[out] grad1d  Variable to store gradient matrix
467 
468   @return An error code: 0 - success, otherwise - failure
469 
470   @ref Backend
471 **/
472 int CeedBasisGetGrad1D(CeedBasis basis, CeedScalar **grad1d) {
473   if (!basis->tensorbasis)
474     // LCOV_EXCL_START
475     return CeedError(basis->ceed, 1, "CeedBasis is not a tensor product basis.");
476   // LCOV_EXCL_STOP
477 
478   *grad1d = basis->grad1d;
479 
480   return 0;
481 }
482 
483 /**
484   @brief Get backend data of a CeedBasis
485 
486   @param basis      CeedBasis
487   @param[out] data  Variable to store data
488 
489   @return An error code: 0 - success, otherwise - failure
490 
491   @ref Backend
492 **/
493 int CeedBasisGetData(CeedBasis basis, void **data) {
494   *data = basis->data;
495   return 0;
496 }
497 
498 /**
499   @brief Set backend data of a CeedBasis
500 
501   @param[out] basis  CeedBasis
502   @param data        Data to set
503 
504   @return An error code: 0 - success, otherwise - failure
505 
506   @ref Backend
507 **/
508 int CeedBasisSetData(CeedBasis basis, void **data) {
509   basis->data = *data;
510   return 0;
511 }
512 
513 /**
514   @brief Get dimension for given CeedElemTopology
515 
516   @param topo      CeedElemTopology
517   @param[out] dim  Variable to store dimension of topology
518 
519   @return An error code: 0 - success, otherwise - failure
520 
521   @ref Backend
522 **/
523 int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) {
524   *dim = (CeedInt) topo >> 16;
525   return 0;
526 }
527 
528 /**
529   @brief Get CeedTensorContract of a CeedBasis
530 
531   @param basis          CeedBasis
532   @param[out] contract  Variable to store CeedTensorContract
533 
534   @return An error code: 0 - success, otherwise - failure
535 
536   @ref Backend
537 **/
538 int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) {
539   *contract = basis->contract;
540   return 0;
541 }
542 
543 /**
544   @brief Set CeedTensorContract of a CeedBasis
545 
546   @param[out] basis     CeedBasis
547   @param contract       CeedTensorContract to set
548 
549   @return An error code: 0 - success, otherwise - failure
550 
551   @ref Backend
552 **/
553 int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract *contract) {
554   basis->contract = *contract;
555   return 0;
556 }
557 
558 /**
559   @brief Return a reference implementation of matrix multiplication C = A B.
560            Note, this is a reference implementation for CPU CeedScalar pointers
561            that is not intended for high performance.
562 
563   @param ceed         A Ceed context for error handling
564   @param[in] matA     Row-major matrix A
565   @param[in] matB     Row-major matrix B
566   @param[out] matC    Row-major output matrix C
567   @param m            Number of rows of C
568   @param n            Number of columns of C
569   @param kk           Number of columns of A/rows of B
570 
571   @return An error code: 0 - success, otherwise - failure
572 
573   @ref Utility
574 **/
575 int CeedMatrixMultiply(Ceed ceed, const CeedScalar *matA,
576                        const CeedScalar *matB, CeedScalar *matC, CeedInt m,
577                        CeedInt n, CeedInt kk) {
578   for (CeedInt i=0; i<m; i++)
579     for (CeedInt j=0; j<n; j++) {
580       CeedScalar sum = 0;
581       for (CeedInt k=0; k<kk; k++)
582         sum += matA[k+i*kk]*matB[j+k*n];
583       matC[j+i*n] = sum;
584     }
585   return 0;
586 }
587 
588 /// @}
589 
590 /// ----------------------------------------------------------------------------
591 /// CeedBasis Public API
592 /// ----------------------------------------------------------------------------
593 /// @addtogroup CeedBasisUser
594 /// @{
595 
596 /**
597   @brief Create a tensor-product basis for H^1 discretizations
598 
599   @param ceed        A Ceed object where the CeedBasis will be created
600   @param dim         Topological dimension
601   @param ncomp       Number of field components (1 for scalar fields)
602   @param P1d         Number of nodes in one dimension
603   @param Q1d         Number of quadrature points in one dimension
604   @param interp1d    Row-major (Q1d * P1d) matrix expressing the values of nodal
605                        basis functions at quadrature points
606   @param grad1d      Row-major (Q1d * P1d) matrix expressing derivatives of nodal
607                        basis functions at quadrature points
608   @param qref1d      Array of length Q1d holding the locations of quadrature points
609                        on the 1D reference element [-1, 1]
610   @param qweight1d   Array of length Q1d holding the quadrature weights on the
611                        reference element
612   @param[out] basis  Address of the variable where the newly created
613                        CeedBasis will be stored.
614 
615   @return An error code: 0 - success, otherwise - failure
616 
617   @ref User
618 **/
619 int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt ncomp, CeedInt P1d,
620                             CeedInt Q1d, const CeedScalar *interp1d,
621                             const CeedScalar *grad1d, const CeedScalar *qref1d,
622                             const CeedScalar *qweight1d, CeedBasis *basis) {
623   int ierr;
624 
625   if (dim<1)
626     // LCOV_EXCL_START
627     return CeedError(ceed, 1, "Basis dimension must be a positive value");
628   // LCOV_EXCL_STOP
629 
630   if (!ceed->BasisCreateTensorH1) {
631     Ceed delegate;
632     ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr);
633 
634     if (!delegate)
635       // LCOV_EXCL_START
636       return CeedError(ceed, 1, "Backend does not support BasisCreateTensorH1");
637     // LCOV_EXCL_STOP
638 
639     ierr = CeedBasisCreateTensorH1(delegate, dim, ncomp, P1d,
640                                    Q1d, interp1d, grad1d, qref1d,
641                                    qweight1d, basis); CeedChk(ierr);
642     return 0;
643   }
644   ierr = CeedCalloc(1,basis); CeedChk(ierr);
645   (*basis)->ceed = ceed;
646   ceed->refcount++;
647   (*basis)->refcount = 1;
648   (*basis)->tensorbasis = 1;
649   (*basis)->dim = dim;
650   (*basis)->ncomp = ncomp;
651   (*basis)->P1d = P1d;
652   (*basis)->Q1d = Q1d;
653   (*basis)->P = CeedIntPow(P1d, dim);
654   (*basis)->Q = CeedIntPow(Q1d, dim);
655   ierr = CeedMalloc(Q1d,&(*basis)->qref1d); CeedChk(ierr);
656   ierr = CeedMalloc(Q1d,&(*basis)->qweight1d); CeedChk(ierr);
657   memcpy((*basis)->qref1d, qref1d, Q1d*sizeof(qref1d[0]));
658   memcpy((*basis)->qweight1d, qweight1d, Q1d*sizeof(qweight1d[0]));
659   ierr = CeedMalloc(Q1d*P1d,&(*basis)->interp1d); CeedChk(ierr);
660   ierr = CeedMalloc(Q1d*P1d,&(*basis)->grad1d); CeedChk(ierr);
661   memcpy((*basis)->interp1d, interp1d, Q1d*P1d*sizeof(interp1d[0]));
662   memcpy((*basis)->grad1d, grad1d, Q1d*P1d*sizeof(grad1d[0]));
663   ierr = ceed->BasisCreateTensorH1(dim, P1d, Q1d, interp1d, grad1d, qref1d,
664                                    qweight1d, *basis); CeedChk(ierr);
665   return 0;
666 }
667 
668 /**
669   @brief Create a tensor-product Lagrange basis
670 
671   @param ceed        A Ceed object where the CeedBasis will be created
672   @param dim         Topological dimension of element
673   @param ncomp       Number of field components (1 for scalar fields)
674   @param P           Number of Gauss-Lobatto nodes in one dimension.  The
675                        polynomial degree of the resulting Q_k element is k=P-1.
676   @param Q           Number of quadrature points in one dimension.
677   @param qmode       Distribution of the Q quadrature points (affects order of
678                        accuracy for the quadrature)
679   @param[out] basis  Address of the variable where the newly created
680                        CeedBasis will be stored.
681 
682   @return An error code: 0 - success, otherwise - failure
683 
684   @ref User
685 **/
686 int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt ncomp,
687                                     CeedInt P, CeedInt Q, CeedQuadMode qmode,
688                                     CeedBasis *basis) {
689   // Allocate
690   int ierr, i, j, k;
691   CeedScalar c1, c2, c3, c4, dx, *nodes, *interp1d, *grad1d, *qref1d, *qweight1d;
692 
693   if (dim<1)
694     // LCOV_EXCL_START
695     return CeedError(ceed, 1, "Basis dimension must be a positive value");
696   // LCOV_EXCL_STOP
697 
698   ierr = CeedCalloc(P*Q, &interp1d); CeedChk(ierr);
699   ierr = CeedCalloc(P*Q, &grad1d); CeedChk(ierr);
700   ierr = CeedCalloc(P, &nodes); CeedChk(ierr);
701   ierr = CeedCalloc(Q, &qref1d); CeedChk(ierr);
702   ierr = CeedCalloc(Q, &qweight1d); CeedChk(ierr);
703   // Get Nodes and Weights
704   ierr = CeedLobattoQuadrature(P, nodes, NULL); CeedChk(ierr);
705   switch (qmode) {
706   case CEED_GAUSS:
707     ierr = CeedGaussQuadrature(Q, qref1d, qweight1d); CeedChk(ierr);
708     break;
709   case CEED_GAUSS_LOBATTO:
710     ierr = CeedLobattoQuadrature(Q, qref1d, qweight1d); CeedChk(ierr);
711     break;
712   }
713   // Build B, D matrix
714   // Fornberg, 1998
715   for (i = 0; i  < Q; i++) {
716     c1 = 1.0;
717     c3 = nodes[0] - qref1d[i];
718     interp1d[i*P+0] = 1.0;
719     for (j = 1; j < P; j++) {
720       c2 = 1.0;
721       c4 = c3;
722       c3 = nodes[j] - qref1d[i];
723       for (k = 0; k < j; k++) {
724         dx = nodes[j] - nodes[k];
725         c2 *= dx;
726         if (k == j - 1) {
727           grad1d[i*P + j] = c1*(interp1d[i*P + k] - c4*grad1d[i*P + k]) / c2;
728           interp1d[i*P + j] = - c1*c4*interp1d[i*P + k] / c2;
729         }
730         grad1d[i*P + k] = (c3*grad1d[i*P + k] - interp1d[i*P + k]) / dx;
731         interp1d[i*P + k] = c3*interp1d[i*P + k] / dx;
732       }
733       c1 = c2;
734     }
735   }
736   //  // Pass to CeedBasisCreateTensorH1
737   ierr = CeedBasisCreateTensorH1(ceed, dim, ncomp, P, Q, interp1d, grad1d, qref1d,
738                                  qweight1d, basis); CeedChk(ierr);
739   ierr = CeedFree(&interp1d); CeedChk(ierr);
740   ierr = CeedFree(&grad1d); CeedChk(ierr);
741   ierr = CeedFree(&nodes); CeedChk(ierr);
742   ierr = CeedFree(&qref1d); CeedChk(ierr);
743   ierr = CeedFree(&qweight1d); CeedChk(ierr);
744   return 0;
745 }
746 
747 /**
748   @brief Create a non tensor-product basis for H^1 discretizations
749 
750   @param ceed        A Ceed object where the CeedBasis will be created
751   @param topo        Topology of element, e.g. hypercube, simplex, ect
752   @param ncomp       Number of field components (1 for scalar fields)
753   @param nnodes      Total number of nodes
754   @param nqpts       Total number of quadrature points
755   @param interp      Row-major (nqpts * nnodes) matrix expressing the values of
756                        nodal basis functions at quadrature points
757   @param grad        Row-major (nqpts * dim * nnodes) matrix expressing
758                        derivatives of nodal basis functions at quadrature points
759   @param qref        Array of length nqpts holding the locations of quadrature
760                        points on the reference element [-1, 1]
761   @param qweight     Array of length nqpts holding the quadrature weights on the
762                        reference element
763   @param[out] basis  Address of the variable where the newly created
764                        CeedBasis will be stored.
765 
766   @return An error code: 0 - success, otherwise - failure
767 
768   @ref User
769 **/
770 int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt ncomp,
771                       CeedInt nnodes, CeedInt nqpts, const CeedScalar *interp,
772                       const CeedScalar *grad, const CeedScalar *qref,
773                       const CeedScalar *qweight, CeedBasis *basis) {
774   int ierr;
775   CeedInt P = nnodes, Q = nqpts, dim = 0;
776 
777   if (!ceed->BasisCreateH1) {
778     Ceed delegate;
779     ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr);
780 
781     if (!delegate)
782       // LCOV_EXCL_START
783       return CeedError(ceed, 1, "Backend does not support BasisCreateH1");
784     // LCOV_EXCL_STOP
785 
786     ierr = CeedBasisCreateH1(delegate, topo, ncomp, nnodes,
787                              nqpts, interp, grad, qref,
788                              qweight, basis); CeedChk(ierr);
789     return 0;
790   }
791 
792   ierr = CeedCalloc(1,basis); CeedChk(ierr);
793 
794   ierr = CeedBasisGetTopologyDimension(topo, &dim); CeedChk(ierr);
795 
796   (*basis)->ceed = ceed;
797   ceed->refcount++;
798   (*basis)->refcount = 1;
799   (*basis)->tensorbasis = 0;
800   (*basis)->dim = dim;
801   (*basis)->ncomp = ncomp;
802   (*basis)->P = P;
803   (*basis)->Q = Q;
804   ierr = CeedMalloc(Q*dim,&(*basis)->qref1d); CeedChk(ierr);
805   ierr = CeedMalloc(Q,&(*basis)->qweight1d); CeedChk(ierr);
806   memcpy((*basis)->qref1d, qref, Q*dim*sizeof(qref[0]));
807   memcpy((*basis)->qweight1d, qweight, Q*sizeof(qweight[0]));
808   ierr = CeedMalloc(Q*P, &(*basis)->interp); CeedChk(ierr);
809   ierr = CeedMalloc(dim*Q*P, &(*basis)->grad); CeedChk(ierr);
810   memcpy((*basis)->interp, interp, Q*P*sizeof(interp[0]));
811   memcpy((*basis)->grad, grad, dim*Q*P*sizeof(grad[0]));
812   ierr = ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, qref,
813                              qweight, *basis); CeedChk(ierr);
814   return 0;
815 }
816 
817 /**
818   @brief View a CeedBasis
819 
820   @param basis   CeedBasis to view
821   @param stream  Stream to view to, e.g., stdout
822 
823   @return An error code: 0 - success, otherwise - failure
824 
825   @ref User
826 **/
827 int CeedBasisView(CeedBasis basis, FILE *stream) {
828   int ierr;
829 
830   if (basis->tensorbasis) {
831     fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P1d,
832             basis->Q1d);
833     ierr = CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q1d, basis->qref1d,
834                           stream); CeedChk(ierr);
835     ierr = CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q1d,
836                           basis->qweight1d, stream); CeedChk(ierr);
837     ierr = CeedScalarView("interp1d", "\t% 12.8f", basis->Q1d, basis->P1d,
838                           basis->interp1d, stream); CeedChk(ierr);
839     ierr = CeedScalarView("grad1d", "\t% 12.8f", basis->Q1d, basis->P1d,
840                           basis->grad1d, stream); CeedChk(ierr);
841   } else {
842     fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P,
843             basis->Q);
844     ierr = CeedScalarView("qref", "\t% 12.8f", 1, basis->Q*basis->dim,
845                           basis->qref1d,
846                           stream); CeedChk(ierr);
847     ierr = CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->qweight1d,
848                           stream); CeedChk(ierr);
849     ierr = CeedScalarView("interp", "\t% 12.8f", basis->Q, basis->P,
850                           basis->interp, stream); CeedChk(ierr);
851     ierr = CeedScalarView("grad", "\t% 12.8f", basis->dim*basis->Q, basis->P,
852                           basis->grad, stream); CeedChk(ierr);
853   }
854   return 0;
855 }
856 
857 /**
858   @brief Get total number of nodes (in dim dimensions) of a CeedBasis
859 
860   @param basis   CeedBasis
861   @param[out] P  Variable to store number of nodes
862 
863   @return An error code: 0 - success, otherwise - failure
864 
865   @ref Utility
866 **/
867 int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) {
868   *P = basis->P;
869   return 0;
870 }
871 
872 /**
873   @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis
874 
875   @param basis   CeedBasis
876   @param[out] Q  Variable to store number of quadrature points
877 
878   @return An error code: 0 - success, otherwise - failure
879 
880   @ref Utility
881 **/
882 int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) {
883   *Q = basis->Q;
884   return 0;
885 }
886 
887 /**
888   @brief Apply basis evaluation from nodes to quadrature points or vice versa
889 
890   @param basis   CeedBasis to evaluate
891   @param nelem   The number of elements to apply the basis evaluation to;
892                    the backend will specify the ordering in
893                    ElemRestrictionCreateBlocked
894   @param tmode   \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature
895                    points, \ref CEED_TRANSPOSE to apply the transpose, mapping
896                    from quadrature points to nodes
897   @param emode   \ref CEED_EVAL_NONE to use values directly,
898                    \ref CEED_EVAL_INTERP to use interpolated values,
899                    \ref CEED_EVAL_GRAD to use gradients,
900                    \ref CEED_EVAL_WEIGHT to use quadrature weights.
901   @param[in] u   Input CeedVector
902   @param[out] v  Output CeedVector
903 
904   @return An error code: 0 - success, otherwise - failure
905 
906   @ref User
907 **/
908 int CeedBasisApply(CeedBasis basis, CeedInt nelem, CeedTransposeMode tmode,
909                    CeedEvalMode emode, CeedVector u, CeedVector v) {
910   int ierr;
911   CeedInt ulength = 0, vlength, nnodes, nqpt;
912   if (!basis->Apply)
913     // LCOV_EXCL_START
914     return CeedError(basis->ceed, 1, "Backend does not support BasisApply");
915   // LCOV_EXCL_STOP
916 
917   // Check compatibility of topological and geometrical dimensions
918   ierr = CeedBasisGetNumNodes(basis, &nnodes); CeedChk(ierr);
919   ierr = CeedBasisGetNumQuadraturePoints(basis, &nqpt); CeedChk(ierr);
920   ierr = CeedVectorGetLength(v, &vlength); CeedChk(ierr);
921 
922   if (u) {
923     ierr = CeedVectorGetLength(u, &ulength); CeedChk(ierr);
924   }
925 
926   if ((tmode == CEED_TRANSPOSE && (vlength%nnodes != 0 || ulength%nqpt != 0)) ||
927       (tmode == CEED_NOTRANSPOSE && (ulength%nnodes != 0 || vlength%nqpt != 0)))
928     return CeedError(basis->ceed, 1, "Length of input/output vectors "
929                      "incompatible with basis dimensions");
930 
931   ierr = basis->Apply(basis, nelem, tmode, emode, u, v); CeedChk(ierr);
932   return 0;
933 }
934 
935 /**
936   @brief Destroy a CeedBasis
937 
938   @param basis CeedBasis to destroy
939 
940   @return An error code: 0 - success, otherwise - failure
941 
942   @ref User
943 **/
944 int CeedBasisDestroy(CeedBasis *basis) {
945   int ierr;
946 
947   if (!*basis || --(*basis)->refcount > 0)
948     return 0;
949   if ((*basis)->Destroy) {
950     ierr = (*basis)->Destroy(*basis); CeedChk(ierr);
951   }
952   ierr = CeedFree(&(*basis)->interp); CeedChk(ierr);
953   ierr = CeedFree(&(*basis)->interp1d); CeedChk(ierr);
954   ierr = CeedFree(&(*basis)->grad); CeedChk(ierr);
955   ierr = CeedFree(&(*basis)->grad1d); CeedChk(ierr);
956   ierr = CeedFree(&(*basis)->qref1d); CeedChk(ierr);
957   ierr = CeedFree(&(*basis)->qweight1d); CeedChk(ierr);
958   ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr);
959   ierr = CeedFree(basis); CeedChk(ierr);
960   return 0;
961 }
962 
963 /**
964   @brief Construct a Gauss-Legendre quadrature
965 
966   @param Q               Number of quadrature points (integrates polynomials of
967                            degree 2*Q-1 exactly)
968   @param[out] qref1d     Array of length Q to hold the abscissa on [-1, 1]
969   @param[out] qweight1d  Array of length Q to hold the weights
970 
971   @return An error code: 0 - success, otherwise - failure
972 
973   @ref Utility
974 **/
975 int CeedGaussQuadrature(CeedInt Q, CeedScalar *qref1d, CeedScalar *qweight1d) {
976   // Allocate
977   CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0*atan(1.0);
978   // Build qref1d, qweight1d
979   for (int i = 0; i <= Q/2; i++) {
980     // Guess
981     xi = cos(PI*(CeedScalar)(2*i+1)/((CeedScalar)(2*Q)));
982     // Pn(xi)
983     P0 = 1.0;
984     P1 = xi;
985     P2 = 0.0;
986     for (int j = 2; j <= Q; j++) {
987       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
988       P0 = P1;
989       P1 = P2;
990     }
991     // First Newton Step
992     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
993     xi = xi-P2/dP2;
994     // Newton to convergence
995     for (int k=0; k<100 && fabs(P2)>10*CEED_EPSILON; k++) {
996       P0 = 1.0;
997       P1 = xi;
998       for (int j = 2; j <= Q; j++) {
999         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1000         P0 = P1;
1001         P1 = P2;
1002       }
1003       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1004       xi = xi-P2/dP2;
1005     }
1006     // Save xi, wi
1007     wi = 2.0/((1.0-xi*xi)*dP2*dP2);
1008     qweight1d[i] = wi;
1009     qweight1d[Q-1-i] = wi;
1010     qref1d[i] = -xi;
1011     qref1d[Q-1-i]= xi;
1012   }
1013   return 0;
1014 }
1015 
1016 /**
1017   @brief Construct a Gauss-Legendre-Lobatto quadrature
1018 
1019   @param Q               Number of quadrature points (integrates polynomials of
1020                            degree 2*Q-3 exactly)
1021   @param[out] qref1d     Array of length Q to hold the abscissa on [-1, 1]
1022   @param[out] qweight1d  Array of length Q to hold the weights
1023 
1024   @return An error code: 0 - success, otherwise - failure
1025 
1026   @ref Utility
1027 **/
1028 int CeedLobattoQuadrature(CeedInt Q, CeedScalar *qref1d,
1029                           CeedScalar *qweight1d) {
1030   // Allocate
1031   CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0*atan(1.0);
1032   // Build qref1d, qweight1d
1033   // Set endpoints
1034   if (Q < 2)
1035     return CeedError(NULL, 1,
1036                      "Cannot create Lobatto quadrature with Q=%d < 2 points", Q);
1037   wi = 2.0/((CeedScalar)(Q*(Q-1)));
1038   if (qweight1d) {
1039     qweight1d[0] = wi;
1040     qweight1d[Q-1] = wi;
1041   }
1042   qref1d[0] = -1.0;
1043   qref1d[Q-1] = 1.0;
1044   // Interior
1045   for (int i = 1; i <= (Q-1)/2; i++) {
1046     // Guess
1047     xi = cos(PI*(CeedScalar)(i)/(CeedScalar)(Q-1));
1048     // Pn(xi)
1049     P0 = 1.0;
1050     P1 = xi;
1051     P2 = 0.0;
1052     for (int j = 2; j < Q; j++) {
1053       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1054       P0 = P1;
1055       P1 = P2;
1056     }
1057     // First Newton step
1058     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1059     d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
1060     xi = xi-dP2/d2P2;
1061     // Newton to convergence
1062     for (int k=0; k<100 && fabs(dP2)>10*CEED_EPSILON; k++) {
1063       P0 = 1.0;
1064       P1 = xi;
1065       for (int j = 2; j < Q; j++) {
1066         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1067         P0 = P1;
1068         P1 = P2;
1069       }
1070       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1071       d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
1072       xi = xi-dP2/d2P2;
1073     }
1074     // Save xi, wi
1075     wi = 2.0/(((CeedScalar)(Q*(Q-1)))*P2*P2);
1076     if (qweight1d) {
1077       qweight1d[i] = wi;
1078       qweight1d[Q-1-i] = wi;
1079     }
1080     qref1d[i] = -xi;
1081     qref1d[Q-1-i]= xi;
1082   }
1083   return 0;
1084 }
1085 
1086 /**
1087   @brief Return QR Factorization of a matrix
1088 
1089   @param ceed         A Ceed context for error handling
1090   @param[in,out] mat  Row-major matrix to be factorized in place
1091   @param[in,out] tau  Vector of length m of scaling factors
1092   @param m            Number of rows
1093   @param n            Number of columns
1094 
1095   @return An error code: 0 - success, otherwise - failure
1096 
1097   @ref Utility
1098 **/
1099 int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau,
1100                         CeedInt m, CeedInt n) {
1101   CeedScalar v[m];
1102 
1103   // Check m >= n
1104   if (n > m)
1105     // LCOV_EXCL_START
1106     return CeedError(ceed, 1, "Cannot compute QR factorization with n > m");
1107   // LCOV_EXCL_STOP
1108 
1109   for (CeedInt i=0; i<n; i++) {
1110     // Calculate Householder vector, magnitude
1111     CeedScalar sigma = 0.0;
1112     v[i] = mat[i+n*i];
1113     for (CeedInt j=i+1; j<m; j++) {
1114       v[j] = mat[i+n*j];
1115       sigma += v[j] * v[j];
1116     }
1117     CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m]
1118     CeedScalar Rii = -copysign(norm, v[i]);
1119     v[i] -= Rii;
1120     // norm of v[i:m] after modification above and scaling below
1121     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
1122     //   tau = 2 / (norm*norm)
1123     tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma);
1124 
1125     for (CeedInt j=i+1; j<m; j++)
1126       v[j] /= v[i];
1127 
1128     // Apply Householder reflector to lower right panel
1129     CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1);
1130     // Save v
1131     mat[i+n*i] = Rii;
1132     for (CeedInt j=i+1; j<m; j++)
1133       mat[i+n*j] = v[j];
1134   }
1135 
1136   return 0;
1137 }
1138 
1139 /**
1140   @brief Return symmetric Schur decomposition of the symmetric matrix mat via
1141            symmetric QR factorization
1142 
1143   @param ceed         A Ceed context for error handling
1144   @param[in,out] mat  Row-major matrix to be factorized in place
1145   @param[out] lambda  Vector of length n of eigenvalues
1146   @param n            Number of rows/columns
1147 
1148   @return An error code: 0 - success, otherwise - failure
1149 
1150   @ref Utility
1151 **/
1152 int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat,
1153                                     CeedScalar *lambda, CeedInt n) {
1154   // Check bounds for clang-tidy
1155   if (n<2)
1156     // LCOV_EXCL_START
1157     return CeedError(ceed, 1,
1158                      "Cannot compute symmetric Schur decomposition of scalars");
1159   // LCOV_EXCL_STOP
1160 
1161   CeedScalar v[n-1], tau[n-1], matT[n*n];
1162 
1163   // Copy mat to matT and set mat to I
1164   memcpy(matT, mat, n*n*sizeof(mat[0]));
1165   for (CeedInt i=0; i<n; i++)
1166     for (CeedInt j=0; j<n; j++)
1167       mat[j+n*i] = (i==j) ? 1 : 0;
1168 
1169   // Reduce to tridiagonal
1170   for (CeedInt i=0; i<n-1; i++) {
1171     // Calculate Householder vector, magnitude
1172     CeedScalar sigma = 0.0;
1173     v[i] = matT[i+n*(i+1)];
1174     for (CeedInt j=i+1; j<n-1; j++) {
1175       v[j] = matT[i+n*(j+1)];
1176       sigma += v[j] * v[j];
1177     }
1178     CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:n-1]
1179     CeedScalar Rii = -copysign(norm, v[i]);
1180     v[i] -= Rii;
1181     // norm of v[i:m] after modification above and scaling below
1182     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
1183     //   tau = 2 / (norm*norm)
1184     if (sigma > 10*CEED_EPSILON)
1185       tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma);
1186     else
1187       tau[i] = 0;
1188 
1189     for (CeedInt j=i+1; j<n-1; j++)
1190       v[j] /= v[i];
1191 
1192     // Update sub and super diagonal
1193     matT[i+n*(i+1)] = Rii;
1194     matT[(i+1)+n*i] = Rii;
1195     for (CeedInt j=i+2; j<n; j++) {
1196       matT[i+n*j] = 0; matT[j+n*i] = 0;
1197     }
1198     // Apply symmetric Householder reflector to lower right panel
1199     CeedHouseholderReflect(&matT[(i+1)+n*(i+1)], &v[i], tau[i],
1200                            n-(i+1), n-(i+1), n, 1);
1201     CeedHouseholderReflect(&matT[(i+1)+n*(i+1)], &v[i], tau[i],
1202                            n-(i+1), n-(i+1), 1, n);
1203     // Save v
1204     for (CeedInt j=i+1; j<n-1; j++) {
1205       matT[i+n*(j+1)] = v[j];
1206     }
1207   }
1208   // Backwards accumulation of Q
1209   for (CeedInt i=n-2; i>=0; i--) {
1210     v[i] = 1;
1211     for (CeedInt j=i+1; j<n-1; j++) {
1212       v[j] = matT[i+n*(j+1)];
1213       matT[i+n*(j+1)] = 0;
1214     }
1215     CeedHouseholderReflect(&mat[(i+1)+n*(i+1)], &v[i], tau[i],
1216                            n-(i+1), n-(i+1), n, 1);
1217   }
1218 
1219   // Reduce sub and super diagonal
1220   CeedInt p = 0, q = 0, itr = 0, maxitr = n*n*n;
1221   CeedScalar tol = 10*CEED_EPSILON;
1222 
1223   while (q < n && itr < maxitr) {
1224     // Update p, q, size of reduced portions of diagonal
1225     p = 0; q = 0;
1226     for (CeedInt i=n-2; i>=0; i--) {
1227       if (fabs(matT[i+n*(i+1)]) < tol)
1228         q += 1;
1229       else
1230         break;
1231     }
1232     for (CeedInt i=0; i<n-1-q; i++) {
1233       if (fabs(matT[i+n*(i+1)]) < tol)
1234         p += 1;
1235       else
1236         break;
1237     }
1238     if (q == n-1) break; // Finished reducing
1239 
1240     // Reduce tridiagonal portion
1241     CeedScalar tnn = matT[(n-1-q)+n*(n-1-q)],
1242                tnnm1 = matT[(n-2-q)+n*(n-1-q)];
1243     CeedScalar d = (matT[(n-2-q)+n*(n-2-q)] - tnn)/2;
1244     CeedScalar mu = tnn - tnnm1*tnnm1 /
1245                     (d + copysign(sqrt(d*d + tnnm1*tnnm1), d));
1246     CeedScalar x = matT[p+n*p] - mu;
1247     CeedScalar z = matT[p+n*(p+1)];
1248     for (CeedInt k=p; k<n-1-q; k++) {
1249       // Compute Givens rotation
1250       CeedScalar c = 1, s = 0;
1251       if (fabs(z) > tol) {
1252         if (fabs(z) > fabs(x)) {
1253           CeedScalar tau = -x/z;
1254           s = 1/sqrt(1+tau*tau), c = s*tau;
1255         } else {
1256           CeedScalar tau = -z/x;
1257           c = 1/sqrt(1+tau*tau), s = c*tau;
1258         }
1259       }
1260 
1261       // Apply Givens rotation to T
1262       CeedGivensRotation(matT, c, s, CEED_NOTRANSPOSE, k, k+1, n, n);
1263       CeedGivensRotation(matT, c, s, CEED_TRANSPOSE, k, k+1, n, n);
1264 
1265       // Apply Givens rotation to Q
1266       CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k+1, n, n);
1267 
1268       // Update x, z
1269       if (k < n-q-2) {
1270         x = matT[k+n*(k+1)];
1271         z = matT[k+n*(k+2)];
1272       }
1273     }
1274     itr++;
1275   }
1276   // Save eigenvalues
1277   for (CeedInt i=0; i<n; i++)
1278     lambda[i] = matT[i+n*i];
1279 
1280   // Check convergence
1281   if (itr == maxitr && q < n-1)
1282     // LCOV_EXCL_START
1283     return CeedError(ceed, 1, "Symmetric QR failed to converge");
1284   // LCOV_EXCL_STOP
1285 
1286   return 0;
1287 }
1288 
1289 /**
1290   @brief Return Simultaneous Diagonalization of two matrices. This solves the
1291            generalized eigenvalue problem A x = lambda B x, where A and B
1292            are symmetric and B is positive definite. We generate the matrix X
1293            and vector Lambda such that X^T A X = Lambda and X^T B X = I. This
1294            is equivalent to the LAPACK routine 'sygv' with TYPE = 1.
1295 
1296   @param ceed         A Ceed context for error handling
1297   @param[in] matA     Row-major matrix to be factorized with eigenvalues
1298   @param[in] matB     Row-major matrix to be factorized to identity
1299   @param[out] x       Row-major orthogonal matrix
1300   @param[out] lambda  Vector of length n of generalized eigenvalues
1301   @param n            Number of rows/columns
1302 
1303   @return An error code: 0 - success, otherwise - failure
1304 
1305   @ref Utility
1306 **/
1307 int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *matA,
1308                                     CeedScalar *matB, CeedScalar *x,
1309                                     CeedScalar *lambda, CeedInt n) {
1310   int ierr;
1311   CeedScalar matC[n*n], matG[n*n], vecD[n];
1312 
1313   // Compute B = G D G^T
1314   memcpy(matG, matB, n*n*sizeof(matB[0]));
1315   ierr = CeedSymmetricSchurDecomposition(ceed, matG, vecD, n); CeedChk(ierr);
1316   for (CeedInt i=0; i<n; i++)
1317     vecD[i] = sqrt(vecD[i]);
1318 
1319   // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T
1320   //           = D^-1/2 G^T A G D^-1/2
1321   for (CeedInt i=0; i<n; i++)
1322     for (CeedInt j=0; j<n; j++)
1323       matC[j+i*n] = matG[i+j*n] / vecD[i];
1324   ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)matC,
1325                             (const CeedScalar *)matA, x, n, n, n);
1326   CeedChk(ierr);
1327   for (CeedInt i=0; i<n; i++)
1328     for (CeedInt j=0; j<n; j++)
1329       matG[j+i*n] = matG[j+i*n] / vecD[j];
1330   ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)x,
1331                             (const CeedScalar *)matG, matC, n, n, n);
1332   CeedChk(ierr);
1333 
1334   // Compute Q^T C Q = lambda
1335   ierr = CeedSymmetricSchurDecomposition(ceed, matC, lambda, n); CeedChk(ierr);
1336 
1337   // Set x = (G D^1/2)^-T Q
1338   //       = G D^-1/2 Q
1339   ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)matG,
1340                             (const CeedScalar *)matC, x, n, n, n);
1341   CeedChk(ierr);
1342 
1343   return 0;
1344 }
1345 
1346 /// @}
1347