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