xref: /libCEED/rust/libceed-sys/c-src/interface/ceed-basis.c (revision 2247a93f19c75a71369fb56e8f7b7446bda8ddb1)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3d7b241e6Sjeremylt //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5d7b241e6Sjeremylt //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7d7b241e6Sjeremylt 
83d576824SJeremy L Thompson #include <ceed-impl.h>
949aac155SJeremy L Thompson #include <ceed.h>
102b730f8bSJeremy L Thompson #include <ceed/backend.h>
11d7b241e6Sjeremylt #include <math.h>
123d576824SJeremy L Thompson #include <stdbool.h>
13d7b241e6Sjeremylt #include <stdio.h>
14d7b241e6Sjeremylt #include <string.h>
15d7b241e6Sjeremylt 
167a982d89SJeremy L. Thompson /// @file
177a982d89SJeremy L. Thompson /// Implementation of CeedBasis interfaces
187a982d89SJeremy L. Thompson 
19d7b241e6Sjeremylt /// @cond DOXYGEN_SKIP
20356036faSJeremy L Thompson static struct CeedBasis_private ceed_basis_none;
21d7b241e6Sjeremylt /// @endcond
22d7b241e6Sjeremylt 
237a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser
247a982d89SJeremy L. Thompson /// @{
257a982d89SJeremy L. Thompson 
26ca94c3ddSJeremy L Thompson /// Argument for @ref CeedOperatorSetField() indicating that the field does not require a `CeedBasis`
27356036faSJeremy L Thompson const CeedBasis CEED_BASIS_NONE = &ceed_basis_none;
28356036faSJeremy L Thompson 
297a982d89SJeremy L. Thompson /// @}
307a982d89SJeremy L. Thompson 
317a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
327a982d89SJeremy L. Thompson /// CeedBasis Library Internal Functions
337a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
347a982d89SJeremy L. Thompson /// @addtogroup CeedBasisDeveloper
357a982d89SJeremy L. Thompson /// @{
367a982d89SJeremy L. Thompson 
377a982d89SJeremy L. Thompson /**
383778dbaaSJeremy L Thompson   @brief Compute Chebyshev polynomial values at a point
393778dbaaSJeremy L Thompson 
403778dbaaSJeremy L Thompson   @param[in]  x           Coordinate to evaluate Chebyshev polynomials at
41ca94c3ddSJeremy L Thompson   @param[in]  n           Number of Chebyshev polynomials to evaluate, `n >= 2`
423778dbaaSJeremy L Thompson   @param[out] chebyshev_x Array of Chebyshev polynomial values
433778dbaaSJeremy L Thompson 
443778dbaaSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
453778dbaaSJeremy L Thompson 
463778dbaaSJeremy L Thompson   @ref Developer
473778dbaaSJeremy L Thompson **/
483778dbaaSJeremy L Thompson static int CeedChebyshevPolynomialsAtPoint(CeedScalar x, CeedInt n, CeedScalar *chebyshev_x) {
493778dbaaSJeremy L Thompson   chebyshev_x[0] = 1.0;
503778dbaaSJeremy L Thompson   chebyshev_x[1] = 2 * x;
513778dbaaSJeremy L Thompson   for (CeedInt i = 2; i < n; i++) chebyshev_x[i] = 2 * x * chebyshev_x[i - 1] - chebyshev_x[i - 2];
523778dbaaSJeremy L Thompson   return CEED_ERROR_SUCCESS;
533778dbaaSJeremy L Thompson }
543778dbaaSJeremy L Thompson 
553778dbaaSJeremy L Thompson /**
563778dbaaSJeremy L Thompson   @brief Compute values of the derivative of Chebyshev polynomials at a point
573778dbaaSJeremy L Thompson 
583778dbaaSJeremy L Thompson   @param[in]  x            Coordinate to evaluate derivative of Chebyshev polynomials at
59ca94c3ddSJeremy L Thompson   @param[in]  n            Number of Chebyshev polynomials to evaluate, `n >= 2`
606cec60aaSJed Brown   @param[out] chebyshev_dx Array of Chebyshev polynomial derivative values
613778dbaaSJeremy L Thompson 
623778dbaaSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
633778dbaaSJeremy L Thompson 
643778dbaaSJeremy L Thompson   @ref Developer
653778dbaaSJeremy L Thompson **/
663778dbaaSJeremy L Thompson static int CeedChebyshevDerivativeAtPoint(CeedScalar x, CeedInt n, CeedScalar *chebyshev_dx) {
673778dbaaSJeremy L Thompson   CeedScalar chebyshev_x[3];
683778dbaaSJeremy L Thompson 
693778dbaaSJeremy L Thompson   chebyshev_x[1]  = 1.0;
703778dbaaSJeremy L Thompson   chebyshev_x[2]  = 2 * x;
713778dbaaSJeremy L Thompson   chebyshev_dx[0] = 0.0;
723778dbaaSJeremy L Thompson   chebyshev_dx[1] = 2.0;
733778dbaaSJeremy L Thompson   for (CeedInt i = 2; i < n; i++) {
743778dbaaSJeremy L Thompson     chebyshev_x[0]  = chebyshev_x[1];
753778dbaaSJeremy L Thompson     chebyshev_x[1]  = chebyshev_x[2];
763778dbaaSJeremy L Thompson     chebyshev_x[2]  = 2 * x * chebyshev_x[1] - chebyshev_x[0];
773778dbaaSJeremy L Thompson     chebyshev_dx[i] = 2 * x * chebyshev_dx[i - 1] + 2 * chebyshev_x[1] - chebyshev_dx[i - 2];
783778dbaaSJeremy L Thompson   }
793778dbaaSJeremy L Thompson   return CEED_ERROR_SUCCESS;
803778dbaaSJeremy L Thompson }
813778dbaaSJeremy L Thompson 
823778dbaaSJeremy L Thompson /**
83ca94c3ddSJeremy L Thompson   @brief Compute Householder reflection.
847a982d89SJeremy L. Thompson 
85ca94c3ddSJeremy L Thompson   Computes \f$A = (I - b v v^T) A\f$, where \f$A\f$ is an \f$m \times n\f$ matrix indexed as `A[i*row + j*col]`.
867a982d89SJeremy L. Thompson 
877a982d89SJeremy L. Thompson   @param[in,out] A   Matrix to apply Householder reflection to, in place
88ea61e9acSJeremy L Thompson   @param[in]     v   Householder vector
89ea61e9acSJeremy L Thompson   @param[in]     b   Scaling factor
90ca94c3ddSJeremy L Thompson   @param[in]     m   Number of rows in `A`
91ca94c3ddSJeremy L Thompson   @param[in]     n   Number of columns in `A`
92ea61e9acSJeremy L Thompson   @param[in]     row Row stride
93ea61e9acSJeremy L Thompson   @param[in]     col Col stride
947a982d89SJeremy L. Thompson 
957a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
967a982d89SJeremy L. Thompson 
977a982d89SJeremy L. Thompson   @ref Developer
987a982d89SJeremy L. Thompson **/
992b730f8bSJeremy L Thompson static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v, CeedScalar b, CeedInt m, CeedInt n, CeedInt row, CeedInt col) {
1007a982d89SJeremy L. Thompson   for (CeedInt j = 0; j < n; j++) {
1017a982d89SJeremy L. Thompson     CeedScalar w = A[0 * row + j * col];
1021c66c397SJeremy L Thompson 
1032b730f8bSJeremy L Thompson     for (CeedInt i = 1; i < m; i++) w += v[i] * A[i * row + j * col];
1047a982d89SJeremy L. Thompson     A[0 * row + j * col] -= b * w;
1052b730f8bSJeremy L Thompson     for (CeedInt i = 1; i < m; i++) A[i * row + j * col] -= b * w * v[i];
1067a982d89SJeremy L. Thompson   }
107e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1087a982d89SJeremy L. Thompson }
1097a982d89SJeremy L. Thompson 
1107a982d89SJeremy L. Thompson /**
1117a982d89SJeremy L. Thompson   @brief Compute Givens rotation
1127a982d89SJeremy L. Thompson 
113ca94c3ddSJeremy L Thompson   Computes \f$A = G A\f$ (or \f$G^T A\f$ in transpose mode), where \f$A\f$ is an \f$m \times n\f$ matrix indexed as `A[i*n + j*m]`.
1147a982d89SJeremy L. Thompson 
1157a982d89SJeremy L. Thompson   @param[in,out] A      Row major matrix to apply Givens rotation to, in place
116ea61e9acSJeremy L Thompson   @param[in]     c      Cosine factor
117ea61e9acSJeremy L Thompson   @param[in]     s      Sine factor
118ca94c3ddSJeremy L Thompson   @param[in]     t_mode @ref CEED_NOTRANSPOSE to rotate the basis counter-clockwise, which has the effect of rotating columns of `A` clockwise;
1194cc79fe7SJed Brown                           @ref CEED_TRANSPOSE for the opposite rotation
120ea61e9acSJeremy L Thompson   @param[in]     i      First row/column to apply rotation
121ea61e9acSJeremy L Thompson   @param[in]     k      Second row/column to apply rotation
122ca94c3ddSJeremy L Thompson   @param[in]     m      Number of rows in `A`
123ca94c3ddSJeremy L Thompson   @param[in]     n      Number of columns in `A`
1247a982d89SJeremy L. Thompson 
1257a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
1267a982d89SJeremy L. Thompson 
1277a982d89SJeremy L. Thompson   @ref Developer
1287a982d89SJeremy L. Thompson **/
1292b730f8bSJeremy L Thompson static int CeedGivensRotation(CeedScalar *A, CeedScalar c, CeedScalar s, CeedTransposeMode t_mode, CeedInt i, CeedInt k, CeedInt m, CeedInt n) {
130d1d35e2fSjeremylt   CeedInt stride_j = 1, stride_ik = m, num_its = n;
1311c66c397SJeremy L Thompson 
132d1d35e2fSjeremylt   if (t_mode == CEED_NOTRANSPOSE) {
1332b730f8bSJeremy L Thompson     stride_j  = n;
1342b730f8bSJeremy L Thompson     stride_ik = 1;
1352b730f8bSJeremy L Thompson     num_its   = m;
1367a982d89SJeremy L. Thompson   }
1377a982d89SJeremy L. Thompson 
1387a982d89SJeremy L. Thompson   // Apply rotation
139d1d35e2fSjeremylt   for (CeedInt j = 0; j < num_its; j++) {
140d1d35e2fSjeremylt     CeedScalar tau1 = A[i * stride_ik + j * stride_j], tau2 = A[k * stride_ik + j * stride_j];
1411c66c397SJeremy L Thompson 
142d1d35e2fSjeremylt     A[i * stride_ik + j * stride_j] = c * tau1 - s * tau2;
143d1d35e2fSjeremylt     A[k * stride_ik + j * stride_j] = s * tau1 + c * tau2;
1447a982d89SJeremy L. Thompson   }
145e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1467a982d89SJeremy L. Thompson }
1477a982d89SJeremy L. Thompson 
1487a982d89SJeremy L. Thompson /**
149ca94c3ddSJeremy L Thompson   @brief View an array stored in a `CeedBasis`
1507a982d89SJeremy L. Thompson 
1510a0da059Sjeremylt   @param[in] name   Name of array
152d1d35e2fSjeremylt   @param[in] fp_fmt Printing format
1530a0da059Sjeremylt   @param[in] m      Number of rows in array
1540a0da059Sjeremylt   @param[in] n      Number of columns in array
1550a0da059Sjeremylt   @param[in] a      Array to be viewed
156ca94c3ddSJeremy L Thompson   @param[in] stream Stream to view to, e.g., `stdout`
1577a982d89SJeremy L. Thompson 
1587a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
1597a982d89SJeremy L. Thompson 
1607a982d89SJeremy L. Thompson   @ref Developer
1617a982d89SJeremy L. Thompson **/
1622b730f8bSJeremy L Thompson static int CeedScalarView(const char *name, const char *fp_fmt, CeedInt m, CeedInt n, const CeedScalar *a, FILE *stream) {
163edf04919SJeremy L Thompson   if (m > 1) {
164edf04919SJeremy L Thompson     fprintf(stream, "  %s:\n", name);
165edf04919SJeremy L Thompson   } else {
166edf04919SJeremy L Thompson     char padded_name[12];
167edf04919SJeremy L Thompson 
168edf04919SJeremy L Thompson     snprintf(padded_name, 11, "%s:", name);
169edf04919SJeremy L Thompson     fprintf(stream, "  %-10s", padded_name);
170edf04919SJeremy L Thompson   }
17192ae7e47SJeremy L Thompson   for (CeedInt i = 0; i < m; i++) {
172edf04919SJeremy L Thompson     if (m > 1) fprintf(stream, "    [%" CeedInt_FMT "]", i);
1732b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < n; j++) fprintf(stream, fp_fmt, fabs(a[i * n + j]) > 1E-14 ? a[i * n + j] : 0);
1747a982d89SJeremy L. Thompson     fputs("\n", stream);
1757a982d89SJeremy L. Thompson   }
176e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1777a982d89SJeremy L. Thompson }
1787a982d89SJeremy L. Thompson 
179a76a04e7SJeremy L Thompson /**
180ea61e9acSJeremy L Thompson   @brief Create the interpolation and gradient matrices for projection from the nodes of `basis_from` to the nodes of `basis_to`.
181ba59ac12SSebastian Grimberg 
18215ad3917SSebastian Grimberg   The interpolation is given by `interp_project = interp_to^+ * interp_from`, where the pseudoinverse `interp_to^+` is given by QR factorization.
183ca94c3ddSJeremy L Thompson   The gradient is given by `grad_project = interp_to^+ * grad_from`, and is only computed for \f$H^1\f$ spaces otherwise it should not be used.
18415ad3917SSebastian Grimberg 
185ba59ac12SSebastian Grimberg   Note: `basis_from` and `basis_to` must have compatible quadrature spaces.
186a76a04e7SJeremy L Thompson 
187ca94c3ddSJeremy L Thompson   @param[in]  basis_from     `CeedBasis` to project from
188ca94c3ddSJeremy L Thompson   @param[in]  basis_to       `CeedBasis` to project to
189ca94c3ddSJeremy L Thompson   @param[out] interp_project Address of the variable where the newly created interpolation matrix will be stored
190ca94c3ddSJeremy L Thompson   @param[out] grad_project   Address of the variable where the newly created gradient matrix will be stored
191a76a04e7SJeremy L Thompson 
192a76a04e7SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
193a76a04e7SJeremy L Thompson 
194a76a04e7SJeremy L Thompson   @ref Developer
195a76a04e7SJeremy L Thompson **/
1962b730f8bSJeremy L Thompson static int CeedBasisCreateProjectionMatrices(CeedBasis basis_from, CeedBasis basis_to, CeedScalar **interp_project, CeedScalar **grad_project) {
197a76a04e7SJeremy L Thompson   Ceed    ceed;
1981c66c397SJeremy L Thompson   bool    is_tensor_to, is_tensor_from;
1991c66c397SJeremy L Thompson   CeedInt Q, Q_to, Q_from, P_to, P_from;
2001c66c397SJeremy L Thompson 
2012b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetCeed(basis_to, &ceed));
202a76a04e7SJeremy L Thompson 
203a76a04e7SJeremy L Thompson   // Check for compatible quadrature spaces
2042b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumQuadraturePoints(basis_to, &Q_to));
2052b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumQuadraturePoints(basis_from, &Q_from));
2066574a04fSJeremy L Thompson   CeedCheck(Q_to == Q_from, ceed, CEED_ERROR_DIMENSION, "Bases must have compatible quadrature spaces");
2071c66c397SJeremy L Thompson   Q = Q_to;
208a76a04e7SJeremy L Thompson 
20914556e63SJeremy L Thompson   // Check for matching tensor or non-tensor
2102b730f8bSJeremy L Thompson   CeedCall(CeedBasisIsTensor(basis_to, &is_tensor_to));
2112b730f8bSJeremy L Thompson   CeedCall(CeedBasisIsTensor(basis_from, &is_tensor_from));
2126574a04fSJeremy L Thompson   CeedCheck(is_tensor_to == is_tensor_from, ceed, CEED_ERROR_MINOR, "Bases must both be tensor or non-tensor");
2136574a04fSJeremy L Thompson   if (is_tensor_to) {
2142b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_to));
2152b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_from));
2162b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumQuadraturePoints1D(basis_from, &Q));
2176574a04fSJeremy L Thompson   } else {
2182b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis_to, &P_to));
2192b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis_from, &P_from));
220a76a04e7SJeremy L Thompson   }
221a76a04e7SJeremy L Thompson 
22215ad3917SSebastian Grimberg   // Check for matching FE space
22315ad3917SSebastian Grimberg   CeedFESpace fe_space_to, fe_space_from;
22415ad3917SSebastian Grimberg   CeedCall(CeedBasisGetFESpace(basis_to, &fe_space_to));
22515ad3917SSebastian Grimberg   CeedCall(CeedBasisGetFESpace(basis_from, &fe_space_from));
2266574a04fSJeremy L Thompson   CeedCheck(fe_space_to == fe_space_from, ceed, CEED_ERROR_MINOR, "Bases must both be the same FE space type");
22715ad3917SSebastian Grimberg 
22814556e63SJeremy L Thompson   // Get source matrices
22915ad3917SSebastian Grimberg   CeedInt           dim, q_comp = 1;
230*2247a93fSRezgar Shakeri   CeedScalar       *interp_to_inv, *interp_from;
2311c66c397SJeremy L Thompson   const CeedScalar *interp_to_source = NULL, *interp_from_source = NULL, *grad_from_source = NULL;
2321c66c397SJeremy L Thompson 
2332b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetDimension(basis_to, &dim));
234a76a04e7SJeremy L Thompson   if (is_tensor_to) {
2352b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetInterp1D(basis_to, &interp_to_source));
2362b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetInterp1D(basis_from, &interp_from_source));
237a76a04e7SJeremy L Thompson   } else {
23815ad3917SSebastian Grimberg     CeedCall(CeedBasisGetNumQuadratureComponents(basis_from, CEED_EVAL_INTERP, &q_comp));
2392b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetInterp(basis_to, &interp_to_source));
2402b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetInterp(basis_from, &interp_from_source));
24115ad3917SSebastian Grimberg   }
24215ad3917SSebastian Grimberg   CeedCall(CeedMalloc(Q * P_from * q_comp, &interp_from));
24315ad3917SSebastian Grimberg   CeedCall(CeedCalloc(P_to * P_from, interp_project));
24415ad3917SSebastian Grimberg 
24515ad3917SSebastian Grimberg   // `grad_project = interp_to^+ * grad_from` is computed for the H^1 space case so the
246de05fbb2SSebastian Grimberg   // projection basis will have a gradient operation (allocated even if not H^1 for the
247de05fbb2SSebastian Grimberg   // basis construction later on)
24815ad3917SSebastian Grimberg   if (fe_space_to == CEED_FE_SPACE_H1) {
24915ad3917SSebastian Grimberg     if (is_tensor_to) {
25015ad3917SSebastian Grimberg       CeedCall(CeedBasisGetGrad1D(basis_from, &grad_from_source));
25115ad3917SSebastian Grimberg     } else {
2522b730f8bSJeremy L Thompson       CeedCall(CeedBasisGetGrad(basis_from, &grad_from_source));
253a76a04e7SJeremy L Thompson     }
254de05fbb2SSebastian Grimberg   }
25515ad3917SSebastian Grimberg   CeedCall(CeedCalloc(P_to * P_from * (is_tensor_to ? 1 : dim), grad_project));
25615ad3917SSebastian Grimberg 
257*2247a93fSRezgar Shakeri   // Compute interp_to^+, pseudoinverse of interp_to
258*2247a93fSRezgar Shakeri   CeedCall(CeedCalloc(Q * q_comp * P_to, &interp_to_inv));
259*2247a93fSRezgar Shakeri   CeedCall(CeedMatrixPseudoinverse(ceed, (CeedScalar *)interp_to_source, Q * q_comp, P_to, interp_to_inv));
26014556e63SJeremy L Thompson   // Build matrices
26115ad3917SSebastian Grimberg   CeedInt     num_matrices = 1 + (fe_space_to == CEED_FE_SPACE_H1) * (is_tensor_to ? 1 : dim);
26214556e63SJeremy L Thompson   CeedScalar *input_from[num_matrices], *output_project[num_matrices];
2631c66c397SJeremy L Thompson 
26414556e63SJeremy L Thompson   input_from[0]     = (CeedScalar *)interp_from_source;
26514556e63SJeremy L Thompson   output_project[0] = *interp_project;
26614556e63SJeremy L Thompson   for (CeedInt m = 1; m < num_matrices; m++) {
26714556e63SJeremy L Thompson     input_from[m]     = (CeedScalar *)&grad_from_source[(m - 1) * Q * P_from];
26802af4036SJeremy L Thompson     output_project[m] = &((*grad_project)[(m - 1) * P_to * P_from]);
26914556e63SJeremy L Thompson   }
27014556e63SJeremy L Thompson   for (CeedInt m = 0; m < num_matrices; m++) {
271*2247a93fSRezgar Shakeri     // output_project = interp_to^+ * interp_from
27215ad3917SSebastian Grimberg     memcpy(interp_from, input_from[m], Q * P_from * q_comp * sizeof(input_from[m][0]));
273*2247a93fSRezgar Shakeri     CeedCall(CeedMatrixMatrixMultiply(ceed, interp_to_inv, input_from[m], output_project[m], P_to, P_from, Q * q_comp));
274*2247a93fSRezgar Shakeri     // Round zero to machine precision
275*2247a93fSRezgar Shakeri     for (CeedInt i = 0; i < P_to * P_from; i++) {
276*2247a93fSRezgar Shakeri       if (fabs(output_project[m][i]) < 10 * CEED_EPSILON) output_project[m][i] = 0.0;
277a76a04e7SJeremy L Thompson     }
27814556e63SJeremy L Thompson   }
27914556e63SJeremy L Thompson 
28014556e63SJeremy L Thompson   // Cleanup
281*2247a93fSRezgar Shakeri   CeedCall(CeedFree(&interp_to_inv));
2822b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_from));
283a76a04e7SJeremy L Thompson   return CEED_ERROR_SUCCESS;
284a76a04e7SJeremy L Thompson }
285a76a04e7SJeremy L Thompson 
2867a982d89SJeremy L. Thompson /// @}
2877a982d89SJeremy L. Thompson 
2887a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
2897a982d89SJeremy L. Thompson /// Ceed Backend API
2907a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
2917a982d89SJeremy L. Thompson /// @addtogroup CeedBasisBackend
2927a982d89SJeremy L. Thompson /// @{
2937a982d89SJeremy L. Thompson 
2947a982d89SJeremy L. Thompson /**
295ca94c3ddSJeremy L Thompson   @brief Return collocated gradient matrix
2967a982d89SJeremy L. Thompson 
297ca94c3ddSJeremy L Thompson   @param[in]  basis         `CeedBasis`
298ca94c3ddSJeremy L Thompson   @param[out] collo_grad_1d Row-major (`Q_1d * Q_1d`) matrix expressing derivatives of basis functions at quadrature points
2997a982d89SJeremy L. Thompson 
3007a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3017a982d89SJeremy L. Thompson 
3027a982d89SJeremy L. Thompson   @ref Backend
3037a982d89SJeremy L. Thompson **/
304d1d35e2fSjeremylt int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *collo_grad_1d) {
3057a982d89SJeremy L. Thompson   Ceed        ceed;
306*2247a93fSRezgar Shakeri   CeedInt     P_1d, Q_1d;
307*2247a93fSRezgar Shakeri   CeedScalar *interp_1d_pinv;
308ea61e9acSJeremy 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.
309*2247a93fSRezgar Shakeri   CeedCall(CeedBasisGetCeed(basis, &ceed));
310*2247a93fSRezgar Shakeri   CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d));
311*2247a93fSRezgar Shakeri   CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
3127a982d89SJeremy L. Thompson 
313*2247a93fSRezgar Shakeri   // Compute interp_1d^+, pseudoinverse of interp_1d
314*2247a93fSRezgar Shakeri   CeedCall(CeedCalloc(P_1d * Q_1d, &interp_1d_pinv));
3157a982d89SJeremy L. Thompson 
316*2247a93fSRezgar Shakeri   CeedCall(CeedMatrixPseudoinverse(ceed, basis->interp_1d, Q_1d, P_1d, interp_1d_pinv));
317*2247a93fSRezgar Shakeri   CeedCall(CeedMatrixMatrixMultiply(ceed, basis->grad_1d, (const CeedScalar *)interp_1d_pinv, collo_grad_1d, Q_1d, Q_1d, P_1d));
3187a982d89SJeremy L. Thompson 
319*2247a93fSRezgar Shakeri   CeedCall(CeedFree(&interp_1d_pinv));
320e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3217a982d89SJeremy L. Thompson }
3227a982d89SJeremy L. Thompson 
3237a982d89SJeremy L. Thompson /**
324ca94c3ddSJeremy L Thompson   @brief Get tensor status for given `CeedBasis`
3257a982d89SJeremy L. Thompson 
326ca94c3ddSJeremy L Thompson   @param[in]  basis     `CeedBasis`
327d1d35e2fSjeremylt   @param[out] is_tensor Variable to store tensor status
3287a982d89SJeremy L. Thompson 
3297a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3307a982d89SJeremy L. Thompson 
3317a982d89SJeremy L. Thompson   @ref Backend
3327a982d89SJeremy L. Thompson **/
333d1d35e2fSjeremylt int CeedBasisIsTensor(CeedBasis basis, bool *is_tensor) {
3346402da51SJeremy L Thompson   *is_tensor = basis->is_tensor_basis;
335e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3367a982d89SJeremy L. Thompson }
3377a982d89SJeremy L. Thompson 
3387a982d89SJeremy L. Thompson /**
339ca94c3ddSJeremy L Thompson   @brief Get backend data of a `CeedBasis`
3407a982d89SJeremy L. Thompson 
341ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
3427a982d89SJeremy L. Thompson   @param[out] data  Variable to store data
3437a982d89SJeremy L. Thompson 
3447a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3457a982d89SJeremy L. Thompson 
3467a982d89SJeremy L. Thompson   @ref Backend
3477a982d89SJeremy L. Thompson **/
348777ff853SJeremy L Thompson int CeedBasisGetData(CeedBasis basis, void *data) {
349777ff853SJeremy L Thompson   *(void **)data = basis->data;
350e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3517a982d89SJeremy L. Thompson }
3527a982d89SJeremy L. Thompson 
3537a982d89SJeremy L. Thompson /**
354ca94c3ddSJeremy L Thompson   @brief Set backend data of a `CeedBasis`
3557a982d89SJeremy L. Thompson 
356ca94c3ddSJeremy L Thompson   @param[in,out] basis  `CeedBasis`
357ea61e9acSJeremy L Thompson   @param[in]     data   Data to set
3587a982d89SJeremy L. Thompson 
3597a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3607a982d89SJeremy L. Thompson 
3617a982d89SJeremy L. Thompson   @ref Backend
3627a982d89SJeremy L. Thompson **/
363777ff853SJeremy L Thompson int CeedBasisSetData(CeedBasis basis, void *data) {
364777ff853SJeremy L Thompson   basis->data = data;
365e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3667a982d89SJeremy L. Thompson }
3677a982d89SJeremy L. Thompson 
3687a982d89SJeremy L. Thompson /**
369ca94c3ddSJeremy L Thompson   @brief Increment the reference counter for a `CeedBasis`
37034359f16Sjeremylt 
371ca94c3ddSJeremy L Thompson   @param[in,out] basis `CeedBasis` to increment the reference counter
37234359f16Sjeremylt 
37334359f16Sjeremylt   @return An error code: 0 - success, otherwise - failure
37434359f16Sjeremylt 
37534359f16Sjeremylt   @ref Backend
37634359f16Sjeremylt **/
3779560d06aSjeremylt int CeedBasisReference(CeedBasis basis) {
37834359f16Sjeremylt   basis->ref_count++;
37934359f16Sjeremylt   return CEED_ERROR_SUCCESS;
38034359f16Sjeremylt }
38134359f16Sjeremylt 
38234359f16Sjeremylt /**
383ca94c3ddSJeremy L Thompson   @brief Get number of Q-vector components for given `CeedBasis`
384c4e3f59bSSebastian Grimberg 
385ca94c3ddSJeremy L Thompson   @param[in]  basis     `CeedBasis`
386ca94c3ddSJeremy L Thompson   @param[in]  eval_mode @ref CEED_EVAL_INTERP to use interpolated values,
387ca94c3ddSJeremy L Thompson                           @ref CEED_EVAL_GRAD to use gradients,
388ca94c3ddSJeremy L Thompson                           @ref CEED_EVAL_DIV to use divergence,
389ca94c3ddSJeremy L Thompson                           @ref CEED_EVAL_CURL to use curl
390c4e3f59bSSebastian Grimberg   @param[out] q_comp    Variable to store number of Q-vector components of basis
391c4e3f59bSSebastian Grimberg 
392c4e3f59bSSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
393c4e3f59bSSebastian Grimberg 
394c4e3f59bSSebastian Grimberg   @ref Backend
395c4e3f59bSSebastian Grimberg **/
396c4e3f59bSSebastian Grimberg int CeedBasisGetNumQuadratureComponents(CeedBasis basis, CeedEvalMode eval_mode, CeedInt *q_comp) {
397c4e3f59bSSebastian Grimberg   switch (eval_mode) {
398c4e3f59bSSebastian Grimberg     case CEED_EVAL_INTERP:
399c4e3f59bSSebastian Grimberg       *q_comp = (basis->fe_space == CEED_FE_SPACE_H1) ? 1 : basis->dim;
400c4e3f59bSSebastian Grimberg       break;
401c4e3f59bSSebastian Grimberg     case CEED_EVAL_GRAD:
402c4e3f59bSSebastian Grimberg       *q_comp = basis->dim;
403c4e3f59bSSebastian Grimberg       break;
404c4e3f59bSSebastian Grimberg     case CEED_EVAL_DIV:
405c4e3f59bSSebastian Grimberg       *q_comp = 1;
406c4e3f59bSSebastian Grimberg       break;
407c4e3f59bSSebastian Grimberg     case CEED_EVAL_CURL:
408c4e3f59bSSebastian Grimberg       *q_comp = (basis->dim < 3) ? 1 : basis->dim;
409c4e3f59bSSebastian Grimberg       break;
410c4e3f59bSSebastian Grimberg     case CEED_EVAL_NONE:
411c4e3f59bSSebastian Grimberg     case CEED_EVAL_WEIGHT:
412352a5e7cSSebastian Grimberg       *q_comp = 1;
413c4e3f59bSSebastian Grimberg       break;
414c4e3f59bSSebastian Grimberg   }
415c4e3f59bSSebastian Grimberg   return CEED_ERROR_SUCCESS;
416c4e3f59bSSebastian Grimberg }
417c4e3f59bSSebastian Grimberg 
418c4e3f59bSSebastian Grimberg /**
419ca94c3ddSJeremy L Thompson   @brief Estimate number of FLOPs required to apply `CeedBasis` in `t_mode` and `eval_mode`
4206e15d496SJeremy L Thompson 
421ca94c3ddSJeremy L Thompson   @param[in]  basis     `CeedBasis` to estimate FLOPs for
422ea61e9acSJeremy L Thompson   @param[in]  t_mode    Apply basis or transpose
423ca94c3ddSJeremy L Thompson   @param[in]  eval_mode @ref CeedEvalMode
424ea61e9acSJeremy L Thompson   @param[out] flops     Address of variable to hold FLOPs estimate
4256e15d496SJeremy L Thompson 
4266e15d496SJeremy L Thompson   @ref Backend
4276e15d496SJeremy L Thompson **/
4282b730f8bSJeremy L Thompson int CeedBasisGetFlopsEstimate(CeedBasis basis, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedSize *flops) {
4296e15d496SJeremy L Thompson   bool is_tensor;
4306e15d496SJeremy L Thompson 
4312b730f8bSJeremy L Thompson   CeedCall(CeedBasisIsTensor(basis, &is_tensor));
4326e15d496SJeremy L Thompson   if (is_tensor) {
4336e15d496SJeremy L Thompson     CeedInt dim, num_comp, P_1d, Q_1d;
4341c66c397SJeremy L Thompson 
4352b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetDimension(basis, &dim));
4362b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
4372b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d));
4382b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
4396e15d496SJeremy L Thompson     if (t_mode == CEED_TRANSPOSE) {
4402b730f8bSJeremy L Thompson       P_1d = Q_1d;
4412b730f8bSJeremy L Thompson       Q_1d = P_1d;
4426e15d496SJeremy L Thompson     }
4436e15d496SJeremy L Thompson     CeedInt tensor_flops = 0, pre = num_comp * CeedIntPow(P_1d, dim - 1), post = 1;
4446e15d496SJeremy L Thompson     for (CeedInt d = 0; d < dim; d++) {
4456e15d496SJeremy L Thompson       tensor_flops += 2 * pre * P_1d * post * Q_1d;
4466e15d496SJeremy L Thompson       pre /= P_1d;
4476e15d496SJeremy L Thompson       post *= Q_1d;
4486e15d496SJeremy L Thompson     }
4496e15d496SJeremy L Thompson     switch (eval_mode) {
4502b730f8bSJeremy L Thompson       case CEED_EVAL_NONE:
4512b730f8bSJeremy L Thompson         *flops = 0;
4522b730f8bSJeremy L Thompson         break;
4532b730f8bSJeremy L Thompson       case CEED_EVAL_INTERP:
4542b730f8bSJeremy L Thompson         *flops = tensor_flops;
4552b730f8bSJeremy L Thompson         break;
4562b730f8bSJeremy L Thompson       case CEED_EVAL_GRAD:
4572b730f8bSJeremy L Thompson         *flops = tensor_flops * 2;
4582b730f8bSJeremy L Thompson         break;
4596e15d496SJeremy L Thompson       case CEED_EVAL_DIV:
4606e15d496SJeremy L Thompson       case CEED_EVAL_CURL:
4616574a04fSJeremy L Thompson         // LCOV_EXCL_START
4626574a04fSJeremy L Thompson         return CeedError(basis->ceed, CEED_ERROR_INCOMPATIBLE, "Tensor basis evaluation for %s not supported", CeedEvalModes[eval_mode]);
4632b730f8bSJeremy L Thompson         break;
4646e15d496SJeremy L Thompson       // LCOV_EXCL_STOP
4652b730f8bSJeremy L Thompson       case CEED_EVAL_WEIGHT:
4662b730f8bSJeremy L Thompson         *flops = dim * CeedIntPow(Q_1d, dim);
4672b730f8bSJeremy L Thompson         break;
4686e15d496SJeremy L Thompson     }
4696e15d496SJeremy L Thompson   } else {
470c4e3f59bSSebastian Grimberg     CeedInt dim, num_comp, q_comp, num_nodes, num_qpts;
4711c66c397SJeremy L Thompson 
4722b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetDimension(basis, &dim));
4732b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
474c4e3f59bSSebastian Grimberg     CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp));
4752b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis, &num_nodes));
4762b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts));
4776e15d496SJeremy L Thompson     switch (eval_mode) {
4782b730f8bSJeremy L Thompson       case CEED_EVAL_NONE:
4792b730f8bSJeremy L Thompson         *flops = 0;
4802b730f8bSJeremy L Thompson         break;
4812b730f8bSJeremy L Thompson       case CEED_EVAL_INTERP:
4822b730f8bSJeremy L Thompson       case CEED_EVAL_GRAD:
4832b730f8bSJeremy L Thompson       case CEED_EVAL_DIV:
4842b730f8bSJeremy L Thompson       case CEED_EVAL_CURL:
485c4e3f59bSSebastian Grimberg         *flops = num_nodes * num_qpts * num_comp * q_comp;
4862b730f8bSJeremy L Thompson         break;
4872b730f8bSJeremy L Thompson       case CEED_EVAL_WEIGHT:
4882b730f8bSJeremy L Thompson         *flops = 0;
4892b730f8bSJeremy L Thompson         break;
4906e15d496SJeremy L Thompson     }
4916e15d496SJeremy L Thompson   }
4926e15d496SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4936e15d496SJeremy L Thompson }
4946e15d496SJeremy L Thompson 
4956e15d496SJeremy L Thompson /**
496ca94c3ddSJeremy L Thompson   @brief Get `CeedFESpace` for a `CeedBasis`
497c4e3f59bSSebastian Grimberg 
498ca94c3ddSJeremy L Thompson   @param[in]  basis    `CeedBasis`
499ca94c3ddSJeremy L Thompson   @param[out] fe_space Variable to store `CeedFESpace`
500c4e3f59bSSebastian Grimberg 
501c4e3f59bSSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
502c4e3f59bSSebastian Grimberg 
503c4e3f59bSSebastian Grimberg   @ref Backend
504c4e3f59bSSebastian Grimberg **/
505c4e3f59bSSebastian Grimberg int CeedBasisGetFESpace(CeedBasis basis, CeedFESpace *fe_space) {
506c4e3f59bSSebastian Grimberg   *fe_space = basis->fe_space;
507c4e3f59bSSebastian Grimberg   return CEED_ERROR_SUCCESS;
508c4e3f59bSSebastian Grimberg }
509c4e3f59bSSebastian Grimberg 
510c4e3f59bSSebastian Grimberg /**
511ca94c3ddSJeremy L Thompson   @brief Get dimension for given `CeedElemTopology`
5127a982d89SJeremy L. Thompson 
513ca94c3ddSJeremy L Thompson   @param[in]  topo `CeedElemTopology`
5147a982d89SJeremy L. Thompson   @param[out] dim  Variable to store dimension of topology
5157a982d89SJeremy L. Thompson 
5167a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
5177a982d89SJeremy L. Thompson 
5187a982d89SJeremy L. Thompson   @ref Backend
5197a982d89SJeremy L. Thompson **/
5207a982d89SJeremy L. Thompson int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) {
5217a982d89SJeremy L. Thompson   *dim = (CeedInt)topo >> 16;
522e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5237a982d89SJeremy L. Thompson }
5247a982d89SJeremy L. Thompson 
5257a982d89SJeremy L. Thompson /**
526ca94c3ddSJeremy L Thompson   @brief Get `CeedTensorContract` of a `CeedBasis`
5277a982d89SJeremy L. Thompson 
528ca94c3ddSJeremy L Thompson   @param[in]  basis     `CeedBasis`
529ca94c3ddSJeremy L Thompson   @param[out] contract  Variable to store `CeedTensorContract`
5307a982d89SJeremy L. Thompson 
5317a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
5327a982d89SJeremy L. Thompson 
5337a982d89SJeremy L. Thompson   @ref Backend
5347a982d89SJeremy L. Thompson **/
5357a982d89SJeremy L. Thompson int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) {
5367a982d89SJeremy L. Thompson   *contract = basis->contract;
537e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5387a982d89SJeremy L. Thompson }
5397a982d89SJeremy L. Thompson 
5407a982d89SJeremy L. Thompson /**
541ca94c3ddSJeremy L Thompson   @brief Set `CeedTensorContract` of a `CeedBasis`
5427a982d89SJeremy L. Thompson 
543ca94c3ddSJeremy L Thompson   @param[in,out] basis    `CeedBasis`
544ca94c3ddSJeremy L Thompson   @param[in]     contract `CeedTensorContract` to set
5457a982d89SJeremy L. Thompson 
5467a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
5477a982d89SJeremy L. Thompson 
5487a982d89SJeremy L. Thompson   @ref Backend
5497a982d89SJeremy L. Thompson **/
55034359f16Sjeremylt int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract contract) {
55134359f16Sjeremylt   basis->contract = contract;
5522b730f8bSJeremy L Thompson   CeedCall(CeedTensorContractReference(contract));
553e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5547a982d89SJeremy L. Thompson }
5557a982d89SJeremy L. Thompson 
5567a982d89SJeremy L. Thompson /**
557ca94c3ddSJeremy L Thompson   @brief Return a reference implementation of matrix multiplication \f$C = A B\f$.
558ba59ac12SSebastian Grimberg 
559ca94c3ddSJeremy L Thompson   Note: This is a reference implementation for CPU `CeedScalar` pointers that is not intended for high performance.
5607a982d89SJeremy L. Thompson 
561ca94c3ddSJeremy L Thompson   @param[in]  ceed  `Ceed` context for error handling
562ca94c3ddSJeremy L Thompson   @param[in]  mat_A Row-major matrix `A`
563ca94c3ddSJeremy L Thompson   @param[in]  mat_B Row-major matrix `B`
564ca94c3ddSJeremy L Thompson   @param[out] mat_C Row-major output matrix `C`
565ca94c3ddSJeremy L Thompson   @param[in]  m     Number of rows of `C`
566ca94c3ddSJeremy L Thompson   @param[in]  n     Number of columns of `C`
567ca94c3ddSJeremy L Thompson   @param[in]  kk    Number of columns of `A`/rows of `B`
5687a982d89SJeremy L. Thompson 
5697a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
5707a982d89SJeremy L. Thompson 
5717a982d89SJeremy L. Thompson   @ref Utility
5727a982d89SJeremy L. Thompson **/
5732b730f8bSJeremy L Thompson int CeedMatrixMatrixMultiply(Ceed ceed, const CeedScalar *mat_A, const CeedScalar *mat_B, CeedScalar *mat_C, CeedInt m, CeedInt n, CeedInt kk) {
5742b730f8bSJeremy L Thompson   for (CeedInt i = 0; i < m; i++) {
5757a982d89SJeremy L. Thompson     for (CeedInt j = 0; j < n; j++) {
5767a982d89SJeremy L. Thompson       CeedScalar sum = 0;
5771c66c397SJeremy L Thompson 
5782b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < kk; k++) sum += mat_A[k + i * kk] * mat_B[j + k * n];
579d1d35e2fSjeremylt       mat_C[j + i * n] = sum;
5807a982d89SJeremy L. Thompson     }
5812b730f8bSJeremy L Thompson   }
582e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5837a982d89SJeremy L. Thompson }
5847a982d89SJeremy L. Thompson 
585ba59ac12SSebastian Grimberg /**
586ba59ac12SSebastian Grimberg   @brief Return QR Factorization of a matrix
587ba59ac12SSebastian Grimberg 
588ca94c3ddSJeremy L Thompson   @param[in]     ceed `Ceed` context for error handling
589ba59ac12SSebastian Grimberg   @param[in,out] mat  Row-major matrix to be factorized in place
590ca94c3ddSJeremy L Thompson   @param[in,out] tau  Vector of length `m` of scaling factors
591ba59ac12SSebastian Grimberg   @param[in]     m    Number of rows
592ba59ac12SSebastian Grimberg   @param[in]     n    Number of columns
593ba59ac12SSebastian Grimberg 
594ba59ac12SSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
595ba59ac12SSebastian Grimberg 
596ba59ac12SSebastian Grimberg   @ref Utility
597ba59ac12SSebastian Grimberg **/
598ba59ac12SSebastian Grimberg int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau, CeedInt m, CeedInt n) {
599ba59ac12SSebastian Grimberg   CeedScalar v[m];
600ba59ac12SSebastian Grimberg 
601ba59ac12SSebastian Grimberg   // Check matrix shape
6026574a04fSJeremy L Thompson   CeedCheck(n <= m, ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute QR factorization with n > m");
603ba59ac12SSebastian Grimberg 
604ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) {
6051c66c397SJeremy L Thompson     CeedScalar sigma = 0.0;
6061c66c397SJeremy L Thompson 
607ba59ac12SSebastian Grimberg     if (i >= m - 1) {  // last row of matrix, no reflection needed
608ba59ac12SSebastian Grimberg       tau[i] = 0.;
609ba59ac12SSebastian Grimberg       break;
610ba59ac12SSebastian Grimberg     }
611ba59ac12SSebastian Grimberg     // Calculate Householder vector, magnitude
612ba59ac12SSebastian Grimberg     v[i] = mat[i + n * i];
613ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < m; j++) {
614ba59ac12SSebastian Grimberg       v[j] = mat[i + n * j];
615ba59ac12SSebastian Grimberg       sigma += v[j] * v[j];
616ba59ac12SSebastian Grimberg     }
6171c66c397SJeremy L Thompson     const CeedScalar norm = sqrt(v[i] * v[i] + sigma);  // norm of v[i:m]
6181c66c397SJeremy L Thompson     const CeedScalar R_ii = -copysign(norm, v[i]);
6191c66c397SJeremy L Thompson 
620ba59ac12SSebastian Grimberg     v[i] -= R_ii;
621ba59ac12SSebastian Grimberg     // norm of v[i:m] after modification above and scaling below
622ba59ac12SSebastian Grimberg     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
623ba59ac12SSebastian Grimberg     //   tau = 2 / (norm*norm)
624ba59ac12SSebastian Grimberg     tau[i] = 2 * v[i] * v[i] / (v[i] * v[i] + sigma);
625ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < m; j++) v[j] /= v[i];
626ba59ac12SSebastian Grimberg 
627ba59ac12SSebastian Grimberg     // Apply Householder reflector to lower right panel
628ba59ac12SSebastian Grimberg     CeedHouseholderReflect(&mat[i * n + i + 1], &v[i], tau[i], m - i, n - i - 1, n, 1);
629ba59ac12SSebastian Grimberg     // Save v
630ba59ac12SSebastian Grimberg     mat[i + n * i] = R_ii;
631ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < m; j++) mat[i + n * j] = v[j];
632ba59ac12SSebastian Grimberg   }
633ba59ac12SSebastian Grimberg   return CEED_ERROR_SUCCESS;
634ba59ac12SSebastian Grimberg }
635ba59ac12SSebastian Grimberg 
636ba59ac12SSebastian Grimberg /**
637ba59ac12SSebastian Grimberg   @brief Apply Householder Q matrix
638ba59ac12SSebastian Grimberg 
639ca94c3ddSJeremy 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$.
640ba59ac12SSebastian Grimberg 
641ba59ac12SSebastian Grimberg   @param[in,out] mat_A  Matrix to apply Householder Q to, in place
642ba59ac12SSebastian Grimberg   @param[in]     mat_Q  Householder Q matrix
643ba59ac12SSebastian Grimberg   @param[in]     tau    Householder scaling factors
644ba59ac12SSebastian Grimberg   @param[in]     t_mode Transpose mode for application
645ca94c3ddSJeremy L Thompson   @param[in]     m      Number of rows in `A`
646ca94c3ddSJeremy L Thompson   @param[in]     n      Number of columns in `A`
647ca94c3ddSJeremy L Thompson   @param[in]     k      Number of elementary reflectors in Q, `k < m`
648ca94c3ddSJeremy L Thompson   @param[in]     row    Row stride in `A`
649ca94c3ddSJeremy L Thompson   @param[in]     col    Col stride in `A`
650ba59ac12SSebastian Grimberg 
651ba59ac12SSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
652ba59ac12SSebastian Grimberg 
653c4e3f59bSSebastian Grimberg   @ref Utility
654ba59ac12SSebastian Grimberg **/
655ba59ac12SSebastian Grimberg int CeedHouseholderApplyQ(CeedScalar *mat_A, const CeedScalar *mat_Q, const CeedScalar *tau, CeedTransposeMode t_mode, CeedInt m, CeedInt n,
656ba59ac12SSebastian Grimberg                           CeedInt k, CeedInt row, CeedInt col) {
657ba59ac12SSebastian Grimberg   CeedScalar *v;
6581c66c397SJeremy L Thompson 
659ba59ac12SSebastian Grimberg   CeedCall(CeedMalloc(m, &v));
660ba59ac12SSebastian Grimberg   for (CeedInt ii = 0; ii < k; ii++) {
661ba59ac12SSebastian Grimberg     CeedInt i = t_mode == CEED_TRANSPOSE ? ii : k - 1 - ii;
662ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < m; j++) v[j] = mat_Q[j * k + i];
663ba59ac12SSebastian Grimberg     // Apply Householder reflector (I - tau v v^T) collo_grad_1d^T
664ba59ac12SSebastian Grimberg     CeedCall(CeedHouseholderReflect(&mat_A[i * row], &v[i], tau[i], m - i, n, row, col));
665ba59ac12SSebastian Grimberg   }
666ba59ac12SSebastian Grimberg   CeedCall(CeedFree(&v));
667ba59ac12SSebastian Grimberg   return CEED_ERROR_SUCCESS;
668ba59ac12SSebastian Grimberg }
669ba59ac12SSebastian Grimberg 
670ba59ac12SSebastian Grimberg /**
671*2247a93fSRezgar Shakeri   @brief Return pseudoinverse of a matrix
672*2247a93fSRezgar Shakeri 
673*2247a93fSRezgar Shakeri   @param[in]     ceed      Ceed context for error handling
674*2247a93fSRezgar Shakeri   @param[in]     mat       Row-major matrix to compute pseudoinverse of
675*2247a93fSRezgar Shakeri   @param[in]     m         Number of rows
676*2247a93fSRezgar Shakeri   @param[in]     n         Number of columns
677*2247a93fSRezgar Shakeri   @param[out]    mat_pinv  Row-major pseudoinverse matrix
678*2247a93fSRezgar Shakeri 
679*2247a93fSRezgar Shakeri   @return An error code: 0 - success, otherwise - failure
680*2247a93fSRezgar Shakeri 
681*2247a93fSRezgar Shakeri   @ref Utility
682*2247a93fSRezgar Shakeri **/
683*2247a93fSRezgar Shakeri int CeedMatrixPseudoinverse(Ceed ceed, CeedScalar *mat, CeedInt m, CeedInt n, CeedScalar *mat_pinv) {
684*2247a93fSRezgar Shakeri   CeedScalar *tau, *I, *mat_copy;
685*2247a93fSRezgar Shakeri 
686*2247a93fSRezgar Shakeri   CeedCall(CeedCalloc(m, &tau));
687*2247a93fSRezgar Shakeri   CeedCall(CeedCalloc(m * m, &I));
688*2247a93fSRezgar Shakeri   CeedCall(CeedCalloc(m * n, &mat_copy));
689*2247a93fSRezgar Shakeri   memcpy(mat_copy, mat, m * n * sizeof mat[0]);
690*2247a93fSRezgar Shakeri 
691*2247a93fSRezgar Shakeri   // QR Factorization, mat = Q R
692*2247a93fSRezgar Shakeri   CeedCall(CeedQRFactorization(ceed, mat_copy, tau, m, n));
693*2247a93fSRezgar Shakeri 
694*2247a93fSRezgar Shakeri   // -- Apply Q^T, I = Q^T * I
695*2247a93fSRezgar Shakeri   for (CeedInt i = 0; i < m; i++) I[i * m + i] = 1.0;
696*2247a93fSRezgar Shakeri   CeedCall(CeedHouseholderApplyQ(I, mat_copy, tau, CEED_TRANSPOSE, m, m, n, m, 1));
697*2247a93fSRezgar Shakeri   // -- Apply R_inv, mat_pinv = R_inv * Q^T
698*2247a93fSRezgar Shakeri   for (CeedInt j = 0; j < m; j++) {  // Column j
699*2247a93fSRezgar Shakeri     mat_pinv[j + m * (n - 1)] = I[j + m * (n - 1)] / mat_copy[n * n - 1];
700*2247a93fSRezgar Shakeri     for (CeedInt i = n - 2; i >= 0; i--) {  // Row i
701*2247a93fSRezgar Shakeri       mat_pinv[j + m * i] = I[j + m * i];
702*2247a93fSRezgar Shakeri       for (CeedInt k = i + 1; k < n; k++) mat_pinv[j + m * i] -= mat_copy[k + n * i] * mat_pinv[j + m * k];
703*2247a93fSRezgar Shakeri       mat_pinv[j + m * i] /= mat_copy[i + n * i];
704*2247a93fSRezgar Shakeri     }
705*2247a93fSRezgar Shakeri   }
706*2247a93fSRezgar Shakeri 
707*2247a93fSRezgar Shakeri   // Cleanup
708*2247a93fSRezgar Shakeri   CeedCall(CeedFree(&I));
709*2247a93fSRezgar Shakeri   CeedCall(CeedFree(&tau));
710*2247a93fSRezgar Shakeri   CeedCall(CeedFree(&mat_copy));
711*2247a93fSRezgar Shakeri   return CEED_ERROR_SUCCESS;
712*2247a93fSRezgar Shakeri }
713*2247a93fSRezgar Shakeri 
714*2247a93fSRezgar Shakeri /**
715ba59ac12SSebastian Grimberg   @brief Return symmetric Schur decomposition of the symmetric matrix mat via symmetric QR factorization
716ba59ac12SSebastian Grimberg 
717ca94c3ddSJeremy L Thompson   @param[in]     ceed   `Ceed` context for error handling
718ba59ac12SSebastian Grimberg   @param[in,out] mat    Row-major matrix to be factorized in place
719ba59ac12SSebastian Grimberg   @param[out]    lambda Vector of length n of eigenvalues
720ba59ac12SSebastian Grimberg   @param[in]     n      Number of rows/columns
721ba59ac12SSebastian Grimberg 
722ba59ac12SSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
723ba59ac12SSebastian Grimberg 
724ba59ac12SSebastian Grimberg   @ref Utility
725ba59ac12SSebastian Grimberg **/
7262c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOff
7272c2ea1dbSJeremy L Thompson int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat, CeedScalar *lambda, CeedInt n) {
728ba59ac12SSebastian Grimberg   // Check bounds for clang-tidy
7296574a04fSJeremy L Thompson   CeedCheck(n > 1, ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute symmetric Schur decomposition of scalars");
730ba59ac12SSebastian Grimberg 
731ba59ac12SSebastian Grimberg   CeedScalar v[n - 1], tau[n - 1], mat_T[n * n];
732ba59ac12SSebastian Grimberg 
733ba59ac12SSebastian Grimberg   // Copy mat to mat_T and set mat to I
734ba59ac12SSebastian Grimberg   memcpy(mat_T, mat, n * n * sizeof(mat[0]));
735ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) {
736ba59ac12SSebastian Grimberg     for (CeedInt j = 0; j < n; j++) mat[j + n * i] = (i == j) ? 1 : 0;
737ba59ac12SSebastian Grimberg   }
738ba59ac12SSebastian Grimberg 
739ba59ac12SSebastian Grimberg   // Reduce to tridiagonal
740ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n - 1; i++) {
741ba59ac12SSebastian Grimberg     // Calculate Householder vector, magnitude
742ba59ac12SSebastian Grimberg     CeedScalar sigma = 0.0;
7431c66c397SJeremy L Thompson 
744ba59ac12SSebastian Grimberg     v[i] = mat_T[i + n * (i + 1)];
745ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < n - 1; j++) {
746ba59ac12SSebastian Grimberg       v[j] = mat_T[i + n * (j + 1)];
747ba59ac12SSebastian Grimberg       sigma += v[j] * v[j];
748ba59ac12SSebastian Grimberg     }
7491c66c397SJeremy L Thompson     const CeedScalar norm = sqrt(v[i] * v[i] + sigma);  // norm of v[i:n-1]
7501c66c397SJeremy L Thompson     const CeedScalar R_ii = -copysign(norm, v[i]);
7511c66c397SJeremy L Thompson 
752ba59ac12SSebastian Grimberg     v[i] -= R_ii;
753ba59ac12SSebastian Grimberg     // norm of v[i:m] after modification above and scaling below
754ba59ac12SSebastian Grimberg     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
755ba59ac12SSebastian Grimberg     //   tau = 2 / (norm*norm)
756ba59ac12SSebastian Grimberg     tau[i] = i == n - 2 ? 2 : 2 * v[i] * v[i] / (v[i] * v[i] + sigma);
757ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < n - 1; j++) v[j] /= v[i];
758ba59ac12SSebastian Grimberg 
759ba59ac12SSebastian Grimberg     // Update sub and super diagonal
760ba59ac12SSebastian Grimberg     for (CeedInt j = i + 2; j < n; j++) {
761ba59ac12SSebastian Grimberg       mat_T[i + n * j] = 0;
762ba59ac12SSebastian Grimberg       mat_T[j + n * i] = 0;
763ba59ac12SSebastian Grimberg     }
764ba59ac12SSebastian Grimberg     // Apply symmetric Householder reflector to lower right panel
765ba59ac12SSebastian Grimberg     CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1);
766ba59ac12SSebastian Grimberg     CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), 1, n);
767ba59ac12SSebastian Grimberg 
768ba59ac12SSebastian Grimberg     // Save v
769ba59ac12SSebastian Grimberg     mat_T[i + n * (i + 1)] = R_ii;
770ba59ac12SSebastian Grimberg     mat_T[(i + 1) + n * i] = R_ii;
771ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < n - 1; j++) {
772ba59ac12SSebastian Grimberg       mat_T[i + n * (j + 1)] = v[j];
773ba59ac12SSebastian Grimberg     }
774ba59ac12SSebastian Grimberg   }
775ba59ac12SSebastian Grimberg   // Backwards accumulation of Q
776ba59ac12SSebastian Grimberg   for (CeedInt i = n - 2; i >= 0; i--) {
777ba59ac12SSebastian Grimberg     if (tau[i] > 0.0) {
778ba59ac12SSebastian Grimberg       v[i] = 1;
779ba59ac12SSebastian Grimberg       for (CeedInt j = i + 1; j < n - 1; j++) {
780ba59ac12SSebastian Grimberg         v[j]                   = mat_T[i + n * (j + 1)];
781ba59ac12SSebastian Grimberg         mat_T[i + n * (j + 1)] = 0;
782ba59ac12SSebastian Grimberg       }
783ba59ac12SSebastian Grimberg       CeedHouseholderReflect(&mat[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1);
784ba59ac12SSebastian Grimberg     }
785ba59ac12SSebastian Grimberg   }
786ba59ac12SSebastian Grimberg 
787ba59ac12SSebastian Grimberg   // Reduce sub and super diagonal
788ba59ac12SSebastian Grimberg   CeedInt    p = 0, q = 0, itr = 0, max_itr = n * n * n * n;
789ba59ac12SSebastian Grimberg   CeedScalar tol = CEED_EPSILON;
790ba59ac12SSebastian Grimberg 
791ba59ac12SSebastian Grimberg   while (itr < max_itr) {
792ba59ac12SSebastian Grimberg     // Update p, q, size of reduced portions of diagonal
793ba59ac12SSebastian Grimberg     p = 0;
794ba59ac12SSebastian Grimberg     q = 0;
795ba59ac12SSebastian Grimberg     for (CeedInt i = n - 2; i >= 0; i--) {
796ba59ac12SSebastian Grimberg       if (fabs(mat_T[i + n * (i + 1)]) < tol) q += 1;
797ba59ac12SSebastian Grimberg       else break;
798ba59ac12SSebastian Grimberg     }
799ba59ac12SSebastian Grimberg     for (CeedInt i = 0; i < n - q - 1; i++) {
800ba59ac12SSebastian Grimberg       if (fabs(mat_T[i + n * (i + 1)]) < tol) p += 1;
801ba59ac12SSebastian Grimberg       else break;
802ba59ac12SSebastian Grimberg     }
803ba59ac12SSebastian Grimberg     if (q == n - 1) break;  // Finished reducing
804ba59ac12SSebastian Grimberg 
805ba59ac12SSebastian Grimberg     // Reduce tridiagonal portion
806ba59ac12SSebastian Grimberg     CeedScalar t_nn = mat_T[(n - 1 - q) + n * (n - 1 - q)], t_nnm1 = mat_T[(n - 2 - q) + n * (n - 1 - q)];
807ba59ac12SSebastian Grimberg     CeedScalar d  = (mat_T[(n - 2 - q) + n * (n - 2 - q)] - t_nn) / 2;
808ba59ac12SSebastian Grimberg     CeedScalar mu = t_nn - t_nnm1 * t_nnm1 / (d + copysign(sqrt(d * d + t_nnm1 * t_nnm1), d));
809ba59ac12SSebastian Grimberg     CeedScalar x  = mat_T[p + n * p] - mu;
810ba59ac12SSebastian Grimberg     CeedScalar z  = mat_T[p + n * (p + 1)];
8111c66c397SJeremy L Thompson 
812ba59ac12SSebastian Grimberg     for (CeedInt k = p; k < n - q - 1; k++) {
813ba59ac12SSebastian Grimberg       // Compute Givens rotation
814ba59ac12SSebastian Grimberg       CeedScalar c = 1, s = 0;
8151c66c397SJeremy L Thompson 
816ba59ac12SSebastian Grimberg       if (fabs(z) > tol) {
817ba59ac12SSebastian Grimberg         if (fabs(z) > fabs(x)) {
8181c66c397SJeremy L Thompson           const CeedScalar tau = -x / z;
8191c66c397SJeremy L Thompson 
8201c66c397SJeremy L Thompson           s = 1 / sqrt(1 + tau * tau);
8211c66c397SJeremy L Thompson           c = s * tau;
822ba59ac12SSebastian Grimberg         } else {
8231c66c397SJeremy L Thompson           const CeedScalar tau = -z / x;
8241c66c397SJeremy L Thompson 
8251c66c397SJeremy L Thompson           c = 1 / sqrt(1 + tau * tau);
8261c66c397SJeremy L Thompson           s = c * tau;
827ba59ac12SSebastian Grimberg         }
828ba59ac12SSebastian Grimberg       }
829ba59ac12SSebastian Grimberg 
830ba59ac12SSebastian Grimberg       // Apply Givens rotation to T
831ba59ac12SSebastian Grimberg       CeedGivensRotation(mat_T, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n);
832ba59ac12SSebastian Grimberg       CeedGivensRotation(mat_T, c, s, CEED_TRANSPOSE, k, k + 1, n, n);
833ba59ac12SSebastian Grimberg 
834ba59ac12SSebastian Grimberg       // Apply Givens rotation to Q
835ba59ac12SSebastian Grimberg       CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n);
836ba59ac12SSebastian Grimberg 
837ba59ac12SSebastian Grimberg       // Update x, z
838ba59ac12SSebastian Grimberg       if (k < n - q - 2) {
839ba59ac12SSebastian Grimberg         x = mat_T[k + n * (k + 1)];
840ba59ac12SSebastian Grimberg         z = mat_T[k + n * (k + 2)];
841ba59ac12SSebastian Grimberg       }
842ba59ac12SSebastian Grimberg     }
843ba59ac12SSebastian Grimberg     itr++;
844ba59ac12SSebastian Grimberg   }
845ba59ac12SSebastian Grimberg 
846ba59ac12SSebastian Grimberg   // Save eigenvalues
847ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) lambda[i] = mat_T[i + n * i];
848ba59ac12SSebastian Grimberg 
849ba59ac12SSebastian Grimberg   // Check convergence
8506574a04fSJeremy L Thompson   CeedCheck(itr < max_itr || q > n, ceed, CEED_ERROR_MINOR, "Symmetric QR failed to converge");
851ba59ac12SSebastian Grimberg   return CEED_ERROR_SUCCESS;
852ba59ac12SSebastian Grimberg }
8532c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOn
854ba59ac12SSebastian Grimberg 
855ba59ac12SSebastian Grimberg /**
856ba59ac12SSebastian Grimberg   @brief Return Simultaneous Diagonalization of two matrices.
857ba59ac12SSebastian Grimberg 
858ca94c3ddSJeremy L Thompson   This solves the generalized eigenvalue problem `A x = lambda B x`, where `A` and `B` are symmetric and `B` is positive definite.
859ca94c3ddSJeremy L Thompson   We generate the matrix `X` and vector `Lambda` such that `X^T A X = Lambda` and `X^T B X = I`.
860ca94c3ddSJeremy L Thompson   This is equivalent to the LAPACK routine 'sygv' with `TYPE = 1`.
861ba59ac12SSebastian Grimberg 
862ca94c3ddSJeremy L Thompson   @param[in]  ceed   `Ceed` context for error handling
863ba59ac12SSebastian Grimberg   @param[in]  mat_A  Row-major matrix to be factorized with eigenvalues
864ba59ac12SSebastian Grimberg   @param[in]  mat_B  Row-major matrix to be factorized to identity
865ba59ac12SSebastian Grimberg   @param[out] mat_X  Row-major orthogonal matrix
866ca94c3ddSJeremy L Thompson   @param[out] lambda Vector of length `n` of generalized eigenvalues
867ba59ac12SSebastian Grimberg   @param[in]  n      Number of rows/columns
868ba59ac12SSebastian Grimberg 
869ba59ac12SSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
870ba59ac12SSebastian Grimberg 
871ba59ac12SSebastian Grimberg   @ref Utility
872ba59ac12SSebastian Grimberg **/
8732c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOff
8742c2ea1dbSJeremy L Thompson int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *mat_A, CeedScalar *mat_B, CeedScalar *mat_X, CeedScalar *lambda, CeedInt n) {
875ba59ac12SSebastian Grimberg   CeedScalar *mat_C, *mat_G, *vec_D;
8761c66c397SJeremy L Thompson 
877ba59ac12SSebastian Grimberg   CeedCall(CeedCalloc(n * n, &mat_C));
878ba59ac12SSebastian Grimberg   CeedCall(CeedCalloc(n * n, &mat_G));
879ba59ac12SSebastian Grimberg   CeedCall(CeedCalloc(n, &vec_D));
880ba59ac12SSebastian Grimberg 
881ba59ac12SSebastian Grimberg   // Compute B = G D G^T
882ba59ac12SSebastian Grimberg   memcpy(mat_G, mat_B, n * n * sizeof(mat_B[0]));
883ba59ac12SSebastian Grimberg   CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_G, vec_D, n));
884ba59ac12SSebastian Grimberg 
885ba59ac12SSebastian Grimberg   // Sort eigenvalues
886ba59ac12SSebastian Grimberg   for (CeedInt i = n - 1; i >= 0; i--) {
887ba59ac12SSebastian Grimberg     for (CeedInt j = 0; j < i; j++) {
888ba59ac12SSebastian Grimberg       if (fabs(vec_D[j]) > fabs(vec_D[j + 1])) {
8891c66c397SJeremy L Thompson         CeedScalarSwap(vec_D[j], vec_D[j + 1]);
8901c66c397SJeremy L Thompson         for (CeedInt k = 0; k < n; k++) CeedScalarSwap(mat_G[k * n + j], mat_G[k * n + j + 1]);
891ba59ac12SSebastian Grimberg       }
892ba59ac12SSebastian Grimberg     }
893ba59ac12SSebastian Grimberg   }
894ba59ac12SSebastian Grimberg 
895ba59ac12SSebastian Grimberg   // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T
896ba59ac12SSebastian Grimberg   //           = D^-1/2 G^T A G D^-1/2
897ba59ac12SSebastian Grimberg   // -- D = D^-1/2
898ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) vec_D[i] = 1. / sqrt(vec_D[i]);
899ba59ac12SSebastian Grimberg   // -- G = G D^-1/2
900ba59ac12SSebastian Grimberg   // -- C = D^-1/2 G^T
901ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) {
902ba59ac12SSebastian Grimberg     for (CeedInt j = 0; j < n; j++) {
903ba59ac12SSebastian Grimberg       mat_G[i * n + j] *= vec_D[j];
904ba59ac12SSebastian Grimberg       mat_C[j * n + i] = mat_G[i * n + j];
905ba59ac12SSebastian Grimberg     }
906ba59ac12SSebastian Grimberg   }
907ba59ac12SSebastian Grimberg   // -- X = (D^-1/2 G^T) A
908ba59ac12SSebastian Grimberg   CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_C, (const CeedScalar *)mat_A, mat_X, n, n, n));
909ba59ac12SSebastian Grimberg   // -- C = (D^-1/2 G^T A) (G D^-1/2)
910ba59ac12SSebastian Grimberg   CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_X, (const CeedScalar *)mat_G, mat_C, n, n, n));
911ba59ac12SSebastian Grimberg 
912ba59ac12SSebastian Grimberg   // Compute Q^T C Q = lambda
913ba59ac12SSebastian Grimberg   CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_C, lambda, n));
914ba59ac12SSebastian Grimberg 
915ba59ac12SSebastian Grimberg   // Sort eigenvalues
916ba59ac12SSebastian Grimberg   for (CeedInt i = n - 1; i >= 0; i--) {
917ba59ac12SSebastian Grimberg     for (CeedInt j = 0; j < i; j++) {
918ba59ac12SSebastian Grimberg       if (fabs(lambda[j]) > fabs(lambda[j + 1])) {
9191c66c397SJeremy L Thompson         CeedScalarSwap(lambda[j], lambda[j + 1]);
9201c66c397SJeremy L Thompson         for (CeedInt k = 0; k < n; k++) CeedScalarSwap(mat_C[k * n + j], mat_C[k * n + j + 1]);
921ba59ac12SSebastian Grimberg       }
922ba59ac12SSebastian Grimberg     }
923ba59ac12SSebastian Grimberg   }
924ba59ac12SSebastian Grimberg 
925ba59ac12SSebastian Grimberg   // Set X = (G D^1/2)^-T Q
926ba59ac12SSebastian Grimberg   //       = G D^-1/2 Q
927ba59ac12SSebastian Grimberg   CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_G, (const CeedScalar *)mat_C, mat_X, n, n, n));
928ba59ac12SSebastian Grimberg 
929ba59ac12SSebastian Grimberg   // Cleanup
930ba59ac12SSebastian Grimberg   CeedCall(CeedFree(&mat_C));
931ba59ac12SSebastian Grimberg   CeedCall(CeedFree(&mat_G));
932ba59ac12SSebastian Grimberg   CeedCall(CeedFree(&vec_D));
933ba59ac12SSebastian Grimberg   return CEED_ERROR_SUCCESS;
934ba59ac12SSebastian Grimberg }
9352c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOn
936ba59ac12SSebastian Grimberg 
9377a982d89SJeremy L. Thompson /// @}
9387a982d89SJeremy L. Thompson 
9397a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
9407a982d89SJeremy L. Thompson /// CeedBasis Public API
9417a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
9427a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser
943d7b241e6Sjeremylt /// @{
944d7b241e6Sjeremylt 
945b11c1e72Sjeremylt /**
946ca94c3ddSJeremy L Thompson   @brief Create a tensor-product basis for \f$H^1\f$ discretizations
947b11c1e72Sjeremylt 
948ca94c3ddSJeremy L Thompson   @param[in]  ceed        `Ceed` object used to create the `CeedBasis`
949ea61e9acSJeremy L Thompson   @param[in]  dim         Topological dimension
950ea61e9acSJeremy L Thompson   @param[in]  num_comp    Number of field components (1 for scalar fields)
951ea61e9acSJeremy L Thompson   @param[in]  P_1d        Number of nodes in one dimension
952ea61e9acSJeremy L Thompson   @param[in]  Q_1d        Number of quadrature points in one dimension
953ca94c3ddSJeremy L Thompson   @param[in]  interp_1d   Row-major (`Q_1d * P_1d`) matrix expressing the values of nodal basis functions at quadrature points
954ca94c3ddSJeremy L Thompson   @param[in]  grad_1d     Row-major (`Q_1d * P_1d`) matrix expressing derivatives of nodal basis functions at quadrature points
955ca94c3ddSJeremy 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]`
956ca94c3ddSJeremy L Thompson   @param[in]  q_weight_1d Array of length `Q_1d` holding the quadrature weights on the reference element
957ca94c3ddSJeremy L Thompson   @param[out] basis       Address of the variable where the newly created `CeedBasis` will be stored
958b11c1e72Sjeremylt 
959b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
960dfdf5a53Sjeremylt 
9617a982d89SJeremy L. Thompson   @ref User
962b11c1e72Sjeremylt **/
9632b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d,
9642b730f8bSJeremy L Thompson                             const CeedScalar *grad_1d, const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis *basis) {
9655fe0d4faSjeremylt   if (!ceed->BasisCreateTensorH1) {
9665fe0d4faSjeremylt     Ceed delegate;
9676574a04fSJeremy L Thompson 
9682b730f8bSJeremy L Thompson     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
9696574a04fSJeremy L Thompson     CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisCreateTensorH1");
9702b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateTensorH1(delegate, dim, num_comp, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis));
971e15f9bd0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
9725fe0d4faSjeremylt   }
973e15f9bd0SJeremy L Thompson 
974ca94c3ddSJeremy L Thompson   CeedCheck(dim > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis dimension must be a positive value");
975ca94c3ddSJeremy L Thompson   CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 component");
976ca94c3ddSJeremy L Thompson   CeedCheck(P_1d > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 node");
977ca94c3ddSJeremy L Thompson   CeedCheck(Q_1d > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 quadrature point");
978227444bfSJeremy L Thompson 
9792b730f8bSJeremy L Thompson   CeedElemTopology topo = dim == 1 ? CEED_TOPOLOGY_LINE : dim == 2 ? CEED_TOPOLOGY_QUAD : CEED_TOPOLOGY_HEX;
980e15f9bd0SJeremy L Thompson 
9812b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(1, basis));
982db002c03SJeremy L Thompson   CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed));
983d1d35e2fSjeremylt   (*basis)->ref_count       = 1;
9846402da51SJeremy L Thompson   (*basis)->is_tensor_basis = true;
985d7b241e6Sjeremylt   (*basis)->dim             = dim;
986d99fa3c5SJeremy L Thompson   (*basis)->topo            = topo;
987d1d35e2fSjeremylt   (*basis)->num_comp        = num_comp;
988d1d35e2fSjeremylt   (*basis)->P_1d            = P_1d;
989d1d35e2fSjeremylt   (*basis)->Q_1d            = Q_1d;
990d1d35e2fSjeremylt   (*basis)->P               = CeedIntPow(P_1d, dim);
991d1d35e2fSjeremylt   (*basis)->Q               = CeedIntPow(Q_1d, dim);
992c4e3f59bSSebastian Grimberg   (*basis)->fe_space        = CEED_FE_SPACE_H1;
9932b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d, &(*basis)->q_ref_1d));
9942b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d, &(*basis)->q_weight_1d));
995ff3a0f91SJeremy L Thompson   if (q_ref_1d) memcpy((*basis)->q_ref_1d, q_ref_1d, Q_1d * sizeof(q_ref_1d[0]));
9962b730f8bSJeremy L Thompson   if (q_weight_1d) memcpy((*basis)->q_weight_1d, q_weight_1d, Q_1d * sizeof(q_weight_1d[0]));
9972b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->interp_1d));
9982b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->grad_1d));
9992b730f8bSJeremy L Thompson   if (interp_1d) memcpy((*basis)->interp_1d, interp_1d, Q_1d * P_1d * sizeof(interp_1d[0]));
1000ff3a0f91SJeremy L Thompson   if (grad_1d) memcpy((*basis)->grad_1d, grad_1d, Q_1d * P_1d * sizeof(grad_1d[0]));
10012b730f8bSJeremy L Thompson   CeedCall(ceed->BasisCreateTensorH1(dim, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, *basis));
1002e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1003d7b241e6Sjeremylt }
1004d7b241e6Sjeremylt 
1005b11c1e72Sjeremylt /**
1006ca94c3ddSJeremy L Thompson   @brief Create a tensor-product \f$H^1\f$ Lagrange basis
1007b11c1e72Sjeremylt 
1008ca94c3ddSJeremy L Thompson   @param[in]  ceed      `Ceed` object used to create the `CeedBasis`
1009ea61e9acSJeremy L Thompson   @param[in]  dim       Topological dimension of element
1010ea61e9acSJeremy L Thompson   @param[in]  num_comp  Number of field components (1 for scalar fields)
1011ea61e9acSJeremy L Thompson   @param[in]  P         Number of Gauss-Lobatto nodes in one dimension.
1012ca94c3ddSJeremy L Thompson                           The polynomial degree of the resulting `Q_k` element is `k = P - 1`.
1013ea61e9acSJeremy L Thompson   @param[in]  Q         Number of quadrature points in one dimension.
1014ca94c3ddSJeremy L Thompson   @param[in]  quad_mode Distribution of the `Q` quadrature points (affects order of accuracy for the quadrature)
1015ca94c3ddSJeremy L Thompson   @param[out] basis     Address of the variable where the newly created `CeedBasis` will be stored
1016b11c1e72Sjeremylt 
1017b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1018dfdf5a53Sjeremylt 
10197a982d89SJeremy L. Thompson   @ref User
1020b11c1e72Sjeremylt **/
10212b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P, CeedInt Q, CeedQuadMode quad_mode, CeedBasis *basis) {
1022d7b241e6Sjeremylt   // Allocate
1023c8c3fa7dSJeremy L Thompson   int        ierr = CEED_ERROR_SUCCESS;
10242b730f8bSJeremy L Thompson   CeedScalar c1, c2, c3, c4, dx, *nodes, *interp_1d, *grad_1d, *q_ref_1d, *q_weight_1d;
10254d537eeaSYohann 
1026ca94c3ddSJeremy L Thompson   CeedCheck(dim > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis dimension must be a positive value");
1027ca94c3ddSJeremy L Thompson   CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 component");
1028ca94c3ddSJeremy L Thompson   CeedCheck(P > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 node");
1029ca94c3ddSJeremy L Thompson   CeedCheck(Q > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 quadrature point");
1030227444bfSJeremy L Thompson 
1031e15f9bd0SJeremy L Thompson   // Get Nodes and Weights
10322b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(P * Q, &interp_1d));
10332b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(P * Q, &grad_1d));
10342b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(P, &nodes));
10352b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q, &q_ref_1d));
10362b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q, &q_weight_1d));
10372b730f8bSJeremy L Thompson   if (CeedLobattoQuadrature(P, nodes, NULL) != CEED_ERROR_SUCCESS) goto cleanup;
1038d1d35e2fSjeremylt   switch (quad_mode) {
1039d7b241e6Sjeremylt     case CEED_GAUSS:
1040d1d35e2fSjeremylt       ierr = CeedGaussQuadrature(Q, q_ref_1d, q_weight_1d);
1041d7b241e6Sjeremylt       break;
1042d7b241e6Sjeremylt     case CEED_GAUSS_LOBATTO:
1043d1d35e2fSjeremylt       ierr = CeedLobattoQuadrature(Q, q_ref_1d, q_weight_1d);
1044d7b241e6Sjeremylt       break;
1045d7b241e6Sjeremylt   }
10462b730f8bSJeremy L Thompson   if (ierr != CEED_ERROR_SUCCESS) goto cleanup;
1047e15f9bd0SJeremy L Thompson 
1048d7b241e6Sjeremylt   // Build B, D matrix
1049d7b241e6Sjeremylt   // Fornberg, 1998
1050c8c3fa7dSJeremy L Thompson   for (CeedInt i = 0; i < Q; i++) {
1051d7b241e6Sjeremylt     c1                   = 1.0;
1052d1d35e2fSjeremylt     c3                   = nodes[0] - q_ref_1d[i];
1053d1d35e2fSjeremylt     interp_1d[i * P + 0] = 1.0;
1054c8c3fa7dSJeremy L Thompson     for (CeedInt j = 1; j < P; j++) {
1055d7b241e6Sjeremylt       c2 = 1.0;
1056d7b241e6Sjeremylt       c4 = c3;
1057d1d35e2fSjeremylt       c3 = nodes[j] - q_ref_1d[i];
1058c8c3fa7dSJeremy L Thompson       for (CeedInt k = 0; k < j; k++) {
1059d7b241e6Sjeremylt         dx = nodes[j] - nodes[k];
1060d7b241e6Sjeremylt         c2 *= dx;
1061d7b241e6Sjeremylt         if (k == j - 1) {
1062d1d35e2fSjeremylt           grad_1d[i * P + j]   = c1 * (interp_1d[i * P + k] - c4 * grad_1d[i * P + k]) / c2;
1063d1d35e2fSjeremylt           interp_1d[i * P + j] = -c1 * c4 * interp_1d[i * P + k] / c2;
1064d7b241e6Sjeremylt         }
1065d1d35e2fSjeremylt         grad_1d[i * P + k]   = (c3 * grad_1d[i * P + k] - interp_1d[i * P + k]) / dx;
1066d1d35e2fSjeremylt         interp_1d[i * P + k] = c3 * interp_1d[i * P + k] / dx;
1067d7b241e6Sjeremylt       }
1068d7b241e6Sjeremylt       c1 = c2;
1069d7b241e6Sjeremylt     }
1070d7b241e6Sjeremylt   }
10719ac7b42eSJeremy L Thompson   // Pass to CeedBasisCreateTensorH1
10722b730f8bSJeremy L Thompson   CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P, Q, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis));
1073e15f9bd0SJeremy L Thompson cleanup:
10742b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_1d));
10752b730f8bSJeremy L Thompson   CeedCall(CeedFree(&grad_1d));
10762b730f8bSJeremy L Thompson   CeedCall(CeedFree(&nodes));
10772b730f8bSJeremy L Thompson   CeedCall(CeedFree(&q_ref_1d));
10782b730f8bSJeremy L Thompson   CeedCall(CeedFree(&q_weight_1d));
1079e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1080d7b241e6Sjeremylt }
1081d7b241e6Sjeremylt 
1082b11c1e72Sjeremylt /**
1083ca94c3ddSJeremy L Thompson   @brief Create a non tensor-product basis for \f$H^1\f$ discretizations
1084a8de75f0Sjeremylt 
1085ca94c3ddSJeremy L Thompson   @param[in]  ceed      `Ceed` object used to create the `CeedBasis`
1086ea61e9acSJeremy L Thompson   @param[in]  topo      Topology of element, e.g. hypercube, simplex, ect
1087ea61e9acSJeremy L Thompson   @param[in]  num_comp  Number of field components (1 for scalar fields)
1088ea61e9acSJeremy L Thompson   @param[in]  num_nodes Total number of nodes
1089ea61e9acSJeremy L Thompson   @param[in]  num_qpts  Total number of quadrature points
1090ca94c3ddSJeremy L Thompson   @param[in]  interp    Row-major (`num_qpts * num_nodes`) matrix expressing the values of nodal basis functions at quadrature points
1091ca94c3ddSJeremy L Thompson   @param[in]  grad      Row-major (`dim * num_qpts * num_nodes`) matrix expressing derivatives of nodal basis functions at quadrature points
1092ca94c3ddSJeremy L Thompson   @param[in]  q_ref     Array of length `num_qpts` * dim holding the locations of quadrature points on the reference element
1093ca94c3ddSJeremy L Thompson   @param[in]  q_weight  Array of length `num_qpts` holding the quadrature weights on the reference element
1094ca94c3ddSJeremy L Thompson   @param[out] basis     Address of the variable where the newly created `CeedBasis` will be stored
1095a8de75f0Sjeremylt 
1096a8de75f0Sjeremylt   @return An error code: 0 - success, otherwise - failure
1097a8de75f0Sjeremylt 
10987a982d89SJeremy L. Thompson   @ref User
1099a8de75f0Sjeremylt **/
11002b730f8bSJeremy L Thompson int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
11012b730f8bSJeremy L Thompson                       const CeedScalar *grad, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) {
1102d1d35e2fSjeremylt   CeedInt P = num_nodes, Q = num_qpts, dim = 0;
1103a8de75f0Sjeremylt 
11045fe0d4faSjeremylt   if (!ceed->BasisCreateH1) {
11055fe0d4faSjeremylt     Ceed delegate;
11066574a04fSJeremy L Thompson 
11072b730f8bSJeremy L Thompson     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
11086574a04fSJeremy L Thompson     CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisCreateH1");
11092b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateH1(delegate, topo, num_comp, num_nodes, num_qpts, interp, grad, q_ref, q_weight, basis));
1110e15f9bd0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
11115fe0d4faSjeremylt   }
11125fe0d4faSjeremylt 
1113ca94c3ddSJeremy L Thompson   CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 component");
1114ca94c3ddSJeremy L Thompson   CeedCheck(num_nodes > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 node");
1115ca94c3ddSJeremy L Thompson   CeedCheck(num_qpts > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 quadrature point");
1116227444bfSJeremy L Thompson 
11172b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetTopologyDimension(topo, &dim));
1118a8de75f0Sjeremylt 
1119db002c03SJeremy L Thompson   CeedCall(CeedCalloc(1, basis));
1120db002c03SJeremy L Thompson   CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed));
1121d1d35e2fSjeremylt   (*basis)->ref_count       = 1;
11226402da51SJeremy L Thompson   (*basis)->is_tensor_basis = false;
1123a8de75f0Sjeremylt   (*basis)->dim             = dim;
1124d99fa3c5SJeremy L Thompson   (*basis)->topo            = topo;
1125d1d35e2fSjeremylt   (*basis)->num_comp        = num_comp;
1126a8de75f0Sjeremylt   (*basis)->P               = P;
1127a8de75f0Sjeremylt   (*basis)->Q               = Q;
1128c4e3f59bSSebastian Grimberg   (*basis)->fe_space        = CEED_FE_SPACE_H1;
11292b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q * dim, &(*basis)->q_ref_1d));
11302b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q, &(*basis)->q_weight_1d));
1131ff3a0f91SJeremy L Thompson   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0]));
1132ff3a0f91SJeremy L Thompson   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0]));
11332b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q * P, &(*basis)->interp));
11342b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(dim * Q * P, &(*basis)->grad));
1135ff3a0f91SJeremy L Thompson   if (interp) memcpy((*basis)->interp, interp, Q * P * sizeof(interp[0]));
1136ff3a0f91SJeremy L Thompson   if (grad) memcpy((*basis)->grad, grad, dim * Q * P * sizeof(grad[0]));
11372b730f8bSJeremy L Thompson   CeedCall(ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, q_ref, q_weight, *basis));
1138e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1139a8de75f0Sjeremylt }
1140a8de75f0Sjeremylt 
1141a8de75f0Sjeremylt /**
1142859c15bbSJames Wright   @brief Create a non tensor-product basis for \f$H(\mathrm{div})\f$ discretizations
114350c301a5SRezgar Shakeri 
1144ca94c3ddSJeremy L Thompson   @param[in]  ceed      `Ceed` object used to create the `CeedBasis`
1145ea61e9acSJeremy 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
1146ea61e9acSJeremy L Thompson   @param[in]  num_comp  Number of components (usually 1 for vectors in H(div) bases)
1147ca94c3ddSJeremy L Thompson   @param[in]  num_nodes Total number of nodes (DoFs per element)
1148ea61e9acSJeremy L Thompson   @param[in]  num_qpts  Total number of quadrature points
1149ca94c3ddSJeremy L Thompson   @param[in]  interp    Row-major (`dim * num_qpts * num_nodes`) matrix expressing the values of basis functions at quadrature points
1150ca94c3ddSJeremy L Thompson   @param[in]  div       Row-major (`num_qpts * num_nodes`) matrix expressing divergence of basis functions at quadrature points
1151ca94c3ddSJeremy L Thompson   @param[in]  q_ref     Array of length `num_qpts` * dim holding the locations of quadrature points on the reference element
1152ca94c3ddSJeremy L Thompson   @param[in]  q_weight  Array of length `num_qpts` holding the quadrature weights on the reference element
1153ca94c3ddSJeremy L Thompson   @param[out] basis     Address of the variable where the newly created `CeedBasis` will be stored
115450c301a5SRezgar Shakeri 
115550c301a5SRezgar Shakeri   @return An error code: 0 - success, otherwise - failure
115650c301a5SRezgar Shakeri 
115750c301a5SRezgar Shakeri   @ref User
115850c301a5SRezgar Shakeri **/
11592b730f8bSJeremy L Thompson int CeedBasisCreateHdiv(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
11602b730f8bSJeremy L Thompson                         const CeedScalar *div, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) {
116150c301a5SRezgar Shakeri   CeedInt Q = num_qpts, P = num_nodes, dim = 0;
1162c4e3f59bSSebastian Grimberg 
116350c301a5SRezgar Shakeri   if (!ceed->BasisCreateHdiv) {
116450c301a5SRezgar Shakeri     Ceed delegate;
11656574a04fSJeremy L Thompson 
11662b730f8bSJeremy L Thompson     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
11676574a04fSJeremy L Thompson     CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateHdiv");
11682b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateHdiv(delegate, topo, num_comp, num_nodes, num_qpts, interp, div, q_ref, q_weight, basis));
116950c301a5SRezgar Shakeri     return CEED_ERROR_SUCCESS;
117050c301a5SRezgar Shakeri   }
117150c301a5SRezgar Shakeri 
1172ca94c3ddSJeremy L Thompson   CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 component");
1173ca94c3ddSJeremy L Thompson   CeedCheck(num_nodes > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 node");
1174ca94c3ddSJeremy L Thompson   CeedCheck(num_qpts > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 quadrature point");
1175227444bfSJeremy L Thompson 
1176c4e3f59bSSebastian Grimberg   CeedCall(CeedBasisGetTopologyDimension(topo, &dim));
1177c4e3f59bSSebastian Grimberg 
1178db002c03SJeremy L Thompson   CeedCall(CeedCalloc(1, basis));
1179db002c03SJeremy L Thompson   CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed));
118050c301a5SRezgar Shakeri   (*basis)->ref_count       = 1;
11816402da51SJeremy L Thompson   (*basis)->is_tensor_basis = false;
118250c301a5SRezgar Shakeri   (*basis)->dim             = dim;
118350c301a5SRezgar Shakeri   (*basis)->topo            = topo;
118450c301a5SRezgar Shakeri   (*basis)->num_comp        = num_comp;
118550c301a5SRezgar Shakeri   (*basis)->P               = P;
118650c301a5SRezgar Shakeri   (*basis)->Q               = Q;
1187c4e3f59bSSebastian Grimberg   (*basis)->fe_space        = CEED_FE_SPACE_HDIV;
11882b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q * dim, &(*basis)->q_ref_1d));
11892b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q, &(*basis)->q_weight_1d));
119050c301a5SRezgar Shakeri   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0]));
119150c301a5SRezgar Shakeri   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0]));
11922b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(dim * Q * P, &(*basis)->interp));
11932b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q * P, &(*basis)->div));
119450c301a5SRezgar Shakeri   if (interp) memcpy((*basis)->interp, interp, dim * Q * P * sizeof(interp[0]));
119550c301a5SRezgar Shakeri   if (div) memcpy((*basis)->div, div, Q * P * sizeof(div[0]));
11962b730f8bSJeremy L Thompson   CeedCall(ceed->BasisCreateHdiv(topo, dim, P, Q, interp, div, q_ref, q_weight, *basis));
119750c301a5SRezgar Shakeri   return CEED_ERROR_SUCCESS;
119850c301a5SRezgar Shakeri }
119950c301a5SRezgar Shakeri 
120050c301a5SRezgar Shakeri /**
12014385fb7fSSebastian Grimberg   @brief Create a non tensor-product basis for \f$H(\mathrm{curl})\f$ discretizations
1202c4e3f59bSSebastian Grimberg 
1203ca94c3ddSJeremy L Thompson   @param[in]  ceed      `Ceed` object used to create the `CeedBasis`
1204c4e3f59bSSebastian Grimberg   @param[in]  topo      Topology of element (`CEED_TOPOLOGY_QUAD`, `CEED_TOPOLOGY_PRISM`, etc.), dimension of which is used in some array sizes below
1205ca94c3ddSJeremy L Thompson   @param[in]  num_comp  Number of components (usually 1 for vectors in \f$H(\mathrm{curl})\f$ bases)
1206ca94c3ddSJeremy L Thompson   @param[in]  num_nodes Total number of nodes (DoFs per element)
1207c4e3f59bSSebastian Grimberg   @param[in]  num_qpts  Total number of quadrature points
1208ca94c3ddSJeremy L Thompson   @param[in]  interp    Row-major (`dim * num_qpts * num_nodes`) matrix expressing the values of basis functions at quadrature points
1209ca94c3ddSJeremy 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
1210ca94c3ddSJeremy L Thompson   @param[in]  q_ref     Array of length `num_qpts * dim` holding the locations of quadrature points on the reference element
1211ca94c3ddSJeremy L Thompson   @param[in]  q_weight  Array of length `num_qpts` holding the quadrature weights on the reference element
1212ca94c3ddSJeremy L Thompson   @param[out] basis     Address of the variable where the newly created `CeedBasis` will be stored
1213c4e3f59bSSebastian Grimberg 
1214c4e3f59bSSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
1215c4e3f59bSSebastian Grimberg 
1216c4e3f59bSSebastian Grimberg   @ref User
1217c4e3f59bSSebastian Grimberg **/
1218c4e3f59bSSebastian Grimberg int CeedBasisCreateHcurl(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
1219c4e3f59bSSebastian Grimberg                          const CeedScalar *curl, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) {
1220c4e3f59bSSebastian Grimberg   CeedInt Q = num_qpts, P = num_nodes, dim = 0, curl_comp = 0;
1221c4e3f59bSSebastian Grimberg 
1222d075f50bSSebastian Grimberg   if (!ceed->BasisCreateHcurl) {
1223c4e3f59bSSebastian Grimberg     Ceed delegate;
12246574a04fSJeremy L Thompson 
1225c4e3f59bSSebastian Grimberg     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
12266574a04fSJeremy L Thompson     CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateHcurl");
1227c4e3f59bSSebastian Grimberg     CeedCall(CeedBasisCreateHcurl(delegate, topo, num_comp, num_nodes, num_qpts, interp, curl, q_ref, q_weight, basis));
1228c4e3f59bSSebastian Grimberg     return CEED_ERROR_SUCCESS;
1229c4e3f59bSSebastian Grimberg   }
1230c4e3f59bSSebastian Grimberg 
1231ca94c3ddSJeremy L Thompson   CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 component");
1232ca94c3ddSJeremy L Thompson   CeedCheck(num_nodes > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 node");
1233ca94c3ddSJeremy L Thompson   CeedCheck(num_qpts > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 quadrature point");
1234c4e3f59bSSebastian Grimberg 
1235c4e3f59bSSebastian Grimberg   CeedCall(CeedBasisGetTopologyDimension(topo, &dim));
1236c4e3f59bSSebastian Grimberg   curl_comp = (dim < 3) ? 1 : dim;
1237c4e3f59bSSebastian Grimberg 
1238db002c03SJeremy L Thompson   CeedCall(CeedCalloc(1, basis));
1239db002c03SJeremy L Thompson   CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed));
1240c4e3f59bSSebastian Grimberg   (*basis)->ref_count       = 1;
12416402da51SJeremy L Thompson   (*basis)->is_tensor_basis = false;
1242c4e3f59bSSebastian Grimberg   (*basis)->dim             = dim;
1243c4e3f59bSSebastian Grimberg   (*basis)->topo            = topo;
1244c4e3f59bSSebastian Grimberg   (*basis)->num_comp        = num_comp;
1245c4e3f59bSSebastian Grimberg   (*basis)->P               = P;
1246c4e3f59bSSebastian Grimberg   (*basis)->Q               = Q;
1247c4e3f59bSSebastian Grimberg   (*basis)->fe_space        = CEED_FE_SPACE_HCURL;
1248c4e3f59bSSebastian Grimberg   CeedCall(CeedMalloc(Q * dim, &(*basis)->q_ref_1d));
1249c4e3f59bSSebastian Grimberg   CeedCall(CeedMalloc(Q, &(*basis)->q_weight_1d));
1250c4e3f59bSSebastian Grimberg   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0]));
1251c4e3f59bSSebastian Grimberg   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0]));
1252c4e3f59bSSebastian Grimberg   CeedCall(CeedMalloc(dim * Q * P, &(*basis)->interp));
1253c4e3f59bSSebastian Grimberg   CeedCall(CeedMalloc(curl_comp * Q * P, &(*basis)->curl));
1254c4e3f59bSSebastian Grimberg   if (interp) memcpy((*basis)->interp, interp, dim * Q * P * sizeof(interp[0]));
1255c4e3f59bSSebastian Grimberg   if (curl) memcpy((*basis)->curl, curl, curl_comp * Q * P * sizeof(curl[0]));
1256c4e3f59bSSebastian Grimberg   CeedCall(ceed->BasisCreateHcurl(topo, dim, P, Q, interp, curl, q_ref, q_weight, *basis));
1257c4e3f59bSSebastian Grimberg   return CEED_ERROR_SUCCESS;
1258c4e3f59bSSebastian Grimberg }
1259c4e3f59bSSebastian Grimberg 
1260c4e3f59bSSebastian Grimberg /**
1261ca94c3ddSJeremy L Thompson   @brief Create a `CeedBasis` for projection from the nodes of `basis_from` to the nodes of `basis_to`.
1262ba59ac12SSebastian Grimberg 
1263ca94c3ddSJeremy L Thompson   Only @ref CEED_EVAL_INTERP will be valid for the new basis, `basis_project`.
1264ca94c3ddSJeremy L Thompson   For \f$H^1\f$ spaces, @ref CEED_EVAL_GRAD will also be valid.
1265ca94c3ddSJeremy L Thompson   The interpolation is given by `interp_project = interp_to^+ * interp_from`, where the pseudoinverse `interp_to^+` is given by QR factorization.
1266ca94c3ddSJeremy L Thompson   The gradient (for the \f$H^1\f$ case) is given by `grad_project = interp_to^+ * grad_from`.
126715ad3917SSebastian Grimberg 
126815ad3917SSebastian Grimberg   Note: `basis_from` and `basis_to` must have compatible quadrature spaces.
126915ad3917SSebastian Grimberg 
12709fd66db6SSebastian Grimberg   Note: `basis_project` will have the same number of components as `basis_from`, regardless of the number of components that `basis_to` has.
12719fd66db6SSebastian Grimberg         If `basis_from` has 3 components and `basis_to` has 5 components, then `basis_project` will have 3 components.
1272f113e5dcSJeremy L Thompson 
1273ca94c3ddSJeremy L Thompson   @param[in]  basis_from    `CeedBasis` to prolong from
1274ca94c3ddSJeremy L Thompson   @param[in]  basis_to      `CeedBasis` to prolong to
1275ca94c3ddSJeremy L Thompson   @param[out] basis_project Address of the variable where the newly created `CeedBasis` will be stored
1276f113e5dcSJeremy L Thompson 
1277f113e5dcSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1278f113e5dcSJeremy L Thompson 
1279f113e5dcSJeremy L Thompson   @ref User
1280f113e5dcSJeremy L Thompson **/
12812b730f8bSJeremy L Thompson int CeedBasisCreateProjection(CeedBasis basis_from, CeedBasis basis_to, CeedBasis *basis_project) {
1282f113e5dcSJeremy L Thompson   Ceed        ceed;
12831c66c397SJeremy L Thompson   bool        is_tensor;
12841c66c397SJeremy L Thompson   CeedInt     dim, num_comp;
12851c66c397SJeremy L Thompson   CeedScalar *q_ref, *q_weight, *interp_project, *grad_project;
12861c66c397SJeremy L Thompson 
12872b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetCeed(basis_to, &ceed));
1288f113e5dcSJeremy L Thompson 
1289ecc88aebSJeremy L Thompson   // Create projection matrix
12902b730f8bSJeremy L Thompson   CeedCall(CeedBasisCreateProjectionMatrices(basis_from, basis_to, &interp_project, &grad_project));
1291f113e5dcSJeremy L Thompson 
1292f113e5dcSJeremy L Thompson   // Build basis
12932b730f8bSJeremy L Thompson   CeedCall(CeedBasisIsTensor(basis_to, &is_tensor));
12942b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetDimension(basis_to, &dim));
12952b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumComponents(basis_from, &num_comp));
1296f113e5dcSJeremy L Thompson   if (is_tensor) {
1297f113e5dcSJeremy L Thompson     CeedInt P_1d_to, P_1d_from;
12981c66c397SJeremy L Thompson 
12992b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_1d_from));
13002b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_1d_to));
13012b730f8bSJeremy L Thompson     CeedCall(CeedCalloc(P_1d_to, &q_ref));
13022b730f8bSJeremy L Thompson     CeedCall(CeedCalloc(P_1d_to, &q_weight));
13032b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P_1d_from, P_1d_to, interp_project, grad_project, q_ref, q_weight, basis_project));
1304f113e5dcSJeremy L Thompson   } else {
1305de05fbb2SSebastian Grimberg     // Even if basis_to and basis_from are not H1, the resulting basis is H1 for interpolation to work
1306f113e5dcSJeremy L Thompson     CeedInt          num_nodes_to, num_nodes_from;
13071c66c397SJeremy L Thompson     CeedElemTopology topo;
13081c66c397SJeremy L Thompson 
13091c66c397SJeremy L Thompson     CeedCall(CeedBasisGetTopology(basis_to, &topo));
13102b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis_from, &num_nodes_from));
13112b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis_to, &num_nodes_to));
13122b730f8bSJeremy L Thompson     CeedCall(CeedCalloc(num_nodes_to * dim, &q_ref));
13132b730f8bSJeremy L Thompson     CeedCall(CeedCalloc(num_nodes_to, &q_weight));
13142b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateH1(ceed, topo, num_comp, num_nodes_from, num_nodes_to, interp_project, grad_project, q_ref, q_weight, basis_project));
1315f113e5dcSJeremy L Thompson   }
1316f113e5dcSJeremy L Thompson 
1317f113e5dcSJeremy L Thompson   // Cleanup
13182b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_project));
13192b730f8bSJeremy L Thompson   CeedCall(CeedFree(&grad_project));
13202b730f8bSJeremy L Thompson   CeedCall(CeedFree(&q_ref));
13212b730f8bSJeremy L Thompson   CeedCall(CeedFree(&q_weight));
1322f113e5dcSJeremy L Thompson   return CEED_ERROR_SUCCESS;
1323f113e5dcSJeremy L Thompson }
1324f113e5dcSJeremy L Thompson 
1325f113e5dcSJeremy L Thompson /**
1326ca94c3ddSJeremy L Thompson   @brief Copy the pointer to a `CeedBasis`.
13279560d06aSjeremylt 
1328ca94c3ddSJeremy 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`.
1329ca94c3ddSJeremy L Thompson         This `CeedBasis` will be destroyed if `*basis_copy` is the only reference to this `CeedBasis`.
1330ea61e9acSJeremy L Thompson 
1331ca94c3ddSJeremy L Thompson   @param[in]     basis      `CeedBasis` to copy reference to
1332ea61e9acSJeremy L Thompson   @param[in,out] basis_copy Variable to store copied reference
13339560d06aSjeremylt 
13349560d06aSjeremylt   @return An error code: 0 - success, otherwise - failure
13359560d06aSjeremylt 
13369560d06aSjeremylt   @ref User
13379560d06aSjeremylt **/
13389560d06aSjeremylt int CeedBasisReferenceCopy(CeedBasis basis, CeedBasis *basis_copy) {
1339356036faSJeremy L Thompson   if (basis != CEED_BASIS_NONE) CeedCall(CeedBasisReference(basis));
13402b730f8bSJeremy L Thompson   CeedCall(CeedBasisDestroy(basis_copy));
13419560d06aSjeremylt   *basis_copy = basis;
13429560d06aSjeremylt   return CEED_ERROR_SUCCESS;
13439560d06aSjeremylt }
13449560d06aSjeremylt 
13459560d06aSjeremylt /**
1346ca94c3ddSJeremy L Thompson   @brief View a `CeedBasis`
13477a982d89SJeremy L. Thompson 
1348ca94c3ddSJeremy L Thompson   @param[in] basis  `CeedBasis` to view
1349ca94c3ddSJeremy L Thompson   @param[in] stream Stream to view to, e.g., `stdout`
13507a982d89SJeremy L. Thompson 
13517a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
13527a982d89SJeremy L. Thompson 
13537a982d89SJeremy L. Thompson   @ref User
13547a982d89SJeremy L. Thompson **/
13557a982d89SJeremy L. Thompson int CeedBasisView(CeedBasis basis, FILE *stream) {
13561c66c397SJeremy L Thompson   CeedInt          q_comp   = 0;
135750c301a5SRezgar Shakeri   CeedElemTopology topo     = basis->topo;
1358c4e3f59bSSebastian Grimberg   CeedFESpace      fe_space = basis->fe_space;
13592b730f8bSJeremy L Thompson 
136050c301a5SRezgar Shakeri   // Print FE space and element topology of the basis
1361edf04919SJeremy L Thompson   fprintf(stream, "CeedBasis in a %s on a %s element\n", CeedFESpaces[fe_space], CeedElemTopologies[topo]);
13626402da51SJeremy L Thompson   if (basis->is_tensor_basis) {
1363edf04919SJeremy L Thompson     fprintf(stream, "  P: %" CeedInt_FMT "\n  Q: %" CeedInt_FMT "\n", basis->P_1d, basis->Q_1d);
136450c301a5SRezgar Shakeri   } else {
1365edf04919SJeremy L Thompson     fprintf(stream, "  P: %" CeedInt_FMT "\n  Q: %" CeedInt_FMT "\n", basis->P, basis->Q);
136650c301a5SRezgar Shakeri   }
1367edf04919SJeremy L Thompson   fprintf(stream, "  dimension: %" CeedInt_FMT "\n  field components: %" CeedInt_FMT "\n", basis->dim, basis->num_comp);
1368ea61e9acSJeremy L Thompson   // Print quadrature data, interpolation/gradient/divergence/curl of the basis
13696402da51SJeremy L Thompson   if (basis->is_tensor_basis) {  // tensor basis
13702b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_ref_1d, stream));
13712b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_weight_1d, stream));
13722b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("interp1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, basis->interp_1d, stream));
13732b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("grad1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, basis->grad_1d, stream));
137450c301a5SRezgar Shakeri   } else {  // non-tensor basis
13752b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("qref", "\t% 12.8f", 1, basis->Q * basis->dim, basis->q_ref_1d, stream));
13762b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->q_weight_1d, stream));
1377c4e3f59bSSebastian Grimberg     CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp));
1378c4e3f59bSSebastian Grimberg     CeedCall(CeedScalarView("interp", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->interp, stream));
137950c301a5SRezgar Shakeri     if (basis->grad) {
1380c4e3f59bSSebastian Grimberg       CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp));
1381c4e3f59bSSebastian Grimberg       CeedCall(CeedScalarView("grad", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->grad, stream));
13827a982d89SJeremy L. Thompson     }
138350c301a5SRezgar Shakeri     if (basis->div) {
1384c4e3f59bSSebastian Grimberg       CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_DIV, &q_comp));
1385c4e3f59bSSebastian Grimberg       CeedCall(CeedScalarView("div", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->div, stream));
1386c4e3f59bSSebastian Grimberg     }
1387c4e3f59bSSebastian Grimberg     if (basis->curl) {
1388c4e3f59bSSebastian Grimberg       CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_CURL, &q_comp));
1389c4e3f59bSSebastian Grimberg       CeedCall(CeedScalarView("curl", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->curl, stream));
139050c301a5SRezgar Shakeri     }
139150c301a5SRezgar Shakeri   }
1392e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
13937a982d89SJeremy L. Thompson }
13947a982d89SJeremy L. Thompson 
13957a982d89SJeremy L. Thompson /**
13967a982d89SJeremy L. Thompson   @brief Apply basis evaluation from nodes to quadrature points or vice versa
13977a982d89SJeremy L. Thompson 
1398ca94c3ddSJeremy L Thompson   @param[in]  basis     `CeedBasis` to evaluate
1399ea61e9acSJeremy L Thompson   @param[in]  num_elem  The number of elements to apply the basis evaluation to;
1400ca94c3ddSJeremy L Thompson                           the backend will specify the ordering in @ref CeedElemRestrictionCreate()
1401ca94c3ddSJeremy L Thompson   @param[in]  t_mode    @ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature points;
1402ca94c3ddSJeremy L Thompson                           @ref CEED_TRANSPOSE to apply the transpose, mapping from quadrature points to nodes
1403ca94c3ddSJeremy L Thompson   @param[in]  eval_mode @ref CEED_EVAL_NONE to use values directly,
1404ca94c3ddSJeremy L Thompson                           @ref CEED_EVAL_INTERP to use interpolated values,
1405ca94c3ddSJeremy L Thompson                           @ref CEED_EVAL_GRAD to use gradients,
1406ca94c3ddSJeremy L Thompson                           @ref CEED_EVAL_DIV to use divergence,
1407ca94c3ddSJeremy L Thompson                           @ref CEED_EVAL_CURL to use curl,
1408ca94c3ddSJeremy L Thompson                           @ref CEED_EVAL_WEIGHT to use quadrature weights
1409ca94c3ddSJeremy L Thompson   @param[in]  u         Input `CeedVector`
1410ca94c3ddSJeremy L Thompson   @param[out] v         Output `CeedVector`
14117a982d89SJeremy L. Thompson 
14127a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
14137a982d89SJeremy L. Thompson 
14147a982d89SJeremy L. Thompson   @ref User
14157a982d89SJeremy L. Thompson **/
14162b730f8bSJeremy L Thompson int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
1417c4e3f59bSSebastian Grimberg   CeedInt  dim, num_comp, q_comp, num_nodes, num_qpts;
14181c66c397SJeremy L Thompson   CeedSize u_length = 0, v_length;
14191c66c397SJeremy L Thompson 
14202b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetDimension(basis, &dim));
14212b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
1422c4e3f59bSSebastian Grimberg   CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp));
14232b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumNodes(basis, &num_nodes));
14242b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts));
14252b730f8bSJeremy L Thompson   CeedCall(CeedVectorGetLength(v, &v_length));
1426c8c3fa7dSJeremy L Thompson   if (u) CeedCall(CeedVectorGetLength(u, &u_length));
14277a982d89SJeremy L. Thompson 
1428ca94c3ddSJeremy L Thompson   CeedCheck(basis->Apply, basis->ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support CeedBasisApply");
1429e15f9bd0SJeremy L Thompson 
1430e15f9bd0SJeremy L Thompson   // Check compatibility of topological and geometrical dimensions
14316574a04fSJeremy L Thompson   CeedCheck((t_mode == CEED_TRANSPOSE && v_length % num_nodes == 0 && u_length % num_qpts == 0) ||
14326574a04fSJeremy L Thompson                 (t_mode == CEED_NOTRANSPOSE && u_length % num_nodes == 0 && v_length % num_qpts == 0),
14336574a04fSJeremy L Thompson             basis->ceed, CEED_ERROR_DIMENSION, "Length of input/output vectors incompatible with basis dimensions");
14347a982d89SJeremy L. Thompson 
1435e15f9bd0SJeremy L Thompson   // Check vector lengths to prevent out of bounds issues
14366574a04fSJeremy L Thompson   bool good_dims = true;
1437d1d35e2fSjeremylt   switch (eval_mode) {
1438e15f9bd0SJeremy L Thompson     case CEED_EVAL_NONE:
14392b730f8bSJeremy L Thompson     case CEED_EVAL_INTERP:
14402b730f8bSJeremy L Thompson     case CEED_EVAL_GRAD:
1441c4e3f59bSSebastian Grimberg     case CEED_EVAL_DIV:
1442c4e3f59bSSebastian Grimberg     case CEED_EVAL_CURL:
14436574a04fSJeremy L Thompson       good_dims =
14446574a04fSJeremy L Thompson           ((t_mode == CEED_TRANSPOSE && u_length >= num_elem * num_comp * num_qpts * q_comp && v_length >= num_elem * num_comp * num_nodes) ||
14456574a04fSJeremy L Thompson            (t_mode == CEED_NOTRANSPOSE && v_length >= num_elem * num_qpts * num_comp * q_comp && u_length >= num_elem * num_comp * num_nodes));
1446e15f9bd0SJeremy L Thompson       break;
1447e15f9bd0SJeremy L Thompson     case CEED_EVAL_WEIGHT:
14486574a04fSJeremy L Thompson       good_dims = v_length >= num_elem * num_qpts;
1449e15f9bd0SJeremy L Thompson       break;
1450e15f9bd0SJeremy L Thompson   }
14516574a04fSJeremy L Thompson   CeedCheck(good_dims, basis->ceed, CEED_ERROR_DIMENSION, "Input/output vectors too short for basis and evaluation mode");
1452e15f9bd0SJeremy L Thompson 
14532b730f8bSJeremy L Thompson   CeedCall(basis->Apply(basis, num_elem, t_mode, eval_mode, u, v));
1454e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
14557a982d89SJeremy L. Thompson }
14567a982d89SJeremy L. Thompson 
14577a982d89SJeremy L. Thompson /**
1458c8c3fa7dSJeremy L Thompson   @brief Apply basis evaluation from nodes to arbitrary points
1459c8c3fa7dSJeremy L Thompson 
1460ca94c3ddSJeremy L Thompson   @param[in]  basis      `CeedBasis` to evaluate
1461c8c3fa7dSJeremy L Thompson   @param[in]  num_points The number of points to apply the basis evaluation to
1462ca94c3ddSJeremy L Thompson   @param[in]  t_mode     @ref CEED_NOTRANSPOSE to evaluate from nodes to points;
1463ca94c3ddSJeremy L Thompson                            @ref CEED_TRANSPOSE to apply the transpose, mapping from points to nodes
1464ca94c3ddSJeremy L Thompson   @param[in]  eval_mode  @ref CEED_EVAL_INTERP to use interpolated values,
1465ca94c3ddSJeremy L Thompson                            @ref CEED_EVAL_GRAD to use gradients,
1466ca94c3ddSJeremy L Thompson                            @ref CEED_EVAL_WEIGHT to use quadrature weights
1467ca94c3ddSJeremy L Thompson   @param[in]  x_ref      `CeedVector` holding reference coordinates of each point
1468ca94c3ddSJeremy L Thompson   @param[in]  u          Input `CeedVector`, of length `num_nodes * num_comp` for @ref CEED_NOTRANSPOSE
1469ca94c3ddSJeremy L Thompson   @param[out] v          Output `CeedVector`, of length `num_points * num_q_comp` for @ref CEED_NOTRANSPOSE with @ref CEED_EVAL_INTERP
1470c8c3fa7dSJeremy L Thompson 
1471c8c3fa7dSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1472c8c3fa7dSJeremy L Thompson 
1473c8c3fa7dSJeremy L Thompson   @ref User
1474c8c3fa7dSJeremy L Thompson **/
1475c8c3fa7dSJeremy L Thompson int CeedBasisApplyAtPoints(CeedBasis basis, CeedInt num_points, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u,
1476c8c3fa7dSJeremy L Thompson                            CeedVector v) {
1477c8c3fa7dSJeremy L Thompson   CeedInt  dim, num_comp, num_q_comp, num_nodes, P_1d = 1, Q_1d = 1;
14781c66c397SJeremy L Thompson   CeedSize x_length = 0, u_length = 0, v_length;
1479c8c3fa7dSJeremy L Thompson 
1480c8c3fa7dSJeremy L Thompson   CeedCall(CeedBasisGetDimension(basis, &dim));
1481c8c3fa7dSJeremy L Thompson   CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d));
1482c8c3fa7dSJeremy L Thompson   CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
1483c8c3fa7dSJeremy L Thompson   CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
1484c8c3fa7dSJeremy L Thompson   CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &num_q_comp));
1485c8c3fa7dSJeremy L Thompson   CeedCall(CeedBasisGetNumNodes(basis, &num_nodes));
1486c8c3fa7dSJeremy L Thompson   CeedCall(CeedVectorGetLength(v, &v_length));
1487953190f4SJeremy L Thompson   if (x_ref != CEED_VECTOR_NONE) CeedCall(CeedVectorGetLength(x_ref, &x_length));
1488953190f4SJeremy L Thompson   if (u != CEED_VECTOR_NONE) CeedCall(CeedVectorGetLength(u, &u_length));
1489c8c3fa7dSJeremy L Thompson 
1490c8c3fa7dSJeremy L Thompson   // Check compatibility of topological and geometrical dimensions
1491953190f4SJeremy L Thompson   CeedCheck((t_mode == CEED_TRANSPOSE && v_length % num_nodes == 0) || (t_mode == CEED_NOTRANSPOSE && u_length % num_nodes == 0) ||
1492953190f4SJeremy L Thompson                 (eval_mode == CEED_EVAL_WEIGHT),
1493953190f4SJeremy L Thompson             basis->ceed, CEED_ERROR_DIMENSION, "Length of input/output vectors incompatible with basis dimensions and number of points");
1494c8c3fa7dSJeremy L Thompson 
1495c8c3fa7dSJeremy L Thompson   // Check compatibility coordinates vector
1496953190f4SJeremy L Thompson   CeedCheck((x_length >= num_points * dim) || (eval_mode == CEED_EVAL_WEIGHT), basis->ceed, CEED_ERROR_DIMENSION,
1497c8c3fa7dSJeremy L Thompson             "Length of reference coordinate vector incompatible with basis dimension and number of points");
1498c8c3fa7dSJeremy L Thompson 
1499953190f4SJeremy L Thompson   // Check CEED_EVAL_WEIGHT only on CEED_NOTRANSPOSE
1500953190f4SJeremy L Thompson   CeedCheck(eval_mode != CEED_EVAL_WEIGHT || t_mode == CEED_NOTRANSPOSE, basis->ceed, CEED_ERROR_UNSUPPORTED,
1501953190f4SJeremy L Thompson             "CEED_EVAL_WEIGHT only supported with CEED_NOTRANSPOSE");
1502953190f4SJeremy L Thompson 
1503c8c3fa7dSJeremy L Thompson   // Check vector lengths to prevent out of bounds issues
1504c8c3fa7dSJeremy L Thompson   bool good_dims = false;
1505c8c3fa7dSJeremy L Thompson   switch (eval_mode) {
1506c8c3fa7dSJeremy L Thompson     case CEED_EVAL_INTERP:
1507c8c3fa7dSJeremy L Thompson       good_dims = ((t_mode == CEED_TRANSPOSE && (u_length >= num_points * num_q_comp || v_length >= num_nodes * num_comp)) ||
1508c8c3fa7dSJeremy L Thompson                    (t_mode == CEED_NOTRANSPOSE && (v_length >= num_points * num_q_comp || u_length >= num_nodes * num_comp)));
1509c8c3fa7dSJeremy L Thompson       break;
1510c8c3fa7dSJeremy L Thompson     case CEED_EVAL_GRAD:
1511edfbf3a6SJeremy L Thompson       good_dims = ((t_mode == CEED_TRANSPOSE && (u_length >= num_points * num_q_comp * dim || v_length >= num_nodes * num_comp)) ||
1512edfbf3a6SJeremy L Thompson                    (t_mode == CEED_NOTRANSPOSE && (v_length >= num_points * num_q_comp * dim || u_length >= num_nodes * num_comp)));
1513edfbf3a6SJeremy L Thompson       break;
1514c8c3fa7dSJeremy L Thompson     case CEED_EVAL_WEIGHT:
1515953190f4SJeremy L Thompson       good_dims = t_mode == CEED_NOTRANSPOSE && (v_length >= num_points);
1516953190f4SJeremy L Thompson       break;
1517953190f4SJeremy L Thompson     case CEED_EVAL_NONE:
1518c8c3fa7dSJeremy L Thompson     case CEED_EVAL_DIV:
1519c8c3fa7dSJeremy L Thompson     case CEED_EVAL_CURL:
1520c8c3fa7dSJeremy L Thompson       // LCOV_EXCL_START
1521c8c3fa7dSJeremy L Thompson       return CeedError(basis->ceed, CEED_ERROR_UNSUPPORTED, "Evaluation at arbitrary points not supported for %s", CeedEvalModes[eval_mode]);
1522c8c3fa7dSJeremy L Thompson       // LCOV_EXCL_STOP
1523c8c3fa7dSJeremy L Thompson   }
1524c8c3fa7dSJeremy L Thompson   CeedCheck(good_dims, basis->ceed, CEED_ERROR_DIMENSION, "Input/output vectors too short for basis and evaluation mode");
1525c8c3fa7dSJeremy L Thompson 
1526c8c3fa7dSJeremy L Thompson   // Backend method
1527c8c3fa7dSJeremy L Thompson   if (basis->ApplyAtPoints) {
1528c8c3fa7dSJeremy L Thompson     CeedCall(basis->ApplyAtPoints(basis, num_points, t_mode, eval_mode, x_ref, u, v));
1529c8c3fa7dSJeremy L Thompson     return CEED_ERROR_SUCCESS;
1530c8c3fa7dSJeremy L Thompson   }
1531c8c3fa7dSJeremy L Thompson 
1532c8c3fa7dSJeremy L Thompson   // Default implementation
1533c8c3fa7dSJeremy L Thompson   CeedCheck(basis->is_tensor_basis, basis->ceed, CEED_ERROR_UNSUPPORTED, "Evaluation at arbitrary points only supported for tensor product bases");
1534953190f4SJeremy L Thompson   if (eval_mode == CEED_EVAL_WEIGHT) {
1535953190f4SJeremy L Thompson     CeedCall(CeedVectorSetValue(v, 1.0));
1536953190f4SJeremy L Thompson     return CEED_ERROR_SUCCESS;
1537953190f4SJeremy L Thompson   }
1538c8c3fa7dSJeremy L Thompson   if (!basis->basis_chebyshev) {
1539c8c3fa7dSJeremy L Thompson     // Build matrix mapping from quadrature point values to Chebyshev coefficients
1540*2247a93fSRezgar Shakeri     CeedScalar       *C, *chebyshev_coeffs_1d_inv;
1541c8c3fa7dSJeremy L Thompson     const CeedScalar *q_ref_1d;
1542c8c3fa7dSJeremy L Thompson 
1543c8c3fa7dSJeremy L Thompson     // Build coefficient matrix
1544c8c3fa7dSJeremy L Thompson     // -- Note: Clang-tidy needs this check because it does not understand the is_tensor_basis check above
1545ca94c3ddSJeremy L Thompson     CeedCheck(P_1d > 0 && Q_1d > 0, basis->ceed, CEED_ERROR_INCOMPATIBLE, "CeedBasis dimensions are malformed");
1546c8c3fa7dSJeremy L Thompson     CeedCall(CeedCalloc(Q_1d * Q_1d, &C));
1547c8c3fa7dSJeremy L Thompson     CeedCall(CeedBasisGetQRef(basis, &q_ref_1d));
15483778dbaaSJeremy L Thompson     for (CeedInt i = 0; i < Q_1d; i++) CeedCall(CeedChebyshevPolynomialsAtPoint(q_ref_1d[i], Q_1d, &C[i * Q_1d]));
1549c8c3fa7dSJeremy L Thompson 
1550*2247a93fSRezgar Shakeri     // Compute C^+, pseudoinverse of coefficient matrix
1551*2247a93fSRezgar Shakeri     CeedCall(CeedCalloc(Q_1d * Q_1d, &chebyshev_coeffs_1d_inv));
1552*2247a93fSRezgar Shakeri     CeedCall(CeedMatrixPseudoinverse(basis->ceed, C, Q_1d, Q_1d, chebyshev_coeffs_1d_inv));
1553c8c3fa7dSJeremy L Thompson 
1554c8c3fa7dSJeremy L Thompson     // Build basis mapping from nodes to Chebyshev coefficients
1555c8c3fa7dSJeremy L Thompson     CeedScalar       *chebyshev_interp_1d, *chebyshev_grad_1d, *chebyshev_q_weight_1d;
1556c8c3fa7dSJeremy L Thompson     const CeedScalar *interp_1d;
1557c8c3fa7dSJeremy L Thompson 
155871a83b88SJeremy L Thompson     CeedCall(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
155971a83b88SJeremy L Thompson     CeedCall(CeedCalloc(P_1d * Q_1d, &chebyshev_grad_1d));
1560c8c3fa7dSJeremy L Thompson     CeedCall(CeedCalloc(Q_1d, &chebyshev_q_weight_1d));
1561c8c3fa7dSJeremy L Thompson     CeedCall(CeedBasisGetInterp1D(basis, &interp_1d));
1562*2247a93fSRezgar Shakeri     CeedCall(CeedMatrixMatrixMultiply(basis->ceed, chebyshev_coeffs_1d_inv, interp_1d, chebyshev_interp_1d, Q_1d, P_1d, Q_1d));
1563c8c3fa7dSJeremy L Thompson 
1564c8c3fa7dSJeremy L Thompson     CeedCall(CeedVectorCreate(basis->ceed, num_comp * CeedIntPow(Q_1d, dim), &basis->vec_chebyshev));
156571a83b88SJeremy L Thompson     CeedCall(CeedBasisCreateTensorH1(basis->ceed, dim, num_comp, P_1d, Q_1d, chebyshev_interp_1d, chebyshev_grad_1d, q_ref_1d, chebyshev_q_weight_1d,
1566c8c3fa7dSJeremy L Thompson                                      &basis->basis_chebyshev));
1567c8c3fa7dSJeremy L Thompson 
1568c8c3fa7dSJeremy L Thompson     // Cleanup
1569c8c3fa7dSJeremy L Thompson     CeedCall(CeedFree(&C));
1570*2247a93fSRezgar Shakeri     CeedCall(CeedFree(&chebyshev_coeffs_1d_inv));
1571c8c3fa7dSJeremy L Thompson     CeedCall(CeedFree(&chebyshev_interp_1d));
1572c8c3fa7dSJeremy L Thompson     CeedCall(CeedFree(&chebyshev_grad_1d));
1573c8c3fa7dSJeremy L Thompson     CeedCall(CeedFree(&chebyshev_q_weight_1d));
1574c8c3fa7dSJeremy L Thompson   }
1575c8c3fa7dSJeremy L Thompson 
1576c8c3fa7dSJeremy L Thompson   // Create TensorContract object if needed, such as a basis from the GPU backends
1577c8c3fa7dSJeremy L Thompson   if (!basis->contract) {
1578c8c3fa7dSJeremy L Thompson     Ceed      ceed_ref;
1579585a562dSJeremy L Thompson     CeedBasis basis_ref = NULL;
1580c8c3fa7dSJeremy L Thompson 
1581c8c3fa7dSJeremy L Thompson     CeedCall(CeedInit("/cpu/self", &ceed_ref));
1582c8c3fa7dSJeremy L Thompson     // Only need matching tensor contraction dimensions, any type of basis will work
158371a83b88SJeremy L Thompson     CeedCall(CeedBasisCreateTensorH1Lagrange(ceed_ref, dim, num_comp, P_1d, Q_1d, CEED_GAUSS, &basis_ref));
1584585a562dSJeremy L Thompson     // Note - clang-tidy doesn't know basis_ref->contract must be valid here
1585585a562dSJeremy L Thompson     CeedCheck(basis_ref && basis_ref->contract, basis->ceed, CEED_ERROR_UNSUPPORTED,
15861c66c397SJeremy L Thompson               "Reference CPU ceed failed to create a tensor contraction object");
1587585a562dSJeremy L Thompson     CeedCall(CeedTensorContractReferenceCopy(basis_ref->contract, &basis->contract));
1588c8c3fa7dSJeremy L Thompson     CeedCall(CeedBasisDestroy(&basis_ref));
1589c8c3fa7dSJeremy L Thompson     CeedCall(CeedDestroy(&ceed_ref));
1590c8c3fa7dSJeremy L Thompson   }
1591c8c3fa7dSJeremy L Thompson 
1592c8c3fa7dSJeremy L Thompson   // Basis evaluation
1593c8c3fa7dSJeremy L Thompson   switch (t_mode) {
1594c8c3fa7dSJeremy L Thompson     case CEED_NOTRANSPOSE: {
1595c8c3fa7dSJeremy L Thompson       // Nodes to arbitrary points
1596c8c3fa7dSJeremy L Thompson       CeedScalar       *v_array;
1597c8c3fa7dSJeremy L Thompson       const CeedScalar *chebyshev_coeffs, *x_array_read;
1598c8c3fa7dSJeremy L Thompson 
1599c8c3fa7dSJeremy L Thompson       // -- Interpolate to Chebyshev coefficients
1600c8c3fa7dSJeremy L Thompson       CeedCall(CeedBasisApply(basis->basis_chebyshev, 1, CEED_NOTRANSPOSE, CEED_EVAL_INTERP, u, basis->vec_chebyshev));
1601c8c3fa7dSJeremy L Thompson 
1602c8c3fa7dSJeremy L Thompson       // -- Evaluate Chebyshev polynomials at arbitrary points
1603c8c3fa7dSJeremy L Thompson       CeedCall(CeedVectorGetArrayRead(basis->vec_chebyshev, CEED_MEM_HOST, &chebyshev_coeffs));
1604c8c3fa7dSJeremy L Thompson       CeedCall(CeedVectorGetArrayRead(x_ref, CEED_MEM_HOST, &x_array_read));
1605c8c3fa7dSJeremy L Thompson       CeedCall(CeedVectorGetArrayWrite(v, CEED_MEM_HOST, &v_array));
1606edfbf3a6SJeremy L Thompson       switch (eval_mode) {
1607edfbf3a6SJeremy L Thompson         case CEED_EVAL_INTERP: {
1608c8c3fa7dSJeremy L Thompson           CeedScalar tmp[2][num_comp * CeedIntPow(Q_1d, dim)], chebyshev_x[Q_1d];
1609c8c3fa7dSJeremy L Thompson 
1610c8c3fa7dSJeremy L Thompson           // ---- Values at point
1611c8c3fa7dSJeremy L Thompson           for (CeedInt p = 0; p < num_points; p++) {
1612c8c3fa7dSJeremy L Thompson             CeedInt pre = num_comp * CeedIntPow(Q_1d, dim - 1), post = 1;
1613c8c3fa7dSJeremy L Thompson 
161453ef2869SZach Atkins             for (CeedInt d = 0; d < dim; d++) {
16153778dbaaSJeremy L Thompson               // ------ Tensor contract with current Chebyshev polynomial values
16169c34f28eSJeremy L Thompson               CeedCall(CeedChebyshevPolynomialsAtPoint(x_array_read[d * num_points + p], Q_1d, chebyshev_x));
1617c8c3fa7dSJeremy L Thompson               CeedCall(CeedTensorContractApply(basis->contract, pre, Q_1d, post, 1, chebyshev_x, t_mode, false,
16184608bdaaSJeremy L Thompson                                                d == 0 ? chebyshev_coeffs : tmp[d % 2], tmp[(d + 1) % 2]));
1619c8c3fa7dSJeremy L Thompson               pre /= Q_1d;
1620c8c3fa7dSJeremy L Thompson               post *= 1;
1621c8c3fa7dSJeremy L Thompson             }
16224608bdaaSJeremy L Thompson             for (CeedInt c = 0; c < num_comp; c++) v_array[c * num_points + p] = tmp[dim % 2][c];
1623c8c3fa7dSJeremy L Thompson           }
1624edfbf3a6SJeremy L Thompson           break;
1625edfbf3a6SJeremy L Thompson         }
1626edfbf3a6SJeremy L Thompson         case CEED_EVAL_GRAD: {
1627edfbf3a6SJeremy L Thompson           CeedScalar tmp[2][num_comp * CeedIntPow(Q_1d, dim)], chebyshev_x[Q_1d];
1628edfbf3a6SJeremy L Thompson 
1629edfbf3a6SJeremy L Thompson           // ---- Values at point
1630edfbf3a6SJeremy L Thompson           for (CeedInt p = 0; p < num_points; p++) {
1631edfbf3a6SJeremy L Thompson             // Dim**2 contractions, apply grad when pass == dim
163253ef2869SZach Atkins             for (CeedInt pass = 0; pass < dim; pass++) {
1633edfbf3a6SJeremy L Thompson               CeedInt pre = num_comp * CeedIntPow(Q_1d, dim - 1), post = 1;
1634edfbf3a6SJeremy L Thompson 
163553ef2869SZach Atkins               for (CeedInt d = 0; d < dim; d++) {
16363778dbaaSJeremy L Thompson                 // ------ Tensor contract with current Chebyshev polynomial values
16379c34f28eSJeremy L Thompson                 if (pass == d) CeedCall(CeedChebyshevDerivativeAtPoint(x_array_read[d * num_points + p], Q_1d, chebyshev_x));
16389c34f28eSJeremy L Thompson                 else CeedCall(CeedChebyshevPolynomialsAtPoint(x_array_read[d * num_points + p], Q_1d, chebyshev_x));
1639edfbf3a6SJeremy L Thompson                 CeedCall(CeedTensorContractApply(basis->contract, pre, Q_1d, post, 1, chebyshev_x, t_mode, false,
16404608bdaaSJeremy L Thompson                                                  d == 0 ? chebyshev_coeffs : tmp[d % 2], tmp[(d + 1) % 2]));
1641edfbf3a6SJeremy L Thompson                 pre /= Q_1d;
1642edfbf3a6SJeremy L Thompson                 post *= 1;
1643edfbf3a6SJeremy L Thompson               }
16444608bdaaSJeremy L Thompson               for (CeedInt c = 0; c < num_comp; c++) v_array[(pass * num_comp + c) * num_points + p] = tmp[dim % 2][c];
1645edfbf3a6SJeremy L Thompson             }
1646edfbf3a6SJeremy L Thompson           }
1647edfbf3a6SJeremy L Thompson           break;
1648edfbf3a6SJeremy L Thompson         }
1649edfbf3a6SJeremy L Thompson         default:
1650953190f4SJeremy L Thompson           // Nothing to do, excluded above
1651edfbf3a6SJeremy L Thompson           break;
1652c8c3fa7dSJeremy L Thompson       }
1653c8c3fa7dSJeremy L Thompson       CeedCall(CeedVectorRestoreArrayRead(basis->vec_chebyshev, &chebyshev_coeffs));
1654c8c3fa7dSJeremy L Thompson       CeedCall(CeedVectorRestoreArrayRead(x_ref, &x_array_read));
1655c8c3fa7dSJeremy L Thompson       CeedCall(CeedVectorRestoreArray(v, &v_array));
1656c8c3fa7dSJeremy L Thompson       break;
1657c8c3fa7dSJeremy L Thompson     }
16582a94f45fSJeremy L Thompson     case CEED_TRANSPOSE: {
16593778dbaaSJeremy L Thompson       // Note: No switch on e_mode here because only CEED_EVAL_INTERP is supported at this time
16602a94f45fSJeremy L Thompson       // Arbitrary points to nodes
16612a94f45fSJeremy L Thompson       CeedScalar       *chebyshev_coeffs;
16622a94f45fSJeremy L Thompson       const CeedScalar *u_array, *x_array_read;
16632a94f45fSJeremy L Thompson 
16641c66c397SJeremy L Thompson       // -- Transpose of evaluation of Chebyshev polynomials at arbitrary points
16652a94f45fSJeremy L Thompson       CeedCall(CeedVectorGetArrayWrite(basis->vec_chebyshev, CEED_MEM_HOST, &chebyshev_coeffs));
16662a94f45fSJeremy L Thompson       CeedCall(CeedVectorGetArrayRead(x_ref, CEED_MEM_HOST, &x_array_read));
16672a94f45fSJeremy L Thompson       CeedCall(CeedVectorGetArrayRead(u, CEED_MEM_HOST, &u_array));
1668038a8942SZach Atkins 
1669038a8942SZach Atkins       switch (eval_mode) {
1670038a8942SZach Atkins         case CEED_EVAL_INTERP: {
16712a94f45fSJeremy L Thompson           CeedScalar tmp[2][num_comp * CeedIntPow(Q_1d, dim)], chebyshev_x[Q_1d];
16722a94f45fSJeremy L Thompson 
16732a94f45fSJeremy L Thompson           // ---- Values at point
16742a94f45fSJeremy L Thompson           for (CeedInt p = 0; p < num_points; p++) {
16752a94f45fSJeremy L Thompson             CeedInt pre = num_comp * 1, post = 1;
16762a94f45fSJeremy L Thompson 
16774608bdaaSJeremy L Thompson             for (CeedInt c = 0; c < num_comp; c++) tmp[0][c] = u_array[c * num_points + p];
167853ef2869SZach Atkins             for (CeedInt d = 0; d < dim; d++) {
16793778dbaaSJeremy L Thompson               // ------ Tensor contract with current Chebyshev polynomial values
16809c34f28eSJeremy L Thompson               CeedCall(CeedChebyshevPolynomialsAtPoint(x_array_read[d * num_points + p], Q_1d, chebyshev_x));
16814608bdaaSJeremy L Thompson               CeedCall(CeedTensorContractApply(basis->contract, pre, 1, post, Q_1d, chebyshev_x, t_mode, p > 0 && d == (dim - 1), tmp[d % 2],
16824608bdaaSJeremy L Thompson                                                d == (dim - 1) ? chebyshev_coeffs : tmp[(d + 1) % 2]));
16832a94f45fSJeremy L Thompson               pre /= 1;
16842a94f45fSJeremy L Thompson               post *= Q_1d;
16852a94f45fSJeremy L Thompson             }
16862a94f45fSJeremy L Thompson           }
1687038a8942SZach Atkins           break;
1688038a8942SZach Atkins         }
1689038a8942SZach Atkins         case CEED_EVAL_GRAD: {
1690038a8942SZach Atkins           CeedScalar tmp[2][num_comp * CeedIntPow(Q_1d, dim)], chebyshev_x[Q_1d];
1691038a8942SZach Atkins 
1692038a8942SZach Atkins           // ---- Values at point
1693038a8942SZach Atkins           for (CeedInt p = 0; p < num_points; p++) {
1694038a8942SZach Atkins             // Dim**2 contractions, apply grad when pass == dim
1695038a8942SZach Atkins             for (CeedInt pass = 0; pass < dim; pass++) {
1696038a8942SZach Atkins               CeedInt pre = num_comp * 1, post = 1;
1697038a8942SZach Atkins 
16984608bdaaSJeremy L Thompson               for (CeedInt c = 0; c < num_comp; c++) tmp[0][c] = u_array[(pass * num_comp + c) * num_points + p];
1699038a8942SZach Atkins               for (CeedInt d = 0; d < dim; d++) {
1700038a8942SZach Atkins                 // ------ Tensor contract with current Chebyshev polynomial values
17019c34f28eSJeremy L Thompson                 if (pass == d) CeedCall(CeedChebyshevDerivativeAtPoint(x_array_read[d * num_points + p], Q_1d, chebyshev_x));
17029c34f28eSJeremy L Thompson                 else CeedCall(CeedChebyshevPolynomialsAtPoint(x_array_read[d * num_points + p], Q_1d, chebyshev_x));
17034608bdaaSJeremy L Thompson                 CeedCall(CeedTensorContractApply(basis->contract, pre, 1, post, Q_1d, chebyshev_x, t_mode,
17044608bdaaSJeremy L Thompson                                                  (p > 0 || (p == 0 && pass > 0)) && d == (dim - 1), tmp[d % 2],
17054608bdaaSJeremy L Thompson                                                  d == (dim - 1) ? chebyshev_coeffs : tmp[(d + 1) % 2]));
1706038a8942SZach Atkins                 pre /= 1;
1707038a8942SZach Atkins                 post *= Q_1d;
1708038a8942SZach Atkins               }
1709038a8942SZach Atkins             }
1710038a8942SZach Atkins           }
1711038a8942SZach Atkins           break;
1712038a8942SZach Atkins         }
1713038a8942SZach Atkins         default:
1714038a8942SZach Atkins           // Nothing to do, excluded above
1715038a8942SZach Atkins           break;
17162a94f45fSJeremy L Thompson       }
17172a94f45fSJeremy L Thompson       CeedCall(CeedVectorRestoreArray(basis->vec_chebyshev, &chebyshev_coeffs));
17182a94f45fSJeremy L Thompson       CeedCall(CeedVectorRestoreArrayRead(x_ref, &x_array_read));
17192a94f45fSJeremy L Thompson       CeedCall(CeedVectorRestoreArrayRead(u, &u_array));
17202a94f45fSJeremy L Thompson 
17212a94f45fSJeremy L Thompson       // -- Interpolate transpose from Chebyshev coefficients
17222a94f45fSJeremy L Thompson       CeedCall(CeedBasisApply(basis->basis_chebyshev, 1, CEED_TRANSPOSE, CEED_EVAL_INTERP, basis->vec_chebyshev, v));
17232a94f45fSJeremy L Thompson       break;
17242a94f45fSJeremy L Thompson     }
1725c8c3fa7dSJeremy L Thompson   }
1726c8c3fa7dSJeremy L Thompson   return CEED_ERROR_SUCCESS;
1727c8c3fa7dSJeremy L Thompson }
1728c8c3fa7dSJeremy L Thompson 
1729c8c3fa7dSJeremy L Thompson /**
1730ca94c3ddSJeremy L Thompson   @brief Get `Ceed` associated with a `CeedBasis`
1731b7c9bbdaSJeremy L Thompson 
1732ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
1733ca94c3ddSJeremy L Thompson   @param[out] ceed  Variable to store `Ceed`
1734b7c9bbdaSJeremy L Thompson 
1735b7c9bbdaSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1736b7c9bbdaSJeremy L Thompson 
1737b7c9bbdaSJeremy L Thompson   @ref Advanced
1738b7c9bbdaSJeremy L Thompson **/
1739b7c9bbdaSJeremy L Thompson int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) {
1740b7c9bbdaSJeremy L Thompson   *ceed = basis->ceed;
1741b7c9bbdaSJeremy L Thompson   return CEED_ERROR_SUCCESS;
1742b7c9bbdaSJeremy L Thompson }
1743b7c9bbdaSJeremy L Thompson 
1744b7c9bbdaSJeremy L Thompson /**
1745ca94c3ddSJeremy L Thompson   @brief Get dimension for given `CeedBasis`
17469d007619Sjeremylt 
1747ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
17489d007619Sjeremylt   @param[out] dim   Variable to store dimension of basis
17499d007619Sjeremylt 
17509d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
17519d007619Sjeremylt 
1752b7c9bbdaSJeremy L Thompson   @ref Advanced
17539d007619Sjeremylt **/
17549d007619Sjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) {
17559d007619Sjeremylt   *dim = basis->dim;
1756e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
17579d007619Sjeremylt }
17589d007619Sjeremylt 
17599d007619Sjeremylt /**
1760ca94c3ddSJeremy L Thompson   @brief Get topology for given `CeedBasis`
1761d99fa3c5SJeremy L Thompson 
1762ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
1763d99fa3c5SJeremy L Thompson   @param[out] topo  Variable to store topology of basis
1764d99fa3c5SJeremy L Thompson 
1765d99fa3c5SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1766d99fa3c5SJeremy L Thompson 
1767b7c9bbdaSJeremy L Thompson   @ref Advanced
1768d99fa3c5SJeremy L Thompson **/
1769d99fa3c5SJeremy L Thompson int CeedBasisGetTopology(CeedBasis basis, CeedElemTopology *topo) {
1770d99fa3c5SJeremy L Thompson   *topo = basis->topo;
1771e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1772d99fa3c5SJeremy L Thompson }
1773d99fa3c5SJeremy L Thompson 
1774d99fa3c5SJeremy L Thompson /**
1775ca94c3ddSJeremy L Thompson   @brief Get number of components for given `CeedBasis`
17769d007619Sjeremylt 
1777ca94c3ddSJeremy L Thompson   @param[in]  basis    `CeedBasis`
1778ca94c3ddSJeremy L Thompson   @param[out] num_comp Variable to store number of components
17799d007619Sjeremylt 
17809d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
17819d007619Sjeremylt 
1782b7c9bbdaSJeremy L Thompson   @ref Advanced
17839d007619Sjeremylt **/
1784d1d35e2fSjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *num_comp) {
1785d1d35e2fSjeremylt   *num_comp = basis->num_comp;
1786e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
17879d007619Sjeremylt }
17889d007619Sjeremylt 
17899d007619Sjeremylt /**
1790ca94c3ddSJeremy L Thompson   @brief Get total number of nodes (in `dim` dimensions) of a `CeedBasis`
17919d007619Sjeremylt 
1792ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
17939d007619Sjeremylt   @param[out] P     Variable to store number of nodes
17949d007619Sjeremylt 
17959d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
17969d007619Sjeremylt 
17979d007619Sjeremylt   @ref Utility
17989d007619Sjeremylt **/
17999d007619Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) {
18009d007619Sjeremylt   *P = basis->P;
1801e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
18029d007619Sjeremylt }
18039d007619Sjeremylt 
18049d007619Sjeremylt /**
1805ca94c3ddSJeremy L Thompson   @brief Get total number of nodes (in 1 dimension) of a `CeedBasis`
18069d007619Sjeremylt 
1807ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
1808d1d35e2fSjeremylt   @param[out] P_1d  Variable to store number of nodes
18099d007619Sjeremylt 
18109d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
18119d007619Sjeremylt 
1812b7c9bbdaSJeremy L Thompson   @ref Advanced
18139d007619Sjeremylt **/
1814d1d35e2fSjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P_1d) {
1815ca94c3ddSJeremy L Thompson   CeedCheck(basis->is_tensor_basis, basis->ceed, CEED_ERROR_MINOR, "Cannot supply P_1d for non-tensor CeedBasis");
1816d1d35e2fSjeremylt   *P_1d = basis->P_1d;
1817e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
18189d007619Sjeremylt }
18199d007619Sjeremylt 
18209d007619Sjeremylt /**
1821ca94c3ddSJeremy L Thompson   @brief Get total number of quadrature points (in `dim` dimensions) of a `CeedBasis`
18229d007619Sjeremylt 
1823ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
18249d007619Sjeremylt   @param[out] Q     Variable to store number of quadrature points
18259d007619Sjeremylt 
18269d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
18279d007619Sjeremylt 
18289d007619Sjeremylt   @ref Utility
18299d007619Sjeremylt **/
18309d007619Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) {
18319d007619Sjeremylt   *Q = basis->Q;
1832e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
18339d007619Sjeremylt }
18349d007619Sjeremylt 
18359d007619Sjeremylt /**
1836ca94c3ddSJeremy L Thompson   @brief Get total number of quadrature points (in 1 dimension) of a `CeedBasis`
18379d007619Sjeremylt 
1838ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
1839d1d35e2fSjeremylt   @param[out] Q_1d  Variable to store number of quadrature points
18409d007619Sjeremylt 
18419d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
18429d007619Sjeremylt 
1843b7c9bbdaSJeremy L Thompson   @ref Advanced
18449d007619Sjeremylt **/
1845d1d35e2fSjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q_1d) {
1846ca94c3ddSJeremy L Thompson   CeedCheck(basis->is_tensor_basis, basis->ceed, CEED_ERROR_MINOR, "Cannot supply Q_1d for non-tensor CeedBasis");
1847d1d35e2fSjeremylt   *Q_1d = basis->Q_1d;
1848e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
18499d007619Sjeremylt }
18509d007619Sjeremylt 
18519d007619Sjeremylt /**
1852ca94c3ddSJeremy L Thompson   @brief Get reference coordinates of quadrature points (in `dim` dimensions) of a `CeedBasis`
18539d007619Sjeremylt 
1854ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
1855d1d35e2fSjeremylt   @param[out] q_ref Variable to store reference coordinates of quadrature points
18569d007619Sjeremylt 
18579d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
18589d007619Sjeremylt 
1859b7c9bbdaSJeremy L Thompson   @ref Advanced
18609d007619Sjeremylt **/
1861d1d35e2fSjeremylt int CeedBasisGetQRef(CeedBasis basis, const CeedScalar **q_ref) {
1862d1d35e2fSjeremylt   *q_ref = basis->q_ref_1d;
1863e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
18649d007619Sjeremylt }
18659d007619Sjeremylt 
18669d007619Sjeremylt /**
1867ca94c3ddSJeremy L Thompson   @brief Get quadrature weights of quadrature points (in `dim` dimensions) of a `CeedBasis`
18689d007619Sjeremylt 
1869ca94c3ddSJeremy L Thompson   @param[in]  basis    `CeedBasis`
1870d1d35e2fSjeremylt   @param[out] q_weight Variable to store quadrature weights
18719d007619Sjeremylt 
18729d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
18739d007619Sjeremylt 
1874b7c9bbdaSJeremy L Thompson   @ref Advanced
18759d007619Sjeremylt **/
1876d1d35e2fSjeremylt int CeedBasisGetQWeights(CeedBasis basis, const CeedScalar **q_weight) {
1877d1d35e2fSjeremylt   *q_weight = basis->q_weight_1d;
1878e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
18799d007619Sjeremylt }
18809d007619Sjeremylt 
18819d007619Sjeremylt /**
1882ca94c3ddSJeremy L Thompson   @brief Get interpolation matrix of a `CeedBasis`
18839d007619Sjeremylt 
1884ca94c3ddSJeremy L Thompson   @param[in]  basis  `CeedBasis`
18859d007619Sjeremylt   @param[out] interp Variable to store interpolation matrix
18869d007619Sjeremylt 
18879d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
18889d007619Sjeremylt 
1889b7c9bbdaSJeremy L Thompson   @ref Advanced
18909d007619Sjeremylt **/
18916c58de82SJeremy L Thompson int CeedBasisGetInterp(CeedBasis basis, const CeedScalar **interp) {
18926402da51SJeremy L Thompson   if (!basis->interp && basis->is_tensor_basis) {
18939d007619Sjeremylt     // Allocate
18942b730f8bSJeremy L Thompson     CeedCall(CeedMalloc(basis->Q * basis->P, &basis->interp));
18959d007619Sjeremylt 
18969d007619Sjeremylt     // Initialize
18972b730f8bSJeremy L Thompson     for (CeedInt i = 0; i < basis->Q * basis->P; i++) basis->interp[i] = 1.0;
18989d007619Sjeremylt 
18999d007619Sjeremylt     // Calculate
19002b730f8bSJeremy L Thompson     for (CeedInt d = 0; d < basis->dim; d++) {
19012b730f8bSJeremy L Thompson       for (CeedInt qpt = 0; qpt < basis->Q; qpt++) {
19029d007619Sjeremylt         for (CeedInt node = 0; node < basis->P; node++) {
1903d1d35e2fSjeremylt           CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d;
1904d1d35e2fSjeremylt           CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d;
19051c66c397SJeremy L Thompson 
1906d1d35e2fSjeremylt           basis->interp[qpt * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p];
19079d007619Sjeremylt         }
19089d007619Sjeremylt       }
19092b730f8bSJeremy L Thompson     }
19102b730f8bSJeremy L Thompson   }
19119d007619Sjeremylt   *interp = basis->interp;
1912e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
19139d007619Sjeremylt }
19149d007619Sjeremylt 
19159d007619Sjeremylt /**
1916ca94c3ddSJeremy L Thompson   @brief Get 1D interpolation matrix of a tensor product `CeedBasis`
19179d007619Sjeremylt 
1918ca94c3ddSJeremy L Thompson   @param[in]  basis     `CeedBasis`
1919d1d35e2fSjeremylt   @param[out] interp_1d Variable to store interpolation matrix
19209d007619Sjeremylt 
19219d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
19229d007619Sjeremylt 
19239d007619Sjeremylt   @ref Backend
19249d007619Sjeremylt **/
1925d1d35e2fSjeremylt int CeedBasisGetInterp1D(CeedBasis basis, const CeedScalar **interp_1d) {
1926ca94c3ddSJeremy L Thompson   CeedCheck(basis->is_tensor_basis, basis->ceed, CEED_ERROR_MINOR, "CeedBasis is not a tensor product CeedBasis");
1927d1d35e2fSjeremylt   *interp_1d = basis->interp_1d;
1928e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
19299d007619Sjeremylt }
19309d007619Sjeremylt 
19319d007619Sjeremylt /**
1932ca94c3ddSJeremy L Thompson   @brief Get gradient matrix of a `CeedBasis`
19339d007619Sjeremylt 
1934ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
19359d007619Sjeremylt   @param[out] grad  Variable to store gradient matrix
19369d007619Sjeremylt 
19379d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
19389d007619Sjeremylt 
1939b7c9bbdaSJeremy L Thompson   @ref Advanced
19409d007619Sjeremylt **/
19416c58de82SJeremy L Thompson int CeedBasisGetGrad(CeedBasis basis, const CeedScalar **grad) {
19426402da51SJeremy L Thompson   if (!basis->grad && basis->is_tensor_basis) {
19439d007619Sjeremylt     // Allocate
19442b730f8bSJeremy L Thompson     CeedCall(CeedMalloc(basis->dim * basis->Q * basis->P, &basis->grad));
19459d007619Sjeremylt 
19469d007619Sjeremylt     // Initialize
19472b730f8bSJeremy L Thompson     for (CeedInt i = 0; i < basis->dim * basis->Q * basis->P; i++) basis->grad[i] = 1.0;
19489d007619Sjeremylt 
19499d007619Sjeremylt     // Calculate
19502b730f8bSJeremy L Thompson     for (CeedInt d = 0; d < basis->dim; d++) {
19512b730f8bSJeremy L Thompson       for (CeedInt i = 0; i < basis->dim; i++) {
19522b730f8bSJeremy L Thompson         for (CeedInt qpt = 0; qpt < basis->Q; qpt++) {
19539d007619Sjeremylt           for (CeedInt node = 0; node < basis->P; node++) {
1954d1d35e2fSjeremylt             CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d;
1955d1d35e2fSjeremylt             CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d;
19561c66c397SJeremy L Thompson 
19572b730f8bSJeremy L Thompson             if (i == d) basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->grad_1d[q * basis->P_1d + p];
19582b730f8bSJeremy L Thompson             else basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p];
19592b730f8bSJeremy L Thompson           }
19602b730f8bSJeremy L Thompson         }
19612b730f8bSJeremy L Thompson       }
19629d007619Sjeremylt     }
19639d007619Sjeremylt   }
19649d007619Sjeremylt   *grad = basis->grad;
1965e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
19669d007619Sjeremylt }
19679d007619Sjeremylt 
19689d007619Sjeremylt /**
1969ca94c3ddSJeremy L Thompson   @brief Get 1D gradient matrix of a tensor product `CeedBasis`
19709d007619Sjeremylt 
1971ca94c3ddSJeremy L Thompson   @param[in]  basis   `CeedBasis`
1972d1d35e2fSjeremylt   @param[out] grad_1d Variable to store gradient matrix
19739d007619Sjeremylt 
19749d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
19759d007619Sjeremylt 
1976b7c9bbdaSJeremy L Thompson   @ref Advanced
19779d007619Sjeremylt **/
1978d1d35e2fSjeremylt int CeedBasisGetGrad1D(CeedBasis basis, const CeedScalar **grad_1d) {
1979ca94c3ddSJeremy L Thompson   CeedCheck(basis->is_tensor_basis, basis->ceed, CEED_ERROR_MINOR, "CeedBasis is not a tensor product CeedBasis");
1980d1d35e2fSjeremylt   *grad_1d = basis->grad_1d;
1981e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
19829d007619Sjeremylt }
19839d007619Sjeremylt 
19849d007619Sjeremylt /**
1985ca94c3ddSJeremy L Thompson   @brief Get divergence matrix of a `CeedBasis`
198650c301a5SRezgar Shakeri 
1987ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
198850c301a5SRezgar Shakeri   @param[out] div   Variable to store divergence matrix
198950c301a5SRezgar Shakeri 
199050c301a5SRezgar Shakeri   @return An error code: 0 - success, otherwise - failure
199150c301a5SRezgar Shakeri 
199250c301a5SRezgar Shakeri   @ref Advanced
199350c301a5SRezgar Shakeri **/
199450c301a5SRezgar Shakeri int CeedBasisGetDiv(CeedBasis basis, const CeedScalar **div) {
1995ca94c3ddSJeremy L Thompson   CeedCheck(basis->div, basis->ceed, CEED_ERROR_MINOR, "CeedBasis does not have divergence matrix");
199650c301a5SRezgar Shakeri   *div = basis->div;
199750c301a5SRezgar Shakeri   return CEED_ERROR_SUCCESS;
199850c301a5SRezgar Shakeri }
199950c301a5SRezgar Shakeri 
200050c301a5SRezgar Shakeri /**
2001ca94c3ddSJeremy L Thompson   @brief Get curl matrix of a `CeedBasis`
2002c4e3f59bSSebastian Grimberg 
2003ca94c3ddSJeremy L Thompson   @param[in]  basis `CeedBasis`
2004c4e3f59bSSebastian Grimberg   @param[out] curl  Variable to store curl matrix
2005c4e3f59bSSebastian Grimberg 
2006c4e3f59bSSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
2007c4e3f59bSSebastian Grimberg 
2008c4e3f59bSSebastian Grimberg   @ref Advanced
2009c4e3f59bSSebastian Grimberg **/
2010c4e3f59bSSebastian Grimberg int CeedBasisGetCurl(CeedBasis basis, const CeedScalar **curl) {
2011ca94c3ddSJeremy L Thompson   CeedCheck(basis->curl, basis->ceed, CEED_ERROR_MINOR, "CeedBasis does not have curl matrix");
2012c4e3f59bSSebastian Grimberg   *curl = basis->curl;
2013c4e3f59bSSebastian Grimberg   return CEED_ERROR_SUCCESS;
2014c4e3f59bSSebastian Grimberg }
2015c4e3f59bSSebastian Grimberg 
2016c4e3f59bSSebastian Grimberg /**
2017ca94c3ddSJeremy L Thompson   @brief Destroy a @ref  CeedBasis
20187a982d89SJeremy L. Thompson 
2019ca94c3ddSJeremy L Thompson   @param[in,out] basis `CeedBasis` to destroy
20207a982d89SJeremy L. Thompson 
20217a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
20227a982d89SJeremy L. Thompson 
20237a982d89SJeremy L. Thompson   @ref User
20247a982d89SJeremy L. Thompson **/
20257a982d89SJeremy L. Thompson int CeedBasisDestroy(CeedBasis *basis) {
2026356036faSJeremy L Thompson   if (!*basis || *basis == CEED_BASIS_NONE || --(*basis)->ref_count > 0) {
2027ad6481ceSJeremy L Thompson     *basis = NULL;
2028ad6481ceSJeremy L Thompson     return CEED_ERROR_SUCCESS;
2029ad6481ceSJeremy L Thompson   }
20302b730f8bSJeremy L Thompson   if ((*basis)->Destroy) CeedCall((*basis)->Destroy(*basis));
20319831d45aSJeremy L Thompson   CeedCall(CeedTensorContractDestroy(&(*basis)->contract));
2032c4e3f59bSSebastian Grimberg   CeedCall(CeedFree(&(*basis)->q_ref_1d));
2033c4e3f59bSSebastian Grimberg   CeedCall(CeedFree(&(*basis)->q_weight_1d));
20342b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->interp));
20352b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->interp_1d));
20362b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->grad));
20372b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->grad_1d));
2038c4e3f59bSSebastian Grimberg   CeedCall(CeedFree(&(*basis)->div));
2039c4e3f59bSSebastian Grimberg   CeedCall(CeedFree(&(*basis)->curl));
2040c8c3fa7dSJeremy L Thompson   CeedCall(CeedVectorDestroy(&(*basis)->vec_chebyshev));
2041c8c3fa7dSJeremy L Thompson   CeedCall(CeedBasisDestroy(&(*basis)->basis_chebyshev));
20422b730f8bSJeremy L Thompson   CeedCall(CeedDestroy(&(*basis)->ceed));
20432b730f8bSJeremy L Thompson   CeedCall(CeedFree(basis));
2044e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
20457a982d89SJeremy L. Thompson }
20467a982d89SJeremy L. Thompson 
20477a982d89SJeremy L. Thompson /**
2048b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre quadrature
2049b11c1e72Sjeremylt 
2050ca94c3ddSJeremy L Thompson   @param[in]  Q           Number of quadrature points (integrates polynomials of degree `2*Q-1` exactly)
2051ca94c3ddSJeremy L Thompson   @param[out] q_ref_1d    Array of length `Q` to hold the abscissa on `[-1, 1]`
2052ca94c3ddSJeremy L Thompson   @param[out] q_weight_1d Array of length `Q` to hold the weights
2053b11c1e72Sjeremylt 
2054b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
2055dfdf5a53Sjeremylt 
2056dfdf5a53Sjeremylt   @ref Utility
2057b11c1e72Sjeremylt **/
20582b730f8bSJeremy L Thompson int CeedGaussQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) {
2059d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0 * atan(1.0);
20601c66c397SJeremy L Thompson 
2061d1d35e2fSjeremylt   // Build q_ref_1d, q_weight_1d
206292ae7e47SJeremy L Thompson   for (CeedInt i = 0; i <= Q / 2; i++) {
2063d7b241e6Sjeremylt     // Guess
2064d7b241e6Sjeremylt     xi = cos(PI * (CeedScalar)(2 * i + 1) / ((CeedScalar)(2 * Q)));
2065d7b241e6Sjeremylt     // Pn(xi)
2066d7b241e6Sjeremylt     P0 = 1.0;
2067d7b241e6Sjeremylt     P1 = xi;
2068d7b241e6Sjeremylt     P2 = 0.0;
206992ae7e47SJeremy L Thompson     for (CeedInt j = 2; j <= Q; j++) {
2070d7b241e6Sjeremylt       P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
2071d7b241e6Sjeremylt       P0 = P1;
2072d7b241e6Sjeremylt       P1 = P2;
2073d7b241e6Sjeremylt     }
2074d7b241e6Sjeremylt     // First Newton Step
2075d7b241e6Sjeremylt     dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
2076d7b241e6Sjeremylt     xi  = xi - P2 / dP2;
2077d7b241e6Sjeremylt     // Newton to convergence
207892ae7e47SJeremy L Thompson     for (CeedInt k = 0; k < 100 && fabs(P2) > 10 * CEED_EPSILON; k++) {
2079d7b241e6Sjeremylt       P0 = 1.0;
2080d7b241e6Sjeremylt       P1 = xi;
208192ae7e47SJeremy L Thompson       for (CeedInt j = 2; j <= Q; j++) {
2082d7b241e6Sjeremylt         P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
2083d7b241e6Sjeremylt         P0 = P1;
2084d7b241e6Sjeremylt         P1 = P2;
2085d7b241e6Sjeremylt       }
2086d7b241e6Sjeremylt       dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
2087d7b241e6Sjeremylt       xi  = xi - P2 / dP2;
2088d7b241e6Sjeremylt     }
2089d7b241e6Sjeremylt     // Save xi, wi
2090d7b241e6Sjeremylt     wi                     = 2.0 / ((1.0 - xi * xi) * dP2 * dP2);
2091d1d35e2fSjeremylt     q_weight_1d[i]         = wi;
2092d1d35e2fSjeremylt     q_weight_1d[Q - 1 - i] = wi;
2093d1d35e2fSjeremylt     q_ref_1d[i]            = -xi;
2094d1d35e2fSjeremylt     q_ref_1d[Q - 1 - i]    = xi;
2095d7b241e6Sjeremylt   }
2096e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2097d7b241e6Sjeremylt }
2098d7b241e6Sjeremylt 
2099b11c1e72Sjeremylt /**
2100b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre-Lobatto quadrature
2101b11c1e72Sjeremylt 
2102ca94c3ddSJeremy L Thompson   @param[in]  Q           Number of quadrature points (integrates polynomials of degree `2*Q-3` exactly)
2103ca94c3ddSJeremy L Thompson   @param[out] q_ref_1d    Array of length `Q` to hold the abscissa on `[-1, 1]`
2104ca94c3ddSJeremy L Thompson   @param[out] q_weight_1d Array of length `Q` to hold the weights
2105b11c1e72Sjeremylt 
2106b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
2107dfdf5a53Sjeremylt 
2108dfdf5a53Sjeremylt   @ref Utility
2109b11c1e72Sjeremylt **/
21102b730f8bSJeremy L Thompson int CeedLobattoQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) {
2111d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0 * atan(1.0);
21121c66c397SJeremy L Thompson 
2113d1d35e2fSjeremylt   // Build q_ref_1d, q_weight_1d
2114d7b241e6Sjeremylt   // Set endpoints
21156574a04fSJeremy L Thompson   CeedCheck(Q > 1, NULL, CEED_ERROR_DIMENSION, "Cannot create Lobatto quadrature with Q=%" CeedInt_FMT " < 2 points", Q);
2116d7b241e6Sjeremylt   wi = 2.0 / ((CeedScalar)(Q * (Q - 1)));
2117d1d35e2fSjeremylt   if (q_weight_1d) {
2118d1d35e2fSjeremylt     q_weight_1d[0]     = wi;
2119d1d35e2fSjeremylt     q_weight_1d[Q - 1] = wi;
2120d7b241e6Sjeremylt   }
2121d1d35e2fSjeremylt   q_ref_1d[0]     = -1.0;
2122d1d35e2fSjeremylt   q_ref_1d[Q - 1] = 1.0;
2123d7b241e6Sjeremylt   // Interior
212492ae7e47SJeremy L Thompson   for (CeedInt i = 1; i <= (Q - 1) / 2; i++) {
2125d7b241e6Sjeremylt     // Guess
2126d7b241e6Sjeremylt     xi = cos(PI * (CeedScalar)(i) / (CeedScalar)(Q - 1));
2127d7b241e6Sjeremylt     // Pn(xi)
2128d7b241e6Sjeremylt     P0 = 1.0;
2129d7b241e6Sjeremylt     P1 = xi;
2130d7b241e6Sjeremylt     P2 = 0.0;
213192ae7e47SJeremy L Thompson     for (CeedInt j = 2; j < Q; j++) {
2132d7b241e6Sjeremylt       P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
2133d7b241e6Sjeremylt       P0 = P1;
2134d7b241e6Sjeremylt       P1 = P2;
2135d7b241e6Sjeremylt     }
2136d7b241e6Sjeremylt     // First Newton step
2137d7b241e6Sjeremylt     dP2  = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
2138d7b241e6Sjeremylt     d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi);
2139d7b241e6Sjeremylt     xi   = xi - dP2 / d2P2;
2140d7b241e6Sjeremylt     // Newton to convergence
214192ae7e47SJeremy L Thompson     for (CeedInt k = 0; k < 100 && fabs(dP2) > 10 * CEED_EPSILON; k++) {
2142d7b241e6Sjeremylt       P0 = 1.0;
2143d7b241e6Sjeremylt       P1 = xi;
214492ae7e47SJeremy L Thompson       for (CeedInt j = 2; j < Q; j++) {
2145d7b241e6Sjeremylt         P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
2146d7b241e6Sjeremylt         P0 = P1;
2147d7b241e6Sjeremylt         P1 = P2;
2148d7b241e6Sjeremylt       }
2149d7b241e6Sjeremylt       dP2  = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
2150d7b241e6Sjeremylt       d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi);
2151d7b241e6Sjeremylt       xi   = xi - dP2 / d2P2;
2152d7b241e6Sjeremylt     }
2153d7b241e6Sjeremylt     // Save xi, wi
2154d7b241e6Sjeremylt     wi = 2.0 / (((CeedScalar)(Q * (Q - 1))) * P2 * P2);
2155d1d35e2fSjeremylt     if (q_weight_1d) {
2156d1d35e2fSjeremylt       q_weight_1d[i]         = wi;
2157d1d35e2fSjeremylt       q_weight_1d[Q - 1 - i] = wi;
2158d7b241e6Sjeremylt     }
2159d1d35e2fSjeremylt     q_ref_1d[i]         = -xi;
2160d1d35e2fSjeremylt     q_ref_1d[Q - 1 - i] = xi;
2161d7b241e6Sjeremylt   }
2162e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2163d7b241e6Sjeremylt }
2164d7b241e6Sjeremylt 
2165d7b241e6Sjeremylt /// @}
2166