xref: /libCEED/rust/libceed-sys/c-src/interface/ceed-basis.c (revision 227444bf53606b684b95466f0a013ac5dc8b8cba)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3d7b241e6Sjeremylt //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5d7b241e6Sjeremylt //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7d7b241e6Sjeremylt 
8ec3da8bcSJed Brown #include <ceed/ceed.h>
9ec3da8bcSJed Brown #include <ceed/backend.h>
103d576824SJeremy L Thompson #include <ceed-impl.h>
11d7b241e6Sjeremylt #include <math.h>
123d576824SJeremy L Thompson #include <stdbool.h>
13d7b241e6Sjeremylt #include <stdio.h>
14d7b241e6Sjeremylt #include <string.h>
15d7b241e6Sjeremylt 
167a982d89SJeremy L. Thompson /// @file
177a982d89SJeremy L. Thompson /// Implementation of CeedBasis interfaces
187a982d89SJeremy L. Thompson 
19d7b241e6Sjeremylt /// @cond DOXYGEN_SKIP
20783c99b3SValeria Barra static struct CeedBasis_private ceed_basis_collocated;
21d7b241e6Sjeremylt /// @endcond
22d7b241e6Sjeremylt 
237a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser
247a982d89SJeremy L. Thompson /// @{
257a982d89SJeremy L. Thompson 
267a982d89SJeremy L. Thompson /// Indicate that the quadrature points are collocated with the nodes
277a982d89SJeremy L. Thompson const CeedBasis CEED_BASIS_COLLOCATED = &ceed_basis_collocated;
287a982d89SJeremy L. Thompson 
297a982d89SJeremy L. Thompson /// @}
307a982d89SJeremy L. Thompson 
317a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
327a982d89SJeremy L. Thompson /// CeedBasis Library Internal Functions
337a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
347a982d89SJeremy L. Thompson /// @addtogroup CeedBasisDeveloper
357a982d89SJeremy L. Thompson /// @{
367a982d89SJeremy L. Thompson 
377a982d89SJeremy L. Thompson /**
387a982d89SJeremy L. Thompson   @brief Compute Householder reflection
397a982d89SJeremy L. Thompson 
407a982d89SJeremy L. Thompson     Computes A = (I - b v v^T) A
417a982d89SJeremy L. Thompson     where A is an mxn matrix indexed as A[i*row + j*col]
427a982d89SJeremy L. Thompson 
437a982d89SJeremy L. Thompson   @param[in,out] A  Matrix to apply Householder reflection to, in place
447a982d89SJeremy L. Thompson   @param v          Householder vector
457a982d89SJeremy L. Thompson   @param b          Scaling factor
467a982d89SJeremy L. Thompson   @param m          Number of rows in A
477a982d89SJeremy L. Thompson   @param n          Number of columns in A
487a982d89SJeremy L. Thompson   @param row        Row stride
497a982d89SJeremy L. Thompson   @param col        Col stride
507a982d89SJeremy L. Thompson 
517a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
527a982d89SJeremy L. Thompson 
537a982d89SJeremy L. Thompson   @ref Developer
547a982d89SJeremy L. Thompson **/
557a982d89SJeremy L. Thompson static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v,
567a982d89SJeremy L. Thompson                                   CeedScalar b, CeedInt m, CeedInt n,
577a982d89SJeremy L. Thompson                                   CeedInt row, CeedInt col) {
587a982d89SJeremy L. Thompson   for (CeedInt j=0; j<n; j++) {
597a982d89SJeremy L. Thompson     CeedScalar w = A[0*row + j*col];
607a982d89SJeremy L. Thompson     for (CeedInt i=1; i<m; i++)
617a982d89SJeremy L. Thompson       w += v[i] * A[i*row + j*col];
627a982d89SJeremy L. Thompson     A[0*row + j*col] -= b * w;
637a982d89SJeremy L. Thompson     for (CeedInt i=1; i<m; i++)
647a982d89SJeremy L. Thompson       A[i*row + j*col] -= b * w * v[i];
657a982d89SJeremy L. Thompson   }
66e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
677a982d89SJeremy L. Thompson }
687a982d89SJeremy L. Thompson 
697a982d89SJeremy L. Thompson /**
707a982d89SJeremy L. Thompson   @brief Apply Householder Q matrix
717a982d89SJeremy L. Thompson 
727a982d89SJeremy L. Thompson     Compute A = Q A where Q is mxm and A is mxn.
737a982d89SJeremy L. Thompson 
747a982d89SJeremy L. Thompson   @param[in,out] A  Matrix to apply Householder Q to, in place
757a982d89SJeremy L. Thompson   @param Q          Householder Q matrix
767a982d89SJeremy L. Thompson   @param tau        Householder scaling factors
77d1d35e2fSjeremylt   @param t_mode     Transpose mode for application
787a982d89SJeremy L. Thompson   @param m          Number of rows in A
797a982d89SJeremy L. Thompson   @param n          Number of columns in A
807a982d89SJeremy L. Thompson   @param k          Number of elementary reflectors in Q, k<m
817a982d89SJeremy L. Thompson   @param row        Row stride in A
827a982d89SJeremy L. Thompson   @param col        Col stride in A
837a982d89SJeremy L. Thompson 
847a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
857a982d89SJeremy L. Thompson 
867a982d89SJeremy L. Thompson   @ref Developer
877a982d89SJeremy L. Thompson **/
88d99fa3c5SJeremy L Thompson int CeedHouseholderApplyQ(CeedScalar *A, const CeedScalar *Q,
89d1d35e2fSjeremylt                           const CeedScalar *tau, CeedTransposeMode t_mode,
907a982d89SJeremy L. Thompson                           CeedInt m, CeedInt n, CeedInt k,
917a982d89SJeremy L. Thompson                           CeedInt row, CeedInt col) {
92e15f9bd0SJeremy L Thompson   int ierr;
9378464608Sjeremylt   CeedScalar *v;
9478464608Sjeremylt   ierr = CeedMalloc(m, &v); CeedChk(ierr);
957a982d89SJeremy L. Thompson   for (CeedInt ii=0; ii<k; ii++) {
96d1d35e2fSjeremylt     CeedInt i = t_mode == CEED_TRANSPOSE ? ii : k-1-ii;
977a982d89SJeremy L. Thompson     for (CeedInt j=i+1; j<m; j++)
987a982d89SJeremy L. Thompson       v[j] = Q[j*k+i];
99d1d35e2fSjeremylt     // Apply Householder reflector (I - tau v v^T) collo_grad_1d^T
100e15f9bd0SJeremy L Thompson     ierr = CeedHouseholderReflect(&A[i*row], &v[i], tau[i], m-i, n, row, col);
101e15f9bd0SJeremy L Thompson     CeedChk(ierr);
1027a982d89SJeremy L. Thompson   }
10378464608Sjeremylt   ierr = CeedFree(&v); CeedChk(ierr);
104e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1057a982d89SJeremy L. Thompson }
1067a982d89SJeremy L. Thompson 
1077a982d89SJeremy L. Thompson /**
1087a982d89SJeremy L. Thompson   @brief Compute Givens rotation
1097a982d89SJeremy L. Thompson 
1107a982d89SJeremy L. Thompson     Computes A = G A (or G^T A in transpose mode)
1117a982d89SJeremy L. Thompson     where A is an mxn matrix indexed as A[i*n + j*m]
1127a982d89SJeremy L. Thompson 
1137a982d89SJeremy L. Thompson   @param[in,out] A  Row major matrix to apply Givens rotation to, in place
1147a982d89SJeremy L. Thompson   @param c          Cosine factor
1157a982d89SJeremy L. Thompson   @param s          Sine factor
116d1d35e2fSjeremylt   @param t_mode     @ref CEED_NOTRANSPOSE to rotate the basis counter-clockwise,
1174c4400c7SValeria Barra                     which has the effect of rotating columns of A clockwise;
1184cc79fe7SJed Brown                     @ref CEED_TRANSPOSE for the opposite rotation
1197a982d89SJeremy L. Thompson   @param i          First row/column to apply rotation
1207a982d89SJeremy L. Thompson   @param k          Second row/column to apply rotation
1217a982d89SJeremy L. Thompson   @param m          Number of rows in A
1227a982d89SJeremy L. Thompson   @param n          Number of columns in A
1237a982d89SJeremy L. Thompson 
1247a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
1257a982d89SJeremy L. Thompson 
1267a982d89SJeremy L. Thompson   @ref Developer
1277a982d89SJeremy L. Thompson **/
1287a982d89SJeremy L. Thompson static int CeedGivensRotation(CeedScalar *A, CeedScalar c, CeedScalar s,
129d1d35e2fSjeremylt                               CeedTransposeMode t_mode, CeedInt i, CeedInt k,
1307a982d89SJeremy L. Thompson                               CeedInt m, CeedInt n) {
131d1d35e2fSjeremylt   CeedInt stride_j = 1, stride_ik = m, num_its = n;
132d1d35e2fSjeremylt   if (t_mode == CEED_NOTRANSPOSE) {
133d1d35e2fSjeremylt     stride_j = n; stride_ik = 1; num_its = m;
1347a982d89SJeremy L. Thompson   }
1357a982d89SJeremy L. Thompson 
1367a982d89SJeremy L. Thompson   // Apply rotation
137d1d35e2fSjeremylt   for (CeedInt j=0; j<num_its; j++) {
138d1d35e2fSjeremylt     CeedScalar tau1 = A[i*stride_ik+j*stride_j], tau2 = A[k*stride_ik+j*stride_j];
139d1d35e2fSjeremylt     A[i*stride_ik+j*stride_j] = c*tau1 - s*tau2;
140d1d35e2fSjeremylt     A[k*stride_ik+j*stride_j] = s*tau1 + c*tau2;
1417a982d89SJeremy L. Thompson   }
142e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1437a982d89SJeremy L. Thompson }
1447a982d89SJeremy L. Thompson 
1457a982d89SJeremy L. Thompson /**
1467a982d89SJeremy L. Thompson   @brief View an array stored in a CeedBasis
1477a982d89SJeremy L. Thompson 
1480a0da059Sjeremylt   @param[in] name      Name of array
149d1d35e2fSjeremylt   @param[in] fp_fmt    Printing format
1500a0da059Sjeremylt   @param[in] m         Number of rows in array
1510a0da059Sjeremylt   @param[in] n         Number of columns in array
1520a0da059Sjeremylt   @param[in] a         Array to be viewed
1530a0da059Sjeremylt   @param[in] stream    Stream to view to, e.g., stdout
1547a982d89SJeremy L. Thompson 
1557a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
1567a982d89SJeremy L. Thompson 
1577a982d89SJeremy L. Thompson   @ref Developer
1587a982d89SJeremy L. Thompson **/
159d1d35e2fSjeremylt static int CeedScalarView(const char *name, const char *fp_fmt, CeedInt m,
1607a982d89SJeremy L. Thompson                           CeedInt n, const CeedScalar *a, FILE *stream) {
1617a982d89SJeremy L. Thompson   for (int i=0; i<m; i++) {
1627a982d89SJeremy L. Thompson     if (m > 1)
1637a982d89SJeremy L. Thompson       fprintf(stream, "%12s[%d]:", name, i);
1647a982d89SJeremy L. Thompson     else
1657a982d89SJeremy L. Thompson       fprintf(stream, "%12s:", name);
1667a982d89SJeremy L. Thompson     for (int j=0; j<n; j++)
167d1d35e2fSjeremylt       fprintf(stream, fp_fmt, fabs(a[i*n+j]) > 1E-14 ? a[i*n+j] : 0);
1687a982d89SJeremy L. Thompson     fputs("\n", stream);
1697a982d89SJeremy L. Thompson   }
170e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1717a982d89SJeremy L. Thompson }
1727a982d89SJeremy L. Thompson 
1737a982d89SJeremy L. Thompson /// @}
1747a982d89SJeremy L. Thompson 
1757a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
1767a982d89SJeremy L. Thompson /// Ceed Backend API
1777a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
1787a982d89SJeremy L. Thompson /// @addtogroup CeedBasisBackend
1797a982d89SJeremy L. Thompson /// @{
1807a982d89SJeremy L. Thompson 
1817a982d89SJeremy L. Thompson /**
1827a982d89SJeremy L. Thompson   @brief Return collocated grad matrix
1837a982d89SJeremy L. Thompson 
1847a982d89SJeremy L. Thompson   @param basis               CeedBasis
185d1d35e2fSjeremylt   @param[out] collo_grad_1d  Row-major (Q_1d * Q_1d) matrix expressing derivatives of
1867a982d89SJeremy L. Thompson                                basis functions at quadrature points
1877a982d89SJeremy L. Thompson 
1887a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
1897a982d89SJeremy L. Thompson 
1907a982d89SJeremy L. Thompson   @ref Backend
1917a982d89SJeremy L. Thompson **/
192d1d35e2fSjeremylt int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *collo_grad_1d) {
1937a982d89SJeremy L. Thompson   int i, j, k;
1947a982d89SJeremy L. Thompson   Ceed ceed;
195d1d35e2fSjeremylt   CeedInt ierr, P_1d=(basis)->P_1d, Q_1d=(basis)->Q_1d;
19678464608Sjeremylt   CeedScalar *interp_1d, *grad_1d, *tau;
1977a982d89SJeremy L. Thompson 
198d1d35e2fSjeremylt   ierr = CeedMalloc(Q_1d*P_1d, &interp_1d); CeedChk(ierr);
199d1d35e2fSjeremylt   ierr = CeedMalloc(Q_1d*P_1d, &grad_1d); CeedChk(ierr);
20078464608Sjeremylt   ierr = CeedMalloc(Q_1d, &tau); CeedChk(ierr);
201d1d35e2fSjeremylt   memcpy(interp_1d, (basis)->interp_1d, Q_1d*P_1d*sizeof(basis)->interp_1d[0]);
202d1d35e2fSjeremylt   memcpy(grad_1d, (basis)->grad_1d, Q_1d*P_1d*sizeof(basis)->interp_1d[0]);
2037a982d89SJeremy L. Thompson 
204d1d35e2fSjeremylt   // QR Factorization, interp_1d = Q R
2057a982d89SJeremy L. Thompson   ierr = CeedBasisGetCeed(basis, &ceed); CeedChk(ierr);
206d1d35e2fSjeremylt   ierr = CeedQRFactorization(ceed, interp_1d, tau, Q_1d, P_1d); CeedChk(ierr);
207e15f9bd0SJeremy L Thompson   // Note: This function is for backend use, so all errors are terminal
208e15f9bd0SJeremy L Thompson   //   and we do not need to clean up memory on failure.
2097a982d89SJeremy L. Thompson 
210d1d35e2fSjeremylt   // Apply Rinv, collo_grad_1d = grad_1d Rinv
211d1d35e2fSjeremylt   for (i=0; i<Q_1d; i++) { // Row i
212d1d35e2fSjeremylt     collo_grad_1d[Q_1d*i] = grad_1d[P_1d*i]/interp_1d[0];
213d1d35e2fSjeremylt     for (j=1; j<P_1d; j++) { // Column j
214d1d35e2fSjeremylt       collo_grad_1d[j+Q_1d*i] = grad_1d[j+P_1d*i];
2157a982d89SJeremy L. Thompson       for (k=0; k<j; k++)
216d1d35e2fSjeremylt         collo_grad_1d[j+Q_1d*i] -= interp_1d[j+P_1d*k]*collo_grad_1d[k+Q_1d*i];
217d1d35e2fSjeremylt       collo_grad_1d[j+Q_1d*i] /= interp_1d[j+P_1d*j];
2187a982d89SJeremy L. Thompson     }
219d1d35e2fSjeremylt     for (j=P_1d; j<Q_1d; j++)
220d1d35e2fSjeremylt       collo_grad_1d[j+Q_1d*i] = 0;
2217a982d89SJeremy L. Thompson   }
2227a982d89SJeremy L. Thompson 
223673160d7Sjeremylt   // Apply Qtranspose, collo_grad = collo_grad Q_transpose
224d1d35e2fSjeremylt   ierr = CeedHouseholderApplyQ(collo_grad_1d, interp_1d, tau, CEED_NOTRANSPOSE,
225d1d35e2fSjeremylt                                Q_1d, Q_1d, P_1d, 1, Q_1d); CeedChk(ierr);
2267a982d89SJeremy L. Thompson 
227d1d35e2fSjeremylt   ierr = CeedFree(&interp_1d); CeedChk(ierr);
228d1d35e2fSjeremylt   ierr = CeedFree(&grad_1d); CeedChk(ierr);
22978464608Sjeremylt   ierr = CeedFree(&tau); CeedChk(ierr);
230e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2317a982d89SJeremy L. Thompson }
2327a982d89SJeremy L. Thompson 
2337a982d89SJeremy L. Thompson /**
2347a982d89SJeremy L. Thompson   @brief Get tensor status for given CeedBasis
2357a982d89SJeremy L. Thompson 
2367a982d89SJeremy L. Thompson   @param basis           CeedBasis
237d1d35e2fSjeremylt   @param[out] is_tensor  Variable to store tensor status
2387a982d89SJeremy L. Thompson 
2397a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
2407a982d89SJeremy L. Thompson 
2417a982d89SJeremy L. Thompson   @ref Backend
2427a982d89SJeremy L. Thompson **/
243d1d35e2fSjeremylt int CeedBasisIsTensor(CeedBasis basis, bool *is_tensor) {
244d1d35e2fSjeremylt   *is_tensor = basis->tensor_basis;
245e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2467a982d89SJeremy L. Thompson }
2477a982d89SJeremy L. Thompson 
2487a982d89SJeremy L. Thompson /**
2497a982d89SJeremy L. Thompson   @brief Get backend data of a CeedBasis
2507a982d89SJeremy L. Thompson 
2517a982d89SJeremy L. Thompson   @param basis      CeedBasis
2527a982d89SJeremy L. Thompson   @param[out] data  Variable to store data
2537a982d89SJeremy L. Thompson 
2547a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
2557a982d89SJeremy L. Thompson 
2567a982d89SJeremy L. Thompson   @ref Backend
2577a982d89SJeremy L. Thompson **/
258777ff853SJeremy L Thompson int CeedBasisGetData(CeedBasis basis, void *data) {
259777ff853SJeremy L Thompson   *(void **)data = basis->data;
260e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2617a982d89SJeremy L. Thompson }
2627a982d89SJeremy L. Thompson 
2637a982d89SJeremy L. Thompson /**
2647a982d89SJeremy L. Thompson   @brief Set backend data of a CeedBasis
2657a982d89SJeremy L. Thompson 
2667a982d89SJeremy L. Thompson   @param[out] basis  CeedBasis
2677a982d89SJeremy L. Thompson   @param data        Data to set
2687a982d89SJeremy L. Thompson 
2697a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
2707a982d89SJeremy L. Thompson 
2717a982d89SJeremy L. Thompson   @ref Backend
2727a982d89SJeremy L. Thompson **/
273777ff853SJeremy L Thompson int CeedBasisSetData(CeedBasis basis, void *data) {
274777ff853SJeremy L Thompson   basis->data = data;
275e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2767a982d89SJeremy L. Thompson }
2777a982d89SJeremy L. Thompson 
2787a982d89SJeremy L. Thompson /**
27934359f16Sjeremylt   @brief Increment the reference counter for a CeedBasis
28034359f16Sjeremylt 
28134359f16Sjeremylt   @param basis  Basis to increment the reference counter
28234359f16Sjeremylt 
28334359f16Sjeremylt   @return An error code: 0 - success, otherwise - failure
28434359f16Sjeremylt 
28534359f16Sjeremylt   @ref Backend
28634359f16Sjeremylt **/
2879560d06aSjeremylt int CeedBasisReference(CeedBasis basis) {
28834359f16Sjeremylt   basis->ref_count++;
28934359f16Sjeremylt   return CEED_ERROR_SUCCESS;
29034359f16Sjeremylt }
29134359f16Sjeremylt 
29234359f16Sjeremylt /**
2936e15d496SJeremy L Thompson   @brief Estimate number of FLOPs required to apply CeedBasis in t_mode and eval_mode
2946e15d496SJeremy L Thompson 
2956e15d496SJeremy L Thompson   @param basis     Basis to estimate FLOPs for
2966e15d496SJeremy L Thompson   @param t_mode    Apply basis or transpose
2976e15d496SJeremy L Thompson   @param eval_mode Basis evaluation mode
2986e15d496SJeremy L Thompson   @param flops     Address of variable to hold FLOPs estimate
2996e15d496SJeremy L Thompson 
3006e15d496SJeremy L Thompson   @ref Backend
3016e15d496SJeremy L Thompson **/
3026e15d496SJeremy L Thompson int CeedBasisGetFlopsEstimate(CeedBasis basis, CeedTransposeMode t_mode,
3039d36ca50SJeremy L Thompson                               CeedEvalMode eval_mode, CeedSize *flops) {
3046e15d496SJeremy L Thompson   int ierr;
3056e15d496SJeremy L Thompson   bool is_tensor;
3066e15d496SJeremy L Thompson 
3076e15d496SJeremy L Thompson   ierr = CeedBasisIsTensor(basis, &is_tensor); CeedChk(ierr);
3086e15d496SJeremy L Thompson   if (is_tensor) {
3096e15d496SJeremy L Thompson     CeedInt dim, num_comp, P_1d, Q_1d;
3106e15d496SJeremy L Thompson     ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr);
3116e15d496SJeremy L Thompson     ierr = CeedBasisGetNumComponents(basis, &num_comp); CeedChk(ierr);
3126e15d496SJeremy L Thompson     ierr = CeedBasisGetNumNodes1D(basis, &P_1d);  CeedChk(ierr);
3136e15d496SJeremy L Thompson     ierr = CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d);  CeedChk(ierr);
3146e15d496SJeremy L Thompson     if (t_mode == CEED_TRANSPOSE) {
3156e15d496SJeremy L Thompson       P_1d = Q_1d; Q_1d = P_1d;
3166e15d496SJeremy L Thompson     }
3176e15d496SJeremy L Thompson     CeedInt tensor_flops = 0, pre = num_comp * CeedIntPow(P_1d, dim-1), post = 1;
3186e15d496SJeremy L Thompson     for (CeedInt d = 0; d < dim; d++) {
3196e15d496SJeremy L Thompson       tensor_flops += 2 * pre * P_1d * post * Q_1d;
3206e15d496SJeremy L Thompson       pre /= P_1d;
3216e15d496SJeremy L Thompson       post *= Q_1d;
3226e15d496SJeremy L Thompson     }
3236e15d496SJeremy L Thompson     switch (eval_mode) {
3246e15d496SJeremy L Thompson     case CEED_EVAL_NONE:   *flops = 0; break;
3256e15d496SJeremy L Thompson     case CEED_EVAL_INTERP: *flops = tensor_flops; break;
3266e15d496SJeremy L Thompson     case CEED_EVAL_GRAD:   *flops = tensor_flops * 2; break;
3276e15d496SJeremy L Thompson     case CEED_EVAL_DIV:
3286e15d496SJeremy L Thompson       // LCOV_EXCL_START
3296e15d496SJeremy L Thompson       return CeedError(basis->ceed, CEED_ERROR_INCOMPATIBLE,
3306e15d496SJeremy L Thompson                        "Tensor CEED_EVAL_DIV not supported"); break;
3316e15d496SJeremy L Thompson     case CEED_EVAL_CURL:
3326e15d496SJeremy L Thompson       return CeedError(basis->ceed, CEED_ERROR_INCOMPATIBLE,
3336e15d496SJeremy L Thompson                        "Tensor CEED_EVAL_CURL not supported"); break;
3346e15d496SJeremy L Thompson     // LCOV_EXCL_STOP
3356e15d496SJeremy L Thompson     case CEED_EVAL_WEIGHT: *flops = dim * CeedIntPow(Q_1d, dim); break;
3366e15d496SJeremy L Thompson     }
3376e15d496SJeremy L Thompson   } else {
3386e15d496SJeremy L Thompson     CeedInt dim, num_comp, num_nodes, num_qpts, Q_comp;
3396e15d496SJeremy L Thompson     ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr);
3406e15d496SJeremy L Thompson     ierr = CeedBasisGetNumComponents(basis, &num_comp); CeedChk(ierr);
3416e15d496SJeremy L Thompson     ierr = CeedBasisGetNumNodes(basis, &num_nodes); CeedChk(ierr);
3426e15d496SJeremy L Thompson     ierr = CeedBasisGetNumQuadraturePoints(basis, &num_qpts); CeedChk(ierr);
3436e15d496SJeremy L Thompson     ierr = CeedBasisGetNumQuadratureComponents(basis, &Q_comp); CeedChk(ierr);
3446e15d496SJeremy L Thompson     switch (eval_mode) {
3456e15d496SJeremy L Thompson     case CEED_EVAL_NONE:   *flops = 0; break;
3466e15d496SJeremy L Thompson     case CEED_EVAL_INTERP: *flops = num_nodes * num_qpts * num_comp; break;
3476e15d496SJeremy L Thompson     case CEED_EVAL_GRAD:   *flops = num_nodes * num_qpts * num_comp * dim; break;
3486e15d496SJeremy L Thompson     case CEED_EVAL_DIV:    *flops = num_nodes * num_qpts; break;
3496e15d496SJeremy L Thompson     case CEED_EVAL_CURL:   *flops = num_nodes * num_qpts * dim; break;
3506e15d496SJeremy L Thompson     case CEED_EVAL_WEIGHT: *flops = 0; break;
3516e15d496SJeremy L Thompson     }
3526e15d496SJeremy L Thompson   }
3536e15d496SJeremy L Thompson 
3546e15d496SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3556e15d496SJeremy L Thompson }
3566e15d496SJeremy L Thompson 
3576e15d496SJeremy L Thompson /**
3587a982d89SJeremy L. Thompson   @brief Get dimension for given CeedElemTopology
3597a982d89SJeremy L. Thompson 
3607a982d89SJeremy L. Thompson   @param topo      CeedElemTopology
3617a982d89SJeremy L. Thompson   @param[out] dim  Variable to store dimension of topology
3627a982d89SJeremy L. Thompson 
3637a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3647a982d89SJeremy L. Thompson 
3657a982d89SJeremy L. Thompson   @ref Backend
3667a982d89SJeremy L. Thompson **/
3677a982d89SJeremy L. Thompson int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) {
3687a982d89SJeremy L. Thompson   *dim = (CeedInt) topo >> 16;
369e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3707a982d89SJeremy L. Thompson }
3717a982d89SJeremy L. Thompson 
3727a982d89SJeremy L. Thompson /**
3737a982d89SJeremy L. Thompson   @brief Get CeedTensorContract of a CeedBasis
3747a982d89SJeremy L. Thompson 
3757a982d89SJeremy L. Thompson   @param basis          CeedBasis
3767a982d89SJeremy L. Thompson   @param[out] contract  Variable to store CeedTensorContract
3777a982d89SJeremy L. Thompson 
3787a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3797a982d89SJeremy L. Thompson 
3807a982d89SJeremy L. Thompson   @ref Backend
3817a982d89SJeremy L. Thompson **/
3827a982d89SJeremy L. Thompson int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) {
3837a982d89SJeremy L. Thompson   *contract = basis->contract;
384e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3857a982d89SJeremy L. Thompson }
3867a982d89SJeremy L. Thompson 
3877a982d89SJeremy L. Thompson /**
3887a982d89SJeremy L. Thompson   @brief Set CeedTensorContract of a CeedBasis
3897a982d89SJeremy L. Thompson 
3907a982d89SJeremy L. Thompson   @param[out] basis  CeedBasis
3917a982d89SJeremy L. Thompson   @param contract    CeedTensorContract to set
3927a982d89SJeremy L. Thompson 
3937a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3947a982d89SJeremy L. Thompson 
3957a982d89SJeremy L. Thompson   @ref Backend
3967a982d89SJeremy L. Thompson **/
39734359f16Sjeremylt int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract contract) {
39834359f16Sjeremylt   int ierr;
39934359f16Sjeremylt   basis->contract = contract;
4009560d06aSjeremylt   ierr = CeedTensorContractReference(contract); CeedChk(ierr);
401e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4027a982d89SJeremy L. Thompson }
4037a982d89SJeremy L. Thompson 
4047a982d89SJeremy L. Thompson /**
4057a982d89SJeremy L. Thompson   @brief Return a reference implementation of matrix multiplication C = A B.
4067a982d89SJeremy L. Thompson            Note, this is a reference implementation for CPU CeedScalar pointers
4077a982d89SJeremy L. Thompson            that is not intended for high performance.
4087a982d89SJeremy L. Thompson 
4097a982d89SJeremy L. Thompson   @param ceed        A Ceed context for error handling
410d1d35e2fSjeremylt   @param[in] mat_A   Row-major matrix A
411d1d35e2fSjeremylt   @param[in] mat_B   Row-major matrix B
412d1d35e2fSjeremylt   @param[out] mat_C  Row-major output matrix C
4137a982d89SJeremy L. Thompson   @param m           Number of rows of C
4147a982d89SJeremy L. Thompson   @param n           Number of columns of C
4157a982d89SJeremy L. Thompson   @param kk          Number of columns of A/rows of B
4167a982d89SJeremy L. Thompson 
4177a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
4187a982d89SJeremy L. Thompson 
4197a982d89SJeremy L. Thompson   @ref Utility
4207a982d89SJeremy L. Thompson **/
421d1d35e2fSjeremylt int CeedMatrixMultiply(Ceed ceed, const CeedScalar *mat_A,
422d1d35e2fSjeremylt                        const CeedScalar *mat_B, CeedScalar *mat_C, CeedInt m,
4237a982d89SJeremy L. Thompson                        CeedInt n, CeedInt kk) {
4247a982d89SJeremy L. Thompson   for (CeedInt i=0; i<m; i++)
4257a982d89SJeremy L. Thompson     for (CeedInt j=0; j<n; j++) {
4267a982d89SJeremy L. Thompson       CeedScalar sum = 0;
4277a982d89SJeremy L. Thompson       for (CeedInt k=0; k<kk; k++)
428d1d35e2fSjeremylt         sum += mat_A[k+i*kk]*mat_B[j+k*n];
429d1d35e2fSjeremylt       mat_C[j+i*n] = sum;
4307a982d89SJeremy L. Thompson     }
431e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4327a982d89SJeremy L. Thompson }
4337a982d89SJeremy L. Thompson 
4347a982d89SJeremy L. Thompson /// @}
4357a982d89SJeremy L. Thompson 
4367a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
4377a982d89SJeremy L. Thompson /// CeedBasis Public API
4387a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
4397a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser
440d7b241e6Sjeremylt /// @{
441d7b241e6Sjeremylt 
442b11c1e72Sjeremylt /**
44395bb1877Svaleriabarra   @brief Create a tensor-product basis for H^1 discretizations
444b11c1e72Sjeremylt 
445b11c1e72Sjeremylt   @param ceed        A Ceed object where the CeedBasis will be created
446b11c1e72Sjeremylt   @param dim         Topological dimension
447d1d35e2fSjeremylt   @param num_comp    Number of field components (1 for scalar fields)
448d1d35e2fSjeremylt   @param P_1d        Number of nodes in one dimension
449d1d35e2fSjeremylt   @param Q_1d        Number of quadrature points in one dimension
450d1d35e2fSjeremylt   @param interp_1d   Row-major (Q_1d * P_1d) matrix expressing the values of nodal
451b11c1e72Sjeremylt                        basis functions at quadrature points
452d1d35e2fSjeremylt   @param grad_1d     Row-major (Q_1d * P_1d) matrix expressing derivatives of nodal
453b11c1e72Sjeremylt                        basis functions at quadrature points
454d1d35e2fSjeremylt   @param q_ref_1d    Array of length Q_1d holding the locations of quadrature points
455b11c1e72Sjeremylt                        on the 1D reference element [-1, 1]
456d1d35e2fSjeremylt   @param q_weight_1d Array of length Q_1d holding the quadrature weights on the
457b11c1e72Sjeremylt                        reference element
458b11c1e72Sjeremylt   @param[out] basis  Address of the variable where the newly created
459b11c1e72Sjeremylt                        CeedBasis will be stored.
460b11c1e72Sjeremylt 
461b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
462dfdf5a53Sjeremylt 
4637a982d89SJeremy L. Thompson   @ref User
464b11c1e72Sjeremylt **/
465d1d35e2fSjeremylt int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt num_comp,
466d1d35e2fSjeremylt                             CeedInt P_1d, CeedInt Q_1d,
467d1d35e2fSjeremylt                             const CeedScalar *interp_1d,
468d1d35e2fSjeremylt                             const CeedScalar *grad_1d, const CeedScalar *q_ref_1d,
469d1d35e2fSjeremylt                             const CeedScalar *q_weight_1d, CeedBasis *basis) {
470d7b241e6Sjeremylt   int ierr;
471d7b241e6Sjeremylt 
4725fe0d4faSjeremylt   if (!ceed->BasisCreateTensorH1) {
4735fe0d4faSjeremylt     Ceed delegate;
474aefd8378Sjeremylt     ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr);
4755fe0d4faSjeremylt 
4765fe0d4faSjeremylt     if (!delegate)
477c042f62fSJeremy L Thompson       // LCOV_EXCL_START
478e15f9bd0SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_UNSUPPORTED,
479e15f9bd0SJeremy L Thompson                        "Backend does not support BasisCreateTensorH1");
480c042f62fSJeremy L Thompson     // LCOV_EXCL_STOP
4815fe0d4faSjeremylt 
482d1d35e2fSjeremylt     ierr = CeedBasisCreateTensorH1(delegate, dim, num_comp, P_1d,
483d1d35e2fSjeremylt                                    Q_1d, interp_1d, grad_1d, q_ref_1d,
484d1d35e2fSjeremylt                                    q_weight_1d, basis); CeedChk(ierr);
485e15f9bd0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
4865fe0d4faSjeremylt   }
487e15f9bd0SJeremy L Thompson 
488e15f9bd0SJeremy L Thompson   if (dim < 1)
489e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
490e15f9bd0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
491e15f9bd0SJeremy L Thompson                      "Basis dimension must be a positive value");
492e15f9bd0SJeremy L Thompson   // LCOV_EXCL_STOP
493*227444bfSJeremy L Thompson 
494*227444bfSJeremy L Thompson   if (num_comp < 1)
495*227444bfSJeremy L Thompson     // LCOV_EXCL_START
496*227444bfSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
497*227444bfSJeremy L Thompson                      "Basis must have at least 1 component");
498*227444bfSJeremy L Thompson   // LCOV_EXCL_STOP
499*227444bfSJeremy L Thompson 
500*227444bfSJeremy L Thompson   if (P_1d < 1)
501*227444bfSJeremy L Thompson     // LCOV_EXCL_START
502*227444bfSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
503*227444bfSJeremy L Thompson                      "Basis must have at least 1 node");
504*227444bfSJeremy L Thompson   // LCOV_EXCL_STOP
505*227444bfSJeremy L Thompson 
506*227444bfSJeremy L Thompson   if (Q_1d < 1)
507*227444bfSJeremy L Thompson     // LCOV_EXCL_START
508*227444bfSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
509*227444bfSJeremy L Thompson                      "Basis must have at least 1 quadrature point");
510*227444bfSJeremy L Thompson   // LCOV_EXCL_STOP
511*227444bfSJeremy L Thompson 
51250c301a5SRezgar Shakeri   CeedElemTopology topo = dim == 1 ? CEED_TOPOLOGY_LINE
51350c301a5SRezgar Shakeri                           : dim == 2 ? CEED_TOPOLOGY_QUAD
51450c301a5SRezgar Shakeri                           : CEED_TOPOLOGY_HEX;
515e15f9bd0SJeremy L Thompson 
516d7b241e6Sjeremylt   ierr = CeedCalloc(1, basis); CeedChk(ierr);
517d7b241e6Sjeremylt   (*basis)->ceed = ceed;
5189560d06aSjeremylt   ierr = CeedReference(ceed); CeedChk(ierr);
519d1d35e2fSjeremylt   (*basis)->ref_count = 1;
520d1d35e2fSjeremylt   (*basis)->tensor_basis = 1;
521d7b241e6Sjeremylt   (*basis)->dim = dim;
522d99fa3c5SJeremy L Thompson   (*basis)->topo = topo;
523d1d35e2fSjeremylt   (*basis)->num_comp = num_comp;
524d1d35e2fSjeremylt   (*basis)->P_1d = P_1d;
525d1d35e2fSjeremylt   (*basis)->Q_1d = Q_1d;
526d1d35e2fSjeremylt   (*basis)->P = CeedIntPow(P_1d, dim);
527d1d35e2fSjeremylt   (*basis)->Q = CeedIntPow(Q_1d, dim);
52850c301a5SRezgar Shakeri   (*basis)->Q_comp = 1;
52950c301a5SRezgar Shakeri   (*basis)->basis_space = 1; // 1 for H^1 space
530ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(Q_1d, &(*basis)->q_ref_1d); CeedChk(ierr);
531ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(Q_1d, &(*basis)->q_weight_1d); CeedChk(ierr);
532ff3a0f91SJeremy L Thompson   if (q_ref_1d) memcpy((*basis)->q_ref_1d, q_ref_1d, Q_1d*sizeof(q_ref_1d[0]));
533ff3a0f91SJeremy L Thompson   if (q_weight_1d) memcpy((*basis)->q_weight_1d, q_weight_1d,
534ff3a0f91SJeremy L Thompson                             Q_1d*sizeof(q_weight_1d[0]));
535ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(Q_1d*P_1d, &(*basis)->interp_1d); CeedChk(ierr);
536ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(Q_1d*P_1d, &(*basis)->grad_1d); CeedChk(ierr);
537ff3a0f91SJeremy L Thompson   if (interp_1d) memcpy((*basis)->interp_1d, interp_1d,
538ff3a0f91SJeremy L Thompson                           Q_1d*P_1d*sizeof(interp_1d[0]));
539ff3a0f91SJeremy L Thompson   if (grad_1d) memcpy((*basis)->grad_1d, grad_1d, Q_1d*P_1d*sizeof(grad_1d[0]));
540d1d35e2fSjeremylt   ierr = ceed->BasisCreateTensorH1(dim, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d,
541d1d35e2fSjeremylt                                    q_weight_1d, *basis); CeedChk(ierr);
542e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
543d7b241e6Sjeremylt }
544d7b241e6Sjeremylt 
545b11c1e72Sjeremylt /**
54695bb1877Svaleriabarra   @brief Create a tensor-product Lagrange basis
547b11c1e72Sjeremylt 
548b11c1e72Sjeremylt   @param ceed        A Ceed object where the CeedBasis will be created
549b11c1e72Sjeremylt   @param dim         Topological dimension of element
550d1d35e2fSjeremylt   @param num_comp      Number of field components (1 for scalar fields)
551b11c1e72Sjeremylt   @param P           Number of Gauss-Lobatto nodes in one dimension.  The
552b11c1e72Sjeremylt                        polynomial degree of the resulting Q_k element is k=P-1.
553b11c1e72Sjeremylt   @param Q           Number of quadrature points in one dimension.
554d1d35e2fSjeremylt   @param quad_mode   Distribution of the Q quadrature points (affects order of
555b11c1e72Sjeremylt                        accuracy for the quadrature)
556b11c1e72Sjeremylt   @param[out] basis  Address of the variable where the newly created
557b11c1e72Sjeremylt                        CeedBasis will be stored.
558b11c1e72Sjeremylt 
559b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
560dfdf5a53Sjeremylt 
5617a982d89SJeremy L. Thompson   @ref User
562b11c1e72Sjeremylt **/
563d1d35e2fSjeremylt int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt num_comp,
564d1d35e2fSjeremylt                                     CeedInt P, CeedInt Q, CeedQuadMode quad_mode,
565692c2638Sjeremylt                                     CeedBasis *basis) {
566d7b241e6Sjeremylt   // Allocate
567e15f9bd0SJeremy L Thompson   int ierr, ierr2, i, j, k;
568d1d35e2fSjeremylt   CeedScalar c1, c2, c3, c4, dx, *nodes, *interp_1d, *grad_1d, *q_ref_1d,
569d1d35e2fSjeremylt              *q_weight_1d;
5704d537eeaSYohann 
5714d537eeaSYohann   if (dim < 1)
572c042f62fSJeremy L Thompson     // LCOV_EXCL_START
573e15f9bd0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
574e15f9bd0SJeremy L Thompson                      "Basis dimension must be a positive value");
575c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
5764d537eeaSYohann 
577*227444bfSJeremy L Thompson   if (num_comp < 1)
578*227444bfSJeremy L Thompson     // LCOV_EXCL_START
579*227444bfSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
580*227444bfSJeremy L Thompson                      "Basis must have at least 1 component");
581*227444bfSJeremy L Thompson   // LCOV_EXCL_STOP
582*227444bfSJeremy L Thompson 
583*227444bfSJeremy L Thompson   if (P < 1)
584*227444bfSJeremy L Thompson     // LCOV_EXCL_START
585*227444bfSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
586*227444bfSJeremy L Thompson                      "Basis must have at least 1 node");
587*227444bfSJeremy L Thompson   // LCOV_EXCL_STOP
588*227444bfSJeremy L Thompson 
589*227444bfSJeremy L Thompson   if (Q < 1)
590*227444bfSJeremy L Thompson     // LCOV_EXCL_START
591*227444bfSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
592*227444bfSJeremy L Thompson                      "Basis must have at least 1 quadrature point");
593*227444bfSJeremy L Thompson   // LCOV_EXCL_STOP
594*227444bfSJeremy L Thompson 
595e15f9bd0SJeremy L Thompson   // Get Nodes and Weights
596d1d35e2fSjeremylt   ierr = CeedCalloc(P*Q, &interp_1d); CeedChk(ierr);
597d1d35e2fSjeremylt   ierr = CeedCalloc(P*Q, &grad_1d); CeedChk(ierr);
598d7b241e6Sjeremylt   ierr = CeedCalloc(P, &nodes); CeedChk(ierr);
599d1d35e2fSjeremylt   ierr = CeedCalloc(Q, &q_ref_1d); CeedChk(ierr);
600d1d35e2fSjeremylt   ierr = CeedCalloc(Q, &q_weight_1d); CeedChk(ierr);
601e15f9bd0SJeremy L Thompson   ierr = CeedLobattoQuadrature(P, nodes, NULL);
602e15f9bd0SJeremy L Thompson   if (ierr) { goto cleanup; } CeedChk(ierr);
603d1d35e2fSjeremylt   switch (quad_mode) {
604d7b241e6Sjeremylt   case CEED_GAUSS:
605d1d35e2fSjeremylt     ierr = CeedGaussQuadrature(Q, q_ref_1d, q_weight_1d);
606d7b241e6Sjeremylt     break;
607d7b241e6Sjeremylt   case CEED_GAUSS_LOBATTO:
608d1d35e2fSjeremylt     ierr = CeedLobattoQuadrature(Q, q_ref_1d, q_weight_1d);
609d7b241e6Sjeremylt     break;
610d7b241e6Sjeremylt   }
611e15f9bd0SJeremy L Thompson   if (ierr) { goto cleanup; } CeedChk(ierr);
612e15f9bd0SJeremy L Thompson 
613d7b241e6Sjeremylt   // Build B, D matrix
614d7b241e6Sjeremylt   // Fornberg, 1998
615d7b241e6Sjeremylt   for (i = 0; i  < Q; i++) {
616d7b241e6Sjeremylt     c1 = 1.0;
617d1d35e2fSjeremylt     c3 = nodes[0] - q_ref_1d[i];
618d1d35e2fSjeremylt     interp_1d[i*P+0] = 1.0;
619d7b241e6Sjeremylt     for (j = 1; j < P; j++) {
620d7b241e6Sjeremylt       c2 = 1.0;
621d7b241e6Sjeremylt       c4 = c3;
622d1d35e2fSjeremylt       c3 = nodes[j] - q_ref_1d[i];
623d7b241e6Sjeremylt       for (k = 0; k < j; k++) {
624d7b241e6Sjeremylt         dx = nodes[j] - nodes[k];
625d7b241e6Sjeremylt         c2 *= dx;
626d7b241e6Sjeremylt         if (k == j - 1) {
627d1d35e2fSjeremylt           grad_1d[i*P + j] = c1*(interp_1d[i*P + k] - c4*grad_1d[i*P + k]) / c2;
628d1d35e2fSjeremylt           interp_1d[i*P + j] = - c1*c4*interp_1d[i*P + k] / c2;
629d7b241e6Sjeremylt         }
630d1d35e2fSjeremylt         grad_1d[i*P + k] = (c3*grad_1d[i*P + k] - interp_1d[i*P + k]) / dx;
631d1d35e2fSjeremylt         interp_1d[i*P + k] = c3*interp_1d[i*P + k] / dx;
632d7b241e6Sjeremylt       }
633d7b241e6Sjeremylt       c1 = c2;
634d7b241e6Sjeremylt     }
635d7b241e6Sjeremylt   }
6369ac7b42eSJeremy L Thompson   // Pass to CeedBasisCreateTensorH1
637d1d35e2fSjeremylt   ierr = CeedBasisCreateTensorH1(ceed, dim, num_comp, P, Q, interp_1d, grad_1d,
6389ac7b42eSJeremy L Thompson                                  q_ref_1d, q_weight_1d, basis); CeedChk(ierr);
639e15f9bd0SJeremy L Thompson cleanup:
640d1d35e2fSjeremylt   ierr2 = CeedFree(&interp_1d); CeedChk(ierr2);
641d1d35e2fSjeremylt   ierr2 = CeedFree(&grad_1d); CeedChk(ierr2);
642e15f9bd0SJeremy L Thompson   ierr2 = CeedFree(&nodes); CeedChk(ierr2);
643d1d35e2fSjeremylt   ierr2 = CeedFree(&q_ref_1d); CeedChk(ierr2);
644d1d35e2fSjeremylt   ierr2 = CeedFree(&q_weight_1d); CeedChk(ierr2);
645e15f9bd0SJeremy L Thompson   CeedChk(ierr);
646e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
647d7b241e6Sjeremylt }
648d7b241e6Sjeremylt 
649b11c1e72Sjeremylt /**
65095bb1877Svaleriabarra   @brief Create a non tensor-product basis for H^1 discretizations
651a8de75f0Sjeremylt 
652a8de75f0Sjeremylt   @param ceed        A Ceed object where the CeedBasis will be created
653a8de75f0Sjeremylt   @param topo        Topology of element, e.g. hypercube, simplex, ect
654d1d35e2fSjeremylt   @param num_comp    Number of field components (1 for scalar fields)
655d1d35e2fSjeremylt   @param num_nodes   Total number of nodes
656d1d35e2fSjeremylt   @param num_qpts    Total number of quadrature points
657d1d35e2fSjeremylt   @param interp      Row-major (num_qpts * num_nodes) matrix expressing the values of
6588795c945Sjeremylt                        nodal basis functions at quadrature points
659d1d35e2fSjeremylt   @param grad        Row-major (num_qpts * dim * num_nodes) matrix expressing
6608795c945Sjeremylt                        derivatives of nodal basis functions at quadrature points
661d1d35e2fSjeremylt   @param q_ref       Array of length num_qpts holding the locations of quadrature
66250c301a5SRezgar Shakeri                        points on the reference element
663d1d35e2fSjeremylt   @param q_weight    Array of length num_qpts holding the quadrature weights on the
664a8de75f0Sjeremylt                        reference element
665a8de75f0Sjeremylt   @param[out] basis  Address of the variable where the newly created
666a8de75f0Sjeremylt                        CeedBasis will be stored.
667a8de75f0Sjeremylt 
668a8de75f0Sjeremylt   @return An error code: 0 - success, otherwise - failure
669a8de75f0Sjeremylt 
6707a982d89SJeremy L. Thompson   @ref User
671a8de75f0Sjeremylt **/
672d1d35e2fSjeremylt int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt num_comp,
673d1d35e2fSjeremylt                       CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
674d1d35e2fSjeremylt                       const CeedScalar *grad, const CeedScalar *q_ref,
675d1d35e2fSjeremylt                       const CeedScalar *q_weight, CeedBasis *basis) {
676a8de75f0Sjeremylt   int ierr;
677d1d35e2fSjeremylt   CeedInt P = num_nodes, Q = num_qpts, dim = 0;
678a8de75f0Sjeremylt 
6795fe0d4faSjeremylt   if (!ceed->BasisCreateH1) {
6805fe0d4faSjeremylt     Ceed delegate;
681aefd8378Sjeremylt     ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr);
6825fe0d4faSjeremylt 
6835fe0d4faSjeremylt     if (!delegate)
684c042f62fSJeremy L Thompson       // LCOV_EXCL_START
685e15f9bd0SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_UNSUPPORTED,
686e15f9bd0SJeremy L Thompson                        "Backend does not support BasisCreateH1");
687c042f62fSJeremy L Thompson     // LCOV_EXCL_STOP
6885fe0d4faSjeremylt 
689d1d35e2fSjeremylt     ierr = CeedBasisCreateH1(delegate, topo, num_comp, num_nodes,
690d1d35e2fSjeremylt                              num_qpts, interp, grad, q_ref,
691d1d35e2fSjeremylt                              q_weight, basis); CeedChk(ierr);
692e15f9bd0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
6935fe0d4faSjeremylt   }
6945fe0d4faSjeremylt 
695*227444bfSJeremy L Thompson   if (num_comp < 1)
696*227444bfSJeremy L Thompson     // LCOV_EXCL_START
697*227444bfSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
698*227444bfSJeremy L Thompson                      "Basis must have at least 1 component");
699*227444bfSJeremy L Thompson   // LCOV_EXCL_STOP
700*227444bfSJeremy L Thompson 
701*227444bfSJeremy L Thompson   if (num_nodes < 1)
702*227444bfSJeremy L Thompson     // LCOV_EXCL_START
703*227444bfSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
704*227444bfSJeremy L Thompson                      "Basis must have at least 1 node");
705*227444bfSJeremy L Thompson   // LCOV_EXCL_STOP
706*227444bfSJeremy L Thompson 
707*227444bfSJeremy L Thompson   if (num_qpts < 1)
708*227444bfSJeremy L Thompson     // LCOV_EXCL_START
709*227444bfSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
710*227444bfSJeremy L Thompson                      "Basis must have at least 1 quadrature point");
711*227444bfSJeremy L Thompson   // LCOV_EXCL_STOP
712*227444bfSJeremy L Thompson 
713a8de75f0Sjeremylt   ierr = CeedCalloc(1, basis); CeedChk(ierr);
714a8de75f0Sjeremylt 
715a8de75f0Sjeremylt   ierr = CeedBasisGetTopologyDimension(topo, &dim); CeedChk(ierr);
716a8de75f0Sjeremylt 
717a8de75f0Sjeremylt   (*basis)->ceed = ceed;
7189560d06aSjeremylt   ierr = CeedReference(ceed); CeedChk(ierr);
719d1d35e2fSjeremylt   (*basis)->ref_count = 1;
720d1d35e2fSjeremylt   (*basis)->tensor_basis = 0;
721a8de75f0Sjeremylt   (*basis)->dim = dim;
722d99fa3c5SJeremy L Thompson   (*basis)->topo = topo;
723d1d35e2fSjeremylt   (*basis)->num_comp = num_comp;
724a8de75f0Sjeremylt   (*basis)->P = P;
725a8de75f0Sjeremylt   (*basis)->Q = Q;
72650c301a5SRezgar Shakeri   (*basis)->Q_comp = 1;
72750c301a5SRezgar Shakeri   (*basis)->basis_space = 1; // 1 for H^1 space
728ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(Q*dim, &(*basis)->q_ref_1d); CeedChk(ierr);
729ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(Q, &(*basis)->q_weight_1d); CeedChk(ierr);
730ff3a0f91SJeremy L Thompson   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q*dim*sizeof(q_ref[0]));
731ff3a0f91SJeremy L Thompson   if(q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q*sizeof(q_weight[0]));
732ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(Q*P, &(*basis)->interp); CeedChk(ierr);
733ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(dim*Q*P, &(*basis)->grad); CeedChk(ierr);
734ff3a0f91SJeremy L Thompson   if(interp) memcpy((*basis)->interp, interp, Q*P*sizeof(interp[0]));
735ff3a0f91SJeremy L Thompson   if(grad) memcpy((*basis)->grad, grad, dim*Q*P*sizeof(grad[0]));
736d1d35e2fSjeremylt   ierr = ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, q_ref,
737d1d35e2fSjeremylt                              q_weight, *basis); CeedChk(ierr);
738e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
739a8de75f0Sjeremylt }
740a8de75f0Sjeremylt 
741a8de75f0Sjeremylt /**
74250c301a5SRezgar Shakeri   @brief Create a non tensor-product basis for H(div) discretizations
74350c301a5SRezgar Shakeri 
74450c301a5SRezgar Shakeri   @param ceed        A Ceed object where the CeedBasis will be created
74550c301a5SRezgar Shakeri   @param topo        Topology of element (`CEED_TOPOLOGY_QUAD`, `CEED_TOPOLOGY_PRISM`, etc.),
74650c301a5SRezgar Shakeri                      dimension of which is used in some array sizes below
74750c301a5SRezgar Shakeri   @param num_comp    Number of components (usually 1 for vectors in H(div) bases)
74850c301a5SRezgar Shakeri   @param num_nodes   Total number of nodes (dofs per element)
74950c301a5SRezgar Shakeri   @param num_qpts    Total number of quadrature points
75050c301a5SRezgar Shakeri   @param interp      Row-major (dim*num_qpts * num_nodes) matrix expressing the values of
75150c301a5SRezgar Shakeri                        nodal basis functions at quadrature points
75250c301a5SRezgar Shakeri   @param div        Row-major (num_qpts * num_nodes) matrix expressing
75350c301a5SRezgar Shakeri                        divergence of nodal basis functions at quadrature points
75450c301a5SRezgar Shakeri   @param q_ref       Array of length num_qpts holding the locations of quadrature
75550c301a5SRezgar Shakeri                        points on the reference element
75650c301a5SRezgar Shakeri   @param q_weight    Array of length num_qpts holding the quadrature weights on the
75750c301a5SRezgar Shakeri                        reference element
75850c301a5SRezgar Shakeri   @param[out] basis  Address of the variable where the newly created
75950c301a5SRezgar Shakeri                        CeedBasis will be stored.
76050c301a5SRezgar Shakeri 
76150c301a5SRezgar Shakeri   @return An error code: 0 - success, otherwise - failure
76250c301a5SRezgar Shakeri 
76350c301a5SRezgar Shakeri   @ref User
76450c301a5SRezgar Shakeri **/
76550c301a5SRezgar Shakeri int CeedBasisCreateHdiv(Ceed ceed, CeedElemTopology topo, CeedInt num_comp,
76650c301a5SRezgar Shakeri                         CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
76750c301a5SRezgar Shakeri                         const CeedScalar *div, const CeedScalar *q_ref,
76850c301a5SRezgar Shakeri                         const CeedScalar *q_weight, CeedBasis *basis) {
76950c301a5SRezgar Shakeri   int ierr;
77050c301a5SRezgar Shakeri   CeedInt Q = num_qpts, P = num_nodes, dim = 0;
77150c301a5SRezgar Shakeri   ierr = CeedBasisGetTopologyDimension(topo, &dim); CeedChk(ierr);
77250c301a5SRezgar Shakeri   if (!ceed->BasisCreateHdiv) {
77350c301a5SRezgar Shakeri     Ceed delegate;
77450c301a5SRezgar Shakeri     ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr);
77550c301a5SRezgar Shakeri 
77650c301a5SRezgar Shakeri     if (!delegate)
77750c301a5SRezgar Shakeri       // LCOV_EXCL_START
77850c301a5SRezgar Shakeri       return CeedError(ceed, CEED_ERROR_UNSUPPORTED,
77950c301a5SRezgar Shakeri                        "Backend does not implement BasisCreateHdiv");
78050c301a5SRezgar Shakeri     // LCOV_EXCL_STOP
78150c301a5SRezgar Shakeri 
78250c301a5SRezgar Shakeri     ierr = CeedBasisCreateHdiv(delegate, topo, num_comp, num_nodes,
78350c301a5SRezgar Shakeri                                num_qpts, interp, div, q_ref,
78450c301a5SRezgar Shakeri                                q_weight, basis); CeedChk(ierr);
78550c301a5SRezgar Shakeri     return CEED_ERROR_SUCCESS;
78650c301a5SRezgar Shakeri   }
78750c301a5SRezgar Shakeri 
788*227444bfSJeremy L Thompson   if (num_comp < 1)
789*227444bfSJeremy L Thompson     // LCOV_EXCL_START
790*227444bfSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
791*227444bfSJeremy L Thompson                      "Basis must have at least 1 component");
792*227444bfSJeremy L Thompson   // LCOV_EXCL_STOP
793*227444bfSJeremy L Thompson 
794*227444bfSJeremy L Thompson   if (num_nodes < 1)
795*227444bfSJeremy L Thompson     // LCOV_EXCL_START
796*227444bfSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
797*227444bfSJeremy L Thompson                      "Basis must have at least 1 node");
798*227444bfSJeremy L Thompson   // LCOV_EXCL_STOP
799*227444bfSJeremy L Thompson 
800*227444bfSJeremy L Thompson   if (num_qpts < 1)
801*227444bfSJeremy L Thompson     // LCOV_EXCL_START
802*227444bfSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
803*227444bfSJeremy L Thompson                      "Basis must have at least 1 quadrature point");
804*227444bfSJeremy L Thompson   // LCOV_EXCL_STOP
805*227444bfSJeremy L Thompson 
80650c301a5SRezgar Shakeri   ierr = CeedCalloc(1, basis); CeedChk(ierr);
80750c301a5SRezgar Shakeri 
80850c301a5SRezgar Shakeri   (*basis)->ceed = ceed;
80950c301a5SRezgar Shakeri   ierr = CeedReference(ceed); CeedChk(ierr);
81050c301a5SRezgar Shakeri   (*basis)->ref_count = 1;
81150c301a5SRezgar Shakeri   (*basis)->tensor_basis = 0;
81250c301a5SRezgar Shakeri   (*basis)->dim = dim;
81350c301a5SRezgar Shakeri   (*basis)->topo = topo;
81450c301a5SRezgar Shakeri   (*basis)->num_comp = num_comp;
81550c301a5SRezgar Shakeri   (*basis)->P = P;
81650c301a5SRezgar Shakeri   (*basis)->Q = Q;
81750c301a5SRezgar Shakeri   (*basis)->Q_comp = dim;
81850c301a5SRezgar Shakeri   (*basis)->basis_space = 2; // 2 for H(div) space
81950c301a5SRezgar Shakeri   ierr = CeedMalloc(Q*dim, &(*basis)->q_ref_1d); CeedChk(ierr);
82050c301a5SRezgar Shakeri   ierr = CeedMalloc(Q, &(*basis)->q_weight_1d); CeedChk(ierr);
82150c301a5SRezgar Shakeri   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q*dim*sizeof(q_ref[0]));
82250c301a5SRezgar Shakeri   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q*sizeof(q_weight[0]));
82350c301a5SRezgar Shakeri   ierr = CeedMalloc(dim*Q*P, &(*basis)->interp); CeedChk(ierr);
82450c301a5SRezgar Shakeri   ierr = CeedMalloc(Q*P, &(*basis)->div); CeedChk(ierr);
82550c301a5SRezgar Shakeri   if (interp) memcpy((*basis)->interp, interp, dim*Q*P*sizeof(interp[0]));
82650c301a5SRezgar Shakeri   if (div) memcpy((*basis)->div, div, Q*P*sizeof(div[0]));
82750c301a5SRezgar Shakeri   ierr = ceed->BasisCreateHdiv(topo, dim, P, Q, interp, div, q_ref,
82850c301a5SRezgar Shakeri                                q_weight, *basis); CeedChk(ierr);
82950c301a5SRezgar Shakeri   return CEED_ERROR_SUCCESS;
83050c301a5SRezgar Shakeri }
83150c301a5SRezgar Shakeri 
83250c301a5SRezgar Shakeri /**
8339560d06aSjeremylt   @brief Copy the pointer to a CeedBasis. Both pointers should
8349560d06aSjeremylt            be destroyed with `CeedBasisDestroy()`;
8359560d06aSjeremylt            Note: If `*basis_copy` is non-NULL, then it is assumed that
8369560d06aSjeremylt            `*basis_copy` is a pointer to a CeedBasis. This CeedBasis
8379560d06aSjeremylt            will be destroyed if `*basis_copy` is the only
8389560d06aSjeremylt            reference to this CeedBasis.
8399560d06aSjeremylt 
8409560d06aSjeremylt   @param basis            CeedBasis to copy reference to
8419560d06aSjeremylt   @param[out] basis_copy  Variable to store copied reference
8429560d06aSjeremylt 
8439560d06aSjeremylt   @return An error code: 0 - success, otherwise - failure
8449560d06aSjeremylt 
8459560d06aSjeremylt   @ref User
8469560d06aSjeremylt **/
8479560d06aSjeremylt int CeedBasisReferenceCopy(CeedBasis basis, CeedBasis *basis_copy) {
8489560d06aSjeremylt   int ierr;
8499560d06aSjeremylt 
8509560d06aSjeremylt   ierr = CeedBasisReference(basis); CeedChk(ierr);
8519560d06aSjeremylt   ierr = CeedBasisDestroy(basis_copy); CeedChk(ierr);
8529560d06aSjeremylt   *basis_copy = basis;
8539560d06aSjeremylt   return CEED_ERROR_SUCCESS;
8549560d06aSjeremylt }
8559560d06aSjeremylt 
8569560d06aSjeremylt /**
8577a982d89SJeremy L. Thompson   @brief View a CeedBasis
8587a982d89SJeremy L. Thompson 
8597a982d89SJeremy L. Thompson   @param basis   CeedBasis to view
8607a982d89SJeremy L. Thompson   @param stream  Stream to view to, e.g., stdout
8617a982d89SJeremy L. Thompson 
8627a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
8637a982d89SJeremy L. Thompson 
8647a982d89SJeremy L. Thompson   @ref User
8657a982d89SJeremy L. Thompson **/
8667a982d89SJeremy L. Thompson int CeedBasisView(CeedBasis basis, FILE *stream) {
8677a982d89SJeremy L. Thompson   int ierr;
86850c301a5SRezgar Shakeri   CeedFESpace FE_space = basis->basis_space;
86950c301a5SRezgar Shakeri   CeedElemTopology topo = basis->topo;
87050c301a5SRezgar Shakeri   // Print FE space and element topology of the basis
871d1d35e2fSjeremylt   if (basis->tensor_basis) {
87250c301a5SRezgar Shakeri     fprintf(stream, "CeedBasis (%s on a %s element): dim=%d P=%d Q=%d\n",
87350c301a5SRezgar Shakeri             CeedFESpaces[FE_space], CeedElemTopologies[topo],
87450c301a5SRezgar Shakeri             basis->dim, basis->P_1d, basis->Q_1d);
87550c301a5SRezgar Shakeri   } else {
87650c301a5SRezgar Shakeri     fprintf(stream, "CeedBasis (%s on a %s element): dim=%d P=%d Q=%d\n",
87750c301a5SRezgar Shakeri             CeedFESpaces[FE_space], CeedElemTopologies[topo],
87850c301a5SRezgar Shakeri             basis->dim, basis->P, basis->Q);
87950c301a5SRezgar Shakeri   }
88050c301a5SRezgar Shakeri   // Print quadrature data, interpolation/gradient/divergene/curl of the basis
88150c301a5SRezgar Shakeri   if (basis->tensor_basis) { // tensor basis
882d1d35e2fSjeremylt     ierr = CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_ref_1d,
8837a982d89SJeremy L. Thompson                           stream); CeedChk(ierr);
884d1d35e2fSjeremylt     ierr = CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q_1d,
885d1d35e2fSjeremylt                           basis->q_weight_1d, stream); CeedChk(ierr);
886d1d35e2fSjeremylt     ierr = CeedScalarView("interp1d", "\t% 12.8f", basis->Q_1d, basis->P_1d,
887d1d35e2fSjeremylt                           basis->interp_1d, stream); CeedChk(ierr);
888d1d35e2fSjeremylt     ierr = CeedScalarView("grad1d", "\t% 12.8f", basis->Q_1d, basis->P_1d,
889d1d35e2fSjeremylt                           basis->grad_1d, stream); CeedChk(ierr);
89050c301a5SRezgar Shakeri   } else { // non-tensor basis
8917a982d89SJeremy L. Thompson     ierr = CeedScalarView("qref", "\t% 12.8f", 1, basis->Q*basis->dim,
892d1d35e2fSjeremylt                           basis->q_ref_1d,
8937a982d89SJeremy L. Thompson                           stream); CeedChk(ierr);
894d1d35e2fSjeremylt     ierr = CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->q_weight_1d,
8957a982d89SJeremy L. Thompson                           stream); CeedChk(ierr);
89650c301a5SRezgar Shakeri     ierr = CeedScalarView("interp", "\t% 12.8f", basis->Q_comp*basis->Q, basis->P,
8977a982d89SJeremy L. Thompson                           basis->interp, stream); CeedChk(ierr);
89850c301a5SRezgar Shakeri     if (basis->grad) {
8997a982d89SJeremy L. Thompson       ierr = CeedScalarView("grad", "\t% 12.8f", basis->dim*basis->Q, basis->P,
9007a982d89SJeremy L. Thompson                             basis->grad, stream); CeedChk(ierr);
9017a982d89SJeremy L. Thompson     }
90250c301a5SRezgar Shakeri     if (basis->div) {
90350c301a5SRezgar Shakeri       ierr = CeedScalarView("div", "\t% 12.8f", basis->Q, basis->P,
90450c301a5SRezgar Shakeri                             basis->div, stream); CeedChk(ierr);
90550c301a5SRezgar Shakeri     }
90650c301a5SRezgar Shakeri   }
907e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9087a982d89SJeremy L. Thompson }
9097a982d89SJeremy L. Thompson 
9107a982d89SJeremy L. Thompson /**
9117a982d89SJeremy L. Thompson   @brief Apply basis evaluation from nodes to quadrature points or vice versa
9127a982d89SJeremy L. Thompson 
9137a982d89SJeremy L. Thompson   @param basis     CeedBasis to evaluate
914d1d35e2fSjeremylt   @param num_elem  The number of elements to apply the basis evaluation to;
9157a982d89SJeremy L. Thompson                      the backend will specify the ordering in
9164cc79fe7SJed Brown                      CeedElemRestrictionCreateBlocked()
917d1d35e2fSjeremylt   @param t_mode    \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature
9187a982d89SJeremy L. Thompson                      points, \ref CEED_TRANSPOSE to apply the transpose, mapping
9197a982d89SJeremy L. Thompson                      from quadrature points to nodes
920d1d35e2fSjeremylt   @param eval_mode \ref CEED_EVAL_NONE to use values directly,
9217a982d89SJeremy L. Thompson                      \ref CEED_EVAL_INTERP to use interpolated values,
9227a982d89SJeremy L. Thompson                      \ref CEED_EVAL_GRAD to use gradients,
9237a982d89SJeremy L. Thompson                      \ref CEED_EVAL_WEIGHT to use quadrature weights.
9247a982d89SJeremy L. Thompson   @param[in] u     Input CeedVector
9257a982d89SJeremy L. Thompson   @param[out] v    Output CeedVector
9267a982d89SJeremy L. Thompson 
9277a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
9287a982d89SJeremy L. Thompson 
9297a982d89SJeremy L. Thompson   @ref User
9307a982d89SJeremy L. Thompson **/
931d1d35e2fSjeremylt int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode,
932d1d35e2fSjeremylt                    CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
9337a982d89SJeremy L. Thompson   int ierr;
9341f9221feSJeremy L Thompson   CeedSize u_length = 0, v_length;
9351f9221feSJeremy L Thompson   CeedInt dim, num_comp, num_nodes, num_qpts;
936e15f9bd0SJeremy L Thompson   ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr);
937d1d35e2fSjeremylt   ierr = CeedBasisGetNumComponents(basis, &num_comp); CeedChk(ierr);
938d1d35e2fSjeremylt   ierr = CeedBasisGetNumNodes(basis, &num_nodes); CeedChk(ierr);
939d1d35e2fSjeremylt   ierr = CeedBasisGetNumQuadraturePoints(basis, &num_qpts); CeedChk(ierr);
940d1d35e2fSjeremylt   ierr = CeedVectorGetLength(v, &v_length); CeedChk(ierr);
9417a982d89SJeremy L. Thompson   if (u) {
942d1d35e2fSjeremylt     ierr = CeedVectorGetLength(u, &u_length); CeedChk(ierr);
9437a982d89SJeremy L. Thompson   }
9447a982d89SJeremy L. Thompson 
945e15f9bd0SJeremy L Thompson   if (!basis->Apply)
946e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
947e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_UNSUPPORTED,
948e15f9bd0SJeremy L Thompson                      "Backend does not support BasisApply");
949e15f9bd0SJeremy L Thompson   // LCOV_EXCL_STOP
950e15f9bd0SJeremy L Thompson 
951e15f9bd0SJeremy L Thompson   // Check compatibility of topological and geometrical dimensions
952d1d35e2fSjeremylt   if ((t_mode == CEED_TRANSPOSE && (v_length%num_nodes != 0 ||
953d1d35e2fSjeremylt                                     u_length%num_qpts != 0)) ||
954d1d35e2fSjeremylt       (t_mode == CEED_NOTRANSPOSE && (u_length%num_nodes != 0 ||
955d1d35e2fSjeremylt                                       v_length%num_qpts != 0)))
9568229195eSjeremylt     // LCOV_EXCL_START
957e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_DIMENSION,
958e15f9bd0SJeremy L Thompson                      "Length of input/output vectors "
9597a982d89SJeremy L. Thompson                      "incompatible with basis dimensions");
9608229195eSjeremylt   // LCOV_EXCL_STOP
9617a982d89SJeremy L. Thompson 
962e15f9bd0SJeremy L Thompson   // Check vector lengths to prevent out of bounds issues
963d1d35e2fSjeremylt   bool bad_dims = false;
964d1d35e2fSjeremylt   switch (eval_mode) {
965e15f9bd0SJeremy L Thompson   case CEED_EVAL_NONE:
966d1d35e2fSjeremylt   case CEED_EVAL_INTERP: bad_dims =
967d1d35e2fSjeremylt       ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts ||
968d1d35e2fSjeremylt                                      v_length < num_elem*num_comp*num_nodes)) ||
969d1d35e2fSjeremylt        (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp ||
970d1d35e2fSjeremylt                                        u_length < num_elem*num_comp*num_nodes)));
971e15f9bd0SJeremy L Thompson     break;
972d1d35e2fSjeremylt   case CEED_EVAL_GRAD: bad_dims =
973d1d35e2fSjeremylt       ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts*dim ||
974d1d35e2fSjeremylt                                      v_length < num_elem*num_comp*num_nodes)) ||
975d1d35e2fSjeremylt        (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp*dim ||
976d1d35e2fSjeremylt                                        u_length < num_elem*num_comp*num_nodes)));
977e15f9bd0SJeremy L Thompson     break;
978e15f9bd0SJeremy L Thompson   case CEED_EVAL_WEIGHT:
979d1d35e2fSjeremylt     bad_dims = v_length < num_elem*num_qpts;
980e15f9bd0SJeremy L Thompson     break;
981e15f9bd0SJeremy L Thompson   // LCOV_EXCL_START
982d1d35e2fSjeremylt   case CEED_EVAL_DIV: bad_dims =
983d1d35e2fSjeremylt       ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts ||
984d1d35e2fSjeremylt                                      v_length < num_elem*num_comp*num_nodes)) ||
985d1d35e2fSjeremylt        (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp ||
986d1d35e2fSjeremylt                                        u_length < num_elem*num_comp*num_nodes)));
987e15f9bd0SJeremy L Thompson     break;
988d1d35e2fSjeremylt   case CEED_EVAL_CURL: bad_dims =
989d1d35e2fSjeremylt       ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts ||
990d1d35e2fSjeremylt                                      v_length < num_elem*num_comp*num_nodes)) ||
991d1d35e2fSjeremylt        (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp ||
992d1d35e2fSjeremylt                                        u_length < num_elem*num_comp*num_nodes)));
993e15f9bd0SJeremy L Thompson     break;
994e15f9bd0SJeremy L Thompson     // LCOV_EXCL_STOP
995e15f9bd0SJeremy L Thompson   }
996d1d35e2fSjeremylt   if (bad_dims)
997e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
998e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_DIMENSION,
999d1d35e2fSjeremylt                      "Input/output vectors too short for basis and evaluation mode");
1000e15f9bd0SJeremy L Thompson   // LCOV_EXCL_STOP
1001e15f9bd0SJeremy L Thompson 
1002d1d35e2fSjeremylt   ierr = basis->Apply(basis, num_elem, t_mode, eval_mode, u, v); CeedChk(ierr);
1003e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
10047a982d89SJeremy L. Thompson }
10057a982d89SJeremy L. Thompson 
10067a982d89SJeremy L. Thompson /**
1007b7c9bbdaSJeremy L Thompson   @brief Get Ceed associated with a CeedBasis
1008b7c9bbdaSJeremy L Thompson 
1009b7c9bbdaSJeremy L Thompson   @param basis      CeedBasis
1010b7c9bbdaSJeremy L Thompson   @param[out] ceed  Variable to store Ceed
1011b7c9bbdaSJeremy L Thompson 
1012b7c9bbdaSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1013b7c9bbdaSJeremy L Thompson 
1014b7c9bbdaSJeremy L Thompson   @ref Advanced
1015b7c9bbdaSJeremy L Thompson **/
1016b7c9bbdaSJeremy L Thompson int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) {
1017b7c9bbdaSJeremy L Thompson   *ceed = basis->ceed;
1018b7c9bbdaSJeremy L Thompson   return CEED_ERROR_SUCCESS;
1019b7c9bbdaSJeremy L Thompson }
1020b7c9bbdaSJeremy L Thompson 
1021b7c9bbdaSJeremy L Thompson /**
10229d007619Sjeremylt   @brief Get dimension for given CeedBasis
10239d007619Sjeremylt 
10249d007619Sjeremylt   @param basis     CeedBasis
10259d007619Sjeremylt   @param[out] dim  Variable to store dimension of basis
10269d007619Sjeremylt 
10279d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
10289d007619Sjeremylt 
1029b7c9bbdaSJeremy L Thompson   @ref Advanced
10309d007619Sjeremylt **/
10319d007619Sjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) {
10329d007619Sjeremylt   *dim = basis->dim;
1033e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
10349d007619Sjeremylt }
10359d007619Sjeremylt 
10369d007619Sjeremylt /**
1037d99fa3c5SJeremy L Thompson   @brief Get topology for given CeedBasis
1038d99fa3c5SJeremy L Thompson 
1039d99fa3c5SJeremy L Thompson   @param basis      CeedBasis
1040d99fa3c5SJeremy L Thompson   @param[out] topo  Variable to store topology of basis
1041d99fa3c5SJeremy L Thompson 
1042d99fa3c5SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1043d99fa3c5SJeremy L Thompson 
1044b7c9bbdaSJeremy L Thompson   @ref Advanced
1045d99fa3c5SJeremy L Thompson **/
1046d99fa3c5SJeremy L Thompson int CeedBasisGetTopology(CeedBasis basis, CeedElemTopology *topo) {
1047d99fa3c5SJeremy L Thompson   *topo = basis->topo;
1048e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1049d99fa3c5SJeremy L Thompson }
1050d99fa3c5SJeremy L Thompson 
1051d99fa3c5SJeremy L Thompson /**
105250c301a5SRezgar Shakeri   @brief Get number of Q-vector components for given CeedBasis
105350c301a5SRezgar Shakeri 
105450c301a5SRezgar Shakeri   @param basis          CeedBasis
105550c301a5SRezgar Shakeri   @param[out] Q_comp  Variable to store number of Q-vector components of basis
105650c301a5SRezgar Shakeri 
105750c301a5SRezgar Shakeri   @return An error code: 0 - success, otherwise - failure
105850c301a5SRezgar Shakeri 
105950c301a5SRezgar Shakeri   @ref Advanced
106050c301a5SRezgar Shakeri **/
106150c301a5SRezgar Shakeri int CeedBasisGetNumQuadratureComponents(CeedBasis basis, CeedInt *Q_comp) {
106250c301a5SRezgar Shakeri   *Q_comp = basis->Q_comp;
106350c301a5SRezgar Shakeri   return CEED_ERROR_SUCCESS;
106450c301a5SRezgar Shakeri }
106550c301a5SRezgar Shakeri 
106650c301a5SRezgar Shakeri /**
10679d007619Sjeremylt   @brief Get number of components for given CeedBasis
10689d007619Sjeremylt 
10699d007619Sjeremylt   @param basis          CeedBasis
1070d1d35e2fSjeremylt   @param[out] num_comp  Variable to store number of components of basis
10719d007619Sjeremylt 
10729d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
10739d007619Sjeremylt 
1074b7c9bbdaSJeremy L Thompson   @ref Advanced
10759d007619Sjeremylt **/
1076d1d35e2fSjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *num_comp) {
1077d1d35e2fSjeremylt   *num_comp = basis->num_comp;
1078e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
10799d007619Sjeremylt }
10809d007619Sjeremylt 
10819d007619Sjeremylt /**
10829d007619Sjeremylt   @brief Get total number of nodes (in dim dimensions) of a CeedBasis
10839d007619Sjeremylt 
10849d007619Sjeremylt   @param basis   CeedBasis
10859d007619Sjeremylt   @param[out] P  Variable to store number of nodes
10869d007619Sjeremylt 
10879d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
10889d007619Sjeremylt 
10899d007619Sjeremylt   @ref Utility
10909d007619Sjeremylt **/
10919d007619Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) {
10929d007619Sjeremylt   *P = basis->P;
1093e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
10949d007619Sjeremylt }
10959d007619Sjeremylt 
10969d007619Sjeremylt /**
10979d007619Sjeremylt   @brief Get total number of nodes (in 1 dimension) of a CeedBasis
10989d007619Sjeremylt 
10999d007619Sjeremylt   @param basis     CeedBasis
1100d1d35e2fSjeremylt   @param[out] P_1d  Variable to store number of nodes
11019d007619Sjeremylt 
11029d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
11039d007619Sjeremylt 
1104b7c9bbdaSJeremy L Thompson   @ref Advanced
11059d007619Sjeremylt **/
1106d1d35e2fSjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P_1d) {
1107d1d35e2fSjeremylt   if (!basis->tensor_basis)
11089d007619Sjeremylt     // LCOV_EXCL_START
1109e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR,
1110d1d35e2fSjeremylt                      "Cannot supply P_1d for non-tensor basis");
11119d007619Sjeremylt   // LCOV_EXCL_STOP
11129d007619Sjeremylt 
1113d1d35e2fSjeremylt   *P_1d = basis->P_1d;
1114e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
11159d007619Sjeremylt }
11169d007619Sjeremylt 
11179d007619Sjeremylt /**
11189d007619Sjeremylt   @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis
11199d007619Sjeremylt 
11209d007619Sjeremylt   @param basis   CeedBasis
11219d007619Sjeremylt   @param[out] Q  Variable to store number of quadrature points
11229d007619Sjeremylt 
11239d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
11249d007619Sjeremylt 
11259d007619Sjeremylt   @ref Utility
11269d007619Sjeremylt **/
11279d007619Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) {
11289d007619Sjeremylt   *Q = basis->Q;
1129e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
11309d007619Sjeremylt }
11319d007619Sjeremylt 
11329d007619Sjeremylt /**
11339d007619Sjeremylt   @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis
11349d007619Sjeremylt 
11359d007619Sjeremylt   @param basis      CeedBasis
1136d1d35e2fSjeremylt   @param[out] Q_1d  Variable to store number of quadrature points
11379d007619Sjeremylt 
11389d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
11399d007619Sjeremylt 
1140b7c9bbdaSJeremy L Thompson   @ref Advanced
11419d007619Sjeremylt **/
1142d1d35e2fSjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q_1d) {
1143d1d35e2fSjeremylt   if (!basis->tensor_basis)
11449d007619Sjeremylt     // LCOV_EXCL_START
1145e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR,
1146d1d35e2fSjeremylt                      "Cannot supply Q_1d for non-tensor basis");
11479d007619Sjeremylt   // LCOV_EXCL_STOP
11489d007619Sjeremylt 
1149d1d35e2fSjeremylt   *Q_1d = basis->Q_1d;
1150e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
11519d007619Sjeremylt }
11529d007619Sjeremylt 
11539d007619Sjeremylt /**
11549d007619Sjeremylt   @brief Get reference coordinates of quadrature points (in dim dimensions)
11559d007619Sjeremylt          of a CeedBasis
11569d007619Sjeremylt 
11579d007619Sjeremylt   @param basis       CeedBasis
1158d1d35e2fSjeremylt   @param[out] q_ref  Variable to store reference coordinates of quadrature points
11599d007619Sjeremylt 
11609d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
11619d007619Sjeremylt 
1162b7c9bbdaSJeremy L Thompson   @ref Advanced
11639d007619Sjeremylt **/
1164d1d35e2fSjeremylt int CeedBasisGetQRef(CeedBasis basis, const CeedScalar **q_ref) {
1165d1d35e2fSjeremylt   *q_ref = basis->q_ref_1d;
1166e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
11679d007619Sjeremylt }
11689d007619Sjeremylt 
11699d007619Sjeremylt /**
11709d007619Sjeremylt   @brief Get quadrature weights of quadrature points (in dim dimensions)
11719d007619Sjeremylt          of a CeedBasis
11729d007619Sjeremylt 
11739d007619Sjeremylt   @param basis          CeedBasis
1174d1d35e2fSjeremylt   @param[out] q_weight  Variable to store quadrature weights
11759d007619Sjeremylt 
11769d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
11779d007619Sjeremylt 
1178b7c9bbdaSJeremy L Thompson   @ref Advanced
11799d007619Sjeremylt **/
1180d1d35e2fSjeremylt int CeedBasisGetQWeights(CeedBasis basis, const CeedScalar **q_weight) {
1181d1d35e2fSjeremylt   *q_weight = basis->q_weight_1d;
1182e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
11839d007619Sjeremylt }
11849d007619Sjeremylt 
11859d007619Sjeremylt /**
11869d007619Sjeremylt   @brief Get interpolation matrix of a CeedBasis
11879d007619Sjeremylt 
11889d007619Sjeremylt   @param basis        CeedBasis
11899d007619Sjeremylt   @param[out] interp  Variable to store interpolation matrix
11909d007619Sjeremylt 
11919d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
11929d007619Sjeremylt 
1193b7c9bbdaSJeremy L Thompson   @ref Advanced
11949d007619Sjeremylt **/
11956c58de82SJeremy L Thompson int CeedBasisGetInterp(CeedBasis basis, const CeedScalar **interp) {
1196d1d35e2fSjeremylt   if (!basis->interp && basis->tensor_basis) {
11979d007619Sjeremylt     // Allocate
11989d007619Sjeremylt     int ierr;
11999d007619Sjeremylt     ierr = CeedMalloc(basis->Q*basis->P, &basis->interp); CeedChk(ierr);
12009d007619Sjeremylt 
12019d007619Sjeremylt     // Initialize
12029d007619Sjeremylt     for (CeedInt i=0; i<basis->Q*basis->P; i++)
12039d007619Sjeremylt       basis->interp[i] = 1.0;
12049d007619Sjeremylt 
12059d007619Sjeremylt     // Calculate
12069d007619Sjeremylt     for (CeedInt d=0; d<basis->dim; d++)
12079d007619Sjeremylt       for (CeedInt qpt=0; qpt<basis->Q; qpt++)
12089d007619Sjeremylt         for (CeedInt node=0; node<basis->P; node++) {
1209d1d35e2fSjeremylt           CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d;
1210d1d35e2fSjeremylt           CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d;
1211d1d35e2fSjeremylt           basis->interp[qpt*(basis->P)+node] *= basis->interp_1d[q*basis->P_1d+p];
12129d007619Sjeremylt         }
12139d007619Sjeremylt   }
12149d007619Sjeremylt   *interp = basis->interp;
1215e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
12169d007619Sjeremylt }
12179d007619Sjeremylt 
12189d007619Sjeremylt /**
12199d007619Sjeremylt   @brief Get 1D interpolation matrix of a tensor product CeedBasis
12209d007619Sjeremylt 
12219d007619Sjeremylt   @param basis           CeedBasis
1222d1d35e2fSjeremylt   @param[out] interp_1d  Variable to store interpolation matrix
12239d007619Sjeremylt 
12249d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
12259d007619Sjeremylt 
12269d007619Sjeremylt   @ref Backend
12279d007619Sjeremylt **/
1228d1d35e2fSjeremylt int CeedBasisGetInterp1D(CeedBasis basis, const CeedScalar **interp_1d) {
1229d1d35e2fSjeremylt   if (!basis->tensor_basis)
12309d007619Sjeremylt     // LCOV_EXCL_START
1231e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR,
1232e15f9bd0SJeremy L Thompson                      "CeedBasis is not a tensor product basis.");
12339d007619Sjeremylt   // LCOV_EXCL_STOP
12349d007619Sjeremylt 
1235d1d35e2fSjeremylt   *interp_1d = basis->interp_1d;
1236e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
12379d007619Sjeremylt }
12389d007619Sjeremylt 
12399d007619Sjeremylt /**
12409d007619Sjeremylt   @brief Get gradient matrix of a CeedBasis
12419d007619Sjeremylt 
12429d007619Sjeremylt   @param basis      CeedBasis
12439d007619Sjeremylt   @param[out] grad  Variable to store gradient matrix
12449d007619Sjeremylt 
12459d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
12469d007619Sjeremylt 
1247b7c9bbdaSJeremy L Thompson   @ref Advanced
12489d007619Sjeremylt **/
12496c58de82SJeremy L Thompson int CeedBasisGetGrad(CeedBasis basis, const CeedScalar **grad) {
1250d1d35e2fSjeremylt   if (!basis->grad && basis->tensor_basis) {
12519d007619Sjeremylt     // Allocate
12529d007619Sjeremylt     int ierr;
12539d007619Sjeremylt     ierr = CeedMalloc(basis->dim*basis->Q*basis->P, &basis->grad);
12549d007619Sjeremylt     CeedChk(ierr);
12559d007619Sjeremylt 
12569d007619Sjeremylt     // Initialize
12579d007619Sjeremylt     for (CeedInt i=0; i<basis->dim*basis->Q*basis->P; i++)
12589d007619Sjeremylt       basis->grad[i] = 1.0;
12599d007619Sjeremylt 
12609d007619Sjeremylt     // Calculate
12619d007619Sjeremylt     for (CeedInt d=0; d<basis->dim; d++)
12629d007619Sjeremylt       for (CeedInt i=0; i<basis->dim; i++)
12639d007619Sjeremylt         for (CeedInt qpt=0; qpt<basis->Q; qpt++)
12649d007619Sjeremylt           for (CeedInt node=0; node<basis->P; node++) {
1265d1d35e2fSjeremylt             CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d;
1266d1d35e2fSjeremylt             CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d;
12679d007619Sjeremylt             if (i == d)
12689d007619Sjeremylt               basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *=
1269d1d35e2fSjeremylt                 basis->grad_1d[q*basis->P_1d+p];
12709d007619Sjeremylt             else
12719d007619Sjeremylt               basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *=
1272d1d35e2fSjeremylt                 basis->interp_1d[q*basis->P_1d+p];
12739d007619Sjeremylt           }
12749d007619Sjeremylt   }
12759d007619Sjeremylt   *grad = basis->grad;
1276e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
12779d007619Sjeremylt }
12789d007619Sjeremylt 
12799d007619Sjeremylt /**
12809d007619Sjeremylt   @brief Get 1D gradient matrix of a tensor product CeedBasis
12819d007619Sjeremylt 
12829d007619Sjeremylt   @param basis         CeedBasis
1283d1d35e2fSjeremylt   @param[out] grad_1d  Variable to store gradient matrix
12849d007619Sjeremylt 
12859d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
12869d007619Sjeremylt 
1287b7c9bbdaSJeremy L Thompson   @ref Advanced
12889d007619Sjeremylt **/
1289d1d35e2fSjeremylt int CeedBasisGetGrad1D(CeedBasis basis, const CeedScalar **grad_1d) {
1290d1d35e2fSjeremylt   if (!basis->tensor_basis)
12919d007619Sjeremylt     // LCOV_EXCL_START
1292e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR,
1293e15f9bd0SJeremy L Thompson                      "CeedBasis is not a tensor product basis.");
12949d007619Sjeremylt   // LCOV_EXCL_STOP
12959d007619Sjeremylt 
1296d1d35e2fSjeremylt   *grad_1d = basis->grad_1d;
1297e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
12989d007619Sjeremylt }
12999d007619Sjeremylt 
13009d007619Sjeremylt /**
130150c301a5SRezgar Shakeri   @brief Get divergence matrix of a CeedBasis
130250c301a5SRezgar Shakeri 
130350c301a5SRezgar Shakeri   @param basis     CeedBasis
130450c301a5SRezgar Shakeri   @param[out] div  Variable to store divergence matrix
130550c301a5SRezgar Shakeri 
130650c301a5SRezgar Shakeri   @return An error code: 0 - success, otherwise - failure
130750c301a5SRezgar Shakeri 
130850c301a5SRezgar Shakeri   @ref Advanced
130950c301a5SRezgar Shakeri **/
131050c301a5SRezgar Shakeri int CeedBasisGetDiv(CeedBasis basis, const CeedScalar **div) {
131150c301a5SRezgar Shakeri   if (!basis->div)
131250c301a5SRezgar Shakeri     // LCOV_EXCL_START
131350c301a5SRezgar Shakeri     return CeedError(basis->ceed, CEED_ERROR_MINOR,
131450c301a5SRezgar Shakeri                      "CeedBasis does not have divergence matrix.");
131550c301a5SRezgar Shakeri   // LCOV_EXCL_STOP
131650c301a5SRezgar Shakeri 
131750c301a5SRezgar Shakeri   *div = basis->div;
131850c301a5SRezgar Shakeri   return CEED_ERROR_SUCCESS;
131950c301a5SRezgar Shakeri }
132050c301a5SRezgar Shakeri 
132150c301a5SRezgar Shakeri /**
13227a982d89SJeremy L. Thompson   @brief Destroy a CeedBasis
13237a982d89SJeremy L. Thompson 
13247a982d89SJeremy L. Thompson   @param basis CeedBasis to destroy
13257a982d89SJeremy L. Thompson 
13267a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
13277a982d89SJeremy L. Thompson 
13287a982d89SJeremy L. Thompson   @ref User
13297a982d89SJeremy L. Thompson **/
13307a982d89SJeremy L. Thompson int CeedBasisDestroy(CeedBasis *basis) {
13317a982d89SJeremy L. Thompson   int ierr;
13327a982d89SJeremy L. Thompson 
1333d1d35e2fSjeremylt   if (!*basis || --(*basis)->ref_count > 0) return CEED_ERROR_SUCCESS;
13347a982d89SJeremy L. Thompson   if ((*basis)->Destroy) {
13357a982d89SJeremy L. Thompson     ierr = (*basis)->Destroy(*basis); CeedChk(ierr);
13367a982d89SJeremy L. Thompson   }
133734359f16Sjeremylt   if ((*basis)->contract) {
133834359f16Sjeremylt     ierr = CeedTensorContractDestroy(&(*basis)->contract); CeedChk(ierr);
133934359f16Sjeremylt   }
13407a982d89SJeremy L. Thompson   ierr = CeedFree(&(*basis)->interp); CeedChk(ierr);
1341d1d35e2fSjeremylt   ierr = CeedFree(&(*basis)->interp_1d); CeedChk(ierr);
13427a982d89SJeremy L. Thompson   ierr = CeedFree(&(*basis)->grad); CeedChk(ierr);
134350c301a5SRezgar Shakeri   ierr = CeedFree(&(*basis)->div); CeedChk(ierr);
1344d1d35e2fSjeremylt   ierr = CeedFree(&(*basis)->grad_1d); CeedChk(ierr);
1345d1d35e2fSjeremylt   ierr = CeedFree(&(*basis)->q_ref_1d); CeedChk(ierr);
1346d1d35e2fSjeremylt   ierr = CeedFree(&(*basis)->q_weight_1d); CeedChk(ierr);
13477a982d89SJeremy L. Thompson   ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr);
13487a982d89SJeremy L. Thompson   ierr = CeedFree(basis); CeedChk(ierr);
1349e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
13507a982d89SJeremy L. Thompson }
13517a982d89SJeremy L. Thompson 
13527a982d89SJeremy L. Thompson /**
1353b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre quadrature
1354b11c1e72Sjeremylt 
1355b11c1e72Sjeremylt   @param Q                 Number of quadrature points (integrates polynomials of
1356b11c1e72Sjeremylt                              degree 2*Q-1 exactly)
1357d1d35e2fSjeremylt   @param[out] q_ref_1d     Array of length Q to hold the abscissa on [-1, 1]
1358d1d35e2fSjeremylt   @param[out] q_weight_1d  Array of length Q to hold the weights
1359b11c1e72Sjeremylt 
1360b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1361dfdf5a53Sjeremylt 
1362dfdf5a53Sjeremylt   @ref Utility
1363b11c1e72Sjeremylt **/
1364d1d35e2fSjeremylt int CeedGaussQuadrature(CeedInt Q, CeedScalar *q_ref_1d,
1365d1d35e2fSjeremylt                         CeedScalar *q_weight_1d) {
1366d7b241e6Sjeremylt   // Allocate
1367d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0*atan(1.0);
1368d1d35e2fSjeremylt   // Build q_ref_1d, q_weight_1d
1369d7b241e6Sjeremylt   for (int i = 0; i <= Q/2; i++) {
1370d7b241e6Sjeremylt     // Guess
1371d7b241e6Sjeremylt     xi = cos(PI*(CeedScalar)(2*i+1)/((CeedScalar)(2*Q)));
1372d7b241e6Sjeremylt     // Pn(xi)
1373d7b241e6Sjeremylt     P0 = 1.0;
1374d7b241e6Sjeremylt     P1 = xi;
1375d7b241e6Sjeremylt     P2 = 0.0;
1376d7b241e6Sjeremylt     for (int j = 2; j <= Q; j++) {
1377d7b241e6Sjeremylt       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1378d7b241e6Sjeremylt       P0 = P1;
1379d7b241e6Sjeremylt       P1 = P2;
1380d7b241e6Sjeremylt     }
1381d7b241e6Sjeremylt     // First Newton Step
1382d7b241e6Sjeremylt     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1383d7b241e6Sjeremylt     xi = xi-P2/dP2;
1384d7b241e6Sjeremylt     // Newton to convergence
13850e4d4210Sjeremylt     for (int k=0; k<100 && fabs(P2)>10*CEED_EPSILON; k++) {
1386d7b241e6Sjeremylt       P0 = 1.0;
1387d7b241e6Sjeremylt       P1 = xi;
1388d7b241e6Sjeremylt       for (int j = 2; j <= Q; j++) {
1389d7b241e6Sjeremylt         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1390d7b241e6Sjeremylt         P0 = P1;
1391d7b241e6Sjeremylt         P1 = P2;
1392d7b241e6Sjeremylt       }
1393d7b241e6Sjeremylt       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1394d7b241e6Sjeremylt       xi = xi-P2/dP2;
1395d7b241e6Sjeremylt     }
1396d7b241e6Sjeremylt     // Save xi, wi
1397d7b241e6Sjeremylt     wi = 2.0/((1.0-xi*xi)*dP2*dP2);
1398d1d35e2fSjeremylt     q_weight_1d[i] = wi;
1399d1d35e2fSjeremylt     q_weight_1d[Q-1-i] = wi;
1400d1d35e2fSjeremylt     q_ref_1d[i] = -xi;
1401d1d35e2fSjeremylt     q_ref_1d[Q-1-i]= xi;
1402d7b241e6Sjeremylt   }
1403e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1404d7b241e6Sjeremylt }
1405d7b241e6Sjeremylt 
1406b11c1e72Sjeremylt /**
1407b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre-Lobatto quadrature
1408b11c1e72Sjeremylt 
1409b11c1e72Sjeremylt   @param Q                 Number of quadrature points (integrates polynomials of
1410b11c1e72Sjeremylt                              degree 2*Q-3 exactly)
1411d1d35e2fSjeremylt   @param[out] q_ref_1d     Array of length Q to hold the abscissa on [-1, 1]
1412d1d35e2fSjeremylt   @param[out] q_weight_1d  Array of length Q to hold the weights
1413b11c1e72Sjeremylt 
1414b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1415dfdf5a53Sjeremylt 
1416dfdf5a53Sjeremylt   @ref Utility
1417b11c1e72Sjeremylt **/
1418d1d35e2fSjeremylt int CeedLobattoQuadrature(CeedInt Q, CeedScalar *q_ref_1d,
1419d1d35e2fSjeremylt                           CeedScalar *q_weight_1d) {
1420d7b241e6Sjeremylt   // Allocate
1421d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0*atan(1.0);
1422d1d35e2fSjeremylt   // Build q_ref_1d, q_weight_1d
1423d7b241e6Sjeremylt   // Set endpoints
142430a100c3SJed Brown   if (Q < 2)
1425b0d62198Sjeremylt     // LCOV_EXCL_START
1426e15f9bd0SJeremy L Thompson     return CeedError(NULL, CEED_ERROR_DIMENSION,
14277ed52d01Sjeremylt                      "Cannot create Lobatto quadrature with Q=%d < 2 points", Q);
1428b0d62198Sjeremylt   // LCOV_EXCL_STOP
1429d7b241e6Sjeremylt   wi = 2.0/((CeedScalar)(Q*(Q-1)));
1430d1d35e2fSjeremylt   if (q_weight_1d) {
1431d1d35e2fSjeremylt     q_weight_1d[0] = wi;
1432d1d35e2fSjeremylt     q_weight_1d[Q-1] = wi;
1433d7b241e6Sjeremylt   }
1434d1d35e2fSjeremylt   q_ref_1d[0] = -1.0;
1435d1d35e2fSjeremylt   q_ref_1d[Q-1] = 1.0;
1436d7b241e6Sjeremylt   // Interior
1437d7b241e6Sjeremylt   for (int i = 1; i <= (Q-1)/2; i++) {
1438d7b241e6Sjeremylt     // Guess
1439d7b241e6Sjeremylt     xi = cos(PI*(CeedScalar)(i)/(CeedScalar)(Q-1));
1440d7b241e6Sjeremylt     // Pn(xi)
1441d7b241e6Sjeremylt     P0 = 1.0;
1442d7b241e6Sjeremylt     P1 = xi;
1443d7b241e6Sjeremylt     P2 = 0.0;
1444d7b241e6Sjeremylt     for (int j = 2; j < Q; j++) {
1445d7b241e6Sjeremylt       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1446d7b241e6Sjeremylt       P0 = P1;
1447d7b241e6Sjeremylt       P1 = P2;
1448d7b241e6Sjeremylt     }
1449d7b241e6Sjeremylt     // First Newton step
1450d7b241e6Sjeremylt     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1451d7b241e6Sjeremylt     d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
1452d7b241e6Sjeremylt     xi = xi-dP2/d2P2;
1453d7b241e6Sjeremylt     // Newton to convergence
14540e4d4210Sjeremylt     for (int k=0; k<100 && fabs(dP2)>10*CEED_EPSILON; k++) {
1455d7b241e6Sjeremylt       P0 = 1.0;
1456d7b241e6Sjeremylt       P1 = xi;
1457d7b241e6Sjeremylt       for (int j = 2; j < Q; j++) {
1458d7b241e6Sjeremylt         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1459d7b241e6Sjeremylt         P0 = P1;
1460d7b241e6Sjeremylt         P1 = P2;
1461d7b241e6Sjeremylt       }
1462d7b241e6Sjeremylt       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1463d7b241e6Sjeremylt       d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
1464d7b241e6Sjeremylt       xi = xi-dP2/d2P2;
1465d7b241e6Sjeremylt     }
1466d7b241e6Sjeremylt     // Save xi, wi
1467d7b241e6Sjeremylt     wi = 2.0/(((CeedScalar)(Q*(Q-1)))*P2*P2);
1468d1d35e2fSjeremylt     if (q_weight_1d) {
1469d1d35e2fSjeremylt       q_weight_1d[i] = wi;
1470d1d35e2fSjeremylt       q_weight_1d[Q-1-i] = wi;
1471d7b241e6Sjeremylt     }
1472d1d35e2fSjeremylt     q_ref_1d[i] = -xi;
1473d1d35e2fSjeremylt     q_ref_1d[Q-1-i]= xi;
1474d7b241e6Sjeremylt   }
1475e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1476d7b241e6Sjeremylt }
1477d7b241e6Sjeremylt 
1478dfdf5a53Sjeremylt /**
147995bb1877Svaleriabarra   @brief Return QR Factorization of a matrix
1480b11c1e72Sjeremylt 
148177645d7bSjeremylt   @param ceed         A Ceed context for error handling
148252bfb9bbSJeremy L Thompson   @param[in,out] mat  Row-major matrix to be factorized in place
148352bfb9bbSJeremy L Thompson   @param[in,out] tau  Vector of length m of scaling factors
1484b11c1e72Sjeremylt   @param m            Number of rows
1485b11c1e72Sjeremylt   @param n            Number of columns
1486b11c1e72Sjeremylt 
1487b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1488dfdf5a53Sjeremylt 
1489dfdf5a53Sjeremylt   @ref Utility
1490b11c1e72Sjeremylt **/
1491a7bd39daSjeremylt int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau,
1492d7b241e6Sjeremylt                         CeedInt m, CeedInt n) {
1493d7b241e6Sjeremylt   CeedScalar v[m];
1494d7b241e6Sjeremylt 
1495a7bd39daSjeremylt   // Check m >= n
1496a7bd39daSjeremylt   if (n > m)
1497c042f62fSJeremy L Thompson     // LCOV_EXCL_START
1498e15f9bd0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_UNSUPPORTED,
1499e15f9bd0SJeremy L Thompson                      "Cannot compute QR factorization with n > m");
1500c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
1501a7bd39daSjeremylt 
150252bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++) {
1503bde37e8cSJed Brown     if (i >= m-1) { // last row of matrix, no reflection needed
1504bde37e8cSJed Brown       tau[i] = 0.;
1505bde37e8cSJed Brown       break;
1506bde37e8cSJed Brown     }
1507d7b241e6Sjeremylt     // Calculate Householder vector, magnitude
1508d7b241e6Sjeremylt     CeedScalar sigma = 0.0;
1509d7b241e6Sjeremylt     v[i] = mat[i+n*i];
151052bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<m; j++) {
1511d7b241e6Sjeremylt       v[j] = mat[i+n*j];
1512d7b241e6Sjeremylt       sigma += v[j] * v[j];
1513d7b241e6Sjeremylt     }
1514d7b241e6Sjeremylt     CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m]
1515673160d7Sjeremylt     CeedScalar R_ii = -copysign(norm, v[i]);
1516673160d7Sjeremylt     v[i] -= R_ii;
1517d7b241e6Sjeremylt     // norm of v[i:m] after modification above and scaling below
1518d7b241e6Sjeremylt     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
1519d7b241e6Sjeremylt     //   tau = 2 / (norm*norm)
1520d7b241e6Sjeremylt     tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma);
15211d102b48SJeremy L Thompson     for (CeedInt j=i+1; j<m; j++)
15221d102b48SJeremy L Thompson       v[j] /= v[i];
1523d7b241e6Sjeremylt 
1524d7b241e6Sjeremylt     // Apply Householder reflector to lower right panel
1525d7b241e6Sjeremylt     CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1);
1526d7b241e6Sjeremylt     // Save v
1527673160d7Sjeremylt     mat[i+n*i] = R_ii;
15281d102b48SJeremy L Thompson     for (CeedInt j=i+1; j<m; j++)
1529d7b241e6Sjeremylt       mat[i+n*j] = v[j];
1530d7b241e6Sjeremylt   }
1531e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1532d7b241e6Sjeremylt }
1533d7b241e6Sjeremylt 
1534b11c1e72Sjeremylt /**
153552bfb9bbSJeremy L Thompson   @brief Return symmetric Schur decomposition of the symmetric matrix mat via
153652bfb9bbSJeremy L Thompson            symmetric QR factorization
153752bfb9bbSJeremy L Thompson 
153877645d7bSjeremylt   @param ceed         A Ceed context for error handling
153952bfb9bbSJeremy L Thompson   @param[in,out] mat  Row-major matrix to be factorized in place
1540460bf743SValeria Barra   @param[out] lambda  Vector of length n of eigenvalues
154152bfb9bbSJeremy L Thompson   @param n            Number of rows/columns
154252bfb9bbSJeremy L Thompson 
154352bfb9bbSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
154452bfb9bbSJeremy L Thompson 
154552bfb9bbSJeremy L Thompson   @ref Utility
154652bfb9bbSJeremy L Thompson **/
154703d18186Sjeremylt CeedPragmaOptimizeOff
154852bfb9bbSJeremy L Thompson int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat,
154952bfb9bbSJeremy L Thompson                                     CeedScalar *lambda, CeedInt n) {
155052bfb9bbSJeremy L Thompson   // Check bounds for clang-tidy
155152bfb9bbSJeremy L Thompson   if (n<2)
1552c042f62fSJeremy L Thompson     // LCOV_EXCL_START
1553e15f9bd0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_UNSUPPORTED,
1554c042f62fSJeremy L Thompson                      "Cannot compute symmetric Schur decomposition of scalars");
1555c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
155652bfb9bbSJeremy L Thompson 
1557673160d7Sjeremylt   CeedScalar v[n-1], tau[n-1], mat_T[n*n];
155852bfb9bbSJeremy L Thompson 
1559673160d7Sjeremylt   // Copy mat to mat_T and set mat to I
1560673160d7Sjeremylt   memcpy(mat_T, mat, n*n*sizeof(mat[0]));
156152bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
156252bfb9bbSJeremy L Thompson     for (CeedInt j=0; j<n; j++)
156352bfb9bbSJeremy L Thompson       mat[j+n*i] = (i==j) ? 1 : 0;
156452bfb9bbSJeremy L Thompson 
156552bfb9bbSJeremy L Thompson   // Reduce to tridiagonal
156652bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n-1; i++) {
156752bfb9bbSJeremy L Thompson     // Calculate Householder vector, magnitude
156852bfb9bbSJeremy L Thompson     CeedScalar sigma = 0.0;
1569673160d7Sjeremylt     v[i] = mat_T[i+n*(i+1)];
157052bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) {
1571673160d7Sjeremylt       v[j] = mat_T[i+n*(j+1)];
157252bfb9bbSJeremy L Thompson       sigma += v[j] * v[j];
157352bfb9bbSJeremy L Thompson     }
157452bfb9bbSJeremy L Thompson     CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:n-1]
1575673160d7Sjeremylt     CeedScalar R_ii = -copysign(norm, v[i]);
1576673160d7Sjeremylt     v[i] -= R_ii;
157752bfb9bbSJeremy L Thompson     // norm of v[i:m] after modification above and scaling below
157852bfb9bbSJeremy L Thompson     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
157952bfb9bbSJeremy L Thompson     //   tau = 2 / (norm*norm)
158080a9ef05SNatalie Beams     tau[i] = i == n - 2 ? 2 : 2 * v[i]*v[i] / (v[i]*v[i] + sigma);
1581fb551037Sjeremylt     for (CeedInt j=i+1; j<n-1; j++)
1582fb551037Sjeremylt       v[j] /= v[i];
158352bfb9bbSJeremy L Thompson 
158452bfb9bbSJeremy L Thompson     // Update sub and super diagonal
158552bfb9bbSJeremy L Thompson     for (CeedInt j=i+2; j<n; j++) {
1586673160d7Sjeremylt       mat_T[i+n*j] = 0; mat_T[j+n*i] = 0;
158752bfb9bbSJeremy L Thompson     }
158852bfb9bbSJeremy L Thompson     // Apply symmetric Householder reflector to lower right panel
1589673160d7Sjeremylt     CeedHouseholderReflect(&mat_T[(i+1)+n*(i+1)], &v[i], tau[i],
159052bfb9bbSJeremy L Thompson                            n-(i+1), n-(i+1), n, 1);
1591673160d7Sjeremylt     CeedHouseholderReflect(&mat_T[(i+1)+n*(i+1)], &v[i], tau[i],
159252bfb9bbSJeremy L Thompson                            n-(i+1), n-(i+1), 1, n);
1593673160d7Sjeremylt 
159452bfb9bbSJeremy L Thompson     // Save v
1595673160d7Sjeremylt     mat_T[i+n*(i+1)] = R_ii;
1596673160d7Sjeremylt     mat_T[(i+1)+n*i] = R_ii;
159752bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) {
1598673160d7Sjeremylt       mat_T[i+n*(j+1)] = v[j];
159952bfb9bbSJeremy L Thompson     }
160052bfb9bbSJeremy L Thompson   }
160152bfb9bbSJeremy L Thompson   // Backwards accumulation of Q
160252bfb9bbSJeremy L Thompson   for (CeedInt i=n-2; i>=0; i--) {
160385cf89eaSjeremylt     if (tau[i] > 0.0) {
160452bfb9bbSJeremy L Thompson       v[i] = 1;
160552bfb9bbSJeremy L Thompson       for (CeedInt j=i+1; j<n-1; j++) {
1606673160d7Sjeremylt         v[j] = mat_T[i+n*(j+1)];
1607673160d7Sjeremylt         mat_T[i+n*(j+1)] = 0;
160852bfb9bbSJeremy L Thompson       }
160952bfb9bbSJeremy L Thompson       CeedHouseholderReflect(&mat[(i+1)+n*(i+1)], &v[i], tau[i],
161052bfb9bbSJeremy L Thompson                              n-(i+1), n-(i+1), n, 1);
161152bfb9bbSJeremy L Thompson     }
161285cf89eaSjeremylt   }
161352bfb9bbSJeremy L Thompson 
161452bfb9bbSJeremy L Thompson   // Reduce sub and super diagonal
1615673160d7Sjeremylt   CeedInt p = 0, q = 0, itr = 0, max_itr = n*n*n*n;
1616673160d7Sjeremylt   CeedScalar tol = CEED_EPSILON;
161752bfb9bbSJeremy L Thompson 
1618673160d7Sjeremylt   while (itr < max_itr) {
161952bfb9bbSJeremy L Thompson     // Update p, q, size of reduced portions of diagonal
162052bfb9bbSJeremy L Thompson     p = 0; q = 0;
162152bfb9bbSJeremy L Thompson     for (CeedInt i=n-2; i>=0; i--) {
1622673160d7Sjeremylt       if (fabs(mat_T[i+n*(i+1)]) < tol)
162352bfb9bbSJeremy L Thompson         q += 1;
162452bfb9bbSJeremy L Thompson       else
162552bfb9bbSJeremy L Thompson         break;
162652bfb9bbSJeremy L Thompson     }
1627673160d7Sjeremylt     for (CeedInt i=0; i<n-q-1; i++) {
1628673160d7Sjeremylt       if (fabs(mat_T[i+n*(i+1)]) < tol)
162952bfb9bbSJeremy L Thompson         p += 1;
163052bfb9bbSJeremy L Thompson       else
163152bfb9bbSJeremy L Thompson         break;
163252bfb9bbSJeremy L Thompson     }
163352bfb9bbSJeremy L Thompson     if (q == n-1) break; // Finished reducing
163452bfb9bbSJeremy L Thompson 
163552bfb9bbSJeremy L Thompson     // Reduce tridiagonal portion
1636673160d7Sjeremylt     CeedScalar t_nn = mat_T[(n-1-q)+n*(n-1-q)],
1637673160d7Sjeremylt                t_nnm1 = mat_T[(n-2-q)+n*(n-1-q)];
1638673160d7Sjeremylt     CeedScalar d = (mat_T[(n-2-q)+n*(n-2-q)] - t_nn)/2;
1639673160d7Sjeremylt     CeedScalar mu = t_nn - t_nnm1*t_nnm1 /
1640673160d7Sjeremylt                     (d + copysign(sqrt(d*d + t_nnm1*t_nnm1), d));
1641673160d7Sjeremylt     CeedScalar x = mat_T[p+n*p] - mu;
1642673160d7Sjeremylt     CeedScalar z = mat_T[p+n*(p+1)];
1643673160d7Sjeremylt     for (CeedInt k=p; k<n-q-1; k++) {
164452bfb9bbSJeremy L Thompson       // Compute Givens rotation
164552bfb9bbSJeremy L Thompson       CeedScalar c = 1, s = 0;
164652bfb9bbSJeremy L Thompson       if (fabs(z) > tol) {
164752bfb9bbSJeremy L Thompson         if (fabs(z) > fabs(x)) {
164852bfb9bbSJeremy L Thompson           CeedScalar tau = -x/z;
164952bfb9bbSJeremy L Thompson           s = 1/sqrt(1+tau*tau), c = s*tau;
165052bfb9bbSJeremy L Thompson         } else {
165152bfb9bbSJeremy L Thompson           CeedScalar tau = -z/x;
165252bfb9bbSJeremy L Thompson           c = 1/sqrt(1+tau*tau), s = c*tau;
165352bfb9bbSJeremy L Thompson         }
165452bfb9bbSJeremy L Thompson       }
165552bfb9bbSJeremy L Thompson 
165652bfb9bbSJeremy L Thompson       // Apply Givens rotation to T
1657673160d7Sjeremylt       CeedGivensRotation(mat_T, c, s, CEED_NOTRANSPOSE, k, k+1, n, n);
1658673160d7Sjeremylt       CeedGivensRotation(mat_T, c, s, CEED_TRANSPOSE, k, k+1, n, n);
165952bfb9bbSJeremy L Thompson 
166052bfb9bbSJeremy L Thompson       // Apply Givens rotation to Q
166152bfb9bbSJeremy L Thompson       CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k+1, n, n);
166252bfb9bbSJeremy L Thompson 
166352bfb9bbSJeremy L Thompson       // Update x, z
166452bfb9bbSJeremy L Thompson       if (k < n-q-2) {
1665673160d7Sjeremylt         x = mat_T[k+n*(k+1)];
1666673160d7Sjeremylt         z = mat_T[k+n*(k+2)];
166752bfb9bbSJeremy L Thompson       }
166852bfb9bbSJeremy L Thompson     }
166952bfb9bbSJeremy L Thompson     itr++;
167052bfb9bbSJeremy L Thompson   }
1671673160d7Sjeremylt 
167252bfb9bbSJeremy L Thompson   // Save eigenvalues
167352bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
1674673160d7Sjeremylt     lambda[i] = mat_T[i+n*i];
167552bfb9bbSJeremy L Thompson 
167652bfb9bbSJeremy L Thompson   // Check convergence
1677673160d7Sjeremylt   if (itr == max_itr && q < n-1)
1678c042f62fSJeremy L Thompson     // LCOV_EXCL_START
1679e15f9bd0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_MINOR,
1680e15f9bd0SJeremy L Thompson                      "Symmetric QR failed to converge");
1681c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
1682e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
168352bfb9bbSJeremy L Thompson }
168403d18186Sjeremylt CeedPragmaOptimizeOn
168552bfb9bbSJeremy L Thompson 
168652bfb9bbSJeremy L Thompson /**
168752bfb9bbSJeremy L Thompson   @brief Return Simultaneous Diagonalization of two matrices. This solves the
168852bfb9bbSJeremy L Thompson            generalized eigenvalue problem A x = lambda B x, where A and B
168952bfb9bbSJeremy L Thompson            are symmetric and B is positive definite. We generate the matrix X
169052bfb9bbSJeremy L Thompson            and vector Lambda such that X^T A X = Lambda and X^T B X = I. This
169152bfb9bbSJeremy L Thompson            is equivalent to the LAPACK routine 'sygv' with TYPE = 1.
169252bfb9bbSJeremy L Thompson 
169377645d7bSjeremylt   @param ceed         A Ceed context for error handling
1694d1d35e2fSjeremylt   @param[in] mat_A    Row-major matrix to be factorized with eigenvalues
1695d1d35e2fSjeremylt   @param[in] mat_B    Row-major matrix to be factorized to identity
1696d3331725Sjeremylt   @param[out] mat_X   Row-major orthogonal matrix
1697460bf743SValeria Barra   @param[out] lambda  Vector of length n of generalized eigenvalues
169852bfb9bbSJeremy L Thompson   @param n            Number of rows/columns
169952bfb9bbSJeremy L Thompson 
170052bfb9bbSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
170152bfb9bbSJeremy L Thompson 
170252bfb9bbSJeremy L Thompson   @ref Utility
170352bfb9bbSJeremy L Thompson **/
170403d18186Sjeremylt CeedPragmaOptimizeOff
1705d1d35e2fSjeremylt int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *mat_A,
1706d3331725Sjeremylt                                     CeedScalar *mat_B, CeedScalar *mat_X,
170752bfb9bbSJeremy L Thompson                                     CeedScalar *lambda, CeedInt n) {
170852bfb9bbSJeremy L Thompson   int ierr;
1709d3331725Sjeremylt   CeedScalar *mat_C, *mat_G, *vec_D;
171078464608Sjeremylt   ierr = CeedCalloc(n*n, &mat_C); CeedChk(ierr);
171178464608Sjeremylt   ierr = CeedCalloc(n*n, &mat_G); CeedChk(ierr);
1712d3331725Sjeremylt   ierr = CeedCalloc(n, &vec_D); CeedChk(ierr);
171352bfb9bbSJeremy L Thompson 
171452bfb9bbSJeremy L Thompson   // Compute B = G D G^T
171578464608Sjeremylt   memcpy(mat_G, mat_B, n*n*sizeof(mat_B[0]));
1716d3331725Sjeremylt   ierr = CeedSymmetricSchurDecomposition(ceed, mat_G, vec_D, n); CeedChk(ierr);
171752bfb9bbSJeremy L Thompson 
171885cf89eaSjeremylt   // Sort eigenvalues
171985cf89eaSjeremylt   for (CeedInt i=n-1; i>=0; i--)
172085cf89eaSjeremylt     for (CeedInt j=0; j<i; j++) {
172185cf89eaSjeremylt       if (fabs(vec_D[j]) > fabs(vec_D[j+1])) {
172285cf89eaSjeremylt         CeedScalar temp;
172385cf89eaSjeremylt         temp = vec_D[j]; vec_D[j] = vec_D[j+1]; vec_D[j+1] = temp;
172485cf89eaSjeremylt         for (CeedInt k=0; k<n; k++) {
172585cf89eaSjeremylt           temp = mat_G[k*n+j]; mat_G[k*n+j] = mat_G[k*n+j+1]; mat_G[k*n+j+1] = temp;
172685cf89eaSjeremylt         }
172785cf89eaSjeremylt       }
172885cf89eaSjeremylt     }
172985cf89eaSjeremylt 
1730fb551037Sjeremylt   // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T
1731fb551037Sjeremylt   //           = D^-1/2 G^T A G D^-1/2
1732d3331725Sjeremylt   // -- D = D^-1/2
173352bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
1734d3331725Sjeremylt     vec_D[i] = 1./sqrt(vec_D[i]);
1735d3331725Sjeremylt   // -- G = G D^-1/2
1736d3331725Sjeremylt   // -- C = D^-1/2 G^T
1737d3331725Sjeremylt   for (CeedInt i=0; i<n; i++)
1738d3331725Sjeremylt     for (CeedInt j=0; j<n; j++) {
1739673160d7Sjeremylt       mat_G[i*n+j] *= vec_D[j];
1740673160d7Sjeremylt       mat_C[j*n+i]  = mat_G[i*n+j];
1741d3331725Sjeremylt     }
1742673160d7Sjeremylt   // -- X = (D^-1/2 G^T) A
1743d1d35e2fSjeremylt   ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)mat_C,
1744d3331725Sjeremylt                             (const CeedScalar *)mat_A, mat_X, n, n, n);
17459289e5bfSjeremylt   CeedChk(ierr);
1746673160d7Sjeremylt   // -- C = (D^-1/2 G^T A) (G D^-1/2)
1747d3331725Sjeremylt   ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)mat_X,
174878464608Sjeremylt                             (const CeedScalar *)mat_G, mat_C, n, n, n);
17499289e5bfSjeremylt   CeedChk(ierr);
175052bfb9bbSJeremy L Thompson 
175152bfb9bbSJeremy L Thompson   // Compute Q^T C Q = lambda
1752d1d35e2fSjeremylt   ierr = CeedSymmetricSchurDecomposition(ceed, mat_C, lambda, n); CeedChk(ierr);
175352bfb9bbSJeremy L Thompson 
175485cf89eaSjeremylt   // Sort eigenvalues
175585cf89eaSjeremylt   for (CeedInt i=n-1; i>=0; i--)
175685cf89eaSjeremylt     for (CeedInt j=0; j<i; j++) {
175785cf89eaSjeremylt       if (fabs(lambda[j]) > fabs(lambda[j+1])) {
175885cf89eaSjeremylt         CeedScalar temp;
175985cf89eaSjeremylt         temp = lambda[j]; lambda[j] = lambda[j+1]; lambda[j+1] = temp;
176085cf89eaSjeremylt         for (CeedInt k=0; k<n; k++) {
176185cf89eaSjeremylt           temp = mat_C[k*n+j]; mat_C[k*n+j] = mat_C[k*n+j+1]; mat_C[k*n+j+1] = temp;
176285cf89eaSjeremylt         }
176385cf89eaSjeremylt       }
176485cf89eaSjeremylt     }
176585cf89eaSjeremylt 
1766d3331725Sjeremylt   // Set X = (G D^1/2)^-T Q
1767fb551037Sjeremylt   //       = G D^-1/2 Q
176878464608Sjeremylt   ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)mat_G,
1769d3331725Sjeremylt                             (const CeedScalar *)mat_C, mat_X, n, n, n);
17709289e5bfSjeremylt   CeedChk(ierr);
177178464608Sjeremylt 
177278464608Sjeremylt   // Cleanup
177378464608Sjeremylt   ierr = CeedFree(&mat_C); CeedChk(ierr);
177478464608Sjeremylt   ierr = CeedFree(&mat_G); CeedChk(ierr);
1775d3331725Sjeremylt   ierr = CeedFree(&vec_D); CeedChk(ierr);
1776e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
177752bfb9bbSJeremy L Thompson }
177803d18186Sjeremylt CeedPragmaOptimizeOn
177952bfb9bbSJeremy L Thompson 
1780d7b241e6Sjeremylt /// @}
1781