xref: /libCEED/interface/ceed-basis.c (revision 288c044332e33f37503f09b6484fec9d0a55fba1)
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     return CeedError(ceed, 1, "Basis dimension must be a positive value");
65 
66   if (!ceed->BasisCreateTensorH1) {
67     Ceed delegate;
68     ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr);
69 
70     if (!delegate)
71       return CeedError(ceed, 1, "Backend does not support BasisCreateTensorH1");
72 
73     ierr = CeedBasisCreateTensorH1(delegate, dim, ncomp, P1d,
74                                    Q1d, interp1d, grad1d, qref1d,
75                                    qweight1d, basis); CeedChk(ierr);
76     return 0;
77   }
78   ierr = CeedCalloc(1,basis); CeedChk(ierr);
79   (*basis)->ceed = ceed;
80   ceed->refcount++;
81   (*basis)->refcount = 1;
82   (*basis)->tensorbasis = 1;
83   (*basis)->dim = dim;
84   (*basis)->ncomp = ncomp;
85   (*basis)->P1d = P1d;
86   (*basis)->Q1d = Q1d;
87   (*basis)->P = CeedIntPow(P1d, dim);
88   (*basis)->Q = CeedIntPow(Q1d, dim);
89   ierr = CeedMalloc(Q1d,&(*basis)->qref1d); CeedChk(ierr);
90   ierr = CeedMalloc(Q1d,&(*basis)->qweight1d); CeedChk(ierr);
91   memcpy((*basis)->qref1d, qref1d, Q1d*sizeof(qref1d[0]));
92   memcpy((*basis)->qweight1d, qweight1d, Q1d*sizeof(qweight1d[0]));
93   ierr = CeedMalloc(Q1d*P1d,&(*basis)->interp1d); CeedChk(ierr);
94   ierr = CeedMalloc(Q1d*P1d,&(*basis)->grad1d); CeedChk(ierr);
95   memcpy((*basis)->interp1d, interp1d, Q1d*P1d*sizeof(interp1d[0]));
96   memcpy((*basis)->grad1d, grad1d, Q1d*P1d*sizeof(grad1d[0]));
97   ierr = ceed->BasisCreateTensorH1(dim, P1d, Q1d, interp1d, grad1d, qref1d,
98                                    qweight1d, *basis); CeedChk(ierr);
99   return 0;
100 }
101 
102 /**
103   @brief Create a tensor product Lagrange basis
104 
105   @param ceed       A Ceed object where the CeedBasis will be created
106   @param dim        Topological dimension of element
107   @param ncomp      Number of field components
108   @param P          Number of Gauss-Lobatto nodes in one dimension.  The
109                       polynomial degree of the resulting Q_k element is k=P-1.
110   @param Q          Number of quadrature points in one dimension.
111   @param qmode      Distribution of the Q quadrature points (affects order of
112                       accuracy for the quadrature)
113   @param[out] basis Address of the variable where the newly created
114                       CeedBasis will be stored.
115 
116   @return An error code: 0 - success, otherwise - failure
117 
118   @ref Basic
119 **/
120 int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt ncomp,
121                                     CeedInt P, CeedInt Q,
122                                     CeedQuadMode qmode, CeedBasis *basis) {
123   // Allocate
124   int ierr, i, j, k;
125   CeedScalar c1, c2, c3, c4, dx, *nodes, *interp1d, *grad1d, *qref1d, *qweight1d;
126 
127   if (dim<1)
128     return CeedError(ceed, 1, "Basis dimension must be a positive value");
129 
130   ierr = CeedCalloc(P*Q, &interp1d); CeedChk(ierr);
131   ierr = CeedCalloc(P*Q, &grad1d); CeedChk(ierr);
132   ierr = CeedCalloc(P, &nodes); CeedChk(ierr);
133   ierr = CeedCalloc(Q, &qref1d); CeedChk(ierr);
134   ierr = CeedCalloc(Q, &qweight1d); CeedChk(ierr);
135   // Get Nodes and Weights
136   ierr = CeedLobattoQuadrature(P, nodes, NULL); CeedChk(ierr);
137   switch (qmode) {
138   case CEED_GAUSS:
139     ierr = CeedGaussQuadrature(Q, qref1d, qweight1d); CeedChk(ierr);
140     break;
141   case CEED_GAUSS_LOBATTO:
142     ierr = CeedLobattoQuadrature(Q, qref1d, qweight1d); CeedChk(ierr);
143     break;
144   }
145   // Build B, D matrix
146   // Fornberg, 1998
147   for (i = 0; i  < Q; i++) {
148     c1 = 1.0;
149     c3 = nodes[0] - qref1d[i];
150     interp1d[i*P+0] = 1.0;
151     for (j = 1; j < P; j++) {
152       c2 = 1.0;
153       c4 = c3;
154       c3 = nodes[j] - qref1d[i];
155       for (k = 0; k < j; k++) {
156         dx = nodes[j] - nodes[k];
157         c2 *= dx;
158         if (k == j - 1) {
159           grad1d[i*P + j] = c1*(interp1d[i*P + k] - c4*grad1d[i*P + k]) / c2;
160           interp1d[i*P + j] = - c1*c4*interp1d[i*P + k] / c2;
161         }
162         grad1d[i*P + k] = (c3*grad1d[i*P + k] - interp1d[i*P + k]) / dx;
163         interp1d[i*P + k] = c3*interp1d[i*P + k] / dx;
164       }
165       c1 = c2;
166     }
167   }
168   //  // Pass to CeedBasisCreateTensorH1
169   ierr = CeedBasisCreateTensorH1(ceed, dim, ncomp, P, Q, interp1d, grad1d, qref1d,
170                                  qweight1d, basis); CeedChk(ierr);
171   ierr = CeedFree(&interp1d); CeedChk(ierr);
172   ierr = CeedFree(&grad1d); CeedChk(ierr);
173   ierr = CeedFree(&nodes); CeedChk(ierr);
174   ierr = CeedFree(&qref1d); CeedChk(ierr);
175   ierr = CeedFree(&qweight1d); CeedChk(ierr);
176   return 0;
177 }
178 
179 /**
180   @brief Create a non tensor product basis for H^1 discretizations
181 
182   @param ceed       A Ceed object where the CeedBasis will be created
183   @param topo       Topology of element, e.g. hypercube, simplex, ect
184   @param ncomp      Number of field components (1 for scalar fields)
185   @param nnodes       Total number of nodes
186   @param nqpts      Total number of quadrature points
187   @param interp     Row-major nqpts × nnodes matrix expressing the values of
188                       nodal basis functions at quadrature points
189   @param grad       Row-major (nqpts x dim) × nnodes matrix expressing
190                       derivatives of nodal basis functions at quadrature points
191   @param qref       Array of length nqpts holding the locations of quadrature
192                       points on the reference element [-1, 1]
193   @param qweight    Array of length nqpts holding the quadrature weights on the
194                       reference element
195   @param[out] basis Address of the variable where the newly created
196                       CeedBasis will be stored.
197 
198   @return An error code: 0 - success, otherwise - failure
199 
200   @ref Basic
201 **/
202 int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt ncomp,
203                       CeedInt nnodes, CeedInt nqpts,
204                       const CeedScalar *interp,
205                       const CeedScalar *grad, const CeedScalar *qref,
206                       const CeedScalar *qweight, CeedBasis *basis) {
207   int ierr;
208   CeedInt P = nnodes, Q = nqpts, dim = 0;
209 
210   if (!ceed->BasisCreateH1) {
211     Ceed delegate;
212     ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr);
213 
214     if (!delegate)
215       return CeedError(ceed, 1, "Backend does not support BasisCreateH1");
216 
217     ierr = CeedBasisCreateH1(delegate, topo, ncomp, nnodes,
218                              nqpts, interp, grad, qref,
219                              qweight, basis); CeedChk(ierr);
220     return 0;
221   }
222 
223   ierr = CeedCalloc(1,basis); CeedChk(ierr);
224 
225   ierr = CeedBasisGetTopologyDimension(topo, &dim); CeedChk(ierr);
226 
227   (*basis)->ceed = ceed;
228   ceed->refcount++;
229   (*basis)->refcount = 1;
230   (*basis)->tensorbasis = 0;
231   (*basis)->dim = dim;
232   (*basis)->ncomp = ncomp;
233   (*basis)->P = P;
234   (*basis)->Q = Q;
235   ierr = CeedMalloc(Q*dim,&(*basis)->qref1d); CeedChk(ierr);
236   ierr = CeedMalloc(Q,&(*basis)->qweight1d); CeedChk(ierr);
237   memcpy((*basis)->qref1d, qref, Q*dim*sizeof(qref[0]));
238   memcpy((*basis)->qweight1d, qweight, Q*sizeof(qweight[0]));
239   ierr = CeedMalloc(Q*P,&(*basis)->interp1d); CeedChk(ierr);
240   ierr = CeedMalloc(dim*Q*P,&(*basis)->grad1d); CeedChk(ierr);
241   memcpy((*basis)->interp1d, interp, Q*P*sizeof(interp[0]));
242   memcpy((*basis)->grad1d, grad, dim*Q*P*sizeof(grad[0]));
243   ierr = ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, qref,
244                              qweight, *basis); CeedChk(ierr);
245   return 0;
246 }
247 
248 /**
249   @brief Construct a Gauss-Legendre quadrature
250 
251   @param Q              Number of quadrature points (integrates polynomials of
252                           degree 2*Q-1 exactly)
253   @param[out] qref1d    Array of length Q to hold the abscissa on [-1, 1]
254   @param[out] qweight1d Array of length Q to hold the weights
255 
256   @return An error code: 0 - success, otherwise - failure
257 
258   @ref Utility
259 **/
260 int CeedGaussQuadrature(CeedInt Q, CeedScalar *qref1d, CeedScalar *qweight1d) {
261   // Allocate
262   CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0*atan(1.0);
263   // Build qref1d, qweight1d
264   for (int i = 0; i <= Q/2; i++) {
265     // Guess
266     xi = cos(PI*(CeedScalar)(2*i+1)/((CeedScalar)(2*Q)));
267     // Pn(xi)
268     P0 = 1.0;
269     P1 = xi;
270     P2 = 0.0;
271     for (int j = 2; j <= Q; j++) {
272       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
273       P0 = P1;
274       P1 = P2;
275     }
276     // First Newton Step
277     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
278     xi = xi-P2/dP2;
279     // Newton to convergence
280     for (int k=0; k<100 && fabs(P2)>1e-15; k++) {
281       P0 = 1.0;
282       P1 = xi;
283       for (int j = 2; j <= Q; j++) {
284         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
285         P0 = P1;
286         P1 = P2;
287       }
288       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
289       xi = xi-P2/dP2;
290     }
291     // Save xi, wi
292     wi = 2.0/((1.0-xi*xi)*dP2*dP2);
293     qweight1d[i] = wi;
294     qweight1d[Q-1-i] = wi;
295     qref1d[i] = -xi;
296     qref1d[Q-1-i]= xi;
297   }
298   return 0;
299 }
300 
301 /**
302   @brief Construct a Gauss-Legendre-Lobatto quadrature
303 
304   @param Q              Number of quadrature points (integrates polynomials of
305                           degree 2*Q-3 exactly)
306   @param[out] qref1d    Array of length Q to hold the abscissa on [-1, 1]
307   @param[out] qweight1d Array of length Q to hold the weights
308 
309   @return An error code: 0 - success, otherwise - failure
310 
311   @ref Utility
312 **/
313 int CeedLobattoQuadrature(CeedInt Q, CeedScalar *qref1d,
314                           CeedScalar *qweight1d) {
315   // Allocate
316   CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0*atan(1.0);
317   // Build qref1d, qweight1d
318   // Set endpoints
319   wi = 2.0/((CeedScalar)(Q*(Q-1)));
320   if (qweight1d) {
321     qweight1d[0] = wi;
322     qweight1d[Q-1] = wi;
323   }
324   qref1d[0] = -1.0;
325   qref1d[Q-1] = 1.0;
326   // Interior
327   for (int i = 1; i <= (Q-1)/2; i++) {
328     // Guess
329     xi = cos(PI*(CeedScalar)(i)/(CeedScalar)(Q-1));
330     // Pn(xi)
331     P0 = 1.0;
332     P1 = xi;
333     P2 = 0.0;
334     for (int j = 2; j < Q; j++) {
335       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
336       P0 = P1;
337       P1 = P2;
338     }
339     // First Newton step
340     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
341     d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
342     xi = xi-dP2/d2P2;
343     // Newton to convergence
344     for (int k=0; k<100 && fabs(dP2)>1e-15; k++) {
345       P0 = 1.0;
346       P1 = xi;
347       for (int j = 2; j < Q; j++) {
348         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
349         P0 = P1;
350         P1 = P2;
351       }
352       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
353       d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
354       xi = xi-dP2/d2P2;
355     }
356     // Save xi, wi
357     wi = 2.0/(((CeedScalar)(Q*(Q-1)))*P2*P2);
358     if (qweight1d) {
359       qweight1d[i] = wi;
360       qweight1d[Q-1-i] = wi;
361     }
362     qref1d[i] = -xi;
363     qref1d[Q-1-i]= xi;
364   }
365   return 0;
366 }
367 
368 /**
369   @brief View an array stored in a CeedBasis
370 
371   @param name      Name of array
372   @param fpformat  Printing format
373   @param m         Number of rows in array
374   @param n         Number of columns in array
375   @param a         Array to be viewed
376   @param stream    Stream to view to, e.g., stdout
377 
378   @return An error code: 0 - success, otherwise - failure
379 
380   @ref Utility
381 **/
382 static int CeedScalarView(const char *name, const char *fpformat, CeedInt m,
383                           CeedInt n, const CeedScalar *a, FILE *stream) {
384   for (int i=0; i<m; i++) {
385     if (m > 1) fprintf(stream, "%12s[%d]:", name, i);
386     else fprintf(stream, "%12s:", name);
387     for (int j=0; j<n; j++) {
388       fprintf(stream, fpformat, fabs(a[i*n+j]) > 1E-14 ? a[i*n+j] : 0);
389     }
390     fputs("\n", stream);
391   }
392   return 0;
393 }
394 
395 /**
396   @brief View a CeedBasis
397 
398   @param basis  CeedBasis to view
399   @param stream Stream to view to, e.g., stdout
400 
401   @return An error code: 0 - success, otherwise - failure
402 
403   @ref Utility
404 **/
405 int CeedBasisView(CeedBasis basis, FILE *stream) {
406   int ierr;
407 
408   if (basis->tensorbasis) {
409     fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P1d,
410             basis->Q1d);
411     ierr = CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q1d, basis->qref1d,
412                           stream); CeedChk(ierr);
413     ierr = CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q1d,
414                           basis->qweight1d, stream); CeedChk(ierr);
415     ierr = CeedScalarView("interp1d", "\t% 12.8f", basis->Q1d, basis->P1d,
416                           basis->interp1d, stream); CeedChk(ierr);
417     ierr = CeedScalarView("grad1d", "\t% 12.8f", basis->Q1d, basis->P1d,
418                           basis->grad1d, stream); CeedChk(ierr);
419   } else {
420     fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P,
421             basis->Q);
422     ierr = CeedScalarView("qref", "\t% 12.8f", 1, basis->Q*basis->dim,
423                           basis->qref1d,
424                           stream); CeedChk(ierr);
425     ierr = CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->qweight1d,
426                           stream); CeedChk(ierr);
427     ierr = CeedScalarView("interp", "\t% 12.8f", basis->Q, basis->P,
428                           basis->interp1d, stream); CeedChk(ierr);
429     ierr = CeedScalarView("grad", "\t% 12.8f", basis->dim*basis->Q, basis->P,
430                           basis->grad1d, stream); CeedChk(ierr);
431   }
432   return 0;
433 }
434 
435 /**
436   @brief Compute Householder Reflection
437 
438     Computes A = (I - b v v^T) A
439     where A is an mxn matrix indexed as A[i*row + j*col]
440 
441   @param[out] A  Matrix to apply Householder reflection to, in place
442   @param v       Householder vector
443   @param b       Scaling factor
444   @param m       Number of rows in A
445   @param n       Number of columns in A
446   @param row     Col stride
447   @param col     Row stride
448 
449   @return An error code: 0 - success, otherwise - failure
450 
451   @ref Developer
452 **/
453 static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v,
454                                   CeedScalar b, CeedInt m, CeedInt n,
455                                   CeedInt row, CeedInt col) {
456   for (CeedInt j=0; j<n; j++) {
457     CeedScalar w = A[0*row + j*col];
458     for (CeedInt i=1; i<m; i++) w += v[i] * A[i*row + j*col];
459     A[0*row + j*col] -= b * w;
460     for (CeedInt i=1; i<m; i++) A[i*row + j*col] -= b * w * v[i];
461   }
462   return 0;
463 }
464 
465 /**
466   @brief Apply Householder Q matrix
467 
468     Compute A = Q A where Q is mxk and A is mxn. k<m
469 
470   @param[out] A  Matrix to apply Householder Q to, in place
471   @param Q       Householder Q matrix
472   @param tau     Householder scaling factors
473   @param tmode   Transpose mode for application
474   @param m       Number of rows in A
475   @param n       Number of columns in A
476   @param k       Index of row targeted
477   @param row     Col stride
478   @param col     Row stride
479 
480   @return An error code: 0 - success, otherwise - failure
481 
482   @ref Developer
483 **/
484 static int CeedHouseholderApplyQ(CeedScalar *A, const CeedScalar *Q,
485                                  const CeedScalar *tau, CeedTransposeMode tmode,
486                                  CeedInt m, CeedInt n, CeedInt k,
487                                  CeedInt row, CeedInt col) {
488   CeedScalar v[m];
489   for (CeedInt ii=0; ii<k; ii++) {
490     CeedInt i = tmode == CEED_TRANSPOSE ? ii : k-1-ii;
491     for (CeedInt j=i+1; j<m; j++) {
492       v[j] = Q[j*k+i];
493     }
494     // Apply Householder reflector (I - tau v v^T) colograd1d^T
495     CeedHouseholderReflect(&A[i*row], &v[i], tau[i], m-i, n, row, col);
496   }
497   return 0;
498 }
499 
500 /**
501   @brief Return QR Factorization of matrix
502 
503   @param ceed      A Ceed object currently in use
504   @param[out] mat  Row-major matrix to be factorized in place
505   @param[out] tau  Vector of length m of scaling factors
506   @param m         Number of rows
507   @param n         Number of columns
508 
509   @return An error code: 0 - success, otherwise - failure
510 
511   @ref Utility
512 **/
513 int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau,
514                         CeedInt m, CeedInt n) {
515   CeedInt i, j;
516   CeedScalar v[m];
517 
518   // Check m >= n
519   if (n > m)
520     return CeedError(ceed, 1, "Cannot compute QR factorization with n > m");
521 
522   for (i=0; i<n; i++) {
523     // Calculate Householder vector, magnitude
524     CeedScalar sigma = 0.0;
525     v[i] = mat[i+n*i];
526     for (j=i+1; j<m; j++) {
527       v[j] = mat[i+n*j];
528       sigma += v[j] * v[j];
529     }
530     CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m]
531     CeedScalar Rii = -copysign(norm, v[i]);
532     v[i] -= Rii;
533     // norm of v[i:m] after modification above and scaling below
534     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
535     //   tau = 2 / (norm*norm)
536     tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma);
537     for (j=i+1; j<m; j++) v[j] /= v[i];
538 
539     // Apply Householder reflector to lower right panel
540     CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1);
541     // Save v
542     mat[i+n*i] = Rii;
543     for (j=i+1; j<m; j++) {
544       mat[i+n*j] = v[j];
545     }
546   }
547 
548   return 0;
549 }
550 
551 /**
552   @brief Return collocated grad matrix
553 
554   @param basis           CeedBasis
555   @param[out] colograd1d Row-major Q1d × Q1d matrix expressing derivatives of
556                            basis functions at quadrature points
557 
558   @return An error code: 0 - success, otherwise - failure
559 
560   @ref Advanced
561 **/
562 int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *colograd1d) {
563   int i, j, k;
564   Ceed ceed;
565   CeedInt ierr, P1d=(basis)->P1d, Q1d=(basis)->Q1d;
566   CeedScalar *interp1d, *grad1d, tau[Q1d];
567 
568   ierr = CeedMalloc(Q1d*P1d, &interp1d); CeedChk(ierr);
569   ierr = CeedMalloc(Q1d*P1d, &grad1d); CeedChk(ierr);
570   memcpy(interp1d, (basis)->interp1d, Q1d*P1d*sizeof(basis)->interp1d[0]);
571   memcpy(grad1d, (basis)->grad1d, Q1d*P1d*sizeof(basis)->interp1d[0]);
572 
573   // QR Factorization, interp1d = Q R
574   ierr = CeedBasisGetCeed(basis, &ceed); CeedChk(ierr);
575   ierr = CeedQRFactorization(ceed, interp1d, tau, Q1d, P1d); CeedChk(ierr);
576 
577   // Apply Rinv, colograd1d = grad1d Rinv
578   for (i=0; i<Q1d; i++) { // Row i
579     colograd1d[Q1d*i] = grad1d[P1d*i]/interp1d[0];
580     for (j=1; j<P1d; j++) { // Column j
581       colograd1d[j+Q1d*i] = grad1d[j+P1d*i];
582       for (k=0; k<j; k++) {
583         colograd1d[j+Q1d*i] -= interp1d[j+P1d*k]*colograd1d[k+Q1d*i];
584       }
585       colograd1d[j+Q1d*i] /= interp1d[j+P1d*j];
586     }
587     for (j=P1d; j<Q1d; j++) {
588       colograd1d[j+Q1d*i] = 0;
589     }
590   }
591 
592   // Apply Qtranspose, colograd = colograd Qtranspose
593   CeedHouseholderApplyQ(colograd1d, interp1d, tau, CEED_NOTRANSPOSE,
594                         Q1d, Q1d, P1d, 1, Q1d);
595 
596   ierr = CeedFree(&interp1d); CeedChk(ierr);
597   ierr = CeedFree(&grad1d); CeedChk(ierr);
598 
599   return 0;
600 }
601 
602 /**
603   @brief Apply basis evaluation from nodes to quadrature points or vice-versa
604 
605   @param basis  CeedBasis to evaluate
606   @param nelem  The number of elements to apply the basis evaluation to;
607                   the backend will specify the ordering in
608                   ElemRestrictionCreateBlocked
609   @param tmode  \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature
610                   points, \ref CEED_TRANSPOSE to apply the transpose, mapping
611                   from quadrature points to nodes
612   @param emode  \ref CEED_EVAL_INTERP to obtain interpolated values,
613                   \ref CEED_EVAL_GRAD to obtain gradients.
614   @param[in] u  Input array
615   @param[out] v Output array
616 
617   @return An error code: 0 - success, otherwise - failure
618 
619   @ref Advanced
620 **/
621 int CeedBasisApply(CeedBasis basis, CeedInt nelem, CeedTransposeMode tmode,
622                    CeedEvalMode emode, CeedVector u, CeedVector v) {
623   int ierr;
624   CeedInt ulength = 0, vlength, nnodes, nqpt;
625   if (!basis->Apply) return CeedError(basis->ceed, 1,
626                                         "Backend does not support BasisApply");
627   // check compatibility of topological and geometrical dimensions
628   ierr = CeedBasisGetNumNodes(basis, &nnodes); CeedChk(ierr);
629   ierr = CeedBasisGetNumQuadraturePoints(basis, &nqpt); CeedChk(ierr);
630   ierr = CeedVectorGetLength(v, &vlength); CeedChk(ierr);
631 
632   if (u) {
633     ierr = CeedVectorGetLength(u, &ulength); CeedChk(ierr);
634   }
635 
636   if ((tmode == CEED_TRANSPOSE   && (vlength % nnodes != 0
637                                      || ulength % nqpt != 0))
638       ||
639       (tmode == CEED_NOTRANSPOSE && (ulength % nnodes != 0 || vlength % nqpt != 0)))
640     return CeedError(basis->ceed, 1,
641                      "Length of input/output vectors incompatible with basis dimensions");
642 
643   ierr = basis->Apply(basis, nelem, tmode, emode, u, v); CeedChk(ierr);
644   return 0;
645 }
646 
647 /**
648   @brief Get Ceed associated with a CeedBasis
649 
650   @param basis      CeedBasis
651   @param[out] ceed  Variable to store Ceed
652 
653   @return An error code: 0 - success, otherwise - failure
654 
655   @ref Advanced
656 **/
657 int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) {
658   *ceed = basis->ceed;
659 
660   return 0;
661 };
662 
663 /**
664   @brief Get dimension for given CeedBasis
665 
666   @param basis     CeedBasis
667   @param[out] dim  Variable to store dimension of basis
668 
669   @return An error code: 0 - success, otherwise - failure
670 
671   @ref Advanced
672 **/
673 int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) {
674   *dim = basis->dim;
675 
676   return 0;
677 };
678 
679 /**
680   @brief Get tensor status for given CeedBasis
681 
682   @param basis        CeedBasis
683   @param[out] tensor  Variable to store tensor status
684 
685   @return An error code: 0 - success, otherwise - failure
686 
687   @ref Advanced
688 **/
689 int CeedBasisGetTensorStatus(CeedBasis basis, bool *tensor) {
690   *tensor = basis->tensorbasis;
691 
692   return 0;
693 };
694 
695 /**
696   @brief Get number of components for given CeedBasis
697 
698   @param basis        CeedBasis
699   @param[out] numcomp Variable to store number of components of basis
700 
701   @return An error code: 0 - success, otherwise - failure
702 
703   @ref Advanced
704 **/
705 int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *numcomp) {
706   *numcomp = basis->ncomp;
707 
708   return 0;
709 };
710 
711 /**
712   @brief Get total number of nodes (in 1 dimension) of a CeedBasis
713 
714   @param basis     CeedBasis
715   @param[out] P1d  Variable to store number of nodes
716 
717   @return An error code: 0 - success, otherwise - failure
718 
719   @ref Advanced
720 **/
721 int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P1d) {
722   if (!basis->tensorbasis) return CeedError(basis->ceed, 1,
723                                     "Cannot supply P1d for non-tensor basis");
724   *P1d = basis->P1d;
725   return 0;
726 }
727 
728 /**
729   @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis
730 
731   @param basis     CeedBasis
732   @param[out] Q1d  Variable to store number of quadrature points
733 
734   @return An error code: 0 - success, otherwise - failure
735 
736   @ref Advanced
737 **/
738 int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q1d) {
739   if (!basis->tensorbasis) return CeedError(basis->ceed, 1,
740                                     "Cannot supply Q1d for non-tensor basis");
741   *Q1d = basis->Q1d;
742   return 0;
743 }
744 
745 /**
746   @brief Get total number of nodes (in dim dimensions) of a CeedBasis
747 
748   @param basis   CeedBasis
749   @param[out] P  Variable to store number of nodes
750 
751   @return An error code: 0 - success, otherwise - failure
752 
753   @ref Utility
754 **/
755 int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) {
756   *P = basis->P;
757   return 0;
758 }
759 
760 /**
761   @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis
762 
763   @param basis   CeedBasis
764   @param[out] Q  Variable to store number of quadrature points
765 
766   @return An error code: 0 - success, otherwise - failure
767 
768   @ref Utility
769 **/
770 int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) {
771   *Q = basis->Q;
772   return 0;
773 }
774 
775 /**
776   @brief Get reference coordinates of quadrature points (in dim dimensions)
777          of a CeedBasis
778 
779   @param basis      CeedBasis
780   @param[out] qref  Variable to store reference coordinates of quadrature points
781 
782   @return An error code: 0 - success, otherwise - failure
783 
784   @ref Advanced
785 **/
786 int CeedBasisGetQRef(CeedBasis basis, CeedScalar* *qref) {
787   *qref = basis->qref1d;
788   return 0;
789 }
790 
791 /**
792   @brief Get quadrature weights of quadrature points (in dim dimensions)
793          of a CeedBasis
794 
795   @param basis         CeedBasis
796   @param[out] qweight  Variable to store quadrature weights
797 
798   @return An error code: 0 - success, otherwise - failure
799 
800   @ref Advanced
801 **/
802 int CeedBasisGetQWeights(CeedBasis basis, CeedScalar* *qweight) {
803   *qweight = basis->qweight1d;
804   return 0;
805 }
806 
807 /**
808   @brief Get interpolation matrix of a CeedBasis
809 
810   @param basis       CeedBasis
811   @param[out] interp Variable to store interpolation matrix
812 
813   @return An error code: 0 - success, otherwise - failure
814 
815   @ref Advanced
816 **/
817 int CeedBasisGetInterp(CeedBasis basis, CeedScalar* *interp) {
818   *interp = basis->interp1d;
819   return 0;
820 }
821 
822 /**
823   @brief Get gradient matrix of a CeedBasis
824 
825   @param basis      CeedBasis
826   @param[out] grad  Variable to store gradient matrix
827 
828   @return An error code: 0 - success, otherwise - failure
829 
830   @ref Advanced
831 **/
832 int CeedBasisGetGrad(CeedBasis basis, CeedScalar* *grad) {
833   *grad = basis->grad1d;
834   return 0;
835 }
836 
837 /**
838   @brief Get backend data of a CeedBasis
839 
840   @param basis      CeedBasis
841   @param[out] data  Variable to store data
842 
843   @return An error code: 0 - success, otherwise - failure
844 
845   @ref Advanced
846 **/
847 int CeedBasisGetData(CeedBasis basis, void* *data) {
848   *data = basis->data;
849   return 0;
850 }
851 
852 /**
853   @brief Set backend data of a CeedBasis
854 
855   @param[out] basis CeedBasis
856   @param data       Data to set
857 
858   @return An error code: 0 - success, otherwise - failure
859 
860   @ref Advanced
861 **/
862 int CeedBasisSetData(CeedBasis basis, void* *data) {
863   basis->data = *data;
864   return 0;
865 }
866 
867 /**
868   @brief Get CeedTensorContract of a CeedBasis
869 
870   @param basis          CeedBasis
871   @param[out] contract  Variable to store CeedTensorContract
872 
873   @return An error code: 0 - success, otherwise - failure
874 
875   @ref Advanced
876 **/
877 int CeedBasisGetTensorContract(CeedBasis basis,
878                                CeedTensorContract *contract) {
879   *contract = basis->contract;
880   return 0;
881 }
882 
883 /**
884   @brief Set CeedTensorContract of a CeedBasis
885 
886   @param[out] basis     CeedBasis
887   @param contract       CeedTensorContract to set
888 
889   @return An error code: 0 - success, otherwise - failure
890 
891   @ref Advanced
892 **/
893 int CeedBasisSetTensorContract(CeedBasis basis,
894                                CeedTensorContract *contract) {
895   basis->contract = *contract;
896   return 0;
897 }
898 
899 /**
900   @brief Get dimension for given CeedElemTopology
901 
902   @param topo      CeedElemTopology
903   @param[out] dim  Variable to store dimension of topology
904 
905   @return An error code: 0 - success, otherwise - failure
906 
907   @ref Advanced
908 **/
909 int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) {
910   *dim = (CeedInt) topo >> 16;
911 
912   return 0;
913 };
914 
915 /**
916   @brief Destroy a CeedBasis
917 
918   @param basis CeedBasis to destroy
919 
920   @return An error code: 0 - success, otherwise - failure
921 
922   @ref Basic
923 **/
924 int CeedBasisDestroy(CeedBasis *basis) {
925   int ierr;
926 
927   if (!*basis || --(*basis)->refcount > 0) return 0;
928   if ((*basis)->Destroy) {
929     ierr = (*basis)->Destroy(*basis); CeedChk(ierr);
930   }
931   ierr = CeedFree(&(*basis)->interp1d); CeedChk(ierr);
932   ierr = CeedFree(&(*basis)->grad1d); CeedChk(ierr);
933   ierr = CeedFree(&(*basis)->qref1d); CeedChk(ierr);
934   ierr = CeedFree(&(*basis)->qweight1d); CeedChk(ierr);
935   ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr);
936   ierr = CeedFree(basis); CeedChk(ierr);
937   return 0;
938 }
939 
940 /// @cond DOXYGEN_SKIP
941 // Indicate that the quadrature points are collocated with the nodes
942 CeedBasis CEED_BASIS_COLLOCATED = &ceed_basis_collocated;
943 /// @endcond
944 /// @}
945