xref: /libCEED/rust/libceed-sys/c-src/interface/ceed-basis.c (revision 1d102b48beba01a4e67ef28409a80cc546ae2651)
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)
64c042f62fSJeremy L Thompson     // LCOV_EXCL_START
654d537eeaSYohann     return CeedError(ceed, 1, "Basis dimension must be a positive value");
66c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
674d537eeaSYohann 
685fe0d4faSjeremylt   if (!ceed->BasisCreateTensorH1) {
695fe0d4faSjeremylt     Ceed delegate;
70aefd8378Sjeremylt     ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr);
715fe0d4faSjeremylt 
725fe0d4faSjeremylt     if (!delegate)
73c042f62fSJeremy L Thompson       // LCOV_EXCL_START
74d7b241e6Sjeremylt       return CeedError(ceed, 1, "Backend does not support BasisCreateTensorH1");
75c042f62fSJeremy L Thompson     // LCOV_EXCL_STOP
765fe0d4faSjeremylt 
775fe0d4faSjeremylt     ierr = CeedBasisCreateTensorH1(delegate, dim, ncomp, P1d,
785fe0d4faSjeremylt                                    Q1d, interp1d, grad1d, qref1d,
795fe0d4faSjeremylt                                    qweight1d, basis); CeedChk(ierr);
805fe0d4faSjeremylt     return 0;
815fe0d4faSjeremylt   }
82d7b241e6Sjeremylt   ierr = CeedCalloc(1,basis); CeedChk(ierr);
83d7b241e6Sjeremylt   (*basis)->ceed = ceed;
84d7b241e6Sjeremylt   ceed->refcount++;
85d7b241e6Sjeremylt   (*basis)->refcount = 1;
86a8de75f0Sjeremylt   (*basis)->tensorbasis = 1;
87d7b241e6Sjeremylt   (*basis)->dim = dim;
88d7b241e6Sjeremylt   (*basis)->ncomp = ncomp;
89d7b241e6Sjeremylt   (*basis)->P1d = P1d;
90d7b241e6Sjeremylt   (*basis)->Q1d = Q1d;
91a8de75f0Sjeremylt   (*basis)->P = CeedIntPow(P1d, dim);
92a8de75f0Sjeremylt   (*basis)->Q = CeedIntPow(Q1d, dim);
93d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d,&(*basis)->qref1d); CeedChk(ierr);
94d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d,&(*basis)->qweight1d); CeedChk(ierr);
95d7b241e6Sjeremylt   memcpy((*basis)->qref1d, qref1d, Q1d*sizeof(qref1d[0]));
96d7b241e6Sjeremylt   memcpy((*basis)->qweight1d, qweight1d, Q1d*sizeof(qweight1d[0]));
97d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d*P1d,&(*basis)->interp1d); CeedChk(ierr);
98d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d*P1d,&(*basis)->grad1d); CeedChk(ierr);
99d7b241e6Sjeremylt   memcpy((*basis)->interp1d, interp1d, Q1d*P1d*sizeof(interp1d[0]));
10009486605Sjeremylt   memcpy((*basis)->grad1d, grad1d, Q1d*P1d*sizeof(grad1d[0]));
101667bc5fcSjeremylt   ierr = ceed->BasisCreateTensorH1(dim, P1d, Q1d, interp1d, grad1d, qref1d,
102d7b241e6Sjeremylt                                    qweight1d, *basis); CeedChk(ierr);
103d7b241e6Sjeremylt   return 0;
104d7b241e6Sjeremylt }
105d7b241e6Sjeremylt 
106b11c1e72Sjeremylt /**
107b11c1e72Sjeremylt   @brief Create a tensor product Lagrange basis
108b11c1e72Sjeremylt 
109b11c1e72Sjeremylt   @param ceed       A Ceed object where the CeedBasis will be created
110b11c1e72Sjeremylt   @param dim        Topological dimension of element
111b11c1e72Sjeremylt   @param ncomp      Number of field components
112b11c1e72Sjeremylt   @param P          Number of Gauss-Lobatto nodes in one dimension.  The
113b11c1e72Sjeremylt                       polynomial degree of the resulting Q_k element is k=P-1.
114b11c1e72Sjeremylt   @param Q          Number of quadrature points in one dimension.
115b11c1e72Sjeremylt   @param qmode      Distribution of the Q quadrature points (affects order of
116b11c1e72Sjeremylt                       accuracy for the quadrature)
117b11c1e72Sjeremylt   @param[out] basis Address of the variable where the newly created
118b11c1e72Sjeremylt                       CeedBasis will be stored.
119b11c1e72Sjeremylt 
120b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
121dfdf5a53Sjeremylt 
122dfdf5a53Sjeremylt   @ref Basic
123b11c1e72Sjeremylt **/
124d7b241e6Sjeremylt int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt ncomp,
125d7b241e6Sjeremylt                                     CeedInt P, CeedInt Q,
126d7b241e6Sjeremylt                                     CeedQuadMode qmode, CeedBasis *basis) {
127d7b241e6Sjeremylt   // Allocate
128d7b241e6Sjeremylt   int ierr, i, j, k;
129d7b241e6Sjeremylt   CeedScalar c1, c2, c3, c4, dx, *nodes, *interp1d, *grad1d, *qref1d, *qweight1d;
1304d537eeaSYohann 
1314d537eeaSYohann   if (dim<1)
132c042f62fSJeremy L Thompson     // LCOV_EXCL_START
1334d537eeaSYohann     return CeedError(ceed, 1, "Basis dimension must be a positive value");
134c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
1354d537eeaSYohann 
136d7b241e6Sjeremylt   ierr = CeedCalloc(P*Q, &interp1d); CeedChk(ierr);
137d7b241e6Sjeremylt   ierr = CeedCalloc(P*Q, &grad1d); CeedChk(ierr);
138d7b241e6Sjeremylt   ierr = CeedCalloc(P, &nodes); CeedChk(ierr);
139d7b241e6Sjeremylt   ierr = CeedCalloc(Q, &qref1d); CeedChk(ierr);
140d7b241e6Sjeremylt   ierr = CeedCalloc(Q, &qweight1d); CeedChk(ierr);
141d7b241e6Sjeremylt   // Get Nodes and Weights
142d7b241e6Sjeremylt   ierr = CeedLobattoQuadrature(P, nodes, NULL); CeedChk(ierr);
143d7b241e6Sjeremylt   switch (qmode) {
144d7b241e6Sjeremylt   case CEED_GAUSS:
145d7b241e6Sjeremylt     ierr = CeedGaussQuadrature(Q, qref1d, qweight1d); CeedChk(ierr);
146d7b241e6Sjeremylt     break;
147d7b241e6Sjeremylt   case CEED_GAUSS_LOBATTO:
148d7b241e6Sjeremylt     ierr = CeedLobattoQuadrature(Q, qref1d, qweight1d); CeedChk(ierr);
149d7b241e6Sjeremylt     break;
150d7b241e6Sjeremylt   }
151d7b241e6Sjeremylt   // Build B, D matrix
152d7b241e6Sjeremylt   // Fornberg, 1998
153d7b241e6Sjeremylt   for (i = 0; i  < Q; i++) {
154d7b241e6Sjeremylt     c1 = 1.0;
155d7b241e6Sjeremylt     c3 = nodes[0] - qref1d[i];
156d7b241e6Sjeremylt     interp1d[i*P+0] = 1.0;
157d7b241e6Sjeremylt     for (j = 1; j < P; j++) {
158d7b241e6Sjeremylt       c2 = 1.0;
159d7b241e6Sjeremylt       c4 = c3;
160d7b241e6Sjeremylt       c3 = nodes[j] - qref1d[i];
161d7b241e6Sjeremylt       for (k = 0; k < j; k++) {
162d7b241e6Sjeremylt         dx = nodes[j] - nodes[k];
163d7b241e6Sjeremylt         c2 *= dx;
164d7b241e6Sjeremylt         if (k == j - 1) {
165d7b241e6Sjeremylt           grad1d[i*P + j] = c1*(interp1d[i*P + k] - c4*grad1d[i*P + k]) / c2;
166d7b241e6Sjeremylt           interp1d[i*P + j] = - c1*c4*interp1d[i*P + k] / c2;
167d7b241e6Sjeremylt         }
168d7b241e6Sjeremylt         grad1d[i*P + k] = (c3*grad1d[i*P + k] - interp1d[i*P + k]) / dx;
169d7b241e6Sjeremylt         interp1d[i*P + k] = c3*interp1d[i*P + k] / dx;
170d7b241e6Sjeremylt       }
171d7b241e6Sjeremylt       c1 = c2;
172d7b241e6Sjeremylt     }
173d7b241e6Sjeremylt   }
174d7b241e6Sjeremylt   //  // Pass to CeedBasisCreateTensorH1
175d7b241e6Sjeremylt   ierr = CeedBasisCreateTensorH1(ceed, dim, ncomp, P, Q, interp1d, grad1d, qref1d,
176d7b241e6Sjeremylt                                  qweight1d, basis); CeedChk(ierr);
177d7b241e6Sjeremylt   ierr = CeedFree(&interp1d); CeedChk(ierr);
178d7b241e6Sjeremylt   ierr = CeedFree(&grad1d); CeedChk(ierr);
179d7b241e6Sjeremylt   ierr = CeedFree(&nodes); CeedChk(ierr);
180d7b241e6Sjeremylt   ierr = CeedFree(&qref1d); CeedChk(ierr);
181d7b241e6Sjeremylt   ierr = CeedFree(&qweight1d); CeedChk(ierr);
182d7b241e6Sjeremylt   return 0;
183d7b241e6Sjeremylt }
184d7b241e6Sjeremylt 
185b11c1e72Sjeremylt /**
186a8de75f0Sjeremylt   @brief Create a non tensor product basis for H^1 discretizations
187a8de75f0Sjeremylt 
188a8de75f0Sjeremylt   @param ceed       A Ceed object where the CeedBasis will be created
189a8de75f0Sjeremylt   @param topo       Topology of element, e.g. hypercube, simplex, ect
190a8de75f0Sjeremylt   @param ncomp      Number of field components (1 for scalar fields)
1918795c945Sjeremylt   @param nnodes       Total number of nodes
192a8de75f0Sjeremylt   @param nqpts      Total number of quadrature points
1938795c945Sjeremylt   @param interp     Row-major nqpts × nnodes matrix expressing the values of
1948795c945Sjeremylt                       nodal basis functions at quadrature points
1958795c945Sjeremylt   @param grad       Row-major (nqpts x dim) × nnodes matrix expressing
1968795c945Sjeremylt                       derivatives of nodal basis functions at quadrature points
1978795c945Sjeremylt   @param qref       Array of length nqpts holding the locations of quadrature
1988795c945Sjeremylt                       points on the reference element [-1, 1]
199a8de75f0Sjeremylt   @param qweight    Array of length nqpts holding the quadrature weights on the
200a8de75f0Sjeremylt                       reference element
201a8de75f0Sjeremylt   @param[out] basis Address of the variable where the newly created
202a8de75f0Sjeremylt                       CeedBasis will be stored.
203a8de75f0Sjeremylt 
204a8de75f0Sjeremylt   @return An error code: 0 - success, otherwise - failure
205a8de75f0Sjeremylt 
206a8de75f0Sjeremylt   @ref Basic
207a8de75f0Sjeremylt **/
208a8de75f0Sjeremylt int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt ncomp,
2098795c945Sjeremylt                       CeedInt nnodes, CeedInt nqpts,
210a8de75f0Sjeremylt                       const CeedScalar *interp,
211a8de75f0Sjeremylt                       const CeedScalar *grad, const CeedScalar *qref,
212a8de75f0Sjeremylt                       const CeedScalar *qweight, CeedBasis *basis) {
213a8de75f0Sjeremylt   int ierr;
2148795c945Sjeremylt   CeedInt P = nnodes, Q = nqpts, dim = 0;
215a8de75f0Sjeremylt 
2165fe0d4faSjeremylt   if (!ceed->BasisCreateH1) {
2175fe0d4faSjeremylt     Ceed delegate;
218aefd8378Sjeremylt     ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr);
2195fe0d4faSjeremylt 
2205fe0d4faSjeremylt     if (!delegate)
221c042f62fSJeremy L Thompson       // LCOV_EXCL_START
222a8de75f0Sjeremylt       return CeedError(ceed, 1, "Backend does not support BasisCreateH1");
223c042f62fSJeremy L Thompson     // LCOV_EXCL_STOP
2245fe0d4faSjeremylt 
2258795c945Sjeremylt     ierr = CeedBasisCreateH1(delegate, topo, ncomp, nnodes,
2265fe0d4faSjeremylt                              nqpts, interp, grad, qref,
2275fe0d4faSjeremylt                              qweight, basis); CeedChk(ierr);
2285fe0d4faSjeremylt     return 0;
2295fe0d4faSjeremylt   }
2305fe0d4faSjeremylt 
231a8de75f0Sjeremylt   ierr = CeedCalloc(1,basis); CeedChk(ierr);
232a8de75f0Sjeremylt 
233a8de75f0Sjeremylt   ierr = CeedBasisGetTopologyDimension(topo, &dim); CeedChk(ierr);
234a8de75f0Sjeremylt 
235a8de75f0Sjeremylt   (*basis)->ceed = ceed;
236a8de75f0Sjeremylt   ceed->refcount++;
237a8de75f0Sjeremylt   (*basis)->refcount = 1;
238a8de75f0Sjeremylt   (*basis)->tensorbasis = 0;
239a8de75f0Sjeremylt   (*basis)->dim = dim;
240a8de75f0Sjeremylt   (*basis)->ncomp = ncomp;
241a8de75f0Sjeremylt   (*basis)->P = P;
242a8de75f0Sjeremylt   (*basis)->Q = Q;
243a8de75f0Sjeremylt   ierr = CeedMalloc(Q*dim,&(*basis)->qref1d); CeedChk(ierr);
244a8de75f0Sjeremylt   ierr = CeedMalloc(Q,&(*basis)->qweight1d); CeedChk(ierr);
245a8de75f0Sjeremylt   memcpy((*basis)->qref1d, qref, Q*dim*sizeof(qref[0]));
246a8de75f0Sjeremylt   memcpy((*basis)->qweight1d, qweight, Q*sizeof(qweight[0]));
247a8de75f0Sjeremylt   ierr = CeedMalloc(Q*P,&(*basis)->interp1d); CeedChk(ierr);
248a8de75f0Sjeremylt   ierr = CeedMalloc(dim*Q*P,&(*basis)->grad1d); CeedChk(ierr);
249a8de75f0Sjeremylt   memcpy((*basis)->interp1d, interp, Q*P*sizeof(interp[0]));
250a8de75f0Sjeremylt   memcpy((*basis)->grad1d, grad, dim*Q*P*sizeof(grad[0]));
251667bc5fcSjeremylt   ierr = ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, qref,
252a8de75f0Sjeremylt                              qweight, *basis); CeedChk(ierr);
253a8de75f0Sjeremylt   return 0;
254a8de75f0Sjeremylt }
255a8de75f0Sjeremylt 
256a8de75f0Sjeremylt /**
257b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre quadrature
258b11c1e72Sjeremylt 
259b11c1e72Sjeremylt   @param Q              Number of quadrature points (integrates polynomials of
260b11c1e72Sjeremylt                           degree 2*Q-1 exactly)
261b11c1e72Sjeremylt   @param[out] qref1d    Array of length Q to hold the abscissa on [-1, 1]
262b11c1e72Sjeremylt   @param[out] qweight1d Array of length Q to hold the weights
263b11c1e72Sjeremylt 
264b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
265dfdf5a53Sjeremylt 
266dfdf5a53Sjeremylt   @ref Utility
267b11c1e72Sjeremylt **/
268d7b241e6Sjeremylt int CeedGaussQuadrature(CeedInt Q, CeedScalar *qref1d, CeedScalar *qweight1d) {
269d7b241e6Sjeremylt   // Allocate
270d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0*atan(1.0);
271d7b241e6Sjeremylt   // Build qref1d, qweight1d
272d7b241e6Sjeremylt   for (int i = 0; i <= Q/2; i++) {
273d7b241e6Sjeremylt     // Guess
274d7b241e6Sjeremylt     xi = cos(PI*(CeedScalar)(2*i+1)/((CeedScalar)(2*Q)));
275d7b241e6Sjeremylt     // Pn(xi)
276d7b241e6Sjeremylt     P0 = 1.0;
277d7b241e6Sjeremylt     P1 = xi;
278d7b241e6Sjeremylt     P2 = 0.0;
279d7b241e6Sjeremylt     for (int j = 2; j <= Q; j++) {
280d7b241e6Sjeremylt       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
281d7b241e6Sjeremylt       P0 = P1;
282d7b241e6Sjeremylt       P1 = P2;
283d7b241e6Sjeremylt     }
284d7b241e6Sjeremylt     // First Newton Step
285d7b241e6Sjeremylt     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
286d7b241e6Sjeremylt     xi = xi-P2/dP2;
287d7b241e6Sjeremylt     // Newton to convergence
288d7b241e6Sjeremylt     for (int k=0; k<100 && fabs(P2)>1e-15; k++) {
289d7b241e6Sjeremylt       P0 = 1.0;
290d7b241e6Sjeremylt       P1 = xi;
291d7b241e6Sjeremylt       for (int j = 2; j <= Q; j++) {
292d7b241e6Sjeremylt         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
293d7b241e6Sjeremylt         P0 = P1;
294d7b241e6Sjeremylt         P1 = P2;
295d7b241e6Sjeremylt       }
296d7b241e6Sjeremylt       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
297d7b241e6Sjeremylt       xi = xi-P2/dP2;
298d7b241e6Sjeremylt     }
299d7b241e6Sjeremylt     // Save xi, wi
300d7b241e6Sjeremylt     wi = 2.0/((1.0-xi*xi)*dP2*dP2);
301d7b241e6Sjeremylt     qweight1d[i] = wi;
302d7b241e6Sjeremylt     qweight1d[Q-1-i] = wi;
303d7b241e6Sjeremylt     qref1d[i] = -xi;
304d7b241e6Sjeremylt     qref1d[Q-1-i]= xi;
305d7b241e6Sjeremylt   }
306d7b241e6Sjeremylt   return 0;
307d7b241e6Sjeremylt }
308d7b241e6Sjeremylt 
309b11c1e72Sjeremylt /**
310b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre-Lobatto quadrature
311b11c1e72Sjeremylt 
312b11c1e72Sjeremylt   @param Q              Number of quadrature points (integrates polynomials of
313b11c1e72Sjeremylt                           degree 2*Q-3 exactly)
314b11c1e72Sjeremylt   @param[out] qref1d    Array of length Q to hold the abscissa on [-1, 1]
315b11c1e72Sjeremylt   @param[out] qweight1d Array of length Q to hold the weights
316b11c1e72Sjeremylt 
317b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
318dfdf5a53Sjeremylt 
319dfdf5a53Sjeremylt   @ref Utility
320b11c1e72Sjeremylt **/
321d7b241e6Sjeremylt int CeedLobattoQuadrature(CeedInt Q, CeedScalar *qref1d,
322d7b241e6Sjeremylt                           CeedScalar *qweight1d) {
323d7b241e6Sjeremylt   // Allocate
324d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0*atan(1.0);
325d7b241e6Sjeremylt   // Build qref1d, qweight1d
326d7b241e6Sjeremylt   // Set endpoints
327d7b241e6Sjeremylt   wi = 2.0/((CeedScalar)(Q*(Q-1)));
328d7b241e6Sjeremylt   if (qweight1d) {
329d7b241e6Sjeremylt     qweight1d[0] = wi;
330d7b241e6Sjeremylt     qweight1d[Q-1] = wi;
331d7b241e6Sjeremylt   }
332d7b241e6Sjeremylt   qref1d[0] = -1.0;
333d7b241e6Sjeremylt   qref1d[Q-1] = 1.0;
334d7b241e6Sjeremylt   // Interior
335d7b241e6Sjeremylt   for (int i = 1; i <= (Q-1)/2; i++) {
336d7b241e6Sjeremylt     // Guess
337d7b241e6Sjeremylt     xi = cos(PI*(CeedScalar)(i)/(CeedScalar)(Q-1));
338d7b241e6Sjeremylt     // Pn(xi)
339d7b241e6Sjeremylt     P0 = 1.0;
340d7b241e6Sjeremylt     P1 = xi;
341d7b241e6Sjeremylt     P2 = 0.0;
342d7b241e6Sjeremylt     for (int j = 2; j < Q; j++) {
343d7b241e6Sjeremylt       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
344d7b241e6Sjeremylt       P0 = P1;
345d7b241e6Sjeremylt       P1 = P2;
346d7b241e6Sjeremylt     }
347d7b241e6Sjeremylt     // First Newton step
348d7b241e6Sjeremylt     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
349d7b241e6Sjeremylt     d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
350d7b241e6Sjeremylt     xi = xi-dP2/d2P2;
351d7b241e6Sjeremylt     // Newton to convergence
352d7b241e6Sjeremylt     for (int k=0; k<100 && fabs(dP2)>1e-15; k++) {
353d7b241e6Sjeremylt       P0 = 1.0;
354d7b241e6Sjeremylt       P1 = xi;
355d7b241e6Sjeremylt       for (int j = 2; j < Q; j++) {
356d7b241e6Sjeremylt         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
357d7b241e6Sjeremylt         P0 = P1;
358d7b241e6Sjeremylt         P1 = P2;
359d7b241e6Sjeremylt       }
360d7b241e6Sjeremylt       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
361d7b241e6Sjeremylt       d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
362d7b241e6Sjeremylt       xi = xi-dP2/d2P2;
363d7b241e6Sjeremylt     }
364d7b241e6Sjeremylt     // Save xi, wi
365d7b241e6Sjeremylt     wi = 2.0/(((CeedScalar)(Q*(Q-1)))*P2*P2);
366d7b241e6Sjeremylt     if (qweight1d) {
367d7b241e6Sjeremylt       qweight1d[i] = wi;
368d7b241e6Sjeremylt       qweight1d[Q-1-i] = wi;
369d7b241e6Sjeremylt     }
370d7b241e6Sjeremylt     qref1d[i] = -xi;
371d7b241e6Sjeremylt     qref1d[Q-1-i]= xi;
372d7b241e6Sjeremylt   }
373d7b241e6Sjeremylt   return 0;
374d7b241e6Sjeremylt }
375d7b241e6Sjeremylt 
376dfdf5a53Sjeremylt /**
377dfdf5a53Sjeremylt   @brief View an array stored in a CeedBasis
378dfdf5a53Sjeremylt 
379dfdf5a53Sjeremylt   @param name      Name of array
380dfdf5a53Sjeremylt   @param fpformat  Printing format
381dfdf5a53Sjeremylt   @param m         Number of rows in array
382dfdf5a53Sjeremylt   @param n         Number of columns in array
383dfdf5a53Sjeremylt   @param a         Array to be viewed
384dfdf5a53Sjeremylt   @param stream    Stream to view to, e.g., stdout
385dfdf5a53Sjeremylt 
386dfdf5a53Sjeremylt   @return An error code: 0 - success, otherwise - failure
387dfdf5a53Sjeremylt 
388dfdf5a53Sjeremylt   @ref Utility
389dfdf5a53Sjeremylt **/
390d7b241e6Sjeremylt static int CeedScalarView(const char *name, const char *fpformat, CeedInt m,
391d7b241e6Sjeremylt                           CeedInt n, const CeedScalar *a, FILE *stream) {
392d7b241e6Sjeremylt   for (int i=0; i<m; i++) {
393*1d102b48SJeremy L Thompson     if (m > 1)
394*1d102b48SJeremy L Thompson       fprintf(stream, "%12s[%d]:", name, i);
395*1d102b48SJeremy L Thompson     else
396*1d102b48SJeremy L Thompson       fprintf(stream, "%12s:", name);
397*1d102b48SJeremy L Thompson     for (int j=0; j<n; j++)
398d7b241e6Sjeremylt       fprintf(stream, fpformat, fabs(a[i*n+j]) > 1E-14 ? a[i*n+j] : 0);
399d7b241e6Sjeremylt     fputs("\n", stream);
400d7b241e6Sjeremylt   }
401d7b241e6Sjeremylt   return 0;
402d7b241e6Sjeremylt }
403d7b241e6Sjeremylt 
404b11c1e72Sjeremylt /**
405b11c1e72Sjeremylt   @brief View a CeedBasis
406b11c1e72Sjeremylt 
407b11c1e72Sjeremylt   @param basis  CeedBasis to view
408b11c1e72Sjeremylt   @param stream Stream to view to, e.g., stdout
409b11c1e72Sjeremylt 
410b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
411dfdf5a53Sjeremylt 
412dfdf5a53Sjeremylt   @ref Utility
413b11c1e72Sjeremylt **/
414d7b241e6Sjeremylt int CeedBasisView(CeedBasis basis, FILE *stream) {
415d7b241e6Sjeremylt   int ierr;
416d7b241e6Sjeremylt 
417a8de75f0Sjeremylt   if (basis->tensorbasis) {
418d7b241e6Sjeremylt     fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P1d,
419d7b241e6Sjeremylt             basis->Q1d);
420d7b241e6Sjeremylt     ierr = CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q1d, basis->qref1d,
421d7b241e6Sjeremylt                           stream); CeedChk(ierr);
4228795c945Sjeremylt     ierr = CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q1d,
4238795c945Sjeremylt                           basis->qweight1d, stream); CeedChk(ierr);
424d7b241e6Sjeremylt     ierr = CeedScalarView("interp1d", "\t% 12.8f", basis->Q1d, basis->P1d,
425d7b241e6Sjeremylt                           basis->interp1d, stream); CeedChk(ierr);
426d7b241e6Sjeremylt     ierr = CeedScalarView("grad1d", "\t% 12.8f", basis->Q1d, basis->P1d,
427d7b241e6Sjeremylt                           basis->grad1d, stream); CeedChk(ierr);
428a8de75f0Sjeremylt   } else {
429a8de75f0Sjeremylt     fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P,
430a8de75f0Sjeremylt             basis->Q);
431a8de75f0Sjeremylt     ierr = CeedScalarView("qref", "\t% 12.8f", 1, basis->Q*basis->dim,
432a8de75f0Sjeremylt                           basis->qref1d,
433a8de75f0Sjeremylt                           stream); CeedChk(ierr);
434a8de75f0Sjeremylt     ierr = CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->qweight1d,
435a8de75f0Sjeremylt                           stream); CeedChk(ierr);
436a8de75f0Sjeremylt     ierr = CeedScalarView("interp", "\t% 12.8f", basis->Q, basis->P,
437a8de75f0Sjeremylt                           basis->interp1d, stream); CeedChk(ierr);
438a8de75f0Sjeremylt     ierr = CeedScalarView("grad", "\t% 12.8f", basis->dim*basis->Q, basis->P,
439a8de75f0Sjeremylt                           basis->grad1d, stream); CeedChk(ierr);
440a8de75f0Sjeremylt   }
441d7b241e6Sjeremylt   return 0;
442d7b241e6Sjeremylt }
443d7b241e6Sjeremylt 
444dfdf5a53Sjeremylt /**
44552bfb9bbSJeremy L Thompson   @brief Compute Householder reflection
446dfdf5a53Sjeremylt 
447dfdf5a53Sjeremylt     Computes A = (I - b v v^T) A
448dfdf5a53Sjeremylt     where A is an mxn matrix indexed as A[i*row + j*col]
449dfdf5a53Sjeremylt 
45052bfb9bbSJeremy L Thompson   @param[in,out] A  Matrix to apply Householder reflection to, in place
451dfdf5a53Sjeremylt   @param v          Householder vector
452dfdf5a53Sjeremylt   @param b          Scaling factor
453dfdf5a53Sjeremylt   @param m          Number of rows in A
454dfdf5a53Sjeremylt   @param n          Number of columns in A
45552bfb9bbSJeremy L Thompson   @param row        Row stride
45652bfb9bbSJeremy L Thompson   @param col        Col stride
457dfdf5a53Sjeremylt 
458dfdf5a53Sjeremylt   @return An error code: 0 - success, otherwise - failure
459dfdf5a53Sjeremylt 
460dfdf5a53Sjeremylt   @ref Developer
461dfdf5a53Sjeremylt **/
462d7b241e6Sjeremylt static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v,
463d7b241e6Sjeremylt                                   CeedScalar b, CeedInt m, CeedInt n,
464d7b241e6Sjeremylt                                   CeedInt row, CeedInt col) {
465d7b241e6Sjeremylt   for (CeedInt j=0; j<n; j++) {
466d7b241e6Sjeremylt     CeedScalar w = A[0*row + j*col];
467*1d102b48SJeremy L Thompson     for (CeedInt i=1; i<m; i++)
468*1d102b48SJeremy L Thompson       w += v[i] * A[i*row + j*col];
469d7b241e6Sjeremylt     A[0*row + j*col] -= b * w;
470*1d102b48SJeremy L Thompson     for (CeedInt i=1; i<m; i++)
471*1d102b48SJeremy L Thompson       A[i*row + j*col] -= b * w * v[i];
472d7b241e6Sjeremylt   }
473d7b241e6Sjeremylt   return 0;
474d7b241e6Sjeremylt }
475d7b241e6Sjeremylt 
476dfdf5a53Sjeremylt /**
477dfdf5a53Sjeremylt   @brief Apply Householder Q matrix
478dfdf5a53Sjeremylt 
47952bfb9bbSJeremy L Thompson     Compute A = Q A where Q is mxm and A is mxn.
480dfdf5a53Sjeremylt 
48152bfb9bbSJeremy L Thompson   @param[in,out] A  Matrix to apply Householder Q to, in place
482dfdf5a53Sjeremylt   @param Q          Householder Q matrix
483dfdf5a53Sjeremylt   @param tau        Householder scaling factors
484dfdf5a53Sjeremylt   @param tmode      Transpose mode for application
485dfdf5a53Sjeremylt   @param m          Number of rows in A
486dfdf5a53Sjeremylt   @param n          Number of columns in A
48752bfb9bbSJeremy L Thompson   @param k          Number of elementary reflectors in Q, k<m
48852bfb9bbSJeremy L Thompson   @param row        Row stride in A
48952bfb9bbSJeremy L Thompson   @param col        Col stride in A
490dfdf5a53Sjeremylt 
491dfdf5a53Sjeremylt   @return An error code: 0 - success, otherwise - failure
492dfdf5a53Sjeremylt 
493dfdf5a53Sjeremylt   @ref Developer
494dfdf5a53Sjeremylt **/
495d7b241e6Sjeremylt static int CeedHouseholderApplyQ(CeedScalar *A, const CeedScalar *Q,
496d7b241e6Sjeremylt                                  const CeedScalar *tau, CeedTransposeMode tmode,
497d7b241e6Sjeremylt                                  CeedInt m, CeedInt n, CeedInt k,
498d7b241e6Sjeremylt                                  CeedInt row, CeedInt col) {
499d7b241e6Sjeremylt   CeedScalar v[m];
500d7b241e6Sjeremylt   for (CeedInt ii=0; ii<k; ii++) {
501d7b241e6Sjeremylt     CeedInt i = tmode == CEED_TRANSPOSE ? ii : k-1-ii;
50252bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<m; j++)
503d7b241e6Sjeremylt       v[j] = Q[j*k+i];
504d7b241e6Sjeremylt     // Apply Householder reflector (I - tau v v^T) colograd1d^T
505d7b241e6Sjeremylt     CeedHouseholderReflect(&A[i*row], &v[i], tau[i], m-i, n, row, col);
506d7b241e6Sjeremylt   }
507d7b241e6Sjeremylt   return 0;
508d7b241e6Sjeremylt }
509d7b241e6Sjeremylt 
510b11c1e72Sjeremylt /**
51152bfb9bbSJeremy L Thompson   @brief Compute Givens rotation
51252bfb9bbSJeremy L Thompson 
51352bfb9bbSJeremy L Thompson     Computes A = G A (or G^T A in transpose mode)
51452bfb9bbSJeremy L Thompson     where A is an mxn matrix indexed as A[i*n + j*m]
51552bfb9bbSJeremy L Thompson 
51652bfb9bbSJeremy L Thompson   @param[in,out] A  Row major matrix to apply Givens rotation to, in place
51752bfb9bbSJeremy L Thompson   @param c          Cosine factor
51852bfb9bbSJeremy L Thompson   @param s          Sine factor
51952bfb9bbSJeremy L Thompson   @param i          First row/column to apply rotation
52052bfb9bbSJeremy L Thompson   @param k          Second row/column to apply rotation
52152bfb9bbSJeremy L Thompson   @param m          Number of rows in A
52252bfb9bbSJeremy L Thompson   @param n          Number of columns in A
52352bfb9bbSJeremy L Thompson 
52452bfb9bbSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
52552bfb9bbSJeremy L Thompson 
52652bfb9bbSJeremy L Thompson   @ref Developer
52752bfb9bbSJeremy L Thompson **/
52852bfb9bbSJeremy L Thompson static int CeedGivensRotation(CeedScalar *A, CeedScalar c, CeedScalar s,
52952bfb9bbSJeremy L Thompson                               CeedTransposeMode tmode, CeedInt i, CeedInt k,
53052bfb9bbSJeremy L Thompson                               CeedInt m, CeedInt n) {
53152bfb9bbSJeremy L Thompson   CeedInt stridej = 1, strideik = m, numits = n;
53252bfb9bbSJeremy L Thompson   if (tmode == CEED_NOTRANSPOSE) {
53352bfb9bbSJeremy L Thompson     stridej = n; strideik = 1; numits = m;
53452bfb9bbSJeremy L Thompson   }
53552bfb9bbSJeremy L Thompson 
53652bfb9bbSJeremy L Thompson   // Apply rotation
53752bfb9bbSJeremy L Thompson   for (CeedInt j=0; j<numits; j++) {
53852bfb9bbSJeremy L Thompson     CeedScalar tau1 = A[i*strideik+j*stridej], tau2 = A[k*strideik+j*stridej];
53952bfb9bbSJeremy L Thompson     A[i*strideik+j*stridej] = c*tau1 - s*tau2;
54052bfb9bbSJeremy L Thompson     A[k*strideik+j*stridej] = s*tau1 + c*tau2;
54152bfb9bbSJeremy L Thompson   }
54252bfb9bbSJeremy L Thompson 
54352bfb9bbSJeremy L Thompson   return 0;
54452bfb9bbSJeremy L Thompson }
54552bfb9bbSJeremy L Thompson 
54652bfb9bbSJeremy L Thompson /**
547b11c1e72Sjeremylt   @brief Return QR Factorization of matrix
548b11c1e72Sjeremylt 
549288c0443SJeremy L Thompson   @param ceed         A Ceed object currently in use
55052bfb9bbSJeremy L Thompson   @param[in,out] mat  Row-major matrix to be factorized in place
55152bfb9bbSJeremy L Thompson   @param[in,out] tau  Vector of length m of scaling factors
552b11c1e72Sjeremylt   @param m            Number of rows
553b11c1e72Sjeremylt   @param n            Number of columns
554b11c1e72Sjeremylt 
555b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
556dfdf5a53Sjeremylt 
557dfdf5a53Sjeremylt   @ref Utility
558b11c1e72Sjeremylt **/
559a7bd39daSjeremylt int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau,
560d7b241e6Sjeremylt                         CeedInt m, CeedInt n) {
561d7b241e6Sjeremylt   CeedScalar v[m];
562d7b241e6Sjeremylt 
563a7bd39daSjeremylt   // Check m >= n
564a7bd39daSjeremylt   if (n > m)
565c042f62fSJeremy L Thompson     // LCOV_EXCL_START
566a7bd39daSjeremylt     return CeedError(ceed, 1, "Cannot compute QR factorization with n > m");
567c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
568a7bd39daSjeremylt 
56952bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++) {
570d7b241e6Sjeremylt     // Calculate Householder vector, magnitude
571d7b241e6Sjeremylt     CeedScalar sigma = 0.0;
572d7b241e6Sjeremylt     v[i] = mat[i+n*i];
57352bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<m; j++) {
574d7b241e6Sjeremylt       v[j] = mat[i+n*j];
575d7b241e6Sjeremylt       sigma += v[j] * v[j];
576d7b241e6Sjeremylt     }
577d7b241e6Sjeremylt     CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m]
578d7b241e6Sjeremylt     CeedScalar Rii = -copysign(norm, v[i]);
579d7b241e6Sjeremylt     v[i] -= Rii;
580d7b241e6Sjeremylt     // norm of v[i:m] after modification above and scaling below
581d7b241e6Sjeremylt     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
582d7b241e6Sjeremylt     //   tau = 2 / (norm*norm)
583d7b241e6Sjeremylt     tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma);
584*1d102b48SJeremy L Thompson     for (CeedInt j=i+1; j<m; j++)
585*1d102b48SJeremy L Thompson       v[j] /= v[i];
586d7b241e6Sjeremylt 
587d7b241e6Sjeremylt     // Apply Householder reflector to lower right panel
588d7b241e6Sjeremylt     CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1);
589d7b241e6Sjeremylt     // Save v
590d7b241e6Sjeremylt     mat[i+n*i] = Rii;
591*1d102b48SJeremy L Thompson     for (CeedInt j=i+1; j<m; j++)
592d7b241e6Sjeremylt       mat[i+n*j] = v[j];
593d7b241e6Sjeremylt   }
594d7b241e6Sjeremylt 
595d7b241e6Sjeremylt   return 0;
596d7b241e6Sjeremylt }
597d7b241e6Sjeremylt 
598b11c1e72Sjeremylt /**
59952bfb9bbSJeremy L Thompson   @brief Return symmetric Schur decomposition of the symmetric matrix mat via
60052bfb9bbSJeremy L Thompson            symmetric QR factorization
60152bfb9bbSJeremy L Thompson 
60252bfb9bbSJeremy L Thompson   @param[in,out] mat  Row-major matrix to be factorized in place
60352bfb9bbSJeremy L Thompson   @param[out] lambda  Vector of length m of eigenvalues
60452bfb9bbSJeremy L Thompson   @param n            Number of rows/columns
60552bfb9bbSJeremy L Thompson 
60652bfb9bbSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
60752bfb9bbSJeremy L Thompson 
60852bfb9bbSJeremy L Thompson   @ref Utility
60952bfb9bbSJeremy L Thompson **/
61052bfb9bbSJeremy L Thompson int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat,
61152bfb9bbSJeremy L Thompson                                     CeedScalar *lambda, CeedInt n) {
61252bfb9bbSJeremy L Thompson   // Check bounds for clang-tidy
61352bfb9bbSJeremy L Thompson   if (n<2)
614c042f62fSJeremy L Thompson     // LCOV_EXCL_START
615c042f62fSJeremy L Thompson     return CeedError(ceed, 1,
616c042f62fSJeremy L Thompson                      "Cannot compute symmetric Schur decomposition of scalars");
617c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
61852bfb9bbSJeremy L Thompson 
61952bfb9bbSJeremy L Thompson   CeedScalar v[n-1], tau[n-1], matT[n*n];
62052bfb9bbSJeremy L Thompson 
62152bfb9bbSJeremy L Thompson   // Copy mat to matT and set mat to I
62252bfb9bbSJeremy L Thompson   memcpy(matT, mat, n*n*sizeof(mat[0]));
62352bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
62452bfb9bbSJeremy L Thompson     for (CeedInt j=0; j<n; j++)
62552bfb9bbSJeremy L Thompson       mat[j+n*i] = (i==j) ? 1 : 0;
62652bfb9bbSJeremy L Thompson 
62752bfb9bbSJeremy L Thompson   // Reduce to tridiagonal
62852bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n-1; i++) {
62952bfb9bbSJeremy L Thompson     // Calculate Householder vector, magnitude
63052bfb9bbSJeremy L Thompson     CeedScalar sigma = 0.0;
63152bfb9bbSJeremy L Thompson     v[i] = matT[i+n*(i+1)];
63252bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) {
63352bfb9bbSJeremy L Thompson       v[j] = matT[i+n*(j+1)];
63452bfb9bbSJeremy L Thompson       sigma += v[j] * v[j];
63552bfb9bbSJeremy L Thompson     }
63652bfb9bbSJeremy L Thompson     CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:n-1]
63752bfb9bbSJeremy L Thompson     CeedScalar Rii = -copysign(norm, v[i]);
63852bfb9bbSJeremy L Thompson     v[i] -= Rii;
63952bfb9bbSJeremy L Thompson     // norm of v[i:m] after modification above and scaling below
64052bfb9bbSJeremy L Thompson     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
64152bfb9bbSJeremy L Thompson     //   tau = 2 / (norm*norm)
64252bfb9bbSJeremy L Thompson     tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma);
64352bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) v[j] /= v[i];
64452bfb9bbSJeremy L Thompson 
64552bfb9bbSJeremy L Thompson     // Update sub and super diagonal
64652bfb9bbSJeremy L Thompson     matT[i+n*(i+1)] = Rii;
64752bfb9bbSJeremy L Thompson     matT[(i+1)+n*i] = Rii;
64852bfb9bbSJeremy L Thompson     for (CeedInt j=i+2; j<n; j++) {
64952bfb9bbSJeremy L Thompson       matT[i+n*j] = 0; matT[j+n*i] = 0;
65052bfb9bbSJeremy L Thompson     }
65152bfb9bbSJeremy L Thompson     // Apply symmetric Householder reflector to lower right panel
65252bfb9bbSJeremy L Thompson     CeedHouseholderReflect(&matT[(i+1)+n*(i+1)], &v[i], tau[i],
65352bfb9bbSJeremy L Thompson                            n-(i+1), n-(i+1), n, 1);
65452bfb9bbSJeremy L Thompson     CeedHouseholderReflect(&matT[(i+1)+n*(i+1)], &v[i], tau[i],
65552bfb9bbSJeremy L Thompson                            n-(i+1), n-(i+1), 1, n);
65652bfb9bbSJeremy L Thompson     // Save v
65752bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) {
65852bfb9bbSJeremy L Thompson       matT[i+n*(j+1)] = v[j];
65952bfb9bbSJeremy L Thompson     }
66052bfb9bbSJeremy L Thompson   }
66152bfb9bbSJeremy L Thompson   // Backwards accumulation of Q
66252bfb9bbSJeremy L Thompson   for (CeedInt i=n-2; i>=0; i--) {
66352bfb9bbSJeremy L Thompson     v[i] = 1;
66452bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) {
66552bfb9bbSJeremy L Thompson       v[j] = matT[i+n*(j+1)];
66652bfb9bbSJeremy L Thompson       matT[i+n*(j+1)] = 0;
66752bfb9bbSJeremy L Thompson     }
66852bfb9bbSJeremy L Thompson     CeedHouseholderReflect(&mat[(i+1)+n*(i+1)], &v[i], tau[i],
66952bfb9bbSJeremy L Thompson                            n-(i+1), n-(i+1), n, 1);
67052bfb9bbSJeremy L Thompson   }
67152bfb9bbSJeremy L Thompson 
67252bfb9bbSJeremy L Thompson   // Reduce sub and super diagonal
67352bfb9bbSJeremy L Thompson   CeedInt p = 0, q = 0, itr = 0, maxitr = n*n*n;
67452bfb9bbSJeremy L Thompson   CeedScalar tol = 1e-15;
67552bfb9bbSJeremy L Thompson 
67652bfb9bbSJeremy L Thompson   while (q < n && itr < maxitr) {
67752bfb9bbSJeremy L Thompson     // Update p, q, size of reduced portions of diagonal
67852bfb9bbSJeremy L Thompson     p = 0; q = 0;
67952bfb9bbSJeremy L Thompson     for (CeedInt i=n-2; i>=0; i--) {
68052bfb9bbSJeremy L Thompson       if (fabs(matT[i+n*(i+1)]) < tol)
68152bfb9bbSJeremy L Thompson         q += 1;
68252bfb9bbSJeremy L Thompson       else
68352bfb9bbSJeremy L Thompson          break;
68452bfb9bbSJeremy L Thompson     }
68552bfb9bbSJeremy L Thompson     for (CeedInt i=0; i<n-1-q; i++) {
68652bfb9bbSJeremy L Thompson       if (fabs(matT[i+n*(i+1)]) < tol)
68752bfb9bbSJeremy L Thompson         p += 1;
68852bfb9bbSJeremy L Thompson       else
68952bfb9bbSJeremy L Thompson         break;
69052bfb9bbSJeremy L Thompson     }
69152bfb9bbSJeremy L Thompson     if (q == n-1) break; // Finished reducing
69252bfb9bbSJeremy L Thompson 
69352bfb9bbSJeremy L Thompson     // Reduce tridiagonal portion
69452bfb9bbSJeremy L Thompson     CeedScalar tnn = matT[(n-1-q)+n*(n-1-q)],
69552bfb9bbSJeremy L Thompson                tnnm1 = matT[(n-2-q)+n*(n-1-q)];
69652bfb9bbSJeremy L Thompson     CeedScalar d = (matT[(n-2-q)+n*(n-2-q)] - tnn)/2;
69752bfb9bbSJeremy L Thompson     CeedScalar mu = tnn - tnnm1*tnnm1 /
69852bfb9bbSJeremy L Thompson                       (d + copysign(sqrt(d*d + tnnm1*tnnm1), d));
69952bfb9bbSJeremy L Thompson     CeedScalar x = matT[p+n*p] - mu;
70052bfb9bbSJeremy L Thompson     CeedScalar z = matT[p+n*(p+1)];
70152bfb9bbSJeremy L Thompson     for (CeedInt k=p; k<n-1-q; k++) {
70252bfb9bbSJeremy L Thompson       // Compute Givens rotation
70352bfb9bbSJeremy L Thompson       CeedScalar c = 1, s = 0;
70452bfb9bbSJeremy L Thompson       if (fabs(z) > tol) {
70552bfb9bbSJeremy L Thompson         if (fabs(z) > fabs(x)) {
70652bfb9bbSJeremy L Thompson           CeedScalar tau = -x/z;
70752bfb9bbSJeremy L Thompson           s = 1/sqrt(1+tau*tau), c = s*tau;
70852bfb9bbSJeremy L Thompson         } else {
70952bfb9bbSJeremy L Thompson           CeedScalar tau = -z/x;
71052bfb9bbSJeremy L Thompson           c = 1/sqrt(1+tau*tau), s = c*tau;
71152bfb9bbSJeremy L Thompson         }
71252bfb9bbSJeremy L Thompson       }
71352bfb9bbSJeremy L Thompson 
71452bfb9bbSJeremy L Thompson       // Apply Givens rotation to T
71552bfb9bbSJeremy L Thompson       CeedGivensRotation(matT, c, s, CEED_NOTRANSPOSE, k, k+1, n, n);
71652bfb9bbSJeremy L Thompson       CeedGivensRotation(matT, c, s, CEED_TRANSPOSE, k, k+1, n, n);
71752bfb9bbSJeremy L Thompson 
71852bfb9bbSJeremy L Thompson       // Apply Givens rotation to Q
71952bfb9bbSJeremy L Thompson       CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k+1, n, n);
72052bfb9bbSJeremy L Thompson 
72152bfb9bbSJeremy L Thompson       // Update x, z
72252bfb9bbSJeremy L Thompson       if (k < n-q-2) {
72352bfb9bbSJeremy L Thompson         x = matT[k+n*(k+1)];
72452bfb9bbSJeremy L Thompson         z = matT[k+n*(k+2)];
72552bfb9bbSJeremy L Thompson       }
72652bfb9bbSJeremy L Thompson     }
72752bfb9bbSJeremy L Thompson     itr++;
72852bfb9bbSJeremy L Thompson   }
72952bfb9bbSJeremy L Thompson   // Save eigenvalues
73052bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
73152bfb9bbSJeremy L Thompson     lambda[i] = matT[i+n*i];
73252bfb9bbSJeremy L Thompson 
73352bfb9bbSJeremy L Thompson   // Check convergence
73452bfb9bbSJeremy L Thompson   if (itr == maxitr && q < n-1)
735c042f62fSJeremy L Thompson     // LCOV_EXCL_START
73652bfb9bbSJeremy L Thompson     return CeedError(ceed, 1, "Symmetric QR failed to converge");
737c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
73852bfb9bbSJeremy L Thompson 
73952bfb9bbSJeremy L Thompson   return 0;
74052bfb9bbSJeremy L Thompson }
74152bfb9bbSJeremy L Thompson 
74252bfb9bbSJeremy L Thompson /**
74352bfb9bbSJeremy L Thompson   @brief Return C = A B
74452bfb9bbSJeremy L Thompson 
74552bfb9bbSJeremy L Thompson   @param[in] matA     Row-major matrix A
74652bfb9bbSJeremy L Thompson   @param[in] matB     Row-major matrix B
74752bfb9bbSJeremy L Thompson   @param[out] matC    Row-major output matrix C
74852bfb9bbSJeremy L Thompson   @param m            Number of rows of C
74952bfb9bbSJeremy L Thompson   @param n            Number of columns of C
75052bfb9bbSJeremy L Thompson   @param kk           Number of columns of A/rows of B
75152bfb9bbSJeremy L Thompson 
75252bfb9bbSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
75352bfb9bbSJeremy L Thompson 
75452bfb9bbSJeremy L Thompson   @ref Utility
75552bfb9bbSJeremy L Thompson **/
75652bfb9bbSJeremy L Thompson static int CeedMatrixMultiply(Ceed ceed, CeedScalar *matA, CeedScalar *matB,
75752bfb9bbSJeremy L Thompson                               CeedScalar *matC, CeedInt m, CeedInt n,
75852bfb9bbSJeremy L Thompson                               CeedInt kk) {
75952bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<m; i++)
76052bfb9bbSJeremy L Thompson     for (CeedInt j=0; j<n; j++) {
76152bfb9bbSJeremy L Thompson       CeedScalar sum = 0;
76252bfb9bbSJeremy L Thompson       for (CeedInt k=0; k<kk; k++)
76352bfb9bbSJeremy L Thompson         sum += matA[k+i*kk]*matB[j+k*n];
76452bfb9bbSJeremy L Thompson       matC[j+i*n] = sum;
76552bfb9bbSJeremy L Thompson     }
76652bfb9bbSJeremy L Thompson   return 0;
76752bfb9bbSJeremy L Thompson }
76852bfb9bbSJeremy L Thompson 
76952bfb9bbSJeremy L Thompson /**
77052bfb9bbSJeremy L Thompson   @brief Return Simultaneous Diagonalization of two matrices. This solves the
77152bfb9bbSJeremy L Thompson            generalized eigenvalue problem A x = lambda B x, where A and B
77252bfb9bbSJeremy L Thompson            are symmetric and B is positive definite. We generate the matrix X
77352bfb9bbSJeremy L Thompson            and vector Lambda such that X^T A X = Lambda and X^T B X = I. This
77452bfb9bbSJeremy L Thompson            is equivalent to the LAPACK routine 'sygv' with TYPE = 1.
77552bfb9bbSJeremy L Thompson 
77652bfb9bbSJeremy L Thompson   @param[in] matA     Row-major matrix to be factorized with eigenvalues
77752bfb9bbSJeremy L Thompson   @param[in] matB     Row-major matrix to be factorized to identity
77852bfb9bbSJeremy L Thompson   @param[out] x       Row-major orthogonal matrix
77952bfb9bbSJeremy L Thompson   @param[out] lambda  Vector of length m of generalized eigenvalues
78052bfb9bbSJeremy L Thompson   @param n            Number of rows/columns
78152bfb9bbSJeremy L Thompson 
78252bfb9bbSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
78352bfb9bbSJeremy L Thompson 
78452bfb9bbSJeremy L Thompson   @ref Utility
78552bfb9bbSJeremy L Thompson **/
78652bfb9bbSJeremy L Thompson int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *matA,
78752bfb9bbSJeremy L Thompson                                     CeedScalar *matB, CeedScalar *x,
78852bfb9bbSJeremy L Thompson                                     CeedScalar *lambda, CeedInt n) {
78952bfb9bbSJeremy L Thompson   int ierr;
79052bfb9bbSJeremy L Thompson   CeedScalar matC[n*n], matG[n*n], vecD[n];
79152bfb9bbSJeremy L Thompson 
79252bfb9bbSJeremy L Thompson   // Compute B = G D G^T
79352bfb9bbSJeremy L Thompson   memcpy(matG, matB, n*n*sizeof(matB[0]));
79452bfb9bbSJeremy L Thompson   ierr = CeedSymmetricSchurDecomposition(ceed, matG, vecD, n); CeedChk(ierr);
79552bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++) vecD[i] = sqrt(vecD[i]);
79652bfb9bbSJeremy L Thompson 
79752bfb9bbSJeremy L Thompson   // Compute C = (G D^-1/2)^-1 A (G D^-1/2)^-T
79852bfb9bbSJeremy L Thompson   //           = D^1/2 G^T A D^1/2 G
79952bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
80052bfb9bbSJeremy L Thompson     for (CeedInt j=0; j<n; j++)
80152bfb9bbSJeremy L Thompson       matC[j+i*n] = vecD[i] * matG[i+j*n];
80252bfb9bbSJeremy L Thompson   CeedMatrixMultiply(ceed, matC, matA, x, n, n, n);
80352bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
80452bfb9bbSJeremy L Thompson     for (CeedInt j=0; j<n; j++)
80552bfb9bbSJeremy L Thompson       matG[j+i*n] = vecD[i] * matG[j+i*n];
80652bfb9bbSJeremy L Thompson   CeedMatrixMultiply(ceed, x, matG, matC, n, n, n);
80752bfb9bbSJeremy L Thompson 
80852bfb9bbSJeremy L Thompson   // Compute Q^T C Q = lambda
80952bfb9bbSJeremy L Thompson   ierr = CeedSymmetricSchurDecomposition(ceed, matC, lambda, n); CeedChk(ierr);
81052bfb9bbSJeremy L Thompson 
81152bfb9bbSJeremy L Thompson   // Set x = (G D^-1/2)^-T Q
81252bfb9bbSJeremy L Thompson   //       = D^1/2 G Q
81352bfb9bbSJeremy L Thompson   CeedMatrixMultiply(ceed, matG, matC, x, n, n, n);
81452bfb9bbSJeremy L Thompson 
81552bfb9bbSJeremy L Thompson   return 0;
81652bfb9bbSJeremy L Thompson }
81752bfb9bbSJeremy L Thompson 
81852bfb9bbSJeremy L Thompson /**
819783c99b3SValeria Barra   @brief Return collocated grad matrix
820b11c1e72Sjeremylt 
821b11c1e72Sjeremylt   @param basis           CeedBasis
822b11c1e72Sjeremylt   @param[out] colograd1d Row-major Q1d × Q1d matrix expressing derivatives of
823b11c1e72Sjeremylt                            basis functions at quadrature points
824b11c1e72Sjeremylt 
825b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
826dfdf5a53Sjeremylt 
827dfdf5a53Sjeremylt   @ref Advanced
828b11c1e72Sjeremylt **/
829783c99b3SValeria Barra int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *colograd1d) {
830d7b241e6Sjeremylt   int i, j, k;
831a7bd39daSjeremylt   Ceed ceed;
832d7b241e6Sjeremylt   CeedInt ierr, P1d=(basis)->P1d, Q1d=(basis)->Q1d;
833d7b241e6Sjeremylt   CeedScalar *interp1d, *grad1d, tau[Q1d];
834d7b241e6Sjeremylt 
835d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d*P1d, &interp1d); CeedChk(ierr);
836d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d*P1d, &grad1d); CeedChk(ierr);
837d7b241e6Sjeremylt   memcpy(interp1d, (basis)->interp1d, Q1d*P1d*sizeof(basis)->interp1d[0]);
838d7b241e6Sjeremylt   memcpy(grad1d, (basis)->grad1d, Q1d*P1d*sizeof(basis)->interp1d[0]);
839d7b241e6Sjeremylt 
840d7b241e6Sjeremylt   // QR Factorization, interp1d = Q R
841a7bd39daSjeremylt   ierr = CeedBasisGetCeed(basis, &ceed); CeedChk(ierr);
842a7bd39daSjeremylt   ierr = CeedQRFactorization(ceed, interp1d, tau, Q1d, P1d); CeedChk(ierr);
843d7b241e6Sjeremylt 
844d7b241e6Sjeremylt   // Apply Rinv, colograd1d = grad1d Rinv
845d7b241e6Sjeremylt   for (i=0; i<Q1d; i++) { // Row i
846d7b241e6Sjeremylt     colograd1d[Q1d*i] = grad1d[P1d*i]/interp1d[0];
847d7b241e6Sjeremylt     for (j=1; j<P1d; j++) { // Column j
848d7b241e6Sjeremylt       colograd1d[j+Q1d*i] = grad1d[j+P1d*i];
849*1d102b48SJeremy L Thompson       for (k=0; k<j; k++)
850d7b241e6Sjeremylt         colograd1d[j+Q1d*i] -= interp1d[j+P1d*k]*colograd1d[k+Q1d*i];
851d7b241e6Sjeremylt       colograd1d[j+Q1d*i] /= interp1d[j+P1d*j];
852d7b241e6Sjeremylt     }
853*1d102b48SJeremy L Thompson     for (j=P1d; j<Q1d; j++)
854d7b241e6Sjeremylt       colograd1d[j+Q1d*i] = 0;
855d7b241e6Sjeremylt   }
856d7b241e6Sjeremylt 
857d7b241e6Sjeremylt   // Apply Qtranspose, colograd = colograd Qtranspose
858d7b241e6Sjeremylt   CeedHouseholderApplyQ(colograd1d, interp1d, tau, CEED_NOTRANSPOSE,
859d7b241e6Sjeremylt                         Q1d, Q1d, P1d, 1, Q1d);
860d7b241e6Sjeremylt 
861d7b241e6Sjeremylt   ierr = CeedFree(&interp1d); CeedChk(ierr);
862d7b241e6Sjeremylt   ierr = CeedFree(&grad1d); CeedChk(ierr);
863d7b241e6Sjeremylt 
864d7b241e6Sjeremylt   return 0;
865d7b241e6Sjeremylt }
866d7b241e6Sjeremylt 
867b11c1e72Sjeremylt /**
868b11c1e72Sjeremylt   @brief Apply basis evaluation from nodes to quadrature points or vice-versa
869b11c1e72Sjeremylt 
870b11c1e72Sjeremylt   @param basis  CeedBasis to evaluate
871b11c1e72Sjeremylt   @param nelem  The number of elements to apply the basis evaluation to;
872b11c1e72Sjeremylt                   the backend will specify the ordering in
873b11c1e72Sjeremylt                   ElemRestrictionCreateBlocked
874b11c1e72Sjeremylt   @param tmode  \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature
875b11c1e72Sjeremylt                   points, \ref CEED_TRANSPOSE to apply the transpose, mapping
876b11c1e72Sjeremylt                   from quadrature points to nodes
877b11c1e72Sjeremylt   @param emode  \ref CEED_EVAL_INTERP to obtain interpolated values,
878b11c1e72Sjeremylt                   \ref CEED_EVAL_GRAD to obtain gradients.
879b11c1e72Sjeremylt   @param[in] u  Input array
880b11c1e72Sjeremylt   @param[out] v Output array
881b11c1e72Sjeremylt 
882b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
883dfdf5a53Sjeremylt 
884dfdf5a53Sjeremylt   @ref Advanced
885b11c1e72Sjeremylt **/
886d7b241e6Sjeremylt int CeedBasisApply(CeedBasis basis, CeedInt nelem, CeedTransposeMode tmode,
887aedaa0e5Sjeremylt                    CeedEvalMode emode, CeedVector u, CeedVector v) {
888d7b241e6Sjeremylt   int ierr;
8898795c945Sjeremylt   CeedInt ulength = 0, vlength, nnodes, nqpt;
890c042f62fSJeremy L Thompson   if (!basis->Apply)
891c042f62fSJeremy L Thompson     // LCOV_EXCL_START
892c042f62fSJeremy L Thompson     return CeedError(basis->ceed, 1, "Backend does not support BasisApply");
893c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
894c042f62fSJeremy L Thompson 
895c042f62fSJeremy L Thompson   // Check compatibility of topological and geometrical dimensions
8968795c945Sjeremylt   ierr = CeedBasisGetNumNodes(basis, &nnodes); CeedChk(ierr);
897b502e64cSValeria Barra   ierr = CeedBasisGetNumQuadraturePoints(basis, &nqpt); CeedChk(ierr);
898b502e64cSValeria Barra   ierr = CeedVectorGetLength(v, &vlength); CeedChk(ierr);
899b502e64cSValeria Barra 
900b502e64cSValeria Barra   if (u) {
901b502e64cSValeria Barra     ierr = CeedVectorGetLength(u, &ulength); CeedChk(ierr);
902b502e64cSValeria Barra   }
903b502e64cSValeria Barra 
904*1d102b48SJeremy L Thompson   if ((tmode == CEED_TRANSPOSE && (vlength%nnodes != 0 || ulength%nqpt != 0)) ||
9058795c945Sjeremylt       (tmode == CEED_NOTRANSPOSE && (ulength%nnodes != 0 || vlength%nqpt != 0)))
906*1d102b48SJeremy L Thompson     return CeedError(basis->ceed, 1, "Length of input/output vectors "
907*1d102b48SJeremy L Thompson                      "incompatible with basis dimensions");
908b502e64cSValeria Barra 
909d7b241e6Sjeremylt   ierr = basis->Apply(basis, nelem, tmode, emode, u, v); CeedChk(ierr);
910d7b241e6Sjeremylt   return 0;
911d7b241e6Sjeremylt }
912d7b241e6Sjeremylt 
913b11c1e72Sjeremylt /**
9144ce2993fSjeremylt   @brief Get Ceed associated with a CeedBasis
915b11c1e72Sjeremylt 
916b11c1e72Sjeremylt   @param basis      CeedBasis
9174ce2993fSjeremylt   @param[out] ceed  Variable to store Ceed
9184ce2993fSjeremylt 
9194ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
9204ce2993fSjeremylt 
92123617272Sjeremylt   @ref Advanced
9224ce2993fSjeremylt **/
9234ce2993fSjeremylt int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) {
9244ce2993fSjeremylt   *ceed = basis->ceed;
9254ce2993fSjeremylt   return 0;
9264ce2993fSjeremylt };
9274ce2993fSjeremylt 
9284ce2993fSjeremylt /**
9294ce2993fSjeremylt   @brief Get dimension for given CeedBasis
9304ce2993fSjeremylt 
9314ce2993fSjeremylt   @param basis     CeedBasis
9324ce2993fSjeremylt   @param[out] dim  Variable to store dimension of basis
9334ce2993fSjeremylt 
9344ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
9354ce2993fSjeremylt 
93623617272Sjeremylt   @ref Advanced
9374ce2993fSjeremylt **/
9384ce2993fSjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) {
9394ce2993fSjeremylt   *dim = basis->dim;
9404ce2993fSjeremylt   return 0;
9414ce2993fSjeremylt };
9424ce2993fSjeremylt 
9434ce2993fSjeremylt /**
9444ce2993fSjeremylt   @brief Get tensor status for given CeedBasis
9454ce2993fSjeremylt 
9464ce2993fSjeremylt   @param basis        CeedBasis
9474ce2993fSjeremylt   @param[out] tensor  Variable to store tensor status
9484ce2993fSjeremylt 
9494ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
9504ce2993fSjeremylt 
95123617272Sjeremylt   @ref Advanced
9524ce2993fSjeremylt **/
9534ce2993fSjeremylt int CeedBasisGetTensorStatus(CeedBasis basis, bool *tensor) {
9544ce2993fSjeremylt   *tensor = basis->tensorbasis;
9554ce2993fSjeremylt   return 0;
9564ce2993fSjeremylt };
9574ce2993fSjeremylt 
9584ce2993fSjeremylt /**
9594ce2993fSjeremylt   @brief Get number of components for given CeedBasis
9604ce2993fSjeremylt 
9614ce2993fSjeremylt   @param basis        CeedBasis
962288c0443SJeremy L Thompson   @param[out] numcomp Variable to store number of components of basis
9634ce2993fSjeremylt 
9644ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
9654ce2993fSjeremylt 
96623617272Sjeremylt   @ref Advanced
9674ce2993fSjeremylt **/
9684ce2993fSjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *numcomp) {
9694ce2993fSjeremylt   *numcomp = basis->ncomp;
9704ce2993fSjeremylt   return 0;
9714ce2993fSjeremylt };
9724ce2993fSjeremylt 
9734ce2993fSjeremylt /**
9744ce2993fSjeremylt   @brief Get total number of nodes (in 1 dimension) of a CeedBasis
9754ce2993fSjeremylt 
9764ce2993fSjeremylt   @param basis     CeedBasis
9774ce2993fSjeremylt   @param[out] P1d  Variable to store number of nodes
9784ce2993fSjeremylt 
9794ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
9804ce2993fSjeremylt 
98123617272Sjeremylt   @ref Advanced
9824ce2993fSjeremylt **/
9834ce2993fSjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P1d) {
984c042f62fSJeremy L Thompson   if (!basis->tensorbasis)
985c042f62fSJeremy L Thompson     // LCOV_EXCL_START
986c042f62fSJeremy L Thompson     return CeedError(basis->ceed, 1, "Cannot supply P1d for non-tensor basis");
987c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
988c042f62fSJeremy L Thompson 
9894ce2993fSjeremylt   *P1d = basis->P1d;
9904ce2993fSjeremylt   return 0;
9914ce2993fSjeremylt }
9924ce2993fSjeremylt 
9934ce2993fSjeremylt /**
9944ce2993fSjeremylt   @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis
9954ce2993fSjeremylt 
9964ce2993fSjeremylt   @param basis     CeedBasis
9974ce2993fSjeremylt   @param[out] Q1d  Variable to store number of quadrature points
9984ce2993fSjeremylt 
9994ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
10004ce2993fSjeremylt 
100123617272Sjeremylt   @ref Advanced
10024ce2993fSjeremylt **/
10034ce2993fSjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q1d) {
1004c042f62fSJeremy L Thompson   if (!basis->tensorbasis)
1005c042f62fSJeremy L Thompson     // LCOV_EXCL_START
1006c042f62fSJeremy L Thompson     return CeedError(basis->ceed, 1, "Cannot supply Q1d for non-tensor basis");
1007c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
1008c042f62fSJeremy L Thompson 
10094ce2993fSjeremylt   *Q1d = basis->Q1d;
10104ce2993fSjeremylt   return 0;
10114ce2993fSjeremylt }
10124ce2993fSjeremylt 
10134ce2993fSjeremylt /**
10144ce2993fSjeremylt   @brief Get total number of nodes (in dim dimensions) of a CeedBasis
10154ce2993fSjeremylt 
10164ce2993fSjeremylt   @param basis   CeedBasis
10174ce2993fSjeremylt   @param[out] P  Variable to store number of nodes
1018b11c1e72Sjeremylt 
1019b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1020dfdf5a53Sjeremylt 
1021dfdf5a53Sjeremylt   @ref Utility
1022b11c1e72Sjeremylt **/
1023d7b241e6Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) {
1024a8de75f0Sjeremylt   *P = basis->P;
1025d7b241e6Sjeremylt   return 0;
1026d7b241e6Sjeremylt }
1027d7b241e6Sjeremylt 
1028b11c1e72Sjeremylt /**
10294ce2993fSjeremylt   @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis
1030b11c1e72Sjeremylt 
1031b11c1e72Sjeremylt   @param basis   CeedBasis
10324ce2993fSjeremylt   @param[out] Q  Variable to store number of quadrature points
1033b11c1e72Sjeremylt 
1034b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1035dfdf5a53Sjeremylt 
1036dfdf5a53Sjeremylt   @ref Utility
1037b11c1e72Sjeremylt **/
1038d7b241e6Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) {
1039a8de75f0Sjeremylt   *Q = basis->Q;
1040d7b241e6Sjeremylt   return 0;
1041d7b241e6Sjeremylt }
1042d7b241e6Sjeremylt 
1043b11c1e72Sjeremylt /**
10448c91a0c9SJeremy L Thompson   @brief Get reference coordinates of quadrature points (in dim dimensions)
10454ce2993fSjeremylt          of a CeedBasis
10464ce2993fSjeremylt 
10474ce2993fSjeremylt   @param basis      CeedBasis
10488c91a0c9SJeremy L Thompson   @param[out] qref  Variable to store reference coordinates of quadrature points
10494ce2993fSjeremylt 
10504ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
10514ce2993fSjeremylt 
105223617272Sjeremylt   @ref Advanced
10534ce2993fSjeremylt **/
10544ce2993fSjeremylt int CeedBasisGetQRef(CeedBasis basis, CeedScalar* *qref) {
10554ce2993fSjeremylt   *qref = basis->qref1d;
10564ce2993fSjeremylt   return 0;
10574ce2993fSjeremylt }
10584ce2993fSjeremylt 
10594ce2993fSjeremylt /**
10604ce2993fSjeremylt   @brief Get quadrature weights of quadrature points (in dim dimensions)
10614ce2993fSjeremylt          of a CeedBasis
10624ce2993fSjeremylt 
10634ce2993fSjeremylt   @param basis         CeedBasis
10644ce2993fSjeremylt   @param[out] qweight  Variable to store quadrature weights
10654ce2993fSjeremylt 
10664ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
10674ce2993fSjeremylt 
106823617272Sjeremylt   @ref Advanced
10694ce2993fSjeremylt **/
10704ce2993fSjeremylt int CeedBasisGetQWeights(CeedBasis basis, CeedScalar* *qweight) {
10714ce2993fSjeremylt   *qweight = basis->qweight1d;
10724ce2993fSjeremylt   return 0;
10734ce2993fSjeremylt }
10744ce2993fSjeremylt 
10754ce2993fSjeremylt /**
10764ce2993fSjeremylt   @brief Get interpolation matrix of a CeedBasis
10774ce2993fSjeremylt 
10784ce2993fSjeremylt   @param basis       CeedBasis
1079288c0443SJeremy L Thompson   @param[out] interp Variable to store interpolation matrix
10804ce2993fSjeremylt 
10814ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
10824ce2993fSjeremylt 
108323617272Sjeremylt   @ref Advanced
10844ce2993fSjeremylt **/
10854ce2993fSjeremylt int CeedBasisGetInterp(CeedBasis basis, CeedScalar* *interp) {
10864ce2993fSjeremylt   *interp = basis->interp1d;
10874ce2993fSjeremylt   return 0;
10884ce2993fSjeremylt }
10894ce2993fSjeremylt 
10904ce2993fSjeremylt /**
10914ce2993fSjeremylt   @brief Get gradient matrix of a CeedBasis
10924ce2993fSjeremylt 
10934ce2993fSjeremylt   @param basis      CeedBasis
1094288c0443SJeremy L Thompson   @param[out] grad  Variable to store gradient matrix
10954ce2993fSjeremylt 
10964ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
10974ce2993fSjeremylt 
109823617272Sjeremylt   @ref Advanced
10994ce2993fSjeremylt **/
11004ce2993fSjeremylt int CeedBasisGetGrad(CeedBasis basis, CeedScalar* *grad) {
11014ce2993fSjeremylt   *grad = basis->grad1d;
11024ce2993fSjeremylt   return 0;
11034ce2993fSjeremylt }
11044ce2993fSjeremylt 
11054ce2993fSjeremylt /**
11064ce2993fSjeremylt   @brief Get backend data of a CeedBasis
11074ce2993fSjeremylt 
11084ce2993fSjeremylt   @param basis      CeedBasis
11094ce2993fSjeremylt   @param[out] data  Variable to store data
11104ce2993fSjeremylt 
11114ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
11124ce2993fSjeremylt 
111323617272Sjeremylt   @ref Advanced
11144ce2993fSjeremylt **/
11154ce2993fSjeremylt int CeedBasisGetData(CeedBasis basis, void* *data) {
11164ce2993fSjeremylt   *data = basis->data;
11174ce2993fSjeremylt   return 0;
11184ce2993fSjeremylt }
11194ce2993fSjeremylt 
11204ce2993fSjeremylt /**
1121fe2413ffSjeremylt   @brief Set backend data of a CeedBasis
1122fe2413ffSjeremylt 
1123fe2413ffSjeremylt   @param[out] basis CeedBasis
1124fe2413ffSjeremylt   @param data       Data to set
1125fe2413ffSjeremylt 
1126fe2413ffSjeremylt   @return An error code: 0 - success, otherwise - failure
1127fe2413ffSjeremylt 
1128fe2413ffSjeremylt   @ref Advanced
1129fe2413ffSjeremylt **/
1130fe2413ffSjeremylt int CeedBasisSetData(CeedBasis basis, void* *data) {
1131fe2413ffSjeremylt   basis->data = *data;
1132fe2413ffSjeremylt   return 0;
1133fe2413ffSjeremylt }
1134fe2413ffSjeremylt 
1135fe2413ffSjeremylt /**
11362f86a920SJeremy L Thompson   @brief Get CeedTensorContract of a CeedBasis
11372f86a920SJeremy L Thompson 
11382f86a920SJeremy L Thompson   @param basis          CeedBasis
11392f86a920SJeremy L Thompson   @param[out] contract  Variable to store CeedTensorContract
11402f86a920SJeremy L Thompson 
11412f86a920SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
11422f86a920SJeremy L Thompson 
11432f86a920SJeremy L Thompson   @ref Advanced
11442f86a920SJeremy L Thompson **/
1145*1d102b48SJeremy L Thompson int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) {
11462f86a920SJeremy L Thompson   *contract = basis->contract;
11472f86a920SJeremy L Thompson   return 0;
11482f86a920SJeremy L Thompson }
11492f86a920SJeremy L Thompson 
11502f86a920SJeremy L Thompson /**
11512f86a920SJeremy L Thompson   @brief Set CeedTensorContract of a CeedBasis
11522f86a920SJeremy L Thompson 
11532f86a920SJeremy L Thompson   @param[out] basis     CeedBasis
11542f86a920SJeremy L Thompson   @param contract       CeedTensorContract to set
11552f86a920SJeremy L Thompson 
11562f86a920SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
11572f86a920SJeremy L Thompson 
11582f86a920SJeremy L Thompson   @ref Advanced
11592f86a920SJeremy L Thompson **/
1160*1d102b48SJeremy L Thompson int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract *contract) {
11612f86a920SJeremy L Thompson   basis->contract = *contract;
11622f86a920SJeremy L Thompson   return 0;
11632f86a920SJeremy L Thompson }
11642f86a920SJeremy L Thompson 
11652f86a920SJeremy L Thompson /**
1166a8de75f0Sjeremylt   @brief Get dimension for given CeedElemTopology
1167a8de75f0Sjeremylt 
1168a8de75f0Sjeremylt   @param topo      CeedElemTopology
11694ce2993fSjeremylt   @param[out] dim  Variable to store dimension of topology
1170a8de75f0Sjeremylt 
1171a8de75f0Sjeremylt   @return An error code: 0 - success, otherwise - failure
1172a8de75f0Sjeremylt 
117323617272Sjeremylt   @ref Advanced
1174a8de75f0Sjeremylt **/
1175a8de75f0Sjeremylt int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) {
1176a8de75f0Sjeremylt   *dim = (CeedInt) topo >> 16;
1177a8de75f0Sjeremylt   return 0;
1178a8de75f0Sjeremylt };
1179a8de75f0Sjeremylt 
1180a8de75f0Sjeremylt /**
1181b11c1e72Sjeremylt   @brief Destroy a CeedBasis
1182b11c1e72Sjeremylt 
1183b11c1e72Sjeremylt   @param basis CeedBasis to destroy
1184b11c1e72Sjeremylt 
1185b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1186dfdf5a53Sjeremylt 
1187dfdf5a53Sjeremylt   @ref Basic
1188b11c1e72Sjeremylt **/
1189d7b241e6Sjeremylt int CeedBasisDestroy(CeedBasis *basis) {
1190d7b241e6Sjeremylt   int ierr;
1191d7b241e6Sjeremylt 
1192*1d102b48SJeremy L Thompson   if (!*basis || --(*basis)->refcount > 0)
1193*1d102b48SJeremy L Thompson     return 0;
1194d7b241e6Sjeremylt   if ((*basis)->Destroy) {
1195d7b241e6Sjeremylt     ierr = (*basis)->Destroy(*basis); CeedChk(ierr);
1196d7b241e6Sjeremylt   }
1197d7b241e6Sjeremylt   ierr = CeedFree(&(*basis)->interp1d); CeedChk(ierr);
1198d7b241e6Sjeremylt   ierr = CeedFree(&(*basis)->grad1d); CeedChk(ierr);
1199d7b241e6Sjeremylt   ierr = CeedFree(&(*basis)->qref1d); CeedChk(ierr);
1200d7b241e6Sjeremylt   ierr = CeedFree(&(*basis)->qweight1d); CeedChk(ierr);
1201d7b241e6Sjeremylt   ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr);
1202d7b241e6Sjeremylt   ierr = CeedFree(basis); CeedChk(ierr);
1203d7b241e6Sjeremylt   return 0;
1204d7b241e6Sjeremylt }
1205d7b241e6Sjeremylt 
120633e6becaSjeremylt /// @cond DOXYGEN_SKIP
12078795c945Sjeremylt // Indicate that the quadrature points are collocated with the nodes
1208783c99b3SValeria Barra CeedBasis CEED_BASIS_COLLOCATED = &ceed_basis_collocated;
120933e6becaSjeremylt /// @endcond
1210d7b241e6Sjeremylt /// @}
1211