xref: /libCEED/rust/libceed-sys/c-src/interface/ceed-basis.c (revision 0948660574688e7f772302e3ec183ecf15d81fe4)
1d7b241e6Sjeremylt // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2d7b241e6Sjeremylt // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3d7b241e6Sjeremylt // reserved. See files LICENSE and NOTICE for details.
4d7b241e6Sjeremylt //
5d7b241e6Sjeremylt // This file is part of CEED, a collection of benchmarks, miniapps, software
6d7b241e6Sjeremylt // libraries and APIs for efficient high-order finite element and spectral
7d7b241e6Sjeremylt // element discretizations for exascale applications. For more information and
8d7b241e6Sjeremylt // source code availability see http://github.com/ceed.
9d7b241e6Sjeremylt //
10d7b241e6Sjeremylt // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11d7b241e6Sjeremylt // a collaborative effort of two U.S. Department of Energy organizations (Office
12d7b241e6Sjeremylt // of Science and the National Nuclear Security Administration) responsible for
13d7b241e6Sjeremylt // the planning and preparation of a capable exascale ecosystem, including
14d7b241e6Sjeremylt // software, applications, hardware, advanced system engineering and early
15d7b241e6Sjeremylt // testbed platforms, in support of the nation's exascale computing imperative.
16d7b241e6Sjeremylt 
17d7b241e6Sjeremylt #include <ceed-impl.h>
18d7b241e6Sjeremylt #include <math.h>
19d7b241e6Sjeremylt #include <stdio.h>
20d7b241e6Sjeremylt #include <stdlib.h>
21d7b241e6Sjeremylt #include <string.h>
22d7b241e6Sjeremylt 
23d7b241e6Sjeremylt /// @cond DOXYGEN_SKIP
24d7b241e6Sjeremylt static struct CeedBasis_private ceed_basis_colocated;
25d7b241e6Sjeremylt /// @endcond
26d7b241e6Sjeremylt 
27d7b241e6Sjeremylt /// @file
28d7b241e6Sjeremylt /// Implementation of public CeedBasis interfaces
29d7b241e6Sjeremylt ///
30dfdf5a53Sjeremylt /// @addtogroup CeedBasis
31d7b241e6Sjeremylt /// @{
32d7b241e6Sjeremylt 
33b11c1e72Sjeremylt /**
34b11c1e72Sjeremylt   @brief Create a tensor product basis for H^1 discretizations
35b11c1e72Sjeremylt 
36b11c1e72Sjeremylt   @param ceed       A Ceed object where the CeedBasis will be created
37b11c1e72Sjeremylt   @param dim        Topological dimension
38b11c1e72Sjeremylt   @param ncomp      Number of field components (1 for scalar fields)
39b11c1e72Sjeremylt   @param P1d        Number of nodes in one dimension
40b11c1e72Sjeremylt   @param Q1d        Number of quadrature points in one dimension
41b11c1e72Sjeremylt   @param interp1d   Row-major Q1d × P1d matrix expressing the values of nodal
42b11c1e72Sjeremylt                       basis functions at quadrature points
43b11c1e72Sjeremylt   @param grad1d     Row-major Q1d × P1d matrix expressing derivatives of nodal
44b11c1e72Sjeremylt                       basis functions at quadrature points
45b11c1e72Sjeremylt   @param qref1d     Array of length Q1d holding the locations of quadrature points
46b11c1e72Sjeremylt                       on the 1D reference element [-1, 1]
47b11c1e72Sjeremylt   @param qweight1d  Array of length Q1d holding the quadrature weights on the
48b11c1e72Sjeremylt                       reference element
49b11c1e72Sjeremylt   @param[out] basis Address of the variable where the newly created
50b11c1e72Sjeremylt                       CeedBasis will be stored.
51b11c1e72Sjeremylt 
52b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
53dfdf5a53Sjeremylt 
54dfdf5a53Sjeremylt   @ref Basic
55b11c1e72Sjeremylt **/
56d7b241e6Sjeremylt int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt ncomp, CeedInt P1d,
57d7b241e6Sjeremylt                             CeedInt Q1d, const CeedScalar *interp1d,
58d7b241e6Sjeremylt                             const CeedScalar *grad1d, const CeedScalar *qref1d,
59d7b241e6Sjeremylt                             const CeedScalar *qweight1d, CeedBasis *basis) {
60d7b241e6Sjeremylt   int ierr;
61d7b241e6Sjeremylt 
62d7b241e6Sjeremylt   if (!ceed->BasisCreateTensorH1)
63d7b241e6Sjeremylt     return CeedError(ceed, 1, "Backend does not support BasisCreateTensorH1");
64d7b241e6Sjeremylt   ierr = CeedCalloc(1,basis); CeedChk(ierr);
65d7b241e6Sjeremylt   (*basis)->ceed = ceed;
66d7b241e6Sjeremylt   ceed->refcount++;
67d7b241e6Sjeremylt   (*basis)->refcount = 1;
68d7b241e6Sjeremylt   (*basis)->dim = dim;
69d7b241e6Sjeremylt   (*basis)->ncomp = ncomp;
70d7b241e6Sjeremylt   (*basis)->P1d = P1d;
71d7b241e6Sjeremylt   (*basis)->Q1d = Q1d;
72d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d,&(*basis)->qref1d); CeedChk(ierr);
73d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d,&(*basis)->qweight1d); CeedChk(ierr);
74d7b241e6Sjeremylt   memcpy((*basis)->qref1d, qref1d, Q1d*sizeof(qref1d[0]));
75d7b241e6Sjeremylt   memcpy((*basis)->qweight1d, qweight1d, Q1d*sizeof(qweight1d[0]));
76d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d*P1d,&(*basis)->interp1d); CeedChk(ierr);
77d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d*P1d,&(*basis)->grad1d); CeedChk(ierr);
78d7b241e6Sjeremylt   memcpy((*basis)->interp1d, interp1d, Q1d*P1d*sizeof(interp1d[0]));
79*09486605Sjeremylt   memcpy((*basis)->grad1d, grad1d, Q1d*P1d*sizeof(grad1d[0]));
80d7b241e6Sjeremylt   ierr = ceed->BasisCreateTensorH1(ceed, dim, P1d, Q1d, interp1d, grad1d, qref1d,
81d7b241e6Sjeremylt                                    qweight1d, *basis); CeedChk(ierr);
82d7b241e6Sjeremylt   return 0;
83d7b241e6Sjeremylt }
84d7b241e6Sjeremylt 
85b11c1e72Sjeremylt /**
86b11c1e72Sjeremylt   @brief Create a tensor product Lagrange basis
87b11c1e72Sjeremylt 
88b11c1e72Sjeremylt   @param ceed       A Ceed object where the CeedBasis will be created
89b11c1e72Sjeremylt   @param dim        Topological dimension of element
90b11c1e72Sjeremylt   @param ncomp      Number of field components
91b11c1e72Sjeremylt   @param P          Number of Gauss-Lobatto nodes in one dimension.  The
92b11c1e72Sjeremylt                       polynomial degree of the resulting Q_k element is k=P-1.
93b11c1e72Sjeremylt   @param Q          Number of quadrature points in one dimension.
94b11c1e72Sjeremylt   @param qmode      Distribution of the Q quadrature points (affects order of
95b11c1e72Sjeremylt                       accuracy for the quadrature)
96b11c1e72Sjeremylt   @param[out] basis Address of the variable where the newly created
97b11c1e72Sjeremylt                       CeedBasis will be stored.
98b11c1e72Sjeremylt 
99b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
100dfdf5a53Sjeremylt 
101dfdf5a53Sjeremylt   @ref Basic
102b11c1e72Sjeremylt **/
103d7b241e6Sjeremylt int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt ncomp,
104d7b241e6Sjeremylt                                     CeedInt P, CeedInt Q,
105d7b241e6Sjeremylt                                     CeedQuadMode qmode, CeedBasis *basis) {
106d7b241e6Sjeremylt   // Allocate
107d7b241e6Sjeremylt   int ierr, i, j, k;
108d7b241e6Sjeremylt   CeedScalar c1, c2, c3, c4, dx, *nodes, *interp1d, *grad1d, *qref1d, *qweight1d;
109d7b241e6Sjeremylt   ierr = CeedCalloc(P*Q, &interp1d); CeedChk(ierr);
110d7b241e6Sjeremylt   ierr = CeedCalloc(P*Q, &grad1d); CeedChk(ierr);
111d7b241e6Sjeremylt   ierr = CeedCalloc(P, &nodes); CeedChk(ierr);
112d7b241e6Sjeremylt   ierr = CeedCalloc(Q, &qref1d); CeedChk(ierr);
113d7b241e6Sjeremylt   ierr = CeedCalloc(Q, &qweight1d); CeedChk(ierr);
114d7b241e6Sjeremylt   // Get Nodes and Weights
115d7b241e6Sjeremylt   ierr = CeedLobattoQuadrature(P, nodes, NULL); CeedChk(ierr);
116d7b241e6Sjeremylt   switch (qmode) {
117d7b241e6Sjeremylt   case CEED_GAUSS:
118d7b241e6Sjeremylt     ierr = CeedGaussQuadrature(Q, qref1d, qweight1d); CeedChk(ierr);
119d7b241e6Sjeremylt     break;
120d7b241e6Sjeremylt   case CEED_GAUSS_LOBATTO:
121d7b241e6Sjeremylt     ierr = CeedLobattoQuadrature(Q, qref1d, qweight1d); CeedChk(ierr);
122d7b241e6Sjeremylt     break;
123d7b241e6Sjeremylt   }
124d7b241e6Sjeremylt   // Build B, D matrix
125d7b241e6Sjeremylt   // Fornberg, 1998
126d7b241e6Sjeremylt   for (i = 0; i  < Q; i++) {
127d7b241e6Sjeremylt     c1 = 1.0;
128d7b241e6Sjeremylt     c3 = nodes[0] - qref1d[i];
129d7b241e6Sjeremylt     interp1d[i*P+0] = 1.0;
130d7b241e6Sjeremylt     for (j = 1; j < P; j++) {
131d7b241e6Sjeremylt       c2 = 1.0;
132d7b241e6Sjeremylt       c4 = c3;
133d7b241e6Sjeremylt       c3 = nodes[j] - qref1d[i];
134d7b241e6Sjeremylt       for (k = 0; k < j; k++) {
135d7b241e6Sjeremylt         dx = nodes[j] - nodes[k];
136d7b241e6Sjeremylt         c2 *= dx;
137d7b241e6Sjeremylt         if (k == j - 1) {
138d7b241e6Sjeremylt           grad1d[i*P + j] = c1*(interp1d[i*P + k] - c4*grad1d[i*P + k]) / c2;
139d7b241e6Sjeremylt           interp1d[i*P + j] = - c1*c4*interp1d[i*P + k] / c2;
140d7b241e6Sjeremylt         }
141d7b241e6Sjeremylt         grad1d[i*P + k] = (c3*grad1d[i*P + k] - interp1d[i*P + k]) / dx;
142d7b241e6Sjeremylt         interp1d[i*P + k] = c3*interp1d[i*P + k] / dx;
143d7b241e6Sjeremylt       }
144d7b241e6Sjeremylt       c1 = c2;
145d7b241e6Sjeremylt     }
146d7b241e6Sjeremylt   }
147d7b241e6Sjeremylt   //  // Pass to CeedBasisCreateTensorH1
148d7b241e6Sjeremylt   ierr = CeedBasisCreateTensorH1(ceed, dim, ncomp, P, Q, interp1d, grad1d, qref1d,
149d7b241e6Sjeremylt                                  qweight1d, basis); CeedChk(ierr);
150d7b241e6Sjeremylt   ierr = CeedFree(&interp1d); CeedChk(ierr);
151d7b241e6Sjeremylt   ierr = CeedFree(&grad1d); CeedChk(ierr);
152d7b241e6Sjeremylt   ierr = CeedFree(&nodes); CeedChk(ierr);
153d7b241e6Sjeremylt   ierr = CeedFree(&qref1d); CeedChk(ierr);
154d7b241e6Sjeremylt   ierr = CeedFree(&qweight1d); CeedChk(ierr);
155d7b241e6Sjeremylt   return 0;
156d7b241e6Sjeremylt }
157d7b241e6Sjeremylt 
158b11c1e72Sjeremylt /**
159b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre quadrature
160b11c1e72Sjeremylt 
161b11c1e72Sjeremylt   @param Q              Number of quadrature points (integrates polynomials of
162b11c1e72Sjeremylt                           degree 2*Q-1 exactly)
163b11c1e72Sjeremylt   @param[out] qref1d    Array of length Q to hold the abscissa on [-1, 1]
164b11c1e72Sjeremylt   @param[out] qweight1d Array of length Q to hold the weights
165b11c1e72Sjeremylt 
166b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
167dfdf5a53Sjeremylt 
168dfdf5a53Sjeremylt   @ref Utility
169b11c1e72Sjeremylt **/
170d7b241e6Sjeremylt int CeedGaussQuadrature(CeedInt Q, CeedScalar *qref1d, CeedScalar *qweight1d) {
171d7b241e6Sjeremylt   // Allocate
172d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0*atan(1.0);
173d7b241e6Sjeremylt   // Build qref1d, qweight1d
174d7b241e6Sjeremylt   for (int i = 0; i <= Q/2; i++) {
175d7b241e6Sjeremylt     // Guess
176d7b241e6Sjeremylt     xi = cos(PI*(CeedScalar)(2*i+1)/((CeedScalar)(2*Q)));
177d7b241e6Sjeremylt     // Pn(xi)
178d7b241e6Sjeremylt     P0 = 1.0;
179d7b241e6Sjeremylt     P1 = xi;
180d7b241e6Sjeremylt     P2 = 0.0;
181d7b241e6Sjeremylt     for (int j = 2; j <= Q; j++) {
182d7b241e6Sjeremylt       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
183d7b241e6Sjeremylt       P0 = P1;
184d7b241e6Sjeremylt       P1 = P2;
185d7b241e6Sjeremylt     }
186d7b241e6Sjeremylt     // First Newton Step
187d7b241e6Sjeremylt     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
188d7b241e6Sjeremylt     xi = xi-P2/dP2;
189d7b241e6Sjeremylt     // Newton to convergence
190d7b241e6Sjeremylt     for (int k=0; k<100 && fabs(P2)>1e-15; k++) {
191d7b241e6Sjeremylt       P0 = 1.0;
192d7b241e6Sjeremylt       P1 = xi;
193d7b241e6Sjeremylt       for (int j = 2; j <= Q; j++) {
194d7b241e6Sjeremylt         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
195d7b241e6Sjeremylt         P0 = P1;
196d7b241e6Sjeremylt         P1 = P2;
197d7b241e6Sjeremylt       }
198d7b241e6Sjeremylt       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
199d7b241e6Sjeremylt       xi = xi-P2/dP2;
200d7b241e6Sjeremylt     }
201d7b241e6Sjeremylt     // Save xi, wi
202d7b241e6Sjeremylt     wi = 2.0/((1.0-xi*xi)*dP2*dP2);
203d7b241e6Sjeremylt     qweight1d[i] = wi;
204d7b241e6Sjeremylt     qweight1d[Q-1-i] = wi;
205d7b241e6Sjeremylt     qref1d[i] = -xi;
206d7b241e6Sjeremylt     qref1d[Q-1-i]= xi;
207d7b241e6Sjeremylt   }
208d7b241e6Sjeremylt   return 0;
209d7b241e6Sjeremylt }
210d7b241e6Sjeremylt 
211b11c1e72Sjeremylt /**
212b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre-Lobatto quadrature
213b11c1e72Sjeremylt 
214b11c1e72Sjeremylt   @param Q              Number of quadrature points (integrates polynomials of
215b11c1e72Sjeremylt                           degree 2*Q-3 exactly)
216b11c1e72Sjeremylt   @param[out] qref1d    Array of length Q to hold the abscissa on [-1, 1]
217b11c1e72Sjeremylt   @param[out] qweight1d Array of length Q to hold the weights
218b11c1e72Sjeremylt 
219b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
220dfdf5a53Sjeremylt 
221dfdf5a53Sjeremylt   @ref Utility
222b11c1e72Sjeremylt **/
223d7b241e6Sjeremylt int CeedLobattoQuadrature(CeedInt Q, CeedScalar *qref1d,
224d7b241e6Sjeremylt                           CeedScalar *qweight1d) {
225d7b241e6Sjeremylt   // Allocate
226d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0*atan(1.0);
227d7b241e6Sjeremylt   // Build qref1d, qweight1d
228d7b241e6Sjeremylt   // Set endpoints
229d7b241e6Sjeremylt   wi = 2.0/((CeedScalar)(Q*(Q-1)));
230d7b241e6Sjeremylt   if (qweight1d) {
231d7b241e6Sjeremylt     qweight1d[0] = wi;
232d7b241e6Sjeremylt     qweight1d[Q-1] = wi;
233d7b241e6Sjeremylt   }
234d7b241e6Sjeremylt   qref1d[0] = -1.0;
235d7b241e6Sjeremylt   qref1d[Q-1] = 1.0;
236d7b241e6Sjeremylt   // Interior
237d7b241e6Sjeremylt   for (int i = 1; i <= (Q-1)/2; i++) {
238d7b241e6Sjeremylt     // Guess
239d7b241e6Sjeremylt     xi = cos(PI*(CeedScalar)(i)/(CeedScalar)(Q-1));
240d7b241e6Sjeremylt     // Pn(xi)
241d7b241e6Sjeremylt     P0 = 1.0;
242d7b241e6Sjeremylt     P1 = xi;
243d7b241e6Sjeremylt     P2 = 0.0;
244d7b241e6Sjeremylt     for (int j = 2; j < Q; j++) {
245d7b241e6Sjeremylt       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
246d7b241e6Sjeremylt       P0 = P1;
247d7b241e6Sjeremylt       P1 = P2;
248d7b241e6Sjeremylt     }
249d7b241e6Sjeremylt     // First Newton step
250d7b241e6Sjeremylt     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
251d7b241e6Sjeremylt     d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
252d7b241e6Sjeremylt     xi = xi-dP2/d2P2;
253d7b241e6Sjeremylt     // Newton to convergence
254d7b241e6Sjeremylt     for (int k=0; k<100 && fabs(dP2)>1e-15; k++) {
255d7b241e6Sjeremylt       P0 = 1.0;
256d7b241e6Sjeremylt       P1 = xi;
257d7b241e6Sjeremylt       for (int j = 2; j < Q; j++) {
258d7b241e6Sjeremylt         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
259d7b241e6Sjeremylt         P0 = P1;
260d7b241e6Sjeremylt         P1 = P2;
261d7b241e6Sjeremylt       }
262d7b241e6Sjeremylt       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
263d7b241e6Sjeremylt       d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
264d7b241e6Sjeremylt       xi = xi-dP2/d2P2;
265d7b241e6Sjeremylt     }
266d7b241e6Sjeremylt     // Save xi, wi
267d7b241e6Sjeremylt     wi = 2.0/(((CeedScalar)(Q*(Q-1)))*P2*P2);
268d7b241e6Sjeremylt     if (qweight1d) {
269d7b241e6Sjeremylt       qweight1d[i] = wi;
270d7b241e6Sjeremylt       qweight1d[Q-1-i] = wi;
271d7b241e6Sjeremylt     }
272d7b241e6Sjeremylt     qref1d[i] = -xi;
273d7b241e6Sjeremylt     qref1d[Q-1-i]= xi;
274d7b241e6Sjeremylt   }
275d7b241e6Sjeremylt   return 0;
276d7b241e6Sjeremylt }
277d7b241e6Sjeremylt 
278dfdf5a53Sjeremylt /**
279dfdf5a53Sjeremylt   @brief View an array stored in a CeedBasis
280dfdf5a53Sjeremylt 
281dfdf5a53Sjeremylt   @param name      Name of array
282dfdf5a53Sjeremylt   @param fpformat  Printing format
283dfdf5a53Sjeremylt   @param m         Number of rows in array
284dfdf5a53Sjeremylt   @param n         Number of columns in array
285dfdf5a53Sjeremylt   @param a         Array to be viewed
286dfdf5a53Sjeremylt   @param stream    Stream to view to, e.g., stdout
287dfdf5a53Sjeremylt 
288dfdf5a53Sjeremylt   @return An error code: 0 - success, otherwise - failure
289dfdf5a53Sjeremylt 
290dfdf5a53Sjeremylt   @ref Utility
291dfdf5a53Sjeremylt **/
292d7b241e6Sjeremylt static int CeedScalarView(const char *name, const char *fpformat, CeedInt m,
293d7b241e6Sjeremylt                           CeedInt n, const CeedScalar *a, FILE *stream) {
294d7b241e6Sjeremylt   for (int i=0; i<m; i++) {
295d7b241e6Sjeremylt     if (m > 1) fprintf(stream, "%12s[%d]:", name, i);
296d7b241e6Sjeremylt     else fprintf(stream, "%12s:", name);
297d7b241e6Sjeremylt     for (int j=0; j<n; j++) {
298d7b241e6Sjeremylt       fprintf(stream, fpformat, fabs(a[i*n+j]) > 1E-14 ? a[i*n+j] : 0);
299d7b241e6Sjeremylt     }
300d7b241e6Sjeremylt     fputs("\n", stream);
301d7b241e6Sjeremylt   }
302d7b241e6Sjeremylt   return 0;
303d7b241e6Sjeremylt }
304d7b241e6Sjeremylt 
305b11c1e72Sjeremylt /**
306b11c1e72Sjeremylt   @brief View a CeedBasis
307b11c1e72Sjeremylt 
308b11c1e72Sjeremylt   @param basis  CeedBasis to view
309b11c1e72Sjeremylt   @param stream Stream to view to, e.g., stdout
310b11c1e72Sjeremylt 
311b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
312dfdf5a53Sjeremylt 
313dfdf5a53Sjeremylt   @ref Utility
314b11c1e72Sjeremylt **/
315d7b241e6Sjeremylt int CeedBasisView(CeedBasis basis, FILE *stream) {
316d7b241e6Sjeremylt   int ierr;
317d7b241e6Sjeremylt 
318d7b241e6Sjeremylt   fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P1d,
319d7b241e6Sjeremylt           basis->Q1d);
320d7b241e6Sjeremylt   ierr = CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q1d, basis->qref1d,
321d7b241e6Sjeremylt                         stream); CeedChk(ierr);
322d7b241e6Sjeremylt   ierr = CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q1d, basis->qweight1d,
323d7b241e6Sjeremylt                         stream); CeedChk(ierr);
324d7b241e6Sjeremylt   ierr = CeedScalarView("interp1d", "\t% 12.8f", basis->Q1d, basis->P1d,
325d7b241e6Sjeremylt                         basis->interp1d, stream); CeedChk(ierr);
326d7b241e6Sjeremylt   ierr = CeedScalarView("grad1d", "\t% 12.8f", basis->Q1d, basis->P1d,
327d7b241e6Sjeremylt                         basis->grad1d, stream); CeedChk(ierr);
328d7b241e6Sjeremylt   return 0;
329d7b241e6Sjeremylt }
330d7b241e6Sjeremylt 
331dfdf5a53Sjeremylt /**
332dfdf5a53Sjeremylt   @brief Compute Householder Reflection
333dfdf5a53Sjeremylt 
334dfdf5a53Sjeremylt     Computes A = (I - b v v^T) A
335dfdf5a53Sjeremylt     where A is an mxn matrix indexed as A[i*row + j*col]
336dfdf5a53Sjeremylt 
337dfdf5a53Sjeremylt   @param[out] A  Matrix to apply Householder reflection to, in place
338dfdf5a53Sjeremylt   @param v       Householder vector
339dfdf5a53Sjeremylt   @param b       Scaling factor
340dfdf5a53Sjeremylt   @param m       Number of rows in A
341dfdf5a53Sjeremylt   @param n       Number of columns in A
342dfdf5a53Sjeremylt   @param row     Col stride
343dfdf5a53Sjeremylt   @param col     Row stride
344dfdf5a53Sjeremylt 
345dfdf5a53Sjeremylt   @return An error code: 0 - success, otherwise - failure
346dfdf5a53Sjeremylt 
347dfdf5a53Sjeremylt   @ref Developer
348dfdf5a53Sjeremylt **/
349d7b241e6Sjeremylt static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v,
350d7b241e6Sjeremylt                                   CeedScalar b, CeedInt m, CeedInt n,
351d7b241e6Sjeremylt                                   CeedInt row, CeedInt col) {
352d7b241e6Sjeremylt   for (CeedInt j=0; j<n; j++) {
353d7b241e6Sjeremylt     CeedScalar w = A[0*row + j*col];
354d7b241e6Sjeremylt     for (CeedInt i=1; i<m; i++) w += v[i] * A[i*row + j*col];
355d7b241e6Sjeremylt     A[0*row + j*col] -= b * w;
356d7b241e6Sjeremylt     for (CeedInt i=1; i<m; i++) A[i*row + j*col] -= b * w * v[i];
357d7b241e6Sjeremylt   }
358d7b241e6Sjeremylt   return 0;
359d7b241e6Sjeremylt }
360d7b241e6Sjeremylt 
361dfdf5a53Sjeremylt /**
362dfdf5a53Sjeremylt   @brief Apply Householder Q matrix
363dfdf5a53Sjeremylt 
364dfdf5a53Sjeremylt     Compute A = Q A where Q is mxk and A is mxn. k<m
365dfdf5a53Sjeremylt 
366dfdf5a53Sjeremylt   @param[out] A  Matrix to apply Householder Q to, in place
367dfdf5a53Sjeremylt   @param Q       Householder Q matrix
368dfdf5a53Sjeremylt   @param tau     Householder scaling factors
369dfdf5a53Sjeremylt   @param tmode   Transpose mode for application
370dfdf5a53Sjeremylt   @param m       Number of rows in A
371dfdf5a53Sjeremylt   @param n       Number of columns in A
372dfdf5a53Sjeremylt   @param k       Index of row targeted
373dfdf5a53Sjeremylt   @param row     Col stride
374dfdf5a53Sjeremylt   @param col     Row stride
375dfdf5a53Sjeremylt 
376dfdf5a53Sjeremylt   @return An error code: 0 - success, otherwise - failure
377dfdf5a53Sjeremylt 
378dfdf5a53Sjeremylt   @ref Developer
379dfdf5a53Sjeremylt **/
380d7b241e6Sjeremylt static int CeedHouseholderApplyQ(CeedScalar *A, const CeedScalar *Q,
381d7b241e6Sjeremylt                                  const CeedScalar *tau, CeedTransposeMode tmode,
382d7b241e6Sjeremylt                                  CeedInt m, CeedInt n, CeedInt k,
383d7b241e6Sjeremylt                                  CeedInt row, CeedInt col) {
384d7b241e6Sjeremylt   CeedScalar v[m];
385d7b241e6Sjeremylt   for (CeedInt ii=0; ii<k; ii++) {
386d7b241e6Sjeremylt     CeedInt i = tmode == CEED_TRANSPOSE ? ii : k-1-ii;
387d7b241e6Sjeremylt     for (CeedInt j=i+1; j<m; j++) {
388d7b241e6Sjeremylt       v[j] = Q[j*k+i];
389d7b241e6Sjeremylt     }
390d7b241e6Sjeremylt     // Apply Householder reflector (I - tau v v^T) colograd1d^T
391d7b241e6Sjeremylt     CeedHouseholderReflect(&A[i*row], &v[i], tau[i], m-i, n, row, col);
392d7b241e6Sjeremylt   }
393d7b241e6Sjeremylt   return 0;
394d7b241e6Sjeremylt }
395d7b241e6Sjeremylt 
396b11c1e72Sjeremylt /**
397b11c1e72Sjeremylt   @brief Return QR Factorization of matrix
398b11c1e72Sjeremylt 
399b11c1e72Sjeremylt   @param[out] mat  Row-major matrix to be factorized in place
400b11c1e72Sjeremylt   @param[out] tau  Vector of length m of scaling fators
401b11c1e72Sjeremylt   @param m         Number of rows
402b11c1e72Sjeremylt   @param n         Number of columns
403b11c1e72Sjeremylt 
404b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
405dfdf5a53Sjeremylt 
406dfdf5a53Sjeremylt   @ref Utility
407b11c1e72Sjeremylt **/
408d7b241e6Sjeremylt int CeedQRFactorization(CeedScalar *mat, CeedScalar *tau,
409d7b241e6Sjeremylt                         CeedInt m, CeedInt n) {
410d7b241e6Sjeremylt   CeedInt i, j;
411d7b241e6Sjeremylt   CeedScalar v[m];
412d7b241e6Sjeremylt 
413d7b241e6Sjeremylt   for (i=0; i<n; i++) {
414d7b241e6Sjeremylt     // Calculate Householder vector, magnitude
415d7b241e6Sjeremylt     CeedScalar sigma = 0.0;
416d7b241e6Sjeremylt     v[i] = mat[i+n*i];
417d7b241e6Sjeremylt     for (j=i+1; j<m; j++) {
418d7b241e6Sjeremylt       v[j] = mat[i+n*j];
419d7b241e6Sjeremylt       sigma += v[j] * v[j];
420d7b241e6Sjeremylt     }
421d7b241e6Sjeremylt     CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m]
422d7b241e6Sjeremylt     CeedScalar Rii = -copysign(norm, v[i]);
423d7b241e6Sjeremylt     v[i] -= Rii;
424d7b241e6Sjeremylt     // norm of v[i:m] after modification above and scaling below
425d7b241e6Sjeremylt     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
426d7b241e6Sjeremylt     //   tau = 2 / (norm*norm)
427d7b241e6Sjeremylt     tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma);
428d7b241e6Sjeremylt     for (j=i+1; j<m; j++) v[j] /= v[i];
429d7b241e6Sjeremylt 
430d7b241e6Sjeremylt     // Apply Householder reflector to lower right panel
431d7b241e6Sjeremylt     CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1);
432d7b241e6Sjeremylt     // Save v
433d7b241e6Sjeremylt     mat[i+n*i] = Rii;
434d7b241e6Sjeremylt     for (j=i+1; j<m; j++) {
435d7b241e6Sjeremylt       mat[i+n*j] = v[j];
436d7b241e6Sjeremylt     }
437d7b241e6Sjeremylt   }
438d7b241e6Sjeremylt 
439d7b241e6Sjeremylt   return 0;
440d7b241e6Sjeremylt }
441d7b241e6Sjeremylt 
442b11c1e72Sjeremylt /**
443b11c1e72Sjeremylt   @brief Return colocated grad matrix
444b11c1e72Sjeremylt 
445b11c1e72Sjeremylt   @param basis           CeedBasis
446b11c1e72Sjeremylt   @param[out] colograd1d Row-major Q1d × Q1d matrix expressing derivatives of
447b11c1e72Sjeremylt                            basis functions at quadrature points
448b11c1e72Sjeremylt 
449b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
450dfdf5a53Sjeremylt 
451dfdf5a53Sjeremylt   @ref Advanced
452b11c1e72Sjeremylt **/
453d7b241e6Sjeremylt int CeedBasisGetColocatedGrad(CeedBasis basis, CeedScalar *colograd1d) {
454d7b241e6Sjeremylt   int i, j, k;
455d7b241e6Sjeremylt   CeedInt ierr, P1d=(basis)->P1d, Q1d=(basis)->Q1d;
456d7b241e6Sjeremylt   CeedScalar *interp1d, *grad1d, tau[Q1d];
457d7b241e6Sjeremylt 
458d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d*P1d, &interp1d); CeedChk(ierr);
459d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d*P1d, &grad1d); CeedChk(ierr);
460d7b241e6Sjeremylt   memcpy(interp1d, (basis)->interp1d, Q1d*P1d*sizeof(basis)->interp1d[0]);
461d7b241e6Sjeremylt   memcpy(grad1d, (basis)->grad1d, Q1d*P1d*sizeof(basis)->interp1d[0]);
462d7b241e6Sjeremylt 
463d7b241e6Sjeremylt   // QR Factorization, interp1d = Q R
464d7b241e6Sjeremylt   ierr = CeedQRFactorization(interp1d, tau, Q1d, P1d); CeedChk(ierr);
465d7b241e6Sjeremylt 
466d7b241e6Sjeremylt   // Apply Rinv, colograd1d = grad1d Rinv
467d7b241e6Sjeremylt   for (i=0; i<Q1d; i++) { // Row i
468d7b241e6Sjeremylt     colograd1d[Q1d*i] = grad1d[P1d*i]/interp1d[0];
469d7b241e6Sjeremylt     for (j=1; j<P1d; j++) { // Column j
470d7b241e6Sjeremylt       colograd1d[j+Q1d*i] = grad1d[j+P1d*i];
471d7b241e6Sjeremylt       for (k=0; k<j; k++) {
472d7b241e6Sjeremylt         colograd1d[j+Q1d*i] -= interp1d[j+P1d*k]*colograd1d[k+Q1d*i];
473d7b241e6Sjeremylt       }
474d7b241e6Sjeremylt       colograd1d[j+Q1d*i] /= interp1d[j+P1d*j];
475d7b241e6Sjeremylt     }
476d7b241e6Sjeremylt     for (j=P1d; j<Q1d; j++) {
477d7b241e6Sjeremylt       colograd1d[j+Q1d*i] = 0;
478d7b241e6Sjeremylt     }
479d7b241e6Sjeremylt   }
480d7b241e6Sjeremylt 
481d7b241e6Sjeremylt   // Apply Qtranspose, colograd = colograd Qtranspose
482d7b241e6Sjeremylt   CeedHouseholderApplyQ(colograd1d, interp1d, tau, CEED_NOTRANSPOSE,
483d7b241e6Sjeremylt                         Q1d, Q1d, P1d, 1, Q1d);
484d7b241e6Sjeremylt 
485d7b241e6Sjeremylt   ierr = CeedFree(&interp1d); CeedChk(ierr);
486d7b241e6Sjeremylt   ierr = CeedFree(&grad1d); CeedChk(ierr);
487d7b241e6Sjeremylt 
488d7b241e6Sjeremylt   return 0;
489d7b241e6Sjeremylt }
490d7b241e6Sjeremylt 
491b11c1e72Sjeremylt /**
492b11c1e72Sjeremylt   @brief Apply basis evaluation from nodes to quadrature points or vice-versa
493b11c1e72Sjeremylt 
494b11c1e72Sjeremylt   @param basis  CeedBasis to evaluate
495b11c1e72Sjeremylt   @param nelem  The number of elements to apply the basis evaluation to;
496b11c1e72Sjeremylt                   the backend will specify the ordering in
497b11c1e72Sjeremylt                   ElemRestrictionCreateBlocked
498b11c1e72Sjeremylt   @param tmode  \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature
499b11c1e72Sjeremylt                   points, \ref CEED_TRANSPOSE to apply the transpose, mapping
500b11c1e72Sjeremylt                   from quadrature points to nodes
501b11c1e72Sjeremylt   @param emode  \ref CEED_EVAL_INTERP to obtain interpolated values,
502b11c1e72Sjeremylt                   \ref CEED_EVAL_GRAD to obtain gradients.
503b11c1e72Sjeremylt   @param[in] u  Input array
504b11c1e72Sjeremylt   @param[out] v Output array
505b11c1e72Sjeremylt 
506b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
507dfdf5a53Sjeremylt 
508dfdf5a53Sjeremylt   @ref Advanced
509b11c1e72Sjeremylt **/
510d7b241e6Sjeremylt int CeedBasisApply(CeedBasis basis, CeedInt nelem, CeedTransposeMode tmode,
511d7b241e6Sjeremylt                    CeedEvalMode emode, const CeedScalar *u, CeedScalar *v) {
512d7b241e6Sjeremylt   int ierr;
513d7b241e6Sjeremylt   if (!basis->Apply) return CeedError(basis->ceed, 1,
514d7b241e6Sjeremylt                                         "Backend does not support BasisApply");
515d7b241e6Sjeremylt   ierr = basis->Apply(basis, nelem, tmode, emode, u, v); CeedChk(ierr);
516d7b241e6Sjeremylt   return 0;
517d7b241e6Sjeremylt }
518d7b241e6Sjeremylt 
519b11c1e72Sjeremylt /**
520b11c1e72Sjeremylt   @brief Get total number of nodes (in dim dimensions)
521b11c1e72Sjeremylt 
522b11c1e72Sjeremylt   @param basis   CeedBasis
523b11c1e72Sjeremylt   @param[out] P  Number of nodes
524b11c1e72Sjeremylt 
525b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
526dfdf5a53Sjeremylt 
527dfdf5a53Sjeremylt   @ref Utility
528b11c1e72Sjeremylt **/
529d7b241e6Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) {
530b5cf12eeSjeremylt   *P = CeedIntPow(basis->P1d, basis->dim);
531d7b241e6Sjeremylt   return 0;
532d7b241e6Sjeremylt }
533d7b241e6Sjeremylt 
534b11c1e72Sjeremylt /**
535b11c1e72Sjeremylt   @brief Get total number of quadrature points (in dim dimensions)
536b11c1e72Sjeremylt 
537b11c1e72Sjeremylt   @param basis   CeedBasis
538b11c1e72Sjeremylt   @param[out] Q  Number of quadrature points
539b11c1e72Sjeremylt 
540b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
541dfdf5a53Sjeremylt 
542dfdf5a53Sjeremylt   @ref Utility
543b11c1e72Sjeremylt **/
544d7b241e6Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) {
545b5cf12eeSjeremylt   *Q = CeedIntPow(basis->Q1d, basis->dim);
546d7b241e6Sjeremylt   return 0;
547d7b241e6Sjeremylt }
548d7b241e6Sjeremylt 
549b11c1e72Sjeremylt /**
550b11c1e72Sjeremylt   @brief Destroy a CeedBasis
551b11c1e72Sjeremylt 
552b11c1e72Sjeremylt   @param basis CeedBasis to destroy
553b11c1e72Sjeremylt 
554b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
555dfdf5a53Sjeremylt 
556dfdf5a53Sjeremylt   @ref Basic
557b11c1e72Sjeremylt **/
558d7b241e6Sjeremylt int CeedBasisDestroy(CeedBasis *basis) {
559d7b241e6Sjeremylt   int ierr;
560d7b241e6Sjeremylt 
561d7b241e6Sjeremylt   if (!*basis || --(*basis)->refcount > 0) return 0;
562d7b241e6Sjeremylt   if ((*basis)->Destroy) {
563d7b241e6Sjeremylt     ierr = (*basis)->Destroy(*basis); CeedChk(ierr);
564d7b241e6Sjeremylt   }
565d7b241e6Sjeremylt   ierr = CeedFree(&(*basis)->interp1d); CeedChk(ierr);
566d7b241e6Sjeremylt   ierr = CeedFree(&(*basis)->grad1d); CeedChk(ierr);
567d7b241e6Sjeremylt   ierr = CeedFree(&(*basis)->qref1d); CeedChk(ierr);
568d7b241e6Sjeremylt   ierr = CeedFree(&(*basis)->qweight1d); CeedChk(ierr);
569d7b241e6Sjeremylt   ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr);
570d7b241e6Sjeremylt   ierr = CeedFree(basis); CeedChk(ierr);
571d7b241e6Sjeremylt   return 0;
572d7b241e6Sjeremylt }
573d7b241e6Sjeremylt 
57433e6becaSjeremylt /// @cond DOXYGEN_SKIP
57533e6becaSjeremylt // Indicate that the quadrature points are colocated with the dofs
576d7b241e6Sjeremylt CeedBasis CEED_BASIS_COLOCATED = &ceed_basis_colocated;
57733e6becaSjeremylt /// @endcond
578d7b241e6Sjeremylt /// @}
579