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