xref: /libCEED/rust/libceed-sys/c-src/interface/ceed-basis.c (revision fe2413ff2a5f6e0d0808f191532abb277c4b8bc7)
1d7b241e6Sjeremylt // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2d7b241e6Sjeremylt // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3d7b241e6Sjeremylt // reserved. See files LICENSE and NOTICE for details.
4d7b241e6Sjeremylt //
5d7b241e6Sjeremylt // This file is part of CEED, a collection of benchmarks, miniapps, software
6d7b241e6Sjeremylt // libraries and APIs for efficient high-order finite element and spectral
7d7b241e6Sjeremylt // element discretizations for exascale applications. For more information and
8d7b241e6Sjeremylt // source code availability see http://github.com/ceed.
9d7b241e6Sjeremylt //
10d7b241e6Sjeremylt // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11d7b241e6Sjeremylt // a collaborative effort of two U.S. Department of Energy organizations (Office
12d7b241e6Sjeremylt // of Science and the National Nuclear Security Administration) responsible for
13d7b241e6Sjeremylt // the planning and preparation of a capable exascale ecosystem, including
14d7b241e6Sjeremylt // software, applications, hardware, advanced system engineering and early
15d7b241e6Sjeremylt // testbed platforms, in support of the nation's exascale computing imperative.
16d7b241e6Sjeremylt 
17d7b241e6Sjeremylt #include <ceed-impl.h>
18d863ab9bSjeremylt #include <ceed-backend.h>
19d7b241e6Sjeremylt #include <math.h>
20d7b241e6Sjeremylt #include <stdio.h>
21d7b241e6Sjeremylt #include <stdlib.h>
22d7b241e6Sjeremylt #include <string.h>
23d7b241e6Sjeremylt 
24d7b241e6Sjeremylt /// @cond DOXYGEN_SKIP
25783c99b3SValeria Barra static struct CeedBasis_private ceed_basis_collocated;
26d7b241e6Sjeremylt /// @endcond
27d7b241e6Sjeremylt 
28d7b241e6Sjeremylt /// @file
29d7b241e6Sjeremylt /// Implementation of public CeedBasis interfaces
30d7b241e6Sjeremylt ///
31dfdf5a53Sjeremylt /// @addtogroup CeedBasis
32d7b241e6Sjeremylt /// @{
33d7b241e6Sjeremylt 
34b11c1e72Sjeremylt /**
35b11c1e72Sjeremylt   @brief Create a tensor product basis for H^1 discretizations
36b11c1e72Sjeremylt 
37b11c1e72Sjeremylt   @param ceed       A Ceed object where the CeedBasis will be created
38b11c1e72Sjeremylt   @param dim        Topological dimension
39b11c1e72Sjeremylt   @param ncomp      Number of field components (1 for scalar fields)
40b11c1e72Sjeremylt   @param P1d        Number of nodes in one dimension
41b11c1e72Sjeremylt   @param Q1d        Number of quadrature points in one dimension
42b11c1e72Sjeremylt   @param interp1d   Row-major Q1d × P1d matrix expressing the values of nodal
43b11c1e72Sjeremylt                       basis functions at quadrature points
44b11c1e72Sjeremylt   @param grad1d     Row-major Q1d × P1d matrix expressing derivatives of nodal
45b11c1e72Sjeremylt                       basis functions at quadrature points
46b11c1e72Sjeremylt   @param qref1d     Array of length Q1d holding the locations of quadrature points
47b11c1e72Sjeremylt                       on the 1D reference element [-1, 1]
48b11c1e72Sjeremylt   @param qweight1d  Array of length Q1d holding the quadrature weights on the
49b11c1e72Sjeremylt                       reference element
50b11c1e72Sjeremylt   @param[out] basis Address of the variable where the newly created
51b11c1e72Sjeremylt                       CeedBasis will be stored.
52b11c1e72Sjeremylt 
53b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
54dfdf5a53Sjeremylt 
55dfdf5a53Sjeremylt   @ref Basic
56b11c1e72Sjeremylt **/
57d7b241e6Sjeremylt int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt ncomp, CeedInt P1d,
58d7b241e6Sjeremylt                             CeedInt Q1d, const CeedScalar *interp1d,
59d7b241e6Sjeremylt                             const CeedScalar *grad1d, const CeedScalar *qref1d,
60d7b241e6Sjeremylt                             const CeedScalar *qweight1d, CeedBasis *basis) {
61d7b241e6Sjeremylt   int ierr;
62d7b241e6Sjeremylt 
635fe0d4faSjeremylt   if (!ceed->BasisCreateTensorH1) {
645fe0d4faSjeremylt     Ceed delegate;
655fe0d4faSjeremylt     ierr = CeedGetDelegate(ceed, &delegate); CeedChk(ierr);
665fe0d4faSjeremylt 
675fe0d4faSjeremylt     if (!delegate)
68d7b241e6Sjeremylt       return CeedError(ceed, 1, "Backend does not support BasisCreateTensorH1");
695fe0d4faSjeremylt 
705fe0d4faSjeremylt     ierr = CeedBasisCreateTensorH1(delegate, dim, ncomp, P1d,
715fe0d4faSjeremylt                             Q1d, interp1d, grad1d, qref1d,
725fe0d4faSjeremylt                             qweight1d, basis); CeedChk(ierr);
735fe0d4faSjeremylt     return 0;
745fe0d4faSjeremylt   }
75d7b241e6Sjeremylt   ierr = CeedCalloc(1,basis); CeedChk(ierr);
76d7b241e6Sjeremylt   (*basis)->ceed = ceed;
77d7b241e6Sjeremylt   ceed->refcount++;
78d7b241e6Sjeremylt   (*basis)->refcount = 1;
79a8de75f0Sjeremylt   (*basis)->tensorbasis = 1;
80d7b241e6Sjeremylt   (*basis)->dim = dim;
81d7b241e6Sjeremylt   (*basis)->ncomp = ncomp;
82d7b241e6Sjeremylt   (*basis)->P1d = P1d;
83d7b241e6Sjeremylt   (*basis)->Q1d = Q1d;
84a8de75f0Sjeremylt   (*basis)->P = CeedIntPow(P1d, dim);
85a8de75f0Sjeremylt   (*basis)->Q = CeedIntPow(Q1d, dim);
86d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d,&(*basis)->qref1d); CeedChk(ierr);
87d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d,&(*basis)->qweight1d); CeedChk(ierr);
88d7b241e6Sjeremylt   memcpy((*basis)->qref1d, qref1d, Q1d*sizeof(qref1d[0]));
89d7b241e6Sjeremylt   memcpy((*basis)->qweight1d, qweight1d, Q1d*sizeof(qweight1d[0]));
90d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d*P1d,&(*basis)->interp1d); CeedChk(ierr);
91d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d*P1d,&(*basis)->grad1d); CeedChk(ierr);
92d7b241e6Sjeremylt   memcpy((*basis)->interp1d, interp1d, Q1d*P1d*sizeof(interp1d[0]));
9309486605Sjeremylt   memcpy((*basis)->grad1d, grad1d, Q1d*P1d*sizeof(grad1d[0]));
94667bc5fcSjeremylt   ierr = ceed->BasisCreateTensorH1(dim, P1d, Q1d, interp1d, grad1d, qref1d,
95d7b241e6Sjeremylt                                    qweight1d, *basis); CeedChk(ierr);
96d7b241e6Sjeremylt   return 0;
97d7b241e6Sjeremylt }
98d7b241e6Sjeremylt 
99b11c1e72Sjeremylt /**
100b11c1e72Sjeremylt   @brief Create a tensor product Lagrange basis
101b11c1e72Sjeremylt 
102b11c1e72Sjeremylt   @param ceed       A Ceed object where the CeedBasis will be created
103b11c1e72Sjeremylt   @param dim        Topological dimension of element
104b11c1e72Sjeremylt   @param ncomp      Number of field components
105b11c1e72Sjeremylt   @param P          Number of Gauss-Lobatto nodes in one dimension.  The
106b11c1e72Sjeremylt                       polynomial degree of the resulting Q_k element is k=P-1.
107b11c1e72Sjeremylt   @param Q          Number of quadrature points in one dimension.
108b11c1e72Sjeremylt   @param qmode      Distribution of the Q quadrature points (affects order of
109b11c1e72Sjeremylt                       accuracy for the quadrature)
110b11c1e72Sjeremylt   @param[out] basis Address of the variable where the newly created
111b11c1e72Sjeremylt                       CeedBasis will be stored.
112b11c1e72Sjeremylt 
113b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
114dfdf5a53Sjeremylt 
115dfdf5a53Sjeremylt   @ref Basic
116b11c1e72Sjeremylt **/
117d7b241e6Sjeremylt int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt ncomp,
118d7b241e6Sjeremylt                                     CeedInt P, CeedInt Q,
119d7b241e6Sjeremylt                                     CeedQuadMode qmode, CeedBasis *basis) {
120d7b241e6Sjeremylt   // Allocate
121d7b241e6Sjeremylt   int ierr, i, j, k;
122d7b241e6Sjeremylt   CeedScalar c1, c2, c3, c4, dx, *nodes, *interp1d, *grad1d, *qref1d, *qweight1d;
123d7b241e6Sjeremylt   ierr = CeedCalloc(P*Q, &interp1d); CeedChk(ierr);
124d7b241e6Sjeremylt   ierr = CeedCalloc(P*Q, &grad1d); CeedChk(ierr);
125d7b241e6Sjeremylt   ierr = CeedCalloc(P, &nodes); CeedChk(ierr);
126d7b241e6Sjeremylt   ierr = CeedCalloc(Q, &qref1d); CeedChk(ierr);
127d7b241e6Sjeremylt   ierr = CeedCalloc(Q, &qweight1d); CeedChk(ierr);
128d7b241e6Sjeremylt   // Get Nodes and Weights
129d7b241e6Sjeremylt   ierr = CeedLobattoQuadrature(P, nodes, NULL); CeedChk(ierr);
130d7b241e6Sjeremylt   switch (qmode) {
131d7b241e6Sjeremylt   case CEED_GAUSS:
132d7b241e6Sjeremylt     ierr = CeedGaussQuadrature(Q, qref1d, qweight1d); CeedChk(ierr);
133d7b241e6Sjeremylt     break;
134d7b241e6Sjeremylt   case CEED_GAUSS_LOBATTO:
135d7b241e6Sjeremylt     ierr = CeedLobattoQuadrature(Q, qref1d, qweight1d); CeedChk(ierr);
136d7b241e6Sjeremylt     break;
137d7b241e6Sjeremylt   }
138d7b241e6Sjeremylt   // Build B, D matrix
139d7b241e6Sjeremylt   // Fornberg, 1998
140d7b241e6Sjeremylt   for (i = 0; i  < Q; i++) {
141d7b241e6Sjeremylt     c1 = 1.0;
142d7b241e6Sjeremylt     c3 = nodes[0] - qref1d[i];
143d7b241e6Sjeremylt     interp1d[i*P+0] = 1.0;
144d7b241e6Sjeremylt     for (j = 1; j < P; j++) {
145d7b241e6Sjeremylt       c2 = 1.0;
146d7b241e6Sjeremylt       c4 = c3;
147d7b241e6Sjeremylt       c3 = nodes[j] - qref1d[i];
148d7b241e6Sjeremylt       for (k = 0; k < j; k++) {
149d7b241e6Sjeremylt         dx = nodes[j] - nodes[k];
150d7b241e6Sjeremylt         c2 *= dx;
151d7b241e6Sjeremylt         if (k == j - 1) {
152d7b241e6Sjeremylt           grad1d[i*P + j] = c1*(interp1d[i*P + k] - c4*grad1d[i*P + k]) / c2;
153d7b241e6Sjeremylt           interp1d[i*P + j] = - c1*c4*interp1d[i*P + k] / c2;
154d7b241e6Sjeremylt         }
155d7b241e6Sjeremylt         grad1d[i*P + k] = (c3*grad1d[i*P + k] - interp1d[i*P + k]) / dx;
156d7b241e6Sjeremylt         interp1d[i*P + k] = c3*interp1d[i*P + k] / dx;
157d7b241e6Sjeremylt       }
158d7b241e6Sjeremylt       c1 = c2;
159d7b241e6Sjeremylt     }
160d7b241e6Sjeremylt   }
161d7b241e6Sjeremylt   //  // Pass to CeedBasisCreateTensorH1
162d7b241e6Sjeremylt   ierr = CeedBasisCreateTensorH1(ceed, dim, ncomp, P, Q, interp1d, grad1d, qref1d,
163d7b241e6Sjeremylt                                  qweight1d, basis); CeedChk(ierr);
164d7b241e6Sjeremylt   ierr = CeedFree(&interp1d); CeedChk(ierr);
165d7b241e6Sjeremylt   ierr = CeedFree(&grad1d); CeedChk(ierr);
166d7b241e6Sjeremylt   ierr = CeedFree(&nodes); CeedChk(ierr);
167d7b241e6Sjeremylt   ierr = CeedFree(&qref1d); CeedChk(ierr);
168d7b241e6Sjeremylt   ierr = CeedFree(&qweight1d); CeedChk(ierr);
169d7b241e6Sjeremylt   return 0;
170d7b241e6Sjeremylt }
171d7b241e6Sjeremylt 
172b11c1e72Sjeremylt /**
173a8de75f0Sjeremylt   @brief Create a non tensor product basis for H^1 discretizations
174a8de75f0Sjeremylt 
175a8de75f0Sjeremylt   @param ceed       A Ceed object where the CeedBasis will be created
176a8de75f0Sjeremylt   @param topo       Topology of element, e.g. hypercube, simplex, ect
177a8de75f0Sjeremylt   @param ncomp      Number of field components (1 for scalar fields)
178a8de75f0Sjeremylt   @param ndof       Total number of nodes
179a8de75f0Sjeremylt   @param nqpts      Total number of quadrature points
180a8de75f0Sjeremylt   @param interp     Row-major nqpts × ndof matrix expressing the values of nodal
181a8de75f0Sjeremylt                       basis functions at quadrature points
182a8de75f0Sjeremylt   @param grad       Row-major (nqpts x dim) × ndof matrix expressing derivatives
183a8de75f0Sjeremylt                       of nodal basis functions at quadrature points
184a8de75f0Sjeremylt   @param qref       Array of length nqpts holding the locations of quadrature points
185a8de75f0Sjeremylt                       on the reference element [-1, 1]
186a8de75f0Sjeremylt   @param qweight    Array of length nqpts holding the quadrature weights on the
187a8de75f0Sjeremylt                       reference element
188a8de75f0Sjeremylt   @param[out] basis Address of the variable where the newly created
189a8de75f0Sjeremylt                       CeedBasis will be stored.
190a8de75f0Sjeremylt 
191a8de75f0Sjeremylt   @return An error code: 0 - success, otherwise - failure
192a8de75f0Sjeremylt 
193a8de75f0Sjeremylt   @ref Basic
194a8de75f0Sjeremylt **/
195a8de75f0Sjeremylt int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt ncomp,
196a8de75f0Sjeremylt                       CeedInt ndof, CeedInt nqpts,
197a8de75f0Sjeremylt                       const CeedScalar *interp,
198a8de75f0Sjeremylt                       const CeedScalar *grad, const CeedScalar *qref,
199a8de75f0Sjeremylt                       const CeedScalar *qweight, CeedBasis *basis) {
200a8de75f0Sjeremylt   int ierr;
201a8de75f0Sjeremylt   CeedInt P = ndof, Q = nqpts, dim = 0;
202a8de75f0Sjeremylt 
2035fe0d4faSjeremylt   if (!ceed->BasisCreateH1) {
2045fe0d4faSjeremylt     Ceed delegate;
2055fe0d4faSjeremylt     ierr = CeedGetDelegate(ceed, &delegate); CeedChk(ierr);
2065fe0d4faSjeremylt 
2075fe0d4faSjeremylt     if (!delegate)
208a8de75f0Sjeremylt       return CeedError(ceed, 1, "Backend does not support BasisCreateH1");
2095fe0d4faSjeremylt 
2105fe0d4faSjeremylt     ierr = CeedBasisCreateH1(delegate, topo, ncomp, ndof,
2115fe0d4faSjeremylt                             nqpts, interp, grad, qref,
2125fe0d4faSjeremylt                             qweight, basis); CeedChk(ierr);
2135fe0d4faSjeremylt     return 0;
2145fe0d4faSjeremylt   }
2155fe0d4faSjeremylt 
216a8de75f0Sjeremylt   ierr = CeedCalloc(1,basis); CeedChk(ierr);
217a8de75f0Sjeremylt 
218a8de75f0Sjeremylt   ierr = CeedBasisGetTopologyDimension(topo, &dim); CeedChk(ierr);
219a8de75f0Sjeremylt 
220a8de75f0Sjeremylt   (*basis)->ceed = ceed;
221a8de75f0Sjeremylt   ceed->refcount++;
222a8de75f0Sjeremylt   (*basis)->refcount = 1;
223a8de75f0Sjeremylt   (*basis)->tensorbasis = 0;
224a8de75f0Sjeremylt   (*basis)->dim = dim;
225a8de75f0Sjeremylt   (*basis)->ncomp = ncomp;
226a8de75f0Sjeremylt   (*basis)->P = P;
227a8de75f0Sjeremylt   (*basis)->Q = Q;
228a8de75f0Sjeremylt   ierr = CeedMalloc(Q*dim,&(*basis)->qref1d); CeedChk(ierr);
229a8de75f0Sjeremylt   ierr = CeedMalloc(Q,&(*basis)->qweight1d); CeedChk(ierr);
230a8de75f0Sjeremylt   memcpy((*basis)->qref1d, qref, Q*dim*sizeof(qref[0]));
231a8de75f0Sjeremylt   memcpy((*basis)->qweight1d, qweight, Q*sizeof(qweight[0]));
232a8de75f0Sjeremylt   ierr = CeedMalloc(Q*P,&(*basis)->interp1d); CeedChk(ierr);
233a8de75f0Sjeremylt   ierr = CeedMalloc(dim*Q*P,&(*basis)->grad1d); CeedChk(ierr);
234a8de75f0Sjeremylt   memcpy((*basis)->interp1d, interp, Q*P*sizeof(interp[0]));
235a8de75f0Sjeremylt   memcpy((*basis)->grad1d, grad, dim*Q*P*sizeof(grad[0]));
236667bc5fcSjeremylt   ierr = ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, qref,
237a8de75f0Sjeremylt                              qweight, *basis); CeedChk(ierr);
238a8de75f0Sjeremylt   return 0;
239a8de75f0Sjeremylt }
240a8de75f0Sjeremylt 
241a8de75f0Sjeremylt /**
242b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre quadrature
243b11c1e72Sjeremylt 
244b11c1e72Sjeremylt   @param Q              Number of quadrature points (integrates polynomials of
245b11c1e72Sjeremylt                           degree 2*Q-1 exactly)
246b11c1e72Sjeremylt   @param[out] qref1d    Array of length Q to hold the abscissa on [-1, 1]
247b11c1e72Sjeremylt   @param[out] qweight1d Array of length Q to hold the weights
248b11c1e72Sjeremylt 
249b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
250dfdf5a53Sjeremylt 
251dfdf5a53Sjeremylt   @ref Utility
252b11c1e72Sjeremylt **/
253d7b241e6Sjeremylt int CeedGaussQuadrature(CeedInt Q, CeedScalar *qref1d, CeedScalar *qweight1d) {
254d7b241e6Sjeremylt   // Allocate
255d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0*atan(1.0);
256d7b241e6Sjeremylt   // Build qref1d, qweight1d
257d7b241e6Sjeremylt   for (int i = 0; i <= Q/2; i++) {
258d7b241e6Sjeremylt     // Guess
259d7b241e6Sjeremylt     xi = cos(PI*(CeedScalar)(2*i+1)/((CeedScalar)(2*Q)));
260d7b241e6Sjeremylt     // Pn(xi)
261d7b241e6Sjeremylt     P0 = 1.0;
262d7b241e6Sjeremylt     P1 = xi;
263d7b241e6Sjeremylt     P2 = 0.0;
264d7b241e6Sjeremylt     for (int j = 2; j <= Q; j++) {
265d7b241e6Sjeremylt       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
266d7b241e6Sjeremylt       P0 = P1;
267d7b241e6Sjeremylt       P1 = P2;
268d7b241e6Sjeremylt     }
269d7b241e6Sjeremylt     // First Newton Step
270d7b241e6Sjeremylt     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
271d7b241e6Sjeremylt     xi = xi-P2/dP2;
272d7b241e6Sjeremylt     // Newton to convergence
273d7b241e6Sjeremylt     for (int k=0; k<100 && fabs(P2)>1e-15; k++) {
274d7b241e6Sjeremylt       P0 = 1.0;
275d7b241e6Sjeremylt       P1 = xi;
276d7b241e6Sjeremylt       for (int j = 2; j <= Q; j++) {
277d7b241e6Sjeremylt         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
278d7b241e6Sjeremylt         P0 = P1;
279d7b241e6Sjeremylt         P1 = P2;
280d7b241e6Sjeremylt       }
281d7b241e6Sjeremylt       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
282d7b241e6Sjeremylt       xi = xi-P2/dP2;
283d7b241e6Sjeremylt     }
284d7b241e6Sjeremylt     // Save xi, wi
285d7b241e6Sjeremylt     wi = 2.0/((1.0-xi*xi)*dP2*dP2);
286d7b241e6Sjeremylt     qweight1d[i] = wi;
287d7b241e6Sjeremylt     qweight1d[Q-1-i] = wi;
288d7b241e6Sjeremylt     qref1d[i] = -xi;
289d7b241e6Sjeremylt     qref1d[Q-1-i]= xi;
290d7b241e6Sjeremylt   }
291d7b241e6Sjeremylt   return 0;
292d7b241e6Sjeremylt }
293d7b241e6Sjeremylt 
294b11c1e72Sjeremylt /**
295b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre-Lobatto quadrature
296b11c1e72Sjeremylt 
297b11c1e72Sjeremylt   @param Q              Number of quadrature points (integrates polynomials of
298b11c1e72Sjeremylt                           degree 2*Q-3 exactly)
299b11c1e72Sjeremylt   @param[out] qref1d    Array of length Q to hold the abscissa on [-1, 1]
300b11c1e72Sjeremylt   @param[out] qweight1d Array of length Q to hold the weights
301b11c1e72Sjeremylt 
302b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
303dfdf5a53Sjeremylt 
304dfdf5a53Sjeremylt   @ref Utility
305b11c1e72Sjeremylt **/
306d7b241e6Sjeremylt int CeedLobattoQuadrature(CeedInt Q, CeedScalar *qref1d,
307d7b241e6Sjeremylt                           CeedScalar *qweight1d) {
308d7b241e6Sjeremylt   // Allocate
309d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0*atan(1.0);
310d7b241e6Sjeremylt   // Build qref1d, qweight1d
311d7b241e6Sjeremylt   // Set endpoints
312d7b241e6Sjeremylt   wi = 2.0/((CeedScalar)(Q*(Q-1)));
313d7b241e6Sjeremylt   if (qweight1d) {
314d7b241e6Sjeremylt     qweight1d[0] = wi;
315d7b241e6Sjeremylt     qweight1d[Q-1] = wi;
316d7b241e6Sjeremylt   }
317d7b241e6Sjeremylt   qref1d[0] = -1.0;
318d7b241e6Sjeremylt   qref1d[Q-1] = 1.0;
319d7b241e6Sjeremylt   // Interior
320d7b241e6Sjeremylt   for (int i = 1; i <= (Q-1)/2; i++) {
321d7b241e6Sjeremylt     // Guess
322d7b241e6Sjeremylt     xi = cos(PI*(CeedScalar)(i)/(CeedScalar)(Q-1));
323d7b241e6Sjeremylt     // Pn(xi)
324d7b241e6Sjeremylt     P0 = 1.0;
325d7b241e6Sjeremylt     P1 = xi;
326d7b241e6Sjeremylt     P2 = 0.0;
327d7b241e6Sjeremylt     for (int j = 2; j < Q; j++) {
328d7b241e6Sjeremylt       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
329d7b241e6Sjeremylt       P0 = P1;
330d7b241e6Sjeremylt       P1 = P2;
331d7b241e6Sjeremylt     }
332d7b241e6Sjeremylt     // First Newton step
333d7b241e6Sjeremylt     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
334d7b241e6Sjeremylt     d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
335d7b241e6Sjeremylt     xi = xi-dP2/d2P2;
336d7b241e6Sjeremylt     // Newton to convergence
337d7b241e6Sjeremylt     for (int k=0; k<100 && fabs(dP2)>1e-15; k++) {
338d7b241e6Sjeremylt       P0 = 1.0;
339d7b241e6Sjeremylt       P1 = xi;
340d7b241e6Sjeremylt       for (int j = 2; j < Q; j++) {
341d7b241e6Sjeremylt         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
342d7b241e6Sjeremylt         P0 = P1;
343d7b241e6Sjeremylt         P1 = P2;
344d7b241e6Sjeremylt       }
345d7b241e6Sjeremylt       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
346d7b241e6Sjeremylt       d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
347d7b241e6Sjeremylt       xi = xi-dP2/d2P2;
348d7b241e6Sjeremylt     }
349d7b241e6Sjeremylt     // Save xi, wi
350d7b241e6Sjeremylt     wi = 2.0/(((CeedScalar)(Q*(Q-1)))*P2*P2);
351d7b241e6Sjeremylt     if (qweight1d) {
352d7b241e6Sjeremylt       qweight1d[i] = wi;
353d7b241e6Sjeremylt       qweight1d[Q-1-i] = wi;
354d7b241e6Sjeremylt     }
355d7b241e6Sjeremylt     qref1d[i] = -xi;
356d7b241e6Sjeremylt     qref1d[Q-1-i]= xi;
357d7b241e6Sjeremylt   }
358d7b241e6Sjeremylt   return 0;
359d7b241e6Sjeremylt }
360d7b241e6Sjeremylt 
361dfdf5a53Sjeremylt /**
362dfdf5a53Sjeremylt   @brief View an array stored in a CeedBasis
363dfdf5a53Sjeremylt 
364dfdf5a53Sjeremylt   @param name      Name of array
365dfdf5a53Sjeremylt   @param fpformat  Printing format
366dfdf5a53Sjeremylt   @param m         Number of rows in array
367dfdf5a53Sjeremylt   @param n         Number of columns in array
368dfdf5a53Sjeremylt   @param a         Array to be viewed
369dfdf5a53Sjeremylt   @param stream    Stream to view to, e.g., stdout
370dfdf5a53Sjeremylt 
371dfdf5a53Sjeremylt   @return An error code: 0 - success, otherwise - failure
372dfdf5a53Sjeremylt 
373dfdf5a53Sjeremylt   @ref Utility
374dfdf5a53Sjeremylt **/
375d7b241e6Sjeremylt static int CeedScalarView(const char *name, const char *fpformat, CeedInt m,
376d7b241e6Sjeremylt                           CeedInt n, const CeedScalar *a, FILE *stream) {
377d7b241e6Sjeremylt   for (int i=0; i<m; i++) {
378d7b241e6Sjeremylt     if (m > 1) fprintf(stream, "%12s[%d]:", name, i);
379d7b241e6Sjeremylt     else fprintf(stream, "%12s:", name);
380d7b241e6Sjeremylt     for (int j=0; j<n; j++) {
381d7b241e6Sjeremylt       fprintf(stream, fpformat, fabs(a[i*n+j]) > 1E-14 ? a[i*n+j] : 0);
382d7b241e6Sjeremylt     }
383d7b241e6Sjeremylt     fputs("\n", stream);
384d7b241e6Sjeremylt   }
385d7b241e6Sjeremylt   return 0;
386d7b241e6Sjeremylt }
387d7b241e6Sjeremylt 
388b11c1e72Sjeremylt /**
389b11c1e72Sjeremylt   @brief View a CeedBasis
390b11c1e72Sjeremylt 
391b11c1e72Sjeremylt   @param basis  CeedBasis to view
392b11c1e72Sjeremylt   @param stream Stream to view to, e.g., stdout
393b11c1e72Sjeremylt 
394b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
395dfdf5a53Sjeremylt 
396dfdf5a53Sjeremylt   @ref Utility
397b11c1e72Sjeremylt **/
398d7b241e6Sjeremylt int CeedBasisView(CeedBasis basis, FILE *stream) {
399d7b241e6Sjeremylt   int ierr;
400d7b241e6Sjeremylt 
401a8de75f0Sjeremylt   if (basis->tensorbasis) {
402d7b241e6Sjeremylt     fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P1d,
403d7b241e6Sjeremylt             basis->Q1d);
404d7b241e6Sjeremylt     ierr = CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q1d, basis->qref1d,
405d7b241e6Sjeremylt                           stream); CeedChk(ierr);
406d7b241e6Sjeremylt     ierr = CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q1d, basis->qweight1d,
407d7b241e6Sjeremylt                           stream); CeedChk(ierr);
408d7b241e6Sjeremylt     ierr = CeedScalarView("interp1d", "\t% 12.8f", basis->Q1d, basis->P1d,
409d7b241e6Sjeremylt                           basis->interp1d, stream); CeedChk(ierr);
410d7b241e6Sjeremylt     ierr = CeedScalarView("grad1d", "\t% 12.8f", basis->Q1d, basis->P1d,
411d7b241e6Sjeremylt                           basis->grad1d, stream); CeedChk(ierr);
412a8de75f0Sjeremylt   } else {
413a8de75f0Sjeremylt     fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P,
414a8de75f0Sjeremylt             basis->Q);
415a8de75f0Sjeremylt     ierr = CeedScalarView("qref", "\t% 12.8f", 1, basis->Q*basis->dim,
416a8de75f0Sjeremylt                           basis->qref1d,
417a8de75f0Sjeremylt                           stream); CeedChk(ierr);
418a8de75f0Sjeremylt     ierr = CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->qweight1d,
419a8de75f0Sjeremylt                           stream); CeedChk(ierr);
420a8de75f0Sjeremylt     ierr = CeedScalarView("interp", "\t% 12.8f", basis->Q, basis->P,
421a8de75f0Sjeremylt                           basis->interp1d, stream); CeedChk(ierr);
422a8de75f0Sjeremylt     ierr = CeedScalarView("grad", "\t% 12.8f", basis->dim*basis->Q, basis->P,
423a8de75f0Sjeremylt                           basis->grad1d, stream); CeedChk(ierr);
424a8de75f0Sjeremylt   }
425d7b241e6Sjeremylt   return 0;
426d7b241e6Sjeremylt }
427d7b241e6Sjeremylt 
428dfdf5a53Sjeremylt /**
429dfdf5a53Sjeremylt   @brief Compute Householder Reflection
430dfdf5a53Sjeremylt 
431dfdf5a53Sjeremylt     Computes A = (I - b v v^T) A
432dfdf5a53Sjeremylt     where A is an mxn matrix indexed as A[i*row + j*col]
433dfdf5a53Sjeremylt 
434dfdf5a53Sjeremylt   @param[out] A  Matrix to apply Householder reflection to, in place
435dfdf5a53Sjeremylt   @param v       Householder vector
436dfdf5a53Sjeremylt   @param b       Scaling factor
437dfdf5a53Sjeremylt   @param m       Number of rows in A
438dfdf5a53Sjeremylt   @param n       Number of columns in A
439dfdf5a53Sjeremylt   @param row     Col stride
440dfdf5a53Sjeremylt   @param col     Row stride
441dfdf5a53Sjeremylt 
442dfdf5a53Sjeremylt   @return An error code: 0 - success, otherwise - failure
443dfdf5a53Sjeremylt 
444dfdf5a53Sjeremylt   @ref Developer
445dfdf5a53Sjeremylt **/
446d7b241e6Sjeremylt static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v,
447d7b241e6Sjeremylt                                   CeedScalar b, CeedInt m, CeedInt n,
448d7b241e6Sjeremylt                                   CeedInt row, CeedInt col) {
449d7b241e6Sjeremylt   for (CeedInt j=0; j<n; j++) {
450d7b241e6Sjeremylt     CeedScalar w = A[0*row + j*col];
451d7b241e6Sjeremylt     for (CeedInt i=1; i<m; i++) w += v[i] * A[i*row + j*col];
452d7b241e6Sjeremylt     A[0*row + j*col] -= b * w;
453d7b241e6Sjeremylt     for (CeedInt i=1; i<m; i++) A[i*row + j*col] -= b * w * v[i];
454d7b241e6Sjeremylt   }
455d7b241e6Sjeremylt   return 0;
456d7b241e6Sjeremylt }
457d7b241e6Sjeremylt 
458dfdf5a53Sjeremylt /**
459dfdf5a53Sjeremylt   @brief Apply Householder Q matrix
460dfdf5a53Sjeremylt 
461dfdf5a53Sjeremylt     Compute A = Q A where Q is mxk and A is mxn. k<m
462dfdf5a53Sjeremylt 
463dfdf5a53Sjeremylt   @param[out] A  Matrix to apply Householder Q to, in place
464dfdf5a53Sjeremylt   @param Q       Householder Q matrix
465dfdf5a53Sjeremylt   @param tau     Householder scaling factors
466dfdf5a53Sjeremylt   @param tmode   Transpose mode for application
467dfdf5a53Sjeremylt   @param m       Number of rows in A
468dfdf5a53Sjeremylt   @param n       Number of columns in A
469dfdf5a53Sjeremylt   @param k       Index of row targeted
470dfdf5a53Sjeremylt   @param row     Col stride
471dfdf5a53Sjeremylt   @param col     Row stride
472dfdf5a53Sjeremylt 
473dfdf5a53Sjeremylt   @return An error code: 0 - success, otherwise - failure
474dfdf5a53Sjeremylt 
475dfdf5a53Sjeremylt   @ref Developer
476dfdf5a53Sjeremylt **/
477d7b241e6Sjeremylt static int CeedHouseholderApplyQ(CeedScalar *A, const CeedScalar *Q,
478d7b241e6Sjeremylt                                  const CeedScalar *tau, CeedTransposeMode tmode,
479d7b241e6Sjeremylt                                  CeedInt m, CeedInt n, CeedInt k,
480d7b241e6Sjeremylt                                  CeedInt row, CeedInt col) {
481d7b241e6Sjeremylt   CeedScalar v[m];
482d7b241e6Sjeremylt   for (CeedInt ii=0; ii<k; ii++) {
483d7b241e6Sjeremylt     CeedInt i = tmode == CEED_TRANSPOSE ? ii : k-1-ii;
484d7b241e6Sjeremylt     for (CeedInt j=i+1; j<m; j++) {
485d7b241e6Sjeremylt       v[j] = Q[j*k+i];
486d7b241e6Sjeremylt     }
487d7b241e6Sjeremylt     // Apply Householder reflector (I - tau v v^T) colograd1d^T
488d7b241e6Sjeremylt     CeedHouseholderReflect(&A[i*row], &v[i], tau[i], m-i, n, row, col);
489d7b241e6Sjeremylt   }
490d7b241e6Sjeremylt   return 0;
491d7b241e6Sjeremylt }
492d7b241e6Sjeremylt 
493b11c1e72Sjeremylt /**
494b11c1e72Sjeremylt   @brief Return QR Factorization of matrix
495b11c1e72Sjeremylt 
496b11c1e72Sjeremylt   @param[out] mat  Row-major matrix to be factorized in place
497b11c1e72Sjeremylt   @param[out] tau  Vector of length m of scaling fators
498b11c1e72Sjeremylt   @param m         Number of rows
499b11c1e72Sjeremylt   @param n         Number of columns
500b11c1e72Sjeremylt 
501b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
502dfdf5a53Sjeremylt 
503dfdf5a53Sjeremylt   @ref Utility
504b11c1e72Sjeremylt **/
505d7b241e6Sjeremylt int CeedQRFactorization(CeedScalar *mat, CeedScalar *tau,
506d7b241e6Sjeremylt                         CeedInt m, CeedInt n) {
507d7b241e6Sjeremylt   CeedInt i, j;
508d7b241e6Sjeremylt   CeedScalar v[m];
509d7b241e6Sjeremylt 
510d7b241e6Sjeremylt   for (i=0; i<n; i++) {
511d7b241e6Sjeremylt     // Calculate Householder vector, magnitude
512d7b241e6Sjeremylt     CeedScalar sigma = 0.0;
513d7b241e6Sjeremylt     v[i] = mat[i+n*i];
514d7b241e6Sjeremylt     for (j=i+1; j<m; j++) {
515d7b241e6Sjeremylt       v[j] = mat[i+n*j];
516d7b241e6Sjeremylt       sigma += v[j] * v[j];
517d7b241e6Sjeremylt     }
518d7b241e6Sjeremylt     CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m]
519d7b241e6Sjeremylt     CeedScalar Rii = -copysign(norm, v[i]);
520d7b241e6Sjeremylt     v[i] -= Rii;
521d7b241e6Sjeremylt     // norm of v[i:m] after modification above and scaling below
522d7b241e6Sjeremylt     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
523d7b241e6Sjeremylt     //   tau = 2 / (norm*norm)
524d7b241e6Sjeremylt     tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma);
525d7b241e6Sjeremylt     for (j=i+1; j<m; j++) v[j] /= v[i];
526d7b241e6Sjeremylt 
527d7b241e6Sjeremylt     // Apply Householder reflector to lower right panel
528d7b241e6Sjeremylt     CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1);
529d7b241e6Sjeremylt     // Save v
530d7b241e6Sjeremylt     mat[i+n*i] = Rii;
531d7b241e6Sjeremylt     for (j=i+1; j<m; j++) {
532d7b241e6Sjeremylt       mat[i+n*j] = v[j];
533d7b241e6Sjeremylt     }
534d7b241e6Sjeremylt   }
535d7b241e6Sjeremylt 
536d7b241e6Sjeremylt   return 0;
537d7b241e6Sjeremylt }
538d7b241e6Sjeremylt 
539b11c1e72Sjeremylt /**
540783c99b3SValeria Barra   @brief Return collocated grad matrix
541b11c1e72Sjeremylt 
542b11c1e72Sjeremylt   @param basis           CeedBasis
543b11c1e72Sjeremylt   @param[out] colograd1d Row-major Q1d × Q1d matrix expressing derivatives of
544b11c1e72Sjeremylt                            basis functions at quadrature points
545b11c1e72Sjeremylt 
546b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
547dfdf5a53Sjeremylt 
548dfdf5a53Sjeremylt   @ref Advanced
549b11c1e72Sjeremylt **/
550783c99b3SValeria Barra int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *colograd1d) {
551d7b241e6Sjeremylt   int i, j, k;
552d7b241e6Sjeremylt   CeedInt ierr, P1d=(basis)->P1d, Q1d=(basis)->Q1d;
553d7b241e6Sjeremylt   CeedScalar *interp1d, *grad1d, tau[Q1d];
554d7b241e6Sjeremylt 
555d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d*P1d, &interp1d); CeedChk(ierr);
556d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d*P1d, &grad1d); CeedChk(ierr);
557d7b241e6Sjeremylt   memcpy(interp1d, (basis)->interp1d, Q1d*P1d*sizeof(basis)->interp1d[0]);
558d7b241e6Sjeremylt   memcpy(grad1d, (basis)->grad1d, Q1d*P1d*sizeof(basis)->interp1d[0]);
559d7b241e6Sjeremylt 
560d7b241e6Sjeremylt   // QR Factorization, interp1d = Q R
561d7b241e6Sjeremylt   ierr = CeedQRFactorization(interp1d, tau, Q1d, P1d); CeedChk(ierr);
562d7b241e6Sjeremylt 
563d7b241e6Sjeremylt   // Apply Rinv, colograd1d = grad1d Rinv
564d7b241e6Sjeremylt   for (i=0; i<Q1d; i++) { // Row i
565d7b241e6Sjeremylt     colograd1d[Q1d*i] = grad1d[P1d*i]/interp1d[0];
566d7b241e6Sjeremylt     for (j=1; j<P1d; j++) { // Column j
567d7b241e6Sjeremylt       colograd1d[j+Q1d*i] = grad1d[j+P1d*i];
568d7b241e6Sjeremylt       for (k=0; k<j; k++) {
569d7b241e6Sjeremylt         colograd1d[j+Q1d*i] -= interp1d[j+P1d*k]*colograd1d[k+Q1d*i];
570d7b241e6Sjeremylt       }
571d7b241e6Sjeremylt       colograd1d[j+Q1d*i] /= interp1d[j+P1d*j];
572d7b241e6Sjeremylt     }
573d7b241e6Sjeremylt     for (j=P1d; j<Q1d; j++) {
574d7b241e6Sjeremylt       colograd1d[j+Q1d*i] = 0;
575d7b241e6Sjeremylt     }
576d7b241e6Sjeremylt   }
577d7b241e6Sjeremylt 
578d7b241e6Sjeremylt   // Apply Qtranspose, colograd = colograd Qtranspose
579d7b241e6Sjeremylt   CeedHouseholderApplyQ(colograd1d, interp1d, tau, CEED_NOTRANSPOSE,
580d7b241e6Sjeremylt                         Q1d, Q1d, P1d, 1, Q1d);
581d7b241e6Sjeremylt 
582d7b241e6Sjeremylt   ierr = CeedFree(&interp1d); CeedChk(ierr);
583d7b241e6Sjeremylt   ierr = CeedFree(&grad1d); CeedChk(ierr);
584d7b241e6Sjeremylt 
585d7b241e6Sjeremylt   return 0;
586d7b241e6Sjeremylt }
587d7b241e6Sjeremylt 
588b11c1e72Sjeremylt /**
589b11c1e72Sjeremylt   @brief Apply basis evaluation from nodes to quadrature points or vice-versa
590b11c1e72Sjeremylt 
591b11c1e72Sjeremylt   @param basis  CeedBasis to evaluate
592b11c1e72Sjeremylt   @param nelem  The number of elements to apply the basis evaluation to;
593b11c1e72Sjeremylt                   the backend will specify the ordering in
594b11c1e72Sjeremylt                   ElemRestrictionCreateBlocked
595b11c1e72Sjeremylt   @param tmode  \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature
596b11c1e72Sjeremylt                   points, \ref CEED_TRANSPOSE to apply the transpose, mapping
597b11c1e72Sjeremylt                   from quadrature points to nodes
598b11c1e72Sjeremylt   @param emode  \ref CEED_EVAL_INTERP to obtain interpolated values,
599b11c1e72Sjeremylt                   \ref CEED_EVAL_GRAD to obtain gradients.
600b11c1e72Sjeremylt   @param[in] u  Input array
601b11c1e72Sjeremylt   @param[out] v Output array
602b11c1e72Sjeremylt 
603b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
604dfdf5a53Sjeremylt 
605dfdf5a53Sjeremylt   @ref Advanced
606b11c1e72Sjeremylt **/
607d7b241e6Sjeremylt int CeedBasisApply(CeedBasis basis, CeedInt nelem, CeedTransposeMode tmode,
608d7b241e6Sjeremylt                    CeedEvalMode emode, const CeedScalar *u, CeedScalar *v) {
609d7b241e6Sjeremylt   int ierr;
610d7b241e6Sjeremylt   if (!basis->Apply) return CeedError(basis->ceed, 1,
611d7b241e6Sjeremylt                                         "Backend does not support BasisApply");
612d7b241e6Sjeremylt   ierr = basis->Apply(basis, nelem, tmode, emode, u, v); CeedChk(ierr);
613d7b241e6Sjeremylt   return 0;
614d7b241e6Sjeremylt }
615d7b241e6Sjeremylt 
616b11c1e72Sjeremylt /**
6174ce2993fSjeremylt   @brief Get Ceed associated with a CeedBasis
618b11c1e72Sjeremylt 
619b11c1e72Sjeremylt   @param basis      CeedBasis
6204ce2993fSjeremylt   @param[out] ceed  Variable to store Ceed
6214ce2993fSjeremylt 
6224ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
6234ce2993fSjeremylt 
62423617272Sjeremylt   @ref Advanced
6254ce2993fSjeremylt **/
6264ce2993fSjeremylt int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) {
6274ce2993fSjeremylt   *ceed = basis->ceed;
6284ce2993fSjeremylt 
6294ce2993fSjeremylt   return 0;
6304ce2993fSjeremylt };
6314ce2993fSjeremylt 
6324ce2993fSjeremylt /**
6334ce2993fSjeremylt   @brief Get dimension for given CeedBasis
6344ce2993fSjeremylt 
6354ce2993fSjeremylt   @param basis     CeedBasis
6364ce2993fSjeremylt   @param[out] dim  Variable to store dimension of basis
6374ce2993fSjeremylt 
6384ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
6394ce2993fSjeremylt 
64023617272Sjeremylt   @ref Advanced
6414ce2993fSjeremylt **/
6424ce2993fSjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) {
6434ce2993fSjeremylt   *dim = basis->dim;
6444ce2993fSjeremylt 
6454ce2993fSjeremylt   return 0;
6464ce2993fSjeremylt };
6474ce2993fSjeremylt 
6484ce2993fSjeremylt /**
6494ce2993fSjeremylt   @brief Get tensor status for given CeedBasis
6504ce2993fSjeremylt 
6514ce2993fSjeremylt   @param basis        CeedBasis
6524ce2993fSjeremylt   @param[out] tensor  Variable to store tensor status
6534ce2993fSjeremylt 
6544ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
6554ce2993fSjeremylt 
65623617272Sjeremylt   @ref Advanced
6574ce2993fSjeremylt **/
6584ce2993fSjeremylt int CeedBasisGetTensorStatus(CeedBasis basis, bool *tensor) {
6594ce2993fSjeremylt   *tensor = basis->tensorbasis;
6604ce2993fSjeremylt 
6614ce2993fSjeremylt   return 0;
6624ce2993fSjeremylt };
6634ce2993fSjeremylt 
6644ce2993fSjeremylt /**
6654ce2993fSjeremylt   @brief Get number of components for given CeedBasis
6664ce2993fSjeremylt 
6674ce2993fSjeremylt   @param basis     CeedBasis
6684ce2993fSjeremylt   @param[out] dim  Variable to store number of components of basis
6694ce2993fSjeremylt 
6704ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
6714ce2993fSjeremylt 
67223617272Sjeremylt   @ref Advanced
6734ce2993fSjeremylt **/
6744ce2993fSjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *numcomp) {
6754ce2993fSjeremylt   *numcomp = basis->ncomp;
6764ce2993fSjeremylt 
6774ce2993fSjeremylt   return 0;
6784ce2993fSjeremylt };
6794ce2993fSjeremylt 
6804ce2993fSjeremylt /**
6814ce2993fSjeremylt   @brief Get total number of nodes (in 1 dimension) of a CeedBasis
6824ce2993fSjeremylt 
6834ce2993fSjeremylt   @param basis     CeedBasis
6844ce2993fSjeremylt   @param[out] P1d  Variable to store number of nodes
6854ce2993fSjeremylt 
6864ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
6874ce2993fSjeremylt 
68823617272Sjeremylt   @ref Advanced
6894ce2993fSjeremylt **/
6904ce2993fSjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P1d) {
6914ce2993fSjeremylt   if (!basis->tensorbasis) return CeedError(basis->ceed, 1,
6924ce2993fSjeremylt                                       "Cannot supply P1d for non-tensor basis");
6934ce2993fSjeremylt   *P1d = basis->P1d;
6944ce2993fSjeremylt   return 0;
6954ce2993fSjeremylt }
6964ce2993fSjeremylt 
6974ce2993fSjeremylt /**
6984ce2993fSjeremylt   @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis
6994ce2993fSjeremylt 
7004ce2993fSjeremylt   @param basis     CeedBasis
7014ce2993fSjeremylt   @param[out] Q1d  Variable to store number of quadrature points
7024ce2993fSjeremylt 
7034ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
7044ce2993fSjeremylt 
70523617272Sjeremylt   @ref Advanced
7064ce2993fSjeremylt **/
7074ce2993fSjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q1d) {
7084ce2993fSjeremylt   if (!basis->tensorbasis) return CeedError(basis->ceed, 1,
7094ce2993fSjeremylt                                       "Cannot supply Q1d for non-tensor basis");
7104ce2993fSjeremylt   *Q1d = basis->Q1d;
7114ce2993fSjeremylt   return 0;
7124ce2993fSjeremylt }
7134ce2993fSjeremylt 
7144ce2993fSjeremylt /**
7154ce2993fSjeremylt   @brief Get total number of nodes (in dim dimensions) of a CeedBasis
7164ce2993fSjeremylt 
7174ce2993fSjeremylt   @param basis   CeedBasis
7184ce2993fSjeremylt   @param[out] P  Variable to store number of nodes
719b11c1e72Sjeremylt 
720b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
721dfdf5a53Sjeremylt 
722dfdf5a53Sjeremylt   @ref Utility
723b11c1e72Sjeremylt **/
724d7b241e6Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) {
725a8de75f0Sjeremylt   *P = basis->P;
726d7b241e6Sjeremylt   return 0;
727d7b241e6Sjeremylt }
728d7b241e6Sjeremylt 
729b11c1e72Sjeremylt /**
7304ce2993fSjeremylt   @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis
731b11c1e72Sjeremylt 
732b11c1e72Sjeremylt   @param basis   CeedBasis
7334ce2993fSjeremylt   @param[out] Q  Variable to store number of quadrature points
734b11c1e72Sjeremylt 
735b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
736dfdf5a53Sjeremylt 
737dfdf5a53Sjeremylt   @ref Utility
738b11c1e72Sjeremylt **/
739d7b241e6Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) {
740a8de75f0Sjeremylt   *Q = basis->Q;
741d7b241e6Sjeremylt   return 0;
742d7b241e6Sjeremylt }
743d7b241e6Sjeremylt 
744b11c1e72Sjeremylt /**
7454ce2993fSjeremylt   @brief Get refrence coordinates of quadrature points (in dim dimensions)
7464ce2993fSjeremylt          of a CeedBasis
7474ce2993fSjeremylt 
7484ce2993fSjeremylt   @param basis      CeedBasis
7494ce2993fSjeremylt   @param[out] qref  Variable to store refrence coordinates of quadrature points
7504ce2993fSjeremylt 
7514ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
7524ce2993fSjeremylt 
75323617272Sjeremylt   @ref Advanced
7544ce2993fSjeremylt **/
7554ce2993fSjeremylt int CeedBasisGetQRef(CeedBasis basis, CeedScalar* *qref) {
7564ce2993fSjeremylt   *qref = basis->qref1d;
7574ce2993fSjeremylt   return 0;
7584ce2993fSjeremylt }
7594ce2993fSjeremylt 
7604ce2993fSjeremylt /**
7614ce2993fSjeremylt   @brief Get quadrature weights of quadrature points (in dim dimensions)
7624ce2993fSjeremylt          of a CeedBasis
7634ce2993fSjeremylt 
7644ce2993fSjeremylt   @param basis         CeedBasis
7654ce2993fSjeremylt   @param[out] qweight  Variable to store quadrature weights
7664ce2993fSjeremylt 
7674ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
7684ce2993fSjeremylt 
76923617272Sjeremylt   @ref Advanced
7704ce2993fSjeremylt **/
7714ce2993fSjeremylt int CeedBasisGetQWeights(CeedBasis basis, CeedScalar* *qweight) {
7724ce2993fSjeremylt   *qweight = basis->qweight1d;
7734ce2993fSjeremylt   return 0;
7744ce2993fSjeremylt }
7754ce2993fSjeremylt 
7764ce2993fSjeremylt /**
7774ce2993fSjeremylt   @brief Get interpolation matrix of a CeedBasis
7784ce2993fSjeremylt 
7794ce2993fSjeremylt   @param basis      CeedBasis
7804ce2993fSjeremylt   @param[out] qref  Variable to store interpolation matrix
7814ce2993fSjeremylt 
7824ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
7834ce2993fSjeremylt 
78423617272Sjeremylt   @ref Advanced
7854ce2993fSjeremylt **/
7864ce2993fSjeremylt int CeedBasisGetInterp(CeedBasis basis, CeedScalar* *interp) {
7874ce2993fSjeremylt   *interp = basis->interp1d;
7884ce2993fSjeremylt   return 0;
7894ce2993fSjeremylt }
7904ce2993fSjeremylt 
7914ce2993fSjeremylt /**
7924ce2993fSjeremylt   @brief Get gradient matrix of a CeedBasis
7934ce2993fSjeremylt 
7944ce2993fSjeremylt   @param basis      CeedBasis
7954ce2993fSjeremylt   @param[out] qref  Variable to store gradient matrix
7964ce2993fSjeremylt 
7974ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
7984ce2993fSjeremylt 
79923617272Sjeremylt   @ref Advanced
8004ce2993fSjeremylt **/
8014ce2993fSjeremylt int CeedBasisGetGrad(CeedBasis basis, CeedScalar* *grad) {
8024ce2993fSjeremylt   *grad = basis->grad1d;
8034ce2993fSjeremylt   return 0;
8044ce2993fSjeremylt }
8054ce2993fSjeremylt 
8064ce2993fSjeremylt /**
8074ce2993fSjeremylt   @brief Get backend data of a CeedBasis
8084ce2993fSjeremylt 
8094ce2993fSjeremylt   @param basis      CeedBasis
8104ce2993fSjeremylt   @param[out] data  Variable to store data
8114ce2993fSjeremylt 
8124ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
8134ce2993fSjeremylt 
81423617272Sjeremylt   @ref Advanced
8154ce2993fSjeremylt **/
8164ce2993fSjeremylt int CeedBasisGetData(CeedBasis basis, void* *data) {
8174ce2993fSjeremylt   *data = basis->data;
8184ce2993fSjeremylt   return 0;
8194ce2993fSjeremylt }
8204ce2993fSjeremylt 
8214ce2993fSjeremylt /**
822*fe2413ffSjeremylt   @brief Set backend data of a CeedBasis
823*fe2413ffSjeremylt 
824*fe2413ffSjeremylt   @param[out] basis CeedBasis
825*fe2413ffSjeremylt   @param data       Data to set
826*fe2413ffSjeremylt 
827*fe2413ffSjeremylt   @return An error code: 0 - success, otherwise - failure
828*fe2413ffSjeremylt 
829*fe2413ffSjeremylt   @ref Advanced
830*fe2413ffSjeremylt **/
831*fe2413ffSjeremylt int CeedBasisSetData(CeedBasis basis, void* *data) {
832*fe2413ffSjeremylt   basis->data = *data;
833*fe2413ffSjeremylt   return 0;
834*fe2413ffSjeremylt }
835*fe2413ffSjeremylt 
836*fe2413ffSjeremylt /**
837a8de75f0Sjeremylt   @brief Get dimension for given CeedElemTopology
838a8de75f0Sjeremylt 
839a8de75f0Sjeremylt   @param topo      CeedElemTopology
8404ce2993fSjeremylt   @param[out] dim  Variable to store dimension of topology
841a8de75f0Sjeremylt 
842a8de75f0Sjeremylt   @return An error code: 0 - success, otherwise - failure
843a8de75f0Sjeremylt 
84423617272Sjeremylt   @ref Advanced
845a8de75f0Sjeremylt **/
846a8de75f0Sjeremylt int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) {
847a8de75f0Sjeremylt   *dim = (CeedInt) topo >> 16;
848a8de75f0Sjeremylt 
849a8de75f0Sjeremylt   return 0;
850a8de75f0Sjeremylt };
851a8de75f0Sjeremylt 
852a8de75f0Sjeremylt /**
853b11c1e72Sjeremylt   @brief Destroy a CeedBasis
854b11c1e72Sjeremylt 
855b11c1e72Sjeremylt   @param basis CeedBasis to destroy
856b11c1e72Sjeremylt 
857b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
858dfdf5a53Sjeremylt 
859dfdf5a53Sjeremylt   @ref Basic
860b11c1e72Sjeremylt **/
861d7b241e6Sjeremylt int CeedBasisDestroy(CeedBasis *basis) {
862d7b241e6Sjeremylt   int ierr;
863d7b241e6Sjeremylt 
864d7b241e6Sjeremylt   if (!*basis || --(*basis)->refcount > 0) return 0;
865d7b241e6Sjeremylt   if ((*basis)->Destroy) {
866d7b241e6Sjeremylt     ierr = (*basis)->Destroy(*basis); CeedChk(ierr);
867d7b241e6Sjeremylt   }
868d7b241e6Sjeremylt   ierr = CeedFree(&(*basis)->interp1d); CeedChk(ierr);
869d7b241e6Sjeremylt   ierr = CeedFree(&(*basis)->grad1d); CeedChk(ierr);
870d7b241e6Sjeremylt   ierr = CeedFree(&(*basis)->qref1d); CeedChk(ierr);
871d7b241e6Sjeremylt   ierr = CeedFree(&(*basis)->qweight1d); CeedChk(ierr);
872d7b241e6Sjeremylt   ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr);
873d7b241e6Sjeremylt   ierr = CeedFree(basis); CeedChk(ierr);
874d7b241e6Sjeremylt   return 0;
875d7b241e6Sjeremylt }
876d7b241e6Sjeremylt 
87733e6becaSjeremylt /// @cond DOXYGEN_SKIP
878783c99b3SValeria Barra // Indicate that the quadrature points are collocated with the dofs
879783c99b3SValeria Barra CeedBasis CEED_BASIS_COLLOCATED = &ceed_basis_collocated;
88033e6becaSjeremylt /// @endcond
881d7b241e6Sjeremylt /// @}
882