xref: /libCEED/interface/ceed-basis.c (revision 7a982d89c751e293e39d23a7c44a161cef1fcd06)
1d7b241e6Sjeremylt // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2d7b241e6Sjeremylt // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3d7b241e6Sjeremylt // reserved. See files LICENSE and NOTICE for details.
4d7b241e6Sjeremylt //
5d7b241e6Sjeremylt // This file is part of CEED, a collection of benchmarks, miniapps, software
6d7b241e6Sjeremylt // libraries and APIs for efficient high-order finite element and spectral
7d7b241e6Sjeremylt // element discretizations for exascale applications. For more information and
8d7b241e6Sjeremylt // source code availability see http://github.com/ceed.
9d7b241e6Sjeremylt //
10d7b241e6Sjeremylt // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11d7b241e6Sjeremylt // a collaborative effort of two U.S. Department of Energy organizations (Office
12d7b241e6Sjeremylt // of Science and the National Nuclear Security Administration) responsible for
13d7b241e6Sjeremylt // the planning and preparation of a capable exascale ecosystem, including
14d7b241e6Sjeremylt // software, applications, hardware, advanced system engineering and early
15d7b241e6Sjeremylt // testbed platforms, in support of the nation's exascale computing imperative.
16d7b241e6Sjeremylt 
17d7b241e6Sjeremylt #include <ceed-impl.h>
18d863ab9bSjeremylt #include <ceed-backend.h>
19d7b241e6Sjeremylt #include <math.h>
20d7b241e6Sjeremylt #include <stdio.h>
21d7b241e6Sjeremylt #include <stdlib.h>
22d7b241e6Sjeremylt #include <string.h>
23d7b241e6Sjeremylt 
24*7a982d89SJeremy L. Thompson /// @file
25*7a982d89SJeremy L. Thompson /// Implementation of CeedBasis interfaces
26*7a982d89SJeremy L. Thompson 
27d7b241e6Sjeremylt /// @cond DOXYGEN_SKIP
28783c99b3SValeria Barra static struct CeedBasis_private ceed_basis_collocated;
29d7b241e6Sjeremylt /// @endcond
30d7b241e6Sjeremylt 
31*7a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser
32*7a982d89SJeremy L. Thompson /// @{
33*7a982d89SJeremy L. Thompson 
34*7a982d89SJeremy L. Thompson /// Indicate that the quadrature points are collocated with the nodes
35*7a982d89SJeremy L. Thompson const CeedBasis CEED_BASIS_COLLOCATED = &ceed_basis_collocated;
36*7a982d89SJeremy L. Thompson 
37*7a982d89SJeremy L. Thompson /// @}
38*7a982d89SJeremy L. Thompson 
39*7a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
40*7a982d89SJeremy L. Thompson /// CeedBasis Library Internal Functions
41*7a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
42*7a982d89SJeremy L. Thompson /// @addtogroup CeedBasisDeveloper
43*7a982d89SJeremy L. Thompson /// @{
44*7a982d89SJeremy L. Thompson 
45*7a982d89SJeremy L. Thompson /**
46*7a982d89SJeremy L. Thompson   @brief Compute Householder reflection
47*7a982d89SJeremy L. Thompson 
48*7a982d89SJeremy L. Thompson     Computes A = (I - b v v^T) A
49*7a982d89SJeremy L. Thompson     where A is an mxn matrix indexed as A[i*row + j*col]
50*7a982d89SJeremy L. Thompson 
51*7a982d89SJeremy L. Thompson   @param[in,out] A  Matrix to apply Householder reflection to, in place
52*7a982d89SJeremy L. Thompson   @param v          Householder vector
53*7a982d89SJeremy L. Thompson   @param b          Scaling factor
54*7a982d89SJeremy L. Thompson   @param m          Number of rows in A
55*7a982d89SJeremy L. Thompson   @param n          Number of columns in A
56*7a982d89SJeremy L. Thompson   @param row        Row stride
57*7a982d89SJeremy L. Thompson   @param col        Col stride
58*7a982d89SJeremy L. Thompson 
59*7a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
60*7a982d89SJeremy L. Thompson 
61*7a982d89SJeremy L. Thompson   @ref Developer
62*7a982d89SJeremy L. Thompson **/
63*7a982d89SJeremy L. Thompson static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v,
64*7a982d89SJeremy L. Thompson                                   CeedScalar b, CeedInt m, CeedInt n,
65*7a982d89SJeremy L. Thompson                                   CeedInt row, CeedInt col) {
66*7a982d89SJeremy L. Thompson   for (CeedInt j=0; j<n; j++) {
67*7a982d89SJeremy L. Thompson     CeedScalar w = A[0*row + j*col];
68*7a982d89SJeremy L. Thompson     for (CeedInt i=1; i<m; i++)
69*7a982d89SJeremy L. Thompson       w += v[i] * A[i*row + j*col];
70*7a982d89SJeremy L. Thompson     A[0*row + j*col] -= b * w;
71*7a982d89SJeremy L. Thompson     for (CeedInt i=1; i<m; i++)
72*7a982d89SJeremy L. Thompson       A[i*row + j*col] -= b * w * v[i];
73*7a982d89SJeremy L. Thompson   }
74*7a982d89SJeremy L. Thompson   return 0;
75*7a982d89SJeremy L. Thompson }
76*7a982d89SJeremy L. Thompson 
77*7a982d89SJeremy L. Thompson /**
78*7a982d89SJeremy L. Thompson   @brief Apply Householder Q matrix
79*7a982d89SJeremy L. Thompson 
80*7a982d89SJeremy L. Thompson     Compute A = Q A where Q is mxm and A is mxn.
81*7a982d89SJeremy L. Thompson 
82*7a982d89SJeremy L. Thompson   @param[in,out] A  Matrix to apply Householder Q to, in place
83*7a982d89SJeremy L. Thompson   @param Q          Householder Q matrix
84*7a982d89SJeremy L. Thompson   @param tau        Householder scaling factors
85*7a982d89SJeremy L. Thompson   @param tmode      Transpose mode for application
86*7a982d89SJeremy L. Thompson   @param m          Number of rows in A
87*7a982d89SJeremy L. Thompson   @param n          Number of columns in A
88*7a982d89SJeremy L. Thompson   @param k          Number of elementary reflectors in Q, k<m
89*7a982d89SJeremy L. Thompson   @param row        Row stride in A
90*7a982d89SJeremy L. Thompson   @param col        Col stride in A
91*7a982d89SJeremy L. Thompson 
92*7a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
93*7a982d89SJeremy L. Thompson 
94*7a982d89SJeremy L. Thompson   @ref Developer
95*7a982d89SJeremy L. Thompson **/
96*7a982d89SJeremy L. Thompson static int CeedHouseholderApplyQ(CeedScalar *A, const CeedScalar *Q,
97*7a982d89SJeremy L. Thompson                                  const CeedScalar *tau, CeedTransposeMode tmode,
98*7a982d89SJeremy L. Thompson                                  CeedInt m, CeedInt n, CeedInt k,
99*7a982d89SJeremy L. Thompson                                  CeedInt row, CeedInt col) {
100*7a982d89SJeremy L. Thompson   CeedScalar v[m];
101*7a982d89SJeremy L. Thompson   for (CeedInt ii=0; ii<k; ii++) {
102*7a982d89SJeremy L. Thompson     CeedInt i = tmode == CEED_TRANSPOSE ? ii : k-1-ii;
103*7a982d89SJeremy L. Thompson     for (CeedInt j=i+1; j<m; j++)
104*7a982d89SJeremy L. Thompson       v[j] = Q[j*k+i];
105*7a982d89SJeremy L. Thompson     // Apply Householder reflector (I - tau v v^T) collograd1d^T
106*7a982d89SJeremy L. Thompson     CeedHouseholderReflect(&A[i*row], &v[i], tau[i], m-i, n, row, col);
107*7a982d89SJeremy L. Thompson   }
108*7a982d89SJeremy L. Thompson   return 0;
109*7a982d89SJeremy L. Thompson }
110*7a982d89SJeremy L. Thompson 
111*7a982d89SJeremy L. Thompson /**
112*7a982d89SJeremy L. Thompson   @brief Compute Givens rotation
113*7a982d89SJeremy L. Thompson 
114*7a982d89SJeremy L. Thompson     Computes A = G A (or G^T A in transpose mode)
115*7a982d89SJeremy L. Thompson     where A is an mxn matrix indexed as A[i*n + j*m]
116*7a982d89SJeremy L. Thompson 
117*7a982d89SJeremy L. Thompson   @param[in,out] A  Row major matrix to apply Givens rotation to, in place
118*7a982d89SJeremy L. Thompson   @param c          Cosine factor
119*7a982d89SJeremy L. Thompson   @param s          Sine factor
120*7a982d89SJeremy L. Thompson   @param i          First row/column to apply rotation
121*7a982d89SJeremy L. Thompson   @param k          Second row/column to apply rotation
122*7a982d89SJeremy L. Thompson   @param m          Number of rows in A
123*7a982d89SJeremy L. Thompson   @param n          Number of columns in A
124*7a982d89SJeremy L. Thompson 
125*7a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
126*7a982d89SJeremy L. Thompson 
127*7a982d89SJeremy L. Thompson   @ref Developer
128*7a982d89SJeremy L. Thompson **/
129*7a982d89SJeremy L. Thompson static int CeedGivensRotation(CeedScalar *A, CeedScalar c, CeedScalar s,
130*7a982d89SJeremy L. Thompson                               CeedTransposeMode tmode, CeedInt i, CeedInt k,
131*7a982d89SJeremy L. Thompson                               CeedInt m, CeedInt n) {
132*7a982d89SJeremy L. Thompson   CeedInt stridej = 1, strideik = m, numits = n;
133*7a982d89SJeremy L. Thompson   if (tmode == CEED_NOTRANSPOSE) {
134*7a982d89SJeremy L. Thompson     stridej = n; strideik = 1; numits = m;
135*7a982d89SJeremy L. Thompson   }
136*7a982d89SJeremy L. Thompson 
137*7a982d89SJeremy L. Thompson   // Apply rotation
138*7a982d89SJeremy L. Thompson   for (CeedInt j=0; j<numits; j++) {
139*7a982d89SJeremy L. Thompson     CeedScalar tau1 = A[i*strideik+j*stridej], tau2 = A[k*strideik+j*stridej];
140*7a982d89SJeremy L. Thompson     A[i*strideik+j*stridej] = c*tau1 - s*tau2;
141*7a982d89SJeremy L. Thompson     A[k*strideik+j*stridej] = s*tau1 + c*tau2;
142*7a982d89SJeremy L. Thompson   }
143*7a982d89SJeremy L. Thompson 
144*7a982d89SJeremy L. Thompson   return 0;
145*7a982d89SJeremy L. Thompson }
146*7a982d89SJeremy L. Thompson 
147*7a982d89SJeremy L. Thompson /**
148*7a982d89SJeremy L. Thompson   @brief View an array stored in a CeedBasis
149*7a982d89SJeremy L. Thompson 
150*7a982d89SJeremy L. Thompson   @param name      Name of array
151*7a982d89SJeremy L. Thompson   @param fpformat  Printing format
152*7a982d89SJeremy L. Thompson   @param m         Number of rows in array
153*7a982d89SJeremy L. Thompson   @param n         Number of columns in array
154*7a982d89SJeremy L. Thompson   @param a         Array to be viewed
155*7a982d89SJeremy L. Thompson   @param stream    Stream to view to, e.g., stdout
156*7a982d89SJeremy L. Thompson 
157*7a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
158*7a982d89SJeremy L. Thompson 
159*7a982d89SJeremy L. Thompson   @ref Developer
160*7a982d89SJeremy L. Thompson **/
161*7a982d89SJeremy L. Thompson static int CeedScalarView(const char *name, const char *fpformat, CeedInt m,
162*7a982d89SJeremy L. Thompson                           CeedInt n, const CeedScalar *a, FILE *stream) {
163*7a982d89SJeremy L. Thompson   for (int i=0; i<m; i++) {
164*7a982d89SJeremy L. Thompson     if (m > 1)
165*7a982d89SJeremy L. Thompson       fprintf(stream, "%12s[%d]:", name, i);
166*7a982d89SJeremy L. Thompson     else
167*7a982d89SJeremy L. Thompson       fprintf(stream, "%12s:", name);
168*7a982d89SJeremy L. Thompson     for (int j=0; j<n; j++)
169*7a982d89SJeremy L. Thompson       fprintf(stream, fpformat, fabs(a[i*n+j]) > 1E-14 ? a[i*n+j] : 0);
170*7a982d89SJeremy L. Thompson     fputs("\n", stream);
171*7a982d89SJeremy L. Thompson   }
172*7a982d89SJeremy L. Thompson   return 0;
173*7a982d89SJeremy L. Thompson }
174*7a982d89SJeremy L. Thompson 
175*7a982d89SJeremy L. Thompson /// @}
176*7a982d89SJeremy L. Thompson 
177*7a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
178*7a982d89SJeremy L. Thompson /// Ceed Backend API
179*7a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
180*7a982d89SJeremy L. Thompson /// @addtogroup CeedBasisBackend
181*7a982d89SJeremy L. Thompson /// @{
182*7a982d89SJeremy L. Thompson 
183*7a982d89SJeremy L. Thompson /**
184*7a982d89SJeremy L. Thompson   @brief Return collocated grad matrix
185*7a982d89SJeremy L. Thompson 
186*7a982d89SJeremy L. Thompson   @param basis             CeedBasis
187*7a982d89SJeremy L. Thompson   @param[out] collograd1d  Row-major (Q1d * Q1d) matrix expressing derivatives of
188*7a982d89SJeremy L. Thompson                             basis functions at quadrature points
189*7a982d89SJeremy L. Thompson 
190*7a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
191*7a982d89SJeremy L. Thompson 
192*7a982d89SJeremy L. Thompson   @ref Backend
193*7a982d89SJeremy L. Thompson **/
194*7a982d89SJeremy L. Thompson int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *collograd1d) {
195*7a982d89SJeremy L. Thompson   int i, j, k;
196*7a982d89SJeremy L. Thompson   Ceed ceed;
197*7a982d89SJeremy L. Thompson   CeedInt ierr, P1d=(basis)->P1d, Q1d=(basis)->Q1d;
198*7a982d89SJeremy L. Thompson   CeedScalar *interp1d, *grad1d, tau[Q1d];
199*7a982d89SJeremy L. Thompson 
200*7a982d89SJeremy L. Thompson   ierr = CeedMalloc(Q1d*P1d, &interp1d); CeedChk(ierr);
201*7a982d89SJeremy L. Thompson   ierr = CeedMalloc(Q1d*P1d, &grad1d); CeedChk(ierr);
202*7a982d89SJeremy L. Thompson   memcpy(interp1d, (basis)->interp1d, Q1d*P1d*sizeof(basis)->interp1d[0]);
203*7a982d89SJeremy L. Thompson   memcpy(grad1d, (basis)->grad1d, Q1d*P1d*sizeof(basis)->interp1d[0]);
204*7a982d89SJeremy L. Thompson 
205*7a982d89SJeremy L. Thompson   // QR Factorization, interp1d = Q R
206*7a982d89SJeremy L. Thompson   ierr = CeedBasisGetCeed(basis, &ceed); CeedChk(ierr);
207*7a982d89SJeremy L. Thompson   ierr = CeedQRFactorization(ceed, interp1d, tau, Q1d, P1d); CeedChk(ierr);
208*7a982d89SJeremy L. Thompson 
209*7a982d89SJeremy L. Thompson   // Apply Rinv, collograd1d = grad1d Rinv
210*7a982d89SJeremy L. Thompson   for (i=0; i<Q1d; i++) { // Row i
211*7a982d89SJeremy L. Thompson     collograd1d[Q1d*i] = grad1d[P1d*i]/interp1d[0];
212*7a982d89SJeremy L. Thompson     for (j=1; j<P1d; j++) { // Column j
213*7a982d89SJeremy L. Thompson       collograd1d[j+Q1d*i] = grad1d[j+P1d*i];
214*7a982d89SJeremy L. Thompson       for (k=0; k<j; k++)
215*7a982d89SJeremy L. Thompson         collograd1d[j+Q1d*i] -= interp1d[j+P1d*k]*collograd1d[k+Q1d*i];
216*7a982d89SJeremy L. Thompson       collograd1d[j+Q1d*i] /= interp1d[j+P1d*j];
217*7a982d89SJeremy L. Thompson     }
218*7a982d89SJeremy L. Thompson     for (j=P1d; j<Q1d; j++)
219*7a982d89SJeremy L. Thompson       collograd1d[j+Q1d*i] = 0;
220*7a982d89SJeremy L. Thompson   }
221*7a982d89SJeremy L. Thompson 
222*7a982d89SJeremy L. Thompson   // Apply Qtranspose, collograd = collograd Qtranspose
223*7a982d89SJeremy L. Thompson   CeedHouseholderApplyQ(collograd1d, interp1d, tau, CEED_NOTRANSPOSE,
224*7a982d89SJeremy L. Thompson                         Q1d, Q1d, P1d, 1, Q1d);
225*7a982d89SJeremy L. Thompson 
226*7a982d89SJeremy L. Thompson   ierr = CeedFree(&interp1d); CeedChk(ierr);
227*7a982d89SJeremy L. Thompson   ierr = CeedFree(&grad1d); CeedChk(ierr);
228*7a982d89SJeremy L. Thompson 
229*7a982d89SJeremy L. Thompson   return 0;
230*7a982d89SJeremy L. Thompson }
231*7a982d89SJeremy L. Thompson 
232*7a982d89SJeremy L. Thompson /**
233*7a982d89SJeremy L. Thompson   @brief Get Ceed associated with a CeedBasis
234*7a982d89SJeremy L. Thompson 
235*7a982d89SJeremy L. Thompson   @param basis      CeedBasis
236*7a982d89SJeremy L. Thompson   @param[out] ceed  Variable to store Ceed
237*7a982d89SJeremy L. Thompson 
238*7a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
239*7a982d89SJeremy L. Thompson 
240*7a982d89SJeremy L. Thompson   @ref Backend
241*7a982d89SJeremy L. Thompson **/
242*7a982d89SJeremy L. Thompson int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) {
243*7a982d89SJeremy L. Thompson   *ceed = basis->ceed;
244*7a982d89SJeremy L. Thompson   return 0;
245*7a982d89SJeremy L. Thompson }
246*7a982d89SJeremy L. Thompson 
247*7a982d89SJeremy L. Thompson /**
248*7a982d89SJeremy L. Thompson   @brief Get tensor status for given CeedBasis
249*7a982d89SJeremy L. Thompson 
250*7a982d89SJeremy L. Thompson   @param basis        CeedBasis
251*7a982d89SJeremy L. Thompson   @param[out] tensor  Variable to store tensor status
252*7a982d89SJeremy L. Thompson 
253*7a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
254*7a982d89SJeremy L. Thompson 
255*7a982d89SJeremy L. Thompson   @ref Backend
256*7a982d89SJeremy L. Thompson **/
257*7a982d89SJeremy L. Thompson int CeedBasisGetTensorStatus(CeedBasis basis, bool *tensor) {
258*7a982d89SJeremy L. Thompson   *tensor = basis->tensorbasis;
259*7a982d89SJeremy L. Thompson   return 0;
260*7a982d89SJeremy L. Thompson }
261*7a982d89SJeremy L. Thompson 
262*7a982d89SJeremy L. Thompson /**
263*7a982d89SJeremy L. Thompson   @brief Get dimension for given CeedBasis
264*7a982d89SJeremy L. Thompson 
265*7a982d89SJeremy L. Thompson   @param basis     CeedBasis
266*7a982d89SJeremy L. Thompson   @param[out] dim  Variable to store dimension of basis
267*7a982d89SJeremy L. Thompson 
268*7a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
269*7a982d89SJeremy L. Thompson 
270*7a982d89SJeremy L. Thompson   @ref Backend
271*7a982d89SJeremy L. Thompson **/
272*7a982d89SJeremy L. Thompson int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) {
273*7a982d89SJeremy L. Thompson   *dim = basis->dim;
274*7a982d89SJeremy L. Thompson   return 0;
275*7a982d89SJeremy L. Thompson }
276*7a982d89SJeremy L. Thompson 
277*7a982d89SJeremy L. Thompson /**
278*7a982d89SJeremy L. Thompson   @brief Get number of components for given CeedBasis
279*7a982d89SJeremy L. Thompson 
280*7a982d89SJeremy L. Thompson   @param basis         CeedBasis
281*7a982d89SJeremy L. Thompson   @param[out] numcomp  Variable to store number of components of basis
282*7a982d89SJeremy L. Thompson 
283*7a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
284*7a982d89SJeremy L. Thompson 
285*7a982d89SJeremy L. Thompson   @ref Backend
286*7a982d89SJeremy L. Thompson **/
287*7a982d89SJeremy L. Thompson int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *numcomp) {
288*7a982d89SJeremy L. Thompson   *numcomp = basis->ncomp;
289*7a982d89SJeremy L. Thompson   return 0;
290*7a982d89SJeremy L. Thompson }
291*7a982d89SJeremy L. Thompson 
292*7a982d89SJeremy L. Thompson /**
293*7a982d89SJeremy L. Thompson   @brief Get total number of nodes (in 1 dimension) of a CeedBasis
294*7a982d89SJeremy L. Thompson 
295*7a982d89SJeremy L. Thompson   @param basis     CeedBasis
296*7a982d89SJeremy L. Thompson   @param[out] P1d  Variable to store number of nodes
297*7a982d89SJeremy L. Thompson 
298*7a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
299*7a982d89SJeremy L. Thompson 
300*7a982d89SJeremy L. Thompson   @ref Backend
301*7a982d89SJeremy L. Thompson **/
302*7a982d89SJeremy L. Thompson int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P1d) {
303*7a982d89SJeremy L. Thompson   if (!basis->tensorbasis)
304*7a982d89SJeremy L. Thompson     // LCOV_EXCL_START
305*7a982d89SJeremy L. Thompson     return CeedError(basis->ceed, 1, "Cannot supply P1d for non-tensor basis");
306*7a982d89SJeremy L. Thompson   // LCOV_EXCL_STOP
307*7a982d89SJeremy L. Thompson 
308*7a982d89SJeremy L. Thompson   *P1d = basis->P1d;
309*7a982d89SJeremy L. Thompson   return 0;
310*7a982d89SJeremy L. Thompson }
311*7a982d89SJeremy L. Thompson 
312*7a982d89SJeremy L. Thompson /**
313*7a982d89SJeremy L. Thompson   @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis
314*7a982d89SJeremy L. Thompson 
315*7a982d89SJeremy L. Thompson   @param basis     CeedBasis
316*7a982d89SJeremy L. Thompson   @param[out] Q1d  Variable to store number of quadrature points
317*7a982d89SJeremy L. Thompson 
318*7a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
319*7a982d89SJeremy L. Thompson 
320*7a982d89SJeremy L. Thompson   @ref Backend
321*7a982d89SJeremy L. Thompson **/
322*7a982d89SJeremy L. Thompson int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q1d) {
323*7a982d89SJeremy L. Thompson   if (!basis->tensorbasis)
324*7a982d89SJeremy L. Thompson     // LCOV_EXCL_START
325*7a982d89SJeremy L. Thompson     return CeedError(basis->ceed, 1, "Cannot supply Q1d for non-tensor basis");
326*7a982d89SJeremy L. Thompson   // LCOV_EXCL_STOP
327*7a982d89SJeremy L. Thompson 
328*7a982d89SJeremy L. Thompson   *Q1d = basis->Q1d;
329*7a982d89SJeremy L. Thompson   return 0;
330*7a982d89SJeremy L. Thompson }
331*7a982d89SJeremy L. Thompson 
332*7a982d89SJeremy L. Thompson /**
333*7a982d89SJeremy L. Thompson   @brief Get reference coordinates of quadrature points (in dim dimensions)
334*7a982d89SJeremy L. Thompson          of a CeedBasis
335*7a982d89SJeremy L. Thompson 
336*7a982d89SJeremy L. Thompson   @param basis      CeedBasis
337*7a982d89SJeremy L. Thompson   @param[out] qref  Variable to store reference coordinates of quadrature points
338*7a982d89SJeremy L. Thompson 
339*7a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
340*7a982d89SJeremy L. Thompson 
341*7a982d89SJeremy L. Thompson   @ref Backend
342*7a982d89SJeremy L. Thompson **/
343*7a982d89SJeremy L. Thompson int CeedBasisGetQRef(CeedBasis basis, CeedScalar **qref) {
344*7a982d89SJeremy L. Thompson   *qref = basis->qref1d;
345*7a982d89SJeremy L. Thompson   return 0;
346*7a982d89SJeremy L. Thompson }
347*7a982d89SJeremy L. Thompson 
348*7a982d89SJeremy L. Thompson /**
349*7a982d89SJeremy L. Thompson   @brief Get quadrature weights of quadrature points (in dim dimensions)
350*7a982d89SJeremy L. Thompson          of a CeedBasis
351*7a982d89SJeremy L. Thompson 
352*7a982d89SJeremy L. Thompson   @param basis         CeedBasis
353*7a982d89SJeremy L. Thompson   @param[out] qweight  Variable to store quadrature weights
354*7a982d89SJeremy L. Thompson 
355*7a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
356*7a982d89SJeremy L. Thompson 
357*7a982d89SJeremy L. Thompson   @ref Backend
358*7a982d89SJeremy L. Thompson **/
359*7a982d89SJeremy L. Thompson int CeedBasisGetQWeights(CeedBasis basis, CeedScalar **qweight) {
360*7a982d89SJeremy L. Thompson   *qweight = basis->qweight1d;
361*7a982d89SJeremy L. Thompson   return 0;
362*7a982d89SJeremy L. Thompson }
363*7a982d89SJeremy L. Thompson 
364*7a982d89SJeremy L. Thompson /**
365*7a982d89SJeremy L. Thompson   @brief Get interpolation matrix of a CeedBasis
366*7a982d89SJeremy L. Thompson 
367*7a982d89SJeremy L. Thompson   @param basis        CeedBasis
368*7a982d89SJeremy L. Thompson   @param[out] interp  Variable to store interpolation matrix
369*7a982d89SJeremy L. Thompson 
370*7a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
371*7a982d89SJeremy L. Thompson 
372*7a982d89SJeremy L. Thompson   @ref Backend
373*7a982d89SJeremy L. Thompson **/
374*7a982d89SJeremy L. Thompson int CeedBasisGetInterp(CeedBasis basis, CeedScalar **interp) {
375*7a982d89SJeremy L. Thompson   if (!basis->interp && basis->tensorbasis) {
376*7a982d89SJeremy L. Thompson     // Allocate
377*7a982d89SJeremy L. Thompson     int ierr;
378*7a982d89SJeremy L. Thompson     ierr = CeedMalloc(basis->Q*basis->P, &basis->interp); CeedChk(ierr);
379*7a982d89SJeremy L. Thompson 
380*7a982d89SJeremy L. Thompson     // Initialize
381*7a982d89SJeremy L. Thompson     for (CeedInt i=0; i<basis->Q*basis->P; i++)
382*7a982d89SJeremy L. Thompson       basis->interp[i] = 1.0;
383*7a982d89SJeremy L. Thompson 
384*7a982d89SJeremy L. Thompson     // Calculate
385*7a982d89SJeremy L. Thompson     for (CeedInt d=0; d<basis->dim; d++)
386*7a982d89SJeremy L. Thompson       for (CeedInt qpt=0; qpt<basis->Q; qpt++)
387*7a982d89SJeremy L. Thompson         for (CeedInt node=0; node<basis->P; node++) {
388*7a982d89SJeremy L. Thompson           CeedInt p = (node / CeedIntPow(basis->P1d, d)) % basis->P1d;
389*7a982d89SJeremy L. Thompson           CeedInt q = (qpt / CeedIntPow(basis->Q1d, d)) % basis->Q1d;
390*7a982d89SJeremy L. Thompson           basis->interp[qpt*(basis->P)+node] *= basis->interp1d[q*basis->P1d+p];
391*7a982d89SJeremy L. Thompson         }
392*7a982d89SJeremy L. Thompson   }
393*7a982d89SJeremy L. Thompson 
394*7a982d89SJeremy L. Thompson   *interp = basis->interp;
395*7a982d89SJeremy L. Thompson 
396*7a982d89SJeremy L. Thompson   return 0;
397*7a982d89SJeremy L. Thompson }
398*7a982d89SJeremy L. Thompson 
399*7a982d89SJeremy L. Thompson /**
400*7a982d89SJeremy L. Thompson   @brief Get 1D interpolation matrix of a tensor product CeedBasis
401*7a982d89SJeremy L. Thompson 
402*7a982d89SJeremy L. Thompson   @param basis          CeedBasis
403*7a982d89SJeremy L. Thompson   @param[out] interp1d  Variable to store interpolation matrix
404*7a982d89SJeremy L. Thompson 
405*7a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
406*7a982d89SJeremy L. Thompson 
407*7a982d89SJeremy L. Thompson   @ref Backend
408*7a982d89SJeremy L. Thompson **/
409*7a982d89SJeremy L. Thompson int CeedBasisGetInterp1D(CeedBasis basis, CeedScalar **interp1d) {
410*7a982d89SJeremy L. Thompson   if (!basis->tensorbasis)
411*7a982d89SJeremy L. Thompson     // LCOV_EXCL_START
412*7a982d89SJeremy L. Thompson     return CeedError(basis->ceed, 1, "CeedBasis is not a tensor product basis.");
413*7a982d89SJeremy L. Thompson   // LCOV_EXCL_STOP
414*7a982d89SJeremy L. Thompson 
415*7a982d89SJeremy L. Thompson   *interp1d = basis->interp1d;
416*7a982d89SJeremy L. Thompson 
417*7a982d89SJeremy L. Thompson   return 0;
418*7a982d89SJeremy L. Thompson }
419*7a982d89SJeremy L. Thompson 
420*7a982d89SJeremy L. Thompson /**
421*7a982d89SJeremy L. Thompson   @brief Get gradient matrix of a CeedBasis
422*7a982d89SJeremy L. Thompson 
423*7a982d89SJeremy L. Thompson   @param basis      CeedBasis
424*7a982d89SJeremy L. Thompson   @param[out] grad  Variable to store gradient matrix
425*7a982d89SJeremy L. Thompson 
426*7a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
427*7a982d89SJeremy L. Thompson 
428*7a982d89SJeremy L. Thompson   @ref Backend
429*7a982d89SJeremy L. Thompson **/
430*7a982d89SJeremy L. Thompson int CeedBasisGetGrad(CeedBasis basis, CeedScalar **grad) {
431*7a982d89SJeremy L. Thompson   if (!basis->grad && basis->tensorbasis) {
432*7a982d89SJeremy L. Thompson     // Allocate
433*7a982d89SJeremy L. Thompson     int ierr;
434*7a982d89SJeremy L. Thompson     ierr = CeedMalloc(basis->dim*basis->Q*basis->P, &basis->grad);
435*7a982d89SJeremy L. Thompson     CeedChk(ierr);
436*7a982d89SJeremy L. Thompson 
437*7a982d89SJeremy L. Thompson     // Initialize
438*7a982d89SJeremy L. Thompson     for (CeedInt i=0; i<basis->dim*basis->Q*basis->P; i++)
439*7a982d89SJeremy L. Thompson       basis->grad[i] = 1.0;
440*7a982d89SJeremy L. Thompson 
441*7a982d89SJeremy L. Thompson     // Calculate
442*7a982d89SJeremy L. Thompson     for (CeedInt d=0; d<basis->dim; d++)
443*7a982d89SJeremy L. Thompson       for (CeedInt i=0; i<basis->dim; i++)
444*7a982d89SJeremy L. Thompson         for (CeedInt qpt=0; qpt<basis->Q; qpt++)
445*7a982d89SJeremy L. Thompson           for (CeedInt node=0; node<basis->P; node++) {
446*7a982d89SJeremy L. Thompson             CeedInt p = (node / CeedIntPow(basis->P1d, d)) % basis->P1d;
447*7a982d89SJeremy L. Thompson             CeedInt q = (qpt / CeedIntPow(basis->Q1d, d)) % basis->Q1d;
448*7a982d89SJeremy L. Thompson             if (i == d)
449*7a982d89SJeremy L. Thompson               basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *=
450*7a982d89SJeremy L. Thompson                 basis->grad1d[q*basis->P1d+p];
451*7a982d89SJeremy L. Thompson             else
452*7a982d89SJeremy L. Thompson               basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *=
453*7a982d89SJeremy L. Thompson                 basis->interp1d[q*basis->P1d+p];
454*7a982d89SJeremy L. Thompson           }
455*7a982d89SJeremy L. Thompson   }
456*7a982d89SJeremy L. Thompson 
457*7a982d89SJeremy L. Thompson   *grad = basis->grad;
458*7a982d89SJeremy L. Thompson 
459*7a982d89SJeremy L. Thompson   return 0;
460*7a982d89SJeremy L. Thompson }
461*7a982d89SJeremy L. Thompson 
462*7a982d89SJeremy L. Thompson /**
463*7a982d89SJeremy L. Thompson   @brief Get 1D gradient matrix of a tensor product CeedBasis
464*7a982d89SJeremy L. Thompson 
465*7a982d89SJeremy L. Thompson   @param basis        CeedBasis
466*7a982d89SJeremy L. Thompson   @param[out] grad1d  Variable to store gradient matrix
467*7a982d89SJeremy L. Thompson 
468*7a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
469*7a982d89SJeremy L. Thompson 
470*7a982d89SJeremy L. Thompson   @ref Backend
471*7a982d89SJeremy L. Thompson **/
472*7a982d89SJeremy L. Thompson int CeedBasisGetGrad1D(CeedBasis basis, CeedScalar **grad1d) {
473*7a982d89SJeremy L. Thompson   if (!basis->tensorbasis)
474*7a982d89SJeremy L. Thompson     // LCOV_EXCL_START
475*7a982d89SJeremy L. Thompson     return CeedError(basis->ceed, 1, "CeedBasis is not a tensor product basis.");
476*7a982d89SJeremy L. Thompson   // LCOV_EXCL_STOP
477*7a982d89SJeremy L. Thompson 
478*7a982d89SJeremy L. Thompson   *grad1d = basis->grad1d;
479*7a982d89SJeremy L. Thompson 
480*7a982d89SJeremy L. Thompson   return 0;
481*7a982d89SJeremy L. Thompson }
482*7a982d89SJeremy L. Thompson 
483*7a982d89SJeremy L. Thompson /**
484*7a982d89SJeremy L. Thompson   @brief Get backend data of a CeedBasis
485*7a982d89SJeremy L. Thompson 
486*7a982d89SJeremy L. Thompson   @param basis      CeedBasis
487*7a982d89SJeremy L. Thompson   @param[out] data  Variable to store data
488*7a982d89SJeremy L. Thompson 
489*7a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
490*7a982d89SJeremy L. Thompson 
491*7a982d89SJeremy L. Thompson   @ref Backend
492*7a982d89SJeremy L. Thompson **/
493*7a982d89SJeremy L. Thompson int CeedBasisGetData(CeedBasis basis, void **data) {
494*7a982d89SJeremy L. Thompson   *data = basis->data;
495*7a982d89SJeremy L. Thompson   return 0;
496*7a982d89SJeremy L. Thompson }
497*7a982d89SJeremy L. Thompson 
498*7a982d89SJeremy L. Thompson /**
499*7a982d89SJeremy L. Thompson   @brief Set backend data of a CeedBasis
500*7a982d89SJeremy L. Thompson 
501*7a982d89SJeremy L. Thompson   @param[out] basis  CeedBasis
502*7a982d89SJeremy L. Thompson   @param data        Data to set
503*7a982d89SJeremy L. Thompson 
504*7a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
505*7a982d89SJeremy L. Thompson 
506*7a982d89SJeremy L. Thompson   @ref Backend
507*7a982d89SJeremy L. Thompson **/
508*7a982d89SJeremy L. Thompson int CeedBasisSetData(CeedBasis basis, void **data) {
509*7a982d89SJeremy L. Thompson   basis->data = *data;
510*7a982d89SJeremy L. Thompson   return 0;
511*7a982d89SJeremy L. Thompson }
512*7a982d89SJeremy L. Thompson 
513*7a982d89SJeremy L. Thompson /**
514*7a982d89SJeremy L. Thompson   @brief Get dimension for given CeedElemTopology
515*7a982d89SJeremy L. Thompson 
516*7a982d89SJeremy L. Thompson   @param topo      CeedElemTopology
517*7a982d89SJeremy L. Thompson   @param[out] dim  Variable to store dimension of topology
518*7a982d89SJeremy L. Thompson 
519*7a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
520*7a982d89SJeremy L. Thompson 
521*7a982d89SJeremy L. Thompson   @ref Backend
522*7a982d89SJeremy L. Thompson **/
523*7a982d89SJeremy L. Thompson int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) {
524*7a982d89SJeremy L. Thompson   *dim = (CeedInt) topo >> 16;
525*7a982d89SJeremy L. Thompson   return 0;
526*7a982d89SJeremy L. Thompson }
527*7a982d89SJeremy L. Thompson 
528*7a982d89SJeremy L. Thompson /**
529*7a982d89SJeremy L. Thompson   @brief Get CeedTensorContract of a CeedBasis
530*7a982d89SJeremy L. Thompson 
531*7a982d89SJeremy L. Thompson   @param basis          CeedBasis
532*7a982d89SJeremy L. Thompson   @param[out] contract  Variable to store CeedTensorContract
533*7a982d89SJeremy L. Thompson 
534*7a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
535*7a982d89SJeremy L. Thompson 
536*7a982d89SJeremy L. Thompson   @ref Backend
537*7a982d89SJeremy L. Thompson **/
538*7a982d89SJeremy L. Thompson int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) {
539*7a982d89SJeremy L. Thompson   *contract = basis->contract;
540*7a982d89SJeremy L. Thompson   return 0;
541*7a982d89SJeremy L. Thompson }
542*7a982d89SJeremy L. Thompson 
543*7a982d89SJeremy L. Thompson /**
544*7a982d89SJeremy L. Thompson   @brief Set CeedTensorContract of a CeedBasis
545*7a982d89SJeremy L. Thompson 
546*7a982d89SJeremy L. Thompson   @param[out] basis     CeedBasis
547*7a982d89SJeremy L. Thompson   @param contract       CeedTensorContract to set
548*7a982d89SJeremy L. Thompson 
549*7a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
550*7a982d89SJeremy L. Thompson 
551*7a982d89SJeremy L. Thompson   @ref Backend
552*7a982d89SJeremy L. Thompson **/
553*7a982d89SJeremy L. Thompson int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract *contract) {
554*7a982d89SJeremy L. Thompson   basis->contract = *contract;
555*7a982d89SJeremy L. Thompson   return 0;
556*7a982d89SJeremy L. Thompson }
557*7a982d89SJeremy L. Thompson 
558*7a982d89SJeremy L. Thompson /**
559*7a982d89SJeremy L. Thompson   @brief Return a reference implementation of matrix multiplication C = A B.
560*7a982d89SJeremy L. Thompson            Note, this is a reference implementation for CPU CeedScalar pointers
561*7a982d89SJeremy L. Thompson            that is not intended for high performance.
562*7a982d89SJeremy L. Thompson 
563*7a982d89SJeremy L. Thompson   @param ceed         A Ceed context for error handling
564*7a982d89SJeremy L. Thompson   @param[in] matA     Row-major matrix A
565*7a982d89SJeremy L. Thompson   @param[in] matB     Row-major matrix B
566*7a982d89SJeremy L. Thompson   @param[out] matC    Row-major output matrix C
567*7a982d89SJeremy L. Thompson   @param m            Number of rows of C
568*7a982d89SJeremy L. Thompson   @param n            Number of columns of C
569*7a982d89SJeremy L. Thompson   @param kk           Number of columns of A/rows of B
570*7a982d89SJeremy L. Thompson 
571*7a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
572*7a982d89SJeremy L. Thompson 
573*7a982d89SJeremy L. Thompson   @ref Utility
574*7a982d89SJeremy L. Thompson **/
575*7a982d89SJeremy L. Thompson int CeedMatrixMultiply(Ceed ceed, const CeedScalar *matA,
576*7a982d89SJeremy L. Thompson                        const CeedScalar *matB, CeedScalar *matC, CeedInt m,
577*7a982d89SJeremy L. Thompson                        CeedInt n, CeedInt kk) {
578*7a982d89SJeremy L. Thompson   for (CeedInt i=0; i<m; i++)
579*7a982d89SJeremy L. Thompson     for (CeedInt j=0; j<n; j++) {
580*7a982d89SJeremy L. Thompson       CeedScalar sum = 0;
581*7a982d89SJeremy L. Thompson       for (CeedInt k=0; k<kk; k++)
582*7a982d89SJeremy L. Thompson         sum += matA[k+i*kk]*matB[j+k*n];
583*7a982d89SJeremy L. Thompson       matC[j+i*n] = sum;
584*7a982d89SJeremy L. Thompson     }
585*7a982d89SJeremy L. Thompson   return 0;
586*7a982d89SJeremy L. Thompson }
587*7a982d89SJeremy L. Thompson 
588*7a982d89SJeremy L. Thompson /// @}
589*7a982d89SJeremy L. Thompson 
590*7a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
591*7a982d89SJeremy L. Thompson /// CeedBasis Public API
592*7a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
593*7a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser
594d7b241e6Sjeremylt /// @{
595d7b241e6Sjeremylt 
596b11c1e72Sjeremylt /**
59795bb1877Svaleriabarra   @brief Create a tensor-product basis for H^1 discretizations
598b11c1e72Sjeremylt 
599b11c1e72Sjeremylt   @param ceed        A Ceed object where the CeedBasis will be created
600b11c1e72Sjeremylt   @param dim         Topological dimension
601b11c1e72Sjeremylt   @param ncomp       Number of field components (1 for scalar fields)
602b11c1e72Sjeremylt   @param P1d         Number of nodes in one dimension
603b11c1e72Sjeremylt   @param Q1d         Number of quadrature points in one dimension
60495bb1877Svaleriabarra   @param interp1d    Row-major (Q1d * P1d) matrix expressing the values of nodal
605b11c1e72Sjeremylt                        basis functions at quadrature points
60695bb1877Svaleriabarra   @param grad1d      Row-major (Q1d * P1d) matrix expressing derivatives of nodal
607b11c1e72Sjeremylt                        basis functions at quadrature points
608b11c1e72Sjeremylt   @param qref1d      Array of length Q1d holding the locations of quadrature points
609b11c1e72Sjeremylt                        on the 1D reference element [-1, 1]
610b11c1e72Sjeremylt   @param qweight1d   Array of length Q1d holding the quadrature weights on the
611b11c1e72Sjeremylt                        reference element
612b11c1e72Sjeremylt   @param[out] basis  Address of the variable where the newly created
613b11c1e72Sjeremylt                        CeedBasis will be stored.
614b11c1e72Sjeremylt 
615b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
616dfdf5a53Sjeremylt 
617*7a982d89SJeremy L. Thompson   @ref User
618b11c1e72Sjeremylt **/
619d7b241e6Sjeremylt int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt ncomp, CeedInt P1d,
620d7b241e6Sjeremylt                             CeedInt Q1d, const CeedScalar *interp1d,
621d7b241e6Sjeremylt                             const CeedScalar *grad1d, const CeedScalar *qref1d,
622d7b241e6Sjeremylt                             const CeedScalar *qweight1d, CeedBasis *basis) {
623d7b241e6Sjeremylt   int ierr;
624d7b241e6Sjeremylt 
6254d537eeaSYohann   if (dim<1)
626c042f62fSJeremy L Thompson     // LCOV_EXCL_START
6274d537eeaSYohann     return CeedError(ceed, 1, "Basis dimension must be a positive value");
628c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
6294d537eeaSYohann 
6305fe0d4faSjeremylt   if (!ceed->BasisCreateTensorH1) {
6315fe0d4faSjeremylt     Ceed delegate;
632aefd8378Sjeremylt     ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr);
6335fe0d4faSjeremylt 
6345fe0d4faSjeremylt     if (!delegate)
635c042f62fSJeremy L Thompson       // LCOV_EXCL_START
636d7b241e6Sjeremylt       return CeedError(ceed, 1, "Backend does not support BasisCreateTensorH1");
637c042f62fSJeremy L Thompson     // LCOV_EXCL_STOP
6385fe0d4faSjeremylt 
6395fe0d4faSjeremylt     ierr = CeedBasisCreateTensorH1(delegate, dim, ncomp, P1d,
6405fe0d4faSjeremylt                                    Q1d, interp1d, grad1d, qref1d,
6415fe0d4faSjeremylt                                    qweight1d, basis); CeedChk(ierr);
6425fe0d4faSjeremylt     return 0;
6435fe0d4faSjeremylt   }
644d7b241e6Sjeremylt   ierr = CeedCalloc(1,basis); CeedChk(ierr);
645d7b241e6Sjeremylt   (*basis)->ceed = ceed;
646d7b241e6Sjeremylt   ceed->refcount++;
647d7b241e6Sjeremylt   (*basis)->refcount = 1;
648a8de75f0Sjeremylt   (*basis)->tensorbasis = 1;
649d7b241e6Sjeremylt   (*basis)->dim = dim;
650d7b241e6Sjeremylt   (*basis)->ncomp = ncomp;
651d7b241e6Sjeremylt   (*basis)->P1d = P1d;
652d7b241e6Sjeremylt   (*basis)->Q1d = Q1d;
653a8de75f0Sjeremylt   (*basis)->P = CeedIntPow(P1d, dim);
654a8de75f0Sjeremylt   (*basis)->Q = CeedIntPow(Q1d, dim);
655d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d,&(*basis)->qref1d); CeedChk(ierr);
656d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d,&(*basis)->qweight1d); CeedChk(ierr);
657d7b241e6Sjeremylt   memcpy((*basis)->qref1d, qref1d, Q1d*sizeof(qref1d[0]));
658d7b241e6Sjeremylt   memcpy((*basis)->qweight1d, qweight1d, Q1d*sizeof(qweight1d[0]));
659d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d*P1d,&(*basis)->interp1d); CeedChk(ierr);
660d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d*P1d,&(*basis)->grad1d); CeedChk(ierr);
661d7b241e6Sjeremylt   memcpy((*basis)->interp1d, interp1d, Q1d*P1d*sizeof(interp1d[0]));
66209486605Sjeremylt   memcpy((*basis)->grad1d, grad1d, Q1d*P1d*sizeof(grad1d[0]));
663667bc5fcSjeremylt   ierr = ceed->BasisCreateTensorH1(dim, P1d, Q1d, interp1d, grad1d, qref1d,
664d7b241e6Sjeremylt                                    qweight1d, *basis); CeedChk(ierr);
665d7b241e6Sjeremylt   return 0;
666d7b241e6Sjeremylt }
667d7b241e6Sjeremylt 
668b11c1e72Sjeremylt /**
66995bb1877Svaleriabarra   @brief Create a tensor-product Lagrange basis
670b11c1e72Sjeremylt 
671b11c1e72Sjeremylt   @param ceed        A Ceed object where the CeedBasis will be created
672b11c1e72Sjeremylt   @param dim         Topological dimension of element
67395bb1877Svaleriabarra   @param ncomp       Number of field components (1 for scalar fields)
674b11c1e72Sjeremylt   @param P           Number of Gauss-Lobatto nodes in one dimension.  The
675b11c1e72Sjeremylt                        polynomial degree of the resulting Q_k element is k=P-1.
676b11c1e72Sjeremylt   @param Q           Number of quadrature points in one dimension.
677b11c1e72Sjeremylt   @param qmode       Distribution of the Q quadrature points (affects order of
678b11c1e72Sjeremylt                        accuracy for the quadrature)
679b11c1e72Sjeremylt   @param[out] basis  Address of the variable where the newly created
680b11c1e72Sjeremylt                        CeedBasis will be stored.
681b11c1e72Sjeremylt 
682b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
683dfdf5a53Sjeremylt 
684*7a982d89SJeremy L. Thompson   @ref User
685b11c1e72Sjeremylt **/
686d7b241e6Sjeremylt int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt ncomp,
687692c2638Sjeremylt                                     CeedInt P, CeedInt Q, CeedQuadMode qmode,
688692c2638Sjeremylt                                     CeedBasis *basis) {
689d7b241e6Sjeremylt   // Allocate
690d7b241e6Sjeremylt   int ierr, i, j, k;
691d7b241e6Sjeremylt   CeedScalar c1, c2, c3, c4, dx, *nodes, *interp1d, *grad1d, *qref1d, *qweight1d;
6924d537eeaSYohann 
6934d537eeaSYohann   if (dim<1)
694c042f62fSJeremy L Thompson     // LCOV_EXCL_START
6954d537eeaSYohann     return CeedError(ceed, 1, "Basis dimension must be a positive value");
696c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
6974d537eeaSYohann 
698d7b241e6Sjeremylt   ierr = CeedCalloc(P*Q, &interp1d); CeedChk(ierr);
699d7b241e6Sjeremylt   ierr = CeedCalloc(P*Q, &grad1d); CeedChk(ierr);
700d7b241e6Sjeremylt   ierr = CeedCalloc(P, &nodes); CeedChk(ierr);
701d7b241e6Sjeremylt   ierr = CeedCalloc(Q, &qref1d); CeedChk(ierr);
702d7b241e6Sjeremylt   ierr = CeedCalloc(Q, &qweight1d); CeedChk(ierr);
703d7b241e6Sjeremylt   // Get Nodes and Weights
704d7b241e6Sjeremylt   ierr = CeedLobattoQuadrature(P, nodes, NULL); CeedChk(ierr);
705d7b241e6Sjeremylt   switch (qmode) {
706d7b241e6Sjeremylt   case CEED_GAUSS:
707d7b241e6Sjeremylt     ierr = CeedGaussQuadrature(Q, qref1d, qweight1d); CeedChk(ierr);
708d7b241e6Sjeremylt     break;
709d7b241e6Sjeremylt   case CEED_GAUSS_LOBATTO:
710d7b241e6Sjeremylt     ierr = CeedLobattoQuadrature(Q, qref1d, qweight1d); CeedChk(ierr);
711d7b241e6Sjeremylt     break;
712d7b241e6Sjeremylt   }
713d7b241e6Sjeremylt   // Build B, D matrix
714d7b241e6Sjeremylt   // Fornberg, 1998
715d7b241e6Sjeremylt   for (i = 0; i  < Q; i++) {
716d7b241e6Sjeremylt     c1 = 1.0;
717d7b241e6Sjeremylt     c3 = nodes[0] - qref1d[i];
718d7b241e6Sjeremylt     interp1d[i*P+0] = 1.0;
719d7b241e6Sjeremylt     for (j = 1; j < P; j++) {
720d7b241e6Sjeremylt       c2 = 1.0;
721d7b241e6Sjeremylt       c4 = c3;
722d7b241e6Sjeremylt       c3 = nodes[j] - qref1d[i];
723d7b241e6Sjeremylt       for (k = 0; k < j; k++) {
724d7b241e6Sjeremylt         dx = nodes[j] - nodes[k];
725d7b241e6Sjeremylt         c2 *= dx;
726d7b241e6Sjeremylt         if (k == j - 1) {
727d7b241e6Sjeremylt           grad1d[i*P + j] = c1*(interp1d[i*P + k] - c4*grad1d[i*P + k]) / c2;
728d7b241e6Sjeremylt           interp1d[i*P + j] = - c1*c4*interp1d[i*P + k] / c2;
729d7b241e6Sjeremylt         }
730d7b241e6Sjeremylt         grad1d[i*P + k] = (c3*grad1d[i*P + k] - interp1d[i*P + k]) / dx;
731d7b241e6Sjeremylt         interp1d[i*P + k] = c3*interp1d[i*P + k] / dx;
732d7b241e6Sjeremylt       }
733d7b241e6Sjeremylt       c1 = c2;
734d7b241e6Sjeremylt     }
735d7b241e6Sjeremylt   }
736d7b241e6Sjeremylt   //  // Pass to CeedBasisCreateTensorH1
737d7b241e6Sjeremylt   ierr = CeedBasisCreateTensorH1(ceed, dim, ncomp, P, Q, interp1d, grad1d, qref1d,
738d7b241e6Sjeremylt                                  qweight1d, basis); CeedChk(ierr);
739d7b241e6Sjeremylt   ierr = CeedFree(&interp1d); CeedChk(ierr);
740d7b241e6Sjeremylt   ierr = CeedFree(&grad1d); CeedChk(ierr);
741d7b241e6Sjeremylt   ierr = CeedFree(&nodes); CeedChk(ierr);
742d7b241e6Sjeremylt   ierr = CeedFree(&qref1d); CeedChk(ierr);
743d7b241e6Sjeremylt   ierr = CeedFree(&qweight1d); CeedChk(ierr);
744d7b241e6Sjeremylt   return 0;
745d7b241e6Sjeremylt }
746d7b241e6Sjeremylt 
747b11c1e72Sjeremylt /**
74895bb1877Svaleriabarra   @brief Create a non tensor-product basis for H^1 discretizations
749a8de75f0Sjeremylt 
750a8de75f0Sjeremylt   @param ceed        A Ceed object where the CeedBasis will be created
751a8de75f0Sjeremylt   @param topo        Topology of element, e.g. hypercube, simplex, ect
752a8de75f0Sjeremylt   @param ncomp       Number of field components (1 for scalar fields)
7538795c945Sjeremylt   @param nnodes      Total number of nodes
754a8de75f0Sjeremylt   @param nqpts       Total number of quadrature points
75595bb1877Svaleriabarra   @param interp      Row-major (nqpts * nnodes) matrix expressing the values of
7568795c945Sjeremylt                        nodal basis functions at quadrature points
75795bb1877Svaleriabarra   @param grad        Row-major (nqpts * dim * nnodes) matrix expressing
7588795c945Sjeremylt                        derivatives of nodal basis functions at quadrature points
7598795c945Sjeremylt   @param qref        Array of length nqpts holding the locations of quadrature
7608795c945Sjeremylt                        points on the reference element [-1, 1]
761a8de75f0Sjeremylt   @param qweight     Array of length nqpts holding the quadrature weights on the
762a8de75f0Sjeremylt                        reference element
763a8de75f0Sjeremylt   @param[out] basis  Address of the variable where the newly created
764a8de75f0Sjeremylt                        CeedBasis will be stored.
765a8de75f0Sjeremylt 
766a8de75f0Sjeremylt   @return An error code: 0 - success, otherwise - failure
767a8de75f0Sjeremylt 
768*7a982d89SJeremy L. Thompson   @ref User
769a8de75f0Sjeremylt **/
770a8de75f0Sjeremylt int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt ncomp,
771692c2638Sjeremylt                       CeedInt nnodes, CeedInt nqpts, const CeedScalar *interp,
772a8de75f0Sjeremylt                       const CeedScalar *grad, const CeedScalar *qref,
773a8de75f0Sjeremylt                       const CeedScalar *qweight, CeedBasis *basis) {
774a8de75f0Sjeremylt   int ierr;
7758795c945Sjeremylt   CeedInt P = nnodes, Q = nqpts, dim = 0;
776a8de75f0Sjeremylt 
7775fe0d4faSjeremylt   if (!ceed->BasisCreateH1) {
7785fe0d4faSjeremylt     Ceed delegate;
779aefd8378Sjeremylt     ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr);
7805fe0d4faSjeremylt 
7815fe0d4faSjeremylt     if (!delegate)
782c042f62fSJeremy L Thompson       // LCOV_EXCL_START
783a8de75f0Sjeremylt       return CeedError(ceed, 1, "Backend does not support BasisCreateH1");
784c042f62fSJeremy L Thompson     // LCOV_EXCL_STOP
7855fe0d4faSjeremylt 
7868795c945Sjeremylt     ierr = CeedBasisCreateH1(delegate, topo, ncomp, nnodes,
7875fe0d4faSjeremylt                              nqpts, interp, grad, qref,
7885fe0d4faSjeremylt                              qweight, basis); CeedChk(ierr);
7895fe0d4faSjeremylt     return 0;
7905fe0d4faSjeremylt   }
7915fe0d4faSjeremylt 
792a8de75f0Sjeremylt   ierr = CeedCalloc(1,basis); CeedChk(ierr);
793a8de75f0Sjeremylt 
794a8de75f0Sjeremylt   ierr = CeedBasisGetTopologyDimension(topo, &dim); CeedChk(ierr);
795a8de75f0Sjeremylt 
796a8de75f0Sjeremylt   (*basis)->ceed = ceed;
797a8de75f0Sjeremylt   ceed->refcount++;
798a8de75f0Sjeremylt   (*basis)->refcount = 1;
799a8de75f0Sjeremylt   (*basis)->tensorbasis = 0;
800a8de75f0Sjeremylt   (*basis)->dim = dim;
801a8de75f0Sjeremylt   (*basis)->ncomp = ncomp;
802a8de75f0Sjeremylt   (*basis)->P = P;
803a8de75f0Sjeremylt   (*basis)->Q = Q;
804a8de75f0Sjeremylt   ierr = CeedMalloc(Q*dim,&(*basis)->qref1d); CeedChk(ierr);
805a8de75f0Sjeremylt   ierr = CeedMalloc(Q,&(*basis)->qweight1d); CeedChk(ierr);
806a8de75f0Sjeremylt   memcpy((*basis)->qref1d, qref, Q*dim*sizeof(qref[0]));
807a8de75f0Sjeremylt   memcpy((*basis)->qweight1d, qweight, Q*sizeof(qweight[0]));
80800f91b2bSjeremylt   ierr = CeedMalloc(Q*P, &(*basis)->interp); CeedChk(ierr);
80900f91b2bSjeremylt   ierr = CeedMalloc(dim*Q*P, &(*basis)->grad); CeedChk(ierr);
81000f91b2bSjeremylt   memcpy((*basis)->interp, interp, Q*P*sizeof(interp[0]));
81100f91b2bSjeremylt   memcpy((*basis)->grad, grad, dim*Q*P*sizeof(grad[0]));
812667bc5fcSjeremylt   ierr = ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, qref,
813a8de75f0Sjeremylt                              qweight, *basis); CeedChk(ierr);
814a8de75f0Sjeremylt   return 0;
815a8de75f0Sjeremylt }
816a8de75f0Sjeremylt 
817a8de75f0Sjeremylt /**
818*7a982d89SJeremy L. Thompson   @brief View a CeedBasis
819*7a982d89SJeremy L. Thompson 
820*7a982d89SJeremy L. Thompson   @param basis   CeedBasis to view
821*7a982d89SJeremy L. Thompson   @param stream  Stream to view to, e.g., stdout
822*7a982d89SJeremy L. Thompson 
823*7a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
824*7a982d89SJeremy L. Thompson 
825*7a982d89SJeremy L. Thompson   @ref User
826*7a982d89SJeremy L. Thompson **/
827*7a982d89SJeremy L. Thompson int CeedBasisView(CeedBasis basis, FILE *stream) {
828*7a982d89SJeremy L. Thompson   int ierr;
829*7a982d89SJeremy L. Thompson 
830*7a982d89SJeremy L. Thompson   if (basis->tensorbasis) {
831*7a982d89SJeremy L. Thompson     fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P1d,
832*7a982d89SJeremy L. Thompson             basis->Q1d);
833*7a982d89SJeremy L. Thompson     ierr = CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q1d, basis->qref1d,
834*7a982d89SJeremy L. Thompson                           stream); CeedChk(ierr);
835*7a982d89SJeremy L. Thompson     ierr = CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q1d,
836*7a982d89SJeremy L. Thompson                           basis->qweight1d, stream); CeedChk(ierr);
837*7a982d89SJeremy L. Thompson     ierr = CeedScalarView("interp1d", "\t% 12.8f", basis->Q1d, basis->P1d,
838*7a982d89SJeremy L. Thompson                           basis->interp1d, stream); CeedChk(ierr);
839*7a982d89SJeremy L. Thompson     ierr = CeedScalarView("grad1d", "\t% 12.8f", basis->Q1d, basis->P1d,
840*7a982d89SJeremy L. Thompson                           basis->grad1d, stream); CeedChk(ierr);
841*7a982d89SJeremy L. Thompson   } else {
842*7a982d89SJeremy L. Thompson     fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P,
843*7a982d89SJeremy L. Thompson             basis->Q);
844*7a982d89SJeremy L. Thompson     ierr = CeedScalarView("qref", "\t% 12.8f", 1, basis->Q*basis->dim,
845*7a982d89SJeremy L. Thompson                           basis->qref1d,
846*7a982d89SJeremy L. Thompson                           stream); CeedChk(ierr);
847*7a982d89SJeremy L. Thompson     ierr = CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->qweight1d,
848*7a982d89SJeremy L. Thompson                           stream); CeedChk(ierr);
849*7a982d89SJeremy L. Thompson     ierr = CeedScalarView("interp", "\t% 12.8f", basis->Q, basis->P,
850*7a982d89SJeremy L. Thompson                           basis->interp, stream); CeedChk(ierr);
851*7a982d89SJeremy L. Thompson     ierr = CeedScalarView("grad", "\t% 12.8f", basis->dim*basis->Q, basis->P,
852*7a982d89SJeremy L. Thompson                           basis->grad, stream); CeedChk(ierr);
853*7a982d89SJeremy L. Thompson   }
854*7a982d89SJeremy L. Thompson   return 0;
855*7a982d89SJeremy L. Thompson }
856*7a982d89SJeremy L. Thompson 
857*7a982d89SJeremy L. Thompson /**
858*7a982d89SJeremy L. Thompson   @brief Get total number of nodes (in dim dimensions) of a CeedBasis
859*7a982d89SJeremy L. Thompson 
860*7a982d89SJeremy L. Thompson   @param basis   CeedBasis
861*7a982d89SJeremy L. Thompson   @param[out] P  Variable to store number of nodes
862*7a982d89SJeremy L. Thompson 
863*7a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
864*7a982d89SJeremy L. Thompson 
865*7a982d89SJeremy L. Thompson   @ref Utility
866*7a982d89SJeremy L. Thompson **/
867*7a982d89SJeremy L. Thompson int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) {
868*7a982d89SJeremy L. Thompson   *P = basis->P;
869*7a982d89SJeremy L. Thompson   return 0;
870*7a982d89SJeremy L. Thompson }
871*7a982d89SJeremy L. Thompson 
872*7a982d89SJeremy L. Thompson /**
873*7a982d89SJeremy L. Thompson   @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis
874*7a982d89SJeremy L. Thompson 
875*7a982d89SJeremy L. Thompson   @param basis   CeedBasis
876*7a982d89SJeremy L. Thompson   @param[out] Q  Variable to store number of quadrature points
877*7a982d89SJeremy L. Thompson 
878*7a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
879*7a982d89SJeremy L. Thompson 
880*7a982d89SJeremy L. Thompson   @ref Utility
881*7a982d89SJeremy L. Thompson **/
882*7a982d89SJeremy L. Thompson int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) {
883*7a982d89SJeremy L. Thompson   *Q = basis->Q;
884*7a982d89SJeremy L. Thompson   return 0;
885*7a982d89SJeremy L. Thompson }
886*7a982d89SJeremy L. Thompson 
887*7a982d89SJeremy L. Thompson /**
888*7a982d89SJeremy L. Thompson   @brief Apply basis evaluation from nodes to quadrature points or vice versa
889*7a982d89SJeremy L. Thompson 
890*7a982d89SJeremy L. Thompson   @param basis   CeedBasis to evaluate
891*7a982d89SJeremy L. Thompson   @param nelem   The number of elements to apply the basis evaluation to;
892*7a982d89SJeremy L. Thompson                    the backend will specify the ordering in
893*7a982d89SJeremy L. Thompson                    ElemRestrictionCreateBlocked
894*7a982d89SJeremy L. Thompson   @param tmode   \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature
895*7a982d89SJeremy L. Thompson                    points, \ref CEED_TRANSPOSE to apply the transpose, mapping
896*7a982d89SJeremy L. Thompson                    from quadrature points to nodes
897*7a982d89SJeremy L. Thompson   @param emode   \ref CEED_EVAL_NONE to use values directly,
898*7a982d89SJeremy L. Thompson                    \ref CEED_EVAL_INTERP to use interpolated values,
899*7a982d89SJeremy L. Thompson                    \ref CEED_EVAL_GRAD to use gradients,
900*7a982d89SJeremy L. Thompson                    \ref CEED_EVAL_WEIGHT to use quadrature weights.
901*7a982d89SJeremy L. Thompson   @param[in] u   Input CeedVector
902*7a982d89SJeremy L. Thompson   @param[out] v  Output CeedVector
903*7a982d89SJeremy L. Thompson 
904*7a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
905*7a982d89SJeremy L. Thompson 
906*7a982d89SJeremy L. Thompson   @ref User
907*7a982d89SJeremy L. Thompson **/
908*7a982d89SJeremy L. Thompson int CeedBasisApply(CeedBasis basis, CeedInt nelem, CeedTransposeMode tmode,
909*7a982d89SJeremy L. Thompson                    CeedEvalMode emode, CeedVector u, CeedVector v) {
910*7a982d89SJeremy L. Thompson   int ierr;
911*7a982d89SJeremy L. Thompson   CeedInt ulength = 0, vlength, nnodes, nqpt;
912*7a982d89SJeremy L. Thompson   if (!basis->Apply)
913*7a982d89SJeremy L. Thompson     // LCOV_EXCL_START
914*7a982d89SJeremy L. Thompson     return CeedError(basis->ceed, 1, "Backend does not support BasisApply");
915*7a982d89SJeremy L. Thompson   // LCOV_EXCL_STOP
916*7a982d89SJeremy L. Thompson 
917*7a982d89SJeremy L. Thompson   // Check compatibility of topological and geometrical dimensions
918*7a982d89SJeremy L. Thompson   ierr = CeedBasisGetNumNodes(basis, &nnodes); CeedChk(ierr);
919*7a982d89SJeremy L. Thompson   ierr = CeedBasisGetNumQuadraturePoints(basis, &nqpt); CeedChk(ierr);
920*7a982d89SJeremy L. Thompson   ierr = CeedVectorGetLength(v, &vlength); CeedChk(ierr);
921*7a982d89SJeremy L. Thompson 
922*7a982d89SJeremy L. Thompson   if (u) {
923*7a982d89SJeremy L. Thompson     ierr = CeedVectorGetLength(u, &ulength); CeedChk(ierr);
924*7a982d89SJeremy L. Thompson   }
925*7a982d89SJeremy L. Thompson 
926*7a982d89SJeremy L. Thompson   if ((tmode == CEED_TRANSPOSE && (vlength%nnodes != 0 || ulength%nqpt != 0)) ||
927*7a982d89SJeremy L. Thompson       (tmode == CEED_NOTRANSPOSE && (ulength%nnodes != 0 || vlength%nqpt != 0)))
928*7a982d89SJeremy L. Thompson     return CeedError(basis->ceed, 1, "Length of input/output vectors "
929*7a982d89SJeremy L. Thompson                      "incompatible with basis dimensions");
930*7a982d89SJeremy L. Thompson 
931*7a982d89SJeremy L. Thompson   ierr = basis->Apply(basis, nelem, tmode, emode, u, v); CeedChk(ierr);
932*7a982d89SJeremy L. Thompson   return 0;
933*7a982d89SJeremy L. Thompson }
934*7a982d89SJeremy L. Thompson 
935*7a982d89SJeremy L. Thompson /**
936*7a982d89SJeremy L. Thompson   @brief Destroy a CeedBasis
937*7a982d89SJeremy L. Thompson 
938*7a982d89SJeremy L. Thompson   @param basis CeedBasis to destroy
939*7a982d89SJeremy L. Thompson 
940*7a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
941*7a982d89SJeremy L. Thompson 
942*7a982d89SJeremy L. Thompson   @ref User
943*7a982d89SJeremy L. Thompson **/
944*7a982d89SJeremy L. Thompson int CeedBasisDestroy(CeedBasis *basis) {
945*7a982d89SJeremy L. Thompson   int ierr;
946*7a982d89SJeremy L. Thompson 
947*7a982d89SJeremy L. Thompson   if (!*basis || --(*basis)->refcount > 0)
948*7a982d89SJeremy L. Thompson     return 0;
949*7a982d89SJeremy L. Thompson   if ((*basis)->Destroy) {
950*7a982d89SJeremy L. Thompson     ierr = (*basis)->Destroy(*basis); CeedChk(ierr);
951*7a982d89SJeremy L. Thompson   }
952*7a982d89SJeremy L. Thompson   ierr = CeedFree(&(*basis)->interp); CeedChk(ierr);
953*7a982d89SJeremy L. Thompson   ierr = CeedFree(&(*basis)->interp1d); CeedChk(ierr);
954*7a982d89SJeremy L. Thompson   ierr = CeedFree(&(*basis)->grad); CeedChk(ierr);
955*7a982d89SJeremy L. Thompson   ierr = CeedFree(&(*basis)->grad1d); CeedChk(ierr);
956*7a982d89SJeremy L. Thompson   ierr = CeedFree(&(*basis)->qref1d); CeedChk(ierr);
957*7a982d89SJeremy L. Thompson   ierr = CeedFree(&(*basis)->qweight1d); CeedChk(ierr);
958*7a982d89SJeremy L. Thompson   ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr);
959*7a982d89SJeremy L. Thompson   ierr = CeedFree(basis); CeedChk(ierr);
960*7a982d89SJeremy L. Thompson   return 0;
961*7a982d89SJeremy L. Thompson }
962*7a982d89SJeremy L. Thompson 
963*7a982d89SJeremy L. Thompson /**
964b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre quadrature
965b11c1e72Sjeremylt 
966b11c1e72Sjeremylt   @param Q               Number of quadrature points (integrates polynomials of
967b11c1e72Sjeremylt                            degree 2*Q-1 exactly)
968b11c1e72Sjeremylt   @param[out] qref1d     Array of length Q to hold the abscissa on [-1, 1]
969b11c1e72Sjeremylt   @param[out] qweight1d  Array of length Q to hold the weights
970b11c1e72Sjeremylt 
971b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
972dfdf5a53Sjeremylt 
973dfdf5a53Sjeremylt   @ref Utility
974b11c1e72Sjeremylt **/
975d7b241e6Sjeremylt int CeedGaussQuadrature(CeedInt Q, CeedScalar *qref1d, CeedScalar *qweight1d) {
976d7b241e6Sjeremylt   // Allocate
977d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0*atan(1.0);
978d7b241e6Sjeremylt   // Build qref1d, qweight1d
979d7b241e6Sjeremylt   for (int i = 0; i <= Q/2; i++) {
980d7b241e6Sjeremylt     // Guess
981d7b241e6Sjeremylt     xi = cos(PI*(CeedScalar)(2*i+1)/((CeedScalar)(2*Q)));
982d7b241e6Sjeremylt     // Pn(xi)
983d7b241e6Sjeremylt     P0 = 1.0;
984d7b241e6Sjeremylt     P1 = xi;
985d7b241e6Sjeremylt     P2 = 0.0;
986d7b241e6Sjeremylt     for (int j = 2; j <= Q; j++) {
987d7b241e6Sjeremylt       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
988d7b241e6Sjeremylt       P0 = P1;
989d7b241e6Sjeremylt       P1 = P2;
990d7b241e6Sjeremylt     }
991d7b241e6Sjeremylt     // First Newton Step
992d7b241e6Sjeremylt     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
993d7b241e6Sjeremylt     xi = xi-P2/dP2;
994d7b241e6Sjeremylt     // Newton to convergence
9950e4d4210Sjeremylt     for (int k=0; k<100 && fabs(P2)>10*CEED_EPSILON; k++) {
996d7b241e6Sjeremylt       P0 = 1.0;
997d7b241e6Sjeremylt       P1 = xi;
998d7b241e6Sjeremylt       for (int j = 2; j <= Q; j++) {
999d7b241e6Sjeremylt         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1000d7b241e6Sjeremylt         P0 = P1;
1001d7b241e6Sjeremylt         P1 = P2;
1002d7b241e6Sjeremylt       }
1003d7b241e6Sjeremylt       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1004d7b241e6Sjeremylt       xi = xi-P2/dP2;
1005d7b241e6Sjeremylt     }
1006d7b241e6Sjeremylt     // Save xi, wi
1007d7b241e6Sjeremylt     wi = 2.0/((1.0-xi*xi)*dP2*dP2);
1008d7b241e6Sjeremylt     qweight1d[i] = wi;
1009d7b241e6Sjeremylt     qweight1d[Q-1-i] = wi;
1010d7b241e6Sjeremylt     qref1d[i] = -xi;
1011d7b241e6Sjeremylt     qref1d[Q-1-i]= xi;
1012d7b241e6Sjeremylt   }
1013d7b241e6Sjeremylt   return 0;
1014d7b241e6Sjeremylt }
1015d7b241e6Sjeremylt 
1016b11c1e72Sjeremylt /**
1017b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre-Lobatto quadrature
1018b11c1e72Sjeremylt 
1019b11c1e72Sjeremylt   @param Q               Number of quadrature points (integrates polynomials of
1020b11c1e72Sjeremylt                            degree 2*Q-3 exactly)
1021b11c1e72Sjeremylt   @param[out] qref1d     Array of length Q to hold the abscissa on [-1, 1]
1022b11c1e72Sjeremylt   @param[out] qweight1d  Array of length Q to hold the weights
1023b11c1e72Sjeremylt 
1024b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1025dfdf5a53Sjeremylt 
1026dfdf5a53Sjeremylt   @ref Utility
1027b11c1e72Sjeremylt **/
1028d7b241e6Sjeremylt int CeedLobattoQuadrature(CeedInt Q, CeedScalar *qref1d,
1029d7b241e6Sjeremylt                           CeedScalar *qweight1d) {
1030d7b241e6Sjeremylt   // Allocate
1031d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0*atan(1.0);
1032d7b241e6Sjeremylt   // Build qref1d, qweight1d
1033d7b241e6Sjeremylt   // Set endpoints
103430a100c3SJed Brown   if (Q < 2)
10357ed52d01Sjeremylt     return CeedError(NULL, 1,
10367ed52d01Sjeremylt                      "Cannot create Lobatto quadrature with Q=%d < 2 points", Q);
1037d7b241e6Sjeremylt   wi = 2.0/((CeedScalar)(Q*(Q-1)));
1038d7b241e6Sjeremylt   if (qweight1d) {
1039d7b241e6Sjeremylt     qweight1d[0] = wi;
1040d7b241e6Sjeremylt     qweight1d[Q-1] = wi;
1041d7b241e6Sjeremylt   }
1042d7b241e6Sjeremylt   qref1d[0] = -1.0;
1043d7b241e6Sjeremylt   qref1d[Q-1] = 1.0;
1044d7b241e6Sjeremylt   // Interior
1045d7b241e6Sjeremylt   for (int i = 1; i <= (Q-1)/2; i++) {
1046d7b241e6Sjeremylt     // Guess
1047d7b241e6Sjeremylt     xi = cos(PI*(CeedScalar)(i)/(CeedScalar)(Q-1));
1048d7b241e6Sjeremylt     // Pn(xi)
1049d7b241e6Sjeremylt     P0 = 1.0;
1050d7b241e6Sjeremylt     P1 = xi;
1051d7b241e6Sjeremylt     P2 = 0.0;
1052d7b241e6Sjeremylt     for (int j = 2; j < Q; j++) {
1053d7b241e6Sjeremylt       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1054d7b241e6Sjeremylt       P0 = P1;
1055d7b241e6Sjeremylt       P1 = P2;
1056d7b241e6Sjeremylt     }
1057d7b241e6Sjeremylt     // First Newton step
1058d7b241e6Sjeremylt     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1059d7b241e6Sjeremylt     d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
1060d7b241e6Sjeremylt     xi = xi-dP2/d2P2;
1061d7b241e6Sjeremylt     // Newton to convergence
10620e4d4210Sjeremylt     for (int k=0; k<100 && fabs(dP2)>10*CEED_EPSILON; k++) {
1063d7b241e6Sjeremylt       P0 = 1.0;
1064d7b241e6Sjeremylt       P1 = xi;
1065d7b241e6Sjeremylt       for (int j = 2; j < Q; j++) {
1066d7b241e6Sjeremylt         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1067d7b241e6Sjeremylt         P0 = P1;
1068d7b241e6Sjeremylt         P1 = P2;
1069d7b241e6Sjeremylt       }
1070d7b241e6Sjeremylt       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1071d7b241e6Sjeremylt       d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
1072d7b241e6Sjeremylt       xi = xi-dP2/d2P2;
1073d7b241e6Sjeremylt     }
1074d7b241e6Sjeremylt     // Save xi, wi
1075d7b241e6Sjeremylt     wi = 2.0/(((CeedScalar)(Q*(Q-1)))*P2*P2);
1076d7b241e6Sjeremylt     if (qweight1d) {
1077d7b241e6Sjeremylt       qweight1d[i] = wi;
1078d7b241e6Sjeremylt       qweight1d[Q-1-i] = wi;
1079d7b241e6Sjeremylt     }
1080d7b241e6Sjeremylt     qref1d[i] = -xi;
1081d7b241e6Sjeremylt     qref1d[Q-1-i]= xi;
1082d7b241e6Sjeremylt   }
1083d7b241e6Sjeremylt   return 0;
1084d7b241e6Sjeremylt }
1085d7b241e6Sjeremylt 
1086dfdf5a53Sjeremylt /**
108795bb1877Svaleriabarra   @brief Return QR Factorization of a matrix
1088b11c1e72Sjeremylt 
108977645d7bSjeremylt   @param ceed         A Ceed context for error handling
109052bfb9bbSJeremy L Thompson   @param[in,out] mat  Row-major matrix to be factorized in place
109152bfb9bbSJeremy L Thompson   @param[in,out] tau  Vector of length m of scaling factors
1092b11c1e72Sjeremylt   @param m            Number of rows
1093b11c1e72Sjeremylt   @param n            Number of columns
1094b11c1e72Sjeremylt 
1095b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1096dfdf5a53Sjeremylt 
1097dfdf5a53Sjeremylt   @ref Utility
1098b11c1e72Sjeremylt **/
1099a7bd39daSjeremylt int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau,
1100d7b241e6Sjeremylt                         CeedInt m, CeedInt n) {
1101d7b241e6Sjeremylt   CeedScalar v[m];
1102d7b241e6Sjeremylt 
1103a7bd39daSjeremylt   // Check m >= n
1104a7bd39daSjeremylt   if (n > m)
1105c042f62fSJeremy L Thompson     // LCOV_EXCL_START
1106a7bd39daSjeremylt     return CeedError(ceed, 1, "Cannot compute QR factorization with n > m");
1107c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
1108a7bd39daSjeremylt 
110952bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++) {
1110d7b241e6Sjeremylt     // Calculate Householder vector, magnitude
1111d7b241e6Sjeremylt     CeedScalar sigma = 0.0;
1112d7b241e6Sjeremylt     v[i] = mat[i+n*i];
111352bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<m; j++) {
1114d7b241e6Sjeremylt       v[j] = mat[i+n*j];
1115d7b241e6Sjeremylt       sigma += v[j] * v[j];
1116d7b241e6Sjeremylt     }
1117d7b241e6Sjeremylt     CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m]
1118d7b241e6Sjeremylt     CeedScalar Rii = -copysign(norm, v[i]);
1119d7b241e6Sjeremylt     v[i] -= Rii;
1120d7b241e6Sjeremylt     // norm of v[i:m] after modification above and scaling below
1121d7b241e6Sjeremylt     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
1122d7b241e6Sjeremylt     //   tau = 2 / (norm*norm)
1123d7b241e6Sjeremylt     tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma);
1124fb551037Sjeremylt 
11251d102b48SJeremy L Thompson     for (CeedInt j=i+1; j<m; j++)
11261d102b48SJeremy L Thompson       v[j] /= v[i];
1127d7b241e6Sjeremylt 
1128d7b241e6Sjeremylt     // Apply Householder reflector to lower right panel
1129d7b241e6Sjeremylt     CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1);
1130d7b241e6Sjeremylt     // Save v
1131d7b241e6Sjeremylt     mat[i+n*i] = Rii;
11321d102b48SJeremy L Thompson     for (CeedInt j=i+1; j<m; j++)
1133d7b241e6Sjeremylt       mat[i+n*j] = v[j];
1134d7b241e6Sjeremylt   }
1135d7b241e6Sjeremylt 
1136d7b241e6Sjeremylt   return 0;
1137d7b241e6Sjeremylt }
1138d7b241e6Sjeremylt 
1139b11c1e72Sjeremylt /**
114052bfb9bbSJeremy L Thompson   @brief Return symmetric Schur decomposition of the symmetric matrix mat via
114152bfb9bbSJeremy L Thompson            symmetric QR factorization
114252bfb9bbSJeremy L Thompson 
114377645d7bSjeremylt   @param ceed         A Ceed context for error handling
114452bfb9bbSJeremy L Thompson   @param[in,out] mat  Row-major matrix to be factorized in place
1145460bf743SValeria Barra   @param[out] lambda  Vector of length n of eigenvalues
114652bfb9bbSJeremy L Thompson   @param n            Number of rows/columns
114752bfb9bbSJeremy L Thompson 
114852bfb9bbSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
114952bfb9bbSJeremy L Thompson 
115052bfb9bbSJeremy L Thompson   @ref Utility
115152bfb9bbSJeremy L Thompson **/
115252bfb9bbSJeremy L Thompson int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat,
115352bfb9bbSJeremy L Thompson                                     CeedScalar *lambda, CeedInt n) {
115452bfb9bbSJeremy L Thompson   // Check bounds for clang-tidy
115552bfb9bbSJeremy L Thompson   if (n<2)
1156c042f62fSJeremy L Thompson     // LCOV_EXCL_START
1157c042f62fSJeremy L Thompson     return CeedError(ceed, 1,
1158c042f62fSJeremy L Thompson                      "Cannot compute symmetric Schur decomposition of scalars");
1159c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
116052bfb9bbSJeremy L Thompson 
116152bfb9bbSJeremy L Thompson   CeedScalar v[n-1], tau[n-1], matT[n*n];
116252bfb9bbSJeremy L Thompson 
116352bfb9bbSJeremy L Thompson   // Copy mat to matT and set mat to I
116452bfb9bbSJeremy L Thompson   memcpy(matT, mat, n*n*sizeof(mat[0]));
116552bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
116652bfb9bbSJeremy L Thompson     for (CeedInt j=0; j<n; j++)
116752bfb9bbSJeremy L Thompson       mat[j+n*i] = (i==j) ? 1 : 0;
116852bfb9bbSJeremy L Thompson 
116952bfb9bbSJeremy L Thompson   // Reduce to tridiagonal
117052bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n-1; i++) {
117152bfb9bbSJeremy L Thompson     // Calculate Householder vector, magnitude
117252bfb9bbSJeremy L Thompson     CeedScalar sigma = 0.0;
117352bfb9bbSJeremy L Thompson     v[i] = matT[i+n*(i+1)];
117452bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) {
117552bfb9bbSJeremy L Thompson       v[j] = matT[i+n*(j+1)];
117652bfb9bbSJeremy L Thompson       sigma += v[j] * v[j];
117752bfb9bbSJeremy L Thompson     }
117852bfb9bbSJeremy L Thompson     CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:n-1]
117952bfb9bbSJeremy L Thompson     CeedScalar Rii = -copysign(norm, v[i]);
118052bfb9bbSJeremy L Thompson     v[i] -= Rii;
118152bfb9bbSJeremy L Thompson     // norm of v[i:m] after modification above and scaling below
118252bfb9bbSJeremy L Thompson     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
118352bfb9bbSJeremy L Thompson     //   tau = 2 / (norm*norm)
11840e4d4210Sjeremylt     if (sigma > 10*CEED_EPSILON)
118552bfb9bbSJeremy L Thompson       tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma);
1186fb551037Sjeremylt     else
1187fb551037Sjeremylt       tau[i] = 0;
1188fb551037Sjeremylt 
1189fb551037Sjeremylt     for (CeedInt j=i+1; j<n-1; j++)
1190fb551037Sjeremylt       v[j] /= v[i];
119152bfb9bbSJeremy L Thompson 
119252bfb9bbSJeremy L Thompson     // Update sub and super diagonal
119352bfb9bbSJeremy L Thompson     matT[i+n*(i+1)] = Rii;
119452bfb9bbSJeremy L Thompson     matT[(i+1)+n*i] = Rii;
119552bfb9bbSJeremy L Thompson     for (CeedInt j=i+2; j<n; j++) {
119652bfb9bbSJeremy L Thompson       matT[i+n*j] = 0; matT[j+n*i] = 0;
119752bfb9bbSJeremy L Thompson     }
119852bfb9bbSJeremy L Thompson     // Apply symmetric Householder reflector to lower right panel
119952bfb9bbSJeremy L Thompson     CeedHouseholderReflect(&matT[(i+1)+n*(i+1)], &v[i], tau[i],
120052bfb9bbSJeremy L Thompson                            n-(i+1), n-(i+1), n, 1);
120152bfb9bbSJeremy L Thompson     CeedHouseholderReflect(&matT[(i+1)+n*(i+1)], &v[i], tau[i],
120252bfb9bbSJeremy L Thompson                            n-(i+1), n-(i+1), 1, n);
120352bfb9bbSJeremy L Thompson     // Save v
120452bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) {
120552bfb9bbSJeremy L Thompson       matT[i+n*(j+1)] = v[j];
120652bfb9bbSJeremy L Thompson     }
120752bfb9bbSJeremy L Thompson   }
120852bfb9bbSJeremy L Thompson   // Backwards accumulation of Q
120952bfb9bbSJeremy L Thompson   for (CeedInt i=n-2; i>=0; i--) {
121052bfb9bbSJeremy L Thompson     v[i] = 1;
121152bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) {
121252bfb9bbSJeremy L Thompson       v[j] = matT[i+n*(j+1)];
121352bfb9bbSJeremy L Thompson       matT[i+n*(j+1)] = 0;
121452bfb9bbSJeremy L Thompson     }
121552bfb9bbSJeremy L Thompson     CeedHouseholderReflect(&mat[(i+1)+n*(i+1)], &v[i], tau[i],
121652bfb9bbSJeremy L Thompson                            n-(i+1), n-(i+1), n, 1);
121752bfb9bbSJeremy L Thompson   }
121852bfb9bbSJeremy L Thompson 
121952bfb9bbSJeremy L Thompson   // Reduce sub and super diagonal
122052bfb9bbSJeremy L Thompson   CeedInt p = 0, q = 0, itr = 0, maxitr = n*n*n;
12210e4d4210Sjeremylt   CeedScalar tol = 10*CEED_EPSILON;
122252bfb9bbSJeremy L Thompson 
122352bfb9bbSJeremy L Thompson   while (q < n && itr < maxitr) {
122452bfb9bbSJeremy L Thompson     // Update p, q, size of reduced portions of diagonal
122552bfb9bbSJeremy L Thompson     p = 0; q = 0;
122652bfb9bbSJeremy L Thompson     for (CeedInt i=n-2; i>=0; i--) {
122752bfb9bbSJeremy L Thompson       if (fabs(matT[i+n*(i+1)]) < tol)
122852bfb9bbSJeremy L Thompson         q += 1;
122952bfb9bbSJeremy L Thompson       else
123052bfb9bbSJeremy L Thompson         break;
123152bfb9bbSJeremy L Thompson     }
123252bfb9bbSJeremy L Thompson     for (CeedInt i=0; i<n-1-q; i++) {
123352bfb9bbSJeremy L Thompson       if (fabs(matT[i+n*(i+1)]) < tol)
123452bfb9bbSJeremy L Thompson         p += 1;
123552bfb9bbSJeremy L Thompson       else
123652bfb9bbSJeremy L Thompson         break;
123752bfb9bbSJeremy L Thompson     }
123852bfb9bbSJeremy L Thompson     if (q == n-1) break; // Finished reducing
123952bfb9bbSJeremy L Thompson 
124052bfb9bbSJeremy L Thompson     // Reduce tridiagonal portion
124152bfb9bbSJeremy L Thompson     CeedScalar tnn = matT[(n-1-q)+n*(n-1-q)],
124252bfb9bbSJeremy L Thompson                tnnm1 = matT[(n-2-q)+n*(n-1-q)];
124352bfb9bbSJeremy L Thompson     CeedScalar d = (matT[(n-2-q)+n*(n-2-q)] - tnn)/2;
124452bfb9bbSJeremy L Thompson     CeedScalar mu = tnn - tnnm1*tnnm1 /
124552bfb9bbSJeremy L Thompson                     (d + copysign(sqrt(d*d + tnnm1*tnnm1), d));
124652bfb9bbSJeremy L Thompson     CeedScalar x = matT[p+n*p] - mu;
124752bfb9bbSJeremy L Thompson     CeedScalar z = matT[p+n*(p+1)];
124852bfb9bbSJeremy L Thompson     for (CeedInt k=p; k<n-1-q; k++) {
124952bfb9bbSJeremy L Thompson       // Compute Givens rotation
125052bfb9bbSJeremy L Thompson       CeedScalar c = 1, s = 0;
125152bfb9bbSJeremy L Thompson       if (fabs(z) > tol) {
125252bfb9bbSJeremy L Thompson         if (fabs(z) > fabs(x)) {
125352bfb9bbSJeremy L Thompson           CeedScalar tau = -x/z;
125452bfb9bbSJeremy L Thompson           s = 1/sqrt(1+tau*tau), c = s*tau;
125552bfb9bbSJeremy L Thompson         } else {
125652bfb9bbSJeremy L Thompson           CeedScalar tau = -z/x;
125752bfb9bbSJeremy L Thompson           c = 1/sqrt(1+tau*tau), s = c*tau;
125852bfb9bbSJeremy L Thompson         }
125952bfb9bbSJeremy L Thompson       }
126052bfb9bbSJeremy L Thompson 
126152bfb9bbSJeremy L Thompson       // Apply Givens rotation to T
126252bfb9bbSJeremy L Thompson       CeedGivensRotation(matT, c, s, CEED_NOTRANSPOSE, k, k+1, n, n);
126352bfb9bbSJeremy L Thompson       CeedGivensRotation(matT, c, s, CEED_TRANSPOSE, k, k+1, n, n);
126452bfb9bbSJeremy L Thompson 
126552bfb9bbSJeremy L Thompson       // Apply Givens rotation to Q
126652bfb9bbSJeremy L Thompson       CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k+1, n, n);
126752bfb9bbSJeremy L Thompson 
126852bfb9bbSJeremy L Thompson       // Update x, z
126952bfb9bbSJeremy L Thompson       if (k < n-q-2) {
127052bfb9bbSJeremy L Thompson         x = matT[k+n*(k+1)];
127152bfb9bbSJeremy L Thompson         z = matT[k+n*(k+2)];
127252bfb9bbSJeremy L Thompson       }
127352bfb9bbSJeremy L Thompson     }
127452bfb9bbSJeremy L Thompson     itr++;
127552bfb9bbSJeremy L Thompson   }
127652bfb9bbSJeremy L Thompson   // Save eigenvalues
127752bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
127852bfb9bbSJeremy L Thompson     lambda[i] = matT[i+n*i];
127952bfb9bbSJeremy L Thompson 
128052bfb9bbSJeremy L Thompson   // Check convergence
128152bfb9bbSJeremy L Thompson   if (itr == maxitr && q < n-1)
1282c042f62fSJeremy L Thompson     // LCOV_EXCL_START
128352bfb9bbSJeremy L Thompson     return CeedError(ceed, 1, "Symmetric QR failed to converge");
1284c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
128552bfb9bbSJeremy L Thompson 
128652bfb9bbSJeremy L Thompson   return 0;
128752bfb9bbSJeremy L Thompson }
128852bfb9bbSJeremy L Thompson 
128952bfb9bbSJeremy L Thompson /**
129052bfb9bbSJeremy L Thompson   @brief Return Simultaneous Diagonalization of two matrices. This solves the
129152bfb9bbSJeremy L Thompson            generalized eigenvalue problem A x = lambda B x, where A and B
129252bfb9bbSJeremy L Thompson            are symmetric and B is positive definite. We generate the matrix X
129352bfb9bbSJeremy L Thompson            and vector Lambda such that X^T A X = Lambda and X^T B X = I. This
129452bfb9bbSJeremy L Thompson            is equivalent to the LAPACK routine 'sygv' with TYPE = 1.
129552bfb9bbSJeremy L Thompson 
129677645d7bSjeremylt   @param ceed         A Ceed context for error handling
129752bfb9bbSJeremy L Thompson   @param[in] matA     Row-major matrix to be factorized with eigenvalues
129852bfb9bbSJeremy L Thompson   @param[in] matB     Row-major matrix to be factorized to identity
129952bfb9bbSJeremy L Thompson   @param[out] x       Row-major orthogonal matrix
1300460bf743SValeria Barra   @param[out] lambda  Vector of length n of generalized eigenvalues
130152bfb9bbSJeremy L Thompson   @param n            Number of rows/columns
130252bfb9bbSJeremy L Thompson 
130352bfb9bbSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
130452bfb9bbSJeremy L Thompson 
130552bfb9bbSJeremy L Thompson   @ref Utility
130652bfb9bbSJeremy L Thompson **/
130752bfb9bbSJeremy L Thompson int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *matA,
130852bfb9bbSJeremy L Thompson                                     CeedScalar *matB, CeedScalar *x,
130952bfb9bbSJeremy L Thompson                                     CeedScalar *lambda, CeedInt n) {
131052bfb9bbSJeremy L Thompson   int ierr;
131152bfb9bbSJeremy L Thompson   CeedScalar matC[n*n], matG[n*n], vecD[n];
131252bfb9bbSJeremy L Thompson 
131352bfb9bbSJeremy L Thompson   // Compute B = G D G^T
131452bfb9bbSJeremy L Thompson   memcpy(matG, matB, n*n*sizeof(matB[0]));
131552bfb9bbSJeremy L Thompson   ierr = CeedSymmetricSchurDecomposition(ceed, matG, vecD, n); CeedChk(ierr);
1316fb551037Sjeremylt   for (CeedInt i=0; i<n; i++)
1317fb551037Sjeremylt     vecD[i] = sqrt(vecD[i]);
131852bfb9bbSJeremy L Thompson 
1319fb551037Sjeremylt   // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T
1320fb551037Sjeremylt   //           = D^-1/2 G^T A G D^-1/2
132152bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
132252bfb9bbSJeremy L Thompson     for (CeedInt j=0; j<n; j++)
1323fb551037Sjeremylt       matC[j+i*n] = matG[i+j*n] / vecD[i];
13249289e5bfSjeremylt   ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)matC,
13259289e5bfSjeremylt                             (const CeedScalar *)matA, x, n, n, n);
13269289e5bfSjeremylt   CeedChk(ierr);
132752bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
132852bfb9bbSJeremy L Thompson     for (CeedInt j=0; j<n; j++)
1329fb551037Sjeremylt       matG[j+i*n] = matG[j+i*n] / vecD[j];
13309289e5bfSjeremylt   ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)x,
13319289e5bfSjeremylt                             (const CeedScalar *)matG, matC, n, n, n);
13329289e5bfSjeremylt   CeedChk(ierr);
133352bfb9bbSJeremy L Thompson 
133452bfb9bbSJeremy L Thompson   // Compute Q^T C Q = lambda
133552bfb9bbSJeremy L Thompson   ierr = CeedSymmetricSchurDecomposition(ceed, matC, lambda, n); CeedChk(ierr);
133652bfb9bbSJeremy L Thompson 
1337fb551037Sjeremylt   // Set x = (G D^1/2)^-T Q
1338fb551037Sjeremylt   //       = G D^-1/2 Q
13399289e5bfSjeremylt   ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)matG,
13409289e5bfSjeremylt                             (const CeedScalar *)matC, x, n, n, n);
13419289e5bfSjeremylt   CeedChk(ierr);
134252bfb9bbSJeremy L Thompson 
134352bfb9bbSJeremy L Thompson   return 0;
134452bfb9bbSJeremy L Thompson }
134552bfb9bbSJeremy L Thompson 
1346d7b241e6Sjeremylt /// @}
1347