xref: /libCEED/interface/ceed-basis.c (revision 19a04db8f2245ce5b2896f00805cba58aafa343d)
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, 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 
83d576824SJeremy L Thompson #include <ceed-impl.h>
949aac155SJeremy L Thompson #include <ceed.h>
102b730f8bSJeremy L Thompson #include <ceed/backend.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
20356036faSJeremy L Thompson static struct CeedBasis_private ceed_basis_none;
21d7b241e6Sjeremylt /// @endcond
22d7b241e6Sjeremylt 
237a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser
247a982d89SJeremy L. Thompson /// @{
257a982d89SJeremy L. Thompson 
26ca94c3ddSJeremy L Thompson /// Argument for @ref CeedOperatorSetField() indicating that the field does not require a `CeedBasis`
27356036faSJeremy L Thompson const CeedBasis CEED_BASIS_NONE = &ceed_basis_none;
28356036faSJeremy 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 /**
383778dbaaSJeremy L Thompson   @brief Compute Chebyshev polynomial values at a point
393778dbaaSJeremy L Thompson 
403778dbaaSJeremy L Thompson   @param[in]  x           Coordinate to evaluate Chebyshev polynomials at
41ca94c3ddSJeremy L Thompson   @param[in]  n           Number of Chebyshev polynomials to evaluate, `n >= 2`
423778dbaaSJeremy L Thompson   @param[out] chebyshev_x Array of Chebyshev polynomial values
433778dbaaSJeremy L Thompson 
443778dbaaSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
453778dbaaSJeremy L Thompson 
463778dbaaSJeremy L Thompson   @ref Developer
473778dbaaSJeremy L Thompson **/
483778dbaaSJeremy L Thompson static int CeedChebyshevPolynomialsAtPoint(CeedScalar x, CeedInt n, CeedScalar *chebyshev_x) {
493778dbaaSJeremy L Thompson   chebyshev_x[0] = 1.0;
503778dbaaSJeremy L Thompson   chebyshev_x[1] = 2 * x;
513778dbaaSJeremy L Thompson   for (CeedInt i = 2; i < n; i++) chebyshev_x[i] = 2 * x * chebyshev_x[i - 1] - chebyshev_x[i - 2];
523778dbaaSJeremy L Thompson   return CEED_ERROR_SUCCESS;
533778dbaaSJeremy L Thompson }
543778dbaaSJeremy L Thompson 
553778dbaaSJeremy L Thompson /**
563778dbaaSJeremy L Thompson   @brief Compute values of the derivative of Chebyshev polynomials at a point
573778dbaaSJeremy L Thompson 
583778dbaaSJeremy L Thompson   @param[in]  x            Coordinate to evaluate derivative of Chebyshev polynomials at
59ca94c3ddSJeremy L Thompson   @param[in]  n            Number of Chebyshev polynomials to evaluate, `n >= 2`
606cec60aaSJed Brown   @param[out] chebyshev_dx Array of Chebyshev polynomial derivative values
613778dbaaSJeremy L Thompson 
623778dbaaSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
633778dbaaSJeremy L Thompson 
643778dbaaSJeremy L Thompson   @ref Developer
653778dbaaSJeremy L Thompson **/
663778dbaaSJeremy L Thompson static int CeedChebyshevDerivativeAtPoint(CeedScalar x, CeedInt n, CeedScalar *chebyshev_dx) {
673778dbaaSJeremy L Thompson   CeedScalar chebyshev_x[3];
683778dbaaSJeremy L Thompson 
693778dbaaSJeremy L Thompson   chebyshev_x[1]  = 1.0;
703778dbaaSJeremy L Thompson   chebyshev_x[2]  = 2 * x;
713778dbaaSJeremy L Thompson   chebyshev_dx[0] = 0.0;
723778dbaaSJeremy L Thompson   chebyshev_dx[1] = 2.0;
733778dbaaSJeremy L Thompson   for (CeedInt i = 2; i < n; i++) {
743778dbaaSJeremy L Thompson     chebyshev_x[0]  = chebyshev_x[1];
753778dbaaSJeremy L Thompson     chebyshev_x[1]  = chebyshev_x[2];
763778dbaaSJeremy L Thompson     chebyshev_x[2]  = 2 * x * chebyshev_x[1] - chebyshev_x[0];
773778dbaaSJeremy L Thompson     chebyshev_dx[i] = 2 * x * chebyshev_dx[i - 1] + 2 * chebyshev_x[1] - chebyshev_dx[i - 2];
783778dbaaSJeremy L Thompson   }
793778dbaaSJeremy L Thompson   return CEED_ERROR_SUCCESS;
803778dbaaSJeremy L Thompson }
813778dbaaSJeremy L Thompson 
823778dbaaSJeremy L Thompson /**
83ca94c3ddSJeremy L Thompson   @brief Compute Householder reflection.
847a982d89SJeremy L. Thompson 
85ca94c3ddSJeremy L Thompson   Computes \f$A = (I - b v v^T) A\f$, where \f$A\f$ is an \f$m \times n\f$ matrix indexed as `A[i*row + j*col]`.
867a982d89SJeremy L. Thompson 
877a982d89SJeremy L. Thompson   @param[in,out] A   Matrix to apply Householder reflection to, in place
88ea61e9acSJeremy L Thompson   @param[in]     v   Householder vector
89ea61e9acSJeremy L Thompson   @param[in]     b   Scaling factor
90ca94c3ddSJeremy L Thompson   @param[in]     m   Number of rows in `A`
91ca94c3ddSJeremy L Thompson   @param[in]     n   Number of columns in `A`
92ea61e9acSJeremy L Thompson   @param[in]     row Row stride
93ea61e9acSJeremy L Thompson   @param[in]     col Col stride
947a982d89SJeremy L. Thompson 
957a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
967a982d89SJeremy L. Thompson 
977a982d89SJeremy L. Thompson   @ref Developer
987a982d89SJeremy L. Thompson **/
992b730f8bSJeremy L Thompson static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v, CeedScalar b, CeedInt m, CeedInt n, CeedInt row, CeedInt col) {
1007a982d89SJeremy L. Thompson   for (CeedInt j = 0; j < n; j++) {
1017a982d89SJeremy L. Thompson     CeedScalar w = A[0 * row + j * col];
1021c66c397SJeremy L Thompson 
1032b730f8bSJeremy L Thompson     for (CeedInt i = 1; i < m; i++) w += v[i] * A[i * row + j * col];
1047a982d89SJeremy L. Thompson     A[0 * row + j * col] -= b * w;
1052b730f8bSJeremy L Thompson     for (CeedInt i = 1; i < m; i++) A[i * row + j * col] -= b * w * v[i];
1067a982d89SJeremy L. Thompson   }
107e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1087a982d89SJeremy L. Thompson }
1097a982d89SJeremy L. Thompson 
1107a982d89SJeremy L. Thompson /**
1117a982d89SJeremy L. Thompson   @brief Compute Givens rotation
1127a982d89SJeremy L. Thompson 
113ca94c3ddSJeremy L Thompson   Computes \f$A = G A\f$ (or \f$G^T A\f$ in transpose mode), where \f$A\f$ is an \f$m \times n\f$ matrix indexed as `A[i*n + j*m]`.
1147a982d89SJeremy L. Thompson 
1157a982d89SJeremy L. Thompson   @param[in,out] A      Row major matrix to apply Givens rotation to, in place
116ea61e9acSJeremy L Thompson   @param[in]     c      Cosine factor
117ea61e9acSJeremy L Thompson   @param[in]     s      Sine factor
118ca94c3ddSJeremy L Thompson   @param[in]     t_mode @ref CEED_NOTRANSPOSE to rotate the basis counter-clockwise, which has the effect of rotating columns of `A` clockwise;
1194cc79fe7SJed Brown                           @ref CEED_TRANSPOSE for the opposite rotation
120ea61e9acSJeremy L Thompson   @param[in]     i      First row/column to apply rotation
121ea61e9acSJeremy L Thompson   @param[in]     k      Second row/column to apply rotation
122ca94c3ddSJeremy L Thompson   @param[in]     m      Number of rows in `A`
123ca94c3ddSJeremy L Thompson   @param[in]     n      Number of columns in `A`
1247a982d89SJeremy L. Thompson 
1257a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
1267a982d89SJeremy L. Thompson 
1277a982d89SJeremy L. Thompson   @ref Developer
1287a982d89SJeremy L. Thompson **/
1292b730f8bSJeremy L Thompson static int CeedGivensRotation(CeedScalar *A, CeedScalar c, CeedScalar s, CeedTransposeMode t_mode, CeedInt i, CeedInt k, CeedInt m, CeedInt n) {
130d1d35e2fSjeremylt   CeedInt stride_j = 1, stride_ik = m, num_its = n;
1311c66c397SJeremy L Thompson 
132d1d35e2fSjeremylt   if (t_mode == CEED_NOTRANSPOSE) {
1332b730f8bSJeremy L Thompson     stride_j  = n;
1342b730f8bSJeremy L Thompson     stride_ik = 1;
1352b730f8bSJeremy L Thompson     num_its   = m;
1367a982d89SJeremy L. Thompson   }
1377a982d89SJeremy L. Thompson 
1387a982d89SJeremy L. Thompson   // Apply rotation
139d1d35e2fSjeremylt   for (CeedInt j = 0; j < num_its; j++) {
140d1d35e2fSjeremylt     CeedScalar tau1 = A[i * stride_ik + j * stride_j], tau2 = A[k * stride_ik + j * stride_j];
1411c66c397SJeremy L Thompson 
142d1d35e2fSjeremylt     A[i * stride_ik + j * stride_j] = c * tau1 - s * tau2;
143d1d35e2fSjeremylt     A[k * stride_ik + j * stride_j] = s * tau1 + c * tau2;
1447a982d89SJeremy L. Thompson   }
145e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1467a982d89SJeremy L. Thompson }
1477a982d89SJeremy L. Thompson 
1487a982d89SJeremy L. Thompson /**
149ca94c3ddSJeremy L Thompson   @brief View an array stored in a `CeedBasis`
1507a982d89SJeremy L. Thompson 
1510a0da059Sjeremylt   @param[in] name   Name of array
152d1d35e2fSjeremylt   @param[in] fp_fmt Printing format
1530a0da059Sjeremylt   @param[in] m      Number of rows in array
1540a0da059Sjeremylt   @param[in] n      Number of columns in array
1550a0da059Sjeremylt   @param[in] a      Array to be viewed
156ca94c3ddSJeremy L Thompson   @param[in] stream Stream to view to, e.g., `stdout`
1577a982d89SJeremy L. Thompson 
1587a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
1597a982d89SJeremy L. Thompson 
1607a982d89SJeremy L. Thompson   @ref Developer
1617a982d89SJeremy L. Thompson **/
1622b730f8bSJeremy L Thompson static int CeedScalarView(const char *name, const char *fp_fmt, CeedInt m, CeedInt n, const CeedScalar *a, FILE *stream) {
163edf04919SJeremy L Thompson   if (m > 1) {
164edf04919SJeremy L Thompson     fprintf(stream, "  %s:\n", name);
165edf04919SJeremy L Thompson   } else {
166edf04919SJeremy L Thompson     char padded_name[12];
167edf04919SJeremy L Thompson 
168edf04919SJeremy L Thompson     snprintf(padded_name, 11, "%s:", name);
169edf04919SJeremy L Thompson     fprintf(stream, "  %-10s", padded_name);
170edf04919SJeremy L Thompson   }
17192ae7e47SJeremy L Thompson   for (CeedInt i = 0; i < m; i++) {
172edf04919SJeremy L Thompson     if (m > 1) fprintf(stream, "    [%" CeedInt_FMT "]", i);
1732b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < n; j++) fprintf(stream, fp_fmt, fabs(a[i * n + j]) > 1E-14 ? a[i * n + j] : 0);
1747a982d89SJeremy L. Thompson     fputs("\n", stream);
1757a982d89SJeremy L. Thompson   }
176e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1777a982d89SJeremy L. Thompson }
1787a982d89SJeremy L. Thompson 
179a76a04e7SJeremy L Thompson /**
180ea61e9acSJeremy L Thompson   @brief Create the interpolation and gradient matrices for projection from the nodes of `basis_from` to the nodes of `basis_to`.
181ba59ac12SSebastian Grimberg 
18215ad3917SSebastian Grimberg   The interpolation is given by `interp_project = interp_to^+ * interp_from`, where the pseudoinverse `interp_to^+` is given by QR factorization.
183ca94c3ddSJeremy L Thompson   The gradient is given by `grad_project = interp_to^+ * grad_from`, and is only computed for \f$H^1\f$ spaces otherwise it should not be used.
18415ad3917SSebastian Grimberg 
185ba59ac12SSebastian Grimberg   Note: `basis_from` and `basis_to` must have compatible quadrature spaces.
186a76a04e7SJeremy L Thompson 
187ca94c3ddSJeremy L Thompson   @param[in]  basis_from     `CeedBasis` to project from
188ca94c3ddSJeremy L Thompson   @param[in]  basis_to       `CeedBasis` to project to
189ca94c3ddSJeremy L Thompson   @param[out] interp_project Address of the variable where the newly created interpolation matrix will be stored
190ca94c3ddSJeremy L Thompson   @param[out] grad_project   Address of the variable where the newly created gradient matrix will be stored
191a76a04e7SJeremy L Thompson 
192a76a04e7SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
193a76a04e7SJeremy L Thompson 
194a76a04e7SJeremy L Thompson   @ref Developer
195a76a04e7SJeremy L Thompson **/
1962b730f8bSJeremy L Thompson static int CeedBasisCreateProjectionMatrices(CeedBasis basis_from, CeedBasis basis_to, CeedScalar **interp_project, CeedScalar **grad_project) {
197a76a04e7SJeremy L Thompson   Ceed    ceed;
198e104ad11SJames Wright   bool    are_both_tensor;
1991c66c397SJeremy L Thompson   CeedInt Q, Q_to, Q_from, P_to, P_from;
2001c66c397SJeremy L Thompson 
2012b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetCeed(basis_to, &ceed));
202a76a04e7SJeremy L Thompson 
203a76a04e7SJeremy L Thompson   // Check for compatible quadrature spaces
2042b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumQuadraturePoints(basis_to, &Q_to));
2052b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumQuadraturePoints(basis_from, &Q_from));
2063f08121cSJeremy L Thompson   CeedCheck(Q_to == Q_from, ceed, CEED_ERROR_DIMENSION,
2073f08121cSJeremy L Thompson             "Bases must have compatible quadrature spaces."
20823622755SJeremy L Thompson             " 'basis_from' has %" CeedInt_FMT " points and 'basis_to' has %" CeedInt_FMT,
2093f08121cSJeremy L Thompson             Q_from, Q_to);
2101c66c397SJeremy L Thompson   Q = Q_to;
211a76a04e7SJeremy L Thompson 
21214556e63SJeremy L Thompson   // Check for matching tensor or non-tensor
213e104ad11SJames Wright   {
214e104ad11SJames Wright     bool is_tensor_to, is_tensor_from;
215e104ad11SJames Wright 
2162b730f8bSJeremy L Thompson     CeedCall(CeedBasisIsTensor(basis_to, &is_tensor_to));
2172b730f8bSJeremy L Thompson     CeedCall(CeedBasisIsTensor(basis_from, &is_tensor_from));
218e104ad11SJames Wright     are_both_tensor = is_tensor_to && is_tensor_from;
219e104ad11SJames Wright   }
220e104ad11SJames Wright   if (are_both_tensor) {
2212b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_to));
2222b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_from));
2232b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumQuadraturePoints1D(basis_from, &Q));
2246574a04fSJeremy L Thompson   } else {
2252b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis_to, &P_to));
2262b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis_from, &P_from));
227a76a04e7SJeremy L Thompson   }
228a76a04e7SJeremy L Thompson 
22915ad3917SSebastian Grimberg   // Check for matching FE space
23015ad3917SSebastian Grimberg   CeedFESpace fe_space_to, fe_space_from;
2313f08121cSJeremy L Thompson 
23215ad3917SSebastian Grimberg   CeedCall(CeedBasisGetFESpace(basis_to, &fe_space_to));
23315ad3917SSebastian Grimberg   CeedCall(CeedBasisGetFESpace(basis_from, &fe_space_from));
2343f08121cSJeremy L Thompson   CeedCheck(fe_space_to == fe_space_from, ceed, CEED_ERROR_MINOR,
2353f08121cSJeremy L Thompson             "Bases must both be the same FE space type."
2363f08121cSJeremy L Thompson             " 'basis_from' is a %s and 'basis_to' is a %s",
2373f08121cSJeremy L Thompson             CeedFESpaces[fe_space_from], CeedFESpaces[fe_space_to]);
23815ad3917SSebastian Grimberg 
23914556e63SJeremy L Thompson   // Get source matrices
24015ad3917SSebastian Grimberg   CeedInt           dim, q_comp = 1;
2412247a93fSRezgar Shakeri   CeedScalar       *interp_to_inv, *interp_from;
2421c66c397SJeremy L Thompson   const CeedScalar *interp_to_source = NULL, *interp_from_source = NULL, *grad_from_source = NULL;
2431c66c397SJeremy L Thompson 
244b3ed00e5SJames Wright   CeedCall(CeedBasisGetDimension(basis_from, &dim));
245e104ad11SJames Wright   if (are_both_tensor) {
2462b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetInterp1D(basis_to, &interp_to_source));
2472b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetInterp1D(basis_from, &interp_from_source));
248a76a04e7SJeremy L Thompson   } else {
24915ad3917SSebastian Grimberg     CeedCall(CeedBasisGetNumQuadratureComponents(basis_from, CEED_EVAL_INTERP, &q_comp));
2502b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetInterp(basis_to, &interp_to_source));
2512b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetInterp(basis_from, &interp_from_source));
25215ad3917SSebastian Grimberg   }
25315ad3917SSebastian Grimberg   CeedCall(CeedMalloc(Q * P_from * q_comp, &interp_from));
25415ad3917SSebastian Grimberg   CeedCall(CeedCalloc(P_to * P_from, interp_project));
25515ad3917SSebastian Grimberg 
25615ad3917SSebastian Grimberg   // `grad_project = interp_to^+ * grad_from` is computed for the H^1 space case so the
257de05fbb2SSebastian Grimberg   // projection basis will have a gradient operation (allocated even if not H^1 for the
258de05fbb2SSebastian Grimberg   // basis construction later on)
25915ad3917SSebastian Grimberg   if (fe_space_to == CEED_FE_SPACE_H1) {
260e104ad11SJames Wright     if (are_both_tensor) {
26115ad3917SSebastian Grimberg       CeedCall(CeedBasisGetGrad1D(basis_from, &grad_from_source));
26215ad3917SSebastian Grimberg     } else {
2632b730f8bSJeremy L Thompson       CeedCall(CeedBasisGetGrad(basis_from, &grad_from_source));
264a76a04e7SJeremy L Thompson     }
265de05fbb2SSebastian Grimberg   }
266e104ad11SJames Wright   CeedCall(CeedCalloc(P_to * P_from * (are_both_tensor ? 1 : dim), grad_project));
26715ad3917SSebastian Grimberg 
2682247a93fSRezgar Shakeri   // Compute interp_to^+, pseudoinverse of interp_to
2692247a93fSRezgar Shakeri   CeedCall(CeedCalloc(Q * q_comp * P_to, &interp_to_inv));
2701203703bSJeremy L Thompson   CeedCall(CeedMatrixPseudoinverse(ceed, interp_to_source, Q * q_comp, P_to, interp_to_inv));
27114556e63SJeremy L Thompson   // Build matrices
272e104ad11SJames Wright   CeedInt     num_matrices = 1 + (fe_space_to == CEED_FE_SPACE_H1) * (are_both_tensor ? 1 : dim);
27314556e63SJeremy L Thompson   CeedScalar *input_from[num_matrices], *output_project[num_matrices];
2741c66c397SJeremy L Thompson 
27514556e63SJeremy L Thompson   input_from[0]     = (CeedScalar *)interp_from_source;
27614556e63SJeremy L Thompson   output_project[0] = *interp_project;
27714556e63SJeremy L Thompson   for (CeedInt m = 1; m < num_matrices; m++) {
27814556e63SJeremy L Thompson     input_from[m]     = (CeedScalar *)&grad_from_source[(m - 1) * Q * P_from];
27902af4036SJeremy L Thompson     output_project[m] = &((*grad_project)[(m - 1) * P_to * P_from]);
28014556e63SJeremy L Thompson   }
28114556e63SJeremy L Thompson   for (CeedInt m = 0; m < num_matrices; m++) {
2822247a93fSRezgar Shakeri     // output_project = interp_to^+ * interp_from
28315ad3917SSebastian Grimberg     memcpy(interp_from, input_from[m], Q * P_from * q_comp * sizeof(input_from[m][0]));
2842247a93fSRezgar Shakeri     CeedCall(CeedMatrixMatrixMultiply(ceed, interp_to_inv, input_from[m], output_project[m], P_to, P_from, Q * q_comp));
2852247a93fSRezgar Shakeri     // Round zero to machine precision
2862247a93fSRezgar Shakeri     for (CeedInt i = 0; i < P_to * P_from; i++) {
2872247a93fSRezgar Shakeri       if (fabs(output_project[m][i]) < 10 * CEED_EPSILON) output_project[m][i] = 0.0;
288a76a04e7SJeremy L Thompson     }
28914556e63SJeremy L Thompson   }
29014556e63SJeremy L Thompson 
29114556e63SJeremy L Thompson   // Cleanup
2922247a93fSRezgar Shakeri   CeedCall(CeedFree(&interp_to_inv));
2932b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_from));
294a76a04e7SJeremy L Thompson   return CEED_ERROR_SUCCESS;
295a76a04e7SJeremy L Thompson }
296a76a04e7SJeremy L Thompson 
2970b31fde2SJeremy L Thompson /**
2980b31fde2SJeremy L Thompson   @brief Check input vector dimensions for CeedBasisApply[Add]AtPoints
2990b31fde2SJeremy L Thompson 
3000b31fde2SJeremy L Thompson   @param[in]  basis      `CeedBasis` to evaluate
3010b31fde2SJeremy L Thompson   @param[in]  num_elem   The number of elements to apply the basis evaluation to;
3020b31fde2SJeremy L Thompson                           the backend will specify the ordering in @ref CeedElemRestrictionCreate()
3030b31fde2SJeremy L Thompson   @param[in]  num_points Array of the number of points to apply the basis evaluation to in each element, size `num_elem`
3040b31fde2SJeremy L Thompson   @param[in]  t_mode     @ref CEED_NOTRANSPOSE to evaluate from nodes to points;
3050b31fde2SJeremy L Thompson                            @ref CEED_TRANSPOSE to apply the transpose, mapping from points to nodes
3060b31fde2SJeremy L Thompson   @param[in]  eval_mode  @ref CEED_EVAL_INTERP to use interpolated values,
3070b31fde2SJeremy L Thompson                            @ref CEED_EVAL_GRAD to use gradients,
3080b31fde2SJeremy L Thompson                            @ref CEED_EVAL_WEIGHT to use quadrature weights
3090b31fde2SJeremy L Thompson   @param[in]  x_ref      `CeedVector` holding reference coordinates of each point
3100b31fde2SJeremy L Thompson   @param[in]  u          Input `CeedVector`, of length `num_nodes * num_comp` for @ref CEED_NOTRANSPOSE
3110b31fde2SJeremy L Thompson   @param[out] v          Output `CeedVector`, of length `num_points * num_q_comp` for @ref CEED_NOTRANSPOSE with @ref CEED_EVAL_INTERP
3120b31fde2SJeremy L Thompson 
3130b31fde2SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
3140b31fde2SJeremy L Thompson 
3150b31fde2SJeremy L Thompson   @ref Developer
3160b31fde2SJeremy L Thompson **/
3170b31fde2SJeremy L Thompson static int CeedBasisApplyAtPointsCheckDims(CeedBasis basis, CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode,
3180b31fde2SJeremy L Thompson                                            CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) {
3190b31fde2SJeremy L Thompson   CeedInt  dim, num_comp, num_q_comp, num_nodes, P_1d = 1, Q_1d = 1, total_num_points = 0;
3200b31fde2SJeremy L Thompson   CeedSize x_length = 0, u_length = 0, v_length;
3210b31fde2SJeremy L Thompson   Ceed     ceed;
3220b31fde2SJeremy L Thompson 
3230b31fde2SJeremy L Thompson   CeedCall(CeedBasisGetCeed(basis, &ceed));
3240b31fde2SJeremy L Thompson   CeedCall(CeedBasisGetDimension(basis, &dim));
3250b31fde2SJeremy L Thompson   CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d));
3260b31fde2SJeremy L Thompson   CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
3270b31fde2SJeremy L Thompson   CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
3280b31fde2SJeremy L Thompson   CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &num_q_comp));
3290b31fde2SJeremy L Thompson   CeedCall(CeedBasisGetNumNodes(basis, &num_nodes));
3300b31fde2SJeremy L Thompson   CeedCall(CeedVectorGetLength(v, &v_length));
3310b31fde2SJeremy L Thompson   if (x_ref != CEED_VECTOR_NONE) CeedCall(CeedVectorGetLength(x_ref, &x_length));
3320b31fde2SJeremy L Thompson   if (u != CEED_VECTOR_NONE) CeedCall(CeedVectorGetLength(u, &u_length));
3330b31fde2SJeremy L Thompson 
3340b31fde2SJeremy L Thompson   // Check compatibility coordinates vector
3350b31fde2SJeremy L Thompson   for (CeedInt i = 0; i < num_elem; i++) total_num_points += num_points[i];
336*19a04db8SJeremy L Thompson   CeedCheck((x_length >= (CeedSize)total_num_points * (CeedSize)dim) || (eval_mode == CEED_EVAL_WEIGHT), ceed, CEED_ERROR_DIMENSION,
3370b31fde2SJeremy L Thompson             "Length of reference coordinate vector incompatible with basis dimension and number of points."
3380b31fde2SJeremy L Thompson             " Found reference coordinate vector of length %" CeedSize_FMT ", not of length %" CeedSize_FMT ".",
339*19a04db8SJeremy L Thompson             x_length, (CeedSize)total_num_points * (CeedSize)dim);
3400b31fde2SJeremy L Thompson 
3410b31fde2SJeremy L Thompson   // Check CEED_EVAL_WEIGHT only on CEED_NOTRANSPOSE
3420b31fde2SJeremy L Thompson   CeedCheck(eval_mode != CEED_EVAL_WEIGHT || t_mode == CEED_NOTRANSPOSE, ceed, CEED_ERROR_UNSUPPORTED,
3430b31fde2SJeremy L Thompson             "CEED_EVAL_WEIGHT only supported with CEED_NOTRANSPOSE");
3440b31fde2SJeremy L Thompson 
3450b31fde2SJeremy L Thompson   // Check vector lengths to prevent out of bounds issues
3460b31fde2SJeremy L Thompson   bool has_good_dims = true;
3470b31fde2SJeremy L Thompson   switch (eval_mode) {
3480b31fde2SJeremy L Thompson     case CEED_EVAL_INTERP:
349*19a04db8SJeremy L Thompson       has_good_dims = ((t_mode == CEED_TRANSPOSE && (u_length >= (CeedSize)total_num_points * (CeedSize)num_q_comp ||
350*19a04db8SJeremy L Thompson                                                      v_length >= (CeedSize)num_elem * (CeedSize)num_nodes * (CeedSize)num_comp)) ||
351*19a04db8SJeremy L Thompson                        (t_mode == CEED_NOTRANSPOSE && (v_length >= (CeedSize)total_num_points * (CeedSize)num_q_comp ||
352*19a04db8SJeremy L Thompson                                                        u_length >= (CeedSize)num_elem * (CeedSize)num_nodes * (CeedSize)num_comp)));
3530b31fde2SJeremy L Thompson       break;
3540b31fde2SJeremy L Thompson     case CEED_EVAL_GRAD:
355*19a04db8SJeremy L Thompson       has_good_dims = ((t_mode == CEED_TRANSPOSE && (u_length >= (CeedSize)total_num_points * (CeedSize)num_q_comp * (CeedSize)dim ||
356*19a04db8SJeremy L Thompson                                                      v_length >= (CeedSize)num_elem * (CeedSize)num_nodes * (CeedSize)num_comp)) ||
357*19a04db8SJeremy L Thompson                        (t_mode == CEED_NOTRANSPOSE && (v_length >= (CeedSize)total_num_points * (CeedSize)num_q_comp * (CeedSize)dim ||
358*19a04db8SJeremy L Thompson                                                        u_length >= (CeedSize)num_elem * (CeedSize)num_nodes * (CeedSize)num_comp)));
3590b31fde2SJeremy L Thompson       break;
3600b31fde2SJeremy L Thompson     case CEED_EVAL_WEIGHT:
3610b31fde2SJeremy L Thompson       has_good_dims = t_mode == CEED_NOTRANSPOSE && (v_length >= total_num_points);
3620b31fde2SJeremy L Thompson       break;
3630b31fde2SJeremy L Thompson       // LCOV_EXCL_START
3640b31fde2SJeremy L Thompson     case CEED_EVAL_NONE:
3650b31fde2SJeremy L Thompson     case CEED_EVAL_DIV:
3660b31fde2SJeremy L Thompson     case CEED_EVAL_CURL:
3670b31fde2SJeremy L Thompson       return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Evaluation at arbitrary points not supported for %s", CeedEvalModes[eval_mode]);
3680b31fde2SJeremy L Thompson       // LCOV_EXCL_STOP
3690b31fde2SJeremy L Thompson   }
3700b31fde2SJeremy L Thompson   CeedCheck(has_good_dims, ceed, CEED_ERROR_DIMENSION, "Input/output vectors too short for basis and evaluation mode");
3710b31fde2SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3720b31fde2SJeremy L Thompson }
3730b31fde2SJeremy L Thompson 
3740b31fde2SJeremy L Thompson /**
3750b31fde2SJeremy L Thompson   @brief Default implimentation to apply basis evaluation from nodes to arbitrary points
3760b31fde2SJeremy L Thompson 
3770b31fde2SJeremy L Thompson   @param[in]  basis      `CeedBasis` to evaluate
3780b31fde2SJeremy L Thompson   @param[in]  apply_add  Sum result into target vector or overwrite
3790b31fde2SJeremy L Thompson   @param[in]  num_elem   The number of elements to apply the basis evaluation to;
3800b31fde2SJeremy L Thompson                           the backend will specify the ordering in @ref CeedElemRestrictionCreate()
3810b31fde2SJeremy L Thompson   @param[in]  num_points Array of the number of points to apply the basis evaluation to in each element, size `num_elem`
3820b31fde2SJeremy L Thompson   @param[in]  t_mode     @ref CEED_NOTRANSPOSE to evaluate from nodes to points;
3830b31fde2SJeremy L Thompson                            @ref CEED_TRANSPOSE to apply the transpose, mapping from points to nodes
3840b31fde2SJeremy L Thompson   @param[in]  eval_mode  @ref CEED_EVAL_INTERP to use interpolated values,
3850b31fde2SJeremy L Thompson                            @ref CEED_EVAL_GRAD to use gradients,
3860b31fde2SJeremy L Thompson                            @ref CEED_EVAL_WEIGHT to use quadrature weights
3870b31fde2SJeremy L Thompson   @param[in]  x_ref      `CeedVector` holding reference coordinates of each point
3880b31fde2SJeremy L Thompson   @param[in]  u          Input `CeedVector`, of length `num_nodes * num_comp` for @ref CEED_NOTRANSPOSE
3890b31fde2SJeremy L Thompson   @param[out] v          Output `CeedVector`, of length `num_points * num_q_comp` for @ref CEED_NOTRANSPOSE with @ref CEED_EVAL_INTERP
3900b31fde2SJeremy L Thompson 
3910b31fde2SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
3920b31fde2SJeremy L Thompson 
3930b31fde2SJeremy L Thompson   @ref Developer
3940b31fde2SJeremy L Thompson **/
3950b31fde2SJeremy L Thompson static int CeedBasisApplyAtPoints_Core(CeedBasis basis, bool apply_add, CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode,
3960b31fde2SJeremy L Thompson                                        CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) {
3970b31fde2SJeremy L Thompson   CeedInt dim, num_comp, P_1d = 1, Q_1d = 1, total_num_points = num_points[0];
3980b31fde2SJeremy L Thompson   Ceed    ceed;
3990b31fde2SJeremy L Thompson 
4000b31fde2SJeremy L Thompson   CeedCall(CeedBasisGetCeed(basis, &ceed));
4010b31fde2SJeremy L Thompson   CeedCall(CeedBasisGetDimension(basis, &dim));
4020b31fde2SJeremy L Thompson   // Inserting check because clang-tidy doesn't understand this cannot occur
4030b31fde2SJeremy L Thompson   CeedCheck(dim > 0, ceed, CEED_ERROR_UNSUPPORTED, "Malformed CeedBasis, dim > 0 is required");
4040b31fde2SJeremy L Thompson   CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d));
4050b31fde2SJeremy L Thompson   CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
4060b31fde2SJeremy L Thompson   CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
4070b31fde2SJeremy L Thompson 
4080b31fde2SJeremy L Thompson   // Default implementation
4090b31fde2SJeremy L Thompson   {
4100b31fde2SJeremy L Thompson     bool is_tensor_basis;
4110b31fde2SJeremy L Thompson 
4120b31fde2SJeremy L Thompson     CeedCall(CeedBasisIsTensor(basis, &is_tensor_basis));
4130b31fde2SJeremy L Thompson     CeedCheck(is_tensor_basis, ceed, CEED_ERROR_UNSUPPORTED, "Evaluation at arbitrary points only supported for tensor product bases");
4140b31fde2SJeremy L Thompson   }
4150b31fde2SJeremy L Thompson   CeedCheck(num_elem == 1, ceed, CEED_ERROR_UNSUPPORTED, "Evaluation at arbitrary  points only supported for a single element at a time");
4160b31fde2SJeremy L Thompson   if (eval_mode == CEED_EVAL_WEIGHT) {
4170b31fde2SJeremy L Thompson     CeedCall(CeedVectorSetValue(v, 1.0));
4180b31fde2SJeremy L Thompson     return CEED_ERROR_SUCCESS;
4190b31fde2SJeremy L Thompson   }
4200b31fde2SJeremy L Thompson   if (!basis->basis_chebyshev) {
4210b31fde2SJeremy L Thompson     // Build basis mapping from nodes to Chebyshev coefficients
4220b31fde2SJeremy L Thompson     CeedScalar       *chebyshev_interp_1d, *chebyshev_grad_1d, *chebyshev_q_weight_1d;
4230b31fde2SJeremy L Thompson     const CeedScalar *q_ref_1d;
4240b31fde2SJeremy L Thompson 
4250b31fde2SJeremy L Thompson     CeedCall(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
4260b31fde2SJeremy L Thompson     CeedCall(CeedCalloc(P_1d * Q_1d, &chebyshev_grad_1d));
4270b31fde2SJeremy L Thompson     CeedCall(CeedCalloc(Q_1d, &chebyshev_q_weight_1d));
4280b31fde2SJeremy L Thompson     CeedCall(CeedBasisGetQRef(basis, &q_ref_1d));
4290b31fde2SJeremy L Thompson     CeedCall(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d));
4300b31fde2SJeremy L Thompson 
4310b31fde2SJeremy L Thompson     CeedCall(CeedVectorCreate(ceed, num_comp * CeedIntPow(Q_1d, dim), &basis->vec_chebyshev));
4320b31fde2SJeremy L Thompson     CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P_1d, Q_1d, chebyshev_interp_1d, chebyshev_grad_1d, q_ref_1d, chebyshev_q_weight_1d,
4330b31fde2SJeremy L Thompson                                      &basis->basis_chebyshev));
4340b31fde2SJeremy L Thompson 
4350b31fde2SJeremy L Thompson     // Cleanup
4360b31fde2SJeremy L Thompson     CeedCall(CeedFree(&chebyshev_interp_1d));
4370b31fde2SJeremy L Thompson     CeedCall(CeedFree(&chebyshev_grad_1d));
4380b31fde2SJeremy L Thompson     CeedCall(CeedFree(&chebyshev_q_weight_1d));
4390b31fde2SJeremy L Thompson   }
4400b31fde2SJeremy L Thompson 
4410b31fde2SJeremy L Thompson   // Create TensorContract object if needed, such as a basis from the GPU backends
4420b31fde2SJeremy L Thompson   if (!basis->contract) {
4430b31fde2SJeremy L Thompson     Ceed      ceed_ref;
4440b31fde2SJeremy L Thompson     CeedBasis basis_ref = NULL;
4450b31fde2SJeremy L Thompson 
4460b31fde2SJeremy L Thompson     CeedCall(CeedInit("/cpu/self", &ceed_ref));
4470b31fde2SJeremy L Thompson     // Only need matching tensor contraction dimensions, any type of basis will work
4480b31fde2SJeremy L Thompson     CeedCall(CeedBasisCreateTensorH1Lagrange(ceed_ref, dim, num_comp, P_1d, Q_1d, CEED_GAUSS, &basis_ref));
4490b31fde2SJeremy L Thompson     // Note - clang-tidy doesn't know basis_ref->contract must be valid here
4500b31fde2SJeremy L Thompson     CeedCheck(basis_ref && basis_ref->contract, ceed, CEED_ERROR_UNSUPPORTED, "Reference CPU ceed failed to create a tensor contraction object");
4510b31fde2SJeremy L Thompson     CeedCall(CeedTensorContractReferenceCopy(basis_ref->contract, &basis->contract));
4520b31fde2SJeremy L Thompson     CeedCall(CeedBasisDestroy(&basis_ref));
4530b31fde2SJeremy L Thompson     CeedCall(CeedDestroy(&ceed_ref));
4540b31fde2SJeremy L Thompson   }
4550b31fde2SJeremy L Thompson 
4560b31fde2SJeremy L Thompson   // Basis evaluation
4570b31fde2SJeremy L Thompson   switch (t_mode) {
4580b31fde2SJeremy L Thompson     case CEED_NOTRANSPOSE: {
4590b31fde2SJeremy L Thompson       // Nodes to arbitrary points
4600b31fde2SJeremy L Thompson       CeedScalar       *v_array;
4610b31fde2SJeremy L Thompson       const CeedScalar *chebyshev_coeffs, *x_array_read;
4620b31fde2SJeremy L Thompson 
4630b31fde2SJeremy L Thompson       // -- Interpolate to Chebyshev coefficients
4640b31fde2SJeremy L Thompson       CeedCall(CeedBasisApply(basis->basis_chebyshev, 1, CEED_NOTRANSPOSE, CEED_EVAL_INTERP, u, basis->vec_chebyshev));
4650b31fde2SJeremy L Thompson 
4660b31fde2SJeremy L Thompson       // -- Evaluate Chebyshev polynomials at arbitrary points
4670b31fde2SJeremy L Thompson       CeedCall(CeedVectorGetArrayRead(basis->vec_chebyshev, CEED_MEM_HOST, &chebyshev_coeffs));
4680b31fde2SJeremy L Thompson       CeedCall(CeedVectorGetArrayRead(x_ref, CEED_MEM_HOST, &x_array_read));
4690b31fde2SJeremy L Thompson       CeedCall(CeedVectorGetArrayWrite(v, CEED_MEM_HOST, &v_array));
4700b31fde2SJeremy L Thompson       switch (eval_mode) {
4710b31fde2SJeremy L Thompson         case CEED_EVAL_INTERP: {
4720b31fde2SJeremy L Thompson           CeedScalar tmp[2][num_comp * CeedIntPow(Q_1d, dim)], chebyshev_x[Q_1d];
4730b31fde2SJeremy L Thompson 
4740b31fde2SJeremy L Thompson           // ---- Values at point
4750b31fde2SJeremy L Thompson           for (CeedInt p = 0; p < total_num_points; p++) {
4760b31fde2SJeremy L Thompson             CeedInt pre = num_comp * CeedIntPow(Q_1d, dim - 1), post = 1;
4770b31fde2SJeremy L Thompson 
4780b31fde2SJeremy L Thompson             for (CeedInt d = 0; d < dim; d++) {
4790b31fde2SJeremy L Thompson               // ------ Tensor contract with current Chebyshev polynomial values
4800b31fde2SJeremy L Thompson               CeedCall(CeedChebyshevPolynomialsAtPoint(x_array_read[d * total_num_points + p], Q_1d, chebyshev_x));
4810b31fde2SJeremy L Thompson               CeedCall(CeedTensorContractApply(basis->contract, pre, Q_1d, post, 1, chebyshev_x, t_mode, false,
4820b31fde2SJeremy L Thompson                                                d == 0 ? chebyshev_coeffs : tmp[d % 2], tmp[(d + 1) % 2]));
4830b31fde2SJeremy L Thompson               pre /= Q_1d;
4840b31fde2SJeremy L Thompson               post *= 1;
4850b31fde2SJeremy L Thompson             }
4860b31fde2SJeremy L Thompson             for (CeedInt c = 0; c < num_comp; c++) v_array[c * total_num_points + p] = tmp[dim % 2][c];
4870b31fde2SJeremy L Thompson           }
4880b31fde2SJeremy L Thompson           break;
4890b31fde2SJeremy L Thompson         }
4900b31fde2SJeremy L Thompson         case CEED_EVAL_GRAD: {
4910b31fde2SJeremy L Thompson           CeedScalar tmp[2][num_comp * CeedIntPow(Q_1d, dim)], chebyshev_x[Q_1d];
4920b31fde2SJeremy L Thompson 
4930b31fde2SJeremy L Thompson           // ---- Values at point
4940b31fde2SJeremy L Thompson           for (CeedInt p = 0; p < total_num_points; p++) {
4950b31fde2SJeremy L Thompson             // Dim**2 contractions, apply grad when pass == dim
4960b31fde2SJeremy L Thompson             for (CeedInt pass = 0; pass < dim; pass++) {
4970b31fde2SJeremy L Thompson               CeedInt pre = num_comp * CeedIntPow(Q_1d, dim - 1), post = 1;
4980b31fde2SJeremy L Thompson 
4990b31fde2SJeremy L Thompson               for (CeedInt d = 0; d < dim; d++) {
5000b31fde2SJeremy L Thompson                 // ------ Tensor contract with current Chebyshev polynomial values
5010b31fde2SJeremy L Thompson                 if (pass == d) CeedCall(CeedChebyshevDerivativeAtPoint(x_array_read[d * total_num_points + p], Q_1d, chebyshev_x));
5020b31fde2SJeremy L Thompson                 else CeedCall(CeedChebyshevPolynomialsAtPoint(x_array_read[d * total_num_points + p], Q_1d, chebyshev_x));
5030b31fde2SJeremy L Thompson                 CeedCall(CeedTensorContractApply(basis->contract, pre, Q_1d, post, 1, chebyshev_x, t_mode, false,
5040b31fde2SJeremy L Thompson                                                  d == 0 ? chebyshev_coeffs : tmp[d % 2], tmp[(d + 1) % 2]));
5050b31fde2SJeremy L Thompson                 pre /= Q_1d;
5060b31fde2SJeremy L Thompson                 post *= 1;
5070b31fde2SJeremy L Thompson               }
5080b31fde2SJeremy L Thompson               for (CeedInt c = 0; c < num_comp; c++) v_array[(pass * num_comp + c) * total_num_points + p] = tmp[dim % 2][c];
5090b31fde2SJeremy L Thompson             }
5100b31fde2SJeremy L Thompson           }
5110b31fde2SJeremy L Thompson           break;
5120b31fde2SJeremy L Thompson         }
5130b31fde2SJeremy L Thompson         default:
5140b31fde2SJeremy L Thompson           // Nothing to do, excluded above
5150b31fde2SJeremy L Thompson           break;
5160b31fde2SJeremy L Thompson       }
5170b31fde2SJeremy L Thompson       CeedCall(CeedVectorRestoreArrayRead(basis->vec_chebyshev, &chebyshev_coeffs));
5180b31fde2SJeremy L Thompson       CeedCall(CeedVectorRestoreArrayRead(x_ref, &x_array_read));
5190b31fde2SJeremy L Thompson       CeedCall(CeedVectorRestoreArray(v, &v_array));
5200b31fde2SJeremy L Thompson       break;
5210b31fde2SJeremy L Thompson     }
5220b31fde2SJeremy L Thompson     case CEED_TRANSPOSE: {
5230b31fde2SJeremy L Thompson       // Note: No switch on e_mode here because only CEED_EVAL_INTERP is supported at this time
5240b31fde2SJeremy L Thompson       // Arbitrary points to nodes
5250b31fde2SJeremy L Thompson       CeedScalar       *chebyshev_coeffs;
5260b31fde2SJeremy L Thompson       const CeedScalar *u_array, *x_array_read;
5270b31fde2SJeremy L Thompson 
5280b31fde2SJeremy L Thompson       // -- Transpose of evaluation of Chebyshev polynomials at arbitrary points
5290b31fde2SJeremy L Thompson       CeedCall(CeedVectorGetArrayWrite(basis->vec_chebyshev, CEED_MEM_HOST, &chebyshev_coeffs));
5300b31fde2SJeremy L Thompson       CeedCall(CeedVectorGetArrayRead(x_ref, CEED_MEM_HOST, &x_array_read));
5310b31fde2SJeremy L Thompson       CeedCall(CeedVectorGetArrayRead(u, CEED_MEM_HOST, &u_array));
5320b31fde2SJeremy L Thompson 
5330b31fde2SJeremy L Thompson       switch (eval_mode) {
5340b31fde2SJeremy L Thompson         case CEED_EVAL_INTERP: {
5350b31fde2SJeremy L Thompson           CeedScalar tmp[2][num_comp * CeedIntPow(Q_1d, dim)], chebyshev_x[Q_1d];
5360b31fde2SJeremy L Thompson 
5370b31fde2SJeremy L Thompson           // ---- Values at point
5380b31fde2SJeremy L Thompson           for (CeedInt p = 0; p < total_num_points; p++) {
5390b31fde2SJeremy L Thompson             CeedInt pre = num_comp * 1, post = 1;
5400b31fde2SJeremy L Thompson 
5410b31fde2SJeremy L Thompson             for (CeedInt c = 0; c < num_comp; c++) tmp[0][c] = u_array[c * total_num_points + p];
5420b31fde2SJeremy L Thompson             for (CeedInt d = 0; d < dim; d++) {
5430b31fde2SJeremy L Thompson               // ------ Tensor contract with current Chebyshev polynomial values
5440b31fde2SJeremy L Thompson               CeedCall(CeedChebyshevPolynomialsAtPoint(x_array_read[d * total_num_points + p], Q_1d, chebyshev_x));
5450b31fde2SJeremy L Thompson               CeedCall(CeedTensorContractApply(basis->contract, pre, 1, post, Q_1d, chebyshev_x, t_mode, p > 0 && d == (dim - 1), tmp[d % 2],
5460b31fde2SJeremy L Thompson                                                d == (dim - 1) ? chebyshev_coeffs : tmp[(d + 1) % 2]));
5470b31fde2SJeremy L Thompson               pre /= 1;
5480b31fde2SJeremy L Thompson               post *= Q_1d;
5490b31fde2SJeremy L Thompson             }
5500b31fde2SJeremy L Thompson           }
5510b31fde2SJeremy L Thompson           break;
5520b31fde2SJeremy L Thompson         }
5530b31fde2SJeremy L Thompson         case CEED_EVAL_GRAD: {
5540b31fde2SJeremy L Thompson           CeedScalar tmp[2][num_comp * CeedIntPow(Q_1d, dim)], chebyshev_x[Q_1d];
5550b31fde2SJeremy L Thompson 
5560b31fde2SJeremy L Thompson           // ---- Values at point
5570b31fde2SJeremy L Thompson           for (CeedInt p = 0; p < total_num_points; p++) {
5580b31fde2SJeremy L Thompson             // Dim**2 contractions, apply grad when pass == dim
5590b31fde2SJeremy L Thompson             for (CeedInt pass = 0; pass < dim; pass++) {
5600b31fde2SJeremy L Thompson               CeedInt pre = num_comp * 1, post = 1;
5610b31fde2SJeremy L Thompson 
5620b31fde2SJeremy L Thompson               for (CeedInt c = 0; c < num_comp; c++) tmp[0][c] = u_array[(pass * num_comp + c) * total_num_points + p];
5630b31fde2SJeremy L Thompson               for (CeedInt d = 0; d < dim; d++) {
5640b31fde2SJeremy L Thompson                 // ------ Tensor contract with current Chebyshev polynomial values
5650b31fde2SJeremy L Thompson                 if (pass == d) CeedCall(CeedChebyshevDerivativeAtPoint(x_array_read[d * total_num_points + p], Q_1d, chebyshev_x));
5660b31fde2SJeremy L Thompson                 else CeedCall(CeedChebyshevPolynomialsAtPoint(x_array_read[d * total_num_points + p], Q_1d, chebyshev_x));
5670b31fde2SJeremy L Thompson                 CeedCall(CeedTensorContractApply(basis->contract, pre, 1, post, Q_1d, chebyshev_x, t_mode,
5680b31fde2SJeremy L Thompson                                                  (p > 0 || (p == 0 && pass > 0)) && d == (dim - 1), tmp[d % 2],
5690b31fde2SJeremy L Thompson                                                  d == (dim - 1) ? chebyshev_coeffs : tmp[(d + 1) % 2]));
5700b31fde2SJeremy L Thompson                 pre /= 1;
5710b31fde2SJeremy L Thompson                 post *= Q_1d;
5720b31fde2SJeremy L Thompson               }
5730b31fde2SJeremy L Thompson             }
5740b31fde2SJeremy L Thompson           }
5750b31fde2SJeremy L Thompson           break;
5760b31fde2SJeremy L Thompson         }
5770b31fde2SJeremy L Thompson         default:
5780b31fde2SJeremy L Thompson           // Nothing to do, excluded above
5790b31fde2SJeremy L Thompson           break;
5800b31fde2SJeremy L Thompson       }
5810b31fde2SJeremy L Thompson       CeedCall(CeedVectorRestoreArray(basis->vec_chebyshev, &chebyshev_coeffs));
5820b31fde2SJeremy L Thompson       CeedCall(CeedVectorRestoreArrayRead(x_ref, &x_array_read));
5830b31fde2SJeremy L Thompson       CeedCall(CeedVectorRestoreArrayRead(u, &u_array));
5840b31fde2SJeremy L Thompson 
5850b31fde2SJeremy L Thompson       // -- Interpolate transpose from Chebyshev coefficients
5860b31fde2SJeremy L Thompson       if (apply_add) CeedCall(CeedBasisApplyAdd(basis->basis_chebyshev, 1, CEED_TRANSPOSE, CEED_EVAL_INTERP, basis->vec_chebyshev, v));
5870b31fde2SJeremy L Thompson       else CeedCall(CeedBasisApply(basis->basis_chebyshev, 1, CEED_TRANSPOSE, CEED_EVAL_INTERP, basis->vec_chebyshev, v));
5880b31fde2SJeremy L Thompson       break;
5890b31fde2SJeremy L Thompson     }
5900b31fde2SJeremy L Thompson   }
5910b31fde2SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5920b31fde2SJeremy L Thompson }
5930b31fde2SJeremy L Thompson 
5947a982d89SJeremy L. Thompson /// @}
5957a982d89SJeremy L. Thompson 
5967a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
5977a982d89SJeremy L. Thompson /// Ceed Backend API
5987a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
5997a982d89SJeremy L. Thompson /// @addtogroup CeedBasisBackend
6007a982d89SJeremy L. Thompson /// @{
6017a982d89SJeremy L. Thompson 
6027a982d89SJeremy L. Thompson /**
603ca94c3ddSJeremy L Thompson   @brief Return collocated gradient matrix
6047a982d89SJeremy L. Thompson 
605ca94c3ddSJeremy L Thompson   @param[in]  basis         `CeedBasis`
606ca94c3ddSJeremy L Thompson   @param[out] collo_grad_1d Row-major (`Q_1d * Q_1d`) matrix expressing derivatives of basis functions at quadrature points
6077a982d89SJeremy L. Thompson 
6087a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
6097a982d89SJeremy L. Thompson 
6107a982d89SJeremy L. Thompson   @ref Backend
6117a982d89SJeremy L. Thompson **/
612d1d35e2fSjeremylt int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *collo_grad_1d) {
6137a982d89SJeremy L. Thompson   Ceed              ceed;
6142247a93fSRezgar Shakeri   CeedInt           P_1d, Q_1d;
6152247a93fSRezgar Shakeri   CeedScalar       *interp_1d_pinv;
6161203703bSJeremy L Thompson   const CeedScalar *grad_1d, *interp_1d;
6171203703bSJeremy L Thompson 
618ea61e9acSJeremy L Thompson   // Note: This function is for backend use, so all errors are terminal and we do not need to clean up memory on failure.
6192247a93fSRezgar Shakeri   CeedCall(CeedBasisGetCeed(basis, &ceed));
6202247a93fSRezgar Shakeri   CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d));
6212247a93fSRezgar Shakeri   CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
6227a982d89SJeremy L. Thompson 
6232247a93fSRezgar Shakeri   // Compute interp_1d^+, pseudoinverse of interp_1d
6242247a93fSRezgar Shakeri   CeedCall(CeedCalloc(P_1d * Q_1d, &interp_1d_pinv));
6251203703bSJeremy L Thompson   CeedCall(CeedBasisGetInterp1D(basis, &interp_1d));
6261203703bSJeremy L Thompson   CeedCall(CeedMatrixPseudoinverse(ceed, interp_1d, Q_1d, P_1d, interp_1d_pinv));
6271203703bSJeremy L Thompson   CeedCall(CeedBasisGetGrad1D(basis, &grad_1d));
6281203703bSJeremy L Thompson   CeedCall(CeedMatrixMatrixMultiply(ceed, grad_1d, (const CeedScalar *)interp_1d_pinv, collo_grad_1d, Q_1d, Q_1d, P_1d));
6297a982d89SJeremy L. Thompson 
6302247a93fSRezgar Shakeri   CeedCall(CeedFree(&interp_1d_pinv));
631e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
6327a982d89SJeremy L. Thompson }
6337a982d89SJeremy L. Thompson 
6347a982d89SJeremy L. Thompson /**
635b0cc4569SJeremy L Thompson   @brief Return 1D interpolation matrix to Chebyshev polynomial coefficients on quadrature space
636b0cc4569SJeremy L Thompson 
637b0cc4569SJeremy L Thompson   @param[in]  basis               `CeedBasis`
638b0cc4569SJeremy L Thompson   @param[out] chebyshev_interp_1d Row-major (`P_1d * Q_1d`) matrix interpolating from basis nodes to Chebyshev polynomial coefficients
639b0cc4569SJeremy L Thompson 
640b0cc4569SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
641b0cc4569SJeremy L Thompson 
642b0cc4569SJeremy L Thompson   @ref Backend
643b0cc4569SJeremy L Thompson **/
644b0cc4569SJeremy L Thompson int CeedBasisGetChebyshevInterp1D(CeedBasis basis, CeedScalar *chebyshev_interp_1d) {
645b0cc4569SJeremy L Thompson   CeedInt           P_1d, Q_1d;
646b0cc4569SJeremy L Thompson   CeedScalar       *C, *chebyshev_coeffs_1d_inv;
647b0cc4569SJeremy L Thompson   const CeedScalar *interp_1d, *q_ref_1d;
648b0cc4569SJeremy L Thompson   Ceed              ceed;
649b0cc4569SJeremy L Thompson 
650b0cc4569SJeremy L Thompson   CeedCall(CeedBasisGetCeed(basis, &ceed));
651b0cc4569SJeremy L Thompson   CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d));
652b0cc4569SJeremy L Thompson   CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
653b0cc4569SJeremy L Thompson 
654b0cc4569SJeremy L Thompson   // Build coefficient matrix
655bd83cbc5SJeremy L Thompson   // -- Note: Clang-tidy needs this check
656bd83cbc5SJeremy L Thompson   CeedCheck((P_1d > 0) && (Q_1d > 0), ceed, CEED_ERROR_INCOMPATIBLE, "CeedBasis dimensions are malformed");
657b0cc4569SJeremy L Thompson   CeedCall(CeedCalloc(Q_1d * Q_1d, &C));
658b0cc4569SJeremy L Thompson   CeedCall(CeedBasisGetQRef(basis, &q_ref_1d));
659b0cc4569SJeremy L Thompson   for (CeedInt i = 0; i < Q_1d; i++) CeedCall(CeedChebyshevPolynomialsAtPoint(q_ref_1d[i], Q_1d, &C[i * Q_1d]));
660b0cc4569SJeremy L Thompson 
661b0cc4569SJeremy L Thompson   // Compute C^+, pseudoinverse of coefficient matrix
662b0cc4569SJeremy L Thompson   CeedCall(CeedCalloc(Q_1d * Q_1d, &chebyshev_coeffs_1d_inv));
663b0cc4569SJeremy L Thompson   CeedCall(CeedMatrixPseudoinverse(ceed, C, Q_1d, Q_1d, chebyshev_coeffs_1d_inv));
664b0cc4569SJeremy L Thompson 
665b0cc4569SJeremy L Thompson   // Build mapping from nodes to Chebyshev coefficients
666b0cc4569SJeremy L Thompson   CeedCall(CeedBasisGetInterp1D(basis, &interp_1d));
667b0cc4569SJeremy L Thompson   CeedCall(CeedMatrixMatrixMultiply(ceed, chebyshev_coeffs_1d_inv, interp_1d, chebyshev_interp_1d, Q_1d, P_1d, Q_1d));
668b0cc4569SJeremy L Thompson 
669b0cc4569SJeremy L Thompson   // Cleanup
670b0cc4569SJeremy L Thompson   CeedCall(CeedFree(&C));
671b0cc4569SJeremy L Thompson   CeedCall(CeedFree(&chebyshev_coeffs_1d_inv));
672b0cc4569SJeremy L Thompson   return CEED_ERROR_SUCCESS;
673b0cc4569SJeremy L Thompson }
674b0cc4569SJeremy L Thompson 
675b0cc4569SJeremy L Thompson /**
676ca94c3ddSJeremy L Thompson   @brief Get tensor status for given `CeedBasis`
6777a982d89SJeremy L. Thompson 
678ca94c3ddSJeremy L Thompson   @param[in]  basis     `CeedBasis`
679d1d35e2fSjeremylt   @param[out] is_tensor Variable to store tensor status
6807a982d89SJeremy L. Thompson 
6817a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
6827a982d89SJeremy L. Thompson 
6837a982d89SJeremy L. Thompson   @ref Backend
6847a982d89SJeremy L. Thompson **/
685d1d35e2fSjeremylt int CeedBasisIsTensor(CeedBasis basis, bool *is_tensor) {
6866402da51SJeremy L Thompson   *is_tensor = basis->is_tensor_basis;
687e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
6887a982d89SJeremy L. Thompson }
6897a982d89SJeremy L. Thompson 
6907a982d89SJeremy L. Thompson /**
691ca94c3ddSJeremy L Thompson   @brief Get backend data of a `CeedBasis`
6927a982d89SJeremy L. Thompson 
693ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
6947a982d89SJeremy L. Thompson   @param[out] data  Variable to store data
6957a982d89SJeremy L. Thompson 
6967a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
6977a982d89SJeremy L. Thompson 
6987a982d89SJeremy L. Thompson   @ref Backend
6997a982d89SJeremy L. Thompson **/
700777ff853SJeremy L Thompson int CeedBasisGetData(CeedBasis basis, void *data) {
701777ff853SJeremy L Thompson   *(void **)data = basis->data;
702e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
7037a982d89SJeremy L. Thompson }
7047a982d89SJeremy L. Thompson 
7057a982d89SJeremy L. Thompson /**
706ca94c3ddSJeremy L Thompson   @brief Set backend data of a `CeedBasis`
7077a982d89SJeremy L. Thompson 
708ca94c3ddSJeremy L Thompson   @param[in,out] basis  `CeedBasis`
709ea61e9acSJeremy L Thompson   @param[in]     data   Data to set
7107a982d89SJeremy L. Thompson 
7117a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
7127a982d89SJeremy L. Thompson 
7137a982d89SJeremy L. Thompson   @ref Backend
7147a982d89SJeremy L. Thompson **/
715777ff853SJeremy L Thompson int CeedBasisSetData(CeedBasis basis, void *data) {
716777ff853SJeremy L Thompson   basis->data = data;
717e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
7187a982d89SJeremy L. Thompson }
7197a982d89SJeremy L. Thompson 
7207a982d89SJeremy L. Thompson /**
721ca94c3ddSJeremy L Thompson   @brief Increment the reference counter for a `CeedBasis`
72234359f16Sjeremylt 
723ca94c3ddSJeremy L Thompson   @param[in,out] basis `CeedBasis` to increment the reference counter
72434359f16Sjeremylt 
72534359f16Sjeremylt   @return An error code: 0 - success, otherwise - failure
72634359f16Sjeremylt 
72734359f16Sjeremylt   @ref Backend
72834359f16Sjeremylt **/
7299560d06aSjeremylt int CeedBasisReference(CeedBasis basis) {
73034359f16Sjeremylt   basis->ref_count++;
73134359f16Sjeremylt   return CEED_ERROR_SUCCESS;
73234359f16Sjeremylt }
73334359f16Sjeremylt 
73434359f16Sjeremylt /**
735ca94c3ddSJeremy L Thompson   @brief Get number of Q-vector components for given `CeedBasis`
736c4e3f59bSSebastian Grimberg 
737ca94c3ddSJeremy L Thompson   @param[in]  basis     `CeedBasis`
738ca94c3ddSJeremy L Thompson   @param[in]  eval_mode @ref CEED_EVAL_INTERP to use interpolated values,
739ca94c3ddSJeremy L Thompson                           @ref CEED_EVAL_GRAD to use gradients,
740ca94c3ddSJeremy L Thompson                           @ref CEED_EVAL_DIV to use divergence,
741ca94c3ddSJeremy L Thompson                           @ref CEED_EVAL_CURL to use curl
742c4e3f59bSSebastian Grimberg   @param[out] q_comp    Variable to store number of Q-vector components of basis
743c4e3f59bSSebastian Grimberg 
744c4e3f59bSSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
745c4e3f59bSSebastian Grimberg 
746c4e3f59bSSebastian Grimberg   @ref Backend
747c4e3f59bSSebastian Grimberg **/
748c4e3f59bSSebastian Grimberg int CeedBasisGetNumQuadratureComponents(CeedBasis basis, CeedEvalMode eval_mode, CeedInt *q_comp) {
7491203703bSJeremy L Thompson   CeedInt dim;
7501203703bSJeremy L Thompson 
7511203703bSJeremy L Thompson   CeedCall(CeedBasisGetDimension(basis, &dim));
752c4e3f59bSSebastian Grimberg   switch (eval_mode) {
7531203703bSJeremy L Thompson     case CEED_EVAL_INTERP: {
7541203703bSJeremy L Thompson       CeedFESpace fe_space;
7551203703bSJeremy L Thompson 
7561203703bSJeremy L Thompson       CeedCall(CeedBasisGetFESpace(basis, &fe_space));
7571203703bSJeremy L Thompson       *q_comp = (fe_space == CEED_FE_SPACE_H1) ? 1 : dim;
7581203703bSJeremy L Thompson     } break;
759c4e3f59bSSebastian Grimberg     case CEED_EVAL_GRAD:
7601203703bSJeremy L Thompson       *q_comp = dim;
761c4e3f59bSSebastian Grimberg       break;
762c4e3f59bSSebastian Grimberg     case CEED_EVAL_DIV:
763c4e3f59bSSebastian Grimberg       *q_comp = 1;
764c4e3f59bSSebastian Grimberg       break;
765c4e3f59bSSebastian Grimberg     case CEED_EVAL_CURL:
7661203703bSJeremy L Thompson       *q_comp = (dim < 3) ? 1 : dim;
767c4e3f59bSSebastian Grimberg       break;
768c4e3f59bSSebastian Grimberg     case CEED_EVAL_NONE:
769c4e3f59bSSebastian Grimberg     case CEED_EVAL_WEIGHT:
770352a5e7cSSebastian Grimberg       *q_comp = 1;
771c4e3f59bSSebastian Grimberg       break;
772c4e3f59bSSebastian Grimberg   }
773c4e3f59bSSebastian Grimberg   return CEED_ERROR_SUCCESS;
774c4e3f59bSSebastian Grimberg }
775c4e3f59bSSebastian Grimberg 
776c4e3f59bSSebastian Grimberg /**
777ca94c3ddSJeremy L Thompson   @brief Estimate number of FLOPs required to apply `CeedBasis` in `t_mode` and `eval_mode`
7786e15d496SJeremy L Thompson 
779ca94c3ddSJeremy L Thompson   @param[in]  basis     `CeedBasis` to estimate FLOPs for
780ea61e9acSJeremy L Thompson   @param[in]  t_mode    Apply basis or transpose
781ca94c3ddSJeremy L Thompson   @param[in]  eval_mode @ref CeedEvalMode
782ea61e9acSJeremy L Thompson   @param[out] flops     Address of variable to hold FLOPs estimate
7836e15d496SJeremy L Thompson 
7846e15d496SJeremy L Thompson   @ref Backend
7856e15d496SJeremy L Thompson **/
7862b730f8bSJeremy L Thompson int CeedBasisGetFlopsEstimate(CeedBasis basis, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedSize *flops) {
7876e15d496SJeremy L Thompson   bool is_tensor;
7886e15d496SJeremy L Thompson 
7892b730f8bSJeremy L Thompson   CeedCall(CeedBasisIsTensor(basis, &is_tensor));
7906e15d496SJeremy L Thompson   if (is_tensor) {
7916e15d496SJeremy L Thompson     CeedInt dim, num_comp, P_1d, Q_1d;
7921c66c397SJeremy L Thompson 
7932b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetDimension(basis, &dim));
7942b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
7952b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d));
7962b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
7976e15d496SJeremy L Thompson     if (t_mode == CEED_TRANSPOSE) {
7982b730f8bSJeremy L Thompson       P_1d = Q_1d;
7992b730f8bSJeremy L Thompson       Q_1d = P_1d;
8006e15d496SJeremy L Thompson     }
8016e15d496SJeremy L Thompson     CeedInt tensor_flops = 0, pre = num_comp * CeedIntPow(P_1d, dim - 1), post = 1;
8026e15d496SJeremy L Thompson     for (CeedInt d = 0; d < dim; d++) {
8036e15d496SJeremy L Thompson       tensor_flops += 2 * pre * P_1d * post * Q_1d;
8046e15d496SJeremy L Thompson       pre /= P_1d;
8056e15d496SJeremy L Thompson       post *= Q_1d;
8066e15d496SJeremy L Thompson     }
8076e15d496SJeremy L Thompson     switch (eval_mode) {
8082b730f8bSJeremy L Thompson       case CEED_EVAL_NONE:
8092b730f8bSJeremy L Thompson         *flops = 0;
8102b730f8bSJeremy L Thompson         break;
8112b730f8bSJeremy L Thompson       case CEED_EVAL_INTERP:
8122b730f8bSJeremy L Thompson         *flops = tensor_flops;
8132b730f8bSJeremy L Thompson         break;
8142b730f8bSJeremy L Thompson       case CEED_EVAL_GRAD:
8152b730f8bSJeremy L Thompson         *flops = tensor_flops * 2;
8162b730f8bSJeremy L Thompson         break;
8176e15d496SJeremy L Thompson       case CEED_EVAL_DIV:
8181203703bSJeremy L Thompson       case CEED_EVAL_CURL: {
8196574a04fSJeremy L Thompson         // LCOV_EXCL_START
8206e536b99SJeremy L Thompson         return CeedError(CeedBasisReturnCeed(basis), CEED_ERROR_INCOMPATIBLE, "Tensor basis evaluation for %s not supported",
8216e536b99SJeremy L Thompson                          CeedEvalModes[eval_mode]);
8222b730f8bSJeremy L Thompson         break;
8236e15d496SJeremy L Thompson         // LCOV_EXCL_STOP
8241203703bSJeremy L Thompson       }
8252b730f8bSJeremy L Thompson       case CEED_EVAL_WEIGHT:
8262b730f8bSJeremy L Thompson         *flops = dim * CeedIntPow(Q_1d, dim);
8272b730f8bSJeremy L Thompson         break;
8286e15d496SJeremy L Thompson     }
8296e15d496SJeremy L Thompson   } else {
830c4e3f59bSSebastian Grimberg     CeedInt dim, num_comp, q_comp, num_nodes, num_qpts;
8311c66c397SJeremy L Thompson 
8322b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetDimension(basis, &dim));
8332b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
834c4e3f59bSSebastian Grimberg     CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp));
8352b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis, &num_nodes));
8362b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts));
8376e15d496SJeremy L Thompson     switch (eval_mode) {
8382b730f8bSJeremy L Thompson       case CEED_EVAL_NONE:
8392b730f8bSJeremy L Thompson         *flops = 0;
8402b730f8bSJeremy L Thompson         break;
8412b730f8bSJeremy L Thompson       case CEED_EVAL_INTERP:
8422b730f8bSJeremy L Thompson       case CEED_EVAL_GRAD:
8432b730f8bSJeremy L Thompson       case CEED_EVAL_DIV:
8442b730f8bSJeremy L Thompson       case CEED_EVAL_CURL:
845c4e3f59bSSebastian Grimberg         *flops = num_nodes * num_qpts * num_comp * q_comp;
8462b730f8bSJeremy L Thompson         break;
8472b730f8bSJeremy L Thompson       case CEED_EVAL_WEIGHT:
8482b730f8bSJeremy L Thompson         *flops = 0;
8492b730f8bSJeremy L Thompson         break;
8506e15d496SJeremy L Thompson     }
8516e15d496SJeremy L Thompson   }
8526e15d496SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8536e15d496SJeremy L Thompson }
8546e15d496SJeremy L Thompson 
8556e15d496SJeremy L Thompson /**
856ca94c3ddSJeremy L Thompson   @brief Get `CeedFESpace` for a `CeedBasis`
857c4e3f59bSSebastian Grimberg 
858ca94c3ddSJeremy L Thompson   @param[in]  basis    `CeedBasis`
859ca94c3ddSJeremy L Thompson   @param[out] fe_space Variable to store `CeedFESpace`
860c4e3f59bSSebastian Grimberg 
861c4e3f59bSSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
862c4e3f59bSSebastian Grimberg 
863c4e3f59bSSebastian Grimberg   @ref Backend
864c4e3f59bSSebastian Grimberg **/
865c4e3f59bSSebastian Grimberg int CeedBasisGetFESpace(CeedBasis basis, CeedFESpace *fe_space) {
866c4e3f59bSSebastian Grimberg   *fe_space = basis->fe_space;
867c4e3f59bSSebastian Grimberg   return CEED_ERROR_SUCCESS;
868c4e3f59bSSebastian Grimberg }
869c4e3f59bSSebastian Grimberg 
870c4e3f59bSSebastian Grimberg /**
871ca94c3ddSJeremy L Thompson   @brief Get dimension for given `CeedElemTopology`
8727a982d89SJeremy L. Thompson 
873ca94c3ddSJeremy L Thompson   @param[in]  topo `CeedElemTopology`
8747a982d89SJeremy L. Thompson   @param[out] dim  Variable to store dimension of topology
8757a982d89SJeremy L. Thompson 
8767a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
8777a982d89SJeremy L. Thompson 
8787a982d89SJeremy L. Thompson   @ref Backend
8797a982d89SJeremy L. Thompson **/
8807a982d89SJeremy L. Thompson int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) {
8817a982d89SJeremy L. Thompson   *dim = (CeedInt)topo >> 16;
882e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8837a982d89SJeremy L. Thompson }
8847a982d89SJeremy L. Thompson 
8857a982d89SJeremy L. Thompson /**
886ca94c3ddSJeremy L Thompson   @brief Get `CeedTensorContract` of a `CeedBasis`
8877a982d89SJeremy L. Thompson 
888ca94c3ddSJeremy L Thompson   @param[in]  basis     `CeedBasis`
889ca94c3ddSJeremy L Thompson   @param[out] contract  Variable to store `CeedTensorContract`
8907a982d89SJeremy L. Thompson 
8917a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
8927a982d89SJeremy L. Thompson 
8937a982d89SJeremy L. Thompson   @ref Backend
8947a982d89SJeremy L. Thompson **/
8957a982d89SJeremy L. Thompson int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) {
8967a982d89SJeremy L. Thompson   *contract = basis->contract;
897e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8987a982d89SJeremy L. Thompson }
8997a982d89SJeremy L. Thompson 
9007a982d89SJeremy L. Thompson /**
901ca94c3ddSJeremy L Thompson   @brief Set `CeedTensorContract` of a `CeedBasis`
9027a982d89SJeremy L. Thompson 
903ca94c3ddSJeremy L Thompson   @param[in,out] basis    `CeedBasis`
904ca94c3ddSJeremy L Thompson   @param[in]     contract `CeedTensorContract` to set
9057a982d89SJeremy L. Thompson 
9067a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
9077a982d89SJeremy L. Thompson 
9087a982d89SJeremy L. Thompson   @ref Backend
9097a982d89SJeremy L. Thompson **/
91034359f16Sjeremylt int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract contract) {
91134359f16Sjeremylt   basis->contract = contract;
9122b730f8bSJeremy L Thompson   CeedCall(CeedTensorContractReference(contract));
913e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9147a982d89SJeremy L. Thompson }
9157a982d89SJeremy L. Thompson 
9167a982d89SJeremy L. Thompson /**
917ca94c3ddSJeremy L Thompson   @brief Return a reference implementation of matrix multiplication \f$C = A B\f$.
918ba59ac12SSebastian Grimberg 
919ca94c3ddSJeremy L Thompson   Note: This is a reference implementation for CPU `CeedScalar` pointers that is not intended for high performance.
9207a982d89SJeremy L. Thompson 
921ca94c3ddSJeremy L Thompson   @param[in]  ceed  `Ceed` context for error handling
922ca94c3ddSJeremy L Thompson   @param[in]  mat_A Row-major matrix `A`
923ca94c3ddSJeremy L Thompson   @param[in]  mat_B Row-major matrix `B`
924ca94c3ddSJeremy L Thompson   @param[out] mat_C Row-major output matrix `C`
925ca94c3ddSJeremy L Thompson   @param[in]  m     Number of rows of `C`
926ca94c3ddSJeremy L Thompson   @param[in]  n     Number of columns of `C`
927ca94c3ddSJeremy L Thompson   @param[in]  kk    Number of columns of `A`/rows of `B`
9287a982d89SJeremy L. Thompson 
9297a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
9307a982d89SJeremy L. Thompson 
9317a982d89SJeremy L. Thompson   @ref Utility
9327a982d89SJeremy L. Thompson **/
9332b730f8bSJeremy L Thompson int CeedMatrixMatrixMultiply(Ceed ceed, const CeedScalar *mat_A, const CeedScalar *mat_B, CeedScalar *mat_C, CeedInt m, CeedInt n, CeedInt kk) {
9342b730f8bSJeremy L Thompson   for (CeedInt i = 0; i < m; i++) {
9357a982d89SJeremy L. Thompson     for (CeedInt j = 0; j < n; j++) {
9367a982d89SJeremy L. Thompson       CeedScalar sum = 0;
9371c66c397SJeremy L Thompson 
9382b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < kk; k++) sum += mat_A[k + i * kk] * mat_B[j + k * n];
939d1d35e2fSjeremylt       mat_C[j + i * n] = sum;
9407a982d89SJeremy L. Thompson     }
9412b730f8bSJeremy L Thompson   }
942e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9437a982d89SJeremy L. Thompson }
9447a982d89SJeremy L. Thompson 
945ba59ac12SSebastian Grimberg /**
946ba59ac12SSebastian Grimberg   @brief Return QR Factorization of a matrix
947ba59ac12SSebastian Grimberg 
948ca94c3ddSJeremy L Thompson   @param[in]     ceed `Ceed` context for error handling
949ba59ac12SSebastian Grimberg   @param[in,out] mat  Row-major matrix to be factorized in place
950ca94c3ddSJeremy L Thompson   @param[in,out] tau  Vector of length `m` of scaling factors
951ba59ac12SSebastian Grimberg   @param[in]     m    Number of rows
952ba59ac12SSebastian Grimberg   @param[in]     n    Number of columns
953ba59ac12SSebastian Grimberg 
954ba59ac12SSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
955ba59ac12SSebastian Grimberg 
956ba59ac12SSebastian Grimberg   @ref Utility
957ba59ac12SSebastian Grimberg **/
958ba59ac12SSebastian Grimberg int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau, CeedInt m, CeedInt n) {
959ba59ac12SSebastian Grimberg   CeedScalar v[m];
960ba59ac12SSebastian Grimberg 
961ba59ac12SSebastian Grimberg   // Check matrix shape
9626574a04fSJeremy L Thompson   CeedCheck(n <= m, ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute QR factorization with n > m");
963ba59ac12SSebastian Grimberg 
964ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) {
9651c66c397SJeremy L Thompson     CeedScalar sigma = 0.0;
9661c66c397SJeremy L Thompson 
967ba59ac12SSebastian Grimberg     if (i >= m - 1) {  // last row of matrix, no reflection needed
968ba59ac12SSebastian Grimberg       tau[i] = 0.;
969ba59ac12SSebastian Grimberg       break;
970ba59ac12SSebastian Grimberg     }
971ba59ac12SSebastian Grimberg     // Calculate Householder vector, magnitude
972ba59ac12SSebastian Grimberg     v[i] = mat[i + n * i];
973ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < m; j++) {
974ba59ac12SSebastian Grimberg       v[j] = mat[i + n * j];
975ba59ac12SSebastian Grimberg       sigma += v[j] * v[j];
976ba59ac12SSebastian Grimberg     }
9771c66c397SJeremy L Thompson     const CeedScalar norm = sqrt(v[i] * v[i] + sigma);  // norm of v[i:m]
9781c66c397SJeremy L Thompson     const CeedScalar R_ii = -copysign(norm, v[i]);
9791c66c397SJeremy L Thompson 
980ba59ac12SSebastian Grimberg     v[i] -= R_ii;
981ba59ac12SSebastian Grimberg     // norm of v[i:m] after modification above and scaling below
982ba59ac12SSebastian Grimberg     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
983ba59ac12SSebastian Grimberg     //   tau = 2 / (norm*norm)
984ba59ac12SSebastian Grimberg     tau[i] = 2 * v[i] * v[i] / (v[i] * v[i] + sigma);
985ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < m; j++) v[j] /= v[i];
986ba59ac12SSebastian Grimberg 
987ba59ac12SSebastian Grimberg     // Apply Householder reflector to lower right panel
988ba59ac12SSebastian Grimberg     CeedHouseholderReflect(&mat[i * n + i + 1], &v[i], tau[i], m - i, n - i - 1, n, 1);
989ba59ac12SSebastian Grimberg     // Save v
990ba59ac12SSebastian Grimberg     mat[i + n * i] = R_ii;
991ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < m; j++) mat[i + n * j] = v[j];
992ba59ac12SSebastian Grimberg   }
993ba59ac12SSebastian Grimberg   return CEED_ERROR_SUCCESS;
994ba59ac12SSebastian Grimberg }
995ba59ac12SSebastian Grimberg 
996ba59ac12SSebastian Grimberg /**
997ba59ac12SSebastian Grimberg   @brief Apply Householder Q matrix
998ba59ac12SSebastian Grimberg 
999ca94c3ddSJeremy L Thompson   Compute `mat_A = mat_Q mat_A`, where `mat_Q` is \f$m \times m\f$ and `mat_A` is \f$m \times n\f$.
1000ba59ac12SSebastian Grimberg 
1001ba59ac12SSebastian Grimberg   @param[in,out] mat_A  Matrix to apply Householder Q to, in place
1002ba59ac12SSebastian Grimberg   @param[in]     mat_Q  Householder Q matrix
1003ba59ac12SSebastian Grimberg   @param[in]     tau    Householder scaling factors
1004ba59ac12SSebastian Grimberg   @param[in]     t_mode Transpose mode for application
1005ca94c3ddSJeremy L Thompson   @param[in]     m      Number of rows in `A`
1006ca94c3ddSJeremy L Thompson   @param[in]     n      Number of columns in `A`
1007ca94c3ddSJeremy L Thompson   @param[in]     k      Number of elementary reflectors in Q, `k < m`
1008ca94c3ddSJeremy L Thompson   @param[in]     row    Row stride in `A`
1009ca94c3ddSJeremy L Thompson   @param[in]     col    Col stride in `A`
1010ba59ac12SSebastian Grimberg 
1011ba59ac12SSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
1012ba59ac12SSebastian Grimberg 
1013c4e3f59bSSebastian Grimberg   @ref Utility
1014ba59ac12SSebastian Grimberg **/
1015ba59ac12SSebastian Grimberg int CeedHouseholderApplyQ(CeedScalar *mat_A, const CeedScalar *mat_Q, const CeedScalar *tau, CeedTransposeMode t_mode, CeedInt m, CeedInt n,
1016ba59ac12SSebastian Grimberg                           CeedInt k, CeedInt row, CeedInt col) {
1017ba59ac12SSebastian Grimberg   CeedScalar *v;
10181c66c397SJeremy L Thompson 
1019ba59ac12SSebastian Grimberg   CeedCall(CeedMalloc(m, &v));
1020ba59ac12SSebastian Grimberg   for (CeedInt ii = 0; ii < k; ii++) {
1021ba59ac12SSebastian Grimberg     CeedInt i = t_mode == CEED_TRANSPOSE ? ii : k - 1 - ii;
1022ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < m; j++) v[j] = mat_Q[j * k + i];
1023ba59ac12SSebastian Grimberg     // Apply Householder reflector (I - tau v v^T) collo_grad_1d^T
1024ba59ac12SSebastian Grimberg     CeedCall(CeedHouseholderReflect(&mat_A[i * row], &v[i], tau[i], m - i, n, row, col));
1025ba59ac12SSebastian Grimberg   }
1026ba59ac12SSebastian Grimberg   CeedCall(CeedFree(&v));
1027ba59ac12SSebastian Grimberg   return CEED_ERROR_SUCCESS;
1028ba59ac12SSebastian Grimberg }
1029ba59ac12SSebastian Grimberg 
1030ba59ac12SSebastian Grimberg /**
10312247a93fSRezgar Shakeri   @brief Return pseudoinverse of a matrix
10322247a93fSRezgar Shakeri 
10332247a93fSRezgar Shakeri   @param[in]     ceed      Ceed context for error handling
10342247a93fSRezgar Shakeri   @param[in]     mat       Row-major matrix to compute pseudoinverse of
10352247a93fSRezgar Shakeri   @param[in]     m         Number of rows
10362247a93fSRezgar Shakeri   @param[in]     n         Number of columns
10372247a93fSRezgar Shakeri   @param[out]    mat_pinv  Row-major pseudoinverse matrix
10382247a93fSRezgar Shakeri 
10392247a93fSRezgar Shakeri   @return An error code: 0 - success, otherwise - failure
10402247a93fSRezgar Shakeri 
10412247a93fSRezgar Shakeri   @ref Utility
10422247a93fSRezgar Shakeri **/
10431203703bSJeremy L Thompson int CeedMatrixPseudoinverse(Ceed ceed, const CeedScalar *mat, CeedInt m, CeedInt n, CeedScalar *mat_pinv) {
10442247a93fSRezgar Shakeri   CeedScalar *tau, *I, *mat_copy;
10452247a93fSRezgar Shakeri 
10462247a93fSRezgar Shakeri   CeedCall(CeedCalloc(m, &tau));
10472247a93fSRezgar Shakeri   CeedCall(CeedCalloc(m * m, &I));
10482247a93fSRezgar Shakeri   CeedCall(CeedCalloc(m * n, &mat_copy));
10492247a93fSRezgar Shakeri   memcpy(mat_copy, mat, m * n * sizeof mat[0]);
10502247a93fSRezgar Shakeri 
10512247a93fSRezgar Shakeri   // QR Factorization, mat = Q R
10522247a93fSRezgar Shakeri   CeedCall(CeedQRFactorization(ceed, mat_copy, tau, m, n));
10532247a93fSRezgar Shakeri 
10542247a93fSRezgar Shakeri   // -- Apply Q^T, I = Q^T * I
10552247a93fSRezgar Shakeri   for (CeedInt i = 0; i < m; i++) I[i * m + i] = 1.0;
10562247a93fSRezgar Shakeri   CeedCall(CeedHouseholderApplyQ(I, mat_copy, tau, CEED_TRANSPOSE, m, m, n, m, 1));
10572247a93fSRezgar Shakeri   // -- Apply R_inv, mat_pinv = R_inv * Q^T
10582247a93fSRezgar Shakeri   for (CeedInt j = 0; j < m; j++) {  // Column j
10592247a93fSRezgar Shakeri     mat_pinv[j + m * (n - 1)] = I[j + m * (n - 1)] / mat_copy[n * n - 1];
10602247a93fSRezgar Shakeri     for (CeedInt i = n - 2; i >= 0; i--) {  // Row i
10612247a93fSRezgar Shakeri       mat_pinv[j + m * i] = I[j + m * i];
10622247a93fSRezgar Shakeri       for (CeedInt k = i + 1; k < n; k++) mat_pinv[j + m * i] -= mat_copy[k + n * i] * mat_pinv[j + m * k];
10632247a93fSRezgar Shakeri       mat_pinv[j + m * i] /= mat_copy[i + n * i];
10642247a93fSRezgar Shakeri     }
10652247a93fSRezgar Shakeri   }
10662247a93fSRezgar Shakeri 
10672247a93fSRezgar Shakeri   // Cleanup
10682247a93fSRezgar Shakeri   CeedCall(CeedFree(&I));
10692247a93fSRezgar Shakeri   CeedCall(CeedFree(&tau));
10702247a93fSRezgar Shakeri   CeedCall(CeedFree(&mat_copy));
10712247a93fSRezgar Shakeri   return CEED_ERROR_SUCCESS;
10722247a93fSRezgar Shakeri }
10732247a93fSRezgar Shakeri 
10742247a93fSRezgar Shakeri /**
1075ba59ac12SSebastian Grimberg   @brief Return symmetric Schur decomposition of the symmetric matrix mat via symmetric QR factorization
1076ba59ac12SSebastian Grimberg 
1077ca94c3ddSJeremy L Thompson   @param[in]     ceed   `Ceed` context for error handling
1078ba59ac12SSebastian Grimberg   @param[in,out] mat    Row-major matrix to be factorized in place
1079ba59ac12SSebastian Grimberg   @param[out]    lambda Vector of length n of eigenvalues
1080ba59ac12SSebastian Grimberg   @param[in]     n      Number of rows/columns
1081ba59ac12SSebastian Grimberg 
1082ba59ac12SSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
1083ba59ac12SSebastian Grimberg 
1084ba59ac12SSebastian Grimberg   @ref Utility
1085ba59ac12SSebastian Grimberg **/
10862c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOff
10872c2ea1dbSJeremy L Thompson int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat, CeedScalar *lambda, CeedInt n) {
1088ba59ac12SSebastian Grimberg   // Check bounds for clang-tidy
10896574a04fSJeremy L Thompson   CeedCheck(n > 1, ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute symmetric Schur decomposition of scalars");
1090ba59ac12SSebastian Grimberg 
1091ba59ac12SSebastian Grimberg   CeedScalar v[n - 1], tau[n - 1], mat_T[n * n];
1092ba59ac12SSebastian Grimberg 
1093ba59ac12SSebastian Grimberg   // Copy mat to mat_T and set mat to I
1094ba59ac12SSebastian Grimberg   memcpy(mat_T, mat, n * n * sizeof(mat[0]));
1095ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) {
1096ba59ac12SSebastian Grimberg     for (CeedInt j = 0; j < n; j++) mat[j + n * i] = (i == j) ? 1 : 0;
1097ba59ac12SSebastian Grimberg   }
1098ba59ac12SSebastian Grimberg 
1099ba59ac12SSebastian Grimberg   // Reduce to tridiagonal
1100ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n - 1; i++) {
1101ba59ac12SSebastian Grimberg     // Calculate Householder vector, magnitude
1102ba59ac12SSebastian Grimberg     CeedScalar sigma = 0.0;
11031c66c397SJeremy L Thompson 
1104ba59ac12SSebastian Grimberg     v[i] = mat_T[i + n * (i + 1)];
1105ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < n - 1; j++) {
1106ba59ac12SSebastian Grimberg       v[j] = mat_T[i + n * (j + 1)];
1107ba59ac12SSebastian Grimberg       sigma += v[j] * v[j];
1108ba59ac12SSebastian Grimberg     }
11091c66c397SJeremy L Thompson     const CeedScalar norm = sqrt(v[i] * v[i] + sigma);  // norm of v[i:n-1]
11101c66c397SJeremy L Thompson     const CeedScalar R_ii = -copysign(norm, v[i]);
11111c66c397SJeremy L Thompson 
1112ba59ac12SSebastian Grimberg     v[i] -= R_ii;
1113ba59ac12SSebastian Grimberg     // norm of v[i:m] after modification above and scaling below
1114ba59ac12SSebastian Grimberg     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
1115ba59ac12SSebastian Grimberg     //   tau = 2 / (norm*norm)
1116ba59ac12SSebastian Grimberg     tau[i] = i == n - 2 ? 2 : 2 * v[i] * v[i] / (v[i] * v[i] + sigma);
1117ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < n - 1; j++) v[j] /= v[i];
1118ba59ac12SSebastian Grimberg 
1119ba59ac12SSebastian Grimberg     // Update sub and super diagonal
1120ba59ac12SSebastian Grimberg     for (CeedInt j = i + 2; j < n; j++) {
1121ba59ac12SSebastian Grimberg       mat_T[i + n * j] = 0;
1122ba59ac12SSebastian Grimberg       mat_T[j + n * i] = 0;
1123ba59ac12SSebastian Grimberg     }
1124ba59ac12SSebastian Grimberg     // Apply symmetric Householder reflector to lower right panel
1125ba59ac12SSebastian Grimberg     CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1);
1126ba59ac12SSebastian Grimberg     CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), 1, n);
1127ba59ac12SSebastian Grimberg 
1128ba59ac12SSebastian Grimberg     // Save v
1129ba59ac12SSebastian Grimberg     mat_T[i + n * (i + 1)] = R_ii;
1130ba59ac12SSebastian Grimberg     mat_T[(i + 1) + n * i] = R_ii;
1131ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < n - 1; j++) {
1132ba59ac12SSebastian Grimberg       mat_T[i + n * (j + 1)] = v[j];
1133ba59ac12SSebastian Grimberg     }
1134ba59ac12SSebastian Grimberg   }
1135ba59ac12SSebastian Grimberg   // Backwards accumulation of Q
1136ba59ac12SSebastian Grimberg   for (CeedInt i = n - 2; i >= 0; i--) {
1137ba59ac12SSebastian Grimberg     if (tau[i] > 0.0) {
1138ba59ac12SSebastian Grimberg       v[i] = 1;
1139ba59ac12SSebastian Grimberg       for (CeedInt j = i + 1; j < n - 1; j++) {
1140ba59ac12SSebastian Grimberg         v[j]                   = mat_T[i + n * (j + 1)];
1141ba59ac12SSebastian Grimberg         mat_T[i + n * (j + 1)] = 0;
1142ba59ac12SSebastian Grimberg       }
1143ba59ac12SSebastian Grimberg       CeedHouseholderReflect(&mat[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1);
1144ba59ac12SSebastian Grimberg     }
1145ba59ac12SSebastian Grimberg   }
1146ba59ac12SSebastian Grimberg 
1147ba59ac12SSebastian Grimberg   // Reduce sub and super diagonal
1148ba59ac12SSebastian Grimberg   CeedInt    p = 0, q = 0, itr = 0, max_itr = n * n * n * n;
1149ba59ac12SSebastian Grimberg   CeedScalar tol = CEED_EPSILON;
1150ba59ac12SSebastian Grimberg 
1151ba59ac12SSebastian Grimberg   while (itr < max_itr) {
1152ba59ac12SSebastian Grimberg     // Update p, q, size of reduced portions of diagonal
1153ba59ac12SSebastian Grimberg     p = 0;
1154ba59ac12SSebastian Grimberg     q = 0;
1155ba59ac12SSebastian Grimberg     for (CeedInt i = n - 2; i >= 0; i--) {
1156ba59ac12SSebastian Grimberg       if (fabs(mat_T[i + n * (i + 1)]) < tol) q += 1;
1157ba59ac12SSebastian Grimberg       else break;
1158ba59ac12SSebastian Grimberg     }
1159ba59ac12SSebastian Grimberg     for (CeedInt i = 0; i < n - q - 1; i++) {
1160ba59ac12SSebastian Grimberg       if (fabs(mat_T[i + n * (i + 1)]) < tol) p += 1;
1161ba59ac12SSebastian Grimberg       else break;
1162ba59ac12SSebastian Grimberg     }
1163ba59ac12SSebastian Grimberg     if (q == n - 1) break;  // Finished reducing
1164ba59ac12SSebastian Grimberg 
1165ba59ac12SSebastian Grimberg     // Reduce tridiagonal portion
1166ba59ac12SSebastian Grimberg     CeedScalar t_nn = mat_T[(n - 1 - q) + n * (n - 1 - q)], t_nnm1 = mat_T[(n - 2 - q) + n * (n - 1 - q)];
1167ba59ac12SSebastian Grimberg     CeedScalar d  = (mat_T[(n - 2 - q) + n * (n - 2 - q)] - t_nn) / 2;
1168ba59ac12SSebastian Grimberg     CeedScalar mu = t_nn - t_nnm1 * t_nnm1 / (d + copysign(sqrt(d * d + t_nnm1 * t_nnm1), d));
1169ba59ac12SSebastian Grimberg     CeedScalar x  = mat_T[p + n * p] - mu;
1170ba59ac12SSebastian Grimberg     CeedScalar z  = mat_T[p + n * (p + 1)];
11711c66c397SJeremy L Thompson 
1172ba59ac12SSebastian Grimberg     for (CeedInt k = p; k < n - q - 1; k++) {
1173ba59ac12SSebastian Grimberg       // Compute Givens rotation
1174ba59ac12SSebastian Grimberg       CeedScalar c = 1, s = 0;
11751c66c397SJeremy L Thompson 
1176ba59ac12SSebastian Grimberg       if (fabs(z) > tol) {
1177ba59ac12SSebastian Grimberg         if (fabs(z) > fabs(x)) {
11781c66c397SJeremy L Thompson           const CeedScalar tau = -x / z;
11791c66c397SJeremy L Thompson 
11801c66c397SJeremy L Thompson           s = 1 / sqrt(1 + tau * tau);
11811c66c397SJeremy L Thompson           c = s * tau;
1182ba59ac12SSebastian Grimberg         } else {
11831c66c397SJeremy L Thompson           const CeedScalar tau = -z / x;
11841c66c397SJeremy L Thompson 
11851c66c397SJeremy L Thompson           c = 1 / sqrt(1 + tau * tau);
11861c66c397SJeremy L Thompson           s = c * tau;
1187ba59ac12SSebastian Grimberg         }
1188ba59ac12SSebastian Grimberg       }
1189ba59ac12SSebastian Grimberg 
1190ba59ac12SSebastian Grimberg       // Apply Givens rotation to T
1191ba59ac12SSebastian Grimberg       CeedGivensRotation(mat_T, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n);
1192ba59ac12SSebastian Grimberg       CeedGivensRotation(mat_T, c, s, CEED_TRANSPOSE, k, k + 1, n, n);
1193ba59ac12SSebastian Grimberg 
1194ba59ac12SSebastian Grimberg       // Apply Givens rotation to Q
1195ba59ac12SSebastian Grimberg       CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n);
1196ba59ac12SSebastian Grimberg 
1197ba59ac12SSebastian Grimberg       // Update x, z
1198ba59ac12SSebastian Grimberg       if (k < n - q - 2) {
1199ba59ac12SSebastian Grimberg         x = mat_T[k + n * (k + 1)];
1200ba59ac12SSebastian Grimberg         z = mat_T[k + n * (k + 2)];
1201ba59ac12SSebastian Grimberg       }
1202ba59ac12SSebastian Grimberg     }
1203ba59ac12SSebastian Grimberg     itr++;
1204ba59ac12SSebastian Grimberg   }
1205ba59ac12SSebastian Grimberg 
1206ba59ac12SSebastian Grimberg   // Save eigenvalues
1207ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) lambda[i] = mat_T[i + n * i];
1208ba59ac12SSebastian Grimberg 
1209ba59ac12SSebastian Grimberg   // Check convergence
12106574a04fSJeremy L Thompson   CeedCheck(itr < max_itr || q > n, ceed, CEED_ERROR_MINOR, "Symmetric QR failed to converge");
1211ba59ac12SSebastian Grimberg   return CEED_ERROR_SUCCESS;
1212ba59ac12SSebastian Grimberg }
12132c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOn
1214ba59ac12SSebastian Grimberg 
1215ba59ac12SSebastian Grimberg /**
1216ba59ac12SSebastian Grimberg   @brief Return Simultaneous Diagonalization of two matrices.
1217ba59ac12SSebastian Grimberg 
1218ca94c3ddSJeremy L Thompson   This solves the generalized eigenvalue problem `A x = lambda B x`, where `A` and `B` are symmetric and `B` is positive definite.
1219ca94c3ddSJeremy L Thompson   We generate the matrix `X` and vector `Lambda` such that `X^T A X = Lambda` and `X^T B X = I`.
1220ca94c3ddSJeremy L Thompson   This is equivalent to the LAPACK routine 'sygv' with `TYPE = 1`.
1221ba59ac12SSebastian Grimberg 
1222ca94c3ddSJeremy L Thompson   @param[in]  ceed   `Ceed` context for error handling
1223ba59ac12SSebastian Grimberg   @param[in]  mat_A  Row-major matrix to be factorized with eigenvalues
1224ba59ac12SSebastian Grimberg   @param[in]  mat_B  Row-major matrix to be factorized to identity
1225ba59ac12SSebastian Grimberg   @param[out] mat_X  Row-major orthogonal matrix
1226ca94c3ddSJeremy L Thompson   @param[out] lambda Vector of length `n` of generalized eigenvalues
1227ba59ac12SSebastian Grimberg   @param[in]  n      Number of rows/columns
1228ba59ac12SSebastian Grimberg 
1229ba59ac12SSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
1230ba59ac12SSebastian Grimberg 
1231ba59ac12SSebastian Grimberg   @ref Utility
1232ba59ac12SSebastian Grimberg **/
12332c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOff
12342c2ea1dbSJeremy L Thompson int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *mat_A, CeedScalar *mat_B, CeedScalar *mat_X, CeedScalar *lambda, CeedInt n) {
1235ba59ac12SSebastian Grimberg   CeedScalar *mat_C, *mat_G, *vec_D;
12361c66c397SJeremy L Thompson 
1237ba59ac12SSebastian Grimberg   CeedCall(CeedCalloc(n * n, &mat_C));
1238ba59ac12SSebastian Grimberg   CeedCall(CeedCalloc(n * n, &mat_G));
1239ba59ac12SSebastian Grimberg   CeedCall(CeedCalloc(n, &vec_D));
1240ba59ac12SSebastian Grimberg 
1241ba59ac12SSebastian Grimberg   // Compute B = G D G^T
1242ba59ac12SSebastian Grimberg   memcpy(mat_G, mat_B, n * n * sizeof(mat_B[0]));
1243ba59ac12SSebastian Grimberg   CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_G, vec_D, n));
1244ba59ac12SSebastian Grimberg 
1245ba59ac12SSebastian Grimberg   // Sort eigenvalues
1246ba59ac12SSebastian Grimberg   for (CeedInt i = n - 1; i >= 0; i--) {
1247ba59ac12SSebastian Grimberg     for (CeedInt j = 0; j < i; j++) {
1248ba59ac12SSebastian Grimberg       if (fabs(vec_D[j]) > fabs(vec_D[j + 1])) {
12491c66c397SJeremy L Thompson         CeedScalarSwap(vec_D[j], vec_D[j + 1]);
12501c66c397SJeremy L Thompson         for (CeedInt k = 0; k < n; k++) CeedScalarSwap(mat_G[k * n + j], mat_G[k * n + j + 1]);
1251ba59ac12SSebastian Grimberg       }
1252ba59ac12SSebastian Grimberg     }
1253ba59ac12SSebastian Grimberg   }
1254ba59ac12SSebastian Grimberg 
1255ba59ac12SSebastian Grimberg   // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T
1256ba59ac12SSebastian Grimberg   //           = D^-1/2 G^T A G D^-1/2
1257ba59ac12SSebastian Grimberg   // -- D = D^-1/2
1258ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) vec_D[i] = 1. / sqrt(vec_D[i]);
1259ba59ac12SSebastian Grimberg   // -- G = G D^-1/2
1260ba59ac12SSebastian Grimberg   // -- C = D^-1/2 G^T
1261ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) {
1262ba59ac12SSebastian Grimberg     for (CeedInt j = 0; j < n; j++) {
1263ba59ac12SSebastian Grimberg       mat_G[i * n + j] *= vec_D[j];
1264ba59ac12SSebastian Grimberg       mat_C[j * n + i] = mat_G[i * n + j];
1265ba59ac12SSebastian Grimberg     }
1266ba59ac12SSebastian Grimberg   }
1267ba59ac12SSebastian Grimberg   // -- X = (D^-1/2 G^T) A
1268ba59ac12SSebastian Grimberg   CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_C, (const CeedScalar *)mat_A, mat_X, n, n, n));
1269ba59ac12SSebastian Grimberg   // -- C = (D^-1/2 G^T A) (G D^-1/2)
1270ba59ac12SSebastian Grimberg   CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_X, (const CeedScalar *)mat_G, mat_C, n, n, n));
1271ba59ac12SSebastian Grimberg 
1272ba59ac12SSebastian Grimberg   // Compute Q^T C Q = lambda
1273ba59ac12SSebastian Grimberg   CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_C, lambda, n));
1274ba59ac12SSebastian Grimberg 
1275ba59ac12SSebastian Grimberg   // Sort eigenvalues
1276ba59ac12SSebastian Grimberg   for (CeedInt i = n - 1; i >= 0; i--) {
1277ba59ac12SSebastian Grimberg     for (CeedInt j = 0; j < i; j++) {
1278ba59ac12SSebastian Grimberg       if (fabs(lambda[j]) > fabs(lambda[j + 1])) {
12791c66c397SJeremy L Thompson         CeedScalarSwap(lambda[j], lambda[j + 1]);
12801c66c397SJeremy L Thompson         for (CeedInt k = 0; k < n; k++) CeedScalarSwap(mat_C[k * n + j], mat_C[k * n + j + 1]);
1281ba59ac12SSebastian Grimberg       }
1282ba59ac12SSebastian Grimberg     }
1283ba59ac12SSebastian Grimberg   }
1284ba59ac12SSebastian Grimberg 
1285ba59ac12SSebastian Grimberg   // Set X = (G D^1/2)^-T Q
1286ba59ac12SSebastian Grimberg   //       = G D^-1/2 Q
1287ba59ac12SSebastian Grimberg   CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_G, (const CeedScalar *)mat_C, mat_X, n, n, n));
1288ba59ac12SSebastian Grimberg 
1289ba59ac12SSebastian Grimberg   // Cleanup
1290ba59ac12SSebastian Grimberg   CeedCall(CeedFree(&mat_C));
1291ba59ac12SSebastian Grimberg   CeedCall(CeedFree(&mat_G));
1292ba59ac12SSebastian Grimberg   CeedCall(CeedFree(&vec_D));
1293ba59ac12SSebastian Grimberg   return CEED_ERROR_SUCCESS;
1294ba59ac12SSebastian Grimberg }
12952c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOn
1296ba59ac12SSebastian Grimberg 
12977a982d89SJeremy L. Thompson /// @}
12987a982d89SJeremy L. Thompson 
12997a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
13007a982d89SJeremy L. Thompson /// CeedBasis Public API
13017a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
13027a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser
1303d7b241e6Sjeremylt /// @{
1304d7b241e6Sjeremylt 
1305b11c1e72Sjeremylt /**
1306ca94c3ddSJeremy L Thompson   @brief Create a tensor-product basis for \f$H^1\f$ discretizations
1307b11c1e72Sjeremylt 
1308ca94c3ddSJeremy L Thompson   @param[in]  ceed        `Ceed` object used to create the `CeedBasis`
1309ea61e9acSJeremy L Thompson   @param[in]  dim         Topological dimension
1310ea61e9acSJeremy L Thompson   @param[in]  num_comp    Number of field components (1 for scalar fields)
1311ea61e9acSJeremy L Thompson   @param[in]  P_1d        Number of nodes in one dimension
1312ea61e9acSJeremy L Thompson   @param[in]  Q_1d        Number of quadrature points in one dimension
1313ca94c3ddSJeremy L Thompson   @param[in]  interp_1d   Row-major (`Q_1d * P_1d`) matrix expressing the values of nodal basis functions at quadrature points
1314ca94c3ddSJeremy L Thompson   @param[in]  grad_1d     Row-major (`Q_1d * P_1d`) matrix expressing derivatives of nodal basis functions at quadrature points
1315ca94c3ddSJeremy L Thompson   @param[in]  q_ref_1d    Array of length `Q_1d` holding the locations of quadrature points on the 1D reference element `[-1, 1]`
1316ca94c3ddSJeremy L Thompson   @param[in]  q_weight_1d Array of length `Q_1d` holding the quadrature weights on the reference element
1317ca94c3ddSJeremy L Thompson   @param[out] basis       Address of the variable where the newly created `CeedBasis` will be stored
1318b11c1e72Sjeremylt 
1319b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1320dfdf5a53Sjeremylt 
13217a982d89SJeremy L. Thompson   @ref User
1322b11c1e72Sjeremylt **/
13232b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d,
13242b730f8bSJeremy L Thompson                             const CeedScalar *grad_1d, const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis *basis) {
13255fe0d4faSjeremylt   if (!ceed->BasisCreateTensorH1) {
13265fe0d4faSjeremylt     Ceed delegate;
13276574a04fSJeremy L Thompson 
13282b730f8bSJeremy L Thompson     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
13291ef3a2a9SJeremy L Thompson     CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateTensorH1");
13302b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateTensorH1(delegate, dim, num_comp, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis));
1331e15f9bd0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
13325fe0d4faSjeremylt   }
1333e15f9bd0SJeremy L Thompson 
1334ca94c3ddSJeremy L Thompson   CeedCheck(dim > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis dimension must be a positive value");
1335ca94c3ddSJeremy L Thompson   CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 component");
1336ca94c3ddSJeremy L Thompson   CeedCheck(P_1d > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 node");
1337ca94c3ddSJeremy L Thompson   CeedCheck(Q_1d > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 quadrature point");
1338227444bfSJeremy L Thompson 
13392b730f8bSJeremy L Thompson   CeedElemTopology topo = dim == 1 ? CEED_TOPOLOGY_LINE : dim == 2 ? CEED_TOPOLOGY_QUAD : CEED_TOPOLOGY_HEX;
1340e15f9bd0SJeremy L Thompson 
13412b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(1, basis));
1342db002c03SJeremy L Thompson   CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed));
1343d1d35e2fSjeremylt   (*basis)->ref_count       = 1;
13446402da51SJeremy L Thompson   (*basis)->is_tensor_basis = true;
1345d7b241e6Sjeremylt   (*basis)->dim             = dim;
1346d99fa3c5SJeremy L Thompson   (*basis)->topo            = topo;
1347d1d35e2fSjeremylt   (*basis)->num_comp        = num_comp;
1348d1d35e2fSjeremylt   (*basis)->P_1d            = P_1d;
1349d1d35e2fSjeremylt   (*basis)->Q_1d            = Q_1d;
1350d1d35e2fSjeremylt   (*basis)->P               = CeedIntPow(P_1d, dim);
1351d1d35e2fSjeremylt   (*basis)->Q               = CeedIntPow(Q_1d, dim);
1352c4e3f59bSSebastian Grimberg   (*basis)->fe_space        = CEED_FE_SPACE_H1;
13532b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d, &(*basis)->q_ref_1d));
13542b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d, &(*basis)->q_weight_1d));
1355ff3a0f91SJeremy L Thompson   if (q_ref_1d) memcpy((*basis)->q_ref_1d, q_ref_1d, Q_1d * sizeof(q_ref_1d[0]));
13562b730f8bSJeremy L Thompson   if (q_weight_1d) memcpy((*basis)->q_weight_1d, q_weight_1d, Q_1d * sizeof(q_weight_1d[0]));
13572b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->interp_1d));
13582b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->grad_1d));
13592b730f8bSJeremy L Thompson   if (interp_1d) memcpy((*basis)->interp_1d, interp_1d, Q_1d * P_1d * sizeof(interp_1d[0]));
1360ff3a0f91SJeremy L Thompson   if (grad_1d) memcpy((*basis)->grad_1d, grad_1d, Q_1d * P_1d * sizeof(grad_1d[0]));
13612b730f8bSJeremy L Thompson   CeedCall(ceed->BasisCreateTensorH1(dim, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, *basis));
1362e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1363d7b241e6Sjeremylt }
1364d7b241e6Sjeremylt 
1365b11c1e72Sjeremylt /**
1366ca94c3ddSJeremy L Thompson   @brief Create a tensor-product \f$H^1\f$ Lagrange basis
1367b11c1e72Sjeremylt 
1368ca94c3ddSJeremy L Thompson   @param[in]  ceed      `Ceed` object used to create the `CeedBasis`
1369ea61e9acSJeremy L Thompson   @param[in]  dim       Topological dimension of element
1370ea61e9acSJeremy L Thompson   @param[in]  num_comp  Number of field components (1 for scalar fields)
1371ea61e9acSJeremy L Thompson   @param[in]  P         Number of Gauss-Lobatto nodes in one dimension.
1372ca94c3ddSJeremy L Thompson                           The polynomial degree of the resulting `Q_k` element is `k = P - 1`.
1373ea61e9acSJeremy L Thompson   @param[in]  Q         Number of quadrature points in one dimension.
1374ca94c3ddSJeremy L Thompson   @param[in]  quad_mode Distribution of the `Q` quadrature points (affects order of accuracy for the quadrature)
1375ca94c3ddSJeremy L Thompson   @param[out] basis     Address of the variable where the newly created `CeedBasis` will be stored
1376b11c1e72Sjeremylt 
1377b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1378dfdf5a53Sjeremylt 
13797a982d89SJeremy L. Thompson   @ref User
1380b11c1e72Sjeremylt **/
13812b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P, CeedInt Q, CeedQuadMode quad_mode, CeedBasis *basis) {
1382d7b241e6Sjeremylt   // Allocate
1383c8c3fa7dSJeremy L Thompson   int        ierr = CEED_ERROR_SUCCESS;
13842b730f8bSJeremy L Thompson   CeedScalar c1, c2, c3, c4, dx, *nodes, *interp_1d, *grad_1d, *q_ref_1d, *q_weight_1d;
13854d537eeaSYohann 
1386ca94c3ddSJeremy L Thompson   CeedCheck(dim > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis dimension must be a positive value");
1387ca94c3ddSJeremy L Thompson   CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 component");
1388ca94c3ddSJeremy L Thompson   CeedCheck(P > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 node");
1389ca94c3ddSJeremy L Thompson   CeedCheck(Q > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 quadrature point");
1390227444bfSJeremy L Thompson 
1391e15f9bd0SJeremy L Thompson   // Get Nodes and Weights
13922b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(P * Q, &interp_1d));
13932b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(P * Q, &grad_1d));
13942b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(P, &nodes));
13952b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q, &q_ref_1d));
13962b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q, &q_weight_1d));
13972b730f8bSJeremy L Thompson   if (CeedLobattoQuadrature(P, nodes, NULL) != CEED_ERROR_SUCCESS) goto cleanup;
1398d1d35e2fSjeremylt   switch (quad_mode) {
1399d7b241e6Sjeremylt     case CEED_GAUSS:
1400d1d35e2fSjeremylt       ierr = CeedGaussQuadrature(Q, q_ref_1d, q_weight_1d);
1401d7b241e6Sjeremylt       break;
1402d7b241e6Sjeremylt     case CEED_GAUSS_LOBATTO:
1403d1d35e2fSjeremylt       ierr = CeedLobattoQuadrature(Q, q_ref_1d, q_weight_1d);
1404d7b241e6Sjeremylt       break;
1405d7b241e6Sjeremylt   }
14062b730f8bSJeremy L Thompson   if (ierr != CEED_ERROR_SUCCESS) goto cleanup;
1407e15f9bd0SJeremy L Thompson 
1408d7b241e6Sjeremylt   // Build B, D matrix
1409d7b241e6Sjeremylt   // Fornberg, 1998
1410c8c3fa7dSJeremy L Thompson   for (CeedInt i = 0; i < Q; i++) {
1411d7b241e6Sjeremylt     c1                   = 1.0;
1412d1d35e2fSjeremylt     c3                   = nodes[0] - q_ref_1d[i];
1413d1d35e2fSjeremylt     interp_1d[i * P + 0] = 1.0;
1414c8c3fa7dSJeremy L Thompson     for (CeedInt j = 1; j < P; j++) {
1415d7b241e6Sjeremylt       c2 = 1.0;
1416d7b241e6Sjeremylt       c4 = c3;
1417d1d35e2fSjeremylt       c3 = nodes[j] - q_ref_1d[i];
1418c8c3fa7dSJeremy L Thompson       for (CeedInt k = 0; k < j; k++) {
1419d7b241e6Sjeremylt         dx = nodes[j] - nodes[k];
1420d7b241e6Sjeremylt         c2 *= dx;
1421d7b241e6Sjeremylt         if (k == j - 1) {
1422d1d35e2fSjeremylt           grad_1d[i * P + j]   = c1 * (interp_1d[i * P + k] - c4 * grad_1d[i * P + k]) / c2;
1423d1d35e2fSjeremylt           interp_1d[i * P + j] = -c1 * c4 * interp_1d[i * P + k] / c2;
1424d7b241e6Sjeremylt         }
1425d1d35e2fSjeremylt         grad_1d[i * P + k]   = (c3 * grad_1d[i * P + k] - interp_1d[i * P + k]) / dx;
1426d1d35e2fSjeremylt         interp_1d[i * P + k] = c3 * interp_1d[i * P + k] / dx;
1427d7b241e6Sjeremylt       }
1428d7b241e6Sjeremylt       c1 = c2;
1429d7b241e6Sjeremylt     }
1430d7b241e6Sjeremylt   }
14319ac7b42eSJeremy L Thompson   // Pass to CeedBasisCreateTensorH1
14322b730f8bSJeremy L Thompson   CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P, Q, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis));
1433e15f9bd0SJeremy L Thompson cleanup:
14342b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_1d));
14352b730f8bSJeremy L Thompson   CeedCall(CeedFree(&grad_1d));
14362b730f8bSJeremy L Thompson   CeedCall(CeedFree(&nodes));
14372b730f8bSJeremy L Thompson   CeedCall(CeedFree(&q_ref_1d));
14382b730f8bSJeremy L Thompson   CeedCall(CeedFree(&q_weight_1d));
1439e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1440d7b241e6Sjeremylt }
1441d7b241e6Sjeremylt 
1442b11c1e72Sjeremylt /**
1443ca94c3ddSJeremy L Thompson   @brief Create a non tensor-product basis for \f$H^1\f$ discretizations
1444a8de75f0Sjeremylt 
1445ca94c3ddSJeremy L Thompson   @param[in]  ceed      `Ceed` object used to create the `CeedBasis`
1446e00f3be8SJames Wright   @param[in]  topo      Topology of element, e.g. hypercube, simplex, etc
1447ea61e9acSJeremy L Thompson   @param[in]  num_comp  Number of field components (1 for scalar fields)
1448ea61e9acSJeremy L Thompson   @param[in]  num_nodes Total number of nodes
1449ea61e9acSJeremy L Thompson   @param[in]  num_qpts  Total number of quadrature points
1450ca94c3ddSJeremy L Thompson   @param[in]  interp    Row-major (`num_qpts * num_nodes`) matrix expressing the values of nodal basis functions at quadrature points
1451ca94c3ddSJeremy L Thompson   @param[in]  grad      Row-major (`dim * num_qpts * num_nodes`) matrix expressing derivatives of nodal basis functions at quadrature points
1452ca94c3ddSJeremy L Thompson   @param[in]  q_ref     Array of length `num_qpts` * dim holding the locations of quadrature points on the reference element
1453ca94c3ddSJeremy L Thompson   @param[in]  q_weight  Array of length `num_qpts` holding the quadrature weights on the reference element
1454ca94c3ddSJeremy L Thompson   @param[out] basis     Address of the variable where the newly created `CeedBasis` will be stored
1455a8de75f0Sjeremylt 
1456a8de75f0Sjeremylt   @return An error code: 0 - success, otherwise - failure
1457a8de75f0Sjeremylt 
14587a982d89SJeremy L. Thompson   @ref User
1459a8de75f0Sjeremylt **/
14602b730f8bSJeremy L Thompson int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
14612b730f8bSJeremy L Thompson                       const CeedScalar *grad, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) {
1462d1d35e2fSjeremylt   CeedInt P = num_nodes, Q = num_qpts, dim = 0;
1463a8de75f0Sjeremylt 
14645fe0d4faSjeremylt   if (!ceed->BasisCreateH1) {
14655fe0d4faSjeremylt     Ceed delegate;
14666574a04fSJeremy L Thompson 
14672b730f8bSJeremy L Thompson     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
14681ef3a2a9SJeremy L Thompson     CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateH1");
14692b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateH1(delegate, topo, num_comp, num_nodes, num_qpts, interp, grad, q_ref, q_weight, basis));
1470e15f9bd0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
14715fe0d4faSjeremylt   }
14725fe0d4faSjeremylt 
1473ca94c3ddSJeremy L Thompson   CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 component");
1474ca94c3ddSJeremy L Thompson   CeedCheck(num_nodes > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 node");
1475ca94c3ddSJeremy L Thompson   CeedCheck(num_qpts > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 quadrature point");
1476227444bfSJeremy L Thompson 
14772b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetTopologyDimension(topo, &dim));
1478a8de75f0Sjeremylt 
1479db002c03SJeremy L Thompson   CeedCall(CeedCalloc(1, basis));
1480db002c03SJeremy L Thompson   CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed));
1481d1d35e2fSjeremylt   (*basis)->ref_count       = 1;
14826402da51SJeremy L Thompson   (*basis)->is_tensor_basis = false;
1483a8de75f0Sjeremylt   (*basis)->dim             = dim;
1484d99fa3c5SJeremy L Thompson   (*basis)->topo            = topo;
1485d1d35e2fSjeremylt   (*basis)->num_comp        = num_comp;
1486a8de75f0Sjeremylt   (*basis)->P               = P;
1487a8de75f0Sjeremylt   (*basis)->Q               = Q;
1488c4e3f59bSSebastian Grimberg   (*basis)->fe_space        = CEED_FE_SPACE_H1;
14892b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q * dim, &(*basis)->q_ref_1d));
14902b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q, &(*basis)->q_weight_1d));
1491ff3a0f91SJeremy L Thompson   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0]));
1492ff3a0f91SJeremy L Thompson   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0]));
14932b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q * P, &(*basis)->interp));
14942b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(dim * Q * P, &(*basis)->grad));
1495ff3a0f91SJeremy L Thompson   if (interp) memcpy((*basis)->interp, interp, Q * P * sizeof(interp[0]));
1496ff3a0f91SJeremy L Thompson   if (grad) memcpy((*basis)->grad, grad, dim * Q * P * sizeof(grad[0]));
14972b730f8bSJeremy L Thompson   CeedCall(ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, q_ref, q_weight, *basis));
1498e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1499a8de75f0Sjeremylt }
1500a8de75f0Sjeremylt 
1501a8de75f0Sjeremylt /**
1502859c15bbSJames Wright   @brief Create a non tensor-product basis for \f$H(\mathrm{div})\f$ discretizations
150350c301a5SRezgar Shakeri 
1504ca94c3ddSJeremy L Thompson   @param[in]  ceed      `Ceed` object used to create the `CeedBasis`
1505ea61e9acSJeremy L Thompson   @param[in]  topo      Topology of element (`CEED_TOPOLOGY_QUAD`, `CEED_TOPOLOGY_PRISM`, etc.), dimension of which is used in some array sizes below
1506ea61e9acSJeremy L Thompson   @param[in]  num_comp  Number of components (usually 1 for vectors in H(div) bases)
1507ca94c3ddSJeremy L Thompson   @param[in]  num_nodes Total number of nodes (DoFs per element)
1508ea61e9acSJeremy L Thompson   @param[in]  num_qpts  Total number of quadrature points
1509ca94c3ddSJeremy L Thompson   @param[in]  interp    Row-major (`dim * num_qpts * num_nodes`) matrix expressing the values of basis functions at quadrature points
1510ca94c3ddSJeremy L Thompson   @param[in]  div       Row-major (`num_qpts * num_nodes`) matrix expressing divergence of basis functions at quadrature points
1511ca94c3ddSJeremy L Thompson   @param[in]  q_ref     Array of length `num_qpts` * dim holding the locations of quadrature points on the reference element
1512ca94c3ddSJeremy L Thompson   @param[in]  q_weight  Array of length `num_qpts` holding the quadrature weights on the reference element
1513ca94c3ddSJeremy L Thompson   @param[out] basis     Address of the variable where the newly created `CeedBasis` will be stored
151450c301a5SRezgar Shakeri 
151550c301a5SRezgar Shakeri   @return An error code: 0 - success, otherwise - failure
151650c301a5SRezgar Shakeri 
151750c301a5SRezgar Shakeri   @ref User
151850c301a5SRezgar Shakeri **/
15192b730f8bSJeremy L Thompson int CeedBasisCreateHdiv(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
15202b730f8bSJeremy L Thompson                         const CeedScalar *div, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) {
152150c301a5SRezgar Shakeri   CeedInt Q = num_qpts, P = num_nodes, dim = 0;
1522c4e3f59bSSebastian Grimberg 
152350c301a5SRezgar Shakeri   if (!ceed->BasisCreateHdiv) {
152450c301a5SRezgar Shakeri     Ceed delegate;
15256574a04fSJeremy L Thompson 
15262b730f8bSJeremy L Thompson     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
15276574a04fSJeremy L Thompson     CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateHdiv");
15282b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateHdiv(delegate, topo, num_comp, num_nodes, num_qpts, interp, div, q_ref, q_weight, basis));
152950c301a5SRezgar Shakeri     return CEED_ERROR_SUCCESS;
153050c301a5SRezgar Shakeri   }
153150c301a5SRezgar Shakeri 
1532ca94c3ddSJeremy L Thompson   CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 component");
1533ca94c3ddSJeremy L Thompson   CeedCheck(num_nodes > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 node");
1534ca94c3ddSJeremy L Thompson   CeedCheck(num_qpts > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 quadrature point");
1535227444bfSJeremy L Thompson 
1536c4e3f59bSSebastian Grimberg   CeedCall(CeedBasisGetTopologyDimension(topo, &dim));
1537c4e3f59bSSebastian Grimberg 
1538db002c03SJeremy L Thompson   CeedCall(CeedCalloc(1, basis));
1539db002c03SJeremy L Thompson   CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed));
154050c301a5SRezgar Shakeri   (*basis)->ref_count       = 1;
15416402da51SJeremy L Thompson   (*basis)->is_tensor_basis = false;
154250c301a5SRezgar Shakeri   (*basis)->dim             = dim;
154350c301a5SRezgar Shakeri   (*basis)->topo            = topo;
154450c301a5SRezgar Shakeri   (*basis)->num_comp        = num_comp;
154550c301a5SRezgar Shakeri   (*basis)->P               = P;
154650c301a5SRezgar Shakeri   (*basis)->Q               = Q;
1547c4e3f59bSSebastian Grimberg   (*basis)->fe_space        = CEED_FE_SPACE_HDIV;
15482b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q * dim, &(*basis)->q_ref_1d));
15492b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q, &(*basis)->q_weight_1d));
155050c301a5SRezgar Shakeri   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0]));
155150c301a5SRezgar Shakeri   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0]));
15522b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(dim * Q * P, &(*basis)->interp));
15532b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q * P, &(*basis)->div));
155450c301a5SRezgar Shakeri   if (interp) memcpy((*basis)->interp, interp, dim * Q * P * sizeof(interp[0]));
155550c301a5SRezgar Shakeri   if (div) memcpy((*basis)->div, div, Q * P * sizeof(div[0]));
15562b730f8bSJeremy L Thompson   CeedCall(ceed->BasisCreateHdiv(topo, dim, P, Q, interp, div, q_ref, q_weight, *basis));
155750c301a5SRezgar Shakeri   return CEED_ERROR_SUCCESS;
155850c301a5SRezgar Shakeri }
155950c301a5SRezgar Shakeri 
156050c301a5SRezgar Shakeri /**
15614385fb7fSSebastian Grimberg   @brief Create a non tensor-product basis for \f$H(\mathrm{curl})\f$ discretizations
1562c4e3f59bSSebastian Grimberg 
1563ca94c3ddSJeremy L Thompson   @param[in]  ceed      `Ceed` object used to create the `CeedBasis`
1564c4e3f59bSSebastian Grimberg   @param[in]  topo      Topology of element (`CEED_TOPOLOGY_QUAD`, `CEED_TOPOLOGY_PRISM`, etc.), dimension of which is used in some array sizes below
1565ca94c3ddSJeremy L Thompson   @param[in]  num_comp  Number of components (usually 1 for vectors in \f$H(\mathrm{curl})\f$ bases)
1566ca94c3ddSJeremy L Thompson   @param[in]  num_nodes Total number of nodes (DoFs per element)
1567c4e3f59bSSebastian Grimberg   @param[in]  num_qpts  Total number of quadrature points
1568ca94c3ddSJeremy L Thompson   @param[in]  interp    Row-major (`dim * num_qpts * num_nodes`) matrix expressing the values of basis functions at quadrature points
1569ca94c3ddSJeremy L Thompson   @param[in]  curl      Row-major (`curl_comp * num_qpts * num_nodes`, `curl_comp = 1` if `dim < 3` otherwise `curl_comp = dim`) matrix expressing curl of basis functions at quadrature points
1570ca94c3ddSJeremy L Thompson   @param[in]  q_ref     Array of length `num_qpts * dim` holding the locations of quadrature points on the reference element
1571ca94c3ddSJeremy L Thompson   @param[in]  q_weight  Array of length `num_qpts` holding the quadrature weights on the reference element
1572ca94c3ddSJeremy L Thompson   @param[out] basis     Address of the variable where the newly created `CeedBasis` will be stored
1573c4e3f59bSSebastian Grimberg 
1574c4e3f59bSSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
1575c4e3f59bSSebastian Grimberg 
1576c4e3f59bSSebastian Grimberg   @ref User
1577c4e3f59bSSebastian Grimberg **/
1578c4e3f59bSSebastian Grimberg int CeedBasisCreateHcurl(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
1579c4e3f59bSSebastian Grimberg                          const CeedScalar *curl, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) {
1580c4e3f59bSSebastian Grimberg   CeedInt Q = num_qpts, P = num_nodes, dim = 0, curl_comp = 0;
1581c4e3f59bSSebastian Grimberg 
1582d075f50bSSebastian Grimberg   if (!ceed->BasisCreateHcurl) {
1583c4e3f59bSSebastian Grimberg     Ceed delegate;
15846574a04fSJeremy L Thompson 
1585c4e3f59bSSebastian Grimberg     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
15866574a04fSJeremy L Thompson     CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateHcurl");
1587c4e3f59bSSebastian Grimberg     CeedCall(CeedBasisCreateHcurl(delegate, topo, num_comp, num_nodes, num_qpts, interp, curl, q_ref, q_weight, basis));
1588c4e3f59bSSebastian Grimberg     return CEED_ERROR_SUCCESS;
1589c4e3f59bSSebastian Grimberg   }
1590c4e3f59bSSebastian Grimberg 
1591ca94c3ddSJeremy L Thompson   CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 component");
1592ca94c3ddSJeremy L Thompson   CeedCheck(num_nodes > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 node");
1593ca94c3ddSJeremy L Thompson   CeedCheck(num_qpts > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 quadrature point");
1594c4e3f59bSSebastian Grimberg 
1595c4e3f59bSSebastian Grimberg   CeedCall(CeedBasisGetTopologyDimension(topo, &dim));
1596c4e3f59bSSebastian Grimberg   curl_comp = (dim < 3) ? 1 : dim;
1597c4e3f59bSSebastian Grimberg 
1598db002c03SJeremy L Thompson   CeedCall(CeedCalloc(1, basis));
1599db002c03SJeremy L Thompson   CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed));
1600c4e3f59bSSebastian Grimberg   (*basis)->ref_count       = 1;
16016402da51SJeremy L Thompson   (*basis)->is_tensor_basis = false;
1602c4e3f59bSSebastian Grimberg   (*basis)->dim             = dim;
1603c4e3f59bSSebastian Grimberg   (*basis)->topo            = topo;
1604c4e3f59bSSebastian Grimberg   (*basis)->num_comp        = num_comp;
1605c4e3f59bSSebastian Grimberg   (*basis)->P               = P;
1606c4e3f59bSSebastian Grimberg   (*basis)->Q               = Q;
1607c4e3f59bSSebastian Grimberg   (*basis)->fe_space        = CEED_FE_SPACE_HCURL;
1608c4e3f59bSSebastian Grimberg   CeedCall(CeedMalloc(Q * dim, &(*basis)->q_ref_1d));
1609c4e3f59bSSebastian Grimberg   CeedCall(CeedMalloc(Q, &(*basis)->q_weight_1d));
1610c4e3f59bSSebastian Grimberg   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0]));
1611c4e3f59bSSebastian Grimberg   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0]));
1612c4e3f59bSSebastian Grimberg   CeedCall(CeedMalloc(dim * Q * P, &(*basis)->interp));
1613c4e3f59bSSebastian Grimberg   CeedCall(CeedMalloc(curl_comp * Q * P, &(*basis)->curl));
1614c4e3f59bSSebastian Grimberg   if (interp) memcpy((*basis)->interp, interp, dim * Q * P * sizeof(interp[0]));
1615c4e3f59bSSebastian Grimberg   if (curl) memcpy((*basis)->curl, curl, curl_comp * Q * P * sizeof(curl[0]));
1616c4e3f59bSSebastian Grimberg   CeedCall(ceed->BasisCreateHcurl(topo, dim, P, Q, interp, curl, q_ref, q_weight, *basis));
1617c4e3f59bSSebastian Grimberg   return CEED_ERROR_SUCCESS;
1618c4e3f59bSSebastian Grimberg }
1619c4e3f59bSSebastian Grimberg 
1620c4e3f59bSSebastian Grimberg /**
1621ca94c3ddSJeremy L Thompson   @brief Create a `CeedBasis` for projection from the nodes of `basis_from` to the nodes of `basis_to`.
1622ba59ac12SSebastian Grimberg 
1623ca94c3ddSJeremy L Thompson   Only @ref CEED_EVAL_INTERP will be valid for the new basis, `basis_project`.
1624ca94c3ddSJeremy L Thompson   For \f$H^1\f$ spaces, @ref CEED_EVAL_GRAD will also be valid.
1625ca94c3ddSJeremy L Thompson   The interpolation is given by `interp_project = interp_to^+ * interp_from`, where the pseudoinverse `interp_to^+` is given by QR factorization.
1626ca94c3ddSJeremy L Thompson   The gradient (for the \f$H^1\f$ case) is given by `grad_project = interp_to^+ * grad_from`.
162715ad3917SSebastian Grimberg 
162815ad3917SSebastian Grimberg   Note: `basis_from` and `basis_to` must have compatible quadrature spaces.
162915ad3917SSebastian Grimberg 
16309fd66db6SSebastian Grimberg   Note: `basis_project` will have the same number of components as `basis_from`, regardless of the number of components that `basis_to` has.
16319fd66db6SSebastian Grimberg         If `basis_from` has 3 components and `basis_to` has 5 components, then `basis_project` will have 3 components.
1632f113e5dcSJeremy L Thompson 
1633e104ad11SJames Wright   Note: If either `basis_from` or `basis_to` are non-tensor, then `basis_project` will also be non-tensor
1634e104ad11SJames Wright 
1635ca94c3ddSJeremy L Thompson   @param[in]  basis_from    `CeedBasis` to prolong from
1636ca94c3ddSJeremy L Thompson   @param[in]  basis_to      `CeedBasis` to prolong to
1637ca94c3ddSJeremy L Thompson   @param[out] basis_project Address of the variable where the newly created `CeedBasis` will be stored
1638f113e5dcSJeremy L Thompson 
1639f113e5dcSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1640f113e5dcSJeremy L Thompson 
1641f113e5dcSJeremy L Thompson   @ref User
1642f113e5dcSJeremy L Thompson **/
16432b730f8bSJeremy L Thompson int CeedBasisCreateProjection(CeedBasis basis_from, CeedBasis basis_to, CeedBasis *basis_project) {
1644f113e5dcSJeremy L Thompson   Ceed        ceed;
1645e104ad11SJames Wright   bool        create_tensor;
16461c66c397SJeremy L Thompson   CeedInt     dim, num_comp;
1647097cc795SJames Wright   CeedScalar *interp_project, *grad_project;
16481c66c397SJeremy L Thompson 
16492b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetCeed(basis_to, &ceed));
1650f113e5dcSJeremy L Thompson 
1651ecc88aebSJeremy L Thompson   // Create projection matrix
16522b730f8bSJeremy L Thompson   CeedCall(CeedBasisCreateProjectionMatrices(basis_from, basis_to, &interp_project, &grad_project));
1653f113e5dcSJeremy L Thompson 
1654f113e5dcSJeremy L Thompson   // Build basis
1655e104ad11SJames Wright   {
1656e104ad11SJames Wright     bool is_tensor_to, is_tensor_from;
1657e104ad11SJames Wright 
1658e104ad11SJames Wright     CeedCall(CeedBasisIsTensor(basis_to, &is_tensor_to));
1659e104ad11SJames Wright     CeedCall(CeedBasisIsTensor(basis_from, &is_tensor_from));
1660e104ad11SJames Wright     create_tensor = is_tensor_from && is_tensor_to;
1661e104ad11SJames Wright   }
16622b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetDimension(basis_to, &dim));
16632b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumComponents(basis_from, &num_comp));
1664e104ad11SJames Wright   if (create_tensor) {
1665f113e5dcSJeremy L Thompson     CeedInt P_1d_to, P_1d_from;
16661c66c397SJeremy L Thompson 
16672b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_1d_from));
16682b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_1d_to));
1669097cc795SJames Wright     CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P_1d_from, P_1d_to, interp_project, grad_project, NULL, NULL, basis_project));
1670f113e5dcSJeremy L Thompson   } else {
1671de05fbb2SSebastian Grimberg     // Even if basis_to and basis_from are not H1, the resulting basis is H1 for interpolation to work
1672f113e5dcSJeremy L Thompson     CeedInt          num_nodes_to, num_nodes_from;
16731c66c397SJeremy L Thompson     CeedElemTopology topo;
16741c66c397SJeremy L Thompson 
1675e00f3be8SJames Wright     CeedCall(CeedBasisGetTopology(basis_from, &topo));
16762b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis_from, &num_nodes_from));
16772b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis_to, &num_nodes_to));
1678097cc795SJames Wright     CeedCall(CeedBasisCreateH1(ceed, topo, num_comp, num_nodes_from, num_nodes_to, interp_project, grad_project, NULL, NULL, basis_project));
1679f113e5dcSJeremy L Thompson   }
1680f113e5dcSJeremy L Thompson 
1681f113e5dcSJeremy L Thompson   // Cleanup
16822b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_project));
16832b730f8bSJeremy L Thompson   CeedCall(CeedFree(&grad_project));
1684f113e5dcSJeremy L Thompson   return CEED_ERROR_SUCCESS;
1685f113e5dcSJeremy L Thompson }
1686f113e5dcSJeremy L Thompson 
1687f113e5dcSJeremy L Thompson /**
1688ca94c3ddSJeremy L Thompson   @brief Copy the pointer to a `CeedBasis`.
16899560d06aSjeremylt 
1690ca94c3ddSJeremy L Thompson   Note: If the value of `*basis_copy` passed into this function is non-`NULL`, then it is assumed that `*basis_copy` is a pointer to a `CeedBasis`.
1691ca94c3ddSJeremy L Thompson         This `CeedBasis` will be destroyed if `*basis_copy` is the only reference to this `CeedBasis`.
1692ea61e9acSJeremy L Thompson 
1693ca94c3ddSJeremy L Thompson   @param[in]     basis      `CeedBasis` to copy reference to
1694ea61e9acSJeremy L Thompson   @param[in,out] basis_copy Variable to store copied reference
16959560d06aSjeremylt 
16969560d06aSjeremylt   @return An error code: 0 - success, otherwise - failure
16979560d06aSjeremylt 
16989560d06aSjeremylt   @ref User
16999560d06aSjeremylt **/
17009560d06aSjeremylt int CeedBasisReferenceCopy(CeedBasis basis, CeedBasis *basis_copy) {
1701356036faSJeremy L Thompson   if (basis != CEED_BASIS_NONE) CeedCall(CeedBasisReference(basis));
17022b730f8bSJeremy L Thompson   CeedCall(CeedBasisDestroy(basis_copy));
17039560d06aSjeremylt   *basis_copy = basis;
17049560d06aSjeremylt   return CEED_ERROR_SUCCESS;
17059560d06aSjeremylt }
17069560d06aSjeremylt 
17079560d06aSjeremylt /**
1708ca94c3ddSJeremy L Thompson   @brief View a `CeedBasis`
17097a982d89SJeremy L. Thompson 
1710ca94c3ddSJeremy L Thompson   @param[in] basis  `CeedBasis` to view
1711ca94c3ddSJeremy L Thompson   @param[in] stream Stream to view to, e.g., `stdout`
17127a982d89SJeremy L. Thompson 
17137a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
17147a982d89SJeremy L. Thompson 
17157a982d89SJeremy L. Thompson   @ref User
17167a982d89SJeremy L. Thompson **/
17177a982d89SJeremy L. Thompson int CeedBasisView(CeedBasis basis, FILE *stream) {
17181203703bSJeremy L Thompson   bool             is_tensor_basis;
17191203703bSJeremy L Thompson   CeedElemTopology topo;
17201203703bSJeremy L Thompson   CeedFESpace      fe_space;
17211203703bSJeremy L Thompson 
17221203703bSJeremy L Thompson   // Basis data
17231203703bSJeremy L Thompson   CeedCall(CeedBasisIsTensor(basis, &is_tensor_basis));
17241203703bSJeremy L Thompson   CeedCall(CeedBasisGetTopology(basis, &topo));
17251203703bSJeremy L Thompson   CeedCall(CeedBasisGetFESpace(basis, &fe_space));
17262b730f8bSJeremy L Thompson 
172750c301a5SRezgar Shakeri   // Print FE space and element topology of the basis
1728edf04919SJeremy L Thompson   fprintf(stream, "CeedBasis in a %s on a %s element\n", CeedFESpaces[fe_space], CeedElemTopologies[topo]);
17291203703bSJeremy L Thompson   if (is_tensor_basis) {
1730edf04919SJeremy L Thompson     fprintf(stream, "  P: %" CeedInt_FMT "\n  Q: %" CeedInt_FMT "\n", basis->P_1d, basis->Q_1d);
173150c301a5SRezgar Shakeri   } else {
1732edf04919SJeremy L Thompson     fprintf(stream, "  P: %" CeedInt_FMT "\n  Q: %" CeedInt_FMT "\n", basis->P, basis->Q);
173350c301a5SRezgar Shakeri   }
1734edf04919SJeremy L Thompson   fprintf(stream, "  dimension: %" CeedInt_FMT "\n  field components: %" CeedInt_FMT "\n", basis->dim, basis->num_comp);
1735ea61e9acSJeremy L Thompson   // Print quadrature data, interpolation/gradient/divergence/curl of the basis
17361203703bSJeremy L Thompson   if (is_tensor_basis) {  // tensor basis
17371203703bSJeremy L Thompson     CeedInt           P_1d, Q_1d;
17381203703bSJeremy L Thompson     const CeedScalar *q_ref_1d, *q_weight_1d, *interp_1d, *grad_1d;
17391203703bSJeremy L Thompson 
17401203703bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d));
17411203703bSJeremy L Thompson     CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
17421203703bSJeremy L Thompson     CeedCall(CeedBasisGetQRef(basis, &q_ref_1d));
17431203703bSJeremy L Thompson     CeedCall(CeedBasisGetQWeights(basis, &q_weight_1d));
17441203703bSJeremy L Thompson     CeedCall(CeedBasisGetInterp1D(basis, &interp_1d));
17451203703bSJeremy L Thompson     CeedCall(CeedBasisGetGrad1D(basis, &grad_1d));
17461203703bSJeremy L Thompson 
17471203703bSJeremy L Thompson     CeedCall(CeedScalarView("qref1d", "\t% 12.8f", 1, Q_1d, q_ref_1d, stream));
17481203703bSJeremy L Thompson     CeedCall(CeedScalarView("qweight1d", "\t% 12.8f", 1, Q_1d, q_weight_1d, stream));
17491203703bSJeremy L Thompson     CeedCall(CeedScalarView("interp1d", "\t% 12.8f", Q_1d, P_1d, interp_1d, stream));
17501203703bSJeremy L Thompson     CeedCall(CeedScalarView("grad1d", "\t% 12.8f", Q_1d, P_1d, grad_1d, stream));
175150c301a5SRezgar Shakeri   } else {  // non-tensor basis
17521203703bSJeremy L Thompson     CeedInt           P, Q, dim, q_comp;
17531203703bSJeremy L Thompson     const CeedScalar *q_ref, *q_weight, *interp, *grad, *div, *curl;
17541203703bSJeremy L Thompson 
17551203703bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis, &P));
17561203703bSJeremy L Thompson     CeedCall(CeedBasisGetNumQuadraturePoints(basis, &Q));
17571203703bSJeremy L Thompson     CeedCall(CeedBasisGetDimension(basis, &dim));
17581203703bSJeremy L Thompson     CeedCall(CeedBasisGetQRef(basis, &q_ref));
17591203703bSJeremy L Thompson     CeedCall(CeedBasisGetQWeights(basis, &q_weight));
17601203703bSJeremy L Thompson     CeedCall(CeedBasisGetInterp(basis, &interp));
17611203703bSJeremy L Thompson     CeedCall(CeedBasisGetGrad(basis, &grad));
17621203703bSJeremy L Thompson     CeedCall(CeedBasisGetDiv(basis, &div));
17631203703bSJeremy L Thompson     CeedCall(CeedBasisGetCurl(basis, &curl));
17641203703bSJeremy L Thompson 
17651203703bSJeremy L Thompson     CeedCall(CeedScalarView("qref", "\t% 12.8f", 1, Q * dim, q_ref, stream));
17661203703bSJeremy L Thompson     CeedCall(CeedScalarView("qweight", "\t% 12.8f", 1, Q, q_weight, stream));
1767c4e3f59bSSebastian Grimberg     CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp));
17681203703bSJeremy L Thompson     CeedCall(CeedScalarView("interp", "\t% 12.8f", q_comp * Q, P, interp, stream));
17691203703bSJeremy L Thompson     if (grad) {
1770c4e3f59bSSebastian Grimberg       CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp));
17711203703bSJeremy L Thompson       CeedCall(CeedScalarView("grad", "\t% 12.8f", q_comp * Q, P, grad, stream));
17727a982d89SJeremy L. Thompson     }
17731203703bSJeremy L Thompson     if (div) {
1774c4e3f59bSSebastian Grimberg       CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_DIV, &q_comp));
17751203703bSJeremy L Thompson       CeedCall(CeedScalarView("div", "\t% 12.8f", q_comp * Q, P, div, stream));
1776c4e3f59bSSebastian Grimberg     }
17771203703bSJeremy L Thompson     if (curl) {
1778c4e3f59bSSebastian Grimberg       CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_CURL, &q_comp));
17791203703bSJeremy L Thompson       CeedCall(CeedScalarView("curl", "\t% 12.8f", q_comp * Q, P, curl, stream));
178050c301a5SRezgar Shakeri     }
178150c301a5SRezgar Shakeri   }
1782e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
17837a982d89SJeremy L. Thompson }
17847a982d89SJeremy L. Thompson 
17857a982d89SJeremy L. Thompson /**
1786db2becc9SJeremy L Thompson   @brief Check input vector dimensions for CeedBasisApply[Add]
17877a982d89SJeremy L. Thompson 
1788ca94c3ddSJeremy L Thompson   @param[in]  basis     `CeedBasis` to evaluate
1789ea61e9acSJeremy L Thompson   @param[in]  num_elem  The number of elements to apply the basis evaluation to;
1790ca94c3ddSJeremy L Thompson                           the backend will specify the ordering in @ref CeedElemRestrictionCreate()
1791ca94c3ddSJeremy L Thompson   @param[in]  t_mode    @ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature points;
1792ca94c3ddSJeremy L Thompson                           @ref CEED_TRANSPOSE to apply the transpose, mapping from quadrature points to nodes
1793ca94c3ddSJeremy L Thompson   @param[in]  eval_mode @ref CEED_EVAL_NONE to use values directly,
1794ca94c3ddSJeremy L Thompson                           @ref CEED_EVAL_INTERP to use interpolated values,
1795ca94c3ddSJeremy L Thompson                           @ref CEED_EVAL_GRAD to use gradients,
1796ca94c3ddSJeremy L Thompson                           @ref CEED_EVAL_DIV to use divergence,
1797ca94c3ddSJeremy L Thompson                           @ref CEED_EVAL_CURL to use curl,
1798ca94c3ddSJeremy L Thompson                           @ref CEED_EVAL_WEIGHT to use quadrature weights
1799ca94c3ddSJeremy L Thompson   @param[in]  u         Input `CeedVector`
1800ca94c3ddSJeremy L Thompson   @param[out] v         Output `CeedVector`
18017a982d89SJeremy L. Thompson 
18027a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
18037a982d89SJeremy L. Thompson 
1804db2becc9SJeremy L Thompson   @ref Developer
18057a982d89SJeremy L. Thompson **/
1806db2becc9SJeremy L Thompson static int CeedBasisApplyCheckDims(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
1807c4e3f59bSSebastian Grimberg   CeedInt  dim, num_comp, q_comp, num_nodes, num_qpts;
18081c66c397SJeremy L Thompson   CeedSize u_length = 0, v_length;
18091203703bSJeremy L Thompson   Ceed     ceed;
18101c66c397SJeremy L Thompson 
18111203703bSJeremy L Thompson   CeedCall(CeedBasisGetCeed(basis, &ceed));
18122b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetDimension(basis, &dim));
18132b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
1814c4e3f59bSSebastian Grimberg   CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp));
18152b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumNodes(basis, &num_nodes));
18162b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts));
18172b730f8bSJeremy L Thompson   CeedCall(CeedVectorGetLength(v, &v_length));
1818c8c3fa7dSJeremy L Thompson   if (u) CeedCall(CeedVectorGetLength(u, &u_length));
18197a982d89SJeremy L. Thompson 
1820e15f9bd0SJeremy L Thompson   // Check vector lengths to prevent out of bounds issues
182199e754f0SJeremy L Thompson   bool has_good_dims = true;
1822d1d35e2fSjeremylt   switch (eval_mode) {
1823e15f9bd0SJeremy L Thompson     case CEED_EVAL_NONE:
18242b730f8bSJeremy L Thompson     case CEED_EVAL_INTERP:
18252b730f8bSJeremy L Thompson     case CEED_EVAL_GRAD:
1826c4e3f59bSSebastian Grimberg     case CEED_EVAL_DIV:
1827c4e3f59bSSebastian Grimberg     case CEED_EVAL_CURL:
1828*19a04db8SJeremy L Thompson       has_good_dims = ((t_mode == CEED_TRANSPOSE && u_length >= (CeedSize)num_elem * (CeedSize)num_comp * (CeedSize)num_qpts * (CeedSize)q_comp &&
1829*19a04db8SJeremy L Thompson                         v_length >= (CeedSize)num_elem * (CeedSize)num_comp * (CeedSize)num_nodes) ||
1830*19a04db8SJeremy L Thompson                        (t_mode == CEED_NOTRANSPOSE && v_length >= (CeedSize)num_elem * (CeedSize)num_qpts * (CeedSize)num_comp * (CeedSize)q_comp &&
1831*19a04db8SJeremy L Thompson                         u_length >= (CeedSize)num_elem * (CeedSize)num_comp * (CeedSize)num_nodes));
1832e15f9bd0SJeremy L Thompson       break;
1833e15f9bd0SJeremy L Thompson     case CEED_EVAL_WEIGHT:
1834*19a04db8SJeremy L Thompson       has_good_dims = v_length >= (CeedSize)num_elem * (CeedSize)num_qpts;
1835e15f9bd0SJeremy L Thompson       break;
1836e15f9bd0SJeremy L Thompson   }
183799e754f0SJeremy L Thompson   CeedCheck(has_good_dims, ceed, CEED_ERROR_DIMENSION, "Input/output vectors too short for basis and evaluation mode");
1838db2becc9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1839db2becc9SJeremy L Thompson }
1840e15f9bd0SJeremy L Thompson 
1841db2becc9SJeremy L Thompson /**
1842db2becc9SJeremy L Thompson   @brief Apply basis evaluation from nodes to quadrature points or vice versa
1843db2becc9SJeremy L Thompson 
1844db2becc9SJeremy L Thompson   @param[in]  basis     `CeedBasis` to evaluate
1845db2becc9SJeremy L Thompson   @param[in]  num_elem  The number of elements to apply the basis evaluation to;
1846db2becc9SJeremy L Thompson                           the backend will specify the ordering in @ref CeedElemRestrictionCreate()
1847db2becc9SJeremy L Thompson   @param[in]  t_mode    @ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature points;
1848db2becc9SJeremy L Thompson                           @ref CEED_TRANSPOSE to apply the transpose, mapping from quadrature points to nodes
1849db2becc9SJeremy L Thompson   @param[in]  eval_mode @ref CEED_EVAL_NONE to use values directly,
1850db2becc9SJeremy L Thompson                           @ref CEED_EVAL_INTERP to use interpolated values,
1851db2becc9SJeremy L Thompson                           @ref CEED_EVAL_GRAD to use gradients,
1852db2becc9SJeremy L Thompson                           @ref CEED_EVAL_DIV to use divergence,
1853db2becc9SJeremy L Thompson                           @ref CEED_EVAL_CURL to use curl,
1854db2becc9SJeremy L Thompson                           @ref CEED_EVAL_WEIGHT to use quadrature weights
1855db2becc9SJeremy L Thompson   @param[in]  u         Input `CeedVector`
1856db2becc9SJeremy L Thompson   @param[out] v         Output `CeedVector`
1857db2becc9SJeremy L Thompson 
1858db2becc9SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1859db2becc9SJeremy L Thompson 
1860db2becc9SJeremy L Thompson   @ref User
1861db2becc9SJeremy L Thompson **/
1862db2becc9SJeremy L Thompson int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
1863db2becc9SJeremy L Thompson   CeedCall(CeedBasisApplyCheckDims(basis, num_elem, t_mode, eval_mode, u, v));
1864db2becc9SJeremy L Thompson   CeedCheck(basis->Apply, CeedBasisReturnCeed(basis), CEED_ERROR_UNSUPPORTED, "Backend does not support CeedBasisApply");
18652b730f8bSJeremy L Thompson   CeedCall(basis->Apply(basis, num_elem, t_mode, eval_mode, u, v));
1866e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
18677a982d89SJeremy L. Thompson }
18687a982d89SJeremy L. Thompson 
18697a982d89SJeremy L. Thompson /**
1870db2becc9SJeremy L Thompson   @brief Apply basis evaluation from quadrature points to nodes and sum into target vector
1871db2becc9SJeremy L Thompson 
1872db2becc9SJeremy L Thompson   @param[in]  basis     `CeedBasis` to evaluate
1873db2becc9SJeremy L Thompson   @param[in]  num_elem  The number of elements to apply the basis evaluation to;
1874db2becc9SJeremy L Thompson                           the backend will specify the ordering in @ref CeedElemRestrictionCreate()
1875db2becc9SJeremy L Thompson   @param[in]  t_mode    @ref CEED_TRANSPOSE to apply the transpose, mapping from quadrature points to nodes;
1876db2becc9SJeremy L Thompson                            @ref CEED_NOTRANSPOSE is not valid for `CeedBasisApplyAdd()`
1877db2becc9SJeremy L Thompson   @param[in]  eval_mode @ref CEED_EVAL_NONE to use values directly,
1878db2becc9SJeremy L Thompson                           @ref CEED_EVAL_INTERP to use interpolated values,
1879db2becc9SJeremy L Thompson                           @ref CEED_EVAL_GRAD to use gradients,
1880db2becc9SJeremy L Thompson                           @ref CEED_EVAL_DIV to use divergence,
1881db2becc9SJeremy L Thompson                           @ref CEED_EVAL_CURL to use curl,
1882db2becc9SJeremy L Thompson                           @ref CEED_EVAL_WEIGHT to use quadrature weights
1883db2becc9SJeremy L Thompson   @param[in]  u         Input `CeedVector`
1884db2becc9SJeremy L Thompson   @param[out] v         Output `CeedVector` to sum into
1885db2becc9SJeremy L Thompson 
1886db2becc9SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1887db2becc9SJeremy L Thompson 
1888db2becc9SJeremy L Thompson   @ref User
1889db2becc9SJeremy L Thompson **/
1890db2becc9SJeremy L Thompson int CeedBasisApplyAdd(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
1891db2becc9SJeremy L Thompson   CeedCheck(t_mode == CEED_TRANSPOSE, CeedBasisReturnCeed(basis), CEED_ERROR_UNSUPPORTED, "CeedBasisApplyAdd only supports CEED_TRANSPOSE");
1892db2becc9SJeremy L Thompson   CeedCall(CeedBasisApplyCheckDims(basis, num_elem, t_mode, eval_mode, u, v));
1893db2becc9SJeremy L Thompson   CeedCheck(basis->ApplyAdd, CeedBasisReturnCeed(basis), CEED_ERROR_UNSUPPORTED, "Backend does not implement CeedBasisApplyAdd");
1894db2becc9SJeremy L Thompson   CeedCall(basis->ApplyAdd(basis, num_elem, t_mode, eval_mode, u, v));
1895db2becc9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1896db2becc9SJeremy L Thompson }
1897db2becc9SJeremy L Thompson 
1898db2becc9SJeremy L Thompson /**
1899db2becc9SJeremy L Thompson   @brief Apply basis evaluation from nodes to arbitrary points
1900db2becc9SJeremy L Thompson 
1901db2becc9SJeremy L Thompson   @param[in]  basis      `CeedBasis` to evaluate
1902db2becc9SJeremy L Thompson   @param[in]  num_elem   The number of elements to apply the basis evaluation to;
1903db2becc9SJeremy L Thompson                           the backend will specify the ordering in @ref CeedElemRestrictionCreate()
1904db2becc9SJeremy L Thompson   @param[in]  num_points Array of the number of points to apply the basis evaluation to in each element, size `num_elem`
1905db2becc9SJeremy L Thompson   @param[in]  t_mode     @ref CEED_NOTRANSPOSE to evaluate from nodes to points;
1906db2becc9SJeremy L Thompson                            @ref CEED_TRANSPOSE to apply the transpose, mapping from points to nodes
1907db2becc9SJeremy L Thompson   @param[in]  eval_mode  @ref CEED_EVAL_INTERP to use interpolated values,
1908db2becc9SJeremy L Thompson                            @ref CEED_EVAL_GRAD to use gradients,
1909db2becc9SJeremy L Thompson                            @ref CEED_EVAL_WEIGHT to use quadrature weights
1910db2becc9SJeremy L Thompson   @param[in]  x_ref      `CeedVector` holding reference coordinates of each point
1911db2becc9SJeremy L Thompson   @param[in]  u          Input `CeedVector`, of length `num_nodes * num_comp` for @ref CEED_NOTRANSPOSE
1912db2becc9SJeremy L Thompson   @param[out] v          Output `CeedVector`, of length `num_points * num_q_comp` for @ref CEED_NOTRANSPOSE with @ref CEED_EVAL_INTERP
1913db2becc9SJeremy L Thompson 
1914db2becc9SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1915db2becc9SJeremy L Thompson 
1916db2becc9SJeremy L Thompson   @ref User
1917db2becc9SJeremy L Thompson **/
1918db2becc9SJeremy L Thompson int CeedBasisApplyAtPoints(CeedBasis basis, CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode, CeedEvalMode eval_mode,
1919db2becc9SJeremy L Thompson                            CeedVector x_ref, CeedVector u, CeedVector v) {
1920db2becc9SJeremy L Thompson   CeedCall(CeedBasisApplyAtPointsCheckDims(basis, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
1921db2becc9SJeremy L Thompson   if (basis->ApplyAtPoints) {
1922db2becc9SJeremy L Thompson     CeedCall(basis->ApplyAtPoints(basis, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
1923db2becc9SJeremy L Thompson   } else {
1924db2becc9SJeremy L Thompson     CeedCall(CeedBasisApplyAtPoints_Core(basis, false, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
1925db2becc9SJeremy L Thompson   }
1926db2becc9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1927db2becc9SJeremy L Thompson }
1928db2becc9SJeremy L Thompson 
1929db2becc9SJeremy L Thompson /**
1930db2becc9SJeremy L Thompson   @brief Apply basis evaluation from nodes to arbitrary points and sum into target vector
1931db2becc9SJeremy L Thompson 
1932db2becc9SJeremy L Thompson   @param[in]  basis      `CeedBasis` to evaluate
1933db2becc9SJeremy L Thompson   @param[in]  num_elem   The number of elements to apply the basis evaluation to;
1934db2becc9SJeremy L Thompson                           the backend will specify the ordering in @ref CeedElemRestrictionCreate()
1935db2becc9SJeremy L Thompson   @param[in]  num_points Array of the number of points to apply the basis evaluation to in each element, size `num_elem`
1936db2becc9SJeremy L Thompson   @param[in]  t_mode     @ref CEED_NOTRANSPOSE to evaluate from nodes to points;
1937db2becc9SJeremy L Thompson                            @ref CEED_NOTRANSPOSE is not valid for `CeedBasisApplyAddAtPoints()`
1938db2becc9SJeremy L Thompson   @param[in]  eval_mode  @ref CEED_EVAL_INTERP to use interpolated values,
1939db2becc9SJeremy L Thompson                            @ref CEED_EVAL_GRAD to use gradients,
1940db2becc9SJeremy L Thompson                            @ref CEED_EVAL_WEIGHT to use quadrature weights
1941db2becc9SJeremy L Thompson   @param[in]  x_ref      `CeedVector` holding reference coordinates of each point
1942db2becc9SJeremy L Thompson   @param[in]  u          Input `CeedVector`, of length `num_nodes * num_comp` for @ref CEED_NOTRANSPOSE
1943db2becc9SJeremy L Thompson   @param[out] v          Output `CeedVector`, of length `num_points * num_q_comp` for @ref CEED_NOTRANSPOSE with @ref CEED_EVAL_INTERP
1944db2becc9SJeremy L Thompson 
1945db2becc9SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1946db2becc9SJeremy L Thompson 
1947db2becc9SJeremy L Thompson   @ref User
1948db2becc9SJeremy L Thompson **/
1949db2becc9SJeremy L Thompson int CeedBasisApplyAddAtPoints(CeedBasis basis, CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode, CeedEvalMode eval_mode,
1950db2becc9SJeremy L Thompson                               CeedVector x_ref, CeedVector u, CeedVector v) {
1951db2becc9SJeremy L Thompson   CeedCheck(t_mode == CEED_TRANSPOSE, CeedBasisReturnCeed(basis), CEED_ERROR_UNSUPPORTED, "CeedBasisApplyAddAtPoints only supports CEED_TRANSPOSE");
1952db2becc9SJeremy L Thompson   CeedCall(CeedBasisApplyAtPointsCheckDims(basis, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
1953db2becc9SJeremy L Thompson   if (basis->ApplyAddAtPoints) {
1954db2becc9SJeremy L Thompson     CeedCall(basis->ApplyAddAtPoints(basis, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
1955db2becc9SJeremy L Thompson   } else {
1956db2becc9SJeremy L Thompson     CeedCall(CeedBasisApplyAtPoints_Core(basis, true, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
1957db2becc9SJeremy L Thompson   }
1958db2becc9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1959db2becc9SJeremy L Thompson }
1960db2becc9SJeremy L Thompson 
1961db2becc9SJeremy L Thompson /**
19626e536b99SJeremy L Thompson   @brief Get the `Ceed` associated with a `CeedBasis`
1963b7c9bbdaSJeremy L Thompson 
1964ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
1965ca94c3ddSJeremy L Thompson   @param[out] ceed  Variable to store `Ceed`
1966b7c9bbdaSJeremy L Thompson 
1967b7c9bbdaSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1968b7c9bbdaSJeremy L Thompson 
1969b7c9bbdaSJeremy L Thompson   @ref Advanced
1970b7c9bbdaSJeremy L Thompson **/
1971b7c9bbdaSJeremy L Thompson int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) {
19726e536b99SJeremy L Thompson   *ceed = CeedBasisReturnCeed(basis);
1973b7c9bbdaSJeremy L Thompson   return CEED_ERROR_SUCCESS;
1974b7c9bbdaSJeremy L Thompson }
1975b7c9bbdaSJeremy L Thompson 
1976b7c9bbdaSJeremy L Thompson /**
19776e536b99SJeremy L Thompson   @brief Return the `Ceed` associated with a `CeedBasis`
19786e536b99SJeremy L Thompson 
19796e536b99SJeremy L Thompson   @param[in]  basis `CeedBasis`
19806e536b99SJeremy L Thompson 
19816e536b99SJeremy L Thompson   @return `Ceed` associated with the `basis`
19826e536b99SJeremy L Thompson 
19836e536b99SJeremy L Thompson   @ref Advanced
19846e536b99SJeremy L Thompson **/
19856e536b99SJeremy L Thompson Ceed CeedBasisReturnCeed(CeedBasis basis) { return basis->ceed; }
19866e536b99SJeremy L Thompson 
19876e536b99SJeremy L Thompson /**
1988ca94c3ddSJeremy L Thompson   @brief Get dimension for given `CeedBasis`
19899d007619Sjeremylt 
1990ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
19919d007619Sjeremylt   @param[out] dim   Variable to store dimension of basis
19929d007619Sjeremylt 
19939d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
19949d007619Sjeremylt 
1995b7c9bbdaSJeremy L Thompson   @ref Advanced
19969d007619Sjeremylt **/
19979d007619Sjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) {
19989d007619Sjeremylt   *dim = basis->dim;
1999e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
20009d007619Sjeremylt }
20019d007619Sjeremylt 
20029d007619Sjeremylt /**
2003ca94c3ddSJeremy L Thompson   @brief Get topology for given `CeedBasis`
2004d99fa3c5SJeremy L Thompson 
2005ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
2006d99fa3c5SJeremy L Thompson   @param[out] topo  Variable to store topology of basis
2007d99fa3c5SJeremy L Thompson 
2008d99fa3c5SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
2009d99fa3c5SJeremy L Thompson 
2010b7c9bbdaSJeremy L Thompson   @ref Advanced
2011d99fa3c5SJeremy L Thompson **/
2012d99fa3c5SJeremy L Thompson int CeedBasisGetTopology(CeedBasis basis, CeedElemTopology *topo) {
2013d99fa3c5SJeremy L Thompson   *topo = basis->topo;
2014e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2015d99fa3c5SJeremy L Thompson }
2016d99fa3c5SJeremy L Thompson 
2017d99fa3c5SJeremy L Thompson /**
2018ca94c3ddSJeremy L Thompson   @brief Get number of components for given `CeedBasis`
20199d007619Sjeremylt 
2020ca94c3ddSJeremy L Thompson   @param[in]  basis    `CeedBasis`
2021ca94c3ddSJeremy L Thompson   @param[out] num_comp Variable to store number of components
20229d007619Sjeremylt 
20239d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
20249d007619Sjeremylt 
2025b7c9bbdaSJeremy L Thompson   @ref Advanced
20269d007619Sjeremylt **/
2027d1d35e2fSjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *num_comp) {
2028d1d35e2fSjeremylt   *num_comp = basis->num_comp;
2029e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
20309d007619Sjeremylt }
20319d007619Sjeremylt 
20329d007619Sjeremylt /**
2033ca94c3ddSJeremy L Thompson   @brief Get total number of nodes (in `dim` dimensions) of a `CeedBasis`
20349d007619Sjeremylt 
2035ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
20369d007619Sjeremylt   @param[out] P     Variable to store number of nodes
20379d007619Sjeremylt 
20389d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
20399d007619Sjeremylt 
20409d007619Sjeremylt   @ref Utility
20419d007619Sjeremylt **/
20429d007619Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) {
20439d007619Sjeremylt   *P = basis->P;
2044e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
20459d007619Sjeremylt }
20469d007619Sjeremylt 
20479d007619Sjeremylt /**
2048ca94c3ddSJeremy L Thompson   @brief Get total number of nodes (in 1 dimension) of a `CeedBasis`
20499d007619Sjeremylt 
2050ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
2051d1d35e2fSjeremylt   @param[out] P_1d  Variable to store number of nodes
20529d007619Sjeremylt 
20539d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
20549d007619Sjeremylt 
2055b7c9bbdaSJeremy L Thompson   @ref Advanced
20569d007619Sjeremylt **/
2057d1d35e2fSjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P_1d) {
20586e536b99SJeremy L Thompson   CeedCheck(basis->is_tensor_basis, CeedBasisReturnCeed(basis), CEED_ERROR_MINOR, "Cannot supply P_1d for non-tensor CeedBasis");
2059d1d35e2fSjeremylt   *P_1d = basis->P_1d;
2060e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
20619d007619Sjeremylt }
20629d007619Sjeremylt 
20639d007619Sjeremylt /**
2064ca94c3ddSJeremy L Thompson   @brief Get total number of quadrature points (in `dim` dimensions) of a `CeedBasis`
20659d007619Sjeremylt 
2066ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
20679d007619Sjeremylt   @param[out] Q     Variable to store number of quadrature points
20689d007619Sjeremylt 
20699d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
20709d007619Sjeremylt 
20719d007619Sjeremylt   @ref Utility
20729d007619Sjeremylt **/
20739d007619Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) {
20749d007619Sjeremylt   *Q = basis->Q;
2075e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
20769d007619Sjeremylt }
20779d007619Sjeremylt 
20789d007619Sjeremylt /**
2079ca94c3ddSJeremy L Thompson   @brief Get total number of quadrature points (in 1 dimension) of a `CeedBasis`
20809d007619Sjeremylt 
2081ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
2082d1d35e2fSjeremylt   @param[out] Q_1d  Variable to store number of quadrature points
20839d007619Sjeremylt 
20849d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
20859d007619Sjeremylt 
2086b7c9bbdaSJeremy L Thompson   @ref Advanced
20879d007619Sjeremylt **/
2088d1d35e2fSjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q_1d) {
20896e536b99SJeremy L Thompson   CeedCheck(basis->is_tensor_basis, CeedBasisReturnCeed(basis), CEED_ERROR_MINOR, "Cannot supply Q_1d for non-tensor CeedBasis");
2090d1d35e2fSjeremylt   *Q_1d = basis->Q_1d;
2091e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
20929d007619Sjeremylt }
20939d007619Sjeremylt 
20949d007619Sjeremylt /**
2095ca94c3ddSJeremy L Thompson   @brief Get reference coordinates of quadrature points (in `dim` dimensions) of a `CeedBasis`
20969d007619Sjeremylt 
2097ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
2098d1d35e2fSjeremylt   @param[out] q_ref Variable to store reference coordinates of quadrature points
20999d007619Sjeremylt 
21009d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
21019d007619Sjeremylt 
2102b7c9bbdaSJeremy L Thompson   @ref Advanced
21039d007619Sjeremylt **/
2104d1d35e2fSjeremylt int CeedBasisGetQRef(CeedBasis basis, const CeedScalar **q_ref) {
2105d1d35e2fSjeremylt   *q_ref = basis->q_ref_1d;
2106e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
21079d007619Sjeremylt }
21089d007619Sjeremylt 
21099d007619Sjeremylt /**
2110ca94c3ddSJeremy L Thompson   @brief Get quadrature weights of quadrature points (in `dim` dimensions) of a `CeedBasis`
21119d007619Sjeremylt 
2112ca94c3ddSJeremy L Thompson   @param[in]  basis    `CeedBasis`
2113d1d35e2fSjeremylt   @param[out] q_weight Variable to store quadrature weights
21149d007619Sjeremylt 
21159d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
21169d007619Sjeremylt 
2117b7c9bbdaSJeremy L Thompson   @ref Advanced
21189d007619Sjeremylt **/
2119d1d35e2fSjeremylt int CeedBasisGetQWeights(CeedBasis basis, const CeedScalar **q_weight) {
2120d1d35e2fSjeremylt   *q_weight = basis->q_weight_1d;
2121e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
21229d007619Sjeremylt }
21239d007619Sjeremylt 
21249d007619Sjeremylt /**
2125ca94c3ddSJeremy L Thompson   @brief Get interpolation matrix of a `CeedBasis`
21269d007619Sjeremylt 
2127ca94c3ddSJeremy L Thompson   @param[in]  basis  `CeedBasis`
21289d007619Sjeremylt   @param[out] interp Variable to store interpolation matrix
21299d007619Sjeremylt 
21309d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
21319d007619Sjeremylt 
2132b7c9bbdaSJeremy L Thompson   @ref Advanced
21339d007619Sjeremylt **/
21346c58de82SJeremy L Thompson int CeedBasisGetInterp(CeedBasis basis, const CeedScalar **interp) {
21356402da51SJeremy L Thompson   if (!basis->interp && basis->is_tensor_basis) {
21369d007619Sjeremylt     // Allocate
21372b730f8bSJeremy L Thompson     CeedCall(CeedMalloc(basis->Q * basis->P, &basis->interp));
21389d007619Sjeremylt 
21399d007619Sjeremylt     // Initialize
21402b730f8bSJeremy L Thompson     for (CeedInt i = 0; i < basis->Q * basis->P; i++) basis->interp[i] = 1.0;
21419d007619Sjeremylt 
21429d007619Sjeremylt     // Calculate
21432b730f8bSJeremy L Thompson     for (CeedInt d = 0; d < basis->dim; d++) {
21442b730f8bSJeremy L Thompson       for (CeedInt qpt = 0; qpt < basis->Q; qpt++) {
21459d007619Sjeremylt         for (CeedInt node = 0; node < basis->P; node++) {
2146d1d35e2fSjeremylt           CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d;
2147d1d35e2fSjeremylt           CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d;
21481c66c397SJeremy L Thompson 
2149d1d35e2fSjeremylt           basis->interp[qpt * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p];
21509d007619Sjeremylt         }
21519d007619Sjeremylt       }
21522b730f8bSJeremy L Thompson     }
21532b730f8bSJeremy L Thompson   }
21549d007619Sjeremylt   *interp = basis->interp;
2155e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
21569d007619Sjeremylt }
21579d007619Sjeremylt 
21589d007619Sjeremylt /**
2159ca94c3ddSJeremy L Thompson   @brief Get 1D interpolation matrix of a tensor product `CeedBasis`
21609d007619Sjeremylt 
2161ca94c3ddSJeremy L Thompson   @param[in]  basis     `CeedBasis`
2162d1d35e2fSjeremylt   @param[out] interp_1d Variable to store interpolation matrix
21639d007619Sjeremylt 
21649d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
21659d007619Sjeremylt 
21669d007619Sjeremylt   @ref Backend
21679d007619Sjeremylt **/
2168d1d35e2fSjeremylt int CeedBasisGetInterp1D(CeedBasis basis, const CeedScalar **interp_1d) {
21691203703bSJeremy L Thompson   bool is_tensor_basis;
21701203703bSJeremy L Thompson 
21711203703bSJeremy L Thompson   CeedCall(CeedBasisIsTensor(basis, &is_tensor_basis));
21726e536b99SJeremy L Thompson   CeedCheck(is_tensor_basis, CeedBasisReturnCeed(basis), CEED_ERROR_MINOR, "CeedBasis is not a tensor product CeedBasis");
2173d1d35e2fSjeremylt   *interp_1d = basis->interp_1d;
2174e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
21759d007619Sjeremylt }
21769d007619Sjeremylt 
21779d007619Sjeremylt /**
2178ca94c3ddSJeremy L Thompson   @brief Get gradient matrix of a `CeedBasis`
21799d007619Sjeremylt 
2180ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
21819d007619Sjeremylt   @param[out] grad  Variable to store gradient matrix
21829d007619Sjeremylt 
21839d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
21849d007619Sjeremylt 
2185b7c9bbdaSJeremy L Thompson   @ref Advanced
21869d007619Sjeremylt **/
21876c58de82SJeremy L Thompson int CeedBasisGetGrad(CeedBasis basis, const CeedScalar **grad) {
21886402da51SJeremy L Thompson   if (!basis->grad && basis->is_tensor_basis) {
21899d007619Sjeremylt     // Allocate
21902b730f8bSJeremy L Thompson     CeedCall(CeedMalloc(basis->dim * basis->Q * basis->P, &basis->grad));
21919d007619Sjeremylt 
21929d007619Sjeremylt     // Initialize
21932b730f8bSJeremy L Thompson     for (CeedInt i = 0; i < basis->dim * basis->Q * basis->P; i++) basis->grad[i] = 1.0;
21949d007619Sjeremylt 
21959d007619Sjeremylt     // Calculate
21962b730f8bSJeremy L Thompson     for (CeedInt d = 0; d < basis->dim; d++) {
21972b730f8bSJeremy L Thompson       for (CeedInt i = 0; i < basis->dim; i++) {
21982b730f8bSJeremy L Thompson         for (CeedInt qpt = 0; qpt < basis->Q; qpt++) {
21999d007619Sjeremylt           for (CeedInt node = 0; node < basis->P; node++) {
2200d1d35e2fSjeremylt             CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d;
2201d1d35e2fSjeremylt             CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d;
22021c66c397SJeremy L Thompson 
22032b730f8bSJeremy L Thompson             if (i == d) basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->grad_1d[q * basis->P_1d + p];
22042b730f8bSJeremy L Thompson             else basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p];
22052b730f8bSJeremy L Thompson           }
22062b730f8bSJeremy L Thompson         }
22072b730f8bSJeremy L Thompson       }
22089d007619Sjeremylt     }
22099d007619Sjeremylt   }
22109d007619Sjeremylt   *grad = basis->grad;
2211e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
22129d007619Sjeremylt }
22139d007619Sjeremylt 
22149d007619Sjeremylt /**
2215ca94c3ddSJeremy L Thompson   @brief Get 1D gradient matrix of a tensor product `CeedBasis`
22169d007619Sjeremylt 
2217ca94c3ddSJeremy L Thompson   @param[in]  basis   `CeedBasis`
2218d1d35e2fSjeremylt   @param[out] grad_1d Variable to store gradient matrix
22199d007619Sjeremylt 
22209d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
22219d007619Sjeremylt 
2222b7c9bbdaSJeremy L Thompson   @ref Advanced
22239d007619Sjeremylt **/
2224d1d35e2fSjeremylt int CeedBasisGetGrad1D(CeedBasis basis, const CeedScalar **grad_1d) {
22251203703bSJeremy L Thompson   bool is_tensor_basis;
22261203703bSJeremy L Thompson 
22271203703bSJeremy L Thompson   CeedCall(CeedBasisIsTensor(basis, &is_tensor_basis));
22286e536b99SJeremy L Thompson   CeedCheck(is_tensor_basis, CeedBasisReturnCeed(basis), CEED_ERROR_MINOR, "CeedBasis is not a tensor product CeedBasis");
2229d1d35e2fSjeremylt   *grad_1d = basis->grad_1d;
2230e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
22319d007619Sjeremylt }
22329d007619Sjeremylt 
22339d007619Sjeremylt /**
2234ca94c3ddSJeremy L Thompson   @brief Get divergence matrix of a `CeedBasis`
223550c301a5SRezgar Shakeri 
2236ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
223750c301a5SRezgar Shakeri   @param[out] div   Variable to store divergence matrix
223850c301a5SRezgar Shakeri 
223950c301a5SRezgar Shakeri   @return An error code: 0 - success, otherwise - failure
224050c301a5SRezgar Shakeri 
224150c301a5SRezgar Shakeri   @ref Advanced
224250c301a5SRezgar Shakeri **/
224350c301a5SRezgar Shakeri int CeedBasisGetDiv(CeedBasis basis, const CeedScalar **div) {
224450c301a5SRezgar Shakeri   *div = basis->div;
224550c301a5SRezgar Shakeri   return CEED_ERROR_SUCCESS;
224650c301a5SRezgar Shakeri }
224750c301a5SRezgar Shakeri 
224850c301a5SRezgar Shakeri /**
2249ca94c3ddSJeremy L Thompson   @brief Get curl matrix of a `CeedBasis`
2250c4e3f59bSSebastian Grimberg 
2251ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
2252c4e3f59bSSebastian Grimberg   @param[out] curl  Variable to store curl matrix
2253c4e3f59bSSebastian Grimberg 
2254c4e3f59bSSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
2255c4e3f59bSSebastian Grimberg 
2256c4e3f59bSSebastian Grimberg   @ref Advanced
2257c4e3f59bSSebastian Grimberg **/
2258c4e3f59bSSebastian Grimberg int CeedBasisGetCurl(CeedBasis basis, const CeedScalar **curl) {
2259c4e3f59bSSebastian Grimberg   *curl = basis->curl;
2260c4e3f59bSSebastian Grimberg   return CEED_ERROR_SUCCESS;
2261c4e3f59bSSebastian Grimberg }
2262c4e3f59bSSebastian Grimberg 
2263c4e3f59bSSebastian Grimberg /**
2264ca94c3ddSJeremy L Thompson   @brief Destroy a @ref  CeedBasis
22657a982d89SJeremy L. Thompson 
2266ca94c3ddSJeremy L Thompson   @param[in,out] basis `CeedBasis` to destroy
22677a982d89SJeremy L. Thompson 
22687a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
22697a982d89SJeremy L. Thompson 
22707a982d89SJeremy L. Thompson   @ref User
22717a982d89SJeremy L. Thompson **/
22727a982d89SJeremy L. Thompson int CeedBasisDestroy(CeedBasis *basis) {
2273356036faSJeremy L Thompson   if (!*basis || *basis == CEED_BASIS_NONE || --(*basis)->ref_count > 0) {
2274ad6481ceSJeremy L Thompson     *basis = NULL;
2275ad6481ceSJeremy L Thompson     return CEED_ERROR_SUCCESS;
2276ad6481ceSJeremy L Thompson   }
22772b730f8bSJeremy L Thompson   if ((*basis)->Destroy) CeedCall((*basis)->Destroy(*basis));
22789831d45aSJeremy L Thompson   CeedCall(CeedTensorContractDestroy(&(*basis)->contract));
2279c4e3f59bSSebastian Grimberg   CeedCall(CeedFree(&(*basis)->q_ref_1d));
2280c4e3f59bSSebastian Grimberg   CeedCall(CeedFree(&(*basis)->q_weight_1d));
22812b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->interp));
22822b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->interp_1d));
22832b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->grad));
22842b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->grad_1d));
2285c4e3f59bSSebastian Grimberg   CeedCall(CeedFree(&(*basis)->div));
2286c4e3f59bSSebastian Grimberg   CeedCall(CeedFree(&(*basis)->curl));
2287c8c3fa7dSJeremy L Thompson   CeedCall(CeedVectorDestroy(&(*basis)->vec_chebyshev));
2288c8c3fa7dSJeremy L Thompson   CeedCall(CeedBasisDestroy(&(*basis)->basis_chebyshev));
22892b730f8bSJeremy L Thompson   CeedCall(CeedDestroy(&(*basis)->ceed));
22902b730f8bSJeremy L Thompson   CeedCall(CeedFree(basis));
2291e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
22927a982d89SJeremy L. Thompson }
22937a982d89SJeremy L. Thompson 
22947a982d89SJeremy L. Thompson /**
2295b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre quadrature
2296b11c1e72Sjeremylt 
2297ca94c3ddSJeremy L Thompson   @param[in]  Q           Number of quadrature points (integrates polynomials of degree `2*Q-1` exactly)
2298ca94c3ddSJeremy L Thompson   @param[out] q_ref_1d    Array of length `Q` to hold the abscissa on `[-1, 1]`
2299ca94c3ddSJeremy L Thompson   @param[out] q_weight_1d Array of length `Q` to hold the weights
2300b11c1e72Sjeremylt 
2301b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
2302dfdf5a53Sjeremylt 
2303dfdf5a53Sjeremylt   @ref Utility
2304b11c1e72Sjeremylt **/
23052b730f8bSJeremy L Thompson int CeedGaussQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) {
2306d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0 * atan(1.0);
23071c66c397SJeremy L Thompson 
2308d1d35e2fSjeremylt   // Build q_ref_1d, q_weight_1d
230992ae7e47SJeremy L Thompson   for (CeedInt i = 0; i <= Q / 2; i++) {
2310d7b241e6Sjeremylt     // Guess
2311d7b241e6Sjeremylt     xi = cos(PI * (CeedScalar)(2 * i + 1) / ((CeedScalar)(2 * Q)));
2312d7b241e6Sjeremylt     // Pn(xi)
2313d7b241e6Sjeremylt     P0 = 1.0;
2314d7b241e6Sjeremylt     P1 = xi;
2315d7b241e6Sjeremylt     P2 = 0.0;
231692ae7e47SJeremy L Thompson     for (CeedInt j = 2; j <= Q; j++) {
2317d7b241e6Sjeremylt       P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
2318d7b241e6Sjeremylt       P0 = P1;
2319d7b241e6Sjeremylt       P1 = P2;
2320d7b241e6Sjeremylt     }
2321d7b241e6Sjeremylt     // First Newton Step
2322d7b241e6Sjeremylt     dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
2323d7b241e6Sjeremylt     xi  = xi - P2 / dP2;
2324d7b241e6Sjeremylt     // Newton to convergence
232592ae7e47SJeremy L Thompson     for (CeedInt k = 0; k < 100 && fabs(P2) > 10 * CEED_EPSILON; k++) {
2326d7b241e6Sjeremylt       P0 = 1.0;
2327d7b241e6Sjeremylt       P1 = xi;
232892ae7e47SJeremy L Thompson       for (CeedInt j = 2; j <= Q; j++) {
2329d7b241e6Sjeremylt         P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
2330d7b241e6Sjeremylt         P0 = P1;
2331d7b241e6Sjeremylt         P1 = P2;
2332d7b241e6Sjeremylt       }
2333d7b241e6Sjeremylt       dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
2334d7b241e6Sjeremylt       xi  = xi - P2 / dP2;
2335d7b241e6Sjeremylt     }
2336d7b241e6Sjeremylt     // Save xi, wi
2337d7b241e6Sjeremylt     wi                     = 2.0 / ((1.0 - xi * xi) * dP2 * dP2);
2338d1d35e2fSjeremylt     q_weight_1d[i]         = wi;
2339d1d35e2fSjeremylt     q_weight_1d[Q - 1 - i] = wi;
2340d1d35e2fSjeremylt     q_ref_1d[i]            = -xi;
2341d1d35e2fSjeremylt     q_ref_1d[Q - 1 - i]    = xi;
2342d7b241e6Sjeremylt   }
2343e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2344d7b241e6Sjeremylt }
2345d7b241e6Sjeremylt 
2346b11c1e72Sjeremylt /**
2347b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre-Lobatto quadrature
2348b11c1e72Sjeremylt 
2349ca94c3ddSJeremy L Thompson   @param[in]  Q           Number of quadrature points (integrates polynomials of degree `2*Q-3` exactly)
2350ca94c3ddSJeremy L Thompson   @param[out] q_ref_1d    Array of length `Q` to hold the abscissa on `[-1, 1]`
2351ca94c3ddSJeremy L Thompson   @param[out] q_weight_1d Array of length `Q` to hold the weights
2352b11c1e72Sjeremylt 
2353b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
2354dfdf5a53Sjeremylt 
2355dfdf5a53Sjeremylt   @ref Utility
2356b11c1e72Sjeremylt **/
23572b730f8bSJeremy L Thompson int CeedLobattoQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) {
2358d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0 * atan(1.0);
23591c66c397SJeremy L Thompson 
2360d1d35e2fSjeremylt   // Build q_ref_1d, q_weight_1d
2361d7b241e6Sjeremylt   // Set endpoints
23626574a04fSJeremy L Thompson   CeedCheck(Q > 1, NULL, CEED_ERROR_DIMENSION, "Cannot create Lobatto quadrature with Q=%" CeedInt_FMT " < 2 points", Q);
2363d7b241e6Sjeremylt   wi = 2.0 / ((CeedScalar)(Q * (Q - 1)));
2364d1d35e2fSjeremylt   if (q_weight_1d) {
2365d1d35e2fSjeremylt     q_weight_1d[0]     = wi;
2366d1d35e2fSjeremylt     q_weight_1d[Q - 1] = wi;
2367d7b241e6Sjeremylt   }
2368d1d35e2fSjeremylt   q_ref_1d[0]     = -1.0;
2369d1d35e2fSjeremylt   q_ref_1d[Q - 1] = 1.0;
2370d7b241e6Sjeremylt   // Interior
237192ae7e47SJeremy L Thompson   for (CeedInt i = 1; i <= (Q - 1) / 2; i++) {
2372d7b241e6Sjeremylt     // Guess
2373d7b241e6Sjeremylt     xi = cos(PI * (CeedScalar)(i) / (CeedScalar)(Q - 1));
2374d7b241e6Sjeremylt     // Pn(xi)
2375d7b241e6Sjeremylt     P0 = 1.0;
2376d7b241e6Sjeremylt     P1 = xi;
2377d7b241e6Sjeremylt     P2 = 0.0;
237892ae7e47SJeremy L Thompson     for (CeedInt j = 2; j < Q; j++) {
2379d7b241e6Sjeremylt       P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
2380d7b241e6Sjeremylt       P0 = P1;
2381d7b241e6Sjeremylt       P1 = P2;
2382d7b241e6Sjeremylt     }
2383d7b241e6Sjeremylt     // First Newton step
2384d7b241e6Sjeremylt     dP2  = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
2385d7b241e6Sjeremylt     d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi);
2386d7b241e6Sjeremylt     xi   = xi - dP2 / d2P2;
2387d7b241e6Sjeremylt     // Newton to convergence
238892ae7e47SJeremy L Thompson     for (CeedInt k = 0; k < 100 && fabs(dP2) > 10 * CEED_EPSILON; k++) {
2389d7b241e6Sjeremylt       P0 = 1.0;
2390d7b241e6Sjeremylt       P1 = xi;
239192ae7e47SJeremy L Thompson       for (CeedInt j = 2; j < Q; j++) {
2392d7b241e6Sjeremylt         P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
2393d7b241e6Sjeremylt         P0 = P1;
2394d7b241e6Sjeremylt         P1 = P2;
2395d7b241e6Sjeremylt       }
2396d7b241e6Sjeremylt       dP2  = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
2397d7b241e6Sjeremylt       d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi);
2398d7b241e6Sjeremylt       xi   = xi - dP2 / d2P2;
2399d7b241e6Sjeremylt     }
2400d7b241e6Sjeremylt     // Save xi, wi
2401d7b241e6Sjeremylt     wi = 2.0 / (((CeedScalar)(Q * (Q - 1))) * P2 * P2);
2402d1d35e2fSjeremylt     if (q_weight_1d) {
2403d1d35e2fSjeremylt       q_weight_1d[i]         = wi;
2404d1d35e2fSjeremylt       q_weight_1d[Q - 1 - i] = wi;
2405d7b241e6Sjeremylt     }
2406d1d35e2fSjeremylt     q_ref_1d[i]         = -xi;
2407d1d35e2fSjeremylt     q_ref_1d[Q - 1 - i] = xi;
2408d7b241e6Sjeremylt   }
2409e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2410d7b241e6Sjeremylt }
2411d7b241e6Sjeremylt 
2412d7b241e6Sjeremylt /// @}
2413