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