xref: /libCEED/rust/libceed-sys/c-src/interface/ceed-basis.c (revision 784646086dc2cd22411d95107c6b787565214f18)
1d7b241e6Sjeremylt // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2d7b241e6Sjeremylt // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3d7b241e6Sjeremylt // reserved. See files LICENSE and NOTICE for details.
4d7b241e6Sjeremylt //
5d7b241e6Sjeremylt // This file is part of CEED, a collection of benchmarks, miniapps, software
6d7b241e6Sjeremylt // libraries and APIs for efficient high-order finite element and spectral
7d7b241e6Sjeremylt // element discretizations for exascale applications. For more information and
8d7b241e6Sjeremylt // source code availability see http://github.com/ceed.
9d7b241e6Sjeremylt //
10d7b241e6Sjeremylt // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11d7b241e6Sjeremylt // a collaborative effort of two U.S. Department of Energy organizations (Office
12d7b241e6Sjeremylt // of Science and the National Nuclear Security Administration) responsible for
13d7b241e6Sjeremylt // the planning and preparation of a capable exascale ecosystem, including
14d7b241e6Sjeremylt // software, applications, hardware, advanced system engineering and early
15d7b241e6Sjeremylt // testbed platforms, in support of the nation's exascale computing imperative.
16d7b241e6Sjeremylt 
17ec3da8bcSJed Brown #include <ceed/ceed.h>
18ec3da8bcSJed Brown #include <ceed/backend.h>
193d576824SJeremy L Thompson #include <ceed-impl.h>
20d7b241e6Sjeremylt #include <math.h>
213d576824SJeremy L Thompson #include <stdbool.h>
22d7b241e6Sjeremylt #include <stdio.h>
23d7b241e6Sjeremylt #include <string.h>
24d7b241e6Sjeremylt 
257a982d89SJeremy L. Thompson /// @file
267a982d89SJeremy L. Thompson /// Implementation of CeedBasis interfaces
277a982d89SJeremy L. Thompson 
28d7b241e6Sjeremylt /// @cond DOXYGEN_SKIP
29783c99b3SValeria Barra static struct CeedBasis_private ceed_basis_collocated;
30d7b241e6Sjeremylt /// @endcond
31d7b241e6Sjeremylt 
327a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser
337a982d89SJeremy L. Thompson /// @{
347a982d89SJeremy L. Thompson 
357a982d89SJeremy L. Thompson /// Indicate that the quadrature points are collocated with the nodes
367a982d89SJeremy L. Thompson const CeedBasis CEED_BASIS_COLLOCATED = &ceed_basis_collocated;
377a982d89SJeremy L. Thompson 
387a982d89SJeremy L. Thompson /// @}
397a982d89SJeremy L. Thompson 
407a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
417a982d89SJeremy L. Thompson /// CeedBasis Library Internal Functions
427a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
437a982d89SJeremy L. Thompson /// @addtogroup CeedBasisDeveloper
447a982d89SJeremy L. Thompson /// @{
457a982d89SJeremy L. Thompson 
467a982d89SJeremy L. Thompson /**
477a982d89SJeremy L. Thompson   @brief Compute Householder reflection
487a982d89SJeremy L. Thompson 
497a982d89SJeremy L. Thompson     Computes A = (I - b v v^T) A
507a982d89SJeremy L. Thompson     where A is an mxn matrix indexed as A[i*row + j*col]
517a982d89SJeremy L. Thompson 
527a982d89SJeremy L. Thompson   @param[in,out] A  Matrix to apply Householder reflection to, in place
537a982d89SJeremy L. Thompson   @param v          Householder vector
547a982d89SJeremy L. Thompson   @param b          Scaling factor
557a982d89SJeremy L. Thompson   @param m          Number of rows in A
567a982d89SJeremy L. Thompson   @param n          Number of columns in A
577a982d89SJeremy L. Thompson   @param row        Row stride
587a982d89SJeremy L. Thompson   @param col        Col stride
597a982d89SJeremy L. Thompson 
607a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
617a982d89SJeremy L. Thompson 
627a982d89SJeremy L. Thompson   @ref Developer
637a982d89SJeremy L. Thompson **/
647a982d89SJeremy L. Thompson static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v,
657a982d89SJeremy L. Thompson                                   CeedScalar b, CeedInt m, CeedInt n,
667a982d89SJeremy L. Thompson                                   CeedInt row, CeedInt col) {
677a982d89SJeremy L. Thompson   for (CeedInt j=0; j<n; j++) {
687a982d89SJeremy L. Thompson     CeedScalar w = A[0*row + j*col];
697a982d89SJeremy L. Thompson     for (CeedInt i=1; i<m; i++)
707a982d89SJeremy L. Thompson       w += v[i] * A[i*row + j*col];
717a982d89SJeremy L. Thompson     A[0*row + j*col] -= b * w;
727a982d89SJeremy L. Thompson     for (CeedInt i=1; i<m; i++)
737a982d89SJeremy L. Thompson       A[i*row + j*col] -= b * w * v[i];
747a982d89SJeremy L. Thompson   }
75e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
767a982d89SJeremy L. Thompson }
777a982d89SJeremy L. Thompson 
787a982d89SJeremy L. Thompson /**
797a982d89SJeremy L. Thompson   @brief Apply Householder Q matrix
807a982d89SJeremy L. Thompson 
817a982d89SJeremy L. Thompson     Compute A = Q A where Q is mxm and A is mxn.
827a982d89SJeremy L. Thompson 
837a982d89SJeremy L. Thompson   @param[in,out] A  Matrix to apply Householder Q to, in place
847a982d89SJeremy L. Thompson   @param Q          Householder Q matrix
857a982d89SJeremy L. Thompson   @param tau        Householder scaling factors
86d1d35e2fSjeremylt   @param t_mode    Transpose mode for application
877a982d89SJeremy L. Thompson   @param m          Number of rows in A
887a982d89SJeremy L. Thompson   @param n          Number of columns in A
897a982d89SJeremy L. Thompson   @param k          Number of elementary reflectors in Q, k<m
907a982d89SJeremy L. Thompson   @param row        Row stride in A
917a982d89SJeremy L. Thompson   @param col        Col stride in A
927a982d89SJeremy L. Thompson 
937a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
947a982d89SJeremy L. Thompson 
957a982d89SJeremy L. Thompson   @ref Developer
967a982d89SJeremy L. Thompson **/
97d99fa3c5SJeremy L Thompson int CeedHouseholderApplyQ(CeedScalar *A, const CeedScalar *Q,
98d1d35e2fSjeremylt                           const CeedScalar *tau, CeedTransposeMode t_mode,
997a982d89SJeremy L. Thompson                           CeedInt m, CeedInt n, CeedInt k,
1007a982d89SJeremy L. Thompson                           CeedInt row, CeedInt col) {
101e15f9bd0SJeremy L Thompson   int ierr;
102*78464608Sjeremylt   CeedScalar *v;
103*78464608Sjeremylt   ierr = CeedMalloc(m, &v); CeedChk(ierr);
1047a982d89SJeremy L. Thompson   for (CeedInt ii=0; ii<k; ii++) {
105d1d35e2fSjeremylt     CeedInt i = t_mode == CEED_TRANSPOSE ? ii : k-1-ii;
1067a982d89SJeremy L. Thompson     for (CeedInt j=i+1; j<m; j++)
1077a982d89SJeremy L. Thompson       v[j] = Q[j*k+i];
108d1d35e2fSjeremylt     // Apply Householder reflector (I - tau v v^T) collo_grad_1d^T
109e15f9bd0SJeremy L Thompson     ierr = CeedHouseholderReflect(&A[i*row], &v[i], tau[i], m-i, n, row, col);
110e15f9bd0SJeremy L Thompson     CeedChk(ierr);
1117a982d89SJeremy L. Thompson   }
112*78464608Sjeremylt   ierr = CeedFree(&v); CeedChk(ierr);
113e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1147a982d89SJeremy L. Thompson }
1157a982d89SJeremy L. Thompson 
1167a982d89SJeremy L. Thompson /**
1177a982d89SJeremy L. Thompson   @brief Compute Givens rotation
1187a982d89SJeremy L. Thompson 
1197a982d89SJeremy L. Thompson     Computes A = G A (or G^T A in transpose mode)
1207a982d89SJeremy L. Thompson     where A is an mxn matrix indexed as A[i*n + j*m]
1217a982d89SJeremy L. Thompson 
1227a982d89SJeremy L. Thompson   @param[in,out] A  Row major matrix to apply Givens rotation to, in place
1237a982d89SJeremy L. Thompson   @param c          Cosine factor
1247a982d89SJeremy L. Thompson   @param s          Sine factor
125d1d35e2fSjeremylt   @param t_mode     @ref CEED_NOTRANSPOSE to rotate the basis counter-clockwise,
1264c4400c7SValeria Barra                     which has the effect of rotating columns of A clockwise;
1274cc79fe7SJed Brown                     @ref CEED_TRANSPOSE for the opposite rotation
1287a982d89SJeremy L. Thompson   @param i          First row/column to apply rotation
1297a982d89SJeremy L. Thompson   @param k          Second row/column to apply rotation
1307a982d89SJeremy L. Thompson   @param m          Number of rows in A
1317a982d89SJeremy L. Thompson   @param n          Number of columns in A
1327a982d89SJeremy L. Thompson 
1337a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
1347a982d89SJeremy L. Thompson 
1357a982d89SJeremy L. Thompson   @ref Developer
1367a982d89SJeremy L. Thompson **/
1377a982d89SJeremy L. Thompson static int CeedGivensRotation(CeedScalar *A, CeedScalar c, CeedScalar s,
138d1d35e2fSjeremylt                               CeedTransposeMode t_mode, CeedInt i, CeedInt k,
1397a982d89SJeremy L. Thompson                               CeedInt m, CeedInt n) {
140d1d35e2fSjeremylt   CeedInt stride_j = 1, stride_ik = m, num_its = n;
141d1d35e2fSjeremylt   if (t_mode == CEED_NOTRANSPOSE) {
142d1d35e2fSjeremylt     stride_j = n; stride_ik = 1; num_its = m;
1437a982d89SJeremy L. Thompson   }
1447a982d89SJeremy L. Thompson 
1457a982d89SJeremy L. Thompson   // Apply rotation
146d1d35e2fSjeremylt   for (CeedInt j=0; j<num_its; j++) {
147d1d35e2fSjeremylt     CeedScalar tau1 = A[i*stride_ik+j*stride_j], tau2 = A[k*stride_ik+j*stride_j];
148d1d35e2fSjeremylt     A[i*stride_ik+j*stride_j] = c*tau1 - s*tau2;
149d1d35e2fSjeremylt     A[k*stride_ik+j*stride_j] = s*tau1 + c*tau2;
1507a982d89SJeremy L. Thompson   }
151e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1527a982d89SJeremy L. Thompson }
1537a982d89SJeremy L. Thompson 
1547a982d89SJeremy L. Thompson /**
1557a982d89SJeremy L. Thompson   @brief View an array stored in a CeedBasis
1567a982d89SJeremy L. Thompson 
1570a0da059Sjeremylt   @param[in] name      Name of array
158d1d35e2fSjeremylt   @param[in] fp_fmt    Printing format
1590a0da059Sjeremylt   @param[in] m         Number of rows in array
1600a0da059Sjeremylt   @param[in] n         Number of columns in array
1610a0da059Sjeremylt   @param[in] a         Array to be viewed
1620a0da059Sjeremylt   @param[in] stream    Stream to view to, e.g., stdout
1637a982d89SJeremy L. Thompson 
1647a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
1657a982d89SJeremy L. Thompson 
1667a982d89SJeremy L. Thompson   @ref Developer
1677a982d89SJeremy L. Thompson **/
168d1d35e2fSjeremylt static int CeedScalarView(const char *name, const char *fp_fmt, CeedInt m,
1697a982d89SJeremy L. Thompson                           CeedInt n, const CeedScalar *a, FILE *stream) {
1707a982d89SJeremy L. Thompson   for (int i=0; i<m; i++) {
1717a982d89SJeremy L. Thompson     if (m > 1)
1727a982d89SJeremy L. Thompson       fprintf(stream, "%12s[%d]:", name, i);
1737a982d89SJeremy L. Thompson     else
1747a982d89SJeremy L. Thompson       fprintf(stream, "%12s:", name);
1757a982d89SJeremy L. Thompson     for (int j=0; j<n; j++)
176d1d35e2fSjeremylt       fprintf(stream, fp_fmt, fabs(a[i*n+j]) > 1E-14 ? a[i*n+j] : 0);
1777a982d89SJeremy L. Thompson     fputs("\n", stream);
1787a982d89SJeremy L. Thompson   }
179e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1807a982d89SJeremy L. Thompson }
1817a982d89SJeremy L. Thompson 
1827a982d89SJeremy L. Thompson /// @}
1837a982d89SJeremy L. Thompson 
1847a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
1857a982d89SJeremy L. Thompson /// Ceed Backend API
1867a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
1877a982d89SJeremy L. Thompson /// @addtogroup CeedBasisBackend
1887a982d89SJeremy L. Thompson /// @{
1897a982d89SJeremy L. Thompson 
1907a982d89SJeremy L. Thompson /**
1917a982d89SJeremy L. Thompson   @brief Return collocated grad matrix
1927a982d89SJeremy L. Thompson 
1937a982d89SJeremy L. Thompson   @param basis               CeedBasis
194d1d35e2fSjeremylt   @param[out] collo_grad_1d  Row-major (Q_1d * Q_1d) matrix expressing derivatives of
1957a982d89SJeremy L. Thompson                                basis functions at quadrature points
1967a982d89SJeremy L. Thompson 
1977a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
1987a982d89SJeremy L. Thompson 
1997a982d89SJeremy L. Thompson   @ref Backend
2007a982d89SJeremy L. Thompson **/
201d1d35e2fSjeremylt int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *collo_grad_1d) {
2027a982d89SJeremy L. Thompson   int i, j, k;
2037a982d89SJeremy L. Thompson   Ceed ceed;
204d1d35e2fSjeremylt   CeedInt ierr, P_1d=(basis)->P_1d, Q_1d=(basis)->Q_1d;
205*78464608Sjeremylt   CeedScalar *interp_1d, *grad_1d, *tau;
2067a982d89SJeremy L. Thompson 
207d1d35e2fSjeremylt   ierr = CeedMalloc(Q_1d*P_1d, &interp_1d); CeedChk(ierr);
208d1d35e2fSjeremylt   ierr = CeedMalloc(Q_1d*P_1d, &grad_1d); CeedChk(ierr);
209*78464608Sjeremylt   ierr = CeedMalloc(Q_1d, &tau); CeedChk(ierr);
210d1d35e2fSjeremylt   memcpy(interp_1d, (basis)->interp_1d, Q_1d*P_1d*sizeof(basis)->interp_1d[0]);
211d1d35e2fSjeremylt   memcpy(grad_1d, (basis)->grad_1d, Q_1d*P_1d*sizeof(basis)->interp_1d[0]);
2127a982d89SJeremy L. Thompson 
213d1d35e2fSjeremylt   // QR Factorization, interp_1d = Q R
2147a982d89SJeremy L. Thompson   ierr = CeedBasisGetCeed(basis, &ceed); CeedChk(ierr);
215d1d35e2fSjeremylt   ierr = CeedQRFactorization(ceed, interp_1d, tau, Q_1d, P_1d); CeedChk(ierr);
216e15f9bd0SJeremy L Thompson   // Note: This function is for backend use, so all errors are terminal
217e15f9bd0SJeremy L Thompson   //   and we do not need to clean up memory on failure.
2187a982d89SJeremy L. Thompson 
219d1d35e2fSjeremylt   // Apply Rinv, collo_grad_1d = grad_1d Rinv
220d1d35e2fSjeremylt   for (i=0; i<Q_1d; i++) { // Row i
221d1d35e2fSjeremylt     collo_grad_1d[Q_1d*i] = grad_1d[P_1d*i]/interp_1d[0];
222d1d35e2fSjeremylt     for (j=1; j<P_1d; j++) { // Column j
223d1d35e2fSjeremylt       collo_grad_1d[j+Q_1d*i] = grad_1d[j+P_1d*i];
2247a982d89SJeremy L. Thompson       for (k=0; k<j; k++)
225d1d35e2fSjeremylt         collo_grad_1d[j+Q_1d*i] -= interp_1d[j+P_1d*k]*collo_grad_1d[k+Q_1d*i];
226d1d35e2fSjeremylt       collo_grad_1d[j+Q_1d*i] /= interp_1d[j+P_1d*j];
2277a982d89SJeremy L. Thompson     }
228d1d35e2fSjeremylt     for (j=P_1d; j<Q_1d; j++)
229d1d35e2fSjeremylt       collo_grad_1d[j+Q_1d*i] = 0;
2307a982d89SJeremy L. Thompson   }
2317a982d89SJeremy L. Thompson 
232d1d35e2fSjeremylt   // Apply Qtranspose, collograd = collo_grad Q_transpose
233d1d35e2fSjeremylt   ierr = CeedHouseholderApplyQ(collo_grad_1d, interp_1d, tau, CEED_NOTRANSPOSE,
234d1d35e2fSjeremylt                                Q_1d, Q_1d, P_1d, 1, Q_1d); CeedChk(ierr);
2357a982d89SJeremy L. Thompson 
236d1d35e2fSjeremylt   ierr = CeedFree(&interp_1d); CeedChk(ierr);
237d1d35e2fSjeremylt   ierr = CeedFree(&grad_1d); CeedChk(ierr);
238*78464608Sjeremylt   ierr = CeedFree(&tau); CeedChk(ierr);
239e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2407a982d89SJeremy L. Thompson }
2417a982d89SJeremy L. Thompson 
2427a982d89SJeremy L. Thompson /**
2437a982d89SJeremy L. Thompson   @brief Get Ceed associated with a CeedBasis
2447a982d89SJeremy L. Thompson 
2457a982d89SJeremy L. Thompson   @param basis      CeedBasis
2467a982d89SJeremy L. Thompson   @param[out] ceed  Variable to store Ceed
2477a982d89SJeremy L. Thompson 
2487a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
2497a982d89SJeremy L. Thompson 
2507a982d89SJeremy L. Thompson   @ref Backend
2517a982d89SJeremy L. Thompson **/
2527a982d89SJeremy L. Thompson int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) {
2537a982d89SJeremy L. Thompson   *ceed = basis->ceed;
254e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2557a982d89SJeremy L. Thompson }
2567a982d89SJeremy L. Thompson 
2577a982d89SJeremy L. Thompson /**
2587a982d89SJeremy L. Thompson   @brief Get tensor status for given CeedBasis
2597a982d89SJeremy L. Thompson 
2607a982d89SJeremy L. Thompson   @param basis           CeedBasis
261d1d35e2fSjeremylt   @param[out] is_tensor  Variable to store tensor status
2627a982d89SJeremy L. Thompson 
2637a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
2647a982d89SJeremy L. Thompson 
2657a982d89SJeremy L. Thompson   @ref Backend
2667a982d89SJeremy L. Thompson **/
267d1d35e2fSjeremylt int CeedBasisIsTensor(CeedBasis basis, bool *is_tensor) {
268d1d35e2fSjeremylt   *is_tensor = basis->tensor_basis;
269e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2707a982d89SJeremy L. Thompson }
2717a982d89SJeremy L. Thompson 
2727a982d89SJeremy L. Thompson /**
2737a982d89SJeremy L. Thompson   @brief Get backend data of a CeedBasis
2747a982d89SJeremy L. Thompson 
2757a982d89SJeremy L. Thompson   @param basis      CeedBasis
2767a982d89SJeremy L. Thompson   @param[out] data  Variable to store data
2777a982d89SJeremy L. Thompson 
2787a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
2797a982d89SJeremy L. Thompson 
2807a982d89SJeremy L. Thompson   @ref Backend
2817a982d89SJeremy L. Thompson **/
282777ff853SJeremy L Thompson int CeedBasisGetData(CeedBasis basis, void *data) {
283777ff853SJeremy L Thompson   *(void **)data = basis->data;
284e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2857a982d89SJeremy L. Thompson }
2867a982d89SJeremy L. Thompson 
2877a982d89SJeremy L. Thompson /**
2887a982d89SJeremy L. Thompson   @brief Set backend data of a CeedBasis
2897a982d89SJeremy L. Thompson 
2907a982d89SJeremy L. Thompson   @param[out] basis  CeedBasis
2917a982d89SJeremy L. Thompson   @param data        Data to set
2927a982d89SJeremy L. Thompson 
2937a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
2947a982d89SJeremy L. Thompson 
2957a982d89SJeremy L. Thompson   @ref Backend
2967a982d89SJeremy L. Thompson **/
297777ff853SJeremy L Thompson int CeedBasisSetData(CeedBasis basis, void *data) {
298777ff853SJeremy L Thompson   basis->data = data;
299e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3007a982d89SJeremy L. Thompson }
3017a982d89SJeremy L. Thompson 
3027a982d89SJeremy L. Thompson /**
3037a982d89SJeremy L. Thompson   @brief Get dimension for given CeedElemTopology
3047a982d89SJeremy L. Thompson 
3057a982d89SJeremy L. Thompson   @param topo      CeedElemTopology
3067a982d89SJeremy L. Thompson   @param[out] dim  Variable to store dimension of topology
3077a982d89SJeremy L. Thompson 
3087a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3097a982d89SJeremy L. Thompson 
3107a982d89SJeremy L. Thompson   @ref Backend
3117a982d89SJeremy L. Thompson **/
3127a982d89SJeremy L. Thompson int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) {
3137a982d89SJeremy L. Thompson   *dim = (CeedInt) topo >> 16;
314e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3157a982d89SJeremy L. Thompson }
3167a982d89SJeremy L. Thompson 
3177a982d89SJeremy L. Thompson /**
3187a982d89SJeremy L. Thompson   @brief Get CeedTensorContract of a CeedBasis
3197a982d89SJeremy L. Thompson 
3207a982d89SJeremy L. Thompson   @param basis          CeedBasis
3217a982d89SJeremy L. Thompson   @param[out] contract  Variable to store CeedTensorContract
3227a982d89SJeremy L. Thompson 
3237a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3247a982d89SJeremy L. Thompson 
3257a982d89SJeremy L. Thompson   @ref Backend
3267a982d89SJeremy L. Thompson **/
3277a982d89SJeremy L. Thompson int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) {
3287a982d89SJeremy L. Thompson   *contract = basis->contract;
329e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3307a982d89SJeremy L. Thompson }
3317a982d89SJeremy L. Thompson 
3327a982d89SJeremy L. Thompson /**
3337a982d89SJeremy L. Thompson   @brief Set CeedTensorContract of a CeedBasis
3347a982d89SJeremy L. Thompson 
3357a982d89SJeremy L. Thompson   @param[out] basis  CeedBasis
3367a982d89SJeremy L. Thompson   @param contract    CeedTensorContract to set
3377a982d89SJeremy L. Thompson 
3387a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3397a982d89SJeremy L. Thompson 
3407a982d89SJeremy L. Thompson   @ref Backend
3417a982d89SJeremy L. Thompson **/
3427a982d89SJeremy L. Thompson int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract *contract) {
3437a982d89SJeremy L. Thompson   basis->contract = *contract;
344e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3457a982d89SJeremy L. Thompson }
3467a982d89SJeremy L. Thompson 
3477a982d89SJeremy L. Thompson /**
3487a982d89SJeremy L. Thompson   @brief Return a reference implementation of matrix multiplication C = A B.
3497a982d89SJeremy L. Thompson            Note, this is a reference implementation for CPU CeedScalar pointers
3507a982d89SJeremy L. Thompson            that is not intended for high performance.
3517a982d89SJeremy L. Thompson 
3527a982d89SJeremy L. Thompson   @param ceed        A Ceed context for error handling
353d1d35e2fSjeremylt   @param[in] mat_A   Row-major matrix A
354d1d35e2fSjeremylt   @param[in] mat_B   Row-major matrix B
355d1d35e2fSjeremylt   @param[out] mat_C  Row-major output matrix C
3567a982d89SJeremy L. Thompson   @param m           Number of rows of C
3577a982d89SJeremy L. Thompson   @param n           Number of columns of C
3587a982d89SJeremy L. Thompson   @param kk          Number of columns of A/rows of B
3597a982d89SJeremy L. Thompson 
3607a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3617a982d89SJeremy L. Thompson 
3627a982d89SJeremy L. Thompson   @ref Utility
3637a982d89SJeremy L. Thompson **/
364d1d35e2fSjeremylt int CeedMatrixMultiply(Ceed ceed, const CeedScalar *mat_A,
365d1d35e2fSjeremylt                        const CeedScalar *mat_B, CeedScalar *mat_C, CeedInt m,
3667a982d89SJeremy L. Thompson                        CeedInt n, CeedInt kk) {
3677a982d89SJeremy L. Thompson   for (CeedInt i=0; i<m; i++)
3687a982d89SJeremy L. Thompson     for (CeedInt j=0; j<n; j++) {
3697a982d89SJeremy L. Thompson       CeedScalar sum = 0;
3707a982d89SJeremy L. Thompson       for (CeedInt k=0; k<kk; k++)
371d1d35e2fSjeremylt         sum += mat_A[k+i*kk]*mat_B[j+k*n];
372d1d35e2fSjeremylt       mat_C[j+i*n] = sum;
3737a982d89SJeremy L. Thompson     }
374e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3757a982d89SJeremy L. Thompson }
3767a982d89SJeremy L. Thompson 
3777a982d89SJeremy L. Thompson /// @}
3787a982d89SJeremy L. Thompson 
3797a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
3807a982d89SJeremy L. Thompson /// CeedBasis Public API
3817a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
3827a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser
383d7b241e6Sjeremylt /// @{
384d7b241e6Sjeremylt 
385b11c1e72Sjeremylt /**
38695bb1877Svaleriabarra   @brief Create a tensor-product basis for H^1 discretizations
387b11c1e72Sjeremylt 
388b11c1e72Sjeremylt   @param ceed        A Ceed object where the CeedBasis will be created
389b11c1e72Sjeremylt   @param dim         Topological dimension
390d1d35e2fSjeremylt   @param num_comp    Number of field components (1 for scalar fields)
391d1d35e2fSjeremylt   @param P_1d        Number of nodes in one dimension
392d1d35e2fSjeremylt   @param Q_1d        Number of quadrature points in one dimension
393d1d35e2fSjeremylt   @param interp_1d   Row-major (Q_1d * P_1d) matrix expressing the values of nodal
394b11c1e72Sjeremylt                        basis functions at quadrature points
395d1d35e2fSjeremylt   @param grad_1d     Row-major (Q_1d * P_1d) matrix expressing derivatives of nodal
396b11c1e72Sjeremylt                        basis functions at quadrature points
397d1d35e2fSjeremylt   @param q_ref_1d    Array of length Q_1d holding the locations of quadrature points
398b11c1e72Sjeremylt                        on the 1D reference element [-1, 1]
399d1d35e2fSjeremylt   @param q_weight_1d Array of length Q_1d holding the quadrature weights on the
400b11c1e72Sjeremylt                        reference element
401b11c1e72Sjeremylt   @param[out] basis  Address of the variable where the newly created
402b11c1e72Sjeremylt                        CeedBasis will be stored.
403b11c1e72Sjeremylt 
404b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
405dfdf5a53Sjeremylt 
4067a982d89SJeremy L. Thompson   @ref User
407b11c1e72Sjeremylt **/
408d1d35e2fSjeremylt int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt num_comp,
409d1d35e2fSjeremylt                             CeedInt P_1d, CeedInt Q_1d,
410d1d35e2fSjeremylt                             const CeedScalar *interp_1d,
411d1d35e2fSjeremylt                             const CeedScalar *grad_1d, const CeedScalar *q_ref_1d,
412d1d35e2fSjeremylt                             const CeedScalar *q_weight_1d, CeedBasis *basis) {
413d7b241e6Sjeremylt   int ierr;
414d7b241e6Sjeremylt 
4155fe0d4faSjeremylt   if (!ceed->BasisCreateTensorH1) {
4165fe0d4faSjeremylt     Ceed delegate;
417aefd8378Sjeremylt     ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr);
4185fe0d4faSjeremylt 
4195fe0d4faSjeremylt     if (!delegate)
420c042f62fSJeremy L Thompson       // LCOV_EXCL_START
421e15f9bd0SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_UNSUPPORTED,
422e15f9bd0SJeremy L Thompson                        "Backend does not support BasisCreateTensorH1");
423c042f62fSJeremy L Thompson     // LCOV_EXCL_STOP
4245fe0d4faSjeremylt 
425d1d35e2fSjeremylt     ierr = CeedBasisCreateTensorH1(delegate, dim, num_comp, P_1d,
426d1d35e2fSjeremylt                                    Q_1d, interp_1d, grad_1d, q_ref_1d,
427d1d35e2fSjeremylt                                    q_weight_1d, basis); CeedChk(ierr);
428e15f9bd0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
4295fe0d4faSjeremylt   }
430e15f9bd0SJeremy L Thompson 
431e15f9bd0SJeremy L Thompson   if (dim<1)
432e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
433e15f9bd0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
434e15f9bd0SJeremy L Thompson                      "Basis dimension must be a positive value");
435e15f9bd0SJeremy L Thompson   // LCOV_EXCL_STOP
436d1d35e2fSjeremylt   CeedElemTopology topo = dim == 1 ? CEED_LINE
437d1d35e2fSjeremylt                           : dim == 2 ? CEED_QUAD
438d1d35e2fSjeremylt                           : CEED_HEX;
439e15f9bd0SJeremy L Thompson 
440d7b241e6Sjeremylt   ierr = CeedCalloc(1, basis); CeedChk(ierr);
441d7b241e6Sjeremylt   (*basis)->ceed = ceed;
442d1d35e2fSjeremylt   ceed->ref_count++;
443d1d35e2fSjeremylt   (*basis)->ref_count = 1;
444d1d35e2fSjeremylt   (*basis)->tensor_basis = 1;
445d7b241e6Sjeremylt   (*basis)->dim = dim;
446d99fa3c5SJeremy L Thompson   (*basis)->topo = topo;
447d1d35e2fSjeremylt   (*basis)->num_comp = num_comp;
448d1d35e2fSjeremylt   (*basis)->P_1d = P_1d;
449d1d35e2fSjeremylt   (*basis)->Q_1d = Q_1d;
450d1d35e2fSjeremylt   (*basis)->P = CeedIntPow(P_1d, dim);
451d1d35e2fSjeremylt   (*basis)->Q = CeedIntPow(Q_1d, dim);
452d1d35e2fSjeremylt   ierr = CeedMalloc(Q_1d,&(*basis)->q_ref_1d); CeedChk(ierr);
453d1d35e2fSjeremylt   ierr = CeedMalloc(Q_1d,&(*basis)->q_weight_1d); CeedChk(ierr);
454d1d35e2fSjeremylt   memcpy((*basis)->q_ref_1d, q_ref_1d, Q_1d*sizeof(q_ref_1d[0]));
455d1d35e2fSjeremylt   memcpy((*basis)->q_weight_1d, q_weight_1d, Q_1d*sizeof(q_weight_1d[0]));
456d1d35e2fSjeremylt   ierr = CeedMalloc(Q_1d*P_1d,&(*basis)->interp_1d); CeedChk(ierr);
457d1d35e2fSjeremylt   ierr = CeedMalloc(Q_1d*P_1d,&(*basis)->grad_1d); CeedChk(ierr);
458d1d35e2fSjeremylt   memcpy((*basis)->interp_1d, interp_1d, Q_1d*P_1d*sizeof(interp_1d[0]));
459d1d35e2fSjeremylt   memcpy((*basis)->grad_1d, grad_1d, Q_1d*P_1d*sizeof(grad_1d[0]));
460d1d35e2fSjeremylt   ierr = ceed->BasisCreateTensorH1(dim, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d,
461d1d35e2fSjeremylt                                    q_weight_1d, *basis); CeedChk(ierr);
462e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
463d7b241e6Sjeremylt }
464d7b241e6Sjeremylt 
465b11c1e72Sjeremylt /**
46695bb1877Svaleriabarra   @brief Create a tensor-product Lagrange basis
467b11c1e72Sjeremylt 
468b11c1e72Sjeremylt   @param ceed        A Ceed object where the CeedBasis will be created
469b11c1e72Sjeremylt   @param dim         Topological dimension of element
470d1d35e2fSjeremylt   @param num_comp      Number of field components (1 for scalar fields)
471b11c1e72Sjeremylt   @param P           Number of Gauss-Lobatto nodes in one dimension.  The
472b11c1e72Sjeremylt                        polynomial degree of the resulting Q_k element is k=P-1.
473b11c1e72Sjeremylt   @param Q           Number of quadrature points in one dimension.
474d1d35e2fSjeremylt   @param quad_mode   Distribution of the Q quadrature points (affects order of
475b11c1e72Sjeremylt                        accuracy for the quadrature)
476b11c1e72Sjeremylt   @param[out] basis  Address of the variable where the newly created
477b11c1e72Sjeremylt                        CeedBasis will be stored.
478b11c1e72Sjeremylt 
479b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
480dfdf5a53Sjeremylt 
4817a982d89SJeremy L. Thompson   @ref User
482b11c1e72Sjeremylt **/
483d1d35e2fSjeremylt int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt num_comp,
484d1d35e2fSjeremylt                                     CeedInt P, CeedInt Q, CeedQuadMode quad_mode,
485692c2638Sjeremylt                                     CeedBasis *basis) {
486d7b241e6Sjeremylt   // Allocate
487e15f9bd0SJeremy L Thompson   int ierr, ierr2, i, j, k;
488d1d35e2fSjeremylt   CeedScalar c1, c2, c3, c4, dx, *nodes, *interp_1d, *grad_1d, *q_ref_1d,
489d1d35e2fSjeremylt              *q_weight_1d;
4904d537eeaSYohann 
4914d537eeaSYohann   if (dim<1)
492c042f62fSJeremy L Thompson     // LCOV_EXCL_START
493e15f9bd0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
494e15f9bd0SJeremy L Thompson                      "Basis dimension must be a positive value");
495c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
4964d537eeaSYohann 
497e15f9bd0SJeremy L Thompson   // Get Nodes and Weights
498d1d35e2fSjeremylt   ierr = CeedCalloc(P*Q, &interp_1d); CeedChk(ierr);
499d1d35e2fSjeremylt   ierr = CeedCalloc(P*Q, &grad_1d); CeedChk(ierr);
500d7b241e6Sjeremylt   ierr = CeedCalloc(P, &nodes); CeedChk(ierr);
501d1d35e2fSjeremylt   ierr = CeedCalloc(Q, &q_ref_1d); CeedChk(ierr);
502d1d35e2fSjeremylt   ierr = CeedCalloc(Q, &q_weight_1d); CeedChk(ierr);
503e15f9bd0SJeremy L Thompson   ierr = CeedLobattoQuadrature(P, nodes, NULL);
504e15f9bd0SJeremy L Thompson   if (ierr) { goto cleanup; } CeedChk(ierr);
505d1d35e2fSjeremylt   switch (quad_mode) {
506d7b241e6Sjeremylt   case CEED_GAUSS:
507d1d35e2fSjeremylt     ierr = CeedGaussQuadrature(Q, q_ref_1d, q_weight_1d);
508d7b241e6Sjeremylt     break;
509d7b241e6Sjeremylt   case CEED_GAUSS_LOBATTO:
510d1d35e2fSjeremylt     ierr = CeedLobattoQuadrature(Q, q_ref_1d, q_weight_1d);
511d7b241e6Sjeremylt     break;
512d7b241e6Sjeremylt   }
513e15f9bd0SJeremy L Thompson   if (ierr) { goto cleanup; } CeedChk(ierr);
514e15f9bd0SJeremy L Thompson 
515d7b241e6Sjeremylt   // Build B, D matrix
516d7b241e6Sjeremylt   // Fornberg, 1998
517d7b241e6Sjeremylt   for (i = 0; i  < Q; i++) {
518d7b241e6Sjeremylt     c1 = 1.0;
519d1d35e2fSjeremylt     c3 = nodes[0] - q_ref_1d[i];
520d1d35e2fSjeremylt     interp_1d[i*P+0] = 1.0;
521d7b241e6Sjeremylt     for (j = 1; j < P; j++) {
522d7b241e6Sjeremylt       c2 = 1.0;
523d7b241e6Sjeremylt       c4 = c3;
524d1d35e2fSjeremylt       c3 = nodes[j] - q_ref_1d[i];
525d7b241e6Sjeremylt       for (k = 0; k < j; k++) {
526d7b241e6Sjeremylt         dx = nodes[j] - nodes[k];
527d7b241e6Sjeremylt         c2 *= dx;
528d7b241e6Sjeremylt         if (k == j - 1) {
529d1d35e2fSjeremylt           grad_1d[i*P + j] = c1*(interp_1d[i*P + k] - c4*grad_1d[i*P + k]) / c2;
530d1d35e2fSjeremylt           interp_1d[i*P + j] = - c1*c4*interp_1d[i*P + k] / c2;
531d7b241e6Sjeremylt         }
532d1d35e2fSjeremylt         grad_1d[i*P + k] = (c3*grad_1d[i*P + k] - interp_1d[i*P + k]) / dx;
533d1d35e2fSjeremylt         interp_1d[i*P + k] = c3*interp_1d[i*P + k] / dx;
534d7b241e6Sjeremylt       }
535d7b241e6Sjeremylt       c1 = c2;
536d7b241e6Sjeremylt     }
537d7b241e6Sjeremylt   }
538d7b241e6Sjeremylt   //  // Pass to CeedBasisCreateTensorH1
539d1d35e2fSjeremylt   ierr = CeedBasisCreateTensorH1(ceed, dim, num_comp, P, Q, interp_1d, grad_1d,
540d1d35e2fSjeremylt                                  q_ref_1d,
541d1d35e2fSjeremylt                                  q_weight_1d, basis); CeedChk(ierr);
542e15f9bd0SJeremy L Thompson cleanup:
543d1d35e2fSjeremylt   ierr2 = CeedFree(&interp_1d); CeedChk(ierr2);
544d1d35e2fSjeremylt   ierr2 = CeedFree(&grad_1d); CeedChk(ierr2);
545e15f9bd0SJeremy L Thompson   ierr2 = CeedFree(&nodes); CeedChk(ierr2);
546d1d35e2fSjeremylt   ierr2 = CeedFree(&q_ref_1d); CeedChk(ierr2);
547d1d35e2fSjeremylt   ierr2 = CeedFree(&q_weight_1d); CeedChk(ierr2);
548e15f9bd0SJeremy L Thompson   CeedChk(ierr);
549e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
550d7b241e6Sjeremylt }
551d7b241e6Sjeremylt 
552b11c1e72Sjeremylt /**
55395bb1877Svaleriabarra   @brief Create a non tensor-product basis for H^1 discretizations
554a8de75f0Sjeremylt 
555a8de75f0Sjeremylt   @param ceed        A Ceed object where the CeedBasis will be created
556a8de75f0Sjeremylt   @param topo        Topology of element, e.g. hypercube, simplex, ect
557d1d35e2fSjeremylt   @param num_comp    Number of field components (1 for scalar fields)
558d1d35e2fSjeremylt   @param num_nodes   Total number of nodes
559d1d35e2fSjeremylt   @param num_qpts    Total number of quadrature points
560d1d35e2fSjeremylt   @param interp      Row-major (num_qpts * num_nodes) matrix expressing the values of
5618795c945Sjeremylt                        nodal basis functions at quadrature points
562d1d35e2fSjeremylt   @param grad        Row-major (num_qpts * dim * num_nodes) matrix expressing
5638795c945Sjeremylt                        derivatives of nodal basis functions at quadrature points
564d1d35e2fSjeremylt   @param q_ref       Array of length num_qpts holding the locations of quadrature
5658795c945Sjeremylt                        points on the reference element [-1, 1]
566d1d35e2fSjeremylt   @param q_weight    Array of length num_qpts holding the quadrature weights on the
567a8de75f0Sjeremylt                        reference element
568a8de75f0Sjeremylt   @param[out] basis  Address of the variable where the newly created
569a8de75f0Sjeremylt                        CeedBasis will be stored.
570a8de75f0Sjeremylt 
571a8de75f0Sjeremylt   @return An error code: 0 - success, otherwise - failure
572a8de75f0Sjeremylt 
5737a982d89SJeremy L. Thompson   @ref User
574a8de75f0Sjeremylt **/
575d1d35e2fSjeremylt int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt num_comp,
576d1d35e2fSjeremylt                       CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
577d1d35e2fSjeremylt                       const CeedScalar *grad, const CeedScalar *q_ref,
578d1d35e2fSjeremylt                       const CeedScalar *q_weight, CeedBasis *basis) {
579a8de75f0Sjeremylt   int ierr;
580d1d35e2fSjeremylt   CeedInt P = num_nodes, Q = num_qpts, dim = 0;
581a8de75f0Sjeremylt 
5825fe0d4faSjeremylt   if (!ceed->BasisCreateH1) {
5835fe0d4faSjeremylt     Ceed delegate;
584aefd8378Sjeremylt     ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr);
5855fe0d4faSjeremylt 
5865fe0d4faSjeremylt     if (!delegate)
587c042f62fSJeremy L Thompson       // LCOV_EXCL_START
588e15f9bd0SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_UNSUPPORTED,
589e15f9bd0SJeremy L Thompson                        "Backend does not support BasisCreateH1");
590c042f62fSJeremy L Thompson     // LCOV_EXCL_STOP
5915fe0d4faSjeremylt 
592d1d35e2fSjeremylt     ierr = CeedBasisCreateH1(delegate, topo, num_comp, num_nodes,
593d1d35e2fSjeremylt                              num_qpts, interp, grad, q_ref,
594d1d35e2fSjeremylt                              q_weight, basis); CeedChk(ierr);
595e15f9bd0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
5965fe0d4faSjeremylt   }
5975fe0d4faSjeremylt 
598a8de75f0Sjeremylt   ierr = CeedCalloc(1,basis); CeedChk(ierr);
599a8de75f0Sjeremylt 
600a8de75f0Sjeremylt   ierr = CeedBasisGetTopologyDimension(topo, &dim); CeedChk(ierr);
601a8de75f0Sjeremylt 
602a8de75f0Sjeremylt   (*basis)->ceed = ceed;
603d1d35e2fSjeremylt   ceed->ref_count++;
604d1d35e2fSjeremylt   (*basis)->ref_count = 1;
605d1d35e2fSjeremylt   (*basis)->tensor_basis = 0;
606a8de75f0Sjeremylt   (*basis)->dim = dim;
607d99fa3c5SJeremy L Thompson   (*basis)->topo = topo;
608d1d35e2fSjeremylt   (*basis)->num_comp = num_comp;
609a8de75f0Sjeremylt   (*basis)->P = P;
610a8de75f0Sjeremylt   (*basis)->Q = Q;
611d1d35e2fSjeremylt   ierr = CeedMalloc(Q*dim,&(*basis)->q_ref_1d); CeedChk(ierr);
612d1d35e2fSjeremylt   ierr = CeedMalloc(Q,&(*basis)->q_weight_1d); CeedChk(ierr);
613d1d35e2fSjeremylt   memcpy((*basis)->q_ref_1d, q_ref, Q*dim*sizeof(q_ref[0]));
614d1d35e2fSjeremylt   memcpy((*basis)->q_weight_1d, q_weight, Q*sizeof(q_weight[0]));
61500f91b2bSjeremylt   ierr = CeedMalloc(Q*P, &(*basis)->interp); CeedChk(ierr);
61600f91b2bSjeremylt   ierr = CeedMalloc(dim*Q*P, &(*basis)->grad); CeedChk(ierr);
61700f91b2bSjeremylt   memcpy((*basis)->interp, interp, Q*P*sizeof(interp[0]));
61800f91b2bSjeremylt   memcpy((*basis)->grad, grad, dim*Q*P*sizeof(grad[0]));
619d1d35e2fSjeremylt   ierr = ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, q_ref,
620d1d35e2fSjeremylt                              q_weight, *basis); CeedChk(ierr);
621e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
622a8de75f0Sjeremylt }
623a8de75f0Sjeremylt 
624a8de75f0Sjeremylt /**
6257a982d89SJeremy L. Thompson   @brief View a CeedBasis
6267a982d89SJeremy L. Thompson 
6277a982d89SJeremy L. Thompson   @param basis   CeedBasis to view
6287a982d89SJeremy L. Thompson   @param stream  Stream to view to, e.g., stdout
6297a982d89SJeremy L. Thompson 
6307a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
6317a982d89SJeremy L. Thompson 
6327a982d89SJeremy L. Thompson   @ref User
6337a982d89SJeremy L. Thompson **/
6347a982d89SJeremy L. Thompson int CeedBasisView(CeedBasis basis, FILE *stream) {
6357a982d89SJeremy L. Thompson   int ierr;
6367a982d89SJeremy L. Thompson 
637d1d35e2fSjeremylt   if (basis->tensor_basis) {
638d1d35e2fSjeremylt     fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P_1d,
639d1d35e2fSjeremylt             basis->Q_1d);
640d1d35e2fSjeremylt     ierr = CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_ref_1d,
6417a982d89SJeremy L. Thompson                           stream); CeedChk(ierr);
642d1d35e2fSjeremylt     ierr = CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q_1d,
643d1d35e2fSjeremylt                           basis->q_weight_1d, stream); CeedChk(ierr);
644d1d35e2fSjeremylt     ierr = CeedScalarView("interp1d", "\t% 12.8f", basis->Q_1d, basis->P_1d,
645d1d35e2fSjeremylt                           basis->interp_1d, stream); CeedChk(ierr);
646d1d35e2fSjeremylt     ierr = CeedScalarView("grad1d", "\t% 12.8f", basis->Q_1d, basis->P_1d,
647d1d35e2fSjeremylt                           basis->grad_1d, stream); CeedChk(ierr);
6487a982d89SJeremy L. Thompson   } else {
6497a982d89SJeremy L. Thompson     fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P,
6507a982d89SJeremy L. Thompson             basis->Q);
6517a982d89SJeremy L. Thompson     ierr = CeedScalarView("qref", "\t% 12.8f", 1, basis->Q*basis->dim,
652d1d35e2fSjeremylt                           basis->q_ref_1d,
6537a982d89SJeremy L. Thompson                           stream); CeedChk(ierr);
654d1d35e2fSjeremylt     ierr = CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->q_weight_1d,
6557a982d89SJeremy L. Thompson                           stream); CeedChk(ierr);
6567a982d89SJeremy L. Thompson     ierr = CeedScalarView("interp", "\t% 12.8f", basis->Q, basis->P,
6577a982d89SJeremy L. Thompson                           basis->interp, stream); CeedChk(ierr);
6587a982d89SJeremy L. Thompson     ierr = CeedScalarView("grad", "\t% 12.8f", basis->dim*basis->Q, basis->P,
6597a982d89SJeremy L. Thompson                           basis->grad, stream); CeedChk(ierr);
6607a982d89SJeremy L. Thompson   }
661e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
6627a982d89SJeremy L. Thompson }
6637a982d89SJeremy L. Thompson 
6647a982d89SJeremy L. Thompson /**
6657a982d89SJeremy L. Thompson   @brief Apply basis evaluation from nodes to quadrature points or vice versa
6667a982d89SJeremy L. Thompson 
6677a982d89SJeremy L. Thompson   @param basis     CeedBasis to evaluate
668d1d35e2fSjeremylt   @param num_elem  The number of elements to apply the basis evaluation to;
6697a982d89SJeremy L. Thompson                      the backend will specify the ordering in
6704cc79fe7SJed Brown                      CeedElemRestrictionCreateBlocked()
671d1d35e2fSjeremylt   @param t_mode    \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature
6727a982d89SJeremy L. Thompson                      points, \ref CEED_TRANSPOSE to apply the transpose, mapping
6737a982d89SJeremy L. Thompson                      from quadrature points to nodes
674d1d35e2fSjeremylt   @param eval_mode \ref CEED_EVAL_NONE to use values directly,
6757a982d89SJeremy L. Thompson                      \ref CEED_EVAL_INTERP to use interpolated values,
6767a982d89SJeremy L. Thompson                      \ref CEED_EVAL_GRAD to use gradients,
6777a982d89SJeremy L. Thompson                      \ref CEED_EVAL_WEIGHT to use quadrature weights.
6787a982d89SJeremy L. Thompson   @param[in] u     Input CeedVector
6797a982d89SJeremy L. Thompson   @param[out] v    Output CeedVector
6807a982d89SJeremy L. Thompson 
6817a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
6827a982d89SJeremy L. Thompson 
6837a982d89SJeremy L. Thompson   @ref User
6847a982d89SJeremy L. Thompson **/
685d1d35e2fSjeremylt int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode,
686d1d35e2fSjeremylt                    CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
6877a982d89SJeremy L. Thompson   int ierr;
688d1d35e2fSjeremylt   CeedInt u_length = 0, v_length, dim, num_comp, num_nodes, num_qpts;
689e15f9bd0SJeremy L Thompson   ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr);
690d1d35e2fSjeremylt   ierr = CeedBasisGetNumComponents(basis, &num_comp); CeedChk(ierr);
691d1d35e2fSjeremylt   ierr = CeedBasisGetNumNodes(basis, &num_nodes); CeedChk(ierr);
692d1d35e2fSjeremylt   ierr = CeedBasisGetNumQuadraturePoints(basis, &num_qpts); CeedChk(ierr);
693d1d35e2fSjeremylt   ierr = CeedVectorGetLength(v, &v_length); CeedChk(ierr);
6947a982d89SJeremy L. Thompson   if (u) {
695d1d35e2fSjeremylt     ierr = CeedVectorGetLength(u, &u_length); CeedChk(ierr);
6967a982d89SJeremy L. Thompson   }
6977a982d89SJeremy L. Thompson 
698e15f9bd0SJeremy L Thompson   if (!basis->Apply)
699e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
700e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_UNSUPPORTED,
701e15f9bd0SJeremy L Thompson                      "Backend does not support BasisApply");
702e15f9bd0SJeremy L Thompson   // LCOV_EXCL_STOP
703e15f9bd0SJeremy L Thompson 
704e15f9bd0SJeremy L Thompson   // Check compatibility of topological and geometrical dimensions
705d1d35e2fSjeremylt   if ((t_mode == CEED_TRANSPOSE && (v_length%num_nodes != 0 ||
706d1d35e2fSjeremylt                                     u_length%num_qpts != 0)) ||
707d1d35e2fSjeremylt       (t_mode == CEED_NOTRANSPOSE && (u_length%num_nodes != 0 ||
708d1d35e2fSjeremylt                                       v_length%num_qpts != 0)))
709e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_DIMENSION,
710e15f9bd0SJeremy L Thompson                      "Length of input/output vectors "
7117a982d89SJeremy L. Thompson                      "incompatible with basis dimensions");
7127a982d89SJeremy L. Thompson 
713e15f9bd0SJeremy L Thompson   // Check vector lengths to prevent out of bounds issues
714d1d35e2fSjeremylt   bool bad_dims = false;
715d1d35e2fSjeremylt   switch (eval_mode) {
716e15f9bd0SJeremy L Thompson   case CEED_EVAL_NONE:
717d1d35e2fSjeremylt   case CEED_EVAL_INTERP: bad_dims =
718d1d35e2fSjeremylt       ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts ||
719d1d35e2fSjeremylt                                      v_length < num_elem*num_comp*num_nodes)) ||
720d1d35e2fSjeremylt        (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp ||
721d1d35e2fSjeremylt                                        u_length < num_elem*num_comp*num_nodes)));
722e15f9bd0SJeremy L Thompson     break;
723d1d35e2fSjeremylt   case CEED_EVAL_GRAD: bad_dims =
724d1d35e2fSjeremylt       ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts*dim ||
725d1d35e2fSjeremylt                                      v_length < num_elem*num_comp*num_nodes)) ||
726d1d35e2fSjeremylt        (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp*dim ||
727d1d35e2fSjeremylt                                        u_length < num_elem*num_comp*num_nodes)));
728e15f9bd0SJeremy L Thompson     break;
729e15f9bd0SJeremy L Thompson   case CEED_EVAL_WEIGHT:
730d1d35e2fSjeremylt     bad_dims = v_length < num_elem*num_qpts;
731e15f9bd0SJeremy L Thompson     break;
732e15f9bd0SJeremy L Thompson   // LCOV_EXCL_START
733d1d35e2fSjeremylt   case CEED_EVAL_DIV: bad_dims =
734d1d35e2fSjeremylt       ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts ||
735d1d35e2fSjeremylt                                      v_length < num_elem*num_comp*num_nodes)) ||
736d1d35e2fSjeremylt        (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp ||
737d1d35e2fSjeremylt                                        u_length < num_elem*num_comp*num_nodes)));
738e15f9bd0SJeremy L Thompson     break;
739d1d35e2fSjeremylt   case CEED_EVAL_CURL: bad_dims =
740d1d35e2fSjeremylt       ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts ||
741d1d35e2fSjeremylt                                      v_length < num_elem*num_comp*num_nodes)) ||
742d1d35e2fSjeremylt        (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp ||
743d1d35e2fSjeremylt                                        u_length < num_elem*num_comp*num_nodes)));
744e15f9bd0SJeremy L Thompson     break;
745e15f9bd0SJeremy L Thompson     // LCOV_EXCL_STOP
746e15f9bd0SJeremy L Thompson   }
747d1d35e2fSjeremylt   if (bad_dims)
748e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
749e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_DIMENSION,
750d1d35e2fSjeremylt                      "Input/output vectors too short for basis and evaluation mode");
751e15f9bd0SJeremy L Thompson   // LCOV_EXCL_STOP
752e15f9bd0SJeremy L Thompson 
753d1d35e2fSjeremylt   ierr = basis->Apply(basis, num_elem, t_mode, eval_mode, u, v); CeedChk(ierr);
754e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
7557a982d89SJeremy L. Thompson }
7567a982d89SJeremy L. Thompson 
7577a982d89SJeremy L. Thompson /**
7589d007619Sjeremylt   @brief Get dimension for given CeedBasis
7599d007619Sjeremylt 
7609d007619Sjeremylt   @param basis     CeedBasis
7619d007619Sjeremylt   @param[out] dim  Variable to store dimension of basis
7629d007619Sjeremylt 
7639d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
7649d007619Sjeremylt 
7659d007619Sjeremylt   @ref Backend
7669d007619Sjeremylt **/
7679d007619Sjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) {
7689d007619Sjeremylt   *dim = basis->dim;
769e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
7709d007619Sjeremylt }
7719d007619Sjeremylt 
7729d007619Sjeremylt /**
773d99fa3c5SJeremy L Thompson   @brief Get topology for given CeedBasis
774d99fa3c5SJeremy L Thompson 
775d99fa3c5SJeremy L Thompson   @param basis      CeedBasis
776d99fa3c5SJeremy L Thompson   @param[out] topo  Variable to store topology of basis
777d99fa3c5SJeremy L Thompson 
778d99fa3c5SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
779d99fa3c5SJeremy L Thompson 
780d99fa3c5SJeremy L Thompson   @ref Backend
781d99fa3c5SJeremy L Thompson **/
782d99fa3c5SJeremy L Thompson int CeedBasisGetTopology(CeedBasis basis, CeedElemTopology *topo) {
783d99fa3c5SJeremy L Thompson   *topo = basis->topo;
784e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
785d99fa3c5SJeremy L Thompson }
786d99fa3c5SJeremy L Thompson 
787d99fa3c5SJeremy L Thompson /**
7889d007619Sjeremylt   @brief Get number of components for given CeedBasis
7899d007619Sjeremylt 
7909d007619Sjeremylt   @param basis          CeedBasis
791d1d35e2fSjeremylt   @param[out] num_comp  Variable to store number of components of basis
7929d007619Sjeremylt 
7939d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
7949d007619Sjeremylt 
7959d007619Sjeremylt   @ref Backend
7969d007619Sjeremylt **/
797d1d35e2fSjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *num_comp) {
798d1d35e2fSjeremylt   *num_comp = basis->num_comp;
799e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8009d007619Sjeremylt }
8019d007619Sjeremylt 
8029d007619Sjeremylt /**
8039d007619Sjeremylt   @brief Get total number of nodes (in dim dimensions) of a CeedBasis
8049d007619Sjeremylt 
8059d007619Sjeremylt   @param basis   CeedBasis
8069d007619Sjeremylt   @param[out] P  Variable to store number of nodes
8079d007619Sjeremylt 
8089d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
8099d007619Sjeremylt 
8109d007619Sjeremylt   @ref Utility
8119d007619Sjeremylt **/
8129d007619Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) {
8139d007619Sjeremylt   *P = basis->P;
814e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8159d007619Sjeremylt }
8169d007619Sjeremylt 
8179d007619Sjeremylt /**
8189d007619Sjeremylt   @brief Get total number of nodes (in 1 dimension) of a CeedBasis
8199d007619Sjeremylt 
8209d007619Sjeremylt   @param basis     CeedBasis
821d1d35e2fSjeremylt   @param[out] P_1d  Variable to store number of nodes
8229d007619Sjeremylt 
8239d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
8249d007619Sjeremylt 
8259d007619Sjeremylt   @ref Backend
8269d007619Sjeremylt **/
827d1d35e2fSjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P_1d) {
828d1d35e2fSjeremylt   if (!basis->tensor_basis)
8299d007619Sjeremylt     // LCOV_EXCL_START
830e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR,
831d1d35e2fSjeremylt                      "Cannot supply P_1d for non-tensor basis");
8329d007619Sjeremylt   // LCOV_EXCL_STOP
8339d007619Sjeremylt 
834d1d35e2fSjeremylt   *P_1d = basis->P_1d;
835e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8369d007619Sjeremylt }
8379d007619Sjeremylt 
8389d007619Sjeremylt /**
8399d007619Sjeremylt   @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis
8409d007619Sjeremylt 
8419d007619Sjeremylt   @param basis   CeedBasis
8429d007619Sjeremylt   @param[out] Q  Variable to store number of quadrature points
8439d007619Sjeremylt 
8449d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
8459d007619Sjeremylt 
8469d007619Sjeremylt   @ref Utility
8479d007619Sjeremylt **/
8489d007619Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) {
8499d007619Sjeremylt   *Q = basis->Q;
850e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8519d007619Sjeremylt }
8529d007619Sjeremylt 
8539d007619Sjeremylt /**
8549d007619Sjeremylt   @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis
8559d007619Sjeremylt 
8569d007619Sjeremylt   @param basis     CeedBasis
857d1d35e2fSjeremylt   @param[out] Q_1d  Variable to store number of quadrature points
8589d007619Sjeremylt 
8599d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
8609d007619Sjeremylt 
8619d007619Sjeremylt   @ref Backend
8629d007619Sjeremylt **/
863d1d35e2fSjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q_1d) {
864d1d35e2fSjeremylt   if (!basis->tensor_basis)
8659d007619Sjeremylt     // LCOV_EXCL_START
866e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR,
867d1d35e2fSjeremylt                      "Cannot supply Q_1d for non-tensor basis");
8689d007619Sjeremylt   // LCOV_EXCL_STOP
8699d007619Sjeremylt 
870d1d35e2fSjeremylt   *Q_1d = basis->Q_1d;
871e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8729d007619Sjeremylt }
8739d007619Sjeremylt 
8749d007619Sjeremylt /**
8759d007619Sjeremylt   @brief Get reference coordinates of quadrature points (in dim dimensions)
8769d007619Sjeremylt          of a CeedBasis
8779d007619Sjeremylt 
8789d007619Sjeremylt   @param basis       CeedBasis
879d1d35e2fSjeremylt   @param[out] q_ref  Variable to store reference coordinates of quadrature points
8809d007619Sjeremylt 
8819d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
8829d007619Sjeremylt 
8839d007619Sjeremylt   @ref Backend
8849d007619Sjeremylt **/
885d1d35e2fSjeremylt int CeedBasisGetQRef(CeedBasis basis, const CeedScalar **q_ref) {
886d1d35e2fSjeremylt   *q_ref = basis->q_ref_1d;
887e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8889d007619Sjeremylt }
8899d007619Sjeremylt 
8909d007619Sjeremylt /**
8919d007619Sjeremylt   @brief Get quadrature weights of quadrature points (in dim dimensions)
8929d007619Sjeremylt          of a CeedBasis
8939d007619Sjeremylt 
8949d007619Sjeremylt   @param basis         CeedBasis
895d1d35e2fSjeremylt   @param[out] q_weight  Variable to store quadrature weights
8969d007619Sjeremylt 
8979d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
8989d007619Sjeremylt 
8999d007619Sjeremylt   @ref Backend
9009d007619Sjeremylt **/
901d1d35e2fSjeremylt int CeedBasisGetQWeights(CeedBasis basis, const CeedScalar **q_weight) {
902d1d35e2fSjeremylt   *q_weight = basis->q_weight_1d;
903e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9049d007619Sjeremylt }
9059d007619Sjeremylt 
9069d007619Sjeremylt /**
9079d007619Sjeremylt   @brief Get interpolation matrix of a CeedBasis
9089d007619Sjeremylt 
9099d007619Sjeremylt   @param basis        CeedBasis
9109d007619Sjeremylt   @param[out] interp  Variable to store interpolation matrix
9119d007619Sjeremylt 
9129d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
9139d007619Sjeremylt 
9149d007619Sjeremylt   @ref Backend
9159d007619Sjeremylt **/
9166c58de82SJeremy L Thompson int CeedBasisGetInterp(CeedBasis basis, const CeedScalar **interp) {
917d1d35e2fSjeremylt   if (!basis->interp && basis->tensor_basis) {
9189d007619Sjeremylt     // Allocate
9199d007619Sjeremylt     int ierr;
9209d007619Sjeremylt     ierr = CeedMalloc(basis->Q*basis->P, &basis->interp); CeedChk(ierr);
9219d007619Sjeremylt 
9229d007619Sjeremylt     // Initialize
9239d007619Sjeremylt     for (CeedInt i=0; i<basis->Q*basis->P; i++)
9249d007619Sjeremylt       basis->interp[i] = 1.0;
9259d007619Sjeremylt 
9269d007619Sjeremylt     // Calculate
9279d007619Sjeremylt     for (CeedInt d=0; d<basis->dim; d++)
9289d007619Sjeremylt       for (CeedInt qpt=0; qpt<basis->Q; qpt++)
9299d007619Sjeremylt         for (CeedInt node=0; node<basis->P; node++) {
930d1d35e2fSjeremylt           CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d;
931d1d35e2fSjeremylt           CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d;
932d1d35e2fSjeremylt           basis->interp[qpt*(basis->P)+node] *= basis->interp_1d[q*basis->P_1d+p];
9339d007619Sjeremylt         }
9349d007619Sjeremylt   }
9359d007619Sjeremylt   *interp = basis->interp;
936e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9379d007619Sjeremylt }
9389d007619Sjeremylt 
9399d007619Sjeremylt /**
9409d007619Sjeremylt   @brief Get 1D interpolation matrix of a tensor product CeedBasis
9419d007619Sjeremylt 
9429d007619Sjeremylt   @param basis           CeedBasis
943d1d35e2fSjeremylt   @param[out] interp_1d  Variable to store interpolation matrix
9449d007619Sjeremylt 
9459d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
9469d007619Sjeremylt 
9479d007619Sjeremylt   @ref Backend
9489d007619Sjeremylt **/
949d1d35e2fSjeremylt int CeedBasisGetInterp1D(CeedBasis basis, const CeedScalar **interp_1d) {
950d1d35e2fSjeremylt   if (!basis->tensor_basis)
9519d007619Sjeremylt     // LCOV_EXCL_START
952e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR,
953e15f9bd0SJeremy L Thompson                      "CeedBasis is not a tensor product basis.");
9549d007619Sjeremylt   // LCOV_EXCL_STOP
9559d007619Sjeremylt 
956d1d35e2fSjeremylt   *interp_1d = basis->interp_1d;
957e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9589d007619Sjeremylt }
9599d007619Sjeremylt 
9609d007619Sjeremylt /**
9619d007619Sjeremylt   @brief Get gradient matrix of a CeedBasis
9629d007619Sjeremylt 
9639d007619Sjeremylt   @param basis      CeedBasis
9649d007619Sjeremylt   @param[out] grad  Variable to store gradient matrix
9659d007619Sjeremylt 
9669d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
9679d007619Sjeremylt 
9689d007619Sjeremylt   @ref Backend
9699d007619Sjeremylt **/
9706c58de82SJeremy L Thompson int CeedBasisGetGrad(CeedBasis basis, const CeedScalar **grad) {
971d1d35e2fSjeremylt   if (!basis->grad && basis->tensor_basis) {
9729d007619Sjeremylt     // Allocate
9739d007619Sjeremylt     int ierr;
9749d007619Sjeremylt     ierr = CeedMalloc(basis->dim*basis->Q*basis->P, &basis->grad);
9759d007619Sjeremylt     CeedChk(ierr);
9769d007619Sjeremylt 
9779d007619Sjeremylt     // Initialize
9789d007619Sjeremylt     for (CeedInt i=0; i<basis->dim*basis->Q*basis->P; i++)
9799d007619Sjeremylt       basis->grad[i] = 1.0;
9809d007619Sjeremylt 
9819d007619Sjeremylt     // Calculate
9829d007619Sjeremylt     for (CeedInt d=0; d<basis->dim; d++)
9839d007619Sjeremylt       for (CeedInt i=0; i<basis->dim; i++)
9849d007619Sjeremylt         for (CeedInt qpt=0; qpt<basis->Q; qpt++)
9859d007619Sjeremylt           for (CeedInt node=0; node<basis->P; node++) {
986d1d35e2fSjeremylt             CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d;
987d1d35e2fSjeremylt             CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d;
9889d007619Sjeremylt             if (i == d)
9899d007619Sjeremylt               basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *=
990d1d35e2fSjeremylt                 basis->grad_1d[q*basis->P_1d+p];
9919d007619Sjeremylt             else
9929d007619Sjeremylt               basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *=
993d1d35e2fSjeremylt                 basis->interp_1d[q*basis->P_1d+p];
9949d007619Sjeremylt           }
9959d007619Sjeremylt   }
9969d007619Sjeremylt   *grad = basis->grad;
997e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9989d007619Sjeremylt }
9999d007619Sjeremylt 
10009d007619Sjeremylt /**
10019d007619Sjeremylt   @brief Get 1D gradient matrix of a tensor product CeedBasis
10029d007619Sjeremylt 
10039d007619Sjeremylt   @param basis         CeedBasis
1004d1d35e2fSjeremylt   @param[out] grad_1d  Variable to store gradient matrix
10059d007619Sjeremylt 
10069d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
10079d007619Sjeremylt 
10089d007619Sjeremylt   @ref Backend
10099d007619Sjeremylt **/
1010d1d35e2fSjeremylt int CeedBasisGetGrad1D(CeedBasis basis, const CeedScalar **grad_1d) {
1011d1d35e2fSjeremylt   if (!basis->tensor_basis)
10129d007619Sjeremylt     // LCOV_EXCL_START
1013e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR,
1014e15f9bd0SJeremy L Thompson                      "CeedBasis is not a tensor product basis.");
10159d007619Sjeremylt   // LCOV_EXCL_STOP
10169d007619Sjeremylt 
1017d1d35e2fSjeremylt   *grad_1d = basis->grad_1d;
1018e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
10199d007619Sjeremylt }
10209d007619Sjeremylt 
10219d007619Sjeremylt /**
10227a982d89SJeremy L. Thompson   @brief Destroy a CeedBasis
10237a982d89SJeremy L. Thompson 
10247a982d89SJeremy L. Thompson   @param basis CeedBasis to destroy
10257a982d89SJeremy L. Thompson 
10267a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
10277a982d89SJeremy L. Thompson 
10287a982d89SJeremy L. Thompson   @ref User
10297a982d89SJeremy L. Thompson **/
10307a982d89SJeremy L. Thompson int CeedBasisDestroy(CeedBasis *basis) {
10317a982d89SJeremy L. Thompson   int ierr;
10327a982d89SJeremy L. Thompson 
1033d1d35e2fSjeremylt   if (!*basis || --(*basis)->ref_count > 0) return CEED_ERROR_SUCCESS;
10347a982d89SJeremy L. Thompson   if ((*basis)->Destroy) {
10357a982d89SJeremy L. Thompson     ierr = (*basis)->Destroy(*basis); CeedChk(ierr);
10367a982d89SJeremy L. Thompson   }
10377a982d89SJeremy L. Thompson   ierr = CeedFree(&(*basis)->interp); CeedChk(ierr);
1038d1d35e2fSjeremylt   ierr = CeedFree(&(*basis)->interp_1d); CeedChk(ierr);
10397a982d89SJeremy L. Thompson   ierr = CeedFree(&(*basis)->grad); CeedChk(ierr);
1040d1d35e2fSjeremylt   ierr = CeedFree(&(*basis)->grad_1d); CeedChk(ierr);
1041d1d35e2fSjeremylt   ierr = CeedFree(&(*basis)->q_ref_1d); CeedChk(ierr);
1042d1d35e2fSjeremylt   ierr = CeedFree(&(*basis)->q_weight_1d); CeedChk(ierr);
10437a982d89SJeremy L. Thompson   ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr);
10447a982d89SJeremy L. Thompson   ierr = CeedFree(basis); CeedChk(ierr);
1045e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
10467a982d89SJeremy L. Thompson }
10477a982d89SJeremy L. Thompson 
10487a982d89SJeremy L. Thompson /**
1049b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre quadrature
1050b11c1e72Sjeremylt 
1051b11c1e72Sjeremylt   @param Q               Number of quadrature points (integrates polynomials of
1052b11c1e72Sjeremylt                            degree 2*Q-1 exactly)
1053d1d35e2fSjeremylt   @param[out] q_ref_1d     Array of length Q to hold the abscissa on [-1, 1]
1054d1d35e2fSjeremylt   @param[out] q_weight_1d  Array of length Q to hold the weights
1055b11c1e72Sjeremylt 
1056b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1057dfdf5a53Sjeremylt 
1058dfdf5a53Sjeremylt   @ref Utility
1059b11c1e72Sjeremylt **/
1060d1d35e2fSjeremylt int CeedGaussQuadrature(CeedInt Q, CeedScalar *q_ref_1d,
1061d1d35e2fSjeremylt                         CeedScalar *q_weight_1d) {
1062d7b241e6Sjeremylt   // Allocate
1063d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0*atan(1.0);
1064d1d35e2fSjeremylt   // Build q_ref_1d, q_weight_1d
1065d7b241e6Sjeremylt   for (int i = 0; i <= Q/2; i++) {
1066d7b241e6Sjeremylt     // Guess
1067d7b241e6Sjeremylt     xi = cos(PI*(CeedScalar)(2*i+1)/((CeedScalar)(2*Q)));
1068d7b241e6Sjeremylt     // Pn(xi)
1069d7b241e6Sjeremylt     P0 = 1.0;
1070d7b241e6Sjeremylt     P1 = xi;
1071d7b241e6Sjeremylt     P2 = 0.0;
1072d7b241e6Sjeremylt     for (int j = 2; j <= Q; j++) {
1073d7b241e6Sjeremylt       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1074d7b241e6Sjeremylt       P0 = P1;
1075d7b241e6Sjeremylt       P1 = P2;
1076d7b241e6Sjeremylt     }
1077d7b241e6Sjeremylt     // First Newton Step
1078d7b241e6Sjeremylt     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1079d7b241e6Sjeremylt     xi = xi-P2/dP2;
1080d7b241e6Sjeremylt     // Newton to convergence
10810e4d4210Sjeremylt     for (int k=0; k<100 && fabs(P2)>10*CEED_EPSILON; k++) {
1082d7b241e6Sjeremylt       P0 = 1.0;
1083d7b241e6Sjeremylt       P1 = xi;
1084d7b241e6Sjeremylt       for (int j = 2; j <= Q; j++) {
1085d7b241e6Sjeremylt         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1086d7b241e6Sjeremylt         P0 = P1;
1087d7b241e6Sjeremylt         P1 = P2;
1088d7b241e6Sjeremylt       }
1089d7b241e6Sjeremylt       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1090d7b241e6Sjeremylt       xi = xi-P2/dP2;
1091d7b241e6Sjeremylt     }
1092d7b241e6Sjeremylt     // Save xi, wi
1093d7b241e6Sjeremylt     wi = 2.0/((1.0-xi*xi)*dP2*dP2);
1094d1d35e2fSjeremylt     q_weight_1d[i] = wi;
1095d1d35e2fSjeremylt     q_weight_1d[Q-1-i] = wi;
1096d1d35e2fSjeremylt     q_ref_1d[i] = -xi;
1097d1d35e2fSjeremylt     q_ref_1d[Q-1-i]= xi;
1098d7b241e6Sjeremylt   }
1099e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1100d7b241e6Sjeremylt }
1101d7b241e6Sjeremylt 
1102b11c1e72Sjeremylt /**
1103b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre-Lobatto quadrature
1104b11c1e72Sjeremylt 
1105b11c1e72Sjeremylt   @param Q               Number of quadrature points (integrates polynomials of
1106b11c1e72Sjeremylt                            degree 2*Q-3 exactly)
1107d1d35e2fSjeremylt   @param[out] q_ref_1d     Array of length Q to hold the abscissa on [-1, 1]
1108d1d35e2fSjeremylt   @param[out] q_weight_1d  Array of length Q to hold the weights
1109b11c1e72Sjeremylt 
1110b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1111dfdf5a53Sjeremylt 
1112dfdf5a53Sjeremylt   @ref Utility
1113b11c1e72Sjeremylt **/
1114d1d35e2fSjeremylt int CeedLobattoQuadrature(CeedInt Q, CeedScalar *q_ref_1d,
1115d1d35e2fSjeremylt                           CeedScalar *q_weight_1d) {
1116d7b241e6Sjeremylt   // Allocate
1117d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0*atan(1.0);
1118d1d35e2fSjeremylt   // Build q_ref_1d, q_weight_1d
1119d7b241e6Sjeremylt   // Set endpoints
112030a100c3SJed Brown   if (Q < 2)
1121b0d62198Sjeremylt     // LCOV_EXCL_START
1122e15f9bd0SJeremy L Thompson     return CeedError(NULL, CEED_ERROR_DIMENSION,
11237ed52d01Sjeremylt                      "Cannot create Lobatto quadrature with Q=%d < 2 points", Q);
1124b0d62198Sjeremylt   // LCOV_EXCL_STOP
1125d7b241e6Sjeremylt   wi = 2.0/((CeedScalar)(Q*(Q-1)));
1126d1d35e2fSjeremylt   if (q_weight_1d) {
1127d1d35e2fSjeremylt     q_weight_1d[0] = wi;
1128d1d35e2fSjeremylt     q_weight_1d[Q-1] = wi;
1129d7b241e6Sjeremylt   }
1130d1d35e2fSjeremylt   q_ref_1d[0] = -1.0;
1131d1d35e2fSjeremylt   q_ref_1d[Q-1] = 1.0;
1132d7b241e6Sjeremylt   // Interior
1133d7b241e6Sjeremylt   for (int i = 1; i <= (Q-1)/2; i++) {
1134d7b241e6Sjeremylt     // Guess
1135d7b241e6Sjeremylt     xi = cos(PI*(CeedScalar)(i)/(CeedScalar)(Q-1));
1136d7b241e6Sjeremylt     // Pn(xi)
1137d7b241e6Sjeremylt     P0 = 1.0;
1138d7b241e6Sjeremylt     P1 = xi;
1139d7b241e6Sjeremylt     P2 = 0.0;
1140d7b241e6Sjeremylt     for (int j = 2; j < Q; j++) {
1141d7b241e6Sjeremylt       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1142d7b241e6Sjeremylt       P0 = P1;
1143d7b241e6Sjeremylt       P1 = P2;
1144d7b241e6Sjeremylt     }
1145d7b241e6Sjeremylt     // First Newton step
1146d7b241e6Sjeremylt     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1147d7b241e6Sjeremylt     d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
1148d7b241e6Sjeremylt     xi = xi-dP2/d2P2;
1149d7b241e6Sjeremylt     // Newton to convergence
11500e4d4210Sjeremylt     for (int k=0; k<100 && fabs(dP2)>10*CEED_EPSILON; k++) {
1151d7b241e6Sjeremylt       P0 = 1.0;
1152d7b241e6Sjeremylt       P1 = xi;
1153d7b241e6Sjeremylt       for (int j = 2; j < Q; j++) {
1154d7b241e6Sjeremylt         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1155d7b241e6Sjeremylt         P0 = P1;
1156d7b241e6Sjeremylt         P1 = P2;
1157d7b241e6Sjeremylt       }
1158d7b241e6Sjeremylt       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1159d7b241e6Sjeremylt       d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
1160d7b241e6Sjeremylt       xi = xi-dP2/d2P2;
1161d7b241e6Sjeremylt     }
1162d7b241e6Sjeremylt     // Save xi, wi
1163d7b241e6Sjeremylt     wi = 2.0/(((CeedScalar)(Q*(Q-1)))*P2*P2);
1164d1d35e2fSjeremylt     if (q_weight_1d) {
1165d1d35e2fSjeremylt       q_weight_1d[i] = wi;
1166d1d35e2fSjeremylt       q_weight_1d[Q-1-i] = wi;
1167d7b241e6Sjeremylt     }
1168d1d35e2fSjeremylt     q_ref_1d[i] = -xi;
1169d1d35e2fSjeremylt     q_ref_1d[Q-1-i]= xi;
1170d7b241e6Sjeremylt   }
1171e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1172d7b241e6Sjeremylt }
1173d7b241e6Sjeremylt 
1174dfdf5a53Sjeremylt /**
117595bb1877Svaleriabarra   @brief Return QR Factorization of a matrix
1176b11c1e72Sjeremylt 
117777645d7bSjeremylt   @param ceed         A Ceed context for error handling
117852bfb9bbSJeremy L Thompson   @param[in,out] mat  Row-major matrix to be factorized in place
117952bfb9bbSJeremy L Thompson   @param[in,out] tau  Vector of length m of scaling factors
1180b11c1e72Sjeremylt   @param m            Number of rows
1181b11c1e72Sjeremylt   @param n            Number of columns
1182b11c1e72Sjeremylt 
1183b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1184dfdf5a53Sjeremylt 
1185dfdf5a53Sjeremylt   @ref Utility
1186b11c1e72Sjeremylt **/
1187a7bd39daSjeremylt int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau,
1188d7b241e6Sjeremylt                         CeedInt m, CeedInt n) {
1189d7b241e6Sjeremylt   CeedScalar v[m];
1190d7b241e6Sjeremylt 
1191a7bd39daSjeremylt   // Check m >= n
1192a7bd39daSjeremylt   if (n > m)
1193c042f62fSJeremy L Thompson     // LCOV_EXCL_START
1194e15f9bd0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_UNSUPPORTED,
1195e15f9bd0SJeremy L Thompson                      "Cannot compute QR factorization with n > m");
1196c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
1197a7bd39daSjeremylt 
119852bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++) {
1199d7b241e6Sjeremylt     // Calculate Householder vector, magnitude
1200d7b241e6Sjeremylt     CeedScalar sigma = 0.0;
1201d7b241e6Sjeremylt     v[i] = mat[i+n*i];
120252bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<m; j++) {
1203d7b241e6Sjeremylt       v[j] = mat[i+n*j];
1204d7b241e6Sjeremylt       sigma += v[j] * v[j];
1205d7b241e6Sjeremylt     }
1206d7b241e6Sjeremylt     CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m]
1207d7b241e6Sjeremylt     CeedScalar Rii = -copysign(norm, v[i]);
1208d7b241e6Sjeremylt     v[i] -= Rii;
1209d7b241e6Sjeremylt     // norm of v[i:m] after modification above and scaling below
1210d7b241e6Sjeremylt     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
1211d7b241e6Sjeremylt     //   tau = 2 / (norm*norm)
1212d7b241e6Sjeremylt     tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma);
1213fb551037Sjeremylt 
12141d102b48SJeremy L Thompson     for (CeedInt j=i+1; j<m; j++)
12151d102b48SJeremy L Thompson       v[j] /= v[i];
1216d7b241e6Sjeremylt 
1217d7b241e6Sjeremylt     // Apply Householder reflector to lower right panel
1218d7b241e6Sjeremylt     CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1);
1219d7b241e6Sjeremylt     // Save v
1220d7b241e6Sjeremylt     mat[i+n*i] = Rii;
12211d102b48SJeremy L Thompson     for (CeedInt j=i+1; j<m; j++)
1222d7b241e6Sjeremylt       mat[i+n*j] = v[j];
1223d7b241e6Sjeremylt   }
1224e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1225d7b241e6Sjeremylt }
1226d7b241e6Sjeremylt 
1227b11c1e72Sjeremylt /**
122852bfb9bbSJeremy L Thompson   @brief Return symmetric Schur decomposition of the symmetric matrix mat via
122952bfb9bbSJeremy L Thompson            symmetric QR factorization
123052bfb9bbSJeremy L Thompson 
123177645d7bSjeremylt   @param ceed         A Ceed context for error handling
123252bfb9bbSJeremy L Thompson   @param[in,out] mat  Row-major matrix to be factorized in place
1233460bf743SValeria Barra   @param[out] lambda  Vector of length n of eigenvalues
123452bfb9bbSJeremy L Thompson   @param n            Number of rows/columns
123552bfb9bbSJeremy L Thompson 
123652bfb9bbSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
123752bfb9bbSJeremy L Thompson 
123852bfb9bbSJeremy L Thompson   @ref Utility
123952bfb9bbSJeremy L Thompson **/
124052bfb9bbSJeremy L Thompson int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat,
124152bfb9bbSJeremy L Thompson                                     CeedScalar *lambda, CeedInt n) {
124252bfb9bbSJeremy L Thompson   // Check bounds for clang-tidy
124352bfb9bbSJeremy L Thompson   if (n<2)
1244c042f62fSJeremy L Thompson     // LCOV_EXCL_START
1245e15f9bd0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_UNSUPPORTED,
1246c042f62fSJeremy L Thompson                      "Cannot compute symmetric Schur decomposition of scalars");
1247c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
124852bfb9bbSJeremy L Thompson 
124952bfb9bbSJeremy L Thompson   CeedScalar v[n-1], tau[n-1], matT[n*n];
125052bfb9bbSJeremy L Thompson 
125152bfb9bbSJeremy L Thompson   // Copy mat to matT and set mat to I
125252bfb9bbSJeremy L Thompson   memcpy(matT, mat, n*n*sizeof(mat[0]));
125352bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
125452bfb9bbSJeremy L Thompson     for (CeedInt j=0; j<n; j++)
125552bfb9bbSJeremy L Thompson       mat[j+n*i] = (i==j) ? 1 : 0;
125652bfb9bbSJeremy L Thompson 
125752bfb9bbSJeremy L Thompson   // Reduce to tridiagonal
125852bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n-1; i++) {
125952bfb9bbSJeremy L Thompson     // Calculate Householder vector, magnitude
126052bfb9bbSJeremy L Thompson     CeedScalar sigma = 0.0;
126152bfb9bbSJeremy L Thompson     v[i] = matT[i+n*(i+1)];
126252bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) {
126352bfb9bbSJeremy L Thompson       v[j] = matT[i+n*(j+1)];
126452bfb9bbSJeremy L Thompson       sigma += v[j] * v[j];
126552bfb9bbSJeremy L Thompson     }
126652bfb9bbSJeremy L Thompson     CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:n-1]
126752bfb9bbSJeremy L Thompson     CeedScalar Rii = -copysign(norm, v[i]);
126852bfb9bbSJeremy L Thompson     v[i] -= Rii;
126952bfb9bbSJeremy L Thompson     // norm of v[i:m] after modification above and scaling below
127052bfb9bbSJeremy L Thompson     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
127152bfb9bbSJeremy L Thompson     //   tau = 2 / (norm*norm)
12720e4d4210Sjeremylt     if (sigma > 10*CEED_EPSILON)
127352bfb9bbSJeremy L Thompson       tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma);
1274fb551037Sjeremylt     else
1275fb551037Sjeremylt       tau[i] = 0;
1276fb551037Sjeremylt 
1277fb551037Sjeremylt     for (CeedInt j=i+1; j<n-1; j++)
1278fb551037Sjeremylt       v[j] /= v[i];
127952bfb9bbSJeremy L Thompson 
128052bfb9bbSJeremy L Thompson     // Update sub and super diagonal
128152bfb9bbSJeremy L Thompson     matT[i+n*(i+1)] = Rii;
128252bfb9bbSJeremy L Thompson     matT[(i+1)+n*i] = Rii;
128352bfb9bbSJeremy L Thompson     for (CeedInt j=i+2; j<n; j++) {
128452bfb9bbSJeremy L Thompson       matT[i+n*j] = 0; matT[j+n*i] = 0;
128552bfb9bbSJeremy L Thompson     }
128652bfb9bbSJeremy L Thompson     // Apply symmetric Householder reflector to lower right panel
128752bfb9bbSJeremy L Thompson     CeedHouseholderReflect(&matT[(i+1)+n*(i+1)], &v[i], tau[i],
128852bfb9bbSJeremy L Thompson                            n-(i+1), n-(i+1), n, 1);
128952bfb9bbSJeremy L Thompson     CeedHouseholderReflect(&matT[(i+1)+n*(i+1)], &v[i], tau[i],
129052bfb9bbSJeremy L Thompson                            n-(i+1), n-(i+1), 1, n);
129152bfb9bbSJeremy L Thompson     // Save v
129252bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) {
129352bfb9bbSJeremy L Thompson       matT[i+n*(j+1)] = v[j];
129452bfb9bbSJeremy L Thompson     }
129552bfb9bbSJeremy L Thompson   }
129652bfb9bbSJeremy L Thompson   // Backwards accumulation of Q
129752bfb9bbSJeremy L Thompson   for (CeedInt i=n-2; i>=0; i--) {
129852bfb9bbSJeremy L Thompson     v[i] = 1;
129952bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) {
130052bfb9bbSJeremy L Thompson       v[j] = matT[i+n*(j+1)];
130152bfb9bbSJeremy L Thompson       matT[i+n*(j+1)] = 0;
130252bfb9bbSJeremy L Thompson     }
130352bfb9bbSJeremy L Thompson     CeedHouseholderReflect(&mat[(i+1)+n*(i+1)], &v[i], tau[i],
130452bfb9bbSJeremy L Thompson                            n-(i+1), n-(i+1), n, 1);
130552bfb9bbSJeremy L Thompson   }
130652bfb9bbSJeremy L Thompson 
130752bfb9bbSJeremy L Thompson   // Reduce sub and super diagonal
130852bfb9bbSJeremy L Thompson   CeedInt p = 0, q = 0, itr = 0, maxitr = n*n*n;
13090e4d4210Sjeremylt   CeedScalar tol = 10*CEED_EPSILON;
131052bfb9bbSJeremy L Thompson 
131152bfb9bbSJeremy L Thompson   while (q < n && itr < maxitr) {
131252bfb9bbSJeremy L Thompson     // Update p, q, size of reduced portions of diagonal
131352bfb9bbSJeremy L Thompson     p = 0; q = 0;
131452bfb9bbSJeremy L Thompson     for (CeedInt i=n-2; i>=0; i--) {
131552bfb9bbSJeremy L Thompson       if (fabs(matT[i+n*(i+1)]) < tol)
131652bfb9bbSJeremy L Thompson         q += 1;
131752bfb9bbSJeremy L Thompson       else
131852bfb9bbSJeremy L Thompson         break;
131952bfb9bbSJeremy L Thompson     }
132052bfb9bbSJeremy L Thompson     for (CeedInt i=0; i<n-1-q; i++) {
132152bfb9bbSJeremy L Thompson       if (fabs(matT[i+n*(i+1)]) < tol)
132252bfb9bbSJeremy L Thompson         p += 1;
132352bfb9bbSJeremy L Thompson       else
132452bfb9bbSJeremy L Thompson         break;
132552bfb9bbSJeremy L Thompson     }
132652bfb9bbSJeremy L Thompson     if (q == n-1) break; // Finished reducing
132752bfb9bbSJeremy L Thompson 
132852bfb9bbSJeremy L Thompson     // Reduce tridiagonal portion
132952bfb9bbSJeremy L Thompson     CeedScalar tnn = matT[(n-1-q)+n*(n-1-q)],
133052bfb9bbSJeremy L Thompson                tnnm1 = matT[(n-2-q)+n*(n-1-q)];
133152bfb9bbSJeremy L Thompson     CeedScalar d = (matT[(n-2-q)+n*(n-2-q)] - tnn)/2;
133252bfb9bbSJeremy L Thompson     CeedScalar mu = tnn - tnnm1*tnnm1 /
133352bfb9bbSJeremy L Thompson                     (d + copysign(sqrt(d*d + tnnm1*tnnm1), d));
133452bfb9bbSJeremy L Thompson     CeedScalar x = matT[p+n*p] - mu;
133552bfb9bbSJeremy L Thompson     CeedScalar z = matT[p+n*(p+1)];
133652bfb9bbSJeremy L Thompson     for (CeedInt k=p; k<n-1-q; k++) {
133752bfb9bbSJeremy L Thompson       // Compute Givens rotation
133852bfb9bbSJeremy L Thompson       CeedScalar c = 1, s = 0;
133952bfb9bbSJeremy L Thompson       if (fabs(z) > tol) {
134052bfb9bbSJeremy L Thompson         if (fabs(z) > fabs(x)) {
134152bfb9bbSJeremy L Thompson           CeedScalar tau = -x/z;
134252bfb9bbSJeremy L Thompson           s = 1/sqrt(1+tau*tau), c = s*tau;
134352bfb9bbSJeremy L Thompson         } else {
134452bfb9bbSJeremy L Thompson           CeedScalar tau = -z/x;
134552bfb9bbSJeremy L Thompson           c = 1/sqrt(1+tau*tau), s = c*tau;
134652bfb9bbSJeremy L Thompson         }
134752bfb9bbSJeremy L Thompson       }
134852bfb9bbSJeremy L Thompson 
134952bfb9bbSJeremy L Thompson       // Apply Givens rotation to T
135052bfb9bbSJeremy L Thompson       CeedGivensRotation(matT, c, s, CEED_NOTRANSPOSE, k, k+1, n, n);
135152bfb9bbSJeremy L Thompson       CeedGivensRotation(matT, c, s, CEED_TRANSPOSE, k, k+1, n, n);
135252bfb9bbSJeremy L Thompson 
135352bfb9bbSJeremy L Thompson       // Apply Givens rotation to Q
135452bfb9bbSJeremy L Thompson       CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k+1, n, n);
135552bfb9bbSJeremy L Thompson 
135652bfb9bbSJeremy L Thompson       // Update x, z
135752bfb9bbSJeremy L Thompson       if (k < n-q-2) {
135852bfb9bbSJeremy L Thompson         x = matT[k+n*(k+1)];
135952bfb9bbSJeremy L Thompson         z = matT[k+n*(k+2)];
136052bfb9bbSJeremy L Thompson       }
136152bfb9bbSJeremy L Thompson     }
136252bfb9bbSJeremy L Thompson     itr++;
136352bfb9bbSJeremy L Thompson   }
136452bfb9bbSJeremy L Thompson   // Save eigenvalues
136552bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
136652bfb9bbSJeremy L Thompson     lambda[i] = matT[i+n*i];
136752bfb9bbSJeremy L Thompson 
136852bfb9bbSJeremy L Thompson   // Check convergence
136952bfb9bbSJeremy L Thompson   if (itr == maxitr && q < n-1)
1370c042f62fSJeremy L Thompson     // LCOV_EXCL_START
1371e15f9bd0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_MINOR,
1372e15f9bd0SJeremy L Thompson                      "Symmetric QR failed to converge");
1373c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
1374e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
137552bfb9bbSJeremy L Thompson }
137652bfb9bbSJeremy L Thompson 
137752bfb9bbSJeremy L Thompson /**
137852bfb9bbSJeremy L Thompson   @brief Return Simultaneous Diagonalization of two matrices. This solves the
137952bfb9bbSJeremy L Thompson            generalized eigenvalue problem A x = lambda B x, where A and B
138052bfb9bbSJeremy L Thompson            are symmetric and B is positive definite. We generate the matrix X
138152bfb9bbSJeremy L Thompson            and vector Lambda such that X^T A X = Lambda and X^T B X = I. This
138252bfb9bbSJeremy L Thompson            is equivalent to the LAPACK routine 'sygv' with TYPE = 1.
138352bfb9bbSJeremy L Thompson 
138477645d7bSjeremylt   @param ceed         A Ceed context for error handling
1385d1d35e2fSjeremylt   @param[in] mat_A     Row-major matrix to be factorized with eigenvalues
1386d1d35e2fSjeremylt   @param[in] mat_B     Row-major matrix to be factorized to identity
138752bfb9bbSJeremy L Thompson   @param[out] x       Row-major orthogonal matrix
1388460bf743SValeria Barra   @param[out] lambda  Vector of length n of generalized eigenvalues
138952bfb9bbSJeremy L Thompson   @param n            Number of rows/columns
139052bfb9bbSJeremy L Thompson 
139152bfb9bbSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
139252bfb9bbSJeremy L Thompson 
139352bfb9bbSJeremy L Thompson   @ref Utility
139452bfb9bbSJeremy L Thompson **/
1395d1d35e2fSjeremylt int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *mat_A,
1396d1d35e2fSjeremylt                                     CeedScalar *mat_B, CeedScalar *x,
139752bfb9bbSJeremy L Thompson                                     CeedScalar *lambda, CeedInt n) {
139852bfb9bbSJeremy L Thompson   int ierr;
1399*78464608Sjeremylt   CeedScalar *mat_C, *mat_G, *vecD;
1400*78464608Sjeremylt   ierr = CeedCalloc(n*n, &mat_C); CeedChk(ierr);
1401*78464608Sjeremylt   ierr = CeedCalloc(n*n, &mat_G); CeedChk(ierr);
1402*78464608Sjeremylt   ierr = CeedCalloc(n, &vecD); CeedChk(ierr);
140352bfb9bbSJeremy L Thompson 
140452bfb9bbSJeremy L Thompson   // Compute B = G D G^T
1405*78464608Sjeremylt   memcpy(mat_G, mat_B, n*n*sizeof(mat_B[0]));
1406*78464608Sjeremylt   ierr = CeedSymmetricSchurDecomposition(ceed, mat_G, vecD, n); CeedChk(ierr);
1407fb551037Sjeremylt   for (CeedInt i=0; i<n; i++)
1408fb551037Sjeremylt     vecD[i] = sqrt(vecD[i]);
140952bfb9bbSJeremy L Thompson 
1410fb551037Sjeremylt   // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T
1411fb551037Sjeremylt   //           = D^-1/2 G^T A G D^-1/2
141252bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
141352bfb9bbSJeremy L Thompson     for (CeedInt j=0; j<n; j++)
1414*78464608Sjeremylt       mat_C[j+i*n] = mat_G[i+j*n] / vecD[i];
1415d1d35e2fSjeremylt   ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)mat_C,
1416d1d35e2fSjeremylt                             (const CeedScalar *)mat_A, x, n, n, n);
14179289e5bfSjeremylt   CeedChk(ierr);
141852bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
141952bfb9bbSJeremy L Thompson     for (CeedInt j=0; j<n; j++)
1420*78464608Sjeremylt       mat_G[j+i*n] = mat_G[j+i*n] / vecD[j];
14219289e5bfSjeremylt   ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)x,
1422*78464608Sjeremylt                             (const CeedScalar *)mat_G, mat_C, n, n, n);
14239289e5bfSjeremylt   CeedChk(ierr);
142452bfb9bbSJeremy L Thompson 
142552bfb9bbSJeremy L Thompson   // Compute Q^T C Q = lambda
1426d1d35e2fSjeremylt   ierr = CeedSymmetricSchurDecomposition(ceed, mat_C, lambda, n); CeedChk(ierr);
142752bfb9bbSJeremy L Thompson 
1428fb551037Sjeremylt   // Set x = (G D^1/2)^-T Q
1429fb551037Sjeremylt   //       = G D^-1/2 Q
1430*78464608Sjeremylt   ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)mat_G,
1431d1d35e2fSjeremylt                             (const CeedScalar *)mat_C, x, n, n, n);
14329289e5bfSjeremylt   CeedChk(ierr);
1433*78464608Sjeremylt 
1434*78464608Sjeremylt   // Cleanup
1435*78464608Sjeremylt   ierr = CeedFree(&mat_C); CeedChk(ierr);
1436*78464608Sjeremylt   ierr = CeedFree(&mat_G); CeedChk(ierr);
1437*78464608Sjeremylt   ierr = CeedFree(&vecD); CeedChk(ierr);
1438e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
143952bfb9bbSJeremy L Thompson }
144052bfb9bbSJeremy L Thompson 
1441d7b241e6Sjeremylt /// @}
1442