xref: /libCEED/interface/ceed-basis.c (revision ff3a0f9146a8db6126e027cdd0b3d17c2d26b6e0)
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;
10278464608Sjeremylt   CeedScalar *v;
10378464608Sjeremylt   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   }
11278464608Sjeremylt   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;
20578464608Sjeremylt   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);
20978464608Sjeremylt   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 
232673160d7Sjeremylt   // Apply Qtranspose, collo_grad = 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);
23878464608Sjeremylt   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 tensor status for given CeedBasis
2447a982d89SJeremy L. Thompson 
2457a982d89SJeremy L. Thompson   @param basis           CeedBasis
246d1d35e2fSjeremylt   @param[out] is_tensor  Variable to store tensor status
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 **/
252d1d35e2fSjeremylt int CeedBasisIsTensor(CeedBasis basis, bool *is_tensor) {
253d1d35e2fSjeremylt   *is_tensor = basis->tensor_basis;
254e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2557a982d89SJeremy L. Thompson }
2567a982d89SJeremy L. Thompson 
2577a982d89SJeremy L. Thompson /**
2587a982d89SJeremy L. Thompson   @brief Get backend data of a CeedBasis
2597a982d89SJeremy L. Thompson 
2607a982d89SJeremy L. Thompson   @param basis      CeedBasis
2617a982d89SJeremy L. Thompson   @param[out] data  Variable to store data
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 **/
267777ff853SJeremy L Thompson int CeedBasisGetData(CeedBasis basis, void *data) {
268777ff853SJeremy L Thompson   *(void **)data = basis->data;
269e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2707a982d89SJeremy L. Thompson }
2717a982d89SJeremy L. Thompson 
2727a982d89SJeremy L. Thompson /**
2737a982d89SJeremy L. Thompson   @brief Set backend data of a CeedBasis
2747a982d89SJeremy L. Thompson 
2757a982d89SJeremy L. Thompson   @param[out] basis  CeedBasis
2767a982d89SJeremy L. Thompson   @param data        Data to set
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 CeedBasisSetData(CeedBasis basis, void *data) {
283777ff853SJeremy L Thompson   basis->data = data;
284e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2857a982d89SJeremy L. Thompson }
2867a982d89SJeremy L. Thompson 
2877a982d89SJeremy L. Thompson /**
28834359f16Sjeremylt   @brief Increment the reference counter for a CeedBasis
28934359f16Sjeremylt 
29034359f16Sjeremylt   @param basis  Basis to increment the reference counter
29134359f16Sjeremylt 
29234359f16Sjeremylt   @return An error code: 0 - success, otherwise - failure
29334359f16Sjeremylt 
29434359f16Sjeremylt   @ref Backend
29534359f16Sjeremylt **/
2969560d06aSjeremylt int CeedBasisReference(CeedBasis basis) {
29734359f16Sjeremylt   basis->ref_count++;
29834359f16Sjeremylt   return CEED_ERROR_SUCCESS;
29934359f16Sjeremylt }
30034359f16Sjeremylt 
30134359f16Sjeremylt /**
3027a982d89SJeremy L. Thompson   @brief Get dimension for given CeedElemTopology
3037a982d89SJeremy L. Thompson 
3047a982d89SJeremy L. Thompson   @param topo      CeedElemTopology
3057a982d89SJeremy L. Thompson   @param[out] dim  Variable to store dimension of topology
3067a982d89SJeremy L. Thompson 
3077a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3087a982d89SJeremy L. Thompson 
3097a982d89SJeremy L. Thompson   @ref Backend
3107a982d89SJeremy L. Thompson **/
3117a982d89SJeremy L. Thompson int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) {
3127a982d89SJeremy L. Thompson   *dim = (CeedInt) topo >> 16;
313e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3147a982d89SJeremy L. Thompson }
3157a982d89SJeremy L. Thompson 
3167a982d89SJeremy L. Thompson /**
3177a982d89SJeremy L. Thompson   @brief Get CeedTensorContract of a CeedBasis
3187a982d89SJeremy L. Thompson 
3197a982d89SJeremy L. Thompson   @param basis          CeedBasis
3207a982d89SJeremy L. Thompson   @param[out] contract  Variable to store CeedTensorContract
3217a982d89SJeremy L. Thompson 
3227a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3237a982d89SJeremy L. Thompson 
3247a982d89SJeremy L. Thompson   @ref Backend
3257a982d89SJeremy L. Thompson **/
3267a982d89SJeremy L. Thompson int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) {
3277a982d89SJeremy L. Thompson   *contract = basis->contract;
328e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3297a982d89SJeremy L. Thompson }
3307a982d89SJeremy L. Thompson 
3317a982d89SJeremy L. Thompson /**
3327a982d89SJeremy L. Thompson   @brief Set CeedTensorContract of a CeedBasis
3337a982d89SJeremy L. Thompson 
3347a982d89SJeremy L. Thompson   @param[out] basis  CeedBasis
3357a982d89SJeremy L. Thompson   @param contract    CeedTensorContract to set
3367a982d89SJeremy L. Thompson 
3377a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3387a982d89SJeremy L. Thompson 
3397a982d89SJeremy L. Thompson   @ref Backend
3407a982d89SJeremy L. Thompson **/
34134359f16Sjeremylt int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract contract) {
34234359f16Sjeremylt   int ierr;
34334359f16Sjeremylt   basis->contract = contract;
3449560d06aSjeremylt   ierr = CeedTensorContractReference(contract); CeedChk(ierr);
345e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3467a982d89SJeremy L. Thompson }
3477a982d89SJeremy L. Thompson 
3487a982d89SJeremy L. Thompson /**
3497a982d89SJeremy L. Thompson   @brief Return a reference implementation of matrix multiplication C = A B.
3507a982d89SJeremy L. Thompson            Note, this is a reference implementation for CPU CeedScalar pointers
3517a982d89SJeremy L. Thompson            that is not intended for high performance.
3527a982d89SJeremy L. Thompson 
3537a982d89SJeremy L. Thompson   @param ceed        A Ceed context for error handling
354d1d35e2fSjeremylt   @param[in] mat_A   Row-major matrix A
355d1d35e2fSjeremylt   @param[in] mat_B   Row-major matrix B
356d1d35e2fSjeremylt   @param[out] mat_C  Row-major output matrix C
3577a982d89SJeremy L. Thompson   @param m           Number of rows of C
3587a982d89SJeremy L. Thompson   @param n           Number of columns of C
3597a982d89SJeremy L. Thompson   @param kk          Number of columns of A/rows of B
3607a982d89SJeremy L. Thompson 
3617a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3627a982d89SJeremy L. Thompson 
3637a982d89SJeremy L. Thompson   @ref Utility
3647a982d89SJeremy L. Thompson **/
365d1d35e2fSjeremylt int CeedMatrixMultiply(Ceed ceed, const CeedScalar *mat_A,
366d1d35e2fSjeremylt                        const CeedScalar *mat_B, CeedScalar *mat_C, CeedInt m,
3677a982d89SJeremy L. Thompson                        CeedInt n, CeedInt kk) {
3687a982d89SJeremy L. Thompson   for (CeedInt i=0; i<m; i++)
3697a982d89SJeremy L. Thompson     for (CeedInt j=0; j<n; j++) {
3707a982d89SJeremy L. Thompson       CeedScalar sum = 0;
3717a982d89SJeremy L. Thompson       for (CeedInt k=0; k<kk; k++)
372d1d35e2fSjeremylt         sum += mat_A[k+i*kk]*mat_B[j+k*n];
373d1d35e2fSjeremylt       mat_C[j+i*n] = sum;
3747a982d89SJeremy L. Thompson     }
375e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3767a982d89SJeremy L. Thompson }
3777a982d89SJeremy L. Thompson 
3787a982d89SJeremy L. Thompson /// @}
3797a982d89SJeremy L. Thompson 
3807a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
3817a982d89SJeremy L. Thompson /// CeedBasis Public API
3827a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
3837a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser
384d7b241e6Sjeremylt /// @{
385d7b241e6Sjeremylt 
386b11c1e72Sjeremylt /**
38795bb1877Svaleriabarra   @brief Create a tensor-product basis for H^1 discretizations
388b11c1e72Sjeremylt 
389b11c1e72Sjeremylt   @param ceed        A Ceed object where the CeedBasis will be created
390b11c1e72Sjeremylt   @param dim         Topological dimension
391d1d35e2fSjeremylt   @param num_comp    Number of field components (1 for scalar fields)
392d1d35e2fSjeremylt   @param P_1d        Number of nodes in one dimension
393d1d35e2fSjeremylt   @param Q_1d        Number of quadrature points in one dimension
394d1d35e2fSjeremylt   @param interp_1d   Row-major (Q_1d * P_1d) matrix expressing the values of nodal
395b11c1e72Sjeremylt                        basis functions at quadrature points
396d1d35e2fSjeremylt   @param grad_1d     Row-major (Q_1d * P_1d) matrix expressing derivatives of nodal
397b11c1e72Sjeremylt                        basis functions at quadrature points
398d1d35e2fSjeremylt   @param q_ref_1d    Array of length Q_1d holding the locations of quadrature points
399b11c1e72Sjeremylt                        on the 1D reference element [-1, 1]
400d1d35e2fSjeremylt   @param q_weight_1d Array of length Q_1d holding the quadrature weights on the
401b11c1e72Sjeremylt                        reference element
402b11c1e72Sjeremylt   @param[out] basis  Address of the variable where the newly created
403b11c1e72Sjeremylt                        CeedBasis will be stored.
404b11c1e72Sjeremylt 
405b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
406dfdf5a53Sjeremylt 
4077a982d89SJeremy L. Thompson   @ref User
408b11c1e72Sjeremylt **/
409d1d35e2fSjeremylt int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt num_comp,
410d1d35e2fSjeremylt                             CeedInt P_1d, CeedInt Q_1d,
411d1d35e2fSjeremylt                             const CeedScalar *interp_1d,
412d1d35e2fSjeremylt                             const CeedScalar *grad_1d, const CeedScalar *q_ref_1d,
413d1d35e2fSjeremylt                             const CeedScalar *q_weight_1d, CeedBasis *basis) {
414d7b241e6Sjeremylt   int ierr;
415d7b241e6Sjeremylt 
4165fe0d4faSjeremylt   if (!ceed->BasisCreateTensorH1) {
4175fe0d4faSjeremylt     Ceed delegate;
418aefd8378Sjeremylt     ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr);
4195fe0d4faSjeremylt 
4205fe0d4faSjeremylt     if (!delegate)
421c042f62fSJeremy L Thompson       // LCOV_EXCL_START
422e15f9bd0SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_UNSUPPORTED,
423e15f9bd0SJeremy L Thompson                        "Backend does not support BasisCreateTensorH1");
424c042f62fSJeremy L Thompson     // LCOV_EXCL_STOP
4255fe0d4faSjeremylt 
426d1d35e2fSjeremylt     ierr = CeedBasisCreateTensorH1(delegate, dim, num_comp, P_1d,
427d1d35e2fSjeremylt                                    Q_1d, interp_1d, grad_1d, q_ref_1d,
428d1d35e2fSjeremylt                                    q_weight_1d, basis); CeedChk(ierr);
429e15f9bd0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
4305fe0d4faSjeremylt   }
431e15f9bd0SJeremy L Thompson 
432e15f9bd0SJeremy L Thompson   if (dim<1)
433e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
434e15f9bd0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
435e15f9bd0SJeremy L Thompson                      "Basis dimension must be a positive value");
436e15f9bd0SJeremy L Thompson   // LCOV_EXCL_STOP
437d1d35e2fSjeremylt   CeedElemTopology topo = dim == 1 ? CEED_LINE
438d1d35e2fSjeremylt                           : dim == 2 ? CEED_QUAD
439d1d35e2fSjeremylt                           : CEED_HEX;
440e15f9bd0SJeremy L Thompson 
441d7b241e6Sjeremylt   ierr = CeedCalloc(1, basis); CeedChk(ierr);
442d7b241e6Sjeremylt   (*basis)->ceed = ceed;
4439560d06aSjeremylt   ierr = CeedReference(ceed); CeedChk(ierr);
444d1d35e2fSjeremylt   (*basis)->ref_count = 1;
445d1d35e2fSjeremylt   (*basis)->tensor_basis = 1;
446d7b241e6Sjeremylt   (*basis)->dim = dim;
447d99fa3c5SJeremy L Thompson   (*basis)->topo = topo;
448d1d35e2fSjeremylt   (*basis)->num_comp = num_comp;
449d1d35e2fSjeremylt   (*basis)->P_1d = P_1d;
450d1d35e2fSjeremylt   (*basis)->Q_1d = Q_1d;
451d1d35e2fSjeremylt   (*basis)->P = CeedIntPow(P_1d, dim);
452d1d35e2fSjeremylt   (*basis)->Q = CeedIntPow(Q_1d, dim);
453*ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(Q_1d, &(*basis)->q_ref_1d); CeedChk(ierr);
454*ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(Q_1d, &(*basis)->q_weight_1d); CeedChk(ierr);
455*ff3a0f91SJeremy L Thompson   if (q_ref_1d) memcpy((*basis)->q_ref_1d, q_ref_1d, Q_1d*sizeof(q_ref_1d[0]));
456*ff3a0f91SJeremy L Thompson   if (q_weight_1d) memcpy((*basis)->q_weight_1d, q_weight_1d,
457*ff3a0f91SJeremy L Thompson                             Q_1d*sizeof(q_weight_1d[0]));
458*ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(Q_1d*P_1d, &(*basis)->interp_1d); CeedChk(ierr);
459*ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(Q_1d*P_1d, &(*basis)->grad_1d); CeedChk(ierr);
460*ff3a0f91SJeremy L Thompson   if (interp_1d) memcpy((*basis)->interp_1d, interp_1d,
461*ff3a0f91SJeremy L Thompson                           Q_1d*P_1d*sizeof(interp_1d[0]));
462*ff3a0f91SJeremy L Thompson   if (grad_1d) memcpy((*basis)->grad_1d, grad_1d, Q_1d*P_1d*sizeof(grad_1d[0]));
463d1d35e2fSjeremylt   ierr = ceed->BasisCreateTensorH1(dim, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d,
464d1d35e2fSjeremylt                                    q_weight_1d, *basis); CeedChk(ierr);
465e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
466d7b241e6Sjeremylt }
467d7b241e6Sjeremylt 
468b11c1e72Sjeremylt /**
46995bb1877Svaleriabarra   @brief Create a tensor-product Lagrange basis
470b11c1e72Sjeremylt 
471b11c1e72Sjeremylt   @param ceed        A Ceed object where the CeedBasis will be created
472b11c1e72Sjeremylt   @param dim         Topological dimension of element
473d1d35e2fSjeremylt   @param num_comp      Number of field components (1 for scalar fields)
474b11c1e72Sjeremylt   @param P           Number of Gauss-Lobatto nodes in one dimension.  The
475b11c1e72Sjeremylt                        polynomial degree of the resulting Q_k element is k=P-1.
476b11c1e72Sjeremylt   @param Q           Number of quadrature points in one dimension.
477d1d35e2fSjeremylt   @param quad_mode   Distribution of the Q quadrature points (affects order of
478b11c1e72Sjeremylt                        accuracy for the quadrature)
479b11c1e72Sjeremylt   @param[out] basis  Address of the variable where the newly created
480b11c1e72Sjeremylt                        CeedBasis will be stored.
481b11c1e72Sjeremylt 
482b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
483dfdf5a53Sjeremylt 
4847a982d89SJeremy L. Thompson   @ref User
485b11c1e72Sjeremylt **/
486d1d35e2fSjeremylt int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt num_comp,
487d1d35e2fSjeremylt                                     CeedInt P, CeedInt Q, CeedQuadMode quad_mode,
488692c2638Sjeremylt                                     CeedBasis *basis) {
489d7b241e6Sjeremylt   // Allocate
490e15f9bd0SJeremy L Thompson   int ierr, ierr2, i, j, k;
491d1d35e2fSjeremylt   CeedScalar c1, c2, c3, c4, dx, *nodes, *interp_1d, *grad_1d, *q_ref_1d,
492d1d35e2fSjeremylt              *q_weight_1d;
4934d537eeaSYohann 
4944d537eeaSYohann   if (dim<1)
495c042f62fSJeremy L Thompson     // LCOV_EXCL_START
496e15f9bd0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
497e15f9bd0SJeremy L Thompson                      "Basis dimension must be a positive value");
498c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
4994d537eeaSYohann 
500e15f9bd0SJeremy L Thompson   // Get Nodes and Weights
501d1d35e2fSjeremylt   ierr = CeedCalloc(P*Q, &interp_1d); CeedChk(ierr);
502d1d35e2fSjeremylt   ierr = CeedCalloc(P*Q, &grad_1d); CeedChk(ierr);
503d7b241e6Sjeremylt   ierr = CeedCalloc(P, &nodes); CeedChk(ierr);
504d1d35e2fSjeremylt   ierr = CeedCalloc(Q, &q_ref_1d); CeedChk(ierr);
505d1d35e2fSjeremylt   ierr = CeedCalloc(Q, &q_weight_1d); CeedChk(ierr);
506e15f9bd0SJeremy L Thompson   ierr = CeedLobattoQuadrature(P, nodes, NULL);
507e15f9bd0SJeremy L Thompson   if (ierr) { goto cleanup; } CeedChk(ierr);
508d1d35e2fSjeremylt   switch (quad_mode) {
509d7b241e6Sjeremylt   case CEED_GAUSS:
510d1d35e2fSjeremylt     ierr = CeedGaussQuadrature(Q, q_ref_1d, q_weight_1d);
511d7b241e6Sjeremylt     break;
512d7b241e6Sjeremylt   case CEED_GAUSS_LOBATTO:
513d1d35e2fSjeremylt     ierr = CeedLobattoQuadrature(Q, q_ref_1d, q_weight_1d);
514d7b241e6Sjeremylt     break;
515d7b241e6Sjeremylt   }
516e15f9bd0SJeremy L Thompson   if (ierr) { goto cleanup; } CeedChk(ierr);
517e15f9bd0SJeremy L Thompson 
518d7b241e6Sjeremylt   // Build B, D matrix
519d7b241e6Sjeremylt   // Fornberg, 1998
520d7b241e6Sjeremylt   for (i = 0; i  < Q; i++) {
521d7b241e6Sjeremylt     c1 = 1.0;
522d1d35e2fSjeremylt     c3 = nodes[0] - q_ref_1d[i];
523d1d35e2fSjeremylt     interp_1d[i*P+0] = 1.0;
524d7b241e6Sjeremylt     for (j = 1; j < P; j++) {
525d7b241e6Sjeremylt       c2 = 1.0;
526d7b241e6Sjeremylt       c4 = c3;
527d1d35e2fSjeremylt       c3 = nodes[j] - q_ref_1d[i];
528d7b241e6Sjeremylt       for (k = 0; k < j; k++) {
529d7b241e6Sjeremylt         dx = nodes[j] - nodes[k];
530d7b241e6Sjeremylt         c2 *= dx;
531d7b241e6Sjeremylt         if (k == j - 1) {
532d1d35e2fSjeremylt           grad_1d[i*P + j] = c1*(interp_1d[i*P + k] - c4*grad_1d[i*P + k]) / c2;
533d1d35e2fSjeremylt           interp_1d[i*P + j] = - c1*c4*interp_1d[i*P + k] / c2;
534d7b241e6Sjeremylt         }
535d1d35e2fSjeremylt         grad_1d[i*P + k] = (c3*grad_1d[i*P + k] - interp_1d[i*P + k]) / dx;
536d1d35e2fSjeremylt         interp_1d[i*P + k] = c3*interp_1d[i*P + k] / dx;
537d7b241e6Sjeremylt       }
538d7b241e6Sjeremylt       c1 = c2;
539d7b241e6Sjeremylt     }
540d7b241e6Sjeremylt   }
5419ac7b42eSJeremy L Thompson   // Pass to CeedBasisCreateTensorH1
542d1d35e2fSjeremylt   ierr = CeedBasisCreateTensorH1(ceed, dim, num_comp, P, Q, interp_1d, grad_1d,
5439ac7b42eSJeremy L Thompson                                  q_ref_1d, q_weight_1d, basis); CeedChk(ierr);
544e15f9bd0SJeremy L Thompson cleanup:
545d1d35e2fSjeremylt   ierr2 = CeedFree(&interp_1d); CeedChk(ierr2);
546d1d35e2fSjeremylt   ierr2 = CeedFree(&grad_1d); CeedChk(ierr2);
547e15f9bd0SJeremy L Thompson   ierr2 = CeedFree(&nodes); CeedChk(ierr2);
548d1d35e2fSjeremylt   ierr2 = CeedFree(&q_ref_1d); CeedChk(ierr2);
549d1d35e2fSjeremylt   ierr2 = CeedFree(&q_weight_1d); CeedChk(ierr2);
550e15f9bd0SJeremy L Thompson   CeedChk(ierr);
551e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
552d7b241e6Sjeremylt }
553d7b241e6Sjeremylt 
554b11c1e72Sjeremylt /**
55595bb1877Svaleriabarra   @brief Create a non tensor-product basis for H^1 discretizations
556a8de75f0Sjeremylt 
557a8de75f0Sjeremylt   @param ceed        A Ceed object where the CeedBasis will be created
558a8de75f0Sjeremylt   @param topo        Topology of element, e.g. hypercube, simplex, ect
559d1d35e2fSjeremylt   @param num_comp    Number of field components (1 for scalar fields)
560d1d35e2fSjeremylt   @param num_nodes   Total number of nodes
561d1d35e2fSjeremylt   @param num_qpts    Total number of quadrature points
562d1d35e2fSjeremylt   @param interp      Row-major (num_qpts * num_nodes) matrix expressing the values of
5638795c945Sjeremylt                        nodal basis functions at quadrature points
564d1d35e2fSjeremylt   @param grad        Row-major (num_qpts * dim * num_nodes) matrix expressing
5658795c945Sjeremylt                        derivatives of nodal basis functions at quadrature points
566d1d35e2fSjeremylt   @param q_ref       Array of length num_qpts holding the locations of quadrature
5678795c945Sjeremylt                        points on the reference element [-1, 1]
568d1d35e2fSjeremylt   @param q_weight    Array of length num_qpts holding the quadrature weights on the
569a8de75f0Sjeremylt                        reference element
570a8de75f0Sjeremylt   @param[out] basis  Address of the variable where the newly created
571a8de75f0Sjeremylt                        CeedBasis will be stored.
572a8de75f0Sjeremylt 
573a8de75f0Sjeremylt   @return An error code: 0 - success, otherwise - failure
574a8de75f0Sjeremylt 
5757a982d89SJeremy L. Thompson   @ref User
576a8de75f0Sjeremylt **/
577d1d35e2fSjeremylt int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt num_comp,
578d1d35e2fSjeremylt                       CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
579d1d35e2fSjeremylt                       const CeedScalar *grad, const CeedScalar *q_ref,
580d1d35e2fSjeremylt                       const CeedScalar *q_weight, CeedBasis *basis) {
581a8de75f0Sjeremylt   int ierr;
582d1d35e2fSjeremylt   CeedInt P = num_nodes, Q = num_qpts, dim = 0;
583a8de75f0Sjeremylt 
5845fe0d4faSjeremylt   if (!ceed->BasisCreateH1) {
5855fe0d4faSjeremylt     Ceed delegate;
586aefd8378Sjeremylt     ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr);
5875fe0d4faSjeremylt 
5885fe0d4faSjeremylt     if (!delegate)
589c042f62fSJeremy L Thompson       // LCOV_EXCL_START
590e15f9bd0SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_UNSUPPORTED,
591e15f9bd0SJeremy L Thompson                        "Backend does not support BasisCreateH1");
592c042f62fSJeremy L Thompson     // LCOV_EXCL_STOP
5935fe0d4faSjeremylt 
594d1d35e2fSjeremylt     ierr = CeedBasisCreateH1(delegate, topo, num_comp, num_nodes,
595d1d35e2fSjeremylt                              num_qpts, interp, grad, q_ref,
596d1d35e2fSjeremylt                              q_weight, basis); CeedChk(ierr);
597e15f9bd0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
5985fe0d4faSjeremylt   }
5995fe0d4faSjeremylt 
600a8de75f0Sjeremylt   ierr = CeedCalloc(1, basis); CeedChk(ierr);
601a8de75f0Sjeremylt 
602a8de75f0Sjeremylt   ierr = CeedBasisGetTopologyDimension(topo, &dim); CeedChk(ierr);
603a8de75f0Sjeremylt 
604a8de75f0Sjeremylt   (*basis)->ceed = ceed;
6059560d06aSjeremylt   ierr = CeedReference(ceed); CeedChk(ierr);
606d1d35e2fSjeremylt   (*basis)->ref_count = 1;
607d1d35e2fSjeremylt   (*basis)->tensor_basis = 0;
608a8de75f0Sjeremylt   (*basis)->dim = dim;
609d99fa3c5SJeremy L Thompson   (*basis)->topo = topo;
610d1d35e2fSjeremylt   (*basis)->num_comp = num_comp;
611a8de75f0Sjeremylt   (*basis)->P = P;
612a8de75f0Sjeremylt   (*basis)->Q = Q;
613*ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(Q*dim, &(*basis)->q_ref_1d); CeedChk(ierr);
614*ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(Q, &(*basis)->q_weight_1d); CeedChk(ierr);
615*ff3a0f91SJeremy L Thompson   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q*dim*sizeof(q_ref[0]));
616*ff3a0f91SJeremy L Thompson   if(q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q*sizeof(q_weight[0]));
617*ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(Q*P, &(*basis)->interp); CeedChk(ierr);
618*ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(dim*Q*P, &(*basis)->grad); CeedChk(ierr);
619*ff3a0f91SJeremy L Thompson   if(interp) memcpy((*basis)->interp, interp, Q*P*sizeof(interp[0]));
620*ff3a0f91SJeremy L Thompson   if(grad) memcpy((*basis)->grad, grad, dim*Q*P*sizeof(grad[0]));
621d1d35e2fSjeremylt   ierr = ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, q_ref,
622d1d35e2fSjeremylt                              q_weight, *basis); CeedChk(ierr);
623e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
624a8de75f0Sjeremylt }
625a8de75f0Sjeremylt 
626a8de75f0Sjeremylt /**
6279560d06aSjeremylt   @brief Copy the pointer to a CeedBasis. Both pointers should
6289560d06aSjeremylt            be destroyed with `CeedBasisDestroy()`;
6299560d06aSjeremylt            Note: If `*basis_copy` is non-NULL, then it is assumed that
6309560d06aSjeremylt            `*basis_copy` is a pointer to a CeedBasis. This CeedBasis
6319560d06aSjeremylt            will be destroyed if `*basis_copy` is the only
6329560d06aSjeremylt            reference to this CeedBasis.
6339560d06aSjeremylt 
6349560d06aSjeremylt   @param basis            CeedBasis to copy reference to
6359560d06aSjeremylt   @param[out] basis_copy  Variable to store copied reference
6369560d06aSjeremylt 
6379560d06aSjeremylt   @return An error code: 0 - success, otherwise - failure
6389560d06aSjeremylt 
6399560d06aSjeremylt   @ref User
6409560d06aSjeremylt **/
6419560d06aSjeremylt int CeedBasisReferenceCopy(CeedBasis basis, CeedBasis *basis_copy) {
6429560d06aSjeremylt   int ierr;
6439560d06aSjeremylt 
6449560d06aSjeremylt   ierr = CeedBasisReference(basis); CeedChk(ierr);
6459560d06aSjeremylt   ierr = CeedBasisDestroy(basis_copy); CeedChk(ierr);
6469560d06aSjeremylt   *basis_copy = basis;
6479560d06aSjeremylt   return CEED_ERROR_SUCCESS;
6489560d06aSjeremylt }
6499560d06aSjeremylt 
6509560d06aSjeremylt /**
6517a982d89SJeremy L. Thompson   @brief View a CeedBasis
6527a982d89SJeremy L. Thompson 
6537a982d89SJeremy L. Thompson   @param basis   CeedBasis to view
6547a982d89SJeremy L. Thompson   @param stream  Stream to view to, e.g., stdout
6557a982d89SJeremy L. Thompson 
6567a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
6577a982d89SJeremy L. Thompson 
6587a982d89SJeremy L. Thompson   @ref User
6597a982d89SJeremy L. Thompson **/
6607a982d89SJeremy L. Thompson int CeedBasisView(CeedBasis basis, FILE *stream) {
6617a982d89SJeremy L. Thompson   int ierr;
6627a982d89SJeremy L. Thompson 
663d1d35e2fSjeremylt   if (basis->tensor_basis) {
664d1d35e2fSjeremylt     fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P_1d,
665d1d35e2fSjeremylt             basis->Q_1d);
666d1d35e2fSjeremylt     ierr = CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_ref_1d,
6677a982d89SJeremy L. Thompson                           stream); CeedChk(ierr);
668d1d35e2fSjeremylt     ierr = CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q_1d,
669d1d35e2fSjeremylt                           basis->q_weight_1d, stream); CeedChk(ierr);
670d1d35e2fSjeremylt     ierr = CeedScalarView("interp1d", "\t% 12.8f", basis->Q_1d, basis->P_1d,
671d1d35e2fSjeremylt                           basis->interp_1d, stream); CeedChk(ierr);
672d1d35e2fSjeremylt     ierr = CeedScalarView("grad1d", "\t% 12.8f", basis->Q_1d, basis->P_1d,
673d1d35e2fSjeremylt                           basis->grad_1d, stream); CeedChk(ierr);
6747a982d89SJeremy L. Thompson   } else {
6757a982d89SJeremy L. Thompson     fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P,
6767a982d89SJeremy L. Thompson             basis->Q);
6777a982d89SJeremy L. Thompson     ierr = CeedScalarView("qref", "\t% 12.8f", 1, basis->Q*basis->dim,
678d1d35e2fSjeremylt                           basis->q_ref_1d,
6797a982d89SJeremy L. Thompson                           stream); CeedChk(ierr);
680d1d35e2fSjeremylt     ierr = CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->q_weight_1d,
6817a982d89SJeremy L. Thompson                           stream); CeedChk(ierr);
6827a982d89SJeremy L. Thompson     ierr = CeedScalarView("interp", "\t% 12.8f", basis->Q, basis->P,
6837a982d89SJeremy L. Thompson                           basis->interp, stream); CeedChk(ierr);
6847a982d89SJeremy L. Thompson     ierr = CeedScalarView("grad", "\t% 12.8f", basis->dim*basis->Q, basis->P,
6857a982d89SJeremy L. Thompson                           basis->grad, stream); CeedChk(ierr);
6867a982d89SJeremy L. Thompson   }
687e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
6887a982d89SJeremy L. Thompson }
6897a982d89SJeremy L. Thompson 
6907a982d89SJeremy L. Thompson /**
6917a982d89SJeremy L. Thompson   @brief Apply basis evaluation from nodes to quadrature points or vice versa
6927a982d89SJeremy L. Thompson 
6937a982d89SJeremy L. Thompson   @param basis     CeedBasis to evaluate
694d1d35e2fSjeremylt   @param num_elem  The number of elements to apply the basis evaluation to;
6957a982d89SJeremy L. Thompson                      the backend will specify the ordering in
6964cc79fe7SJed Brown                      CeedElemRestrictionCreateBlocked()
697d1d35e2fSjeremylt   @param t_mode    \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature
6987a982d89SJeremy L. Thompson                      points, \ref CEED_TRANSPOSE to apply the transpose, mapping
6997a982d89SJeremy L. Thompson                      from quadrature points to nodes
700d1d35e2fSjeremylt   @param eval_mode \ref CEED_EVAL_NONE to use values directly,
7017a982d89SJeremy L. Thompson                      \ref CEED_EVAL_INTERP to use interpolated values,
7027a982d89SJeremy L. Thompson                      \ref CEED_EVAL_GRAD to use gradients,
7037a982d89SJeremy L. Thompson                      \ref CEED_EVAL_WEIGHT to use quadrature weights.
7047a982d89SJeremy L. Thompson   @param[in] u     Input CeedVector
7057a982d89SJeremy L. Thompson   @param[out] v    Output CeedVector
7067a982d89SJeremy L. Thompson 
7077a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
7087a982d89SJeremy L. Thompson 
7097a982d89SJeremy L. Thompson   @ref User
7107a982d89SJeremy L. Thompson **/
711d1d35e2fSjeremylt int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode,
712d1d35e2fSjeremylt                    CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
7137a982d89SJeremy L. Thompson   int ierr;
714d1d35e2fSjeremylt   CeedInt u_length = 0, v_length, dim, num_comp, num_nodes, num_qpts;
715e15f9bd0SJeremy L Thompson   ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr);
716d1d35e2fSjeremylt   ierr = CeedBasisGetNumComponents(basis, &num_comp); CeedChk(ierr);
717d1d35e2fSjeremylt   ierr = CeedBasisGetNumNodes(basis, &num_nodes); CeedChk(ierr);
718d1d35e2fSjeremylt   ierr = CeedBasisGetNumQuadraturePoints(basis, &num_qpts); CeedChk(ierr);
719d1d35e2fSjeremylt   ierr = CeedVectorGetLength(v, &v_length); CeedChk(ierr);
7207a982d89SJeremy L. Thompson   if (u) {
721d1d35e2fSjeremylt     ierr = CeedVectorGetLength(u, &u_length); CeedChk(ierr);
7227a982d89SJeremy L. Thompson   }
7237a982d89SJeremy L. Thompson 
724e15f9bd0SJeremy L Thompson   if (!basis->Apply)
725e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
726e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_UNSUPPORTED,
727e15f9bd0SJeremy L Thompson                      "Backend does not support BasisApply");
728e15f9bd0SJeremy L Thompson   // LCOV_EXCL_STOP
729e15f9bd0SJeremy L Thompson 
730e15f9bd0SJeremy L Thompson   // Check compatibility of topological and geometrical dimensions
731d1d35e2fSjeremylt   if ((t_mode == CEED_TRANSPOSE && (v_length%num_nodes != 0 ||
732d1d35e2fSjeremylt                                     u_length%num_qpts != 0)) ||
733d1d35e2fSjeremylt       (t_mode == CEED_NOTRANSPOSE && (u_length%num_nodes != 0 ||
734d1d35e2fSjeremylt                                       v_length%num_qpts != 0)))
7358229195eSjeremylt     // LCOV_EXCL_START
736e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_DIMENSION,
737e15f9bd0SJeremy L Thompson                      "Length of input/output vectors "
7387a982d89SJeremy L. Thompson                      "incompatible with basis dimensions");
7398229195eSjeremylt   // LCOV_EXCL_STOP
7407a982d89SJeremy L. Thompson 
741e15f9bd0SJeremy L Thompson   // Check vector lengths to prevent out of bounds issues
742d1d35e2fSjeremylt   bool bad_dims = false;
743d1d35e2fSjeremylt   switch (eval_mode) {
744e15f9bd0SJeremy L Thompson   case CEED_EVAL_NONE:
745d1d35e2fSjeremylt   case CEED_EVAL_INTERP: bad_dims =
746d1d35e2fSjeremylt       ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts ||
747d1d35e2fSjeremylt                                      v_length < num_elem*num_comp*num_nodes)) ||
748d1d35e2fSjeremylt        (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp ||
749d1d35e2fSjeremylt                                        u_length < num_elem*num_comp*num_nodes)));
750e15f9bd0SJeremy L Thompson     break;
751d1d35e2fSjeremylt   case CEED_EVAL_GRAD: bad_dims =
752d1d35e2fSjeremylt       ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts*dim ||
753d1d35e2fSjeremylt                                      v_length < num_elem*num_comp*num_nodes)) ||
754d1d35e2fSjeremylt        (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp*dim ||
755d1d35e2fSjeremylt                                        u_length < num_elem*num_comp*num_nodes)));
756e15f9bd0SJeremy L Thompson     break;
757e15f9bd0SJeremy L Thompson   case CEED_EVAL_WEIGHT:
758d1d35e2fSjeremylt     bad_dims = v_length < num_elem*num_qpts;
759e15f9bd0SJeremy L Thompson     break;
760e15f9bd0SJeremy L Thompson   // LCOV_EXCL_START
761d1d35e2fSjeremylt   case CEED_EVAL_DIV: bad_dims =
762d1d35e2fSjeremylt       ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts ||
763d1d35e2fSjeremylt                                      v_length < num_elem*num_comp*num_nodes)) ||
764d1d35e2fSjeremylt        (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp ||
765d1d35e2fSjeremylt                                        u_length < num_elem*num_comp*num_nodes)));
766e15f9bd0SJeremy L Thompson     break;
767d1d35e2fSjeremylt   case CEED_EVAL_CURL: bad_dims =
768d1d35e2fSjeremylt       ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts ||
769d1d35e2fSjeremylt                                      v_length < num_elem*num_comp*num_nodes)) ||
770d1d35e2fSjeremylt        (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp ||
771d1d35e2fSjeremylt                                        u_length < num_elem*num_comp*num_nodes)));
772e15f9bd0SJeremy L Thompson     break;
773e15f9bd0SJeremy L Thompson     // LCOV_EXCL_STOP
774e15f9bd0SJeremy L Thompson   }
775d1d35e2fSjeremylt   if (bad_dims)
776e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
777e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_DIMENSION,
778d1d35e2fSjeremylt                      "Input/output vectors too short for basis and evaluation mode");
779e15f9bd0SJeremy L Thompson   // LCOV_EXCL_STOP
780e15f9bd0SJeremy L Thompson 
781d1d35e2fSjeremylt   ierr = basis->Apply(basis, num_elem, t_mode, eval_mode, u, v); CeedChk(ierr);
782e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
7837a982d89SJeremy L. Thompson }
7847a982d89SJeremy L. Thompson 
7857a982d89SJeremy L. Thompson /**
786b7c9bbdaSJeremy L Thompson   @brief Get Ceed associated with a CeedBasis
787b7c9bbdaSJeremy L Thompson 
788b7c9bbdaSJeremy L Thompson   @param basis      CeedBasis
789b7c9bbdaSJeremy L Thompson   @param[out] ceed  Variable to store Ceed
790b7c9bbdaSJeremy L Thompson 
791b7c9bbdaSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
792b7c9bbdaSJeremy L Thompson 
793b7c9bbdaSJeremy L Thompson   @ref Advanced
794b7c9bbdaSJeremy L Thompson **/
795b7c9bbdaSJeremy L Thompson int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) {
796b7c9bbdaSJeremy L Thompson   *ceed = basis->ceed;
797b7c9bbdaSJeremy L Thompson   return CEED_ERROR_SUCCESS;
798b7c9bbdaSJeremy L Thompson }
799b7c9bbdaSJeremy L Thompson 
800b7c9bbdaSJeremy L Thompson /**
8019d007619Sjeremylt   @brief Get dimension for given CeedBasis
8029d007619Sjeremylt 
8039d007619Sjeremylt   @param basis     CeedBasis
8049d007619Sjeremylt   @param[out] dim  Variable to store dimension of basis
8059d007619Sjeremylt 
8069d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
8079d007619Sjeremylt 
808b7c9bbdaSJeremy L Thompson   @ref Advanced
8099d007619Sjeremylt **/
8109d007619Sjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) {
8119d007619Sjeremylt   *dim = basis->dim;
812e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8139d007619Sjeremylt }
8149d007619Sjeremylt 
8159d007619Sjeremylt /**
816d99fa3c5SJeremy L Thompson   @brief Get topology for given CeedBasis
817d99fa3c5SJeremy L Thompson 
818d99fa3c5SJeremy L Thompson   @param basis      CeedBasis
819d99fa3c5SJeremy L Thompson   @param[out] topo  Variable to store topology of basis
820d99fa3c5SJeremy L Thompson 
821d99fa3c5SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
822d99fa3c5SJeremy L Thompson 
823b7c9bbdaSJeremy L Thompson   @ref Advanced
824d99fa3c5SJeremy L Thompson **/
825d99fa3c5SJeremy L Thompson int CeedBasisGetTopology(CeedBasis basis, CeedElemTopology *topo) {
826d99fa3c5SJeremy L Thompson   *topo = basis->topo;
827e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
828d99fa3c5SJeremy L Thompson }
829d99fa3c5SJeremy L Thompson 
830d99fa3c5SJeremy L Thompson /**
8319d007619Sjeremylt   @brief Get number of components for given CeedBasis
8329d007619Sjeremylt 
8339d007619Sjeremylt   @param basis          CeedBasis
834d1d35e2fSjeremylt   @param[out] num_comp  Variable to store number of components of basis
8359d007619Sjeremylt 
8369d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
8379d007619Sjeremylt 
838b7c9bbdaSJeremy L Thompson   @ref Advanced
8399d007619Sjeremylt **/
840d1d35e2fSjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *num_comp) {
841d1d35e2fSjeremylt   *num_comp = basis->num_comp;
842e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8439d007619Sjeremylt }
8449d007619Sjeremylt 
8459d007619Sjeremylt /**
8469d007619Sjeremylt   @brief Get total number of nodes (in dim dimensions) of a CeedBasis
8479d007619Sjeremylt 
8489d007619Sjeremylt   @param basis   CeedBasis
8499d007619Sjeremylt   @param[out] P  Variable to store number of nodes
8509d007619Sjeremylt 
8519d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
8529d007619Sjeremylt 
8539d007619Sjeremylt   @ref Utility
8549d007619Sjeremylt **/
8559d007619Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) {
8569d007619Sjeremylt   *P = basis->P;
857e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8589d007619Sjeremylt }
8599d007619Sjeremylt 
8609d007619Sjeremylt /**
8619d007619Sjeremylt   @brief Get total number of nodes (in 1 dimension) of a CeedBasis
8629d007619Sjeremylt 
8639d007619Sjeremylt   @param basis     CeedBasis
864d1d35e2fSjeremylt   @param[out] P_1d  Variable to store number of nodes
8659d007619Sjeremylt 
8669d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
8679d007619Sjeremylt 
868b7c9bbdaSJeremy L Thompson   @ref Advanced
8699d007619Sjeremylt **/
870d1d35e2fSjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P_1d) {
871d1d35e2fSjeremylt   if (!basis->tensor_basis)
8729d007619Sjeremylt     // LCOV_EXCL_START
873e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR,
874d1d35e2fSjeremylt                      "Cannot supply P_1d for non-tensor basis");
8759d007619Sjeremylt   // LCOV_EXCL_STOP
8769d007619Sjeremylt 
877d1d35e2fSjeremylt   *P_1d = basis->P_1d;
878e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8799d007619Sjeremylt }
8809d007619Sjeremylt 
8819d007619Sjeremylt /**
8829d007619Sjeremylt   @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis
8839d007619Sjeremylt 
8849d007619Sjeremylt   @param basis   CeedBasis
8859d007619Sjeremylt   @param[out] Q  Variable to store number of quadrature points
8869d007619Sjeremylt 
8879d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
8889d007619Sjeremylt 
8899d007619Sjeremylt   @ref Utility
8909d007619Sjeremylt **/
8919d007619Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) {
8929d007619Sjeremylt   *Q = basis->Q;
893e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8949d007619Sjeremylt }
8959d007619Sjeremylt 
8969d007619Sjeremylt /**
8979d007619Sjeremylt   @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis
8989d007619Sjeremylt 
8999d007619Sjeremylt   @param basis      CeedBasis
900d1d35e2fSjeremylt   @param[out] Q_1d  Variable to store number of quadrature points
9019d007619Sjeremylt 
9029d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
9039d007619Sjeremylt 
904b7c9bbdaSJeremy L Thompson   @ref Advanced
9059d007619Sjeremylt **/
906d1d35e2fSjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q_1d) {
907d1d35e2fSjeremylt   if (!basis->tensor_basis)
9089d007619Sjeremylt     // LCOV_EXCL_START
909e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR,
910d1d35e2fSjeremylt                      "Cannot supply Q_1d for non-tensor basis");
9119d007619Sjeremylt   // LCOV_EXCL_STOP
9129d007619Sjeremylt 
913d1d35e2fSjeremylt   *Q_1d = basis->Q_1d;
914e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9159d007619Sjeremylt }
9169d007619Sjeremylt 
9179d007619Sjeremylt /**
9189d007619Sjeremylt   @brief Get reference coordinates of quadrature points (in dim dimensions)
9199d007619Sjeremylt          of a CeedBasis
9209d007619Sjeremylt 
9219d007619Sjeremylt   @param basis       CeedBasis
922d1d35e2fSjeremylt   @param[out] q_ref  Variable to store reference coordinates of quadrature points
9239d007619Sjeremylt 
9249d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
9259d007619Sjeremylt 
926b7c9bbdaSJeremy L Thompson   @ref Advanced
9279d007619Sjeremylt **/
928d1d35e2fSjeremylt int CeedBasisGetQRef(CeedBasis basis, const CeedScalar **q_ref) {
929d1d35e2fSjeremylt   *q_ref = basis->q_ref_1d;
930e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9319d007619Sjeremylt }
9329d007619Sjeremylt 
9339d007619Sjeremylt /**
9349d007619Sjeremylt   @brief Get quadrature weights of quadrature points (in dim dimensions)
9359d007619Sjeremylt          of a CeedBasis
9369d007619Sjeremylt 
9379d007619Sjeremylt   @param basis          CeedBasis
938d1d35e2fSjeremylt   @param[out] q_weight  Variable to store quadrature weights
9399d007619Sjeremylt 
9409d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
9419d007619Sjeremylt 
942b7c9bbdaSJeremy L Thompson   @ref Advanced
9439d007619Sjeremylt **/
944d1d35e2fSjeremylt int CeedBasisGetQWeights(CeedBasis basis, const CeedScalar **q_weight) {
945d1d35e2fSjeremylt   *q_weight = basis->q_weight_1d;
946e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9479d007619Sjeremylt }
9489d007619Sjeremylt 
9499d007619Sjeremylt /**
9509d007619Sjeremylt   @brief Get interpolation matrix of a CeedBasis
9519d007619Sjeremylt 
9529d007619Sjeremylt   @param basis        CeedBasis
9539d007619Sjeremylt   @param[out] interp  Variable to store interpolation matrix
9549d007619Sjeremylt 
9559d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
9569d007619Sjeremylt 
957b7c9bbdaSJeremy L Thompson   @ref Advanced
9589d007619Sjeremylt **/
9596c58de82SJeremy L Thompson int CeedBasisGetInterp(CeedBasis basis, const CeedScalar **interp) {
960d1d35e2fSjeremylt   if (!basis->interp && basis->tensor_basis) {
9619d007619Sjeremylt     // Allocate
9629d007619Sjeremylt     int ierr;
9639d007619Sjeremylt     ierr = CeedMalloc(basis->Q*basis->P, &basis->interp); CeedChk(ierr);
9649d007619Sjeremylt 
9659d007619Sjeremylt     // Initialize
9669d007619Sjeremylt     for (CeedInt i=0; i<basis->Q*basis->P; i++)
9679d007619Sjeremylt       basis->interp[i] = 1.0;
9689d007619Sjeremylt 
9699d007619Sjeremylt     // Calculate
9709d007619Sjeremylt     for (CeedInt d=0; d<basis->dim; d++)
9719d007619Sjeremylt       for (CeedInt qpt=0; qpt<basis->Q; qpt++)
9729d007619Sjeremylt         for (CeedInt node=0; node<basis->P; node++) {
973d1d35e2fSjeremylt           CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d;
974d1d35e2fSjeremylt           CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d;
975d1d35e2fSjeremylt           basis->interp[qpt*(basis->P)+node] *= basis->interp_1d[q*basis->P_1d+p];
9769d007619Sjeremylt         }
9779d007619Sjeremylt   }
9789d007619Sjeremylt   *interp = basis->interp;
979e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9809d007619Sjeremylt }
9819d007619Sjeremylt 
9829d007619Sjeremylt /**
9839d007619Sjeremylt   @brief Get 1D interpolation matrix of a tensor product CeedBasis
9849d007619Sjeremylt 
9859d007619Sjeremylt   @param basis           CeedBasis
986d1d35e2fSjeremylt   @param[out] interp_1d  Variable to store interpolation matrix
9879d007619Sjeremylt 
9889d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
9899d007619Sjeremylt 
9909d007619Sjeremylt   @ref Backend
9919d007619Sjeremylt **/
992d1d35e2fSjeremylt int CeedBasisGetInterp1D(CeedBasis basis, const CeedScalar **interp_1d) {
993d1d35e2fSjeremylt   if (!basis->tensor_basis)
9949d007619Sjeremylt     // LCOV_EXCL_START
995e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR,
996e15f9bd0SJeremy L Thompson                      "CeedBasis is not a tensor product basis.");
9979d007619Sjeremylt   // LCOV_EXCL_STOP
9989d007619Sjeremylt 
999d1d35e2fSjeremylt   *interp_1d = basis->interp_1d;
1000e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
10019d007619Sjeremylt }
10029d007619Sjeremylt 
10039d007619Sjeremylt /**
10049d007619Sjeremylt   @brief Get gradient matrix of a CeedBasis
10059d007619Sjeremylt 
10069d007619Sjeremylt   @param basis      CeedBasis
10079d007619Sjeremylt   @param[out] grad  Variable to store gradient matrix
10089d007619Sjeremylt 
10099d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
10109d007619Sjeremylt 
1011b7c9bbdaSJeremy L Thompson   @ref Advanced
10129d007619Sjeremylt **/
10136c58de82SJeremy L Thompson int CeedBasisGetGrad(CeedBasis basis, const CeedScalar **grad) {
1014d1d35e2fSjeremylt   if (!basis->grad && basis->tensor_basis) {
10159d007619Sjeremylt     // Allocate
10169d007619Sjeremylt     int ierr;
10179d007619Sjeremylt     ierr = CeedMalloc(basis->dim*basis->Q*basis->P, &basis->grad);
10189d007619Sjeremylt     CeedChk(ierr);
10199d007619Sjeremylt 
10209d007619Sjeremylt     // Initialize
10219d007619Sjeremylt     for (CeedInt i=0; i<basis->dim*basis->Q*basis->P; i++)
10229d007619Sjeremylt       basis->grad[i] = 1.0;
10239d007619Sjeremylt 
10249d007619Sjeremylt     // Calculate
10259d007619Sjeremylt     for (CeedInt d=0; d<basis->dim; d++)
10269d007619Sjeremylt       for (CeedInt i=0; i<basis->dim; i++)
10279d007619Sjeremylt         for (CeedInt qpt=0; qpt<basis->Q; qpt++)
10289d007619Sjeremylt           for (CeedInt node=0; node<basis->P; node++) {
1029d1d35e2fSjeremylt             CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d;
1030d1d35e2fSjeremylt             CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d;
10319d007619Sjeremylt             if (i == d)
10329d007619Sjeremylt               basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *=
1033d1d35e2fSjeremylt                 basis->grad_1d[q*basis->P_1d+p];
10349d007619Sjeremylt             else
10359d007619Sjeremylt               basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *=
1036d1d35e2fSjeremylt                 basis->interp_1d[q*basis->P_1d+p];
10379d007619Sjeremylt           }
10389d007619Sjeremylt   }
10399d007619Sjeremylt   *grad = basis->grad;
1040e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
10419d007619Sjeremylt }
10429d007619Sjeremylt 
10439d007619Sjeremylt /**
10449d007619Sjeremylt   @brief Get 1D gradient matrix of a tensor product CeedBasis
10459d007619Sjeremylt 
10469d007619Sjeremylt   @param basis         CeedBasis
1047d1d35e2fSjeremylt   @param[out] grad_1d  Variable to store gradient matrix
10489d007619Sjeremylt 
10499d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
10509d007619Sjeremylt 
1051b7c9bbdaSJeremy L Thompson   @ref Advanced
10529d007619Sjeremylt **/
1053d1d35e2fSjeremylt int CeedBasisGetGrad1D(CeedBasis basis, const CeedScalar **grad_1d) {
1054d1d35e2fSjeremylt   if (!basis->tensor_basis)
10559d007619Sjeremylt     // LCOV_EXCL_START
1056e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR,
1057e15f9bd0SJeremy L Thompson                      "CeedBasis is not a tensor product basis.");
10589d007619Sjeremylt   // LCOV_EXCL_STOP
10599d007619Sjeremylt 
1060d1d35e2fSjeremylt   *grad_1d = basis->grad_1d;
1061e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
10629d007619Sjeremylt }
10639d007619Sjeremylt 
10649d007619Sjeremylt /**
10657a982d89SJeremy L. Thompson   @brief Destroy a CeedBasis
10667a982d89SJeremy L. Thompson 
10677a982d89SJeremy L. Thompson   @param basis CeedBasis to destroy
10687a982d89SJeremy L. Thompson 
10697a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
10707a982d89SJeremy L. Thompson 
10717a982d89SJeremy L. Thompson   @ref User
10727a982d89SJeremy L. Thompson **/
10737a982d89SJeremy L. Thompson int CeedBasisDestroy(CeedBasis *basis) {
10747a982d89SJeremy L. Thompson   int ierr;
10757a982d89SJeremy L. Thompson 
1076d1d35e2fSjeremylt   if (!*basis || --(*basis)->ref_count > 0) return CEED_ERROR_SUCCESS;
10777a982d89SJeremy L. Thompson   if ((*basis)->Destroy) {
10787a982d89SJeremy L. Thompson     ierr = (*basis)->Destroy(*basis); CeedChk(ierr);
10797a982d89SJeremy L. Thompson   }
108034359f16Sjeremylt   if ((*basis)->contract) {
108134359f16Sjeremylt     ierr = CeedTensorContractDestroy(&(*basis)->contract); CeedChk(ierr);
108234359f16Sjeremylt   }
10837a982d89SJeremy L. Thompson   ierr = CeedFree(&(*basis)->interp); CeedChk(ierr);
1084d1d35e2fSjeremylt   ierr = CeedFree(&(*basis)->interp_1d); CeedChk(ierr);
10857a982d89SJeremy L. Thompson   ierr = CeedFree(&(*basis)->grad); CeedChk(ierr);
1086d1d35e2fSjeremylt   ierr = CeedFree(&(*basis)->grad_1d); CeedChk(ierr);
1087d1d35e2fSjeremylt   ierr = CeedFree(&(*basis)->q_ref_1d); CeedChk(ierr);
1088d1d35e2fSjeremylt   ierr = CeedFree(&(*basis)->q_weight_1d); CeedChk(ierr);
10897a982d89SJeremy L. Thompson   ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr);
10907a982d89SJeremy L. Thompson   ierr = CeedFree(basis); CeedChk(ierr);
1091e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
10927a982d89SJeremy L. Thompson }
10937a982d89SJeremy L. Thompson 
10947a982d89SJeremy L. Thompson /**
1095b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre quadrature
1096b11c1e72Sjeremylt 
1097b11c1e72Sjeremylt   @param Q                 Number of quadrature points (integrates polynomials of
1098b11c1e72Sjeremylt                              degree 2*Q-1 exactly)
1099d1d35e2fSjeremylt   @param[out] q_ref_1d     Array of length Q to hold the abscissa on [-1, 1]
1100d1d35e2fSjeremylt   @param[out] q_weight_1d  Array of length Q to hold the weights
1101b11c1e72Sjeremylt 
1102b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1103dfdf5a53Sjeremylt 
1104dfdf5a53Sjeremylt   @ref Utility
1105b11c1e72Sjeremylt **/
1106d1d35e2fSjeremylt int CeedGaussQuadrature(CeedInt Q, CeedScalar *q_ref_1d,
1107d1d35e2fSjeremylt                         CeedScalar *q_weight_1d) {
1108d7b241e6Sjeremylt   // Allocate
1109d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0*atan(1.0);
1110d1d35e2fSjeremylt   // Build q_ref_1d, q_weight_1d
1111d7b241e6Sjeremylt   for (int i = 0; i <= Q/2; i++) {
1112d7b241e6Sjeremylt     // Guess
1113d7b241e6Sjeremylt     xi = cos(PI*(CeedScalar)(2*i+1)/((CeedScalar)(2*Q)));
1114d7b241e6Sjeremylt     // Pn(xi)
1115d7b241e6Sjeremylt     P0 = 1.0;
1116d7b241e6Sjeremylt     P1 = xi;
1117d7b241e6Sjeremylt     P2 = 0.0;
1118d7b241e6Sjeremylt     for (int j = 2; j <= Q; j++) {
1119d7b241e6Sjeremylt       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1120d7b241e6Sjeremylt       P0 = P1;
1121d7b241e6Sjeremylt       P1 = P2;
1122d7b241e6Sjeremylt     }
1123d7b241e6Sjeremylt     // First Newton Step
1124d7b241e6Sjeremylt     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1125d7b241e6Sjeremylt     xi = xi-P2/dP2;
1126d7b241e6Sjeremylt     // Newton to convergence
11270e4d4210Sjeremylt     for (int k=0; k<100 && fabs(P2)>10*CEED_EPSILON; k++) {
1128d7b241e6Sjeremylt       P0 = 1.0;
1129d7b241e6Sjeremylt       P1 = xi;
1130d7b241e6Sjeremylt       for (int j = 2; j <= Q; j++) {
1131d7b241e6Sjeremylt         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1132d7b241e6Sjeremylt         P0 = P1;
1133d7b241e6Sjeremylt         P1 = P2;
1134d7b241e6Sjeremylt       }
1135d7b241e6Sjeremylt       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1136d7b241e6Sjeremylt       xi = xi-P2/dP2;
1137d7b241e6Sjeremylt     }
1138d7b241e6Sjeremylt     // Save xi, wi
1139d7b241e6Sjeremylt     wi = 2.0/((1.0-xi*xi)*dP2*dP2);
1140d1d35e2fSjeremylt     q_weight_1d[i] = wi;
1141d1d35e2fSjeremylt     q_weight_1d[Q-1-i] = wi;
1142d1d35e2fSjeremylt     q_ref_1d[i] = -xi;
1143d1d35e2fSjeremylt     q_ref_1d[Q-1-i]= xi;
1144d7b241e6Sjeremylt   }
1145e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1146d7b241e6Sjeremylt }
1147d7b241e6Sjeremylt 
1148b11c1e72Sjeremylt /**
1149b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre-Lobatto quadrature
1150b11c1e72Sjeremylt 
1151b11c1e72Sjeremylt   @param Q                 Number of quadrature points (integrates polynomials of
1152b11c1e72Sjeremylt                              degree 2*Q-3 exactly)
1153d1d35e2fSjeremylt   @param[out] q_ref_1d     Array of length Q to hold the abscissa on [-1, 1]
1154d1d35e2fSjeremylt   @param[out] q_weight_1d  Array of length Q to hold the weights
1155b11c1e72Sjeremylt 
1156b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1157dfdf5a53Sjeremylt 
1158dfdf5a53Sjeremylt   @ref Utility
1159b11c1e72Sjeremylt **/
1160d1d35e2fSjeremylt int CeedLobattoQuadrature(CeedInt Q, CeedScalar *q_ref_1d,
1161d1d35e2fSjeremylt                           CeedScalar *q_weight_1d) {
1162d7b241e6Sjeremylt   // Allocate
1163d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0*atan(1.0);
1164d1d35e2fSjeremylt   // Build q_ref_1d, q_weight_1d
1165d7b241e6Sjeremylt   // Set endpoints
116630a100c3SJed Brown   if (Q < 2)
1167b0d62198Sjeremylt     // LCOV_EXCL_START
1168e15f9bd0SJeremy L Thompson     return CeedError(NULL, CEED_ERROR_DIMENSION,
11697ed52d01Sjeremylt                      "Cannot create Lobatto quadrature with Q=%d < 2 points", Q);
1170b0d62198Sjeremylt   // LCOV_EXCL_STOP
1171d7b241e6Sjeremylt   wi = 2.0/((CeedScalar)(Q*(Q-1)));
1172d1d35e2fSjeremylt   if (q_weight_1d) {
1173d1d35e2fSjeremylt     q_weight_1d[0] = wi;
1174d1d35e2fSjeremylt     q_weight_1d[Q-1] = wi;
1175d7b241e6Sjeremylt   }
1176d1d35e2fSjeremylt   q_ref_1d[0] = -1.0;
1177d1d35e2fSjeremylt   q_ref_1d[Q-1] = 1.0;
1178d7b241e6Sjeremylt   // Interior
1179d7b241e6Sjeremylt   for (int i = 1; i <= (Q-1)/2; i++) {
1180d7b241e6Sjeremylt     // Guess
1181d7b241e6Sjeremylt     xi = cos(PI*(CeedScalar)(i)/(CeedScalar)(Q-1));
1182d7b241e6Sjeremylt     // Pn(xi)
1183d7b241e6Sjeremylt     P0 = 1.0;
1184d7b241e6Sjeremylt     P1 = xi;
1185d7b241e6Sjeremylt     P2 = 0.0;
1186d7b241e6Sjeremylt     for (int j = 2; j < Q; j++) {
1187d7b241e6Sjeremylt       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1188d7b241e6Sjeremylt       P0 = P1;
1189d7b241e6Sjeremylt       P1 = P2;
1190d7b241e6Sjeremylt     }
1191d7b241e6Sjeremylt     // First Newton step
1192d7b241e6Sjeremylt     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1193d7b241e6Sjeremylt     d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
1194d7b241e6Sjeremylt     xi = xi-dP2/d2P2;
1195d7b241e6Sjeremylt     // Newton to convergence
11960e4d4210Sjeremylt     for (int k=0; k<100 && fabs(dP2)>10*CEED_EPSILON; k++) {
1197d7b241e6Sjeremylt       P0 = 1.0;
1198d7b241e6Sjeremylt       P1 = xi;
1199d7b241e6Sjeremylt       for (int j = 2; j < Q; j++) {
1200d7b241e6Sjeremylt         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1201d7b241e6Sjeremylt         P0 = P1;
1202d7b241e6Sjeremylt         P1 = P2;
1203d7b241e6Sjeremylt       }
1204d7b241e6Sjeremylt       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1205d7b241e6Sjeremylt       d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
1206d7b241e6Sjeremylt       xi = xi-dP2/d2P2;
1207d7b241e6Sjeremylt     }
1208d7b241e6Sjeremylt     // Save xi, wi
1209d7b241e6Sjeremylt     wi = 2.0/(((CeedScalar)(Q*(Q-1)))*P2*P2);
1210d1d35e2fSjeremylt     if (q_weight_1d) {
1211d1d35e2fSjeremylt       q_weight_1d[i] = wi;
1212d1d35e2fSjeremylt       q_weight_1d[Q-1-i] = wi;
1213d7b241e6Sjeremylt     }
1214d1d35e2fSjeremylt     q_ref_1d[i] = -xi;
1215d1d35e2fSjeremylt     q_ref_1d[Q-1-i]= xi;
1216d7b241e6Sjeremylt   }
1217e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1218d7b241e6Sjeremylt }
1219d7b241e6Sjeremylt 
1220dfdf5a53Sjeremylt /**
122195bb1877Svaleriabarra   @brief Return QR Factorization of a matrix
1222b11c1e72Sjeremylt 
122377645d7bSjeremylt   @param ceed         A Ceed context for error handling
122452bfb9bbSJeremy L Thompson   @param[in,out] mat  Row-major matrix to be factorized in place
122552bfb9bbSJeremy L Thompson   @param[in,out] tau  Vector of length m of scaling factors
1226b11c1e72Sjeremylt   @param m            Number of rows
1227b11c1e72Sjeremylt   @param n            Number of columns
1228b11c1e72Sjeremylt 
1229b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1230dfdf5a53Sjeremylt 
1231dfdf5a53Sjeremylt   @ref Utility
1232b11c1e72Sjeremylt **/
1233a7bd39daSjeremylt int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau,
1234d7b241e6Sjeremylt                         CeedInt m, CeedInt n) {
1235d7b241e6Sjeremylt   CeedScalar v[m];
1236d7b241e6Sjeremylt 
1237a7bd39daSjeremylt   // Check m >= n
1238a7bd39daSjeremylt   if (n > m)
1239c042f62fSJeremy L Thompson     // LCOV_EXCL_START
1240e15f9bd0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_UNSUPPORTED,
1241e15f9bd0SJeremy L Thompson                      "Cannot compute QR factorization with n > m");
1242c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
1243a7bd39daSjeremylt 
124452bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++) {
1245bde37e8cSJed Brown     if (i >= m-1) { // last row of matrix, no reflection needed
1246bde37e8cSJed Brown       tau[i] = 0.;
1247bde37e8cSJed Brown       break;
1248bde37e8cSJed Brown     }
1249d7b241e6Sjeremylt     // Calculate Householder vector, magnitude
1250d7b241e6Sjeremylt     CeedScalar sigma = 0.0;
1251d7b241e6Sjeremylt     v[i] = mat[i+n*i];
125252bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<m; j++) {
1253d7b241e6Sjeremylt       v[j] = mat[i+n*j];
1254d7b241e6Sjeremylt       sigma += v[j] * v[j];
1255d7b241e6Sjeremylt     }
1256d7b241e6Sjeremylt     CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m]
1257673160d7Sjeremylt     CeedScalar R_ii = -copysign(norm, v[i]);
1258673160d7Sjeremylt     v[i] -= R_ii;
1259d7b241e6Sjeremylt     // norm of v[i:m] after modification above and scaling below
1260d7b241e6Sjeremylt     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
1261d7b241e6Sjeremylt     //   tau = 2 / (norm*norm)
1262d7b241e6Sjeremylt     tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma);
12631d102b48SJeremy L Thompson     for (CeedInt j=i+1; j<m; j++)
12641d102b48SJeremy L Thompson       v[j] /= v[i];
1265d7b241e6Sjeremylt 
1266d7b241e6Sjeremylt     // Apply Householder reflector to lower right panel
1267d7b241e6Sjeremylt     CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1);
1268d7b241e6Sjeremylt     // Save v
1269673160d7Sjeremylt     mat[i+n*i] = R_ii;
12701d102b48SJeremy L Thompson     for (CeedInt j=i+1; j<m; j++)
1271d7b241e6Sjeremylt       mat[i+n*j] = v[j];
1272d7b241e6Sjeremylt   }
1273e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1274d7b241e6Sjeremylt }
1275d7b241e6Sjeremylt 
1276b11c1e72Sjeremylt /**
127752bfb9bbSJeremy L Thompson   @brief Return symmetric Schur decomposition of the symmetric matrix mat via
127852bfb9bbSJeremy L Thompson            symmetric QR factorization
127952bfb9bbSJeremy L Thompson 
128077645d7bSjeremylt   @param ceed         A Ceed context for error handling
128152bfb9bbSJeremy L Thompson   @param[in,out] mat  Row-major matrix to be factorized in place
1282460bf743SValeria Barra   @param[out] lambda  Vector of length n of eigenvalues
128352bfb9bbSJeremy L Thompson   @param n            Number of rows/columns
128452bfb9bbSJeremy L Thompson 
128552bfb9bbSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
128652bfb9bbSJeremy L Thompson 
128752bfb9bbSJeremy L Thompson   @ref Utility
128852bfb9bbSJeremy L Thompson **/
128903d18186Sjeremylt CeedPragmaOptimizeOff
129052bfb9bbSJeremy L Thompson int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat,
129152bfb9bbSJeremy L Thompson                                     CeedScalar *lambda, CeedInt n) {
129252bfb9bbSJeremy L Thompson   // Check bounds for clang-tidy
129352bfb9bbSJeremy L Thompson   if (n<2)
1294c042f62fSJeremy L Thompson     // LCOV_EXCL_START
1295e15f9bd0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_UNSUPPORTED,
1296c042f62fSJeremy L Thompson                      "Cannot compute symmetric Schur decomposition of scalars");
1297c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
129852bfb9bbSJeremy L Thompson 
1299673160d7Sjeremylt   CeedScalar v[n-1], tau[n-1], mat_T[n*n];
130052bfb9bbSJeremy L Thompson 
1301673160d7Sjeremylt   // Copy mat to mat_T and set mat to I
1302673160d7Sjeremylt   memcpy(mat_T, mat, n*n*sizeof(mat[0]));
130352bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
130452bfb9bbSJeremy L Thompson     for (CeedInt j=0; j<n; j++)
130552bfb9bbSJeremy L Thompson       mat[j+n*i] = (i==j) ? 1 : 0;
130652bfb9bbSJeremy L Thompson 
130752bfb9bbSJeremy L Thompson   // Reduce to tridiagonal
130852bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n-1; i++) {
130952bfb9bbSJeremy L Thompson     // Calculate Householder vector, magnitude
131052bfb9bbSJeremy L Thompson     CeedScalar sigma = 0.0;
1311673160d7Sjeremylt     v[i] = mat_T[i+n*(i+1)];
131252bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) {
1313673160d7Sjeremylt       v[j] = mat_T[i+n*(j+1)];
131452bfb9bbSJeremy L Thompson       sigma += v[j] * v[j];
131552bfb9bbSJeremy L Thompson     }
131652bfb9bbSJeremy L Thompson     CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:n-1]
1317673160d7Sjeremylt     CeedScalar R_ii = -copysign(norm, v[i]);
1318673160d7Sjeremylt     v[i] -= R_ii;
131952bfb9bbSJeremy L Thompson     // norm of v[i:m] after modification above and scaling below
132052bfb9bbSJeremy L Thompson     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
132152bfb9bbSJeremy L Thompson     //   tau = 2 / (norm*norm)
132280a9ef05SNatalie Beams     tau[i] = i == n - 2 ? 2 : 2 * v[i]*v[i] / (v[i]*v[i] + sigma);
1323fb551037Sjeremylt     for (CeedInt j=i+1; j<n-1; j++)
1324fb551037Sjeremylt       v[j] /= v[i];
132552bfb9bbSJeremy L Thompson 
132652bfb9bbSJeremy L Thompson     // Update sub and super diagonal
132752bfb9bbSJeremy L Thompson     for (CeedInt j=i+2; j<n; j++) {
1328673160d7Sjeremylt       mat_T[i+n*j] = 0; mat_T[j+n*i] = 0;
132952bfb9bbSJeremy L Thompson     }
133052bfb9bbSJeremy L Thompson     // Apply symmetric Householder reflector to lower right panel
1331673160d7Sjeremylt     CeedHouseholderReflect(&mat_T[(i+1)+n*(i+1)], &v[i], tau[i],
133252bfb9bbSJeremy L Thompson                            n-(i+1), n-(i+1), n, 1);
1333673160d7Sjeremylt     CeedHouseholderReflect(&mat_T[(i+1)+n*(i+1)], &v[i], tau[i],
133452bfb9bbSJeremy L Thompson                            n-(i+1), n-(i+1), 1, n);
1335673160d7Sjeremylt 
133652bfb9bbSJeremy L Thompson     // Save v
1337673160d7Sjeremylt     mat_T[i+n*(i+1)] = R_ii;
1338673160d7Sjeremylt     mat_T[(i+1)+n*i] = R_ii;
133952bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) {
1340673160d7Sjeremylt       mat_T[i+n*(j+1)] = v[j];
134152bfb9bbSJeremy L Thompson     }
134252bfb9bbSJeremy L Thompson   }
134352bfb9bbSJeremy L Thompson   // Backwards accumulation of Q
134452bfb9bbSJeremy L Thompson   for (CeedInt i=n-2; i>=0; i--) {
134585cf89eaSjeremylt     if (tau[i] > 0.0) {
134652bfb9bbSJeremy L Thompson       v[i] = 1;
134752bfb9bbSJeremy L Thompson       for (CeedInt j=i+1; j<n-1; j++) {
1348673160d7Sjeremylt         v[j] = mat_T[i+n*(j+1)];
1349673160d7Sjeremylt         mat_T[i+n*(j+1)] = 0;
135052bfb9bbSJeremy L Thompson       }
135152bfb9bbSJeremy L Thompson       CeedHouseholderReflect(&mat[(i+1)+n*(i+1)], &v[i], tau[i],
135252bfb9bbSJeremy L Thompson                              n-(i+1), n-(i+1), n, 1);
135352bfb9bbSJeremy L Thompson     }
135485cf89eaSjeremylt   }
135552bfb9bbSJeremy L Thompson 
135652bfb9bbSJeremy L Thompson   // Reduce sub and super diagonal
1357673160d7Sjeremylt   CeedInt p = 0, q = 0, itr = 0, max_itr = n*n*n*n;
1358673160d7Sjeremylt   CeedScalar tol = CEED_EPSILON;
135952bfb9bbSJeremy L Thompson 
1360673160d7Sjeremylt   while (itr < max_itr) {
136152bfb9bbSJeremy L Thompson     // Update p, q, size of reduced portions of diagonal
136252bfb9bbSJeremy L Thompson     p = 0; q = 0;
136352bfb9bbSJeremy L Thompson     for (CeedInt i=n-2; i>=0; i--) {
1364673160d7Sjeremylt       if (fabs(mat_T[i+n*(i+1)]) < tol)
136552bfb9bbSJeremy L Thompson         q += 1;
136652bfb9bbSJeremy L Thompson       else
136752bfb9bbSJeremy L Thompson         break;
136852bfb9bbSJeremy L Thompson     }
1369673160d7Sjeremylt     for (CeedInt i=0; i<n-q-1; i++) {
1370673160d7Sjeremylt       if (fabs(mat_T[i+n*(i+1)]) < tol)
137152bfb9bbSJeremy L Thompson         p += 1;
137252bfb9bbSJeremy L Thompson       else
137352bfb9bbSJeremy L Thompson         break;
137452bfb9bbSJeremy L Thompson     }
137552bfb9bbSJeremy L Thompson     if (q == n-1) break; // Finished reducing
137652bfb9bbSJeremy L Thompson 
137752bfb9bbSJeremy L Thompson     // Reduce tridiagonal portion
1378673160d7Sjeremylt     CeedScalar t_nn = mat_T[(n-1-q)+n*(n-1-q)],
1379673160d7Sjeremylt                t_nnm1 = mat_T[(n-2-q)+n*(n-1-q)];
1380673160d7Sjeremylt     CeedScalar d = (mat_T[(n-2-q)+n*(n-2-q)] - t_nn)/2;
1381673160d7Sjeremylt     CeedScalar mu = t_nn - t_nnm1*t_nnm1 /
1382673160d7Sjeremylt                     (d + copysign(sqrt(d*d + t_nnm1*t_nnm1), d));
1383673160d7Sjeremylt     CeedScalar x = mat_T[p+n*p] - mu;
1384673160d7Sjeremylt     CeedScalar z = mat_T[p+n*(p+1)];
1385673160d7Sjeremylt     for (CeedInt k=p; k<n-q-1; k++) {
138652bfb9bbSJeremy L Thompson       // Compute Givens rotation
138752bfb9bbSJeremy L Thompson       CeedScalar c = 1, s = 0;
138852bfb9bbSJeremy L Thompson       if (fabs(z) > tol) {
138952bfb9bbSJeremy L Thompson         if (fabs(z) > fabs(x)) {
139052bfb9bbSJeremy L Thompson           CeedScalar tau = -x/z;
139152bfb9bbSJeremy L Thompson           s = 1/sqrt(1+tau*tau), c = s*tau;
139252bfb9bbSJeremy L Thompson         } else {
139352bfb9bbSJeremy L Thompson           CeedScalar tau = -z/x;
139452bfb9bbSJeremy L Thompson           c = 1/sqrt(1+tau*tau), s = c*tau;
139552bfb9bbSJeremy L Thompson         }
139652bfb9bbSJeremy L Thompson       }
139752bfb9bbSJeremy L Thompson 
139852bfb9bbSJeremy L Thompson       // Apply Givens rotation to T
1399673160d7Sjeremylt       CeedGivensRotation(mat_T, c, s, CEED_NOTRANSPOSE, k, k+1, n, n);
1400673160d7Sjeremylt       CeedGivensRotation(mat_T, c, s, CEED_TRANSPOSE, k, k+1, n, n);
140152bfb9bbSJeremy L Thompson 
140252bfb9bbSJeremy L Thompson       // Apply Givens rotation to Q
140352bfb9bbSJeremy L Thompson       CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k+1, n, n);
140452bfb9bbSJeremy L Thompson 
140552bfb9bbSJeremy L Thompson       // Update x, z
140652bfb9bbSJeremy L Thompson       if (k < n-q-2) {
1407673160d7Sjeremylt         x = mat_T[k+n*(k+1)];
1408673160d7Sjeremylt         z = mat_T[k+n*(k+2)];
140952bfb9bbSJeremy L Thompson       }
141052bfb9bbSJeremy L Thompson     }
141152bfb9bbSJeremy L Thompson     itr++;
141252bfb9bbSJeremy L Thompson   }
1413673160d7Sjeremylt 
141452bfb9bbSJeremy L Thompson   // Save eigenvalues
141552bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
1416673160d7Sjeremylt     lambda[i] = mat_T[i+n*i];
141752bfb9bbSJeremy L Thompson 
141852bfb9bbSJeremy L Thompson   // Check convergence
1419673160d7Sjeremylt   if (itr == max_itr && q < n-1)
1420c042f62fSJeremy L Thompson     // LCOV_EXCL_START
1421e15f9bd0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_MINOR,
1422e15f9bd0SJeremy L Thompson                      "Symmetric QR failed to converge");
1423c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
1424e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
142552bfb9bbSJeremy L Thompson }
142603d18186Sjeremylt CeedPragmaOptimizeOn
142752bfb9bbSJeremy L Thompson 
142852bfb9bbSJeremy L Thompson /**
142952bfb9bbSJeremy L Thompson   @brief Return Simultaneous Diagonalization of two matrices. This solves the
143052bfb9bbSJeremy L Thompson            generalized eigenvalue problem A x = lambda B x, where A and B
143152bfb9bbSJeremy L Thompson            are symmetric and B is positive definite. We generate the matrix X
143252bfb9bbSJeremy L Thompson            and vector Lambda such that X^T A X = Lambda and X^T B X = I. This
143352bfb9bbSJeremy L Thompson            is equivalent to the LAPACK routine 'sygv' with TYPE = 1.
143452bfb9bbSJeremy L Thompson 
143577645d7bSjeremylt   @param ceed         A Ceed context for error handling
1436d1d35e2fSjeremylt   @param[in] mat_A    Row-major matrix to be factorized with eigenvalues
1437d1d35e2fSjeremylt   @param[in] mat_B    Row-major matrix to be factorized to identity
1438d3331725Sjeremylt   @param[out] mat_X   Row-major orthogonal matrix
1439460bf743SValeria Barra   @param[out] lambda  Vector of length n of generalized eigenvalues
144052bfb9bbSJeremy L Thompson   @param n            Number of rows/columns
144152bfb9bbSJeremy L Thompson 
144252bfb9bbSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
144352bfb9bbSJeremy L Thompson 
144452bfb9bbSJeremy L Thompson   @ref Utility
144552bfb9bbSJeremy L Thompson **/
144603d18186Sjeremylt CeedPragmaOptimizeOff
1447d1d35e2fSjeremylt int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *mat_A,
1448d3331725Sjeremylt                                     CeedScalar *mat_B, CeedScalar *mat_X,
144952bfb9bbSJeremy L Thompson                                     CeedScalar *lambda, CeedInt n) {
145052bfb9bbSJeremy L Thompson   int ierr;
1451d3331725Sjeremylt   CeedScalar *mat_C, *mat_G, *vec_D;
145278464608Sjeremylt   ierr = CeedCalloc(n*n, &mat_C); CeedChk(ierr);
145378464608Sjeremylt   ierr = CeedCalloc(n*n, &mat_G); CeedChk(ierr);
1454d3331725Sjeremylt   ierr = CeedCalloc(n, &vec_D); CeedChk(ierr);
145552bfb9bbSJeremy L Thompson 
145652bfb9bbSJeremy L Thompson   // Compute B = G D G^T
145778464608Sjeremylt   memcpy(mat_G, mat_B, n*n*sizeof(mat_B[0]));
1458d3331725Sjeremylt   ierr = CeedSymmetricSchurDecomposition(ceed, mat_G, vec_D, n); CeedChk(ierr);
145952bfb9bbSJeremy L Thompson 
146085cf89eaSjeremylt   // Sort eigenvalues
146185cf89eaSjeremylt   for (CeedInt i=n-1; i>=0; i--)
146285cf89eaSjeremylt     for (CeedInt j=0; j<i; j++) {
146385cf89eaSjeremylt       if (fabs(vec_D[j]) > fabs(vec_D[j+1])) {
146485cf89eaSjeremylt         CeedScalar temp;
146585cf89eaSjeremylt         temp = vec_D[j]; vec_D[j] = vec_D[j+1]; vec_D[j+1] = temp;
146685cf89eaSjeremylt         for (CeedInt k=0; k<n; k++) {
146785cf89eaSjeremylt           temp = mat_G[k*n+j]; mat_G[k*n+j] = mat_G[k*n+j+1]; mat_G[k*n+j+1] = temp;
146885cf89eaSjeremylt         }
146985cf89eaSjeremylt       }
147085cf89eaSjeremylt     }
147185cf89eaSjeremylt 
1472fb551037Sjeremylt   // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T
1473fb551037Sjeremylt   //           = D^-1/2 G^T A G D^-1/2
1474d3331725Sjeremylt   // -- D = D^-1/2
147552bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
1476d3331725Sjeremylt     vec_D[i] = 1./sqrt(vec_D[i]);
1477d3331725Sjeremylt   // -- G = G D^-1/2
1478d3331725Sjeremylt   // -- C = D^-1/2 G^T
1479d3331725Sjeremylt   for (CeedInt i=0; i<n; i++)
1480d3331725Sjeremylt     for (CeedInt j=0; j<n; j++) {
1481673160d7Sjeremylt       mat_G[i*n+j] *= vec_D[j];
1482673160d7Sjeremylt       mat_C[j*n+i]  = mat_G[i*n+j];
1483d3331725Sjeremylt     }
1484673160d7Sjeremylt   // -- X = (D^-1/2 G^T) A
1485d1d35e2fSjeremylt   ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)mat_C,
1486d3331725Sjeremylt                             (const CeedScalar *)mat_A, mat_X, n, n, n);
14879289e5bfSjeremylt   CeedChk(ierr);
1488673160d7Sjeremylt   // -- C = (D^-1/2 G^T A) (G D^-1/2)
1489d3331725Sjeremylt   ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)mat_X,
149078464608Sjeremylt                             (const CeedScalar *)mat_G, mat_C, n, n, n);
14919289e5bfSjeremylt   CeedChk(ierr);
149252bfb9bbSJeremy L Thompson 
149352bfb9bbSJeremy L Thompson   // Compute Q^T C Q = lambda
1494d1d35e2fSjeremylt   ierr = CeedSymmetricSchurDecomposition(ceed, mat_C, lambda, n); CeedChk(ierr);
149552bfb9bbSJeremy L Thompson 
149685cf89eaSjeremylt   // Sort eigenvalues
149785cf89eaSjeremylt   for (CeedInt i=n-1; i>=0; i--)
149885cf89eaSjeremylt     for (CeedInt j=0; j<i; j++) {
149985cf89eaSjeremylt       if (fabs(lambda[j]) > fabs(lambda[j+1])) {
150085cf89eaSjeremylt         CeedScalar temp;
150185cf89eaSjeremylt         temp = lambda[j]; lambda[j] = lambda[j+1]; lambda[j+1] = temp;
150285cf89eaSjeremylt         for (CeedInt k=0; k<n; k++) {
150385cf89eaSjeremylt           temp = mat_C[k*n+j]; mat_C[k*n+j] = mat_C[k*n+j+1]; mat_C[k*n+j+1] = temp;
150485cf89eaSjeremylt         }
150585cf89eaSjeremylt       }
150685cf89eaSjeremylt     }
150785cf89eaSjeremylt 
1508d3331725Sjeremylt   // Set X = (G D^1/2)^-T Q
1509fb551037Sjeremylt   //       = G D^-1/2 Q
151078464608Sjeremylt   ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)mat_G,
1511d3331725Sjeremylt                             (const CeedScalar *)mat_C, mat_X, n, n, n);
15129289e5bfSjeremylt   CeedChk(ierr);
151378464608Sjeremylt 
151478464608Sjeremylt   // Cleanup
151578464608Sjeremylt   ierr = CeedFree(&mat_C); CeedChk(ierr);
151678464608Sjeremylt   ierr = CeedFree(&mat_G); CeedChk(ierr);
1517d3331725Sjeremylt   ierr = CeedFree(&vec_D); CeedChk(ierr);
1518e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
151952bfb9bbSJeremy L Thompson }
152003d18186Sjeremylt CeedPragmaOptimizeOn
152152bfb9bbSJeremy L Thompson 
1522d7b241e6Sjeremylt /// @}
1523