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