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