xref: /libCEED/rust/libceed-sys/c-src/interface/ceed-basis.c (revision 752c3701a992135134df075f4ef18abc790b3495)
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 
17d7b241e6Sjeremylt #include <ceed-impl.h>
18d863ab9bSjeremylt #include <ceed-backend.h>
19d7b241e6Sjeremylt #include <math.h>
20d7b241e6Sjeremylt #include <stdio.h>
21d7b241e6Sjeremylt #include <stdlib.h>
22d7b241e6Sjeremylt #include <string.h>
23d7b241e6Sjeremylt 
247a982d89SJeremy L. Thompson /// @file
257a982d89SJeremy L. Thompson /// Implementation of CeedBasis interfaces
267a982d89SJeremy L. Thompson 
27d7b241e6Sjeremylt /// @cond DOXYGEN_SKIP
28783c99b3SValeria Barra static struct CeedBasis_private ceed_basis_collocated;
29d7b241e6Sjeremylt /// @endcond
30d7b241e6Sjeremylt 
317a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser
327a982d89SJeremy L. Thompson /// @{
337a982d89SJeremy L. Thompson 
347a982d89SJeremy L. Thompson /// Indicate that the quadrature points are collocated with the nodes
357a982d89SJeremy L. Thompson const CeedBasis CEED_BASIS_COLLOCATED = &ceed_basis_collocated;
367a982d89SJeremy L. Thompson 
377a982d89SJeremy L. Thompson /// @}
387a982d89SJeremy L. Thompson 
397a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
407a982d89SJeremy L. Thompson /// CeedBasis Library Internal Functions
417a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
427a982d89SJeremy L. Thompson /// @addtogroup CeedBasisDeveloper
437a982d89SJeremy L. Thompson /// @{
447a982d89SJeremy L. Thompson 
457a982d89SJeremy L. Thompson /**
467a982d89SJeremy L. Thompson   @brief Compute Householder reflection
477a982d89SJeremy L. Thompson 
487a982d89SJeremy L. Thompson     Computes A = (I - b v v^T) A
497a982d89SJeremy L. Thompson     where A is an mxn matrix indexed as A[i*row + j*col]
507a982d89SJeremy L. Thompson 
517a982d89SJeremy L. Thompson   @param[in,out] A  Matrix to apply Householder reflection to, in place
527a982d89SJeremy L. Thompson   @param v          Householder vector
537a982d89SJeremy L. Thompson   @param b          Scaling factor
547a982d89SJeremy L. Thompson   @param m          Number of rows in A
557a982d89SJeremy L. Thompson   @param n          Number of columns in A
567a982d89SJeremy L. Thompson   @param row        Row stride
577a982d89SJeremy L. Thompson   @param col        Col stride
587a982d89SJeremy L. Thompson 
597a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
607a982d89SJeremy L. Thompson 
617a982d89SJeremy L. Thompson   @ref Developer
627a982d89SJeremy L. Thompson **/
637a982d89SJeremy L. Thompson static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v,
647a982d89SJeremy L. Thompson                                   CeedScalar b, CeedInt m, CeedInt n,
657a982d89SJeremy L. Thompson                                   CeedInt row, CeedInt col) {
667a982d89SJeremy L. Thompson   for (CeedInt j=0; j<n; j++) {
677a982d89SJeremy L. Thompson     CeedScalar w = A[0*row + j*col];
687a982d89SJeremy L. Thompson     for (CeedInt i=1; i<m; i++)
697a982d89SJeremy L. Thompson       w += v[i] * A[i*row + j*col];
707a982d89SJeremy L. Thompson     A[0*row + j*col] -= b * w;
717a982d89SJeremy L. Thompson     for (CeedInt i=1; i<m; i++)
727a982d89SJeremy L. Thompson       A[i*row + j*col] -= b * w * v[i];
737a982d89SJeremy L. Thompson   }
747a982d89SJeremy L. Thompson   return 0;
757a982d89SJeremy L. Thompson }
767a982d89SJeremy L. Thompson 
777a982d89SJeremy L. Thompson /**
787a982d89SJeremy L. Thompson   @brief Apply Householder Q matrix
797a982d89SJeremy L. Thompson 
807a982d89SJeremy L. Thompson     Compute A = Q A where Q is mxm and A is mxn.
817a982d89SJeremy L. Thompson 
827a982d89SJeremy L. Thompson   @param[in,out] A  Matrix to apply Householder Q to, in place
837a982d89SJeremy L. Thompson   @param Q          Householder Q matrix
847a982d89SJeremy L. Thompson   @param tau        Householder scaling factors
857a982d89SJeremy L. Thompson   @param tmode      Transpose mode for application
867a982d89SJeremy L. Thompson   @param m          Number of rows in A
877a982d89SJeremy L. Thompson   @param n          Number of columns in A
887a982d89SJeremy L. Thompson   @param k          Number of elementary reflectors in Q, k<m
897a982d89SJeremy L. Thompson   @param row        Row stride in A
907a982d89SJeremy L. Thompson   @param col        Col stride in A
917a982d89SJeremy L. Thompson 
927a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
937a982d89SJeremy L. Thompson 
947a982d89SJeremy L. Thompson   @ref Developer
957a982d89SJeremy L. Thompson **/
96d99fa3c5SJeremy L Thompson int CeedHouseholderApplyQ(CeedScalar *A, const CeedScalar *Q,
977a982d89SJeremy L. Thompson                           const CeedScalar *tau, CeedTransposeMode tmode,
987a982d89SJeremy L. Thompson                           CeedInt m, CeedInt n, CeedInt k,
997a982d89SJeremy L. Thompson                           CeedInt row, CeedInt col) {
1007a982d89SJeremy L. Thompson   CeedScalar v[m];
1017a982d89SJeremy L. Thompson   for (CeedInt ii=0; ii<k; ii++) {
1027a982d89SJeremy L. Thompson     CeedInt i = tmode == CEED_TRANSPOSE ? ii : k-1-ii;
1037a982d89SJeremy L. Thompson     for (CeedInt j=i+1; j<m; j++)
1047a982d89SJeremy L. Thompson       v[j] = Q[j*k+i];
1057a982d89SJeremy L. Thompson     // Apply Householder reflector (I - tau v v^T) collograd1d^T
1067a982d89SJeremy L. Thompson     CeedHouseholderReflect(&A[i*row], &v[i], tau[i], m-i, n, row, col);
1077a982d89SJeremy L. Thompson   }
1087a982d89SJeremy L. Thompson   return 0;
1097a982d89SJeremy L. Thompson }
1107a982d89SJeremy L. Thompson 
1117a982d89SJeremy L. Thompson /**
1127a982d89SJeremy L. Thompson   @brief Compute Givens rotation
1137a982d89SJeremy L. Thompson 
1147a982d89SJeremy L. Thompson     Computes A = G A (or G^T A in transpose mode)
1157a982d89SJeremy L. Thompson     where A is an mxn matrix indexed as A[i*n + j*m]
1167a982d89SJeremy L. Thompson 
1177a982d89SJeremy L. Thompson   @param[in,out] A  Row major matrix to apply Givens rotation to, in place
1187a982d89SJeremy L. Thompson   @param c          Cosine factor
1197a982d89SJeremy L. Thompson   @param s          Sine factor
1204cc79fe7SJed Brown   @param tmode      @ref CEED_NOTRANSPOSE to rotate the basis counter-clockwise,
1214c4400c7SValeria Barra                     which has the effect of rotating columns of A clockwise;
1224cc79fe7SJed Brown                     @ref CEED_TRANSPOSE for the opposite rotation
1237a982d89SJeremy L. Thompson   @param i          First row/column to apply rotation
1247a982d89SJeremy L. Thompson   @param k          Second row/column to apply rotation
1257a982d89SJeremy L. Thompson   @param m          Number of rows in A
1267a982d89SJeremy L. Thompson   @param n          Number of columns in A
1277a982d89SJeremy L. Thompson 
1287a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
1297a982d89SJeremy L. Thompson 
1307a982d89SJeremy L. Thompson   @ref Developer
1317a982d89SJeremy L. Thompson **/
1327a982d89SJeremy L. Thompson static int CeedGivensRotation(CeedScalar *A, CeedScalar c, CeedScalar s,
1337a982d89SJeremy L. Thompson                               CeedTransposeMode tmode, CeedInt i, CeedInt k,
1347a982d89SJeremy L. Thompson                               CeedInt m, CeedInt n) {
1357a982d89SJeremy L. Thompson   CeedInt stridej = 1, strideik = m, numits = n;
1367a982d89SJeremy L. Thompson   if (tmode == CEED_NOTRANSPOSE) {
1377a982d89SJeremy L. Thompson     stridej = n; strideik = 1; numits = m;
1387a982d89SJeremy L. Thompson   }
1397a982d89SJeremy L. Thompson 
1407a982d89SJeremy L. Thompson   // Apply rotation
1417a982d89SJeremy L. Thompson   for (CeedInt j=0; j<numits; j++) {
1427a982d89SJeremy L. Thompson     CeedScalar tau1 = A[i*strideik+j*stridej], tau2 = A[k*strideik+j*stridej];
1437a982d89SJeremy L. Thompson     A[i*strideik+j*stridej] = c*tau1 - s*tau2;
1447a982d89SJeremy L. Thompson     A[k*strideik+j*stridej] = s*tau1 + c*tau2;
1457a982d89SJeremy L. Thompson   }
1467a982d89SJeremy L. Thompson 
1477a982d89SJeremy L. Thompson   return 0;
1487a982d89SJeremy L. Thompson }
1497a982d89SJeremy L. Thompson 
1507a982d89SJeremy L. Thompson /**
1517a982d89SJeremy L. Thompson   @brief View an array stored in a CeedBasis
1527a982d89SJeremy L. Thompson 
1530a0da059Sjeremylt   @param[in] name      Name of array
1540a0da059Sjeremylt   @param[in] fpformat  Printing format
1550a0da059Sjeremylt   @param[in] m         Number of rows in array
1560a0da059Sjeremylt   @param[in] n         Number of columns in array
1570a0da059Sjeremylt   @param[in] a         Array to be viewed
1580a0da059Sjeremylt   @param[in] stream    Stream to view to, e.g., stdout
1597a982d89SJeremy L. Thompson 
1607a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
1617a982d89SJeremy L. Thompson 
1627a982d89SJeremy L. Thompson   @ref Developer
1637a982d89SJeremy L. Thompson **/
1647a982d89SJeremy L. Thompson static int CeedScalarView(const char *name, const char *fpformat, CeedInt m,
1657a982d89SJeremy L. Thompson                           CeedInt n, const CeedScalar *a, FILE *stream) {
1667a982d89SJeremy L. Thompson   for (int i=0; i<m; i++) {
1677a982d89SJeremy L. Thompson     if (m > 1)
1687a982d89SJeremy L. Thompson       fprintf(stream, "%12s[%d]:", name, i);
1697a982d89SJeremy L. Thompson     else
1707a982d89SJeremy L. Thompson       fprintf(stream, "%12s:", name);
1717a982d89SJeremy L. Thompson     for (int j=0; j<n; j++)
1727a982d89SJeremy L. Thompson       fprintf(stream, fpformat, fabs(a[i*n+j]) > 1E-14 ? a[i*n+j] : 0);
1737a982d89SJeremy L. Thompson     fputs("\n", stream);
1747a982d89SJeremy L. Thompson   }
1757a982d89SJeremy L. Thompson   return 0;
1767a982d89SJeremy L. Thompson }
1777a982d89SJeremy L. Thompson 
1787a982d89SJeremy L. Thompson /// @}
1797a982d89SJeremy L. Thompson 
1807a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
1817a982d89SJeremy L. Thompson /// Ceed Backend API
1827a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
1837a982d89SJeremy L. Thompson /// @addtogroup CeedBasisBackend
1847a982d89SJeremy L. Thompson /// @{
1857a982d89SJeremy L. Thompson 
1867a982d89SJeremy L. Thompson /**
1877a982d89SJeremy L. Thompson   @brief Return collocated grad matrix
1887a982d89SJeremy L. Thompson 
1897a982d89SJeremy L. Thompson   @param basis             CeedBasis
1907a982d89SJeremy L. Thompson   @param[out] collograd1d  Row-major (Q1d * Q1d) matrix expressing derivatives of
1917a982d89SJeremy L. Thompson                             basis functions at quadrature points
1927a982d89SJeremy L. Thompson 
1937a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
1947a982d89SJeremy L. Thompson 
1957a982d89SJeremy L. Thompson   @ref Backend
1967a982d89SJeremy L. Thompson **/
1977a982d89SJeremy L. Thompson int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *collograd1d) {
1987a982d89SJeremy L. Thompson   int i, j, k;
1997a982d89SJeremy L. Thompson   Ceed ceed;
2007a982d89SJeremy L. Thompson   CeedInt ierr, P1d=(basis)->P1d, Q1d=(basis)->Q1d;
2017a982d89SJeremy L. Thompson   CeedScalar *interp1d, *grad1d, tau[Q1d];
2027a982d89SJeremy L. Thompson 
2037a982d89SJeremy L. Thompson   ierr = CeedMalloc(Q1d*P1d, &interp1d); CeedChk(ierr);
2047a982d89SJeremy L. Thompson   ierr = CeedMalloc(Q1d*P1d, &grad1d); CeedChk(ierr);
2057a982d89SJeremy L. Thompson   memcpy(interp1d, (basis)->interp1d, Q1d*P1d*sizeof(basis)->interp1d[0]);
2067a982d89SJeremy L. Thompson   memcpy(grad1d, (basis)->grad1d, Q1d*P1d*sizeof(basis)->interp1d[0]);
2077a982d89SJeremy L. Thompson 
2087a982d89SJeremy L. Thompson   // QR Factorization, interp1d = Q R
2097a982d89SJeremy L. Thompson   ierr = CeedBasisGetCeed(basis, &ceed); CeedChk(ierr);
2107a982d89SJeremy L. Thompson   ierr = CeedQRFactorization(ceed, interp1d, tau, Q1d, P1d); CeedChk(ierr);
2117a982d89SJeremy L. Thompson 
2127a982d89SJeremy L. Thompson   // Apply Rinv, collograd1d = grad1d Rinv
2137a982d89SJeremy L. Thompson   for (i=0; i<Q1d; i++) { // Row i
2147a982d89SJeremy L. Thompson     collograd1d[Q1d*i] = grad1d[P1d*i]/interp1d[0];
2157a982d89SJeremy L. Thompson     for (j=1; j<P1d; j++) { // Column j
2167a982d89SJeremy L. Thompson       collograd1d[j+Q1d*i] = grad1d[j+P1d*i];
2177a982d89SJeremy L. Thompson       for (k=0; k<j; k++)
2187a982d89SJeremy L. Thompson         collograd1d[j+Q1d*i] -= interp1d[j+P1d*k]*collograd1d[k+Q1d*i];
2197a982d89SJeremy L. Thompson       collograd1d[j+Q1d*i] /= interp1d[j+P1d*j];
2207a982d89SJeremy L. Thompson     }
2217a982d89SJeremy L. Thompson     for (j=P1d; j<Q1d; j++)
2227a982d89SJeremy L. Thompson       collograd1d[j+Q1d*i] = 0;
2237a982d89SJeremy L. Thompson   }
2247a982d89SJeremy L. Thompson 
2257a982d89SJeremy L. Thompson   // Apply Qtranspose, collograd = collograd Qtranspose
2267a982d89SJeremy L. Thompson   CeedHouseholderApplyQ(collograd1d, interp1d, tau, CEED_NOTRANSPOSE,
2277a982d89SJeremy L. Thompson                         Q1d, Q1d, P1d, 1, Q1d);
2287a982d89SJeremy L. Thompson 
2297a982d89SJeremy L. Thompson   ierr = CeedFree(&interp1d); CeedChk(ierr);
2307a982d89SJeremy L. Thompson   ierr = CeedFree(&grad1d); CeedChk(ierr);
2317a982d89SJeremy L. Thompson 
2327a982d89SJeremy L. Thompson   return 0;
2337a982d89SJeremy L. Thompson }
2347a982d89SJeremy L. Thompson 
2357a982d89SJeremy L. Thompson /**
2367a982d89SJeremy L. Thompson   @brief Get Ceed associated with a CeedBasis
2377a982d89SJeremy L. Thompson 
2387a982d89SJeremy L. Thompson   @param basis      CeedBasis
2397a982d89SJeremy L. Thompson   @param[out] ceed  Variable to store Ceed
2407a982d89SJeremy L. Thompson 
2417a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
2427a982d89SJeremy L. Thompson 
2437a982d89SJeremy L. Thompson   @ref Backend
2447a982d89SJeremy L. Thompson **/
2457a982d89SJeremy L. Thompson int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) {
2467a982d89SJeremy L. Thompson   *ceed = basis->ceed;
2477a982d89SJeremy L. Thompson   return 0;
2487a982d89SJeremy L. Thompson }
2497a982d89SJeremy L. Thompson 
2507a982d89SJeremy L. Thompson /**
2517a982d89SJeremy L. Thompson   @brief Get tensor status for given CeedBasis
2527a982d89SJeremy L. Thompson 
2537a982d89SJeremy L. Thompson   @param basis          CeedBasis
254fd364f38SJeremy L Thompson   @param[out] istensor  Variable to store tensor status
2557a982d89SJeremy L. Thompson 
2567a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
2577a982d89SJeremy L. Thompson 
2587a982d89SJeremy L. Thompson   @ref Backend
2597a982d89SJeremy L. Thompson **/
260fd364f38SJeremy L Thompson int CeedBasisIsTensor(CeedBasis basis, bool *istensor) {
261fd364f38SJeremy L Thompson   *istensor = basis->tensorbasis;
2627a982d89SJeremy L. Thompson   return 0;
2637a982d89SJeremy L. Thompson }
2647a982d89SJeremy L. Thompson 
2657a982d89SJeremy L. Thompson /**
2667a982d89SJeremy L. Thompson   @brief Get backend data of a CeedBasis
2677a982d89SJeremy L. Thompson 
2687a982d89SJeremy L. Thompson   @param basis      CeedBasis
2697a982d89SJeremy L. Thompson   @param[out] data  Variable to store data
2707a982d89SJeremy L. Thompson 
2717a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
2727a982d89SJeremy L. Thompson 
2737a982d89SJeremy L. Thompson   @ref Backend
2747a982d89SJeremy L. Thompson **/
2757a982d89SJeremy L. Thompson int CeedBasisGetData(CeedBasis basis, void **data) {
2767a982d89SJeremy L. Thompson   *data = basis->data;
2777a982d89SJeremy L. Thompson   return 0;
2787a982d89SJeremy L. Thompson }
2797a982d89SJeremy L. Thompson 
2807a982d89SJeremy L. Thompson /**
2817a982d89SJeremy L. Thompson   @brief Set backend data of a CeedBasis
2827a982d89SJeremy L. Thompson 
2837a982d89SJeremy L. Thompson   @param[out] basis  CeedBasis
2847a982d89SJeremy L. Thompson   @param data        Data to set
2857a982d89SJeremy L. Thompson 
2867a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
2877a982d89SJeremy L. Thompson 
2887a982d89SJeremy L. Thompson   @ref Backend
2897a982d89SJeremy L. Thompson **/
2907a982d89SJeremy L. Thompson int CeedBasisSetData(CeedBasis basis, void **data) {
2917a982d89SJeremy L. Thompson   basis->data = *data;
2927a982d89SJeremy L. Thompson   return 0;
2937a982d89SJeremy L. Thompson }
2947a982d89SJeremy L. Thompson 
2957a982d89SJeremy L. Thompson /**
2967a982d89SJeremy L. Thompson   @brief Get dimension for given CeedElemTopology
2977a982d89SJeremy L. Thompson 
2987a982d89SJeremy L. Thompson   @param topo      CeedElemTopology
2997a982d89SJeremy L. Thompson   @param[out] dim  Variable to store dimension of topology
3007a982d89SJeremy L. Thompson 
3017a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3027a982d89SJeremy L. Thompson 
3037a982d89SJeremy L. Thompson   @ref Backend
3047a982d89SJeremy L. Thompson **/
3057a982d89SJeremy L. Thompson int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) {
3067a982d89SJeremy L. Thompson   *dim = (CeedInt) topo >> 16;
3077a982d89SJeremy L. Thompson   return 0;
3087a982d89SJeremy L. Thompson }
3097a982d89SJeremy L. Thompson 
3107a982d89SJeremy L. Thompson /**
3117a982d89SJeremy L. Thompson   @brief Get CeedTensorContract of a CeedBasis
3127a982d89SJeremy L. Thompson 
3137a982d89SJeremy L. Thompson   @param basis          CeedBasis
3147a982d89SJeremy L. Thompson   @param[out] contract  Variable to store CeedTensorContract
3157a982d89SJeremy L. Thompson 
3167a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3177a982d89SJeremy L. Thompson 
3187a982d89SJeremy L. Thompson   @ref Backend
3197a982d89SJeremy L. Thompson **/
3207a982d89SJeremy L. Thompson int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) {
3217a982d89SJeremy L. Thompson   *contract = basis->contract;
3227a982d89SJeremy L. Thompson   return 0;
3237a982d89SJeremy L. Thompson }
3247a982d89SJeremy L. Thompson 
3257a982d89SJeremy L. Thompson /**
3267a982d89SJeremy L. Thompson   @brief Set CeedTensorContract of a CeedBasis
3277a982d89SJeremy L. Thompson 
3287a982d89SJeremy L. Thompson   @param[out] basis     CeedBasis
3297a982d89SJeremy L. Thompson   @param contract       CeedTensorContract to set
3307a982d89SJeremy L. Thompson 
3317a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3327a982d89SJeremy L. Thompson 
3337a982d89SJeremy L. Thompson   @ref Backend
3347a982d89SJeremy L. Thompson **/
3357a982d89SJeremy L. Thompson int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract *contract) {
3367a982d89SJeremy L. Thompson   basis->contract = *contract;
3377a982d89SJeremy L. Thompson   return 0;
3387a982d89SJeremy L. Thompson }
3397a982d89SJeremy L. Thompson 
3407a982d89SJeremy L. Thompson /**
3417a982d89SJeremy L. Thompson   @brief Return a reference implementation of matrix multiplication C = A B.
3427a982d89SJeremy L. Thompson            Note, this is a reference implementation for CPU CeedScalar pointers
3437a982d89SJeremy L. Thompson            that is not intended for high performance.
3447a982d89SJeremy L. Thompson 
3457a982d89SJeremy L. Thompson   @param ceed         A Ceed context for error handling
3467a982d89SJeremy L. Thompson   @param[in] matA     Row-major matrix A
3477a982d89SJeremy L. Thompson   @param[in] matB     Row-major matrix B
3487a982d89SJeremy L. Thompson   @param[out] matC    Row-major output matrix C
3497a982d89SJeremy L. Thompson   @param m            Number of rows of C
3507a982d89SJeremy L. Thompson   @param n            Number of columns of C
3517a982d89SJeremy L. Thompson   @param kk           Number of columns of A/rows of B
3527a982d89SJeremy L. Thompson 
3537a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3547a982d89SJeremy L. Thompson 
3557a982d89SJeremy L. Thompson   @ref Utility
3567a982d89SJeremy L. Thompson **/
3577a982d89SJeremy L. Thompson int CeedMatrixMultiply(Ceed ceed, const CeedScalar *matA,
3587a982d89SJeremy L. Thompson                        const CeedScalar *matB, CeedScalar *matC, CeedInt m,
3597a982d89SJeremy L. Thompson                        CeedInt n, CeedInt kk) {
3607a982d89SJeremy L. Thompson   for (CeedInt i=0; i<m; i++)
3617a982d89SJeremy L. Thompson     for (CeedInt j=0; j<n; j++) {
3627a982d89SJeremy L. Thompson       CeedScalar sum = 0;
3637a982d89SJeremy L. Thompson       for (CeedInt k=0; k<kk; k++)
3647a982d89SJeremy L. Thompson         sum += matA[k+i*kk]*matB[j+k*n];
3657a982d89SJeremy L. Thompson       matC[j+i*n] = sum;
3667a982d89SJeremy L. Thompson     }
3677a982d89SJeremy L. Thompson   return 0;
3687a982d89SJeremy L. Thompson }
3697a982d89SJeremy L. Thompson 
3707a982d89SJeremy L. Thompson /// @}
3717a982d89SJeremy L. Thompson 
3727a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
3737a982d89SJeremy L. Thompson /// CeedBasis Public API
3747a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
3757a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser
376d7b241e6Sjeremylt /// @{
377d7b241e6Sjeremylt 
378b11c1e72Sjeremylt /**
37995bb1877Svaleriabarra   @brief Create a tensor-product basis for H^1 discretizations
380b11c1e72Sjeremylt 
381b11c1e72Sjeremylt   @param ceed        A Ceed object where the CeedBasis will be created
382b11c1e72Sjeremylt   @param dim         Topological dimension
383b11c1e72Sjeremylt   @param ncomp       Number of field components (1 for scalar fields)
384b11c1e72Sjeremylt   @param P1d         Number of nodes in one dimension
385b11c1e72Sjeremylt   @param Q1d         Number of quadrature points in one dimension
38695bb1877Svaleriabarra   @param interp1d    Row-major (Q1d * P1d) matrix expressing the values of nodal
387b11c1e72Sjeremylt                        basis functions at quadrature points
38895bb1877Svaleriabarra   @param grad1d      Row-major (Q1d * P1d) matrix expressing derivatives of nodal
389b11c1e72Sjeremylt                        basis functions at quadrature points
390b11c1e72Sjeremylt   @param qref1d      Array of length Q1d holding the locations of quadrature points
391b11c1e72Sjeremylt                        on the 1D reference element [-1, 1]
392b11c1e72Sjeremylt   @param qweight1d   Array of length Q1d holding the quadrature weights on the
393b11c1e72Sjeremylt                        reference element
394b11c1e72Sjeremylt   @param[out] basis  Address of the variable where the newly created
395b11c1e72Sjeremylt                        CeedBasis will be stored.
396b11c1e72Sjeremylt 
397b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
398dfdf5a53Sjeremylt 
3997a982d89SJeremy L. Thompson   @ref User
400b11c1e72Sjeremylt **/
401d7b241e6Sjeremylt int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt ncomp, CeedInt P1d,
402d7b241e6Sjeremylt                             CeedInt Q1d, const CeedScalar *interp1d,
403d7b241e6Sjeremylt                             const CeedScalar *grad1d, const CeedScalar *qref1d,
404d7b241e6Sjeremylt                             const CeedScalar *qweight1d, CeedBasis *basis) {
405d7b241e6Sjeremylt   int ierr;
406d7b241e6Sjeremylt 
4074d537eeaSYohann   if (dim<1)
408c042f62fSJeremy L Thompson     // LCOV_EXCL_START
4094d537eeaSYohann     return CeedError(ceed, 1, "Basis dimension must be a positive value");
410c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
411d99fa3c5SJeremy L Thompson   CeedElemTopology topo = dim == 1 ? CEED_LINE :
412d99fa3c5SJeremy L Thompson                           dim == 2 ? CEED_QUAD :
413d99fa3c5SJeremy L Thompson                           CEED_HEX;
4144d537eeaSYohann 
4155fe0d4faSjeremylt   if (!ceed->BasisCreateTensorH1) {
4165fe0d4faSjeremylt     Ceed delegate;
417aefd8378Sjeremylt     ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr);
4185fe0d4faSjeremylt 
4195fe0d4faSjeremylt     if (!delegate)
420c042f62fSJeremy L Thompson       // LCOV_EXCL_START
421d7b241e6Sjeremylt       return CeedError(ceed, 1, "Backend does not support BasisCreateTensorH1");
422c042f62fSJeremy L Thompson     // LCOV_EXCL_STOP
4235fe0d4faSjeremylt 
4245fe0d4faSjeremylt     ierr = CeedBasisCreateTensorH1(delegate, dim, ncomp, P1d,
4255fe0d4faSjeremylt                                    Q1d, interp1d, grad1d, qref1d,
4265fe0d4faSjeremylt                                    qweight1d, basis); CeedChk(ierr);
4275fe0d4faSjeremylt     return 0;
4285fe0d4faSjeremylt   }
429d7b241e6Sjeremylt   ierr = CeedCalloc(1,basis); CeedChk(ierr);
430d7b241e6Sjeremylt   (*basis)->ceed = ceed;
431d7b241e6Sjeremylt   ceed->refcount++;
432d7b241e6Sjeremylt   (*basis)->refcount = 1;
433a8de75f0Sjeremylt   (*basis)->tensorbasis = 1;
434d7b241e6Sjeremylt   (*basis)->dim = dim;
435d99fa3c5SJeremy L Thompson   (*basis)->topo = topo;
436d7b241e6Sjeremylt   (*basis)->ncomp = ncomp;
437d7b241e6Sjeremylt   (*basis)->P1d = P1d;
438d7b241e6Sjeremylt   (*basis)->Q1d = Q1d;
439a8de75f0Sjeremylt   (*basis)->P = CeedIntPow(P1d, dim);
440a8de75f0Sjeremylt   (*basis)->Q = CeedIntPow(Q1d, dim);
441d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d,&(*basis)->qref1d); CeedChk(ierr);
442d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d,&(*basis)->qweight1d); CeedChk(ierr);
443d7b241e6Sjeremylt   memcpy((*basis)->qref1d, qref1d, Q1d*sizeof(qref1d[0]));
444d7b241e6Sjeremylt   memcpy((*basis)->qweight1d, qweight1d, Q1d*sizeof(qweight1d[0]));
445d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d*P1d,&(*basis)->interp1d); CeedChk(ierr);
446d7b241e6Sjeremylt   ierr = CeedMalloc(Q1d*P1d,&(*basis)->grad1d); CeedChk(ierr);
447d7b241e6Sjeremylt   memcpy((*basis)->interp1d, interp1d, Q1d*P1d*sizeof(interp1d[0]));
44809486605Sjeremylt   memcpy((*basis)->grad1d, grad1d, Q1d*P1d*sizeof(grad1d[0]));
449667bc5fcSjeremylt   ierr = ceed->BasisCreateTensorH1(dim, P1d, Q1d, interp1d, grad1d, qref1d,
450d7b241e6Sjeremylt                                    qweight1d, *basis); CeedChk(ierr);
451d7b241e6Sjeremylt   return 0;
452d7b241e6Sjeremylt }
453d7b241e6Sjeremylt 
454b11c1e72Sjeremylt /**
45595bb1877Svaleriabarra   @brief Create a tensor-product Lagrange basis
456b11c1e72Sjeremylt 
457b11c1e72Sjeremylt   @param ceed        A Ceed object where the CeedBasis will be created
458b11c1e72Sjeremylt   @param dim         Topological dimension of element
45995bb1877Svaleriabarra   @param ncomp       Number of field components (1 for scalar fields)
460b11c1e72Sjeremylt   @param P           Number of Gauss-Lobatto nodes in one dimension.  The
461b11c1e72Sjeremylt                        polynomial degree of the resulting Q_k element is k=P-1.
462b11c1e72Sjeremylt   @param Q           Number of quadrature points in one dimension.
463b11c1e72Sjeremylt   @param qmode       Distribution of the Q quadrature points (affects order of
464b11c1e72Sjeremylt                        accuracy for the quadrature)
465b11c1e72Sjeremylt   @param[out] basis  Address of the variable where the newly created
466b11c1e72Sjeremylt                        CeedBasis will be stored.
467b11c1e72Sjeremylt 
468b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
469dfdf5a53Sjeremylt 
4707a982d89SJeremy L. Thompson   @ref User
471b11c1e72Sjeremylt **/
472d7b241e6Sjeremylt int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt ncomp,
473692c2638Sjeremylt                                     CeedInt P, CeedInt Q, CeedQuadMode qmode,
474692c2638Sjeremylt                                     CeedBasis *basis) {
475d7b241e6Sjeremylt   // Allocate
476d7b241e6Sjeremylt   int ierr, i, j, k;
477d7b241e6Sjeremylt   CeedScalar c1, c2, c3, c4, dx, *nodes, *interp1d, *grad1d, *qref1d, *qweight1d;
4784d537eeaSYohann 
4794d537eeaSYohann   if (dim<1)
480c042f62fSJeremy L Thompson     // LCOV_EXCL_START
4814d537eeaSYohann     return CeedError(ceed, 1, "Basis dimension must be a positive value");
482c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
4834d537eeaSYohann 
484d7b241e6Sjeremylt   ierr = CeedCalloc(P*Q, &interp1d); CeedChk(ierr);
485d7b241e6Sjeremylt   ierr = CeedCalloc(P*Q, &grad1d); CeedChk(ierr);
486d7b241e6Sjeremylt   ierr = CeedCalloc(P, &nodes); CeedChk(ierr);
487d7b241e6Sjeremylt   ierr = CeedCalloc(Q, &qref1d); CeedChk(ierr);
488d7b241e6Sjeremylt   ierr = CeedCalloc(Q, &qweight1d); CeedChk(ierr);
489d7b241e6Sjeremylt   // Get Nodes and Weights
490d7b241e6Sjeremylt   ierr = CeedLobattoQuadrature(P, nodes, NULL); CeedChk(ierr);
491d7b241e6Sjeremylt   switch (qmode) {
492d7b241e6Sjeremylt   case CEED_GAUSS:
493d7b241e6Sjeremylt     ierr = CeedGaussQuadrature(Q, qref1d, qweight1d); CeedChk(ierr);
494d7b241e6Sjeremylt     break;
495d7b241e6Sjeremylt   case CEED_GAUSS_LOBATTO:
496d7b241e6Sjeremylt     ierr = CeedLobattoQuadrature(Q, qref1d, qweight1d); CeedChk(ierr);
497d7b241e6Sjeremylt     break;
498d7b241e6Sjeremylt   }
499d7b241e6Sjeremylt   // Build B, D matrix
500d7b241e6Sjeremylt   // Fornberg, 1998
501d7b241e6Sjeremylt   for (i = 0; i  < Q; i++) {
502d7b241e6Sjeremylt     c1 = 1.0;
503d7b241e6Sjeremylt     c3 = nodes[0] - qref1d[i];
504d7b241e6Sjeremylt     interp1d[i*P+0] = 1.0;
505d7b241e6Sjeremylt     for (j = 1; j < P; j++) {
506d7b241e6Sjeremylt       c2 = 1.0;
507d7b241e6Sjeremylt       c4 = c3;
508d7b241e6Sjeremylt       c3 = nodes[j] - qref1d[i];
509d7b241e6Sjeremylt       for (k = 0; k < j; k++) {
510d7b241e6Sjeremylt         dx = nodes[j] - nodes[k];
511d7b241e6Sjeremylt         c2 *= dx;
512d7b241e6Sjeremylt         if (k == j - 1) {
513d7b241e6Sjeremylt           grad1d[i*P + j] = c1*(interp1d[i*P + k] - c4*grad1d[i*P + k]) / c2;
514d7b241e6Sjeremylt           interp1d[i*P + j] = - c1*c4*interp1d[i*P + k] / c2;
515d7b241e6Sjeremylt         }
516d7b241e6Sjeremylt         grad1d[i*P + k] = (c3*grad1d[i*P + k] - interp1d[i*P + k]) / dx;
517d7b241e6Sjeremylt         interp1d[i*P + k] = c3*interp1d[i*P + k] / dx;
518d7b241e6Sjeremylt       }
519d7b241e6Sjeremylt       c1 = c2;
520d7b241e6Sjeremylt     }
521d7b241e6Sjeremylt   }
522d7b241e6Sjeremylt   //  // Pass to CeedBasisCreateTensorH1
523d7b241e6Sjeremylt   ierr = CeedBasisCreateTensorH1(ceed, dim, ncomp, P, Q, interp1d, grad1d, qref1d,
524d7b241e6Sjeremylt                                  qweight1d, basis); CeedChk(ierr);
525d7b241e6Sjeremylt   ierr = CeedFree(&interp1d); CeedChk(ierr);
526d7b241e6Sjeremylt   ierr = CeedFree(&grad1d); CeedChk(ierr);
527d7b241e6Sjeremylt   ierr = CeedFree(&nodes); CeedChk(ierr);
528d7b241e6Sjeremylt   ierr = CeedFree(&qref1d); CeedChk(ierr);
529d7b241e6Sjeremylt   ierr = CeedFree(&qweight1d); CeedChk(ierr);
530d7b241e6Sjeremylt   return 0;
531d7b241e6Sjeremylt }
532d7b241e6Sjeremylt 
533b11c1e72Sjeremylt /**
53495bb1877Svaleriabarra   @brief Create a non tensor-product basis for H^1 discretizations
535a8de75f0Sjeremylt 
536a8de75f0Sjeremylt   @param ceed        A Ceed object where the CeedBasis will be created
537a8de75f0Sjeremylt   @param topo        Topology of element, e.g. hypercube, simplex, ect
538a8de75f0Sjeremylt   @param ncomp       Number of field components (1 for scalar fields)
5398795c945Sjeremylt   @param nnodes      Total number of nodes
540a8de75f0Sjeremylt   @param nqpts       Total number of quadrature points
54195bb1877Svaleriabarra   @param interp      Row-major (nqpts * nnodes) matrix expressing the values of
5428795c945Sjeremylt                        nodal basis functions at quadrature points
54395bb1877Svaleriabarra   @param grad        Row-major (nqpts * dim * nnodes) matrix expressing
5448795c945Sjeremylt                        derivatives of nodal basis functions at quadrature points
5458795c945Sjeremylt   @param qref        Array of length nqpts holding the locations of quadrature
5468795c945Sjeremylt                        points on the reference element [-1, 1]
547a8de75f0Sjeremylt   @param qweight     Array of length nqpts holding the quadrature weights on the
548a8de75f0Sjeremylt                        reference element
549a8de75f0Sjeremylt   @param[out] basis  Address of the variable where the newly created
550a8de75f0Sjeremylt                        CeedBasis will be stored.
551a8de75f0Sjeremylt 
552a8de75f0Sjeremylt   @return An error code: 0 - success, otherwise - failure
553a8de75f0Sjeremylt 
5547a982d89SJeremy L. Thompson   @ref User
555a8de75f0Sjeremylt **/
556a8de75f0Sjeremylt int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt ncomp,
557692c2638Sjeremylt                       CeedInt nnodes, CeedInt nqpts, const CeedScalar *interp,
558a8de75f0Sjeremylt                       const CeedScalar *grad, const CeedScalar *qref,
559a8de75f0Sjeremylt                       const CeedScalar *qweight, CeedBasis *basis) {
560a8de75f0Sjeremylt   int ierr;
5618795c945Sjeremylt   CeedInt P = nnodes, Q = nqpts, dim = 0;
562a8de75f0Sjeremylt 
5635fe0d4faSjeremylt   if (!ceed->BasisCreateH1) {
5645fe0d4faSjeremylt     Ceed delegate;
565aefd8378Sjeremylt     ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr);
5665fe0d4faSjeremylt 
5675fe0d4faSjeremylt     if (!delegate)
568c042f62fSJeremy L Thompson       // LCOV_EXCL_START
569a8de75f0Sjeremylt       return CeedError(ceed, 1, "Backend does not support BasisCreateH1");
570c042f62fSJeremy L Thompson     // LCOV_EXCL_STOP
5715fe0d4faSjeremylt 
5728795c945Sjeremylt     ierr = CeedBasisCreateH1(delegate, topo, ncomp, nnodes,
5735fe0d4faSjeremylt                              nqpts, interp, grad, qref,
5745fe0d4faSjeremylt                              qweight, basis); CeedChk(ierr);
5755fe0d4faSjeremylt     return 0;
5765fe0d4faSjeremylt   }
5775fe0d4faSjeremylt 
578a8de75f0Sjeremylt   ierr = CeedCalloc(1,basis); CeedChk(ierr);
579a8de75f0Sjeremylt 
580a8de75f0Sjeremylt   ierr = CeedBasisGetTopologyDimension(topo, &dim); CeedChk(ierr);
581a8de75f0Sjeremylt 
582a8de75f0Sjeremylt   (*basis)->ceed = ceed;
583a8de75f0Sjeremylt   ceed->refcount++;
584a8de75f0Sjeremylt   (*basis)->refcount = 1;
585a8de75f0Sjeremylt   (*basis)->tensorbasis = 0;
586a8de75f0Sjeremylt   (*basis)->dim = dim;
587d99fa3c5SJeremy L Thompson   (*basis)->topo = topo;
588a8de75f0Sjeremylt   (*basis)->ncomp = ncomp;
589a8de75f0Sjeremylt   (*basis)->P = P;
590a8de75f0Sjeremylt   (*basis)->Q = Q;
591a8de75f0Sjeremylt   ierr = CeedMalloc(Q*dim,&(*basis)->qref1d); CeedChk(ierr);
592a8de75f0Sjeremylt   ierr = CeedMalloc(Q,&(*basis)->qweight1d); CeedChk(ierr);
593a8de75f0Sjeremylt   memcpy((*basis)->qref1d, qref, Q*dim*sizeof(qref[0]));
594a8de75f0Sjeremylt   memcpy((*basis)->qweight1d, qweight, Q*sizeof(qweight[0]));
59500f91b2bSjeremylt   ierr = CeedMalloc(Q*P, &(*basis)->interp); CeedChk(ierr);
59600f91b2bSjeremylt   ierr = CeedMalloc(dim*Q*P, &(*basis)->grad); CeedChk(ierr);
59700f91b2bSjeremylt   memcpy((*basis)->interp, interp, Q*P*sizeof(interp[0]));
59800f91b2bSjeremylt   memcpy((*basis)->grad, grad, dim*Q*P*sizeof(grad[0]));
599667bc5fcSjeremylt   ierr = ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, qref,
600a8de75f0Sjeremylt                              qweight, *basis); CeedChk(ierr);
601a8de75f0Sjeremylt   return 0;
602a8de75f0Sjeremylt }
603a8de75f0Sjeremylt 
604a8de75f0Sjeremylt /**
6057a982d89SJeremy L. Thompson   @brief View a CeedBasis
6067a982d89SJeremy L. Thompson 
6077a982d89SJeremy L. Thompson   @param basis   CeedBasis to view
6087a982d89SJeremy L. Thompson   @param stream  Stream to view to, e.g., stdout
6097a982d89SJeremy L. Thompson 
6107a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
6117a982d89SJeremy L. Thompson 
6127a982d89SJeremy L. Thompson   @ref User
6137a982d89SJeremy L. Thompson **/
6147a982d89SJeremy L. Thompson int CeedBasisView(CeedBasis basis, FILE *stream) {
6157a982d89SJeremy L. Thompson   int ierr;
6167a982d89SJeremy L. Thompson 
6177a982d89SJeremy L. Thompson   if (basis->tensorbasis) {
6187a982d89SJeremy L. Thompson     fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P1d,
6197a982d89SJeremy L. Thompson             basis->Q1d);
6207a982d89SJeremy L. Thompson     ierr = CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q1d, basis->qref1d,
6217a982d89SJeremy L. Thompson                           stream); CeedChk(ierr);
6227a982d89SJeremy L. Thompson     ierr = CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q1d,
6237a982d89SJeremy L. Thompson                           basis->qweight1d, stream); CeedChk(ierr);
6247a982d89SJeremy L. Thompson     ierr = CeedScalarView("interp1d", "\t% 12.8f", basis->Q1d, basis->P1d,
6257a982d89SJeremy L. Thompson                           basis->interp1d, stream); CeedChk(ierr);
6267a982d89SJeremy L. Thompson     ierr = CeedScalarView("grad1d", "\t% 12.8f", basis->Q1d, basis->P1d,
6277a982d89SJeremy L. Thompson                           basis->grad1d, stream); CeedChk(ierr);
6287a982d89SJeremy L. Thompson   } else {
6297a982d89SJeremy L. Thompson     fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P,
6307a982d89SJeremy L. Thompson             basis->Q);
6317a982d89SJeremy L. Thompson     ierr = CeedScalarView("qref", "\t% 12.8f", 1, basis->Q*basis->dim,
6327a982d89SJeremy L. Thompson                           basis->qref1d,
6337a982d89SJeremy L. Thompson                           stream); CeedChk(ierr);
6347a982d89SJeremy L. Thompson     ierr = CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->qweight1d,
6357a982d89SJeremy L. Thompson                           stream); CeedChk(ierr);
6367a982d89SJeremy L. Thompson     ierr = CeedScalarView("interp", "\t% 12.8f", basis->Q, basis->P,
6377a982d89SJeremy L. Thompson                           basis->interp, stream); CeedChk(ierr);
6387a982d89SJeremy L. Thompson     ierr = CeedScalarView("grad", "\t% 12.8f", basis->dim*basis->Q, basis->P,
6397a982d89SJeremy L. Thompson                           basis->grad, stream); CeedChk(ierr);
6407a982d89SJeremy L. Thompson   }
6417a982d89SJeremy L. Thompson   return 0;
6427a982d89SJeremy L. Thompson }
6437a982d89SJeremy L. Thompson 
6447a982d89SJeremy L. Thompson /**
6457a982d89SJeremy L. Thompson   @brief Apply basis evaluation from nodes to quadrature points or vice versa
6467a982d89SJeremy L. Thompson 
6477a982d89SJeremy L. Thompson   @param basis   CeedBasis to evaluate
6487a982d89SJeremy L. Thompson   @param nelem   The number of elements to apply the basis evaluation to;
6497a982d89SJeremy L. Thompson                    the backend will specify the ordering in
6504cc79fe7SJed Brown                    CeedElemRestrictionCreateBlocked()
6517a982d89SJeremy L. Thompson   @param tmode   \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature
6527a982d89SJeremy L. Thompson                    points, \ref CEED_TRANSPOSE to apply the transpose, mapping
6537a982d89SJeremy L. Thompson                    from quadrature points to nodes
6547a982d89SJeremy L. Thompson   @param emode   \ref CEED_EVAL_NONE to use values directly,
6557a982d89SJeremy L. Thompson                    \ref CEED_EVAL_INTERP to use interpolated values,
6567a982d89SJeremy L. Thompson                    \ref CEED_EVAL_GRAD to use gradients,
6577a982d89SJeremy L. Thompson                    \ref CEED_EVAL_WEIGHT to use quadrature weights.
6587a982d89SJeremy L. Thompson   @param[in] u   Input CeedVector
6597a982d89SJeremy L. Thompson   @param[out] v  Output CeedVector
6607a982d89SJeremy L. Thompson 
6617a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
6627a982d89SJeremy L. Thompson 
6637a982d89SJeremy L. Thompson   @ref User
6647a982d89SJeremy L. Thompson **/
6657a982d89SJeremy L. Thompson int CeedBasisApply(CeedBasis basis, CeedInt nelem, CeedTransposeMode tmode,
6667a982d89SJeremy L. Thompson                    CeedEvalMode emode, CeedVector u, CeedVector v) {
6677a982d89SJeremy L. Thompson   int ierr;
6687a982d89SJeremy L. Thompson   CeedInt ulength = 0, vlength, nnodes, nqpt;
6697a982d89SJeremy L. Thompson   if (!basis->Apply)
6707a982d89SJeremy L. Thompson     // LCOV_EXCL_START
6717a982d89SJeremy L. Thompson     return CeedError(basis->ceed, 1, "Backend does not support BasisApply");
6727a982d89SJeremy L. Thompson   // LCOV_EXCL_STOP
6737a982d89SJeremy L. Thompson 
6747a982d89SJeremy L. Thompson   // Check compatibility of topological and geometrical dimensions
6757a982d89SJeremy L. Thompson   ierr = CeedBasisGetNumNodes(basis, &nnodes); CeedChk(ierr);
6767a982d89SJeremy L. Thompson   ierr = CeedBasisGetNumQuadraturePoints(basis, &nqpt); CeedChk(ierr);
6777a982d89SJeremy L. Thompson   ierr = CeedVectorGetLength(v, &vlength); CeedChk(ierr);
6787a982d89SJeremy L. Thompson 
6797a982d89SJeremy L. Thompson   if (u) {
6807a982d89SJeremy L. Thompson     ierr = CeedVectorGetLength(u, &ulength); CeedChk(ierr);
6817a982d89SJeremy L. Thompson   }
6827a982d89SJeremy L. Thompson 
6837a982d89SJeremy L. Thompson   if ((tmode == CEED_TRANSPOSE && (vlength%nnodes != 0 || ulength%nqpt != 0)) ||
6847a982d89SJeremy L. Thompson       (tmode == CEED_NOTRANSPOSE && (ulength%nnodes != 0 || vlength%nqpt != 0)))
6857a982d89SJeremy L. Thompson     return CeedError(basis->ceed, 1, "Length of input/output vectors "
6867a982d89SJeremy L. Thompson                      "incompatible with basis dimensions");
6877a982d89SJeremy L. Thompson 
6887a982d89SJeremy L. Thompson   ierr = basis->Apply(basis, nelem, tmode, emode, u, v); CeedChk(ierr);
6897a982d89SJeremy L. Thompson   return 0;
6907a982d89SJeremy L. Thompson }
6917a982d89SJeremy L. Thompson 
6927a982d89SJeremy L. Thompson /**
6939d007619Sjeremylt   @brief Get dimension for given CeedBasis
6949d007619Sjeremylt 
6959d007619Sjeremylt   @param basis     CeedBasis
6969d007619Sjeremylt   @param[out] dim  Variable to store dimension of basis
6979d007619Sjeremylt 
6989d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
6999d007619Sjeremylt 
7009d007619Sjeremylt   @ref Backend
7019d007619Sjeremylt **/
7029d007619Sjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) {
7039d007619Sjeremylt   *dim = basis->dim;
7049d007619Sjeremylt   return 0;
7059d007619Sjeremylt }
7069d007619Sjeremylt 
7079d007619Sjeremylt /**
708d99fa3c5SJeremy L Thompson   @brief Get topology for given CeedBasis
709d99fa3c5SJeremy L Thompson 
710d99fa3c5SJeremy L Thompson   @param basis      CeedBasis
711d99fa3c5SJeremy L Thompson   @param[out] topo  Variable to store topology of basis
712d99fa3c5SJeremy L Thompson 
713d99fa3c5SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
714d99fa3c5SJeremy L Thompson 
715d99fa3c5SJeremy L Thompson   @ref Backend
716d99fa3c5SJeremy L Thompson **/
717d99fa3c5SJeremy L Thompson int CeedBasisGetTopology(CeedBasis basis, CeedElemTopology *topo) {
718d99fa3c5SJeremy L Thompson   *topo = basis->topo;
719d99fa3c5SJeremy L Thompson   return 0;
720d99fa3c5SJeremy L Thompson }
721d99fa3c5SJeremy L Thompson 
722d99fa3c5SJeremy L Thompson /**
7239d007619Sjeremylt   @brief Get number of components for given CeedBasis
7249d007619Sjeremylt 
7259d007619Sjeremylt   @param basis         CeedBasis
7269d007619Sjeremylt   @param[out] numcomp  Variable to store number of components of basis
7279d007619Sjeremylt 
7289d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
7299d007619Sjeremylt 
7309d007619Sjeremylt   @ref Backend
7319d007619Sjeremylt **/
7329d007619Sjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *numcomp) {
7339d007619Sjeremylt   *numcomp = basis->ncomp;
7349d007619Sjeremylt   return 0;
7359d007619Sjeremylt }
7369d007619Sjeremylt 
7379d007619Sjeremylt /**
7389d007619Sjeremylt   @brief Get total number of nodes (in dim dimensions) of a CeedBasis
7399d007619Sjeremylt 
7409d007619Sjeremylt   @param basis   CeedBasis
7419d007619Sjeremylt   @param[out] P  Variable to store number of nodes
7429d007619Sjeremylt 
7439d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
7449d007619Sjeremylt 
7459d007619Sjeremylt   @ref Utility
7469d007619Sjeremylt **/
7479d007619Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) {
7489d007619Sjeremylt   *P = basis->P;
7499d007619Sjeremylt   return 0;
7509d007619Sjeremylt }
7519d007619Sjeremylt 
7529d007619Sjeremylt /**
7539d007619Sjeremylt   @brief Get total number of nodes (in 1 dimension) of a CeedBasis
7549d007619Sjeremylt 
7559d007619Sjeremylt   @param basis     CeedBasis
7569d007619Sjeremylt   @param[out] P1d  Variable to store number of nodes
7579d007619Sjeremylt 
7589d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
7599d007619Sjeremylt 
7609d007619Sjeremylt   @ref Backend
7619d007619Sjeremylt **/
7629d007619Sjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P1d) {
7639d007619Sjeremylt   if (!basis->tensorbasis)
7649d007619Sjeremylt     // LCOV_EXCL_START
7659d007619Sjeremylt     return CeedError(basis->ceed, 1, "Cannot supply P1d for non-tensor basis");
7669d007619Sjeremylt   // LCOV_EXCL_STOP
7679d007619Sjeremylt 
7689d007619Sjeremylt   *P1d = basis->P1d;
7699d007619Sjeremylt   return 0;
7709d007619Sjeremylt }
7719d007619Sjeremylt 
7729d007619Sjeremylt /**
7739d007619Sjeremylt   @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis
7749d007619Sjeremylt 
7759d007619Sjeremylt   @param basis   CeedBasis
7769d007619Sjeremylt   @param[out] Q  Variable to store number of quadrature points
7779d007619Sjeremylt 
7789d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
7799d007619Sjeremylt 
7809d007619Sjeremylt   @ref Utility
7819d007619Sjeremylt **/
7829d007619Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) {
7839d007619Sjeremylt   *Q = basis->Q;
7849d007619Sjeremylt   return 0;
7859d007619Sjeremylt }
7869d007619Sjeremylt 
7879d007619Sjeremylt /**
7889d007619Sjeremylt   @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis
7899d007619Sjeremylt 
7909d007619Sjeremylt   @param basis     CeedBasis
7919d007619Sjeremylt   @param[out] Q1d  Variable to store number of quadrature points
7929d007619Sjeremylt 
7939d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
7949d007619Sjeremylt 
7959d007619Sjeremylt   @ref Backend
7969d007619Sjeremylt **/
7979d007619Sjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q1d) {
7989d007619Sjeremylt   if (!basis->tensorbasis)
7999d007619Sjeremylt     // LCOV_EXCL_START
8009d007619Sjeremylt     return CeedError(basis->ceed, 1, "Cannot supply Q1d for non-tensor basis");
8019d007619Sjeremylt   // LCOV_EXCL_STOP
8029d007619Sjeremylt 
8039d007619Sjeremylt   *Q1d = basis->Q1d;
8049d007619Sjeremylt   return 0;
8059d007619Sjeremylt }
8069d007619Sjeremylt 
8079d007619Sjeremylt /**
8089d007619Sjeremylt   @brief Get reference coordinates of quadrature points (in dim dimensions)
8099d007619Sjeremylt          of a CeedBasis
8109d007619Sjeremylt 
8119d007619Sjeremylt   @param basis      CeedBasis
8129d007619Sjeremylt   @param[out] qref  Variable to store reference coordinates of quadrature points
8139d007619Sjeremylt 
8149d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
8159d007619Sjeremylt 
8169d007619Sjeremylt   @ref Backend
8179d007619Sjeremylt **/
8186c58de82SJeremy L Thompson int CeedBasisGetQRef(CeedBasis basis, const CeedScalar **qref) {
8199d007619Sjeremylt   *qref = basis->qref1d;
8209d007619Sjeremylt   return 0;
8219d007619Sjeremylt }
8229d007619Sjeremylt 
8239d007619Sjeremylt /**
8249d007619Sjeremylt   @brief Get quadrature weights of quadrature points (in dim dimensions)
8259d007619Sjeremylt          of a CeedBasis
8269d007619Sjeremylt 
8279d007619Sjeremylt   @param basis         CeedBasis
8289d007619Sjeremylt   @param[out] qweight  Variable to store quadrature weights
8299d007619Sjeremylt 
8309d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
8319d007619Sjeremylt 
8329d007619Sjeremylt   @ref Backend
8339d007619Sjeremylt **/
8346c58de82SJeremy L Thompson int CeedBasisGetQWeights(CeedBasis basis, const CeedScalar **qweight) {
8359d007619Sjeremylt   *qweight = basis->qweight1d;
8369d007619Sjeremylt   return 0;
8379d007619Sjeremylt }
8389d007619Sjeremylt 
8399d007619Sjeremylt /**
8409d007619Sjeremylt   @brief Get interpolation matrix of a CeedBasis
8419d007619Sjeremylt 
8429d007619Sjeremylt   @param basis        CeedBasis
8439d007619Sjeremylt   @param[out] interp  Variable to store interpolation matrix
8449d007619Sjeremylt 
8459d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
8469d007619Sjeremylt 
8479d007619Sjeremylt   @ref Backend
8489d007619Sjeremylt **/
8496c58de82SJeremy L Thompson int CeedBasisGetInterp(CeedBasis basis, const CeedScalar **interp) {
8509d007619Sjeremylt   if (!basis->interp && basis->tensorbasis) {
8519d007619Sjeremylt     // Allocate
8529d007619Sjeremylt     int ierr;
8539d007619Sjeremylt     ierr = CeedMalloc(basis->Q*basis->P, &basis->interp); CeedChk(ierr);
8549d007619Sjeremylt 
8559d007619Sjeremylt     // Initialize
8569d007619Sjeremylt     for (CeedInt i=0; i<basis->Q*basis->P; i++)
8579d007619Sjeremylt       basis->interp[i] = 1.0;
8589d007619Sjeremylt 
8599d007619Sjeremylt     // Calculate
8609d007619Sjeremylt     for (CeedInt d=0; d<basis->dim; d++)
8619d007619Sjeremylt       for (CeedInt qpt=0; qpt<basis->Q; qpt++)
8629d007619Sjeremylt         for (CeedInt node=0; node<basis->P; node++) {
8639d007619Sjeremylt           CeedInt p = (node / CeedIntPow(basis->P1d, d)) % basis->P1d;
8649d007619Sjeremylt           CeedInt q = (qpt / CeedIntPow(basis->Q1d, d)) % basis->Q1d;
8659d007619Sjeremylt           basis->interp[qpt*(basis->P)+node] *= basis->interp1d[q*basis->P1d+p];
8669d007619Sjeremylt         }
8679d007619Sjeremylt   }
8689d007619Sjeremylt 
8699d007619Sjeremylt   *interp = basis->interp;
8709d007619Sjeremylt 
8719d007619Sjeremylt   return 0;
8729d007619Sjeremylt }
8739d007619Sjeremylt 
8749d007619Sjeremylt /**
8759d007619Sjeremylt   @brief Get 1D interpolation matrix of a tensor product CeedBasis
8769d007619Sjeremylt 
8779d007619Sjeremylt   @param basis          CeedBasis
8789d007619Sjeremylt   @param[out] interp1d  Variable to store interpolation matrix
8799d007619Sjeremylt 
8809d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
8819d007619Sjeremylt 
8829d007619Sjeremylt   @ref Backend
8839d007619Sjeremylt **/
8846c58de82SJeremy L Thompson int CeedBasisGetInterp1D(CeedBasis basis, const CeedScalar **interp1d) {
8859d007619Sjeremylt   if (!basis->tensorbasis)
8869d007619Sjeremylt     // LCOV_EXCL_START
8879d007619Sjeremylt     return CeedError(basis->ceed, 1, "CeedBasis is not a tensor product basis.");
8889d007619Sjeremylt   // LCOV_EXCL_STOP
8899d007619Sjeremylt 
8909d007619Sjeremylt   *interp1d = basis->interp1d;
8919d007619Sjeremylt 
8929d007619Sjeremylt   return 0;
8939d007619Sjeremylt }
8949d007619Sjeremylt 
8959d007619Sjeremylt /**
8969d007619Sjeremylt   @brief Get gradient matrix of a CeedBasis
8979d007619Sjeremylt 
8989d007619Sjeremylt   @param basis      CeedBasis
8999d007619Sjeremylt   @param[out] grad  Variable to store gradient matrix
9009d007619Sjeremylt 
9019d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
9029d007619Sjeremylt 
9039d007619Sjeremylt   @ref Backend
9049d007619Sjeremylt **/
9056c58de82SJeremy L Thompson int CeedBasisGetGrad(CeedBasis basis, const CeedScalar **grad) {
9069d007619Sjeremylt   if (!basis->grad && basis->tensorbasis) {
9079d007619Sjeremylt     // Allocate
9089d007619Sjeremylt     int ierr;
9099d007619Sjeremylt     ierr = CeedMalloc(basis->dim*basis->Q*basis->P, &basis->grad);
9109d007619Sjeremylt     CeedChk(ierr);
9119d007619Sjeremylt 
9129d007619Sjeremylt     // Initialize
9139d007619Sjeremylt     for (CeedInt i=0; i<basis->dim*basis->Q*basis->P; i++)
9149d007619Sjeremylt       basis->grad[i] = 1.0;
9159d007619Sjeremylt 
9169d007619Sjeremylt     // Calculate
9179d007619Sjeremylt     for (CeedInt d=0; d<basis->dim; d++)
9189d007619Sjeremylt       for (CeedInt i=0; i<basis->dim; i++)
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             if (i == d)
9249d007619Sjeremylt               basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *=
9259d007619Sjeremylt                 basis->grad1d[q*basis->P1d+p];
9269d007619Sjeremylt             else
9279d007619Sjeremylt               basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *=
9289d007619Sjeremylt                 basis->interp1d[q*basis->P1d+p];
9299d007619Sjeremylt           }
9309d007619Sjeremylt   }
9319d007619Sjeremylt 
9329d007619Sjeremylt   *grad = basis->grad;
9339d007619Sjeremylt 
9349d007619Sjeremylt   return 0;
9359d007619Sjeremylt }
9369d007619Sjeremylt 
9379d007619Sjeremylt /**
9389d007619Sjeremylt   @brief Get 1D gradient matrix of a tensor product CeedBasis
9399d007619Sjeremylt 
9409d007619Sjeremylt   @param basis        CeedBasis
9419d007619Sjeremylt   @param[out] grad1d  Variable to store gradient matrix
9429d007619Sjeremylt 
9439d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
9449d007619Sjeremylt 
9459d007619Sjeremylt   @ref Backend
9469d007619Sjeremylt **/
9476c58de82SJeremy L Thompson int CeedBasisGetGrad1D(CeedBasis basis, const CeedScalar **grad1d) {
9489d007619Sjeremylt   if (!basis->tensorbasis)
9499d007619Sjeremylt     // LCOV_EXCL_START
9509d007619Sjeremylt     return CeedError(basis->ceed, 1, "CeedBasis is not a tensor product basis.");
9519d007619Sjeremylt   // LCOV_EXCL_STOP
9529d007619Sjeremylt 
9539d007619Sjeremylt   *grad1d = basis->grad1d;
9549d007619Sjeremylt 
9559d007619Sjeremylt   return 0;
9569d007619Sjeremylt }
9579d007619Sjeremylt 
9589d007619Sjeremylt /**
9597a982d89SJeremy L. Thompson   @brief Destroy a CeedBasis
9607a982d89SJeremy L. Thompson 
9617a982d89SJeremy L. Thompson   @param basis CeedBasis to destroy
9627a982d89SJeremy L. Thompson 
9637a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
9647a982d89SJeremy L. Thompson 
9657a982d89SJeremy L. Thompson   @ref User
9667a982d89SJeremy L. Thompson **/
9677a982d89SJeremy L. Thompson int CeedBasisDestroy(CeedBasis *basis) {
9687a982d89SJeremy L. Thompson   int ierr;
9697a982d89SJeremy L. Thompson 
970*752c3701SJeremy L Thompson   if (!*basis || --(*basis)->refcount > 0) return 0;
9717a982d89SJeremy L. Thompson   if ((*basis)->Destroy) {
9727a982d89SJeremy L. Thompson     ierr = (*basis)->Destroy(*basis); CeedChk(ierr);
9737a982d89SJeremy L. Thompson   }
9747a982d89SJeremy L. Thompson   ierr = CeedFree(&(*basis)->interp); CeedChk(ierr);
9757a982d89SJeremy L. Thompson   ierr = CeedFree(&(*basis)->interp1d); CeedChk(ierr);
9767a982d89SJeremy L. Thompson   ierr = CeedFree(&(*basis)->grad); CeedChk(ierr);
9777a982d89SJeremy L. Thompson   ierr = CeedFree(&(*basis)->grad1d); CeedChk(ierr);
9787a982d89SJeremy L. Thompson   ierr = CeedFree(&(*basis)->qref1d); CeedChk(ierr);
9797a982d89SJeremy L. Thompson   ierr = CeedFree(&(*basis)->qweight1d); CeedChk(ierr);
9807a982d89SJeremy L. Thompson   ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr);
9817a982d89SJeremy L. Thompson   ierr = CeedFree(basis); CeedChk(ierr);
9827a982d89SJeremy L. Thompson   return 0;
9837a982d89SJeremy L. Thompson }
9847a982d89SJeremy L. Thompson 
9857a982d89SJeremy L. Thompson /**
986b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre quadrature
987b11c1e72Sjeremylt 
988b11c1e72Sjeremylt   @param Q               Number of quadrature points (integrates polynomials of
989b11c1e72Sjeremylt                            degree 2*Q-1 exactly)
990b11c1e72Sjeremylt   @param[out] qref1d     Array of length Q to hold the abscissa on [-1, 1]
991b11c1e72Sjeremylt   @param[out] qweight1d  Array of length Q to hold the weights
992b11c1e72Sjeremylt 
993b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
994dfdf5a53Sjeremylt 
995dfdf5a53Sjeremylt   @ref Utility
996b11c1e72Sjeremylt **/
997d7b241e6Sjeremylt int CeedGaussQuadrature(CeedInt Q, CeedScalar *qref1d, CeedScalar *qweight1d) {
998d7b241e6Sjeremylt   // Allocate
999d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0*atan(1.0);
1000d7b241e6Sjeremylt   // Build qref1d, qweight1d
1001d7b241e6Sjeremylt   for (int i = 0; i <= Q/2; i++) {
1002d7b241e6Sjeremylt     // Guess
1003d7b241e6Sjeremylt     xi = cos(PI*(CeedScalar)(2*i+1)/((CeedScalar)(2*Q)));
1004d7b241e6Sjeremylt     // Pn(xi)
1005d7b241e6Sjeremylt     P0 = 1.0;
1006d7b241e6Sjeremylt     P1 = xi;
1007d7b241e6Sjeremylt     P2 = 0.0;
1008d7b241e6Sjeremylt     for (int j = 2; j <= Q; j++) {
1009d7b241e6Sjeremylt       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1010d7b241e6Sjeremylt       P0 = P1;
1011d7b241e6Sjeremylt       P1 = P2;
1012d7b241e6Sjeremylt     }
1013d7b241e6Sjeremylt     // First Newton Step
1014d7b241e6Sjeremylt     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1015d7b241e6Sjeremylt     xi = xi-P2/dP2;
1016d7b241e6Sjeremylt     // Newton to convergence
10170e4d4210Sjeremylt     for (int k=0; k<100 && fabs(P2)>10*CEED_EPSILON; k++) {
1018d7b241e6Sjeremylt       P0 = 1.0;
1019d7b241e6Sjeremylt       P1 = xi;
1020d7b241e6Sjeremylt       for (int j = 2; j <= Q; j++) {
1021d7b241e6Sjeremylt         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1022d7b241e6Sjeremylt         P0 = P1;
1023d7b241e6Sjeremylt         P1 = P2;
1024d7b241e6Sjeremylt       }
1025d7b241e6Sjeremylt       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1026d7b241e6Sjeremylt       xi = xi-P2/dP2;
1027d7b241e6Sjeremylt     }
1028d7b241e6Sjeremylt     // Save xi, wi
1029d7b241e6Sjeremylt     wi = 2.0/((1.0-xi*xi)*dP2*dP2);
1030d7b241e6Sjeremylt     qweight1d[i] = wi;
1031d7b241e6Sjeremylt     qweight1d[Q-1-i] = wi;
1032d7b241e6Sjeremylt     qref1d[i] = -xi;
1033d7b241e6Sjeremylt     qref1d[Q-1-i]= xi;
1034d7b241e6Sjeremylt   }
1035d7b241e6Sjeremylt   return 0;
1036d7b241e6Sjeremylt }
1037d7b241e6Sjeremylt 
1038b11c1e72Sjeremylt /**
1039b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre-Lobatto quadrature
1040b11c1e72Sjeremylt 
1041b11c1e72Sjeremylt   @param Q               Number of quadrature points (integrates polynomials of
1042b11c1e72Sjeremylt                            degree 2*Q-3 exactly)
1043b11c1e72Sjeremylt   @param[out] qref1d     Array of length Q to hold the abscissa on [-1, 1]
1044b11c1e72Sjeremylt   @param[out] qweight1d  Array of length Q to hold the weights
1045b11c1e72Sjeremylt 
1046b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1047dfdf5a53Sjeremylt 
1048dfdf5a53Sjeremylt   @ref Utility
1049b11c1e72Sjeremylt **/
1050d7b241e6Sjeremylt int CeedLobattoQuadrature(CeedInt Q, CeedScalar *qref1d,
1051d7b241e6Sjeremylt                           CeedScalar *qweight1d) {
1052d7b241e6Sjeremylt   // Allocate
1053d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0*atan(1.0);
1054d7b241e6Sjeremylt   // Build qref1d, qweight1d
1055d7b241e6Sjeremylt   // Set endpoints
105630a100c3SJed Brown   if (Q < 2)
1057b0d62198Sjeremylt     // LCOV_EXCL_START
10587ed52d01Sjeremylt     return CeedError(NULL, 1,
10597ed52d01Sjeremylt                      "Cannot create Lobatto quadrature with Q=%d < 2 points", Q);
1060b0d62198Sjeremylt   // LCOV_EXCL_STOP
1061d7b241e6Sjeremylt   wi = 2.0/((CeedScalar)(Q*(Q-1)));
1062d7b241e6Sjeremylt   if (qweight1d) {
1063d7b241e6Sjeremylt     qweight1d[0] = wi;
1064d7b241e6Sjeremylt     qweight1d[Q-1] = wi;
1065d7b241e6Sjeremylt   }
1066d7b241e6Sjeremylt   qref1d[0] = -1.0;
1067d7b241e6Sjeremylt   qref1d[Q-1] = 1.0;
1068d7b241e6Sjeremylt   // Interior
1069d7b241e6Sjeremylt   for (int i = 1; i <= (Q-1)/2; i++) {
1070d7b241e6Sjeremylt     // Guess
1071d7b241e6Sjeremylt     xi = cos(PI*(CeedScalar)(i)/(CeedScalar)(Q-1));
1072d7b241e6Sjeremylt     // Pn(xi)
1073d7b241e6Sjeremylt     P0 = 1.0;
1074d7b241e6Sjeremylt     P1 = xi;
1075d7b241e6Sjeremylt     P2 = 0.0;
1076d7b241e6Sjeremylt     for (int j = 2; j < Q; j++) {
1077d7b241e6Sjeremylt       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1078d7b241e6Sjeremylt       P0 = P1;
1079d7b241e6Sjeremylt       P1 = P2;
1080d7b241e6Sjeremylt     }
1081d7b241e6Sjeremylt     // First Newton step
1082d7b241e6Sjeremylt     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1083d7b241e6Sjeremylt     d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
1084d7b241e6Sjeremylt     xi = xi-dP2/d2P2;
1085d7b241e6Sjeremylt     // Newton to convergence
10860e4d4210Sjeremylt     for (int k=0; k<100 && fabs(dP2)>10*CEED_EPSILON; k++) {
1087d7b241e6Sjeremylt       P0 = 1.0;
1088d7b241e6Sjeremylt       P1 = xi;
1089d7b241e6Sjeremylt       for (int j = 2; j < Q; j++) {
1090d7b241e6Sjeremylt         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1091d7b241e6Sjeremylt         P0 = P1;
1092d7b241e6Sjeremylt         P1 = P2;
1093d7b241e6Sjeremylt       }
1094d7b241e6Sjeremylt       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1095d7b241e6Sjeremylt       d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
1096d7b241e6Sjeremylt       xi = xi-dP2/d2P2;
1097d7b241e6Sjeremylt     }
1098d7b241e6Sjeremylt     // Save xi, wi
1099d7b241e6Sjeremylt     wi = 2.0/(((CeedScalar)(Q*(Q-1)))*P2*P2);
1100d7b241e6Sjeremylt     if (qweight1d) {
1101d7b241e6Sjeremylt       qweight1d[i] = wi;
1102d7b241e6Sjeremylt       qweight1d[Q-1-i] = wi;
1103d7b241e6Sjeremylt     }
1104d7b241e6Sjeremylt     qref1d[i] = -xi;
1105d7b241e6Sjeremylt     qref1d[Q-1-i]= xi;
1106d7b241e6Sjeremylt   }
1107d7b241e6Sjeremylt   return 0;
1108d7b241e6Sjeremylt }
1109d7b241e6Sjeremylt 
1110dfdf5a53Sjeremylt /**
111195bb1877Svaleriabarra   @brief Return QR Factorization of a matrix
1112b11c1e72Sjeremylt 
111377645d7bSjeremylt   @param ceed         A Ceed context for error handling
111452bfb9bbSJeremy L Thompson   @param[in,out] mat  Row-major matrix to be factorized in place
111552bfb9bbSJeremy L Thompson   @param[in,out] tau  Vector of length m of scaling factors
1116b11c1e72Sjeremylt   @param m            Number of rows
1117b11c1e72Sjeremylt   @param n            Number of columns
1118b11c1e72Sjeremylt 
1119b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1120dfdf5a53Sjeremylt 
1121dfdf5a53Sjeremylt   @ref Utility
1122b11c1e72Sjeremylt **/
1123a7bd39daSjeremylt int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau,
1124d7b241e6Sjeremylt                         CeedInt m, CeedInt n) {
1125d7b241e6Sjeremylt   CeedScalar v[m];
1126d7b241e6Sjeremylt 
1127a7bd39daSjeremylt   // Check m >= n
1128a7bd39daSjeremylt   if (n > m)
1129c042f62fSJeremy L Thompson     // LCOV_EXCL_START
1130a7bd39daSjeremylt     return CeedError(ceed, 1, "Cannot compute QR factorization with n > m");
1131c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
1132a7bd39daSjeremylt 
113352bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++) {
1134d7b241e6Sjeremylt     // Calculate Householder vector, magnitude
1135d7b241e6Sjeremylt     CeedScalar sigma = 0.0;
1136d7b241e6Sjeremylt     v[i] = mat[i+n*i];
113752bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<m; j++) {
1138d7b241e6Sjeremylt       v[j] = mat[i+n*j];
1139d7b241e6Sjeremylt       sigma += v[j] * v[j];
1140d7b241e6Sjeremylt     }
1141d7b241e6Sjeremylt     CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m]
1142d7b241e6Sjeremylt     CeedScalar Rii = -copysign(norm, v[i]);
1143d7b241e6Sjeremylt     v[i] -= Rii;
1144d7b241e6Sjeremylt     // norm of v[i:m] after modification above and scaling below
1145d7b241e6Sjeremylt     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
1146d7b241e6Sjeremylt     //   tau = 2 / (norm*norm)
1147d7b241e6Sjeremylt     tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma);
1148fb551037Sjeremylt 
11491d102b48SJeremy L Thompson     for (CeedInt j=i+1; j<m; j++)
11501d102b48SJeremy L Thompson       v[j] /= v[i];
1151d7b241e6Sjeremylt 
1152d7b241e6Sjeremylt     // Apply Householder reflector to lower right panel
1153d7b241e6Sjeremylt     CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1);
1154d7b241e6Sjeremylt     // Save v
1155d7b241e6Sjeremylt     mat[i+n*i] = Rii;
11561d102b48SJeremy L Thompson     for (CeedInt j=i+1; j<m; j++)
1157d7b241e6Sjeremylt       mat[i+n*j] = v[j];
1158d7b241e6Sjeremylt   }
1159d7b241e6Sjeremylt 
1160d7b241e6Sjeremylt   return 0;
1161d7b241e6Sjeremylt }
1162d7b241e6Sjeremylt 
1163b11c1e72Sjeremylt /**
116452bfb9bbSJeremy L Thompson   @brief Return symmetric Schur decomposition of the symmetric matrix mat via
116552bfb9bbSJeremy L Thompson            symmetric QR factorization
116652bfb9bbSJeremy L Thompson 
116777645d7bSjeremylt   @param ceed         A Ceed context for error handling
116852bfb9bbSJeremy L Thompson   @param[in,out] mat  Row-major matrix to be factorized in place
1169460bf743SValeria Barra   @param[out] lambda  Vector of length n of eigenvalues
117052bfb9bbSJeremy L Thompson   @param n            Number of rows/columns
117152bfb9bbSJeremy L Thompson 
117252bfb9bbSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
117352bfb9bbSJeremy L Thompson 
117452bfb9bbSJeremy L Thompson   @ref Utility
117552bfb9bbSJeremy L Thompson **/
117652bfb9bbSJeremy L Thompson int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat,
117752bfb9bbSJeremy L Thompson                                     CeedScalar *lambda, CeedInt n) {
117852bfb9bbSJeremy L Thompson   // Check bounds for clang-tidy
117952bfb9bbSJeremy L Thompson   if (n<2)
1180c042f62fSJeremy L Thompson     // LCOV_EXCL_START
1181c042f62fSJeremy L Thompson     return CeedError(ceed, 1,
1182c042f62fSJeremy L Thompson                      "Cannot compute symmetric Schur decomposition of scalars");
1183c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
118452bfb9bbSJeremy L Thompson 
118552bfb9bbSJeremy L Thompson   CeedScalar v[n-1], tau[n-1], matT[n*n];
118652bfb9bbSJeremy L Thompson 
118752bfb9bbSJeremy L Thompson   // Copy mat to matT and set mat to I
118852bfb9bbSJeremy L Thompson   memcpy(matT, mat, n*n*sizeof(mat[0]));
118952bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
119052bfb9bbSJeremy L Thompson     for (CeedInt j=0; j<n; j++)
119152bfb9bbSJeremy L Thompson       mat[j+n*i] = (i==j) ? 1 : 0;
119252bfb9bbSJeremy L Thompson 
119352bfb9bbSJeremy L Thompson   // Reduce to tridiagonal
119452bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n-1; i++) {
119552bfb9bbSJeremy L Thompson     // Calculate Householder vector, magnitude
119652bfb9bbSJeremy L Thompson     CeedScalar sigma = 0.0;
119752bfb9bbSJeremy L Thompson     v[i] = matT[i+n*(i+1)];
119852bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) {
119952bfb9bbSJeremy L Thompson       v[j] = matT[i+n*(j+1)];
120052bfb9bbSJeremy L Thompson       sigma += v[j] * v[j];
120152bfb9bbSJeremy L Thompson     }
120252bfb9bbSJeremy L Thompson     CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:n-1]
120352bfb9bbSJeremy L Thompson     CeedScalar Rii = -copysign(norm, v[i]);
120452bfb9bbSJeremy L Thompson     v[i] -= Rii;
120552bfb9bbSJeremy L Thompson     // norm of v[i:m] after modification above and scaling below
120652bfb9bbSJeremy L Thompson     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
120752bfb9bbSJeremy L Thompson     //   tau = 2 / (norm*norm)
12080e4d4210Sjeremylt     if (sigma > 10*CEED_EPSILON)
120952bfb9bbSJeremy L Thompson       tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma);
1210fb551037Sjeremylt     else
1211fb551037Sjeremylt       tau[i] = 0;
1212fb551037Sjeremylt 
1213fb551037Sjeremylt     for (CeedInt j=i+1; j<n-1; j++)
1214fb551037Sjeremylt       v[j] /= v[i];
121552bfb9bbSJeremy L Thompson 
121652bfb9bbSJeremy L Thompson     // Update sub and super diagonal
121752bfb9bbSJeremy L Thompson     matT[i+n*(i+1)] = Rii;
121852bfb9bbSJeremy L Thompson     matT[(i+1)+n*i] = Rii;
121952bfb9bbSJeremy L Thompson     for (CeedInt j=i+2; j<n; j++) {
122052bfb9bbSJeremy L Thompson       matT[i+n*j] = 0; matT[j+n*i] = 0;
122152bfb9bbSJeremy L Thompson     }
122252bfb9bbSJeremy L Thompson     // Apply symmetric Householder reflector to lower right panel
122352bfb9bbSJeremy L Thompson     CeedHouseholderReflect(&matT[(i+1)+n*(i+1)], &v[i], tau[i],
122452bfb9bbSJeremy L Thompson                            n-(i+1), n-(i+1), n, 1);
122552bfb9bbSJeremy L Thompson     CeedHouseholderReflect(&matT[(i+1)+n*(i+1)], &v[i], tau[i],
122652bfb9bbSJeremy L Thompson                            n-(i+1), n-(i+1), 1, n);
122752bfb9bbSJeremy L Thompson     // Save v
122852bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) {
122952bfb9bbSJeremy L Thompson       matT[i+n*(j+1)] = v[j];
123052bfb9bbSJeremy L Thompson     }
123152bfb9bbSJeremy L Thompson   }
123252bfb9bbSJeremy L Thompson   // Backwards accumulation of Q
123352bfb9bbSJeremy L Thompson   for (CeedInt i=n-2; i>=0; i--) {
123452bfb9bbSJeremy L Thompson     v[i] = 1;
123552bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) {
123652bfb9bbSJeremy L Thompson       v[j] = matT[i+n*(j+1)];
123752bfb9bbSJeremy L Thompson       matT[i+n*(j+1)] = 0;
123852bfb9bbSJeremy L Thompson     }
123952bfb9bbSJeremy L Thompson     CeedHouseholderReflect(&mat[(i+1)+n*(i+1)], &v[i], tau[i],
124052bfb9bbSJeremy L Thompson                            n-(i+1), n-(i+1), n, 1);
124152bfb9bbSJeremy L Thompson   }
124252bfb9bbSJeremy L Thompson 
124352bfb9bbSJeremy L Thompson   // Reduce sub and super diagonal
124452bfb9bbSJeremy L Thompson   CeedInt p = 0, q = 0, itr = 0, maxitr = n*n*n;
12450e4d4210Sjeremylt   CeedScalar tol = 10*CEED_EPSILON;
124652bfb9bbSJeremy L Thompson 
124752bfb9bbSJeremy L Thompson   while (q < n && itr < maxitr) {
124852bfb9bbSJeremy L Thompson     // Update p, q, size of reduced portions of diagonal
124952bfb9bbSJeremy L Thompson     p = 0; q = 0;
125052bfb9bbSJeremy L Thompson     for (CeedInt i=n-2; i>=0; i--) {
125152bfb9bbSJeremy L Thompson       if (fabs(matT[i+n*(i+1)]) < tol)
125252bfb9bbSJeremy L Thompson         q += 1;
125352bfb9bbSJeremy L Thompson       else
125452bfb9bbSJeremy L Thompson         break;
125552bfb9bbSJeremy L Thompson     }
125652bfb9bbSJeremy L Thompson     for (CeedInt i=0; i<n-1-q; i++) {
125752bfb9bbSJeremy L Thompson       if (fabs(matT[i+n*(i+1)]) < tol)
125852bfb9bbSJeremy L Thompson         p += 1;
125952bfb9bbSJeremy L Thompson       else
126052bfb9bbSJeremy L Thompson         break;
126152bfb9bbSJeremy L Thompson     }
126252bfb9bbSJeremy L Thompson     if (q == n-1) break; // Finished reducing
126352bfb9bbSJeremy L Thompson 
126452bfb9bbSJeremy L Thompson     // Reduce tridiagonal portion
126552bfb9bbSJeremy L Thompson     CeedScalar tnn = matT[(n-1-q)+n*(n-1-q)],
126652bfb9bbSJeremy L Thompson                tnnm1 = matT[(n-2-q)+n*(n-1-q)];
126752bfb9bbSJeremy L Thompson     CeedScalar d = (matT[(n-2-q)+n*(n-2-q)] - tnn)/2;
126852bfb9bbSJeremy L Thompson     CeedScalar mu = tnn - tnnm1*tnnm1 /
126952bfb9bbSJeremy L Thompson                     (d + copysign(sqrt(d*d + tnnm1*tnnm1), d));
127052bfb9bbSJeremy L Thompson     CeedScalar x = matT[p+n*p] - mu;
127152bfb9bbSJeremy L Thompson     CeedScalar z = matT[p+n*(p+1)];
127252bfb9bbSJeremy L Thompson     for (CeedInt k=p; k<n-1-q; k++) {
127352bfb9bbSJeremy L Thompson       // Compute Givens rotation
127452bfb9bbSJeremy L Thompson       CeedScalar c = 1, s = 0;
127552bfb9bbSJeremy L Thompson       if (fabs(z) > tol) {
127652bfb9bbSJeremy L Thompson         if (fabs(z) > fabs(x)) {
127752bfb9bbSJeremy L Thompson           CeedScalar tau = -x/z;
127852bfb9bbSJeremy L Thompson           s = 1/sqrt(1+tau*tau), c = s*tau;
127952bfb9bbSJeremy L Thompson         } else {
128052bfb9bbSJeremy L Thompson           CeedScalar tau = -z/x;
128152bfb9bbSJeremy L Thompson           c = 1/sqrt(1+tau*tau), s = c*tau;
128252bfb9bbSJeremy L Thompson         }
128352bfb9bbSJeremy L Thompson       }
128452bfb9bbSJeremy L Thompson 
128552bfb9bbSJeremy L Thompson       // Apply Givens rotation to T
128652bfb9bbSJeremy L Thompson       CeedGivensRotation(matT, c, s, CEED_NOTRANSPOSE, k, k+1, n, n);
128752bfb9bbSJeremy L Thompson       CeedGivensRotation(matT, c, s, CEED_TRANSPOSE, k, k+1, n, n);
128852bfb9bbSJeremy L Thompson 
128952bfb9bbSJeremy L Thompson       // Apply Givens rotation to Q
129052bfb9bbSJeremy L Thompson       CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k+1, n, n);
129152bfb9bbSJeremy L Thompson 
129252bfb9bbSJeremy L Thompson       // Update x, z
129352bfb9bbSJeremy L Thompson       if (k < n-q-2) {
129452bfb9bbSJeremy L Thompson         x = matT[k+n*(k+1)];
129552bfb9bbSJeremy L Thompson         z = matT[k+n*(k+2)];
129652bfb9bbSJeremy L Thompson       }
129752bfb9bbSJeremy L Thompson     }
129852bfb9bbSJeremy L Thompson     itr++;
129952bfb9bbSJeremy L Thompson   }
130052bfb9bbSJeremy L Thompson   // Save eigenvalues
130152bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
130252bfb9bbSJeremy L Thompson     lambda[i] = matT[i+n*i];
130352bfb9bbSJeremy L Thompson 
130452bfb9bbSJeremy L Thompson   // Check convergence
130552bfb9bbSJeremy L Thompson   if (itr == maxitr && q < n-1)
1306c042f62fSJeremy L Thompson     // LCOV_EXCL_START
130752bfb9bbSJeremy L Thompson     return CeedError(ceed, 1, "Symmetric QR failed to converge");
1308c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
130952bfb9bbSJeremy L Thompson 
131052bfb9bbSJeremy L Thompson   return 0;
131152bfb9bbSJeremy L Thompson }
131252bfb9bbSJeremy L Thompson 
131352bfb9bbSJeremy L Thompson /**
131452bfb9bbSJeremy L Thompson   @brief Return Simultaneous Diagonalization of two matrices. This solves the
131552bfb9bbSJeremy L Thompson            generalized eigenvalue problem A x = lambda B x, where A and B
131652bfb9bbSJeremy L Thompson            are symmetric and B is positive definite. We generate the matrix X
131752bfb9bbSJeremy L Thompson            and vector Lambda such that X^T A X = Lambda and X^T B X = I. This
131852bfb9bbSJeremy L Thompson            is equivalent to the LAPACK routine 'sygv' with TYPE = 1.
131952bfb9bbSJeremy L Thompson 
132077645d7bSjeremylt   @param ceed         A Ceed context for error handling
132152bfb9bbSJeremy L Thompson   @param[in] matA     Row-major matrix to be factorized with eigenvalues
132252bfb9bbSJeremy L Thompson   @param[in] matB     Row-major matrix to be factorized to identity
132352bfb9bbSJeremy L Thompson   @param[out] x       Row-major orthogonal matrix
1324460bf743SValeria Barra   @param[out] lambda  Vector of length n of generalized eigenvalues
132552bfb9bbSJeremy L Thompson   @param n            Number of rows/columns
132652bfb9bbSJeremy L Thompson 
132752bfb9bbSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
132852bfb9bbSJeremy L Thompson 
132952bfb9bbSJeremy L Thompson   @ref Utility
133052bfb9bbSJeremy L Thompson **/
133152bfb9bbSJeremy L Thompson int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *matA,
133252bfb9bbSJeremy L Thompson                                     CeedScalar *matB, CeedScalar *x,
133352bfb9bbSJeremy L Thompson                                     CeedScalar *lambda, CeedInt n) {
133452bfb9bbSJeremy L Thompson   int ierr;
133552bfb9bbSJeremy L Thompson   CeedScalar matC[n*n], matG[n*n], vecD[n];
133652bfb9bbSJeremy L Thompson 
133752bfb9bbSJeremy L Thompson   // Compute B = G D G^T
133852bfb9bbSJeremy L Thompson   memcpy(matG, matB, n*n*sizeof(matB[0]));
133952bfb9bbSJeremy L Thompson   ierr = CeedSymmetricSchurDecomposition(ceed, matG, vecD, n); CeedChk(ierr);
1340fb551037Sjeremylt   for (CeedInt i=0; i<n; i++)
1341fb551037Sjeremylt     vecD[i] = sqrt(vecD[i]);
134252bfb9bbSJeremy L Thompson 
1343fb551037Sjeremylt   // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T
1344fb551037Sjeremylt   //           = D^-1/2 G^T A G D^-1/2
134552bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
134652bfb9bbSJeremy L Thompson     for (CeedInt j=0; j<n; j++)
1347fb551037Sjeremylt       matC[j+i*n] = matG[i+j*n] / vecD[i];
13489289e5bfSjeremylt   ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)matC,
13499289e5bfSjeremylt                             (const CeedScalar *)matA, x, n, n, n);
13509289e5bfSjeremylt   CeedChk(ierr);
135152bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
135252bfb9bbSJeremy L Thompson     for (CeedInt j=0; j<n; j++)
1353fb551037Sjeremylt       matG[j+i*n] = matG[j+i*n] / vecD[j];
13549289e5bfSjeremylt   ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)x,
13559289e5bfSjeremylt                             (const CeedScalar *)matG, matC, n, n, n);
13569289e5bfSjeremylt   CeedChk(ierr);
135752bfb9bbSJeremy L Thompson 
135852bfb9bbSJeremy L Thompson   // Compute Q^T C Q = lambda
135952bfb9bbSJeremy L Thompson   ierr = CeedSymmetricSchurDecomposition(ceed, matC, lambda, n); CeedChk(ierr);
136052bfb9bbSJeremy L Thompson 
1361fb551037Sjeremylt   // Set x = (G D^1/2)^-T Q
1362fb551037Sjeremylt   //       = G D^-1/2 Q
13639289e5bfSjeremylt   ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)matG,
13649289e5bfSjeremylt                             (const CeedScalar *)matC, x, n, n, n);
13659289e5bfSjeremylt   CeedChk(ierr);
136652bfb9bbSJeremy L Thompson 
136752bfb9bbSJeremy L Thompson   return 0;
136852bfb9bbSJeremy L Thompson }
136952bfb9bbSJeremy L Thompson 
1370d7b241e6Sjeremylt /// @}
1371