xref: /libCEED/rust/libceed-sys/c-src/interface/ceed-basis.c (revision 7ed52d01ec037ec8310375bd921b20a1514ec747)
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 /**
3595bb1877Svaleriabarra   @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
4295bb1877Svaleriabarra   @param interp1d    Row-major (Q1d * P1d) matrix expressing the values of nodal
43b11c1e72Sjeremylt                        basis functions at quadrature points
4495bb1877Svaleriabarra   @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 /**
10795bb1877Svaleriabarra   @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
11195bb1877Svaleriabarra   @param ncomp       Number of field components (1 for scalar fields)
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,
125692c2638Sjeremylt                                     CeedInt P, CeedInt Q, CeedQuadMode qmode,
126692c2638Sjeremylt                                     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 /**
18695bb1877Svaleriabarra   @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
19395bb1877Svaleriabarra   @param interp      Row-major (nqpts * nnodes) matrix expressing the values of
1948795c945Sjeremylt                        nodal basis functions at quadrature points
19595bb1877Svaleriabarra   @param grad        Row-major (nqpts * 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,
209692c2638Sjeremylt                       CeedInt nnodes, CeedInt nqpts, const CeedScalar *interp,
210a8de75f0Sjeremylt                       const CeedScalar *grad, const CeedScalar *qref,
211a8de75f0Sjeremylt                       const CeedScalar *qweight, CeedBasis *basis) {
212a8de75f0Sjeremylt   int ierr;
2138795c945Sjeremylt   CeedInt P = nnodes, Q = nqpts, dim = 0;
214a8de75f0Sjeremylt 
2155fe0d4faSjeremylt   if (!ceed->BasisCreateH1) {
2165fe0d4faSjeremylt     Ceed delegate;
217aefd8378Sjeremylt     ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr);
2185fe0d4faSjeremylt 
2195fe0d4faSjeremylt     if (!delegate)
220c042f62fSJeremy L Thompson       // LCOV_EXCL_START
221a8de75f0Sjeremylt       return CeedError(ceed, 1, "Backend does not support BasisCreateH1");
222c042f62fSJeremy L Thompson     // LCOV_EXCL_STOP
2235fe0d4faSjeremylt 
2248795c945Sjeremylt     ierr = CeedBasisCreateH1(delegate, topo, ncomp, nnodes,
2255fe0d4faSjeremylt                              nqpts, interp, grad, qref,
2265fe0d4faSjeremylt                              qweight, basis); CeedChk(ierr);
2275fe0d4faSjeremylt     return 0;
2285fe0d4faSjeremylt   }
2295fe0d4faSjeremylt 
230a8de75f0Sjeremylt   ierr = CeedCalloc(1,basis); CeedChk(ierr);
231a8de75f0Sjeremylt 
232a8de75f0Sjeremylt   ierr = CeedBasisGetTopologyDimension(topo, &dim); CeedChk(ierr);
233a8de75f0Sjeremylt 
234a8de75f0Sjeremylt   (*basis)->ceed = ceed;
235a8de75f0Sjeremylt   ceed->refcount++;
236a8de75f0Sjeremylt   (*basis)->refcount = 1;
237a8de75f0Sjeremylt   (*basis)->tensorbasis = 0;
238a8de75f0Sjeremylt   (*basis)->dim = dim;
239a8de75f0Sjeremylt   (*basis)->ncomp = ncomp;
240a8de75f0Sjeremylt   (*basis)->P = P;
241a8de75f0Sjeremylt   (*basis)->Q = Q;
242a8de75f0Sjeremylt   ierr = CeedMalloc(Q*dim,&(*basis)->qref1d); CeedChk(ierr);
243a8de75f0Sjeremylt   ierr = CeedMalloc(Q,&(*basis)->qweight1d); CeedChk(ierr);
244a8de75f0Sjeremylt   memcpy((*basis)->qref1d, qref, Q*dim*sizeof(qref[0]));
245a8de75f0Sjeremylt   memcpy((*basis)->qweight1d, qweight, Q*sizeof(qweight[0]));
24600f91b2bSjeremylt   ierr = CeedMalloc(Q*P, &(*basis)->interp); CeedChk(ierr);
24700f91b2bSjeremylt   ierr = CeedMalloc(dim*Q*P, &(*basis)->grad); CeedChk(ierr);
24800f91b2bSjeremylt   memcpy((*basis)->interp, interp, Q*P*sizeof(interp[0]));
24900f91b2bSjeremylt   memcpy((*basis)->grad, grad, dim*Q*P*sizeof(grad[0]));
250667bc5fcSjeremylt   ierr = ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, qref,
251a8de75f0Sjeremylt                              qweight, *basis); CeedChk(ierr);
252a8de75f0Sjeremylt   return 0;
253a8de75f0Sjeremylt }
254a8de75f0Sjeremylt 
255a8de75f0Sjeremylt /**
256b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre quadrature
257b11c1e72Sjeremylt 
258b11c1e72Sjeremylt   @param Q               Number of quadrature points (integrates polynomials of
259b11c1e72Sjeremylt                            degree 2*Q-1 exactly)
260b11c1e72Sjeremylt   @param[out] qref1d     Array of length Q to hold the abscissa on [-1, 1]
261b11c1e72Sjeremylt   @param[out] qweight1d  Array of length Q to hold the weights
262b11c1e72Sjeremylt 
263b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
264dfdf5a53Sjeremylt 
265dfdf5a53Sjeremylt   @ref Utility
266b11c1e72Sjeremylt **/
267d7b241e6Sjeremylt int CeedGaussQuadrature(CeedInt Q, CeedScalar *qref1d, CeedScalar *qweight1d) {
268d7b241e6Sjeremylt   // Allocate
269d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0*atan(1.0);
270d7b241e6Sjeremylt   // Build qref1d, qweight1d
271d7b241e6Sjeremylt   for (int i = 0; i <= Q/2; i++) {
272d7b241e6Sjeremylt     // Guess
273d7b241e6Sjeremylt     xi = cos(PI*(CeedScalar)(2*i+1)/((CeedScalar)(2*Q)));
274d7b241e6Sjeremylt     // Pn(xi)
275d7b241e6Sjeremylt     P0 = 1.0;
276d7b241e6Sjeremylt     P1 = xi;
277d7b241e6Sjeremylt     P2 = 0.0;
278d7b241e6Sjeremylt     for (int j = 2; j <= Q; j++) {
279d7b241e6Sjeremylt       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
280d7b241e6Sjeremylt       P0 = P1;
281d7b241e6Sjeremylt       P1 = P2;
282d7b241e6Sjeremylt     }
283d7b241e6Sjeremylt     // First Newton Step
284d7b241e6Sjeremylt     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
285d7b241e6Sjeremylt     xi = xi-P2/dP2;
286d7b241e6Sjeremylt     // Newton to convergence
2870e4d4210Sjeremylt     for (int k=0; k<100 && fabs(P2)>10*CEED_EPSILON; k++) {
288d7b241e6Sjeremylt       P0 = 1.0;
289d7b241e6Sjeremylt       P1 = xi;
290d7b241e6Sjeremylt       for (int j = 2; j <= Q; j++) {
291d7b241e6Sjeremylt         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
292d7b241e6Sjeremylt         P0 = P1;
293d7b241e6Sjeremylt         P1 = P2;
294d7b241e6Sjeremylt       }
295d7b241e6Sjeremylt       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
296d7b241e6Sjeremylt       xi = xi-P2/dP2;
297d7b241e6Sjeremylt     }
298d7b241e6Sjeremylt     // Save xi, wi
299d7b241e6Sjeremylt     wi = 2.0/((1.0-xi*xi)*dP2*dP2);
300d7b241e6Sjeremylt     qweight1d[i] = wi;
301d7b241e6Sjeremylt     qweight1d[Q-1-i] = wi;
302d7b241e6Sjeremylt     qref1d[i] = -xi;
303d7b241e6Sjeremylt     qref1d[Q-1-i]= xi;
304d7b241e6Sjeremylt   }
305d7b241e6Sjeremylt   return 0;
306d7b241e6Sjeremylt }
307d7b241e6Sjeremylt 
308b11c1e72Sjeremylt /**
309b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre-Lobatto quadrature
310b11c1e72Sjeremylt 
311b11c1e72Sjeremylt   @param Q               Number of quadrature points (integrates polynomials of
312b11c1e72Sjeremylt                            degree 2*Q-3 exactly)
313b11c1e72Sjeremylt   @param[out] qref1d     Array of length Q to hold the abscissa on [-1, 1]
314b11c1e72Sjeremylt   @param[out] qweight1d  Array of length Q to hold the weights
315b11c1e72Sjeremylt 
316b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
317dfdf5a53Sjeremylt 
318dfdf5a53Sjeremylt   @ref Utility
319b11c1e72Sjeremylt **/
320d7b241e6Sjeremylt int CeedLobattoQuadrature(CeedInt Q, CeedScalar *qref1d,
321d7b241e6Sjeremylt                           CeedScalar *qweight1d) {
322d7b241e6Sjeremylt   // Allocate
323d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0*atan(1.0);
324d7b241e6Sjeremylt   // Build qref1d, qweight1d
325d7b241e6Sjeremylt   // Set endpoints
32630a100c3SJed Brown   if (Q < 2)
327*7ed52d01Sjeremylt     return CeedError(NULL, 1,
328*7ed52d01Sjeremylt                      "Cannot create Lobatto quadrature with Q=%d < 2 points", Q);
329d7b241e6Sjeremylt   wi = 2.0/((CeedScalar)(Q*(Q-1)));
330d7b241e6Sjeremylt   if (qweight1d) {
331d7b241e6Sjeremylt     qweight1d[0] = wi;
332d7b241e6Sjeremylt     qweight1d[Q-1] = wi;
333d7b241e6Sjeremylt   }
334d7b241e6Sjeremylt   qref1d[0] = -1.0;
335d7b241e6Sjeremylt   qref1d[Q-1] = 1.0;
336d7b241e6Sjeremylt   // Interior
337d7b241e6Sjeremylt   for (int i = 1; i <= (Q-1)/2; i++) {
338d7b241e6Sjeremylt     // Guess
339d7b241e6Sjeremylt     xi = cos(PI*(CeedScalar)(i)/(CeedScalar)(Q-1));
340d7b241e6Sjeremylt     // Pn(xi)
341d7b241e6Sjeremylt     P0 = 1.0;
342d7b241e6Sjeremylt     P1 = xi;
343d7b241e6Sjeremylt     P2 = 0.0;
344d7b241e6Sjeremylt     for (int j = 2; j < Q; j++) {
345d7b241e6Sjeremylt       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
346d7b241e6Sjeremylt       P0 = P1;
347d7b241e6Sjeremylt       P1 = P2;
348d7b241e6Sjeremylt     }
349d7b241e6Sjeremylt     // First Newton step
350d7b241e6Sjeremylt     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
351d7b241e6Sjeremylt     d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
352d7b241e6Sjeremylt     xi = xi-dP2/d2P2;
353d7b241e6Sjeremylt     // Newton to convergence
3540e4d4210Sjeremylt     for (int k=0; k<100 && fabs(dP2)>10*CEED_EPSILON; k++) {
355d7b241e6Sjeremylt       P0 = 1.0;
356d7b241e6Sjeremylt       P1 = xi;
357d7b241e6Sjeremylt       for (int j = 2; j < Q; j++) {
358d7b241e6Sjeremylt         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
359d7b241e6Sjeremylt         P0 = P1;
360d7b241e6Sjeremylt         P1 = P2;
361d7b241e6Sjeremylt       }
362d7b241e6Sjeremylt       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
363d7b241e6Sjeremylt       d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
364d7b241e6Sjeremylt       xi = xi-dP2/d2P2;
365d7b241e6Sjeremylt     }
366d7b241e6Sjeremylt     // Save xi, wi
367d7b241e6Sjeremylt     wi = 2.0/(((CeedScalar)(Q*(Q-1)))*P2*P2);
368d7b241e6Sjeremylt     if (qweight1d) {
369d7b241e6Sjeremylt       qweight1d[i] = wi;
370d7b241e6Sjeremylt       qweight1d[Q-1-i] = wi;
371d7b241e6Sjeremylt     }
372d7b241e6Sjeremylt     qref1d[i] = -xi;
373d7b241e6Sjeremylt     qref1d[Q-1-i]= xi;
374d7b241e6Sjeremylt   }
375d7b241e6Sjeremylt   return 0;
376d7b241e6Sjeremylt }
377d7b241e6Sjeremylt 
378dfdf5a53Sjeremylt /**
379dfdf5a53Sjeremylt   @brief View an array stored in a CeedBasis
380dfdf5a53Sjeremylt 
381dfdf5a53Sjeremylt   @param name      Name of array
382dfdf5a53Sjeremylt   @param fpformat  Printing format
383dfdf5a53Sjeremylt   @param m         Number of rows in array
384dfdf5a53Sjeremylt   @param n         Number of columns in array
385dfdf5a53Sjeremylt   @param a         Array to be viewed
386dfdf5a53Sjeremylt   @param stream    Stream to view to, e.g., stdout
387dfdf5a53Sjeremylt 
388dfdf5a53Sjeremylt   @return An error code: 0 - success, otherwise - failure
389dfdf5a53Sjeremylt 
390dfdf5a53Sjeremylt   @ref Utility
391dfdf5a53Sjeremylt **/
392d7b241e6Sjeremylt static int CeedScalarView(const char *name, const char *fpformat, CeedInt m,
393d7b241e6Sjeremylt                           CeedInt n, const CeedScalar *a, FILE *stream) {
394d7b241e6Sjeremylt   for (int i=0; i<m; i++) {
3951d102b48SJeremy L Thompson     if (m > 1)
3961d102b48SJeremy L Thompson       fprintf(stream, "%12s[%d]:", name, i);
3971d102b48SJeremy L Thompson     else
3981d102b48SJeremy L Thompson       fprintf(stream, "%12s:", name);
3991d102b48SJeremy L Thompson     for (int j=0; j<n; j++)
400d7b241e6Sjeremylt       fprintf(stream, fpformat, fabs(a[i*n+j]) > 1E-14 ? a[i*n+j] : 0);
401d7b241e6Sjeremylt     fputs("\n", stream);
402d7b241e6Sjeremylt   }
403d7b241e6Sjeremylt   return 0;
404d7b241e6Sjeremylt }
405d7b241e6Sjeremylt 
406b11c1e72Sjeremylt /**
407b11c1e72Sjeremylt   @brief View a CeedBasis
408b11c1e72Sjeremylt 
409b11c1e72Sjeremylt   @param basis   CeedBasis to view
410b11c1e72Sjeremylt   @param stream  Stream to view to, e.g., stdout
411b11c1e72Sjeremylt 
412b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
413dfdf5a53Sjeremylt 
414dfdf5a53Sjeremylt   @ref Utility
415b11c1e72Sjeremylt **/
416d7b241e6Sjeremylt int CeedBasisView(CeedBasis basis, FILE *stream) {
417d7b241e6Sjeremylt   int ierr;
418d7b241e6Sjeremylt 
419a8de75f0Sjeremylt   if (basis->tensorbasis) {
420d7b241e6Sjeremylt     fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P1d,
421d7b241e6Sjeremylt             basis->Q1d);
422d7b241e6Sjeremylt     ierr = CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q1d, basis->qref1d,
423d7b241e6Sjeremylt                           stream); CeedChk(ierr);
4248795c945Sjeremylt     ierr = CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q1d,
4258795c945Sjeremylt                           basis->qweight1d, stream); CeedChk(ierr);
426d7b241e6Sjeremylt     ierr = CeedScalarView("interp1d", "\t% 12.8f", basis->Q1d, basis->P1d,
427d7b241e6Sjeremylt                           basis->interp1d, stream); CeedChk(ierr);
428d7b241e6Sjeremylt     ierr = CeedScalarView("grad1d", "\t% 12.8f", basis->Q1d, basis->P1d,
429d7b241e6Sjeremylt                           basis->grad1d, stream); CeedChk(ierr);
430a8de75f0Sjeremylt   } else {
431a8de75f0Sjeremylt     fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P,
432a8de75f0Sjeremylt             basis->Q);
433a8de75f0Sjeremylt     ierr = CeedScalarView("qref", "\t% 12.8f", 1, basis->Q*basis->dim,
434a8de75f0Sjeremylt                           basis->qref1d,
435a8de75f0Sjeremylt                           stream); CeedChk(ierr);
436a8de75f0Sjeremylt     ierr = CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->qweight1d,
437a8de75f0Sjeremylt                           stream); CeedChk(ierr);
438a8de75f0Sjeremylt     ierr = CeedScalarView("interp", "\t% 12.8f", basis->Q, basis->P,
43900f91b2bSjeremylt                           basis->interp, stream); CeedChk(ierr);
440a8de75f0Sjeremylt     ierr = CeedScalarView("grad", "\t% 12.8f", basis->dim*basis->Q, basis->P,
44100f91b2bSjeremylt                           basis->grad, stream); CeedChk(ierr);
442a8de75f0Sjeremylt   }
443d7b241e6Sjeremylt   return 0;
444d7b241e6Sjeremylt }
445d7b241e6Sjeremylt 
446dfdf5a53Sjeremylt /**
44752bfb9bbSJeremy L Thompson   @brief Compute Householder reflection
448dfdf5a53Sjeremylt 
449dfdf5a53Sjeremylt     Computes A = (I - b v v^T) A
450dfdf5a53Sjeremylt     where A is an mxn matrix indexed as A[i*row + j*col]
451dfdf5a53Sjeremylt 
45252bfb9bbSJeremy L Thompson   @param[in,out] A  Matrix to apply Householder reflection to, in place
453dfdf5a53Sjeremylt   @param v          Householder vector
454dfdf5a53Sjeremylt   @param b          Scaling factor
455dfdf5a53Sjeremylt   @param m          Number of rows in A
456dfdf5a53Sjeremylt   @param n          Number of columns in A
45752bfb9bbSJeremy L Thompson   @param row        Row stride
45852bfb9bbSJeremy L Thompson   @param col        Col stride
459dfdf5a53Sjeremylt 
460dfdf5a53Sjeremylt   @return An error code: 0 - success, otherwise - failure
461dfdf5a53Sjeremylt 
462dfdf5a53Sjeremylt   @ref Developer
463dfdf5a53Sjeremylt **/
464d7b241e6Sjeremylt static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v,
465d7b241e6Sjeremylt                                   CeedScalar b, CeedInt m, CeedInt n,
466d7b241e6Sjeremylt                                   CeedInt row, CeedInt col) {
467d7b241e6Sjeremylt   for (CeedInt j=0; j<n; j++) {
468d7b241e6Sjeremylt     CeedScalar w = A[0*row + j*col];
4691d102b48SJeremy L Thompson     for (CeedInt i=1; i<m; i++)
4701d102b48SJeremy L Thompson       w += v[i] * A[i*row + j*col];
471d7b241e6Sjeremylt     A[0*row + j*col] -= b * w;
4721d102b48SJeremy L Thompson     for (CeedInt i=1; i<m; i++)
4731d102b48SJeremy L Thompson       A[i*row + j*col] -= b * w * v[i];
474d7b241e6Sjeremylt   }
475d7b241e6Sjeremylt   return 0;
476d7b241e6Sjeremylt }
477d7b241e6Sjeremylt 
478dfdf5a53Sjeremylt /**
479dfdf5a53Sjeremylt   @brief Apply Householder Q matrix
480dfdf5a53Sjeremylt 
48152bfb9bbSJeremy L Thompson     Compute A = Q A where Q is mxm and A is mxn.
482dfdf5a53Sjeremylt 
48352bfb9bbSJeremy L Thompson   @param[in,out] A  Matrix to apply Householder Q to, in place
484dfdf5a53Sjeremylt   @param Q          Householder Q matrix
485dfdf5a53Sjeremylt   @param tau        Householder scaling factors
486dfdf5a53Sjeremylt   @param tmode      Transpose mode for application
487dfdf5a53Sjeremylt   @param m          Number of rows in A
488dfdf5a53Sjeremylt   @param n          Number of columns in A
48952bfb9bbSJeremy L Thompson   @param k          Number of elementary reflectors in Q, k<m
49052bfb9bbSJeremy L Thompson   @param row        Row stride in A
49152bfb9bbSJeremy L Thompson   @param col        Col stride in A
492dfdf5a53Sjeremylt 
493dfdf5a53Sjeremylt   @return An error code: 0 - success, otherwise - failure
494dfdf5a53Sjeremylt 
495dfdf5a53Sjeremylt   @ref Developer
496dfdf5a53Sjeremylt **/
497d7b241e6Sjeremylt static int CeedHouseholderApplyQ(CeedScalar *A, const CeedScalar *Q,
498d7b241e6Sjeremylt                                  const CeedScalar *tau, CeedTransposeMode tmode,
499d7b241e6Sjeremylt                                  CeedInt m, CeedInt n, CeedInt k,
500d7b241e6Sjeremylt                                  CeedInt row, CeedInt col) {
501d7b241e6Sjeremylt   CeedScalar v[m];
502d7b241e6Sjeremylt   for (CeedInt ii=0; ii<k; ii++) {
503d7b241e6Sjeremylt     CeedInt i = tmode == CEED_TRANSPOSE ? ii : k-1-ii;
50452bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<m; j++)
505d7b241e6Sjeremylt       v[j] = Q[j*k+i];
506dc7d240cSValeria Barra     // Apply Householder reflector (I - tau v v^T) collograd1d^T
507d7b241e6Sjeremylt     CeedHouseholderReflect(&A[i*row], &v[i], tau[i], m-i, n, row, col);
508d7b241e6Sjeremylt   }
509d7b241e6Sjeremylt   return 0;
510d7b241e6Sjeremylt }
511d7b241e6Sjeremylt 
512b11c1e72Sjeremylt /**
51352bfb9bbSJeremy L Thompson   @brief Compute Givens rotation
51452bfb9bbSJeremy L Thompson 
51552bfb9bbSJeremy L Thompson     Computes A = G A (or G^T A in transpose mode)
51652bfb9bbSJeremy L Thompson     where A is an mxn matrix indexed as A[i*n + j*m]
51752bfb9bbSJeremy L Thompson 
51852bfb9bbSJeremy L Thompson   @param[in,out] A  Row major matrix to apply Givens rotation to, in place
51952bfb9bbSJeremy L Thompson   @param c          Cosine factor
52052bfb9bbSJeremy L Thompson   @param s          Sine factor
52152bfb9bbSJeremy L Thompson   @param i          First row/column to apply rotation
52252bfb9bbSJeremy L Thompson   @param k          Second row/column to apply rotation
52352bfb9bbSJeremy L Thompson   @param m          Number of rows in A
52452bfb9bbSJeremy L Thompson   @param n          Number of columns in A
52552bfb9bbSJeremy L Thompson 
52652bfb9bbSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
52752bfb9bbSJeremy L Thompson 
52852bfb9bbSJeremy L Thompson   @ref Developer
52952bfb9bbSJeremy L Thompson **/
53052bfb9bbSJeremy L Thompson static int CeedGivensRotation(CeedScalar *A, CeedScalar c, CeedScalar s,
53152bfb9bbSJeremy L Thompson                               CeedTransposeMode tmode, CeedInt i, CeedInt k,
53252bfb9bbSJeremy L Thompson                               CeedInt m, CeedInt n) {
53352bfb9bbSJeremy L Thompson   CeedInt stridej = 1, strideik = m, numits = n;
53452bfb9bbSJeremy L Thompson   if (tmode == CEED_NOTRANSPOSE) {
53552bfb9bbSJeremy L Thompson     stridej = n; strideik = 1; numits = m;
53652bfb9bbSJeremy L Thompson   }
53752bfb9bbSJeremy L Thompson 
53852bfb9bbSJeremy L Thompson   // Apply rotation
53952bfb9bbSJeremy L Thompson   for (CeedInt j=0; j<numits; j++) {
54052bfb9bbSJeremy L Thompson     CeedScalar tau1 = A[i*strideik+j*stridej], tau2 = A[k*strideik+j*stridej];
54152bfb9bbSJeremy L Thompson     A[i*strideik+j*stridej] = c*tau1 - s*tau2;
54252bfb9bbSJeremy L Thompson     A[k*strideik+j*stridej] = s*tau1 + c*tau2;
54352bfb9bbSJeremy L Thompson   }
54452bfb9bbSJeremy L Thompson 
54552bfb9bbSJeremy L Thompson   return 0;
54652bfb9bbSJeremy L Thompson }
54752bfb9bbSJeremy L Thompson 
54852bfb9bbSJeremy L Thompson /**
54995bb1877Svaleriabarra   @brief Return QR Factorization of a matrix
550b11c1e72Sjeremylt 
55177645d7bSjeremylt   @param ceed         A Ceed context for error handling
55252bfb9bbSJeremy L Thompson   @param[in,out] mat  Row-major matrix to be factorized in place
55352bfb9bbSJeremy L Thompson   @param[in,out] tau  Vector of length m of scaling factors
554b11c1e72Sjeremylt   @param m            Number of rows
555b11c1e72Sjeremylt   @param n            Number of columns
556b11c1e72Sjeremylt 
557b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
558dfdf5a53Sjeremylt 
559dfdf5a53Sjeremylt   @ref Utility
560b11c1e72Sjeremylt **/
561a7bd39daSjeremylt int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau,
562d7b241e6Sjeremylt                         CeedInt m, CeedInt n) {
563d7b241e6Sjeremylt   CeedScalar v[m];
564d7b241e6Sjeremylt 
565a7bd39daSjeremylt   // Check m >= n
566a7bd39daSjeremylt   if (n > m)
567c042f62fSJeremy L Thompson     // LCOV_EXCL_START
568a7bd39daSjeremylt     return CeedError(ceed, 1, "Cannot compute QR factorization with n > m");
569c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
570a7bd39daSjeremylt 
57152bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++) {
572d7b241e6Sjeremylt     // Calculate Householder vector, magnitude
573d7b241e6Sjeremylt     CeedScalar sigma = 0.0;
574d7b241e6Sjeremylt     v[i] = mat[i+n*i];
57552bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<m; j++) {
576d7b241e6Sjeremylt       v[j] = mat[i+n*j];
577d7b241e6Sjeremylt       sigma += v[j] * v[j];
578d7b241e6Sjeremylt     }
579d7b241e6Sjeremylt     CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m]
580d7b241e6Sjeremylt     CeedScalar Rii = -copysign(norm, v[i]);
581d7b241e6Sjeremylt     v[i] -= Rii;
582d7b241e6Sjeremylt     // norm of v[i:m] after modification above and scaling below
583d7b241e6Sjeremylt     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
584d7b241e6Sjeremylt     //   tau = 2 / (norm*norm)
585d7b241e6Sjeremylt     tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma);
586fb551037Sjeremylt 
5871d102b48SJeremy L Thompson     for (CeedInt j=i+1; j<m; j++)
5881d102b48SJeremy L Thompson       v[j] /= v[i];
589d7b241e6Sjeremylt 
590d7b241e6Sjeremylt     // Apply Householder reflector to lower right panel
591d7b241e6Sjeremylt     CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1);
592d7b241e6Sjeremylt     // Save v
593d7b241e6Sjeremylt     mat[i+n*i] = Rii;
5941d102b48SJeremy L Thompson     for (CeedInt j=i+1; j<m; j++)
595d7b241e6Sjeremylt       mat[i+n*j] = v[j];
596d7b241e6Sjeremylt   }
597d7b241e6Sjeremylt 
598d7b241e6Sjeremylt   return 0;
599d7b241e6Sjeremylt }
600d7b241e6Sjeremylt 
601b11c1e72Sjeremylt /**
60252bfb9bbSJeremy L Thompson   @brief Return symmetric Schur decomposition of the symmetric matrix mat via
60352bfb9bbSJeremy L Thompson            symmetric QR factorization
60452bfb9bbSJeremy L Thompson 
60577645d7bSjeremylt   @param ceed         A Ceed context for error handling
60652bfb9bbSJeremy L Thompson   @param[in,out] mat  Row-major matrix to be factorized in place
607460bf743SValeria Barra   @param[out] lambda  Vector of length n of eigenvalues
60852bfb9bbSJeremy L Thompson   @param n            Number of rows/columns
60952bfb9bbSJeremy L Thompson 
61052bfb9bbSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
61152bfb9bbSJeremy L Thompson 
61252bfb9bbSJeremy L Thompson   @ref Utility
61352bfb9bbSJeremy L Thompson **/
61452bfb9bbSJeremy L Thompson int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat,
61552bfb9bbSJeremy L Thompson                                     CeedScalar *lambda, CeedInt n) {
61652bfb9bbSJeremy L Thompson   // Check bounds for clang-tidy
61752bfb9bbSJeremy L Thompson   if (n<2)
618c042f62fSJeremy L Thompson     // LCOV_EXCL_START
619c042f62fSJeremy L Thompson     return CeedError(ceed, 1,
620c042f62fSJeremy L Thompson                      "Cannot compute symmetric Schur decomposition of scalars");
621c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
62252bfb9bbSJeremy L Thompson 
62352bfb9bbSJeremy L Thompson   CeedScalar v[n-1], tau[n-1], matT[n*n];
62452bfb9bbSJeremy L Thompson 
62552bfb9bbSJeremy L Thompson   // Copy mat to matT and set mat to I
62652bfb9bbSJeremy L Thompson   memcpy(matT, mat, n*n*sizeof(mat[0]));
62752bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
62852bfb9bbSJeremy L Thompson     for (CeedInt j=0; j<n; j++)
62952bfb9bbSJeremy L Thompson       mat[j+n*i] = (i==j) ? 1 : 0;
63052bfb9bbSJeremy L Thompson 
63152bfb9bbSJeremy L Thompson   // Reduce to tridiagonal
63252bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n-1; i++) {
63352bfb9bbSJeremy L Thompson     // Calculate Householder vector, magnitude
63452bfb9bbSJeremy L Thompson     CeedScalar sigma = 0.0;
63552bfb9bbSJeremy L Thompson     v[i] = matT[i+n*(i+1)];
63652bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) {
63752bfb9bbSJeremy L Thompson       v[j] = matT[i+n*(j+1)];
63852bfb9bbSJeremy L Thompson       sigma += v[j] * v[j];
63952bfb9bbSJeremy L Thompson     }
64052bfb9bbSJeremy L Thompson     CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:n-1]
64152bfb9bbSJeremy L Thompson     CeedScalar Rii = -copysign(norm, v[i]);
64252bfb9bbSJeremy L Thompson     v[i] -= Rii;
64352bfb9bbSJeremy L Thompson     // norm of v[i:m] after modification above and scaling below
64452bfb9bbSJeremy L Thompson     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
64552bfb9bbSJeremy L Thompson     //   tau = 2 / (norm*norm)
6460e4d4210Sjeremylt     if (sigma > 10*CEED_EPSILON)
64752bfb9bbSJeremy L Thompson       tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma);
648fb551037Sjeremylt     else
649fb551037Sjeremylt       tau[i] = 0;
650fb551037Sjeremylt 
651fb551037Sjeremylt     for (CeedInt j=i+1; j<n-1; j++)
652fb551037Sjeremylt       v[j] /= v[i];
65352bfb9bbSJeremy L Thompson 
65452bfb9bbSJeremy L Thompson     // Update sub and super diagonal
65552bfb9bbSJeremy L Thompson     matT[i+n*(i+1)] = Rii;
65652bfb9bbSJeremy L Thompson     matT[(i+1)+n*i] = Rii;
65752bfb9bbSJeremy L Thompson     for (CeedInt j=i+2; j<n; j++) {
65852bfb9bbSJeremy L Thompson       matT[i+n*j] = 0; matT[j+n*i] = 0;
65952bfb9bbSJeremy L Thompson     }
66052bfb9bbSJeremy L Thompson     // Apply symmetric Householder reflector to lower right panel
66152bfb9bbSJeremy L Thompson     CeedHouseholderReflect(&matT[(i+1)+n*(i+1)], &v[i], tau[i],
66252bfb9bbSJeremy L Thompson                            n-(i+1), n-(i+1), n, 1);
66352bfb9bbSJeremy L Thompson     CeedHouseholderReflect(&matT[(i+1)+n*(i+1)], &v[i], tau[i],
66452bfb9bbSJeremy L Thompson                            n-(i+1), n-(i+1), 1, n);
66552bfb9bbSJeremy L Thompson     // Save v
66652bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) {
66752bfb9bbSJeremy L Thompson       matT[i+n*(j+1)] = v[j];
66852bfb9bbSJeremy L Thompson     }
66952bfb9bbSJeremy L Thompson   }
67052bfb9bbSJeremy L Thompson   // Backwards accumulation of Q
67152bfb9bbSJeremy L Thompson   for (CeedInt i=n-2; i>=0; i--) {
67252bfb9bbSJeremy L Thompson     v[i] = 1;
67352bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) {
67452bfb9bbSJeremy L Thompson       v[j] = matT[i+n*(j+1)];
67552bfb9bbSJeremy L Thompson       matT[i+n*(j+1)] = 0;
67652bfb9bbSJeremy L Thompson     }
67752bfb9bbSJeremy L Thompson     CeedHouseholderReflect(&mat[(i+1)+n*(i+1)], &v[i], tau[i],
67852bfb9bbSJeremy L Thompson                            n-(i+1), n-(i+1), n, 1);
67952bfb9bbSJeremy L Thompson   }
68052bfb9bbSJeremy L Thompson 
68152bfb9bbSJeremy L Thompson   // Reduce sub and super diagonal
68252bfb9bbSJeremy L Thompson   CeedInt p = 0, q = 0, itr = 0, maxitr = n*n*n;
6830e4d4210Sjeremylt   CeedScalar tol = 10*CEED_EPSILON;
68452bfb9bbSJeremy L Thompson 
68552bfb9bbSJeremy L Thompson   while (q < n && itr < maxitr) {
68652bfb9bbSJeremy L Thompson     // Update p, q, size of reduced portions of diagonal
68752bfb9bbSJeremy L Thompson     p = 0; q = 0;
68852bfb9bbSJeremy L Thompson     for (CeedInt i=n-2; i>=0; i--) {
68952bfb9bbSJeremy L Thompson       if (fabs(matT[i+n*(i+1)]) < tol)
69052bfb9bbSJeremy L Thompson         q += 1;
69152bfb9bbSJeremy L Thompson       else
69252bfb9bbSJeremy L Thompson         break;
69352bfb9bbSJeremy L Thompson     }
69452bfb9bbSJeremy L Thompson     for (CeedInt i=0; i<n-1-q; i++) {
69552bfb9bbSJeremy L Thompson       if (fabs(matT[i+n*(i+1)]) < tol)
69652bfb9bbSJeremy L Thompson         p += 1;
69752bfb9bbSJeremy L Thompson       else
69852bfb9bbSJeremy L Thompson         break;
69952bfb9bbSJeremy L Thompson     }
70052bfb9bbSJeremy L Thompson     if (q == n-1) break; // Finished reducing
70152bfb9bbSJeremy L Thompson 
70252bfb9bbSJeremy L Thompson     // Reduce tridiagonal portion
70352bfb9bbSJeremy L Thompson     CeedScalar tnn = matT[(n-1-q)+n*(n-1-q)],
70452bfb9bbSJeremy L Thompson                tnnm1 = matT[(n-2-q)+n*(n-1-q)];
70552bfb9bbSJeremy L Thompson     CeedScalar d = (matT[(n-2-q)+n*(n-2-q)] - tnn)/2;
70652bfb9bbSJeremy L Thompson     CeedScalar mu = tnn - tnnm1*tnnm1 /
70752bfb9bbSJeremy L Thompson                     (d + copysign(sqrt(d*d + tnnm1*tnnm1), d));
70852bfb9bbSJeremy L Thompson     CeedScalar x = matT[p+n*p] - mu;
70952bfb9bbSJeremy L Thompson     CeedScalar z = matT[p+n*(p+1)];
71052bfb9bbSJeremy L Thompson     for (CeedInt k=p; k<n-1-q; k++) {
71152bfb9bbSJeremy L Thompson       // Compute Givens rotation
71252bfb9bbSJeremy L Thompson       CeedScalar c = 1, s = 0;
71352bfb9bbSJeremy L Thompson       if (fabs(z) > tol) {
71452bfb9bbSJeremy L Thompson         if (fabs(z) > fabs(x)) {
71552bfb9bbSJeremy L Thompson           CeedScalar tau = -x/z;
71652bfb9bbSJeremy L Thompson           s = 1/sqrt(1+tau*tau), c = s*tau;
71752bfb9bbSJeremy L Thompson         } else {
71852bfb9bbSJeremy L Thompson           CeedScalar tau = -z/x;
71952bfb9bbSJeremy L Thompson           c = 1/sqrt(1+tau*tau), s = c*tau;
72052bfb9bbSJeremy L Thompson         }
72152bfb9bbSJeremy L Thompson       }
72252bfb9bbSJeremy L Thompson 
72352bfb9bbSJeremy L Thompson       // Apply Givens rotation to T
72452bfb9bbSJeremy L Thompson       CeedGivensRotation(matT, c, s, CEED_NOTRANSPOSE, k, k+1, n, n);
72552bfb9bbSJeremy L Thompson       CeedGivensRotation(matT, c, s, CEED_TRANSPOSE, k, k+1, n, n);
72652bfb9bbSJeremy L Thompson 
72752bfb9bbSJeremy L Thompson       // Apply Givens rotation to Q
72852bfb9bbSJeremy L Thompson       CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k+1, n, n);
72952bfb9bbSJeremy L Thompson 
73052bfb9bbSJeremy L Thompson       // Update x, z
73152bfb9bbSJeremy L Thompson       if (k < n-q-2) {
73252bfb9bbSJeremy L Thompson         x = matT[k+n*(k+1)];
73352bfb9bbSJeremy L Thompson         z = matT[k+n*(k+2)];
73452bfb9bbSJeremy L Thompson       }
73552bfb9bbSJeremy L Thompson     }
73652bfb9bbSJeremy L Thompson     itr++;
73752bfb9bbSJeremy L Thompson   }
73852bfb9bbSJeremy L Thompson   // Save eigenvalues
73952bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
74052bfb9bbSJeremy L Thompson     lambda[i] = matT[i+n*i];
74152bfb9bbSJeremy L Thompson 
74252bfb9bbSJeremy L Thompson   // Check convergence
74352bfb9bbSJeremy L Thompson   if (itr == maxitr && q < n-1)
744c042f62fSJeremy L Thompson     // LCOV_EXCL_START
74552bfb9bbSJeremy L Thompson     return CeedError(ceed, 1, "Symmetric QR failed to converge");
746c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
74752bfb9bbSJeremy L Thompson 
74852bfb9bbSJeremy L Thompson   return 0;
74952bfb9bbSJeremy L Thompson }
75052bfb9bbSJeremy L Thompson 
75152bfb9bbSJeremy L Thompson /**
7529289e5bfSjeremylt   @brief Return a reference implementation of matrix multiplication C = A B.
7539289e5bfSjeremylt            Note, this is a reference implementation for CPU CeedScalar pointers
7549289e5bfSjeremylt            that is not intended for high performance.
75552bfb9bbSJeremy L Thompson 
75677645d7bSjeremylt   @param ceed         A Ceed context for error handling
75752bfb9bbSJeremy L Thompson   @param[in] matA     Row-major matrix A
75852bfb9bbSJeremy L Thompson   @param[in] matB     Row-major matrix B
75952bfb9bbSJeremy L Thompson   @param[out] matC    Row-major output matrix C
76052bfb9bbSJeremy L Thompson   @param m            Number of rows of C
76152bfb9bbSJeremy L Thompson   @param n            Number of columns of C
76252bfb9bbSJeremy L Thompson   @param kk           Number of columns of A/rows of B
76352bfb9bbSJeremy L Thompson 
76452bfb9bbSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
76552bfb9bbSJeremy L Thompson 
76652bfb9bbSJeremy L Thompson   @ref Utility
76752bfb9bbSJeremy L Thompson **/
7689289e5bfSjeremylt int CeedMatrixMultiply(Ceed ceed, const CeedScalar *matA,
7699289e5bfSjeremylt                        const CeedScalar *matB, CeedScalar *matC, CeedInt m,
7709289e5bfSjeremylt                        CeedInt n, CeedInt kk) {
77152bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<m; i++)
77252bfb9bbSJeremy L Thompson     for (CeedInt j=0; j<n; j++) {
77352bfb9bbSJeremy L Thompson       CeedScalar sum = 0;
77452bfb9bbSJeremy L Thompson       for (CeedInt k=0; k<kk; k++)
77552bfb9bbSJeremy L Thompson         sum += matA[k+i*kk]*matB[j+k*n];
77652bfb9bbSJeremy L Thompson       matC[j+i*n] = sum;
77752bfb9bbSJeremy L Thompson     }
77852bfb9bbSJeremy L Thompson   return 0;
77952bfb9bbSJeremy L Thompson }
78052bfb9bbSJeremy L Thompson 
78152bfb9bbSJeremy L Thompson /**
78252bfb9bbSJeremy L Thompson   @brief Return Simultaneous Diagonalization of two matrices. This solves the
78352bfb9bbSJeremy L Thompson            generalized eigenvalue problem A x = lambda B x, where A and B
78452bfb9bbSJeremy L Thompson            are symmetric and B is positive definite. We generate the matrix X
78552bfb9bbSJeremy L Thompson            and vector Lambda such that X^T A X = Lambda and X^T B X = I. This
78652bfb9bbSJeremy L Thompson            is equivalent to the LAPACK routine 'sygv' with TYPE = 1.
78752bfb9bbSJeremy L Thompson 
78877645d7bSjeremylt   @param ceed         A Ceed context for error handling
78952bfb9bbSJeremy L Thompson   @param[in] matA     Row-major matrix to be factorized with eigenvalues
79052bfb9bbSJeremy L Thompson   @param[in] matB     Row-major matrix to be factorized to identity
79152bfb9bbSJeremy L Thompson   @param[out] x       Row-major orthogonal matrix
792460bf743SValeria Barra   @param[out] lambda  Vector of length n of generalized eigenvalues
79352bfb9bbSJeremy L Thompson   @param n            Number of rows/columns
79452bfb9bbSJeremy L Thompson 
79552bfb9bbSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
79652bfb9bbSJeremy L Thompson 
79752bfb9bbSJeremy L Thompson   @ref Utility
79852bfb9bbSJeremy L Thompson **/
79952bfb9bbSJeremy L Thompson int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *matA,
80052bfb9bbSJeremy L Thompson                                     CeedScalar *matB, CeedScalar *x,
80152bfb9bbSJeremy L Thompson                                     CeedScalar *lambda, CeedInt n) {
80252bfb9bbSJeremy L Thompson   int ierr;
80352bfb9bbSJeremy L Thompson   CeedScalar matC[n*n], matG[n*n], vecD[n];
80452bfb9bbSJeremy L Thompson 
80552bfb9bbSJeremy L Thompson   // Compute B = G D G^T
80652bfb9bbSJeremy L Thompson   memcpy(matG, matB, n*n*sizeof(matB[0]));
80752bfb9bbSJeremy L Thompson   ierr = CeedSymmetricSchurDecomposition(ceed, matG, vecD, n); CeedChk(ierr);
808fb551037Sjeremylt   for (CeedInt i=0; i<n; i++)
809fb551037Sjeremylt     vecD[i] = sqrt(vecD[i]);
81052bfb9bbSJeremy L Thompson 
811fb551037Sjeremylt   // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T
812fb551037Sjeremylt   //           = D^-1/2 G^T A G D^-1/2
81352bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
81452bfb9bbSJeremy L Thompson     for (CeedInt j=0; j<n; j++)
815fb551037Sjeremylt       matC[j+i*n] = matG[i+j*n] / vecD[i];
8169289e5bfSjeremylt   ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)matC,
8179289e5bfSjeremylt                             (const CeedScalar *)matA, x, n, n, n);
8189289e5bfSjeremylt   CeedChk(ierr);
81952bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
82052bfb9bbSJeremy L Thompson     for (CeedInt j=0; j<n; j++)
821fb551037Sjeremylt       matG[j+i*n] = matG[j+i*n] / vecD[j];
8229289e5bfSjeremylt   ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)x,
8239289e5bfSjeremylt                             (const CeedScalar *)matG, matC, n, n, n);
8249289e5bfSjeremylt   CeedChk(ierr);
82552bfb9bbSJeremy L Thompson 
82652bfb9bbSJeremy L Thompson   // Compute Q^T C Q = lambda
82752bfb9bbSJeremy L Thompson   ierr = CeedSymmetricSchurDecomposition(ceed, matC, lambda, n); CeedChk(ierr);
82852bfb9bbSJeremy L Thompson 
829fb551037Sjeremylt   // Set x = (G D^1/2)^-T Q
830fb551037Sjeremylt   //       = G D^-1/2 Q
8319289e5bfSjeremylt   ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)matG,
8329289e5bfSjeremylt                             (const CeedScalar *)matC, x, n, n, n);
8339289e5bfSjeremylt   CeedChk(ierr);
83452bfb9bbSJeremy L Thompson 
83552bfb9bbSJeremy L Thompson   return 0;
83652bfb9bbSJeremy L Thompson }
83752bfb9bbSJeremy L Thompson 
83852bfb9bbSJeremy L Thompson /**
839783c99b3SValeria Barra   @brief Return collocated grad matrix
840b11c1e72Sjeremylt 
841b11c1e72Sjeremylt   @param basis             CeedBasis
84295bb1877Svaleriabarra   @param[out] collograd1d  Row-major (Q1d * Q1d) matrix expressing derivatives of
843b11c1e72Sjeremylt                             basis functions at quadrature points
844b11c1e72Sjeremylt 
845b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
846dfdf5a53Sjeremylt 
847dfdf5a53Sjeremylt   @ref Advanced
848b11c1e72Sjeremylt **/
849dc7d240cSValeria Barra int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *collograd1d) {
850d7b241e6Sjeremylt   int i, j, k;
851a7bd39daSjeremylt   Ceed ceed;
852d7b241e6Sjeremylt   CeedInt ierr, P1d=(basis)->P1d, Q1d=(basis)->Q1d;
853d7b241e6Sjeremylt   CeedScalar *interp1d, *grad1d, tau[Q1d];
854d7b241e6Sjeremylt 
855d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d*P1d, &interp1d); CeedChk(ierr);
856d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d*P1d, &grad1d); CeedChk(ierr);
857d7b241e6Sjeremylt   memcpy(interp1d, (basis)->interp1d, Q1d*P1d*sizeof(basis)->interp1d[0]);
858d7b241e6Sjeremylt   memcpy(grad1d, (basis)->grad1d, Q1d*P1d*sizeof(basis)->interp1d[0]);
859d7b241e6Sjeremylt 
860d7b241e6Sjeremylt   // QR Factorization, interp1d = Q R
861a7bd39daSjeremylt   ierr = CeedBasisGetCeed(basis, &ceed); CeedChk(ierr);
862a7bd39daSjeremylt   ierr = CeedQRFactorization(ceed, interp1d, tau, Q1d, P1d); CeedChk(ierr);
863d7b241e6Sjeremylt 
864dc7d240cSValeria Barra   // Apply Rinv, collograd1d = grad1d Rinv
865d7b241e6Sjeremylt   for (i=0; i<Q1d; i++) { // Row i
866dc7d240cSValeria Barra     collograd1d[Q1d*i] = grad1d[P1d*i]/interp1d[0];
867d7b241e6Sjeremylt     for (j=1; j<P1d; j++) { // Column j
868dc7d240cSValeria Barra       collograd1d[j+Q1d*i] = grad1d[j+P1d*i];
8691d102b48SJeremy L Thompson       for (k=0; k<j; k++)
870dc7d240cSValeria Barra         collograd1d[j+Q1d*i] -= interp1d[j+P1d*k]*collograd1d[k+Q1d*i];
871dc7d240cSValeria Barra       collograd1d[j+Q1d*i] /= interp1d[j+P1d*j];
872d7b241e6Sjeremylt     }
8731d102b48SJeremy L Thompson     for (j=P1d; j<Q1d; j++)
874dc7d240cSValeria Barra       collograd1d[j+Q1d*i] = 0;
875d7b241e6Sjeremylt   }
876d7b241e6Sjeremylt 
877dc7d240cSValeria Barra   // Apply Qtranspose, collograd = collograd Qtranspose
878dc7d240cSValeria Barra   CeedHouseholderApplyQ(collograd1d, interp1d, tau, CEED_NOTRANSPOSE,
879d7b241e6Sjeremylt                         Q1d, Q1d, P1d, 1, Q1d);
880d7b241e6Sjeremylt 
881d7b241e6Sjeremylt   ierr = CeedFree(&interp1d); CeedChk(ierr);
882d7b241e6Sjeremylt   ierr = CeedFree(&grad1d); CeedChk(ierr);
883d7b241e6Sjeremylt 
884d7b241e6Sjeremylt   return 0;
885d7b241e6Sjeremylt }
886d7b241e6Sjeremylt 
887b11c1e72Sjeremylt /**
88895bb1877Svaleriabarra   @brief Apply basis evaluation from nodes to quadrature points or vice versa
889b11c1e72Sjeremylt 
890b11c1e72Sjeremylt   @param basis   CeedBasis to evaluate
891b11c1e72Sjeremylt   @param nelem   The number of elements to apply the basis evaluation to;
892b11c1e72Sjeremylt                    the backend will specify the ordering in
893b11c1e72Sjeremylt                    ElemRestrictionCreateBlocked
894b11c1e72Sjeremylt   @param tmode   \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature
895b11c1e72Sjeremylt                    points, \ref CEED_TRANSPOSE to apply the transpose, mapping
896b11c1e72Sjeremylt                    from quadrature points to nodes
897bb64980dSjeremylt   @param emode   \ref CEED_EVAL_NONE to use values directly,
898bb64980dSjeremylt                    \ref CEED_EVAL_INTERP to use interpolated values,
899bb64980dSjeremylt                    \ref CEED_EVAL_GRAD to use gradients,
900bb64980dSjeremylt                    \ref CEED_EVAL_WEIGHT to use quadrature weights.
90134138859Sjeremylt   @param[in] u   Input CeedVector
90234138859Sjeremylt   @param[out] v  Output CeedVector
903b11c1e72Sjeremylt 
904b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
905dfdf5a53Sjeremylt 
906dfdf5a53Sjeremylt   @ref Advanced
907b11c1e72Sjeremylt **/
908d7b241e6Sjeremylt int CeedBasisApply(CeedBasis basis, CeedInt nelem, CeedTransposeMode tmode,
909aedaa0e5Sjeremylt                    CeedEvalMode emode, CeedVector u, CeedVector v) {
910d7b241e6Sjeremylt   int ierr;
9118795c945Sjeremylt   CeedInt ulength = 0, vlength, nnodes, nqpt;
912c042f62fSJeremy L Thompson   if (!basis->Apply)
913c042f62fSJeremy L Thompson     // LCOV_EXCL_START
914c042f62fSJeremy L Thompson     return CeedError(basis->ceed, 1, "Backend does not support BasisApply");
915c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
916c042f62fSJeremy L Thompson 
917c042f62fSJeremy L Thompson   // Check compatibility of topological and geometrical dimensions
9188795c945Sjeremylt   ierr = CeedBasisGetNumNodes(basis, &nnodes); CeedChk(ierr);
919b502e64cSValeria Barra   ierr = CeedBasisGetNumQuadraturePoints(basis, &nqpt); CeedChk(ierr);
920b502e64cSValeria Barra   ierr = CeedVectorGetLength(v, &vlength); CeedChk(ierr);
921b502e64cSValeria Barra 
922b502e64cSValeria Barra   if (u) {
923b502e64cSValeria Barra     ierr = CeedVectorGetLength(u, &ulength); CeedChk(ierr);
924b502e64cSValeria Barra   }
925b502e64cSValeria Barra 
9261d102b48SJeremy L Thompson   if ((tmode == CEED_TRANSPOSE && (vlength%nnodes != 0 || ulength%nqpt != 0)) ||
9278795c945Sjeremylt       (tmode == CEED_NOTRANSPOSE && (ulength%nnodes != 0 || vlength%nqpt != 0)))
9281d102b48SJeremy L Thompson     return CeedError(basis->ceed, 1, "Length of input/output vectors "
9291d102b48SJeremy L Thompson                      "incompatible with basis dimensions");
930b502e64cSValeria Barra 
931d7b241e6Sjeremylt   ierr = basis->Apply(basis, nelem, tmode, emode, u, v); CeedChk(ierr);
932d7b241e6Sjeremylt   return 0;
933d7b241e6Sjeremylt }
934d7b241e6Sjeremylt 
935b11c1e72Sjeremylt /**
9364ce2993fSjeremylt   @brief Get Ceed associated with a CeedBasis
937b11c1e72Sjeremylt 
938b11c1e72Sjeremylt   @param basis      CeedBasis
9394ce2993fSjeremylt   @param[out] ceed  Variable to store Ceed
9404ce2993fSjeremylt 
9414ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
9424ce2993fSjeremylt 
94323617272Sjeremylt   @ref Advanced
9444ce2993fSjeremylt **/
9454ce2993fSjeremylt int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) {
9464ce2993fSjeremylt   *ceed = basis->ceed;
9474ce2993fSjeremylt   return 0;
9484ce2993fSjeremylt };
9494ce2993fSjeremylt 
9504ce2993fSjeremylt /**
9514ce2993fSjeremylt   @brief Get dimension for given CeedBasis
9524ce2993fSjeremylt 
9534ce2993fSjeremylt   @param basis     CeedBasis
9544ce2993fSjeremylt   @param[out] dim  Variable to store dimension of basis
9554ce2993fSjeremylt 
9564ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
9574ce2993fSjeremylt 
95823617272Sjeremylt   @ref Advanced
9594ce2993fSjeremylt **/
9604ce2993fSjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) {
9614ce2993fSjeremylt   *dim = basis->dim;
9624ce2993fSjeremylt   return 0;
9634ce2993fSjeremylt };
9644ce2993fSjeremylt 
9654ce2993fSjeremylt /**
9664ce2993fSjeremylt   @brief Get tensor status for given CeedBasis
9674ce2993fSjeremylt 
9684ce2993fSjeremylt   @param basis        CeedBasis
9694ce2993fSjeremylt   @param[out] tensor  Variable to store tensor status
9704ce2993fSjeremylt 
9714ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
9724ce2993fSjeremylt 
97323617272Sjeremylt   @ref Advanced
9744ce2993fSjeremylt **/
9754ce2993fSjeremylt int CeedBasisGetTensorStatus(CeedBasis basis, bool *tensor) {
9764ce2993fSjeremylt   *tensor = basis->tensorbasis;
9774ce2993fSjeremylt   return 0;
9784ce2993fSjeremylt };
9794ce2993fSjeremylt 
9804ce2993fSjeremylt /**
9814ce2993fSjeremylt   @brief Get number of components for given CeedBasis
9824ce2993fSjeremylt 
9834ce2993fSjeremylt   @param basis         CeedBasis
984288c0443SJeremy L Thompson   @param[out] numcomp  Variable to store number of components of basis
9854ce2993fSjeremylt 
9864ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
9874ce2993fSjeremylt 
98823617272Sjeremylt   @ref Advanced
9894ce2993fSjeremylt **/
9904ce2993fSjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *numcomp) {
9914ce2993fSjeremylt   *numcomp = basis->ncomp;
9924ce2993fSjeremylt   return 0;
9934ce2993fSjeremylt };
9944ce2993fSjeremylt 
9954ce2993fSjeremylt /**
9964ce2993fSjeremylt   @brief Get total number of nodes (in 1 dimension) of a CeedBasis
9974ce2993fSjeremylt 
9984ce2993fSjeremylt   @param basis     CeedBasis
9994ce2993fSjeremylt   @param[out] P1d  Variable to store number of nodes
10004ce2993fSjeremylt 
10014ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
10024ce2993fSjeremylt 
100323617272Sjeremylt   @ref Advanced
10044ce2993fSjeremylt **/
10054ce2993fSjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P1d) {
1006c042f62fSJeremy L Thompson   if (!basis->tensorbasis)
1007c042f62fSJeremy L Thompson     // LCOV_EXCL_START
1008c042f62fSJeremy L Thompson     return CeedError(basis->ceed, 1, "Cannot supply P1d for non-tensor basis");
1009c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
1010c042f62fSJeremy L Thompson 
10114ce2993fSjeremylt   *P1d = basis->P1d;
10124ce2993fSjeremylt   return 0;
10134ce2993fSjeremylt }
10144ce2993fSjeremylt 
10154ce2993fSjeremylt /**
10164ce2993fSjeremylt   @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis
10174ce2993fSjeremylt 
10184ce2993fSjeremylt   @param basis     CeedBasis
10194ce2993fSjeremylt   @param[out] Q1d  Variable to store number of quadrature points
10204ce2993fSjeremylt 
10214ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
10224ce2993fSjeremylt 
102323617272Sjeremylt   @ref Advanced
10244ce2993fSjeremylt **/
10254ce2993fSjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q1d) {
1026c042f62fSJeremy L Thompson   if (!basis->tensorbasis)
1027c042f62fSJeremy L Thompson     // LCOV_EXCL_START
1028c042f62fSJeremy L Thompson     return CeedError(basis->ceed, 1, "Cannot supply Q1d for non-tensor basis");
1029c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
1030c042f62fSJeremy L Thompson 
10314ce2993fSjeremylt   *Q1d = basis->Q1d;
10324ce2993fSjeremylt   return 0;
10334ce2993fSjeremylt }
10344ce2993fSjeremylt 
10354ce2993fSjeremylt /**
10364ce2993fSjeremylt   @brief Get total number of nodes (in dim dimensions) of a CeedBasis
10374ce2993fSjeremylt 
10384ce2993fSjeremylt   @param basis   CeedBasis
10394ce2993fSjeremylt   @param[out] P  Variable to store number of nodes
1040b11c1e72Sjeremylt 
1041b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1042dfdf5a53Sjeremylt 
1043dfdf5a53Sjeremylt   @ref Utility
1044b11c1e72Sjeremylt **/
1045d7b241e6Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) {
1046a8de75f0Sjeremylt   *P = basis->P;
1047d7b241e6Sjeremylt   return 0;
1048d7b241e6Sjeremylt }
1049d7b241e6Sjeremylt 
1050b11c1e72Sjeremylt /**
10514ce2993fSjeremylt   @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis
1052b11c1e72Sjeremylt 
1053b11c1e72Sjeremylt   @param basis   CeedBasis
10544ce2993fSjeremylt   @param[out] Q  Variable to store number of quadrature points
1055b11c1e72Sjeremylt 
1056b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1057dfdf5a53Sjeremylt 
1058dfdf5a53Sjeremylt   @ref Utility
1059b11c1e72Sjeremylt **/
1060d7b241e6Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) {
1061a8de75f0Sjeremylt   *Q = basis->Q;
1062d7b241e6Sjeremylt   return 0;
1063d7b241e6Sjeremylt }
1064d7b241e6Sjeremylt 
1065b11c1e72Sjeremylt /**
10668c91a0c9SJeremy L Thompson   @brief Get reference coordinates of quadrature points (in dim dimensions)
10674ce2993fSjeremylt          of a CeedBasis
10684ce2993fSjeremylt 
10694ce2993fSjeremylt   @param basis      CeedBasis
10708c91a0c9SJeremy L Thompson   @param[out] qref  Variable to store reference coordinates of quadrature points
10714ce2993fSjeremylt 
10724ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
10734ce2993fSjeremylt 
107423617272Sjeremylt   @ref Advanced
10754ce2993fSjeremylt **/
10764ce2993fSjeremylt int CeedBasisGetQRef(CeedBasis basis, CeedScalar **qref) {
10774ce2993fSjeremylt   *qref = basis->qref1d;
10784ce2993fSjeremylt   return 0;
10794ce2993fSjeremylt }
10804ce2993fSjeremylt 
10814ce2993fSjeremylt /**
10824ce2993fSjeremylt   @brief Get quadrature weights of quadrature points (in dim dimensions)
10834ce2993fSjeremylt          of a CeedBasis
10844ce2993fSjeremylt 
10854ce2993fSjeremylt   @param basis         CeedBasis
10864ce2993fSjeremylt   @param[out] qweight  Variable to store quadrature weights
10874ce2993fSjeremylt 
10884ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
10894ce2993fSjeremylt 
109023617272Sjeremylt   @ref Advanced
10914ce2993fSjeremylt **/
10924ce2993fSjeremylt int CeedBasisGetQWeights(CeedBasis basis, CeedScalar **qweight) {
10934ce2993fSjeremylt   *qweight = basis->qweight1d;
10944ce2993fSjeremylt   return 0;
10954ce2993fSjeremylt }
10964ce2993fSjeremylt 
10974ce2993fSjeremylt /**
10984ce2993fSjeremylt   @brief Get interpolation matrix of a CeedBasis
10994ce2993fSjeremylt 
11004ce2993fSjeremylt   @param basis        CeedBasis
1101288c0443SJeremy L Thompson   @param[out] interp  Variable to store interpolation matrix
11024ce2993fSjeremylt 
11034ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
11044ce2993fSjeremylt 
110523617272Sjeremylt   @ref Advanced
11064ce2993fSjeremylt **/
11074ce2993fSjeremylt int CeedBasisGetInterp(CeedBasis basis, CeedScalar **interp) {
110800f91b2bSjeremylt   if (!basis->interp && basis->tensorbasis) {
110900f91b2bSjeremylt     // Allocate
111000f91b2bSjeremylt     int ierr;
111100f91b2bSjeremylt     ierr = CeedMalloc(basis->Q*basis->P, &basis->interp); CeedChk(ierr);
111200f91b2bSjeremylt 
111300f91b2bSjeremylt     // Initialize
111400f91b2bSjeremylt     for (CeedInt i=0; i<basis->Q*basis->P; i++)
111500f91b2bSjeremylt       basis->interp[i] = 1.0;
111600f91b2bSjeremylt 
111700f91b2bSjeremylt     // Calculate
111800f91b2bSjeremylt     for (CeedInt d=0; d<basis->dim; d++)
111900f91b2bSjeremylt       for (CeedInt qpt=0; qpt<basis->Q; qpt++)
112000f91b2bSjeremylt         for (CeedInt node=0; node<basis->P; node++) {
112100f91b2bSjeremylt           CeedInt p = (node / CeedIntPow(basis->P1d, d)) % basis->P1d;
112200f91b2bSjeremylt           CeedInt q = (qpt / CeedIntPow(basis->Q1d, d)) % basis->Q1d;
112300f91b2bSjeremylt           basis->interp[qpt*(basis->P)+node] *= basis->interp1d[q*basis->P1d+p];
112400f91b2bSjeremylt         }
112500f91b2bSjeremylt   }
112600f91b2bSjeremylt 
112700f91b2bSjeremylt   *interp = basis->interp;
112800f91b2bSjeremylt 
112900f91b2bSjeremylt   return 0;
113000f91b2bSjeremylt }
113100f91b2bSjeremylt 
113200f91b2bSjeremylt /**
113300f91b2bSjeremylt   @brief Get 1D interpolation matrix of a tensor product CeedBasis
113400f91b2bSjeremylt 
113500f91b2bSjeremylt   @param basis          CeedBasis
113600f91b2bSjeremylt   @param[out] interp1d  Variable to store interpolation matrix
113700f91b2bSjeremylt 
113800f91b2bSjeremylt   @return An error code: 0 - success, otherwise - failure
113900f91b2bSjeremylt 
114000f91b2bSjeremylt   @ref Advanced
114100f91b2bSjeremylt **/
114200f91b2bSjeremylt int CeedBasisGetInterp1D(CeedBasis basis, CeedScalar **interp1d) {
114300f91b2bSjeremylt   if (!basis->tensorbasis)
114400f91b2bSjeremylt     // LCOV_EXCL_START
114500f91b2bSjeremylt     return CeedError(basis->ceed, 1, "CeedBasis is not a tensor product basis.");
114600f91b2bSjeremylt   // LCOV_EXCL_STOP
114700f91b2bSjeremylt 
114800f91b2bSjeremylt   *interp1d = basis->interp1d;
114900f91b2bSjeremylt 
11504ce2993fSjeremylt   return 0;
11514ce2993fSjeremylt }
11524ce2993fSjeremylt 
11534ce2993fSjeremylt /**
11544ce2993fSjeremylt   @brief Get gradient matrix of a CeedBasis
11554ce2993fSjeremylt 
11564ce2993fSjeremylt   @param basis      CeedBasis
1157288c0443SJeremy L Thompson   @param[out] grad  Variable to store gradient matrix
11584ce2993fSjeremylt 
11594ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
11604ce2993fSjeremylt 
116123617272Sjeremylt   @ref Advanced
11624ce2993fSjeremylt **/
11634ce2993fSjeremylt int CeedBasisGetGrad(CeedBasis basis, CeedScalar **grad) {
116400f91b2bSjeremylt   if (!basis->grad && basis->tensorbasis) {
116500f91b2bSjeremylt     // Allocate
116600f91b2bSjeremylt     int ierr;
116700f91b2bSjeremylt     ierr = CeedMalloc(basis->dim*basis->Q*basis->P, &basis->grad);
116800f91b2bSjeremylt     CeedChk(ierr);
116900f91b2bSjeremylt 
117000f91b2bSjeremylt     // Initialize
117100f91b2bSjeremylt     for (CeedInt i=0; i<basis->dim*basis->Q*basis->P; i++)
117200f91b2bSjeremylt       basis->grad[i] = 1.0;
117300f91b2bSjeremylt 
117400f91b2bSjeremylt     // Calculate
117500f91b2bSjeremylt     for (CeedInt d=0; d<basis->dim; d++)
117600f91b2bSjeremylt       for (CeedInt i=0; i<basis->dim; i++)
117700f91b2bSjeremylt         for (CeedInt qpt=0; qpt<basis->Q; qpt++)
117800f91b2bSjeremylt           for (CeedInt node=0; node<basis->P; node++) {
117900f91b2bSjeremylt             CeedInt p = (node / CeedIntPow(basis->P1d, d)) % basis->P1d;
118000f91b2bSjeremylt             CeedInt q = (qpt / CeedIntPow(basis->Q1d, d)) % basis->Q1d;
118100f91b2bSjeremylt             if (i == d)
118200f91b2bSjeremylt               basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *=
118300f91b2bSjeremylt                 basis->grad1d[q*basis->P1d+p];
118400f91b2bSjeremylt             else
118500f91b2bSjeremylt               basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *=
118600f91b2bSjeremylt                 basis->interp1d[q*basis->P1d+p];
118700f91b2bSjeremylt           }
118800f91b2bSjeremylt   }
118900f91b2bSjeremylt 
119000f91b2bSjeremylt   *grad = basis->grad;
119100f91b2bSjeremylt 
11924ce2993fSjeremylt   return 0;
11934ce2993fSjeremylt }
11944ce2993fSjeremylt 
11954ce2993fSjeremylt /**
119600f91b2bSjeremylt   @brief Get 1D gradient matrix of a tensor product CeedBasis
1197b7ec98d8SJeremy L Thompson 
1198fb551037Sjeremylt   @param basis        CeedBasis
119900f91b2bSjeremylt   @param[out] grad1d  Variable to store gradient matrix
1200b7ec98d8SJeremy L Thompson 
1201b7ec98d8SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1202b7ec98d8SJeremy L Thompson 
1203b7ec98d8SJeremy L Thompson   @ref Advanced
1204b7ec98d8SJeremy L Thompson **/
120500f91b2bSjeremylt int CeedBasisGetGrad1D(CeedBasis basis, CeedScalar **grad1d) {
120600f91b2bSjeremylt   if (!basis->tensorbasis)
1207b7ec98d8SJeremy L Thompson     // LCOV_EXCL_START
120800f91b2bSjeremylt     return CeedError(basis->ceed, 1, "CeedBasis is not a tensor product basis.");
1209b7ec98d8SJeremy L Thompson   // LCOV_EXCL_STOP
121000f91b2bSjeremylt 
121100f91b2bSjeremylt   *grad1d = basis->grad1d;
121200f91b2bSjeremylt 
1213b7ec98d8SJeremy L Thompson   return 0;
1214b7ec98d8SJeremy L Thompson }
1215b7ec98d8SJeremy L Thompson 
1216b7ec98d8SJeremy L Thompson /**
12174ce2993fSjeremylt   @brief Get backend data of a CeedBasis
12184ce2993fSjeremylt 
12194ce2993fSjeremylt   @param basis      CeedBasis
12204ce2993fSjeremylt   @param[out] data  Variable to store data
12214ce2993fSjeremylt 
12224ce2993fSjeremylt   @return An error code: 0 - success, otherwise - failure
12234ce2993fSjeremylt 
122423617272Sjeremylt   @ref Advanced
12254ce2993fSjeremylt **/
12264ce2993fSjeremylt int CeedBasisGetData(CeedBasis basis, void **data) {
12274ce2993fSjeremylt   *data = basis->data;
12284ce2993fSjeremylt   return 0;
12294ce2993fSjeremylt }
12304ce2993fSjeremylt 
12314ce2993fSjeremylt /**
1232fe2413ffSjeremylt   @brief Set backend data of a CeedBasis
1233fe2413ffSjeremylt 
1234fe2413ffSjeremylt   @param[out] basis  CeedBasis
1235fe2413ffSjeremylt   @param data        Data to set
1236fe2413ffSjeremylt 
1237fe2413ffSjeremylt   @return An error code: 0 - success, otherwise - failure
1238fe2413ffSjeremylt 
1239fe2413ffSjeremylt   @ref Advanced
1240fe2413ffSjeremylt **/
1241fe2413ffSjeremylt int CeedBasisSetData(CeedBasis basis, void **data) {
1242fe2413ffSjeremylt   basis->data = *data;
1243fe2413ffSjeremylt   return 0;
1244fe2413ffSjeremylt }
1245fe2413ffSjeremylt 
1246fe2413ffSjeremylt /**
12472f86a920SJeremy L Thompson   @brief Get CeedTensorContract of a CeedBasis
12482f86a920SJeremy L Thompson 
12492f86a920SJeremy L Thompson   @param basis          CeedBasis
12502f86a920SJeremy L Thompson   @param[out] contract  Variable to store CeedTensorContract
12512f86a920SJeremy L Thompson 
12522f86a920SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
12532f86a920SJeremy L Thompson 
12542f86a920SJeremy L Thompson   @ref Advanced
12552f86a920SJeremy L Thompson **/
12561d102b48SJeremy L Thompson int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) {
12572f86a920SJeremy L Thompson   *contract = basis->contract;
12582f86a920SJeremy L Thompson   return 0;
12592f86a920SJeremy L Thompson }
12602f86a920SJeremy L Thompson 
12612f86a920SJeremy L Thompson /**
12622f86a920SJeremy L Thompson   @brief Set CeedTensorContract of a CeedBasis
12632f86a920SJeremy L Thompson 
12642f86a920SJeremy L Thompson   @param[out] basis     CeedBasis
12652f86a920SJeremy L Thompson   @param contract       CeedTensorContract to set
12662f86a920SJeremy L Thompson 
12672f86a920SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
12682f86a920SJeremy L Thompson 
12692f86a920SJeremy L Thompson   @ref Advanced
12702f86a920SJeremy L Thompson **/
12711d102b48SJeremy L Thompson int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract *contract) {
12722f86a920SJeremy L Thompson   basis->contract = *contract;
12732f86a920SJeremy L Thompson   return 0;
12742f86a920SJeremy L Thompson }
12752f86a920SJeremy L Thompson 
12762f86a920SJeremy L Thompson /**
1277a8de75f0Sjeremylt   @brief Get dimension for given CeedElemTopology
1278a8de75f0Sjeremylt 
1279a8de75f0Sjeremylt   @param topo      CeedElemTopology
12804ce2993fSjeremylt   @param[out] dim  Variable to store dimension of topology
1281a8de75f0Sjeremylt 
1282a8de75f0Sjeremylt   @return An error code: 0 - success, otherwise - failure
1283a8de75f0Sjeremylt 
128423617272Sjeremylt   @ref Advanced
1285a8de75f0Sjeremylt **/
1286a8de75f0Sjeremylt int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) {
1287a8de75f0Sjeremylt   *dim = (CeedInt) topo >> 16;
1288a8de75f0Sjeremylt   return 0;
1289a8de75f0Sjeremylt };
1290a8de75f0Sjeremylt 
1291a8de75f0Sjeremylt /**
1292b11c1e72Sjeremylt   @brief Destroy a CeedBasis
1293b11c1e72Sjeremylt 
1294b11c1e72Sjeremylt   @param basis CeedBasis to destroy
1295b11c1e72Sjeremylt 
1296b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1297dfdf5a53Sjeremylt 
1298dfdf5a53Sjeremylt   @ref Basic
1299b11c1e72Sjeremylt **/
1300d7b241e6Sjeremylt int CeedBasisDestroy(CeedBasis *basis) {
1301d7b241e6Sjeremylt   int ierr;
1302d7b241e6Sjeremylt 
13031d102b48SJeremy L Thompson   if (!*basis || --(*basis)->refcount > 0)
13041d102b48SJeremy L Thompson     return 0;
1305d7b241e6Sjeremylt   if ((*basis)->Destroy) {
1306d7b241e6Sjeremylt     ierr = (*basis)->Destroy(*basis); CeedChk(ierr);
1307d7b241e6Sjeremylt   }
130800f91b2bSjeremylt   ierr = CeedFree(&(*basis)->interp); CeedChk(ierr);
1309d7b241e6Sjeremylt   ierr = CeedFree(&(*basis)->interp1d); CeedChk(ierr);
131000f91b2bSjeremylt   ierr = CeedFree(&(*basis)->grad); CeedChk(ierr);
1311d7b241e6Sjeremylt   ierr = CeedFree(&(*basis)->grad1d); CeedChk(ierr);
1312d7b241e6Sjeremylt   ierr = CeedFree(&(*basis)->qref1d); CeedChk(ierr);
1313d7b241e6Sjeremylt   ierr = CeedFree(&(*basis)->qweight1d); CeedChk(ierr);
1314d7b241e6Sjeremylt   ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr);
1315d7b241e6Sjeremylt   ierr = CeedFree(basis); CeedChk(ierr);
1316d7b241e6Sjeremylt   return 0;
1317d7b241e6Sjeremylt }
1318d7b241e6Sjeremylt 
131933e6becaSjeremylt /// @cond DOXYGEN_SKIP
13208795c945Sjeremylt // Indicate that the quadrature points are collocated with the nodes
1321783c99b3SValeria Barra CeedBasis CEED_BASIS_COLLOCATED = &ceed_basis_collocated;
132233e6becaSjeremylt /// @endcond
1323d7b241e6Sjeremylt /// @}
1324