xref: /libCEED/rust/libceed-sys/c-src/interface/ceed-basis.c (revision 6e15d496119073f7bc35f61df7ac100e9376600e)
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 /**
293*6e15d496SJeremy L Thompson   @brief Estimate number of FLOPs required to apply CeedBasis in t_mode and eval_mode
294*6e15d496SJeremy L Thompson 
295*6e15d496SJeremy L Thompson   @param basis     Basis to estimate FLOPs for
296*6e15d496SJeremy L Thompson   @param t_mode    Apply basis or transpose
297*6e15d496SJeremy L Thompson   @param eval_mode Basis evaluation mode
298*6e15d496SJeremy L Thompson   @param flops     Address of variable to hold FLOPs estimate
299*6e15d496SJeremy L Thompson 
300*6e15d496SJeremy L Thompson   @ref Backend
301*6e15d496SJeremy L Thompson **/
302*6e15d496SJeremy L Thompson int CeedBasisGetFlopsEstimate(CeedBasis basis, CeedTransposeMode t_mode,
303*6e15d496SJeremy L Thompson                               CeedEvalMode eval_mode, CeedInt *flops) {
304*6e15d496SJeremy L Thompson   int ierr;
305*6e15d496SJeremy L Thompson   bool is_tensor;
306*6e15d496SJeremy L Thompson 
307*6e15d496SJeremy L Thompson   ierr = CeedBasisIsTensor(basis, &is_tensor); CeedChk(ierr);
308*6e15d496SJeremy L Thompson   if (is_tensor) {
309*6e15d496SJeremy L Thompson     CeedInt dim, num_comp, P_1d, Q_1d;
310*6e15d496SJeremy L Thompson     ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr);
311*6e15d496SJeremy L Thompson     ierr = CeedBasisGetNumComponents(basis, &num_comp); CeedChk(ierr);
312*6e15d496SJeremy L Thompson     ierr = CeedBasisGetNumNodes1D(basis, &P_1d);  CeedChk(ierr);
313*6e15d496SJeremy L Thompson     ierr = CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d);  CeedChk(ierr);
314*6e15d496SJeremy L Thompson     if (t_mode == CEED_TRANSPOSE) {
315*6e15d496SJeremy L Thompson       P_1d = Q_1d; Q_1d = P_1d;
316*6e15d496SJeremy L Thompson     }
317*6e15d496SJeremy L Thompson     CeedInt tensor_flops = 0, pre = num_comp * CeedIntPow(P_1d, dim-1), post = 1;
318*6e15d496SJeremy L Thompson     for (CeedInt d = 0; d < dim; d++) {
319*6e15d496SJeremy L Thompson       tensor_flops += 2 * pre * P_1d * post * Q_1d;
320*6e15d496SJeremy L Thompson       pre /= P_1d;
321*6e15d496SJeremy L Thompson       post *= Q_1d;
322*6e15d496SJeremy L Thompson     }
323*6e15d496SJeremy L Thompson     switch (eval_mode) {
324*6e15d496SJeremy L Thompson     case CEED_EVAL_NONE:   *flops = 0; break;
325*6e15d496SJeremy L Thompson     case CEED_EVAL_INTERP: *flops = tensor_flops; break;
326*6e15d496SJeremy L Thompson     case CEED_EVAL_GRAD:   *flops = tensor_flops * 2; break;
327*6e15d496SJeremy L Thompson     case CEED_EVAL_DIV:
328*6e15d496SJeremy L Thompson       // LCOV_EXCL_START
329*6e15d496SJeremy L Thompson       return CeedError(basis->ceed, CEED_ERROR_INCOMPATIBLE,
330*6e15d496SJeremy L Thompson                        "Tensor CEED_EVAL_DIV not supported"); break;
331*6e15d496SJeremy L Thompson     case CEED_EVAL_CURL:
332*6e15d496SJeremy L Thompson       return CeedError(basis->ceed, CEED_ERROR_INCOMPATIBLE,
333*6e15d496SJeremy L Thompson                        "Tensor CEED_EVAL_CURL not supported"); break;
334*6e15d496SJeremy L Thompson     // LCOV_EXCL_STOP
335*6e15d496SJeremy L Thompson     case CEED_EVAL_WEIGHT: *flops = dim * CeedIntPow(Q_1d, dim); break;
336*6e15d496SJeremy L Thompson     }
337*6e15d496SJeremy L Thompson   } else {
338*6e15d496SJeremy L Thompson     CeedInt dim, num_comp, num_nodes, num_qpts, Q_comp;
339*6e15d496SJeremy L Thompson     ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr);
340*6e15d496SJeremy L Thompson     ierr = CeedBasisGetNumComponents(basis, &num_comp); CeedChk(ierr);
341*6e15d496SJeremy L Thompson     ierr = CeedBasisGetNumNodes(basis, &num_nodes); CeedChk(ierr);
342*6e15d496SJeremy L Thompson     ierr = CeedBasisGetNumQuadraturePoints(basis, &num_qpts); CeedChk(ierr);
343*6e15d496SJeremy L Thompson     ierr = CeedBasisGetNumQuadratureComponents(basis, &Q_comp); CeedChk(ierr);
344*6e15d496SJeremy L Thompson     switch (eval_mode) {
345*6e15d496SJeremy L Thompson     case CEED_EVAL_NONE:   *flops = 0; break;
346*6e15d496SJeremy L Thompson     case CEED_EVAL_INTERP: *flops = num_nodes * num_qpts * num_comp; break;
347*6e15d496SJeremy L Thompson     case CEED_EVAL_GRAD:   *flops = num_nodes * num_qpts * num_comp * dim; break;
348*6e15d496SJeremy L Thompson     case CEED_EVAL_DIV:    *flops = num_nodes * num_qpts; break;
349*6e15d496SJeremy L Thompson     case CEED_EVAL_CURL:   *flops = num_nodes * num_qpts * dim; break;
350*6e15d496SJeremy L Thompson     case CEED_EVAL_WEIGHT: *flops = 0; break;
351*6e15d496SJeremy L Thompson     }
352*6e15d496SJeremy L Thompson   }
353*6e15d496SJeremy L Thompson 
354*6e15d496SJeremy L Thompson   return CEED_ERROR_SUCCESS;
355*6e15d496SJeremy L Thompson }
356*6e15d496SJeremy L Thompson 
357*6e15d496SJeremy 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
49350c301a5SRezgar Shakeri   CeedElemTopology topo = dim == 1 ? CEED_TOPOLOGY_LINE
49450c301a5SRezgar Shakeri                           : dim == 2 ? CEED_TOPOLOGY_QUAD
49550c301a5SRezgar Shakeri                           : CEED_TOPOLOGY_HEX;
496e15f9bd0SJeremy L Thompson 
497d7b241e6Sjeremylt   ierr = CeedCalloc(1, basis); CeedChk(ierr);
498d7b241e6Sjeremylt   (*basis)->ceed = ceed;
4999560d06aSjeremylt   ierr = CeedReference(ceed); CeedChk(ierr);
500d1d35e2fSjeremylt   (*basis)->ref_count = 1;
501d1d35e2fSjeremylt   (*basis)->tensor_basis = 1;
502d7b241e6Sjeremylt   (*basis)->dim = dim;
503d99fa3c5SJeremy L Thompson   (*basis)->topo = topo;
504d1d35e2fSjeremylt   (*basis)->num_comp = num_comp;
505d1d35e2fSjeremylt   (*basis)->P_1d = P_1d;
506d1d35e2fSjeremylt   (*basis)->Q_1d = Q_1d;
507d1d35e2fSjeremylt   (*basis)->P = CeedIntPow(P_1d, dim);
508d1d35e2fSjeremylt   (*basis)->Q = CeedIntPow(Q_1d, dim);
50950c301a5SRezgar Shakeri   (*basis)->Q_comp = 1;
51050c301a5SRezgar Shakeri   (*basis)->basis_space = 1; // 1 for H^1 space
511ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(Q_1d, &(*basis)->q_ref_1d); CeedChk(ierr);
512ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(Q_1d, &(*basis)->q_weight_1d); CeedChk(ierr);
513ff3a0f91SJeremy L Thompson   if (q_ref_1d) memcpy((*basis)->q_ref_1d, q_ref_1d, Q_1d*sizeof(q_ref_1d[0]));
514ff3a0f91SJeremy L Thompson   if (q_weight_1d) memcpy((*basis)->q_weight_1d, q_weight_1d,
515ff3a0f91SJeremy L Thompson                             Q_1d*sizeof(q_weight_1d[0]));
516ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(Q_1d*P_1d, &(*basis)->interp_1d); CeedChk(ierr);
517ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(Q_1d*P_1d, &(*basis)->grad_1d); CeedChk(ierr);
518ff3a0f91SJeremy L Thompson   if (interp_1d) memcpy((*basis)->interp_1d, interp_1d,
519ff3a0f91SJeremy L Thompson                           Q_1d*P_1d*sizeof(interp_1d[0]));
520ff3a0f91SJeremy L Thompson   if (grad_1d) memcpy((*basis)->grad_1d, grad_1d, Q_1d*P_1d*sizeof(grad_1d[0]));
521d1d35e2fSjeremylt   ierr = ceed->BasisCreateTensorH1(dim, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d,
522d1d35e2fSjeremylt                                    q_weight_1d, *basis); CeedChk(ierr);
523e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
524d7b241e6Sjeremylt }
525d7b241e6Sjeremylt 
526b11c1e72Sjeremylt /**
52795bb1877Svaleriabarra   @brief Create a tensor-product Lagrange basis
528b11c1e72Sjeremylt 
529b11c1e72Sjeremylt   @param ceed        A Ceed object where the CeedBasis will be created
530b11c1e72Sjeremylt   @param dim         Topological dimension of element
531d1d35e2fSjeremylt   @param num_comp      Number of field components (1 for scalar fields)
532b11c1e72Sjeremylt   @param P           Number of Gauss-Lobatto nodes in one dimension.  The
533b11c1e72Sjeremylt                        polynomial degree of the resulting Q_k element is k=P-1.
534b11c1e72Sjeremylt   @param Q           Number of quadrature points in one dimension.
535d1d35e2fSjeremylt   @param quad_mode   Distribution of the Q quadrature points (affects order of
536b11c1e72Sjeremylt                        accuracy for the quadrature)
537b11c1e72Sjeremylt   @param[out] basis  Address of the variable where the newly created
538b11c1e72Sjeremylt                        CeedBasis will be stored.
539b11c1e72Sjeremylt 
540b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
541dfdf5a53Sjeremylt 
5427a982d89SJeremy L. Thompson   @ref User
543b11c1e72Sjeremylt **/
544d1d35e2fSjeremylt int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt num_comp,
545d1d35e2fSjeremylt                                     CeedInt P, CeedInt Q, CeedQuadMode quad_mode,
546692c2638Sjeremylt                                     CeedBasis *basis) {
547d7b241e6Sjeremylt   // Allocate
548e15f9bd0SJeremy L Thompson   int ierr, ierr2, i, j, k;
549d1d35e2fSjeremylt   CeedScalar c1, c2, c3, c4, dx, *nodes, *interp_1d, *grad_1d, *q_ref_1d,
550d1d35e2fSjeremylt              *q_weight_1d;
5514d537eeaSYohann 
5524d537eeaSYohann   if (dim<1)
553c042f62fSJeremy L Thompson     // LCOV_EXCL_START
554e15f9bd0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION,
555e15f9bd0SJeremy L Thompson                      "Basis dimension must be a positive value");
556c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
5574d537eeaSYohann 
558e15f9bd0SJeremy L Thompson   // Get Nodes and Weights
559d1d35e2fSjeremylt   ierr = CeedCalloc(P*Q, &interp_1d); CeedChk(ierr);
560d1d35e2fSjeremylt   ierr = CeedCalloc(P*Q, &grad_1d); CeedChk(ierr);
561d7b241e6Sjeremylt   ierr = CeedCalloc(P, &nodes); CeedChk(ierr);
562d1d35e2fSjeremylt   ierr = CeedCalloc(Q, &q_ref_1d); CeedChk(ierr);
563d1d35e2fSjeremylt   ierr = CeedCalloc(Q, &q_weight_1d); CeedChk(ierr);
564e15f9bd0SJeremy L Thompson   ierr = CeedLobattoQuadrature(P, nodes, NULL);
565e15f9bd0SJeremy L Thompson   if (ierr) { goto cleanup; } CeedChk(ierr);
566d1d35e2fSjeremylt   switch (quad_mode) {
567d7b241e6Sjeremylt   case CEED_GAUSS:
568d1d35e2fSjeremylt     ierr = CeedGaussQuadrature(Q, q_ref_1d, q_weight_1d);
569d7b241e6Sjeremylt     break;
570d7b241e6Sjeremylt   case CEED_GAUSS_LOBATTO:
571d1d35e2fSjeremylt     ierr = CeedLobattoQuadrature(Q, q_ref_1d, q_weight_1d);
572d7b241e6Sjeremylt     break;
573d7b241e6Sjeremylt   }
574e15f9bd0SJeremy L Thompson   if (ierr) { goto cleanup; } CeedChk(ierr);
575e15f9bd0SJeremy L Thompson 
576d7b241e6Sjeremylt   // Build B, D matrix
577d7b241e6Sjeremylt   // Fornberg, 1998
578d7b241e6Sjeremylt   for (i = 0; i  < Q; i++) {
579d7b241e6Sjeremylt     c1 = 1.0;
580d1d35e2fSjeremylt     c3 = nodes[0] - q_ref_1d[i];
581d1d35e2fSjeremylt     interp_1d[i*P+0] = 1.0;
582d7b241e6Sjeremylt     for (j = 1; j < P; j++) {
583d7b241e6Sjeremylt       c2 = 1.0;
584d7b241e6Sjeremylt       c4 = c3;
585d1d35e2fSjeremylt       c3 = nodes[j] - q_ref_1d[i];
586d7b241e6Sjeremylt       for (k = 0; k < j; k++) {
587d7b241e6Sjeremylt         dx = nodes[j] - nodes[k];
588d7b241e6Sjeremylt         c2 *= dx;
589d7b241e6Sjeremylt         if (k == j - 1) {
590d1d35e2fSjeremylt           grad_1d[i*P + j] = c1*(interp_1d[i*P + k] - c4*grad_1d[i*P + k]) / c2;
591d1d35e2fSjeremylt           interp_1d[i*P + j] = - c1*c4*interp_1d[i*P + k] / c2;
592d7b241e6Sjeremylt         }
593d1d35e2fSjeremylt         grad_1d[i*P + k] = (c3*grad_1d[i*P + k] - interp_1d[i*P + k]) / dx;
594d1d35e2fSjeremylt         interp_1d[i*P + k] = c3*interp_1d[i*P + k] / dx;
595d7b241e6Sjeremylt       }
596d7b241e6Sjeremylt       c1 = c2;
597d7b241e6Sjeremylt     }
598d7b241e6Sjeremylt   }
5999ac7b42eSJeremy L Thompson   // Pass to CeedBasisCreateTensorH1
600d1d35e2fSjeremylt   ierr = CeedBasisCreateTensorH1(ceed, dim, num_comp, P, Q, interp_1d, grad_1d,
6019ac7b42eSJeremy L Thompson                                  q_ref_1d, q_weight_1d, basis); CeedChk(ierr);
602e15f9bd0SJeremy L Thompson cleanup:
603d1d35e2fSjeremylt   ierr2 = CeedFree(&interp_1d); CeedChk(ierr2);
604d1d35e2fSjeremylt   ierr2 = CeedFree(&grad_1d); CeedChk(ierr2);
605e15f9bd0SJeremy L Thompson   ierr2 = CeedFree(&nodes); CeedChk(ierr2);
606d1d35e2fSjeremylt   ierr2 = CeedFree(&q_ref_1d); CeedChk(ierr2);
607d1d35e2fSjeremylt   ierr2 = CeedFree(&q_weight_1d); CeedChk(ierr2);
608e15f9bd0SJeremy L Thompson   CeedChk(ierr);
609e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
610d7b241e6Sjeremylt }
611d7b241e6Sjeremylt 
612b11c1e72Sjeremylt /**
61395bb1877Svaleriabarra   @brief Create a non tensor-product basis for H^1 discretizations
614a8de75f0Sjeremylt 
615a8de75f0Sjeremylt   @param ceed        A Ceed object where the CeedBasis will be created
616a8de75f0Sjeremylt   @param topo        Topology of element, e.g. hypercube, simplex, ect
617d1d35e2fSjeremylt   @param num_comp    Number of field components (1 for scalar fields)
618d1d35e2fSjeremylt   @param num_nodes   Total number of nodes
619d1d35e2fSjeremylt   @param num_qpts    Total number of quadrature points
620d1d35e2fSjeremylt   @param interp      Row-major (num_qpts * num_nodes) matrix expressing the values of
6218795c945Sjeremylt                        nodal basis functions at quadrature points
622d1d35e2fSjeremylt   @param grad        Row-major (num_qpts * dim * num_nodes) matrix expressing
6238795c945Sjeremylt                        derivatives of nodal basis functions at quadrature points
624d1d35e2fSjeremylt   @param q_ref       Array of length num_qpts holding the locations of quadrature
62550c301a5SRezgar Shakeri                        points on the reference element
626d1d35e2fSjeremylt   @param q_weight    Array of length num_qpts holding the quadrature weights on the
627a8de75f0Sjeremylt                        reference element
628a8de75f0Sjeremylt   @param[out] basis  Address of the variable where the newly created
629a8de75f0Sjeremylt                        CeedBasis will be stored.
630a8de75f0Sjeremylt 
631a8de75f0Sjeremylt   @return An error code: 0 - success, otherwise - failure
632a8de75f0Sjeremylt 
6337a982d89SJeremy L. Thompson   @ref User
634a8de75f0Sjeremylt **/
635d1d35e2fSjeremylt int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt num_comp,
636d1d35e2fSjeremylt                       CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
637d1d35e2fSjeremylt                       const CeedScalar *grad, const CeedScalar *q_ref,
638d1d35e2fSjeremylt                       const CeedScalar *q_weight, CeedBasis *basis) {
639a8de75f0Sjeremylt   int ierr;
640d1d35e2fSjeremylt   CeedInt P = num_nodes, Q = num_qpts, dim = 0;
641a8de75f0Sjeremylt 
6425fe0d4faSjeremylt   if (!ceed->BasisCreateH1) {
6435fe0d4faSjeremylt     Ceed delegate;
644aefd8378Sjeremylt     ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr);
6455fe0d4faSjeremylt 
6465fe0d4faSjeremylt     if (!delegate)
647c042f62fSJeremy L Thompson       // LCOV_EXCL_START
648e15f9bd0SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_UNSUPPORTED,
649e15f9bd0SJeremy L Thompson                        "Backend does not support BasisCreateH1");
650c042f62fSJeremy L Thompson     // LCOV_EXCL_STOP
6515fe0d4faSjeremylt 
652d1d35e2fSjeremylt     ierr = CeedBasisCreateH1(delegate, topo, num_comp, num_nodes,
653d1d35e2fSjeremylt                              num_qpts, interp, grad, q_ref,
654d1d35e2fSjeremylt                              q_weight, basis); CeedChk(ierr);
655e15f9bd0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
6565fe0d4faSjeremylt   }
6575fe0d4faSjeremylt 
658a8de75f0Sjeremylt   ierr = CeedCalloc(1, basis); CeedChk(ierr);
659a8de75f0Sjeremylt 
660a8de75f0Sjeremylt   ierr = CeedBasisGetTopologyDimension(topo, &dim); CeedChk(ierr);
661a8de75f0Sjeremylt 
662a8de75f0Sjeremylt   (*basis)->ceed = ceed;
6639560d06aSjeremylt   ierr = CeedReference(ceed); CeedChk(ierr);
664d1d35e2fSjeremylt   (*basis)->ref_count = 1;
665d1d35e2fSjeremylt   (*basis)->tensor_basis = 0;
666a8de75f0Sjeremylt   (*basis)->dim = dim;
667d99fa3c5SJeremy L Thompson   (*basis)->topo = topo;
668d1d35e2fSjeremylt   (*basis)->num_comp = num_comp;
669a8de75f0Sjeremylt   (*basis)->P = P;
670a8de75f0Sjeremylt   (*basis)->Q = Q;
67150c301a5SRezgar Shakeri   (*basis)->Q_comp = 1;
67250c301a5SRezgar Shakeri   (*basis)->basis_space = 1; // 1 for H^1 space
673ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(Q*dim, &(*basis)->q_ref_1d); CeedChk(ierr);
674ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(Q, &(*basis)->q_weight_1d); CeedChk(ierr);
675ff3a0f91SJeremy L Thompson   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q*dim*sizeof(q_ref[0]));
676ff3a0f91SJeremy L Thompson   if(q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q*sizeof(q_weight[0]));
677ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(Q*P, &(*basis)->interp); CeedChk(ierr);
678ff3a0f91SJeremy L Thompson   ierr = CeedCalloc(dim*Q*P, &(*basis)->grad); CeedChk(ierr);
679ff3a0f91SJeremy L Thompson   if(interp) memcpy((*basis)->interp, interp, Q*P*sizeof(interp[0]));
680ff3a0f91SJeremy L Thompson   if(grad) memcpy((*basis)->grad, grad, dim*Q*P*sizeof(grad[0]));
681d1d35e2fSjeremylt   ierr = ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, q_ref,
682d1d35e2fSjeremylt                              q_weight, *basis); CeedChk(ierr);
683e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
684a8de75f0Sjeremylt }
685a8de75f0Sjeremylt 
686a8de75f0Sjeremylt /**
68750c301a5SRezgar Shakeri   @brief Create a non tensor-product basis for H(div) discretizations
68850c301a5SRezgar Shakeri 
68950c301a5SRezgar Shakeri   @param ceed        A Ceed object where the CeedBasis will be created
69050c301a5SRezgar Shakeri   @param topo        Topology of element (`CEED_TOPOLOGY_QUAD`, `CEED_TOPOLOGY_PRISM`, etc.),
69150c301a5SRezgar Shakeri                      dimension of which is used in some array sizes below
69250c301a5SRezgar Shakeri   @param num_comp    Number of components (usually 1 for vectors in H(div) bases)
69350c301a5SRezgar Shakeri   @param num_nodes   Total number of nodes (dofs per element)
69450c301a5SRezgar Shakeri   @param num_qpts    Total number of quadrature points
69550c301a5SRezgar Shakeri   @param interp      Row-major (dim*num_qpts * num_nodes) matrix expressing the values of
69650c301a5SRezgar Shakeri                        nodal basis functions at quadrature points
69750c301a5SRezgar Shakeri   @param div        Row-major (num_qpts * num_nodes) matrix expressing
69850c301a5SRezgar Shakeri                        divergence of nodal basis functions at quadrature points
69950c301a5SRezgar Shakeri   @param q_ref       Array of length num_qpts holding the locations of quadrature
70050c301a5SRezgar Shakeri                        points on the reference element
70150c301a5SRezgar Shakeri   @param q_weight    Array of length num_qpts holding the quadrature weights on the
70250c301a5SRezgar Shakeri                        reference element
70350c301a5SRezgar Shakeri   @param[out] basis  Address of the variable where the newly created
70450c301a5SRezgar Shakeri                        CeedBasis will be stored.
70550c301a5SRezgar Shakeri 
70650c301a5SRezgar Shakeri   @return An error code: 0 - success, otherwise - failure
70750c301a5SRezgar Shakeri 
70850c301a5SRezgar Shakeri   @ref User
70950c301a5SRezgar Shakeri **/
71050c301a5SRezgar Shakeri int CeedBasisCreateHdiv(Ceed ceed, CeedElemTopology topo, CeedInt num_comp,
71150c301a5SRezgar Shakeri                         CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
71250c301a5SRezgar Shakeri                         const CeedScalar *div, const CeedScalar *q_ref,
71350c301a5SRezgar Shakeri                         const CeedScalar *q_weight, CeedBasis *basis) {
71450c301a5SRezgar Shakeri   int ierr;
71550c301a5SRezgar Shakeri   CeedInt Q = num_qpts, P = num_nodes, dim = 0;
71650c301a5SRezgar Shakeri   ierr = CeedBasisGetTopologyDimension(topo, &dim); CeedChk(ierr);
71750c301a5SRezgar Shakeri   if (!ceed->BasisCreateHdiv) {
71850c301a5SRezgar Shakeri     Ceed delegate;
71950c301a5SRezgar Shakeri     ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr);
72050c301a5SRezgar Shakeri 
72150c301a5SRezgar Shakeri     if (!delegate)
72250c301a5SRezgar Shakeri       // LCOV_EXCL_START
72350c301a5SRezgar Shakeri       return CeedError(ceed, CEED_ERROR_UNSUPPORTED,
72450c301a5SRezgar Shakeri                        "Backend does not implement BasisCreateHdiv");
72550c301a5SRezgar Shakeri     // LCOV_EXCL_STOP
72650c301a5SRezgar Shakeri 
72750c301a5SRezgar Shakeri     ierr = CeedBasisCreateHdiv(delegate, topo, num_comp, num_nodes,
72850c301a5SRezgar Shakeri                                num_qpts, interp, div, q_ref,
72950c301a5SRezgar Shakeri                                q_weight, basis); CeedChk(ierr);
73050c301a5SRezgar Shakeri     return CEED_ERROR_SUCCESS;
73150c301a5SRezgar Shakeri   }
73250c301a5SRezgar Shakeri 
73350c301a5SRezgar Shakeri   ierr = CeedCalloc(1, basis); CeedChk(ierr);
73450c301a5SRezgar Shakeri 
73550c301a5SRezgar Shakeri   (*basis)->ceed = ceed;
73650c301a5SRezgar Shakeri   ierr = CeedReference(ceed); CeedChk(ierr);
73750c301a5SRezgar Shakeri   (*basis)->ref_count = 1;
73850c301a5SRezgar Shakeri   (*basis)->tensor_basis = 0;
73950c301a5SRezgar Shakeri   (*basis)->dim = dim;
74050c301a5SRezgar Shakeri   (*basis)->topo = topo;
74150c301a5SRezgar Shakeri   (*basis)->num_comp = num_comp;
74250c301a5SRezgar Shakeri   (*basis)->P = P;
74350c301a5SRezgar Shakeri   (*basis)->Q = Q;
74450c301a5SRezgar Shakeri   (*basis)->Q_comp = dim;
74550c301a5SRezgar Shakeri   (*basis)->basis_space = 2; // 2 for H(div) space
74650c301a5SRezgar Shakeri   ierr = CeedMalloc(Q*dim, &(*basis)->q_ref_1d); CeedChk(ierr);
74750c301a5SRezgar Shakeri   ierr = CeedMalloc(Q, &(*basis)->q_weight_1d); CeedChk(ierr);
74850c301a5SRezgar Shakeri   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q*dim*sizeof(q_ref[0]));
74950c301a5SRezgar Shakeri   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q*sizeof(q_weight[0]));
75050c301a5SRezgar Shakeri   ierr = CeedMalloc(dim*Q*P, &(*basis)->interp); CeedChk(ierr);
75150c301a5SRezgar Shakeri   ierr = CeedMalloc(Q*P, &(*basis)->div); CeedChk(ierr);
75250c301a5SRezgar Shakeri   if (interp) memcpy((*basis)->interp, interp, dim*Q*P*sizeof(interp[0]));
75350c301a5SRezgar Shakeri   if (div) memcpy((*basis)->div, div, Q*P*sizeof(div[0]));
75450c301a5SRezgar Shakeri   ierr = ceed->BasisCreateHdiv(topo, dim, P, Q, interp, div, q_ref,
75550c301a5SRezgar Shakeri                                q_weight, *basis); CeedChk(ierr);
75650c301a5SRezgar Shakeri   return CEED_ERROR_SUCCESS;
75750c301a5SRezgar Shakeri }
75850c301a5SRezgar Shakeri 
75950c301a5SRezgar Shakeri /**
7609560d06aSjeremylt   @brief Copy the pointer to a CeedBasis. Both pointers should
7619560d06aSjeremylt            be destroyed with `CeedBasisDestroy()`;
7629560d06aSjeremylt            Note: If `*basis_copy` is non-NULL, then it is assumed that
7639560d06aSjeremylt            `*basis_copy` is a pointer to a CeedBasis. This CeedBasis
7649560d06aSjeremylt            will be destroyed if `*basis_copy` is the only
7659560d06aSjeremylt            reference to this CeedBasis.
7669560d06aSjeremylt 
7679560d06aSjeremylt   @param basis            CeedBasis to copy reference to
7689560d06aSjeremylt   @param[out] basis_copy  Variable to store copied reference
7699560d06aSjeremylt 
7709560d06aSjeremylt   @return An error code: 0 - success, otherwise - failure
7719560d06aSjeremylt 
7729560d06aSjeremylt   @ref User
7739560d06aSjeremylt **/
7749560d06aSjeremylt int CeedBasisReferenceCopy(CeedBasis basis, CeedBasis *basis_copy) {
7759560d06aSjeremylt   int ierr;
7769560d06aSjeremylt 
7779560d06aSjeremylt   ierr = CeedBasisReference(basis); CeedChk(ierr);
7789560d06aSjeremylt   ierr = CeedBasisDestroy(basis_copy); CeedChk(ierr);
7799560d06aSjeremylt   *basis_copy = basis;
7809560d06aSjeremylt   return CEED_ERROR_SUCCESS;
7819560d06aSjeremylt }
7829560d06aSjeremylt 
7839560d06aSjeremylt /**
7847a982d89SJeremy L. Thompson   @brief View a CeedBasis
7857a982d89SJeremy L. Thompson 
7867a982d89SJeremy L. Thompson   @param basis   CeedBasis to view
7877a982d89SJeremy L. Thompson   @param stream  Stream to view to, e.g., stdout
7887a982d89SJeremy L. Thompson 
7897a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
7907a982d89SJeremy L. Thompson 
7917a982d89SJeremy L. Thompson   @ref User
7927a982d89SJeremy L. Thompson **/
7937a982d89SJeremy L. Thompson int CeedBasisView(CeedBasis basis, FILE *stream) {
7947a982d89SJeremy L. Thompson   int ierr;
79550c301a5SRezgar Shakeri   CeedFESpace FE_space = basis->basis_space;
79650c301a5SRezgar Shakeri   CeedElemTopology topo = basis->topo;
79750c301a5SRezgar Shakeri   // Print FE space and element topology of the basis
798d1d35e2fSjeremylt   if (basis->tensor_basis) {
79950c301a5SRezgar Shakeri     fprintf(stream, "CeedBasis (%s on a %s element): dim=%d P=%d Q=%d\n",
80050c301a5SRezgar Shakeri             CeedFESpaces[FE_space], CeedElemTopologies[topo],
80150c301a5SRezgar Shakeri             basis->dim, basis->P_1d, basis->Q_1d);
80250c301a5SRezgar Shakeri   } else {
80350c301a5SRezgar Shakeri     fprintf(stream, "CeedBasis (%s on a %s element): dim=%d P=%d Q=%d\n",
80450c301a5SRezgar Shakeri             CeedFESpaces[FE_space], CeedElemTopologies[topo],
80550c301a5SRezgar Shakeri             basis->dim, basis->P, basis->Q);
80650c301a5SRezgar Shakeri   }
80750c301a5SRezgar Shakeri   // Print quadrature data, interpolation/gradient/divergene/curl of the basis
80850c301a5SRezgar Shakeri   if (basis->tensor_basis) { // tensor basis
809d1d35e2fSjeremylt     ierr = CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_ref_1d,
8107a982d89SJeremy L. Thompson                           stream); CeedChk(ierr);
811d1d35e2fSjeremylt     ierr = CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q_1d,
812d1d35e2fSjeremylt                           basis->q_weight_1d, stream); CeedChk(ierr);
813d1d35e2fSjeremylt     ierr = CeedScalarView("interp1d", "\t% 12.8f", basis->Q_1d, basis->P_1d,
814d1d35e2fSjeremylt                           basis->interp_1d, stream); CeedChk(ierr);
815d1d35e2fSjeremylt     ierr = CeedScalarView("grad1d", "\t% 12.8f", basis->Q_1d, basis->P_1d,
816d1d35e2fSjeremylt                           basis->grad_1d, stream); CeedChk(ierr);
81750c301a5SRezgar Shakeri   } else { // non-tensor basis
8187a982d89SJeremy L. Thompson     ierr = CeedScalarView("qref", "\t% 12.8f", 1, basis->Q*basis->dim,
819d1d35e2fSjeremylt                           basis->q_ref_1d,
8207a982d89SJeremy L. Thompson                           stream); CeedChk(ierr);
821d1d35e2fSjeremylt     ierr = CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->q_weight_1d,
8227a982d89SJeremy L. Thompson                           stream); CeedChk(ierr);
82350c301a5SRezgar Shakeri     ierr = CeedScalarView("interp", "\t% 12.8f", basis->Q_comp*basis->Q, basis->P,
8247a982d89SJeremy L. Thompson                           basis->interp, stream); CeedChk(ierr);
82550c301a5SRezgar Shakeri     if (basis->grad) {
8267a982d89SJeremy L. Thompson       ierr = CeedScalarView("grad", "\t% 12.8f", basis->dim*basis->Q, basis->P,
8277a982d89SJeremy L. Thompson                             basis->grad, stream); CeedChk(ierr);
8287a982d89SJeremy L. Thompson     }
82950c301a5SRezgar Shakeri     if (basis->div) {
83050c301a5SRezgar Shakeri       ierr = CeedScalarView("div", "\t% 12.8f", basis->Q, basis->P,
83150c301a5SRezgar Shakeri                             basis->div, stream); CeedChk(ierr);
83250c301a5SRezgar Shakeri     }
83350c301a5SRezgar Shakeri   }
834e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8357a982d89SJeremy L. Thompson }
8367a982d89SJeremy L. Thompson 
8377a982d89SJeremy L. Thompson /**
8387a982d89SJeremy L. Thompson   @brief Apply basis evaluation from nodes to quadrature points or vice versa
8397a982d89SJeremy L. Thompson 
8407a982d89SJeremy L. Thompson   @param basis     CeedBasis to evaluate
841d1d35e2fSjeremylt   @param num_elem  The number of elements to apply the basis evaluation to;
8427a982d89SJeremy L. Thompson                      the backend will specify the ordering in
8434cc79fe7SJed Brown                      CeedElemRestrictionCreateBlocked()
844d1d35e2fSjeremylt   @param t_mode    \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature
8457a982d89SJeremy L. Thompson                      points, \ref CEED_TRANSPOSE to apply the transpose, mapping
8467a982d89SJeremy L. Thompson                      from quadrature points to nodes
847d1d35e2fSjeremylt   @param eval_mode \ref CEED_EVAL_NONE to use values directly,
8487a982d89SJeremy L. Thompson                      \ref CEED_EVAL_INTERP to use interpolated values,
8497a982d89SJeremy L. Thompson                      \ref CEED_EVAL_GRAD to use gradients,
8507a982d89SJeremy L. Thompson                      \ref CEED_EVAL_WEIGHT to use quadrature weights.
8517a982d89SJeremy L. Thompson   @param[in] u     Input CeedVector
8527a982d89SJeremy L. Thompson   @param[out] v    Output CeedVector
8537a982d89SJeremy L. Thompson 
8547a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
8557a982d89SJeremy L. Thompson 
8567a982d89SJeremy L. Thompson   @ref User
8577a982d89SJeremy L. Thompson **/
858d1d35e2fSjeremylt int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode,
859d1d35e2fSjeremylt                    CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
8607a982d89SJeremy L. Thompson   int ierr;
8611f9221feSJeremy L Thompson   CeedSize u_length = 0, v_length;
8621f9221feSJeremy L Thompson   CeedInt dim, num_comp, num_nodes, num_qpts;
863e15f9bd0SJeremy L Thompson   ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr);
864d1d35e2fSjeremylt   ierr = CeedBasisGetNumComponents(basis, &num_comp); CeedChk(ierr);
865d1d35e2fSjeremylt   ierr = CeedBasisGetNumNodes(basis, &num_nodes); CeedChk(ierr);
866d1d35e2fSjeremylt   ierr = CeedBasisGetNumQuadraturePoints(basis, &num_qpts); CeedChk(ierr);
867d1d35e2fSjeremylt   ierr = CeedVectorGetLength(v, &v_length); CeedChk(ierr);
8687a982d89SJeremy L. Thompson   if (u) {
869d1d35e2fSjeremylt     ierr = CeedVectorGetLength(u, &u_length); CeedChk(ierr);
8707a982d89SJeremy L. Thompson   }
8717a982d89SJeremy L. Thompson 
872e15f9bd0SJeremy L Thompson   if (!basis->Apply)
873e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
874e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_UNSUPPORTED,
875e15f9bd0SJeremy L Thompson                      "Backend does not support BasisApply");
876e15f9bd0SJeremy L Thompson   // LCOV_EXCL_STOP
877e15f9bd0SJeremy L Thompson 
878e15f9bd0SJeremy L Thompson   // Check compatibility of topological and geometrical dimensions
879d1d35e2fSjeremylt   if ((t_mode == CEED_TRANSPOSE && (v_length%num_nodes != 0 ||
880d1d35e2fSjeremylt                                     u_length%num_qpts != 0)) ||
881d1d35e2fSjeremylt       (t_mode == CEED_NOTRANSPOSE && (u_length%num_nodes != 0 ||
882d1d35e2fSjeremylt                                       v_length%num_qpts != 0)))
8838229195eSjeremylt     // LCOV_EXCL_START
884e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_DIMENSION,
885e15f9bd0SJeremy L Thompson                      "Length of input/output vectors "
8867a982d89SJeremy L. Thompson                      "incompatible with basis dimensions");
8878229195eSjeremylt   // LCOV_EXCL_STOP
8887a982d89SJeremy L. Thompson 
889e15f9bd0SJeremy L Thompson   // Check vector lengths to prevent out of bounds issues
890d1d35e2fSjeremylt   bool bad_dims = false;
891d1d35e2fSjeremylt   switch (eval_mode) {
892e15f9bd0SJeremy L Thompson   case CEED_EVAL_NONE:
893d1d35e2fSjeremylt   case CEED_EVAL_INTERP: bad_dims =
894d1d35e2fSjeremylt       ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts ||
895d1d35e2fSjeremylt                                      v_length < num_elem*num_comp*num_nodes)) ||
896d1d35e2fSjeremylt        (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp ||
897d1d35e2fSjeremylt                                        u_length < num_elem*num_comp*num_nodes)));
898e15f9bd0SJeremy L Thompson     break;
899d1d35e2fSjeremylt   case CEED_EVAL_GRAD: bad_dims =
900d1d35e2fSjeremylt       ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts*dim ||
901d1d35e2fSjeremylt                                      v_length < num_elem*num_comp*num_nodes)) ||
902d1d35e2fSjeremylt        (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp*dim ||
903d1d35e2fSjeremylt                                        u_length < num_elem*num_comp*num_nodes)));
904e15f9bd0SJeremy L Thompson     break;
905e15f9bd0SJeremy L Thompson   case CEED_EVAL_WEIGHT:
906d1d35e2fSjeremylt     bad_dims = v_length < num_elem*num_qpts;
907e15f9bd0SJeremy L Thompson     break;
908e15f9bd0SJeremy L Thompson   // LCOV_EXCL_START
909d1d35e2fSjeremylt   case CEED_EVAL_DIV: bad_dims =
910d1d35e2fSjeremylt       ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts ||
911d1d35e2fSjeremylt                                      v_length < num_elem*num_comp*num_nodes)) ||
912d1d35e2fSjeremylt        (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp ||
913d1d35e2fSjeremylt                                        u_length < num_elem*num_comp*num_nodes)));
914e15f9bd0SJeremy L Thompson     break;
915d1d35e2fSjeremylt   case CEED_EVAL_CURL: bad_dims =
916d1d35e2fSjeremylt       ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts ||
917d1d35e2fSjeremylt                                      v_length < num_elem*num_comp*num_nodes)) ||
918d1d35e2fSjeremylt        (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp ||
919d1d35e2fSjeremylt                                        u_length < num_elem*num_comp*num_nodes)));
920e15f9bd0SJeremy L Thompson     break;
921e15f9bd0SJeremy L Thompson     // LCOV_EXCL_STOP
922e15f9bd0SJeremy L Thompson   }
923d1d35e2fSjeremylt   if (bad_dims)
924e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
925e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_DIMENSION,
926d1d35e2fSjeremylt                      "Input/output vectors too short for basis and evaluation mode");
927e15f9bd0SJeremy L Thompson   // LCOV_EXCL_STOP
928e15f9bd0SJeremy L Thompson 
929d1d35e2fSjeremylt   ierr = basis->Apply(basis, num_elem, t_mode, eval_mode, u, v); CeedChk(ierr);
930e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9317a982d89SJeremy L. Thompson }
9327a982d89SJeremy L. Thompson 
9337a982d89SJeremy L. Thompson /**
934b7c9bbdaSJeremy L Thompson   @brief Get Ceed associated with a CeedBasis
935b7c9bbdaSJeremy L Thompson 
936b7c9bbdaSJeremy L Thompson   @param basis      CeedBasis
937b7c9bbdaSJeremy L Thompson   @param[out] ceed  Variable to store Ceed
938b7c9bbdaSJeremy L Thompson 
939b7c9bbdaSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
940b7c9bbdaSJeremy L Thompson 
941b7c9bbdaSJeremy L Thompson   @ref Advanced
942b7c9bbdaSJeremy L Thompson **/
943b7c9bbdaSJeremy L Thompson int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) {
944b7c9bbdaSJeremy L Thompson   *ceed = basis->ceed;
945b7c9bbdaSJeremy L Thompson   return CEED_ERROR_SUCCESS;
946b7c9bbdaSJeremy L Thompson }
947b7c9bbdaSJeremy L Thompson 
948b7c9bbdaSJeremy L Thompson /**
9499d007619Sjeremylt   @brief Get dimension for given CeedBasis
9509d007619Sjeremylt 
9519d007619Sjeremylt   @param basis     CeedBasis
9529d007619Sjeremylt   @param[out] dim  Variable to store dimension of basis
9539d007619Sjeremylt 
9549d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
9559d007619Sjeremylt 
956b7c9bbdaSJeremy L Thompson   @ref Advanced
9579d007619Sjeremylt **/
9589d007619Sjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) {
9599d007619Sjeremylt   *dim = basis->dim;
960e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9619d007619Sjeremylt }
9629d007619Sjeremylt 
9639d007619Sjeremylt /**
964d99fa3c5SJeremy L Thompson   @brief Get topology for given CeedBasis
965d99fa3c5SJeremy L Thompson 
966d99fa3c5SJeremy L Thompson   @param basis      CeedBasis
967d99fa3c5SJeremy L Thompson   @param[out] topo  Variable to store topology of basis
968d99fa3c5SJeremy L Thompson 
969d99fa3c5SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
970d99fa3c5SJeremy L Thompson 
971b7c9bbdaSJeremy L Thompson   @ref Advanced
972d99fa3c5SJeremy L Thompson **/
973d99fa3c5SJeremy L Thompson int CeedBasisGetTopology(CeedBasis basis, CeedElemTopology *topo) {
974d99fa3c5SJeremy L Thompson   *topo = basis->topo;
975e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
976d99fa3c5SJeremy L Thompson }
977d99fa3c5SJeremy L Thompson 
978d99fa3c5SJeremy L Thompson /**
97950c301a5SRezgar Shakeri   @brief Get number of Q-vector components for given CeedBasis
98050c301a5SRezgar Shakeri 
98150c301a5SRezgar Shakeri   @param basis          CeedBasis
98250c301a5SRezgar Shakeri   @param[out] Q_comp  Variable to store number of Q-vector components of basis
98350c301a5SRezgar Shakeri 
98450c301a5SRezgar Shakeri   @return An error code: 0 - success, otherwise - failure
98550c301a5SRezgar Shakeri 
98650c301a5SRezgar Shakeri   @ref Advanced
98750c301a5SRezgar Shakeri **/
98850c301a5SRezgar Shakeri int CeedBasisGetNumQuadratureComponents(CeedBasis basis, CeedInt *Q_comp) {
98950c301a5SRezgar Shakeri   *Q_comp = basis->Q_comp;
99050c301a5SRezgar Shakeri   return CEED_ERROR_SUCCESS;
99150c301a5SRezgar Shakeri }
99250c301a5SRezgar Shakeri 
99350c301a5SRezgar Shakeri /**
9949d007619Sjeremylt   @brief Get number of components for given CeedBasis
9959d007619Sjeremylt 
9969d007619Sjeremylt   @param basis          CeedBasis
997d1d35e2fSjeremylt   @param[out] num_comp  Variable to store number of components of basis
9989d007619Sjeremylt 
9999d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
10009d007619Sjeremylt 
1001b7c9bbdaSJeremy L Thompson   @ref Advanced
10029d007619Sjeremylt **/
1003d1d35e2fSjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *num_comp) {
1004d1d35e2fSjeremylt   *num_comp = basis->num_comp;
1005e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
10069d007619Sjeremylt }
10079d007619Sjeremylt 
10089d007619Sjeremylt /**
10099d007619Sjeremylt   @brief Get total number of nodes (in dim dimensions) of a CeedBasis
10109d007619Sjeremylt 
10119d007619Sjeremylt   @param basis   CeedBasis
10129d007619Sjeremylt   @param[out] P  Variable to store number of nodes
10139d007619Sjeremylt 
10149d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
10159d007619Sjeremylt 
10169d007619Sjeremylt   @ref Utility
10179d007619Sjeremylt **/
10189d007619Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) {
10199d007619Sjeremylt   *P = basis->P;
1020e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
10219d007619Sjeremylt }
10229d007619Sjeremylt 
10239d007619Sjeremylt /**
10249d007619Sjeremylt   @brief Get total number of nodes (in 1 dimension) of a CeedBasis
10259d007619Sjeremylt 
10269d007619Sjeremylt   @param basis     CeedBasis
1027d1d35e2fSjeremylt   @param[out] P_1d  Variable to store number of nodes
10289d007619Sjeremylt 
10299d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
10309d007619Sjeremylt 
1031b7c9bbdaSJeremy L Thompson   @ref Advanced
10329d007619Sjeremylt **/
1033d1d35e2fSjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P_1d) {
1034d1d35e2fSjeremylt   if (!basis->tensor_basis)
10359d007619Sjeremylt     // LCOV_EXCL_START
1036e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR,
1037d1d35e2fSjeremylt                      "Cannot supply P_1d for non-tensor basis");
10389d007619Sjeremylt   // LCOV_EXCL_STOP
10399d007619Sjeremylt 
1040d1d35e2fSjeremylt   *P_1d = basis->P_1d;
1041e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
10429d007619Sjeremylt }
10439d007619Sjeremylt 
10449d007619Sjeremylt /**
10459d007619Sjeremylt   @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis
10469d007619Sjeremylt 
10479d007619Sjeremylt   @param basis   CeedBasis
10489d007619Sjeremylt   @param[out] Q  Variable to store number of quadrature points
10499d007619Sjeremylt 
10509d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
10519d007619Sjeremylt 
10529d007619Sjeremylt   @ref Utility
10539d007619Sjeremylt **/
10549d007619Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) {
10559d007619Sjeremylt   *Q = basis->Q;
1056e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
10579d007619Sjeremylt }
10589d007619Sjeremylt 
10599d007619Sjeremylt /**
10609d007619Sjeremylt   @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis
10619d007619Sjeremylt 
10629d007619Sjeremylt   @param basis      CeedBasis
1063d1d35e2fSjeremylt   @param[out] Q_1d  Variable to store number of quadrature points
10649d007619Sjeremylt 
10659d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
10669d007619Sjeremylt 
1067b7c9bbdaSJeremy L Thompson   @ref Advanced
10689d007619Sjeremylt **/
1069d1d35e2fSjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q_1d) {
1070d1d35e2fSjeremylt   if (!basis->tensor_basis)
10719d007619Sjeremylt     // LCOV_EXCL_START
1072e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR,
1073d1d35e2fSjeremylt                      "Cannot supply Q_1d for non-tensor basis");
10749d007619Sjeremylt   // LCOV_EXCL_STOP
10759d007619Sjeremylt 
1076d1d35e2fSjeremylt   *Q_1d = basis->Q_1d;
1077e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
10789d007619Sjeremylt }
10799d007619Sjeremylt 
10809d007619Sjeremylt /**
10819d007619Sjeremylt   @brief Get reference coordinates of quadrature points (in dim dimensions)
10829d007619Sjeremylt          of a CeedBasis
10839d007619Sjeremylt 
10849d007619Sjeremylt   @param basis       CeedBasis
1085d1d35e2fSjeremylt   @param[out] q_ref  Variable to store reference coordinates of quadrature points
10869d007619Sjeremylt 
10879d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
10889d007619Sjeremylt 
1089b7c9bbdaSJeremy L Thompson   @ref Advanced
10909d007619Sjeremylt **/
1091d1d35e2fSjeremylt int CeedBasisGetQRef(CeedBasis basis, const CeedScalar **q_ref) {
1092d1d35e2fSjeremylt   *q_ref = basis->q_ref_1d;
1093e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
10949d007619Sjeremylt }
10959d007619Sjeremylt 
10969d007619Sjeremylt /**
10979d007619Sjeremylt   @brief Get quadrature weights of quadrature points (in dim dimensions)
10989d007619Sjeremylt          of a CeedBasis
10999d007619Sjeremylt 
11009d007619Sjeremylt   @param basis          CeedBasis
1101d1d35e2fSjeremylt   @param[out] q_weight  Variable to store quadrature weights
11029d007619Sjeremylt 
11039d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
11049d007619Sjeremylt 
1105b7c9bbdaSJeremy L Thompson   @ref Advanced
11069d007619Sjeremylt **/
1107d1d35e2fSjeremylt int CeedBasisGetQWeights(CeedBasis basis, const CeedScalar **q_weight) {
1108d1d35e2fSjeremylt   *q_weight = basis->q_weight_1d;
1109e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
11109d007619Sjeremylt }
11119d007619Sjeremylt 
11129d007619Sjeremylt /**
11139d007619Sjeremylt   @brief Get interpolation matrix of a CeedBasis
11149d007619Sjeremylt 
11159d007619Sjeremylt   @param basis        CeedBasis
11169d007619Sjeremylt   @param[out] interp  Variable to store interpolation matrix
11179d007619Sjeremylt 
11189d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
11199d007619Sjeremylt 
1120b7c9bbdaSJeremy L Thompson   @ref Advanced
11219d007619Sjeremylt **/
11226c58de82SJeremy L Thompson int CeedBasisGetInterp(CeedBasis basis, const CeedScalar **interp) {
1123d1d35e2fSjeremylt   if (!basis->interp && basis->tensor_basis) {
11249d007619Sjeremylt     // Allocate
11259d007619Sjeremylt     int ierr;
11269d007619Sjeremylt     ierr = CeedMalloc(basis->Q*basis->P, &basis->interp); CeedChk(ierr);
11279d007619Sjeremylt 
11289d007619Sjeremylt     // Initialize
11299d007619Sjeremylt     for (CeedInt i=0; i<basis->Q*basis->P; i++)
11309d007619Sjeremylt       basis->interp[i] = 1.0;
11319d007619Sjeremylt 
11329d007619Sjeremylt     // Calculate
11339d007619Sjeremylt     for (CeedInt d=0; d<basis->dim; d++)
11349d007619Sjeremylt       for (CeedInt qpt=0; qpt<basis->Q; qpt++)
11359d007619Sjeremylt         for (CeedInt node=0; node<basis->P; node++) {
1136d1d35e2fSjeremylt           CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d;
1137d1d35e2fSjeremylt           CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d;
1138d1d35e2fSjeremylt           basis->interp[qpt*(basis->P)+node] *= basis->interp_1d[q*basis->P_1d+p];
11399d007619Sjeremylt         }
11409d007619Sjeremylt   }
11419d007619Sjeremylt   *interp = basis->interp;
1142e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
11439d007619Sjeremylt }
11449d007619Sjeremylt 
11459d007619Sjeremylt /**
11469d007619Sjeremylt   @brief Get 1D interpolation matrix of a tensor product CeedBasis
11479d007619Sjeremylt 
11489d007619Sjeremylt   @param basis           CeedBasis
1149d1d35e2fSjeremylt   @param[out] interp_1d  Variable to store interpolation matrix
11509d007619Sjeremylt 
11519d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
11529d007619Sjeremylt 
11539d007619Sjeremylt   @ref Backend
11549d007619Sjeremylt **/
1155d1d35e2fSjeremylt int CeedBasisGetInterp1D(CeedBasis basis, const CeedScalar **interp_1d) {
1156d1d35e2fSjeremylt   if (!basis->tensor_basis)
11579d007619Sjeremylt     // LCOV_EXCL_START
1158e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR,
1159e15f9bd0SJeremy L Thompson                      "CeedBasis is not a tensor product basis.");
11609d007619Sjeremylt   // LCOV_EXCL_STOP
11619d007619Sjeremylt 
1162d1d35e2fSjeremylt   *interp_1d = basis->interp_1d;
1163e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
11649d007619Sjeremylt }
11659d007619Sjeremylt 
11669d007619Sjeremylt /**
11679d007619Sjeremylt   @brief Get gradient matrix of a CeedBasis
11689d007619Sjeremylt 
11699d007619Sjeremylt   @param basis      CeedBasis
11709d007619Sjeremylt   @param[out] grad  Variable to store gradient matrix
11719d007619Sjeremylt 
11729d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
11739d007619Sjeremylt 
1174b7c9bbdaSJeremy L Thompson   @ref Advanced
11759d007619Sjeremylt **/
11766c58de82SJeremy L Thompson int CeedBasisGetGrad(CeedBasis basis, const CeedScalar **grad) {
1177d1d35e2fSjeremylt   if (!basis->grad && basis->tensor_basis) {
11789d007619Sjeremylt     // Allocate
11799d007619Sjeremylt     int ierr;
11809d007619Sjeremylt     ierr = CeedMalloc(basis->dim*basis->Q*basis->P, &basis->grad);
11819d007619Sjeremylt     CeedChk(ierr);
11829d007619Sjeremylt 
11839d007619Sjeremylt     // Initialize
11849d007619Sjeremylt     for (CeedInt i=0; i<basis->dim*basis->Q*basis->P; i++)
11859d007619Sjeremylt       basis->grad[i] = 1.0;
11869d007619Sjeremylt 
11879d007619Sjeremylt     // Calculate
11889d007619Sjeremylt     for (CeedInt d=0; d<basis->dim; d++)
11899d007619Sjeremylt       for (CeedInt i=0; i<basis->dim; i++)
11909d007619Sjeremylt         for (CeedInt qpt=0; qpt<basis->Q; qpt++)
11919d007619Sjeremylt           for (CeedInt node=0; node<basis->P; node++) {
1192d1d35e2fSjeremylt             CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d;
1193d1d35e2fSjeremylt             CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d;
11949d007619Sjeremylt             if (i == d)
11959d007619Sjeremylt               basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *=
1196d1d35e2fSjeremylt                 basis->grad_1d[q*basis->P_1d+p];
11979d007619Sjeremylt             else
11989d007619Sjeremylt               basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *=
1199d1d35e2fSjeremylt                 basis->interp_1d[q*basis->P_1d+p];
12009d007619Sjeremylt           }
12019d007619Sjeremylt   }
12029d007619Sjeremylt   *grad = basis->grad;
1203e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
12049d007619Sjeremylt }
12059d007619Sjeremylt 
12069d007619Sjeremylt /**
12079d007619Sjeremylt   @brief Get 1D gradient matrix of a tensor product CeedBasis
12089d007619Sjeremylt 
12099d007619Sjeremylt   @param basis         CeedBasis
1210d1d35e2fSjeremylt   @param[out] grad_1d  Variable to store gradient matrix
12119d007619Sjeremylt 
12129d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
12139d007619Sjeremylt 
1214b7c9bbdaSJeremy L Thompson   @ref Advanced
12159d007619Sjeremylt **/
1216d1d35e2fSjeremylt int CeedBasisGetGrad1D(CeedBasis basis, const CeedScalar **grad_1d) {
1217d1d35e2fSjeremylt   if (!basis->tensor_basis)
12189d007619Sjeremylt     // LCOV_EXCL_START
1219e15f9bd0SJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR,
1220e15f9bd0SJeremy L Thompson                      "CeedBasis is not a tensor product basis.");
12219d007619Sjeremylt   // LCOV_EXCL_STOP
12229d007619Sjeremylt 
1223d1d35e2fSjeremylt   *grad_1d = basis->grad_1d;
1224e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
12259d007619Sjeremylt }
12269d007619Sjeremylt 
12279d007619Sjeremylt /**
122850c301a5SRezgar Shakeri   @brief Get divergence matrix of a CeedBasis
122950c301a5SRezgar Shakeri 
123050c301a5SRezgar Shakeri   @param basis     CeedBasis
123150c301a5SRezgar Shakeri   @param[out] div  Variable to store divergence matrix
123250c301a5SRezgar Shakeri 
123350c301a5SRezgar Shakeri   @return An error code: 0 - success, otherwise - failure
123450c301a5SRezgar Shakeri 
123550c301a5SRezgar Shakeri   @ref Advanced
123650c301a5SRezgar Shakeri **/
123750c301a5SRezgar Shakeri int CeedBasisGetDiv(CeedBasis basis, const CeedScalar **div) {
123850c301a5SRezgar Shakeri   if (!basis->div)
123950c301a5SRezgar Shakeri     // LCOV_EXCL_START
124050c301a5SRezgar Shakeri     return CeedError(basis->ceed, CEED_ERROR_MINOR,
124150c301a5SRezgar Shakeri                      "CeedBasis does not have divergence matrix.");
124250c301a5SRezgar Shakeri   // LCOV_EXCL_STOP
124350c301a5SRezgar Shakeri 
124450c301a5SRezgar Shakeri   *div = basis->div;
124550c301a5SRezgar Shakeri   return CEED_ERROR_SUCCESS;
124650c301a5SRezgar Shakeri }
124750c301a5SRezgar Shakeri 
124850c301a5SRezgar Shakeri /**
12497a982d89SJeremy L. Thompson   @brief Destroy a CeedBasis
12507a982d89SJeremy L. Thompson 
12517a982d89SJeremy L. Thompson   @param basis CeedBasis to destroy
12527a982d89SJeremy L. Thompson 
12537a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
12547a982d89SJeremy L. Thompson 
12557a982d89SJeremy L. Thompson   @ref User
12567a982d89SJeremy L. Thompson **/
12577a982d89SJeremy L. Thompson int CeedBasisDestroy(CeedBasis *basis) {
12587a982d89SJeremy L. Thompson   int ierr;
12597a982d89SJeremy L. Thompson 
1260d1d35e2fSjeremylt   if (!*basis || --(*basis)->ref_count > 0) return CEED_ERROR_SUCCESS;
12617a982d89SJeremy L. Thompson   if ((*basis)->Destroy) {
12627a982d89SJeremy L. Thompson     ierr = (*basis)->Destroy(*basis); CeedChk(ierr);
12637a982d89SJeremy L. Thompson   }
126434359f16Sjeremylt   if ((*basis)->contract) {
126534359f16Sjeremylt     ierr = CeedTensorContractDestroy(&(*basis)->contract); CeedChk(ierr);
126634359f16Sjeremylt   }
12677a982d89SJeremy L. Thompson   ierr = CeedFree(&(*basis)->interp); CeedChk(ierr);
1268d1d35e2fSjeremylt   ierr = CeedFree(&(*basis)->interp_1d); CeedChk(ierr);
12697a982d89SJeremy L. Thompson   ierr = CeedFree(&(*basis)->grad); CeedChk(ierr);
127050c301a5SRezgar Shakeri   ierr = CeedFree(&(*basis)->div); CeedChk(ierr);
1271d1d35e2fSjeremylt   ierr = CeedFree(&(*basis)->grad_1d); CeedChk(ierr);
1272d1d35e2fSjeremylt   ierr = CeedFree(&(*basis)->q_ref_1d); CeedChk(ierr);
1273d1d35e2fSjeremylt   ierr = CeedFree(&(*basis)->q_weight_1d); CeedChk(ierr);
12747a982d89SJeremy L. Thompson   ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr);
12757a982d89SJeremy L. Thompson   ierr = CeedFree(basis); CeedChk(ierr);
1276e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
12777a982d89SJeremy L. Thompson }
12787a982d89SJeremy L. Thompson 
12797a982d89SJeremy L. Thompson /**
1280b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre quadrature
1281b11c1e72Sjeremylt 
1282b11c1e72Sjeremylt   @param Q                 Number of quadrature points (integrates polynomials of
1283b11c1e72Sjeremylt                              degree 2*Q-1 exactly)
1284d1d35e2fSjeremylt   @param[out] q_ref_1d     Array of length Q to hold the abscissa on [-1, 1]
1285d1d35e2fSjeremylt   @param[out] q_weight_1d  Array of length Q to hold the weights
1286b11c1e72Sjeremylt 
1287b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1288dfdf5a53Sjeremylt 
1289dfdf5a53Sjeremylt   @ref Utility
1290b11c1e72Sjeremylt **/
1291d1d35e2fSjeremylt int CeedGaussQuadrature(CeedInt Q, CeedScalar *q_ref_1d,
1292d1d35e2fSjeremylt                         CeedScalar *q_weight_1d) {
1293d7b241e6Sjeremylt   // Allocate
1294d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0*atan(1.0);
1295d1d35e2fSjeremylt   // Build q_ref_1d, q_weight_1d
1296d7b241e6Sjeremylt   for (int i = 0; i <= Q/2; i++) {
1297d7b241e6Sjeremylt     // Guess
1298d7b241e6Sjeremylt     xi = cos(PI*(CeedScalar)(2*i+1)/((CeedScalar)(2*Q)));
1299d7b241e6Sjeremylt     // Pn(xi)
1300d7b241e6Sjeremylt     P0 = 1.0;
1301d7b241e6Sjeremylt     P1 = xi;
1302d7b241e6Sjeremylt     P2 = 0.0;
1303d7b241e6Sjeremylt     for (int j = 2; j <= Q; j++) {
1304d7b241e6Sjeremylt       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1305d7b241e6Sjeremylt       P0 = P1;
1306d7b241e6Sjeremylt       P1 = P2;
1307d7b241e6Sjeremylt     }
1308d7b241e6Sjeremylt     // First Newton Step
1309d7b241e6Sjeremylt     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1310d7b241e6Sjeremylt     xi = xi-P2/dP2;
1311d7b241e6Sjeremylt     // Newton to convergence
13120e4d4210Sjeremylt     for (int k=0; k<100 && fabs(P2)>10*CEED_EPSILON; k++) {
1313d7b241e6Sjeremylt       P0 = 1.0;
1314d7b241e6Sjeremylt       P1 = xi;
1315d7b241e6Sjeremylt       for (int j = 2; j <= Q; j++) {
1316d7b241e6Sjeremylt         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1317d7b241e6Sjeremylt         P0 = P1;
1318d7b241e6Sjeremylt         P1 = P2;
1319d7b241e6Sjeremylt       }
1320d7b241e6Sjeremylt       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1321d7b241e6Sjeremylt       xi = xi-P2/dP2;
1322d7b241e6Sjeremylt     }
1323d7b241e6Sjeremylt     // Save xi, wi
1324d7b241e6Sjeremylt     wi = 2.0/((1.0-xi*xi)*dP2*dP2);
1325d1d35e2fSjeremylt     q_weight_1d[i] = wi;
1326d1d35e2fSjeremylt     q_weight_1d[Q-1-i] = wi;
1327d1d35e2fSjeremylt     q_ref_1d[i] = -xi;
1328d1d35e2fSjeremylt     q_ref_1d[Q-1-i]= xi;
1329d7b241e6Sjeremylt   }
1330e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1331d7b241e6Sjeremylt }
1332d7b241e6Sjeremylt 
1333b11c1e72Sjeremylt /**
1334b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre-Lobatto quadrature
1335b11c1e72Sjeremylt 
1336b11c1e72Sjeremylt   @param Q                 Number of quadrature points (integrates polynomials of
1337b11c1e72Sjeremylt                              degree 2*Q-3 exactly)
1338d1d35e2fSjeremylt   @param[out] q_ref_1d     Array of length Q to hold the abscissa on [-1, 1]
1339d1d35e2fSjeremylt   @param[out] q_weight_1d  Array of length Q to hold the weights
1340b11c1e72Sjeremylt 
1341b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1342dfdf5a53Sjeremylt 
1343dfdf5a53Sjeremylt   @ref Utility
1344b11c1e72Sjeremylt **/
1345d1d35e2fSjeremylt int CeedLobattoQuadrature(CeedInt Q, CeedScalar *q_ref_1d,
1346d1d35e2fSjeremylt                           CeedScalar *q_weight_1d) {
1347d7b241e6Sjeremylt   // Allocate
1348d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0*atan(1.0);
1349d1d35e2fSjeremylt   // Build q_ref_1d, q_weight_1d
1350d7b241e6Sjeremylt   // Set endpoints
135130a100c3SJed Brown   if (Q < 2)
1352b0d62198Sjeremylt     // LCOV_EXCL_START
1353e15f9bd0SJeremy L Thompson     return CeedError(NULL, CEED_ERROR_DIMENSION,
13547ed52d01Sjeremylt                      "Cannot create Lobatto quadrature with Q=%d < 2 points", Q);
1355b0d62198Sjeremylt   // LCOV_EXCL_STOP
1356d7b241e6Sjeremylt   wi = 2.0/((CeedScalar)(Q*(Q-1)));
1357d1d35e2fSjeremylt   if (q_weight_1d) {
1358d1d35e2fSjeremylt     q_weight_1d[0] = wi;
1359d1d35e2fSjeremylt     q_weight_1d[Q-1] = wi;
1360d7b241e6Sjeremylt   }
1361d1d35e2fSjeremylt   q_ref_1d[0] = -1.0;
1362d1d35e2fSjeremylt   q_ref_1d[Q-1] = 1.0;
1363d7b241e6Sjeremylt   // Interior
1364d7b241e6Sjeremylt   for (int i = 1; i <= (Q-1)/2; i++) {
1365d7b241e6Sjeremylt     // Guess
1366d7b241e6Sjeremylt     xi = cos(PI*(CeedScalar)(i)/(CeedScalar)(Q-1));
1367d7b241e6Sjeremylt     // Pn(xi)
1368d7b241e6Sjeremylt     P0 = 1.0;
1369d7b241e6Sjeremylt     P1 = xi;
1370d7b241e6Sjeremylt     P2 = 0.0;
1371d7b241e6Sjeremylt     for (int j = 2; j < Q; j++) {
1372d7b241e6Sjeremylt       P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1373d7b241e6Sjeremylt       P0 = P1;
1374d7b241e6Sjeremylt       P1 = P2;
1375d7b241e6Sjeremylt     }
1376d7b241e6Sjeremylt     // First Newton step
1377d7b241e6Sjeremylt     dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1378d7b241e6Sjeremylt     d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
1379d7b241e6Sjeremylt     xi = xi-dP2/d2P2;
1380d7b241e6Sjeremylt     // Newton to convergence
13810e4d4210Sjeremylt     for (int k=0; k<100 && fabs(dP2)>10*CEED_EPSILON; k++) {
1382d7b241e6Sjeremylt       P0 = 1.0;
1383d7b241e6Sjeremylt       P1 = xi;
1384d7b241e6Sjeremylt       for (int j = 2; j < Q; j++) {
1385d7b241e6Sjeremylt         P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j));
1386d7b241e6Sjeremylt         P0 = P1;
1387d7b241e6Sjeremylt         P1 = P2;
1388d7b241e6Sjeremylt       }
1389d7b241e6Sjeremylt       dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0);
1390d7b241e6Sjeremylt       d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi);
1391d7b241e6Sjeremylt       xi = xi-dP2/d2P2;
1392d7b241e6Sjeremylt     }
1393d7b241e6Sjeremylt     // Save xi, wi
1394d7b241e6Sjeremylt     wi = 2.0/(((CeedScalar)(Q*(Q-1)))*P2*P2);
1395d1d35e2fSjeremylt     if (q_weight_1d) {
1396d1d35e2fSjeremylt       q_weight_1d[i] = wi;
1397d1d35e2fSjeremylt       q_weight_1d[Q-1-i] = wi;
1398d7b241e6Sjeremylt     }
1399d1d35e2fSjeremylt     q_ref_1d[i] = -xi;
1400d1d35e2fSjeremylt     q_ref_1d[Q-1-i]= xi;
1401d7b241e6Sjeremylt   }
1402e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1403d7b241e6Sjeremylt }
1404d7b241e6Sjeremylt 
1405dfdf5a53Sjeremylt /**
140695bb1877Svaleriabarra   @brief Return QR Factorization of a matrix
1407b11c1e72Sjeremylt 
140877645d7bSjeremylt   @param ceed         A Ceed context for error handling
140952bfb9bbSJeremy L Thompson   @param[in,out] mat  Row-major matrix to be factorized in place
141052bfb9bbSJeremy L Thompson   @param[in,out] tau  Vector of length m of scaling factors
1411b11c1e72Sjeremylt   @param m            Number of rows
1412b11c1e72Sjeremylt   @param n            Number of columns
1413b11c1e72Sjeremylt 
1414b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1415dfdf5a53Sjeremylt 
1416dfdf5a53Sjeremylt   @ref Utility
1417b11c1e72Sjeremylt **/
1418a7bd39daSjeremylt int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau,
1419d7b241e6Sjeremylt                         CeedInt m, CeedInt n) {
1420d7b241e6Sjeremylt   CeedScalar v[m];
1421d7b241e6Sjeremylt 
1422a7bd39daSjeremylt   // Check m >= n
1423a7bd39daSjeremylt   if (n > m)
1424c042f62fSJeremy L Thompson     // LCOV_EXCL_START
1425e15f9bd0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_UNSUPPORTED,
1426e15f9bd0SJeremy L Thompson                      "Cannot compute QR factorization with n > m");
1427c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
1428a7bd39daSjeremylt 
142952bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++) {
1430bde37e8cSJed Brown     if (i >= m-1) { // last row of matrix, no reflection needed
1431bde37e8cSJed Brown       tau[i] = 0.;
1432bde37e8cSJed Brown       break;
1433bde37e8cSJed Brown     }
1434d7b241e6Sjeremylt     // Calculate Householder vector, magnitude
1435d7b241e6Sjeremylt     CeedScalar sigma = 0.0;
1436d7b241e6Sjeremylt     v[i] = mat[i+n*i];
143752bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<m; j++) {
1438d7b241e6Sjeremylt       v[j] = mat[i+n*j];
1439d7b241e6Sjeremylt       sigma += v[j] * v[j];
1440d7b241e6Sjeremylt     }
1441d7b241e6Sjeremylt     CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m]
1442673160d7Sjeremylt     CeedScalar R_ii = -copysign(norm, v[i]);
1443673160d7Sjeremylt     v[i] -= R_ii;
1444d7b241e6Sjeremylt     // norm of v[i:m] after modification above and scaling below
1445d7b241e6Sjeremylt     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
1446d7b241e6Sjeremylt     //   tau = 2 / (norm*norm)
1447d7b241e6Sjeremylt     tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma);
14481d102b48SJeremy L Thompson     for (CeedInt j=i+1; j<m; j++)
14491d102b48SJeremy L Thompson       v[j] /= v[i];
1450d7b241e6Sjeremylt 
1451d7b241e6Sjeremylt     // Apply Householder reflector to lower right panel
1452d7b241e6Sjeremylt     CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1);
1453d7b241e6Sjeremylt     // Save v
1454673160d7Sjeremylt     mat[i+n*i] = R_ii;
14551d102b48SJeremy L Thompson     for (CeedInt j=i+1; j<m; j++)
1456d7b241e6Sjeremylt       mat[i+n*j] = v[j];
1457d7b241e6Sjeremylt   }
1458e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1459d7b241e6Sjeremylt }
1460d7b241e6Sjeremylt 
1461b11c1e72Sjeremylt /**
146252bfb9bbSJeremy L Thompson   @brief Return symmetric Schur decomposition of the symmetric matrix mat via
146352bfb9bbSJeremy L Thompson            symmetric QR factorization
146452bfb9bbSJeremy L Thompson 
146577645d7bSjeremylt   @param ceed         A Ceed context for error handling
146652bfb9bbSJeremy L Thompson   @param[in,out] mat  Row-major matrix to be factorized in place
1467460bf743SValeria Barra   @param[out] lambda  Vector of length n of eigenvalues
146852bfb9bbSJeremy L Thompson   @param n            Number of rows/columns
146952bfb9bbSJeremy L Thompson 
147052bfb9bbSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
147152bfb9bbSJeremy L Thompson 
147252bfb9bbSJeremy L Thompson   @ref Utility
147352bfb9bbSJeremy L Thompson **/
147403d18186Sjeremylt CeedPragmaOptimizeOff
147552bfb9bbSJeremy L Thompson int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat,
147652bfb9bbSJeremy L Thompson                                     CeedScalar *lambda, CeedInt n) {
147752bfb9bbSJeremy L Thompson   // Check bounds for clang-tidy
147852bfb9bbSJeremy L Thompson   if (n<2)
1479c042f62fSJeremy L Thompson     // LCOV_EXCL_START
1480e15f9bd0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_UNSUPPORTED,
1481c042f62fSJeremy L Thompson                      "Cannot compute symmetric Schur decomposition of scalars");
1482c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
148352bfb9bbSJeremy L Thompson 
1484673160d7Sjeremylt   CeedScalar v[n-1], tau[n-1], mat_T[n*n];
148552bfb9bbSJeremy L Thompson 
1486673160d7Sjeremylt   // Copy mat to mat_T and set mat to I
1487673160d7Sjeremylt   memcpy(mat_T, mat, n*n*sizeof(mat[0]));
148852bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
148952bfb9bbSJeremy L Thompson     for (CeedInt j=0; j<n; j++)
149052bfb9bbSJeremy L Thompson       mat[j+n*i] = (i==j) ? 1 : 0;
149152bfb9bbSJeremy L Thompson 
149252bfb9bbSJeremy L Thompson   // Reduce to tridiagonal
149352bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n-1; i++) {
149452bfb9bbSJeremy L Thompson     // Calculate Householder vector, magnitude
149552bfb9bbSJeremy L Thompson     CeedScalar sigma = 0.0;
1496673160d7Sjeremylt     v[i] = mat_T[i+n*(i+1)];
149752bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) {
1498673160d7Sjeremylt       v[j] = mat_T[i+n*(j+1)];
149952bfb9bbSJeremy L Thompson       sigma += v[j] * v[j];
150052bfb9bbSJeremy L Thompson     }
150152bfb9bbSJeremy L Thompson     CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:n-1]
1502673160d7Sjeremylt     CeedScalar R_ii = -copysign(norm, v[i]);
1503673160d7Sjeremylt     v[i] -= R_ii;
150452bfb9bbSJeremy L Thompson     // norm of v[i:m] after modification above and scaling below
150552bfb9bbSJeremy L Thompson     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
150652bfb9bbSJeremy L Thompson     //   tau = 2 / (norm*norm)
150780a9ef05SNatalie Beams     tau[i] = i == n - 2 ? 2 : 2 * v[i]*v[i] / (v[i]*v[i] + sigma);
1508fb551037Sjeremylt     for (CeedInt j=i+1; j<n-1; j++)
1509fb551037Sjeremylt       v[j] /= v[i];
151052bfb9bbSJeremy L Thompson 
151152bfb9bbSJeremy L Thompson     // Update sub and super diagonal
151252bfb9bbSJeremy L Thompson     for (CeedInt j=i+2; j<n; j++) {
1513673160d7Sjeremylt       mat_T[i+n*j] = 0; mat_T[j+n*i] = 0;
151452bfb9bbSJeremy L Thompson     }
151552bfb9bbSJeremy L Thompson     // Apply symmetric Householder reflector to lower right panel
1516673160d7Sjeremylt     CeedHouseholderReflect(&mat_T[(i+1)+n*(i+1)], &v[i], tau[i],
151752bfb9bbSJeremy L Thompson                            n-(i+1), n-(i+1), n, 1);
1518673160d7Sjeremylt     CeedHouseholderReflect(&mat_T[(i+1)+n*(i+1)], &v[i], tau[i],
151952bfb9bbSJeremy L Thompson                            n-(i+1), n-(i+1), 1, n);
1520673160d7Sjeremylt 
152152bfb9bbSJeremy L Thompson     // Save v
1522673160d7Sjeremylt     mat_T[i+n*(i+1)] = R_ii;
1523673160d7Sjeremylt     mat_T[(i+1)+n*i] = R_ii;
152452bfb9bbSJeremy L Thompson     for (CeedInt j=i+1; j<n-1; j++) {
1525673160d7Sjeremylt       mat_T[i+n*(j+1)] = v[j];
152652bfb9bbSJeremy L Thompson     }
152752bfb9bbSJeremy L Thompson   }
152852bfb9bbSJeremy L Thompson   // Backwards accumulation of Q
152952bfb9bbSJeremy L Thompson   for (CeedInt i=n-2; i>=0; i--) {
153085cf89eaSjeremylt     if (tau[i] > 0.0) {
153152bfb9bbSJeremy L Thompson       v[i] = 1;
153252bfb9bbSJeremy L Thompson       for (CeedInt j=i+1; j<n-1; j++) {
1533673160d7Sjeremylt         v[j] = mat_T[i+n*(j+1)];
1534673160d7Sjeremylt         mat_T[i+n*(j+1)] = 0;
153552bfb9bbSJeremy L Thompson       }
153652bfb9bbSJeremy L Thompson       CeedHouseholderReflect(&mat[(i+1)+n*(i+1)], &v[i], tau[i],
153752bfb9bbSJeremy L Thompson                              n-(i+1), n-(i+1), n, 1);
153852bfb9bbSJeremy L Thompson     }
153985cf89eaSjeremylt   }
154052bfb9bbSJeremy L Thompson 
154152bfb9bbSJeremy L Thompson   // Reduce sub and super diagonal
1542673160d7Sjeremylt   CeedInt p = 0, q = 0, itr = 0, max_itr = n*n*n*n;
1543673160d7Sjeremylt   CeedScalar tol = CEED_EPSILON;
154452bfb9bbSJeremy L Thompson 
1545673160d7Sjeremylt   while (itr < max_itr) {
154652bfb9bbSJeremy L Thompson     // Update p, q, size of reduced portions of diagonal
154752bfb9bbSJeremy L Thompson     p = 0; q = 0;
154852bfb9bbSJeremy L Thompson     for (CeedInt i=n-2; i>=0; i--) {
1549673160d7Sjeremylt       if (fabs(mat_T[i+n*(i+1)]) < tol)
155052bfb9bbSJeremy L Thompson         q += 1;
155152bfb9bbSJeremy L Thompson       else
155252bfb9bbSJeremy L Thompson         break;
155352bfb9bbSJeremy L Thompson     }
1554673160d7Sjeremylt     for (CeedInt i=0; i<n-q-1; i++) {
1555673160d7Sjeremylt       if (fabs(mat_T[i+n*(i+1)]) < tol)
155652bfb9bbSJeremy L Thompson         p += 1;
155752bfb9bbSJeremy L Thompson       else
155852bfb9bbSJeremy L Thompson         break;
155952bfb9bbSJeremy L Thompson     }
156052bfb9bbSJeremy L Thompson     if (q == n-1) break; // Finished reducing
156152bfb9bbSJeremy L Thompson 
156252bfb9bbSJeremy L Thompson     // Reduce tridiagonal portion
1563673160d7Sjeremylt     CeedScalar t_nn = mat_T[(n-1-q)+n*(n-1-q)],
1564673160d7Sjeremylt                t_nnm1 = mat_T[(n-2-q)+n*(n-1-q)];
1565673160d7Sjeremylt     CeedScalar d = (mat_T[(n-2-q)+n*(n-2-q)] - t_nn)/2;
1566673160d7Sjeremylt     CeedScalar mu = t_nn - t_nnm1*t_nnm1 /
1567673160d7Sjeremylt                     (d + copysign(sqrt(d*d + t_nnm1*t_nnm1), d));
1568673160d7Sjeremylt     CeedScalar x = mat_T[p+n*p] - mu;
1569673160d7Sjeremylt     CeedScalar z = mat_T[p+n*(p+1)];
1570673160d7Sjeremylt     for (CeedInt k=p; k<n-q-1; k++) {
157152bfb9bbSJeremy L Thompson       // Compute Givens rotation
157252bfb9bbSJeremy L Thompson       CeedScalar c = 1, s = 0;
157352bfb9bbSJeremy L Thompson       if (fabs(z) > tol) {
157452bfb9bbSJeremy L Thompson         if (fabs(z) > fabs(x)) {
157552bfb9bbSJeremy L Thompson           CeedScalar tau = -x/z;
157652bfb9bbSJeremy L Thompson           s = 1/sqrt(1+tau*tau), c = s*tau;
157752bfb9bbSJeremy L Thompson         } else {
157852bfb9bbSJeremy L Thompson           CeedScalar tau = -z/x;
157952bfb9bbSJeremy L Thompson           c = 1/sqrt(1+tau*tau), s = c*tau;
158052bfb9bbSJeremy L Thompson         }
158152bfb9bbSJeremy L Thompson       }
158252bfb9bbSJeremy L Thompson 
158352bfb9bbSJeremy L Thompson       // Apply Givens rotation to T
1584673160d7Sjeremylt       CeedGivensRotation(mat_T, c, s, CEED_NOTRANSPOSE, k, k+1, n, n);
1585673160d7Sjeremylt       CeedGivensRotation(mat_T, c, s, CEED_TRANSPOSE, k, k+1, n, n);
158652bfb9bbSJeremy L Thompson 
158752bfb9bbSJeremy L Thompson       // Apply Givens rotation to Q
158852bfb9bbSJeremy L Thompson       CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k+1, n, n);
158952bfb9bbSJeremy L Thompson 
159052bfb9bbSJeremy L Thompson       // Update x, z
159152bfb9bbSJeremy L Thompson       if (k < n-q-2) {
1592673160d7Sjeremylt         x = mat_T[k+n*(k+1)];
1593673160d7Sjeremylt         z = mat_T[k+n*(k+2)];
159452bfb9bbSJeremy L Thompson       }
159552bfb9bbSJeremy L Thompson     }
159652bfb9bbSJeremy L Thompson     itr++;
159752bfb9bbSJeremy L Thompson   }
1598673160d7Sjeremylt 
159952bfb9bbSJeremy L Thompson   // Save eigenvalues
160052bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
1601673160d7Sjeremylt     lambda[i] = mat_T[i+n*i];
160252bfb9bbSJeremy L Thompson 
160352bfb9bbSJeremy L Thompson   // Check convergence
1604673160d7Sjeremylt   if (itr == max_itr && q < n-1)
1605c042f62fSJeremy L Thompson     // LCOV_EXCL_START
1606e15f9bd0SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_MINOR,
1607e15f9bd0SJeremy L Thompson                      "Symmetric QR failed to converge");
1608c042f62fSJeremy L Thompson   // LCOV_EXCL_STOP
1609e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
161052bfb9bbSJeremy L Thompson }
161103d18186Sjeremylt CeedPragmaOptimizeOn
161252bfb9bbSJeremy L Thompson 
161352bfb9bbSJeremy L Thompson /**
161452bfb9bbSJeremy L Thompson   @brief Return Simultaneous Diagonalization of two matrices. This solves the
161552bfb9bbSJeremy L Thompson            generalized eigenvalue problem A x = lambda B x, where A and B
161652bfb9bbSJeremy L Thompson            are symmetric and B is positive definite. We generate the matrix X
161752bfb9bbSJeremy L Thompson            and vector Lambda such that X^T A X = Lambda and X^T B X = I. This
161852bfb9bbSJeremy L Thompson            is equivalent to the LAPACK routine 'sygv' with TYPE = 1.
161952bfb9bbSJeremy L Thompson 
162077645d7bSjeremylt   @param ceed         A Ceed context for error handling
1621d1d35e2fSjeremylt   @param[in] mat_A    Row-major matrix to be factorized with eigenvalues
1622d1d35e2fSjeremylt   @param[in] mat_B    Row-major matrix to be factorized to identity
1623d3331725Sjeremylt   @param[out] mat_X   Row-major orthogonal matrix
1624460bf743SValeria Barra   @param[out] lambda  Vector of length n of generalized eigenvalues
162552bfb9bbSJeremy L Thompson   @param n            Number of rows/columns
162652bfb9bbSJeremy L Thompson 
162752bfb9bbSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
162852bfb9bbSJeremy L Thompson 
162952bfb9bbSJeremy L Thompson   @ref Utility
163052bfb9bbSJeremy L Thompson **/
163103d18186Sjeremylt CeedPragmaOptimizeOff
1632d1d35e2fSjeremylt int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *mat_A,
1633d3331725Sjeremylt                                     CeedScalar *mat_B, CeedScalar *mat_X,
163452bfb9bbSJeremy L Thompson                                     CeedScalar *lambda, CeedInt n) {
163552bfb9bbSJeremy L Thompson   int ierr;
1636d3331725Sjeremylt   CeedScalar *mat_C, *mat_G, *vec_D;
163778464608Sjeremylt   ierr = CeedCalloc(n*n, &mat_C); CeedChk(ierr);
163878464608Sjeremylt   ierr = CeedCalloc(n*n, &mat_G); CeedChk(ierr);
1639d3331725Sjeremylt   ierr = CeedCalloc(n, &vec_D); CeedChk(ierr);
164052bfb9bbSJeremy L Thompson 
164152bfb9bbSJeremy L Thompson   // Compute B = G D G^T
164278464608Sjeremylt   memcpy(mat_G, mat_B, n*n*sizeof(mat_B[0]));
1643d3331725Sjeremylt   ierr = CeedSymmetricSchurDecomposition(ceed, mat_G, vec_D, n); CeedChk(ierr);
164452bfb9bbSJeremy L Thompson 
164585cf89eaSjeremylt   // Sort eigenvalues
164685cf89eaSjeremylt   for (CeedInt i=n-1; i>=0; i--)
164785cf89eaSjeremylt     for (CeedInt j=0; j<i; j++) {
164885cf89eaSjeremylt       if (fabs(vec_D[j]) > fabs(vec_D[j+1])) {
164985cf89eaSjeremylt         CeedScalar temp;
165085cf89eaSjeremylt         temp = vec_D[j]; vec_D[j] = vec_D[j+1]; vec_D[j+1] = temp;
165185cf89eaSjeremylt         for (CeedInt k=0; k<n; k++) {
165285cf89eaSjeremylt           temp = mat_G[k*n+j]; mat_G[k*n+j] = mat_G[k*n+j+1]; mat_G[k*n+j+1] = temp;
165385cf89eaSjeremylt         }
165485cf89eaSjeremylt       }
165585cf89eaSjeremylt     }
165685cf89eaSjeremylt 
1657fb551037Sjeremylt   // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T
1658fb551037Sjeremylt   //           = D^-1/2 G^T A G D^-1/2
1659d3331725Sjeremylt   // -- D = D^-1/2
166052bfb9bbSJeremy L Thompson   for (CeedInt i=0; i<n; i++)
1661d3331725Sjeremylt     vec_D[i] = 1./sqrt(vec_D[i]);
1662d3331725Sjeremylt   // -- G = G D^-1/2
1663d3331725Sjeremylt   // -- C = D^-1/2 G^T
1664d3331725Sjeremylt   for (CeedInt i=0; i<n; i++)
1665d3331725Sjeremylt     for (CeedInt j=0; j<n; j++) {
1666673160d7Sjeremylt       mat_G[i*n+j] *= vec_D[j];
1667673160d7Sjeremylt       mat_C[j*n+i]  = mat_G[i*n+j];
1668d3331725Sjeremylt     }
1669673160d7Sjeremylt   // -- X = (D^-1/2 G^T) A
1670d1d35e2fSjeremylt   ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)mat_C,
1671d3331725Sjeremylt                             (const CeedScalar *)mat_A, mat_X, n, n, n);
16729289e5bfSjeremylt   CeedChk(ierr);
1673673160d7Sjeremylt   // -- C = (D^-1/2 G^T A) (G D^-1/2)
1674d3331725Sjeremylt   ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)mat_X,
167578464608Sjeremylt                             (const CeedScalar *)mat_G, mat_C, n, n, n);
16769289e5bfSjeremylt   CeedChk(ierr);
167752bfb9bbSJeremy L Thompson 
167852bfb9bbSJeremy L Thompson   // Compute Q^T C Q = lambda
1679d1d35e2fSjeremylt   ierr = CeedSymmetricSchurDecomposition(ceed, mat_C, lambda, n); CeedChk(ierr);
168052bfb9bbSJeremy L Thompson 
168185cf89eaSjeremylt   // Sort eigenvalues
168285cf89eaSjeremylt   for (CeedInt i=n-1; i>=0; i--)
168385cf89eaSjeremylt     for (CeedInt j=0; j<i; j++) {
168485cf89eaSjeremylt       if (fabs(lambda[j]) > fabs(lambda[j+1])) {
168585cf89eaSjeremylt         CeedScalar temp;
168685cf89eaSjeremylt         temp = lambda[j]; lambda[j] = lambda[j+1]; lambda[j+1] = temp;
168785cf89eaSjeremylt         for (CeedInt k=0; k<n; k++) {
168885cf89eaSjeremylt           temp = mat_C[k*n+j]; mat_C[k*n+j] = mat_C[k*n+j+1]; mat_C[k*n+j+1] = temp;
168985cf89eaSjeremylt         }
169085cf89eaSjeremylt       }
169185cf89eaSjeremylt     }
169285cf89eaSjeremylt 
1693d3331725Sjeremylt   // Set X = (G D^1/2)^-T Q
1694fb551037Sjeremylt   //       = G D^-1/2 Q
169578464608Sjeremylt   ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)mat_G,
1696d3331725Sjeremylt                             (const CeedScalar *)mat_C, mat_X, n, n, n);
16979289e5bfSjeremylt   CeedChk(ierr);
169878464608Sjeremylt 
169978464608Sjeremylt   // Cleanup
170078464608Sjeremylt   ierr = CeedFree(&mat_C); CeedChk(ierr);
170178464608Sjeremylt   ierr = CeedFree(&mat_G); CeedChk(ierr);
1702d3331725Sjeremylt   ierr = CeedFree(&vec_D); CeedChk(ierr);
1703e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
170452bfb9bbSJeremy L Thompson }
170503d18186Sjeremylt CeedPragmaOptimizeOn
170652bfb9bbSJeremy L Thompson 
1707d7b241e6Sjeremylt /// @}
1708