xref: /libCEED/interface/ceed-basis.c (revision 3f919cbc23e59285304a4fc152ac9ef29c71e01d)
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) {
197e104ad11SJames Wright   bool    are_both_tensor;
1981c66c397SJeremy L Thompson   CeedInt Q, Q_to, Q_from, P_to, P_from;
1991c66c397SJeremy L Thompson 
200a76a04e7SJeremy L Thompson   // Check for compatible quadrature spaces
2012b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumQuadraturePoints(basis_to, &Q_to));
2022b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumQuadraturePoints(basis_from, &Q_from));
2039bc66399SJeremy L Thompson   CeedCheck(Q_to == Q_from, CeedBasisReturnCeed(basis_to), CEED_ERROR_DIMENSION,
2043f08121cSJeremy L Thompson             "Bases must have compatible quadrature spaces."
20523622755SJeremy L Thompson             " 'basis_from' has %" CeedInt_FMT " points and 'basis_to' has %" CeedInt_FMT,
2063f08121cSJeremy L Thompson             Q_from, Q_to);
2071c66c397SJeremy L Thompson   Q = Q_to;
208a76a04e7SJeremy L Thompson 
20914556e63SJeremy L Thompson   // Check for matching tensor or non-tensor
210e104ad11SJames Wright   {
211e104ad11SJames Wright     bool is_tensor_to, is_tensor_from;
212e104ad11SJames Wright 
2132b730f8bSJeremy L Thompson     CeedCall(CeedBasisIsTensor(basis_to, &is_tensor_to));
2142b730f8bSJeremy L Thompson     CeedCall(CeedBasisIsTensor(basis_from, &is_tensor_from));
215e104ad11SJames Wright     are_both_tensor = is_tensor_to && is_tensor_from;
216e104ad11SJames Wright   }
217e104ad11SJames Wright   if (are_both_tensor) {
2182b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_to));
2192b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_from));
2202b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumQuadraturePoints1D(basis_from, &Q));
2216574a04fSJeremy L Thompson   } else {
2222b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis_to, &P_to));
2232b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis_from, &P_from));
224a76a04e7SJeremy L Thompson   }
225a76a04e7SJeremy L Thompson 
22615ad3917SSebastian Grimberg   // Check for matching FE space
22715ad3917SSebastian Grimberg   CeedFESpace fe_space_to, fe_space_from;
2283f08121cSJeremy L Thompson 
22915ad3917SSebastian Grimberg   CeedCall(CeedBasisGetFESpace(basis_to, &fe_space_to));
23015ad3917SSebastian Grimberg   CeedCall(CeedBasisGetFESpace(basis_from, &fe_space_from));
2319bc66399SJeremy L Thompson   CeedCheck(fe_space_to == fe_space_from, CeedBasisReturnCeed(basis_to), CEED_ERROR_MINOR,
2323f08121cSJeremy L Thompson             "Bases must both be the same FE space type."
2333f08121cSJeremy L Thompson             " 'basis_from' is a %s and 'basis_to' is a %s",
2343f08121cSJeremy L Thompson             CeedFESpaces[fe_space_from], CeedFESpaces[fe_space_to]);
23515ad3917SSebastian Grimberg 
23614556e63SJeremy L Thompson   // Get source matrices
23715ad3917SSebastian Grimberg   CeedInt           dim, q_comp = 1;
2382247a93fSRezgar Shakeri   CeedScalar       *interp_to_inv, *interp_from;
2391c66c397SJeremy L Thompson   const CeedScalar *interp_to_source = NULL, *interp_from_source = NULL, *grad_from_source = NULL;
2401c66c397SJeremy L Thompson 
241b3ed00e5SJames Wright   CeedCall(CeedBasisGetDimension(basis_from, &dim));
242e104ad11SJames Wright   if (are_both_tensor) {
2432b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetInterp1D(basis_to, &interp_to_source));
2442b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetInterp1D(basis_from, &interp_from_source));
245a76a04e7SJeremy L Thompson   } else {
24615ad3917SSebastian Grimberg     CeedCall(CeedBasisGetNumQuadratureComponents(basis_from, CEED_EVAL_INTERP, &q_comp));
2472b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetInterp(basis_to, &interp_to_source));
2482b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetInterp(basis_from, &interp_from_source));
24915ad3917SSebastian Grimberg   }
25015ad3917SSebastian Grimberg   CeedCall(CeedMalloc(Q * P_from * q_comp, &interp_from));
25115ad3917SSebastian Grimberg   CeedCall(CeedCalloc(P_to * P_from, interp_project));
25215ad3917SSebastian Grimberg 
25315ad3917SSebastian Grimberg   // `grad_project = interp_to^+ * grad_from` is computed for the H^1 space case so the
254de05fbb2SSebastian Grimberg   // projection basis will have a gradient operation (allocated even if not H^1 for the
255de05fbb2SSebastian Grimberg   // basis construction later on)
25615ad3917SSebastian Grimberg   if (fe_space_to == CEED_FE_SPACE_H1) {
257e104ad11SJames Wright     if (are_both_tensor) {
25815ad3917SSebastian Grimberg       CeedCall(CeedBasisGetGrad1D(basis_from, &grad_from_source));
25915ad3917SSebastian Grimberg     } else {
2602b730f8bSJeremy L Thompson       CeedCall(CeedBasisGetGrad(basis_from, &grad_from_source));
261a76a04e7SJeremy L Thompson     }
262de05fbb2SSebastian Grimberg   }
263e104ad11SJames Wright   CeedCall(CeedCalloc(P_to * P_from * (are_both_tensor ? 1 : dim), grad_project));
26415ad3917SSebastian Grimberg 
2652247a93fSRezgar Shakeri   // Compute interp_to^+, pseudoinverse of interp_to
2662247a93fSRezgar Shakeri   CeedCall(CeedCalloc(Q * q_comp * P_to, &interp_to_inv));
2679bc66399SJeremy L Thompson   CeedCall(CeedMatrixPseudoinverse(CeedBasisReturnCeed(basis_to), interp_to_source, Q * q_comp, P_to, interp_to_inv));
26814556e63SJeremy L Thompson   // Build matrices
269e104ad11SJames Wright   CeedInt     num_matrices = 1 + (fe_space_to == CEED_FE_SPACE_H1) * (are_both_tensor ? 1 : dim);
27014556e63SJeremy L Thompson   CeedScalar *input_from[num_matrices], *output_project[num_matrices];
2711c66c397SJeremy L Thompson 
27214556e63SJeremy L Thompson   input_from[0]     = (CeedScalar *)interp_from_source;
27314556e63SJeremy L Thompson   output_project[0] = *interp_project;
27414556e63SJeremy L Thompson   for (CeedInt m = 1; m < num_matrices; m++) {
27514556e63SJeremy L Thompson     input_from[m]     = (CeedScalar *)&grad_from_source[(m - 1) * Q * P_from];
27602af4036SJeremy L Thompson     output_project[m] = &((*grad_project)[(m - 1) * P_to * P_from]);
27714556e63SJeremy L Thompson   }
27814556e63SJeremy L Thompson   for (CeedInt m = 0; m < num_matrices; m++) {
2792247a93fSRezgar Shakeri     // output_project = interp_to^+ * interp_from
28015ad3917SSebastian Grimberg     memcpy(interp_from, input_from[m], Q * P_from * q_comp * sizeof(input_from[m][0]));
2819bc66399SJeremy L Thompson     CeedCall(CeedMatrixMatrixMultiply(CeedBasisReturnCeed(basis_to), interp_to_inv, input_from[m], output_project[m], P_to, P_from, Q * q_comp));
2822247a93fSRezgar Shakeri     // Round zero to machine precision
2832247a93fSRezgar Shakeri     for (CeedInt i = 0; i < P_to * P_from; i++) {
2842247a93fSRezgar Shakeri       if (fabs(output_project[m][i]) < 10 * CEED_EPSILON) output_project[m][i] = 0.0;
285a76a04e7SJeremy L Thompson     }
28614556e63SJeremy L Thompson   }
28714556e63SJeremy L Thompson 
28814556e63SJeremy L Thompson   // Cleanup
2892247a93fSRezgar Shakeri   CeedCall(CeedFree(&interp_to_inv));
2902b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_from));
291a76a04e7SJeremy L Thompson   return CEED_ERROR_SUCCESS;
292a76a04e7SJeremy L Thompson }
293a76a04e7SJeremy L Thompson 
2940b31fde2SJeremy L Thompson /**
2950b31fde2SJeremy L Thompson   @brief Check input vector dimensions for CeedBasisApply[Add]AtPoints
2960b31fde2SJeremy L Thompson 
2970b31fde2SJeremy L Thompson   @param[in]  basis      `CeedBasis` to evaluate
2980b31fde2SJeremy L Thompson   @param[in]  num_elem   The number of elements to apply the basis evaluation to;
2990b31fde2SJeremy L Thompson                           the backend will specify the ordering in @ref CeedElemRestrictionCreate()
3000b31fde2SJeremy L Thompson   @param[in]  num_points Array of the number of points to apply the basis evaluation to in each element, size `num_elem`
3010b31fde2SJeremy L Thompson   @param[in]  t_mode     @ref CEED_NOTRANSPOSE to evaluate from nodes to points;
3020b31fde2SJeremy L Thompson                            @ref CEED_TRANSPOSE to apply the transpose, mapping from points to nodes
3030b31fde2SJeremy L Thompson   @param[in]  eval_mode  @ref CEED_EVAL_INTERP to use interpolated values,
3040b31fde2SJeremy L Thompson                            @ref CEED_EVAL_GRAD to use gradients,
3050b31fde2SJeremy L Thompson                            @ref CEED_EVAL_WEIGHT to use quadrature weights
3060b31fde2SJeremy L Thompson   @param[in]  x_ref      `CeedVector` holding reference coordinates of each point
3070b31fde2SJeremy L Thompson   @param[in]  u          Input `CeedVector`, of length `num_nodes * num_comp` for @ref CEED_NOTRANSPOSE
3080b31fde2SJeremy L Thompson   @param[out] v          Output `CeedVector`, of length `num_points * num_q_comp` for @ref CEED_NOTRANSPOSE with @ref CEED_EVAL_INTERP
3090b31fde2SJeremy L Thompson 
3100b31fde2SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
3110b31fde2SJeremy L Thompson 
3120b31fde2SJeremy L Thompson   @ref Developer
3130b31fde2SJeremy L Thompson **/
3140b31fde2SJeremy L Thompson static int CeedBasisApplyAtPointsCheckDims(CeedBasis basis, CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode,
3150b31fde2SJeremy L Thompson                                            CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) {
3160b31fde2SJeremy L Thompson   CeedInt  dim, num_comp, num_q_comp, num_nodes, P_1d = 1, Q_1d = 1, total_num_points = 0;
3170b31fde2SJeremy L Thompson   CeedSize x_length = 0, u_length = 0, v_length;
3180b31fde2SJeremy L Thompson 
3190b31fde2SJeremy L Thompson   CeedCall(CeedBasisGetDimension(basis, &dim));
3200b31fde2SJeremy L Thompson   CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d));
3210b31fde2SJeremy L Thompson   CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
3220b31fde2SJeremy L Thompson   CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
3230b31fde2SJeremy L Thompson   CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &num_q_comp));
3240b31fde2SJeremy L Thompson   CeedCall(CeedBasisGetNumNodes(basis, &num_nodes));
3250b31fde2SJeremy L Thompson   CeedCall(CeedVectorGetLength(v, &v_length));
3260b31fde2SJeremy L Thompson   if (x_ref != CEED_VECTOR_NONE) CeedCall(CeedVectorGetLength(x_ref, &x_length));
3270b31fde2SJeremy L Thompson   if (u != CEED_VECTOR_NONE) CeedCall(CeedVectorGetLength(u, &u_length));
3280b31fde2SJeremy L Thompson 
3290b31fde2SJeremy L Thompson   // Check compatibility coordinates vector
3300b31fde2SJeremy L Thompson   for (CeedInt i = 0; i < num_elem; i++) total_num_points += num_points[i];
3319bc66399SJeremy L Thompson   CeedCheck((x_length >= (CeedSize)total_num_points * (CeedSize)dim) || (eval_mode == CEED_EVAL_WEIGHT), CeedBasisReturnCeed(basis),
3329bc66399SJeremy L Thompson             CEED_ERROR_DIMENSION,
3330b31fde2SJeremy L Thompson             "Length of reference coordinate vector incompatible with basis dimension and number of points."
3340b31fde2SJeremy L Thompson             " Found reference coordinate vector of length %" CeedSize_FMT ", not of length %" CeedSize_FMT ".",
33519a04db8SJeremy L Thompson             x_length, (CeedSize)total_num_points * (CeedSize)dim);
3360b31fde2SJeremy L Thompson 
3370b31fde2SJeremy L Thompson   // Check CEED_EVAL_WEIGHT only on CEED_NOTRANSPOSE
3389bc66399SJeremy L Thompson   CeedCheck(eval_mode != CEED_EVAL_WEIGHT || t_mode == CEED_NOTRANSPOSE, CeedBasisReturnCeed(basis), CEED_ERROR_UNSUPPORTED,
3390b31fde2SJeremy L Thompson             "CEED_EVAL_WEIGHT only supported with CEED_NOTRANSPOSE");
3400b31fde2SJeremy L Thompson 
3410b31fde2SJeremy L Thompson   // Check vector lengths to prevent out of bounds issues
3420b31fde2SJeremy L Thompson   bool has_good_dims = true;
3430b31fde2SJeremy L Thompson   switch (eval_mode) {
3440b31fde2SJeremy L Thompson     case CEED_EVAL_INTERP:
34519a04db8SJeremy L Thompson       has_good_dims = ((t_mode == CEED_TRANSPOSE && (u_length >= (CeedSize)total_num_points * (CeedSize)num_q_comp ||
34619a04db8SJeremy L Thompson                                                      v_length >= (CeedSize)num_elem * (CeedSize)num_nodes * (CeedSize)num_comp)) ||
34719a04db8SJeremy L Thompson                        (t_mode == CEED_NOTRANSPOSE && (v_length >= (CeedSize)total_num_points * (CeedSize)num_q_comp ||
34819a04db8SJeremy L Thompson                                                        u_length >= (CeedSize)num_elem * (CeedSize)num_nodes * (CeedSize)num_comp)));
3490b31fde2SJeremy L Thompson       break;
3500b31fde2SJeremy L Thompson     case CEED_EVAL_GRAD:
35119a04db8SJeremy L Thompson       has_good_dims = ((t_mode == CEED_TRANSPOSE && (u_length >= (CeedSize)total_num_points * (CeedSize)num_q_comp * (CeedSize)dim ||
35219a04db8SJeremy L Thompson                                                      v_length >= (CeedSize)num_elem * (CeedSize)num_nodes * (CeedSize)num_comp)) ||
35319a04db8SJeremy L Thompson                        (t_mode == CEED_NOTRANSPOSE && (v_length >= (CeedSize)total_num_points * (CeedSize)num_q_comp * (CeedSize)dim ||
35419a04db8SJeremy L Thompson                                                        u_length >= (CeedSize)num_elem * (CeedSize)num_nodes * (CeedSize)num_comp)));
3550b31fde2SJeremy L Thompson       break;
3560b31fde2SJeremy L Thompson     case CEED_EVAL_WEIGHT:
3570b31fde2SJeremy L Thompson       has_good_dims = t_mode == CEED_NOTRANSPOSE && (v_length >= total_num_points);
3580b31fde2SJeremy L Thompson       break;
3590b31fde2SJeremy L Thompson       // LCOV_EXCL_START
3600b31fde2SJeremy L Thompson     case CEED_EVAL_NONE:
3610b31fde2SJeremy L Thompson     case CEED_EVAL_DIV:
3620b31fde2SJeremy L Thompson     case CEED_EVAL_CURL:
3639bc66399SJeremy L Thompson       return CeedError(CeedBasisReturnCeed(basis), CEED_ERROR_UNSUPPORTED, "Evaluation at arbitrary points not supported for %s",
3649bc66399SJeremy L Thompson                        CeedEvalModes[eval_mode]);
3650b31fde2SJeremy L Thompson       // LCOV_EXCL_STOP
3660b31fde2SJeremy L Thompson   }
3679bc66399SJeremy L Thompson   CeedCheck(has_good_dims, CeedBasisReturnCeed(basis), CEED_ERROR_DIMENSION, "Input/output vectors too short for basis and evaluation mode");
3680b31fde2SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3690b31fde2SJeremy L Thompson }
3700b31fde2SJeremy L Thompson 
3710b31fde2SJeremy L Thompson /**
3720b31fde2SJeremy L Thompson   @brief Default implimentation to apply basis evaluation from nodes to arbitrary points
3730b31fde2SJeremy L Thompson 
3740b31fde2SJeremy L Thompson   @param[in]  basis      `CeedBasis` to evaluate
3750b31fde2SJeremy L Thompson   @param[in]  apply_add  Sum result into target vector or overwrite
3760b31fde2SJeremy L Thompson   @param[in]  num_elem   The number of elements to apply the basis evaluation to;
3770b31fde2SJeremy L Thompson                           the backend will specify the ordering in @ref CeedElemRestrictionCreate()
3780b31fde2SJeremy L Thompson   @param[in]  num_points Array of the number of points to apply the basis evaluation to in each element, size `num_elem`
3790b31fde2SJeremy L Thompson   @param[in]  t_mode     @ref CEED_NOTRANSPOSE to evaluate from nodes to points;
3800b31fde2SJeremy L Thompson                            @ref CEED_TRANSPOSE to apply the transpose, mapping from points to nodes
3810b31fde2SJeremy L Thompson   @param[in]  eval_mode  @ref CEED_EVAL_INTERP to use interpolated values,
3820b31fde2SJeremy L Thompson                            @ref CEED_EVAL_GRAD to use gradients,
3830b31fde2SJeremy L Thompson                            @ref CEED_EVAL_WEIGHT to use quadrature weights
3840b31fde2SJeremy L Thompson   @param[in]  x_ref      `CeedVector` holding reference coordinates of each point
3850b31fde2SJeremy L Thompson   @param[in]  u          Input `CeedVector`, of length `num_nodes * num_comp` for @ref CEED_NOTRANSPOSE
3860b31fde2SJeremy L Thompson   @param[out] v          Output `CeedVector`, of length `num_points * num_q_comp` for @ref CEED_NOTRANSPOSE with @ref CEED_EVAL_INTERP
3870b31fde2SJeremy L Thompson 
3880b31fde2SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
3890b31fde2SJeremy L Thompson 
3900b31fde2SJeremy L Thompson   @ref Developer
3910b31fde2SJeremy L Thompson **/
3920b31fde2SJeremy L Thompson static int CeedBasisApplyAtPoints_Core(CeedBasis basis, bool apply_add, CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode,
3930b31fde2SJeremy L Thompson                                        CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) {
3940b31fde2SJeremy L Thompson   CeedInt dim, num_comp, P_1d = 1, Q_1d = 1, total_num_points = num_points[0];
3950b31fde2SJeremy L Thompson 
3960b31fde2SJeremy L Thompson   CeedCall(CeedBasisGetDimension(basis, &dim));
3970b31fde2SJeremy L Thompson   // Inserting check because clang-tidy doesn't understand this cannot occur
3989bc66399SJeremy L Thompson   CeedCheck(dim > 0, CeedBasisReturnCeed(basis), CEED_ERROR_UNSUPPORTED, "Malformed CeedBasis, dim > 0 is required");
3990b31fde2SJeremy L Thompson   CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d));
4000b31fde2SJeremy L Thompson   CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
4010b31fde2SJeremy L Thompson   CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
4020b31fde2SJeremy L Thompson 
4030b31fde2SJeremy L Thompson   // Default implementation
4040b31fde2SJeremy L Thompson   {
4050b31fde2SJeremy L Thompson     bool is_tensor_basis;
4060b31fde2SJeremy L Thompson 
4070b31fde2SJeremy L Thompson     CeedCall(CeedBasisIsTensor(basis, &is_tensor_basis));
4089bc66399SJeremy L Thompson     CeedCheck(is_tensor_basis, CeedBasisReturnCeed(basis), CEED_ERROR_UNSUPPORTED,
4099bc66399SJeremy L Thompson               "Evaluation at arbitrary points only supported for tensor product bases");
4100b31fde2SJeremy L Thompson   }
4119bc66399SJeremy L Thompson   CeedCheck(num_elem == 1, CeedBasisReturnCeed(basis), CEED_ERROR_UNSUPPORTED,
4129bc66399SJeremy L Thompson             "Evaluation at arbitrary  points only supported for a single element at a time");
4130b31fde2SJeremy L Thompson   if (eval_mode == CEED_EVAL_WEIGHT) {
4140b31fde2SJeremy L Thompson     CeedCall(CeedVectorSetValue(v, 1.0));
4150b31fde2SJeremy L Thompson     return CEED_ERROR_SUCCESS;
4160b31fde2SJeremy L Thompson   }
4170b31fde2SJeremy L Thompson   if (!basis->basis_chebyshev) {
4180b31fde2SJeremy L Thompson     // Build basis mapping from nodes to Chebyshev coefficients
4190b31fde2SJeremy L Thompson     CeedScalar       *chebyshev_interp_1d, *chebyshev_grad_1d, *chebyshev_q_weight_1d;
4200b31fde2SJeremy L Thompson     const CeedScalar *q_ref_1d;
4219bc66399SJeremy L Thompson     Ceed              ceed;
4220b31fde2SJeremy L Thompson 
4230b31fde2SJeremy L Thompson     CeedCall(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
4240b31fde2SJeremy L Thompson     CeedCall(CeedCalloc(P_1d * Q_1d, &chebyshev_grad_1d));
4250b31fde2SJeremy L Thompson     CeedCall(CeedCalloc(Q_1d, &chebyshev_q_weight_1d));
4260b31fde2SJeremy L Thompson     CeedCall(CeedBasisGetQRef(basis, &q_ref_1d));
4270b31fde2SJeremy L Thompson     CeedCall(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d));
4280b31fde2SJeremy L Thompson 
4299bc66399SJeremy L Thompson     CeedCall(CeedBasisGetCeed(basis, &ceed));
4300b31fde2SJeremy L Thompson     CeedCall(CeedVectorCreate(ceed, num_comp * CeedIntPow(Q_1d, dim), &basis->vec_chebyshev));
4310b31fde2SJeremy 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,
4320b31fde2SJeremy L Thompson                                      &basis->basis_chebyshev));
4330b31fde2SJeremy L Thompson 
4340b31fde2SJeremy L Thompson     // Cleanup
4350b31fde2SJeremy L Thompson     CeedCall(CeedFree(&chebyshev_interp_1d));
4360b31fde2SJeremy L Thompson     CeedCall(CeedFree(&chebyshev_grad_1d));
4370b31fde2SJeremy L Thompson     CeedCall(CeedFree(&chebyshev_q_weight_1d));
4389bc66399SJeremy L Thompson     CeedCall(CeedDestroy(&ceed));
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
4509bc66399SJeremy L Thompson     CeedCheck(basis_ref && basis_ref->contract, CeedBasisReturnCeed(basis), CEED_ERROR_UNSUPPORTED,
4519bc66399SJeremy L Thompson               "Reference CPU ceed failed to create a tensor contraction object");
4520b31fde2SJeremy L Thompson     CeedCall(CeedTensorContractReferenceCopy(basis_ref->contract, &basis->contract));
4530b31fde2SJeremy L Thompson     CeedCall(CeedBasisDestroy(&basis_ref));
4540b31fde2SJeremy L Thompson     CeedCall(CeedDestroy(&ceed_ref));
4550b31fde2SJeremy L Thompson   }
4560b31fde2SJeremy L Thompson 
4570b31fde2SJeremy L Thompson   // Basis evaluation
4580b31fde2SJeremy L Thompson   switch (t_mode) {
4590b31fde2SJeremy L Thompson     case CEED_NOTRANSPOSE: {
4600b31fde2SJeremy L Thompson       // Nodes to arbitrary points
4610b31fde2SJeremy L Thompson       CeedScalar       *v_array;
4620b31fde2SJeremy L Thompson       const CeedScalar *chebyshev_coeffs, *x_array_read;
4630b31fde2SJeremy L Thompson 
4640b31fde2SJeremy L Thompson       // -- Interpolate to Chebyshev coefficients
4650b31fde2SJeremy L Thompson       CeedCall(CeedBasisApply(basis->basis_chebyshev, 1, CEED_NOTRANSPOSE, CEED_EVAL_INTERP, u, basis->vec_chebyshev));
4660b31fde2SJeremy L Thompson 
4670b31fde2SJeremy L Thompson       // -- Evaluate Chebyshev polynomials at arbitrary points
4680b31fde2SJeremy L Thompson       CeedCall(CeedVectorGetArrayRead(basis->vec_chebyshev, CEED_MEM_HOST, &chebyshev_coeffs));
4690b31fde2SJeremy L Thompson       CeedCall(CeedVectorGetArrayRead(x_ref, CEED_MEM_HOST, &x_array_read));
4700b31fde2SJeremy L Thompson       CeedCall(CeedVectorGetArrayWrite(v, CEED_MEM_HOST, &v_array));
4710b31fde2SJeremy L Thompson       switch (eval_mode) {
4720b31fde2SJeremy L Thompson         case CEED_EVAL_INTERP: {
4730b31fde2SJeremy L Thompson           CeedScalar tmp[2][num_comp * CeedIntPow(Q_1d, dim)], chebyshev_x[Q_1d];
4740b31fde2SJeremy L Thompson 
4750b31fde2SJeremy L Thompson           // ---- Values at point
4760b31fde2SJeremy L Thompson           for (CeedInt p = 0; p < total_num_points; p++) {
4770b31fde2SJeremy L Thompson             CeedInt pre = num_comp * CeedIntPow(Q_1d, dim - 1), post = 1;
4780b31fde2SJeremy L Thompson 
4790b31fde2SJeremy L Thompson             for (CeedInt d = 0; d < dim; d++) {
4800b31fde2SJeremy L Thompson               // ------ Tensor contract with current Chebyshev polynomial values
4810b31fde2SJeremy L Thompson               CeedCall(CeedChebyshevPolynomialsAtPoint(x_array_read[d * total_num_points + p], Q_1d, chebyshev_x));
4820b31fde2SJeremy L Thompson               CeedCall(CeedTensorContractApply(basis->contract, pre, Q_1d, post, 1, chebyshev_x, t_mode, false,
4830b31fde2SJeremy L Thompson                                                d == 0 ? chebyshev_coeffs : tmp[d % 2], tmp[(d + 1) % 2]));
4840b31fde2SJeremy L Thompson               pre /= Q_1d;
4850b31fde2SJeremy L Thompson               post *= 1;
4860b31fde2SJeremy L Thompson             }
4870b31fde2SJeremy L Thompson             for (CeedInt c = 0; c < num_comp; c++) v_array[c * total_num_points + p] = tmp[dim % 2][c];
4880b31fde2SJeremy L Thompson           }
4890b31fde2SJeremy L Thompson           break;
4900b31fde2SJeremy L Thompson         }
4910b31fde2SJeremy L Thompson         case CEED_EVAL_GRAD: {
4920b31fde2SJeremy L Thompson           CeedScalar tmp[2][num_comp * CeedIntPow(Q_1d, dim)], chebyshev_x[Q_1d];
4930b31fde2SJeremy L Thompson 
4940b31fde2SJeremy L Thompson           // ---- Values at point
4950b31fde2SJeremy L Thompson           for (CeedInt p = 0; p < total_num_points; p++) {
4960b31fde2SJeremy L Thompson             // Dim**2 contractions, apply grad when pass == dim
4970b31fde2SJeremy L Thompson             for (CeedInt pass = 0; pass < dim; pass++) {
4980b31fde2SJeremy L Thompson               CeedInt pre = num_comp * CeedIntPow(Q_1d, dim - 1), post = 1;
4990b31fde2SJeremy L Thompson 
5000b31fde2SJeremy L Thompson               for (CeedInt d = 0; d < dim; d++) {
5010b31fde2SJeremy L Thompson                 // ------ Tensor contract with current Chebyshev polynomial values
5020b31fde2SJeremy L Thompson                 if (pass == d) CeedCall(CeedChebyshevDerivativeAtPoint(x_array_read[d * total_num_points + p], Q_1d, chebyshev_x));
5030b31fde2SJeremy L Thompson                 else CeedCall(CeedChebyshevPolynomialsAtPoint(x_array_read[d * total_num_points + p], Q_1d, chebyshev_x));
5040b31fde2SJeremy L Thompson                 CeedCall(CeedTensorContractApply(basis->contract, pre, Q_1d, post, 1, chebyshev_x, t_mode, false,
5050b31fde2SJeremy L Thompson                                                  d == 0 ? chebyshev_coeffs : tmp[d % 2], tmp[(d + 1) % 2]));
5060b31fde2SJeremy L Thompson                 pre /= Q_1d;
5070b31fde2SJeremy L Thompson                 post *= 1;
5080b31fde2SJeremy L Thompson               }
5090b31fde2SJeremy L Thompson               for (CeedInt c = 0; c < num_comp; c++) v_array[(pass * num_comp + c) * total_num_points + p] = tmp[dim % 2][c];
5100b31fde2SJeremy L Thompson             }
5110b31fde2SJeremy L Thompson           }
5120b31fde2SJeremy L Thompson           break;
5130b31fde2SJeremy L Thompson         }
5140b31fde2SJeremy L Thompson         default:
5150b31fde2SJeremy L Thompson           // Nothing to do, excluded above
5160b31fde2SJeremy L Thompson           break;
5170b31fde2SJeremy L Thompson       }
5180b31fde2SJeremy L Thompson       CeedCall(CeedVectorRestoreArrayRead(basis->vec_chebyshev, &chebyshev_coeffs));
5190b31fde2SJeremy L Thompson       CeedCall(CeedVectorRestoreArrayRead(x_ref, &x_array_read));
5200b31fde2SJeremy L Thompson       CeedCall(CeedVectorRestoreArray(v, &v_array));
5210b31fde2SJeremy L Thompson       break;
5220b31fde2SJeremy L Thompson     }
5230b31fde2SJeremy L Thompson     case CEED_TRANSPOSE: {
5240b31fde2SJeremy L Thompson       // Note: No switch on e_mode here because only CEED_EVAL_INTERP is supported at this time
5250b31fde2SJeremy L Thompson       // Arbitrary points to nodes
5260b31fde2SJeremy L Thompson       CeedScalar       *chebyshev_coeffs;
5270b31fde2SJeremy L Thompson       const CeedScalar *u_array, *x_array_read;
5280b31fde2SJeremy L Thompson 
5290b31fde2SJeremy L Thompson       // -- Transpose of evaluation of Chebyshev polynomials at arbitrary points
5300b31fde2SJeremy L Thompson       CeedCall(CeedVectorGetArrayWrite(basis->vec_chebyshev, CEED_MEM_HOST, &chebyshev_coeffs));
5310b31fde2SJeremy L Thompson       CeedCall(CeedVectorGetArrayRead(x_ref, CEED_MEM_HOST, &x_array_read));
5320b31fde2SJeremy L Thompson       CeedCall(CeedVectorGetArrayRead(u, CEED_MEM_HOST, &u_array));
5330b31fde2SJeremy L Thompson 
5340b31fde2SJeremy L Thompson       switch (eval_mode) {
5350b31fde2SJeremy L Thompson         case CEED_EVAL_INTERP: {
5360b31fde2SJeremy L Thompson           CeedScalar tmp[2][num_comp * CeedIntPow(Q_1d, dim)], chebyshev_x[Q_1d];
5370b31fde2SJeremy L Thompson 
5380b31fde2SJeremy L Thompson           // ---- Values at point
5390b31fde2SJeremy L Thompson           for (CeedInt p = 0; p < total_num_points; p++) {
5400b31fde2SJeremy L Thompson             CeedInt pre = num_comp * 1, post = 1;
5410b31fde2SJeremy L Thompson 
5420b31fde2SJeremy L Thompson             for (CeedInt c = 0; c < num_comp; c++) tmp[0][c] = u_array[c * total_num_points + p];
5430b31fde2SJeremy L Thompson             for (CeedInt d = 0; d < dim; d++) {
5440b31fde2SJeremy L Thompson               // ------ Tensor contract with current Chebyshev polynomial values
5450b31fde2SJeremy L Thompson               CeedCall(CeedChebyshevPolynomialsAtPoint(x_array_read[d * total_num_points + p], Q_1d, chebyshev_x));
5460b31fde2SJeremy L Thompson               CeedCall(CeedTensorContractApply(basis->contract, pre, 1, post, Q_1d, chebyshev_x, t_mode, p > 0 && d == (dim - 1), tmp[d % 2],
5470b31fde2SJeremy L Thompson                                                d == (dim - 1) ? chebyshev_coeffs : tmp[(d + 1) % 2]));
5480b31fde2SJeremy L Thompson               pre /= 1;
5490b31fde2SJeremy L Thompson               post *= Q_1d;
5500b31fde2SJeremy L Thompson             }
5510b31fde2SJeremy L Thompson           }
5520b31fde2SJeremy L Thompson           break;
5530b31fde2SJeremy L Thompson         }
5540b31fde2SJeremy L Thompson         case CEED_EVAL_GRAD: {
5550b31fde2SJeremy L Thompson           CeedScalar tmp[2][num_comp * CeedIntPow(Q_1d, dim)], chebyshev_x[Q_1d];
5560b31fde2SJeremy L Thompson 
5570b31fde2SJeremy L Thompson           // ---- Values at point
5580b31fde2SJeremy L Thompson           for (CeedInt p = 0; p < total_num_points; p++) {
5590b31fde2SJeremy L Thompson             // Dim**2 contractions, apply grad when pass == dim
5600b31fde2SJeremy L Thompson             for (CeedInt pass = 0; pass < dim; pass++) {
5610b31fde2SJeremy L Thompson               CeedInt pre = num_comp * 1, post = 1;
5620b31fde2SJeremy L Thompson 
5630b31fde2SJeremy L Thompson               for (CeedInt c = 0; c < num_comp; c++) tmp[0][c] = u_array[(pass * num_comp + c) * total_num_points + p];
5640b31fde2SJeremy L Thompson               for (CeedInt d = 0; d < dim; d++) {
5650b31fde2SJeremy L Thompson                 // ------ Tensor contract with current Chebyshev polynomial values
5660b31fde2SJeremy L Thompson                 if (pass == d) CeedCall(CeedChebyshevDerivativeAtPoint(x_array_read[d * total_num_points + p], Q_1d, chebyshev_x));
5670b31fde2SJeremy L Thompson                 else CeedCall(CeedChebyshevPolynomialsAtPoint(x_array_read[d * total_num_points + p], Q_1d, chebyshev_x));
5680b31fde2SJeremy L Thompson                 CeedCall(CeedTensorContractApply(basis->contract, pre, 1, post, Q_1d, chebyshev_x, t_mode,
5690b31fde2SJeremy L Thompson                                                  (p > 0 || (p == 0 && pass > 0)) && d == (dim - 1), tmp[d % 2],
5700b31fde2SJeremy L Thompson                                                  d == (dim - 1) ? chebyshev_coeffs : tmp[(d + 1) % 2]));
5710b31fde2SJeremy L Thompson                 pre /= 1;
5720b31fde2SJeremy L Thompson                 post *= Q_1d;
5730b31fde2SJeremy L Thompson               }
5740b31fde2SJeremy L Thompson             }
5750b31fde2SJeremy L Thompson           }
5760b31fde2SJeremy L Thompson           break;
5770b31fde2SJeremy L Thompson         }
5780b31fde2SJeremy L Thompson         default:
5790b31fde2SJeremy L Thompson           // Nothing to do, excluded above
5800b31fde2SJeremy L Thompson           break;
5810b31fde2SJeremy L Thompson       }
5820b31fde2SJeremy L Thompson       CeedCall(CeedVectorRestoreArray(basis->vec_chebyshev, &chebyshev_coeffs));
5830b31fde2SJeremy L Thompson       CeedCall(CeedVectorRestoreArrayRead(x_ref, &x_array_read));
5840b31fde2SJeremy L Thompson       CeedCall(CeedVectorRestoreArrayRead(u, &u_array));
5850b31fde2SJeremy L Thompson 
5860b31fde2SJeremy L Thompson       // -- Interpolate transpose from Chebyshev coefficients
5870b31fde2SJeremy L Thompson       if (apply_add) CeedCall(CeedBasisApplyAdd(basis->basis_chebyshev, 1, CEED_TRANSPOSE, CEED_EVAL_INTERP, basis->vec_chebyshev, v));
5880b31fde2SJeremy L Thompson       else CeedCall(CeedBasisApply(basis->basis_chebyshev, 1, CEED_TRANSPOSE, CEED_EVAL_INTERP, basis->vec_chebyshev, v));
5890b31fde2SJeremy L Thompson       break;
5900b31fde2SJeremy L Thompson     }
5910b31fde2SJeremy L Thompson   }
5920b31fde2SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5930b31fde2SJeremy L Thompson }
5940b31fde2SJeremy L Thompson 
5957a982d89SJeremy L. Thompson /// @}
5967a982d89SJeremy L. Thompson 
5977a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
5987a982d89SJeremy L. Thompson /// Ceed Backend API
5997a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
6007a982d89SJeremy L. Thompson /// @addtogroup CeedBasisBackend
6017a982d89SJeremy L. Thompson /// @{
6027a982d89SJeremy L. Thompson 
6037a982d89SJeremy L. Thompson /**
604ca94c3ddSJeremy L Thompson   @brief Return collocated gradient matrix
6057a982d89SJeremy L. Thompson 
606ca94c3ddSJeremy L Thompson   @param[in]  basis         `CeedBasis`
607ca94c3ddSJeremy L Thompson   @param[out] collo_grad_1d Row-major (`Q_1d * Q_1d`) matrix expressing derivatives of basis functions at quadrature points
6087a982d89SJeremy L. Thompson 
6097a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
6107a982d89SJeremy L. Thompson 
6117a982d89SJeremy L. Thompson   @ref Backend
6127a982d89SJeremy L. Thompson **/
613d1d35e2fSjeremylt int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *collo_grad_1d) {
6147a982d89SJeremy L. Thompson   Ceed              ceed;
6152247a93fSRezgar Shakeri   CeedInt           P_1d, Q_1d;
6162247a93fSRezgar Shakeri   CeedScalar       *interp_1d_pinv;
6171203703bSJeremy L Thompson   const CeedScalar *grad_1d, *interp_1d;
6181203703bSJeremy L Thompson 
619ea61e9acSJeremy 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.
6202247a93fSRezgar Shakeri   CeedCall(CeedBasisGetCeed(basis, &ceed));
6212247a93fSRezgar Shakeri   CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d));
6222247a93fSRezgar Shakeri   CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
6237a982d89SJeremy L. Thompson 
6242247a93fSRezgar Shakeri   // Compute interp_1d^+, pseudoinverse of interp_1d
6252247a93fSRezgar Shakeri   CeedCall(CeedCalloc(P_1d * Q_1d, &interp_1d_pinv));
6261203703bSJeremy L Thompson   CeedCall(CeedBasisGetInterp1D(basis, &interp_1d));
6271203703bSJeremy L Thompson   CeedCall(CeedMatrixPseudoinverse(ceed, interp_1d, Q_1d, P_1d, interp_1d_pinv));
6281203703bSJeremy L Thompson   CeedCall(CeedBasisGetGrad1D(basis, &grad_1d));
6291203703bSJeremy L Thompson   CeedCall(CeedMatrixMatrixMultiply(ceed, grad_1d, (const CeedScalar *)interp_1d_pinv, collo_grad_1d, Q_1d, Q_1d, P_1d));
6307a982d89SJeremy L. Thompson 
6312247a93fSRezgar Shakeri   CeedCall(CeedFree(&interp_1d_pinv));
6329bc66399SJeremy L Thompson   CeedCall(CeedDestroy(&ceed));
633e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
6347a982d89SJeremy L. Thompson }
6357a982d89SJeremy L. Thompson 
6367a982d89SJeremy L. Thompson /**
637b0cc4569SJeremy L Thompson   @brief Return 1D interpolation matrix to Chebyshev polynomial coefficients on quadrature space
638b0cc4569SJeremy L Thompson 
639b0cc4569SJeremy L Thompson   @param[in]  basis               `CeedBasis`
640b0cc4569SJeremy L Thompson   @param[out] chebyshev_interp_1d Row-major (`P_1d * Q_1d`) matrix interpolating from basis nodes to Chebyshev polynomial coefficients
641b0cc4569SJeremy L Thompson 
642b0cc4569SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
643b0cc4569SJeremy L Thompson 
644b0cc4569SJeremy L Thompson   @ref Backend
645b0cc4569SJeremy L Thompson **/
646b0cc4569SJeremy L Thompson int CeedBasisGetChebyshevInterp1D(CeedBasis basis, CeedScalar *chebyshev_interp_1d) {
647b0cc4569SJeremy L Thompson   CeedInt           P_1d, Q_1d;
648b0cc4569SJeremy L Thompson   CeedScalar       *C, *chebyshev_coeffs_1d_inv;
649b0cc4569SJeremy L Thompson   const CeedScalar *interp_1d, *q_ref_1d;
650b0cc4569SJeremy L Thompson   Ceed              ceed;
651b0cc4569SJeremy L Thompson 
652b0cc4569SJeremy L Thompson   CeedCall(CeedBasisGetCeed(basis, &ceed));
653b0cc4569SJeremy L Thompson   CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d));
654b0cc4569SJeremy L Thompson   CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
655b0cc4569SJeremy L Thompson 
656b0cc4569SJeremy L Thompson   // Build coefficient matrix
657bd83cbc5SJeremy L Thompson   // -- Note: Clang-tidy needs this check
658bd83cbc5SJeremy L Thompson   CeedCheck((P_1d > 0) && (Q_1d > 0), ceed, CEED_ERROR_INCOMPATIBLE, "CeedBasis dimensions are malformed");
659b0cc4569SJeremy L Thompson   CeedCall(CeedCalloc(Q_1d * Q_1d, &C));
660b0cc4569SJeremy L Thompson   CeedCall(CeedBasisGetQRef(basis, &q_ref_1d));
661b0cc4569SJeremy L Thompson   for (CeedInt i = 0; i < Q_1d; i++) CeedCall(CeedChebyshevPolynomialsAtPoint(q_ref_1d[i], Q_1d, &C[i * Q_1d]));
662b0cc4569SJeremy L Thompson 
663b0cc4569SJeremy L Thompson   // Compute C^+, pseudoinverse of coefficient matrix
664b0cc4569SJeremy L Thompson   CeedCall(CeedCalloc(Q_1d * Q_1d, &chebyshev_coeffs_1d_inv));
665b0cc4569SJeremy L Thompson   CeedCall(CeedMatrixPseudoinverse(ceed, C, Q_1d, Q_1d, chebyshev_coeffs_1d_inv));
666b0cc4569SJeremy L Thompson 
667b0cc4569SJeremy L Thompson   // Build mapping from nodes to Chebyshev coefficients
668b0cc4569SJeremy L Thompson   CeedCall(CeedBasisGetInterp1D(basis, &interp_1d));
669b0cc4569SJeremy L Thompson   CeedCall(CeedMatrixMatrixMultiply(ceed, chebyshev_coeffs_1d_inv, interp_1d, chebyshev_interp_1d, Q_1d, P_1d, Q_1d));
670b0cc4569SJeremy L Thompson 
671b0cc4569SJeremy L Thompson   // Cleanup
672b0cc4569SJeremy L Thompson   CeedCall(CeedFree(&C));
673b0cc4569SJeremy L Thompson   CeedCall(CeedFree(&chebyshev_coeffs_1d_inv));
6749bc66399SJeremy L Thompson   CeedCall(CeedDestroy(&ceed));
675b0cc4569SJeremy L Thompson   return CEED_ERROR_SUCCESS;
676b0cc4569SJeremy L Thompson }
677b0cc4569SJeremy L Thompson 
678b0cc4569SJeremy L Thompson /**
679ca94c3ddSJeremy L Thompson   @brief Get tensor status for given `CeedBasis`
6807a982d89SJeremy L. Thompson 
681ca94c3ddSJeremy L Thompson   @param[in]  basis     `CeedBasis`
682d1d35e2fSjeremylt   @param[out] is_tensor Variable to store tensor status
6837a982d89SJeremy L. Thompson 
6847a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
6857a982d89SJeremy L. Thompson 
6867a982d89SJeremy L. Thompson   @ref Backend
6877a982d89SJeremy L. Thompson **/
688d1d35e2fSjeremylt int CeedBasisIsTensor(CeedBasis basis, bool *is_tensor) {
6896402da51SJeremy L Thompson   *is_tensor = basis->is_tensor_basis;
690e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
6917a982d89SJeremy L. Thompson }
6927a982d89SJeremy L. Thompson 
6937a982d89SJeremy L. Thompson /**
694ca94c3ddSJeremy L Thompson   @brief Get backend data of a `CeedBasis`
6957a982d89SJeremy L. Thompson 
696ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
6977a982d89SJeremy L. Thompson   @param[out] data  Variable to store data
6987a982d89SJeremy L. Thompson 
6997a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
7007a982d89SJeremy L. Thompson 
7017a982d89SJeremy L. Thompson   @ref Backend
7027a982d89SJeremy L. Thompson **/
703777ff853SJeremy L Thompson int CeedBasisGetData(CeedBasis basis, void *data) {
704777ff853SJeremy L Thompson   *(void **)data = basis->data;
705e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
7067a982d89SJeremy L. Thompson }
7077a982d89SJeremy L. Thompson 
7087a982d89SJeremy L. Thompson /**
709ca94c3ddSJeremy L Thompson   @brief Set backend data of a `CeedBasis`
7107a982d89SJeremy L. Thompson 
711ca94c3ddSJeremy L Thompson   @param[in,out] basis  `CeedBasis`
712ea61e9acSJeremy L Thompson   @param[in]     data   Data to set
7137a982d89SJeremy L. Thompson 
7147a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
7157a982d89SJeremy L. Thompson 
7167a982d89SJeremy L. Thompson   @ref Backend
7177a982d89SJeremy L. Thompson **/
718777ff853SJeremy L Thompson int CeedBasisSetData(CeedBasis basis, void *data) {
719777ff853SJeremy L Thompson   basis->data = data;
720e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
7217a982d89SJeremy L. Thompson }
7227a982d89SJeremy L. Thompson 
7237a982d89SJeremy L. Thompson /**
724ca94c3ddSJeremy L Thompson   @brief Increment the reference counter for a `CeedBasis`
72534359f16Sjeremylt 
726ca94c3ddSJeremy L Thompson   @param[in,out] basis `CeedBasis` to increment the reference counter
72734359f16Sjeremylt 
72834359f16Sjeremylt   @return An error code: 0 - success, otherwise - failure
72934359f16Sjeremylt 
73034359f16Sjeremylt   @ref Backend
73134359f16Sjeremylt **/
7329560d06aSjeremylt int CeedBasisReference(CeedBasis basis) {
73334359f16Sjeremylt   basis->ref_count++;
73434359f16Sjeremylt   return CEED_ERROR_SUCCESS;
73534359f16Sjeremylt }
73634359f16Sjeremylt 
73734359f16Sjeremylt /**
738ca94c3ddSJeremy L Thompson   @brief Get number of Q-vector components for given `CeedBasis`
739c4e3f59bSSebastian Grimberg 
740ca94c3ddSJeremy L Thompson   @param[in]  basis     `CeedBasis`
741ca94c3ddSJeremy L Thompson   @param[in]  eval_mode @ref CEED_EVAL_INTERP to use interpolated values,
742ca94c3ddSJeremy L Thompson                           @ref CEED_EVAL_GRAD to use gradients,
743ca94c3ddSJeremy L Thompson                           @ref CEED_EVAL_DIV to use divergence,
744ca94c3ddSJeremy L Thompson                           @ref CEED_EVAL_CURL to use curl
745c4e3f59bSSebastian Grimberg   @param[out] q_comp    Variable to store number of Q-vector components of basis
746c4e3f59bSSebastian Grimberg 
747c4e3f59bSSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
748c4e3f59bSSebastian Grimberg 
749c4e3f59bSSebastian Grimberg   @ref Backend
750c4e3f59bSSebastian Grimberg **/
751c4e3f59bSSebastian Grimberg int CeedBasisGetNumQuadratureComponents(CeedBasis basis, CeedEvalMode eval_mode, CeedInt *q_comp) {
7521203703bSJeremy L Thompson   CeedInt dim;
7531203703bSJeremy L Thompson 
7541203703bSJeremy L Thompson   CeedCall(CeedBasisGetDimension(basis, &dim));
755c4e3f59bSSebastian Grimberg   switch (eval_mode) {
7561203703bSJeremy L Thompson     case CEED_EVAL_INTERP: {
7571203703bSJeremy L Thompson       CeedFESpace fe_space;
7581203703bSJeremy L Thompson 
7591203703bSJeremy L Thompson       CeedCall(CeedBasisGetFESpace(basis, &fe_space));
7601203703bSJeremy L Thompson       *q_comp = (fe_space == CEED_FE_SPACE_H1) ? 1 : dim;
7611203703bSJeremy L Thompson     } break;
762c4e3f59bSSebastian Grimberg     case CEED_EVAL_GRAD:
7631203703bSJeremy L Thompson       *q_comp = dim;
764c4e3f59bSSebastian Grimberg       break;
765c4e3f59bSSebastian Grimberg     case CEED_EVAL_DIV:
766c4e3f59bSSebastian Grimberg       *q_comp = 1;
767c4e3f59bSSebastian Grimberg       break;
768c4e3f59bSSebastian Grimberg     case CEED_EVAL_CURL:
7691203703bSJeremy L Thompson       *q_comp = (dim < 3) ? 1 : dim;
770c4e3f59bSSebastian Grimberg       break;
771c4e3f59bSSebastian Grimberg     case CEED_EVAL_NONE:
772c4e3f59bSSebastian Grimberg     case CEED_EVAL_WEIGHT:
773352a5e7cSSebastian Grimberg       *q_comp = 1;
774c4e3f59bSSebastian Grimberg       break;
775c4e3f59bSSebastian Grimberg   }
776c4e3f59bSSebastian Grimberg   return CEED_ERROR_SUCCESS;
777c4e3f59bSSebastian Grimberg }
778c4e3f59bSSebastian Grimberg 
779c4e3f59bSSebastian Grimberg /**
780ca94c3ddSJeremy L Thompson   @brief Estimate number of FLOPs required to apply `CeedBasis` in `t_mode` and `eval_mode`
7816e15d496SJeremy L Thompson 
782ca94c3ddSJeremy L Thompson   @param[in]  basis        `CeedBasis` to estimate FLOPs for
783ea61e9acSJeremy L Thompson   @param[in]  t_mode       Apply basis or transpose
784ca94c3ddSJeremy L Thompson   @param[in]  eval_mode    @ref CeedEvalMode
785*3f919cbcSJeremy L Thompson   @param[in]  is_at_points Evaluate the basis at points or quadrature points
786*3f919cbcSJeremy L Thompson   @param[in]  num_points   Number of points basis is evaluated at
787ea61e9acSJeremy L Thompson   @param[out] flops        Address of variable to hold FLOPs estimate
7886e15d496SJeremy L Thompson 
7896e15d496SJeremy L Thompson   @ref Backend
7906e15d496SJeremy L Thompson **/
791*3f919cbcSJeremy L Thompson int CeedBasisGetFlopsEstimate(CeedBasis basis, CeedTransposeMode t_mode, CeedEvalMode eval_mode, bool is_at_points, CeedInt num_points,
792*3f919cbcSJeremy L Thompson                               CeedSize *flops) {
7936e15d496SJeremy L Thompson   bool is_tensor;
7946e15d496SJeremy L Thompson 
7952b730f8bSJeremy L Thompson   CeedCall(CeedBasisIsTensor(basis, &is_tensor));
796*3f919cbcSJeremy L Thompson   CeedCheck(!is_at_points || is_tensor, CeedBasisReturnCeed(basis), CEED_ERROR_INCOMPATIBLE, "Can only evaluate tensor-product bases at points");
7976e15d496SJeremy L Thompson   if (is_tensor) {
7986e15d496SJeremy L Thompson     CeedInt dim, num_comp, P_1d, Q_1d;
7991c66c397SJeremy L Thompson 
8002b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetDimension(basis, &dim));
8012b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
8022b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d));
8032b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
8046e15d496SJeremy L Thompson     if (t_mode == CEED_TRANSPOSE) {
8052b730f8bSJeremy L Thompson       P_1d = Q_1d;
8062b730f8bSJeremy L Thompson       Q_1d = P_1d;
8076e15d496SJeremy L Thompson     }
8086e15d496SJeremy L Thompson     CeedInt tensor_flops = 0, pre = num_comp * CeedIntPow(P_1d, dim - 1), post = 1;
809*3f919cbcSJeremy L Thompson 
8106e15d496SJeremy L Thompson     for (CeedInt d = 0; d < dim; d++) {
8116e15d496SJeremy L Thompson       tensor_flops += 2 * pre * P_1d * post * Q_1d;
8126e15d496SJeremy L Thompson       pre /= P_1d;
8136e15d496SJeremy L Thompson       post *= Q_1d;
8146e15d496SJeremy L Thompson     }
815*3f919cbcSJeremy L Thompson     if (is_at_points) {
816*3f919cbcSJeremy L Thompson       CeedInt chebyshev_flops = (Q_1d - 2) * 3 + 1, d_chebyshev_flops = (Q_1d - 2) * 8 + 1;
817*3f919cbcSJeremy L Thompson       CeedInt point_tensor_flops = 0, pre = CeedIntPow(Q_1d, dim - 1), post = 1;
818*3f919cbcSJeremy L Thompson 
819*3f919cbcSJeremy L Thompson       for (CeedInt d = 0; d < dim; d++) {
820*3f919cbcSJeremy L Thompson         point_tensor_flops += 2 * pre * Q_1d * post * 1;
821*3f919cbcSJeremy L Thompson         pre /= P_1d;
822*3f919cbcSJeremy L Thompson         post *= Q_1d;
823*3f919cbcSJeremy L Thompson       }
824*3f919cbcSJeremy L Thompson 
825*3f919cbcSJeremy L Thompson       switch (eval_mode) {
826*3f919cbcSJeremy L Thompson         case CEED_EVAL_NONE:
827*3f919cbcSJeremy L Thompson           *flops = 0;
828*3f919cbcSJeremy L Thompson           break;
829*3f919cbcSJeremy L Thompson         case CEED_EVAL_INTERP:
830*3f919cbcSJeremy L Thompson           *flops = tensor_flops + num_points * (dim * chebyshev_flops + point_tensor_flops + (t_mode == CEED_TRANSPOSE ? CeedIntPow(Q_1d, dim) : 0));
831*3f919cbcSJeremy L Thompson           break;
832*3f919cbcSJeremy L Thompson         case CEED_EVAL_GRAD:
833*3f919cbcSJeremy L Thompson           *flops = tensor_flops + num_points * (dim * (d_chebyshev_flops + (dim - 1) * chebyshev_flops + point_tensor_flops +
834*3f919cbcSJeremy L Thompson                                                        (t_mode == CEED_TRANSPOSE ? CeedIntPow(Q_1d, dim) : 0)));
835*3f919cbcSJeremy L Thompson           break;
836*3f919cbcSJeremy L Thompson         case CEED_EVAL_DIV:
837*3f919cbcSJeremy L Thompson         case CEED_EVAL_CURL: {
838*3f919cbcSJeremy L Thompson           // LCOV_EXCL_START
839*3f919cbcSJeremy L Thompson           return CeedError(CeedBasisReturnCeed(basis), CEED_ERROR_INCOMPATIBLE, "Tensor basis evaluation for %s not supported",
840*3f919cbcSJeremy L Thompson                            CeedEvalModes[eval_mode]);
841*3f919cbcSJeremy L Thompson           break;
842*3f919cbcSJeremy L Thompson           // LCOV_EXCL_STOP
843*3f919cbcSJeremy L Thompson         }
844*3f919cbcSJeremy L Thompson         case CEED_EVAL_WEIGHT:
845*3f919cbcSJeremy L Thompson           *flops = num_points;
846*3f919cbcSJeremy L Thompson           break;
847*3f919cbcSJeremy L Thompson       }
848*3f919cbcSJeremy L Thompson     } else {
8496e15d496SJeremy L Thompson       switch (eval_mode) {
8502b730f8bSJeremy L Thompson         case CEED_EVAL_NONE:
8512b730f8bSJeremy L Thompson           *flops = 0;
8522b730f8bSJeremy L Thompson           break;
8532b730f8bSJeremy L Thompson         case CEED_EVAL_INTERP:
8542b730f8bSJeremy L Thompson           *flops = tensor_flops;
8552b730f8bSJeremy L Thompson           break;
8562b730f8bSJeremy L Thompson         case CEED_EVAL_GRAD:
8572b730f8bSJeremy L Thompson           *flops = tensor_flops * 2;
8582b730f8bSJeremy L Thompson           break;
8596e15d496SJeremy L Thompson         case CEED_EVAL_DIV:
8601203703bSJeremy L Thompson         case CEED_EVAL_CURL: {
8616574a04fSJeremy L Thompson           // LCOV_EXCL_START
8626e536b99SJeremy L Thompson           return CeedError(CeedBasisReturnCeed(basis), CEED_ERROR_INCOMPATIBLE, "Tensor basis evaluation for %s not supported",
8636e536b99SJeremy L Thompson                            CeedEvalModes[eval_mode]);
8642b730f8bSJeremy L Thompson           break;
8656e15d496SJeremy L Thompson           // LCOV_EXCL_STOP
8661203703bSJeremy L Thompson         }
8672b730f8bSJeremy L Thompson         case CEED_EVAL_WEIGHT:
8682b730f8bSJeremy L Thompson           *flops = dim * CeedIntPow(Q_1d, dim);
8692b730f8bSJeremy L Thompson           break;
8706e15d496SJeremy L Thompson       }
871*3f919cbcSJeremy L Thompson     }
8726e15d496SJeremy L Thompson   } else {
873c4e3f59bSSebastian Grimberg     CeedInt dim, num_comp, q_comp, num_nodes, num_qpts;
8741c66c397SJeremy L Thompson 
8752b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetDimension(basis, &dim));
8762b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
877c4e3f59bSSebastian Grimberg     CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp));
8782b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis, &num_nodes));
8792b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts));
8806e15d496SJeremy L Thompson     switch (eval_mode) {
8812b730f8bSJeremy L Thompson       case CEED_EVAL_NONE:
8822b730f8bSJeremy L Thompson         *flops = 0;
8832b730f8bSJeremy L Thompson         break;
8842b730f8bSJeremy L Thompson       case CEED_EVAL_INTERP:
8852b730f8bSJeremy L Thompson       case CEED_EVAL_GRAD:
8862b730f8bSJeremy L Thompson       case CEED_EVAL_DIV:
8872b730f8bSJeremy L Thompson       case CEED_EVAL_CURL:
888c4e3f59bSSebastian Grimberg         *flops = num_nodes * num_qpts * num_comp * q_comp;
8892b730f8bSJeremy L Thompson         break;
8902b730f8bSJeremy L Thompson       case CEED_EVAL_WEIGHT:
8912b730f8bSJeremy L Thompson         *flops = 0;
8922b730f8bSJeremy L Thompson         break;
8936e15d496SJeremy L Thompson     }
8946e15d496SJeremy L Thompson   }
8956e15d496SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8966e15d496SJeremy L Thompson }
8976e15d496SJeremy L Thompson 
8986e15d496SJeremy L Thompson /**
899ca94c3ddSJeremy L Thompson   @brief Get `CeedFESpace` for a `CeedBasis`
900c4e3f59bSSebastian Grimberg 
901ca94c3ddSJeremy L Thompson   @param[in]  basis    `CeedBasis`
902ca94c3ddSJeremy L Thompson   @param[out] fe_space Variable to store `CeedFESpace`
903c4e3f59bSSebastian Grimberg 
904c4e3f59bSSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
905c4e3f59bSSebastian Grimberg 
906c4e3f59bSSebastian Grimberg   @ref Backend
907c4e3f59bSSebastian Grimberg **/
908c4e3f59bSSebastian Grimberg int CeedBasisGetFESpace(CeedBasis basis, CeedFESpace *fe_space) {
909c4e3f59bSSebastian Grimberg   *fe_space = basis->fe_space;
910c4e3f59bSSebastian Grimberg   return CEED_ERROR_SUCCESS;
911c4e3f59bSSebastian Grimberg }
912c4e3f59bSSebastian Grimberg 
913c4e3f59bSSebastian Grimberg /**
914ca94c3ddSJeremy L Thompson   @brief Get dimension for given `CeedElemTopology`
9157a982d89SJeremy L. Thompson 
916ca94c3ddSJeremy L Thompson   @param[in]  topo `CeedElemTopology`
9177a982d89SJeremy L. Thompson   @param[out] dim  Variable to store dimension of topology
9187a982d89SJeremy L. Thompson 
9197a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
9207a982d89SJeremy L. Thompson 
9217a982d89SJeremy L. Thompson   @ref Backend
9227a982d89SJeremy L. Thompson **/
9237a982d89SJeremy L. Thompson int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) {
9247a982d89SJeremy L. Thompson   *dim = (CeedInt)topo >> 16;
925e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9267a982d89SJeremy L. Thompson }
9277a982d89SJeremy L. Thompson 
9287a982d89SJeremy L. Thompson /**
929ca94c3ddSJeremy L Thompson   @brief Get `CeedTensorContract` of a `CeedBasis`
9307a982d89SJeremy L. Thompson 
931ca94c3ddSJeremy L Thompson   @param[in]  basis     `CeedBasis`
932ca94c3ddSJeremy L Thompson   @param[out] contract  Variable to store `CeedTensorContract`
9337a982d89SJeremy L. Thompson 
9347a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
9357a982d89SJeremy L. Thompson 
9367a982d89SJeremy L. Thompson   @ref Backend
9377a982d89SJeremy L. Thompson **/
9387a982d89SJeremy L. Thompson int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) {
9397a982d89SJeremy L. Thompson   *contract = basis->contract;
940e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9417a982d89SJeremy L. Thompson }
9427a982d89SJeremy L. Thompson 
9437a982d89SJeremy L. Thompson /**
944ca94c3ddSJeremy L Thompson   @brief Set `CeedTensorContract` of a `CeedBasis`
9457a982d89SJeremy L. Thompson 
946ca94c3ddSJeremy L Thompson   @param[in,out] basis    `CeedBasis`
947ca94c3ddSJeremy L Thompson   @param[in]     contract `CeedTensorContract` to set
9487a982d89SJeremy L. Thompson 
9497a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
9507a982d89SJeremy L. Thompson 
9517a982d89SJeremy L. Thompson   @ref Backend
9527a982d89SJeremy L. Thompson **/
95334359f16Sjeremylt int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract contract) {
95434359f16Sjeremylt   basis->contract = contract;
9552b730f8bSJeremy L Thompson   CeedCall(CeedTensorContractReference(contract));
956e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9577a982d89SJeremy L. Thompson }
9587a982d89SJeremy L. Thompson 
9597a982d89SJeremy L. Thompson /**
960ca94c3ddSJeremy L Thompson   @brief Return a reference implementation of matrix multiplication \f$C = A B\f$.
961ba59ac12SSebastian Grimberg 
962ca94c3ddSJeremy L Thompson   Note: This is a reference implementation for CPU `CeedScalar` pointers that is not intended for high performance.
9637a982d89SJeremy L. Thompson 
964ca94c3ddSJeremy L Thompson   @param[in]  ceed  `Ceed` context for error handling
965ca94c3ddSJeremy L Thompson   @param[in]  mat_A Row-major matrix `A`
966ca94c3ddSJeremy L Thompson   @param[in]  mat_B Row-major matrix `B`
967ca94c3ddSJeremy L Thompson   @param[out] mat_C Row-major output matrix `C`
968ca94c3ddSJeremy L Thompson   @param[in]  m     Number of rows of `C`
969ca94c3ddSJeremy L Thompson   @param[in]  n     Number of columns of `C`
970ca94c3ddSJeremy L Thompson   @param[in]  kk    Number of columns of `A`/rows of `B`
9717a982d89SJeremy L. Thompson 
9727a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
9737a982d89SJeremy L. Thompson 
9747a982d89SJeremy L. Thompson   @ref Utility
9757a982d89SJeremy L. Thompson **/
9762b730f8bSJeremy L Thompson int CeedMatrixMatrixMultiply(Ceed ceed, const CeedScalar *mat_A, const CeedScalar *mat_B, CeedScalar *mat_C, CeedInt m, CeedInt n, CeedInt kk) {
9772b730f8bSJeremy L Thompson   for (CeedInt i = 0; i < m; i++) {
9787a982d89SJeremy L. Thompson     for (CeedInt j = 0; j < n; j++) {
9797a982d89SJeremy L. Thompson       CeedScalar sum = 0;
9801c66c397SJeremy L Thompson 
9812b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < kk; k++) sum += mat_A[k + i * kk] * mat_B[j + k * n];
982d1d35e2fSjeremylt       mat_C[j + i * n] = sum;
9837a982d89SJeremy L. Thompson     }
9842b730f8bSJeremy L Thompson   }
985e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9867a982d89SJeremy L. Thompson }
9877a982d89SJeremy L. Thompson 
988ba59ac12SSebastian Grimberg /**
989ba59ac12SSebastian Grimberg   @brief Return QR Factorization of a matrix
990ba59ac12SSebastian Grimberg 
991ca94c3ddSJeremy L Thompson   @param[in]     ceed `Ceed` context for error handling
992ba59ac12SSebastian Grimberg   @param[in,out] mat  Row-major matrix to be factorized in place
993ca94c3ddSJeremy L Thompson   @param[in,out] tau  Vector of length `m` of scaling factors
994ba59ac12SSebastian Grimberg   @param[in]     m    Number of rows
995ba59ac12SSebastian Grimberg   @param[in]     n    Number of columns
996ba59ac12SSebastian Grimberg 
997ba59ac12SSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
998ba59ac12SSebastian Grimberg 
999ba59ac12SSebastian Grimberg   @ref Utility
1000ba59ac12SSebastian Grimberg **/
1001ba59ac12SSebastian Grimberg int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau, CeedInt m, CeedInt n) {
1002ba59ac12SSebastian Grimberg   CeedScalar v[m];
1003ba59ac12SSebastian Grimberg 
1004ba59ac12SSebastian Grimberg   // Check matrix shape
10056574a04fSJeremy L Thompson   CeedCheck(n <= m, ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute QR factorization with n > m");
1006ba59ac12SSebastian Grimberg 
1007ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) {
10081c66c397SJeremy L Thompson     CeedScalar sigma = 0.0;
10091c66c397SJeremy L Thompson 
1010ba59ac12SSebastian Grimberg     if (i >= m - 1) {  // last row of matrix, no reflection needed
1011ba59ac12SSebastian Grimberg       tau[i] = 0.;
1012ba59ac12SSebastian Grimberg       break;
1013ba59ac12SSebastian Grimberg     }
1014ba59ac12SSebastian Grimberg     // Calculate Householder vector, magnitude
1015ba59ac12SSebastian Grimberg     v[i] = mat[i + n * i];
1016ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < m; j++) {
1017ba59ac12SSebastian Grimberg       v[j] = mat[i + n * j];
1018ba59ac12SSebastian Grimberg       sigma += v[j] * v[j];
1019ba59ac12SSebastian Grimberg     }
10201c66c397SJeremy L Thompson     const CeedScalar norm = sqrt(v[i] * v[i] + sigma);  // norm of v[i:m]
10211c66c397SJeremy L Thompson     const CeedScalar R_ii = -copysign(norm, v[i]);
10221c66c397SJeremy L Thompson 
1023ba59ac12SSebastian Grimberg     v[i] -= R_ii;
1024ba59ac12SSebastian Grimberg     // norm of v[i:m] after modification above and scaling below
1025ba59ac12SSebastian Grimberg     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
1026ba59ac12SSebastian Grimberg     //   tau = 2 / (norm*norm)
1027ba59ac12SSebastian Grimberg     tau[i] = 2 * v[i] * v[i] / (v[i] * v[i] + sigma);
1028ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < m; j++) v[j] /= v[i];
1029ba59ac12SSebastian Grimberg 
1030ba59ac12SSebastian Grimberg     // Apply Householder reflector to lower right panel
1031ba59ac12SSebastian Grimberg     CeedHouseholderReflect(&mat[i * n + i + 1], &v[i], tau[i], m - i, n - i - 1, n, 1);
1032ba59ac12SSebastian Grimberg     // Save v
1033ba59ac12SSebastian Grimberg     mat[i + n * i] = R_ii;
1034ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < m; j++) mat[i + n * j] = v[j];
1035ba59ac12SSebastian Grimberg   }
1036ba59ac12SSebastian Grimberg   return CEED_ERROR_SUCCESS;
1037ba59ac12SSebastian Grimberg }
1038ba59ac12SSebastian Grimberg 
1039ba59ac12SSebastian Grimberg /**
1040ba59ac12SSebastian Grimberg   @brief Apply Householder Q matrix
1041ba59ac12SSebastian Grimberg 
1042ca94c3ddSJeremy 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$.
1043ba59ac12SSebastian Grimberg 
1044ba59ac12SSebastian Grimberg   @param[in,out] mat_A  Matrix to apply Householder Q to, in place
1045ba59ac12SSebastian Grimberg   @param[in]     mat_Q  Householder Q matrix
1046ba59ac12SSebastian Grimberg   @param[in]     tau    Householder scaling factors
1047ba59ac12SSebastian Grimberg   @param[in]     t_mode Transpose mode for application
1048ca94c3ddSJeremy L Thompson   @param[in]     m      Number of rows in `A`
1049ca94c3ddSJeremy L Thompson   @param[in]     n      Number of columns in `A`
1050ca94c3ddSJeremy L Thompson   @param[in]     k      Number of elementary reflectors in Q, `k < m`
1051ca94c3ddSJeremy L Thompson   @param[in]     row    Row stride in `A`
1052ca94c3ddSJeremy L Thompson   @param[in]     col    Col stride in `A`
1053ba59ac12SSebastian Grimberg 
1054ba59ac12SSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
1055ba59ac12SSebastian Grimberg 
1056c4e3f59bSSebastian Grimberg   @ref Utility
1057ba59ac12SSebastian Grimberg **/
1058ba59ac12SSebastian Grimberg int CeedHouseholderApplyQ(CeedScalar *mat_A, const CeedScalar *mat_Q, const CeedScalar *tau, CeedTransposeMode t_mode, CeedInt m, CeedInt n,
1059ba59ac12SSebastian Grimberg                           CeedInt k, CeedInt row, CeedInt col) {
1060ba59ac12SSebastian Grimberg   CeedScalar *v;
10611c66c397SJeremy L Thompson 
1062ba59ac12SSebastian Grimberg   CeedCall(CeedMalloc(m, &v));
1063ba59ac12SSebastian Grimberg   for (CeedInt ii = 0; ii < k; ii++) {
1064ba59ac12SSebastian Grimberg     CeedInt i = t_mode == CEED_TRANSPOSE ? ii : k - 1 - ii;
1065ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < m; j++) v[j] = mat_Q[j * k + i];
1066ba59ac12SSebastian Grimberg     // Apply Householder reflector (I - tau v v^T) collo_grad_1d^T
1067ba59ac12SSebastian Grimberg     CeedCall(CeedHouseholderReflect(&mat_A[i * row], &v[i], tau[i], m - i, n, row, col));
1068ba59ac12SSebastian Grimberg   }
1069ba59ac12SSebastian Grimberg   CeedCall(CeedFree(&v));
1070ba59ac12SSebastian Grimberg   return CEED_ERROR_SUCCESS;
1071ba59ac12SSebastian Grimberg }
1072ba59ac12SSebastian Grimberg 
1073ba59ac12SSebastian Grimberg /**
10742247a93fSRezgar Shakeri   @brief Return pseudoinverse of a matrix
10752247a93fSRezgar Shakeri 
10762247a93fSRezgar Shakeri   @param[in]     ceed      Ceed context for error handling
10772247a93fSRezgar Shakeri   @param[in]     mat       Row-major matrix to compute pseudoinverse of
10782247a93fSRezgar Shakeri   @param[in]     m         Number of rows
10792247a93fSRezgar Shakeri   @param[in]     n         Number of columns
10802247a93fSRezgar Shakeri   @param[out]    mat_pinv  Row-major pseudoinverse matrix
10812247a93fSRezgar Shakeri 
10822247a93fSRezgar Shakeri   @return An error code: 0 - success, otherwise - failure
10832247a93fSRezgar Shakeri 
10842247a93fSRezgar Shakeri   @ref Utility
10852247a93fSRezgar Shakeri **/
10861203703bSJeremy L Thompson int CeedMatrixPseudoinverse(Ceed ceed, const CeedScalar *mat, CeedInt m, CeedInt n, CeedScalar *mat_pinv) {
10872247a93fSRezgar Shakeri   CeedScalar *tau, *I, *mat_copy;
10882247a93fSRezgar Shakeri 
10892247a93fSRezgar Shakeri   CeedCall(CeedCalloc(m, &tau));
10902247a93fSRezgar Shakeri   CeedCall(CeedCalloc(m * m, &I));
10912247a93fSRezgar Shakeri   CeedCall(CeedCalloc(m * n, &mat_copy));
10922247a93fSRezgar Shakeri   memcpy(mat_copy, mat, m * n * sizeof mat[0]);
10932247a93fSRezgar Shakeri 
10942247a93fSRezgar Shakeri   // QR Factorization, mat = Q R
10952247a93fSRezgar Shakeri   CeedCall(CeedQRFactorization(ceed, mat_copy, tau, m, n));
10962247a93fSRezgar Shakeri 
10972247a93fSRezgar Shakeri   // -- Apply Q^T, I = Q^T * I
10982247a93fSRezgar Shakeri   for (CeedInt i = 0; i < m; i++) I[i * m + i] = 1.0;
10992247a93fSRezgar Shakeri   CeedCall(CeedHouseholderApplyQ(I, mat_copy, tau, CEED_TRANSPOSE, m, m, n, m, 1));
11002247a93fSRezgar Shakeri   // -- Apply R_inv, mat_pinv = R_inv * Q^T
11012247a93fSRezgar Shakeri   for (CeedInt j = 0; j < m; j++) {  // Column j
11022247a93fSRezgar Shakeri     mat_pinv[j + m * (n - 1)] = I[j + m * (n - 1)] / mat_copy[n * n - 1];
11032247a93fSRezgar Shakeri     for (CeedInt i = n - 2; i >= 0; i--) {  // Row i
11042247a93fSRezgar Shakeri       mat_pinv[j + m * i] = I[j + m * i];
11052247a93fSRezgar Shakeri       for (CeedInt k = i + 1; k < n; k++) mat_pinv[j + m * i] -= mat_copy[k + n * i] * mat_pinv[j + m * k];
11062247a93fSRezgar Shakeri       mat_pinv[j + m * i] /= mat_copy[i + n * i];
11072247a93fSRezgar Shakeri     }
11082247a93fSRezgar Shakeri   }
11092247a93fSRezgar Shakeri 
11102247a93fSRezgar Shakeri   // Cleanup
11112247a93fSRezgar Shakeri   CeedCall(CeedFree(&I));
11122247a93fSRezgar Shakeri   CeedCall(CeedFree(&tau));
11132247a93fSRezgar Shakeri   CeedCall(CeedFree(&mat_copy));
11142247a93fSRezgar Shakeri   return CEED_ERROR_SUCCESS;
11152247a93fSRezgar Shakeri }
11162247a93fSRezgar Shakeri 
11172247a93fSRezgar Shakeri /**
1118ba59ac12SSebastian Grimberg   @brief Return symmetric Schur decomposition of the symmetric matrix mat via symmetric QR factorization
1119ba59ac12SSebastian Grimberg 
1120ca94c3ddSJeremy L Thompson   @param[in]     ceed   `Ceed` context for error handling
1121ba59ac12SSebastian Grimberg   @param[in,out] mat    Row-major matrix to be factorized in place
1122ba59ac12SSebastian Grimberg   @param[out]    lambda Vector of length n of eigenvalues
1123ba59ac12SSebastian Grimberg   @param[in]     n      Number of rows/columns
1124ba59ac12SSebastian Grimberg 
1125ba59ac12SSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
1126ba59ac12SSebastian Grimberg 
1127ba59ac12SSebastian Grimberg   @ref Utility
1128ba59ac12SSebastian Grimberg **/
11292c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOff
11302c2ea1dbSJeremy L Thompson int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat, CeedScalar *lambda, CeedInt n) {
1131ba59ac12SSebastian Grimberg   // Check bounds for clang-tidy
11326574a04fSJeremy L Thompson   CeedCheck(n > 1, ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute symmetric Schur decomposition of scalars");
1133ba59ac12SSebastian Grimberg 
1134ba59ac12SSebastian Grimberg   CeedScalar v[n - 1], tau[n - 1], mat_T[n * n];
1135ba59ac12SSebastian Grimberg 
1136ba59ac12SSebastian Grimberg   // Copy mat to mat_T and set mat to I
1137ba59ac12SSebastian Grimberg   memcpy(mat_T, mat, n * n * sizeof(mat[0]));
1138ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) {
1139ba59ac12SSebastian Grimberg     for (CeedInt j = 0; j < n; j++) mat[j + n * i] = (i == j) ? 1 : 0;
1140ba59ac12SSebastian Grimberg   }
1141ba59ac12SSebastian Grimberg 
1142ba59ac12SSebastian Grimberg   // Reduce to tridiagonal
1143ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n - 1; i++) {
1144ba59ac12SSebastian Grimberg     // Calculate Householder vector, magnitude
1145ba59ac12SSebastian Grimberg     CeedScalar sigma = 0.0;
11461c66c397SJeremy L Thompson 
1147ba59ac12SSebastian Grimberg     v[i] = mat_T[i + n * (i + 1)];
1148ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < n - 1; j++) {
1149ba59ac12SSebastian Grimberg       v[j] = mat_T[i + n * (j + 1)];
1150ba59ac12SSebastian Grimberg       sigma += v[j] * v[j];
1151ba59ac12SSebastian Grimberg     }
11521c66c397SJeremy L Thompson     const CeedScalar norm = sqrt(v[i] * v[i] + sigma);  // norm of v[i:n-1]
11531c66c397SJeremy L Thompson     const CeedScalar R_ii = -copysign(norm, v[i]);
11541c66c397SJeremy L Thompson 
1155ba59ac12SSebastian Grimberg     v[i] -= R_ii;
1156ba59ac12SSebastian Grimberg     // norm of v[i:m] after modification above and scaling below
1157ba59ac12SSebastian Grimberg     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
1158ba59ac12SSebastian Grimberg     //   tau = 2 / (norm*norm)
1159ba59ac12SSebastian Grimberg     tau[i] = i == n - 2 ? 2 : 2 * v[i] * v[i] / (v[i] * v[i] + sigma);
1160ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < n - 1; j++) v[j] /= v[i];
1161ba59ac12SSebastian Grimberg 
1162ba59ac12SSebastian Grimberg     // Update sub and super diagonal
1163ba59ac12SSebastian Grimberg     for (CeedInt j = i + 2; j < n; j++) {
1164ba59ac12SSebastian Grimberg       mat_T[i + n * j] = 0;
1165ba59ac12SSebastian Grimberg       mat_T[j + n * i] = 0;
1166ba59ac12SSebastian Grimberg     }
1167ba59ac12SSebastian Grimberg     // Apply symmetric Householder reflector to lower right panel
1168ba59ac12SSebastian Grimberg     CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1);
1169ba59ac12SSebastian Grimberg     CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), 1, n);
1170ba59ac12SSebastian Grimberg 
1171ba59ac12SSebastian Grimberg     // Save v
1172ba59ac12SSebastian Grimberg     mat_T[i + n * (i + 1)] = R_ii;
1173ba59ac12SSebastian Grimberg     mat_T[(i + 1) + n * i] = R_ii;
1174ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < n - 1; j++) {
1175ba59ac12SSebastian Grimberg       mat_T[i + n * (j + 1)] = v[j];
1176ba59ac12SSebastian Grimberg     }
1177ba59ac12SSebastian Grimberg   }
1178ba59ac12SSebastian Grimberg   // Backwards accumulation of Q
1179ba59ac12SSebastian Grimberg   for (CeedInt i = n - 2; i >= 0; i--) {
1180ba59ac12SSebastian Grimberg     if (tau[i] > 0.0) {
1181ba59ac12SSebastian Grimberg       v[i] = 1;
1182ba59ac12SSebastian Grimberg       for (CeedInt j = i + 1; j < n - 1; j++) {
1183ba59ac12SSebastian Grimberg         v[j]                   = mat_T[i + n * (j + 1)];
1184ba59ac12SSebastian Grimberg         mat_T[i + n * (j + 1)] = 0;
1185ba59ac12SSebastian Grimberg       }
1186ba59ac12SSebastian Grimberg       CeedHouseholderReflect(&mat[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1);
1187ba59ac12SSebastian Grimberg     }
1188ba59ac12SSebastian Grimberg   }
1189ba59ac12SSebastian Grimberg 
1190ba59ac12SSebastian Grimberg   // Reduce sub and super diagonal
1191ba59ac12SSebastian Grimberg   CeedInt    p = 0, q = 0, itr = 0, max_itr = n * n * n * n;
1192ba59ac12SSebastian Grimberg   CeedScalar tol = CEED_EPSILON;
1193ba59ac12SSebastian Grimberg 
1194ba59ac12SSebastian Grimberg   while (itr < max_itr) {
1195ba59ac12SSebastian Grimberg     // Update p, q, size of reduced portions of diagonal
1196ba59ac12SSebastian Grimberg     p = 0;
1197ba59ac12SSebastian Grimberg     q = 0;
1198ba59ac12SSebastian Grimberg     for (CeedInt i = n - 2; i >= 0; i--) {
1199ba59ac12SSebastian Grimberg       if (fabs(mat_T[i + n * (i + 1)]) < tol) q += 1;
1200ba59ac12SSebastian Grimberg       else break;
1201ba59ac12SSebastian Grimberg     }
1202ba59ac12SSebastian Grimberg     for (CeedInt i = 0; i < n - q - 1; i++) {
1203ba59ac12SSebastian Grimberg       if (fabs(mat_T[i + n * (i + 1)]) < tol) p += 1;
1204ba59ac12SSebastian Grimberg       else break;
1205ba59ac12SSebastian Grimberg     }
1206ba59ac12SSebastian Grimberg     if (q == n - 1) break;  // Finished reducing
1207ba59ac12SSebastian Grimberg 
1208ba59ac12SSebastian Grimberg     // Reduce tridiagonal portion
1209ba59ac12SSebastian Grimberg     CeedScalar t_nn = mat_T[(n - 1 - q) + n * (n - 1 - q)], t_nnm1 = mat_T[(n - 2 - q) + n * (n - 1 - q)];
1210ba59ac12SSebastian Grimberg     CeedScalar d  = (mat_T[(n - 2 - q) + n * (n - 2 - q)] - t_nn) / 2;
1211ba59ac12SSebastian Grimberg     CeedScalar mu = t_nn - t_nnm1 * t_nnm1 / (d + copysign(sqrt(d * d + t_nnm1 * t_nnm1), d));
1212ba59ac12SSebastian Grimberg     CeedScalar x  = mat_T[p + n * p] - mu;
1213ba59ac12SSebastian Grimberg     CeedScalar z  = mat_T[p + n * (p + 1)];
12141c66c397SJeremy L Thompson 
1215ba59ac12SSebastian Grimberg     for (CeedInt k = p; k < n - q - 1; k++) {
1216ba59ac12SSebastian Grimberg       // Compute Givens rotation
1217ba59ac12SSebastian Grimberg       CeedScalar c = 1, s = 0;
12181c66c397SJeremy L Thompson 
1219ba59ac12SSebastian Grimberg       if (fabs(z) > tol) {
1220ba59ac12SSebastian Grimberg         if (fabs(z) > fabs(x)) {
12211c66c397SJeremy L Thompson           const CeedScalar tau = -x / z;
12221c66c397SJeremy L Thompson 
12231c66c397SJeremy L Thompson           s = 1 / sqrt(1 + tau * tau);
12241c66c397SJeremy L Thompson           c = s * tau;
1225ba59ac12SSebastian Grimberg         } else {
12261c66c397SJeremy L Thompson           const CeedScalar tau = -z / x;
12271c66c397SJeremy L Thompson 
12281c66c397SJeremy L Thompson           c = 1 / sqrt(1 + tau * tau);
12291c66c397SJeremy L Thompson           s = c * tau;
1230ba59ac12SSebastian Grimberg         }
1231ba59ac12SSebastian Grimberg       }
1232ba59ac12SSebastian Grimberg 
1233ba59ac12SSebastian Grimberg       // Apply Givens rotation to T
1234ba59ac12SSebastian Grimberg       CeedGivensRotation(mat_T, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n);
1235ba59ac12SSebastian Grimberg       CeedGivensRotation(mat_T, c, s, CEED_TRANSPOSE, k, k + 1, n, n);
1236ba59ac12SSebastian Grimberg 
1237ba59ac12SSebastian Grimberg       // Apply Givens rotation to Q
1238ba59ac12SSebastian Grimberg       CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n);
1239ba59ac12SSebastian Grimberg 
1240ba59ac12SSebastian Grimberg       // Update x, z
1241ba59ac12SSebastian Grimberg       if (k < n - q - 2) {
1242ba59ac12SSebastian Grimberg         x = mat_T[k + n * (k + 1)];
1243ba59ac12SSebastian Grimberg         z = mat_T[k + n * (k + 2)];
1244ba59ac12SSebastian Grimberg       }
1245ba59ac12SSebastian Grimberg     }
1246ba59ac12SSebastian Grimberg     itr++;
1247ba59ac12SSebastian Grimberg   }
1248ba59ac12SSebastian Grimberg 
1249ba59ac12SSebastian Grimberg   // Save eigenvalues
1250ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) lambda[i] = mat_T[i + n * i];
1251ba59ac12SSebastian Grimberg 
1252ba59ac12SSebastian Grimberg   // Check convergence
12536574a04fSJeremy L Thompson   CeedCheck(itr < max_itr || q > n, ceed, CEED_ERROR_MINOR, "Symmetric QR failed to converge");
1254ba59ac12SSebastian Grimberg   return CEED_ERROR_SUCCESS;
1255ba59ac12SSebastian Grimberg }
12562c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOn
1257ba59ac12SSebastian Grimberg 
1258ba59ac12SSebastian Grimberg /**
1259ba59ac12SSebastian Grimberg   @brief Return Simultaneous Diagonalization of two matrices.
1260ba59ac12SSebastian Grimberg 
1261ca94c3ddSJeremy L Thompson   This solves the generalized eigenvalue problem `A x = lambda B x`, where `A` and `B` are symmetric and `B` is positive definite.
1262ca94c3ddSJeremy L Thompson   We generate the matrix `X` and vector `Lambda` such that `X^T A X = Lambda` and `X^T B X = I`.
1263ca94c3ddSJeremy L Thompson   This is equivalent to the LAPACK routine 'sygv' with `TYPE = 1`.
1264ba59ac12SSebastian Grimberg 
1265ca94c3ddSJeremy L Thompson   @param[in]  ceed   `Ceed` context for error handling
1266ba59ac12SSebastian Grimberg   @param[in]  mat_A  Row-major matrix to be factorized with eigenvalues
1267ba59ac12SSebastian Grimberg   @param[in]  mat_B  Row-major matrix to be factorized to identity
1268ba59ac12SSebastian Grimberg   @param[out] mat_X  Row-major orthogonal matrix
1269ca94c3ddSJeremy L Thompson   @param[out] lambda Vector of length `n` of generalized eigenvalues
1270ba59ac12SSebastian Grimberg   @param[in]  n      Number of rows/columns
1271ba59ac12SSebastian Grimberg 
1272ba59ac12SSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
1273ba59ac12SSebastian Grimberg 
1274ba59ac12SSebastian Grimberg   @ref Utility
1275ba59ac12SSebastian Grimberg **/
12762c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOff
12772c2ea1dbSJeremy L Thompson int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *mat_A, CeedScalar *mat_B, CeedScalar *mat_X, CeedScalar *lambda, CeedInt n) {
1278ba59ac12SSebastian Grimberg   CeedScalar *mat_C, *mat_G, *vec_D;
12791c66c397SJeremy L Thompson 
1280ba59ac12SSebastian Grimberg   CeedCall(CeedCalloc(n * n, &mat_C));
1281ba59ac12SSebastian Grimberg   CeedCall(CeedCalloc(n * n, &mat_G));
1282ba59ac12SSebastian Grimberg   CeedCall(CeedCalloc(n, &vec_D));
1283ba59ac12SSebastian Grimberg 
1284ba59ac12SSebastian Grimberg   // Compute B = G D G^T
1285ba59ac12SSebastian Grimberg   memcpy(mat_G, mat_B, n * n * sizeof(mat_B[0]));
1286ba59ac12SSebastian Grimberg   CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_G, vec_D, n));
1287ba59ac12SSebastian Grimberg 
1288ba59ac12SSebastian Grimberg   // Sort eigenvalues
1289ba59ac12SSebastian Grimberg   for (CeedInt i = n - 1; i >= 0; i--) {
1290ba59ac12SSebastian Grimberg     for (CeedInt j = 0; j < i; j++) {
1291ba59ac12SSebastian Grimberg       if (fabs(vec_D[j]) > fabs(vec_D[j + 1])) {
12921c66c397SJeremy L Thompson         CeedScalarSwap(vec_D[j], vec_D[j + 1]);
12931c66c397SJeremy L Thompson         for (CeedInt k = 0; k < n; k++) CeedScalarSwap(mat_G[k * n + j], mat_G[k * n + j + 1]);
1294ba59ac12SSebastian Grimberg       }
1295ba59ac12SSebastian Grimberg     }
1296ba59ac12SSebastian Grimberg   }
1297ba59ac12SSebastian Grimberg 
1298ba59ac12SSebastian Grimberg   // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T
1299ba59ac12SSebastian Grimberg   //           = D^-1/2 G^T A G D^-1/2
1300ba59ac12SSebastian Grimberg   // -- D = D^-1/2
1301ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) vec_D[i] = 1. / sqrt(vec_D[i]);
1302ba59ac12SSebastian Grimberg   // -- G = G D^-1/2
1303ba59ac12SSebastian Grimberg   // -- C = D^-1/2 G^T
1304ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) {
1305ba59ac12SSebastian Grimberg     for (CeedInt j = 0; j < n; j++) {
1306ba59ac12SSebastian Grimberg       mat_G[i * n + j] *= vec_D[j];
1307ba59ac12SSebastian Grimberg       mat_C[j * n + i] = mat_G[i * n + j];
1308ba59ac12SSebastian Grimberg     }
1309ba59ac12SSebastian Grimberg   }
1310ba59ac12SSebastian Grimberg   // -- X = (D^-1/2 G^T) A
1311ba59ac12SSebastian Grimberg   CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_C, (const CeedScalar *)mat_A, mat_X, n, n, n));
1312ba59ac12SSebastian Grimberg   // -- C = (D^-1/2 G^T A) (G D^-1/2)
1313ba59ac12SSebastian Grimberg   CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_X, (const CeedScalar *)mat_G, mat_C, n, n, n));
1314ba59ac12SSebastian Grimberg 
1315ba59ac12SSebastian Grimberg   // Compute Q^T C Q = lambda
1316ba59ac12SSebastian Grimberg   CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_C, lambda, n));
1317ba59ac12SSebastian Grimberg 
1318ba59ac12SSebastian Grimberg   // Sort eigenvalues
1319ba59ac12SSebastian Grimberg   for (CeedInt i = n - 1; i >= 0; i--) {
1320ba59ac12SSebastian Grimberg     for (CeedInt j = 0; j < i; j++) {
1321ba59ac12SSebastian Grimberg       if (fabs(lambda[j]) > fabs(lambda[j + 1])) {
13221c66c397SJeremy L Thompson         CeedScalarSwap(lambda[j], lambda[j + 1]);
13231c66c397SJeremy L Thompson         for (CeedInt k = 0; k < n; k++) CeedScalarSwap(mat_C[k * n + j], mat_C[k * n + j + 1]);
1324ba59ac12SSebastian Grimberg       }
1325ba59ac12SSebastian Grimberg     }
1326ba59ac12SSebastian Grimberg   }
1327ba59ac12SSebastian Grimberg 
1328ba59ac12SSebastian Grimberg   // Set X = (G D^1/2)^-T Q
1329ba59ac12SSebastian Grimberg   //       = G D^-1/2 Q
1330ba59ac12SSebastian Grimberg   CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_G, (const CeedScalar *)mat_C, mat_X, n, n, n));
1331ba59ac12SSebastian Grimberg 
1332ba59ac12SSebastian Grimberg   // Cleanup
1333ba59ac12SSebastian Grimberg   CeedCall(CeedFree(&mat_C));
1334ba59ac12SSebastian Grimberg   CeedCall(CeedFree(&mat_G));
1335ba59ac12SSebastian Grimberg   CeedCall(CeedFree(&vec_D));
1336ba59ac12SSebastian Grimberg   return CEED_ERROR_SUCCESS;
1337ba59ac12SSebastian Grimberg }
13382c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOn
1339ba59ac12SSebastian Grimberg 
13407a982d89SJeremy L. Thompson /// @}
13417a982d89SJeremy L. Thompson 
13427a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
13437a982d89SJeremy L. Thompson /// CeedBasis Public API
13447a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
13457a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser
1346d7b241e6Sjeremylt /// @{
1347d7b241e6Sjeremylt 
1348b11c1e72Sjeremylt /**
1349ca94c3ddSJeremy L Thompson   @brief Create a tensor-product basis for \f$H^1\f$ discretizations
1350b11c1e72Sjeremylt 
1351ca94c3ddSJeremy L Thompson   @param[in]  ceed        `Ceed` object used to create the `CeedBasis`
1352ea61e9acSJeremy L Thompson   @param[in]  dim         Topological dimension
1353ea61e9acSJeremy L Thompson   @param[in]  num_comp    Number of field components (1 for scalar fields)
1354ea61e9acSJeremy L Thompson   @param[in]  P_1d        Number of nodes in one dimension
1355ea61e9acSJeremy L Thompson   @param[in]  Q_1d        Number of quadrature points in one dimension
1356ca94c3ddSJeremy L Thompson   @param[in]  interp_1d   Row-major (`Q_1d * P_1d`) matrix expressing the values of nodal basis functions at quadrature points
1357ca94c3ddSJeremy L Thompson   @param[in]  grad_1d     Row-major (`Q_1d * P_1d`) matrix expressing derivatives of nodal basis functions at quadrature points
1358ca94c3ddSJeremy 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]`
1359ca94c3ddSJeremy L Thompson   @param[in]  q_weight_1d Array of length `Q_1d` holding the quadrature weights on the reference element
1360ca94c3ddSJeremy L Thompson   @param[out] basis       Address of the variable where the newly created `CeedBasis` will be stored
1361b11c1e72Sjeremylt 
1362b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1363dfdf5a53Sjeremylt 
13647a982d89SJeremy L. Thompson   @ref User
1365b11c1e72Sjeremylt **/
13662b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d,
13672b730f8bSJeremy L Thompson                             const CeedScalar *grad_1d, const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis *basis) {
13685fe0d4faSjeremylt   if (!ceed->BasisCreateTensorH1) {
13695fe0d4faSjeremylt     Ceed delegate;
13706574a04fSJeremy L Thompson 
13712b730f8bSJeremy L Thompson     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
13721ef3a2a9SJeremy L Thompson     CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateTensorH1");
13732b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateTensorH1(delegate, dim, num_comp, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis));
13749bc66399SJeremy L Thompson     CeedCall(CeedDestroy(&delegate));
1375e15f9bd0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
13765fe0d4faSjeremylt   }
1377e15f9bd0SJeremy L Thompson 
1378ca94c3ddSJeremy L Thompson   CeedCheck(dim > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis dimension must be a positive value");
1379ca94c3ddSJeremy L Thompson   CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 component");
1380ca94c3ddSJeremy L Thompson   CeedCheck(P_1d > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 node");
1381ca94c3ddSJeremy L Thompson   CeedCheck(Q_1d > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 quadrature point");
1382227444bfSJeremy L Thompson 
13832b730f8bSJeremy L Thompson   CeedElemTopology topo = dim == 1 ? CEED_TOPOLOGY_LINE : dim == 2 ? CEED_TOPOLOGY_QUAD : CEED_TOPOLOGY_HEX;
1384e15f9bd0SJeremy L Thompson 
13852b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(1, basis));
1386db002c03SJeremy L Thompson   CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed));
1387d1d35e2fSjeremylt   (*basis)->ref_count       = 1;
13886402da51SJeremy L Thompson   (*basis)->is_tensor_basis = true;
1389d7b241e6Sjeremylt   (*basis)->dim             = dim;
1390d99fa3c5SJeremy L Thompson   (*basis)->topo            = topo;
1391d1d35e2fSjeremylt   (*basis)->num_comp        = num_comp;
1392d1d35e2fSjeremylt   (*basis)->P_1d            = P_1d;
1393d1d35e2fSjeremylt   (*basis)->Q_1d            = Q_1d;
1394d1d35e2fSjeremylt   (*basis)->P               = CeedIntPow(P_1d, dim);
1395d1d35e2fSjeremylt   (*basis)->Q               = CeedIntPow(Q_1d, dim);
1396c4e3f59bSSebastian Grimberg   (*basis)->fe_space        = CEED_FE_SPACE_H1;
13972b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d, &(*basis)->q_ref_1d));
13982b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d, &(*basis)->q_weight_1d));
1399ff3a0f91SJeremy L Thompson   if (q_ref_1d) memcpy((*basis)->q_ref_1d, q_ref_1d, Q_1d * sizeof(q_ref_1d[0]));
14002b730f8bSJeremy L Thompson   if (q_weight_1d) memcpy((*basis)->q_weight_1d, q_weight_1d, Q_1d * sizeof(q_weight_1d[0]));
14012b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->interp_1d));
14022b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->grad_1d));
14032b730f8bSJeremy L Thompson   if (interp_1d) memcpy((*basis)->interp_1d, interp_1d, Q_1d * P_1d * sizeof(interp_1d[0]));
1404ff3a0f91SJeremy L Thompson   if (grad_1d) memcpy((*basis)->grad_1d, grad_1d, Q_1d * P_1d * sizeof(grad_1d[0]));
14052b730f8bSJeremy L Thompson   CeedCall(ceed->BasisCreateTensorH1(dim, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, *basis));
1406e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1407d7b241e6Sjeremylt }
1408d7b241e6Sjeremylt 
1409b11c1e72Sjeremylt /**
1410ca94c3ddSJeremy L Thompson   @brief Create a tensor-product \f$H^1\f$ Lagrange basis
1411b11c1e72Sjeremylt 
1412ca94c3ddSJeremy L Thompson   @param[in]  ceed      `Ceed` object used to create the `CeedBasis`
1413ea61e9acSJeremy L Thompson   @param[in]  dim       Topological dimension of element
1414ea61e9acSJeremy L Thompson   @param[in]  num_comp  Number of field components (1 for scalar fields)
1415ea61e9acSJeremy L Thompson   @param[in]  P         Number of Gauss-Lobatto nodes in one dimension.
1416ca94c3ddSJeremy L Thompson                           The polynomial degree of the resulting `Q_k` element is `k = P - 1`.
1417ea61e9acSJeremy L Thompson   @param[in]  Q         Number of quadrature points in one dimension.
1418ca94c3ddSJeremy L Thompson   @param[in]  quad_mode Distribution of the `Q` quadrature points (affects order of accuracy for the quadrature)
1419ca94c3ddSJeremy L Thompson   @param[out] basis     Address of the variable where the newly created `CeedBasis` will be stored
1420b11c1e72Sjeremylt 
1421b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1422dfdf5a53Sjeremylt 
14237a982d89SJeremy L. Thompson   @ref User
1424b11c1e72Sjeremylt **/
14252b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P, CeedInt Q, CeedQuadMode quad_mode, CeedBasis *basis) {
1426d7b241e6Sjeremylt   // Allocate
1427c8c3fa7dSJeremy L Thompson   int        ierr = CEED_ERROR_SUCCESS;
14282b730f8bSJeremy L Thompson   CeedScalar c1, c2, c3, c4, dx, *nodes, *interp_1d, *grad_1d, *q_ref_1d, *q_weight_1d;
14294d537eeaSYohann 
1430ca94c3ddSJeremy L Thompson   CeedCheck(dim > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis dimension must be a positive value");
1431ca94c3ddSJeremy L Thompson   CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 component");
1432ca94c3ddSJeremy L Thompson   CeedCheck(P > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 node");
1433ca94c3ddSJeremy L Thompson   CeedCheck(Q > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 quadrature point");
1434227444bfSJeremy L Thompson 
1435e15f9bd0SJeremy L Thompson   // Get Nodes and Weights
14362b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(P * Q, &interp_1d));
14372b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(P * Q, &grad_1d));
14382b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(P, &nodes));
14392b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q, &q_ref_1d));
14402b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q, &q_weight_1d));
14412b730f8bSJeremy L Thompson   if (CeedLobattoQuadrature(P, nodes, NULL) != CEED_ERROR_SUCCESS) goto cleanup;
1442d1d35e2fSjeremylt   switch (quad_mode) {
1443d7b241e6Sjeremylt     case CEED_GAUSS:
1444d1d35e2fSjeremylt       ierr = CeedGaussQuadrature(Q, q_ref_1d, q_weight_1d);
1445d7b241e6Sjeremylt       break;
1446d7b241e6Sjeremylt     case CEED_GAUSS_LOBATTO:
1447d1d35e2fSjeremylt       ierr = CeedLobattoQuadrature(Q, q_ref_1d, q_weight_1d);
1448d7b241e6Sjeremylt       break;
1449d7b241e6Sjeremylt   }
14502b730f8bSJeremy L Thompson   if (ierr != CEED_ERROR_SUCCESS) goto cleanup;
1451e15f9bd0SJeremy L Thompson 
1452d7b241e6Sjeremylt   // Build B, D matrix
1453d7b241e6Sjeremylt   // Fornberg, 1998
1454c8c3fa7dSJeremy L Thompson   for (CeedInt i = 0; i < Q; i++) {
1455d7b241e6Sjeremylt     c1                   = 1.0;
1456d1d35e2fSjeremylt     c3                   = nodes[0] - q_ref_1d[i];
1457d1d35e2fSjeremylt     interp_1d[i * P + 0] = 1.0;
1458c8c3fa7dSJeremy L Thompson     for (CeedInt j = 1; j < P; j++) {
1459d7b241e6Sjeremylt       c2 = 1.0;
1460d7b241e6Sjeremylt       c4 = c3;
1461d1d35e2fSjeremylt       c3 = nodes[j] - q_ref_1d[i];
1462c8c3fa7dSJeremy L Thompson       for (CeedInt k = 0; k < j; k++) {
1463d7b241e6Sjeremylt         dx = nodes[j] - nodes[k];
1464d7b241e6Sjeremylt         c2 *= dx;
1465d7b241e6Sjeremylt         if (k == j - 1) {
1466d1d35e2fSjeremylt           grad_1d[i * P + j]   = c1 * (interp_1d[i * P + k] - c4 * grad_1d[i * P + k]) / c2;
1467d1d35e2fSjeremylt           interp_1d[i * P + j] = -c1 * c4 * interp_1d[i * P + k] / c2;
1468d7b241e6Sjeremylt         }
1469d1d35e2fSjeremylt         grad_1d[i * P + k]   = (c3 * grad_1d[i * P + k] - interp_1d[i * P + k]) / dx;
1470d1d35e2fSjeremylt         interp_1d[i * P + k] = c3 * interp_1d[i * P + k] / dx;
1471d7b241e6Sjeremylt       }
1472d7b241e6Sjeremylt       c1 = c2;
1473d7b241e6Sjeremylt     }
1474d7b241e6Sjeremylt   }
14759ac7b42eSJeremy L Thompson   // Pass to CeedBasisCreateTensorH1
14762b730f8bSJeremy L Thompson   CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P, Q, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis));
1477e15f9bd0SJeremy L Thompson cleanup:
14782b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_1d));
14792b730f8bSJeremy L Thompson   CeedCall(CeedFree(&grad_1d));
14802b730f8bSJeremy L Thompson   CeedCall(CeedFree(&nodes));
14812b730f8bSJeremy L Thompson   CeedCall(CeedFree(&q_ref_1d));
14822b730f8bSJeremy L Thompson   CeedCall(CeedFree(&q_weight_1d));
1483e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1484d7b241e6Sjeremylt }
1485d7b241e6Sjeremylt 
1486b11c1e72Sjeremylt /**
1487ca94c3ddSJeremy L Thompson   @brief Create a non tensor-product basis for \f$H^1\f$ discretizations
1488a8de75f0Sjeremylt 
1489ca94c3ddSJeremy L Thompson   @param[in]  ceed      `Ceed` object used to create the `CeedBasis`
1490e00f3be8SJames Wright   @param[in]  topo      Topology of element, e.g. hypercube, simplex, etc
1491ea61e9acSJeremy L Thompson   @param[in]  num_comp  Number of field components (1 for scalar fields)
1492ea61e9acSJeremy L Thompson   @param[in]  num_nodes Total number of nodes
1493ea61e9acSJeremy L Thompson   @param[in]  num_qpts  Total number of quadrature points
1494ca94c3ddSJeremy L Thompson   @param[in]  interp    Row-major (`num_qpts * num_nodes`) matrix expressing the values of nodal basis functions at quadrature points
1495ca94c3ddSJeremy L Thompson   @param[in]  grad      Row-major (`dim * num_qpts * num_nodes`) matrix expressing derivatives of nodal basis functions at quadrature points
1496ca94c3ddSJeremy L Thompson   @param[in]  q_ref     Array of length `num_qpts` * dim holding the locations of quadrature points on the reference element
1497ca94c3ddSJeremy L Thompson   @param[in]  q_weight  Array of length `num_qpts` holding the quadrature weights on the reference element
1498ca94c3ddSJeremy L Thompson   @param[out] basis     Address of the variable where the newly created `CeedBasis` will be stored
1499a8de75f0Sjeremylt 
1500a8de75f0Sjeremylt   @return An error code: 0 - success, otherwise - failure
1501a8de75f0Sjeremylt 
15027a982d89SJeremy L. Thompson   @ref User
1503a8de75f0Sjeremylt **/
15042b730f8bSJeremy L Thompson int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
15052b730f8bSJeremy L Thompson                       const CeedScalar *grad, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) {
1506d1d35e2fSjeremylt   CeedInt P = num_nodes, Q = num_qpts, dim = 0;
1507a8de75f0Sjeremylt 
15085fe0d4faSjeremylt   if (!ceed->BasisCreateH1) {
15095fe0d4faSjeremylt     Ceed delegate;
15106574a04fSJeremy L Thompson 
15112b730f8bSJeremy L Thompson     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
15121ef3a2a9SJeremy L Thompson     CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateH1");
15132b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateH1(delegate, topo, num_comp, num_nodes, num_qpts, interp, grad, q_ref, q_weight, basis));
15149bc66399SJeremy L Thompson     CeedCall(CeedDestroy(&delegate));
1515e15f9bd0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
15165fe0d4faSjeremylt   }
15175fe0d4faSjeremylt 
1518ca94c3ddSJeremy L Thompson   CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 component");
1519ca94c3ddSJeremy L Thompson   CeedCheck(num_nodes > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 node");
1520ca94c3ddSJeremy L Thompson   CeedCheck(num_qpts > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 quadrature point");
1521227444bfSJeremy L Thompson 
15222b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetTopologyDimension(topo, &dim));
1523a8de75f0Sjeremylt 
1524db002c03SJeremy L Thompson   CeedCall(CeedCalloc(1, basis));
1525db002c03SJeremy L Thompson   CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed));
1526d1d35e2fSjeremylt   (*basis)->ref_count       = 1;
15276402da51SJeremy L Thompson   (*basis)->is_tensor_basis = false;
1528a8de75f0Sjeremylt   (*basis)->dim             = dim;
1529d99fa3c5SJeremy L Thompson   (*basis)->topo            = topo;
1530d1d35e2fSjeremylt   (*basis)->num_comp        = num_comp;
1531a8de75f0Sjeremylt   (*basis)->P               = P;
1532a8de75f0Sjeremylt   (*basis)->Q               = Q;
1533c4e3f59bSSebastian Grimberg   (*basis)->fe_space        = CEED_FE_SPACE_H1;
15342b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q * dim, &(*basis)->q_ref_1d));
15352b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q, &(*basis)->q_weight_1d));
1536ff3a0f91SJeremy L Thompson   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0]));
1537ff3a0f91SJeremy L Thompson   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0]));
15382b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q * P, &(*basis)->interp));
15392b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(dim * Q * P, &(*basis)->grad));
1540ff3a0f91SJeremy L Thompson   if (interp) memcpy((*basis)->interp, interp, Q * P * sizeof(interp[0]));
1541ff3a0f91SJeremy L Thompson   if (grad) memcpy((*basis)->grad, grad, dim * Q * P * sizeof(grad[0]));
15422b730f8bSJeremy L Thompson   CeedCall(ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, q_ref, q_weight, *basis));
1543e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1544a8de75f0Sjeremylt }
1545a8de75f0Sjeremylt 
1546a8de75f0Sjeremylt /**
1547859c15bbSJames Wright   @brief Create a non tensor-product basis for \f$H(\mathrm{div})\f$ discretizations
154850c301a5SRezgar Shakeri 
1549ca94c3ddSJeremy L Thompson   @param[in]  ceed      `Ceed` object used to create the `CeedBasis`
1550ea61e9acSJeremy 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
1551ea61e9acSJeremy L Thompson   @param[in]  num_comp  Number of components (usually 1 for vectors in H(div) bases)
1552ca94c3ddSJeremy L Thompson   @param[in]  num_nodes Total number of nodes (DoFs per element)
1553ea61e9acSJeremy L Thompson   @param[in]  num_qpts  Total number of quadrature points
1554ca94c3ddSJeremy L Thompson   @param[in]  interp    Row-major (`dim * num_qpts * num_nodes`) matrix expressing the values of basis functions at quadrature points
1555ca94c3ddSJeremy L Thompson   @param[in]  div       Row-major (`num_qpts * num_nodes`) matrix expressing divergence of basis functions at quadrature points
1556ca94c3ddSJeremy L Thompson   @param[in]  q_ref     Array of length `num_qpts` * dim holding the locations of quadrature points on the reference element
1557ca94c3ddSJeremy L Thompson   @param[in]  q_weight  Array of length `num_qpts` holding the quadrature weights on the reference element
1558ca94c3ddSJeremy L Thompson   @param[out] basis     Address of the variable where the newly created `CeedBasis` will be stored
155950c301a5SRezgar Shakeri 
156050c301a5SRezgar Shakeri   @return An error code: 0 - success, otherwise - failure
156150c301a5SRezgar Shakeri 
156250c301a5SRezgar Shakeri   @ref User
156350c301a5SRezgar Shakeri **/
15642b730f8bSJeremy L Thompson int CeedBasisCreateHdiv(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
15652b730f8bSJeremy L Thompson                         const CeedScalar *div, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) {
156650c301a5SRezgar Shakeri   CeedInt Q = num_qpts, P = num_nodes, dim = 0;
1567c4e3f59bSSebastian Grimberg 
156850c301a5SRezgar Shakeri   if (!ceed->BasisCreateHdiv) {
156950c301a5SRezgar Shakeri     Ceed delegate;
15706574a04fSJeremy L Thompson 
15712b730f8bSJeremy L Thompson     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
15726574a04fSJeremy L Thompson     CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateHdiv");
15732b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateHdiv(delegate, topo, num_comp, num_nodes, num_qpts, interp, div, q_ref, q_weight, basis));
15749bc66399SJeremy L Thompson     CeedCall(CeedDestroy(&delegate));
157550c301a5SRezgar Shakeri     return CEED_ERROR_SUCCESS;
157650c301a5SRezgar Shakeri   }
157750c301a5SRezgar Shakeri 
1578ca94c3ddSJeremy L Thompson   CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 component");
1579ca94c3ddSJeremy L Thompson   CeedCheck(num_nodes > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 node");
1580ca94c3ddSJeremy L Thompson   CeedCheck(num_qpts > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 quadrature point");
1581227444bfSJeremy L Thompson 
1582c4e3f59bSSebastian Grimberg   CeedCall(CeedBasisGetTopologyDimension(topo, &dim));
1583c4e3f59bSSebastian Grimberg 
1584db002c03SJeremy L Thompson   CeedCall(CeedCalloc(1, basis));
1585db002c03SJeremy L Thompson   CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed));
158650c301a5SRezgar Shakeri   (*basis)->ref_count       = 1;
15876402da51SJeremy L Thompson   (*basis)->is_tensor_basis = false;
158850c301a5SRezgar Shakeri   (*basis)->dim             = dim;
158950c301a5SRezgar Shakeri   (*basis)->topo            = topo;
159050c301a5SRezgar Shakeri   (*basis)->num_comp        = num_comp;
159150c301a5SRezgar Shakeri   (*basis)->P               = P;
159250c301a5SRezgar Shakeri   (*basis)->Q               = Q;
1593c4e3f59bSSebastian Grimberg   (*basis)->fe_space        = CEED_FE_SPACE_HDIV;
15942b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q * dim, &(*basis)->q_ref_1d));
15952b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q, &(*basis)->q_weight_1d));
159650c301a5SRezgar Shakeri   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0]));
159750c301a5SRezgar Shakeri   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0]));
15982b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(dim * Q * P, &(*basis)->interp));
15992b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q * P, &(*basis)->div));
160050c301a5SRezgar Shakeri   if (interp) memcpy((*basis)->interp, interp, dim * Q * P * sizeof(interp[0]));
160150c301a5SRezgar Shakeri   if (div) memcpy((*basis)->div, div, Q * P * sizeof(div[0]));
16022b730f8bSJeremy L Thompson   CeedCall(ceed->BasisCreateHdiv(topo, dim, P, Q, interp, div, q_ref, q_weight, *basis));
160350c301a5SRezgar Shakeri   return CEED_ERROR_SUCCESS;
160450c301a5SRezgar Shakeri }
160550c301a5SRezgar Shakeri 
160650c301a5SRezgar Shakeri /**
16074385fb7fSSebastian Grimberg   @brief Create a non tensor-product basis for \f$H(\mathrm{curl})\f$ discretizations
1608c4e3f59bSSebastian Grimberg 
1609ca94c3ddSJeremy L Thompson   @param[in]  ceed      `Ceed` object used to create the `CeedBasis`
1610c4e3f59bSSebastian Grimberg   @param[in]  topo      Topology of element (`CEED_TOPOLOGY_QUAD`, `CEED_TOPOLOGY_PRISM`, etc.), dimension of which is used in some array sizes below
1611ca94c3ddSJeremy L Thompson   @param[in]  num_comp  Number of components (usually 1 for vectors in \f$H(\mathrm{curl})\f$ bases)
1612ca94c3ddSJeremy L Thompson   @param[in]  num_nodes Total number of nodes (DoFs per element)
1613c4e3f59bSSebastian Grimberg   @param[in]  num_qpts  Total number of quadrature points
1614ca94c3ddSJeremy L Thompson   @param[in]  interp    Row-major (`dim * num_qpts * num_nodes`) matrix expressing the values of basis functions at quadrature points
1615ca94c3ddSJeremy 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
1616ca94c3ddSJeremy L Thompson   @param[in]  q_ref     Array of length `num_qpts * dim` holding the locations of quadrature points on the reference element
1617ca94c3ddSJeremy L Thompson   @param[in]  q_weight  Array of length `num_qpts` holding the quadrature weights on the reference element
1618ca94c3ddSJeremy L Thompson   @param[out] basis     Address of the variable where the newly created `CeedBasis` will be stored
1619c4e3f59bSSebastian Grimberg 
1620c4e3f59bSSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
1621c4e3f59bSSebastian Grimberg 
1622c4e3f59bSSebastian Grimberg   @ref User
1623c4e3f59bSSebastian Grimberg **/
1624c4e3f59bSSebastian Grimberg int CeedBasisCreateHcurl(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
1625c4e3f59bSSebastian Grimberg                          const CeedScalar *curl, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) {
1626c4e3f59bSSebastian Grimberg   CeedInt Q = num_qpts, P = num_nodes, dim = 0, curl_comp = 0;
1627c4e3f59bSSebastian Grimberg 
1628d075f50bSSebastian Grimberg   if (!ceed->BasisCreateHcurl) {
1629c4e3f59bSSebastian Grimberg     Ceed delegate;
16306574a04fSJeremy L Thompson 
1631c4e3f59bSSebastian Grimberg     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
16326574a04fSJeremy L Thompson     CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateHcurl");
1633c4e3f59bSSebastian Grimberg     CeedCall(CeedBasisCreateHcurl(delegate, topo, num_comp, num_nodes, num_qpts, interp, curl, q_ref, q_weight, basis));
16349bc66399SJeremy L Thompson     CeedCall(CeedDestroy(&delegate));
1635c4e3f59bSSebastian Grimberg     return CEED_ERROR_SUCCESS;
1636c4e3f59bSSebastian Grimberg   }
1637c4e3f59bSSebastian Grimberg 
1638ca94c3ddSJeremy L Thompson   CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 component");
1639ca94c3ddSJeremy L Thompson   CeedCheck(num_nodes > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 node");
1640ca94c3ddSJeremy L Thompson   CeedCheck(num_qpts > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 quadrature point");
1641c4e3f59bSSebastian Grimberg 
1642c4e3f59bSSebastian Grimberg   CeedCall(CeedBasisGetTopologyDimension(topo, &dim));
1643c4e3f59bSSebastian Grimberg   curl_comp = (dim < 3) ? 1 : dim;
1644c4e3f59bSSebastian Grimberg 
1645db002c03SJeremy L Thompson   CeedCall(CeedCalloc(1, basis));
1646db002c03SJeremy L Thompson   CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed));
1647c4e3f59bSSebastian Grimberg   (*basis)->ref_count       = 1;
16486402da51SJeremy L Thompson   (*basis)->is_tensor_basis = false;
1649c4e3f59bSSebastian Grimberg   (*basis)->dim             = dim;
1650c4e3f59bSSebastian Grimberg   (*basis)->topo            = topo;
1651c4e3f59bSSebastian Grimberg   (*basis)->num_comp        = num_comp;
1652c4e3f59bSSebastian Grimberg   (*basis)->P               = P;
1653c4e3f59bSSebastian Grimberg   (*basis)->Q               = Q;
1654c4e3f59bSSebastian Grimberg   (*basis)->fe_space        = CEED_FE_SPACE_HCURL;
1655c4e3f59bSSebastian Grimberg   CeedCall(CeedMalloc(Q * dim, &(*basis)->q_ref_1d));
1656c4e3f59bSSebastian Grimberg   CeedCall(CeedMalloc(Q, &(*basis)->q_weight_1d));
1657c4e3f59bSSebastian Grimberg   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0]));
1658c4e3f59bSSebastian Grimberg   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0]));
1659c4e3f59bSSebastian Grimberg   CeedCall(CeedMalloc(dim * Q * P, &(*basis)->interp));
1660c4e3f59bSSebastian Grimberg   CeedCall(CeedMalloc(curl_comp * Q * P, &(*basis)->curl));
1661c4e3f59bSSebastian Grimberg   if (interp) memcpy((*basis)->interp, interp, dim * Q * P * sizeof(interp[0]));
1662c4e3f59bSSebastian Grimberg   if (curl) memcpy((*basis)->curl, curl, curl_comp * Q * P * sizeof(curl[0]));
1663c4e3f59bSSebastian Grimberg   CeedCall(ceed->BasisCreateHcurl(topo, dim, P, Q, interp, curl, q_ref, q_weight, *basis));
1664c4e3f59bSSebastian Grimberg   return CEED_ERROR_SUCCESS;
1665c4e3f59bSSebastian Grimberg }
1666c4e3f59bSSebastian Grimberg 
1667c4e3f59bSSebastian Grimberg /**
1668ca94c3ddSJeremy L Thompson   @brief Create a `CeedBasis` for projection from the nodes of `basis_from` to the nodes of `basis_to`.
1669ba59ac12SSebastian Grimberg 
1670ca94c3ddSJeremy L Thompson   Only @ref CEED_EVAL_INTERP will be valid for the new basis, `basis_project`.
1671ca94c3ddSJeremy L Thompson   For \f$H^1\f$ spaces, @ref CEED_EVAL_GRAD will also be valid.
1672ca94c3ddSJeremy L Thompson   The interpolation is given by `interp_project = interp_to^+ * interp_from`, where the pseudoinverse `interp_to^+` is given by QR factorization.
1673ca94c3ddSJeremy L Thompson   The gradient (for the \f$H^1\f$ case) is given by `grad_project = interp_to^+ * grad_from`.
167415ad3917SSebastian Grimberg 
167515ad3917SSebastian Grimberg   Note: `basis_from` and `basis_to` must have compatible quadrature spaces.
167615ad3917SSebastian Grimberg 
16779fd66db6SSebastian Grimberg   Note: `basis_project` will have the same number of components as `basis_from`, regardless of the number of components that `basis_to` has.
16789fd66db6SSebastian Grimberg         If `basis_from` has 3 components and `basis_to` has 5 components, then `basis_project` will have 3 components.
1679f113e5dcSJeremy L Thompson 
1680e104ad11SJames Wright   Note: If either `basis_from` or `basis_to` are non-tensor, then `basis_project` will also be non-tensor
1681e104ad11SJames Wright 
1682ca94c3ddSJeremy L Thompson   @param[in]  basis_from    `CeedBasis` to prolong from
1683ca94c3ddSJeremy L Thompson   @param[in]  basis_to      `CeedBasis` to prolong to
1684ca94c3ddSJeremy L Thompson   @param[out] basis_project Address of the variable where the newly created `CeedBasis` will be stored
1685f113e5dcSJeremy L Thompson 
1686f113e5dcSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1687f113e5dcSJeremy L Thompson 
1688f113e5dcSJeremy L Thompson   @ref User
1689f113e5dcSJeremy L Thompson **/
16902b730f8bSJeremy L Thompson int CeedBasisCreateProjection(CeedBasis basis_from, CeedBasis basis_to, CeedBasis *basis_project) {
1691f113e5dcSJeremy L Thompson   Ceed        ceed;
1692e104ad11SJames Wright   bool        create_tensor;
16931c66c397SJeremy L Thompson   CeedInt     dim, num_comp;
1694097cc795SJames Wright   CeedScalar *interp_project, *grad_project;
16951c66c397SJeremy L Thompson 
16962b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetCeed(basis_to, &ceed));
1697f113e5dcSJeremy L Thompson 
1698ecc88aebSJeremy L Thompson   // Create projection matrix
16992b730f8bSJeremy L Thompson   CeedCall(CeedBasisCreateProjectionMatrices(basis_from, basis_to, &interp_project, &grad_project));
1700f113e5dcSJeremy L Thompson 
1701f113e5dcSJeremy L Thompson   // Build basis
1702e104ad11SJames Wright   {
1703e104ad11SJames Wright     bool is_tensor_to, is_tensor_from;
1704e104ad11SJames Wright 
1705e104ad11SJames Wright     CeedCall(CeedBasisIsTensor(basis_to, &is_tensor_to));
1706e104ad11SJames Wright     CeedCall(CeedBasisIsTensor(basis_from, &is_tensor_from));
1707e104ad11SJames Wright     create_tensor = is_tensor_from && is_tensor_to;
1708e104ad11SJames Wright   }
17092b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetDimension(basis_to, &dim));
17102b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumComponents(basis_from, &num_comp));
1711e104ad11SJames Wright   if (create_tensor) {
1712f113e5dcSJeremy L Thompson     CeedInt P_1d_to, P_1d_from;
17131c66c397SJeremy L Thompson 
17142b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_1d_from));
17152b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_1d_to));
1716097cc795SJames Wright     CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P_1d_from, P_1d_to, interp_project, grad_project, NULL, NULL, basis_project));
1717f113e5dcSJeremy L Thompson   } else {
1718de05fbb2SSebastian Grimberg     // Even if basis_to and basis_from are not H1, the resulting basis is H1 for interpolation to work
1719f113e5dcSJeremy L Thompson     CeedInt          num_nodes_to, num_nodes_from;
17201c66c397SJeremy L Thompson     CeedElemTopology topo;
17211c66c397SJeremy L Thompson 
1722e00f3be8SJames Wright     CeedCall(CeedBasisGetTopology(basis_from, &topo));
17232b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis_from, &num_nodes_from));
17242b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis_to, &num_nodes_to));
1725097cc795SJames Wright     CeedCall(CeedBasisCreateH1(ceed, topo, num_comp, num_nodes_from, num_nodes_to, interp_project, grad_project, NULL, NULL, basis_project));
1726f113e5dcSJeremy L Thompson   }
1727f113e5dcSJeremy L Thompson 
1728f113e5dcSJeremy L Thompson   // Cleanup
17292b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_project));
17302b730f8bSJeremy L Thompson   CeedCall(CeedFree(&grad_project));
17319bc66399SJeremy L Thompson   CeedCall(CeedDestroy(&ceed));
1732f113e5dcSJeremy L Thompson   return CEED_ERROR_SUCCESS;
1733f113e5dcSJeremy L Thompson }
1734f113e5dcSJeremy L Thompson 
1735f113e5dcSJeremy L Thompson /**
1736ca94c3ddSJeremy L Thompson   @brief Copy the pointer to a `CeedBasis`.
17379560d06aSjeremylt 
1738ca94c3ddSJeremy 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`.
1739ca94c3ddSJeremy L Thompson         This `CeedBasis` will be destroyed if `*basis_copy` is the only reference to this `CeedBasis`.
1740ea61e9acSJeremy L Thompson 
1741ca94c3ddSJeremy L Thompson   @param[in]     basis      `CeedBasis` to copy reference to
1742ea61e9acSJeremy L Thompson   @param[in,out] basis_copy Variable to store copied reference
17439560d06aSjeremylt 
17449560d06aSjeremylt   @return An error code: 0 - success, otherwise - failure
17459560d06aSjeremylt 
17469560d06aSjeremylt   @ref User
17479560d06aSjeremylt **/
17489560d06aSjeremylt int CeedBasisReferenceCopy(CeedBasis basis, CeedBasis *basis_copy) {
1749356036faSJeremy L Thompson   if (basis != CEED_BASIS_NONE) CeedCall(CeedBasisReference(basis));
17502b730f8bSJeremy L Thompson   CeedCall(CeedBasisDestroy(basis_copy));
17519560d06aSjeremylt   *basis_copy = basis;
17529560d06aSjeremylt   return CEED_ERROR_SUCCESS;
17539560d06aSjeremylt }
17549560d06aSjeremylt 
17559560d06aSjeremylt /**
1756ca94c3ddSJeremy L Thompson   @brief View a `CeedBasis`
17577a982d89SJeremy L. Thompson 
1758ca94c3ddSJeremy L Thompson   @param[in] basis  `CeedBasis` to view
1759ca94c3ddSJeremy L Thompson   @param[in] stream Stream to view to, e.g., `stdout`
17607a982d89SJeremy L. Thompson 
17617a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
17627a982d89SJeremy L. Thompson 
17637a982d89SJeremy L. Thompson   @ref User
17647a982d89SJeremy L. Thompson **/
17657a982d89SJeremy L. Thompson int CeedBasisView(CeedBasis basis, FILE *stream) {
17661203703bSJeremy L Thompson   bool             is_tensor_basis;
17671203703bSJeremy L Thompson   CeedElemTopology topo;
17681203703bSJeremy L Thompson   CeedFESpace      fe_space;
17691203703bSJeremy L Thompson 
17701203703bSJeremy L Thompson   // Basis data
17711203703bSJeremy L Thompson   CeedCall(CeedBasisIsTensor(basis, &is_tensor_basis));
17721203703bSJeremy L Thompson   CeedCall(CeedBasisGetTopology(basis, &topo));
17731203703bSJeremy L Thompson   CeedCall(CeedBasisGetFESpace(basis, &fe_space));
17742b730f8bSJeremy L Thompson 
177550c301a5SRezgar Shakeri   // Print FE space and element topology of the basis
1776edf04919SJeremy L Thompson   fprintf(stream, "CeedBasis in a %s on a %s element\n", CeedFESpaces[fe_space], CeedElemTopologies[topo]);
17771203703bSJeremy L Thompson   if (is_tensor_basis) {
1778edf04919SJeremy L Thompson     fprintf(stream, "  P: %" CeedInt_FMT "\n  Q: %" CeedInt_FMT "\n", basis->P_1d, basis->Q_1d);
177950c301a5SRezgar Shakeri   } else {
1780edf04919SJeremy L Thompson     fprintf(stream, "  P: %" CeedInt_FMT "\n  Q: %" CeedInt_FMT "\n", basis->P, basis->Q);
178150c301a5SRezgar Shakeri   }
1782edf04919SJeremy L Thompson   fprintf(stream, "  dimension: %" CeedInt_FMT "\n  field components: %" CeedInt_FMT "\n", basis->dim, basis->num_comp);
1783ea61e9acSJeremy L Thompson   // Print quadrature data, interpolation/gradient/divergence/curl of the basis
17841203703bSJeremy L Thompson   if (is_tensor_basis) {  // tensor basis
17851203703bSJeremy L Thompson     CeedInt           P_1d, Q_1d;
17861203703bSJeremy L Thompson     const CeedScalar *q_ref_1d, *q_weight_1d, *interp_1d, *grad_1d;
17871203703bSJeremy L Thompson 
17881203703bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d));
17891203703bSJeremy L Thompson     CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
17901203703bSJeremy L Thompson     CeedCall(CeedBasisGetQRef(basis, &q_ref_1d));
17911203703bSJeremy L Thompson     CeedCall(CeedBasisGetQWeights(basis, &q_weight_1d));
17921203703bSJeremy L Thompson     CeedCall(CeedBasisGetInterp1D(basis, &interp_1d));
17931203703bSJeremy L Thompson     CeedCall(CeedBasisGetGrad1D(basis, &grad_1d));
17941203703bSJeremy L Thompson 
17951203703bSJeremy L Thompson     CeedCall(CeedScalarView("qref1d", "\t% 12.8f", 1, Q_1d, q_ref_1d, stream));
17961203703bSJeremy L Thompson     CeedCall(CeedScalarView("qweight1d", "\t% 12.8f", 1, Q_1d, q_weight_1d, stream));
17971203703bSJeremy L Thompson     CeedCall(CeedScalarView("interp1d", "\t% 12.8f", Q_1d, P_1d, interp_1d, stream));
17981203703bSJeremy L Thompson     CeedCall(CeedScalarView("grad1d", "\t% 12.8f", Q_1d, P_1d, grad_1d, stream));
179950c301a5SRezgar Shakeri   } else {  // non-tensor basis
18001203703bSJeremy L Thompson     CeedInt           P, Q, dim, q_comp;
18011203703bSJeremy L Thompson     const CeedScalar *q_ref, *q_weight, *interp, *grad, *div, *curl;
18021203703bSJeremy L Thompson 
18031203703bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis, &P));
18041203703bSJeremy L Thompson     CeedCall(CeedBasisGetNumQuadraturePoints(basis, &Q));
18051203703bSJeremy L Thompson     CeedCall(CeedBasisGetDimension(basis, &dim));
18061203703bSJeremy L Thompson     CeedCall(CeedBasisGetQRef(basis, &q_ref));
18071203703bSJeremy L Thompson     CeedCall(CeedBasisGetQWeights(basis, &q_weight));
18081203703bSJeremy L Thompson     CeedCall(CeedBasisGetInterp(basis, &interp));
18091203703bSJeremy L Thompson     CeedCall(CeedBasisGetGrad(basis, &grad));
18101203703bSJeremy L Thompson     CeedCall(CeedBasisGetDiv(basis, &div));
18111203703bSJeremy L Thompson     CeedCall(CeedBasisGetCurl(basis, &curl));
18121203703bSJeremy L Thompson 
18131203703bSJeremy L Thompson     CeedCall(CeedScalarView("qref", "\t% 12.8f", 1, Q * dim, q_ref, stream));
18141203703bSJeremy L Thompson     CeedCall(CeedScalarView("qweight", "\t% 12.8f", 1, Q, q_weight, stream));
1815c4e3f59bSSebastian Grimberg     CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp));
18161203703bSJeremy L Thompson     CeedCall(CeedScalarView("interp", "\t% 12.8f", q_comp * Q, P, interp, stream));
18171203703bSJeremy L Thompson     if (grad) {
1818c4e3f59bSSebastian Grimberg       CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp));
18191203703bSJeremy L Thompson       CeedCall(CeedScalarView("grad", "\t% 12.8f", q_comp * Q, P, grad, stream));
18207a982d89SJeremy L. Thompson     }
18211203703bSJeremy L Thompson     if (div) {
1822c4e3f59bSSebastian Grimberg       CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_DIV, &q_comp));
18231203703bSJeremy L Thompson       CeedCall(CeedScalarView("div", "\t% 12.8f", q_comp * Q, P, div, stream));
1824c4e3f59bSSebastian Grimberg     }
18251203703bSJeremy L Thompson     if (curl) {
1826c4e3f59bSSebastian Grimberg       CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_CURL, &q_comp));
18271203703bSJeremy L Thompson       CeedCall(CeedScalarView("curl", "\t% 12.8f", q_comp * Q, P, curl, stream));
182850c301a5SRezgar Shakeri     }
182950c301a5SRezgar Shakeri   }
1830e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
18317a982d89SJeremy L. Thompson }
18327a982d89SJeremy L. Thompson 
18337a982d89SJeremy L. Thompson /**
1834db2becc9SJeremy L Thompson   @brief Check input vector dimensions for CeedBasisApply[Add]
18357a982d89SJeremy L. Thompson 
1836ca94c3ddSJeremy L Thompson   @param[in]  basis     `CeedBasis` to evaluate
1837ea61e9acSJeremy L Thompson   @param[in]  num_elem  The number of elements to apply the basis evaluation to;
1838ca94c3ddSJeremy L Thompson                           the backend will specify the ordering in @ref CeedElemRestrictionCreate()
1839ca94c3ddSJeremy L Thompson   @param[in]  t_mode    @ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature points;
1840ca94c3ddSJeremy L Thompson                           @ref CEED_TRANSPOSE to apply the transpose, mapping from quadrature points to nodes
1841ca94c3ddSJeremy L Thompson   @param[in]  eval_mode @ref CEED_EVAL_NONE to use values directly,
1842ca94c3ddSJeremy L Thompson                           @ref CEED_EVAL_INTERP to use interpolated values,
1843ca94c3ddSJeremy L Thompson                           @ref CEED_EVAL_GRAD to use gradients,
1844ca94c3ddSJeremy L Thompson                           @ref CEED_EVAL_DIV to use divergence,
1845ca94c3ddSJeremy L Thompson                           @ref CEED_EVAL_CURL to use curl,
1846ca94c3ddSJeremy L Thompson                           @ref CEED_EVAL_WEIGHT to use quadrature weights
1847ca94c3ddSJeremy L Thompson   @param[in]  u         Input `CeedVector`
1848ca94c3ddSJeremy L Thompson   @param[out] v         Output `CeedVector`
18497a982d89SJeremy L. Thompson 
18507a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
18517a982d89SJeremy L. Thompson 
1852db2becc9SJeremy L Thompson   @ref Developer
18537a982d89SJeremy L. Thompson **/
1854db2becc9SJeremy L Thompson static int CeedBasisApplyCheckDims(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
1855c4e3f59bSSebastian Grimberg   CeedInt  dim, num_comp, q_comp, num_nodes, num_qpts;
18561c66c397SJeremy L Thompson   CeedSize u_length = 0, v_length;
18571c66c397SJeremy L Thompson 
18582b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetDimension(basis, &dim));
18592b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
1860c4e3f59bSSebastian Grimberg   CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp));
18612b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumNodes(basis, &num_nodes));
18622b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts));
18632b730f8bSJeremy L Thompson   CeedCall(CeedVectorGetLength(v, &v_length));
1864c8c3fa7dSJeremy L Thompson   if (u) CeedCall(CeedVectorGetLength(u, &u_length));
18657a982d89SJeremy L. Thompson 
1866e15f9bd0SJeremy L Thompson   // Check vector lengths to prevent out of bounds issues
186799e754f0SJeremy L Thompson   bool has_good_dims = true;
1868d1d35e2fSjeremylt   switch (eval_mode) {
1869e15f9bd0SJeremy L Thompson     case CEED_EVAL_NONE:
18702b730f8bSJeremy L Thompson     case CEED_EVAL_INTERP:
18712b730f8bSJeremy L Thompson     case CEED_EVAL_GRAD:
1872c4e3f59bSSebastian Grimberg     case CEED_EVAL_DIV:
1873c4e3f59bSSebastian Grimberg     case CEED_EVAL_CURL:
187419a04db8SJeremy L Thompson       has_good_dims = ((t_mode == CEED_TRANSPOSE && u_length >= (CeedSize)num_elem * (CeedSize)num_comp * (CeedSize)num_qpts * (CeedSize)q_comp &&
187519a04db8SJeremy L Thompson                         v_length >= (CeedSize)num_elem * (CeedSize)num_comp * (CeedSize)num_nodes) ||
187619a04db8SJeremy L Thompson                        (t_mode == CEED_NOTRANSPOSE && v_length >= (CeedSize)num_elem * (CeedSize)num_qpts * (CeedSize)num_comp * (CeedSize)q_comp &&
187719a04db8SJeremy L Thompson                         u_length >= (CeedSize)num_elem * (CeedSize)num_comp * (CeedSize)num_nodes));
1878e15f9bd0SJeremy L Thompson       break;
1879e15f9bd0SJeremy L Thompson     case CEED_EVAL_WEIGHT:
188019a04db8SJeremy L Thompson       has_good_dims = v_length >= (CeedSize)num_elem * (CeedSize)num_qpts;
1881e15f9bd0SJeremy L Thompson       break;
1882e15f9bd0SJeremy L Thompson   }
18839bc66399SJeremy L Thompson   CeedCheck(has_good_dims, CeedBasisReturnCeed(basis), CEED_ERROR_DIMENSION, "Input/output vectors too short for basis and evaluation mode");
1884db2becc9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1885db2becc9SJeremy L Thompson }
1886e15f9bd0SJeremy L Thompson 
1887db2becc9SJeremy L Thompson /**
1888db2becc9SJeremy L Thompson   @brief Apply basis evaluation from nodes to quadrature points or vice versa
1889db2becc9SJeremy L Thompson 
1890db2becc9SJeremy L Thompson   @param[in]  basis     `CeedBasis` to evaluate
1891db2becc9SJeremy L Thompson   @param[in]  num_elem  The number of elements to apply the basis evaluation to;
1892db2becc9SJeremy L Thompson                           the backend will specify the ordering in @ref CeedElemRestrictionCreate()
1893db2becc9SJeremy L Thompson   @param[in]  t_mode    @ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature points;
1894db2becc9SJeremy L Thompson                           @ref CEED_TRANSPOSE to apply the transpose, mapping from quadrature points to nodes
1895db2becc9SJeremy L Thompson   @param[in]  eval_mode @ref CEED_EVAL_NONE to use values directly,
1896db2becc9SJeremy L Thompson                           @ref CEED_EVAL_INTERP to use interpolated values,
1897db2becc9SJeremy L Thompson                           @ref CEED_EVAL_GRAD to use gradients,
1898db2becc9SJeremy L Thompson                           @ref CEED_EVAL_DIV to use divergence,
1899db2becc9SJeremy L Thompson                           @ref CEED_EVAL_CURL to use curl,
1900db2becc9SJeremy L Thompson                           @ref CEED_EVAL_WEIGHT to use quadrature weights
1901db2becc9SJeremy L Thompson   @param[in]  u         Input `CeedVector`
1902db2becc9SJeremy L Thompson   @param[out] v         Output `CeedVector`
1903db2becc9SJeremy L Thompson 
1904db2becc9SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1905db2becc9SJeremy L Thompson 
1906db2becc9SJeremy L Thompson   @ref User
1907db2becc9SJeremy L Thompson **/
1908db2becc9SJeremy L Thompson int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
1909db2becc9SJeremy L Thompson   CeedCall(CeedBasisApplyCheckDims(basis, num_elem, t_mode, eval_mode, u, v));
1910db2becc9SJeremy L Thompson   CeedCheck(basis->Apply, CeedBasisReturnCeed(basis), CEED_ERROR_UNSUPPORTED, "Backend does not support CeedBasisApply");
19112b730f8bSJeremy L Thompson   CeedCall(basis->Apply(basis, num_elem, t_mode, eval_mode, u, v));
1912e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
19137a982d89SJeremy L. Thompson }
19147a982d89SJeremy L. Thompson 
19157a982d89SJeremy L. Thompson /**
1916db2becc9SJeremy L Thompson   @brief Apply basis evaluation from quadrature points to nodes and sum into target vector
1917db2becc9SJeremy L Thompson 
1918db2becc9SJeremy L Thompson   @param[in]  basis     `CeedBasis` to evaluate
1919db2becc9SJeremy L Thompson   @param[in]  num_elem  The number of elements to apply the basis evaluation to;
1920db2becc9SJeremy L Thompson                           the backend will specify the ordering in @ref CeedElemRestrictionCreate()
1921db2becc9SJeremy L Thompson   @param[in]  t_mode    @ref CEED_TRANSPOSE to apply the transpose, mapping from quadrature points to nodes;
1922db2becc9SJeremy L Thompson                            @ref CEED_NOTRANSPOSE is not valid for `CeedBasisApplyAdd()`
1923db2becc9SJeremy L Thompson   @param[in]  eval_mode @ref CEED_EVAL_NONE to use values directly,
1924db2becc9SJeremy L Thompson                           @ref CEED_EVAL_INTERP to use interpolated values,
1925db2becc9SJeremy L Thompson                           @ref CEED_EVAL_GRAD to use gradients,
1926db2becc9SJeremy L Thompson                           @ref CEED_EVAL_DIV to use divergence,
1927db2becc9SJeremy L Thompson                           @ref CEED_EVAL_CURL to use curl,
1928db2becc9SJeremy L Thompson                           @ref CEED_EVAL_WEIGHT to use quadrature weights
1929db2becc9SJeremy L Thompson   @param[in]  u         Input `CeedVector`
1930db2becc9SJeremy L Thompson   @param[out] v         Output `CeedVector` to sum into
1931db2becc9SJeremy L Thompson 
1932db2becc9SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1933db2becc9SJeremy L Thompson 
1934db2becc9SJeremy L Thompson   @ref User
1935db2becc9SJeremy L Thompson **/
1936db2becc9SJeremy L Thompson int CeedBasisApplyAdd(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
1937db2becc9SJeremy L Thompson   CeedCheck(t_mode == CEED_TRANSPOSE, CeedBasisReturnCeed(basis), CEED_ERROR_UNSUPPORTED, "CeedBasisApplyAdd only supports CEED_TRANSPOSE");
1938db2becc9SJeremy L Thompson   CeedCall(CeedBasisApplyCheckDims(basis, num_elem, t_mode, eval_mode, u, v));
1939db2becc9SJeremy L Thompson   CeedCheck(basis->ApplyAdd, CeedBasisReturnCeed(basis), CEED_ERROR_UNSUPPORTED, "Backend does not implement CeedBasisApplyAdd");
1940db2becc9SJeremy L Thompson   CeedCall(basis->ApplyAdd(basis, num_elem, t_mode, eval_mode, u, v));
1941db2becc9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1942db2becc9SJeremy L Thompson }
1943db2becc9SJeremy L Thompson 
1944db2becc9SJeremy L Thompson /**
1945db2becc9SJeremy L Thompson   @brief Apply basis evaluation from nodes to arbitrary points
1946db2becc9SJeremy L Thompson 
1947db2becc9SJeremy L Thompson   @param[in]  basis      `CeedBasis` to evaluate
1948db2becc9SJeremy L Thompson   @param[in]  num_elem   The number of elements to apply the basis evaluation to;
1949db2becc9SJeremy L Thompson                           the backend will specify the ordering in @ref CeedElemRestrictionCreate()
1950db2becc9SJeremy L Thompson   @param[in]  num_points Array of the number of points to apply the basis evaluation to in each element, size `num_elem`
1951db2becc9SJeremy L Thompson   @param[in]  t_mode     @ref CEED_NOTRANSPOSE to evaluate from nodes to points;
1952db2becc9SJeremy L Thompson                            @ref CEED_TRANSPOSE to apply the transpose, mapping from points to nodes
1953db2becc9SJeremy L Thompson   @param[in]  eval_mode  @ref CEED_EVAL_INTERP to use interpolated values,
1954db2becc9SJeremy L Thompson                            @ref CEED_EVAL_GRAD to use gradients,
1955db2becc9SJeremy L Thompson                            @ref CEED_EVAL_WEIGHT to use quadrature weights
1956db2becc9SJeremy L Thompson   @param[in]  x_ref      `CeedVector` holding reference coordinates of each point
1957db2becc9SJeremy L Thompson   @param[in]  u          Input `CeedVector`, of length `num_nodes * num_comp` for @ref CEED_NOTRANSPOSE
1958db2becc9SJeremy L Thompson   @param[out] v          Output `CeedVector`, of length `num_points * num_q_comp` for @ref CEED_NOTRANSPOSE with @ref CEED_EVAL_INTERP
1959db2becc9SJeremy L Thompson 
1960db2becc9SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1961db2becc9SJeremy L Thompson 
1962db2becc9SJeremy L Thompson   @ref User
1963db2becc9SJeremy L Thompson **/
1964db2becc9SJeremy L Thompson int CeedBasisApplyAtPoints(CeedBasis basis, CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode, CeedEvalMode eval_mode,
1965db2becc9SJeremy L Thompson                            CeedVector x_ref, CeedVector u, CeedVector v) {
1966db2becc9SJeremy L Thompson   CeedCall(CeedBasisApplyAtPointsCheckDims(basis, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
1967db2becc9SJeremy L Thompson   if (basis->ApplyAtPoints) {
1968db2becc9SJeremy L Thompson     CeedCall(basis->ApplyAtPoints(basis, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
1969db2becc9SJeremy L Thompson   } else {
1970db2becc9SJeremy L Thompson     CeedCall(CeedBasisApplyAtPoints_Core(basis, false, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
1971db2becc9SJeremy L Thompson   }
1972db2becc9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1973db2becc9SJeremy L Thompson }
1974db2becc9SJeremy L Thompson 
1975db2becc9SJeremy L Thompson /**
1976db2becc9SJeremy L Thompson   @brief Apply basis evaluation from nodes to arbitrary points and sum into target vector
1977db2becc9SJeremy L Thompson 
1978db2becc9SJeremy L Thompson   @param[in]  basis      `CeedBasis` to evaluate
1979db2becc9SJeremy L Thompson   @param[in]  num_elem   The number of elements to apply the basis evaluation to;
1980db2becc9SJeremy L Thompson                           the backend will specify the ordering in @ref CeedElemRestrictionCreate()
1981db2becc9SJeremy L Thompson   @param[in]  num_points Array of the number of points to apply the basis evaluation to in each element, size `num_elem`
1982db2becc9SJeremy L Thompson   @param[in]  t_mode     @ref CEED_NOTRANSPOSE to evaluate from nodes to points;
1983db2becc9SJeremy L Thompson                            @ref CEED_NOTRANSPOSE is not valid for `CeedBasisApplyAddAtPoints()`
1984db2becc9SJeremy L Thompson   @param[in]  eval_mode  @ref CEED_EVAL_INTERP to use interpolated values,
1985db2becc9SJeremy L Thompson                            @ref CEED_EVAL_GRAD to use gradients,
1986db2becc9SJeremy L Thompson                            @ref CEED_EVAL_WEIGHT to use quadrature weights
1987db2becc9SJeremy L Thompson   @param[in]  x_ref      `CeedVector` holding reference coordinates of each point
1988db2becc9SJeremy L Thompson   @param[in]  u          Input `CeedVector`, of length `num_nodes * num_comp` for @ref CEED_NOTRANSPOSE
1989db2becc9SJeremy L Thompson   @param[out] v          Output `CeedVector`, of length `num_points * num_q_comp` for @ref CEED_NOTRANSPOSE with @ref CEED_EVAL_INTERP
1990db2becc9SJeremy L Thompson 
1991db2becc9SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1992db2becc9SJeremy L Thompson 
1993db2becc9SJeremy L Thompson   @ref User
1994db2becc9SJeremy L Thompson **/
1995db2becc9SJeremy L Thompson int CeedBasisApplyAddAtPoints(CeedBasis basis, CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode, CeedEvalMode eval_mode,
1996db2becc9SJeremy L Thompson                               CeedVector x_ref, CeedVector u, CeedVector v) {
1997db2becc9SJeremy L Thompson   CeedCheck(t_mode == CEED_TRANSPOSE, CeedBasisReturnCeed(basis), CEED_ERROR_UNSUPPORTED, "CeedBasisApplyAddAtPoints only supports CEED_TRANSPOSE");
1998db2becc9SJeremy L Thompson   CeedCall(CeedBasisApplyAtPointsCheckDims(basis, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
1999db2becc9SJeremy L Thompson   if (basis->ApplyAddAtPoints) {
2000db2becc9SJeremy L Thompson     CeedCall(basis->ApplyAddAtPoints(basis, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
2001db2becc9SJeremy L Thompson   } else {
2002db2becc9SJeremy L Thompson     CeedCall(CeedBasisApplyAtPoints_Core(basis, true, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
2003db2becc9SJeremy L Thompson   }
2004db2becc9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2005db2becc9SJeremy L Thompson }
2006db2becc9SJeremy L Thompson 
2007db2becc9SJeremy L Thompson /**
20086e536b99SJeremy L Thompson   @brief Get the `Ceed` associated with a `CeedBasis`
2009b7c9bbdaSJeremy L Thompson 
2010ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
2011ca94c3ddSJeremy L Thompson   @param[out] ceed  Variable to store `Ceed`
2012b7c9bbdaSJeremy L Thompson 
2013b7c9bbdaSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
2014b7c9bbdaSJeremy L Thompson 
2015b7c9bbdaSJeremy L Thompson   @ref Advanced
2016b7c9bbdaSJeremy L Thompson **/
2017b7c9bbdaSJeremy L Thompson int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) {
20189bc66399SJeremy L Thompson   *ceed = NULL;
20199bc66399SJeremy L Thompson   CeedCall(CeedReferenceCopy(CeedBasisReturnCeed(basis), ceed));
2020b7c9bbdaSJeremy L Thompson   return CEED_ERROR_SUCCESS;
2021b7c9bbdaSJeremy L Thompson }
2022b7c9bbdaSJeremy L Thompson 
2023b7c9bbdaSJeremy L Thompson /**
20246e536b99SJeremy L Thompson   @brief Return the `Ceed` associated with a `CeedBasis`
20256e536b99SJeremy L Thompson 
20266e536b99SJeremy L Thompson   @param[in]  basis `CeedBasis`
20276e536b99SJeremy L Thompson 
20286e536b99SJeremy L Thompson   @return `Ceed` associated with the `basis`
20296e536b99SJeremy L Thompson 
20306e536b99SJeremy L Thompson   @ref Advanced
20316e536b99SJeremy L Thompson **/
20326e536b99SJeremy L Thompson Ceed CeedBasisReturnCeed(CeedBasis basis) { return basis->ceed; }
20336e536b99SJeremy L Thompson 
20346e536b99SJeremy L Thompson /**
2035ca94c3ddSJeremy L Thompson   @brief Get dimension for given `CeedBasis`
20369d007619Sjeremylt 
2037ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
20389d007619Sjeremylt   @param[out] dim   Variable to store dimension of basis
20399d007619Sjeremylt 
20409d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
20419d007619Sjeremylt 
2042b7c9bbdaSJeremy L Thompson   @ref Advanced
20439d007619Sjeremylt **/
20449d007619Sjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) {
20459d007619Sjeremylt   *dim = basis->dim;
2046e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
20479d007619Sjeremylt }
20489d007619Sjeremylt 
20499d007619Sjeremylt /**
2050ca94c3ddSJeremy L Thompson   @brief Get topology for given `CeedBasis`
2051d99fa3c5SJeremy L Thompson 
2052ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
2053d99fa3c5SJeremy L Thompson   @param[out] topo  Variable to store topology of basis
2054d99fa3c5SJeremy L Thompson 
2055d99fa3c5SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
2056d99fa3c5SJeremy L Thompson 
2057b7c9bbdaSJeremy L Thompson   @ref Advanced
2058d99fa3c5SJeremy L Thompson **/
2059d99fa3c5SJeremy L Thompson int CeedBasisGetTopology(CeedBasis basis, CeedElemTopology *topo) {
2060d99fa3c5SJeremy L Thompson   *topo = basis->topo;
2061e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2062d99fa3c5SJeremy L Thompson }
2063d99fa3c5SJeremy L Thompson 
2064d99fa3c5SJeremy L Thompson /**
2065ca94c3ddSJeremy L Thompson   @brief Get number of components for given `CeedBasis`
20669d007619Sjeremylt 
2067ca94c3ddSJeremy L Thompson   @param[in]  basis    `CeedBasis`
2068ca94c3ddSJeremy L Thompson   @param[out] num_comp Variable to store number of components
20699d007619Sjeremylt 
20709d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
20719d007619Sjeremylt 
2072b7c9bbdaSJeremy L Thompson   @ref Advanced
20739d007619Sjeremylt **/
2074d1d35e2fSjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *num_comp) {
2075d1d35e2fSjeremylt   *num_comp = basis->num_comp;
2076e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
20779d007619Sjeremylt }
20789d007619Sjeremylt 
20799d007619Sjeremylt /**
2080ca94c3ddSJeremy L Thompson   @brief Get total number of nodes (in `dim` dimensions) of a `CeedBasis`
20819d007619Sjeremylt 
2082ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
20839d007619Sjeremylt   @param[out] P     Variable to store number of nodes
20849d007619Sjeremylt 
20859d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
20869d007619Sjeremylt 
20879d007619Sjeremylt   @ref Utility
20889d007619Sjeremylt **/
20899d007619Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) {
20909d007619Sjeremylt   *P = basis->P;
2091e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
20929d007619Sjeremylt }
20939d007619Sjeremylt 
20949d007619Sjeremylt /**
2095ca94c3ddSJeremy L Thompson   @brief Get total number of nodes (in 1 dimension) of a `CeedBasis`
20969d007619Sjeremylt 
2097ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
2098d1d35e2fSjeremylt   @param[out] P_1d  Variable to store number of nodes
20999d007619Sjeremylt 
21009d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
21019d007619Sjeremylt 
2102b7c9bbdaSJeremy L Thompson   @ref Advanced
21039d007619Sjeremylt **/
2104d1d35e2fSjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P_1d) {
21056e536b99SJeremy L Thompson   CeedCheck(basis->is_tensor_basis, CeedBasisReturnCeed(basis), CEED_ERROR_MINOR, "Cannot supply P_1d for non-tensor CeedBasis");
2106d1d35e2fSjeremylt   *P_1d = basis->P_1d;
2107e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
21089d007619Sjeremylt }
21099d007619Sjeremylt 
21109d007619Sjeremylt /**
2111ca94c3ddSJeremy L Thompson   @brief Get total number of quadrature points (in `dim` dimensions) of a `CeedBasis`
21129d007619Sjeremylt 
2113ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
21149d007619Sjeremylt   @param[out] Q     Variable to store number of quadrature points
21159d007619Sjeremylt 
21169d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
21179d007619Sjeremylt 
21189d007619Sjeremylt   @ref Utility
21199d007619Sjeremylt **/
21209d007619Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) {
21219d007619Sjeremylt   *Q = basis->Q;
2122e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
21239d007619Sjeremylt }
21249d007619Sjeremylt 
21259d007619Sjeremylt /**
2126ca94c3ddSJeremy L Thompson   @brief Get total number of quadrature points (in 1 dimension) of a `CeedBasis`
21279d007619Sjeremylt 
2128ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
2129d1d35e2fSjeremylt   @param[out] Q_1d  Variable to store number of quadrature points
21309d007619Sjeremylt 
21319d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
21329d007619Sjeremylt 
2133b7c9bbdaSJeremy L Thompson   @ref Advanced
21349d007619Sjeremylt **/
2135d1d35e2fSjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q_1d) {
21366e536b99SJeremy L Thompson   CeedCheck(basis->is_tensor_basis, CeedBasisReturnCeed(basis), CEED_ERROR_MINOR, "Cannot supply Q_1d for non-tensor CeedBasis");
2137d1d35e2fSjeremylt   *Q_1d = basis->Q_1d;
2138e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
21399d007619Sjeremylt }
21409d007619Sjeremylt 
21419d007619Sjeremylt /**
2142ca94c3ddSJeremy L Thompson   @brief Get reference coordinates of quadrature points (in `dim` dimensions) of a `CeedBasis`
21439d007619Sjeremylt 
2144ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
2145d1d35e2fSjeremylt   @param[out] q_ref Variable to store reference coordinates of quadrature points
21469d007619Sjeremylt 
21479d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
21489d007619Sjeremylt 
2149b7c9bbdaSJeremy L Thompson   @ref Advanced
21509d007619Sjeremylt **/
2151d1d35e2fSjeremylt int CeedBasisGetQRef(CeedBasis basis, const CeedScalar **q_ref) {
2152d1d35e2fSjeremylt   *q_ref = basis->q_ref_1d;
2153e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
21549d007619Sjeremylt }
21559d007619Sjeremylt 
21569d007619Sjeremylt /**
2157ca94c3ddSJeremy L Thompson   @brief Get quadrature weights of quadrature points (in `dim` dimensions) of a `CeedBasis`
21589d007619Sjeremylt 
2159ca94c3ddSJeremy L Thompson   @param[in]  basis    `CeedBasis`
2160d1d35e2fSjeremylt   @param[out] q_weight Variable to store quadrature weights
21619d007619Sjeremylt 
21629d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
21639d007619Sjeremylt 
2164b7c9bbdaSJeremy L Thompson   @ref Advanced
21659d007619Sjeremylt **/
2166d1d35e2fSjeremylt int CeedBasisGetQWeights(CeedBasis basis, const CeedScalar **q_weight) {
2167d1d35e2fSjeremylt   *q_weight = basis->q_weight_1d;
2168e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
21699d007619Sjeremylt }
21709d007619Sjeremylt 
21719d007619Sjeremylt /**
2172ca94c3ddSJeremy L Thompson   @brief Get interpolation matrix of a `CeedBasis`
21739d007619Sjeremylt 
2174ca94c3ddSJeremy L Thompson   @param[in]  basis  `CeedBasis`
21759d007619Sjeremylt   @param[out] interp Variable to store interpolation matrix
21769d007619Sjeremylt 
21779d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
21789d007619Sjeremylt 
2179b7c9bbdaSJeremy L Thompson   @ref Advanced
21809d007619Sjeremylt **/
21816c58de82SJeremy L Thompson int CeedBasisGetInterp(CeedBasis basis, const CeedScalar **interp) {
21826402da51SJeremy L Thompson   if (!basis->interp && basis->is_tensor_basis) {
21839d007619Sjeremylt     // Allocate
21842b730f8bSJeremy L Thompson     CeedCall(CeedMalloc(basis->Q * basis->P, &basis->interp));
21859d007619Sjeremylt 
21869d007619Sjeremylt     // Initialize
21872b730f8bSJeremy L Thompson     for (CeedInt i = 0; i < basis->Q * basis->P; i++) basis->interp[i] = 1.0;
21889d007619Sjeremylt 
21899d007619Sjeremylt     // Calculate
21902b730f8bSJeremy L Thompson     for (CeedInt d = 0; d < basis->dim; d++) {
21912b730f8bSJeremy L Thompson       for (CeedInt qpt = 0; qpt < basis->Q; qpt++) {
21929d007619Sjeremylt         for (CeedInt node = 0; node < basis->P; node++) {
2193d1d35e2fSjeremylt           CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d;
2194d1d35e2fSjeremylt           CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d;
21951c66c397SJeremy L Thompson 
2196d1d35e2fSjeremylt           basis->interp[qpt * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p];
21979d007619Sjeremylt         }
21989d007619Sjeremylt       }
21992b730f8bSJeremy L Thompson     }
22002b730f8bSJeremy L Thompson   }
22019d007619Sjeremylt   *interp = basis->interp;
2202e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
22039d007619Sjeremylt }
22049d007619Sjeremylt 
22059d007619Sjeremylt /**
2206ca94c3ddSJeremy L Thompson   @brief Get 1D interpolation matrix of a tensor product `CeedBasis`
22079d007619Sjeremylt 
2208ca94c3ddSJeremy L Thompson   @param[in]  basis     `CeedBasis`
2209d1d35e2fSjeremylt   @param[out] interp_1d Variable to store interpolation matrix
22109d007619Sjeremylt 
22119d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
22129d007619Sjeremylt 
22139d007619Sjeremylt   @ref Backend
22149d007619Sjeremylt **/
2215d1d35e2fSjeremylt int CeedBasisGetInterp1D(CeedBasis basis, const CeedScalar **interp_1d) {
22161203703bSJeremy L Thompson   bool is_tensor_basis;
22171203703bSJeremy L Thompson 
22181203703bSJeremy L Thompson   CeedCall(CeedBasisIsTensor(basis, &is_tensor_basis));
22196e536b99SJeremy L Thompson   CeedCheck(is_tensor_basis, CeedBasisReturnCeed(basis), CEED_ERROR_MINOR, "CeedBasis is not a tensor product CeedBasis");
2220d1d35e2fSjeremylt   *interp_1d = basis->interp_1d;
2221e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
22229d007619Sjeremylt }
22239d007619Sjeremylt 
22249d007619Sjeremylt /**
2225ca94c3ddSJeremy L Thompson   @brief Get gradient matrix of a `CeedBasis`
22269d007619Sjeremylt 
2227ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
22289d007619Sjeremylt   @param[out] grad  Variable to store gradient matrix
22299d007619Sjeremylt 
22309d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
22319d007619Sjeremylt 
2232b7c9bbdaSJeremy L Thompson   @ref Advanced
22339d007619Sjeremylt **/
22346c58de82SJeremy L Thompson int CeedBasisGetGrad(CeedBasis basis, const CeedScalar **grad) {
22356402da51SJeremy L Thompson   if (!basis->grad && basis->is_tensor_basis) {
22369d007619Sjeremylt     // Allocate
22372b730f8bSJeremy L Thompson     CeedCall(CeedMalloc(basis->dim * basis->Q * basis->P, &basis->grad));
22389d007619Sjeremylt 
22399d007619Sjeremylt     // Initialize
22402b730f8bSJeremy L Thompson     for (CeedInt i = 0; i < basis->dim * basis->Q * basis->P; i++) basis->grad[i] = 1.0;
22419d007619Sjeremylt 
22429d007619Sjeremylt     // Calculate
22432b730f8bSJeremy L Thompson     for (CeedInt d = 0; d < basis->dim; d++) {
22442b730f8bSJeremy L Thompson       for (CeedInt i = 0; i < basis->dim; i++) {
22452b730f8bSJeremy L Thompson         for (CeedInt qpt = 0; qpt < basis->Q; qpt++) {
22469d007619Sjeremylt           for (CeedInt node = 0; node < basis->P; node++) {
2247d1d35e2fSjeremylt             CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d;
2248d1d35e2fSjeremylt             CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d;
22491c66c397SJeremy L Thompson 
22502b730f8bSJeremy L Thompson             if (i == d) basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->grad_1d[q * basis->P_1d + p];
22512b730f8bSJeremy L Thompson             else basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p];
22522b730f8bSJeremy L Thompson           }
22532b730f8bSJeremy L Thompson         }
22542b730f8bSJeremy L Thompson       }
22559d007619Sjeremylt     }
22569d007619Sjeremylt   }
22579d007619Sjeremylt   *grad = basis->grad;
2258e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
22599d007619Sjeremylt }
22609d007619Sjeremylt 
22619d007619Sjeremylt /**
2262ca94c3ddSJeremy L Thompson   @brief Get 1D gradient matrix of a tensor product `CeedBasis`
22639d007619Sjeremylt 
2264ca94c3ddSJeremy L Thompson   @param[in]  basis   `CeedBasis`
2265d1d35e2fSjeremylt   @param[out] grad_1d Variable to store gradient matrix
22669d007619Sjeremylt 
22679d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
22689d007619Sjeremylt 
2269b7c9bbdaSJeremy L Thompson   @ref Advanced
22709d007619Sjeremylt **/
2271d1d35e2fSjeremylt int CeedBasisGetGrad1D(CeedBasis basis, const CeedScalar **grad_1d) {
22721203703bSJeremy L Thompson   bool is_tensor_basis;
22731203703bSJeremy L Thompson 
22741203703bSJeremy L Thompson   CeedCall(CeedBasisIsTensor(basis, &is_tensor_basis));
22756e536b99SJeremy L Thompson   CeedCheck(is_tensor_basis, CeedBasisReturnCeed(basis), CEED_ERROR_MINOR, "CeedBasis is not a tensor product CeedBasis");
2276d1d35e2fSjeremylt   *grad_1d = basis->grad_1d;
2277e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
22789d007619Sjeremylt }
22799d007619Sjeremylt 
22809d007619Sjeremylt /**
2281ca94c3ddSJeremy L Thompson   @brief Get divergence matrix of a `CeedBasis`
228250c301a5SRezgar Shakeri 
2283ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
228450c301a5SRezgar Shakeri   @param[out] div   Variable to store divergence matrix
228550c301a5SRezgar Shakeri 
228650c301a5SRezgar Shakeri   @return An error code: 0 - success, otherwise - failure
228750c301a5SRezgar Shakeri 
228850c301a5SRezgar Shakeri   @ref Advanced
228950c301a5SRezgar Shakeri **/
229050c301a5SRezgar Shakeri int CeedBasisGetDiv(CeedBasis basis, const CeedScalar **div) {
229150c301a5SRezgar Shakeri   *div = basis->div;
229250c301a5SRezgar Shakeri   return CEED_ERROR_SUCCESS;
229350c301a5SRezgar Shakeri }
229450c301a5SRezgar Shakeri 
229550c301a5SRezgar Shakeri /**
2296ca94c3ddSJeremy L Thompson   @brief Get curl matrix of a `CeedBasis`
2297c4e3f59bSSebastian Grimberg 
2298ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
2299c4e3f59bSSebastian Grimberg   @param[out] curl  Variable to store curl matrix
2300c4e3f59bSSebastian Grimberg 
2301c4e3f59bSSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
2302c4e3f59bSSebastian Grimberg 
2303c4e3f59bSSebastian Grimberg   @ref Advanced
2304c4e3f59bSSebastian Grimberg **/
2305c4e3f59bSSebastian Grimberg int CeedBasisGetCurl(CeedBasis basis, const CeedScalar **curl) {
2306c4e3f59bSSebastian Grimberg   *curl = basis->curl;
2307c4e3f59bSSebastian Grimberg   return CEED_ERROR_SUCCESS;
2308c4e3f59bSSebastian Grimberg }
2309c4e3f59bSSebastian Grimberg 
2310c4e3f59bSSebastian Grimberg /**
2311ca94c3ddSJeremy L Thompson   @brief Destroy a @ref  CeedBasis
23127a982d89SJeremy L. Thompson 
2313ca94c3ddSJeremy L Thompson   @param[in,out] basis `CeedBasis` to destroy
23147a982d89SJeremy L. Thompson 
23157a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
23167a982d89SJeremy L. Thompson 
23177a982d89SJeremy L. Thompson   @ref User
23187a982d89SJeremy L. Thompson **/
23197a982d89SJeremy L. Thompson int CeedBasisDestroy(CeedBasis *basis) {
2320356036faSJeremy L Thompson   if (!*basis || *basis == CEED_BASIS_NONE || --(*basis)->ref_count > 0) {
2321ad6481ceSJeremy L Thompson     *basis = NULL;
2322ad6481ceSJeremy L Thompson     return CEED_ERROR_SUCCESS;
2323ad6481ceSJeremy L Thompson   }
23242b730f8bSJeremy L Thompson   if ((*basis)->Destroy) CeedCall((*basis)->Destroy(*basis));
23259831d45aSJeremy L Thompson   CeedCall(CeedTensorContractDestroy(&(*basis)->contract));
2326c4e3f59bSSebastian Grimberg   CeedCall(CeedFree(&(*basis)->q_ref_1d));
2327c4e3f59bSSebastian Grimberg   CeedCall(CeedFree(&(*basis)->q_weight_1d));
23282b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->interp));
23292b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->interp_1d));
23302b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->grad));
23312b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->grad_1d));
2332c4e3f59bSSebastian Grimberg   CeedCall(CeedFree(&(*basis)->div));
2333c4e3f59bSSebastian Grimberg   CeedCall(CeedFree(&(*basis)->curl));
2334c8c3fa7dSJeremy L Thompson   CeedCall(CeedVectorDestroy(&(*basis)->vec_chebyshev));
2335c8c3fa7dSJeremy L Thompson   CeedCall(CeedBasisDestroy(&(*basis)->basis_chebyshev));
23362b730f8bSJeremy L Thompson   CeedCall(CeedDestroy(&(*basis)->ceed));
23372b730f8bSJeremy L Thompson   CeedCall(CeedFree(basis));
2338e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
23397a982d89SJeremy L. Thompson }
23407a982d89SJeremy L. Thompson 
23417a982d89SJeremy L. Thompson /**
2342b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre quadrature
2343b11c1e72Sjeremylt 
2344ca94c3ddSJeremy L Thompson   @param[in]  Q           Number of quadrature points (integrates polynomials of degree `2*Q-1` exactly)
2345ca94c3ddSJeremy L Thompson   @param[out] q_ref_1d    Array of length `Q` to hold the abscissa on `[-1, 1]`
2346ca94c3ddSJeremy L Thompson   @param[out] q_weight_1d Array of length `Q` to hold the weights
2347b11c1e72Sjeremylt 
2348b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
2349dfdf5a53Sjeremylt 
2350dfdf5a53Sjeremylt   @ref Utility
2351b11c1e72Sjeremylt **/
23522b730f8bSJeremy L Thompson int CeedGaussQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) {
2353d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0 * atan(1.0);
23541c66c397SJeremy L Thompson 
2355d1d35e2fSjeremylt   // Build q_ref_1d, q_weight_1d
235692ae7e47SJeremy L Thompson   for (CeedInt i = 0; i <= Q / 2; i++) {
2357d7b241e6Sjeremylt     // Guess
2358d7b241e6Sjeremylt     xi = cos(PI * (CeedScalar)(2 * i + 1) / ((CeedScalar)(2 * Q)));
2359d7b241e6Sjeremylt     // Pn(xi)
2360d7b241e6Sjeremylt     P0 = 1.0;
2361d7b241e6Sjeremylt     P1 = xi;
2362d7b241e6Sjeremylt     P2 = 0.0;
236392ae7e47SJeremy L Thompson     for (CeedInt j = 2; j <= Q; j++) {
2364d7b241e6Sjeremylt       P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
2365d7b241e6Sjeremylt       P0 = P1;
2366d7b241e6Sjeremylt       P1 = P2;
2367d7b241e6Sjeremylt     }
2368d7b241e6Sjeremylt     // First Newton Step
2369d7b241e6Sjeremylt     dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
2370d7b241e6Sjeremylt     xi  = xi - P2 / dP2;
2371d7b241e6Sjeremylt     // Newton to convergence
237292ae7e47SJeremy L Thompson     for (CeedInt k = 0; k < 100 && fabs(P2) > 10 * CEED_EPSILON; k++) {
2373d7b241e6Sjeremylt       P0 = 1.0;
2374d7b241e6Sjeremylt       P1 = xi;
237592ae7e47SJeremy L Thompson       for (CeedInt j = 2; j <= Q; j++) {
2376d7b241e6Sjeremylt         P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
2377d7b241e6Sjeremylt         P0 = P1;
2378d7b241e6Sjeremylt         P1 = P2;
2379d7b241e6Sjeremylt       }
2380d7b241e6Sjeremylt       dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
2381d7b241e6Sjeremylt       xi  = xi - P2 / dP2;
2382d7b241e6Sjeremylt     }
2383d7b241e6Sjeremylt     // Save xi, wi
2384d7b241e6Sjeremylt     wi                     = 2.0 / ((1.0 - xi * xi) * dP2 * dP2);
2385d1d35e2fSjeremylt     q_weight_1d[i]         = wi;
2386d1d35e2fSjeremylt     q_weight_1d[Q - 1 - i] = wi;
2387d1d35e2fSjeremylt     q_ref_1d[i]            = -xi;
2388d1d35e2fSjeremylt     q_ref_1d[Q - 1 - i]    = xi;
2389d7b241e6Sjeremylt   }
2390e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2391d7b241e6Sjeremylt }
2392d7b241e6Sjeremylt 
2393b11c1e72Sjeremylt /**
2394b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre-Lobatto quadrature
2395b11c1e72Sjeremylt 
2396ca94c3ddSJeremy L Thompson   @param[in]  Q           Number of quadrature points (integrates polynomials of degree `2*Q-3` exactly)
2397ca94c3ddSJeremy L Thompson   @param[out] q_ref_1d    Array of length `Q` to hold the abscissa on `[-1, 1]`
2398ca94c3ddSJeremy L Thompson   @param[out] q_weight_1d Array of length `Q` to hold the weights
2399b11c1e72Sjeremylt 
2400b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
2401dfdf5a53Sjeremylt 
2402dfdf5a53Sjeremylt   @ref Utility
2403b11c1e72Sjeremylt **/
24042b730f8bSJeremy L Thompson int CeedLobattoQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) {
2405d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0 * atan(1.0);
24061c66c397SJeremy L Thompson 
2407d1d35e2fSjeremylt   // Build q_ref_1d, q_weight_1d
2408d7b241e6Sjeremylt   // Set endpoints
24096574a04fSJeremy L Thompson   CeedCheck(Q > 1, NULL, CEED_ERROR_DIMENSION, "Cannot create Lobatto quadrature with Q=%" CeedInt_FMT " < 2 points", Q);
2410d7b241e6Sjeremylt   wi = 2.0 / ((CeedScalar)(Q * (Q - 1)));
2411d1d35e2fSjeremylt   if (q_weight_1d) {
2412d1d35e2fSjeremylt     q_weight_1d[0]     = wi;
2413d1d35e2fSjeremylt     q_weight_1d[Q - 1] = wi;
2414d7b241e6Sjeremylt   }
2415d1d35e2fSjeremylt   q_ref_1d[0]     = -1.0;
2416d1d35e2fSjeremylt   q_ref_1d[Q - 1] = 1.0;
2417d7b241e6Sjeremylt   // Interior
241892ae7e47SJeremy L Thompson   for (CeedInt i = 1; i <= (Q - 1) / 2; i++) {
2419d7b241e6Sjeremylt     // Guess
2420d7b241e6Sjeremylt     xi = cos(PI * (CeedScalar)(i) / (CeedScalar)(Q - 1));
2421d7b241e6Sjeremylt     // Pn(xi)
2422d7b241e6Sjeremylt     P0 = 1.0;
2423d7b241e6Sjeremylt     P1 = xi;
2424d7b241e6Sjeremylt     P2 = 0.0;
242592ae7e47SJeremy L Thompson     for (CeedInt j = 2; j < Q; j++) {
2426d7b241e6Sjeremylt       P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
2427d7b241e6Sjeremylt       P0 = P1;
2428d7b241e6Sjeremylt       P1 = P2;
2429d7b241e6Sjeremylt     }
2430d7b241e6Sjeremylt     // First Newton step
2431d7b241e6Sjeremylt     dP2  = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
2432d7b241e6Sjeremylt     d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi);
2433d7b241e6Sjeremylt     xi   = xi - dP2 / d2P2;
2434d7b241e6Sjeremylt     // Newton to convergence
243592ae7e47SJeremy L Thompson     for (CeedInt k = 0; k < 100 && fabs(dP2) > 10 * CEED_EPSILON; k++) {
2436d7b241e6Sjeremylt       P0 = 1.0;
2437d7b241e6Sjeremylt       P1 = xi;
243892ae7e47SJeremy L Thompson       for (CeedInt j = 2; j < Q; j++) {
2439d7b241e6Sjeremylt         P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
2440d7b241e6Sjeremylt         P0 = P1;
2441d7b241e6Sjeremylt         P1 = P2;
2442d7b241e6Sjeremylt       }
2443d7b241e6Sjeremylt       dP2  = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
2444d7b241e6Sjeremylt       d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi);
2445d7b241e6Sjeremylt       xi   = xi - dP2 / d2P2;
2446d7b241e6Sjeremylt     }
2447d7b241e6Sjeremylt     // Save xi, wi
2448d7b241e6Sjeremylt     wi = 2.0 / (((CeedScalar)(Q * (Q - 1))) * P2 * P2);
2449d1d35e2fSjeremylt     if (q_weight_1d) {
2450d1d35e2fSjeremylt       q_weight_1d[i]         = wi;
2451d1d35e2fSjeremylt       q_weight_1d[Q - 1 - i] = wi;
2452d7b241e6Sjeremylt     }
2453d1d35e2fSjeremylt     q_ref_1d[i]         = -xi;
2454d1d35e2fSjeremylt     q_ref_1d[Q - 1 - i] = xi;
2455d7b241e6Sjeremylt   }
2456e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2457d7b241e6Sjeremylt }
2458d7b241e6Sjeremylt 
2459d7b241e6Sjeremylt /// @}
2460