xref: /libCEED/rust/libceed-sys/c-src/interface/ceed-basis.c (revision 33e6beca7b305028e4768f6c701e40497ce598ab)
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 ///
30d7b241e6Sjeremylt /// @defgroup CeedBasis CeedBasis: fully discrete finite element-like objects
31d7b241e6Sjeremylt /// @{
32d7b241e6Sjeremylt 
33d7b241e6Sjeremylt /// Create a tensor product basis for H^1 discretizations
34d7b241e6Sjeremylt ///
35d7b241e6Sjeremylt /// @param ceed   Ceed
36d7b241e6Sjeremylt /// @param dim    Topological dimension
37d7b241e6Sjeremylt /// @param ncomp  Number of field components (1 for scalar fields)
38d7b241e6Sjeremylt /// @param P1d    Number of nodes in one dimension
39d7b241e6Sjeremylt /// @param Q1d    Number of quadrature points in one dimension
40d7b241e6Sjeremylt /// @param interp1d Row-major Q1d × P1d matrix expressing the values of nodal
41d7b241e6Sjeremylt ///               basis functions at quadrature points
42d7b241e6Sjeremylt /// @param grad1d  Row-major Q1d × P1d matrix expressing derivatives of nodal
43d7b241e6Sjeremylt ///               basis functions at quadrature points
44d7b241e6Sjeremylt /// @param qref1d Array of length Q1d holding the locations of quadrature points
45d7b241e6Sjeremylt ///               on the 1D reference element [-1, 1]
46d7b241e6Sjeremylt /// @param qweight1d Array of length Q1d holding the quadrature weights on the
47d7b241e6Sjeremylt ///               reference element
48d7b241e6Sjeremylt /// @param[out] basis New basis
49d7b241e6Sjeremylt ///
50d7b241e6Sjeremylt /// @sa CeedBasisCreateTensorH1Lagrange()
51d7b241e6Sjeremylt int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt ncomp, CeedInt P1d,
52d7b241e6Sjeremylt                             CeedInt Q1d, const CeedScalar *interp1d,
53d7b241e6Sjeremylt                             const CeedScalar *grad1d, const CeedScalar *qref1d,
54d7b241e6Sjeremylt                             const CeedScalar *qweight1d, CeedBasis *basis) {
55d7b241e6Sjeremylt   int ierr;
56d7b241e6Sjeremylt 
57d7b241e6Sjeremylt   if (!ceed->BasisCreateTensorH1)
58d7b241e6Sjeremylt     return CeedError(ceed, 1, "Backend does not support BasisCreateTensorH1");
59d7b241e6Sjeremylt   ierr = CeedCalloc(1,basis); CeedChk(ierr);
60d7b241e6Sjeremylt   (*basis)->ceed = ceed;
61d7b241e6Sjeremylt   ceed->refcount++;
62d7b241e6Sjeremylt   (*basis)->refcount = 1;
63d7b241e6Sjeremylt   (*basis)->dim = dim;
64d7b241e6Sjeremylt   (*basis)->ncomp = ncomp;
65d7b241e6Sjeremylt   (*basis)->P1d = P1d;
66d7b241e6Sjeremylt   (*basis)->Q1d = Q1d;
67d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d,&(*basis)->qref1d); CeedChk(ierr);
68d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d,&(*basis)->qweight1d); CeedChk(ierr);
69d7b241e6Sjeremylt   memcpy((*basis)->qref1d, qref1d, Q1d*sizeof(qref1d[0]));
70d7b241e6Sjeremylt   memcpy((*basis)->qweight1d, qweight1d, Q1d*sizeof(qweight1d[0]));
71d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d*P1d,&(*basis)->interp1d); CeedChk(ierr);
72d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d*P1d,&(*basis)->grad1d); CeedChk(ierr);
73d7b241e6Sjeremylt   memcpy((*basis)->interp1d, interp1d, Q1d*P1d*sizeof(interp1d[0]));
74d7b241e6Sjeremylt   memcpy((*basis)->grad1d, grad1d, Q1d*P1d*sizeof(interp1d[0]));
75d7b241e6Sjeremylt   ierr = ceed->BasisCreateTensorH1(ceed, dim, P1d, Q1d, interp1d, grad1d, qref1d,
76d7b241e6Sjeremylt                                    qweight1d, *basis); CeedChk(ierr);
77d7b241e6Sjeremylt   return 0;
78d7b241e6Sjeremylt }
79d7b241e6Sjeremylt 
80d7b241e6Sjeremylt /// Create a tensor product Lagrange basis
81d7b241e6Sjeremylt ///
82d7b241e6Sjeremylt /// @param ceed Ceed
83d7b241e6Sjeremylt /// @param dim Topological dimension of element
84d7b241e6Sjeremylt /// @param ncomp Number of field components
85d7b241e6Sjeremylt /// @param P Number of Gauss-Lobatto nodes in one dimension.  The polynomial degree
86d7b241e6Sjeremylt ///     of the resulting Q_k element is k=P-1.
87d7b241e6Sjeremylt /// @param Q Number of quadrature points in one dimension.
88d7b241e6Sjeremylt /// @param qmode Distribution of the Q quadrature points (affects order of
89d7b241e6Sjeremylt ///     accuracy for the quadrature)
90d7b241e6Sjeremylt /// @param[out] basis New basis
91d7b241e6Sjeremylt ///
92d7b241e6Sjeremylt /// @sa CeedBasisCreateTensorH1()
93d7b241e6Sjeremylt int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt ncomp,
94d7b241e6Sjeremylt                                     CeedInt P, CeedInt Q,
95d7b241e6Sjeremylt                                     CeedQuadMode qmode, CeedBasis *basis) {
96d7b241e6Sjeremylt   // Allocate
97d7b241e6Sjeremylt   int ierr, i, j, k;
98d7b241e6Sjeremylt   CeedScalar c1, c2, c3, c4, dx, *nodes, *interp1d, *grad1d, *qref1d, *qweight1d;
99d7b241e6Sjeremylt   ierr = CeedCalloc(P*Q, &interp1d); CeedChk(ierr);
100d7b241e6Sjeremylt   ierr = CeedCalloc(P*Q, &grad1d); CeedChk(ierr);
101d7b241e6Sjeremylt   ierr = CeedCalloc(P, &nodes); CeedChk(ierr);
102d7b241e6Sjeremylt   ierr = CeedCalloc(Q, &qref1d); CeedChk(ierr);
103d7b241e6Sjeremylt   ierr = CeedCalloc(Q, &qweight1d); CeedChk(ierr);
104d7b241e6Sjeremylt   // Get Nodes and Weights
105d7b241e6Sjeremylt   ierr = CeedLobattoQuadrature(P, nodes, NULL); CeedChk(ierr);
106d7b241e6Sjeremylt   switch (qmode) {
107d7b241e6Sjeremylt   case CEED_GAUSS:
108d7b241e6Sjeremylt     ierr = CeedGaussQuadrature(Q, qref1d, qweight1d); CeedChk(ierr);
109d7b241e6Sjeremylt     break;
110d7b241e6Sjeremylt   case CEED_GAUSS_LOBATTO:
111d7b241e6Sjeremylt     ierr = CeedLobattoQuadrature(Q, qref1d, qweight1d); CeedChk(ierr);
112d7b241e6Sjeremylt     break;
113d7b241e6Sjeremylt   }
114d7b241e6Sjeremylt   // Build B, D matrix
115d7b241e6Sjeremylt   // Fornberg, 1998
116d7b241e6Sjeremylt   for (i = 0; i  < Q; i++) {
117d7b241e6Sjeremylt     c1 = 1.0;
118d7b241e6Sjeremylt     c3 = nodes[0] - qref1d[i];
119d7b241e6Sjeremylt     interp1d[i*P+0] = 1.0;
120d7b241e6Sjeremylt     for (j = 1; j < P; j++) {
121d7b241e6Sjeremylt       c2 = 1.0;
122d7b241e6Sjeremylt       c4 = c3;
123d7b241e6Sjeremylt       c3 = nodes[j] - qref1d[i];
124d7b241e6Sjeremylt       for (k = 0; k < j; k++) {
125d7b241e6Sjeremylt         dx = nodes[j] - nodes[k];
126d7b241e6Sjeremylt         c2 *= dx;
127d7b241e6Sjeremylt         if (k == j - 1) {
128d7b241e6Sjeremylt           grad1d[i*P + j] = c1*(interp1d[i*P + k] - c4*grad1d[i*P + k]) / c2;
129d7b241e6Sjeremylt           interp1d[i*P + j] = - c1*c4*interp1d[i*P + k] / c2;
130d7b241e6Sjeremylt         }
131d7b241e6Sjeremylt         grad1d[i*P + k] = (c3*grad1d[i*P + k] - interp1d[i*P + k]) / dx;
132d7b241e6Sjeremylt         interp1d[i*P + k] = c3*interp1d[i*P + k] / dx;
133d7b241e6Sjeremylt       }
134d7b241e6Sjeremylt       c1 = c2;
135d7b241e6Sjeremylt     }
136d7b241e6Sjeremylt   }
137d7b241e6Sjeremylt   //  // Pass to CeedBasisCreateTensorH1
138d7b241e6Sjeremylt   ierr = CeedBasisCreateTensorH1(ceed, dim, ncomp, P, Q, interp1d, grad1d, qref1d,
139d7b241e6Sjeremylt                                  qweight1d, basis); CeedChk(ierr);
140d7b241e6Sjeremylt   ierr = CeedFree(&interp1d); CeedChk(ierr);
141d7b241e6Sjeremylt   ierr = CeedFree(&grad1d); CeedChk(ierr);
142d7b241e6Sjeremylt   ierr = CeedFree(&nodes); CeedChk(ierr);
143d7b241e6Sjeremylt   ierr = CeedFree(&qref1d); CeedChk(ierr);
144d7b241e6Sjeremylt   ierr = CeedFree(&qweight1d); CeedChk(ierr);
145d7b241e6Sjeremylt   return 0;
146d7b241e6Sjeremylt }
147d7b241e6Sjeremylt 
148d7b241e6Sjeremylt /// Construct a Gauss-Legendre quadrature
149d7b241e6Sjeremylt ///
150d7b241e6Sjeremylt /// @param Q Number of quadrature points (integrates polynomials of degree 2*Q-1 exactly)
151d7b241e6Sjeremylt /// @param qref1d Array of length Q to hold the abscissa on [-1, 1]
152d7b241e6Sjeremylt /// @param qweight1d Array of length Q to hold the weights
153d7b241e6Sjeremylt /// @sa CeedLobattoQuadrature()
154d7b241e6Sjeremylt int CeedGaussQuadrature(CeedInt Q, CeedScalar *qref1d, CeedScalar *qweight1d) {
155d7b241e6Sjeremylt   // Allocate
156d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0*atan(1.0);
157d7b241e6Sjeremylt   // Build qref1d, qweight1d
158d7b241e6Sjeremylt   for (int i = 0; i <= Q/2; i++) {
159d7b241e6Sjeremylt     // Guess
160d7b241e6Sjeremylt     xi = cos(PI*(CeedScalar)(2*i+1)/((CeedScalar)(2*Q)));
161d7b241e6Sjeremylt     // Pn(xi)
162d7b241e6Sjeremylt     P0 = 1.0;
163d7b241e6Sjeremylt     P1 = xi;
164d7b241e6Sjeremylt     P2 = 0.0;
165d7b241e6Sjeremylt     for (int j = 2; j <= Q; j++) {
166d7b241e6Sjeremylt       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
167d7b241e6Sjeremylt       P0 = P1;
168d7b241e6Sjeremylt       P1 = P2;
169d7b241e6Sjeremylt     }
170d7b241e6Sjeremylt     // First Newton Step
171d7b241e6Sjeremylt     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
172d7b241e6Sjeremylt     xi = xi-P2/dP2;
173d7b241e6Sjeremylt     // Newton to convergence
174d7b241e6Sjeremylt     for (int k=0; k<100 && fabs(P2)>1e-15; k++) {
175d7b241e6Sjeremylt       P0 = 1.0;
176d7b241e6Sjeremylt       P1 = xi;
177d7b241e6Sjeremylt       for (int j = 2; j <= Q; j++) {
178d7b241e6Sjeremylt         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
179d7b241e6Sjeremylt         P0 = P1;
180d7b241e6Sjeremylt         P1 = P2;
181d7b241e6Sjeremylt       }
182d7b241e6Sjeremylt       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
183d7b241e6Sjeremylt       xi = xi-P2/dP2;
184d7b241e6Sjeremylt     }
185d7b241e6Sjeremylt     // Save xi, wi
186d7b241e6Sjeremylt     wi = 2.0/((1.0-xi*xi)*dP2*dP2);
187d7b241e6Sjeremylt     qweight1d[i] = wi;
188d7b241e6Sjeremylt     qweight1d[Q-1-i] = wi;
189d7b241e6Sjeremylt     qref1d[i] = -xi;
190d7b241e6Sjeremylt     qref1d[Q-1-i]= xi;
191d7b241e6Sjeremylt   }
192d7b241e6Sjeremylt   return 0;
193d7b241e6Sjeremylt }
194d7b241e6Sjeremylt 
195d7b241e6Sjeremylt /// Construct a Gauss-Legendre-Lobatto quadrature
196d7b241e6Sjeremylt ///
197d7b241e6Sjeremylt /// @param Q Number of quadrature points (integrates polynomials of degree 2*Q-3 exactly)
198d7b241e6Sjeremylt /// @param qref1d Array of length Q to hold the abscissa on [-1, 1]
199d7b241e6Sjeremylt /// @param qweight1d Array of length Q to hold the weights
200d7b241e6Sjeremylt /// @sa CeedGaussQuadrature()
201d7b241e6Sjeremylt int CeedLobattoQuadrature(CeedInt Q, CeedScalar *qref1d,
202d7b241e6Sjeremylt                           CeedScalar *qweight1d) {
203d7b241e6Sjeremylt   // Allocate
204d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0*atan(1.0);
205d7b241e6Sjeremylt   // Build qref1d, qweight1d
206d7b241e6Sjeremylt   // Set endpoints
207d7b241e6Sjeremylt   wi = 2.0/((CeedScalar)(Q*(Q-1)));
208d7b241e6Sjeremylt   if (qweight1d) {
209d7b241e6Sjeremylt     qweight1d[0] = wi;
210d7b241e6Sjeremylt     qweight1d[Q-1] = wi;
211d7b241e6Sjeremylt   }
212d7b241e6Sjeremylt   qref1d[0] = -1.0;
213d7b241e6Sjeremylt   qref1d[Q-1] = 1.0;
214d7b241e6Sjeremylt   // Interior
215d7b241e6Sjeremylt   for (int i = 1; i <= (Q-1)/2; i++) {
216d7b241e6Sjeremylt     // Guess
217d7b241e6Sjeremylt     xi = cos(PI*(CeedScalar)(i)/(CeedScalar)(Q-1));
218d7b241e6Sjeremylt     // Pn(xi)
219d7b241e6Sjeremylt     P0 = 1.0;
220d7b241e6Sjeremylt     P1 = xi;
221d7b241e6Sjeremylt     P2 = 0.0;
222d7b241e6Sjeremylt     for (int j = 2; j < Q; j++) {
223d7b241e6Sjeremylt       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
224d7b241e6Sjeremylt       P0 = P1;
225d7b241e6Sjeremylt       P1 = P2;
226d7b241e6Sjeremylt     }
227d7b241e6Sjeremylt     // First Newton step
228d7b241e6Sjeremylt     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
229d7b241e6Sjeremylt     d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
230d7b241e6Sjeremylt     xi = xi-dP2/d2P2;
231d7b241e6Sjeremylt     // Newton to convergence
232d7b241e6Sjeremylt     for (int k=0; k<100 && fabs(dP2)>1e-15; k++) {
233d7b241e6Sjeremylt       P0 = 1.0;
234d7b241e6Sjeremylt       P1 = xi;
235d7b241e6Sjeremylt       for (int j = 2; j < Q; j++) {
236d7b241e6Sjeremylt         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
237d7b241e6Sjeremylt         P0 = P1;
238d7b241e6Sjeremylt         P1 = P2;
239d7b241e6Sjeremylt       }
240d7b241e6Sjeremylt       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
241d7b241e6Sjeremylt       d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
242d7b241e6Sjeremylt       xi = xi-dP2/d2P2;
243d7b241e6Sjeremylt     }
244d7b241e6Sjeremylt     // Save xi, wi
245d7b241e6Sjeremylt     wi = 2.0/(((CeedScalar)(Q*(Q-1)))*P2*P2);
246d7b241e6Sjeremylt     if (qweight1d) {
247d7b241e6Sjeremylt       qweight1d[i] = wi;
248d7b241e6Sjeremylt       qweight1d[Q-1-i] = wi;
249d7b241e6Sjeremylt     }
250d7b241e6Sjeremylt     qref1d[i] = -xi;
251d7b241e6Sjeremylt     qref1d[Q-1-i]= xi;
252d7b241e6Sjeremylt   }
253d7b241e6Sjeremylt   return 0;
254d7b241e6Sjeremylt }
255d7b241e6Sjeremylt 
256d7b241e6Sjeremylt static int CeedScalarView(const char *name, const char *fpformat, CeedInt m,
257d7b241e6Sjeremylt                           CeedInt n, const CeedScalar *a, FILE *stream) {
258d7b241e6Sjeremylt   for (int i=0; i<m; i++) {
259d7b241e6Sjeremylt     if (m > 1) fprintf(stream, "%12s[%d]:", name, i);
260d7b241e6Sjeremylt     else fprintf(stream, "%12s:", name);
261d7b241e6Sjeremylt     for (int j=0; j<n; j++) {
262d7b241e6Sjeremylt       fprintf(stream, fpformat, fabs(a[i*n+j]) > 1E-14 ? a[i*n+j] : 0);
263d7b241e6Sjeremylt     }
264d7b241e6Sjeremylt     fputs("\n", stream);
265d7b241e6Sjeremylt   }
266d7b241e6Sjeremylt   return 0;
267d7b241e6Sjeremylt }
268d7b241e6Sjeremylt 
269d7b241e6Sjeremylt /// View a basis
270d7b241e6Sjeremylt ///
271d7b241e6Sjeremylt /// @param basis Basis to view
272d7b241e6Sjeremylt /// @param stream Stream to view to, e.g., stdout
273d7b241e6Sjeremylt int CeedBasisView(CeedBasis basis, FILE *stream) {
274d7b241e6Sjeremylt   int ierr;
275d7b241e6Sjeremylt 
276d7b241e6Sjeremylt   fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P1d,
277d7b241e6Sjeremylt           basis->Q1d);
278d7b241e6Sjeremylt   ierr = CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q1d, basis->qref1d,
279d7b241e6Sjeremylt                         stream); CeedChk(ierr);
280d7b241e6Sjeremylt   ierr = CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q1d, basis->qweight1d,
281d7b241e6Sjeremylt                         stream); CeedChk(ierr);
282d7b241e6Sjeremylt   ierr = CeedScalarView("interp1d", "\t% 12.8f", basis->Q1d, basis->P1d,
283d7b241e6Sjeremylt                         basis->interp1d, stream); CeedChk(ierr);
284d7b241e6Sjeremylt   ierr = CeedScalarView("grad1d", "\t% 12.8f", basis->Q1d, basis->P1d,
285d7b241e6Sjeremylt                         basis->grad1d, stream); CeedChk(ierr);
286d7b241e6Sjeremylt   return 0;
287d7b241e6Sjeremylt }
288d7b241e6Sjeremylt 
289d7b241e6Sjeremylt // Computes A = (I - b v v^T) A
290d7b241e6Sjeremylt // where A is an mxn matrix indexed as A[i*row + j*col]
291d7b241e6Sjeremylt static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v,
292d7b241e6Sjeremylt                                   CeedScalar b, CeedInt m, CeedInt n,
293d7b241e6Sjeremylt                                   CeedInt row, CeedInt col) {
294d7b241e6Sjeremylt   for (CeedInt j=0; j<n; j++) {
295d7b241e6Sjeremylt     CeedScalar w = A[0*row + j*col];
296d7b241e6Sjeremylt     for (CeedInt i=1; i<m; i++) w += v[i] * A[i*row + j*col];
297d7b241e6Sjeremylt     A[0*row + j*col] -= b * w;
298d7b241e6Sjeremylt     for (CeedInt i=1; i<m; i++) A[i*row + j*col] -= b * w * v[i];
299d7b241e6Sjeremylt   }
300d7b241e6Sjeremylt   return 0;
301d7b241e6Sjeremylt }
302d7b241e6Sjeremylt 
303d7b241e6Sjeremylt // Compute A = Q A where Q is mxk and A is mxn. k<m
304d7b241e6Sjeremylt static int CeedHouseholderApplyQ(CeedScalar *A, const CeedScalar *Q,
305d7b241e6Sjeremylt                                  const CeedScalar *tau, CeedTransposeMode tmode,
306d7b241e6Sjeremylt                                  CeedInt m, CeedInt n, CeedInt k,
307d7b241e6Sjeremylt                                  CeedInt row, CeedInt col) {
308d7b241e6Sjeremylt   CeedScalar v[m];
309d7b241e6Sjeremylt   for (CeedInt ii=0; ii<k; ii++) {
310d7b241e6Sjeremylt     CeedInt i = tmode == CEED_TRANSPOSE ? ii : k-1-ii;
311d7b241e6Sjeremylt     for (CeedInt j=i+1; j<m; j++) {
312d7b241e6Sjeremylt       v[j] = Q[j*k+i];
313d7b241e6Sjeremylt     }
314d7b241e6Sjeremylt     // Apply Householder reflector (I - tau v v^T) colograd1d^T
315d7b241e6Sjeremylt     CeedHouseholderReflect(&A[i*row], &v[i], tau[i], m-i, n, row, col);
316d7b241e6Sjeremylt   }
317d7b241e6Sjeremylt   return 0;
318d7b241e6Sjeremylt }
319d7b241e6Sjeremylt 
320d7b241e6Sjeremylt /// Return QR Factorization of matrix
321d7b241e6Sjeremylt /// @param mat        Row-major matrix to be factorized in place
322d7b241e6Sjeremylt /// @param tau        Vector of length m of scaling fators
323d7b241e6Sjeremylt /// @param m          Number of rows
324d7b241e6Sjeremylt /// @param n          Number of columns
325d7b241e6Sjeremylt int CeedQRFactorization(CeedScalar *mat, CeedScalar *tau,
326d7b241e6Sjeremylt                         CeedInt m, CeedInt n) {
327d7b241e6Sjeremylt   CeedInt i, j;
328d7b241e6Sjeremylt   CeedScalar v[m];
329d7b241e6Sjeremylt 
330d7b241e6Sjeremylt   for (i=0; i<n; i++) {
331d7b241e6Sjeremylt     // Calculate Householder vector, magnitude
332d7b241e6Sjeremylt     CeedScalar sigma = 0.0;
333d7b241e6Sjeremylt     v[i] = mat[i+n*i];
334d7b241e6Sjeremylt     for (j=i+1; j<m; j++) {
335d7b241e6Sjeremylt       v[j] = mat[i+n*j];
336d7b241e6Sjeremylt       sigma += v[j] * v[j];
337d7b241e6Sjeremylt     }
338d7b241e6Sjeremylt     CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m]
339d7b241e6Sjeremylt     CeedScalar Rii = -copysign(norm, v[i]);
340d7b241e6Sjeremylt     v[i] -= Rii;
341d7b241e6Sjeremylt     // norm of v[i:m] after modification above and scaling below
342d7b241e6Sjeremylt     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
343d7b241e6Sjeremylt     //   tau = 2 / (norm*norm)
344d7b241e6Sjeremylt     tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma);
345d7b241e6Sjeremylt     for (j=i+1; j<m; j++) v[j] /= v[i];
346d7b241e6Sjeremylt 
347d7b241e6Sjeremylt     // Apply Householder reflector to lower right panel
348d7b241e6Sjeremylt     CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1);
349d7b241e6Sjeremylt     // Save v
350d7b241e6Sjeremylt     mat[i+n*i] = Rii;
351d7b241e6Sjeremylt     for (j=i+1; j<m; j++) {
352d7b241e6Sjeremylt       mat[i+n*j] = v[j];
353d7b241e6Sjeremylt     }
354d7b241e6Sjeremylt   }
355d7b241e6Sjeremylt 
356d7b241e6Sjeremylt   return 0;
357d7b241e6Sjeremylt }
358d7b241e6Sjeremylt 
359d7b241e6Sjeremylt /// Return colocated grad matrix
360d7b241e6Sjeremylt /// @param basis      Basis
361d7b241e6Sjeremylt /// @param colograd1d Row-major Q1d × Q1d matrix expressing derivatives of
362d7b241e6Sjeremylt ///                   basis functions at quadrature points
363d7b241e6Sjeremylt int CeedBasisGetColocatedGrad(CeedBasis basis, CeedScalar *colograd1d) {
364d7b241e6Sjeremylt   int i, j, k;
365d7b241e6Sjeremylt   CeedInt ierr, P1d=(basis)->P1d, Q1d=(basis)->Q1d;
366d7b241e6Sjeremylt   CeedScalar *interp1d, *grad1d, tau[Q1d];
367d7b241e6Sjeremylt 
368d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d*P1d, &interp1d); CeedChk(ierr);
369d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d*P1d, &grad1d); CeedChk(ierr);
370d7b241e6Sjeremylt   memcpy(interp1d, (basis)->interp1d, Q1d*P1d*sizeof(basis)->interp1d[0]);
371d7b241e6Sjeremylt   memcpy(grad1d, (basis)->grad1d, Q1d*P1d*sizeof(basis)->interp1d[0]);
372d7b241e6Sjeremylt 
373d7b241e6Sjeremylt   // QR Factorization, interp1d = Q R
374d7b241e6Sjeremylt   ierr = CeedQRFactorization(interp1d, tau, Q1d, P1d); CeedChk(ierr);
375d7b241e6Sjeremylt 
376d7b241e6Sjeremylt   // Apply Rinv, colograd1d = grad1d Rinv
377d7b241e6Sjeremylt   for (i=0; i<Q1d; i++) { // Row i
378d7b241e6Sjeremylt     colograd1d[Q1d*i] = grad1d[P1d*i]/interp1d[0];
379d7b241e6Sjeremylt     for (j=1; j<P1d; j++) { // Column j
380d7b241e6Sjeremylt       colograd1d[j+Q1d*i] = grad1d[j+P1d*i];
381d7b241e6Sjeremylt       for (k=0; k<j; k++) {
382d7b241e6Sjeremylt         colograd1d[j+Q1d*i] -= interp1d[j+P1d*k]*colograd1d[k+Q1d*i];
383d7b241e6Sjeremylt       }
384d7b241e6Sjeremylt       colograd1d[j+Q1d*i] /= interp1d[j+P1d*j];
385d7b241e6Sjeremylt     }
386d7b241e6Sjeremylt     for (j=P1d; j<Q1d; j++) {
387d7b241e6Sjeremylt       colograd1d[j+Q1d*i] = 0;
388d7b241e6Sjeremylt     }
389d7b241e6Sjeremylt   }
390d7b241e6Sjeremylt 
391d7b241e6Sjeremylt   // Apply Qtranspose, colograd = colograd Qtranspose
392d7b241e6Sjeremylt   CeedHouseholderApplyQ(colograd1d, interp1d, tau, CEED_NOTRANSPOSE,
393d7b241e6Sjeremylt                         Q1d, Q1d, P1d, 1, Q1d);
394d7b241e6Sjeremylt 
395d7b241e6Sjeremylt   ierr = CeedFree(&interp1d); CeedChk(ierr);
396d7b241e6Sjeremylt   ierr = CeedFree(&grad1d); CeedChk(ierr);
397d7b241e6Sjeremylt 
398d7b241e6Sjeremylt   return 0;
399d7b241e6Sjeremylt }
400d7b241e6Sjeremylt 
401d7b241e6Sjeremylt /// Apply basis evaluation from nodes to quadrature points or vice-versa
402d7b241e6Sjeremylt ///
403d7b241e6Sjeremylt /// @param basis Basis to evaluate
404d7b241e6Sjeremylt /// @param nelem the number of elements to apply the basis evaluation to;
405d7b241e6Sjeremylt ///     the backend will specify the ordering in ElemRestrictionCreateBlocked
406d7b241e6Sjeremylt /// @param tmode \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature
407d7b241e6Sjeremylt ///     points, \ref CEED_TRANSPOSE to apply the transpose, mapping from
408d7b241e6Sjeremylt ///     quadrature points to nodes
409d7b241e6Sjeremylt /// @param emode \ref CEED_EVAL_INTERP to obtain interpolated values,
410d7b241e6Sjeremylt ///     \ref CEED_EVAL_GRAD to obtain gradients.
411d7b241e6Sjeremylt /// @param u input vector
412d7b241e6Sjeremylt /// @param v output vector
413d7b241e6Sjeremylt int CeedBasisApply(CeedBasis basis, CeedInt nelem, CeedTransposeMode tmode,
414d7b241e6Sjeremylt                    CeedEvalMode emode, const CeedScalar *u, CeedScalar *v) {
415d7b241e6Sjeremylt   int ierr;
416d7b241e6Sjeremylt   if (!basis->Apply) return CeedError(basis->ceed, 1,
417d7b241e6Sjeremylt                                         "Backend does not support BasisApply");
418d7b241e6Sjeremylt   ierr = basis->Apply(basis, nelem, tmode, emode, u, v); CeedChk(ierr);
419d7b241e6Sjeremylt   return 0;
420d7b241e6Sjeremylt }
421d7b241e6Sjeremylt 
422d7b241e6Sjeremylt /// Get total number of nodes (in dim dimensions)
423d7b241e6Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) {
424d7b241e6Sjeremylt   *P = CeedPowInt(basis->P1d, basis->dim);
425d7b241e6Sjeremylt   return 0;
426d7b241e6Sjeremylt }
427d7b241e6Sjeremylt 
428d7b241e6Sjeremylt /// Get total number of quadrature points (in dim dimensions)
429d7b241e6Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) {
430d7b241e6Sjeremylt   *Q = CeedPowInt(basis->Q1d, basis->dim);
431d7b241e6Sjeremylt   return 0;
432d7b241e6Sjeremylt }
433d7b241e6Sjeremylt 
434d7b241e6Sjeremylt /// Destroy a CeedBasis
435d7b241e6Sjeremylt int CeedBasisDestroy(CeedBasis *basis) {
436d7b241e6Sjeremylt   int ierr;
437d7b241e6Sjeremylt 
438d7b241e6Sjeremylt   if (!*basis || --(*basis)->refcount > 0) return 0;
439d7b241e6Sjeremylt   if ((*basis)->Destroy) {
440d7b241e6Sjeremylt     ierr = (*basis)->Destroy(*basis); CeedChk(ierr);
441d7b241e6Sjeremylt   }
442d7b241e6Sjeremylt   ierr = CeedFree(&(*basis)->interp1d); CeedChk(ierr);
443d7b241e6Sjeremylt   ierr = CeedFree(&(*basis)->grad1d); CeedChk(ierr);
444d7b241e6Sjeremylt   ierr = CeedFree(&(*basis)->qref1d); CeedChk(ierr);
445d7b241e6Sjeremylt   ierr = CeedFree(&(*basis)->qweight1d); CeedChk(ierr);
446d7b241e6Sjeremylt   ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr);
447d7b241e6Sjeremylt   ierr = CeedFree(basis); CeedChk(ierr);
448d7b241e6Sjeremylt   return 0;
449d7b241e6Sjeremylt }
450d7b241e6Sjeremylt 
451*33e6becaSjeremylt /// @cond DOXYGEN_SKIP
452*33e6becaSjeremylt // Indicate that the quadrature points are colocated with the dofs
453d7b241e6Sjeremylt CeedBasis CEED_BASIS_COLOCATED = &ceed_basis_colocated;
454*33e6becaSjeremylt /// @endcond
455d7b241e6Sjeremylt /// @}
456