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