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