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