xref: /libCEED/rust/libceed-sys/c-src/interface/ceed-basis.c (revision 4d537eea83dd2f1011134892241345b6ac537f56)
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 
63*4d537eeaSYohann   if (dim<1)
64*4d537eeaSYohann     return CeedError(ceed, 1, "Basis dimension must be a positive value");
65*4d537eeaSYohann 
665fe0d4faSjeremylt   if (!ceed->BasisCreateTensorH1) {
675fe0d4faSjeremylt     Ceed delegate;
68aefd8378Sjeremylt     ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr);
695fe0d4faSjeremylt 
705fe0d4faSjeremylt     if (!delegate)
71d7b241e6Sjeremylt       return CeedError(ceed, 1, "Backend does not support BasisCreateTensorH1");
725fe0d4faSjeremylt 
735fe0d4faSjeremylt     ierr = CeedBasisCreateTensorH1(delegate, dim, ncomp, P1d,
745fe0d4faSjeremylt                                    Q1d, interp1d, grad1d, qref1d,
755fe0d4faSjeremylt                                    qweight1d, basis); CeedChk(ierr);
765fe0d4faSjeremylt     return 0;
775fe0d4faSjeremylt   }
78d7b241e6Sjeremylt   ierr = CeedCalloc(1,basis); CeedChk(ierr);
79d7b241e6Sjeremylt   (*basis)->ceed = ceed;
80d7b241e6Sjeremylt   ceed->refcount++;
81d7b241e6Sjeremylt   (*basis)->refcount = 1;
82a8de75f0Sjeremylt   (*basis)->tensorbasis = 1;
83d7b241e6Sjeremylt   (*basis)->dim = dim;
84d7b241e6Sjeremylt   (*basis)->ncomp = ncomp;
85d7b241e6Sjeremylt   (*basis)->P1d = P1d;
86d7b241e6Sjeremylt   (*basis)->Q1d = Q1d;
87a8de75f0Sjeremylt   (*basis)->P = CeedIntPow(P1d, dim);
88a8de75f0Sjeremylt   (*basis)->Q = CeedIntPow(Q1d, dim);
89d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d,&(*basis)->qref1d); CeedChk(ierr);
90d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d,&(*basis)->qweight1d); CeedChk(ierr);
91d7b241e6Sjeremylt   memcpy((*basis)->qref1d, qref1d, Q1d*sizeof(qref1d[0]));
92d7b241e6Sjeremylt   memcpy((*basis)->qweight1d, qweight1d, Q1d*sizeof(qweight1d[0]));
93d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d*P1d,&(*basis)->interp1d); CeedChk(ierr);
94d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d*P1d,&(*basis)->grad1d); CeedChk(ierr);
95d7b241e6Sjeremylt   memcpy((*basis)->interp1d, interp1d, Q1d*P1d*sizeof(interp1d[0]));
9609486605Sjeremylt   memcpy((*basis)->grad1d, grad1d, Q1d*P1d*sizeof(grad1d[0]));
97667bc5fcSjeremylt   ierr = ceed->BasisCreateTensorH1(dim, P1d, Q1d, interp1d, grad1d, qref1d,
98d7b241e6Sjeremylt                                    qweight1d, *basis); CeedChk(ierr);
99d7b241e6Sjeremylt   return 0;
100d7b241e6Sjeremylt }
101d7b241e6Sjeremylt 
102b11c1e72Sjeremylt /**
103b11c1e72Sjeremylt   @brief Create a tensor product Lagrange basis
104b11c1e72Sjeremylt 
105b11c1e72Sjeremylt   @param ceed       A Ceed object where the CeedBasis will be created
106b11c1e72Sjeremylt   @param dim        Topological dimension of element
107b11c1e72Sjeremylt   @param ncomp      Number of field components
108b11c1e72Sjeremylt   @param P          Number of Gauss-Lobatto nodes in one dimension.  The
109b11c1e72Sjeremylt                       polynomial degree of the resulting Q_k element is k=P-1.
110b11c1e72Sjeremylt   @param Q          Number of quadrature points in one dimension.
111b11c1e72Sjeremylt   @param qmode      Distribution of the Q quadrature points (affects order of
112b11c1e72Sjeremylt                       accuracy for the quadrature)
113b11c1e72Sjeremylt   @param[out] basis Address of the variable where the newly created
114b11c1e72Sjeremylt                       CeedBasis will be stored.
115b11c1e72Sjeremylt 
116b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
117dfdf5a53Sjeremylt 
118dfdf5a53Sjeremylt   @ref Basic
119b11c1e72Sjeremylt **/
120d7b241e6Sjeremylt int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt ncomp,
121d7b241e6Sjeremylt                                     CeedInt P, CeedInt Q,
122d7b241e6Sjeremylt                                     CeedQuadMode qmode, CeedBasis *basis) {
123d7b241e6Sjeremylt   // Allocate
124d7b241e6Sjeremylt   int ierr, i, j, k;
125d7b241e6Sjeremylt   CeedScalar c1, c2, c3, c4, dx, *nodes, *interp1d, *grad1d, *qref1d, *qweight1d;
126*4d537eeaSYohann 
127*4d537eeaSYohann   if (dim<1)
128*4d537eeaSYohann     return CeedError(ceed, 1, "Basis dimension must be a positive value");
129*4d537eeaSYohann 
130d7b241e6Sjeremylt   ierr = CeedCalloc(P*Q, &interp1d); CeedChk(ierr);
131d7b241e6Sjeremylt   ierr = CeedCalloc(P*Q, &grad1d); CeedChk(ierr);
132d7b241e6Sjeremylt   ierr = CeedCalloc(P, &nodes); CeedChk(ierr);
133d7b241e6Sjeremylt   ierr = CeedCalloc(Q, &qref1d); CeedChk(ierr);
134d7b241e6Sjeremylt   ierr = CeedCalloc(Q, &qweight1d); CeedChk(ierr);
135d7b241e6Sjeremylt   // Get Nodes and Weights
136d7b241e6Sjeremylt   ierr = CeedLobattoQuadrature(P, nodes, NULL); CeedChk(ierr);
137d7b241e6Sjeremylt   switch (qmode) {
138d7b241e6Sjeremylt   case CEED_GAUSS:
139d7b241e6Sjeremylt     ierr = CeedGaussQuadrature(Q, qref1d, qweight1d); CeedChk(ierr);
140d7b241e6Sjeremylt     break;
141d7b241e6Sjeremylt   case CEED_GAUSS_LOBATTO:
142d7b241e6Sjeremylt     ierr = CeedLobattoQuadrature(Q, qref1d, qweight1d); CeedChk(ierr);
143d7b241e6Sjeremylt     break;
144d7b241e6Sjeremylt   }
145d7b241e6Sjeremylt   // Build B, D matrix
146d7b241e6Sjeremylt   // Fornberg, 1998
147d7b241e6Sjeremylt   for (i = 0; i  < Q; i++) {
148d7b241e6Sjeremylt     c1 = 1.0;
149d7b241e6Sjeremylt     c3 = nodes[0] - qref1d[i];
150d7b241e6Sjeremylt     interp1d[i*P+0] = 1.0;
151d7b241e6Sjeremylt     for (j = 1; j < P; j++) {
152d7b241e6Sjeremylt       c2 = 1.0;
153d7b241e6Sjeremylt       c4 = c3;
154d7b241e6Sjeremylt       c3 = nodes[j] - qref1d[i];
155d7b241e6Sjeremylt       for (k = 0; k < j; k++) {
156d7b241e6Sjeremylt         dx = nodes[j] - nodes[k];
157d7b241e6Sjeremylt         c2 *= dx;
158d7b241e6Sjeremylt         if (k == j - 1) {
159d7b241e6Sjeremylt           grad1d[i*P + j] = c1*(interp1d[i*P + k] - c4*grad1d[i*P + k]) / c2;
160d7b241e6Sjeremylt           interp1d[i*P + j] = - c1*c4*interp1d[i*P + k] / c2;
161d7b241e6Sjeremylt         }
162d7b241e6Sjeremylt         grad1d[i*P + k] = (c3*grad1d[i*P + k] - interp1d[i*P + k]) / dx;
163d7b241e6Sjeremylt         interp1d[i*P + k] = c3*interp1d[i*P + k] / dx;
164d7b241e6Sjeremylt       }
165d7b241e6Sjeremylt       c1 = c2;
166d7b241e6Sjeremylt     }
167d7b241e6Sjeremylt   }
168d7b241e6Sjeremylt   //  // Pass to CeedBasisCreateTensorH1
169d7b241e6Sjeremylt   ierr = CeedBasisCreateTensorH1(ceed, dim, ncomp, P, Q, interp1d, grad1d, qref1d,
170d7b241e6Sjeremylt                                  qweight1d, basis); CeedChk(ierr);
171d7b241e6Sjeremylt   ierr = CeedFree(&interp1d); CeedChk(ierr);
172d7b241e6Sjeremylt   ierr = CeedFree(&grad1d); CeedChk(ierr);
173d7b241e6Sjeremylt   ierr = CeedFree(&nodes); CeedChk(ierr);
174d7b241e6Sjeremylt   ierr = CeedFree(&qref1d); CeedChk(ierr);
175d7b241e6Sjeremylt   ierr = CeedFree(&qweight1d); CeedChk(ierr);
176d7b241e6Sjeremylt   return 0;
177d7b241e6Sjeremylt }
178d7b241e6Sjeremylt 
179b11c1e72Sjeremylt /**
180a8de75f0Sjeremylt   @brief Create a non tensor product basis for H^1 discretizations
181a8de75f0Sjeremylt 
182a8de75f0Sjeremylt   @param ceed       A Ceed object where the CeedBasis will be created
183a8de75f0Sjeremylt   @param topo       Topology of element, e.g. hypercube, simplex, ect
184a8de75f0Sjeremylt   @param ncomp      Number of field components (1 for scalar fields)
1858795c945Sjeremylt   @param nnodes       Total number of nodes
186a8de75f0Sjeremylt   @param nqpts      Total number of quadrature points
1878795c945Sjeremylt   @param interp     Row-major nqpts × nnodes matrix expressing the values of
1888795c945Sjeremylt                       nodal basis functions at quadrature points
1898795c945Sjeremylt   @param grad       Row-major (nqpts x dim) × nnodes matrix expressing
1908795c945Sjeremylt                       derivatives of nodal basis functions at quadrature points
1918795c945Sjeremylt   @param qref       Array of length nqpts holding the locations of quadrature
1928795c945Sjeremylt                       points on the reference element [-1, 1]
193a8de75f0Sjeremylt   @param qweight    Array of length nqpts holding the quadrature weights on the
194a8de75f0Sjeremylt                       reference element
195a8de75f0Sjeremylt   @param[out] basis Address of the variable where the newly created
196a8de75f0Sjeremylt                       CeedBasis will be stored.
197a8de75f0Sjeremylt 
198a8de75f0Sjeremylt   @return An error code: 0 - success, otherwise - failure
199a8de75f0Sjeremylt 
200a8de75f0Sjeremylt   @ref Basic
201a8de75f0Sjeremylt **/
202a8de75f0Sjeremylt int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt ncomp,
2038795c945Sjeremylt                       CeedInt nnodes, CeedInt nqpts,
204a8de75f0Sjeremylt                       const CeedScalar *interp,
205a8de75f0Sjeremylt                       const CeedScalar *grad, const CeedScalar *qref,
206a8de75f0Sjeremylt                       const CeedScalar *qweight, CeedBasis *basis) {
207a8de75f0Sjeremylt   int ierr;
2088795c945Sjeremylt   CeedInt P = nnodes, Q = nqpts, dim = 0;
209a8de75f0Sjeremylt 
2105fe0d4faSjeremylt   if (!ceed->BasisCreateH1) {
2115fe0d4faSjeremylt     Ceed delegate;
212aefd8378Sjeremylt     ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr);
2135fe0d4faSjeremylt 
2145fe0d4faSjeremylt     if (!delegate)
215a8de75f0Sjeremylt       return CeedError(ceed, 1, "Backend does not support BasisCreateH1");
2165fe0d4faSjeremylt 
2178795c945Sjeremylt     ierr = CeedBasisCreateH1(delegate, topo, ncomp, nnodes,
2185fe0d4faSjeremylt                              nqpts, interp, grad, qref,
2195fe0d4faSjeremylt                              qweight, basis); CeedChk(ierr);
2205fe0d4faSjeremylt     return 0;
2215fe0d4faSjeremylt   }
2225fe0d4faSjeremylt 
223a8de75f0Sjeremylt   ierr = CeedCalloc(1,basis); CeedChk(ierr);
224a8de75f0Sjeremylt 
225a8de75f0Sjeremylt   ierr = CeedBasisGetTopologyDimension(topo, &dim); CeedChk(ierr);
226a8de75f0Sjeremylt 
227a8de75f0Sjeremylt   (*basis)->ceed = ceed;
228a8de75f0Sjeremylt   ceed->refcount++;
229a8de75f0Sjeremylt   (*basis)->refcount = 1;
230a8de75f0Sjeremylt   (*basis)->tensorbasis = 0;
231a8de75f0Sjeremylt   (*basis)->dim = dim;
232a8de75f0Sjeremylt   (*basis)->ncomp = ncomp;
233a8de75f0Sjeremylt   (*basis)->P = P;
234a8de75f0Sjeremylt   (*basis)->Q = Q;
235a8de75f0Sjeremylt   ierr = CeedMalloc(Q*dim,&(*basis)->qref1d); CeedChk(ierr);
236a8de75f0Sjeremylt   ierr = CeedMalloc(Q,&(*basis)->qweight1d); CeedChk(ierr);
237a8de75f0Sjeremylt   memcpy((*basis)->qref1d, qref, Q*dim*sizeof(qref[0]));
238a8de75f0Sjeremylt   memcpy((*basis)->qweight1d, qweight, Q*sizeof(qweight[0]));
239a8de75f0Sjeremylt   ierr = CeedMalloc(Q*P,&(*basis)->interp1d); CeedChk(ierr);
240a8de75f0Sjeremylt   ierr = CeedMalloc(dim*Q*P,&(*basis)->grad1d); CeedChk(ierr);
241a8de75f0Sjeremylt   memcpy((*basis)->interp1d, interp, Q*P*sizeof(interp[0]));
242a8de75f0Sjeremylt   memcpy((*basis)->grad1d, grad, dim*Q*P*sizeof(grad[0]));
243667bc5fcSjeremylt   ierr = ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, qref,
244a8de75f0Sjeremylt                              qweight, *basis); CeedChk(ierr);
245a8de75f0Sjeremylt   return 0;
246a8de75f0Sjeremylt }
247a8de75f0Sjeremylt 
248a8de75f0Sjeremylt /**
249b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre quadrature
250b11c1e72Sjeremylt 
251b11c1e72Sjeremylt   @param Q              Number of quadrature points (integrates polynomials of
252b11c1e72Sjeremylt                           degree 2*Q-1 exactly)
253b11c1e72Sjeremylt   @param[out] qref1d    Array of length Q to hold the abscissa on [-1, 1]
254b11c1e72Sjeremylt   @param[out] qweight1d Array of length Q to hold the weights
255b11c1e72Sjeremylt 
256b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
257dfdf5a53Sjeremylt 
258dfdf5a53Sjeremylt   @ref Utility
259b11c1e72Sjeremylt **/
260d7b241e6Sjeremylt int CeedGaussQuadrature(CeedInt Q, CeedScalar *qref1d, CeedScalar *qweight1d) {
261d7b241e6Sjeremylt   // Allocate
262d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0*atan(1.0);
263d7b241e6Sjeremylt   // Build qref1d, qweight1d
264d7b241e6Sjeremylt   for (int i = 0; i <= Q/2; i++) {
265d7b241e6Sjeremylt     // Guess
266d7b241e6Sjeremylt     xi = cos(PI*(CeedScalar)(2*i+1)/((CeedScalar)(2*Q)));
267d7b241e6Sjeremylt     // Pn(xi)
268d7b241e6Sjeremylt     P0 = 1.0;
269d7b241e6Sjeremylt     P1 = xi;
270d7b241e6Sjeremylt     P2 = 0.0;
271d7b241e6Sjeremylt     for (int j = 2; j <= Q; j++) {
272d7b241e6Sjeremylt       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
273d7b241e6Sjeremylt       P0 = P1;
274d7b241e6Sjeremylt       P1 = P2;
275d7b241e6Sjeremylt     }
276d7b241e6Sjeremylt     // First Newton Step
277d7b241e6Sjeremylt     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
278d7b241e6Sjeremylt     xi = xi-P2/dP2;
279d7b241e6Sjeremylt     // Newton to convergence
280d7b241e6Sjeremylt     for (int k=0; k<100 && fabs(P2)>1e-15; k++) {
281d7b241e6Sjeremylt       P0 = 1.0;
282d7b241e6Sjeremylt       P1 = xi;
283d7b241e6Sjeremylt       for (int j = 2; j <= Q; j++) {
284d7b241e6Sjeremylt         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
285d7b241e6Sjeremylt         P0 = P1;
286d7b241e6Sjeremylt         P1 = P2;
287d7b241e6Sjeremylt       }
288d7b241e6Sjeremylt       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
289d7b241e6Sjeremylt       xi = xi-P2/dP2;
290d7b241e6Sjeremylt     }
291d7b241e6Sjeremylt     // Save xi, wi
292d7b241e6Sjeremylt     wi = 2.0/((1.0-xi*xi)*dP2*dP2);
293d7b241e6Sjeremylt     qweight1d[i] = wi;
294d7b241e6Sjeremylt     qweight1d[Q-1-i] = wi;
295d7b241e6Sjeremylt     qref1d[i] = -xi;
296d7b241e6Sjeremylt     qref1d[Q-1-i]= xi;
297d7b241e6Sjeremylt   }
298d7b241e6Sjeremylt   return 0;
299d7b241e6Sjeremylt }
300d7b241e6Sjeremylt 
301b11c1e72Sjeremylt /**
302b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre-Lobatto quadrature
303b11c1e72Sjeremylt 
304b11c1e72Sjeremylt   @param Q              Number of quadrature points (integrates polynomials of
305b11c1e72Sjeremylt                           degree 2*Q-3 exactly)
306b11c1e72Sjeremylt   @param[out] qref1d    Array of length Q to hold the abscissa on [-1, 1]
307b11c1e72Sjeremylt   @param[out] qweight1d Array of length Q to hold the weights
308b11c1e72Sjeremylt 
309b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
310dfdf5a53Sjeremylt 
311dfdf5a53Sjeremylt   @ref Utility
312b11c1e72Sjeremylt **/
313d7b241e6Sjeremylt int CeedLobattoQuadrature(CeedInt Q, CeedScalar *qref1d,
314d7b241e6Sjeremylt                           CeedScalar *qweight1d) {
315d7b241e6Sjeremylt   // Allocate
316d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0*atan(1.0);
317d7b241e6Sjeremylt   // Build qref1d, qweight1d
318d7b241e6Sjeremylt   // Set endpoints
319d7b241e6Sjeremylt   wi = 2.0/((CeedScalar)(Q*(Q-1)));
320d7b241e6Sjeremylt   if (qweight1d) {
321d7b241e6Sjeremylt     qweight1d[0] = wi;
322d7b241e6Sjeremylt     qweight1d[Q-1] = wi;
323d7b241e6Sjeremylt   }
324d7b241e6Sjeremylt   qref1d[0] = -1.0;
325d7b241e6Sjeremylt   qref1d[Q-1] = 1.0;
326d7b241e6Sjeremylt   // Interior
327d7b241e6Sjeremylt   for (int i = 1; i <= (Q-1)/2; i++) {
328d7b241e6Sjeremylt     // Guess
329d7b241e6Sjeremylt     xi = cos(PI*(CeedScalar)(i)/(CeedScalar)(Q-1));
330d7b241e6Sjeremylt     // Pn(xi)
331d7b241e6Sjeremylt     P0 = 1.0;
332d7b241e6Sjeremylt     P1 = xi;
333d7b241e6Sjeremylt     P2 = 0.0;
334d7b241e6Sjeremylt     for (int j = 2; j < Q; j++) {
335d7b241e6Sjeremylt       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
336d7b241e6Sjeremylt       P0 = P1;
337d7b241e6Sjeremylt       P1 = P2;
338d7b241e6Sjeremylt     }
339d7b241e6Sjeremylt     // First Newton step
340d7b241e6Sjeremylt     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
341d7b241e6Sjeremylt     d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
342d7b241e6Sjeremylt     xi = xi-dP2/d2P2;
343d7b241e6Sjeremylt     // Newton to convergence
344d7b241e6Sjeremylt     for (int k=0; k<100 && fabs(dP2)>1e-15; k++) {
345d7b241e6Sjeremylt       P0 = 1.0;
346d7b241e6Sjeremylt       P1 = xi;
347d7b241e6Sjeremylt       for (int j = 2; j < Q; j++) {
348d7b241e6Sjeremylt         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
349d7b241e6Sjeremylt         P0 = P1;
350d7b241e6Sjeremylt         P1 = P2;
351d7b241e6Sjeremylt       }
352d7b241e6Sjeremylt       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
353d7b241e6Sjeremylt       d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
354d7b241e6Sjeremylt       xi = xi-dP2/d2P2;
355d7b241e6Sjeremylt     }
356d7b241e6Sjeremylt     // Save xi, wi
357d7b241e6Sjeremylt     wi = 2.0/(((CeedScalar)(Q*(Q-1)))*P2*P2);
358d7b241e6Sjeremylt     if (qweight1d) {
359d7b241e6Sjeremylt       qweight1d[i] = wi;
360d7b241e6Sjeremylt       qweight1d[Q-1-i] = wi;
361d7b241e6Sjeremylt     }
362d7b241e6Sjeremylt     qref1d[i] = -xi;
363d7b241e6Sjeremylt     qref1d[Q-1-i]= xi;
364d7b241e6Sjeremylt   }
365d7b241e6Sjeremylt   return 0;
366d7b241e6Sjeremylt }
367d7b241e6Sjeremylt 
368dfdf5a53Sjeremylt /**
369dfdf5a53Sjeremylt   @brief View an array stored in a CeedBasis
370dfdf5a53Sjeremylt 
371dfdf5a53Sjeremylt   @param name      Name of array
372dfdf5a53Sjeremylt   @param fpformat  Printing format
373dfdf5a53Sjeremylt   @param m         Number of rows in array
374dfdf5a53Sjeremylt   @param n         Number of columns in array
375dfdf5a53Sjeremylt   @param a         Array to be viewed
376dfdf5a53Sjeremylt   @param stream    Stream to view to, e.g., stdout
377dfdf5a53Sjeremylt 
378dfdf5a53Sjeremylt   @return An error code: 0 - success, otherwise - failure
379dfdf5a53Sjeremylt 
380dfdf5a53Sjeremylt   @ref Utility
381dfdf5a53Sjeremylt **/
382d7b241e6Sjeremylt static int CeedScalarView(const char *name, const char *fpformat, CeedInt m,
383d7b241e6Sjeremylt                           CeedInt n, const CeedScalar *a, FILE *stream) {
384d7b241e6Sjeremylt   for (int i=0; i<m; i++) {
385d7b241e6Sjeremylt     if (m > 1) fprintf(stream, "%12s[%d]:", name, i);
386d7b241e6Sjeremylt     else fprintf(stream, "%12s:", name);
387d7b241e6Sjeremylt     for (int j=0; j<n; j++) {
388d7b241e6Sjeremylt       fprintf(stream, fpformat, fabs(a[i*n+j]) > 1E-14 ? a[i*n+j] : 0);
389d7b241e6Sjeremylt     }
390d7b241e6Sjeremylt     fputs("\n", stream);
391d7b241e6Sjeremylt   }
392d7b241e6Sjeremylt   return 0;
393d7b241e6Sjeremylt }
394d7b241e6Sjeremylt 
395b11c1e72Sjeremylt /**
396b11c1e72Sjeremylt   @brief View a CeedBasis
397b11c1e72Sjeremylt 
398b11c1e72Sjeremylt   @param basis  CeedBasis to view
399b11c1e72Sjeremylt   @param stream Stream to view to, e.g., stdout
400b11c1e72Sjeremylt 
401b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
402dfdf5a53Sjeremylt 
403dfdf5a53Sjeremylt   @ref Utility
404b11c1e72Sjeremylt **/
405d7b241e6Sjeremylt int CeedBasisView(CeedBasis basis, FILE *stream) {
406d7b241e6Sjeremylt   int ierr;
407d7b241e6Sjeremylt 
408a8de75f0Sjeremylt   if (basis->tensorbasis) {
409d7b241e6Sjeremylt     fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P1d,
410d7b241e6Sjeremylt             basis->Q1d);
411d7b241e6Sjeremylt     ierr = CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q1d, basis->qref1d,
412d7b241e6Sjeremylt                           stream); CeedChk(ierr);
4138795c945Sjeremylt     ierr = CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q1d,
4148795c945Sjeremylt                           basis->qweight1d, stream); CeedChk(ierr);
415d7b241e6Sjeremylt     ierr = CeedScalarView("interp1d", "\t% 12.8f", basis->Q1d, basis->P1d,
416d7b241e6Sjeremylt                           basis->interp1d, stream); CeedChk(ierr);
417d7b241e6Sjeremylt     ierr = CeedScalarView("grad1d", "\t% 12.8f", basis->Q1d, basis->P1d,
418d7b241e6Sjeremylt                           basis->grad1d, stream); CeedChk(ierr);
419a8de75f0Sjeremylt   } else {
420a8de75f0Sjeremylt     fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P,
421a8de75f0Sjeremylt             basis->Q);
422a8de75f0Sjeremylt     ierr = CeedScalarView("qref", "\t% 12.8f", 1, basis->Q*basis->dim,
423a8de75f0Sjeremylt                           basis->qref1d,
424a8de75f0Sjeremylt                           stream); CeedChk(ierr);
425a8de75f0Sjeremylt     ierr = CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->qweight1d,
426a8de75f0Sjeremylt                           stream); CeedChk(ierr);
427a8de75f0Sjeremylt     ierr = CeedScalarView("interp", "\t% 12.8f", basis->Q, basis->P,
428a8de75f0Sjeremylt                           basis->interp1d, stream); CeedChk(ierr);
429a8de75f0Sjeremylt     ierr = CeedScalarView("grad", "\t% 12.8f", basis->dim*basis->Q, basis->P,
430a8de75f0Sjeremylt                           basis->grad1d, stream); CeedChk(ierr);
431a8de75f0Sjeremylt   }
432d7b241e6Sjeremylt   return 0;
433d7b241e6Sjeremylt }
434d7b241e6Sjeremylt 
435dfdf5a53Sjeremylt /**
436dfdf5a53Sjeremylt   @brief Compute Householder Reflection
437dfdf5a53Sjeremylt 
438dfdf5a53Sjeremylt     Computes A = (I - b v v^T) A
439dfdf5a53Sjeremylt     where A is an mxn matrix indexed as A[i*row + j*col]
440dfdf5a53Sjeremylt 
441dfdf5a53Sjeremylt   @param[out] A  Matrix to apply Householder reflection to, in place
442dfdf5a53Sjeremylt   @param v       Householder vector
443dfdf5a53Sjeremylt   @param b       Scaling factor
444dfdf5a53Sjeremylt   @param m       Number of rows in A
445dfdf5a53Sjeremylt   @param n       Number of columns in A
446dfdf5a53Sjeremylt   @param row     Col stride
447dfdf5a53Sjeremylt   @param col     Row stride
448dfdf5a53Sjeremylt 
449dfdf5a53Sjeremylt   @return An error code: 0 - success, otherwise - failure
450dfdf5a53Sjeremylt 
451dfdf5a53Sjeremylt   @ref Developer
452dfdf5a53Sjeremylt **/
453d7b241e6Sjeremylt static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v,
454d7b241e6Sjeremylt                                   CeedScalar b, CeedInt m, CeedInt n,
455d7b241e6Sjeremylt                                   CeedInt row, CeedInt col) {
456d7b241e6Sjeremylt   for (CeedInt j=0; j<n; j++) {
457d7b241e6Sjeremylt     CeedScalar w = A[0*row + j*col];
458d7b241e6Sjeremylt     for (CeedInt i=1; i<m; i++) w += v[i] * A[i*row + j*col];
459d7b241e6Sjeremylt     A[0*row + j*col] -= b * w;
460d7b241e6Sjeremylt     for (CeedInt i=1; i<m; i++) A[i*row + j*col] -= b * w * v[i];
461d7b241e6Sjeremylt   }
462d7b241e6Sjeremylt   return 0;
463d7b241e6Sjeremylt }
464d7b241e6Sjeremylt 
465dfdf5a53Sjeremylt /**
466dfdf5a53Sjeremylt   @brief Apply Householder Q matrix
467dfdf5a53Sjeremylt 
468dfdf5a53Sjeremylt     Compute A = Q A where Q is mxk and A is mxn. k<m
469dfdf5a53Sjeremylt 
470dfdf5a53Sjeremylt   @param[out] A  Matrix to apply Householder Q to, in place
471dfdf5a53Sjeremylt   @param Q       Householder Q matrix
472dfdf5a53Sjeremylt   @param tau     Householder scaling factors
473dfdf5a53Sjeremylt   @param tmode   Transpose mode for application
474dfdf5a53Sjeremylt   @param m       Number of rows in A
475dfdf5a53Sjeremylt   @param n       Number of columns in A
476dfdf5a53Sjeremylt   @param k       Index of row targeted
477dfdf5a53Sjeremylt   @param row     Col stride
478dfdf5a53Sjeremylt   @param col     Row stride
479dfdf5a53Sjeremylt 
480dfdf5a53Sjeremylt   @return An error code: 0 - success, otherwise - failure
481dfdf5a53Sjeremylt 
482dfdf5a53Sjeremylt   @ref Developer
483dfdf5a53Sjeremylt **/
484d7b241e6Sjeremylt static int CeedHouseholderApplyQ(CeedScalar *A, const CeedScalar *Q,
485d7b241e6Sjeremylt                                  const CeedScalar *tau, CeedTransposeMode tmode,
486d7b241e6Sjeremylt                                  CeedInt m, CeedInt n, CeedInt k,
487d7b241e6Sjeremylt                                  CeedInt row, CeedInt col) {
488d7b241e6Sjeremylt   CeedScalar v[m];
489d7b241e6Sjeremylt   for (CeedInt ii=0; ii<k; ii++) {
490d7b241e6Sjeremylt     CeedInt i = tmode == CEED_TRANSPOSE ? ii : k-1-ii;
491d7b241e6Sjeremylt     for (CeedInt j=i+1; j<m; j++) {
492d7b241e6Sjeremylt       v[j] = Q[j*k+i];
493d7b241e6Sjeremylt     }
494d7b241e6Sjeremylt     // Apply Householder reflector (I - tau v v^T) colograd1d^T
495d7b241e6Sjeremylt     CeedHouseholderReflect(&A[i*row], &v[i], tau[i], m-i, n, row, col);
496d7b241e6Sjeremylt   }
497d7b241e6Sjeremylt   return 0;
498d7b241e6Sjeremylt }
499d7b241e6Sjeremylt 
500b11c1e72Sjeremylt /**
501b11c1e72Sjeremylt   @brief Return QR Factorization of matrix
502b11c1e72Sjeremylt 
503b11c1e72Sjeremylt   @param[out] mat  Row-major matrix to be factorized in place
5048c91a0c9SJeremy L Thompson   @param[out] tau  Vector of length m of scaling factors
505b11c1e72Sjeremylt   @param m         Number of rows
506b11c1e72Sjeremylt   @param n         Number of columns
507b11c1e72Sjeremylt 
508b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
509dfdf5a53Sjeremylt 
510dfdf5a53Sjeremylt   @ref Utility
511b11c1e72Sjeremylt **/
512a7bd39daSjeremylt int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau,
513d7b241e6Sjeremylt                         CeedInt m, CeedInt n) {
514d7b241e6Sjeremylt   CeedInt i, j;
515d7b241e6Sjeremylt   CeedScalar v[m];
516d7b241e6Sjeremylt 
517a7bd39daSjeremylt   // Check m >= n
518a7bd39daSjeremylt   if (n > m)
519a7bd39daSjeremylt     return CeedError(ceed, 1, "Cannot compute QR factorization with n > m");
520a7bd39daSjeremylt 
521d7b241e6Sjeremylt   for (i=0; i<n; i++) {
522d7b241e6Sjeremylt     // Calculate Householder vector, magnitude
523d7b241e6Sjeremylt     CeedScalar sigma = 0.0;
524d7b241e6Sjeremylt     v[i] = mat[i+n*i];
525d7b241e6Sjeremylt     for (j=i+1; j<m; j++) {
526d7b241e6Sjeremylt       v[j] = mat[i+n*j];
527d7b241e6Sjeremylt       sigma += v[j] * v[j];
528d7b241e6Sjeremylt     }
529d7b241e6Sjeremylt     CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m]
530d7b241e6Sjeremylt     CeedScalar Rii = -copysign(norm, v[i]);
531d7b241e6Sjeremylt     v[i] -= Rii;
532d7b241e6Sjeremylt     // norm of v[i:m] after modification above and scaling below
533d7b241e6Sjeremylt     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
534d7b241e6Sjeremylt     //   tau = 2 / (norm*norm)
535d7b241e6Sjeremylt     tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma);
536d7b241e6Sjeremylt     for (j=i+1; j<m; j++) v[j] /= v[i];
537d7b241e6Sjeremylt 
538d7b241e6Sjeremylt     // Apply Householder reflector to lower right panel
539d7b241e6Sjeremylt     CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1);
540d7b241e6Sjeremylt     // Save v
541d7b241e6Sjeremylt     mat[i+n*i] = Rii;
542d7b241e6Sjeremylt     for (j=i+1; j<m; j++) {
543d7b241e6Sjeremylt       mat[i+n*j] = v[j];
544d7b241e6Sjeremylt     }
545d7b241e6Sjeremylt   }
546d7b241e6Sjeremylt 
547d7b241e6Sjeremylt   return 0;
548d7b241e6Sjeremylt }
549d7b241e6Sjeremylt 
550b11c1e72Sjeremylt /**
551783c99b3SValeria Barra   @brief Return collocated grad matrix
552b11c1e72Sjeremylt 
553b11c1e72Sjeremylt   @param basis           CeedBasis
554b11c1e72Sjeremylt   @param[out] colograd1d Row-major Q1d × Q1d matrix expressing derivatives of
555b11c1e72Sjeremylt                            basis functions at quadrature points
556b11c1e72Sjeremylt 
557b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
558dfdf5a53Sjeremylt 
559dfdf5a53Sjeremylt   @ref Advanced
560b11c1e72Sjeremylt **/
561783c99b3SValeria Barra int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *colograd1d) {
562d7b241e6Sjeremylt   int i, j, k;
563a7bd39daSjeremylt   Ceed ceed;
564d7b241e6Sjeremylt   CeedInt ierr, P1d=(basis)->P1d, Q1d=(basis)->Q1d;
565d7b241e6Sjeremylt   CeedScalar *interp1d, *grad1d, tau[Q1d];
566d7b241e6Sjeremylt 
567d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d*P1d, &interp1d); CeedChk(ierr);
568d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d*P1d, &grad1d); CeedChk(ierr);
569d7b241e6Sjeremylt   memcpy(interp1d, (basis)->interp1d, Q1d*P1d*sizeof(basis)->interp1d[0]);
570d7b241e6Sjeremylt   memcpy(grad1d, (basis)->grad1d, Q1d*P1d*sizeof(basis)->interp1d[0]);
571d7b241e6Sjeremylt 
572d7b241e6Sjeremylt   // QR Factorization, interp1d = Q R
573a7bd39daSjeremylt   ierr = CeedBasisGetCeed(basis, &ceed); CeedChk(ierr);
574a7bd39daSjeremylt   ierr = CeedQRFactorization(ceed, interp1d, tau, Q1d, P1d); CeedChk(ierr);
575d7b241e6Sjeremylt 
576d7b241e6Sjeremylt   // Apply Rinv, colograd1d = grad1d Rinv
577d7b241e6Sjeremylt   for (i=0; i<Q1d; i++) { // Row i
578d7b241e6Sjeremylt     colograd1d[Q1d*i] = grad1d[P1d*i]/interp1d[0];
579d7b241e6Sjeremylt     for (j=1; j<P1d; j++) { // Column j
580d7b241e6Sjeremylt       colograd1d[j+Q1d*i] = grad1d[j+P1d*i];
581d7b241e6Sjeremylt       for (k=0; k<j; k++) {
582d7b241e6Sjeremylt         colograd1d[j+Q1d*i] -= interp1d[j+P1d*k]*colograd1d[k+Q1d*i];
583d7b241e6Sjeremylt       }
584d7b241e6Sjeremylt       colograd1d[j+Q1d*i] /= interp1d[j+P1d*j];
585d7b241e6Sjeremylt     }
586d7b241e6Sjeremylt     for (j=P1d; j<Q1d; j++) {
587d7b241e6Sjeremylt       colograd1d[j+Q1d*i] = 0;
588d7b241e6Sjeremylt     }
589d7b241e6Sjeremylt   }
590d7b241e6Sjeremylt 
591d7b241e6Sjeremylt   // Apply Qtranspose, colograd = colograd Qtranspose
592d7b241e6Sjeremylt   CeedHouseholderApplyQ(colograd1d, interp1d, tau, CEED_NOTRANSPOSE,
593d7b241e6Sjeremylt                         Q1d, Q1d, P1d, 1, Q1d);
594d7b241e6Sjeremylt 
595d7b241e6Sjeremylt   ierr = CeedFree(&interp1d); CeedChk(ierr);
596d7b241e6Sjeremylt   ierr = CeedFree(&grad1d); CeedChk(ierr);
597d7b241e6Sjeremylt 
598d7b241e6Sjeremylt   return 0;
599d7b241e6Sjeremylt }
600d7b241e6Sjeremylt 
601b11c1e72Sjeremylt /**
602b11c1e72Sjeremylt   @brief Apply basis evaluation from nodes to quadrature points or vice-versa
603b11c1e72Sjeremylt 
604b11c1e72Sjeremylt   @param basis  CeedBasis to evaluate
605b11c1e72Sjeremylt   @param nelem  The number of elements to apply the basis evaluation to;
606b11c1e72Sjeremylt                   the backend will specify the ordering in
607b11c1e72Sjeremylt                   ElemRestrictionCreateBlocked
608b11c1e72Sjeremylt   @param tmode  \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature
609b11c1e72Sjeremylt                   points, \ref CEED_TRANSPOSE to apply the transpose, mapping
610b11c1e72Sjeremylt                   from quadrature points to nodes
611b11c1e72Sjeremylt   @param emode  \ref CEED_EVAL_INTERP to obtain interpolated values,
612b11c1e72Sjeremylt                   \ref CEED_EVAL_GRAD to obtain gradients.
613b11c1e72Sjeremylt   @param[in] u  Input array
614b11c1e72Sjeremylt   @param[out] v Output array
615b11c1e72Sjeremylt 
616b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
617dfdf5a53Sjeremylt 
618dfdf5a53Sjeremylt   @ref Advanced
619b11c1e72Sjeremylt **/
620d7b241e6Sjeremylt int CeedBasisApply(CeedBasis basis, CeedInt nelem, CeedTransposeMode tmode,
621aedaa0e5Sjeremylt                    CeedEvalMode emode, CeedVector u, CeedVector v) {
622d7b241e6Sjeremylt   int ierr;
6238795c945Sjeremylt   CeedInt ulength = 0, vlength, nnodes, nqpt;
624d7b241e6Sjeremylt   if (!basis->Apply) return CeedError(basis->ceed, 1,
625d7b241e6Sjeremylt                                         "Backend does not support BasisApply");
626b502e64cSValeria Barra   // check compatibility of topological and geometrical dimensions
6278795c945Sjeremylt   ierr = CeedBasisGetNumNodes(basis, &nnodes); CeedChk(ierr);
628b502e64cSValeria Barra   ierr = CeedBasisGetNumQuadraturePoints(basis, &nqpt); CeedChk(ierr);
629b502e64cSValeria Barra   ierr = CeedVectorGetLength(v, &vlength); CeedChk(ierr);
630b502e64cSValeria Barra 
631b502e64cSValeria Barra   if (u) {
632b502e64cSValeria Barra     ierr = CeedVectorGetLength(u, &ulength); CeedChk(ierr);
633b502e64cSValeria Barra   }
634b502e64cSValeria Barra 
635f90c8643Sjeremylt   if ((tmode == CEED_TRANSPOSE   && (vlength % nnodes != 0
636f90c8643Sjeremylt                                      || ulength % nqpt != 0))
637cdf4f918Sjeremylt       ||
6388795c945Sjeremylt       (tmode == CEED_NOTRANSPOSE && (ulength % nnodes != 0 || vlength % nqpt != 0)))
639b502e64cSValeria Barra     return CeedError(basis->ceed, 1,
640b502e64cSValeria Barra                      "Length of input/output vectors incompatible with basis dimensions");
641b502e64cSValeria Barra 
642d7b241e6Sjeremylt   ierr = basis->Apply(basis, nelem, tmode, emode, u, v); CeedChk(ierr);
643d7b241e6Sjeremylt   return 0;
644d7b241e6Sjeremylt }
645d7b241e6Sjeremylt 
646b11c1e72Sjeremylt /**
6474ce2993fSjeremylt   @brief Get Ceed associated with a CeedBasis
648b11c1e72Sjeremylt 
649b11c1e72Sjeremylt   @param basis      CeedBasis
6504ce2993fSjeremylt   @param[out] ceed  Variable to store Ceed
6514ce2993fSjeremylt 
6524ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
6534ce2993fSjeremylt 
65423617272Sjeremylt   @ref Advanced
6554ce2993fSjeremylt **/
6564ce2993fSjeremylt int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) {
6574ce2993fSjeremylt   *ceed = basis->ceed;
6584ce2993fSjeremylt 
6594ce2993fSjeremylt   return 0;
6604ce2993fSjeremylt };
6614ce2993fSjeremylt 
6624ce2993fSjeremylt /**
6634ce2993fSjeremylt   @brief Get dimension for given CeedBasis
6644ce2993fSjeremylt 
6654ce2993fSjeremylt   @param basis     CeedBasis
6664ce2993fSjeremylt   @param[out] dim  Variable to store dimension of basis
6674ce2993fSjeremylt 
6684ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
6694ce2993fSjeremylt 
67023617272Sjeremylt   @ref Advanced
6714ce2993fSjeremylt **/
6724ce2993fSjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) {
6734ce2993fSjeremylt   *dim = basis->dim;
6744ce2993fSjeremylt 
6754ce2993fSjeremylt   return 0;
6764ce2993fSjeremylt };
6774ce2993fSjeremylt 
6784ce2993fSjeremylt /**
6794ce2993fSjeremylt   @brief Get tensor status for given CeedBasis
6804ce2993fSjeremylt 
6814ce2993fSjeremylt   @param basis        CeedBasis
6824ce2993fSjeremylt   @param[out] tensor  Variable to store tensor status
6834ce2993fSjeremylt 
6844ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
6854ce2993fSjeremylt 
68623617272Sjeremylt   @ref Advanced
6874ce2993fSjeremylt **/
6884ce2993fSjeremylt int CeedBasisGetTensorStatus(CeedBasis basis, bool *tensor) {
6894ce2993fSjeremylt   *tensor = basis->tensorbasis;
6904ce2993fSjeremylt 
6914ce2993fSjeremylt   return 0;
6924ce2993fSjeremylt };
6934ce2993fSjeremylt 
6944ce2993fSjeremylt /**
6954ce2993fSjeremylt   @brief Get number of components for given CeedBasis
6964ce2993fSjeremylt 
6974ce2993fSjeremylt   @param basis     CeedBasis
6984ce2993fSjeremylt   @param[out] dim  Variable to store number of components of basis
6994ce2993fSjeremylt 
7004ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
7014ce2993fSjeremylt 
70223617272Sjeremylt   @ref Advanced
7034ce2993fSjeremylt **/
7044ce2993fSjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *numcomp) {
7054ce2993fSjeremylt   *numcomp = basis->ncomp;
7064ce2993fSjeremylt 
7074ce2993fSjeremylt   return 0;
7084ce2993fSjeremylt };
7094ce2993fSjeremylt 
7104ce2993fSjeremylt /**
7114ce2993fSjeremylt   @brief Get total number of nodes (in 1 dimension) of a CeedBasis
7124ce2993fSjeremylt 
7134ce2993fSjeremylt   @param basis     CeedBasis
7144ce2993fSjeremylt   @param[out] P1d  Variable to store number of nodes
7154ce2993fSjeremylt 
7164ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
7174ce2993fSjeremylt 
71823617272Sjeremylt   @ref Advanced
7194ce2993fSjeremylt **/
7204ce2993fSjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P1d) {
7214ce2993fSjeremylt   if (!basis->tensorbasis) return CeedError(basis->ceed, 1,
7224ce2993fSjeremylt                                     "Cannot supply P1d for non-tensor basis");
7234ce2993fSjeremylt   *P1d = basis->P1d;
7244ce2993fSjeremylt   return 0;
7254ce2993fSjeremylt }
7264ce2993fSjeremylt 
7274ce2993fSjeremylt /**
7284ce2993fSjeremylt   @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis
7294ce2993fSjeremylt 
7304ce2993fSjeremylt   @param basis     CeedBasis
7314ce2993fSjeremylt   @param[out] Q1d  Variable to store number of quadrature points
7324ce2993fSjeremylt 
7334ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
7344ce2993fSjeremylt 
73523617272Sjeremylt   @ref Advanced
7364ce2993fSjeremylt **/
7374ce2993fSjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q1d) {
7384ce2993fSjeremylt   if (!basis->tensorbasis) return CeedError(basis->ceed, 1,
7394ce2993fSjeremylt                                     "Cannot supply Q1d for non-tensor basis");
7404ce2993fSjeremylt   *Q1d = basis->Q1d;
7414ce2993fSjeremylt   return 0;
7424ce2993fSjeremylt }
7434ce2993fSjeremylt 
7444ce2993fSjeremylt /**
7454ce2993fSjeremylt   @brief Get total number of nodes (in dim dimensions) of a CeedBasis
7464ce2993fSjeremylt 
7474ce2993fSjeremylt   @param basis   CeedBasis
7484ce2993fSjeremylt   @param[out] P  Variable to store number of nodes
749b11c1e72Sjeremylt 
750b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
751dfdf5a53Sjeremylt 
752dfdf5a53Sjeremylt   @ref Utility
753b11c1e72Sjeremylt **/
754d7b241e6Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) {
755a8de75f0Sjeremylt   *P = basis->P;
756d7b241e6Sjeremylt   return 0;
757d7b241e6Sjeremylt }
758d7b241e6Sjeremylt 
759b11c1e72Sjeremylt /**
7604ce2993fSjeremylt   @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis
761b11c1e72Sjeremylt 
762b11c1e72Sjeremylt   @param basis   CeedBasis
7634ce2993fSjeremylt   @param[out] Q  Variable to store number of quadrature points
764b11c1e72Sjeremylt 
765b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
766dfdf5a53Sjeremylt 
767dfdf5a53Sjeremylt   @ref Utility
768b11c1e72Sjeremylt **/
769d7b241e6Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) {
770a8de75f0Sjeremylt   *Q = basis->Q;
771d7b241e6Sjeremylt   return 0;
772d7b241e6Sjeremylt }
773d7b241e6Sjeremylt 
774b11c1e72Sjeremylt /**
7758c91a0c9SJeremy L Thompson   @brief Get reference coordinates of quadrature points (in dim dimensions)
7764ce2993fSjeremylt          of a CeedBasis
7774ce2993fSjeremylt 
7784ce2993fSjeremylt   @param basis      CeedBasis
7798c91a0c9SJeremy L Thompson   @param[out] qref  Variable to store reference coordinates of quadrature points
7804ce2993fSjeremylt 
7814ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
7824ce2993fSjeremylt 
78323617272Sjeremylt   @ref Advanced
7844ce2993fSjeremylt **/
7854ce2993fSjeremylt int CeedBasisGetQRef(CeedBasis basis, CeedScalar* *qref) {
7864ce2993fSjeremylt   *qref = basis->qref1d;
7874ce2993fSjeremylt   return 0;
7884ce2993fSjeremylt }
7894ce2993fSjeremylt 
7904ce2993fSjeremylt /**
7914ce2993fSjeremylt   @brief Get quadrature weights of quadrature points (in dim dimensions)
7924ce2993fSjeremylt          of a CeedBasis
7934ce2993fSjeremylt 
7944ce2993fSjeremylt   @param basis         CeedBasis
7954ce2993fSjeremylt   @param[out] qweight  Variable to store quadrature weights
7964ce2993fSjeremylt 
7974ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
7984ce2993fSjeremylt 
79923617272Sjeremylt   @ref Advanced
8004ce2993fSjeremylt **/
8014ce2993fSjeremylt int CeedBasisGetQWeights(CeedBasis basis, CeedScalar* *qweight) {
8024ce2993fSjeremylt   *qweight = basis->qweight1d;
8034ce2993fSjeremylt   return 0;
8044ce2993fSjeremylt }
8054ce2993fSjeremylt 
8064ce2993fSjeremylt /**
8074ce2993fSjeremylt   @brief Get interpolation matrix of a CeedBasis
8084ce2993fSjeremylt 
8094ce2993fSjeremylt   @param basis      CeedBasis
8104ce2993fSjeremylt   @param[out] qref  Variable to store interpolation matrix
8114ce2993fSjeremylt 
8124ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
8134ce2993fSjeremylt 
81423617272Sjeremylt   @ref Advanced
8154ce2993fSjeremylt **/
8164ce2993fSjeremylt int CeedBasisGetInterp(CeedBasis basis, CeedScalar* *interp) {
8174ce2993fSjeremylt   *interp = basis->interp1d;
8184ce2993fSjeremylt   return 0;
8194ce2993fSjeremylt }
8204ce2993fSjeremylt 
8214ce2993fSjeremylt /**
8224ce2993fSjeremylt   @brief Get gradient matrix of a CeedBasis
8234ce2993fSjeremylt 
8244ce2993fSjeremylt   @param basis      CeedBasis
8254ce2993fSjeremylt   @param[out] qref  Variable to store gradient matrix
8264ce2993fSjeremylt 
8274ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
8284ce2993fSjeremylt 
82923617272Sjeremylt   @ref Advanced
8304ce2993fSjeremylt **/
8314ce2993fSjeremylt int CeedBasisGetGrad(CeedBasis basis, CeedScalar* *grad) {
8324ce2993fSjeremylt   *grad = basis->grad1d;
8334ce2993fSjeremylt   return 0;
8344ce2993fSjeremylt }
8354ce2993fSjeremylt 
8364ce2993fSjeremylt /**
8374ce2993fSjeremylt   @brief Get backend data of a CeedBasis
8384ce2993fSjeremylt 
8394ce2993fSjeremylt   @param basis      CeedBasis
8404ce2993fSjeremylt   @param[out] data  Variable to store data
8414ce2993fSjeremylt 
8424ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
8434ce2993fSjeremylt 
84423617272Sjeremylt   @ref Advanced
8454ce2993fSjeremylt **/
8464ce2993fSjeremylt int CeedBasisGetData(CeedBasis basis, void* *data) {
8474ce2993fSjeremylt   *data = basis->data;
8484ce2993fSjeremylt   return 0;
8494ce2993fSjeremylt }
8504ce2993fSjeremylt 
8514ce2993fSjeremylt /**
852fe2413ffSjeremylt   @brief Set backend data of a CeedBasis
853fe2413ffSjeremylt 
854fe2413ffSjeremylt   @param[out] basis CeedBasis
855fe2413ffSjeremylt   @param data       Data to set
856fe2413ffSjeremylt 
857fe2413ffSjeremylt   @return An error code: 0 - success, otherwise - failure
858fe2413ffSjeremylt 
859fe2413ffSjeremylt   @ref Advanced
860fe2413ffSjeremylt **/
861fe2413ffSjeremylt int CeedBasisSetData(CeedBasis basis, void* *data) {
862fe2413ffSjeremylt   basis->data = *data;
863fe2413ffSjeremylt   return 0;
864fe2413ffSjeremylt }
865fe2413ffSjeremylt 
866fe2413ffSjeremylt /**
8672f86a920SJeremy L Thompson   @brief Get CeedTensorContract of a CeedBasis
8682f86a920SJeremy L Thompson 
8692f86a920SJeremy L Thompson   @param basis          CeedBasis
8702f86a920SJeremy L Thompson   @param[out] contract  Variable to store CeedTensorContract
8712f86a920SJeremy L Thompson 
8722f86a920SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
8732f86a920SJeremy L Thompson 
8742f86a920SJeremy L Thompson   @ref Advanced
8752f86a920SJeremy L Thompson **/
8762f86a920SJeremy L Thompson int CeedBasisGetTensorContract(CeedBasis basis,
8772f86a920SJeremy L Thompson                                CeedTensorContract *contract) {
8782f86a920SJeremy L Thompson   *contract = basis->contract;
8792f86a920SJeremy L Thompson   return 0;
8802f86a920SJeremy L Thompson }
8812f86a920SJeremy L Thompson 
8822f86a920SJeremy L Thompson /**
8832f86a920SJeremy L Thompson   @brief Set CeedTensorContract of a CeedBasis
8842f86a920SJeremy L Thompson 
8852f86a920SJeremy L Thompson   @param[out] basis     CeedBasis
8862f86a920SJeremy L Thompson   @param contract       CeedTensorContract to set
8872f86a920SJeremy L Thompson 
8882f86a920SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
8892f86a920SJeremy L Thompson 
8902f86a920SJeremy L Thompson   @ref Advanced
8912f86a920SJeremy L Thompson **/
8922f86a920SJeremy L Thompson int CeedBasisSetTensorContract(CeedBasis basis,
8932f86a920SJeremy L Thompson                                CeedTensorContract *contract) {
8942f86a920SJeremy L Thompson   basis->contract = *contract;
8952f86a920SJeremy L Thompson   return 0;
8962f86a920SJeremy L Thompson }
8972f86a920SJeremy L Thompson 
8982f86a920SJeremy L Thompson /**
899a8de75f0Sjeremylt   @brief Get dimension for given CeedElemTopology
900a8de75f0Sjeremylt 
901a8de75f0Sjeremylt   @param topo      CeedElemTopology
9024ce2993fSjeremylt   @param[out] dim  Variable to store dimension of topology
903a8de75f0Sjeremylt 
904a8de75f0Sjeremylt   @return An error code: 0 - success, otherwise - failure
905a8de75f0Sjeremylt 
90623617272Sjeremylt   @ref Advanced
907a8de75f0Sjeremylt **/
908a8de75f0Sjeremylt int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) {
909a8de75f0Sjeremylt   *dim = (CeedInt) topo >> 16;
910a8de75f0Sjeremylt 
911a8de75f0Sjeremylt   return 0;
912a8de75f0Sjeremylt };
913a8de75f0Sjeremylt 
914a8de75f0Sjeremylt /**
915b11c1e72Sjeremylt   @brief Destroy a CeedBasis
916b11c1e72Sjeremylt 
917b11c1e72Sjeremylt   @param basis CeedBasis to destroy
918b11c1e72Sjeremylt 
919b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
920dfdf5a53Sjeremylt 
921dfdf5a53Sjeremylt   @ref Basic
922b11c1e72Sjeremylt **/
923d7b241e6Sjeremylt int CeedBasisDestroy(CeedBasis *basis) {
924d7b241e6Sjeremylt   int ierr;
925d7b241e6Sjeremylt 
926d7b241e6Sjeremylt   if (!*basis || --(*basis)->refcount > 0) return 0;
927d7b241e6Sjeremylt   if ((*basis)->Destroy) {
928d7b241e6Sjeremylt     ierr = (*basis)->Destroy(*basis); CeedChk(ierr);
929d7b241e6Sjeremylt   }
930d7b241e6Sjeremylt   ierr = CeedFree(&(*basis)->interp1d); CeedChk(ierr);
931d7b241e6Sjeremylt   ierr = CeedFree(&(*basis)->grad1d); CeedChk(ierr);
932d7b241e6Sjeremylt   ierr = CeedFree(&(*basis)->qref1d); CeedChk(ierr);
933d7b241e6Sjeremylt   ierr = CeedFree(&(*basis)->qweight1d); CeedChk(ierr);
934d7b241e6Sjeremylt   ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr);
935d7b241e6Sjeremylt   ierr = CeedFree(basis); CeedChk(ierr);
936d7b241e6Sjeremylt   return 0;
937d7b241e6Sjeremylt }
938d7b241e6Sjeremylt 
93933e6becaSjeremylt /// @cond DOXYGEN_SKIP
9408795c945Sjeremylt // Indicate that the quadrature points are collocated with the nodes
941783c99b3SValeria Barra CeedBasis CEED_BASIS_COLLOCATED = &ceed_basis_collocated;
94233e6becaSjeremylt /// @endcond
943d7b241e6Sjeremylt /// @}
944