xref: /libCEED/rust/libceed-sys/c-src/interface/ceed-basis.c (revision 52bfb9bbf17f17edbcd45876cdc8689a879bc683)
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 
634d537eeaSYohann   if (dim<1)
644d537eeaSYohann     return CeedError(ceed, 1, "Basis dimension must be a positive value");
654d537eeaSYohann 
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;
1264d537eeaSYohann 
1274d537eeaSYohann   if (dim<1)
1284d537eeaSYohann     return CeedError(ceed, 1, "Basis dimension must be a positive value");
1294d537eeaSYohann 
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 /**
436*52bfb9bbSJeremy L Thompson   @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 
441*52bfb9bbSJeremy L Thompson   @param[in,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
446*52bfb9bbSJeremy L Thompson   @param row        Row stride
447*52bfb9bbSJeremy L Thompson   @param col        Col 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 
468*52bfb9bbSJeremy L Thompson     Compute A = Q A where Q is mxm and A is mxn.
469dfdf5a53Sjeremylt 
470*52bfb9bbSJeremy L Thompson   @param[in,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
476*52bfb9bbSJeremy L Thompson   @param k          Number of elementary reflectors in Q, k<m
477*52bfb9bbSJeremy L Thompson   @param row        Row stride in A
478*52bfb9bbSJeremy L Thompson   @param col        Col stride in A
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;
491*52bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<m; j++)
492d7b241e6Sjeremylt       v[j] = Q[j*k+i];
493d7b241e6Sjeremylt     // Apply Householder reflector (I - tau v v^T) colograd1d^T
494d7b241e6Sjeremylt     CeedHouseholderReflect(&A[i*row], &v[i], tau[i], m-i, n, row, col);
495d7b241e6Sjeremylt   }
496d7b241e6Sjeremylt   return 0;
497d7b241e6Sjeremylt }
498d7b241e6Sjeremylt 
499b11c1e72Sjeremylt /**
500*52bfb9bbSJeremy L Thompson   @brief Compute Givens rotation
501*52bfb9bbSJeremy L Thompson 
502*52bfb9bbSJeremy L Thompson     Computes A = G A (or G^T A in transpose mode)
503*52bfb9bbSJeremy L Thompson     where A is an mxn matrix indexed as A[i*n + j*m]
504*52bfb9bbSJeremy L Thompson 
505*52bfb9bbSJeremy L Thompson   @param[in,out] A  Row major matrix to apply Givens rotation to, in place
506*52bfb9bbSJeremy L Thompson   @param c          Cosine factor
507*52bfb9bbSJeremy L Thompson   @param s          Sine factor
508*52bfb9bbSJeremy L Thompson   @param i          First row/column to apply rotation
509*52bfb9bbSJeremy L Thompson   @param k          Second row/column to apply rotation
510*52bfb9bbSJeremy L Thompson   @param m          Number of rows in A
511*52bfb9bbSJeremy L Thompson   @param n          Number of columns in A
512*52bfb9bbSJeremy L Thompson 
513*52bfb9bbSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
514*52bfb9bbSJeremy L Thompson 
515*52bfb9bbSJeremy L Thompson   @ref Developer
516*52bfb9bbSJeremy L Thompson **/
517*52bfb9bbSJeremy L Thompson static int CeedGivensRotation(CeedScalar *A, CeedScalar c, CeedScalar s,
518*52bfb9bbSJeremy L Thompson                               CeedTransposeMode tmode, CeedInt i, CeedInt k,
519*52bfb9bbSJeremy L Thompson                               CeedInt m, CeedInt n) {
520*52bfb9bbSJeremy L Thompson   CeedInt stridej = 1, strideik = m, numits = n;
521*52bfb9bbSJeremy L Thompson   if (tmode == CEED_NOTRANSPOSE) {
522*52bfb9bbSJeremy L Thompson     stridej = n; strideik = 1; numits = m;
523*52bfb9bbSJeremy L Thompson   }
524*52bfb9bbSJeremy L Thompson 
525*52bfb9bbSJeremy L Thompson   // Apply rotation
526*52bfb9bbSJeremy L Thompson   for (CeedInt j=0; j<numits; j++) {
527*52bfb9bbSJeremy L Thompson     CeedScalar tau1 = A[i*strideik+j*stridej], tau2 = A[k*strideik+j*stridej];
528*52bfb9bbSJeremy L Thompson     A[i*strideik+j*stridej] = c*tau1 - s*tau2;
529*52bfb9bbSJeremy L Thompson     A[k*strideik+j*stridej] = s*tau1 + c*tau2;
530*52bfb9bbSJeremy L Thompson   }
531*52bfb9bbSJeremy L Thompson 
532*52bfb9bbSJeremy L Thompson   return 0;
533*52bfb9bbSJeremy L Thompson }
534*52bfb9bbSJeremy L Thompson 
535*52bfb9bbSJeremy L Thompson /**
536b11c1e72Sjeremylt   @brief Return QR Factorization of matrix
537b11c1e72Sjeremylt 
538288c0443SJeremy L Thompson   @param ceed         A Ceed object currently in use
539*52bfb9bbSJeremy L Thompson   @param[in,out] mat  Row-major matrix to be factorized in place
540*52bfb9bbSJeremy L Thompson   @param[in,out] tau  Vector of length m of scaling factors
541b11c1e72Sjeremylt   @param m            Number of rows
542b11c1e72Sjeremylt   @param n            Number of columns
543b11c1e72Sjeremylt 
544b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
545dfdf5a53Sjeremylt 
546dfdf5a53Sjeremylt   @ref Utility
547b11c1e72Sjeremylt **/
548a7bd39daSjeremylt int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau,
549d7b241e6Sjeremylt                         CeedInt m, CeedInt n) {
550d7b241e6Sjeremylt   CeedScalar v[m];
551d7b241e6Sjeremylt 
552a7bd39daSjeremylt   // Check m >= n
553a7bd39daSjeremylt   if (n > m)
554a7bd39daSjeremylt     return CeedError(ceed, 1, "Cannot compute QR factorization with n > m");
555a7bd39daSjeremylt 
556*52bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++) {
557d7b241e6Sjeremylt     // Calculate Householder vector, magnitude
558d7b241e6Sjeremylt     CeedScalar sigma = 0.0;
559d7b241e6Sjeremylt     v[i] = mat[i+n*i];
560*52bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<m; j++) {
561d7b241e6Sjeremylt       v[j] = mat[i+n*j];
562d7b241e6Sjeremylt       sigma += v[j] * v[j];
563d7b241e6Sjeremylt     }
564d7b241e6Sjeremylt     CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m]
565d7b241e6Sjeremylt     CeedScalar Rii = -copysign(norm, v[i]);
566d7b241e6Sjeremylt     v[i] -= Rii;
567d7b241e6Sjeremylt     // norm of v[i:m] after modification above and scaling below
568d7b241e6Sjeremylt     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
569d7b241e6Sjeremylt     //   tau = 2 / (norm*norm)
570d7b241e6Sjeremylt     tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma);
571*52bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<m; j++) v[j] /= v[i];
572d7b241e6Sjeremylt 
573d7b241e6Sjeremylt     // Apply Householder reflector to lower right panel
574d7b241e6Sjeremylt     CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1);
575d7b241e6Sjeremylt     // Save v
576d7b241e6Sjeremylt     mat[i+n*i] = Rii;
577*52bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<m; j++) {
578d7b241e6Sjeremylt       mat[i+n*j] = v[j];
579d7b241e6Sjeremylt     }
580d7b241e6Sjeremylt   }
581d7b241e6Sjeremylt 
582d7b241e6Sjeremylt   return 0;
583d7b241e6Sjeremylt }
584d7b241e6Sjeremylt 
585b11c1e72Sjeremylt /**
586*52bfb9bbSJeremy L Thompson   @brief Return symmetric Schur decomposition of the symmetric matrix mat via
587*52bfb9bbSJeremy L Thompson            symmetric QR factorization
588*52bfb9bbSJeremy L Thompson 
589*52bfb9bbSJeremy L Thompson   @param[in,out] mat  Row-major matrix to be factorized in place
590*52bfb9bbSJeremy L Thompson   @param[out] lambda  Vector of length m of eigenvalues
591*52bfb9bbSJeremy L Thompson   @param n            Number of rows/columns
592*52bfb9bbSJeremy L Thompson 
593*52bfb9bbSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
594*52bfb9bbSJeremy L Thompson 
595*52bfb9bbSJeremy L Thompson   @ref Utility
596*52bfb9bbSJeremy L Thompson **/
597*52bfb9bbSJeremy L Thompson int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat,
598*52bfb9bbSJeremy L Thompson                                     CeedScalar *lambda, CeedInt n) {
599*52bfb9bbSJeremy L Thompson   // Check bounds for clang-tidy
600*52bfb9bbSJeremy L Thompson   if (n<2)
601*52bfb9bbSJeremy L Thompson     return CeedError(ceed, 1, "Cannot compute symmetric Schur decomposition of scalars");
602*52bfb9bbSJeremy L Thompson 
603*52bfb9bbSJeremy L Thompson   CeedScalar v[n-1], tau[n-1], matT[n*n];
604*52bfb9bbSJeremy L Thompson 
605*52bfb9bbSJeremy L Thompson   // Copy mat to matT and set mat to I
606*52bfb9bbSJeremy L Thompson   memcpy(matT, mat, n*n*sizeof(mat[0]));
607*52bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
608*52bfb9bbSJeremy L Thompson     for (CeedInt j=0; j<n; j++)
609*52bfb9bbSJeremy L Thompson       mat[j+n*i] = (i==j) ? 1 : 0;
610*52bfb9bbSJeremy L Thompson 
611*52bfb9bbSJeremy L Thompson   // Reduce to tridiagonal
612*52bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n-1; i++) {
613*52bfb9bbSJeremy L Thompson     // Calculate Householder vector, magnitude
614*52bfb9bbSJeremy L Thompson     CeedScalar sigma = 0.0;
615*52bfb9bbSJeremy L Thompson     v[i] = matT[i+n*(i+1)];
616*52bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) {
617*52bfb9bbSJeremy L Thompson       v[j] = matT[i+n*(j+1)];
618*52bfb9bbSJeremy L Thompson       sigma += v[j] * v[j];
619*52bfb9bbSJeremy L Thompson     }
620*52bfb9bbSJeremy L Thompson     CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:n-1]
621*52bfb9bbSJeremy L Thompson     CeedScalar Rii = -copysign(norm, v[i]);
622*52bfb9bbSJeremy L Thompson     v[i] -= Rii;
623*52bfb9bbSJeremy L Thompson     // norm of v[i:m] after modification above and scaling below
624*52bfb9bbSJeremy L Thompson     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
625*52bfb9bbSJeremy L Thompson     //   tau = 2 / (norm*norm)
626*52bfb9bbSJeremy L Thompson     tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma);
627*52bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) v[j] /= v[i];
628*52bfb9bbSJeremy L Thompson 
629*52bfb9bbSJeremy L Thompson     // Update sub and super diagonal
630*52bfb9bbSJeremy L Thompson     matT[i+n*(i+1)] = Rii;
631*52bfb9bbSJeremy L Thompson     matT[(i+1)+n*i] = Rii;
632*52bfb9bbSJeremy L Thompson     for (CeedInt j=i+2; j<n; j++) {
633*52bfb9bbSJeremy L Thompson       matT[i+n*j] = 0; matT[j+n*i] = 0;
634*52bfb9bbSJeremy L Thompson     }
635*52bfb9bbSJeremy L Thompson     // Apply symmetric Householder reflector to lower right panel
636*52bfb9bbSJeremy L Thompson     CeedHouseholderReflect(&matT[(i+1)+n*(i+1)], &v[i], tau[i],
637*52bfb9bbSJeremy L Thompson                            n-(i+1), n-(i+1), n, 1);
638*52bfb9bbSJeremy L Thompson     CeedHouseholderReflect(&matT[(i+1)+n*(i+1)], &v[i], tau[i],
639*52bfb9bbSJeremy L Thompson                            n-(i+1), n-(i+1), 1, n);
640*52bfb9bbSJeremy L Thompson     // Save v
641*52bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) {
642*52bfb9bbSJeremy L Thompson       matT[i+n*(j+1)] = v[j];
643*52bfb9bbSJeremy L Thompson     }
644*52bfb9bbSJeremy L Thompson   }
645*52bfb9bbSJeremy L Thompson   // Backwards accumulation of Q
646*52bfb9bbSJeremy L Thompson   for (CeedInt i=n-2; i>=0; i--) {
647*52bfb9bbSJeremy L Thompson     v[i] = 1;
648*52bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) {
649*52bfb9bbSJeremy L Thompson       v[j] = matT[i+n*(j+1)];
650*52bfb9bbSJeremy L Thompson       matT[i+n*(j+1)] = 0;
651*52bfb9bbSJeremy L Thompson     }
652*52bfb9bbSJeremy L Thompson     CeedHouseholderReflect(&mat[(i+1)+n*(i+1)], &v[i], tau[i],
653*52bfb9bbSJeremy L Thompson                            n-(i+1), n-(i+1), n, 1);
654*52bfb9bbSJeremy L Thompson   }
655*52bfb9bbSJeremy L Thompson 
656*52bfb9bbSJeremy L Thompson   // Reduce sub and super diagonal
657*52bfb9bbSJeremy L Thompson   CeedInt p = 0, q = 0, itr = 0, maxitr = n*n*n;
658*52bfb9bbSJeremy L Thompson   CeedScalar tol = 1e-15;
659*52bfb9bbSJeremy L Thompson 
660*52bfb9bbSJeremy L Thompson   while (q < n && itr < maxitr) {
661*52bfb9bbSJeremy L Thompson     // Update p, q, size of reduced portions of diagonal
662*52bfb9bbSJeremy L Thompson     p = 0; q = 0;
663*52bfb9bbSJeremy L Thompson     for (CeedInt i=n-2; i>=0; i--) {
664*52bfb9bbSJeremy L Thompson       if (fabs(matT[i+n*(i+1)]) < tol)
665*52bfb9bbSJeremy L Thompson         q += 1;
666*52bfb9bbSJeremy L Thompson       else
667*52bfb9bbSJeremy L Thompson          break;
668*52bfb9bbSJeremy L Thompson     }
669*52bfb9bbSJeremy L Thompson     for (CeedInt i=0; i<n-1-q; i++) {
670*52bfb9bbSJeremy L Thompson       if (fabs(matT[i+n*(i+1)]) < tol)
671*52bfb9bbSJeremy L Thompson         p += 1;
672*52bfb9bbSJeremy L Thompson       else
673*52bfb9bbSJeremy L Thompson         break;
674*52bfb9bbSJeremy L Thompson     }
675*52bfb9bbSJeremy L Thompson     if (q == n-1) break; // Finished reducing
676*52bfb9bbSJeremy L Thompson 
677*52bfb9bbSJeremy L Thompson     // Reduce tridiagonal portion
678*52bfb9bbSJeremy L Thompson     CeedScalar tnn = matT[(n-1-q)+n*(n-1-q)],
679*52bfb9bbSJeremy L Thompson                tnnm1 = matT[(n-2-q)+n*(n-1-q)];
680*52bfb9bbSJeremy L Thompson     CeedScalar d = (matT[(n-2-q)+n*(n-2-q)] - tnn)/2;
681*52bfb9bbSJeremy L Thompson     CeedScalar mu = tnn - tnnm1*tnnm1 /
682*52bfb9bbSJeremy L Thompson                       (d + copysign(sqrt(d*d + tnnm1*tnnm1), d));
683*52bfb9bbSJeremy L Thompson     CeedScalar x = matT[p+n*p] - mu;
684*52bfb9bbSJeremy L Thompson     CeedScalar z = matT[p+n*(p+1)];
685*52bfb9bbSJeremy L Thompson     for (CeedInt k=p; k<n-1-q; k++) {
686*52bfb9bbSJeremy L Thompson       // Compute Givens rotation
687*52bfb9bbSJeremy L Thompson       CeedScalar c = 1, s = 0;
688*52bfb9bbSJeremy L Thompson       if (fabs(z) > tol) {
689*52bfb9bbSJeremy L Thompson         if (fabs(z) > fabs(x)) {
690*52bfb9bbSJeremy L Thompson           CeedScalar tau = -x/z;
691*52bfb9bbSJeremy L Thompson           s = 1/sqrt(1+tau*tau), c = s*tau;
692*52bfb9bbSJeremy L Thompson         } else {
693*52bfb9bbSJeremy L Thompson           CeedScalar tau = -z/x;
694*52bfb9bbSJeremy L Thompson           c = 1/sqrt(1+tau*tau), s = c*tau;
695*52bfb9bbSJeremy L Thompson         }
696*52bfb9bbSJeremy L Thompson       }
697*52bfb9bbSJeremy L Thompson 
698*52bfb9bbSJeremy L Thompson       // Apply Givens rotation to T
699*52bfb9bbSJeremy L Thompson       CeedGivensRotation(matT, c, s, CEED_NOTRANSPOSE, k, k+1, n, n);
700*52bfb9bbSJeremy L Thompson       CeedGivensRotation(matT, c, s, CEED_TRANSPOSE, k, k+1, n, n);
701*52bfb9bbSJeremy L Thompson 
702*52bfb9bbSJeremy L Thompson       // Apply Givens rotation to Q
703*52bfb9bbSJeremy L Thompson       CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k+1, n, n);
704*52bfb9bbSJeremy L Thompson 
705*52bfb9bbSJeremy L Thompson       // Update x, z
706*52bfb9bbSJeremy L Thompson       if (k < n-q-2) {
707*52bfb9bbSJeremy L Thompson         x = matT[k+n*(k+1)];
708*52bfb9bbSJeremy L Thompson         z = matT[k+n*(k+2)];
709*52bfb9bbSJeremy L Thompson       }
710*52bfb9bbSJeremy L Thompson     }
711*52bfb9bbSJeremy L Thompson     itr++;
712*52bfb9bbSJeremy L Thompson   }
713*52bfb9bbSJeremy L Thompson   // Save eigenvalues
714*52bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
715*52bfb9bbSJeremy L Thompson     lambda[i] = matT[i+n*i];
716*52bfb9bbSJeremy L Thompson 
717*52bfb9bbSJeremy L Thompson   // Check convergence
718*52bfb9bbSJeremy L Thompson   if (itr == maxitr && q < n-1)
719*52bfb9bbSJeremy L Thompson     return CeedError(ceed, 1, "Symmetric QR failed to converge");
720*52bfb9bbSJeremy L Thompson 
721*52bfb9bbSJeremy L Thompson   return 0;
722*52bfb9bbSJeremy L Thompson }
723*52bfb9bbSJeremy L Thompson 
724*52bfb9bbSJeremy L Thompson /**
725*52bfb9bbSJeremy L Thompson   @brief Return C = A B
726*52bfb9bbSJeremy L Thompson 
727*52bfb9bbSJeremy L Thompson   @param[in] matA     Row-major matrix A
728*52bfb9bbSJeremy L Thompson   @param[in] matB     Row-major matrix B
729*52bfb9bbSJeremy L Thompson   @param[out] matC    Row-major output matrix C
730*52bfb9bbSJeremy L Thompson   @param m            Number of rows of C
731*52bfb9bbSJeremy L Thompson   @param n            Number of columns of C
732*52bfb9bbSJeremy L Thompson   @param kk           Number of columns of A/rows of B
733*52bfb9bbSJeremy L Thompson 
734*52bfb9bbSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
735*52bfb9bbSJeremy L Thompson 
736*52bfb9bbSJeremy L Thompson   @ref Utility
737*52bfb9bbSJeremy L Thompson **/
738*52bfb9bbSJeremy L Thompson static int CeedMatrixMultiply(Ceed ceed, CeedScalar *matA, CeedScalar *matB,
739*52bfb9bbSJeremy L Thompson                               CeedScalar *matC, CeedInt m, CeedInt n,
740*52bfb9bbSJeremy L Thompson                               CeedInt kk) {
741*52bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<m; i++)
742*52bfb9bbSJeremy L Thompson     for (CeedInt j=0; j<n; j++) {
743*52bfb9bbSJeremy L Thompson       CeedScalar sum = 0;
744*52bfb9bbSJeremy L Thompson       for (CeedInt k=0; k<kk; k++)
745*52bfb9bbSJeremy L Thompson         sum += matA[k+i*kk]*matB[j+k*n];
746*52bfb9bbSJeremy L Thompson       matC[j+i*n] = sum;
747*52bfb9bbSJeremy L Thompson     }
748*52bfb9bbSJeremy L Thompson   return 0;
749*52bfb9bbSJeremy L Thompson }
750*52bfb9bbSJeremy L Thompson 
751*52bfb9bbSJeremy L Thompson /**
752*52bfb9bbSJeremy L Thompson   @brief Return Simultaneous Diagonalization of two matrices. This solves the
753*52bfb9bbSJeremy L Thompson            generalized eigenvalue problem A x = lambda B x, where A and B
754*52bfb9bbSJeremy L Thompson            are symmetric and B is positive definite. We generate the matrix X
755*52bfb9bbSJeremy L Thompson            and vector Lambda such that X^T A X = Lambda and X^T B X = I. This
756*52bfb9bbSJeremy L Thompson            is equivalent to the LAPACK routine 'sygv' with TYPE = 1.
757*52bfb9bbSJeremy L Thompson 
758*52bfb9bbSJeremy L Thompson   @param[in] matA     Row-major matrix to be factorized with eigenvalues
759*52bfb9bbSJeremy L Thompson   @param[in] matB     Row-major matrix to be factorized to identity
760*52bfb9bbSJeremy L Thompson   @param[out] x       Row-major orthogonal matrix
761*52bfb9bbSJeremy L Thompson   @param[out] lambda  Vector of length m of generalized eigenvalues
762*52bfb9bbSJeremy L Thompson   @param n            Number of rows/columns
763*52bfb9bbSJeremy L Thompson 
764*52bfb9bbSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
765*52bfb9bbSJeremy L Thompson 
766*52bfb9bbSJeremy L Thompson   @ref Utility
767*52bfb9bbSJeremy L Thompson **/
768*52bfb9bbSJeremy L Thompson int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *matA,
769*52bfb9bbSJeremy L Thompson                                     CeedScalar *matB, CeedScalar *x,
770*52bfb9bbSJeremy L Thompson                                     CeedScalar *lambda, CeedInt n) {
771*52bfb9bbSJeremy L Thompson   int ierr;
772*52bfb9bbSJeremy L Thompson   CeedScalar matC[n*n], matG[n*n], vecD[n];
773*52bfb9bbSJeremy L Thompson 
774*52bfb9bbSJeremy L Thompson   // Compute B = G D G^T
775*52bfb9bbSJeremy L Thompson   memcpy(matG, matB, n*n*sizeof(matB[0]));
776*52bfb9bbSJeremy L Thompson   ierr = CeedSymmetricSchurDecomposition(ceed, matG, vecD, n); CeedChk(ierr);
777*52bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++) vecD[i] = sqrt(vecD[i]);
778*52bfb9bbSJeremy L Thompson 
779*52bfb9bbSJeremy L Thompson   // Compute C = (G D^-1/2)^-1 A (G D^-1/2)^-T
780*52bfb9bbSJeremy L Thompson   //           = D^1/2 G^T A D^1/2 G
781*52bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
782*52bfb9bbSJeremy L Thompson     for (CeedInt j=0; j<n; j++)
783*52bfb9bbSJeremy L Thompson       matC[j+i*n] = vecD[i] * matG[i+j*n];
784*52bfb9bbSJeremy L Thompson   CeedMatrixMultiply(ceed, matC, matA, x, n, n, n);
785*52bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
786*52bfb9bbSJeremy L Thompson     for (CeedInt j=0; j<n; j++)
787*52bfb9bbSJeremy L Thompson       matG[j+i*n] = vecD[i] * matG[j+i*n];
788*52bfb9bbSJeremy L Thompson   CeedMatrixMultiply(ceed, x, matG, matC, n, n, n);
789*52bfb9bbSJeremy L Thompson 
790*52bfb9bbSJeremy L Thompson   // Compute Q^T C Q = lambda
791*52bfb9bbSJeremy L Thompson   ierr = CeedSymmetricSchurDecomposition(ceed, matC, lambda, n); CeedChk(ierr);
792*52bfb9bbSJeremy L Thompson 
793*52bfb9bbSJeremy L Thompson   // Set x = (G D^-1/2)^-T Q
794*52bfb9bbSJeremy L Thompson   //       = D^1/2 G Q
795*52bfb9bbSJeremy L Thompson   CeedMatrixMultiply(ceed, matG, matC, x, n, n, n);
796*52bfb9bbSJeremy L Thompson 
797*52bfb9bbSJeremy L Thompson   return 0;
798*52bfb9bbSJeremy L Thompson }
799*52bfb9bbSJeremy L Thompson 
800*52bfb9bbSJeremy L Thompson /**
801783c99b3SValeria Barra   @brief Return collocated grad matrix
802b11c1e72Sjeremylt 
803b11c1e72Sjeremylt   @param basis           CeedBasis
804b11c1e72Sjeremylt   @param[out] colograd1d Row-major Q1d × Q1d matrix expressing derivatives of
805b11c1e72Sjeremylt                            basis functions at quadrature points
806b11c1e72Sjeremylt 
807b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
808dfdf5a53Sjeremylt 
809dfdf5a53Sjeremylt   @ref Advanced
810b11c1e72Sjeremylt **/
811783c99b3SValeria Barra int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *colograd1d) {
812d7b241e6Sjeremylt   int i, j, k;
813a7bd39daSjeremylt   Ceed ceed;
814d7b241e6Sjeremylt   CeedInt ierr, P1d=(basis)->P1d, Q1d=(basis)->Q1d;
815d7b241e6Sjeremylt   CeedScalar *interp1d, *grad1d, tau[Q1d];
816d7b241e6Sjeremylt 
817d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d*P1d, &interp1d); CeedChk(ierr);
818d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d*P1d, &grad1d); CeedChk(ierr);
819d7b241e6Sjeremylt   memcpy(interp1d, (basis)->interp1d, Q1d*P1d*sizeof(basis)->interp1d[0]);
820d7b241e6Sjeremylt   memcpy(grad1d, (basis)->grad1d, Q1d*P1d*sizeof(basis)->interp1d[0]);
821d7b241e6Sjeremylt 
822d7b241e6Sjeremylt   // QR Factorization, interp1d = Q R
823a7bd39daSjeremylt   ierr = CeedBasisGetCeed(basis, &ceed); CeedChk(ierr);
824a7bd39daSjeremylt   ierr = CeedQRFactorization(ceed, interp1d, tau, Q1d, P1d); CeedChk(ierr);
825d7b241e6Sjeremylt 
826d7b241e6Sjeremylt   // Apply Rinv, colograd1d = grad1d Rinv
827d7b241e6Sjeremylt   for (i=0; i<Q1d; i++) { // Row i
828d7b241e6Sjeremylt     colograd1d[Q1d*i] = grad1d[P1d*i]/interp1d[0];
829d7b241e6Sjeremylt     for (j=1; j<P1d; j++) { // Column j
830d7b241e6Sjeremylt       colograd1d[j+Q1d*i] = grad1d[j+P1d*i];
831d7b241e6Sjeremylt       for (k=0; k<j; k++) {
832d7b241e6Sjeremylt         colograd1d[j+Q1d*i] -= interp1d[j+P1d*k]*colograd1d[k+Q1d*i];
833d7b241e6Sjeremylt       }
834d7b241e6Sjeremylt       colograd1d[j+Q1d*i] /= interp1d[j+P1d*j];
835d7b241e6Sjeremylt     }
836d7b241e6Sjeremylt     for (j=P1d; j<Q1d; j++) {
837d7b241e6Sjeremylt       colograd1d[j+Q1d*i] = 0;
838d7b241e6Sjeremylt     }
839d7b241e6Sjeremylt   }
840d7b241e6Sjeremylt 
841d7b241e6Sjeremylt   // Apply Qtranspose, colograd = colograd Qtranspose
842d7b241e6Sjeremylt   CeedHouseholderApplyQ(colograd1d, interp1d, tau, CEED_NOTRANSPOSE,
843d7b241e6Sjeremylt                         Q1d, Q1d, P1d, 1, Q1d);
844d7b241e6Sjeremylt 
845d7b241e6Sjeremylt   ierr = CeedFree(&interp1d); CeedChk(ierr);
846d7b241e6Sjeremylt   ierr = CeedFree(&grad1d); CeedChk(ierr);
847d7b241e6Sjeremylt 
848d7b241e6Sjeremylt   return 0;
849d7b241e6Sjeremylt }
850d7b241e6Sjeremylt 
851b11c1e72Sjeremylt /**
852b11c1e72Sjeremylt   @brief Apply basis evaluation from nodes to quadrature points or vice-versa
853b11c1e72Sjeremylt 
854b11c1e72Sjeremylt   @param basis  CeedBasis to evaluate
855b11c1e72Sjeremylt   @param nelem  The number of elements to apply the basis evaluation to;
856b11c1e72Sjeremylt                   the backend will specify the ordering in
857b11c1e72Sjeremylt                   ElemRestrictionCreateBlocked
858b11c1e72Sjeremylt   @param tmode  \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature
859b11c1e72Sjeremylt                   points, \ref CEED_TRANSPOSE to apply the transpose, mapping
860b11c1e72Sjeremylt                   from quadrature points to nodes
861b11c1e72Sjeremylt   @param emode  \ref CEED_EVAL_INTERP to obtain interpolated values,
862b11c1e72Sjeremylt                   \ref CEED_EVAL_GRAD to obtain gradients.
863b11c1e72Sjeremylt   @param[in] u  Input array
864b11c1e72Sjeremylt   @param[out] v Output array
865b11c1e72Sjeremylt 
866b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
867dfdf5a53Sjeremylt 
868dfdf5a53Sjeremylt   @ref Advanced
869b11c1e72Sjeremylt **/
870d7b241e6Sjeremylt int CeedBasisApply(CeedBasis basis, CeedInt nelem, CeedTransposeMode tmode,
871aedaa0e5Sjeremylt                    CeedEvalMode emode, CeedVector u, CeedVector v) {
872d7b241e6Sjeremylt   int ierr;
8738795c945Sjeremylt   CeedInt ulength = 0, vlength, nnodes, nqpt;
874d7b241e6Sjeremylt   if (!basis->Apply) return CeedError(basis->ceed, 1,
875d7b241e6Sjeremylt                                         "Backend does not support BasisApply");
876b502e64cSValeria Barra   // check compatibility of topological and geometrical dimensions
8778795c945Sjeremylt   ierr = CeedBasisGetNumNodes(basis, &nnodes); CeedChk(ierr);
878b502e64cSValeria Barra   ierr = CeedBasisGetNumQuadraturePoints(basis, &nqpt); CeedChk(ierr);
879b502e64cSValeria Barra   ierr = CeedVectorGetLength(v, &vlength); CeedChk(ierr);
880b502e64cSValeria Barra 
881b502e64cSValeria Barra   if (u) {
882b502e64cSValeria Barra     ierr = CeedVectorGetLength(u, &ulength); CeedChk(ierr);
883b502e64cSValeria Barra   }
884b502e64cSValeria Barra 
885f90c8643Sjeremylt   if ((tmode == CEED_TRANSPOSE   && (vlength % nnodes != 0
886f90c8643Sjeremylt                                      || ulength % nqpt != 0))
887cdf4f918Sjeremylt       ||
8888795c945Sjeremylt       (tmode == CEED_NOTRANSPOSE && (ulength % nnodes != 0 || vlength % nqpt != 0)))
889b502e64cSValeria Barra     return CeedError(basis->ceed, 1,
890b502e64cSValeria Barra                      "Length of input/output vectors incompatible with basis dimensions");
891b502e64cSValeria Barra 
892d7b241e6Sjeremylt   ierr = basis->Apply(basis, nelem, tmode, emode, u, v); CeedChk(ierr);
893d7b241e6Sjeremylt   return 0;
894d7b241e6Sjeremylt }
895d7b241e6Sjeremylt 
896b11c1e72Sjeremylt /**
8974ce2993fSjeremylt   @brief Get Ceed associated with a CeedBasis
898b11c1e72Sjeremylt 
899b11c1e72Sjeremylt   @param basis      CeedBasis
9004ce2993fSjeremylt   @param[out] ceed  Variable to store Ceed
9014ce2993fSjeremylt 
9024ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
9034ce2993fSjeremylt 
90423617272Sjeremylt   @ref Advanced
9054ce2993fSjeremylt **/
9064ce2993fSjeremylt int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) {
9074ce2993fSjeremylt   *ceed = basis->ceed;
9084ce2993fSjeremylt 
9094ce2993fSjeremylt   return 0;
9104ce2993fSjeremylt };
9114ce2993fSjeremylt 
9124ce2993fSjeremylt /**
9134ce2993fSjeremylt   @brief Get dimension for given CeedBasis
9144ce2993fSjeremylt 
9154ce2993fSjeremylt   @param basis     CeedBasis
9164ce2993fSjeremylt   @param[out] dim  Variable to store dimension of basis
9174ce2993fSjeremylt 
9184ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
9194ce2993fSjeremylt 
92023617272Sjeremylt   @ref Advanced
9214ce2993fSjeremylt **/
9224ce2993fSjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) {
9234ce2993fSjeremylt   *dim = basis->dim;
9244ce2993fSjeremylt 
9254ce2993fSjeremylt   return 0;
9264ce2993fSjeremylt };
9274ce2993fSjeremylt 
9284ce2993fSjeremylt /**
9294ce2993fSjeremylt   @brief Get tensor status for given CeedBasis
9304ce2993fSjeremylt 
9314ce2993fSjeremylt   @param basis        CeedBasis
9324ce2993fSjeremylt   @param[out] tensor  Variable to store tensor status
9334ce2993fSjeremylt 
9344ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
9354ce2993fSjeremylt 
93623617272Sjeremylt   @ref Advanced
9374ce2993fSjeremylt **/
9384ce2993fSjeremylt int CeedBasisGetTensorStatus(CeedBasis basis, bool *tensor) {
9394ce2993fSjeremylt   *tensor = basis->tensorbasis;
9404ce2993fSjeremylt 
9414ce2993fSjeremylt   return 0;
9424ce2993fSjeremylt };
9434ce2993fSjeremylt 
9444ce2993fSjeremylt /**
9454ce2993fSjeremylt   @brief Get number of components for given CeedBasis
9464ce2993fSjeremylt 
9474ce2993fSjeremylt   @param basis        CeedBasis
948288c0443SJeremy L Thompson   @param[out] numcomp Variable to store number of components of basis
9494ce2993fSjeremylt 
9504ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
9514ce2993fSjeremylt 
95223617272Sjeremylt   @ref Advanced
9534ce2993fSjeremylt **/
9544ce2993fSjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *numcomp) {
9554ce2993fSjeremylt   *numcomp = basis->ncomp;
9564ce2993fSjeremylt 
9574ce2993fSjeremylt   return 0;
9584ce2993fSjeremylt };
9594ce2993fSjeremylt 
9604ce2993fSjeremylt /**
9614ce2993fSjeremylt   @brief Get total number of nodes (in 1 dimension) of a CeedBasis
9624ce2993fSjeremylt 
9634ce2993fSjeremylt   @param basis     CeedBasis
9644ce2993fSjeremylt   @param[out] P1d  Variable to store number of nodes
9654ce2993fSjeremylt 
9664ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
9674ce2993fSjeremylt 
96823617272Sjeremylt   @ref Advanced
9694ce2993fSjeremylt **/
9704ce2993fSjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P1d) {
9714ce2993fSjeremylt   if (!basis->tensorbasis) return CeedError(basis->ceed, 1,
9724ce2993fSjeremylt                                     "Cannot supply P1d for non-tensor basis");
9734ce2993fSjeremylt   *P1d = basis->P1d;
9744ce2993fSjeremylt   return 0;
9754ce2993fSjeremylt }
9764ce2993fSjeremylt 
9774ce2993fSjeremylt /**
9784ce2993fSjeremylt   @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis
9794ce2993fSjeremylt 
9804ce2993fSjeremylt   @param basis     CeedBasis
9814ce2993fSjeremylt   @param[out] Q1d  Variable to store number of quadrature points
9824ce2993fSjeremylt 
9834ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
9844ce2993fSjeremylt 
98523617272Sjeremylt   @ref Advanced
9864ce2993fSjeremylt **/
9874ce2993fSjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q1d) {
9884ce2993fSjeremylt   if (!basis->tensorbasis) return CeedError(basis->ceed, 1,
9894ce2993fSjeremylt                                     "Cannot supply Q1d for non-tensor basis");
9904ce2993fSjeremylt   *Q1d = basis->Q1d;
9914ce2993fSjeremylt   return 0;
9924ce2993fSjeremylt }
9934ce2993fSjeremylt 
9944ce2993fSjeremylt /**
9954ce2993fSjeremylt   @brief Get total number of nodes (in dim dimensions) of a CeedBasis
9964ce2993fSjeremylt 
9974ce2993fSjeremylt   @param basis   CeedBasis
9984ce2993fSjeremylt   @param[out] P  Variable to store number of nodes
999b11c1e72Sjeremylt 
1000b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1001dfdf5a53Sjeremylt 
1002dfdf5a53Sjeremylt   @ref Utility
1003b11c1e72Sjeremylt **/
1004d7b241e6Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) {
1005a8de75f0Sjeremylt   *P = basis->P;
1006d7b241e6Sjeremylt   return 0;
1007d7b241e6Sjeremylt }
1008d7b241e6Sjeremylt 
1009b11c1e72Sjeremylt /**
10104ce2993fSjeremylt   @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis
1011b11c1e72Sjeremylt 
1012b11c1e72Sjeremylt   @param basis   CeedBasis
10134ce2993fSjeremylt   @param[out] Q  Variable to store number of quadrature points
1014b11c1e72Sjeremylt 
1015b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1016dfdf5a53Sjeremylt 
1017dfdf5a53Sjeremylt   @ref Utility
1018b11c1e72Sjeremylt **/
1019d7b241e6Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) {
1020a8de75f0Sjeremylt   *Q = basis->Q;
1021d7b241e6Sjeremylt   return 0;
1022d7b241e6Sjeremylt }
1023d7b241e6Sjeremylt 
1024b11c1e72Sjeremylt /**
10258c91a0c9SJeremy L Thompson   @brief Get reference coordinates of quadrature points (in dim dimensions)
10264ce2993fSjeremylt          of a CeedBasis
10274ce2993fSjeremylt 
10284ce2993fSjeremylt   @param basis      CeedBasis
10298c91a0c9SJeremy L Thompson   @param[out] qref  Variable to store reference coordinates of quadrature points
10304ce2993fSjeremylt 
10314ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
10324ce2993fSjeremylt 
103323617272Sjeremylt   @ref Advanced
10344ce2993fSjeremylt **/
10354ce2993fSjeremylt int CeedBasisGetQRef(CeedBasis basis, CeedScalar* *qref) {
10364ce2993fSjeremylt   *qref = basis->qref1d;
10374ce2993fSjeremylt   return 0;
10384ce2993fSjeremylt }
10394ce2993fSjeremylt 
10404ce2993fSjeremylt /**
10414ce2993fSjeremylt   @brief Get quadrature weights of quadrature points (in dim dimensions)
10424ce2993fSjeremylt          of a CeedBasis
10434ce2993fSjeremylt 
10444ce2993fSjeremylt   @param basis         CeedBasis
10454ce2993fSjeremylt   @param[out] qweight  Variable to store quadrature weights
10464ce2993fSjeremylt 
10474ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
10484ce2993fSjeremylt 
104923617272Sjeremylt   @ref Advanced
10504ce2993fSjeremylt **/
10514ce2993fSjeremylt int CeedBasisGetQWeights(CeedBasis basis, CeedScalar* *qweight) {
10524ce2993fSjeremylt   *qweight = basis->qweight1d;
10534ce2993fSjeremylt   return 0;
10544ce2993fSjeremylt }
10554ce2993fSjeremylt 
10564ce2993fSjeremylt /**
10574ce2993fSjeremylt   @brief Get interpolation matrix of a CeedBasis
10584ce2993fSjeremylt 
10594ce2993fSjeremylt   @param basis       CeedBasis
1060288c0443SJeremy L Thompson   @param[out] interp Variable to store interpolation matrix
10614ce2993fSjeremylt 
10624ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
10634ce2993fSjeremylt 
106423617272Sjeremylt   @ref Advanced
10654ce2993fSjeremylt **/
10664ce2993fSjeremylt int CeedBasisGetInterp(CeedBasis basis, CeedScalar* *interp) {
10674ce2993fSjeremylt   *interp = basis->interp1d;
10684ce2993fSjeremylt   return 0;
10694ce2993fSjeremylt }
10704ce2993fSjeremylt 
10714ce2993fSjeremylt /**
10724ce2993fSjeremylt   @brief Get gradient matrix of a CeedBasis
10734ce2993fSjeremylt 
10744ce2993fSjeremylt   @param basis      CeedBasis
1075288c0443SJeremy L Thompson   @param[out] grad  Variable to store gradient matrix
10764ce2993fSjeremylt 
10774ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
10784ce2993fSjeremylt 
107923617272Sjeremylt   @ref Advanced
10804ce2993fSjeremylt **/
10814ce2993fSjeremylt int CeedBasisGetGrad(CeedBasis basis, CeedScalar* *grad) {
10824ce2993fSjeremylt   *grad = basis->grad1d;
10834ce2993fSjeremylt   return 0;
10844ce2993fSjeremylt }
10854ce2993fSjeremylt 
10864ce2993fSjeremylt /**
10874ce2993fSjeremylt   @brief Get backend data of a CeedBasis
10884ce2993fSjeremylt 
10894ce2993fSjeremylt   @param basis      CeedBasis
10904ce2993fSjeremylt   @param[out] data  Variable to store data
10914ce2993fSjeremylt 
10924ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
10934ce2993fSjeremylt 
109423617272Sjeremylt   @ref Advanced
10954ce2993fSjeremylt **/
10964ce2993fSjeremylt int CeedBasisGetData(CeedBasis basis, void* *data) {
10974ce2993fSjeremylt   *data = basis->data;
10984ce2993fSjeremylt   return 0;
10994ce2993fSjeremylt }
11004ce2993fSjeremylt 
11014ce2993fSjeremylt /**
1102fe2413ffSjeremylt   @brief Set backend data of a CeedBasis
1103fe2413ffSjeremylt 
1104fe2413ffSjeremylt   @param[out] basis CeedBasis
1105fe2413ffSjeremylt   @param data       Data to set
1106fe2413ffSjeremylt 
1107fe2413ffSjeremylt   @return An error code: 0 - success, otherwise - failure
1108fe2413ffSjeremylt 
1109fe2413ffSjeremylt   @ref Advanced
1110fe2413ffSjeremylt **/
1111fe2413ffSjeremylt int CeedBasisSetData(CeedBasis basis, void* *data) {
1112fe2413ffSjeremylt   basis->data = *data;
1113fe2413ffSjeremylt   return 0;
1114fe2413ffSjeremylt }
1115fe2413ffSjeremylt 
1116fe2413ffSjeremylt /**
11172f86a920SJeremy L Thompson   @brief Get CeedTensorContract of a CeedBasis
11182f86a920SJeremy L Thompson 
11192f86a920SJeremy L Thompson   @param basis          CeedBasis
11202f86a920SJeremy L Thompson   @param[out] contract  Variable to store CeedTensorContract
11212f86a920SJeremy L Thompson 
11222f86a920SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
11232f86a920SJeremy L Thompson 
11242f86a920SJeremy L Thompson   @ref Advanced
11252f86a920SJeremy L Thompson **/
11262f86a920SJeremy L Thompson int CeedBasisGetTensorContract(CeedBasis basis,
11272f86a920SJeremy L Thompson                                CeedTensorContract *contract) {
11282f86a920SJeremy L Thompson   *contract = basis->contract;
11292f86a920SJeremy L Thompson   return 0;
11302f86a920SJeremy L Thompson }
11312f86a920SJeremy L Thompson 
11322f86a920SJeremy L Thompson /**
11332f86a920SJeremy L Thompson   @brief Set CeedTensorContract of a CeedBasis
11342f86a920SJeremy L Thompson 
11352f86a920SJeremy L Thompson   @param[out] basis     CeedBasis
11362f86a920SJeremy L Thompson   @param contract       CeedTensorContract to set
11372f86a920SJeremy L Thompson 
11382f86a920SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
11392f86a920SJeremy L Thompson 
11402f86a920SJeremy L Thompson   @ref Advanced
11412f86a920SJeremy L Thompson **/
11422f86a920SJeremy L Thompson int CeedBasisSetTensorContract(CeedBasis basis,
11432f86a920SJeremy L Thompson                                CeedTensorContract *contract) {
11442f86a920SJeremy L Thompson   basis->contract = *contract;
11452f86a920SJeremy L Thompson   return 0;
11462f86a920SJeremy L Thompson }
11472f86a920SJeremy L Thompson 
11482f86a920SJeremy L Thompson /**
1149a8de75f0Sjeremylt   @brief Get dimension for given CeedElemTopology
1150a8de75f0Sjeremylt 
1151a8de75f0Sjeremylt   @param topo      CeedElemTopology
11524ce2993fSjeremylt   @param[out] dim  Variable to store dimension of topology
1153a8de75f0Sjeremylt 
1154a8de75f0Sjeremylt   @return An error code: 0 - success, otherwise - failure
1155a8de75f0Sjeremylt 
115623617272Sjeremylt   @ref Advanced
1157a8de75f0Sjeremylt **/
1158a8de75f0Sjeremylt int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) {
1159a8de75f0Sjeremylt   *dim = (CeedInt) topo >> 16;
1160a8de75f0Sjeremylt 
1161a8de75f0Sjeremylt   return 0;
1162a8de75f0Sjeremylt };
1163a8de75f0Sjeremylt 
1164a8de75f0Sjeremylt /**
1165b11c1e72Sjeremylt   @brief Destroy a CeedBasis
1166b11c1e72Sjeremylt 
1167b11c1e72Sjeremylt   @param basis CeedBasis to destroy
1168b11c1e72Sjeremylt 
1169b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1170dfdf5a53Sjeremylt 
1171dfdf5a53Sjeremylt   @ref Basic
1172b11c1e72Sjeremylt **/
1173d7b241e6Sjeremylt int CeedBasisDestroy(CeedBasis *basis) {
1174d7b241e6Sjeremylt   int ierr;
1175d7b241e6Sjeremylt 
1176d7b241e6Sjeremylt   if (!*basis || --(*basis)->refcount > 0) return 0;
1177d7b241e6Sjeremylt   if ((*basis)->Destroy) {
1178d7b241e6Sjeremylt     ierr = (*basis)->Destroy(*basis); CeedChk(ierr);
1179d7b241e6Sjeremylt   }
1180d7b241e6Sjeremylt   ierr = CeedFree(&(*basis)->interp1d); CeedChk(ierr);
1181d7b241e6Sjeremylt   ierr = CeedFree(&(*basis)->grad1d); CeedChk(ierr);
1182d7b241e6Sjeremylt   ierr = CeedFree(&(*basis)->qref1d); CeedChk(ierr);
1183d7b241e6Sjeremylt   ierr = CeedFree(&(*basis)->qweight1d); CeedChk(ierr);
1184d7b241e6Sjeremylt   ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr);
1185d7b241e6Sjeremylt   ierr = CeedFree(basis); CeedChk(ierr);
1186d7b241e6Sjeremylt   return 0;
1187d7b241e6Sjeremylt }
1188d7b241e6Sjeremylt 
118933e6becaSjeremylt /// @cond DOXYGEN_SKIP
11908795c945Sjeremylt // Indicate that the quadrature points are collocated with the nodes
1191783c99b3SValeria Barra CeedBasis CEED_BASIS_COLLOCATED = &ceed_basis_collocated;
119233e6becaSjeremylt /// @endcond
1193d7b241e6Sjeremylt /// @}
1194