xref: /libCEED/interface/ceed-basis.c (revision a76a04e77a9db9309c4310fb27f031a0ea1f3c06)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3d7b241e6Sjeremylt //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5d7b241e6Sjeremylt //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7d7b241e6Sjeremylt 
8ec3da8bcSJed Brown #include <ceed/ceed.h>
9ec3da8bcSJed Brown #include <ceed/backend.h>
103d576824SJeremy L Thompson #include <ceed-impl.h>
11d7b241e6Sjeremylt #include <math.h>
123d576824SJeremy L Thompson #include <stdbool.h>
13d7b241e6Sjeremylt #include <stdio.h>
14d7b241e6Sjeremylt #include <string.h>
15d7b241e6Sjeremylt 
167a982d89SJeremy L. Thompson /// @file
177a982d89SJeremy L. Thompson /// Implementation of CeedBasis interfaces
187a982d89SJeremy L. Thompson 
19d7b241e6Sjeremylt /// @cond DOXYGEN_SKIP
20783c99b3SValeria Barra static struct CeedBasis_private ceed_basis_collocated;
21d7b241e6Sjeremylt /// @endcond
22d7b241e6Sjeremylt 
237a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser
247a982d89SJeremy L. Thompson /// @{
257a982d89SJeremy L. Thompson 
267a982d89SJeremy L. Thompson /// Indicate that the quadrature points are collocated with the nodes
277a982d89SJeremy L. Thompson const CeedBasis CEED_BASIS_COLLOCATED = &ceed_basis_collocated;
287a982d89SJeremy L. Thompson 
297a982d89SJeremy L. Thompson /// @}
307a982d89SJeremy L. Thompson 
317a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
327a982d89SJeremy L. Thompson /// CeedBasis Library Internal Functions
337a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
347a982d89SJeremy L. Thompson /// @addtogroup CeedBasisDeveloper
357a982d89SJeremy L. Thompson /// @{
367a982d89SJeremy L. Thompson 
377a982d89SJeremy L. Thompson /**
387a982d89SJeremy L. Thompson   @brief Compute Householder reflection
397a982d89SJeremy L. Thompson 
407a982d89SJeremy L. Thompson     Computes A = (I - b v v^T) A
417a982d89SJeremy L. Thompson     where A is an mxn matrix indexed as A[i*row + j*col]
427a982d89SJeremy L. Thompson 
437a982d89SJeremy L. Thompson   @param[in,out] A  Matrix to apply Householder reflection to, in place
447a982d89SJeremy L. Thompson   @param v          Householder vector
457a982d89SJeremy L. Thompson   @param b          Scaling factor
467a982d89SJeremy L. Thompson   @param m          Number of rows in A
477a982d89SJeremy L. Thompson   @param n          Number of columns in A
487a982d89SJeremy L. Thompson   @param row        Row stride
497a982d89SJeremy L. Thompson   @param col        Col stride
507a982d89SJeremy L. Thompson 
517a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
527a982d89SJeremy L. Thompson 
537a982d89SJeremy L. Thompson   @ref Developer
547a982d89SJeremy L. Thompson **/
557a982d89SJeremy L. Thompson static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v,
567a982d89SJeremy L. Thompson                                   CeedScalar b, CeedInt m, CeedInt n,
577a982d89SJeremy L. Thompson                                   CeedInt row, CeedInt col) {
587a982d89SJeremy L. Thompson   for (CeedInt j=0; j<n; j++) {
597a982d89SJeremy L. Thompson     CeedScalar w = A[0*row + j*col];
607a982d89SJeremy L. Thompson     for (CeedInt i=1; i<m; i++)
617a982d89SJeremy L. Thompson       w += v[i] * A[i*row + j*col];
627a982d89SJeremy L. Thompson     A[0*row + j*col] -= b * w;
637a982d89SJeremy L. Thompson     for (CeedInt i=1; i<m; i++)
647a982d89SJeremy L. Thompson       A[i*row + j*col] -= b * w * v[i];
657a982d89SJeremy L. Thompson   }
66e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
677a982d89SJeremy L. Thompson }
687a982d89SJeremy L. Thompson 
697a982d89SJeremy L. Thompson /**
707a982d89SJeremy L. Thompson   @brief Apply Householder Q matrix
717a982d89SJeremy L. Thompson 
727a982d89SJeremy L. Thompson     Compute A = Q A where Q is mxm and A is mxn.
737a982d89SJeremy L. Thompson 
747a982d89SJeremy L. Thompson   @param[in,out] A  Matrix to apply Householder Q to, in place
757a982d89SJeremy L. Thompson   @param Q          Householder Q matrix
767a982d89SJeremy L. Thompson   @param tau        Householder scaling factors
77d1d35e2fSjeremylt   @param t_mode     Transpose mode for application
787a982d89SJeremy L. Thompson   @param m          Number of rows in A
797a982d89SJeremy L. Thompson   @param n          Number of columns in A
807a982d89SJeremy L. Thompson   @param k          Number of elementary reflectors in Q, k<m
817a982d89SJeremy L. Thompson   @param row        Row stride in A
827a982d89SJeremy L. Thompson   @param col        Col stride in A
837a982d89SJeremy L. Thompson 
847a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
857a982d89SJeremy L. Thompson 
867a982d89SJeremy L. Thompson   @ref Developer
877a982d89SJeremy L. Thompson **/
88d99fa3c5SJeremy L Thompson int CeedHouseholderApplyQ(CeedScalar *A, const CeedScalar *Q,
89d1d35e2fSjeremylt                           const CeedScalar *tau, CeedTransposeMode t_mode,
907a982d89SJeremy L. Thompson                           CeedInt m, CeedInt n, CeedInt k,
917a982d89SJeremy L. Thompson                           CeedInt row, CeedInt col) {
92e15f9bd0SJeremy L Thompson   int ierr;
9378464608Sjeremylt   CeedScalar *v;
9478464608Sjeremylt   ierr = CeedMalloc(m, &v); CeedChk(ierr);
957a982d89SJeremy L. Thompson   for (CeedInt ii=0; ii<k; ii++) {
96d1d35e2fSjeremylt     CeedInt i = t_mode == CEED_TRANSPOSE ? ii : k-1-ii;
977a982d89SJeremy L. Thompson     for (CeedInt j=i+1; j<m; j++)
987a982d89SJeremy L. Thompson       v[j] = Q[j*k+i];
99d1d35e2fSjeremylt     // Apply Householder reflector (I - tau v v^T) collo_grad_1d^T
100e15f9bd0SJeremy L Thompson     ierr = CeedHouseholderReflect(&A[i*row], &v[i], tau[i], m-i, n, row, col);
101e15f9bd0SJeremy L Thompson     CeedChk(ierr);
1027a982d89SJeremy L. Thompson   }
10378464608Sjeremylt   ierr = CeedFree(&v); CeedChk(ierr);
104e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1057a982d89SJeremy L. Thompson }
1067a982d89SJeremy L. Thompson 
1077a982d89SJeremy L. Thompson /**
1087a982d89SJeremy L. Thompson   @brief Compute Givens rotation
1097a982d89SJeremy L. Thompson 
1107a982d89SJeremy L. Thompson     Computes A = G A (or G^T A in transpose mode)
1117a982d89SJeremy L. Thompson     where A is an mxn matrix indexed as A[i*n + j*m]
1127a982d89SJeremy L. Thompson 
1137a982d89SJeremy L. Thompson   @param[in,out] A  Row major matrix to apply Givens rotation to, in place
1147a982d89SJeremy L. Thompson   @param c          Cosine factor
1157a982d89SJeremy L. Thompson   @param s          Sine factor
116d1d35e2fSjeremylt   @param t_mode     @ref CEED_NOTRANSPOSE to rotate the basis counter-clockwise,
1174c4400c7SValeria Barra                     which has the effect of rotating columns of A clockwise;
1184cc79fe7SJed Brown                     @ref CEED_TRANSPOSE for the opposite rotation
1197a982d89SJeremy L. Thompson   @param i          First row/column to apply rotation
1207a982d89SJeremy L. Thompson   @param k          Second row/column to apply rotation
1217a982d89SJeremy L. Thompson   @param m          Number of rows in A
1227a982d89SJeremy L. Thompson   @param n          Number of columns in A
1237a982d89SJeremy L. Thompson 
1247a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
1257a982d89SJeremy L. Thompson 
1267a982d89SJeremy L. Thompson   @ref Developer
1277a982d89SJeremy L. Thompson **/
1287a982d89SJeremy L. Thompson static int CeedGivensRotation(CeedScalar *A, CeedScalar c, CeedScalar s,
129d1d35e2fSjeremylt                               CeedTransposeMode t_mode, CeedInt i, CeedInt k,
1307a982d89SJeremy L. Thompson                               CeedInt m, CeedInt n) {
131d1d35e2fSjeremylt   CeedInt stride_j = 1, stride_ik = m, num_its = n;
132d1d35e2fSjeremylt   if (t_mode == CEED_NOTRANSPOSE) {
133d1d35e2fSjeremylt     stride_j = n; stride_ik = 1; num_its = m;
1347a982d89SJeremy L. Thompson   }
1357a982d89SJeremy L. Thompson 
1367a982d89SJeremy L. Thompson   // Apply rotation
137d1d35e2fSjeremylt   for (CeedInt j=0; j<num_its; j++) {
138d1d35e2fSjeremylt     CeedScalar tau1 = A[i*stride_ik+j*stride_j], tau2 = A[k*stride_ik+j*stride_j];
139d1d35e2fSjeremylt     A[i*stride_ik+j*stride_j] = c*tau1 - s*tau2;
140d1d35e2fSjeremylt     A[k*stride_ik+j*stride_j] = s*tau1 + c*tau2;
1417a982d89SJeremy L. Thompson   }
142e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1437a982d89SJeremy L. Thompson }
1447a982d89SJeremy L. Thompson 
1457a982d89SJeremy L. Thompson /**
1467a982d89SJeremy L. Thompson   @brief View an array stored in a CeedBasis
1477a982d89SJeremy L. Thompson 
1480a0da059Sjeremylt   @param[in] name      Name of array
149d1d35e2fSjeremylt   @param[in] fp_fmt    Printing format
1500a0da059Sjeremylt   @param[in] m         Number of rows in array
1510a0da059Sjeremylt   @param[in] n         Number of columns in array
1520a0da059Sjeremylt   @param[in] a         Array to be viewed
1530a0da059Sjeremylt   @param[in] stream    Stream to view to, e.g., stdout
1547a982d89SJeremy L. Thompson 
1557a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
1567a982d89SJeremy L. Thompson 
1577a982d89SJeremy L. Thompson   @ref Developer
1587a982d89SJeremy L. Thompson **/
159d1d35e2fSjeremylt static int CeedScalarView(const char *name, const char *fp_fmt, CeedInt m,
1607a982d89SJeremy L. Thompson                           CeedInt n, const CeedScalar *a, FILE *stream) {
16192ae7e47SJeremy L Thompson   for (CeedInt i=0; i<m; i++) {
1627a982d89SJeremy L. Thompson     if (m > 1)
163990fdeb6SJeremy L Thompson       fprintf(stream, "%12s[%" CeedInt_FMT "]:", name, i);
1647a982d89SJeremy L. Thompson     else
1657a982d89SJeremy L. Thompson       fprintf(stream, "%12s:", name);
16692ae7e47SJeremy L Thompson     for (CeedInt j=0; j<n; j++)
167d1d35e2fSjeremylt       fprintf(stream, fp_fmt, fabs(a[i*n+j]) > 1E-14 ? a[i*n+j] : 0);
1687a982d89SJeremy L. Thompson     fputs("\n", stream);
1697a982d89SJeremy L. Thompson   }
170e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1717a982d89SJeremy L. Thompson }
1727a982d89SJeremy L. Thompson 
173*a76a04e7SJeremy L Thompson /**
174*a76a04e7SJeremy L Thompson   @brief Create the interpolation matrix for projection from the nodes of
175*a76a04e7SJeremy L Thompson            `basis_from` to the nodes of `basis_to`. This projection is
176*a76a04e7SJeremy L Thompson            given by `interp_project = interp_to^+ * interp_from`, where
177*a76a04e7SJeremy L Thompson            the pesudoinverse `interp_to^+` is given by QR factorization.
178*a76a04e7SJeremy L Thompson          Note: `basis_from` and `basis_to` must have compatible quadrature
179*a76a04e7SJeremy L Thompson            spaces.
180*a76a04e7SJeremy L Thompson 
181*a76a04e7SJeremy L Thompson   @param[in] basis_from       CeedBasis to project from
182*a76a04e7SJeremy L Thompson   @param[in] basis_to         CeedBasis to project to
183*a76a04e7SJeremy L Thompson   @param[out] interp_project  Address of the variable where the newly created
184*a76a04e7SJeremy L Thompson                                 projection matrix will be stored.
185*a76a04e7SJeremy L Thompson 
186*a76a04e7SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
187*a76a04e7SJeremy L Thompson 
188*a76a04e7SJeremy L Thompson   @ref Developer
189*a76a04e7SJeremy L Thompson **/
190*a76a04e7SJeremy L Thompson static int CeedBasisCreateProjectionMatrix(CeedBasis basis_from,
191*a76a04e7SJeremy L Thompson                                     CeedBasis basis_to,
192*a76a04e7SJeremy L Thompson                                     CeedScalar **interp_project) {
193*a76a04e7SJeremy L Thompson   int ierr;
194*a76a04e7SJeremy L Thompson   Ceed ceed;
195*a76a04e7SJeremy L Thompson   ierr = CeedBasisGetCeed(basis_to, &ceed); CeedChk(ierr);
196*a76a04e7SJeremy L Thompson 
197*a76a04e7SJeremy L Thompson   // Check for compatible quadrature spaces
198*a76a04e7SJeremy L Thompson   CeedInt Q_to, Q_from;
199*a76a04e7SJeremy L Thompson   ierr = CeedBasisGetNumQuadraturePoints(basis_to, &Q_to); CeedChk(ierr);
200*a76a04e7SJeremy L Thompson   ierr = CeedBasisGetNumQuadraturePoints(basis_from, &Q_from); CeedChk(ierr);
201*a76a04e7SJeremy L Thompson   if (Q_to != Q_from)
202*a76a04e7SJeremy L Thompson     // LCOV_EXCL_START
203*a76a04e7SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
204*a76a04e7SJeremy L Thompson                      "Bases must have compatible quadrature spaces");
205*a76a04e7SJeremy L Thompson   // LCOV_EXCL_STOP
206*a76a04e7SJeremy L Thompson 
207*a76a04e7SJeremy L Thompson   // Coarse to fine basis
208*a76a04e7SJeremy L Thompson   CeedInt P_to, P_from, Q = Q_to;
209*a76a04e7SJeremy L Thompson   bool is_tensor_to, is_tensor_from;
210*a76a04e7SJeremy L Thompson   ierr = CeedBasisIsTensor(basis_to, &is_tensor_to); CeedChk(ierr);
211*a76a04e7SJeremy L Thompson   ierr = CeedBasisIsTensor(basis_from, &is_tensor_from); CeedChk(ierr);
212*a76a04e7SJeremy L Thompson   CeedScalar *interp_to, *interp_from, *tau;
213*a76a04e7SJeremy L Thompson   if (is_tensor_to && is_tensor_from) {
214*a76a04e7SJeremy L Thompson     ierr = CeedBasisGetNumNodes1D(basis_to, &P_to); CeedChk(ierr);
215*a76a04e7SJeremy L Thompson     ierr = CeedBasisGetNumNodes1D(basis_from, &P_from); CeedChk(ierr);
216*a76a04e7SJeremy L Thompson     ierr = CeedBasisGetNumQuadraturePoints1D(basis_from, &Q); CeedChk(ierr);
217*a76a04e7SJeremy L Thompson   } else if (!is_tensor_to && !is_tensor_from) {
218*a76a04e7SJeremy L Thompson     ierr = CeedBasisGetNumNodes(basis_to, &P_to); CeedChk(ierr);
219*a76a04e7SJeremy L Thompson     ierr = CeedBasisGetNumNodes(basis_from, &P_from); CeedChk(ierr);
220*a76a04e7SJeremy L Thompson   } else {
221*a76a04e7SJeremy L Thompson     // LCOV_EXCL_START
222*a76a04e7SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_MINOR,
223*a76a04e7SJeremy L Thompson                      "Bases must both be tensor or non-tensor");
224*a76a04e7SJeremy L Thompson     // LCOV_EXCL_STOP
225*a76a04e7SJeremy L Thompson   }
226*a76a04e7SJeremy L Thompson 
227*a76a04e7SJeremy L Thompson   ierr = CeedMalloc(Q * P_from, &interp_from); CeedChk(ierr);
228*a76a04e7SJeremy L Thompson   ierr = CeedMalloc(Q * P_to, &interp_to); CeedChk(ierr);
229*a76a04e7SJeremy L Thompson   ierr = CeedCalloc(P_to * P_from, interp_project); CeedChk(ierr);
230*a76a04e7SJeremy L Thompson   ierr = CeedMalloc(Q, &tau); CeedChk(ierr);
231*a76a04e7SJeremy L Thompson   const CeedScalar *interp_to_source = NULL, *interp_from_source = NULL;
232*a76a04e7SJeremy L Thompson   if (is_tensor_to) {
233*a76a04e7SJeremy L Thompson     ierr = CeedBasisGetInterp1D(basis_to, &interp_to_source); CeedChk(ierr);
234*a76a04e7SJeremy L Thompson     ierr = CeedBasisGetInterp1D(basis_from, &interp_from_source); CeedChk(ierr);
235*a76a04e7SJeremy L Thompson   } else {
236*a76a04e7SJeremy L Thompson     ierr = CeedBasisGetInterp(basis_to, &interp_to_source); CeedChk(ierr);
237*a76a04e7SJeremy L Thompson     ierr = CeedBasisGetInterp(basis_from, &interp_from_source); CeedChk(ierr);
238*a76a04e7SJeremy L Thompson   }
239*a76a04e7SJeremy L Thompson   memcpy(interp_to, interp_to_source, Q * P_to * sizeof(interp_to_source[0]));
240*a76a04e7SJeremy L Thompson   memcpy(interp_from, interp_from_source,
241*a76a04e7SJeremy L Thompson          Q * P_from * sizeof(interp_from_source[0]));
242*a76a04e7SJeremy L Thompson 
243*a76a04e7SJeremy L Thompson   // -- QR Factorization, interp_to = Q R
244*a76a04e7SJeremy L Thompson   ierr = CeedQRFactorization(ceed, interp_to, tau, Q, P_to); CeedChk(ierr);
245*a76a04e7SJeremy L Thompson 
246*a76a04e7SJeremy L Thompson   // -- Apply Qtranspose, interp_to = Qtranspose interp_from
247*a76a04e7SJeremy L Thompson   ierr = CeedHouseholderApplyQ(interp_from, interp_to, tau, CEED_TRANSPOSE,
248*a76a04e7SJeremy L Thompson                                Q, P_from, P_to, P_from, 1); CeedChk(ierr);
249*a76a04e7SJeremy L Thompson 
250*a76a04e7SJeremy L Thompson   // -- Apply Rinv, interp_project = Rinv interp_c
251*a76a04e7SJeremy L Thompson   for (CeedInt j = 0; j < P_from; j++) { // Column j
252*a76a04e7SJeremy L Thompson     (*interp_project)[j + P_from * (P_to - 1)] = interp_from[j + P_from *
253*a76a04e7SJeremy L Thompson         (P_to - 1)] / interp_to[P_to * P_to - 1];
254*a76a04e7SJeremy L Thompson     for (CeedInt i = P_to - 2; i >= 0; i--) { // Row i
255*a76a04e7SJeremy L Thompson       (*interp_project)[j + P_from * i] = interp_from[j + P_from * i];
256*a76a04e7SJeremy L Thompson       for (CeedInt k = i+1; k < P_to; k++) {
257*a76a04e7SJeremy L Thompson         (*interp_project)[j + P_from * i] -= interp_to[k + P_to * i]*
258*a76a04e7SJeremy L Thompson                                              (*interp_project)[j + P_from * k];
259*a76a04e7SJeremy L Thompson       }
260*a76a04e7SJeremy L Thompson       (*interp_project)[j + P_from * i] /= interp_to[i + P_to * i];
261*a76a04e7SJeremy L Thompson     }
262*a76a04e7SJeremy L Thompson   }
263*a76a04e7SJeremy L Thompson   ierr = CeedFree(&tau); CeedChk(ierr);
264*a76a04e7SJeremy L Thompson   ierr = CeedFree(&interp_to); CeedChk(ierr);
265*a76a04e7SJeremy L Thompson   ierr = CeedFree(&interp_from); CeedChk(ierr);
266*a76a04e7SJeremy L Thompson 
267*a76a04e7SJeremy L Thompson   return CEED_ERROR_SUCCESS;
268*a76a04e7SJeremy L Thompson }
269*a76a04e7SJeremy L Thompson 
2707a982d89SJeremy L. Thompson /// @}
2717a982d89SJeremy L. Thompson 
2727a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
2737a982d89SJeremy L. Thompson /// Ceed Backend API
2747a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
2757a982d89SJeremy L. Thompson /// @addtogroup CeedBasisBackend
2767a982d89SJeremy L. Thompson /// @{
2777a982d89SJeremy L. Thompson 
2787a982d89SJeremy L. Thompson /**
2797a982d89SJeremy L. Thompson   @brief Return collocated grad matrix
2807a982d89SJeremy L. Thompson 
2817a982d89SJeremy L. Thompson   @param basis               CeedBasis
282d1d35e2fSjeremylt   @param[out] collo_grad_1d  Row-major (Q_1d * Q_1d) matrix expressing derivatives of
2837a982d89SJeremy L. Thompson                                basis functions at quadrature points
2847a982d89SJeremy L. Thompson 
2857a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
2867a982d89SJeremy L. Thompson 
2877a982d89SJeremy L. Thompson   @ref Backend
2887a982d89SJeremy L. Thompson **/
289d1d35e2fSjeremylt int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *collo_grad_1d) {
2907a982d89SJeremy L. Thompson   int i, j, k;
2917a982d89SJeremy L. Thompson   Ceed ceed;
292d1d35e2fSjeremylt   CeedInt ierr, P_1d=(basis)->P_1d, Q_1d=(basis)->Q_1d;
29378464608Sjeremylt   CeedScalar *interp_1d, *grad_1d, *tau;
2947a982d89SJeremy L. Thompson 
295d1d35e2fSjeremylt   ierr = CeedMalloc(Q_1d*P_1d, &interp_1d); CeedChk(ierr);
296d1d35e2fSjeremylt   ierr = CeedMalloc(Q_1d*P_1d, &grad_1d); CeedChk(ierr);
29778464608Sjeremylt   ierr = CeedMalloc(Q_1d, &tau); CeedChk(ierr);
298d1d35e2fSjeremylt   memcpy(interp_1d, (basis)->interp_1d, Q_1d*P_1d*sizeof(basis)->interp_1d[0]);
299d1d35e2fSjeremylt   memcpy(grad_1d, (basis)->grad_1d, Q_1d*P_1d*sizeof(basis)->interp_1d[0]);
3007a982d89SJeremy L. Thompson 
301d1d35e2fSjeremylt   // QR Factorization, interp_1d = Q R
3027a982d89SJeremy L. Thompson   ierr = CeedBasisGetCeed(basis, &ceed); CeedChk(ierr);
303d1d35e2fSjeremylt   ierr = CeedQRFactorization(ceed, interp_1d, tau, Q_1d, P_1d); CeedChk(ierr);
304e15f9bd0SJeremy L Thompson   // Note: This function is for backend use, so all errors are terminal
305e15f9bd0SJeremy L Thompson   //   and we do not need to clean up memory on failure.
3067a982d89SJeremy L. Thompson 
307d1d35e2fSjeremylt   // Apply Rinv, collo_grad_1d = grad_1d Rinv
308d1d35e2fSjeremylt   for (i=0; i<Q_1d; i++) { // Row i
309d1d35e2fSjeremylt     collo_grad_1d[Q_1d*i] = grad_1d[P_1d*i]/interp_1d[0];
310d1d35e2fSjeremylt     for (j=1; j<P_1d; j++) { // Column j
311d1d35e2fSjeremylt       collo_grad_1d[j+Q_1d*i] = grad_1d[j+P_1d*i];
3127a982d89SJeremy L. Thompson       for (k=0; k<j; k++)
313d1d35e2fSjeremylt         collo_grad_1d[j+Q_1d*i] -= interp_1d[j+P_1d*k]*collo_grad_1d[k+Q_1d*i];
314d1d35e2fSjeremylt       collo_grad_1d[j+Q_1d*i] /= interp_1d[j+P_1d*j];
3157a982d89SJeremy L. Thompson     }
316d1d35e2fSjeremylt     for (j=P_1d; j<Q_1d; j++)
317d1d35e2fSjeremylt       collo_grad_1d[j+Q_1d*i] = 0;
3187a982d89SJeremy L. Thompson   }
3197a982d89SJeremy L. Thompson 
320673160d7Sjeremylt   // Apply Qtranspose, collo_grad = collo_grad Q_transpose
321d1d35e2fSjeremylt   ierr = CeedHouseholderApplyQ(collo_grad_1d, interp_1d, tau, CEED_NOTRANSPOSE,
322d1d35e2fSjeremylt                                Q_1d, Q_1d, P_1d, 1, Q_1d); CeedChk(ierr);
3237a982d89SJeremy L. Thompson 
324d1d35e2fSjeremylt   ierr = CeedFree(&interp_1d); CeedChk(ierr);
325d1d35e2fSjeremylt   ierr = CeedFree(&grad_1d); CeedChk(ierr);
32678464608Sjeremylt   ierr = CeedFree(&tau); CeedChk(ierr);
327e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3287a982d89SJeremy L. Thompson }
3297a982d89SJeremy L. Thompson 
3307a982d89SJeremy L. Thompson /**
3317a982d89SJeremy L. Thompson   @brief Get tensor status for given CeedBasis
3327a982d89SJeremy L. Thompson 
3337a982d89SJeremy L. Thompson   @param basis           CeedBasis
334d1d35e2fSjeremylt   @param[out] is_tensor  Variable to store tensor status
3357a982d89SJeremy L. Thompson 
3367a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3377a982d89SJeremy L. Thompson 
3387a982d89SJeremy L. Thompson   @ref Backend
3397a982d89SJeremy L. Thompson **/
340d1d35e2fSjeremylt int CeedBasisIsTensor(CeedBasis basis, bool *is_tensor) {
341d1d35e2fSjeremylt   *is_tensor = basis->tensor_basis;
342e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3437a982d89SJeremy L. Thompson }
3447a982d89SJeremy L. Thompson 
3457a982d89SJeremy L. Thompson /**
3467a982d89SJeremy L. Thompson   @brief Get backend data of a CeedBasis
3477a982d89SJeremy L. Thompson 
3487a982d89SJeremy L. Thompson   @param basis      CeedBasis
3497a982d89SJeremy L. Thompson   @param[out] data  Variable to store data
3507a982d89SJeremy L. Thompson 
3517a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3527a982d89SJeremy L. Thompson 
3537a982d89SJeremy L. Thompson   @ref Backend
3547a982d89SJeremy L. Thompson **/
355777ff853SJeremy L Thompson int CeedBasisGetData(CeedBasis basis, void *data) {
356777ff853SJeremy L Thompson   *(void **)data = basis->data;
357e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3587a982d89SJeremy L. Thompson }
3597a982d89SJeremy L. Thompson 
3607a982d89SJeremy L. Thompson /**
3617a982d89SJeremy L. Thompson   @brief Set backend data of a CeedBasis
3627a982d89SJeremy L. Thompson 
3637a982d89SJeremy L. Thompson   @param[out] basis  CeedBasis
3647a982d89SJeremy L. Thompson   @param data        Data to set
3657a982d89SJeremy L. Thompson 
3667a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3677a982d89SJeremy L. Thompson 
3687a982d89SJeremy L. Thompson   @ref Backend
3697a982d89SJeremy L. Thompson **/
370777ff853SJeremy L Thompson int CeedBasisSetData(CeedBasis basis, void *data) {
371777ff853SJeremy L Thompson   basis->data = data;
372e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3737a982d89SJeremy L. Thompson }
3747a982d89SJeremy L. Thompson 
3757a982d89SJeremy L. Thompson /**
37634359f16Sjeremylt   @brief Increment the reference counter for a CeedBasis
37734359f16Sjeremylt 
37834359f16Sjeremylt   @param basis  Basis to increment the reference counter
37934359f16Sjeremylt 
38034359f16Sjeremylt   @return An error code: 0 - success, otherwise - failure
38134359f16Sjeremylt 
38234359f16Sjeremylt   @ref Backend
38334359f16Sjeremylt **/
3849560d06aSjeremylt int CeedBasisReference(CeedBasis basis) {
38534359f16Sjeremylt   basis->ref_count++;
38634359f16Sjeremylt   return CEED_ERROR_SUCCESS;
38734359f16Sjeremylt }
38834359f16Sjeremylt 
38934359f16Sjeremylt /**
3906e15d496SJeremy L Thompson   @brief Estimate number of FLOPs required to apply CeedBasis in t_mode and eval_mode
3916e15d496SJeremy L Thompson 
3926e15d496SJeremy L Thompson   @param basis     Basis to estimate FLOPs for
3936e15d496SJeremy L Thompson   @param t_mode    Apply basis or transpose
3946e15d496SJeremy L Thompson   @param eval_mode Basis evaluation mode
3956e15d496SJeremy L Thompson   @param flops     Address of variable to hold FLOPs estimate
3966e15d496SJeremy L Thompson 
3976e15d496SJeremy L Thompson   @ref Backend
3986e15d496SJeremy L Thompson **/
3996e15d496SJeremy L Thompson int CeedBasisGetFlopsEstimate(CeedBasis basis, CeedTransposeMode t_mode,
4009d36ca50SJeremy L Thompson                               CeedEvalMode eval_mode, CeedSize *flops) {
4016e15d496SJeremy L Thompson   int ierr;
4026e15d496SJeremy L Thompson   bool is_tensor;
4036e15d496SJeremy L Thompson 
4046e15d496SJeremy L Thompson   ierr = CeedBasisIsTensor(basis, &is_tensor); CeedChk(ierr);
4056e15d496SJeremy L Thompson   if (is_tensor) {
4066e15d496SJeremy L Thompson     CeedInt dim, num_comp, P_1d, Q_1d;
4076e15d496SJeremy L Thompson     ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr);
4086e15d496SJeremy L Thompson     ierr = CeedBasisGetNumComponents(basis, &num_comp); CeedChk(ierr);
4096e15d496SJeremy L Thompson     ierr = CeedBasisGetNumNodes1D(basis, &P_1d);  CeedChk(ierr);
4106e15d496SJeremy L Thompson     ierr = CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d);  CeedChk(ierr);
4116e15d496SJeremy L Thompson     if (t_mode == CEED_TRANSPOSE) {
4126e15d496SJeremy L Thompson       P_1d = Q_1d; Q_1d = P_1d;
4136e15d496SJeremy L Thompson     }
4146e15d496SJeremy L Thompson     CeedInt tensor_flops = 0, pre = num_comp * CeedIntPow(P_1d, dim-1), post = 1;
4156e15d496SJeremy L Thompson     for (CeedInt d = 0; d < dim; d++) {
4166e15d496SJeremy L Thompson       tensor_flops += 2 * pre * P_1d * post * Q_1d;
4176e15d496SJeremy L Thompson       pre /= P_1d;
4186e15d496SJeremy L Thompson       post *= Q_1d;
4196e15d496SJeremy L Thompson     }
4206e15d496SJeremy L Thompson     switch (eval_mode) {
4216e15d496SJeremy L Thompson     case CEED_EVAL_NONE:   *flops = 0; break;
4226e15d496SJeremy L Thompson     case CEED_EVAL_INTERP: *flops = tensor_flops; break;
4236e15d496SJeremy L Thompson     case CEED_EVAL_GRAD:   *flops = tensor_flops * 2; break;
4246e15d496SJeremy L Thompson     case CEED_EVAL_DIV:
4256e15d496SJeremy L Thompson       // LCOV_EXCL_START
4266e15d496SJeremy L Thompson       return CeedError(basis->ceed, CEED_ERROR_INCOMPATIBLE,
4276e15d496SJeremy L Thompson                        "Tensor CEED_EVAL_DIV not supported"); break;
4286e15d496SJeremy L Thompson     case CEED_EVAL_CURL:
4296e15d496SJeremy L Thompson       return CeedError(basis->ceed, CEED_ERROR_INCOMPATIBLE,
4306e15d496SJeremy L Thompson                        "Tensor CEED_EVAL_CURL not supported"); break;
4316e15d496SJeremy L Thompson     // LCOV_EXCL_STOP
4326e15d496SJeremy L Thompson     case CEED_EVAL_WEIGHT: *flops = dim * CeedIntPow(Q_1d, dim); break;
4336e15d496SJeremy L Thompson     }
4346e15d496SJeremy L Thompson   } else {
4356e15d496SJeremy L Thompson     CeedInt dim, num_comp, num_nodes, num_qpts, Q_comp;
4366e15d496SJeremy L Thompson     ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr);
4376e15d496SJeremy L Thompson     ierr = CeedBasisGetNumComponents(basis, &num_comp); CeedChk(ierr);
4386e15d496SJeremy L Thompson     ierr = CeedBasisGetNumNodes(basis, &num_nodes); CeedChk(ierr);
4396e15d496SJeremy L Thompson     ierr = CeedBasisGetNumQuadraturePoints(basis, &num_qpts); CeedChk(ierr);
4406e15d496SJeremy L Thompson     ierr = CeedBasisGetNumQuadratureComponents(basis, &Q_comp); CeedChk(ierr);
4416e15d496SJeremy L Thompson     switch (eval_mode) {
4426e15d496SJeremy L Thompson     case CEED_EVAL_NONE:   *flops = 0; break;
4436e15d496SJeremy L Thompson     case CEED_EVAL_INTERP: *flops = num_nodes * num_qpts * num_comp; break;
4446e15d496SJeremy L Thompson     case CEED_EVAL_GRAD:   *flops = num_nodes * num_qpts * num_comp * dim; break;
4456e15d496SJeremy L Thompson     case CEED_EVAL_DIV:    *flops = num_nodes * num_qpts; break;
4466e15d496SJeremy L Thompson     case CEED_EVAL_CURL:   *flops = num_nodes * num_qpts * dim; break;
4476e15d496SJeremy L Thompson     case CEED_EVAL_WEIGHT: *flops = 0; break;
4486e15d496SJeremy L Thompson     }
4496e15d496SJeremy L Thompson   }
4506e15d496SJeremy L Thompson 
4516e15d496SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4526e15d496SJeremy L Thompson }
4536e15d496SJeremy L Thompson 
4546e15d496SJeremy L Thompson /**
4557a982d89SJeremy L. Thompson   @brief Get dimension for given CeedElemTopology
4567a982d89SJeremy L. Thompson 
4577a982d89SJeremy L. Thompson   @param topo      CeedElemTopology
4587a982d89SJeremy L. Thompson   @param[out] dim  Variable to store dimension of topology
4597a982d89SJeremy L. Thompson 
4607a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
4617a982d89SJeremy L. Thompson 
4627a982d89SJeremy L. Thompson   @ref Backend
4637a982d89SJeremy L. Thompson **/
4647a982d89SJeremy L. Thompson int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) {
4657a982d89SJeremy L. Thompson   *dim = (CeedInt) topo >> 16;
466e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4677a982d89SJeremy L. Thompson }
4687a982d89SJeremy L. Thompson 
4697a982d89SJeremy L. Thompson /**
4707a982d89SJeremy L. Thompson   @brief Get CeedTensorContract of a CeedBasis
4717a982d89SJeremy L. Thompson 
4727a982d89SJeremy L. Thompson   @param basis          CeedBasis
4737a982d89SJeremy L. Thompson   @param[out] contract  Variable to store CeedTensorContract
4747a982d89SJeremy L. Thompson 
4757a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
4767a982d89SJeremy L. Thompson 
4777a982d89SJeremy L. Thompson   @ref Backend
4787a982d89SJeremy L. Thompson **/
4797a982d89SJeremy L. Thompson int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) {
4807a982d89SJeremy L. Thompson   *contract = basis->contract;
481e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4827a982d89SJeremy L. Thompson }
4837a982d89SJeremy L. Thompson 
4847a982d89SJeremy L. Thompson /**
4857a982d89SJeremy L. Thompson   @brief Set CeedTensorContract of a CeedBasis
4867a982d89SJeremy L. Thompson 
4877a982d89SJeremy L. Thompson   @param[out] basis  CeedBasis
4887a982d89SJeremy L. Thompson   @param contract    CeedTensorContract to set
4897a982d89SJeremy L. Thompson 
4907a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
4917a982d89SJeremy L. Thompson 
4927a982d89SJeremy L. Thompson   @ref Backend
4937a982d89SJeremy L. Thompson **/
49434359f16Sjeremylt int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract contract) {
49534359f16Sjeremylt   int ierr;
49634359f16Sjeremylt   basis->contract = contract;
4979560d06aSjeremylt   ierr = CeedTensorContractReference(contract); CeedChk(ierr);
498e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4997a982d89SJeremy L. Thompson }
5007a982d89SJeremy L. Thompson 
5017a982d89SJeremy L. Thompson /**
5027a982d89SJeremy L. Thompson   @brief Return a reference implementation of matrix multiplication C = A B.
5037a982d89SJeremy L. Thompson            Note, this is a reference implementation for CPU CeedScalar pointers
5047a982d89SJeremy L. Thompson            that is not intended for high performance.
5057a982d89SJeremy L. Thompson 
5067a982d89SJeremy L. Thompson   @param ceed        A Ceed context for error handling
507d1d35e2fSjeremylt   @param[in] mat_A   Row-major matrix A
508d1d35e2fSjeremylt   @param[in] mat_B   Row-major matrix B
509d1d35e2fSjeremylt   @param[out] mat_C  Row-major output matrix C
5107a982d89SJeremy L. Thompson   @param m           Number of rows of C
5117a982d89SJeremy L. Thompson   @param n           Number of columns of C
5127a982d89SJeremy L. Thompson   @param kk          Number of columns of A/rows of B
5137a982d89SJeremy L. Thompson 
5147a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
5157a982d89SJeremy L. Thompson 
5167a982d89SJeremy L. Thompson   @ref Utility
5177a982d89SJeremy L. Thompson **/
518ed9e99e6SJeremy L Thompson int CeedMatrixMatrixMultiply(Ceed ceed, const CeedScalar *mat_A,
519ed9e99e6SJeremy L Thompson                              const CeedScalar *mat_B, CeedScalar *mat_C,
520ed9e99e6SJeremy L Thompson                              CeedInt m, CeedInt n, CeedInt kk) {
5217a982d89SJeremy L. Thompson   for (CeedInt i=0; i<m; i++)
5227a982d89SJeremy L. Thompson     for (CeedInt j=0; j<n; j++) {
5237a982d89SJeremy L. Thompson       CeedScalar sum = 0;
5247a982d89SJeremy L. Thompson       for (CeedInt k=0; k<kk; k++)
525d1d35e2fSjeremylt         sum += mat_A[k+i*kk]*mat_B[j+k*n];
526d1d35e2fSjeremylt       mat_C[j+i*n] = sum;
5277a982d89SJeremy L. Thompson     }
528e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5297a982d89SJeremy L. Thompson }
5307a982d89SJeremy L. Thompson 
5317a982d89SJeremy L. Thompson /// @}
5327a982d89SJeremy L. Thompson 
5337a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
5347a982d89SJeremy L. Thompson /// CeedBasis Public API
5357a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
5367a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser
537d7b241e6Sjeremylt /// @{
538d7b241e6Sjeremylt 
539b11c1e72Sjeremylt /**
54095bb1877Svaleriabarra   @brief Create a tensor-product basis for H^1 discretizations
541b11c1e72Sjeremylt 
542b11c1e72Sjeremylt   @param ceed        A Ceed object where the CeedBasis will be created
543b11c1e72Sjeremylt   @param dim         Topological dimension
544d1d35e2fSjeremylt   @param num_comp    Number of field components (1 for scalar fields)
545d1d35e2fSjeremylt   @param P_1d        Number of nodes in one dimension
546d1d35e2fSjeremylt   @param Q_1d        Number of quadrature points in one dimension
547d1d35e2fSjeremylt   @param interp_1d   Row-major (Q_1d * P_1d) matrix expressing the values of nodal
548b11c1e72Sjeremylt                        basis functions at quadrature points
549d1d35e2fSjeremylt   @param grad_1d     Row-major (Q_1d * P_1d) matrix expressing derivatives of nodal
550b11c1e72Sjeremylt                        basis functions at quadrature points
551d1d35e2fSjeremylt   @param q_ref_1d    Array of length Q_1d holding the locations of quadrature points
552b11c1e72Sjeremylt                        on the 1D reference element [-1, 1]
553d1d35e2fSjeremylt   @param q_weight_1d Array of length Q_1d holding the quadrature weights on the
554b11c1e72Sjeremylt                        reference element
555b11c1e72Sjeremylt   @param[out] basis  Address of the variable where the newly created
556b11c1e72Sjeremylt                        CeedBasis will be stored.
557b11c1e72Sjeremylt 
558b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
559dfdf5a53Sjeremylt 
5607a982d89SJeremy L. Thompson   @ref User
561b11c1e72Sjeremylt **/
562d1d35e2fSjeremylt int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt num_comp,
563d1d35e2fSjeremylt                             CeedInt P_1d, CeedInt Q_1d,
564d1d35e2fSjeremylt                             const CeedScalar *interp_1d,
565d1d35e2fSjeremylt                             const CeedScalar *grad_1d, const CeedScalar *q_ref_1d,
566d1d35e2fSjeremylt                             const CeedScalar *q_weight_1d, CeedBasis *basis) {
567d7b241e6Sjeremylt   int ierr;
568d7b241e6Sjeremylt 
5695fe0d4faSjeremylt   if (!ceed->BasisCreateTensorH1) {
5705fe0d4faSjeremylt     Ceed delegate;
571aefd8378Sjeremylt     ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr);
5725fe0d4faSjeremylt 
5735fe0d4faSjeremylt     if (!delegate)
574c042f62fSJeremy L Thompson       // LCOV_EXCL_START
575e15f9bd0SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_UNSUPPORTED,
576e15f9bd0SJeremy L Thompson                        "Backend does not support BasisCreateTensorH1");
577c042f62fSJeremy L Thompson     // LCOV_EXCL_STOP
5785fe0d4faSjeremylt 
579d1d35e2fSjeremylt     ierr = CeedBasisCreateTensorH1(delegate, dim, num_comp, P_1d,
580d1d35e2fSjeremylt                                    Q_1d, interp_1d, grad_1d, q_ref_1d,
581d1d35e2fSjeremylt                                    q_weight_1d, basis); CeedChk(ierr);
582e15f9bd0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
5835fe0d4faSjeremylt   }
584e15f9bd0SJeremy L Thompson 
585e15f9bd0SJeremy L Thompson   if (dim < 1)
586e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
587e15f9bd0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
588e15f9bd0SJeremy L Thompson                      "Basis dimension must be a positive value");
589e15f9bd0SJeremy L Thompson   // LCOV_EXCL_STOP
590227444bfSJeremy L Thompson 
591227444bfSJeremy L Thompson   if (num_comp < 1)
592227444bfSJeremy L Thompson     // LCOV_EXCL_START
593227444bfSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
594227444bfSJeremy L Thompson                      "Basis must have at least 1 component");
595227444bfSJeremy L Thompson   // LCOV_EXCL_STOP
596227444bfSJeremy L Thompson 
597227444bfSJeremy L Thompson   if (P_1d < 1)
598227444bfSJeremy L Thompson     // LCOV_EXCL_START
599227444bfSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
600227444bfSJeremy L Thompson                      "Basis must have at least 1 node");
601227444bfSJeremy L Thompson   // LCOV_EXCL_STOP
602227444bfSJeremy L Thompson 
603227444bfSJeremy L Thompson   if (Q_1d < 1)
604227444bfSJeremy L Thompson     // LCOV_EXCL_START
605227444bfSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
606227444bfSJeremy L Thompson                      "Basis must have at least 1 quadrature point");
607227444bfSJeremy L Thompson   // LCOV_EXCL_STOP
608227444bfSJeremy L Thompson 
60950c301a5SRezgar Shakeri   CeedElemTopology topo = dim == 1 ? CEED_TOPOLOGY_LINE
61050c301a5SRezgar Shakeri                           : dim == 2 ? CEED_TOPOLOGY_QUAD
61150c301a5SRezgar Shakeri                           : CEED_TOPOLOGY_HEX;
612e15f9bd0SJeremy L Thompson 
613d7b241e6Sjeremylt   ierr = CeedCalloc(1, basis); CeedChk(ierr);
614d7b241e6Sjeremylt   (*basis)->ceed = ceed;
6159560d06aSjeremylt   ierr = CeedReference(ceed); CeedChk(ierr);
616d1d35e2fSjeremylt   (*basis)->ref_count = 1;
617d1d35e2fSjeremylt   (*basis)->tensor_basis = 1;
618d7b241e6Sjeremylt   (*basis)->dim = dim;
619d99fa3c5SJeremy L Thompson   (*basis)->topo = topo;
620d1d35e2fSjeremylt   (*basis)->num_comp = num_comp;
621d1d35e2fSjeremylt   (*basis)->P_1d = P_1d;
622d1d35e2fSjeremylt   (*basis)->Q_1d = Q_1d;
623d1d35e2fSjeremylt   (*basis)->P = CeedIntPow(P_1d, dim);
624d1d35e2fSjeremylt   (*basis)->Q = CeedIntPow(Q_1d, dim);
62550c301a5SRezgar Shakeri   (*basis)->Q_comp = 1;
62650c301a5SRezgar Shakeri   (*basis)->basis_space = 1; // 1 for H^1 space
627ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(Q_1d, &(*basis)->q_ref_1d); CeedChk(ierr);
628ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(Q_1d, &(*basis)->q_weight_1d); CeedChk(ierr);
629ff3a0f91SJeremy L Thompson   if (q_ref_1d) memcpy((*basis)->q_ref_1d, q_ref_1d, Q_1d*sizeof(q_ref_1d[0]));
630ff3a0f91SJeremy L Thompson   if (q_weight_1d) memcpy((*basis)->q_weight_1d, q_weight_1d,
631ff3a0f91SJeremy L Thompson                             Q_1d*sizeof(q_weight_1d[0]));
632ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(Q_1d*P_1d, &(*basis)->interp_1d); CeedChk(ierr);
633ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(Q_1d*P_1d, &(*basis)->grad_1d); CeedChk(ierr);
634ff3a0f91SJeremy L Thompson   if (interp_1d) memcpy((*basis)->interp_1d, interp_1d,
635ff3a0f91SJeremy L Thompson                           Q_1d*P_1d*sizeof(interp_1d[0]));
636ff3a0f91SJeremy L Thompson   if (grad_1d) memcpy((*basis)->grad_1d, grad_1d, Q_1d*P_1d*sizeof(grad_1d[0]));
637d1d35e2fSjeremylt   ierr = ceed->BasisCreateTensorH1(dim, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d,
638d1d35e2fSjeremylt                                    q_weight_1d, *basis); CeedChk(ierr);
639e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
640d7b241e6Sjeremylt }
641d7b241e6Sjeremylt 
642b11c1e72Sjeremylt /**
64395bb1877Svaleriabarra   @brief Create a tensor-product Lagrange basis
644b11c1e72Sjeremylt 
645b11c1e72Sjeremylt   @param ceed        A Ceed object where the CeedBasis will be created
646b11c1e72Sjeremylt   @param dim         Topological dimension of element
647d1d35e2fSjeremylt   @param num_comp      Number of field components (1 for scalar fields)
648b11c1e72Sjeremylt   @param P           Number of Gauss-Lobatto nodes in one dimension.  The
649b11c1e72Sjeremylt                        polynomial degree of the resulting Q_k element is k=P-1.
650b11c1e72Sjeremylt   @param Q           Number of quadrature points in one dimension.
651d1d35e2fSjeremylt   @param quad_mode   Distribution of the Q quadrature points (affects order of
652b11c1e72Sjeremylt                        accuracy for the quadrature)
653b11c1e72Sjeremylt   @param[out] basis  Address of the variable where the newly created
654b11c1e72Sjeremylt                        CeedBasis will be stored.
655b11c1e72Sjeremylt 
656b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
657dfdf5a53Sjeremylt 
6587a982d89SJeremy L. Thompson   @ref User
659b11c1e72Sjeremylt **/
660d1d35e2fSjeremylt int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt num_comp,
661d1d35e2fSjeremylt                                     CeedInt P, CeedInt Q, CeedQuadMode quad_mode,
662692c2638Sjeremylt                                     CeedBasis *basis) {
663d7b241e6Sjeremylt   // Allocate
664e15f9bd0SJeremy L Thompson   int ierr, ierr2, i, j, k;
665d1d35e2fSjeremylt   CeedScalar c1, c2, c3, c4, dx, *nodes, *interp_1d, *grad_1d, *q_ref_1d,
666d1d35e2fSjeremylt              *q_weight_1d;
6674d537eeaSYohann 
6684d537eeaSYohann   if (dim < 1)
669c042f62fSJeremy L Thompson     // LCOV_EXCL_START
670e15f9bd0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
671e15f9bd0SJeremy L Thompson                      "Basis dimension must be a positive value");
672c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
6734d537eeaSYohann 
674227444bfSJeremy L Thompson   if (num_comp < 1)
675227444bfSJeremy L Thompson     // LCOV_EXCL_START
676227444bfSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
677227444bfSJeremy L Thompson                      "Basis must have at least 1 component");
678227444bfSJeremy L Thompson   // LCOV_EXCL_STOP
679227444bfSJeremy L Thompson 
680227444bfSJeremy L Thompson   if (P < 1)
681227444bfSJeremy L Thompson     // LCOV_EXCL_START
682227444bfSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
683227444bfSJeremy L Thompson                      "Basis must have at least 1 node");
684227444bfSJeremy L Thompson   // LCOV_EXCL_STOP
685227444bfSJeremy L Thompson 
686227444bfSJeremy L Thompson   if (Q < 1)
687227444bfSJeremy L Thompson     // LCOV_EXCL_START
688227444bfSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
689227444bfSJeremy L Thompson                      "Basis must have at least 1 quadrature point");
690227444bfSJeremy L Thompson   // LCOV_EXCL_STOP
691227444bfSJeremy L Thompson 
692e15f9bd0SJeremy L Thompson   // Get Nodes and Weights
693d1d35e2fSjeremylt   ierr = CeedCalloc(P*Q, &interp_1d); CeedChk(ierr);
694d1d35e2fSjeremylt   ierr = CeedCalloc(P*Q, &grad_1d); CeedChk(ierr);
695d7b241e6Sjeremylt   ierr = CeedCalloc(P, &nodes); CeedChk(ierr);
696d1d35e2fSjeremylt   ierr = CeedCalloc(Q, &q_ref_1d); CeedChk(ierr);
697d1d35e2fSjeremylt   ierr = CeedCalloc(Q, &q_weight_1d); CeedChk(ierr);
698e15f9bd0SJeremy L Thompson   ierr = CeedLobattoQuadrature(P, nodes, NULL);
699e15f9bd0SJeremy L Thompson   if (ierr) { goto cleanup; } CeedChk(ierr);
700d1d35e2fSjeremylt   switch (quad_mode) {
701d7b241e6Sjeremylt   case CEED_GAUSS:
702d1d35e2fSjeremylt     ierr = CeedGaussQuadrature(Q, q_ref_1d, q_weight_1d);
703d7b241e6Sjeremylt     break;
704d7b241e6Sjeremylt   case CEED_GAUSS_LOBATTO:
705d1d35e2fSjeremylt     ierr = CeedLobattoQuadrature(Q, q_ref_1d, q_weight_1d);
706d7b241e6Sjeremylt     break;
707d7b241e6Sjeremylt   }
708e15f9bd0SJeremy L Thompson   if (ierr) { goto cleanup; } CeedChk(ierr);
709e15f9bd0SJeremy L Thompson 
710d7b241e6Sjeremylt   // Build B, D matrix
711d7b241e6Sjeremylt   // Fornberg, 1998
712d7b241e6Sjeremylt   for (i = 0; i  < Q; i++) {
713d7b241e6Sjeremylt     c1 = 1.0;
714d1d35e2fSjeremylt     c3 = nodes[0] - q_ref_1d[i];
715d1d35e2fSjeremylt     interp_1d[i*P+0] = 1.0;
716d7b241e6Sjeremylt     for (j = 1; j < P; j++) {
717d7b241e6Sjeremylt       c2 = 1.0;
718d7b241e6Sjeremylt       c4 = c3;
719d1d35e2fSjeremylt       c3 = nodes[j] - q_ref_1d[i];
720d7b241e6Sjeremylt       for (k = 0; k < j; k++) {
721d7b241e6Sjeremylt         dx = nodes[j] - nodes[k];
722d7b241e6Sjeremylt         c2 *= dx;
723d7b241e6Sjeremylt         if (k == j - 1) {
724d1d35e2fSjeremylt           grad_1d[i*P + j] = c1*(interp_1d[i*P + k] - c4*grad_1d[i*P + k]) / c2;
725d1d35e2fSjeremylt           interp_1d[i*P + j] = - c1*c4*interp_1d[i*P + k] / c2;
726d7b241e6Sjeremylt         }
727d1d35e2fSjeremylt         grad_1d[i*P + k] = (c3*grad_1d[i*P + k] - interp_1d[i*P + k]) / dx;
728d1d35e2fSjeremylt         interp_1d[i*P + k] = c3*interp_1d[i*P + k] / dx;
729d7b241e6Sjeremylt       }
730d7b241e6Sjeremylt       c1 = c2;
731d7b241e6Sjeremylt     }
732d7b241e6Sjeremylt   }
7339ac7b42eSJeremy L Thompson   // Pass to CeedBasisCreateTensorH1
734d1d35e2fSjeremylt   ierr = CeedBasisCreateTensorH1(ceed, dim, num_comp, P, Q, interp_1d, grad_1d,
7359ac7b42eSJeremy L Thompson                                  q_ref_1d, q_weight_1d, basis); CeedChk(ierr);
736e15f9bd0SJeremy L Thompson cleanup:
737d1d35e2fSjeremylt   ierr2 = CeedFree(&interp_1d); CeedChk(ierr2);
738d1d35e2fSjeremylt   ierr2 = CeedFree(&grad_1d); CeedChk(ierr2);
739e15f9bd0SJeremy L Thompson   ierr2 = CeedFree(&nodes); CeedChk(ierr2);
740d1d35e2fSjeremylt   ierr2 = CeedFree(&q_ref_1d); CeedChk(ierr2);
741d1d35e2fSjeremylt   ierr2 = CeedFree(&q_weight_1d); CeedChk(ierr2);
742e15f9bd0SJeremy L Thompson   CeedChk(ierr);
743e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
744d7b241e6Sjeremylt }
745d7b241e6Sjeremylt 
746b11c1e72Sjeremylt /**
74795bb1877Svaleriabarra   @brief Create a non tensor-product basis for H^1 discretizations
748a8de75f0Sjeremylt 
749a8de75f0Sjeremylt   @param ceed        A Ceed object where the CeedBasis will be created
750a8de75f0Sjeremylt   @param topo        Topology of element, e.g. hypercube, simplex, ect
751d1d35e2fSjeremylt   @param num_comp    Number of field components (1 for scalar fields)
752d1d35e2fSjeremylt   @param num_nodes   Total number of nodes
753d1d35e2fSjeremylt   @param num_qpts    Total number of quadrature points
754d1d35e2fSjeremylt   @param interp      Row-major (num_qpts * num_nodes) matrix expressing the values of
7558795c945Sjeremylt                        nodal basis functions at quadrature points
756d1d35e2fSjeremylt   @param grad        Row-major (num_qpts * dim * num_nodes) matrix expressing
7578795c945Sjeremylt                        derivatives of nodal basis functions at quadrature points
758d1d35e2fSjeremylt   @param q_ref       Array of length num_qpts holding the locations of quadrature
75950c301a5SRezgar Shakeri                        points on the reference element
760d1d35e2fSjeremylt   @param q_weight    Array of length num_qpts holding the quadrature weights on the
761a8de75f0Sjeremylt                        reference element
762a8de75f0Sjeremylt   @param[out] basis  Address of the variable where the newly created
763a8de75f0Sjeremylt                        CeedBasis will be stored.
764a8de75f0Sjeremylt 
765a8de75f0Sjeremylt   @return An error code: 0 - success, otherwise - failure
766a8de75f0Sjeremylt 
7677a982d89SJeremy L. Thompson   @ref User
768a8de75f0Sjeremylt **/
769d1d35e2fSjeremylt int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt num_comp,
770d1d35e2fSjeremylt                       CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
771d1d35e2fSjeremylt                       const CeedScalar *grad, const CeedScalar *q_ref,
772d1d35e2fSjeremylt                       const CeedScalar *q_weight, CeedBasis *basis) {
773a8de75f0Sjeremylt   int ierr;
774d1d35e2fSjeremylt   CeedInt P = num_nodes, Q = num_qpts, dim = 0;
775a8de75f0Sjeremylt 
7765fe0d4faSjeremylt   if (!ceed->BasisCreateH1) {
7775fe0d4faSjeremylt     Ceed delegate;
778aefd8378Sjeremylt     ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr);
7795fe0d4faSjeremylt 
7805fe0d4faSjeremylt     if (!delegate)
781c042f62fSJeremy L Thompson       // LCOV_EXCL_START
782e15f9bd0SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_UNSUPPORTED,
783e15f9bd0SJeremy L Thompson                        "Backend does not support BasisCreateH1");
784c042f62fSJeremy L Thompson     // LCOV_EXCL_STOP
7855fe0d4faSjeremylt 
786d1d35e2fSjeremylt     ierr = CeedBasisCreateH1(delegate, topo, num_comp, num_nodes,
787d1d35e2fSjeremylt                              num_qpts, interp, grad, q_ref,
788d1d35e2fSjeremylt                              q_weight, basis); CeedChk(ierr);
789e15f9bd0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
7905fe0d4faSjeremylt   }
7915fe0d4faSjeremylt 
792227444bfSJeremy L Thompson   if (num_comp < 1)
793227444bfSJeremy L Thompson     // LCOV_EXCL_START
794227444bfSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
795227444bfSJeremy L Thompson                      "Basis must have at least 1 component");
796227444bfSJeremy L Thompson   // LCOV_EXCL_STOP
797227444bfSJeremy L Thompson 
798227444bfSJeremy L Thompson   if (num_nodes < 1)
799227444bfSJeremy L Thompson     // LCOV_EXCL_START
800227444bfSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
801227444bfSJeremy L Thompson                      "Basis must have at least 1 node");
802227444bfSJeremy L Thompson   // LCOV_EXCL_STOP
803227444bfSJeremy L Thompson 
804227444bfSJeremy L Thompson   if (num_qpts < 1)
805227444bfSJeremy L Thompson     // LCOV_EXCL_START
806227444bfSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
807227444bfSJeremy L Thompson                      "Basis must have at least 1 quadrature point");
808227444bfSJeremy L Thompson   // LCOV_EXCL_STOP
809227444bfSJeremy L Thompson 
810a8de75f0Sjeremylt   ierr = CeedCalloc(1, basis); CeedChk(ierr);
811a8de75f0Sjeremylt 
812a8de75f0Sjeremylt   ierr = CeedBasisGetTopologyDimension(topo, &dim); CeedChk(ierr);
813a8de75f0Sjeremylt 
814a8de75f0Sjeremylt   (*basis)->ceed = ceed;
8159560d06aSjeremylt   ierr = CeedReference(ceed); CeedChk(ierr);
816d1d35e2fSjeremylt   (*basis)->ref_count = 1;
817d1d35e2fSjeremylt   (*basis)->tensor_basis = 0;
818a8de75f0Sjeremylt   (*basis)->dim = dim;
819d99fa3c5SJeremy L Thompson   (*basis)->topo = topo;
820d1d35e2fSjeremylt   (*basis)->num_comp = num_comp;
821a8de75f0Sjeremylt   (*basis)->P = P;
822a8de75f0Sjeremylt   (*basis)->Q = Q;
82350c301a5SRezgar Shakeri   (*basis)->Q_comp = 1;
82450c301a5SRezgar Shakeri   (*basis)->basis_space = 1; // 1 for H^1 space
825ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(Q*dim, &(*basis)->q_ref_1d); CeedChk(ierr);
826ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(Q, &(*basis)->q_weight_1d); CeedChk(ierr);
827ff3a0f91SJeremy L Thompson   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q*dim*sizeof(q_ref[0]));
828ff3a0f91SJeremy L Thompson   if(q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q*sizeof(q_weight[0]));
829ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(Q*P, &(*basis)->interp); CeedChk(ierr);
830ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(dim*Q*P, &(*basis)->grad); CeedChk(ierr);
831ff3a0f91SJeremy L Thompson   if(interp) memcpy((*basis)->interp, interp, Q*P*sizeof(interp[0]));
832ff3a0f91SJeremy L Thompson   if(grad) memcpy((*basis)->grad, grad, dim*Q*P*sizeof(grad[0]));
833d1d35e2fSjeremylt   ierr = ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, q_ref,
834d1d35e2fSjeremylt                              q_weight, *basis); CeedChk(ierr);
835e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
836a8de75f0Sjeremylt }
837a8de75f0Sjeremylt 
838a8de75f0Sjeremylt /**
83950c301a5SRezgar Shakeri   @brief Create a non tensor-product basis for H(div) discretizations
84050c301a5SRezgar Shakeri 
84150c301a5SRezgar Shakeri   @param ceed        A Ceed object where the CeedBasis will be created
84250c301a5SRezgar Shakeri   @param topo        Topology of element (`CEED_TOPOLOGY_QUAD`, `CEED_TOPOLOGY_PRISM`, etc.),
84350c301a5SRezgar Shakeri                      dimension of which is used in some array sizes below
84450c301a5SRezgar Shakeri   @param num_comp    Number of components (usually 1 for vectors in H(div) bases)
84550c301a5SRezgar Shakeri   @param num_nodes   Total number of nodes (dofs per element)
84650c301a5SRezgar Shakeri   @param num_qpts    Total number of quadrature points
84750c301a5SRezgar Shakeri   @param interp      Row-major (dim*num_qpts * num_nodes) matrix expressing the values of
84850c301a5SRezgar Shakeri                        nodal basis functions at quadrature points
84950c301a5SRezgar Shakeri   @param div        Row-major (num_qpts * num_nodes) matrix expressing
85050c301a5SRezgar Shakeri                        divergence of nodal basis functions at quadrature points
85150c301a5SRezgar Shakeri   @param q_ref       Array of length num_qpts holding the locations of quadrature
85250c301a5SRezgar Shakeri                        points on the reference element
85350c301a5SRezgar Shakeri   @param q_weight    Array of length num_qpts holding the quadrature weights on the
85450c301a5SRezgar Shakeri                        reference element
85550c301a5SRezgar Shakeri   @param[out] basis  Address of the variable where the newly created
85650c301a5SRezgar Shakeri                        CeedBasis will be stored.
85750c301a5SRezgar Shakeri 
85850c301a5SRezgar Shakeri   @return An error code: 0 - success, otherwise - failure
85950c301a5SRezgar Shakeri 
86050c301a5SRezgar Shakeri   @ref User
86150c301a5SRezgar Shakeri **/
86250c301a5SRezgar Shakeri int CeedBasisCreateHdiv(Ceed ceed, CeedElemTopology topo, CeedInt num_comp,
86350c301a5SRezgar Shakeri                         CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
86450c301a5SRezgar Shakeri                         const CeedScalar *div, const CeedScalar *q_ref,
86550c301a5SRezgar Shakeri                         const CeedScalar *q_weight, CeedBasis *basis) {
86650c301a5SRezgar Shakeri   int ierr;
86750c301a5SRezgar Shakeri   CeedInt Q = num_qpts, P = num_nodes, dim = 0;
86850c301a5SRezgar Shakeri   ierr = CeedBasisGetTopologyDimension(topo, &dim); CeedChk(ierr);
86950c301a5SRezgar Shakeri   if (!ceed->BasisCreateHdiv) {
87050c301a5SRezgar Shakeri     Ceed delegate;
87150c301a5SRezgar Shakeri     ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr);
87250c301a5SRezgar Shakeri 
87350c301a5SRezgar Shakeri     if (!delegate)
87450c301a5SRezgar Shakeri       // LCOV_EXCL_START
87550c301a5SRezgar Shakeri       return CeedError(ceed, CEED_ERROR_UNSUPPORTED,
87650c301a5SRezgar Shakeri                        "Backend does not implement BasisCreateHdiv");
87750c301a5SRezgar Shakeri     // LCOV_EXCL_STOP
87850c301a5SRezgar Shakeri 
87950c301a5SRezgar Shakeri     ierr = CeedBasisCreateHdiv(delegate, topo, num_comp, num_nodes,
88050c301a5SRezgar Shakeri                                num_qpts, interp, div, q_ref,
88150c301a5SRezgar Shakeri                                q_weight, basis); CeedChk(ierr);
88250c301a5SRezgar Shakeri     return CEED_ERROR_SUCCESS;
88350c301a5SRezgar Shakeri   }
88450c301a5SRezgar Shakeri 
885227444bfSJeremy L Thompson   if (num_comp < 1)
886227444bfSJeremy L Thompson     // LCOV_EXCL_START
887227444bfSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
888227444bfSJeremy L Thompson                      "Basis must have at least 1 component");
889227444bfSJeremy L Thompson   // LCOV_EXCL_STOP
890227444bfSJeremy L Thompson 
891227444bfSJeremy L Thompson   if (num_nodes < 1)
892227444bfSJeremy L Thompson     // LCOV_EXCL_START
893227444bfSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
894227444bfSJeremy L Thompson                      "Basis must have at least 1 node");
895227444bfSJeremy L Thompson   // LCOV_EXCL_STOP
896227444bfSJeremy L Thompson 
897227444bfSJeremy L Thompson   if (num_qpts < 1)
898227444bfSJeremy L Thompson     // LCOV_EXCL_START
899227444bfSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
900227444bfSJeremy L Thompson                      "Basis must have at least 1 quadrature point");
901227444bfSJeremy L Thompson   // LCOV_EXCL_STOP
902227444bfSJeremy L Thompson 
90350c301a5SRezgar Shakeri   ierr = CeedCalloc(1, basis); CeedChk(ierr);
90450c301a5SRezgar Shakeri 
90550c301a5SRezgar Shakeri   (*basis)->ceed = ceed;
90650c301a5SRezgar Shakeri   ierr = CeedReference(ceed); CeedChk(ierr);
90750c301a5SRezgar Shakeri   (*basis)->ref_count = 1;
90850c301a5SRezgar Shakeri   (*basis)->tensor_basis = 0;
90950c301a5SRezgar Shakeri   (*basis)->dim = dim;
91050c301a5SRezgar Shakeri   (*basis)->topo = topo;
91150c301a5SRezgar Shakeri   (*basis)->num_comp = num_comp;
91250c301a5SRezgar Shakeri   (*basis)->P = P;
91350c301a5SRezgar Shakeri   (*basis)->Q = Q;
91450c301a5SRezgar Shakeri   (*basis)->Q_comp = dim;
91550c301a5SRezgar Shakeri   (*basis)->basis_space = 2; // 2 for H(div) space
91650c301a5SRezgar Shakeri   ierr = CeedMalloc(Q*dim, &(*basis)->q_ref_1d); CeedChk(ierr);
91750c301a5SRezgar Shakeri   ierr = CeedMalloc(Q, &(*basis)->q_weight_1d); CeedChk(ierr);
91850c301a5SRezgar Shakeri   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q*dim*sizeof(q_ref[0]));
91950c301a5SRezgar Shakeri   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q*sizeof(q_weight[0]));
92050c301a5SRezgar Shakeri   ierr = CeedMalloc(dim*Q*P, &(*basis)->interp); CeedChk(ierr);
92150c301a5SRezgar Shakeri   ierr = CeedMalloc(Q*P, &(*basis)->div); CeedChk(ierr);
92250c301a5SRezgar Shakeri   if (interp) memcpy((*basis)->interp, interp, dim*Q*P*sizeof(interp[0]));
92350c301a5SRezgar Shakeri   if (div) memcpy((*basis)->div, div, Q*P*sizeof(div[0]));
92450c301a5SRezgar Shakeri   ierr = ceed->BasisCreateHdiv(topo, dim, P, Q, interp, div, q_ref,
92550c301a5SRezgar Shakeri                                q_weight, *basis); CeedChk(ierr);
92650c301a5SRezgar Shakeri   return CEED_ERROR_SUCCESS;
92750c301a5SRezgar Shakeri }
92850c301a5SRezgar Shakeri 
92950c301a5SRezgar Shakeri /**
930f113e5dcSJeremy L Thompson   @brief Create a CeedBasis for projection from the nodes of `basis_from`
931f113e5dcSJeremy L Thompson            to the nodes of `basis_to`. Only `CEED_EVAL_INTERP` will be
932f113e5dcSJeremy L Thompson            valid for the new basis, `basis_project`. This projection is
933f113e5dcSJeremy L Thompson            given by `interp_project = interp_to^+ * interp_from`, where
934f113e5dcSJeremy L Thompson            the pesudoinverse `interp_to^+` is given by QR factorization.
935f113e5dcSJeremy L Thompson          Note: `basis_from` and `basis_to` must have compatible quadrature
936f113e5dcSJeremy L Thompson            spaces.
937446e7af4SJeremy L Thompson          Note: `basis_project` will have the same number of components as
938446e7af4SJeremy L Thompson            `basis_from`, regardless of the number of components that
939446e7af4SJeremy L Thompson            `basis_to` has. If `basis_from` has 3 components and `basis_to`
940446e7af4SJeremy L Thompson            has 5 components, then `basis_project` will have 3 components.
941f113e5dcSJeremy L Thompson 
942f113e5dcSJeremy L Thompson   @param[in] basis_from      CeedBasis to prolong from
943446e7af4SJeremy L Thompson   @param[in] basis_to        CeedBasis to prolong to
944f113e5dcSJeremy L Thompson   @param[out] basis_project  Address of the variable where the newly created
945f113e5dcSJeremy L Thompson                                CeedBasis will be stored.
946f113e5dcSJeremy L Thompson 
947f113e5dcSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
948f113e5dcSJeremy L Thompson 
949f113e5dcSJeremy L Thompson   @ref User
950f113e5dcSJeremy L Thompson **/
951446e7af4SJeremy L Thompson int CeedBasisCreateProjection(CeedBasis basis_from, CeedBasis basis_to,
952f113e5dcSJeremy L Thompson                               CeedBasis *basis_project) {
953f113e5dcSJeremy L Thompson   int ierr;
954f113e5dcSJeremy L Thompson   Ceed ceed;
955f113e5dcSJeremy L Thompson   ierr = CeedBasisGetCeed(basis_to, &ceed); CeedChk(ierr);
956f113e5dcSJeremy L Thompson 
957f113e5dcSJeremy L Thompson   // Create projectior matrix
958f113e5dcSJeremy L Thompson   CeedScalar *interp_project;
959446e7af4SJeremy L Thompson   ierr = CeedBasisCreateProjectionMatrix(basis_from, basis_to,
960f113e5dcSJeremy L Thompson                                          &interp_project); CeedChk(ierr);
961f113e5dcSJeremy L Thompson 
962f113e5dcSJeremy L Thompson   // Build basis
963f113e5dcSJeremy L Thompson   bool is_tensor;
964f113e5dcSJeremy L Thompson   CeedInt dim, num_comp;
965f113e5dcSJeremy L Thompson   CeedScalar *q_ref, *q_weight, *grad;
966f113e5dcSJeremy L Thompson   ierr = CeedBasisIsTensor(basis_to, &is_tensor); CeedChk(ierr);
967f113e5dcSJeremy L Thompson   ierr = CeedBasisGetDimension(basis_to, &dim); CeedChk(ierr);
968446e7af4SJeremy L Thompson   ierr = CeedBasisGetNumComponents(basis_from, &num_comp); CeedChk(ierr);
969f113e5dcSJeremy L Thompson   if (is_tensor) {
970f113e5dcSJeremy L Thompson     CeedInt P_1d_to, P_1d_from;
971f113e5dcSJeremy L Thompson     ierr = CeedBasisGetNumNodes1D(basis_from, &P_1d_from); CeedChk(ierr);
972f113e5dcSJeremy L Thompson     ierr = CeedBasisGetNumNodes1D(basis_to, &P_1d_to); CeedChk(ierr);
973f113e5dcSJeremy L Thompson     ierr = CeedCalloc(P_1d_to, &q_ref); CeedChk(ierr);
974f113e5dcSJeremy L Thompson     ierr = CeedCalloc(P_1d_to, &q_weight); CeedChk(ierr);
975f113e5dcSJeremy L Thompson     ierr = CeedCalloc(P_1d_to * P_1d_from * dim, &grad); CeedChk(ierr);
976f113e5dcSJeremy L Thompson     ierr = CeedBasisCreateTensorH1(ceed, dim, num_comp, P_1d_from, P_1d_to,
977f113e5dcSJeremy L Thompson                                    interp_project, grad, q_ref, q_weight, basis_project);
978f113e5dcSJeremy L Thompson     CeedChk(ierr);
979f113e5dcSJeremy L Thompson   } else {
980f113e5dcSJeremy L Thompson     CeedElemTopology topo;
981f113e5dcSJeremy L Thompson     ierr = CeedBasisGetTopology(basis_to, &topo); CeedChk(ierr);
982f113e5dcSJeremy L Thompson     CeedInt num_nodes_to, num_nodes_from;
983f113e5dcSJeremy L Thompson     ierr = CeedBasisGetNumNodes(basis_from, &num_nodes_from); CeedChk(ierr);
984f113e5dcSJeremy L Thompson     ierr = CeedBasisGetNumNodes(basis_to, &num_nodes_to); CeedChk(ierr);
985f113e5dcSJeremy L Thompson     ierr = CeedCalloc(num_nodes_to * dim, &q_ref); CeedChk(ierr);
986f113e5dcSJeremy L Thompson     ierr = CeedCalloc(num_nodes_to, &q_weight); CeedChk(ierr);
987f113e5dcSJeremy L Thompson     ierr = CeedCalloc(num_nodes_to * num_nodes_from * dim, &grad); CeedChk(ierr);
988f113e5dcSJeremy L Thompson     ierr = CeedBasisCreateH1(ceed, topo, num_comp, num_nodes_from, num_nodes_to,
989f113e5dcSJeremy L Thompson                              interp_project, grad, q_ref, q_weight, basis_project);
990f113e5dcSJeremy L Thompson     CeedChk(ierr);
991f113e5dcSJeremy L Thompson   }
992f113e5dcSJeremy L Thompson 
993f113e5dcSJeremy L Thompson   // Cleanup
994f113e5dcSJeremy L Thompson   ierr = CeedFree(&interp_project); CeedChk(ierr);
995f113e5dcSJeremy L Thompson   ierr = CeedFree(&q_ref); CeedChk(ierr);
996f113e5dcSJeremy L Thompson   ierr = CeedFree(&q_weight); CeedChk(ierr);
997f113e5dcSJeremy L Thompson   ierr = CeedFree(&grad); CeedChk(ierr);
998f113e5dcSJeremy L Thompson 
999f113e5dcSJeremy L Thompson   return CEED_ERROR_SUCCESS;
1000f113e5dcSJeremy L Thompson }
1001f113e5dcSJeremy L Thompson 
1002f113e5dcSJeremy L Thompson /**
10039560d06aSjeremylt   @brief Copy the pointer to a CeedBasis. Both pointers should
10049560d06aSjeremylt            be destroyed with `CeedBasisDestroy()`;
10059560d06aSjeremylt            Note: If `*basis_copy` is non-NULL, then it is assumed that
10069560d06aSjeremylt            `*basis_copy` is a pointer to a CeedBasis. This CeedBasis
10079560d06aSjeremylt            will be destroyed if `*basis_copy` is the only
10089560d06aSjeremylt            reference to this CeedBasis.
10099560d06aSjeremylt 
10109560d06aSjeremylt   @param basis            CeedBasis to copy reference to
10119560d06aSjeremylt   @param[out] basis_copy  Variable to store copied reference
10129560d06aSjeremylt 
10139560d06aSjeremylt   @return An error code: 0 - success, otherwise - failure
10149560d06aSjeremylt 
10159560d06aSjeremylt   @ref User
10169560d06aSjeremylt **/
10179560d06aSjeremylt int CeedBasisReferenceCopy(CeedBasis basis, CeedBasis *basis_copy) {
10189560d06aSjeremylt   int ierr;
10199560d06aSjeremylt 
10209560d06aSjeremylt   ierr = CeedBasisReference(basis); CeedChk(ierr);
10219560d06aSjeremylt   ierr = CeedBasisDestroy(basis_copy); CeedChk(ierr);
10229560d06aSjeremylt   *basis_copy = basis;
10239560d06aSjeremylt   return CEED_ERROR_SUCCESS;
10249560d06aSjeremylt }
10259560d06aSjeremylt 
10269560d06aSjeremylt /**
10277a982d89SJeremy L. Thompson   @brief View a CeedBasis
10287a982d89SJeremy L. Thompson 
10297a982d89SJeremy L. Thompson   @param basis   CeedBasis to view
10307a982d89SJeremy L. Thompson   @param stream  Stream to view to, e.g., stdout
10317a982d89SJeremy L. Thompson 
10327a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
10337a982d89SJeremy L. Thompson 
10347a982d89SJeremy L. Thompson   @ref User
10357a982d89SJeremy L. Thompson **/
10367a982d89SJeremy L. Thompson int CeedBasisView(CeedBasis basis, FILE *stream) {
10377a982d89SJeremy L. Thompson   int ierr;
103850c301a5SRezgar Shakeri   CeedFESpace FE_space = basis->basis_space;
103950c301a5SRezgar Shakeri   CeedElemTopology topo = basis->topo;
104050c301a5SRezgar Shakeri   // Print FE space and element topology of the basis
1041d1d35e2fSjeremylt   if (basis->tensor_basis) {
1042990fdeb6SJeremy L Thompson     fprintf(stream, "CeedBasis (%s on a %s element): dim=%" CeedInt_FMT " P=%"
1043990fdeb6SJeremy L Thompson             CeedInt_FMT " Q=%" CeedInt_FMT "\n",
104450c301a5SRezgar Shakeri             CeedFESpaces[FE_space], CeedElemTopologies[topo],
104550c301a5SRezgar Shakeri             basis->dim, basis->P_1d, basis->Q_1d);
104650c301a5SRezgar Shakeri   } else {
1047990fdeb6SJeremy L Thompson     fprintf(stream, "CeedBasis (%s on a %s element): dim=%" CeedInt_FMT " P=%"
1048990fdeb6SJeremy L Thompson             CeedInt_FMT " Q=%" CeedInt_FMT "\n",
104950c301a5SRezgar Shakeri             CeedFESpaces[FE_space], CeedElemTopologies[topo],
105050c301a5SRezgar Shakeri             basis->dim, basis->P, basis->Q);
105150c301a5SRezgar Shakeri   }
105250c301a5SRezgar Shakeri   // Print quadrature data, interpolation/gradient/divergene/curl of the basis
105350c301a5SRezgar Shakeri   if (basis->tensor_basis) { // tensor basis
1054d1d35e2fSjeremylt     ierr = CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_ref_1d,
10557a982d89SJeremy L. Thompson                           stream); CeedChk(ierr);
1056d1d35e2fSjeremylt     ierr = CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q_1d,
1057d1d35e2fSjeremylt                           basis->q_weight_1d, stream); CeedChk(ierr);
1058d1d35e2fSjeremylt     ierr = CeedScalarView("interp1d", "\t% 12.8f", basis->Q_1d, basis->P_1d,
1059d1d35e2fSjeremylt                           basis->interp_1d, stream); CeedChk(ierr);
1060d1d35e2fSjeremylt     ierr = CeedScalarView("grad1d", "\t% 12.8f", basis->Q_1d, basis->P_1d,
1061d1d35e2fSjeremylt                           basis->grad_1d, stream); CeedChk(ierr);
106250c301a5SRezgar Shakeri   } else { // non-tensor basis
10637a982d89SJeremy L. Thompson     ierr = CeedScalarView("qref", "\t% 12.8f", 1, basis->Q*basis->dim,
1064d1d35e2fSjeremylt                           basis->q_ref_1d,
10657a982d89SJeremy L. Thompson                           stream); CeedChk(ierr);
1066d1d35e2fSjeremylt     ierr = CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->q_weight_1d,
10677a982d89SJeremy L. Thompson                           stream); CeedChk(ierr);
106850c301a5SRezgar Shakeri     ierr = CeedScalarView("interp", "\t% 12.8f", basis->Q_comp*basis->Q, basis->P,
10697a982d89SJeremy L. Thompson                           basis->interp, stream); CeedChk(ierr);
107050c301a5SRezgar Shakeri     if (basis->grad) {
10717a982d89SJeremy L. Thompson       ierr = CeedScalarView("grad", "\t% 12.8f", basis->dim*basis->Q, basis->P,
10727a982d89SJeremy L. Thompson                             basis->grad, stream); CeedChk(ierr);
10737a982d89SJeremy L. Thompson     }
107450c301a5SRezgar Shakeri     if (basis->div) {
107550c301a5SRezgar Shakeri       ierr = CeedScalarView("div", "\t% 12.8f", basis->Q, basis->P,
107650c301a5SRezgar Shakeri                             basis->div, stream); CeedChk(ierr);
107750c301a5SRezgar Shakeri     }
107850c301a5SRezgar Shakeri   }
1079e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
10807a982d89SJeremy L. Thompson }
10817a982d89SJeremy L. Thompson 
10827a982d89SJeremy L. Thompson /**
10837a982d89SJeremy L. Thompson   @brief Apply basis evaluation from nodes to quadrature points or vice versa
10847a982d89SJeremy L. Thompson 
10857a982d89SJeremy L. Thompson   @param basis     CeedBasis to evaluate
1086d1d35e2fSjeremylt   @param num_elem  The number of elements to apply the basis evaluation to;
10877a982d89SJeremy L. Thompson                      the backend will specify the ordering in
10884cc79fe7SJed Brown                      CeedElemRestrictionCreateBlocked()
1089d1d35e2fSjeremylt   @param t_mode    \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature
10907a982d89SJeremy L. Thompson                      points, \ref CEED_TRANSPOSE to apply the transpose, mapping
10917a982d89SJeremy L. Thompson                      from quadrature points to nodes
1092d1d35e2fSjeremylt   @param eval_mode \ref CEED_EVAL_NONE to use values directly,
10937a982d89SJeremy L. Thompson                      \ref CEED_EVAL_INTERP to use interpolated values,
10947a982d89SJeremy L. Thompson                      \ref CEED_EVAL_GRAD to use gradients,
10957a982d89SJeremy L. Thompson                      \ref CEED_EVAL_WEIGHT to use quadrature weights.
10967a982d89SJeremy L. Thompson   @param[in] u     Input CeedVector
10977a982d89SJeremy L. Thompson   @param[out] v    Output CeedVector
10987a982d89SJeremy L. Thompson 
10997a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
11007a982d89SJeremy L. Thompson 
11017a982d89SJeremy L. Thompson   @ref User
11027a982d89SJeremy L. Thompson **/
1103d1d35e2fSjeremylt int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode,
1104d1d35e2fSjeremylt                    CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
11057a982d89SJeremy L. Thompson   int ierr;
11061f9221feSJeremy L Thompson   CeedSize u_length = 0, v_length;
11071f9221feSJeremy L Thompson   CeedInt dim, num_comp, num_nodes, num_qpts;
1108e15f9bd0SJeremy L Thompson   ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr);
1109d1d35e2fSjeremylt   ierr = CeedBasisGetNumComponents(basis, &num_comp); CeedChk(ierr);
1110d1d35e2fSjeremylt   ierr = CeedBasisGetNumNodes(basis, &num_nodes); CeedChk(ierr);
1111d1d35e2fSjeremylt   ierr = CeedBasisGetNumQuadraturePoints(basis, &num_qpts); CeedChk(ierr);
1112d1d35e2fSjeremylt   ierr = CeedVectorGetLength(v, &v_length); CeedChk(ierr);
11137a982d89SJeremy L. Thompson   if (u) {
1114d1d35e2fSjeremylt     ierr = CeedVectorGetLength(u, &u_length); CeedChk(ierr);
11157a982d89SJeremy L. Thompson   }
11167a982d89SJeremy L. Thompson 
1117e15f9bd0SJeremy L Thompson   if (!basis->Apply)
1118e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
1119e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_UNSUPPORTED,
1120e15f9bd0SJeremy L Thompson                      "Backend does not support BasisApply");
1121e15f9bd0SJeremy L Thompson   // LCOV_EXCL_STOP
1122e15f9bd0SJeremy L Thompson 
1123e15f9bd0SJeremy L Thompson   // Check compatibility of topological and geometrical dimensions
1124d1d35e2fSjeremylt   if ((t_mode == CEED_TRANSPOSE && (v_length%num_nodes != 0 ||
1125d1d35e2fSjeremylt                                     u_length%num_qpts != 0)) ||
1126d1d35e2fSjeremylt       (t_mode == CEED_NOTRANSPOSE && (u_length%num_nodes != 0 ||
1127d1d35e2fSjeremylt                                       v_length%num_qpts != 0)))
11288229195eSjeremylt     // LCOV_EXCL_START
1129e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_DIMENSION,
1130e15f9bd0SJeremy L Thompson                      "Length of input/output vectors "
11317a982d89SJeremy L. Thompson                      "incompatible with basis dimensions");
11328229195eSjeremylt   // LCOV_EXCL_STOP
11337a982d89SJeremy L. Thompson 
1134e15f9bd0SJeremy L Thompson   // Check vector lengths to prevent out of bounds issues
1135d1d35e2fSjeremylt   bool bad_dims = false;
1136d1d35e2fSjeremylt   switch (eval_mode) {
1137e15f9bd0SJeremy L Thompson   case CEED_EVAL_NONE:
1138d1d35e2fSjeremylt   case CEED_EVAL_INTERP: bad_dims =
1139d1d35e2fSjeremylt       ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts ||
1140d1d35e2fSjeremylt                                      v_length < num_elem*num_comp*num_nodes)) ||
1141d1d35e2fSjeremylt        (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp ||
1142d1d35e2fSjeremylt                                        u_length < num_elem*num_comp*num_nodes)));
1143e15f9bd0SJeremy L Thompson     break;
1144d1d35e2fSjeremylt   case CEED_EVAL_GRAD: bad_dims =
1145d1d35e2fSjeremylt       ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts*dim ||
1146d1d35e2fSjeremylt                                      v_length < num_elem*num_comp*num_nodes)) ||
1147d1d35e2fSjeremylt        (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp*dim ||
1148d1d35e2fSjeremylt                                        u_length < num_elem*num_comp*num_nodes)));
1149e15f9bd0SJeremy L Thompson     break;
1150e15f9bd0SJeremy L Thompson   case CEED_EVAL_WEIGHT:
1151d1d35e2fSjeremylt     bad_dims = v_length < num_elem*num_qpts;
1152e15f9bd0SJeremy L Thompson     break;
1153e15f9bd0SJeremy L Thompson   // LCOV_EXCL_START
1154d1d35e2fSjeremylt   case CEED_EVAL_DIV: bad_dims =
1155d1d35e2fSjeremylt       ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts ||
1156d1d35e2fSjeremylt                                      v_length < num_elem*num_comp*num_nodes)) ||
1157d1d35e2fSjeremylt        (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp ||
1158d1d35e2fSjeremylt                                        u_length < num_elem*num_comp*num_nodes)));
1159e15f9bd0SJeremy L Thompson     break;
1160d1d35e2fSjeremylt   case CEED_EVAL_CURL: bad_dims =
1161d1d35e2fSjeremylt       ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts ||
1162d1d35e2fSjeremylt                                      v_length < num_elem*num_comp*num_nodes)) ||
1163d1d35e2fSjeremylt        (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp ||
1164d1d35e2fSjeremylt                                        u_length < num_elem*num_comp*num_nodes)));
1165e15f9bd0SJeremy L Thompson     break;
1166e15f9bd0SJeremy L Thompson     // LCOV_EXCL_STOP
1167e15f9bd0SJeremy L Thompson   }
1168d1d35e2fSjeremylt   if (bad_dims)
1169e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
1170e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_DIMENSION,
1171d1d35e2fSjeremylt                      "Input/output vectors too short for basis and evaluation mode");
1172e15f9bd0SJeremy L Thompson   // LCOV_EXCL_STOP
1173e15f9bd0SJeremy L Thompson 
1174d1d35e2fSjeremylt   ierr = basis->Apply(basis, num_elem, t_mode, eval_mode, u, v); CeedChk(ierr);
1175e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
11767a982d89SJeremy L. Thompson }
11777a982d89SJeremy L. Thompson 
11787a982d89SJeremy L. Thompson /**
1179b7c9bbdaSJeremy L Thompson   @brief Get Ceed associated with a CeedBasis
1180b7c9bbdaSJeremy L Thompson 
1181b7c9bbdaSJeremy L Thompson   @param basis      CeedBasis
1182b7c9bbdaSJeremy L Thompson   @param[out] ceed  Variable to store Ceed
1183b7c9bbdaSJeremy L Thompson 
1184b7c9bbdaSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1185b7c9bbdaSJeremy L Thompson 
1186b7c9bbdaSJeremy L Thompson   @ref Advanced
1187b7c9bbdaSJeremy L Thompson **/
1188b7c9bbdaSJeremy L Thompson int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) {
1189b7c9bbdaSJeremy L Thompson   *ceed = basis->ceed;
1190b7c9bbdaSJeremy L Thompson   return CEED_ERROR_SUCCESS;
1191b7c9bbdaSJeremy L Thompson }
1192b7c9bbdaSJeremy L Thompson 
1193b7c9bbdaSJeremy L Thompson /**
11949d007619Sjeremylt   @brief Get dimension for given CeedBasis
11959d007619Sjeremylt 
11969d007619Sjeremylt   @param basis     CeedBasis
11979d007619Sjeremylt   @param[out] dim  Variable to store dimension of basis
11989d007619Sjeremylt 
11999d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
12009d007619Sjeremylt 
1201b7c9bbdaSJeremy L Thompson   @ref Advanced
12029d007619Sjeremylt **/
12039d007619Sjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) {
12049d007619Sjeremylt   *dim = basis->dim;
1205e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
12069d007619Sjeremylt }
12079d007619Sjeremylt 
12089d007619Sjeremylt /**
1209d99fa3c5SJeremy L Thompson   @brief Get topology for given CeedBasis
1210d99fa3c5SJeremy L Thompson 
1211d99fa3c5SJeremy L Thompson   @param basis      CeedBasis
1212d99fa3c5SJeremy L Thompson   @param[out] topo  Variable to store topology of basis
1213d99fa3c5SJeremy L Thompson 
1214d99fa3c5SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1215d99fa3c5SJeremy L Thompson 
1216b7c9bbdaSJeremy L Thompson   @ref Advanced
1217d99fa3c5SJeremy L Thompson **/
1218d99fa3c5SJeremy L Thompson int CeedBasisGetTopology(CeedBasis basis, CeedElemTopology *topo) {
1219d99fa3c5SJeremy L Thompson   *topo = basis->topo;
1220e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1221d99fa3c5SJeremy L Thompson }
1222d99fa3c5SJeremy L Thompson 
1223d99fa3c5SJeremy L Thompson /**
122450c301a5SRezgar Shakeri   @brief Get number of Q-vector components for given CeedBasis
122550c301a5SRezgar Shakeri 
122650c301a5SRezgar Shakeri   @param basis          CeedBasis
122750c301a5SRezgar Shakeri   @param[out] Q_comp  Variable to store number of Q-vector components of basis
122850c301a5SRezgar Shakeri 
122950c301a5SRezgar Shakeri   @return An error code: 0 - success, otherwise - failure
123050c301a5SRezgar Shakeri 
123150c301a5SRezgar Shakeri   @ref Advanced
123250c301a5SRezgar Shakeri **/
123350c301a5SRezgar Shakeri int CeedBasisGetNumQuadratureComponents(CeedBasis basis, CeedInt *Q_comp) {
123450c301a5SRezgar Shakeri   *Q_comp = basis->Q_comp;
123550c301a5SRezgar Shakeri   return CEED_ERROR_SUCCESS;
123650c301a5SRezgar Shakeri }
123750c301a5SRezgar Shakeri 
123850c301a5SRezgar Shakeri /**
12399d007619Sjeremylt   @brief Get number of components for given CeedBasis
12409d007619Sjeremylt 
12419d007619Sjeremylt   @param basis          CeedBasis
1242d1d35e2fSjeremylt   @param[out] num_comp  Variable to store number of components of basis
12439d007619Sjeremylt 
12449d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
12459d007619Sjeremylt 
1246b7c9bbdaSJeremy L Thompson   @ref Advanced
12479d007619Sjeremylt **/
1248d1d35e2fSjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *num_comp) {
1249d1d35e2fSjeremylt   *num_comp = basis->num_comp;
1250e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
12519d007619Sjeremylt }
12529d007619Sjeremylt 
12539d007619Sjeremylt /**
12549d007619Sjeremylt   @brief Get total number of nodes (in dim dimensions) of a CeedBasis
12559d007619Sjeremylt 
12569d007619Sjeremylt   @param basis   CeedBasis
12579d007619Sjeremylt   @param[out] P  Variable to store number of nodes
12589d007619Sjeremylt 
12599d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
12609d007619Sjeremylt 
12619d007619Sjeremylt   @ref Utility
12629d007619Sjeremylt **/
12639d007619Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) {
12649d007619Sjeremylt   *P = basis->P;
1265e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
12669d007619Sjeremylt }
12679d007619Sjeremylt 
12689d007619Sjeremylt /**
12699d007619Sjeremylt   @brief Get total number of nodes (in 1 dimension) of a CeedBasis
12709d007619Sjeremylt 
12719d007619Sjeremylt   @param basis     CeedBasis
1272d1d35e2fSjeremylt   @param[out] P_1d  Variable to store number of nodes
12739d007619Sjeremylt 
12749d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
12759d007619Sjeremylt 
1276b7c9bbdaSJeremy L Thompson   @ref Advanced
12779d007619Sjeremylt **/
1278d1d35e2fSjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P_1d) {
1279d1d35e2fSjeremylt   if (!basis->tensor_basis)
12809d007619Sjeremylt     // LCOV_EXCL_START
1281e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR,
1282d1d35e2fSjeremylt                      "Cannot supply P_1d for non-tensor basis");
12839d007619Sjeremylt   // LCOV_EXCL_STOP
12849d007619Sjeremylt 
1285d1d35e2fSjeremylt   *P_1d = basis->P_1d;
1286e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
12879d007619Sjeremylt }
12889d007619Sjeremylt 
12899d007619Sjeremylt /**
12909d007619Sjeremylt   @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis
12919d007619Sjeremylt 
12929d007619Sjeremylt   @param basis   CeedBasis
12939d007619Sjeremylt   @param[out] Q  Variable to store number of quadrature points
12949d007619Sjeremylt 
12959d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
12969d007619Sjeremylt 
12979d007619Sjeremylt   @ref Utility
12989d007619Sjeremylt **/
12999d007619Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) {
13009d007619Sjeremylt   *Q = basis->Q;
1301e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
13029d007619Sjeremylt }
13039d007619Sjeremylt 
13049d007619Sjeremylt /**
13059d007619Sjeremylt   @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis
13069d007619Sjeremylt 
13079d007619Sjeremylt   @param basis      CeedBasis
1308d1d35e2fSjeremylt   @param[out] Q_1d  Variable to store number of quadrature points
13099d007619Sjeremylt 
13109d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
13119d007619Sjeremylt 
1312b7c9bbdaSJeremy L Thompson   @ref Advanced
13139d007619Sjeremylt **/
1314d1d35e2fSjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q_1d) {
1315d1d35e2fSjeremylt   if (!basis->tensor_basis)
13169d007619Sjeremylt     // LCOV_EXCL_START
1317e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR,
1318d1d35e2fSjeremylt                      "Cannot supply Q_1d for non-tensor basis");
13199d007619Sjeremylt   // LCOV_EXCL_STOP
13209d007619Sjeremylt 
1321d1d35e2fSjeremylt   *Q_1d = basis->Q_1d;
1322e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
13239d007619Sjeremylt }
13249d007619Sjeremylt 
13259d007619Sjeremylt /**
13269d007619Sjeremylt   @brief Get reference coordinates of quadrature points (in dim dimensions)
13279d007619Sjeremylt          of a CeedBasis
13289d007619Sjeremylt 
13299d007619Sjeremylt   @param basis       CeedBasis
1330d1d35e2fSjeremylt   @param[out] q_ref  Variable to store reference coordinates of quadrature points
13319d007619Sjeremylt 
13329d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
13339d007619Sjeremylt 
1334b7c9bbdaSJeremy L Thompson   @ref Advanced
13359d007619Sjeremylt **/
1336d1d35e2fSjeremylt int CeedBasisGetQRef(CeedBasis basis, const CeedScalar **q_ref) {
1337d1d35e2fSjeremylt   *q_ref = basis->q_ref_1d;
1338e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
13399d007619Sjeremylt }
13409d007619Sjeremylt 
13419d007619Sjeremylt /**
13429d007619Sjeremylt   @brief Get quadrature weights of quadrature points (in dim dimensions)
13439d007619Sjeremylt          of a CeedBasis
13449d007619Sjeremylt 
13459d007619Sjeremylt   @param basis          CeedBasis
1346d1d35e2fSjeremylt   @param[out] q_weight  Variable to store quadrature weights
13479d007619Sjeremylt 
13489d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
13499d007619Sjeremylt 
1350b7c9bbdaSJeremy L Thompson   @ref Advanced
13519d007619Sjeremylt **/
1352d1d35e2fSjeremylt int CeedBasisGetQWeights(CeedBasis basis, const CeedScalar **q_weight) {
1353d1d35e2fSjeremylt   *q_weight = basis->q_weight_1d;
1354e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
13559d007619Sjeremylt }
13569d007619Sjeremylt 
13579d007619Sjeremylt /**
13589d007619Sjeremylt   @brief Get interpolation matrix of a CeedBasis
13599d007619Sjeremylt 
13609d007619Sjeremylt   @param basis        CeedBasis
13619d007619Sjeremylt   @param[out] interp  Variable to store interpolation matrix
13629d007619Sjeremylt 
13639d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
13649d007619Sjeremylt 
1365b7c9bbdaSJeremy L Thompson   @ref Advanced
13669d007619Sjeremylt **/
13676c58de82SJeremy L Thompson int CeedBasisGetInterp(CeedBasis basis, const CeedScalar **interp) {
1368d1d35e2fSjeremylt   if (!basis->interp && basis->tensor_basis) {
13699d007619Sjeremylt     // Allocate
13709d007619Sjeremylt     int ierr;
13719d007619Sjeremylt     ierr = CeedMalloc(basis->Q*basis->P, &basis->interp); CeedChk(ierr);
13729d007619Sjeremylt 
13739d007619Sjeremylt     // Initialize
13749d007619Sjeremylt     for (CeedInt i=0; i<basis->Q*basis->P; i++)
13759d007619Sjeremylt       basis->interp[i] = 1.0;
13769d007619Sjeremylt 
13779d007619Sjeremylt     // Calculate
13789d007619Sjeremylt     for (CeedInt d=0; d<basis->dim; d++)
13799d007619Sjeremylt       for (CeedInt qpt=0; qpt<basis->Q; qpt++)
13809d007619Sjeremylt         for (CeedInt node=0; node<basis->P; node++) {
1381d1d35e2fSjeremylt           CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d;
1382d1d35e2fSjeremylt           CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d;
1383d1d35e2fSjeremylt           basis->interp[qpt*(basis->P)+node] *= basis->interp_1d[q*basis->P_1d+p];
13849d007619Sjeremylt         }
13859d007619Sjeremylt   }
13869d007619Sjeremylt   *interp = basis->interp;
1387e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
13889d007619Sjeremylt }
13899d007619Sjeremylt 
13909d007619Sjeremylt /**
13919d007619Sjeremylt   @brief Get 1D interpolation matrix of a tensor product CeedBasis
13929d007619Sjeremylt 
13939d007619Sjeremylt   @param basis           CeedBasis
1394d1d35e2fSjeremylt   @param[out] interp_1d  Variable to store interpolation matrix
13959d007619Sjeremylt 
13969d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
13979d007619Sjeremylt 
13989d007619Sjeremylt   @ref Backend
13999d007619Sjeremylt **/
1400d1d35e2fSjeremylt int CeedBasisGetInterp1D(CeedBasis basis, const CeedScalar **interp_1d) {
1401d1d35e2fSjeremylt   if (!basis->tensor_basis)
14029d007619Sjeremylt     // LCOV_EXCL_START
1403e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR,
1404e15f9bd0SJeremy L Thompson                      "CeedBasis is not a tensor product basis.");
14059d007619Sjeremylt   // LCOV_EXCL_STOP
14069d007619Sjeremylt 
1407d1d35e2fSjeremylt   *interp_1d = basis->interp_1d;
1408e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
14099d007619Sjeremylt }
14109d007619Sjeremylt 
14119d007619Sjeremylt /**
14129d007619Sjeremylt   @brief Get gradient matrix of a CeedBasis
14139d007619Sjeremylt 
14149d007619Sjeremylt   @param basis      CeedBasis
14159d007619Sjeremylt   @param[out] grad  Variable to store gradient matrix
14169d007619Sjeremylt 
14179d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
14189d007619Sjeremylt 
1419b7c9bbdaSJeremy L Thompson   @ref Advanced
14209d007619Sjeremylt **/
14216c58de82SJeremy L Thompson int CeedBasisGetGrad(CeedBasis basis, const CeedScalar **grad) {
1422d1d35e2fSjeremylt   if (!basis->grad && basis->tensor_basis) {
14239d007619Sjeremylt     // Allocate
14249d007619Sjeremylt     int ierr;
14259d007619Sjeremylt     ierr = CeedMalloc(basis->dim*basis->Q*basis->P, &basis->grad);
14269d007619Sjeremylt     CeedChk(ierr);
14279d007619Sjeremylt 
14289d007619Sjeremylt     // Initialize
14299d007619Sjeremylt     for (CeedInt i=0; i<basis->dim*basis->Q*basis->P; i++)
14309d007619Sjeremylt       basis->grad[i] = 1.0;
14319d007619Sjeremylt 
14329d007619Sjeremylt     // Calculate
14339d007619Sjeremylt     for (CeedInt d=0; d<basis->dim; d++)
14349d007619Sjeremylt       for (CeedInt i=0; i<basis->dim; i++)
14359d007619Sjeremylt         for (CeedInt qpt=0; qpt<basis->Q; qpt++)
14369d007619Sjeremylt           for (CeedInt node=0; node<basis->P; node++) {
1437d1d35e2fSjeremylt             CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d;
1438d1d35e2fSjeremylt             CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d;
14399d007619Sjeremylt             if (i == d)
14409d007619Sjeremylt               basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *=
1441d1d35e2fSjeremylt                 basis->grad_1d[q*basis->P_1d+p];
14429d007619Sjeremylt             else
14439d007619Sjeremylt               basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *=
1444d1d35e2fSjeremylt                 basis->interp_1d[q*basis->P_1d+p];
14459d007619Sjeremylt           }
14469d007619Sjeremylt   }
14479d007619Sjeremylt   *grad = basis->grad;
1448e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
14499d007619Sjeremylt }
14509d007619Sjeremylt 
14519d007619Sjeremylt /**
14529d007619Sjeremylt   @brief Get 1D gradient matrix of a tensor product CeedBasis
14539d007619Sjeremylt 
14549d007619Sjeremylt   @param basis         CeedBasis
1455d1d35e2fSjeremylt   @param[out] grad_1d  Variable to store gradient matrix
14569d007619Sjeremylt 
14579d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
14589d007619Sjeremylt 
1459b7c9bbdaSJeremy L Thompson   @ref Advanced
14609d007619Sjeremylt **/
1461d1d35e2fSjeremylt int CeedBasisGetGrad1D(CeedBasis basis, const CeedScalar **grad_1d) {
1462d1d35e2fSjeremylt   if (!basis->tensor_basis)
14639d007619Sjeremylt     // LCOV_EXCL_START
1464e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR,
1465e15f9bd0SJeremy L Thompson                      "CeedBasis is not a tensor product basis.");
14669d007619Sjeremylt   // LCOV_EXCL_STOP
14679d007619Sjeremylt 
1468d1d35e2fSjeremylt   *grad_1d = basis->grad_1d;
1469e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
14709d007619Sjeremylt }
14719d007619Sjeremylt 
14729d007619Sjeremylt /**
147350c301a5SRezgar Shakeri   @brief Get divergence matrix of a CeedBasis
147450c301a5SRezgar Shakeri 
147550c301a5SRezgar Shakeri   @param basis     CeedBasis
147650c301a5SRezgar Shakeri   @param[out] div  Variable to store divergence matrix
147750c301a5SRezgar Shakeri 
147850c301a5SRezgar Shakeri   @return An error code: 0 - success, otherwise - failure
147950c301a5SRezgar Shakeri 
148050c301a5SRezgar Shakeri   @ref Advanced
148150c301a5SRezgar Shakeri **/
148250c301a5SRezgar Shakeri int CeedBasisGetDiv(CeedBasis basis, const CeedScalar **div) {
148350c301a5SRezgar Shakeri   if (!basis->div)
148450c301a5SRezgar Shakeri     // LCOV_EXCL_START
148550c301a5SRezgar Shakeri     return CeedError(basis->ceed, CEED_ERROR_MINOR,
148650c301a5SRezgar Shakeri                      "CeedBasis does not have divergence matrix.");
148750c301a5SRezgar Shakeri   // LCOV_EXCL_STOP
148850c301a5SRezgar Shakeri 
148950c301a5SRezgar Shakeri   *div = basis->div;
149050c301a5SRezgar Shakeri   return CEED_ERROR_SUCCESS;
149150c301a5SRezgar Shakeri }
149250c301a5SRezgar Shakeri 
149350c301a5SRezgar Shakeri /**
14947a982d89SJeremy L. Thompson   @brief Destroy a CeedBasis
14957a982d89SJeremy L. Thompson 
14967a982d89SJeremy L. Thompson   @param basis CeedBasis to destroy
14977a982d89SJeremy L. Thompson 
14987a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
14997a982d89SJeremy L. Thompson 
15007a982d89SJeremy L. Thompson   @ref User
15017a982d89SJeremy L. Thompson **/
15027a982d89SJeremy L. Thompson int CeedBasisDestroy(CeedBasis *basis) {
15037a982d89SJeremy L. Thompson   int ierr;
15047a982d89SJeremy L. Thompson 
1505d1d35e2fSjeremylt   if (!*basis || --(*basis)->ref_count > 0) return CEED_ERROR_SUCCESS;
15067a982d89SJeremy L. Thompson   if ((*basis)->Destroy) {
15077a982d89SJeremy L. Thompson     ierr = (*basis)->Destroy(*basis); CeedChk(ierr);
15087a982d89SJeremy L. Thompson   }
150934359f16Sjeremylt   if ((*basis)->contract) {
151034359f16Sjeremylt     ierr = CeedTensorContractDestroy(&(*basis)->contract); CeedChk(ierr);
151134359f16Sjeremylt   }
15127a982d89SJeremy L. Thompson   ierr = CeedFree(&(*basis)->interp); CeedChk(ierr);
1513d1d35e2fSjeremylt   ierr = CeedFree(&(*basis)->interp_1d); CeedChk(ierr);
15147a982d89SJeremy L. Thompson   ierr = CeedFree(&(*basis)->grad); CeedChk(ierr);
151550c301a5SRezgar Shakeri   ierr = CeedFree(&(*basis)->div); CeedChk(ierr);
1516d1d35e2fSjeremylt   ierr = CeedFree(&(*basis)->grad_1d); CeedChk(ierr);
1517d1d35e2fSjeremylt   ierr = CeedFree(&(*basis)->q_ref_1d); CeedChk(ierr);
1518d1d35e2fSjeremylt   ierr = CeedFree(&(*basis)->q_weight_1d); CeedChk(ierr);
15197a982d89SJeremy L. Thompson   ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr);
15207a982d89SJeremy L. Thompson   ierr = CeedFree(basis); CeedChk(ierr);
1521e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
15227a982d89SJeremy L. Thompson }
15237a982d89SJeremy L. Thompson 
15247a982d89SJeremy L. Thompson /**
1525b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre quadrature
1526b11c1e72Sjeremylt 
1527b11c1e72Sjeremylt   @param Q                 Number of quadrature points (integrates polynomials of
1528b11c1e72Sjeremylt                              degree 2*Q-1 exactly)
1529d1d35e2fSjeremylt   @param[out] q_ref_1d     Array of length Q to hold the abscissa on [-1, 1]
1530d1d35e2fSjeremylt   @param[out] q_weight_1d  Array of length Q to hold the weights
1531b11c1e72Sjeremylt 
1532b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1533dfdf5a53Sjeremylt 
1534dfdf5a53Sjeremylt   @ref Utility
1535b11c1e72Sjeremylt **/
1536d1d35e2fSjeremylt int CeedGaussQuadrature(CeedInt Q, CeedScalar *q_ref_1d,
1537d1d35e2fSjeremylt                         CeedScalar *q_weight_1d) {
1538d7b241e6Sjeremylt   // Allocate
1539d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0*atan(1.0);
1540d1d35e2fSjeremylt   // Build q_ref_1d, q_weight_1d
154192ae7e47SJeremy L Thompson   for (CeedInt i = 0; i <= Q/2; i++) {
1542d7b241e6Sjeremylt     // Guess
1543d7b241e6Sjeremylt     xi = cos(PI*(CeedScalar)(2*i+1)/((CeedScalar)(2*Q)));
1544d7b241e6Sjeremylt     // Pn(xi)
1545d7b241e6Sjeremylt     P0 = 1.0;
1546d7b241e6Sjeremylt     P1 = xi;
1547d7b241e6Sjeremylt     P2 = 0.0;
154892ae7e47SJeremy L Thompson     for (CeedInt j = 2; j <= Q; j++) {
1549d7b241e6Sjeremylt       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1550d7b241e6Sjeremylt       P0 = P1;
1551d7b241e6Sjeremylt       P1 = P2;
1552d7b241e6Sjeremylt     }
1553d7b241e6Sjeremylt     // First Newton Step
1554d7b241e6Sjeremylt     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1555d7b241e6Sjeremylt     xi = xi-P2/dP2;
1556d7b241e6Sjeremylt     // Newton to convergence
155792ae7e47SJeremy L Thompson     for (CeedInt k=0; k<100 && fabs(P2)>10*CEED_EPSILON; k++) {
1558d7b241e6Sjeremylt       P0 = 1.0;
1559d7b241e6Sjeremylt       P1 = xi;
156092ae7e47SJeremy L Thompson       for (CeedInt j = 2; j <= Q; j++) {
1561d7b241e6Sjeremylt         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1562d7b241e6Sjeremylt         P0 = P1;
1563d7b241e6Sjeremylt         P1 = P2;
1564d7b241e6Sjeremylt       }
1565d7b241e6Sjeremylt       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1566d7b241e6Sjeremylt       xi = xi-P2/dP2;
1567d7b241e6Sjeremylt     }
1568d7b241e6Sjeremylt     // Save xi, wi
1569d7b241e6Sjeremylt     wi = 2.0/((1.0-xi*xi)*dP2*dP2);
1570d1d35e2fSjeremylt     q_weight_1d[i] = wi;
1571d1d35e2fSjeremylt     q_weight_1d[Q-1-i] = wi;
1572d1d35e2fSjeremylt     q_ref_1d[i] = -xi;
1573d1d35e2fSjeremylt     q_ref_1d[Q-1-i]= xi;
1574d7b241e6Sjeremylt   }
1575e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1576d7b241e6Sjeremylt }
1577d7b241e6Sjeremylt 
1578b11c1e72Sjeremylt /**
1579b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre-Lobatto quadrature
1580b11c1e72Sjeremylt 
1581b11c1e72Sjeremylt   @param Q                 Number of quadrature points (integrates polynomials of
1582b11c1e72Sjeremylt                              degree 2*Q-3 exactly)
1583d1d35e2fSjeremylt   @param[out] q_ref_1d     Array of length Q to hold the abscissa on [-1, 1]
1584d1d35e2fSjeremylt   @param[out] q_weight_1d  Array of length Q to hold the weights
1585b11c1e72Sjeremylt 
1586b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1587dfdf5a53Sjeremylt 
1588dfdf5a53Sjeremylt   @ref Utility
1589b11c1e72Sjeremylt **/
1590d1d35e2fSjeremylt int CeedLobattoQuadrature(CeedInt Q, CeedScalar *q_ref_1d,
1591d1d35e2fSjeremylt                           CeedScalar *q_weight_1d) {
1592d7b241e6Sjeremylt   // Allocate
1593d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0*atan(1.0);
1594d1d35e2fSjeremylt   // Build q_ref_1d, q_weight_1d
1595d7b241e6Sjeremylt   // Set endpoints
159630a100c3SJed Brown   if (Q < 2)
1597b0d62198Sjeremylt     // LCOV_EXCL_START
1598e15f9bd0SJeremy L Thompson     return CeedError(NULL, CEED_ERROR_DIMENSION,
1599990fdeb6SJeremy L Thompson                      "Cannot create Lobatto quadrature with Q=%" CeedInt_FMT " < 2 points", Q);
1600b0d62198Sjeremylt   // LCOV_EXCL_STOP
1601d7b241e6Sjeremylt   wi = 2.0/((CeedScalar)(Q*(Q-1)));
1602d1d35e2fSjeremylt   if (q_weight_1d) {
1603d1d35e2fSjeremylt     q_weight_1d[0] = wi;
1604d1d35e2fSjeremylt     q_weight_1d[Q-1] = wi;
1605d7b241e6Sjeremylt   }
1606d1d35e2fSjeremylt   q_ref_1d[0] = -1.0;
1607d1d35e2fSjeremylt   q_ref_1d[Q-1] = 1.0;
1608d7b241e6Sjeremylt   // Interior
160992ae7e47SJeremy L Thompson   for (CeedInt i = 1; i <= (Q-1)/2; i++) {
1610d7b241e6Sjeremylt     // Guess
1611d7b241e6Sjeremylt     xi = cos(PI*(CeedScalar)(i)/(CeedScalar)(Q-1));
1612d7b241e6Sjeremylt     // Pn(xi)
1613d7b241e6Sjeremylt     P0 = 1.0;
1614d7b241e6Sjeremylt     P1 = xi;
1615d7b241e6Sjeremylt     P2 = 0.0;
161692ae7e47SJeremy L Thompson     for (CeedInt j = 2; j < Q; j++) {
1617d7b241e6Sjeremylt       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1618d7b241e6Sjeremylt       P0 = P1;
1619d7b241e6Sjeremylt       P1 = P2;
1620d7b241e6Sjeremylt     }
1621d7b241e6Sjeremylt     // First Newton step
1622d7b241e6Sjeremylt     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1623d7b241e6Sjeremylt     d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
1624d7b241e6Sjeremylt     xi = xi-dP2/d2P2;
1625d7b241e6Sjeremylt     // Newton to convergence
162692ae7e47SJeremy L Thompson     for (CeedInt k=0; k<100 && fabs(dP2)>10*CEED_EPSILON; k++) {
1627d7b241e6Sjeremylt       P0 = 1.0;
1628d7b241e6Sjeremylt       P1 = xi;
162992ae7e47SJeremy L Thompson       for (CeedInt j = 2; j < Q; j++) {
1630d7b241e6Sjeremylt         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1631d7b241e6Sjeremylt         P0 = P1;
1632d7b241e6Sjeremylt         P1 = P2;
1633d7b241e6Sjeremylt       }
1634d7b241e6Sjeremylt       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1635d7b241e6Sjeremylt       d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
1636d7b241e6Sjeremylt       xi = xi-dP2/d2P2;
1637d7b241e6Sjeremylt     }
1638d7b241e6Sjeremylt     // Save xi, wi
1639d7b241e6Sjeremylt     wi = 2.0/(((CeedScalar)(Q*(Q-1)))*P2*P2);
1640d1d35e2fSjeremylt     if (q_weight_1d) {
1641d1d35e2fSjeremylt       q_weight_1d[i] = wi;
1642d1d35e2fSjeremylt       q_weight_1d[Q-1-i] = wi;
1643d7b241e6Sjeremylt     }
1644d1d35e2fSjeremylt     q_ref_1d[i] = -xi;
1645d1d35e2fSjeremylt     q_ref_1d[Q-1-i]= xi;
1646d7b241e6Sjeremylt   }
1647e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1648d7b241e6Sjeremylt }
1649d7b241e6Sjeremylt 
1650dfdf5a53Sjeremylt /**
165195bb1877Svaleriabarra   @brief Return QR Factorization of a matrix
1652b11c1e72Sjeremylt 
165377645d7bSjeremylt   @param ceed         A Ceed context for error handling
165452bfb9bbSJeremy L Thompson   @param[in,out] mat  Row-major matrix to be factorized in place
165552bfb9bbSJeremy L Thompson   @param[in,out] tau  Vector of length m of scaling factors
1656b11c1e72Sjeremylt   @param m            Number of rows
1657b11c1e72Sjeremylt   @param n            Number of columns
1658b11c1e72Sjeremylt 
1659b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1660dfdf5a53Sjeremylt 
1661dfdf5a53Sjeremylt   @ref Utility
1662b11c1e72Sjeremylt **/
1663a7bd39daSjeremylt int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau,
1664d7b241e6Sjeremylt                         CeedInt m, CeedInt n) {
1665d7b241e6Sjeremylt   CeedScalar v[m];
1666d7b241e6Sjeremylt 
1667a7bd39daSjeremylt   // Check m >= n
1668a7bd39daSjeremylt   if (n > m)
1669c042f62fSJeremy L Thompson     // LCOV_EXCL_START
1670e15f9bd0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_UNSUPPORTED,
1671e15f9bd0SJeremy L Thompson                      "Cannot compute QR factorization with n > m");
1672c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
1673a7bd39daSjeremylt 
167452bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++) {
1675bde37e8cSJed Brown     if (i >= m-1) { // last row of matrix, no reflection needed
1676bde37e8cSJed Brown       tau[i] = 0.;
1677bde37e8cSJed Brown       break;
1678bde37e8cSJed Brown     }
1679d7b241e6Sjeremylt     // Calculate Householder vector, magnitude
1680d7b241e6Sjeremylt     CeedScalar sigma = 0.0;
1681d7b241e6Sjeremylt     v[i] = mat[i+n*i];
168252bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<m; j++) {
1683d7b241e6Sjeremylt       v[j] = mat[i+n*j];
1684d7b241e6Sjeremylt       sigma += v[j] * v[j];
1685d7b241e6Sjeremylt     }
1686d7b241e6Sjeremylt     CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m]
1687673160d7Sjeremylt     CeedScalar R_ii = -copysign(norm, v[i]);
1688673160d7Sjeremylt     v[i] -= R_ii;
1689d7b241e6Sjeremylt     // norm of v[i:m] after modification above and scaling below
1690d7b241e6Sjeremylt     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
1691d7b241e6Sjeremylt     //   tau = 2 / (norm*norm)
1692d7b241e6Sjeremylt     tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma);
16931d102b48SJeremy L Thompson     for (CeedInt j=i+1; j<m; j++)
16941d102b48SJeremy L Thompson       v[j] /= v[i];
1695d7b241e6Sjeremylt 
1696d7b241e6Sjeremylt     // Apply Householder reflector to lower right panel
1697d7b241e6Sjeremylt     CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1);
1698d7b241e6Sjeremylt     // Save v
1699673160d7Sjeremylt     mat[i+n*i] = R_ii;
17001d102b48SJeremy L Thompson     for (CeedInt j=i+1; j<m; j++)
1701d7b241e6Sjeremylt       mat[i+n*j] = v[j];
1702d7b241e6Sjeremylt   }
1703e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1704d7b241e6Sjeremylt }
1705d7b241e6Sjeremylt 
1706b11c1e72Sjeremylt /**
170752bfb9bbSJeremy L Thompson   @brief Return symmetric Schur decomposition of the symmetric matrix mat via
170852bfb9bbSJeremy L Thompson            symmetric QR factorization
170952bfb9bbSJeremy L Thompson 
171077645d7bSjeremylt   @param ceed         A Ceed context for error handling
171152bfb9bbSJeremy L Thompson   @param[in,out] mat  Row-major matrix to be factorized in place
1712460bf743SValeria Barra   @param[out] lambda  Vector of length n of eigenvalues
171352bfb9bbSJeremy L Thompson   @param n            Number of rows/columns
171452bfb9bbSJeremy L Thompson 
171552bfb9bbSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
171652bfb9bbSJeremy L Thompson 
171752bfb9bbSJeremy L Thompson   @ref Utility
171852bfb9bbSJeremy L Thompson **/
171903d18186Sjeremylt CeedPragmaOptimizeOff
172052bfb9bbSJeremy L Thompson int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat,
172152bfb9bbSJeremy L Thompson                                     CeedScalar *lambda, CeedInt n) {
172252bfb9bbSJeremy L Thompson   // Check bounds for clang-tidy
172352bfb9bbSJeremy L Thompson   if (n<2)
1724c042f62fSJeremy L Thompson     // LCOV_EXCL_START
1725e15f9bd0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_UNSUPPORTED,
1726c042f62fSJeremy L Thompson                      "Cannot compute symmetric Schur decomposition of scalars");
1727c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
172852bfb9bbSJeremy L Thompson 
1729673160d7Sjeremylt   CeedScalar v[n-1], tau[n-1], mat_T[n*n];
173052bfb9bbSJeremy L Thompson 
1731673160d7Sjeremylt   // Copy mat to mat_T and set mat to I
1732673160d7Sjeremylt   memcpy(mat_T, mat, n*n*sizeof(mat[0]));
173352bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
173452bfb9bbSJeremy L Thompson     for (CeedInt j=0; j<n; j++)
173552bfb9bbSJeremy L Thompson       mat[j+n*i] = (i==j) ? 1 : 0;
173652bfb9bbSJeremy L Thompson 
173752bfb9bbSJeremy L Thompson   // Reduce to tridiagonal
173852bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n-1; i++) {
173952bfb9bbSJeremy L Thompson     // Calculate Householder vector, magnitude
174052bfb9bbSJeremy L Thompson     CeedScalar sigma = 0.0;
1741673160d7Sjeremylt     v[i] = mat_T[i+n*(i+1)];
174252bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) {
1743673160d7Sjeremylt       v[j] = mat_T[i+n*(j+1)];
174452bfb9bbSJeremy L Thompson       sigma += v[j] * v[j];
174552bfb9bbSJeremy L Thompson     }
174652bfb9bbSJeremy L Thompson     CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:n-1]
1747673160d7Sjeremylt     CeedScalar R_ii = -copysign(norm, v[i]);
1748673160d7Sjeremylt     v[i] -= R_ii;
174952bfb9bbSJeremy L Thompson     // norm of v[i:m] after modification above and scaling below
175052bfb9bbSJeremy L Thompson     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
175152bfb9bbSJeremy L Thompson     //   tau = 2 / (norm*norm)
175280a9ef05SNatalie Beams     tau[i] = i == n - 2 ? 2 : 2 * v[i]*v[i] / (v[i]*v[i] + sigma);
1753fb551037Sjeremylt     for (CeedInt j=i+1; j<n-1; j++)
1754fb551037Sjeremylt       v[j] /= v[i];
175552bfb9bbSJeremy L Thompson 
175652bfb9bbSJeremy L Thompson     // Update sub and super diagonal
175752bfb9bbSJeremy L Thompson     for (CeedInt j=i+2; j<n; j++) {
1758673160d7Sjeremylt       mat_T[i+n*j] = 0; mat_T[j+n*i] = 0;
175952bfb9bbSJeremy L Thompson     }
176052bfb9bbSJeremy L Thompson     // Apply symmetric Householder reflector to lower right panel
1761673160d7Sjeremylt     CeedHouseholderReflect(&mat_T[(i+1)+n*(i+1)], &v[i], tau[i],
176252bfb9bbSJeremy L Thompson                            n-(i+1), n-(i+1), n, 1);
1763673160d7Sjeremylt     CeedHouseholderReflect(&mat_T[(i+1)+n*(i+1)], &v[i], tau[i],
176452bfb9bbSJeremy L Thompson                            n-(i+1), n-(i+1), 1, n);
1765673160d7Sjeremylt 
176652bfb9bbSJeremy L Thompson     // Save v
1767673160d7Sjeremylt     mat_T[i+n*(i+1)] = R_ii;
1768673160d7Sjeremylt     mat_T[(i+1)+n*i] = R_ii;
176952bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) {
1770673160d7Sjeremylt       mat_T[i+n*(j+1)] = v[j];
177152bfb9bbSJeremy L Thompson     }
177252bfb9bbSJeremy L Thompson   }
177352bfb9bbSJeremy L Thompson   // Backwards accumulation of Q
177452bfb9bbSJeremy L Thompson   for (CeedInt i=n-2; i>=0; i--) {
177585cf89eaSjeremylt     if (tau[i] > 0.0) {
177652bfb9bbSJeremy L Thompson       v[i] = 1;
177752bfb9bbSJeremy L Thompson       for (CeedInt j=i+1; j<n-1; j++) {
1778673160d7Sjeremylt         v[j] = mat_T[i+n*(j+1)];
1779673160d7Sjeremylt         mat_T[i+n*(j+1)] = 0;
178052bfb9bbSJeremy L Thompson       }
178152bfb9bbSJeremy L Thompson       CeedHouseholderReflect(&mat[(i+1)+n*(i+1)], &v[i], tau[i],
178252bfb9bbSJeremy L Thompson                              n-(i+1), n-(i+1), n, 1);
178352bfb9bbSJeremy L Thompson     }
178485cf89eaSjeremylt   }
178552bfb9bbSJeremy L Thompson 
178652bfb9bbSJeremy L Thompson   // Reduce sub and super diagonal
1787673160d7Sjeremylt   CeedInt p = 0, q = 0, itr = 0, max_itr = n*n*n*n;
1788673160d7Sjeremylt   CeedScalar tol = CEED_EPSILON;
178952bfb9bbSJeremy L Thompson 
1790673160d7Sjeremylt   while (itr < max_itr) {
179152bfb9bbSJeremy L Thompson     // Update p, q, size of reduced portions of diagonal
179252bfb9bbSJeremy L Thompson     p = 0; q = 0;
179352bfb9bbSJeremy L Thompson     for (CeedInt i=n-2; i>=0; i--) {
1794673160d7Sjeremylt       if (fabs(mat_T[i+n*(i+1)]) < tol)
179552bfb9bbSJeremy L Thompson         q += 1;
179652bfb9bbSJeremy L Thompson       else
179752bfb9bbSJeremy L Thompson         break;
179852bfb9bbSJeremy L Thompson     }
1799673160d7Sjeremylt     for (CeedInt i=0; i<n-q-1; i++) {
1800673160d7Sjeremylt       if (fabs(mat_T[i+n*(i+1)]) < tol)
180152bfb9bbSJeremy L Thompson         p += 1;
180252bfb9bbSJeremy L Thompson       else
180352bfb9bbSJeremy L Thompson         break;
180452bfb9bbSJeremy L Thompson     }
180552bfb9bbSJeremy L Thompson     if (q == n-1) break; // Finished reducing
180652bfb9bbSJeremy L Thompson 
180752bfb9bbSJeremy L Thompson     // Reduce tridiagonal portion
1808673160d7Sjeremylt     CeedScalar t_nn = mat_T[(n-1-q)+n*(n-1-q)],
1809673160d7Sjeremylt                t_nnm1 = mat_T[(n-2-q)+n*(n-1-q)];
1810673160d7Sjeremylt     CeedScalar d = (mat_T[(n-2-q)+n*(n-2-q)] - t_nn)/2;
1811673160d7Sjeremylt     CeedScalar mu = t_nn - t_nnm1*t_nnm1 /
1812673160d7Sjeremylt                     (d + copysign(sqrt(d*d + t_nnm1*t_nnm1), d));
1813673160d7Sjeremylt     CeedScalar x = mat_T[p+n*p] - mu;
1814673160d7Sjeremylt     CeedScalar z = mat_T[p+n*(p+1)];
1815673160d7Sjeremylt     for (CeedInt k=p; k<n-q-1; k++) {
181652bfb9bbSJeremy L Thompson       // Compute Givens rotation
181752bfb9bbSJeremy L Thompson       CeedScalar c = 1, s = 0;
181852bfb9bbSJeremy L Thompson       if (fabs(z) > tol) {
181952bfb9bbSJeremy L Thompson         if (fabs(z) > fabs(x)) {
182052bfb9bbSJeremy L Thompson           CeedScalar tau = -x/z;
182152bfb9bbSJeremy L Thompson           s = 1/sqrt(1+tau*tau), c = s*tau;
182252bfb9bbSJeremy L Thompson         } else {
182352bfb9bbSJeremy L Thompson           CeedScalar tau = -z/x;
182452bfb9bbSJeremy L Thompson           c = 1/sqrt(1+tau*tau), s = c*tau;
182552bfb9bbSJeremy L Thompson         }
182652bfb9bbSJeremy L Thompson       }
182752bfb9bbSJeremy L Thompson 
182852bfb9bbSJeremy L Thompson       // Apply Givens rotation to T
1829673160d7Sjeremylt       CeedGivensRotation(mat_T, c, s, CEED_NOTRANSPOSE, k, k+1, n, n);
1830673160d7Sjeremylt       CeedGivensRotation(mat_T, c, s, CEED_TRANSPOSE, k, k+1, n, n);
183152bfb9bbSJeremy L Thompson 
183252bfb9bbSJeremy L Thompson       // Apply Givens rotation to Q
183352bfb9bbSJeremy L Thompson       CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k+1, n, n);
183452bfb9bbSJeremy L Thompson 
183552bfb9bbSJeremy L Thompson       // Update x, z
183652bfb9bbSJeremy L Thompson       if (k < n-q-2) {
1837673160d7Sjeremylt         x = mat_T[k+n*(k+1)];
1838673160d7Sjeremylt         z = mat_T[k+n*(k+2)];
183952bfb9bbSJeremy L Thompson       }
184052bfb9bbSJeremy L Thompson     }
184152bfb9bbSJeremy L Thompson     itr++;
184252bfb9bbSJeremy L Thompson   }
1843673160d7Sjeremylt 
184452bfb9bbSJeremy L Thompson   // Save eigenvalues
184552bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
1846673160d7Sjeremylt     lambda[i] = mat_T[i+n*i];
184752bfb9bbSJeremy L Thompson 
184852bfb9bbSJeremy L Thompson   // Check convergence
1849673160d7Sjeremylt   if (itr == max_itr && q < n-1)
1850c042f62fSJeremy L Thompson     // LCOV_EXCL_START
1851e15f9bd0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_MINOR,
1852e15f9bd0SJeremy L Thompson                      "Symmetric QR failed to converge");
1853c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
1854e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
185552bfb9bbSJeremy L Thompson }
185603d18186Sjeremylt CeedPragmaOptimizeOn
185752bfb9bbSJeremy L Thompson 
185852bfb9bbSJeremy L Thompson /**
185952bfb9bbSJeremy L Thompson   @brief Return Simultaneous Diagonalization of two matrices. This solves the
186052bfb9bbSJeremy L Thompson            generalized eigenvalue problem A x = lambda B x, where A and B
186152bfb9bbSJeremy L Thompson            are symmetric and B is positive definite. We generate the matrix X
186252bfb9bbSJeremy L Thompson            and vector Lambda such that X^T A X = Lambda and X^T B X = I. This
186352bfb9bbSJeremy L Thompson            is equivalent to the LAPACK routine 'sygv' with TYPE = 1.
186452bfb9bbSJeremy L Thompson 
186577645d7bSjeremylt   @param ceed         A Ceed context for error handling
1866d1d35e2fSjeremylt   @param[in] mat_A    Row-major matrix to be factorized with eigenvalues
1867d1d35e2fSjeremylt   @param[in] mat_B    Row-major matrix to be factorized to identity
1868d3331725Sjeremylt   @param[out] mat_X   Row-major orthogonal matrix
1869460bf743SValeria Barra   @param[out] lambda  Vector of length n of generalized eigenvalues
187052bfb9bbSJeremy L Thompson   @param n            Number of rows/columns
187152bfb9bbSJeremy L Thompson 
187252bfb9bbSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
187352bfb9bbSJeremy L Thompson 
187452bfb9bbSJeremy L Thompson   @ref Utility
187552bfb9bbSJeremy L Thompson **/
187603d18186Sjeremylt CeedPragmaOptimizeOff
1877d1d35e2fSjeremylt int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *mat_A,
1878d3331725Sjeremylt                                     CeedScalar *mat_B, CeedScalar *mat_X,
187952bfb9bbSJeremy L Thompson                                     CeedScalar *lambda, CeedInt n) {
188052bfb9bbSJeremy L Thompson   int ierr;
1881d3331725Sjeremylt   CeedScalar *mat_C, *mat_G, *vec_D;
188278464608Sjeremylt   ierr = CeedCalloc(n*n, &mat_C); CeedChk(ierr);
188378464608Sjeremylt   ierr = CeedCalloc(n*n, &mat_G); CeedChk(ierr);
1884d3331725Sjeremylt   ierr = CeedCalloc(n, &vec_D); CeedChk(ierr);
188552bfb9bbSJeremy L Thompson 
188652bfb9bbSJeremy L Thompson   // Compute B = G D G^T
188778464608Sjeremylt   memcpy(mat_G, mat_B, n*n*sizeof(mat_B[0]));
1888d3331725Sjeremylt   ierr = CeedSymmetricSchurDecomposition(ceed, mat_G, vec_D, n); CeedChk(ierr);
188952bfb9bbSJeremy L Thompson 
189085cf89eaSjeremylt   // Sort eigenvalues
189185cf89eaSjeremylt   for (CeedInt i=n-1; i>=0; i--)
189285cf89eaSjeremylt     for (CeedInt j=0; j<i; j++) {
189385cf89eaSjeremylt       if (fabs(vec_D[j]) > fabs(vec_D[j+1])) {
189485cf89eaSjeremylt         CeedScalar temp;
189585cf89eaSjeremylt         temp = vec_D[j]; vec_D[j] = vec_D[j+1]; vec_D[j+1] = temp;
189685cf89eaSjeremylt         for (CeedInt k=0; k<n; k++) {
189785cf89eaSjeremylt           temp = mat_G[k*n+j]; mat_G[k*n+j] = mat_G[k*n+j+1]; mat_G[k*n+j+1] = temp;
189885cf89eaSjeremylt         }
189985cf89eaSjeremylt       }
190085cf89eaSjeremylt     }
190185cf89eaSjeremylt 
1902fb551037Sjeremylt   // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T
1903fb551037Sjeremylt   //           = D^-1/2 G^T A G D^-1/2
1904d3331725Sjeremylt   // -- D = D^-1/2
190552bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
1906d3331725Sjeremylt     vec_D[i] = 1./sqrt(vec_D[i]);
1907d3331725Sjeremylt   // -- G = G D^-1/2
1908d3331725Sjeremylt   // -- C = D^-1/2 G^T
1909d3331725Sjeremylt   for (CeedInt i=0; i<n; i++)
1910d3331725Sjeremylt     for (CeedInt j=0; j<n; j++) {
1911673160d7Sjeremylt       mat_G[i*n+j] *= vec_D[j];
1912673160d7Sjeremylt       mat_C[j*n+i]  = mat_G[i*n+j];
1913d3331725Sjeremylt     }
1914673160d7Sjeremylt   // -- X = (D^-1/2 G^T) A
1915ed9e99e6SJeremy L Thompson   ierr = CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_C,
1916d3331725Sjeremylt                                   (const CeedScalar *)mat_A, mat_X, n, n, n);
19179289e5bfSjeremylt   CeedChk(ierr);
1918673160d7Sjeremylt   // -- C = (D^-1/2 G^T A) (G D^-1/2)
1919ed9e99e6SJeremy L Thompson   ierr = CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_X,
192078464608Sjeremylt                                   (const CeedScalar *)mat_G, mat_C, n, n, n);
19219289e5bfSjeremylt   CeedChk(ierr);
192252bfb9bbSJeremy L Thompson 
192352bfb9bbSJeremy L Thompson   // Compute Q^T C Q = lambda
1924d1d35e2fSjeremylt   ierr = CeedSymmetricSchurDecomposition(ceed, mat_C, lambda, n); CeedChk(ierr);
192552bfb9bbSJeremy L Thompson 
192685cf89eaSjeremylt   // Sort eigenvalues
192785cf89eaSjeremylt   for (CeedInt i=n-1; i>=0; i--)
192885cf89eaSjeremylt     for (CeedInt j=0; j<i; j++) {
192985cf89eaSjeremylt       if (fabs(lambda[j]) > fabs(lambda[j+1])) {
193085cf89eaSjeremylt         CeedScalar temp;
193185cf89eaSjeremylt         temp = lambda[j]; lambda[j] = lambda[j+1]; lambda[j+1] = temp;
193285cf89eaSjeremylt         for (CeedInt k=0; k<n; k++) {
193385cf89eaSjeremylt           temp = mat_C[k*n+j]; mat_C[k*n+j] = mat_C[k*n+j+1]; mat_C[k*n+j+1] = temp;
193485cf89eaSjeremylt         }
193585cf89eaSjeremylt       }
193685cf89eaSjeremylt     }
193785cf89eaSjeremylt 
1938d3331725Sjeremylt   // Set X = (G D^1/2)^-T Q
1939fb551037Sjeremylt   //       = G D^-1/2 Q
1940ed9e99e6SJeremy L Thompson   ierr = CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_G,
1941d3331725Sjeremylt                                   (const CeedScalar *)mat_C, mat_X, n, n, n);
19429289e5bfSjeremylt   CeedChk(ierr);
194378464608Sjeremylt 
194478464608Sjeremylt   // Cleanup
194578464608Sjeremylt   ierr = CeedFree(&mat_C); CeedChk(ierr);
194678464608Sjeremylt   ierr = CeedFree(&mat_G); CeedChk(ierr);
1947d3331725Sjeremylt   ierr = CeedFree(&vec_D); CeedChk(ierr);
1948e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
194952bfb9bbSJeremy L Thompson }
195003d18186Sjeremylt CeedPragmaOptimizeOn
195152bfb9bbSJeremy L Thompson 
1952d7b241e6Sjeremylt /// @}
1953