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