xref: /libCEED/interface/ceed-basis.c (revision e15f9bd09af0280c89b79924fa9af7dd2e3e30be)
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 
173d576824SJeremy L Thompson #include <ceed.h>
18d863ab9bSjeremylt #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   }
75*e15f9bd0SJeremy 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
867a982d89SJeremy L. Thompson   @param tmode      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,
987a982d89SJeremy L. Thompson                           const CeedScalar *tau, CeedTransposeMode tmode,
997a982d89SJeremy L. Thompson                           CeedInt m, CeedInt n, CeedInt k,
1007a982d89SJeremy L. Thompson                           CeedInt row, CeedInt col) {
101*e15f9bd0SJeremy L Thompson   int ierr;
1027a982d89SJeremy L. Thompson   CeedScalar v[m];
1037a982d89SJeremy L. Thompson   for (CeedInt ii=0; ii<k; ii++) {
1047a982d89SJeremy L. Thompson     CeedInt i = tmode == CEED_TRANSPOSE ? ii : k-1-ii;
1057a982d89SJeremy L. Thompson     for (CeedInt j=i+1; j<m; j++)
1067a982d89SJeremy L. Thompson       v[j] = Q[j*k+i];
1077a982d89SJeremy L. Thompson     // Apply Householder reflector (I - tau v v^T) collograd1d^T
108*e15f9bd0SJeremy L Thompson     ierr = CeedHouseholderReflect(&A[i*row], &v[i], tau[i], m-i, n, row, col);
109*e15f9bd0SJeremy L Thompson     CeedChk(ierr);
1107a982d89SJeremy L. Thompson   }
111*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1127a982d89SJeremy L. Thompson }
1137a982d89SJeremy L. Thompson 
1147a982d89SJeremy L. Thompson /**
1157a982d89SJeremy L. Thompson   @brief Compute Givens rotation
1167a982d89SJeremy L. Thompson 
1177a982d89SJeremy L. Thompson     Computes A = G A (or G^T A in transpose mode)
1187a982d89SJeremy L. Thompson     where A is an mxn matrix indexed as A[i*n + j*m]
1197a982d89SJeremy L. Thompson 
1207a982d89SJeremy L. Thompson   @param[in,out] A  Row major matrix to apply Givens rotation to, in place
1217a982d89SJeremy L. Thompson   @param c          Cosine factor
1227a982d89SJeremy L. Thompson   @param s          Sine factor
1234cc79fe7SJed Brown   @param tmode      @ref CEED_NOTRANSPOSE to rotate the basis counter-clockwise,
1244c4400c7SValeria Barra                     which has the effect of rotating columns of A clockwise;
1254cc79fe7SJed Brown                     @ref CEED_TRANSPOSE for the opposite rotation
1267a982d89SJeremy L. Thompson   @param i          First row/column to apply rotation
1277a982d89SJeremy L. Thompson   @param k          Second row/column to apply rotation
1287a982d89SJeremy L. Thompson   @param m          Number of rows in A
1297a982d89SJeremy L. Thompson   @param n          Number of columns in A
1307a982d89SJeremy L. Thompson 
1317a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
1327a982d89SJeremy L. Thompson 
1337a982d89SJeremy L. Thompson   @ref Developer
1347a982d89SJeremy L. Thompson **/
1357a982d89SJeremy L. Thompson static int CeedGivensRotation(CeedScalar *A, CeedScalar c, CeedScalar s,
1367a982d89SJeremy L. Thompson                               CeedTransposeMode tmode, CeedInt i, CeedInt k,
1377a982d89SJeremy L. Thompson                               CeedInt m, CeedInt n) {
1387a982d89SJeremy L. Thompson   CeedInt stridej = 1, strideik = m, numits = n;
1397a982d89SJeremy L. Thompson   if (tmode == CEED_NOTRANSPOSE) {
1407a982d89SJeremy L. Thompson     stridej = n; strideik = 1; numits = m;
1417a982d89SJeremy L. Thompson   }
1427a982d89SJeremy L. Thompson 
1437a982d89SJeremy L. Thompson   // Apply rotation
1447a982d89SJeremy L. Thompson   for (CeedInt j=0; j<numits; j++) {
1457a982d89SJeremy L. Thompson     CeedScalar tau1 = A[i*strideik+j*stridej], tau2 = A[k*strideik+j*stridej];
1467a982d89SJeremy L. Thompson     A[i*strideik+j*stridej] = c*tau1 - s*tau2;
1477a982d89SJeremy L. Thompson     A[k*strideik+j*stridej] = s*tau1 + c*tau2;
1487a982d89SJeremy L. Thompson   }
149*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1507a982d89SJeremy L. Thompson }
1517a982d89SJeremy L. Thompson 
1527a982d89SJeremy L. Thompson /**
1537a982d89SJeremy L. Thompson   @brief View an array stored in a CeedBasis
1547a982d89SJeremy L. Thompson 
1550a0da059Sjeremylt   @param[in] name      Name of array
1560a0da059Sjeremylt   @param[in] fpformat  Printing format
1570a0da059Sjeremylt   @param[in] m         Number of rows in array
1580a0da059Sjeremylt   @param[in] n         Number of columns in array
1590a0da059Sjeremylt   @param[in] a         Array to be viewed
1600a0da059Sjeremylt   @param[in] stream    Stream to view to, e.g., stdout
1617a982d89SJeremy L. Thompson 
1627a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
1637a982d89SJeremy L. Thompson 
1647a982d89SJeremy L. Thompson   @ref Developer
1657a982d89SJeremy L. Thompson **/
1667a982d89SJeremy L. Thompson static int CeedScalarView(const char *name, const char *fpformat, CeedInt m,
1677a982d89SJeremy L. Thompson                           CeedInt n, const CeedScalar *a, FILE *stream) {
1687a982d89SJeremy L. Thompson   for (int i=0; i<m; i++) {
1697a982d89SJeremy L. Thompson     if (m > 1)
1707a982d89SJeremy L. Thompson       fprintf(stream, "%12s[%d]:", name, i);
1717a982d89SJeremy L. Thompson     else
1727a982d89SJeremy L. Thompson       fprintf(stream, "%12s:", name);
1737a982d89SJeremy L. Thompson     for (int j=0; j<n; j++)
1747a982d89SJeremy L. Thompson       fprintf(stream, fpformat, fabs(a[i*n+j]) > 1E-14 ? a[i*n+j] : 0);
1757a982d89SJeremy L. Thompson     fputs("\n", stream);
1767a982d89SJeremy L. Thompson   }
177*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1787a982d89SJeremy L. Thompson }
1797a982d89SJeremy L. Thompson 
1807a982d89SJeremy L. Thompson /// @}
1817a982d89SJeremy L. Thompson 
1827a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
1837a982d89SJeremy L. Thompson /// Ceed Backend API
1847a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
1857a982d89SJeremy L. Thompson /// @addtogroup CeedBasisBackend
1867a982d89SJeremy L. Thompson /// @{
1877a982d89SJeremy L. Thompson 
1887a982d89SJeremy L. Thompson /**
1897a982d89SJeremy L. Thompson   @brief Return collocated grad matrix
1907a982d89SJeremy L. Thompson 
1917a982d89SJeremy L. Thompson   @param basis             CeedBasis
1927a982d89SJeremy L. Thompson   @param[out] collograd1d  Row-major (Q1d * Q1d) matrix expressing derivatives of
1937a982d89SJeremy L. Thompson                             basis functions at quadrature points
1947a982d89SJeremy L. Thompson 
1957a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
1967a982d89SJeremy L. Thompson 
1977a982d89SJeremy L. Thompson   @ref Backend
1987a982d89SJeremy L. Thompson **/
1997a982d89SJeremy L. Thompson int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *collograd1d) {
2007a982d89SJeremy L. Thompson   int i, j, k;
2017a982d89SJeremy L. Thompson   Ceed ceed;
2027a982d89SJeremy L. Thompson   CeedInt ierr, P1d=(basis)->P1d, Q1d=(basis)->Q1d;
2037a982d89SJeremy L. Thompson   CeedScalar *interp1d, *grad1d, tau[Q1d];
2047a982d89SJeremy L. Thompson 
2057a982d89SJeremy L. Thompson   ierr = CeedMalloc(Q1d*P1d, &interp1d); CeedChk(ierr);
2067a982d89SJeremy L. Thompson   ierr = CeedMalloc(Q1d*P1d, &grad1d); CeedChk(ierr);
2077a982d89SJeremy L. Thompson   memcpy(interp1d, (basis)->interp1d, Q1d*P1d*sizeof(basis)->interp1d[0]);
2087a982d89SJeremy L. Thompson   memcpy(grad1d, (basis)->grad1d, Q1d*P1d*sizeof(basis)->interp1d[0]);
2097a982d89SJeremy L. Thompson 
2107a982d89SJeremy L. Thompson   // QR Factorization, interp1d = Q R
2117a982d89SJeremy L. Thompson   ierr = CeedBasisGetCeed(basis, &ceed); CeedChk(ierr);
2127a982d89SJeremy L. Thompson   ierr = CeedQRFactorization(ceed, interp1d, tau, Q1d, P1d); CeedChk(ierr);
213*e15f9bd0SJeremy L Thompson   // Note: This function is for backend use, so all errors are terminal
214*e15f9bd0SJeremy L Thompson   //   and we do not need to clean up memory on failure.
2157a982d89SJeremy L. Thompson 
2167a982d89SJeremy L. Thompson   // Apply Rinv, collograd1d = grad1d Rinv
2177a982d89SJeremy L. Thompson   for (i=0; i<Q1d; i++) { // Row i
2187a982d89SJeremy L. Thompson     collograd1d[Q1d*i] = grad1d[P1d*i]/interp1d[0];
2197a982d89SJeremy L. Thompson     for (j=1; j<P1d; j++) { // Column j
2207a982d89SJeremy L. Thompson       collograd1d[j+Q1d*i] = grad1d[j+P1d*i];
2217a982d89SJeremy L. Thompson       for (k=0; k<j; k++)
2227a982d89SJeremy L. Thompson         collograd1d[j+Q1d*i] -= interp1d[j+P1d*k]*collograd1d[k+Q1d*i];
2237a982d89SJeremy L. Thompson       collograd1d[j+Q1d*i] /= interp1d[j+P1d*j];
2247a982d89SJeremy L. Thompson     }
2257a982d89SJeremy L. Thompson     for (j=P1d; j<Q1d; j++)
2267a982d89SJeremy L. Thompson       collograd1d[j+Q1d*i] = 0;
2277a982d89SJeremy L. Thompson   }
2287a982d89SJeremy L. Thompson 
2297a982d89SJeremy L. Thompson   // Apply Qtranspose, collograd = collograd Qtranspose
230*e15f9bd0SJeremy L Thompson   ierr = CeedHouseholderApplyQ(collograd1d, interp1d, tau, CEED_NOTRANSPOSE,
231*e15f9bd0SJeremy L Thompson                                Q1d, Q1d, P1d, 1, Q1d); CeedChk(ierr);
2327a982d89SJeremy L. Thompson 
2337a982d89SJeremy L. Thompson   ierr = CeedFree(&interp1d); CeedChk(ierr);
2347a982d89SJeremy L. Thompson   ierr = CeedFree(&grad1d); CeedChk(ierr);
235*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2367a982d89SJeremy L. Thompson }
2377a982d89SJeremy L. Thompson 
2387a982d89SJeremy L. Thompson /**
2397a982d89SJeremy L. Thompson   @brief Get Ceed associated with a CeedBasis
2407a982d89SJeremy L. Thompson 
2417a982d89SJeremy L. Thompson   @param basis      CeedBasis
2427a982d89SJeremy L. Thompson   @param[out] ceed  Variable to store Ceed
2437a982d89SJeremy L. Thompson 
2447a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
2457a982d89SJeremy L. Thompson 
2467a982d89SJeremy L. Thompson   @ref Backend
2477a982d89SJeremy L. Thompson **/
2487a982d89SJeremy L. Thompson int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) {
2497a982d89SJeremy L. Thompson   *ceed = basis->ceed;
250*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2517a982d89SJeremy L. Thompson }
2527a982d89SJeremy L. Thompson 
2537a982d89SJeremy L. Thompson /**
2547a982d89SJeremy L. Thompson   @brief Get tensor status for given CeedBasis
2557a982d89SJeremy L. Thompson 
2567a982d89SJeremy L. Thompson   @param basis          CeedBasis
257fd364f38SJeremy L Thompson   @param[out] istensor  Variable to store tensor status
2587a982d89SJeremy L. Thompson 
2597a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
2607a982d89SJeremy L. Thompson 
2617a982d89SJeremy L. Thompson   @ref Backend
2627a982d89SJeremy L. Thompson **/
263fd364f38SJeremy L Thompson int CeedBasisIsTensor(CeedBasis basis, bool *istensor) {
264fd364f38SJeremy L Thompson   *istensor = basis->tensorbasis;
265*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2667a982d89SJeremy L. Thompson }
2677a982d89SJeremy L. Thompson 
2687a982d89SJeremy L. Thompson /**
2697a982d89SJeremy L. Thompson   @brief Get backend data of a CeedBasis
2707a982d89SJeremy L. Thompson 
2717a982d89SJeremy L. Thompson   @param basis      CeedBasis
2727a982d89SJeremy L. Thompson   @param[out] data  Variable to store data
2737a982d89SJeremy L. Thompson 
2747a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
2757a982d89SJeremy L. Thompson 
2767a982d89SJeremy L. Thompson   @ref Backend
2777a982d89SJeremy L. Thompson **/
278777ff853SJeremy L Thompson int CeedBasisGetData(CeedBasis basis, void *data) {
279777ff853SJeremy L Thompson   *(void **)data = basis->data;
280*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2817a982d89SJeremy L. Thompson }
2827a982d89SJeremy L. Thompson 
2837a982d89SJeremy L. Thompson /**
2847a982d89SJeremy L. Thompson   @brief Set backend data of a CeedBasis
2857a982d89SJeremy L. Thompson 
2867a982d89SJeremy L. Thompson   @param[out] basis  CeedBasis
2877a982d89SJeremy L. Thompson   @param data        Data to set
2887a982d89SJeremy L. Thompson 
2897a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
2907a982d89SJeremy L. Thompson 
2917a982d89SJeremy L. Thompson   @ref Backend
2927a982d89SJeremy L. Thompson **/
293777ff853SJeremy L Thompson int CeedBasisSetData(CeedBasis basis, void *data) {
294777ff853SJeremy L Thompson   basis->data = data;
295*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2967a982d89SJeremy L. Thompson }
2977a982d89SJeremy L. Thompson 
2987a982d89SJeremy L. Thompson /**
2997a982d89SJeremy L. Thompson   @brief Get dimension for given CeedElemTopology
3007a982d89SJeremy L. Thompson 
3017a982d89SJeremy L. Thompson   @param topo      CeedElemTopology
3027a982d89SJeremy L. Thompson   @param[out] dim  Variable to store dimension of topology
3037a982d89SJeremy L. Thompson 
3047a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3057a982d89SJeremy L. Thompson 
3067a982d89SJeremy L. Thompson   @ref Backend
3077a982d89SJeremy L. Thompson **/
3087a982d89SJeremy L. Thompson int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) {
3097a982d89SJeremy L. Thompson   *dim = (CeedInt) topo >> 16;
310*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3117a982d89SJeremy L. Thompson }
3127a982d89SJeremy L. Thompson 
3137a982d89SJeremy L. Thompson /**
3147a982d89SJeremy L. Thompson   @brief Get CeedTensorContract of a CeedBasis
3157a982d89SJeremy L. Thompson 
3167a982d89SJeremy L. Thompson   @param basis          CeedBasis
3177a982d89SJeremy L. Thompson   @param[out] contract  Variable to store CeedTensorContract
3187a982d89SJeremy L. Thompson 
3197a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3207a982d89SJeremy L. Thompson 
3217a982d89SJeremy L. Thompson   @ref Backend
3227a982d89SJeremy L. Thompson **/
3237a982d89SJeremy L. Thompson int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) {
3247a982d89SJeremy L. Thompson   *contract = basis->contract;
325*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3267a982d89SJeremy L. Thompson }
3277a982d89SJeremy L. Thompson 
3287a982d89SJeremy L. Thompson /**
3297a982d89SJeremy L. Thompson   @brief Set CeedTensorContract of a CeedBasis
3307a982d89SJeremy L. Thompson 
3317a982d89SJeremy L. Thompson   @param[out] basis     CeedBasis
3327a982d89SJeremy L. Thompson   @param contract       CeedTensorContract to set
3337a982d89SJeremy L. Thompson 
3347a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3357a982d89SJeremy L. Thompson 
3367a982d89SJeremy L. Thompson   @ref Backend
3377a982d89SJeremy L. Thompson **/
3387a982d89SJeremy L. Thompson int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract *contract) {
3397a982d89SJeremy L. Thompson   basis->contract = *contract;
340*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3417a982d89SJeremy L. Thompson }
3427a982d89SJeremy L. Thompson 
3437a982d89SJeremy L. Thompson /**
3447a982d89SJeremy L. Thompson   @brief Return a reference implementation of matrix multiplication C = A B.
3457a982d89SJeremy L. Thompson            Note, this is a reference implementation for CPU CeedScalar pointers
3467a982d89SJeremy L. Thompson            that is not intended for high performance.
3477a982d89SJeremy L. Thompson 
3487a982d89SJeremy L. Thompson   @param ceed         A Ceed context for error handling
3497a982d89SJeremy L. Thompson   @param[in] matA     Row-major matrix A
3507a982d89SJeremy L. Thompson   @param[in] matB     Row-major matrix B
3517a982d89SJeremy L. Thompson   @param[out] matC    Row-major output matrix C
3527a982d89SJeremy L. Thompson   @param m            Number of rows of C
3537a982d89SJeremy L. Thompson   @param n            Number of columns of C
3547a982d89SJeremy L. Thompson   @param kk           Number of columns of A/rows of B
3557a982d89SJeremy L. Thompson 
3567a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3577a982d89SJeremy L. Thompson 
3587a982d89SJeremy L. Thompson   @ref Utility
3597a982d89SJeremy L. Thompson **/
3607a982d89SJeremy L. Thompson int CeedMatrixMultiply(Ceed ceed, const CeedScalar *matA,
3617a982d89SJeremy L. Thompson                        const CeedScalar *matB, CeedScalar *matC, CeedInt m,
3627a982d89SJeremy L. Thompson                        CeedInt n, CeedInt kk) {
3637a982d89SJeremy L. Thompson   for (CeedInt i=0; i<m; i++)
3647a982d89SJeremy L. Thompson     for (CeedInt j=0; j<n; j++) {
3657a982d89SJeremy L. Thompson       CeedScalar sum = 0;
3667a982d89SJeremy L. Thompson       for (CeedInt k=0; k<kk; k++)
3677a982d89SJeremy L. Thompson         sum += matA[k+i*kk]*matB[j+k*n];
3687a982d89SJeremy L. Thompson       matC[j+i*n] = sum;
3697a982d89SJeremy L. Thompson     }
370*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3717a982d89SJeremy L. Thompson }
3727a982d89SJeremy L. Thompson 
3737a982d89SJeremy L. Thompson /// @}
3747a982d89SJeremy L. Thompson 
3757a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
3767a982d89SJeremy L. Thompson /// CeedBasis Public API
3777a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
3787a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser
379d7b241e6Sjeremylt /// @{
380d7b241e6Sjeremylt 
381b11c1e72Sjeremylt /**
38295bb1877Svaleriabarra   @brief Create a tensor-product basis for H^1 discretizations
383b11c1e72Sjeremylt 
384b11c1e72Sjeremylt   @param ceed        A Ceed object where the CeedBasis will be created
385b11c1e72Sjeremylt   @param dim         Topological dimension
386b11c1e72Sjeremylt   @param ncomp       Number of field components (1 for scalar fields)
387b11c1e72Sjeremylt   @param P1d         Number of nodes in one dimension
388b11c1e72Sjeremylt   @param Q1d         Number of quadrature points in one dimension
38995bb1877Svaleriabarra   @param interp1d    Row-major (Q1d * P1d) matrix expressing the values of nodal
390b11c1e72Sjeremylt                        basis functions at quadrature points
39195bb1877Svaleriabarra   @param grad1d      Row-major (Q1d * P1d) matrix expressing derivatives of nodal
392b11c1e72Sjeremylt                        basis functions at quadrature points
393b11c1e72Sjeremylt   @param qref1d      Array of length Q1d holding the locations of quadrature points
394b11c1e72Sjeremylt                        on the 1D reference element [-1, 1]
395b11c1e72Sjeremylt   @param qweight1d   Array of length Q1d holding the quadrature weights on the
396b11c1e72Sjeremylt                        reference element
397b11c1e72Sjeremylt   @param[out] basis  Address of the variable where the newly created
398b11c1e72Sjeremylt                        CeedBasis will be stored.
399b11c1e72Sjeremylt 
400b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
401dfdf5a53Sjeremylt 
4027a982d89SJeremy L. Thompson   @ref User
403b11c1e72Sjeremylt **/
404d7b241e6Sjeremylt int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt ncomp, CeedInt P1d,
405d7b241e6Sjeremylt                             CeedInt Q1d, const CeedScalar *interp1d,
406d7b241e6Sjeremylt                             const CeedScalar *grad1d, const CeedScalar *qref1d,
407d7b241e6Sjeremylt                             const CeedScalar *qweight1d, CeedBasis *basis) {
408d7b241e6Sjeremylt   int ierr;
409d7b241e6Sjeremylt 
4105fe0d4faSjeremylt   if (!ceed->BasisCreateTensorH1) {
4115fe0d4faSjeremylt     Ceed delegate;
412aefd8378Sjeremylt     ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr);
4135fe0d4faSjeremylt 
4145fe0d4faSjeremylt     if (!delegate)
415c042f62fSJeremy L Thompson       // LCOV_EXCL_START
416*e15f9bd0SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_UNSUPPORTED,
417*e15f9bd0SJeremy L Thompson                        "Backend does not support BasisCreateTensorH1");
418c042f62fSJeremy L Thompson     // LCOV_EXCL_STOP
4195fe0d4faSjeremylt 
4205fe0d4faSjeremylt     ierr = CeedBasisCreateTensorH1(delegate, dim, ncomp, P1d,
4215fe0d4faSjeremylt                                    Q1d, interp1d, grad1d, qref1d,
4225fe0d4faSjeremylt                                    qweight1d, basis); CeedChk(ierr);
423*e15f9bd0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
4245fe0d4faSjeremylt   }
425*e15f9bd0SJeremy L Thompson 
426*e15f9bd0SJeremy L Thompson   if (dim<1)
427*e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
428*e15f9bd0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
429*e15f9bd0SJeremy L Thompson                      "Basis dimension must be a positive value");
430*e15f9bd0SJeremy L Thompson   // LCOV_EXCL_STOP
431*e15f9bd0SJeremy L Thompson   CeedElemTopology topo = dim == 1 ? CEED_LINE :
432*e15f9bd0SJeremy L Thompson                           dim == 2 ? CEED_QUAD :
433*e15f9bd0SJeremy L Thompson                           CEED_HEX;
434*e15f9bd0SJeremy L Thompson 
435d7b241e6Sjeremylt   ierr = CeedCalloc(1, basis); CeedChk(ierr);
436d7b241e6Sjeremylt   (*basis)->ceed = ceed;
437d7b241e6Sjeremylt   ceed->refcount++;
438d7b241e6Sjeremylt   (*basis)->refcount = 1;
439a8de75f0Sjeremylt   (*basis)->tensorbasis = 1;
440d7b241e6Sjeremylt   (*basis)->dim = dim;
441d99fa3c5SJeremy L Thompson   (*basis)->topo = topo;
442d7b241e6Sjeremylt   (*basis)->ncomp = ncomp;
443d7b241e6Sjeremylt   (*basis)->P1d = P1d;
444d7b241e6Sjeremylt   (*basis)->Q1d = Q1d;
445a8de75f0Sjeremylt   (*basis)->P = CeedIntPow(P1d, dim);
446a8de75f0Sjeremylt   (*basis)->Q = CeedIntPow(Q1d, dim);
447d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d,&(*basis)->qref1d); CeedChk(ierr);
448d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d,&(*basis)->qweight1d); CeedChk(ierr);
449d7b241e6Sjeremylt   memcpy((*basis)->qref1d, qref1d, Q1d*sizeof(qref1d[0]));
450d7b241e6Sjeremylt   memcpy((*basis)->qweight1d, qweight1d, Q1d*sizeof(qweight1d[0]));
451d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d*P1d,&(*basis)->interp1d); CeedChk(ierr);
452d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d*P1d,&(*basis)->grad1d); CeedChk(ierr);
453d7b241e6Sjeremylt   memcpy((*basis)->interp1d, interp1d, Q1d*P1d*sizeof(interp1d[0]));
45409486605Sjeremylt   memcpy((*basis)->grad1d, grad1d, Q1d*P1d*sizeof(grad1d[0]));
455667bc5fcSjeremylt   ierr = ceed->BasisCreateTensorH1(dim, P1d, Q1d, interp1d, grad1d, qref1d,
456d7b241e6Sjeremylt                                    qweight1d, *basis); CeedChk(ierr);
457*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
458d7b241e6Sjeremylt }
459d7b241e6Sjeremylt 
460b11c1e72Sjeremylt /**
46195bb1877Svaleriabarra   @brief Create a tensor-product Lagrange basis
462b11c1e72Sjeremylt 
463b11c1e72Sjeremylt   @param ceed        A Ceed object where the CeedBasis will be created
464b11c1e72Sjeremylt   @param dim         Topological dimension of element
46595bb1877Svaleriabarra   @param ncomp       Number of field components (1 for scalar fields)
466b11c1e72Sjeremylt   @param P           Number of Gauss-Lobatto nodes in one dimension.  The
467b11c1e72Sjeremylt                        polynomial degree of the resulting Q_k element is k=P-1.
468b11c1e72Sjeremylt   @param Q           Number of quadrature points in one dimension.
469b11c1e72Sjeremylt   @param qmode       Distribution of the Q quadrature points (affects order of
470b11c1e72Sjeremylt                        accuracy for the quadrature)
471b11c1e72Sjeremylt   @param[out] basis  Address of the variable where the newly created
472b11c1e72Sjeremylt                        CeedBasis will be stored.
473b11c1e72Sjeremylt 
474b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
475dfdf5a53Sjeremylt 
4767a982d89SJeremy L. Thompson   @ref User
477b11c1e72Sjeremylt **/
478d7b241e6Sjeremylt int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt ncomp,
479692c2638Sjeremylt                                     CeedInt P, CeedInt Q, CeedQuadMode qmode,
480692c2638Sjeremylt                                     CeedBasis *basis) {
481d7b241e6Sjeremylt   // Allocate
482*e15f9bd0SJeremy L Thompson   int ierr, ierr2, i, j, k;
483d7b241e6Sjeremylt   CeedScalar c1, c2, c3, c4, dx, *nodes, *interp1d, *grad1d, *qref1d, *qweight1d;
4844d537eeaSYohann 
4854d537eeaSYohann   if (dim<1)
486c042f62fSJeremy L Thompson     // LCOV_EXCL_START
487*e15f9bd0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
488*e15f9bd0SJeremy L Thompson                      "Basis dimension must be a positive value");
489c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
4904d537eeaSYohann 
491*e15f9bd0SJeremy L Thompson   // Get Nodes and Weights
492d7b241e6Sjeremylt   ierr = CeedCalloc(P*Q, &interp1d); CeedChk(ierr);
493d7b241e6Sjeremylt   ierr = CeedCalloc(P*Q, &grad1d); CeedChk(ierr);
494d7b241e6Sjeremylt   ierr = CeedCalloc(P, &nodes); CeedChk(ierr);
495d7b241e6Sjeremylt   ierr = CeedCalloc(Q, &qref1d); CeedChk(ierr);
496d7b241e6Sjeremylt   ierr = CeedCalloc(Q, &qweight1d); CeedChk(ierr);
497*e15f9bd0SJeremy L Thompson   ierr = CeedLobattoQuadrature(P, nodes, NULL);
498*e15f9bd0SJeremy L Thompson   if (ierr) { goto cleanup; } CeedChk(ierr);
499d7b241e6Sjeremylt   switch (qmode) {
500d7b241e6Sjeremylt   case CEED_GAUSS:
501*e15f9bd0SJeremy L Thompson     ierr = CeedGaussQuadrature(Q, qref1d, qweight1d);
502d7b241e6Sjeremylt     break;
503d7b241e6Sjeremylt   case CEED_GAUSS_LOBATTO:
504*e15f9bd0SJeremy L Thompson     ierr = CeedLobattoQuadrature(Q, qref1d, qweight1d);
505d7b241e6Sjeremylt     break;
506d7b241e6Sjeremylt   }
507*e15f9bd0SJeremy L Thompson   if (ierr) { goto cleanup; } CeedChk(ierr);
508*e15f9bd0SJeremy L Thompson 
509d7b241e6Sjeremylt   // Build B, D matrix
510d7b241e6Sjeremylt   // Fornberg, 1998
511d7b241e6Sjeremylt   for (i = 0; i  < Q; i++) {
512d7b241e6Sjeremylt     c1 = 1.0;
513d7b241e6Sjeremylt     c3 = nodes[0] - qref1d[i];
514d7b241e6Sjeremylt     interp1d[i*P+0] = 1.0;
515d7b241e6Sjeremylt     for (j = 1; j < P; j++) {
516d7b241e6Sjeremylt       c2 = 1.0;
517d7b241e6Sjeremylt       c4 = c3;
518d7b241e6Sjeremylt       c3 = nodes[j] - qref1d[i];
519d7b241e6Sjeremylt       for (k = 0; k < j; k++) {
520d7b241e6Sjeremylt         dx = nodes[j] - nodes[k];
521d7b241e6Sjeremylt         c2 *= dx;
522d7b241e6Sjeremylt         if (k == j - 1) {
523d7b241e6Sjeremylt           grad1d[i*P + j] = c1*(interp1d[i*P + k] - c4*grad1d[i*P + k]) / c2;
524d7b241e6Sjeremylt           interp1d[i*P + j] = - c1*c4*interp1d[i*P + k] / c2;
525d7b241e6Sjeremylt         }
526d7b241e6Sjeremylt         grad1d[i*P + k] = (c3*grad1d[i*P + k] - interp1d[i*P + k]) / dx;
527d7b241e6Sjeremylt         interp1d[i*P + k] = c3*interp1d[i*P + k] / dx;
528d7b241e6Sjeremylt       }
529d7b241e6Sjeremylt       c1 = c2;
530d7b241e6Sjeremylt     }
531d7b241e6Sjeremylt   }
532d7b241e6Sjeremylt   //  // Pass to CeedBasisCreateTensorH1
533d7b241e6Sjeremylt   ierr = CeedBasisCreateTensorH1(ceed, dim, ncomp, P, Q, interp1d, grad1d, qref1d,
534d7b241e6Sjeremylt                                  qweight1d, basis); CeedChk(ierr);
535*e15f9bd0SJeremy L Thompson cleanup:
536*e15f9bd0SJeremy L Thompson   ierr2 = CeedFree(&interp1d); CeedChk(ierr2);
537*e15f9bd0SJeremy L Thompson   ierr2 = CeedFree(&grad1d); CeedChk(ierr2);
538*e15f9bd0SJeremy L Thompson   ierr2 = CeedFree(&nodes); CeedChk(ierr2);
539*e15f9bd0SJeremy L Thompson   ierr2 = CeedFree(&qref1d); CeedChk(ierr2);
540*e15f9bd0SJeremy L Thompson   ierr2 = CeedFree(&qweight1d); CeedChk(ierr2);
541*e15f9bd0SJeremy L Thompson   CeedChk(ierr);
542*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
543d7b241e6Sjeremylt }
544d7b241e6Sjeremylt 
545b11c1e72Sjeremylt /**
54695bb1877Svaleriabarra   @brief Create a non tensor-product basis for H^1 discretizations
547a8de75f0Sjeremylt 
548a8de75f0Sjeremylt   @param ceed        A Ceed object where the CeedBasis will be created
549a8de75f0Sjeremylt   @param topo        Topology of element, e.g. hypercube, simplex, ect
550a8de75f0Sjeremylt   @param ncomp       Number of field components (1 for scalar fields)
5518795c945Sjeremylt   @param nnodes      Total number of nodes
552a8de75f0Sjeremylt   @param nqpts       Total number of quadrature points
55395bb1877Svaleriabarra   @param interp      Row-major (nqpts * nnodes) matrix expressing the values of
5548795c945Sjeremylt                        nodal basis functions at quadrature points
55595bb1877Svaleriabarra   @param grad        Row-major (nqpts * dim * nnodes) matrix expressing
5568795c945Sjeremylt                        derivatives of nodal basis functions at quadrature points
5578795c945Sjeremylt   @param qref        Array of length nqpts holding the locations of quadrature
5588795c945Sjeremylt                        points on the reference element [-1, 1]
559a8de75f0Sjeremylt   @param qweight     Array of length nqpts holding the quadrature weights on the
560a8de75f0Sjeremylt                        reference element
561a8de75f0Sjeremylt   @param[out] basis  Address of the variable where the newly created
562a8de75f0Sjeremylt                        CeedBasis will be stored.
563a8de75f0Sjeremylt 
564a8de75f0Sjeremylt   @return An error code: 0 - success, otherwise - failure
565a8de75f0Sjeremylt 
5667a982d89SJeremy L. Thompson   @ref User
567a8de75f0Sjeremylt **/
568a8de75f0Sjeremylt int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt ncomp,
569692c2638Sjeremylt                       CeedInt nnodes, CeedInt nqpts, const CeedScalar *interp,
570a8de75f0Sjeremylt                       const CeedScalar *grad, const CeedScalar *qref,
571a8de75f0Sjeremylt                       const CeedScalar *qweight, CeedBasis *basis) {
572a8de75f0Sjeremylt   int ierr;
5738795c945Sjeremylt   CeedInt P = nnodes, Q = nqpts, dim = 0;
574a8de75f0Sjeremylt 
5755fe0d4faSjeremylt   if (!ceed->BasisCreateH1) {
5765fe0d4faSjeremylt     Ceed delegate;
577aefd8378Sjeremylt     ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr);
5785fe0d4faSjeremylt 
5795fe0d4faSjeremylt     if (!delegate)
580c042f62fSJeremy L Thompson       // LCOV_EXCL_START
581*e15f9bd0SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_UNSUPPORTED,
582*e15f9bd0SJeremy L Thompson                        "Backend does not support BasisCreateH1");
583c042f62fSJeremy L Thompson     // LCOV_EXCL_STOP
5845fe0d4faSjeremylt 
5858795c945Sjeremylt     ierr = CeedBasisCreateH1(delegate, topo, ncomp, nnodes,
5865fe0d4faSjeremylt                              nqpts, interp, grad, qref,
5875fe0d4faSjeremylt                              qweight, basis); CeedChk(ierr);
588*e15f9bd0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
5895fe0d4faSjeremylt   }
5905fe0d4faSjeremylt 
591a8de75f0Sjeremylt   ierr = CeedCalloc(1,basis); CeedChk(ierr);
592a8de75f0Sjeremylt 
593a8de75f0Sjeremylt   ierr = CeedBasisGetTopologyDimension(topo, &dim); CeedChk(ierr);
594a8de75f0Sjeremylt 
595a8de75f0Sjeremylt   (*basis)->ceed = ceed;
596a8de75f0Sjeremylt   ceed->refcount++;
597a8de75f0Sjeremylt   (*basis)->refcount = 1;
598a8de75f0Sjeremylt   (*basis)->tensorbasis = 0;
599a8de75f0Sjeremylt   (*basis)->dim = dim;
600d99fa3c5SJeremy L Thompson   (*basis)->topo = topo;
601a8de75f0Sjeremylt   (*basis)->ncomp = ncomp;
602a8de75f0Sjeremylt   (*basis)->P = P;
603a8de75f0Sjeremylt   (*basis)->Q = Q;
604a8de75f0Sjeremylt   ierr = CeedMalloc(Q*dim,&(*basis)->qref1d); CeedChk(ierr);
605a8de75f0Sjeremylt   ierr = CeedMalloc(Q,&(*basis)->qweight1d); CeedChk(ierr);
606a8de75f0Sjeremylt   memcpy((*basis)->qref1d, qref, Q*dim*sizeof(qref[0]));
607a8de75f0Sjeremylt   memcpy((*basis)->qweight1d, qweight, Q*sizeof(qweight[0]));
60800f91b2bSjeremylt   ierr = CeedMalloc(Q*P, &(*basis)->interp); CeedChk(ierr);
60900f91b2bSjeremylt   ierr = CeedMalloc(dim*Q*P, &(*basis)->grad); CeedChk(ierr);
61000f91b2bSjeremylt   memcpy((*basis)->interp, interp, Q*P*sizeof(interp[0]));
61100f91b2bSjeremylt   memcpy((*basis)->grad, grad, dim*Q*P*sizeof(grad[0]));
612667bc5fcSjeremylt   ierr = ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, qref,
613a8de75f0Sjeremylt                              qweight, *basis); CeedChk(ierr);
614*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
615a8de75f0Sjeremylt }
616a8de75f0Sjeremylt 
617a8de75f0Sjeremylt /**
6187a982d89SJeremy L. Thompson   @brief View a CeedBasis
6197a982d89SJeremy L. Thompson 
6207a982d89SJeremy L. Thompson   @param basis   CeedBasis to view
6217a982d89SJeremy L. Thompson   @param stream  Stream to view to, e.g., stdout
6227a982d89SJeremy L. Thompson 
6237a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
6247a982d89SJeremy L. Thompson 
6257a982d89SJeremy L. Thompson   @ref User
6267a982d89SJeremy L. Thompson **/
6277a982d89SJeremy L. Thompson int CeedBasisView(CeedBasis basis, FILE *stream) {
6287a982d89SJeremy L. Thompson   int ierr;
6297a982d89SJeremy L. Thompson 
6307a982d89SJeremy L. Thompson   if (basis->tensorbasis) {
6317a982d89SJeremy L. Thompson     fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P1d,
6327a982d89SJeremy L. Thompson             basis->Q1d);
6337a982d89SJeremy L. Thompson     ierr = CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q1d, basis->qref1d,
6347a982d89SJeremy L. Thompson                           stream); CeedChk(ierr);
6357a982d89SJeremy L. Thompson     ierr = CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q1d,
6367a982d89SJeremy L. Thompson                           basis->qweight1d, stream); CeedChk(ierr);
6377a982d89SJeremy L. Thompson     ierr = CeedScalarView("interp1d", "\t% 12.8f", basis->Q1d, basis->P1d,
6387a982d89SJeremy L. Thompson                           basis->interp1d, stream); CeedChk(ierr);
6397a982d89SJeremy L. Thompson     ierr = CeedScalarView("grad1d", "\t% 12.8f", basis->Q1d, basis->P1d,
6407a982d89SJeremy L. Thompson                           basis->grad1d, stream); CeedChk(ierr);
6417a982d89SJeremy L. Thompson   } else {
6427a982d89SJeremy L. Thompson     fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P,
6437a982d89SJeremy L. Thompson             basis->Q);
6447a982d89SJeremy L. Thompson     ierr = CeedScalarView("qref", "\t% 12.8f", 1, basis->Q*basis->dim,
6457a982d89SJeremy L. Thompson                           basis->qref1d,
6467a982d89SJeremy L. Thompson                           stream); CeedChk(ierr);
6477a982d89SJeremy L. Thompson     ierr = CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->qweight1d,
6487a982d89SJeremy L. Thompson                           stream); CeedChk(ierr);
6497a982d89SJeremy L. Thompson     ierr = CeedScalarView("interp", "\t% 12.8f", basis->Q, basis->P,
6507a982d89SJeremy L. Thompson                           basis->interp, stream); CeedChk(ierr);
6517a982d89SJeremy L. Thompson     ierr = CeedScalarView("grad", "\t% 12.8f", basis->dim*basis->Q, basis->P,
6527a982d89SJeremy L. Thompson                           basis->grad, stream); CeedChk(ierr);
6537a982d89SJeremy L. Thompson   }
654*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
6557a982d89SJeremy L. Thompson }
6567a982d89SJeremy L. Thompson 
6577a982d89SJeremy L. Thompson /**
6587a982d89SJeremy L. Thompson   @brief Apply basis evaluation from nodes to quadrature points or vice versa
6597a982d89SJeremy L. Thompson 
6607a982d89SJeremy L. Thompson   @param basis   CeedBasis to evaluate
6617a982d89SJeremy L. Thompson   @param nelem   The number of elements to apply the basis evaluation to;
6627a982d89SJeremy L. Thompson                    the backend will specify the ordering in
6634cc79fe7SJed Brown                    CeedElemRestrictionCreateBlocked()
6647a982d89SJeremy L. Thompson   @param tmode   \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature
6657a982d89SJeremy L. Thompson                    points, \ref CEED_TRANSPOSE to apply the transpose, mapping
6667a982d89SJeremy L. Thompson                    from quadrature points to nodes
6677a982d89SJeremy L. Thompson   @param emode   \ref CEED_EVAL_NONE to use values directly,
6687a982d89SJeremy L. Thompson                    \ref CEED_EVAL_INTERP to use interpolated values,
6697a982d89SJeremy L. Thompson                    \ref CEED_EVAL_GRAD to use gradients,
6707a982d89SJeremy L. Thompson                    \ref CEED_EVAL_WEIGHT to use quadrature weights.
6717a982d89SJeremy L. Thompson   @param[in] u   Input CeedVector
6727a982d89SJeremy L. Thompson   @param[out] v  Output CeedVector
6737a982d89SJeremy L. Thompson 
6747a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
6757a982d89SJeremy L. Thompson 
6767a982d89SJeremy L. Thompson   @ref User
6777a982d89SJeremy L. Thompson **/
6787a982d89SJeremy L. Thompson int CeedBasisApply(CeedBasis basis, CeedInt nelem, CeedTransposeMode tmode,
6797a982d89SJeremy L. Thompson                    CeedEvalMode emode, CeedVector u, CeedVector v) {
6807a982d89SJeremy L. Thompson   int ierr;
681*e15f9bd0SJeremy L Thompson   CeedInt ulength = 0, vlength, dim, ncomp, nnodes, nqpts;
682*e15f9bd0SJeremy L Thompson   ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr);
683*e15f9bd0SJeremy L Thompson   ierr = CeedBasisGetNumComponents(basis, &ncomp); CeedChk(ierr);
6847a982d89SJeremy L. Thompson   ierr = CeedBasisGetNumNodes(basis, &nnodes); CeedChk(ierr);
685*e15f9bd0SJeremy L Thompson   ierr = CeedBasisGetNumQuadraturePoints(basis, &nqpts); CeedChk(ierr);
6867a982d89SJeremy L. Thompson   ierr = CeedVectorGetLength(v, &vlength); CeedChk(ierr);
6877a982d89SJeremy L. Thompson   if (u) {
6887a982d89SJeremy L. Thompson     ierr = CeedVectorGetLength(u, &ulength); CeedChk(ierr);
6897a982d89SJeremy L. Thompson   }
6907a982d89SJeremy L. Thompson 
691*e15f9bd0SJeremy L Thompson   if (!basis->Apply)
692*e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
693*e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_UNSUPPORTED,
694*e15f9bd0SJeremy L Thompson                      "Backend does not support BasisApply");
695*e15f9bd0SJeremy L Thompson   // LCOV_EXCL_STOP
696*e15f9bd0SJeremy L Thompson 
697*e15f9bd0SJeremy L Thompson   // Check compatibility of topological and geometrical dimensions
698*e15f9bd0SJeremy L Thompson   if ((tmode == CEED_TRANSPOSE && (vlength%nnodes != 0 || ulength%nqpts != 0)) ||
699*e15f9bd0SJeremy L Thompson       (tmode == CEED_NOTRANSPOSE && (ulength%nnodes != 0 || vlength%nqpts != 0)))
700*e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_DIMENSION,
701*e15f9bd0SJeremy L Thompson                      "Length of input/output vectors "
7027a982d89SJeremy L. Thompson                      "incompatible with basis dimensions");
7037a982d89SJeremy L. Thompson 
704*e15f9bd0SJeremy L Thompson   // Check vector lengths to prevent out of bounds issues
705*e15f9bd0SJeremy L Thompson   bool baddims = false;
706*e15f9bd0SJeremy L Thompson   switch (emode) {
707*e15f9bd0SJeremy L Thompson   case CEED_EVAL_NONE:
708*e15f9bd0SJeremy L Thompson   case CEED_EVAL_INTERP: baddims =
709*e15f9bd0SJeremy L Thompson       ((tmode == CEED_TRANSPOSE && (ulength < nelem*ncomp*nqpts
710*e15f9bd0SJeremy L Thompson                                     || vlength < nelem*ncomp*nnodes)) ||
711*e15f9bd0SJeremy L Thompson        (tmode == CEED_NOTRANSPOSE && (vlength < nelem*nqpts*ncomp
712*e15f9bd0SJeremy L Thompson                                       || ulength < nelem*ncomp*nnodes)));
713*e15f9bd0SJeremy L Thompson     break;
714*e15f9bd0SJeremy L Thompson   case CEED_EVAL_GRAD: baddims =
715*e15f9bd0SJeremy L Thompson       ((tmode == CEED_TRANSPOSE && (ulength < nelem*ncomp*nqpts*dim
716*e15f9bd0SJeremy L Thompson                                     || vlength < nelem*ncomp*nnodes)) ||
717*e15f9bd0SJeremy L Thompson        (tmode == CEED_NOTRANSPOSE && (vlength < nelem*nqpts*ncomp*dim
718*e15f9bd0SJeremy L Thompson                                       || ulength < nelem*ncomp*nnodes)));
719*e15f9bd0SJeremy L Thompson     break;
720*e15f9bd0SJeremy L Thompson   case CEED_EVAL_WEIGHT:
721*e15f9bd0SJeremy L Thompson     baddims = vlength < nelem*nqpts;
722*e15f9bd0SJeremy L Thompson     break;
723*e15f9bd0SJeremy L Thompson   // LCOV_EXCL_START
724*e15f9bd0SJeremy L Thompson   case CEED_EVAL_DIV: baddims =
725*e15f9bd0SJeremy L Thompson       ((tmode == CEED_TRANSPOSE && (ulength < nelem*ncomp*nqpts
726*e15f9bd0SJeremy L Thompson                                     || vlength < nelem*ncomp*nnodes)) ||
727*e15f9bd0SJeremy L Thompson        (tmode == CEED_NOTRANSPOSE && (vlength < nelem*nqpts*ncomp
728*e15f9bd0SJeremy L Thompson                                       || ulength < nelem*ncomp*nnodes)));
729*e15f9bd0SJeremy L Thompson     break;
730*e15f9bd0SJeremy L Thompson   case CEED_EVAL_CURL: baddims =
731*e15f9bd0SJeremy L Thompson       ((tmode == CEED_TRANSPOSE && (ulength < nelem*ncomp*nqpts
732*e15f9bd0SJeremy L Thompson                                     || vlength < nelem*ncomp*nnodes)) ||
733*e15f9bd0SJeremy L Thompson        (tmode == CEED_NOTRANSPOSE && (vlength < nelem*nqpts*ncomp
734*e15f9bd0SJeremy L Thompson                                       || ulength < nelem*ncomp*nnodes)));
735*e15f9bd0SJeremy L Thompson     break;
736*e15f9bd0SJeremy L Thompson     // LCOV_EXCL_STOP
737*e15f9bd0SJeremy L Thompson   }
738*e15f9bd0SJeremy L Thompson   if (baddims)
739*e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
740*e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_DIMENSION,
741*e15f9bd0SJeremy L Thompson                      "Input/output vectors too short for basis and evalualtion mode");
742*e15f9bd0SJeremy L Thompson   // LCOV_EXCL_STOP
743*e15f9bd0SJeremy L Thompson 
7447a982d89SJeremy L. Thompson   ierr = basis->Apply(basis, nelem, tmode, emode, u, v); CeedChk(ierr);
745*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
7467a982d89SJeremy L. Thompson }
7477a982d89SJeremy L. Thompson 
7487a982d89SJeremy L. Thompson /**
7499d007619Sjeremylt   @brief Get dimension for given CeedBasis
7509d007619Sjeremylt 
7519d007619Sjeremylt   @param basis     CeedBasis
7529d007619Sjeremylt   @param[out] dim  Variable to store dimension of basis
7539d007619Sjeremylt 
7549d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
7559d007619Sjeremylt 
7569d007619Sjeremylt   @ref Backend
7579d007619Sjeremylt **/
7589d007619Sjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) {
7599d007619Sjeremylt   *dim = basis->dim;
760*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
7619d007619Sjeremylt }
7629d007619Sjeremylt 
7639d007619Sjeremylt /**
764d99fa3c5SJeremy L Thompson   @brief Get topology for given CeedBasis
765d99fa3c5SJeremy L Thompson 
766d99fa3c5SJeremy L Thompson   @param basis      CeedBasis
767d99fa3c5SJeremy L Thompson   @param[out] topo  Variable to store topology of basis
768d99fa3c5SJeremy L Thompson 
769d99fa3c5SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
770d99fa3c5SJeremy L Thompson 
771d99fa3c5SJeremy L Thompson   @ref Backend
772d99fa3c5SJeremy L Thompson **/
773d99fa3c5SJeremy L Thompson int CeedBasisGetTopology(CeedBasis basis, CeedElemTopology *topo) {
774d99fa3c5SJeremy L Thompson   *topo = basis->topo;
775*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
776d99fa3c5SJeremy L Thompson }
777d99fa3c5SJeremy L Thompson 
778d99fa3c5SJeremy L Thompson /**
7799d007619Sjeremylt   @brief Get number of components for given CeedBasis
7809d007619Sjeremylt 
7819d007619Sjeremylt   @param basis         CeedBasis
7829d007619Sjeremylt   @param[out] numcomp  Variable to store number of components of basis
7839d007619Sjeremylt 
7849d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
7859d007619Sjeremylt 
7869d007619Sjeremylt   @ref Backend
7879d007619Sjeremylt **/
7889d007619Sjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *numcomp) {
7899d007619Sjeremylt   *numcomp = basis->ncomp;
790*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
7919d007619Sjeremylt }
7929d007619Sjeremylt 
7939d007619Sjeremylt /**
7949d007619Sjeremylt   @brief Get total number of nodes (in dim dimensions) of a CeedBasis
7959d007619Sjeremylt 
7969d007619Sjeremylt   @param basis   CeedBasis
7979d007619Sjeremylt   @param[out] P  Variable to store number of nodes
7989d007619Sjeremylt 
7999d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
8009d007619Sjeremylt 
8019d007619Sjeremylt   @ref Utility
8029d007619Sjeremylt **/
8039d007619Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) {
8049d007619Sjeremylt   *P = basis->P;
805*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8069d007619Sjeremylt }
8079d007619Sjeremylt 
8089d007619Sjeremylt /**
8099d007619Sjeremylt   @brief Get total number of nodes (in 1 dimension) of a CeedBasis
8109d007619Sjeremylt 
8119d007619Sjeremylt   @param basis     CeedBasis
8129d007619Sjeremylt   @param[out] P1d  Variable to store number of nodes
8139d007619Sjeremylt 
8149d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
8159d007619Sjeremylt 
8169d007619Sjeremylt   @ref Backend
8179d007619Sjeremylt **/
8189d007619Sjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P1d) {
8199d007619Sjeremylt   if (!basis->tensorbasis)
8209d007619Sjeremylt     // LCOV_EXCL_START
821*e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR,
822*e15f9bd0SJeremy L Thompson                      "Cannot supply P1d for non-tensor basis");
8239d007619Sjeremylt   // LCOV_EXCL_STOP
8249d007619Sjeremylt 
8259d007619Sjeremylt   *P1d = basis->P1d;
826*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8279d007619Sjeremylt }
8289d007619Sjeremylt 
8299d007619Sjeremylt /**
8309d007619Sjeremylt   @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis
8319d007619Sjeremylt 
8329d007619Sjeremylt   @param basis   CeedBasis
8339d007619Sjeremylt   @param[out] Q  Variable to store number of quadrature points
8349d007619Sjeremylt 
8359d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
8369d007619Sjeremylt 
8379d007619Sjeremylt   @ref Utility
8389d007619Sjeremylt **/
8399d007619Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) {
8409d007619Sjeremylt   *Q = basis->Q;
841*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8429d007619Sjeremylt }
8439d007619Sjeremylt 
8449d007619Sjeremylt /**
8459d007619Sjeremylt   @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis
8469d007619Sjeremylt 
8479d007619Sjeremylt   @param basis     CeedBasis
8489d007619Sjeremylt   @param[out] Q1d  Variable to store number of quadrature points
8499d007619Sjeremylt 
8509d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
8519d007619Sjeremylt 
8529d007619Sjeremylt   @ref Backend
8539d007619Sjeremylt **/
8549d007619Sjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q1d) {
8559d007619Sjeremylt   if (!basis->tensorbasis)
8569d007619Sjeremylt     // LCOV_EXCL_START
857*e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR,
858*e15f9bd0SJeremy L Thompson                      "Cannot supply Q1d for non-tensor basis");
8599d007619Sjeremylt   // LCOV_EXCL_STOP
8609d007619Sjeremylt 
8619d007619Sjeremylt   *Q1d = basis->Q1d;
862*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8639d007619Sjeremylt }
8649d007619Sjeremylt 
8659d007619Sjeremylt /**
8669d007619Sjeremylt   @brief Get reference coordinates of quadrature points (in dim dimensions)
8679d007619Sjeremylt          of a CeedBasis
8689d007619Sjeremylt 
8699d007619Sjeremylt   @param basis      CeedBasis
8709d007619Sjeremylt   @param[out] qref  Variable to store reference coordinates of quadrature points
8719d007619Sjeremylt 
8729d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
8739d007619Sjeremylt 
8749d007619Sjeremylt   @ref Backend
8759d007619Sjeremylt **/
8766c58de82SJeremy L Thompson int CeedBasisGetQRef(CeedBasis basis, const CeedScalar **qref) {
8779d007619Sjeremylt   *qref = basis->qref1d;
878*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8799d007619Sjeremylt }
8809d007619Sjeremylt 
8819d007619Sjeremylt /**
8829d007619Sjeremylt   @brief Get quadrature weights of quadrature points (in dim dimensions)
8839d007619Sjeremylt          of a CeedBasis
8849d007619Sjeremylt 
8859d007619Sjeremylt   @param basis         CeedBasis
8869d007619Sjeremylt   @param[out] qweight  Variable to store quadrature weights
8879d007619Sjeremylt 
8889d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
8899d007619Sjeremylt 
8909d007619Sjeremylt   @ref Backend
8919d007619Sjeremylt **/
8926c58de82SJeremy L Thompson int CeedBasisGetQWeights(CeedBasis basis, const CeedScalar **qweight) {
8939d007619Sjeremylt   *qweight = basis->qweight1d;
894*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8959d007619Sjeremylt }
8969d007619Sjeremylt 
8979d007619Sjeremylt /**
8989d007619Sjeremylt   @brief Get interpolation matrix of a CeedBasis
8999d007619Sjeremylt 
9009d007619Sjeremylt   @param basis        CeedBasis
9019d007619Sjeremylt   @param[out] interp  Variable to store interpolation matrix
9029d007619Sjeremylt 
9039d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
9049d007619Sjeremylt 
9059d007619Sjeremylt   @ref Backend
9069d007619Sjeremylt **/
9076c58de82SJeremy L Thompson int CeedBasisGetInterp(CeedBasis basis, const CeedScalar **interp) {
9089d007619Sjeremylt   if (!basis->interp && basis->tensorbasis) {
9099d007619Sjeremylt     // Allocate
9109d007619Sjeremylt     int ierr;
9119d007619Sjeremylt     ierr = CeedMalloc(basis->Q*basis->P, &basis->interp); CeedChk(ierr);
9129d007619Sjeremylt 
9139d007619Sjeremylt     // Initialize
9149d007619Sjeremylt     for (CeedInt i=0; i<basis->Q*basis->P; i++)
9159d007619Sjeremylt       basis->interp[i] = 1.0;
9169d007619Sjeremylt 
9179d007619Sjeremylt     // Calculate
9189d007619Sjeremylt     for (CeedInt d=0; d<basis->dim; d++)
9199d007619Sjeremylt       for (CeedInt qpt=0; qpt<basis->Q; qpt++)
9209d007619Sjeremylt         for (CeedInt node=0; node<basis->P; node++) {
9219d007619Sjeremylt           CeedInt p = (node / CeedIntPow(basis->P1d, d)) % basis->P1d;
9229d007619Sjeremylt           CeedInt q = (qpt / CeedIntPow(basis->Q1d, d)) % basis->Q1d;
9239d007619Sjeremylt           basis->interp[qpt*(basis->P)+node] *= basis->interp1d[q*basis->P1d+p];
9249d007619Sjeremylt         }
9259d007619Sjeremylt   }
9269d007619Sjeremylt   *interp = basis->interp;
927*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9289d007619Sjeremylt }
9299d007619Sjeremylt 
9309d007619Sjeremylt /**
9319d007619Sjeremylt   @brief Get 1D interpolation matrix of a tensor product CeedBasis
9329d007619Sjeremylt 
9339d007619Sjeremylt   @param basis          CeedBasis
9349d007619Sjeremylt   @param[out] interp1d  Variable to store interpolation matrix
9359d007619Sjeremylt 
9369d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
9379d007619Sjeremylt 
9389d007619Sjeremylt   @ref Backend
9399d007619Sjeremylt **/
9406c58de82SJeremy L Thompson int CeedBasisGetInterp1D(CeedBasis basis, const CeedScalar **interp1d) {
9419d007619Sjeremylt   if (!basis->tensorbasis)
9429d007619Sjeremylt     // LCOV_EXCL_START
943*e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR,
944*e15f9bd0SJeremy L Thompson                      "CeedBasis is not a tensor product basis.");
9459d007619Sjeremylt   // LCOV_EXCL_STOP
9469d007619Sjeremylt 
9479d007619Sjeremylt   *interp1d = basis->interp1d;
948*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9499d007619Sjeremylt }
9509d007619Sjeremylt 
9519d007619Sjeremylt /**
9529d007619Sjeremylt   @brief Get gradient matrix of a CeedBasis
9539d007619Sjeremylt 
9549d007619Sjeremylt   @param basis      CeedBasis
9559d007619Sjeremylt   @param[out] grad  Variable to store gradient matrix
9569d007619Sjeremylt 
9579d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
9589d007619Sjeremylt 
9599d007619Sjeremylt   @ref Backend
9609d007619Sjeremylt **/
9616c58de82SJeremy L Thompson int CeedBasisGetGrad(CeedBasis basis, const CeedScalar **grad) {
9629d007619Sjeremylt   if (!basis->grad && basis->tensorbasis) {
9639d007619Sjeremylt     // Allocate
9649d007619Sjeremylt     int ierr;
9659d007619Sjeremylt     ierr = CeedMalloc(basis->dim*basis->Q*basis->P, &basis->grad);
9669d007619Sjeremylt     CeedChk(ierr);
9679d007619Sjeremylt 
9689d007619Sjeremylt     // Initialize
9699d007619Sjeremylt     for (CeedInt i=0; i<basis->dim*basis->Q*basis->P; i++)
9709d007619Sjeremylt       basis->grad[i] = 1.0;
9719d007619Sjeremylt 
9729d007619Sjeremylt     // Calculate
9739d007619Sjeremylt     for (CeedInt d=0; d<basis->dim; d++)
9749d007619Sjeremylt       for (CeedInt i=0; i<basis->dim; i++)
9759d007619Sjeremylt         for (CeedInt qpt=0; qpt<basis->Q; qpt++)
9769d007619Sjeremylt           for (CeedInt node=0; node<basis->P; node++) {
9779d007619Sjeremylt             CeedInt p = (node / CeedIntPow(basis->P1d, d)) % basis->P1d;
9789d007619Sjeremylt             CeedInt q = (qpt / CeedIntPow(basis->Q1d, d)) % basis->Q1d;
9799d007619Sjeremylt             if (i == d)
9809d007619Sjeremylt               basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *=
9819d007619Sjeremylt                 basis->grad1d[q*basis->P1d+p];
9829d007619Sjeremylt             else
9839d007619Sjeremylt               basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *=
9849d007619Sjeremylt                 basis->interp1d[q*basis->P1d+p];
9859d007619Sjeremylt           }
9869d007619Sjeremylt   }
9879d007619Sjeremylt   *grad = basis->grad;
988*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9899d007619Sjeremylt }
9909d007619Sjeremylt 
9919d007619Sjeremylt /**
9929d007619Sjeremylt   @brief Get 1D gradient matrix of a tensor product CeedBasis
9939d007619Sjeremylt 
9949d007619Sjeremylt   @param basis        CeedBasis
9959d007619Sjeremylt   @param[out] grad1d  Variable to store gradient matrix
9969d007619Sjeremylt 
9979d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
9989d007619Sjeremylt 
9999d007619Sjeremylt   @ref Backend
10009d007619Sjeremylt **/
10016c58de82SJeremy L Thompson int CeedBasisGetGrad1D(CeedBasis basis, const CeedScalar **grad1d) {
10029d007619Sjeremylt   if (!basis->tensorbasis)
10039d007619Sjeremylt     // LCOV_EXCL_START
1004*e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR,
1005*e15f9bd0SJeremy L Thompson                      "CeedBasis is not a tensor product basis.");
10069d007619Sjeremylt   // LCOV_EXCL_STOP
10079d007619Sjeremylt 
10089d007619Sjeremylt   *grad1d = basis->grad1d;
1009*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
10109d007619Sjeremylt }
10119d007619Sjeremylt 
10129d007619Sjeremylt /**
10137a982d89SJeremy L. Thompson   @brief Destroy a CeedBasis
10147a982d89SJeremy L. Thompson 
10157a982d89SJeremy L. Thompson   @param basis CeedBasis to destroy
10167a982d89SJeremy L. Thompson 
10177a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
10187a982d89SJeremy L. Thompson 
10197a982d89SJeremy L. Thompson   @ref User
10207a982d89SJeremy L. Thompson **/
10217a982d89SJeremy L. Thompson int CeedBasisDestroy(CeedBasis *basis) {
10227a982d89SJeremy L. Thompson   int ierr;
10237a982d89SJeremy L. Thompson 
1024*e15f9bd0SJeremy L Thompson   if (!*basis || --(*basis)->refcount > 0) return CEED_ERROR_SUCCESS;
10257a982d89SJeremy L. Thompson   if ((*basis)->Destroy) {
10267a982d89SJeremy L. Thompson     ierr = (*basis)->Destroy(*basis); CeedChk(ierr);
10277a982d89SJeremy L. Thompson   }
10287a982d89SJeremy L. Thompson   ierr = CeedFree(&(*basis)->interp); CeedChk(ierr);
10297a982d89SJeremy L. Thompson   ierr = CeedFree(&(*basis)->interp1d); CeedChk(ierr);
10307a982d89SJeremy L. Thompson   ierr = CeedFree(&(*basis)->grad); CeedChk(ierr);
10317a982d89SJeremy L. Thompson   ierr = CeedFree(&(*basis)->grad1d); CeedChk(ierr);
10327a982d89SJeremy L. Thompson   ierr = CeedFree(&(*basis)->qref1d); CeedChk(ierr);
10337a982d89SJeremy L. Thompson   ierr = CeedFree(&(*basis)->qweight1d); CeedChk(ierr);
10347a982d89SJeremy L. Thompson   ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr);
10357a982d89SJeremy L. Thompson   ierr = CeedFree(basis); CeedChk(ierr);
1036*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
10377a982d89SJeremy L. Thompson }
10387a982d89SJeremy L. Thompson 
10397a982d89SJeremy L. Thompson /**
1040b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre quadrature
1041b11c1e72Sjeremylt 
1042b11c1e72Sjeremylt   @param Q               Number of quadrature points (integrates polynomials of
1043b11c1e72Sjeremylt                            degree 2*Q-1 exactly)
1044b11c1e72Sjeremylt   @param[out] qref1d     Array of length Q to hold the abscissa on [-1, 1]
1045b11c1e72Sjeremylt   @param[out] qweight1d  Array of length Q to hold the weights
1046b11c1e72Sjeremylt 
1047b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1048dfdf5a53Sjeremylt 
1049dfdf5a53Sjeremylt   @ref Utility
1050b11c1e72Sjeremylt **/
1051d7b241e6Sjeremylt int CeedGaussQuadrature(CeedInt Q, CeedScalar *qref1d, CeedScalar *qweight1d) {
1052d7b241e6Sjeremylt   // Allocate
1053d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0*atan(1.0);
1054d7b241e6Sjeremylt   // Build qref1d, qweight1d
1055d7b241e6Sjeremylt   for (int i = 0; i <= Q/2; i++) {
1056d7b241e6Sjeremylt     // Guess
1057d7b241e6Sjeremylt     xi = cos(PI*(CeedScalar)(2*i+1)/((CeedScalar)(2*Q)));
1058d7b241e6Sjeremylt     // Pn(xi)
1059d7b241e6Sjeremylt     P0 = 1.0;
1060d7b241e6Sjeremylt     P1 = xi;
1061d7b241e6Sjeremylt     P2 = 0.0;
1062d7b241e6Sjeremylt     for (int j = 2; j <= Q; j++) {
1063d7b241e6Sjeremylt       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1064d7b241e6Sjeremylt       P0 = P1;
1065d7b241e6Sjeremylt       P1 = P2;
1066d7b241e6Sjeremylt     }
1067d7b241e6Sjeremylt     // First Newton Step
1068d7b241e6Sjeremylt     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1069d7b241e6Sjeremylt     xi = xi-P2/dP2;
1070d7b241e6Sjeremylt     // Newton to convergence
10710e4d4210Sjeremylt     for (int k=0; k<100 && fabs(P2)>10*CEED_EPSILON; k++) {
1072d7b241e6Sjeremylt       P0 = 1.0;
1073d7b241e6Sjeremylt       P1 = xi;
1074d7b241e6Sjeremylt       for (int j = 2; j <= Q; j++) {
1075d7b241e6Sjeremylt         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1076d7b241e6Sjeremylt         P0 = P1;
1077d7b241e6Sjeremylt         P1 = P2;
1078d7b241e6Sjeremylt       }
1079d7b241e6Sjeremylt       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1080d7b241e6Sjeremylt       xi = xi-P2/dP2;
1081d7b241e6Sjeremylt     }
1082d7b241e6Sjeremylt     // Save xi, wi
1083d7b241e6Sjeremylt     wi = 2.0/((1.0-xi*xi)*dP2*dP2);
1084d7b241e6Sjeremylt     qweight1d[i] = wi;
1085d7b241e6Sjeremylt     qweight1d[Q-1-i] = wi;
1086d7b241e6Sjeremylt     qref1d[i] = -xi;
1087d7b241e6Sjeremylt     qref1d[Q-1-i]= xi;
1088d7b241e6Sjeremylt   }
1089*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1090d7b241e6Sjeremylt }
1091d7b241e6Sjeremylt 
1092b11c1e72Sjeremylt /**
1093b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre-Lobatto quadrature
1094b11c1e72Sjeremylt 
1095b11c1e72Sjeremylt   @param Q               Number of quadrature points (integrates polynomials of
1096b11c1e72Sjeremylt                            degree 2*Q-3 exactly)
1097b11c1e72Sjeremylt   @param[out] qref1d     Array of length Q to hold the abscissa on [-1, 1]
1098b11c1e72Sjeremylt   @param[out] qweight1d  Array of length Q to hold the weights
1099b11c1e72Sjeremylt 
1100b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1101dfdf5a53Sjeremylt 
1102dfdf5a53Sjeremylt   @ref Utility
1103b11c1e72Sjeremylt **/
1104d7b241e6Sjeremylt int CeedLobattoQuadrature(CeedInt Q, CeedScalar *qref1d,
1105d7b241e6Sjeremylt                           CeedScalar *qweight1d) {
1106d7b241e6Sjeremylt   // Allocate
1107d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0*atan(1.0);
1108d7b241e6Sjeremylt   // Build qref1d, qweight1d
1109d7b241e6Sjeremylt   // Set endpoints
111030a100c3SJed Brown   if (Q < 2)
1111b0d62198Sjeremylt     // LCOV_EXCL_START
1112*e15f9bd0SJeremy L Thompson     return CeedError(NULL, CEED_ERROR_DIMENSION,
11137ed52d01Sjeremylt                      "Cannot create Lobatto quadrature with Q=%d < 2 points", Q);
1114b0d62198Sjeremylt   // LCOV_EXCL_STOP
1115d7b241e6Sjeremylt   wi = 2.0/((CeedScalar)(Q*(Q-1)));
1116d7b241e6Sjeremylt   if (qweight1d) {
1117d7b241e6Sjeremylt     qweight1d[0] = wi;
1118d7b241e6Sjeremylt     qweight1d[Q-1] = wi;
1119d7b241e6Sjeremylt   }
1120d7b241e6Sjeremylt   qref1d[0] = -1.0;
1121d7b241e6Sjeremylt   qref1d[Q-1] = 1.0;
1122d7b241e6Sjeremylt   // Interior
1123d7b241e6Sjeremylt   for (int i = 1; i <= (Q-1)/2; i++) {
1124d7b241e6Sjeremylt     // Guess
1125d7b241e6Sjeremylt     xi = cos(PI*(CeedScalar)(i)/(CeedScalar)(Q-1));
1126d7b241e6Sjeremylt     // Pn(xi)
1127d7b241e6Sjeremylt     P0 = 1.0;
1128d7b241e6Sjeremylt     P1 = xi;
1129d7b241e6Sjeremylt     P2 = 0.0;
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     // First Newton step
1136d7b241e6Sjeremylt     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1137d7b241e6Sjeremylt     d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
1138d7b241e6Sjeremylt     xi = xi-dP2/d2P2;
1139d7b241e6Sjeremylt     // Newton to convergence
11400e4d4210Sjeremylt     for (int k=0; k<100 && fabs(dP2)>10*CEED_EPSILON; k++) {
1141d7b241e6Sjeremylt       P0 = 1.0;
1142d7b241e6Sjeremylt       P1 = xi;
1143d7b241e6Sjeremylt       for (int j = 2; j < Q; j++) {
1144d7b241e6Sjeremylt         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1145d7b241e6Sjeremylt         P0 = P1;
1146d7b241e6Sjeremylt         P1 = P2;
1147d7b241e6Sjeremylt       }
1148d7b241e6Sjeremylt       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1149d7b241e6Sjeremylt       d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
1150d7b241e6Sjeremylt       xi = xi-dP2/d2P2;
1151d7b241e6Sjeremylt     }
1152d7b241e6Sjeremylt     // Save xi, wi
1153d7b241e6Sjeremylt     wi = 2.0/(((CeedScalar)(Q*(Q-1)))*P2*P2);
1154d7b241e6Sjeremylt     if (qweight1d) {
1155d7b241e6Sjeremylt       qweight1d[i] = wi;
1156d7b241e6Sjeremylt       qweight1d[Q-1-i] = wi;
1157d7b241e6Sjeremylt     }
1158d7b241e6Sjeremylt     qref1d[i] = -xi;
1159d7b241e6Sjeremylt     qref1d[Q-1-i]= xi;
1160d7b241e6Sjeremylt   }
1161*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1162d7b241e6Sjeremylt }
1163d7b241e6Sjeremylt 
1164dfdf5a53Sjeremylt /**
116595bb1877Svaleriabarra   @brief Return QR Factorization of a matrix
1166b11c1e72Sjeremylt 
116777645d7bSjeremylt   @param ceed         A Ceed context for error handling
116852bfb9bbSJeremy L Thompson   @param[in,out] mat  Row-major matrix to be factorized in place
116952bfb9bbSJeremy L Thompson   @param[in,out] tau  Vector of length m of scaling factors
1170b11c1e72Sjeremylt   @param m            Number of rows
1171b11c1e72Sjeremylt   @param n            Number of columns
1172b11c1e72Sjeremylt 
1173b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1174dfdf5a53Sjeremylt 
1175dfdf5a53Sjeremylt   @ref Utility
1176b11c1e72Sjeremylt **/
1177a7bd39daSjeremylt int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau,
1178d7b241e6Sjeremylt                         CeedInt m, CeedInt n) {
1179d7b241e6Sjeremylt   CeedScalar v[m];
1180d7b241e6Sjeremylt 
1181a7bd39daSjeremylt   // Check m >= n
1182a7bd39daSjeremylt   if (n > m)
1183c042f62fSJeremy L Thompson     // LCOV_EXCL_START
1184*e15f9bd0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_UNSUPPORTED,
1185*e15f9bd0SJeremy L Thompson                      "Cannot compute QR factorization with n > m");
1186c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
1187a7bd39daSjeremylt 
118852bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++) {
1189d7b241e6Sjeremylt     // Calculate Householder vector, magnitude
1190d7b241e6Sjeremylt     CeedScalar sigma = 0.0;
1191d7b241e6Sjeremylt     v[i] = mat[i+n*i];
119252bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<m; j++) {
1193d7b241e6Sjeremylt       v[j] = mat[i+n*j];
1194d7b241e6Sjeremylt       sigma += v[j] * v[j];
1195d7b241e6Sjeremylt     }
1196d7b241e6Sjeremylt     CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m]
1197d7b241e6Sjeremylt     CeedScalar Rii = -copysign(norm, v[i]);
1198d7b241e6Sjeremylt     v[i] -= Rii;
1199d7b241e6Sjeremylt     // norm of v[i:m] after modification above and scaling below
1200d7b241e6Sjeremylt     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
1201d7b241e6Sjeremylt     //   tau = 2 / (norm*norm)
1202d7b241e6Sjeremylt     tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma);
1203fb551037Sjeremylt 
12041d102b48SJeremy L Thompson     for (CeedInt j=i+1; j<m; j++)
12051d102b48SJeremy L Thompson       v[j] /= v[i];
1206d7b241e6Sjeremylt 
1207d7b241e6Sjeremylt     // Apply Householder reflector to lower right panel
1208d7b241e6Sjeremylt     CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1);
1209d7b241e6Sjeremylt     // Save v
1210d7b241e6Sjeremylt     mat[i+n*i] = Rii;
12111d102b48SJeremy L Thompson     for (CeedInt j=i+1; j<m; j++)
1212d7b241e6Sjeremylt       mat[i+n*j] = v[j];
1213d7b241e6Sjeremylt   }
1214*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1215d7b241e6Sjeremylt }
1216d7b241e6Sjeremylt 
1217b11c1e72Sjeremylt /**
121852bfb9bbSJeremy L Thompson   @brief Return symmetric Schur decomposition of the symmetric matrix mat via
121952bfb9bbSJeremy L Thompson            symmetric QR factorization
122052bfb9bbSJeremy L Thompson 
122177645d7bSjeremylt   @param ceed         A Ceed context for error handling
122252bfb9bbSJeremy L Thompson   @param[in,out] mat  Row-major matrix to be factorized in place
1223460bf743SValeria Barra   @param[out] lambda  Vector of length n of eigenvalues
122452bfb9bbSJeremy L Thompson   @param n            Number of rows/columns
122552bfb9bbSJeremy L Thompson 
122652bfb9bbSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
122752bfb9bbSJeremy L Thompson 
122852bfb9bbSJeremy L Thompson   @ref Utility
122952bfb9bbSJeremy L Thompson **/
123052bfb9bbSJeremy L Thompson int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat,
123152bfb9bbSJeremy L Thompson                                     CeedScalar *lambda, CeedInt n) {
123252bfb9bbSJeremy L Thompson   // Check bounds for clang-tidy
123352bfb9bbSJeremy L Thompson   if (n<2)
1234c042f62fSJeremy L Thompson     // LCOV_EXCL_START
1235*e15f9bd0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_UNSUPPORTED,
1236c042f62fSJeremy L Thompson                      "Cannot compute symmetric Schur decomposition of scalars");
1237c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
123852bfb9bbSJeremy L Thompson 
123952bfb9bbSJeremy L Thompson   CeedScalar v[n-1], tau[n-1], matT[n*n];
124052bfb9bbSJeremy L Thompson 
124152bfb9bbSJeremy L Thompson   // Copy mat to matT and set mat to I
124252bfb9bbSJeremy L Thompson   memcpy(matT, mat, n*n*sizeof(mat[0]));
124352bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
124452bfb9bbSJeremy L Thompson     for (CeedInt j=0; j<n; j++)
124552bfb9bbSJeremy L Thompson       mat[j+n*i] = (i==j) ? 1 : 0;
124652bfb9bbSJeremy L Thompson 
124752bfb9bbSJeremy L Thompson   // Reduce to tridiagonal
124852bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n-1; i++) {
124952bfb9bbSJeremy L Thompson     // Calculate Householder vector, magnitude
125052bfb9bbSJeremy L Thompson     CeedScalar sigma = 0.0;
125152bfb9bbSJeremy L Thompson     v[i] = matT[i+n*(i+1)];
125252bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) {
125352bfb9bbSJeremy L Thompson       v[j] = matT[i+n*(j+1)];
125452bfb9bbSJeremy L Thompson       sigma += v[j] * v[j];
125552bfb9bbSJeremy L Thompson     }
125652bfb9bbSJeremy L Thompson     CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:n-1]
125752bfb9bbSJeremy L Thompson     CeedScalar Rii = -copysign(norm, v[i]);
125852bfb9bbSJeremy L Thompson     v[i] -= Rii;
125952bfb9bbSJeremy L Thompson     // norm of v[i:m] after modification above and scaling below
126052bfb9bbSJeremy L Thompson     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
126152bfb9bbSJeremy L Thompson     //   tau = 2 / (norm*norm)
12620e4d4210Sjeremylt     if (sigma > 10*CEED_EPSILON)
126352bfb9bbSJeremy L Thompson       tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma);
1264fb551037Sjeremylt     else
1265fb551037Sjeremylt       tau[i] = 0;
1266fb551037Sjeremylt 
1267fb551037Sjeremylt     for (CeedInt j=i+1; j<n-1; j++)
1268fb551037Sjeremylt       v[j] /= v[i];
126952bfb9bbSJeremy L Thompson 
127052bfb9bbSJeremy L Thompson     // Update sub and super diagonal
127152bfb9bbSJeremy L Thompson     matT[i+n*(i+1)] = Rii;
127252bfb9bbSJeremy L Thompson     matT[(i+1)+n*i] = Rii;
127352bfb9bbSJeremy L Thompson     for (CeedInt j=i+2; j<n; j++) {
127452bfb9bbSJeremy L Thompson       matT[i+n*j] = 0; matT[j+n*i] = 0;
127552bfb9bbSJeremy L Thompson     }
127652bfb9bbSJeremy L Thompson     // Apply symmetric Householder reflector to lower right panel
127752bfb9bbSJeremy L Thompson     CeedHouseholderReflect(&matT[(i+1)+n*(i+1)], &v[i], tau[i],
127852bfb9bbSJeremy L Thompson                            n-(i+1), n-(i+1), n, 1);
127952bfb9bbSJeremy L Thompson     CeedHouseholderReflect(&matT[(i+1)+n*(i+1)], &v[i], tau[i],
128052bfb9bbSJeremy L Thompson                            n-(i+1), n-(i+1), 1, n);
128152bfb9bbSJeremy L Thompson     // Save v
128252bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) {
128352bfb9bbSJeremy L Thompson       matT[i+n*(j+1)] = v[j];
128452bfb9bbSJeremy L Thompson     }
128552bfb9bbSJeremy L Thompson   }
128652bfb9bbSJeremy L Thompson   // Backwards accumulation of Q
128752bfb9bbSJeremy L Thompson   for (CeedInt i=n-2; i>=0; i--) {
128852bfb9bbSJeremy L Thompson     v[i] = 1;
128952bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) {
129052bfb9bbSJeremy L Thompson       v[j] = matT[i+n*(j+1)];
129152bfb9bbSJeremy L Thompson       matT[i+n*(j+1)] = 0;
129252bfb9bbSJeremy L Thompson     }
129352bfb9bbSJeremy L Thompson     CeedHouseholderReflect(&mat[(i+1)+n*(i+1)], &v[i], tau[i],
129452bfb9bbSJeremy L Thompson                            n-(i+1), n-(i+1), n, 1);
129552bfb9bbSJeremy L Thompson   }
129652bfb9bbSJeremy L Thompson 
129752bfb9bbSJeremy L Thompson   // Reduce sub and super diagonal
129852bfb9bbSJeremy L Thompson   CeedInt p = 0, q = 0, itr = 0, maxitr = n*n*n;
12990e4d4210Sjeremylt   CeedScalar tol = 10*CEED_EPSILON;
130052bfb9bbSJeremy L Thompson 
130152bfb9bbSJeremy L Thompson   while (q < n && itr < maxitr) {
130252bfb9bbSJeremy L Thompson     // Update p, q, size of reduced portions of diagonal
130352bfb9bbSJeremy L Thompson     p = 0; q = 0;
130452bfb9bbSJeremy L Thompson     for (CeedInt i=n-2; i>=0; i--) {
130552bfb9bbSJeremy L Thompson       if (fabs(matT[i+n*(i+1)]) < tol)
130652bfb9bbSJeremy L Thompson         q += 1;
130752bfb9bbSJeremy L Thompson       else
130852bfb9bbSJeremy L Thompson         break;
130952bfb9bbSJeremy L Thompson     }
131052bfb9bbSJeremy L Thompson     for (CeedInt i=0; i<n-1-q; i++) {
131152bfb9bbSJeremy L Thompson       if (fabs(matT[i+n*(i+1)]) < tol)
131252bfb9bbSJeremy L Thompson         p += 1;
131352bfb9bbSJeremy L Thompson       else
131452bfb9bbSJeremy L Thompson         break;
131552bfb9bbSJeremy L Thompson     }
131652bfb9bbSJeremy L Thompson     if (q == n-1) break; // Finished reducing
131752bfb9bbSJeremy L Thompson 
131852bfb9bbSJeremy L Thompson     // Reduce tridiagonal portion
131952bfb9bbSJeremy L Thompson     CeedScalar tnn = matT[(n-1-q)+n*(n-1-q)],
132052bfb9bbSJeremy L Thompson                tnnm1 = matT[(n-2-q)+n*(n-1-q)];
132152bfb9bbSJeremy L Thompson     CeedScalar d = (matT[(n-2-q)+n*(n-2-q)] - tnn)/2;
132252bfb9bbSJeremy L Thompson     CeedScalar mu = tnn - tnnm1*tnnm1 /
132352bfb9bbSJeremy L Thompson                     (d + copysign(sqrt(d*d + tnnm1*tnnm1), d));
132452bfb9bbSJeremy L Thompson     CeedScalar x = matT[p+n*p] - mu;
132552bfb9bbSJeremy L Thompson     CeedScalar z = matT[p+n*(p+1)];
132652bfb9bbSJeremy L Thompson     for (CeedInt k=p; k<n-1-q; k++) {
132752bfb9bbSJeremy L Thompson       // Compute Givens rotation
132852bfb9bbSJeremy L Thompson       CeedScalar c = 1, s = 0;
132952bfb9bbSJeremy L Thompson       if (fabs(z) > tol) {
133052bfb9bbSJeremy L Thompson         if (fabs(z) > fabs(x)) {
133152bfb9bbSJeremy L Thompson           CeedScalar tau = -x/z;
133252bfb9bbSJeremy L Thompson           s = 1/sqrt(1+tau*tau), c = s*tau;
133352bfb9bbSJeremy L Thompson         } else {
133452bfb9bbSJeremy L Thompson           CeedScalar tau = -z/x;
133552bfb9bbSJeremy L Thompson           c = 1/sqrt(1+tau*tau), s = c*tau;
133652bfb9bbSJeremy L Thompson         }
133752bfb9bbSJeremy L Thompson       }
133852bfb9bbSJeremy L Thompson 
133952bfb9bbSJeremy L Thompson       // Apply Givens rotation to T
134052bfb9bbSJeremy L Thompson       CeedGivensRotation(matT, c, s, CEED_NOTRANSPOSE, k, k+1, n, n);
134152bfb9bbSJeremy L Thompson       CeedGivensRotation(matT, c, s, CEED_TRANSPOSE, k, k+1, n, n);
134252bfb9bbSJeremy L Thompson 
134352bfb9bbSJeremy L Thompson       // Apply Givens rotation to Q
134452bfb9bbSJeremy L Thompson       CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k+1, n, n);
134552bfb9bbSJeremy L Thompson 
134652bfb9bbSJeremy L Thompson       // Update x, z
134752bfb9bbSJeremy L Thompson       if (k < n-q-2) {
134852bfb9bbSJeremy L Thompson         x = matT[k+n*(k+1)];
134952bfb9bbSJeremy L Thompson         z = matT[k+n*(k+2)];
135052bfb9bbSJeremy L Thompson       }
135152bfb9bbSJeremy L Thompson     }
135252bfb9bbSJeremy L Thompson     itr++;
135352bfb9bbSJeremy L Thompson   }
135452bfb9bbSJeremy L Thompson   // Save eigenvalues
135552bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
135652bfb9bbSJeremy L Thompson     lambda[i] = matT[i+n*i];
135752bfb9bbSJeremy L Thompson 
135852bfb9bbSJeremy L Thompson   // Check convergence
135952bfb9bbSJeremy L Thompson   if (itr == maxitr && q < n-1)
1360c042f62fSJeremy L Thompson     // LCOV_EXCL_START
1361*e15f9bd0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_MINOR,
1362*e15f9bd0SJeremy L Thompson                      "Symmetric QR failed to converge");
1363c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
1364*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
136552bfb9bbSJeremy L Thompson }
136652bfb9bbSJeremy L Thompson 
136752bfb9bbSJeremy L Thompson /**
136852bfb9bbSJeremy L Thompson   @brief Return Simultaneous Diagonalization of two matrices. This solves the
136952bfb9bbSJeremy L Thompson            generalized eigenvalue problem A x = lambda B x, where A and B
137052bfb9bbSJeremy L Thompson            are symmetric and B is positive definite. We generate the matrix X
137152bfb9bbSJeremy L Thompson            and vector Lambda such that X^T A X = Lambda and X^T B X = I. This
137252bfb9bbSJeremy L Thompson            is equivalent to the LAPACK routine 'sygv' with TYPE = 1.
137352bfb9bbSJeremy L Thompson 
137477645d7bSjeremylt   @param ceed         A Ceed context for error handling
137552bfb9bbSJeremy L Thompson   @param[in] matA     Row-major matrix to be factorized with eigenvalues
137652bfb9bbSJeremy L Thompson   @param[in] matB     Row-major matrix to be factorized to identity
137752bfb9bbSJeremy L Thompson   @param[out] x       Row-major orthogonal matrix
1378460bf743SValeria Barra   @param[out] lambda  Vector of length n of generalized eigenvalues
137952bfb9bbSJeremy L Thompson   @param n            Number of rows/columns
138052bfb9bbSJeremy L Thompson 
138152bfb9bbSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
138252bfb9bbSJeremy L Thompson 
138352bfb9bbSJeremy L Thompson   @ref Utility
138452bfb9bbSJeremy L Thompson **/
138552bfb9bbSJeremy L Thompson int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *matA,
138652bfb9bbSJeremy L Thompson                                     CeedScalar *matB, CeedScalar *x,
138752bfb9bbSJeremy L Thompson                                     CeedScalar *lambda, CeedInt n) {
138852bfb9bbSJeremy L Thompson   int ierr;
138952bfb9bbSJeremy L Thompson   CeedScalar matC[n*n], matG[n*n], vecD[n];
139052bfb9bbSJeremy L Thompson 
139152bfb9bbSJeremy L Thompson   // Compute B = G D G^T
139252bfb9bbSJeremy L Thompson   memcpy(matG, matB, n*n*sizeof(matB[0]));
139352bfb9bbSJeremy L Thompson   ierr = CeedSymmetricSchurDecomposition(ceed, matG, vecD, n); CeedChk(ierr);
1394fb551037Sjeremylt   for (CeedInt i=0; i<n; i++)
1395fb551037Sjeremylt     vecD[i] = sqrt(vecD[i]);
139652bfb9bbSJeremy L Thompson 
1397fb551037Sjeremylt   // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T
1398fb551037Sjeremylt   //           = D^-1/2 G^T A G D^-1/2
139952bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
140052bfb9bbSJeremy L Thompson     for (CeedInt j=0; j<n; j++)
1401fb551037Sjeremylt       matC[j+i*n] = matG[i+j*n] / vecD[i];
14029289e5bfSjeremylt   ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)matC,
14039289e5bfSjeremylt                             (const CeedScalar *)matA, x, n, n, n);
14049289e5bfSjeremylt   CeedChk(ierr);
140552bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
140652bfb9bbSJeremy L Thompson     for (CeedInt j=0; j<n; j++)
1407fb551037Sjeremylt       matG[j+i*n] = matG[j+i*n] / vecD[j];
14089289e5bfSjeremylt   ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)x,
14099289e5bfSjeremylt                             (const CeedScalar *)matG, matC, n, n, n);
14109289e5bfSjeremylt   CeedChk(ierr);
141152bfb9bbSJeremy L Thompson 
141252bfb9bbSJeremy L Thompson   // Compute Q^T C Q = lambda
141352bfb9bbSJeremy L Thompson   ierr = CeedSymmetricSchurDecomposition(ceed, matC, lambda, n); CeedChk(ierr);
141452bfb9bbSJeremy L Thompson 
1415fb551037Sjeremylt   // Set x = (G D^1/2)^-T Q
1416fb551037Sjeremylt   //       = G D^-1/2 Q
14179289e5bfSjeremylt   ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)matG,
14189289e5bfSjeremylt                             (const CeedScalar *)matC, x, n, n, n);
14199289e5bfSjeremylt   CeedChk(ierr);
1420*e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
142152bfb9bbSJeremy L Thompson }
142252bfb9bbSJeremy L Thompson 
1423d7b241e6Sjeremylt /// @}
1424