xref: /libCEED/rust/libceed-sys/c-src/interface/ceed-basis.c (revision 3778dbaabce4c4854d6e52be1d445d2636534c17)
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
20783c99b3SValeria Barra static struct CeedBasis_private ceed_basis_collocated;
21d7b241e6Sjeremylt /// @endcond
22d7b241e6Sjeremylt 
237a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser
247a982d89SJeremy L. Thompson /// @{
257a982d89SJeremy L. Thompson 
267a982d89SJeremy L. Thompson /// Indicate that the quadrature points are collocated with the nodes
277a982d89SJeremy L. Thompson const CeedBasis CEED_BASIS_COLLOCATED = &ceed_basis_collocated;
287a982d89SJeremy L. Thompson 
297a982d89SJeremy L. Thompson /// @}
307a982d89SJeremy L. Thompson 
317a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
327a982d89SJeremy L. Thompson /// CeedBasis Library Internal Functions
337a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
347a982d89SJeremy L. Thompson /// @addtogroup CeedBasisDeveloper
357a982d89SJeremy L. Thompson /// @{
367a982d89SJeremy L. Thompson 
377a982d89SJeremy L. Thompson /**
38*3778dbaaSJeremy L Thompson   @brief Compute Chebyshev polynomial values at a point
39*3778dbaaSJeremy L Thompson 
40*3778dbaaSJeremy L Thompson   @param[in]  x           Coordinate to evaluate Chebyshev polynomials at
41*3778dbaaSJeremy L Thompson   @param[in]  n           Number of Chebyshev polynomials to evaluate, n >= 2
42*3778dbaaSJeremy L Thompson   @param[out] chebyshev_x Array of Chebyshev polynomial values
43*3778dbaaSJeremy L Thompson 
44*3778dbaaSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
45*3778dbaaSJeremy L Thompson 
46*3778dbaaSJeremy L Thompson   @ref Developer
47*3778dbaaSJeremy L Thompson **/
48*3778dbaaSJeremy L Thompson static int CeedChebyshevPolynomialsAtPoint(CeedScalar x, CeedInt n, CeedScalar *chebyshev_x) {
49*3778dbaaSJeremy L Thompson   chebyshev_x[0] = 1.0;
50*3778dbaaSJeremy L Thompson   chebyshev_x[1] = 2 * x;
51*3778dbaaSJeremy L Thompson   for (CeedInt i = 2; i < n; i++) chebyshev_x[i] = 2 * x * chebyshev_x[i - 1] - chebyshev_x[i - 2];
52*3778dbaaSJeremy L Thompson 
53*3778dbaaSJeremy L Thompson   return CEED_ERROR_SUCCESS;
54*3778dbaaSJeremy L Thompson }
55*3778dbaaSJeremy L Thompson 
56*3778dbaaSJeremy L Thompson /**
57*3778dbaaSJeremy L Thompson   @brief Compute values of the derivative of Chebyshev polynomials at a point
58*3778dbaaSJeremy L Thompson 
59*3778dbaaSJeremy L Thompson   @param[in]  x           Coordinate to evaluate derivative of Chebyshev polynomials at
60*3778dbaaSJeremy L Thompson   @param[in]  n           Number of Chebyshev polynomials to evaluate, n >= 2
61*3778dbaaSJeremy L Thompson   @param[out] chebyshev_x Array of Chebyshev polynomial derivative values
62*3778dbaaSJeremy L Thompson 
63*3778dbaaSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
64*3778dbaaSJeremy L Thompson 
65*3778dbaaSJeremy L Thompson   @ref Developer
66*3778dbaaSJeremy L Thompson **/
67*3778dbaaSJeremy L Thompson static int CeedChebyshevDerivativeAtPoint(CeedScalar x, CeedInt n, CeedScalar *chebyshev_dx) {
68*3778dbaaSJeremy L Thompson   CeedScalar chebyshev_x[3];
69*3778dbaaSJeremy L Thompson 
70*3778dbaaSJeremy L Thompson   chebyshev_x[1]  = 1.0;
71*3778dbaaSJeremy L Thompson   chebyshev_x[2]  = 2 * x;
72*3778dbaaSJeremy L Thompson   chebyshev_dx[0] = 0.0;
73*3778dbaaSJeremy L Thompson   chebyshev_dx[1] = 2.0;
74*3778dbaaSJeremy L Thompson   for (CeedInt i = 2; i < n; i++) {
75*3778dbaaSJeremy L Thompson     chebyshev_x[0]  = chebyshev_x[1];
76*3778dbaaSJeremy L Thompson     chebyshev_x[1]  = chebyshev_x[2];
77*3778dbaaSJeremy L Thompson     chebyshev_x[2]  = 2 * x * chebyshev_x[1] - chebyshev_x[0];
78*3778dbaaSJeremy L Thompson     chebyshev_dx[i] = 2 * x * chebyshev_dx[i - 1] + 2 * chebyshev_x[1] - chebyshev_dx[i - 2];
79*3778dbaaSJeremy L Thompson   }
80*3778dbaaSJeremy L Thompson 
81*3778dbaaSJeremy L Thompson   return CEED_ERROR_SUCCESS;
82*3778dbaaSJeremy L Thompson }
83*3778dbaaSJeremy L Thompson 
84*3778dbaaSJeremy L Thompson /**
857a982d89SJeremy L. Thompson   @brief Compute Householder reflection
867a982d89SJeremy L. Thompson 
87ea61e9acSJeremy L Thompson   Computes A = (I - b v v^T) A, where A is an mxn matrix indexed as A[i*row + j*col]
887a982d89SJeremy L. Thompson 
897a982d89SJeremy L. Thompson   @param[in,out] A   Matrix to apply Householder reflection to, in place
90ea61e9acSJeremy L Thompson   @param[in]     v   Householder vector
91ea61e9acSJeremy L Thompson   @param[in]     b   Scaling factor
92ea61e9acSJeremy L Thompson   @param[in]     m   Number of rows in A
93ea61e9acSJeremy L Thompson   @param[in]     n   Number of columns in A
94ea61e9acSJeremy L Thompson   @param[in]     row Row stride
95ea61e9acSJeremy L Thompson   @param[in]     col Col stride
967a982d89SJeremy L. Thompson 
977a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
987a982d89SJeremy L. Thompson 
997a982d89SJeremy L. Thompson   @ref Developer
1007a982d89SJeremy L. Thompson **/
1012b730f8bSJeremy L Thompson static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v, CeedScalar b, CeedInt m, CeedInt n, CeedInt row, CeedInt col) {
1027a982d89SJeremy L. Thompson   for (CeedInt j = 0; j < n; j++) {
1037a982d89SJeremy L. Thompson     CeedScalar w = A[0 * row + j * col];
1042b730f8bSJeremy L Thompson     for (CeedInt i = 1; i < m; i++) w += v[i] * A[i * row + j * col];
1057a982d89SJeremy L. Thompson     A[0 * row + j * col] -= b * w;
1062b730f8bSJeremy L Thompson     for (CeedInt i = 1; i < m; i++) A[i * row + j * col] -= b * w * v[i];
1077a982d89SJeremy L. Thompson   }
108e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1097a982d89SJeremy L. Thompson }
1107a982d89SJeremy L. Thompson 
1117a982d89SJeremy L. Thompson /**
1127a982d89SJeremy L. Thompson   @brief Compute Givens rotation
1137a982d89SJeremy L. Thompson 
114ea61e9acSJeremy L Thompson   Computes A = G A (or G^T A in transpose mode), where A is an mxn matrix indexed as A[i*n + j*m]
1157a982d89SJeremy L. Thompson 
1167a982d89SJeremy L. Thompson   @param[in,out] A      Row major matrix to apply Givens rotation to, in place
117ea61e9acSJeremy L Thompson   @param[in]     c      Cosine factor
118ea61e9acSJeremy L Thompson   @param[in]     s      Sine factor
119ea61e9acSJeremy 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;
1204cc79fe7SJed Brown                           @ref CEED_TRANSPOSE for the opposite rotation
121ea61e9acSJeremy L Thompson   @param[in]     i      First row/column to apply rotation
122ea61e9acSJeremy L Thompson   @param[in]     k      Second row/column to apply rotation
123ea61e9acSJeremy L Thompson   @param[in]     m      Number of rows in A
124ea61e9acSJeremy L Thompson   @param[in]     n      Number of columns in A
1257a982d89SJeremy L. Thompson 
1267a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
1277a982d89SJeremy L. Thompson 
1287a982d89SJeremy L. Thompson   @ref Developer
1297a982d89SJeremy L. Thompson **/
1302b730f8bSJeremy L Thompson static int CeedGivensRotation(CeedScalar *A, CeedScalar c, CeedScalar s, CeedTransposeMode t_mode, CeedInt i, CeedInt k, CeedInt m, CeedInt n) {
131d1d35e2fSjeremylt   CeedInt stride_j = 1, stride_ik = m, num_its = n;
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];
141d1d35e2fSjeremylt     A[i * stride_ik + j * stride_j] = c * tau1 - s * tau2;
142d1d35e2fSjeremylt     A[k * stride_ik + j * stride_j] = s * tau1 + c * tau2;
1437a982d89SJeremy L. Thompson   }
144e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1457a982d89SJeremy L. Thompson }
1467a982d89SJeremy L. Thompson 
1477a982d89SJeremy L. Thompson /**
1487a982d89SJeremy L. Thompson   @brief View an array stored in a CeedBasis
1497a982d89SJeremy L. Thompson 
1500a0da059Sjeremylt   @param[in] name   Name of array
151d1d35e2fSjeremylt   @param[in] fp_fmt Printing format
1520a0da059Sjeremylt   @param[in] m      Number of rows in array
1530a0da059Sjeremylt   @param[in] n      Number of columns in array
1540a0da059Sjeremylt   @param[in] a      Array to be viewed
1550a0da059Sjeremylt   @param[in] stream Stream to view to, e.g., stdout
1567a982d89SJeremy L. Thompson 
1577a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
1587a982d89SJeremy L. Thompson 
1597a982d89SJeremy L. Thompson   @ref Developer
1607a982d89SJeremy L. Thompson **/
1612b730f8bSJeremy L Thompson static int CeedScalarView(const char *name, const char *fp_fmt, CeedInt m, CeedInt n, const CeedScalar *a, FILE *stream) {
16292ae7e47SJeremy L Thompson   for (CeedInt i = 0; i < m; i++) {
1632b730f8bSJeremy L Thompson     if (m > 1) fprintf(stream, "%12s[%" CeedInt_FMT "]:", name, i);
1642b730f8bSJeremy L Thompson     else fprintf(stream, "%12s:", name);
1652b730f8bSJeremy 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);
1667a982d89SJeremy L. Thompson     fputs("\n", stream);
1677a982d89SJeremy L. Thompson   }
168e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1697a982d89SJeremy L. Thompson }
1707a982d89SJeremy L. Thompson 
171a76a04e7SJeremy L Thompson /**
172ea61e9acSJeremy L Thompson   @brief Create the interpolation and gradient matrices for projection from the nodes of `basis_from` to the nodes of `basis_to`.
173ba59ac12SSebastian Grimberg 
17415ad3917SSebastian Grimberg   The interpolation is given by `interp_project = interp_to^+ * interp_from`, where the pseudoinverse `interp_to^+` is given by QR factorization.
17515ad3917SSebastian Grimberg   The gradient is given by `grad_project = interp_to^+ * grad_from`, and is only computed for H^1 spaces otherwise it should not be used.
17615ad3917SSebastian Grimberg 
177ba59ac12SSebastian Grimberg   Note: `basis_from` and `basis_to` must have compatible quadrature spaces.
178a76a04e7SJeremy L Thompson 
179a76a04e7SJeremy L Thompson   @param[in]  basis_from     CeedBasis to project from
180a76a04e7SJeremy L Thompson   @param[in]  basis_to       CeedBasis to project to
181ea61e9acSJeremy L Thompson   @param[out] interp_project Address of the variable where the newly created interpolation matrix will be stored.
182ea61e9acSJeremy L Thompson   @param[out] grad_project   Address of the variable where the newly created gradient matrix will be stored.
183a76a04e7SJeremy L Thompson 
184a76a04e7SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
185a76a04e7SJeremy L Thompson 
186a76a04e7SJeremy L Thompson   @ref Developer
187a76a04e7SJeremy L Thompson **/
1882b730f8bSJeremy L Thompson static int CeedBasisCreateProjectionMatrices(CeedBasis basis_from, CeedBasis basis_to, CeedScalar **interp_project, CeedScalar **grad_project) {
189a76a04e7SJeremy L Thompson   Ceed ceed;
1902b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetCeed(basis_to, &ceed));
191a76a04e7SJeremy L Thompson 
192a76a04e7SJeremy L Thompson   // Check for compatible quadrature spaces
193a76a04e7SJeremy L Thompson   CeedInt Q_to, Q_from;
1942b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumQuadraturePoints(basis_to, &Q_to));
1952b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumQuadraturePoints(basis_from, &Q_from));
1966574a04fSJeremy L Thompson   CeedCheck(Q_to == Q_from, ceed, CEED_ERROR_DIMENSION, "Bases must have compatible quadrature spaces");
197a76a04e7SJeremy L Thompson 
19814556e63SJeremy L Thompson   // Check for matching tensor or non-tensor
199a76a04e7SJeremy L Thompson   CeedInt P_to, P_from, Q = Q_to;
200a76a04e7SJeremy L Thompson   bool    is_tensor_to, is_tensor_from;
2012b730f8bSJeremy L Thompson   CeedCall(CeedBasisIsTensor(basis_to, &is_tensor_to));
2022b730f8bSJeremy L Thompson   CeedCall(CeedBasisIsTensor(basis_from, &is_tensor_from));
2036574a04fSJeremy L Thompson   CeedCheck(is_tensor_to == is_tensor_from, ceed, CEED_ERROR_MINOR, "Bases must both be tensor or non-tensor");
2046574a04fSJeremy L Thompson   if (is_tensor_to) {
2052b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_to));
2062b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_from));
2072b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumQuadraturePoints1D(basis_from, &Q));
2086574a04fSJeremy L Thompson   } else {
2092b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis_to, &P_to));
2102b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis_from, &P_from));
211a76a04e7SJeremy L Thompson   }
212a76a04e7SJeremy L Thompson 
21315ad3917SSebastian Grimberg   // Check for matching FE space
21415ad3917SSebastian Grimberg   CeedFESpace fe_space_to, fe_space_from;
21515ad3917SSebastian Grimberg   CeedCall(CeedBasisGetFESpace(basis_to, &fe_space_to));
21615ad3917SSebastian Grimberg   CeedCall(CeedBasisGetFESpace(basis_from, &fe_space_from));
2176574a04fSJeremy L Thompson   CeedCheck(fe_space_to == fe_space_from, ceed, CEED_ERROR_MINOR, "Bases must both be the same FE space type");
21815ad3917SSebastian Grimberg 
21914556e63SJeremy L Thompson   // Get source matrices
22015ad3917SSebastian Grimberg   CeedInt           dim, q_comp = 1;
22115ad3917SSebastian Grimberg   const CeedScalar *interp_to_source = NULL, *interp_from_source = NULL;
22214556e63SJeremy L Thompson   CeedScalar       *interp_to, *interp_from, *tau;
2232b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetDimension(basis_to, &dim));
224a76a04e7SJeremy L Thompson   if (is_tensor_to) {
2252b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetInterp1D(basis_to, &interp_to_source));
2262b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetInterp1D(basis_from, &interp_from_source));
227a76a04e7SJeremy L Thompson   } else {
22815ad3917SSebastian Grimberg     CeedCall(CeedBasisGetNumQuadratureComponents(basis_from, CEED_EVAL_INTERP, &q_comp));
2292b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetInterp(basis_to, &interp_to_source));
2302b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetInterp(basis_from, &interp_from_source));
23115ad3917SSebastian Grimberg   }
23215ad3917SSebastian Grimberg   CeedCall(CeedMalloc(Q * P_from * q_comp, &interp_from));
23315ad3917SSebastian Grimberg   CeedCall(CeedMalloc(Q * P_to * q_comp, &interp_to));
23415ad3917SSebastian Grimberg   CeedCall(CeedCalloc(P_to * P_from, interp_project));
23515ad3917SSebastian Grimberg   CeedCall(CeedMalloc(Q * q_comp, &tau));
23615ad3917SSebastian Grimberg 
23715ad3917SSebastian Grimberg   // `grad_project = interp_to^+ * grad_from` is computed for the H^1 space case so the
238de05fbb2SSebastian Grimberg   // projection basis will have a gradient operation (allocated even if not H^1 for the
239de05fbb2SSebastian Grimberg   // basis construction later on)
24015ad3917SSebastian Grimberg   const CeedScalar *grad_from_source = NULL;
24115ad3917SSebastian Grimberg   if (fe_space_to == CEED_FE_SPACE_H1) {
24215ad3917SSebastian Grimberg     if (is_tensor_to) {
24315ad3917SSebastian Grimberg       CeedCall(CeedBasisGetGrad1D(basis_from, &grad_from_source));
24415ad3917SSebastian Grimberg     } else {
2452b730f8bSJeremy L Thompson       CeedCall(CeedBasisGetGrad(basis_from, &grad_from_source));
246a76a04e7SJeremy L Thompson     }
247de05fbb2SSebastian Grimberg   }
24815ad3917SSebastian Grimberg   CeedCall(CeedCalloc(P_to * P_from * (is_tensor_to ? 1 : dim), grad_project));
24915ad3917SSebastian Grimberg 
25015ad3917SSebastian Grimberg   // QR Factorization, interp_to = Q R
25115ad3917SSebastian Grimberg   memcpy(interp_to, interp_to_source, Q * P_to * q_comp * sizeof(interp_to_source[0]));
25215ad3917SSebastian Grimberg   CeedCall(CeedQRFactorization(ceed, interp_to, tau, Q * q_comp, P_to));
253a76a04e7SJeremy L Thompson 
25414556e63SJeremy L Thompson   // Build matrices
25515ad3917SSebastian Grimberg   CeedInt     num_matrices = 1 + (fe_space_to == CEED_FE_SPACE_H1) * (is_tensor_to ? 1 : dim);
25614556e63SJeremy L Thompson   CeedScalar *input_from[num_matrices], *output_project[num_matrices];
25714556e63SJeremy L Thompson   input_from[0]     = (CeedScalar *)interp_from_source;
25814556e63SJeremy L Thompson   output_project[0] = *interp_project;
25914556e63SJeremy L Thompson   for (CeedInt m = 1; m < num_matrices; m++) {
26014556e63SJeremy L Thompson     input_from[m]     = (CeedScalar *)&grad_from_source[(m - 1) * Q * P_from];
26102af4036SJeremy L Thompson     output_project[m] = &((*grad_project)[(m - 1) * P_to * P_from]);
26214556e63SJeremy L Thompson   }
26314556e63SJeremy L Thompson   for (CeedInt m = 0; m < num_matrices; m++) {
26415ad3917SSebastian Grimberg     // Apply Q^T, interp_from = Q^T interp_from
26515ad3917SSebastian Grimberg     memcpy(interp_from, input_from[m], Q * P_from * q_comp * sizeof(input_from[m][0]));
26615ad3917SSebastian Grimberg     CeedCall(CeedHouseholderApplyQ(interp_from, interp_to, tau, CEED_TRANSPOSE, Q * q_comp, P_from, P_to, P_from, 1));
267a76a04e7SJeremy L Thompson 
26815ad3917SSebastian Grimberg     // Apply Rinv, output_project = Rinv interp_from
269a76a04e7SJeremy L Thompson     for (CeedInt j = 0; j < P_from; j++) {  // Column j
2702b730f8bSJeremy L Thompson       output_project[m][j + P_from * (P_to - 1)] = interp_from[j + P_from * (P_to - 1)] / interp_to[P_to * P_to - 1];
271a76a04e7SJeremy L Thompson       for (CeedInt i = P_to - 2; i >= 0; i--) {  // Row i
27214556e63SJeremy L Thompson         output_project[m][j + P_from * i] = interp_from[j + P_from * i];
273a76a04e7SJeremy L Thompson         for (CeedInt k = i + 1; k < P_to; k++) {
2742b730f8bSJeremy L Thompson           output_project[m][j + P_from * i] -= interp_to[k + P_to * i] * output_project[m][j + P_from * k];
275a76a04e7SJeremy L Thompson         }
27614556e63SJeremy L Thompson         output_project[m][j + P_from * i] /= interp_to[i + P_to * i];
277a76a04e7SJeremy L Thompson       }
278a76a04e7SJeremy L Thompson     }
27914556e63SJeremy L Thompson   }
28014556e63SJeremy L Thompson 
28114556e63SJeremy L Thompson   // Cleanup
2822b730f8bSJeremy L Thompson   CeedCall(CeedFree(&tau));
2832b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_to));
2842b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_from));
285a76a04e7SJeremy L Thompson 
286a76a04e7SJeremy L Thompson   return CEED_ERROR_SUCCESS;
287a76a04e7SJeremy L Thompson }
288a76a04e7SJeremy L Thompson 
2897a982d89SJeremy L. Thompson /// @}
2907a982d89SJeremy L. Thompson 
2917a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
2927a982d89SJeremy L. Thompson /// Ceed Backend API
2937a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
2947a982d89SJeremy L. Thompson /// @addtogroup CeedBasisBackend
2957a982d89SJeremy L. Thompson /// @{
2967a982d89SJeremy L. Thompson 
2977a982d89SJeremy L. Thompson /**
2987a982d89SJeremy L. Thompson   @brief Return collocated grad matrix
2997a982d89SJeremy L. Thompson 
300ea61e9acSJeremy L Thompson   @param[in]  basis         CeedBasis
301ea61e9acSJeremy L Thompson   @param[out] collo_grad_1d Row-major (Q_1d * Q_1d) matrix expressing derivatives of basis functions at quadrature points
3027a982d89SJeremy L. Thompson 
3037a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3047a982d89SJeremy L. Thompson 
3057a982d89SJeremy L. Thompson   @ref Backend
3067a982d89SJeremy L. Thompson **/
307d1d35e2fSjeremylt int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *collo_grad_1d) {
3087a982d89SJeremy L. Thompson   Ceed        ceed;
3092b730f8bSJeremy L Thompson   CeedInt     P_1d = (basis)->P_1d, Q_1d = (basis)->Q_1d;
31078464608Sjeremylt   CeedScalar *interp_1d, *grad_1d, *tau;
3117a982d89SJeremy L. Thompson 
3122b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q_1d * P_1d, &interp_1d));
3132b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q_1d * P_1d, &grad_1d));
3142b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q_1d, &tau));
315d1d35e2fSjeremylt   memcpy(interp_1d, (basis)->interp_1d, Q_1d * P_1d * sizeof(basis)->interp_1d[0]);
316d1d35e2fSjeremylt   memcpy(grad_1d, (basis)->grad_1d, Q_1d * P_1d * sizeof(basis)->interp_1d[0]);
3177a982d89SJeremy L. Thompson 
318d1d35e2fSjeremylt   // QR Factorization, interp_1d = Q R
3192b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetCeed(basis, &ceed));
3202b730f8bSJeremy L Thompson   CeedCall(CeedQRFactorization(ceed, interp_1d, tau, Q_1d, P_1d));
321ea61e9acSJeremy 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.
3227a982d89SJeremy L. Thompson 
323c8c3fa7dSJeremy L Thompson   // Apply R_inv, collo_grad_1d = grad_1d R_inv
324c8c3fa7dSJeremy L Thompson   for (CeedInt i = 0; i < Q_1d; i++) {  // Row i
325d1d35e2fSjeremylt     collo_grad_1d[Q_1d * i] = grad_1d[P_1d * i] / interp_1d[0];
326c8c3fa7dSJeremy L Thompson     for (CeedInt j = 1; j < P_1d; j++) {  // Column j
327d1d35e2fSjeremylt       collo_grad_1d[j + Q_1d * i] = grad_1d[j + P_1d * i];
328c8c3fa7dSJeremy L Thompson       for (CeedInt k = 0; k < j; k++) collo_grad_1d[j + Q_1d * i] -= interp_1d[j + P_1d * k] * collo_grad_1d[k + Q_1d * i];
329d1d35e2fSjeremylt       collo_grad_1d[j + Q_1d * i] /= interp_1d[j + P_1d * j];
3307a982d89SJeremy L. Thompson     }
331c8c3fa7dSJeremy L Thompson     for (CeedInt j = P_1d; j < Q_1d; j++) collo_grad_1d[j + Q_1d * i] = 0;
3327a982d89SJeremy L. Thompson   }
3337a982d89SJeremy L. Thompson 
33415ad3917SSebastian Grimberg   // Apply Q^T, collo_grad_1d = collo_grad_1d Q^T
3352b730f8bSJeremy L Thompson   CeedCall(CeedHouseholderApplyQ(collo_grad_1d, interp_1d, tau, CEED_NOTRANSPOSE, Q_1d, Q_1d, P_1d, 1, Q_1d));
3367a982d89SJeremy L. Thompson 
3372b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_1d));
3382b730f8bSJeremy L Thompson   CeedCall(CeedFree(&grad_1d));
3392b730f8bSJeremy L Thompson   CeedCall(CeedFree(&tau));
340e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3417a982d89SJeremy L. Thompson }
3427a982d89SJeremy L. Thompson 
3437a982d89SJeremy L. Thompson /**
3447a982d89SJeremy L. Thompson   @brief Get tensor status for given CeedBasis
3457a982d89SJeremy L. Thompson 
346ea61e9acSJeremy L Thompson   @param[in]  basis     CeedBasis
347d1d35e2fSjeremylt   @param[out] is_tensor Variable to store tensor status
3487a982d89SJeremy L. Thompson 
3497a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3507a982d89SJeremy L. Thompson 
3517a982d89SJeremy L. Thompson   @ref Backend
3527a982d89SJeremy L. Thompson **/
353d1d35e2fSjeremylt int CeedBasisIsTensor(CeedBasis basis, bool *is_tensor) {
3546402da51SJeremy L Thompson   *is_tensor = basis->is_tensor_basis;
355e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3567a982d89SJeremy L. Thompson }
3577a982d89SJeremy L. Thompson 
3587a982d89SJeremy L. Thompson /**
3597a982d89SJeremy L. Thompson   @brief Get backend data of a CeedBasis
3607a982d89SJeremy L. Thompson 
361ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
3627a982d89SJeremy L. Thompson   @param[out] data  Variable to store data
3637a982d89SJeremy L. Thompson 
3647a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3657a982d89SJeremy L. Thompson 
3667a982d89SJeremy L. Thompson   @ref Backend
3677a982d89SJeremy L. Thompson **/
368777ff853SJeremy L Thompson int CeedBasisGetData(CeedBasis basis, void *data) {
369777ff853SJeremy L Thompson   *(void **)data = basis->data;
370e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3717a982d89SJeremy L. Thompson }
3727a982d89SJeremy L. Thompson 
3737a982d89SJeremy L. Thompson /**
3747a982d89SJeremy L. Thompson   @brief Set backend data of a CeedBasis
3757a982d89SJeremy L. Thompson 
376ea61e9acSJeremy L Thompson   @param[in,out] basis  CeedBasis
377ea61e9acSJeremy L Thompson   @param[in]     data   Data to set
3787a982d89SJeremy L. Thompson 
3797a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3807a982d89SJeremy L. Thompson 
3817a982d89SJeremy L. Thompson   @ref Backend
3827a982d89SJeremy L. Thompson **/
383777ff853SJeremy L Thompson int CeedBasisSetData(CeedBasis basis, void *data) {
384777ff853SJeremy L Thompson   basis->data = data;
385e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3867a982d89SJeremy L. Thompson }
3877a982d89SJeremy L. Thompson 
3887a982d89SJeremy L. Thompson /**
38934359f16Sjeremylt   @brief Increment the reference counter for a CeedBasis
39034359f16Sjeremylt 
391ea61e9acSJeremy L Thompson   @param[in,out] basis Basis to increment the reference counter
39234359f16Sjeremylt 
39334359f16Sjeremylt   @return An error code: 0 - success, otherwise - failure
39434359f16Sjeremylt 
39534359f16Sjeremylt   @ref Backend
39634359f16Sjeremylt **/
3979560d06aSjeremylt int CeedBasisReference(CeedBasis basis) {
39834359f16Sjeremylt   basis->ref_count++;
39934359f16Sjeremylt   return CEED_ERROR_SUCCESS;
40034359f16Sjeremylt }
40134359f16Sjeremylt 
40234359f16Sjeremylt /**
403c4e3f59bSSebastian Grimberg   @brief Get number of Q-vector components for given CeedBasis
404c4e3f59bSSebastian Grimberg 
405c4e3f59bSSebastian Grimberg   @param[in]  basis  CeedBasis
406c4e3f59bSSebastian Grimberg   @param[in]  eval_mode \ref CEED_EVAL_INTERP to use interpolated values,
407c4e3f59bSSebastian Grimberg                           \ref CEED_EVAL_GRAD to use gradients,
408c4e3f59bSSebastian Grimberg                           \ref CEED_EVAL_DIV to use divergence,
409c4e3f59bSSebastian Grimberg                           \ref CEED_EVAL_CURL to use curl.
410c4e3f59bSSebastian Grimberg   @param[out] q_comp Variable to store number of Q-vector components of basis
411c4e3f59bSSebastian Grimberg 
412c4e3f59bSSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
413c4e3f59bSSebastian Grimberg 
414c4e3f59bSSebastian Grimberg   @ref Backend
415c4e3f59bSSebastian Grimberg **/
416c4e3f59bSSebastian Grimberg int CeedBasisGetNumQuadratureComponents(CeedBasis basis, CeedEvalMode eval_mode, CeedInt *q_comp) {
417c4e3f59bSSebastian Grimberg   switch (eval_mode) {
418c4e3f59bSSebastian Grimberg     case CEED_EVAL_INTERP:
419c4e3f59bSSebastian Grimberg       *q_comp = (basis->fe_space == CEED_FE_SPACE_H1) ? 1 : basis->dim;
420c4e3f59bSSebastian Grimberg       break;
421c4e3f59bSSebastian Grimberg     case CEED_EVAL_GRAD:
422c4e3f59bSSebastian Grimberg       *q_comp = basis->dim;
423c4e3f59bSSebastian Grimberg       break;
424c4e3f59bSSebastian Grimberg     case CEED_EVAL_DIV:
425c4e3f59bSSebastian Grimberg       *q_comp = 1;
426c4e3f59bSSebastian Grimberg       break;
427c4e3f59bSSebastian Grimberg     case CEED_EVAL_CURL:
428c4e3f59bSSebastian Grimberg       *q_comp = (basis->dim < 3) ? 1 : basis->dim;
429c4e3f59bSSebastian Grimberg       break;
430c4e3f59bSSebastian Grimberg     case CEED_EVAL_NONE:
431c4e3f59bSSebastian Grimberg     case CEED_EVAL_WEIGHT:
432352a5e7cSSebastian Grimberg       *q_comp = 1;
433c4e3f59bSSebastian Grimberg       break;
434c4e3f59bSSebastian Grimberg   }
435c4e3f59bSSebastian Grimberg   return CEED_ERROR_SUCCESS;
436c4e3f59bSSebastian Grimberg }
437c4e3f59bSSebastian Grimberg 
438c4e3f59bSSebastian Grimberg /**
4396e15d496SJeremy L Thompson   @brief Estimate number of FLOPs required to apply CeedBasis in t_mode and eval_mode
4406e15d496SJeremy L Thompson 
441ea61e9acSJeremy L Thompson   @param[in]  basis     Basis to estimate FLOPs for
442ea61e9acSJeremy L Thompson   @param[in]  t_mode    Apply basis or transpose
443ea61e9acSJeremy L Thompson   @param[in]  eval_mode Basis evaluation mode
444ea61e9acSJeremy L Thompson   @param[out] flops     Address of variable to hold FLOPs estimate
4456e15d496SJeremy L Thompson 
4466e15d496SJeremy L Thompson   @ref Backend
4476e15d496SJeremy L Thompson **/
4482b730f8bSJeremy L Thompson int CeedBasisGetFlopsEstimate(CeedBasis basis, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedSize *flops) {
4496e15d496SJeremy L Thompson   bool is_tensor;
4506e15d496SJeremy L Thompson 
4512b730f8bSJeremy L Thompson   CeedCall(CeedBasisIsTensor(basis, &is_tensor));
4526e15d496SJeremy L Thompson   if (is_tensor) {
4536e15d496SJeremy L Thompson     CeedInt dim, num_comp, P_1d, Q_1d;
4542b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetDimension(basis, &dim));
4552b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
4562b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d));
4572b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
4586e15d496SJeremy L Thompson     if (t_mode == CEED_TRANSPOSE) {
4592b730f8bSJeremy L Thompson       P_1d = Q_1d;
4602b730f8bSJeremy L Thompson       Q_1d = P_1d;
4616e15d496SJeremy L Thompson     }
4626e15d496SJeremy L Thompson     CeedInt tensor_flops = 0, pre = num_comp * CeedIntPow(P_1d, dim - 1), post = 1;
4636e15d496SJeremy L Thompson     for (CeedInt d = 0; d < dim; d++) {
4646e15d496SJeremy L Thompson       tensor_flops += 2 * pre * P_1d * post * Q_1d;
4656e15d496SJeremy L Thompson       pre /= P_1d;
4666e15d496SJeremy L Thompson       post *= Q_1d;
4676e15d496SJeremy L Thompson     }
4686e15d496SJeremy L Thompson     switch (eval_mode) {
4692b730f8bSJeremy L Thompson       case CEED_EVAL_NONE:
4702b730f8bSJeremy L Thompson         *flops = 0;
4712b730f8bSJeremy L Thompson         break;
4722b730f8bSJeremy L Thompson       case CEED_EVAL_INTERP:
4732b730f8bSJeremy L Thompson         *flops = tensor_flops;
4742b730f8bSJeremy L Thompson         break;
4752b730f8bSJeremy L Thompson       case CEED_EVAL_GRAD:
4762b730f8bSJeremy L Thompson         *flops = tensor_flops * 2;
4772b730f8bSJeremy L Thompson         break;
4786e15d496SJeremy L Thompson       case CEED_EVAL_DIV:
4796e15d496SJeremy L Thompson       case CEED_EVAL_CURL:
4806574a04fSJeremy L Thompson         // LCOV_EXCL_START
4816574a04fSJeremy L Thompson         return CeedError(basis->ceed, CEED_ERROR_INCOMPATIBLE, "Tensor basis evaluation for %s not supported", CeedEvalModes[eval_mode]);
4822b730f8bSJeremy L Thompson         break;
4836e15d496SJeremy L Thompson       // LCOV_EXCL_STOP
4842b730f8bSJeremy L Thompson       case CEED_EVAL_WEIGHT:
4852b730f8bSJeremy L Thompson         *flops = dim * CeedIntPow(Q_1d, dim);
4862b730f8bSJeremy L Thompson         break;
4876e15d496SJeremy L Thompson     }
4886e15d496SJeremy L Thompson   } else {
489c4e3f59bSSebastian Grimberg     CeedInt dim, num_comp, q_comp, num_nodes, num_qpts;
4902b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetDimension(basis, &dim));
4912b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
492c4e3f59bSSebastian Grimberg     CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp));
4932b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis, &num_nodes));
4942b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts));
4956e15d496SJeremy L Thompson     switch (eval_mode) {
4962b730f8bSJeremy L Thompson       case CEED_EVAL_NONE:
4972b730f8bSJeremy L Thompson         *flops = 0;
4982b730f8bSJeremy L Thompson         break;
4992b730f8bSJeremy L Thompson       case CEED_EVAL_INTERP:
5002b730f8bSJeremy L Thompson       case CEED_EVAL_GRAD:
5012b730f8bSJeremy L Thompson       case CEED_EVAL_DIV:
5022b730f8bSJeremy L Thompson       case CEED_EVAL_CURL:
503c4e3f59bSSebastian Grimberg         *flops = num_nodes * num_qpts * num_comp * q_comp;
5042b730f8bSJeremy L Thompson         break;
5052b730f8bSJeremy L Thompson       case CEED_EVAL_WEIGHT:
5062b730f8bSJeremy L Thompson         *flops = 0;
5072b730f8bSJeremy L Thompson         break;
5086e15d496SJeremy L Thompson     }
5096e15d496SJeremy L Thompson   }
5106e15d496SJeremy L Thompson 
5116e15d496SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5126e15d496SJeremy L Thompson }
5136e15d496SJeremy L Thompson 
5146e15d496SJeremy L Thompson /**
515c4e3f59bSSebastian Grimberg   @brief Get CeedFESpace for a CeedBasis
516c4e3f59bSSebastian Grimberg 
517c4e3f59bSSebastian Grimberg   @param[in]  basis     CeedBasis
518c4e3f59bSSebastian Grimberg   @param[out] fe_space  Variable to store CeedFESpace
519c4e3f59bSSebastian Grimberg 
520c4e3f59bSSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
521c4e3f59bSSebastian Grimberg 
522c4e3f59bSSebastian Grimberg   @ref Backend
523c4e3f59bSSebastian Grimberg **/
524c4e3f59bSSebastian Grimberg int CeedBasisGetFESpace(CeedBasis basis, CeedFESpace *fe_space) {
525c4e3f59bSSebastian Grimberg   *fe_space = basis->fe_space;
526c4e3f59bSSebastian Grimberg   return CEED_ERROR_SUCCESS;
527c4e3f59bSSebastian Grimberg }
528c4e3f59bSSebastian Grimberg 
529c4e3f59bSSebastian Grimberg /**
5307a982d89SJeremy L. Thompson   @brief Get dimension for given CeedElemTopology
5317a982d89SJeremy L. Thompson 
532ea61e9acSJeremy L Thompson   @param[in]  topo CeedElemTopology
5337a982d89SJeremy L. Thompson   @param[out] dim  Variable to store dimension of topology
5347a982d89SJeremy L. Thompson 
5357a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
5367a982d89SJeremy L. Thompson 
5377a982d89SJeremy L. Thompson   @ref Backend
5387a982d89SJeremy L. Thompson **/
5397a982d89SJeremy L. Thompson int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) {
5407a982d89SJeremy L. Thompson   *dim = (CeedInt)topo >> 16;
541e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5427a982d89SJeremy L. Thompson }
5437a982d89SJeremy L. Thompson 
5447a982d89SJeremy L. Thompson /**
5457a982d89SJeremy L. Thompson   @brief Get CeedTensorContract of a CeedBasis
5467a982d89SJeremy L. Thompson 
547ea61e9acSJeremy L Thompson   @param[in]  basis     CeedBasis
5487a982d89SJeremy L. Thompson   @param[out] contract  Variable to store CeedTensorContract
5497a982d89SJeremy L. Thompson 
5507a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
5517a982d89SJeremy L. Thompson 
5527a982d89SJeremy L. Thompson   @ref Backend
5537a982d89SJeremy L. Thompson **/
5547a982d89SJeremy L. Thompson int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) {
5557a982d89SJeremy L. Thompson   *contract = basis->contract;
556e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5577a982d89SJeremy L. Thompson }
5587a982d89SJeremy L. Thompson 
5597a982d89SJeremy L. Thompson /**
5607a982d89SJeremy L. Thompson   @brief Set CeedTensorContract of a CeedBasis
5617a982d89SJeremy L. Thompson 
562ea61e9acSJeremy L Thompson   @param[in,out] basis    CeedBasis
563ea61e9acSJeremy L Thompson   @param[in]     contract CeedTensorContract to set
5647a982d89SJeremy L. Thompson 
5657a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
5667a982d89SJeremy L. Thompson 
5677a982d89SJeremy L. Thompson   @ref Backend
5687a982d89SJeremy L. Thompson **/
56934359f16Sjeremylt int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract contract) {
57034359f16Sjeremylt   basis->contract = contract;
5712b730f8bSJeremy L Thompson   CeedCall(CeedTensorContractReference(contract));
572e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5737a982d89SJeremy L. Thompson }
5747a982d89SJeremy L. Thompson 
5757a982d89SJeremy L. Thompson /**
5767a982d89SJeremy L. Thompson   @brief Return a reference implementation of matrix multiplication C = A B.
577ba59ac12SSebastian Grimberg 
578ba59ac12SSebastian Grimberg   Note: This is a reference implementation for CPU CeedScalar pointers that is not intended for high performance.
5797a982d89SJeremy L. Thompson 
580ea61e9acSJeremy L Thompson   @param[in]  ceed  Ceed context for error handling
581d1d35e2fSjeremylt   @param[in]  mat_A Row-major matrix A
582d1d35e2fSjeremylt   @param[in]  mat_B Row-major matrix B
583d1d35e2fSjeremylt   @param[out] mat_C Row-major output matrix C
584ea61e9acSJeremy L Thompson   @param[in]  m     Number of rows of C
585ea61e9acSJeremy L Thompson   @param[in]  n     Number of columns of C
586ea61e9acSJeremy L Thompson   @param[in]  kk    Number of columns of A/rows of B
5877a982d89SJeremy L. Thompson 
5887a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
5897a982d89SJeremy L. Thompson 
5907a982d89SJeremy L. Thompson   @ref Utility
5917a982d89SJeremy L. Thompson **/
5922b730f8bSJeremy L Thompson int CeedMatrixMatrixMultiply(Ceed ceed, const CeedScalar *mat_A, const CeedScalar *mat_B, CeedScalar *mat_C, CeedInt m, CeedInt n, CeedInt kk) {
5932b730f8bSJeremy L Thompson   for (CeedInt i = 0; i < m; i++) {
5947a982d89SJeremy L. Thompson     for (CeedInt j = 0; j < n; j++) {
5957a982d89SJeremy L. Thompson       CeedScalar sum = 0;
5962b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < kk; k++) sum += mat_A[k + i * kk] * mat_B[j + k * n];
597d1d35e2fSjeremylt       mat_C[j + i * n] = sum;
5987a982d89SJeremy L. Thompson     }
5992b730f8bSJeremy L Thompson   }
600e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
6017a982d89SJeremy L. Thompson }
6027a982d89SJeremy L. Thompson 
603ba59ac12SSebastian Grimberg /**
604ba59ac12SSebastian Grimberg   @brief Return QR Factorization of a matrix
605ba59ac12SSebastian Grimberg 
606ba59ac12SSebastian Grimberg   @param[in]     ceed Ceed context for error handling
607ba59ac12SSebastian Grimberg   @param[in,out] mat  Row-major matrix to be factorized in place
608ba59ac12SSebastian Grimberg   @param[in,out] tau  Vector of length m of scaling factors
609ba59ac12SSebastian Grimberg   @param[in]     m    Number of rows
610ba59ac12SSebastian Grimberg   @param[in]     n    Number of columns
611ba59ac12SSebastian Grimberg 
612ba59ac12SSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
613ba59ac12SSebastian Grimberg 
614ba59ac12SSebastian Grimberg   @ref Utility
615ba59ac12SSebastian Grimberg **/
616ba59ac12SSebastian Grimberg int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau, CeedInt m, CeedInt n) {
617ba59ac12SSebastian Grimberg   CeedScalar v[m];
618ba59ac12SSebastian Grimberg 
619ba59ac12SSebastian Grimberg   // Check matrix shape
6206574a04fSJeremy L Thompson   CeedCheck(n <= m, ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute QR factorization with n > m");
621ba59ac12SSebastian Grimberg 
622ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) {
623ba59ac12SSebastian Grimberg     if (i >= m - 1) {  // last row of matrix, no reflection needed
624ba59ac12SSebastian Grimberg       tau[i] = 0.;
625ba59ac12SSebastian Grimberg       break;
626ba59ac12SSebastian Grimberg     }
627ba59ac12SSebastian Grimberg     // Calculate Householder vector, magnitude
628ba59ac12SSebastian Grimberg     CeedScalar sigma = 0.0;
629ba59ac12SSebastian Grimberg     v[i]             = mat[i + n * i];
630ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < m; j++) {
631ba59ac12SSebastian Grimberg       v[j] = mat[i + n * j];
632ba59ac12SSebastian Grimberg       sigma += v[j] * v[j];
633ba59ac12SSebastian Grimberg     }
634ba59ac12SSebastian Grimberg     CeedScalar norm = sqrt(v[i] * v[i] + sigma);  // norm of v[i:m]
635ba59ac12SSebastian Grimberg     CeedScalar R_ii = -copysign(norm, v[i]);
636ba59ac12SSebastian Grimberg     v[i] -= R_ii;
637ba59ac12SSebastian Grimberg     // norm of v[i:m] after modification above and scaling below
638ba59ac12SSebastian Grimberg     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
639ba59ac12SSebastian Grimberg     //   tau = 2 / (norm*norm)
640ba59ac12SSebastian Grimberg     tau[i] = 2 * v[i] * v[i] / (v[i] * v[i] + sigma);
641ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < m; j++) v[j] /= v[i];
642ba59ac12SSebastian Grimberg 
643ba59ac12SSebastian Grimberg     // Apply Householder reflector to lower right panel
644ba59ac12SSebastian Grimberg     CeedHouseholderReflect(&mat[i * n + i + 1], &v[i], tau[i], m - i, n - i - 1, n, 1);
645ba59ac12SSebastian Grimberg     // Save v
646ba59ac12SSebastian Grimberg     mat[i + n * i] = R_ii;
647ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < m; j++) mat[i + n * j] = v[j];
648ba59ac12SSebastian Grimberg   }
649ba59ac12SSebastian Grimberg   return CEED_ERROR_SUCCESS;
650ba59ac12SSebastian Grimberg }
651ba59ac12SSebastian Grimberg 
652ba59ac12SSebastian Grimberg /**
653ba59ac12SSebastian Grimberg   @brief Apply Householder Q matrix
654ba59ac12SSebastian Grimberg 
655ba59ac12SSebastian Grimberg   Compute mat_A = mat_Q mat_A, where mat_Q is mxm and mat_A is mxn.
656ba59ac12SSebastian Grimberg 
657ba59ac12SSebastian Grimberg   @param[in,out] mat_A  Matrix to apply Householder Q to, in place
658ba59ac12SSebastian Grimberg   @param[in]     mat_Q  Householder Q matrix
659ba59ac12SSebastian Grimberg   @param[in]     tau    Householder scaling factors
660ba59ac12SSebastian Grimberg   @param[in]     t_mode Transpose mode for application
661ba59ac12SSebastian Grimberg   @param[in]     m      Number of rows in A
662ba59ac12SSebastian Grimberg   @param[in]     n      Number of columns in A
663ba59ac12SSebastian Grimberg   @param[in]     k      Number of elementary reflectors in Q, k<m
664ba59ac12SSebastian Grimberg   @param[in]     row    Row stride in A
665ba59ac12SSebastian Grimberg   @param[in]     col    Col stride in A
666ba59ac12SSebastian Grimberg 
667ba59ac12SSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
668ba59ac12SSebastian Grimberg 
669c4e3f59bSSebastian Grimberg   @ref Utility
670ba59ac12SSebastian Grimberg **/
671ba59ac12SSebastian Grimberg int CeedHouseholderApplyQ(CeedScalar *mat_A, const CeedScalar *mat_Q, const CeedScalar *tau, CeedTransposeMode t_mode, CeedInt m, CeedInt n,
672ba59ac12SSebastian Grimberg                           CeedInt k, CeedInt row, CeedInt col) {
673ba59ac12SSebastian Grimberg   CeedScalar *v;
674ba59ac12SSebastian Grimberg   CeedCall(CeedMalloc(m, &v));
675ba59ac12SSebastian Grimberg   for (CeedInt ii = 0; ii < k; ii++) {
676ba59ac12SSebastian Grimberg     CeedInt i = t_mode == CEED_TRANSPOSE ? ii : k - 1 - ii;
677ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < m; j++) v[j] = mat_Q[j * k + i];
678ba59ac12SSebastian Grimberg     // Apply Householder reflector (I - tau v v^T) collo_grad_1d^T
679ba59ac12SSebastian Grimberg     CeedCall(CeedHouseholderReflect(&mat_A[i * row], &v[i], tau[i], m - i, n, row, col));
680ba59ac12SSebastian Grimberg   }
681ba59ac12SSebastian Grimberg   CeedCall(CeedFree(&v));
682ba59ac12SSebastian Grimberg   return CEED_ERROR_SUCCESS;
683ba59ac12SSebastian Grimberg }
684ba59ac12SSebastian Grimberg 
685ba59ac12SSebastian Grimberg /**
686ba59ac12SSebastian Grimberg   @brief Return symmetric Schur decomposition of the symmetric matrix mat via symmetric QR factorization
687ba59ac12SSebastian Grimberg 
688ba59ac12SSebastian Grimberg   @param[in]     ceed   Ceed context for error handling
689ba59ac12SSebastian Grimberg   @param[in,out] mat    Row-major matrix to be factorized in place
690ba59ac12SSebastian Grimberg   @param[out]    lambda Vector of length n of eigenvalues
691ba59ac12SSebastian Grimberg   @param[in]     n      Number of rows/columns
692ba59ac12SSebastian Grimberg 
693ba59ac12SSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
694ba59ac12SSebastian Grimberg 
695ba59ac12SSebastian Grimberg   @ref Utility
696ba59ac12SSebastian Grimberg **/
6972c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOff
6982c2ea1dbSJeremy L Thompson int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat, CeedScalar *lambda, CeedInt n) {
699ba59ac12SSebastian Grimberg   // Check bounds for clang-tidy
7006574a04fSJeremy L Thompson   CeedCheck(n > 1, ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute symmetric Schur decomposition of scalars");
701ba59ac12SSebastian Grimberg 
702ba59ac12SSebastian Grimberg   CeedScalar v[n - 1], tau[n - 1], mat_T[n * n];
703ba59ac12SSebastian Grimberg 
704ba59ac12SSebastian Grimberg   // Copy mat to mat_T and set mat to I
705ba59ac12SSebastian Grimberg   memcpy(mat_T, mat, n * n * sizeof(mat[0]));
706ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) {
707ba59ac12SSebastian Grimberg     for (CeedInt j = 0; j < n; j++) mat[j + n * i] = (i == j) ? 1 : 0;
708ba59ac12SSebastian Grimberg   }
709ba59ac12SSebastian Grimberg 
710ba59ac12SSebastian Grimberg   // Reduce to tridiagonal
711ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n - 1; i++) {
712ba59ac12SSebastian Grimberg     // Calculate Householder vector, magnitude
713ba59ac12SSebastian Grimberg     CeedScalar sigma = 0.0;
714ba59ac12SSebastian Grimberg     v[i]             = mat_T[i + n * (i + 1)];
715ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < n - 1; j++) {
716ba59ac12SSebastian Grimberg       v[j] = mat_T[i + n * (j + 1)];
717ba59ac12SSebastian Grimberg       sigma += v[j] * v[j];
718ba59ac12SSebastian Grimberg     }
719ba59ac12SSebastian Grimberg     CeedScalar norm = sqrt(v[i] * v[i] + sigma);  // norm of v[i:n-1]
720ba59ac12SSebastian Grimberg     CeedScalar R_ii = -copysign(norm, v[i]);
721ba59ac12SSebastian Grimberg     v[i] -= R_ii;
722ba59ac12SSebastian Grimberg     // norm of v[i:m] after modification above and scaling below
723ba59ac12SSebastian Grimberg     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
724ba59ac12SSebastian Grimberg     //   tau = 2 / (norm*norm)
725ba59ac12SSebastian Grimberg     tau[i] = i == n - 2 ? 2 : 2 * v[i] * v[i] / (v[i] * v[i] + sigma);
726ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < n - 1; j++) v[j] /= v[i];
727ba59ac12SSebastian Grimberg 
728ba59ac12SSebastian Grimberg     // Update sub and super diagonal
729ba59ac12SSebastian Grimberg     for (CeedInt j = i + 2; j < n; j++) {
730ba59ac12SSebastian Grimberg       mat_T[i + n * j] = 0;
731ba59ac12SSebastian Grimberg       mat_T[j + n * i] = 0;
732ba59ac12SSebastian Grimberg     }
733ba59ac12SSebastian Grimberg     // Apply symmetric Householder reflector to lower right panel
734ba59ac12SSebastian Grimberg     CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1);
735ba59ac12SSebastian Grimberg     CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), 1, n);
736ba59ac12SSebastian Grimberg 
737ba59ac12SSebastian Grimberg     // Save v
738ba59ac12SSebastian Grimberg     mat_T[i + n * (i + 1)] = R_ii;
739ba59ac12SSebastian Grimberg     mat_T[(i + 1) + n * i] = R_ii;
740ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < n - 1; j++) {
741ba59ac12SSebastian Grimberg       mat_T[i + n * (j + 1)] = v[j];
742ba59ac12SSebastian Grimberg     }
743ba59ac12SSebastian Grimberg   }
744ba59ac12SSebastian Grimberg   // Backwards accumulation of Q
745ba59ac12SSebastian Grimberg   for (CeedInt i = n - 2; i >= 0; i--) {
746ba59ac12SSebastian Grimberg     if (tau[i] > 0.0) {
747ba59ac12SSebastian Grimberg       v[i] = 1;
748ba59ac12SSebastian Grimberg       for (CeedInt j = i + 1; j < n - 1; j++) {
749ba59ac12SSebastian Grimberg         v[j]                   = mat_T[i + n * (j + 1)];
750ba59ac12SSebastian Grimberg         mat_T[i + n * (j + 1)] = 0;
751ba59ac12SSebastian Grimberg       }
752ba59ac12SSebastian Grimberg       CeedHouseholderReflect(&mat[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1);
753ba59ac12SSebastian Grimberg     }
754ba59ac12SSebastian Grimberg   }
755ba59ac12SSebastian Grimberg 
756ba59ac12SSebastian Grimberg   // Reduce sub and super diagonal
757ba59ac12SSebastian Grimberg   CeedInt    p = 0, q = 0, itr = 0, max_itr = n * n * n * n;
758ba59ac12SSebastian Grimberg   CeedScalar tol = CEED_EPSILON;
759ba59ac12SSebastian Grimberg 
760ba59ac12SSebastian Grimberg   while (itr < max_itr) {
761ba59ac12SSebastian Grimberg     // Update p, q, size of reduced portions of diagonal
762ba59ac12SSebastian Grimberg     p = 0;
763ba59ac12SSebastian Grimberg     q = 0;
764ba59ac12SSebastian Grimberg     for (CeedInt i = n - 2; i >= 0; i--) {
765ba59ac12SSebastian Grimberg       if (fabs(mat_T[i + n * (i + 1)]) < tol) q += 1;
766ba59ac12SSebastian Grimberg       else break;
767ba59ac12SSebastian Grimberg     }
768ba59ac12SSebastian Grimberg     for (CeedInt i = 0; i < n - q - 1; i++) {
769ba59ac12SSebastian Grimberg       if (fabs(mat_T[i + n * (i + 1)]) < tol) p += 1;
770ba59ac12SSebastian Grimberg       else break;
771ba59ac12SSebastian Grimberg     }
772ba59ac12SSebastian Grimberg     if (q == n - 1) break;  // Finished reducing
773ba59ac12SSebastian Grimberg 
774ba59ac12SSebastian Grimberg     // Reduce tridiagonal portion
775ba59ac12SSebastian Grimberg     CeedScalar t_nn = mat_T[(n - 1 - q) + n * (n - 1 - q)], t_nnm1 = mat_T[(n - 2 - q) + n * (n - 1 - q)];
776ba59ac12SSebastian Grimberg     CeedScalar d  = (mat_T[(n - 2 - q) + n * (n - 2 - q)] - t_nn) / 2;
777ba59ac12SSebastian Grimberg     CeedScalar mu = t_nn - t_nnm1 * t_nnm1 / (d + copysign(sqrt(d * d + t_nnm1 * t_nnm1), d));
778ba59ac12SSebastian Grimberg     CeedScalar x  = mat_T[p + n * p] - mu;
779ba59ac12SSebastian Grimberg     CeedScalar z  = mat_T[p + n * (p + 1)];
780ba59ac12SSebastian Grimberg     for (CeedInt k = p; k < n - q - 1; k++) {
781ba59ac12SSebastian Grimberg       // Compute Givens rotation
782ba59ac12SSebastian Grimberg       CeedScalar c = 1, s = 0;
783ba59ac12SSebastian Grimberg       if (fabs(z) > tol) {
784ba59ac12SSebastian Grimberg         if (fabs(z) > fabs(x)) {
785ba59ac12SSebastian Grimberg           CeedScalar tau = -x / z;
786ba59ac12SSebastian Grimberg           s = 1 / sqrt(1 + tau * tau), c = s * tau;
787ba59ac12SSebastian Grimberg         } else {
788ba59ac12SSebastian Grimberg           CeedScalar tau = -z / x;
789ba59ac12SSebastian Grimberg           c = 1 / sqrt(1 + tau * tau), s = c * tau;
790ba59ac12SSebastian Grimberg         }
791ba59ac12SSebastian Grimberg       }
792ba59ac12SSebastian Grimberg 
793ba59ac12SSebastian Grimberg       // Apply Givens rotation to T
794ba59ac12SSebastian Grimberg       CeedGivensRotation(mat_T, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n);
795ba59ac12SSebastian Grimberg       CeedGivensRotation(mat_T, c, s, CEED_TRANSPOSE, k, k + 1, n, n);
796ba59ac12SSebastian Grimberg 
797ba59ac12SSebastian Grimberg       // Apply Givens rotation to Q
798ba59ac12SSebastian Grimberg       CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n);
799ba59ac12SSebastian Grimberg 
800ba59ac12SSebastian Grimberg       // Update x, z
801ba59ac12SSebastian Grimberg       if (k < n - q - 2) {
802ba59ac12SSebastian Grimberg         x = mat_T[k + n * (k + 1)];
803ba59ac12SSebastian Grimberg         z = mat_T[k + n * (k + 2)];
804ba59ac12SSebastian Grimberg       }
805ba59ac12SSebastian Grimberg     }
806ba59ac12SSebastian Grimberg     itr++;
807ba59ac12SSebastian Grimberg   }
808ba59ac12SSebastian Grimberg 
809ba59ac12SSebastian Grimberg   // Save eigenvalues
810ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) lambda[i] = mat_T[i + n * i];
811ba59ac12SSebastian Grimberg 
812ba59ac12SSebastian Grimberg   // Check convergence
8136574a04fSJeremy L Thompson   CeedCheck(itr < max_itr || q > n, ceed, CEED_ERROR_MINOR, "Symmetric QR failed to converge");
814ba59ac12SSebastian Grimberg   return CEED_ERROR_SUCCESS;
815ba59ac12SSebastian Grimberg }
8162c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOn
817ba59ac12SSebastian Grimberg 
818ba59ac12SSebastian Grimberg /**
819ba59ac12SSebastian Grimberg   @brief Return Simultaneous Diagonalization of two matrices.
820ba59ac12SSebastian Grimberg 
821ba59ac12SSebastian Grimberg   This solves the generalized eigenvalue problem A x = lambda B x, where A and B are symmetric and B is positive definite.
822ba59ac12SSebastian Grimberg   We generate the matrix X and vector Lambda such that X^T A X = Lambda and X^T B X = I.
823ba59ac12SSebastian Grimberg   This is equivalent to the LAPACK routine 'sygv' with TYPE = 1.
824ba59ac12SSebastian Grimberg 
825ba59ac12SSebastian Grimberg   @param[in]  ceed   Ceed context for error handling
826ba59ac12SSebastian Grimberg   @param[in]  mat_A  Row-major matrix to be factorized with eigenvalues
827ba59ac12SSebastian Grimberg   @param[in]  mat_B  Row-major matrix to be factorized to identity
828ba59ac12SSebastian Grimberg   @param[out] mat_X  Row-major orthogonal matrix
829ba59ac12SSebastian Grimberg   @param[out] lambda Vector of length n of generalized eigenvalues
830ba59ac12SSebastian Grimberg   @param[in]  n      Number of rows/columns
831ba59ac12SSebastian Grimberg 
832ba59ac12SSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
833ba59ac12SSebastian Grimberg 
834ba59ac12SSebastian Grimberg   @ref Utility
835ba59ac12SSebastian Grimberg **/
8362c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOff
8372c2ea1dbSJeremy L Thompson int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *mat_A, CeedScalar *mat_B, CeedScalar *mat_X, CeedScalar *lambda, CeedInt n) {
838ba59ac12SSebastian Grimberg   CeedScalar *mat_C, *mat_G, *vec_D;
839ba59ac12SSebastian Grimberg   CeedCall(CeedCalloc(n * n, &mat_C));
840ba59ac12SSebastian Grimberg   CeedCall(CeedCalloc(n * n, &mat_G));
841ba59ac12SSebastian Grimberg   CeedCall(CeedCalloc(n, &vec_D));
842ba59ac12SSebastian Grimberg 
843ba59ac12SSebastian Grimberg   // Compute B = G D G^T
844ba59ac12SSebastian Grimberg   memcpy(mat_G, mat_B, n * n * sizeof(mat_B[0]));
845ba59ac12SSebastian Grimberg   CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_G, vec_D, n));
846ba59ac12SSebastian Grimberg 
847ba59ac12SSebastian Grimberg   // Sort eigenvalues
848ba59ac12SSebastian Grimberg   for (CeedInt i = n - 1; i >= 0; i--) {
849ba59ac12SSebastian Grimberg     for (CeedInt j = 0; j < i; j++) {
850ba59ac12SSebastian Grimberg       if (fabs(vec_D[j]) > fabs(vec_D[j + 1])) {
851ba59ac12SSebastian Grimberg         CeedScalar temp;
852ba59ac12SSebastian Grimberg         temp         = vec_D[j];
853ba59ac12SSebastian Grimberg         vec_D[j]     = vec_D[j + 1];
854ba59ac12SSebastian Grimberg         vec_D[j + 1] = temp;
855ba59ac12SSebastian Grimberg         for (CeedInt k = 0; k < n; k++) {
856ba59ac12SSebastian Grimberg           temp                 = mat_G[k * n + j];
857ba59ac12SSebastian Grimberg           mat_G[k * n + j]     = mat_G[k * n + j + 1];
858ba59ac12SSebastian Grimberg           mat_G[k * n + j + 1] = temp;
859ba59ac12SSebastian Grimberg         }
860ba59ac12SSebastian Grimberg       }
861ba59ac12SSebastian Grimberg     }
862ba59ac12SSebastian Grimberg   }
863ba59ac12SSebastian Grimberg 
864ba59ac12SSebastian Grimberg   // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T
865ba59ac12SSebastian Grimberg   //           = D^-1/2 G^T A G D^-1/2
866ba59ac12SSebastian Grimberg   // -- D = D^-1/2
867ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) vec_D[i] = 1. / sqrt(vec_D[i]);
868ba59ac12SSebastian Grimberg   // -- G = G D^-1/2
869ba59ac12SSebastian Grimberg   // -- C = D^-1/2 G^T
870ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) {
871ba59ac12SSebastian Grimberg     for (CeedInt j = 0; j < n; j++) {
872ba59ac12SSebastian Grimberg       mat_G[i * n + j] *= vec_D[j];
873ba59ac12SSebastian Grimberg       mat_C[j * n + i] = mat_G[i * n + j];
874ba59ac12SSebastian Grimberg     }
875ba59ac12SSebastian Grimberg   }
876ba59ac12SSebastian Grimberg   // -- X = (D^-1/2 G^T) A
877ba59ac12SSebastian Grimberg   CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_C, (const CeedScalar *)mat_A, mat_X, n, n, n));
878ba59ac12SSebastian Grimberg   // -- C = (D^-1/2 G^T A) (G D^-1/2)
879ba59ac12SSebastian Grimberg   CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_X, (const CeedScalar *)mat_G, mat_C, n, n, n));
880ba59ac12SSebastian Grimberg 
881ba59ac12SSebastian Grimberg   // Compute Q^T C Q = lambda
882ba59ac12SSebastian Grimberg   CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_C, lambda, n));
883ba59ac12SSebastian Grimberg 
884ba59ac12SSebastian Grimberg   // Sort eigenvalues
885ba59ac12SSebastian Grimberg   for (CeedInt i = n - 1; i >= 0; i--) {
886ba59ac12SSebastian Grimberg     for (CeedInt j = 0; j < i; j++) {
887ba59ac12SSebastian Grimberg       if (fabs(lambda[j]) > fabs(lambda[j + 1])) {
888ba59ac12SSebastian Grimberg         CeedScalar temp;
889ba59ac12SSebastian Grimberg         temp          = lambda[j];
890ba59ac12SSebastian Grimberg         lambda[j]     = lambda[j + 1];
891ba59ac12SSebastian Grimberg         lambda[j + 1] = temp;
892ba59ac12SSebastian Grimberg         for (CeedInt k = 0; k < n; k++) {
893ba59ac12SSebastian Grimberg           temp                 = mat_C[k * n + j];
894ba59ac12SSebastian Grimberg           mat_C[k * n + j]     = mat_C[k * n + j + 1];
895ba59ac12SSebastian Grimberg           mat_C[k * n + j + 1] = temp;
896ba59ac12SSebastian Grimberg         }
897ba59ac12SSebastian Grimberg       }
898ba59ac12SSebastian Grimberg     }
899ba59ac12SSebastian Grimberg   }
900ba59ac12SSebastian Grimberg 
901ba59ac12SSebastian Grimberg   // Set X = (G D^1/2)^-T Q
902ba59ac12SSebastian Grimberg   //       = G D^-1/2 Q
903ba59ac12SSebastian Grimberg   CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_G, (const CeedScalar *)mat_C, mat_X, n, n, n));
904ba59ac12SSebastian Grimberg 
905ba59ac12SSebastian Grimberg   // Cleanup
906ba59ac12SSebastian Grimberg   CeedCall(CeedFree(&mat_C));
907ba59ac12SSebastian Grimberg   CeedCall(CeedFree(&mat_G));
908ba59ac12SSebastian Grimberg   CeedCall(CeedFree(&vec_D));
909ba59ac12SSebastian Grimberg   return CEED_ERROR_SUCCESS;
910ba59ac12SSebastian Grimberg }
9112c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOn
912ba59ac12SSebastian Grimberg 
9137a982d89SJeremy L. Thompson /// @}
9147a982d89SJeremy L. Thompson 
9157a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
9167a982d89SJeremy L. Thompson /// CeedBasis Public API
9177a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
9187a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser
919d7b241e6Sjeremylt /// @{
920d7b241e6Sjeremylt 
921b11c1e72Sjeremylt /**
922ba59ac12SSebastian Grimberg   @brief Create a tensor-product basis for H^1 discretizations
923b11c1e72Sjeremylt 
924ea61e9acSJeremy L Thompson   @param[in]  ceed        Ceed object where the CeedBasis will be created
925ea61e9acSJeremy L Thompson   @param[in]  dim         Topological dimension
926ea61e9acSJeremy L Thompson   @param[in]  num_comp    Number of field components (1 for scalar fields)
927ea61e9acSJeremy L Thompson   @param[in]  P_1d        Number of nodes in one dimension
928ea61e9acSJeremy L Thompson   @param[in]  Q_1d        Number of quadrature points in one dimension
929ea61e9acSJeremy L Thompson   @param[in]  interp_1d   Row-major (Q_1d * P_1d) matrix expressing the values of nodal basis functions at quadrature points
930ea61e9acSJeremy L Thompson   @param[in]  grad_1d     Row-major (Q_1d * P_1d) matrix expressing derivatives of nodal basis functions at quadrature points
931ea61e9acSJeremy 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]
932ea61e9acSJeremy L Thompson   @param[in]  q_weight_1d Array of length Q_1d holding the quadrature weights on the reference element
933ea61e9acSJeremy L Thompson   @param[out] basis       Address of the variable where the newly created CeedBasis will be stored.
934b11c1e72Sjeremylt 
935b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
936dfdf5a53Sjeremylt 
9377a982d89SJeremy L. Thompson   @ref User
938b11c1e72Sjeremylt **/
9392b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d,
9402b730f8bSJeremy L Thompson                             const CeedScalar *grad_1d, const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis *basis) {
9415fe0d4faSjeremylt   if (!ceed->BasisCreateTensorH1) {
9425fe0d4faSjeremylt     Ceed delegate;
9436574a04fSJeremy L Thompson 
9442b730f8bSJeremy L Thompson     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
9456574a04fSJeremy L Thompson     CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisCreateTensorH1");
9462b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateTensorH1(delegate, dim, num_comp, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis));
947e15f9bd0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
9485fe0d4faSjeremylt   }
949e15f9bd0SJeremy L Thompson 
9506574a04fSJeremy L Thompson   CeedCheck(dim > 0, ceed, CEED_ERROR_DIMENSION, "Basis dimension must be a positive value");
9516574a04fSJeremy L Thompson   CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component");
9526574a04fSJeremy L Thompson   CeedCheck(P_1d > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node");
9536574a04fSJeremy L Thompson   CeedCheck(Q_1d > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point");
954227444bfSJeremy L Thompson 
9552b730f8bSJeremy L Thompson   CeedElemTopology topo = dim == 1 ? CEED_TOPOLOGY_LINE : dim == 2 ? CEED_TOPOLOGY_QUAD : CEED_TOPOLOGY_HEX;
956e15f9bd0SJeremy L Thompson 
9572b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(1, basis));
958db002c03SJeremy L Thompson   CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed));
959d1d35e2fSjeremylt   (*basis)->ref_count       = 1;
9606402da51SJeremy L Thompson   (*basis)->is_tensor_basis = true;
961d7b241e6Sjeremylt   (*basis)->dim             = dim;
962d99fa3c5SJeremy L Thompson   (*basis)->topo            = topo;
963d1d35e2fSjeremylt   (*basis)->num_comp        = num_comp;
964d1d35e2fSjeremylt   (*basis)->P_1d            = P_1d;
965d1d35e2fSjeremylt   (*basis)->Q_1d            = Q_1d;
966d1d35e2fSjeremylt   (*basis)->P               = CeedIntPow(P_1d, dim);
967d1d35e2fSjeremylt   (*basis)->Q               = CeedIntPow(Q_1d, dim);
968c4e3f59bSSebastian Grimberg   (*basis)->fe_space        = CEED_FE_SPACE_H1;
9692b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d, &(*basis)->q_ref_1d));
9702b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d, &(*basis)->q_weight_1d));
971ff3a0f91SJeremy L Thompson   if (q_ref_1d) memcpy((*basis)->q_ref_1d, q_ref_1d, Q_1d * sizeof(q_ref_1d[0]));
9722b730f8bSJeremy L Thompson   if (q_weight_1d) memcpy((*basis)->q_weight_1d, q_weight_1d, Q_1d * sizeof(q_weight_1d[0]));
9732b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->interp_1d));
9742b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->grad_1d));
9752b730f8bSJeremy L Thompson   if (interp_1d) memcpy((*basis)->interp_1d, interp_1d, Q_1d * P_1d * sizeof(interp_1d[0]));
976ff3a0f91SJeremy L Thompson   if (grad_1d) memcpy((*basis)->grad_1d, grad_1d, Q_1d * P_1d * sizeof(grad_1d[0]));
9772b730f8bSJeremy L Thompson   CeedCall(ceed->BasisCreateTensorH1(dim, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, *basis));
978e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
979d7b241e6Sjeremylt }
980d7b241e6Sjeremylt 
981b11c1e72Sjeremylt /**
98295bb1877Svaleriabarra   @brief Create a tensor-product Lagrange basis
983b11c1e72Sjeremylt 
984ea61e9acSJeremy L Thompson   @param[in]  ceed      Ceed object where the CeedBasis will be created
985ea61e9acSJeremy L Thompson   @param[in]  dim       Topological dimension of element
986ea61e9acSJeremy L Thompson   @param[in]  num_comp  Number of field components (1 for scalar fields)
987ea61e9acSJeremy L Thompson   @param[in]  P         Number of Gauss-Lobatto nodes in one dimension.
988ea61e9acSJeremy L Thompson                           The polynomial degree of the resulting Q_k element is k=P-1.
989ea61e9acSJeremy L Thompson   @param[in]  Q         Number of quadrature points in one dimension.
990ea61e9acSJeremy L Thompson   @param[in]  quad_mode Distribution of the Q quadrature points (affects order of accuracy for the quadrature)
991ea61e9acSJeremy L Thompson   @param[out] basis     Address of the variable where the newly created CeedBasis will be stored.
992b11c1e72Sjeremylt 
993b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
994dfdf5a53Sjeremylt 
9957a982d89SJeremy L. Thompson   @ref User
996b11c1e72Sjeremylt **/
9972b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P, CeedInt Q, CeedQuadMode quad_mode, CeedBasis *basis) {
998d7b241e6Sjeremylt   // Allocate
999c8c3fa7dSJeremy L Thompson   int        ierr = CEED_ERROR_SUCCESS;
10002b730f8bSJeremy L Thompson   CeedScalar c1, c2, c3, c4, dx, *nodes, *interp_1d, *grad_1d, *q_ref_1d, *q_weight_1d;
10014d537eeaSYohann 
10026574a04fSJeremy L Thompson   CeedCheck(dim > 0, ceed, CEED_ERROR_DIMENSION, "Basis dimension must be a positive value");
10036574a04fSJeremy L Thompson   CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component");
10046574a04fSJeremy L Thompson   CeedCheck(P > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node");
10056574a04fSJeremy L Thompson   CeedCheck(Q > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point");
1006227444bfSJeremy L Thompson 
1007e15f9bd0SJeremy L Thompson   // Get Nodes and Weights
10082b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(P * Q, &interp_1d));
10092b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(P * Q, &grad_1d));
10102b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(P, &nodes));
10112b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q, &q_ref_1d));
10122b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q, &q_weight_1d));
10132b730f8bSJeremy L Thompson   if (CeedLobattoQuadrature(P, nodes, NULL) != CEED_ERROR_SUCCESS) goto cleanup;
1014d1d35e2fSjeremylt   switch (quad_mode) {
1015d7b241e6Sjeremylt     case CEED_GAUSS:
1016d1d35e2fSjeremylt       ierr = CeedGaussQuadrature(Q, q_ref_1d, q_weight_1d);
1017d7b241e6Sjeremylt       break;
1018d7b241e6Sjeremylt     case CEED_GAUSS_LOBATTO:
1019d1d35e2fSjeremylt       ierr = CeedLobattoQuadrature(Q, q_ref_1d, q_weight_1d);
1020d7b241e6Sjeremylt       break;
1021d7b241e6Sjeremylt   }
10222b730f8bSJeremy L Thompson   if (ierr != CEED_ERROR_SUCCESS) goto cleanup;
1023e15f9bd0SJeremy L Thompson 
1024d7b241e6Sjeremylt   // Build B, D matrix
1025d7b241e6Sjeremylt   // Fornberg, 1998
1026c8c3fa7dSJeremy L Thompson   for (CeedInt i = 0; i < Q; i++) {
1027d7b241e6Sjeremylt     c1                   = 1.0;
1028d1d35e2fSjeremylt     c3                   = nodes[0] - q_ref_1d[i];
1029d1d35e2fSjeremylt     interp_1d[i * P + 0] = 1.0;
1030c8c3fa7dSJeremy L Thompson     for (CeedInt j = 1; j < P; j++) {
1031d7b241e6Sjeremylt       c2 = 1.0;
1032d7b241e6Sjeremylt       c4 = c3;
1033d1d35e2fSjeremylt       c3 = nodes[j] - q_ref_1d[i];
1034c8c3fa7dSJeremy L Thompson       for (CeedInt k = 0; k < j; k++) {
1035d7b241e6Sjeremylt         dx = nodes[j] - nodes[k];
1036d7b241e6Sjeremylt         c2 *= dx;
1037d7b241e6Sjeremylt         if (k == j - 1) {
1038d1d35e2fSjeremylt           grad_1d[i * P + j]   = c1 * (interp_1d[i * P + k] - c4 * grad_1d[i * P + k]) / c2;
1039d1d35e2fSjeremylt           interp_1d[i * P + j] = -c1 * c4 * interp_1d[i * P + k] / c2;
1040d7b241e6Sjeremylt         }
1041d1d35e2fSjeremylt         grad_1d[i * P + k]   = (c3 * grad_1d[i * P + k] - interp_1d[i * P + k]) / dx;
1042d1d35e2fSjeremylt         interp_1d[i * P + k] = c3 * interp_1d[i * P + k] / dx;
1043d7b241e6Sjeremylt       }
1044d7b241e6Sjeremylt       c1 = c2;
1045d7b241e6Sjeremylt     }
1046d7b241e6Sjeremylt   }
10479ac7b42eSJeremy L Thompson   // Pass to CeedBasisCreateTensorH1
10482b730f8bSJeremy L Thompson   CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P, Q, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis));
1049e15f9bd0SJeremy L Thompson cleanup:
10502b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_1d));
10512b730f8bSJeremy L Thompson   CeedCall(CeedFree(&grad_1d));
10522b730f8bSJeremy L Thompson   CeedCall(CeedFree(&nodes));
10532b730f8bSJeremy L Thompson   CeedCall(CeedFree(&q_ref_1d));
10542b730f8bSJeremy L Thompson   CeedCall(CeedFree(&q_weight_1d));
1055e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1056d7b241e6Sjeremylt }
1057d7b241e6Sjeremylt 
1058b11c1e72Sjeremylt /**
1059ba59ac12SSebastian Grimberg   @brief Create a non tensor-product basis for H^1 discretizations
1060a8de75f0Sjeremylt 
1061ea61e9acSJeremy L Thompson   @param[in]  ceed      Ceed object where the CeedBasis will be created
1062ea61e9acSJeremy L Thompson   @param[in]  topo      Topology of element, e.g. hypercube, simplex, ect
1063ea61e9acSJeremy L Thompson   @param[in]  num_comp  Number of field components (1 for scalar fields)
1064ea61e9acSJeremy L Thompson   @param[in]  num_nodes Total number of nodes
1065ea61e9acSJeremy L Thompson   @param[in]  num_qpts  Total number of quadrature points
1066ea61e9acSJeremy L Thompson   @param[in]  interp    Row-major (num_qpts * num_nodes) matrix expressing the values of nodal basis functions at quadrature points
1067c4e3f59bSSebastian Grimberg   @param[in]  grad      Row-major (dim * num_qpts * num_nodes) matrix expressing derivatives of nodal basis functions at quadrature points
10689fe083eeSJeremy L Thompson   @param[in]  q_ref     Array of length num_qpts * dim holding the locations of quadrature points on the reference element
1069ea61e9acSJeremy L Thompson   @param[in]  q_weight  Array of length num_qpts holding the quadrature weights on the reference element
1070ea61e9acSJeremy L Thompson   @param[out] basis     Address of the variable where the newly created CeedBasis will be stored.
1071a8de75f0Sjeremylt 
1072a8de75f0Sjeremylt   @return An error code: 0 - success, otherwise - failure
1073a8de75f0Sjeremylt 
10747a982d89SJeremy L. Thompson   @ref User
1075a8de75f0Sjeremylt **/
10762b730f8bSJeremy L Thompson int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
10772b730f8bSJeremy L Thompson                       const CeedScalar *grad, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) {
1078d1d35e2fSjeremylt   CeedInt P = num_nodes, Q = num_qpts, dim = 0;
1079a8de75f0Sjeremylt 
10805fe0d4faSjeremylt   if (!ceed->BasisCreateH1) {
10815fe0d4faSjeremylt     Ceed delegate;
10826574a04fSJeremy L Thompson 
10832b730f8bSJeremy L Thompson     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
10846574a04fSJeremy L Thompson     CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisCreateH1");
10852b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateH1(delegate, topo, num_comp, num_nodes, num_qpts, interp, grad, q_ref, q_weight, basis));
1086e15f9bd0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
10875fe0d4faSjeremylt   }
10885fe0d4faSjeremylt 
10896574a04fSJeremy L Thompson   CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component");
10906574a04fSJeremy L Thompson   CeedCheck(num_nodes > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node");
10916574a04fSJeremy L Thompson   CeedCheck(num_qpts > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point");
1092227444bfSJeremy L Thompson 
10932b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetTopologyDimension(topo, &dim));
1094a8de75f0Sjeremylt 
1095db002c03SJeremy L Thompson   CeedCall(CeedCalloc(1, basis));
1096db002c03SJeremy L Thompson   CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed));
1097d1d35e2fSjeremylt   (*basis)->ref_count       = 1;
10986402da51SJeremy L Thompson   (*basis)->is_tensor_basis = false;
1099a8de75f0Sjeremylt   (*basis)->dim             = dim;
1100d99fa3c5SJeremy L Thompson   (*basis)->topo            = topo;
1101d1d35e2fSjeremylt   (*basis)->num_comp        = num_comp;
1102a8de75f0Sjeremylt   (*basis)->P               = P;
1103a8de75f0Sjeremylt   (*basis)->Q               = Q;
1104c4e3f59bSSebastian Grimberg   (*basis)->fe_space        = CEED_FE_SPACE_H1;
11052b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q * dim, &(*basis)->q_ref_1d));
11062b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q, &(*basis)->q_weight_1d));
1107ff3a0f91SJeremy L Thompson   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0]));
1108ff3a0f91SJeremy L Thompson   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0]));
11092b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q * P, &(*basis)->interp));
11102b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(dim * Q * P, &(*basis)->grad));
1111ff3a0f91SJeremy L Thompson   if (interp) memcpy((*basis)->interp, interp, Q * P * sizeof(interp[0]));
1112ff3a0f91SJeremy L Thompson   if (grad) memcpy((*basis)->grad, grad, dim * Q * P * sizeof(grad[0]));
11132b730f8bSJeremy L Thompson   CeedCall(ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, q_ref, q_weight, *basis));
1114e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1115a8de75f0Sjeremylt }
1116a8de75f0Sjeremylt 
1117a8de75f0Sjeremylt /**
1118859c15bbSJames Wright   @brief Create a non tensor-product basis for \f$H(\mathrm{div})\f$ discretizations
111950c301a5SRezgar Shakeri 
1120ea61e9acSJeremy L Thompson   @param[in]  ceed      Ceed object where the CeedBasis will be created
1121ea61e9acSJeremy 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
1122ea61e9acSJeremy L Thompson   @param[in]  num_comp  Number of components (usually 1 for vectors in H(div) bases)
1123ea61e9acSJeremy L Thompson   @param[in]  num_nodes Total number of nodes (dofs per element)
1124ea61e9acSJeremy L Thompson   @param[in]  num_qpts  Total number of quadrature points
1125c4e3f59bSSebastian Grimberg   @param[in]  interp    Row-major (dim * num_qpts * num_nodes) matrix expressing the values of basis functions at quadrature points
1126c4e3f59bSSebastian Grimberg   @param[in]  div       Row-major (num_qpts * num_nodes) matrix expressing divergence of basis functions at quadrature points
11279fe083eeSJeremy L Thompson   @param[in]  q_ref     Array of length num_qpts * dim holding the locations of quadrature points on the reference element
1128ea61e9acSJeremy L Thompson   @param[in]  q_weight  Array of length num_qpts holding the quadrature weights on the reference element
1129ea61e9acSJeremy L Thompson   @param[out] basis     Address of the variable where the newly created CeedBasis will be stored.
113050c301a5SRezgar Shakeri 
113150c301a5SRezgar Shakeri   @return An error code: 0 - success, otherwise - failure
113250c301a5SRezgar Shakeri 
113350c301a5SRezgar Shakeri   @ref User
113450c301a5SRezgar Shakeri **/
11352b730f8bSJeremy L Thompson int CeedBasisCreateHdiv(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
11362b730f8bSJeremy L Thompson                         const CeedScalar *div, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) {
113750c301a5SRezgar Shakeri   CeedInt Q = num_qpts, P = num_nodes, dim = 0;
1138c4e3f59bSSebastian Grimberg 
113950c301a5SRezgar Shakeri   if (!ceed->BasisCreateHdiv) {
114050c301a5SRezgar Shakeri     Ceed delegate;
11416574a04fSJeremy L Thompson 
11422b730f8bSJeremy L Thompson     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
11436574a04fSJeremy L Thompson     CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateHdiv");
11442b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateHdiv(delegate, topo, num_comp, num_nodes, num_qpts, interp, div, q_ref, q_weight, basis));
114550c301a5SRezgar Shakeri     return CEED_ERROR_SUCCESS;
114650c301a5SRezgar Shakeri   }
114750c301a5SRezgar Shakeri 
11486574a04fSJeremy L Thompson   CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component");
11496574a04fSJeremy L Thompson   CeedCheck(num_nodes > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node");
11506574a04fSJeremy L Thompson   CeedCheck(num_qpts > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point");
1151227444bfSJeremy L Thompson 
1152c4e3f59bSSebastian Grimberg   CeedCall(CeedBasisGetTopologyDimension(topo, &dim));
1153c4e3f59bSSebastian Grimberg 
1154db002c03SJeremy L Thompson   CeedCall(CeedCalloc(1, basis));
1155db002c03SJeremy L Thompson   CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed));
115650c301a5SRezgar Shakeri   (*basis)->ref_count       = 1;
11576402da51SJeremy L Thompson   (*basis)->is_tensor_basis = false;
115850c301a5SRezgar Shakeri   (*basis)->dim             = dim;
115950c301a5SRezgar Shakeri   (*basis)->topo            = topo;
116050c301a5SRezgar Shakeri   (*basis)->num_comp        = num_comp;
116150c301a5SRezgar Shakeri   (*basis)->P               = P;
116250c301a5SRezgar Shakeri   (*basis)->Q               = Q;
1163c4e3f59bSSebastian Grimberg   (*basis)->fe_space        = CEED_FE_SPACE_HDIV;
11642b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q * dim, &(*basis)->q_ref_1d));
11652b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q, &(*basis)->q_weight_1d));
116650c301a5SRezgar Shakeri   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0]));
116750c301a5SRezgar Shakeri   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0]));
11682b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(dim * Q * P, &(*basis)->interp));
11692b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q * P, &(*basis)->div));
117050c301a5SRezgar Shakeri   if (interp) memcpy((*basis)->interp, interp, dim * Q * P * sizeof(interp[0]));
117150c301a5SRezgar Shakeri   if (div) memcpy((*basis)->div, div, Q * P * sizeof(div[0]));
11722b730f8bSJeremy L Thompson   CeedCall(ceed->BasisCreateHdiv(topo, dim, P, Q, interp, div, q_ref, q_weight, *basis));
117350c301a5SRezgar Shakeri   return CEED_ERROR_SUCCESS;
117450c301a5SRezgar Shakeri }
117550c301a5SRezgar Shakeri 
117650c301a5SRezgar Shakeri /**
11774385fb7fSSebastian Grimberg   @brief Create a non tensor-product basis for \f$H(\mathrm{curl})\f$ discretizations
1178c4e3f59bSSebastian Grimberg 
1179c4e3f59bSSebastian Grimberg   @param[in]  ceed      Ceed object where the CeedBasis will be created
1180c4e3f59bSSebastian Grimberg   @param[in]  topo      Topology of element (`CEED_TOPOLOGY_QUAD`, `CEED_TOPOLOGY_PRISM`, etc.), dimension of which is used in some array sizes below
1181c4e3f59bSSebastian Grimberg   @param[in]  num_comp  Number of components (usually 1 for vectors in H(curl) bases)
1182c4e3f59bSSebastian Grimberg   @param[in]  num_nodes Total number of nodes (dofs per element)
1183c4e3f59bSSebastian Grimberg   @param[in]  num_qpts  Total number of quadrature points
1184c4e3f59bSSebastian Grimberg   @param[in]  interp    Row-major (dim * num_qpts * num_nodes) matrix expressing the values of basis functions at quadrature points
1185c4e3f59bSSebastian Grimberg   @param[in]  curl      Row-major (curl_comp * num_qpts * num_nodes, curl_comp = 1 if dim < 3 else dim) matrix expressing curl of basis functions at
1186c4e3f59bSSebastian Grimberg quadrature points
1187c4e3f59bSSebastian Grimberg   @param[in]  q_ref     Array of length num_qpts * dim holding the locations of quadrature points on the reference element
1188c4e3f59bSSebastian Grimberg   @param[in]  q_weight  Array of length num_qpts holding the quadrature weights on the reference element
1189c4e3f59bSSebastian Grimberg   @param[out] basis     Address of the variable where the newly created CeedBasis will be stored.
1190c4e3f59bSSebastian Grimberg 
1191c4e3f59bSSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
1192c4e3f59bSSebastian Grimberg 
1193c4e3f59bSSebastian Grimberg   @ref User
1194c4e3f59bSSebastian Grimberg **/
1195c4e3f59bSSebastian Grimberg int CeedBasisCreateHcurl(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
1196c4e3f59bSSebastian Grimberg                          const CeedScalar *curl, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) {
1197c4e3f59bSSebastian Grimberg   CeedInt Q = num_qpts, P = num_nodes, dim = 0, curl_comp = 0;
1198c4e3f59bSSebastian Grimberg 
1199c4e3f59bSSebastian Grimberg   if (!ceed->BasisCreateHdiv) {
1200c4e3f59bSSebastian Grimberg     Ceed delegate;
12016574a04fSJeremy L Thompson 
1202c4e3f59bSSebastian Grimberg     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
12036574a04fSJeremy L Thompson     CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateHcurl");
1204c4e3f59bSSebastian Grimberg     CeedCall(CeedBasisCreateHcurl(delegate, topo, num_comp, num_nodes, num_qpts, interp, curl, q_ref, q_weight, basis));
1205c4e3f59bSSebastian Grimberg     return CEED_ERROR_SUCCESS;
1206c4e3f59bSSebastian Grimberg   }
1207c4e3f59bSSebastian Grimberg 
12086574a04fSJeremy L Thompson   CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component");
12096574a04fSJeremy L Thompson   CeedCheck(num_nodes > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node");
12106574a04fSJeremy L Thompson   CeedCheck(num_qpts > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point");
1211c4e3f59bSSebastian Grimberg 
1212c4e3f59bSSebastian Grimberg   CeedCall(CeedBasisGetTopologyDimension(topo, &dim));
1213c4e3f59bSSebastian Grimberg   curl_comp = (dim < 3) ? 1 : dim;
1214c4e3f59bSSebastian Grimberg 
1215db002c03SJeremy L Thompson   CeedCall(CeedCalloc(1, basis));
1216db002c03SJeremy L Thompson   CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed));
1217c4e3f59bSSebastian Grimberg   (*basis)->ref_count       = 1;
12186402da51SJeremy L Thompson   (*basis)->is_tensor_basis = false;
1219c4e3f59bSSebastian Grimberg   (*basis)->dim             = dim;
1220c4e3f59bSSebastian Grimberg   (*basis)->topo            = topo;
1221c4e3f59bSSebastian Grimberg   (*basis)->num_comp        = num_comp;
1222c4e3f59bSSebastian Grimberg   (*basis)->P               = P;
1223c4e3f59bSSebastian Grimberg   (*basis)->Q               = Q;
1224c4e3f59bSSebastian Grimberg   (*basis)->fe_space        = CEED_FE_SPACE_HCURL;
1225c4e3f59bSSebastian Grimberg   CeedCall(CeedMalloc(Q * dim, &(*basis)->q_ref_1d));
1226c4e3f59bSSebastian Grimberg   CeedCall(CeedMalloc(Q, &(*basis)->q_weight_1d));
1227c4e3f59bSSebastian Grimberg   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0]));
1228c4e3f59bSSebastian Grimberg   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0]));
1229c4e3f59bSSebastian Grimberg   CeedCall(CeedMalloc(dim * Q * P, &(*basis)->interp));
1230c4e3f59bSSebastian Grimberg   CeedCall(CeedMalloc(curl_comp * Q * P, &(*basis)->curl));
1231c4e3f59bSSebastian Grimberg   if (interp) memcpy((*basis)->interp, interp, dim * Q * P * sizeof(interp[0]));
1232c4e3f59bSSebastian Grimberg   if (curl) memcpy((*basis)->curl, curl, curl_comp * Q * P * sizeof(curl[0]));
1233c4e3f59bSSebastian Grimberg   CeedCall(ceed->BasisCreateHcurl(topo, dim, P, Q, interp, curl, q_ref, q_weight, *basis));
1234c4e3f59bSSebastian Grimberg   return CEED_ERROR_SUCCESS;
1235c4e3f59bSSebastian Grimberg }
1236c4e3f59bSSebastian Grimberg 
1237c4e3f59bSSebastian Grimberg /**
1238ea61e9acSJeremy L Thompson   @brief Create a CeedBasis for projection from the nodes of `basis_from` to the nodes of `basis_to`.
1239ba59ac12SSebastian Grimberg 
12409fd66db6SSebastian Grimberg   Only `CEED_EVAL_INTERP` will be valid for the new basis, `basis_project`.
12419fd66db6SSebastian Grimberg   For H^1 spaces, `CEED_EVAL_GRAD` will also be valid.
1242de05fbb2SSebastian Grimberg   The interpolation is given by `interp_project = interp_to^+ * interp_from`, where the pseudoinverse `interp_to^+` is given by QR
12439fd66db6SSebastian Grimberg factorization.
12449fd66db6SSebastian Grimberg   The gradient (for the H^1 case) is given by `grad_project = interp_to^+ * grad_from`.
124515ad3917SSebastian Grimberg 
124615ad3917SSebastian Grimberg   Note: `basis_from` and `basis_to` must have compatible quadrature spaces.
124715ad3917SSebastian Grimberg 
12489fd66db6SSebastian Grimberg   Note: `basis_project` will have the same number of components as `basis_from`, regardless of the number of components that `basis_to` has.
12499fd66db6SSebastian Grimberg         If `basis_from` has 3 components and `basis_to` has 5 components, then `basis_project` will have 3 components.
1250f113e5dcSJeremy L Thompson 
1251f113e5dcSJeremy L Thompson   @param[in]  basis_from    CeedBasis to prolong from
1252446e7af4SJeremy L Thompson   @param[in]  basis_to      CeedBasis to prolong to
1253ea61e9acSJeremy L Thompson   @param[out] basis_project Address of the variable where the newly created CeedBasis will be stored.
1254f113e5dcSJeremy L Thompson 
1255f113e5dcSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1256f113e5dcSJeremy L Thompson 
1257f113e5dcSJeremy L Thompson   @ref User
1258f113e5dcSJeremy L Thompson **/
12592b730f8bSJeremy L Thompson int CeedBasisCreateProjection(CeedBasis basis_from, CeedBasis basis_to, CeedBasis *basis_project) {
1260f113e5dcSJeremy L Thompson   Ceed ceed;
12612b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetCeed(basis_to, &ceed));
1262f113e5dcSJeremy L Thompson 
1263ecc88aebSJeremy L Thompson   // Create projection matrix
126414556e63SJeremy L Thompson   CeedScalar *interp_project, *grad_project;
12652b730f8bSJeremy L Thompson   CeedCall(CeedBasisCreateProjectionMatrices(basis_from, basis_to, &interp_project, &grad_project));
1266f113e5dcSJeremy L Thompson 
1267f113e5dcSJeremy L Thompson   // Build basis
1268f113e5dcSJeremy L Thompson   bool        is_tensor;
1269f113e5dcSJeremy L Thompson   CeedInt     dim, num_comp;
127014556e63SJeremy L Thompson   CeedScalar *q_ref, *q_weight;
12712b730f8bSJeremy L Thompson   CeedCall(CeedBasisIsTensor(basis_to, &is_tensor));
12722b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetDimension(basis_to, &dim));
12732b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumComponents(basis_from, &num_comp));
1274f113e5dcSJeremy L Thompson   if (is_tensor) {
1275f113e5dcSJeremy L Thompson     CeedInt P_1d_to, P_1d_from;
12762b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_1d_from));
12772b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_1d_to));
12782b730f8bSJeremy L Thompson     CeedCall(CeedCalloc(P_1d_to, &q_ref));
12792b730f8bSJeremy L Thompson     CeedCall(CeedCalloc(P_1d_to, &q_weight));
12802b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P_1d_from, P_1d_to, interp_project, grad_project, q_ref, q_weight, basis_project));
1281f113e5dcSJeremy L Thompson   } else {
1282de05fbb2SSebastian Grimberg     // Even if basis_to and basis_from are not H1, the resulting basis is H1 for interpolation to work
1283f113e5dcSJeremy L Thompson     CeedElemTopology topo;
12842b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetTopology(basis_to, &topo));
1285f113e5dcSJeremy L Thompson     CeedInt num_nodes_to, num_nodes_from;
12862b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis_from, &num_nodes_from));
12872b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis_to, &num_nodes_to));
12882b730f8bSJeremy L Thompson     CeedCall(CeedCalloc(num_nodes_to * dim, &q_ref));
12892b730f8bSJeremy L Thompson     CeedCall(CeedCalloc(num_nodes_to, &q_weight));
12902b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateH1(ceed, topo, num_comp, num_nodes_from, num_nodes_to, interp_project, grad_project, q_ref, q_weight, basis_project));
1291f113e5dcSJeremy L Thompson   }
1292f113e5dcSJeremy L Thompson 
1293f113e5dcSJeremy L Thompson   // Cleanup
12942b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_project));
12952b730f8bSJeremy L Thompson   CeedCall(CeedFree(&grad_project));
12962b730f8bSJeremy L Thompson   CeedCall(CeedFree(&q_ref));
12972b730f8bSJeremy L Thompson   CeedCall(CeedFree(&q_weight));
1298f113e5dcSJeremy L Thompson 
1299f113e5dcSJeremy L Thompson   return CEED_ERROR_SUCCESS;
1300f113e5dcSJeremy L Thompson }
1301f113e5dcSJeremy L Thompson 
1302f113e5dcSJeremy L Thompson /**
1303ea61e9acSJeremy L Thompson   @brief Copy the pointer to a CeedBasis.
13049560d06aSjeremylt 
1305512bb800SJeremy 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.
1306512bb800SJeremy L Thompson         This CeedBasis will be destroyed if `basis_copy` is the only reference to this CeedBasis.
1307ea61e9acSJeremy L Thompson 
1308ea61e9acSJeremy L Thompson   @param[in]     basis      CeedBasis to copy reference to
1309ea61e9acSJeremy L Thompson   @param[in,out] basis_copy Variable to store copied reference
13109560d06aSjeremylt 
13119560d06aSjeremylt   @return An error code: 0 - success, otherwise - failure
13129560d06aSjeremylt 
13139560d06aSjeremylt   @ref User
13149560d06aSjeremylt **/
13159560d06aSjeremylt int CeedBasisReferenceCopy(CeedBasis basis, CeedBasis *basis_copy) {
1316393ac2cdSJeremy L Thompson   if (basis != CEED_BASIS_COLLOCATED) CeedCall(CeedBasisReference(basis));
13172b730f8bSJeremy L Thompson   CeedCall(CeedBasisDestroy(basis_copy));
13189560d06aSjeremylt   *basis_copy = basis;
13199560d06aSjeremylt   return CEED_ERROR_SUCCESS;
13209560d06aSjeremylt }
13219560d06aSjeremylt 
13229560d06aSjeremylt /**
13237a982d89SJeremy L. Thompson   @brief View a CeedBasis
13247a982d89SJeremy L. Thompson 
1325ea61e9acSJeremy L Thompson   @param[in] basis  CeedBasis to view
1326ea61e9acSJeremy L Thompson   @param[in] stream Stream to view to, e.g., stdout
13277a982d89SJeremy L. Thompson 
13287a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
13297a982d89SJeremy L. Thompson 
13307a982d89SJeremy L. Thompson   @ref User
13317a982d89SJeremy L. Thompson **/
13327a982d89SJeremy L. Thompson int CeedBasisView(CeedBasis basis, FILE *stream) {
133350c301a5SRezgar Shakeri   CeedElemTopology topo     = basis->topo;
1334c4e3f59bSSebastian Grimberg   CeedFESpace      fe_space = basis->fe_space;
1335c4e3f59bSSebastian Grimberg   CeedInt          q_comp   = 0;
13362b730f8bSJeremy L Thompson 
133750c301a5SRezgar Shakeri   // Print FE space and element topology of the basis
13386402da51SJeremy L Thompson   if (basis->is_tensor_basis) {
1339c4e3f59bSSebastian Grimberg     fprintf(stream, "CeedBasis (%s on a %s element): dim=%" CeedInt_FMT " P=%" CeedInt_FMT " Q=%" CeedInt_FMT "\n", CeedFESpaces[fe_space],
13402b730f8bSJeremy L Thompson             CeedElemTopologies[topo], basis->dim, basis->P_1d, basis->Q_1d);
134150c301a5SRezgar Shakeri   } else {
1342c4e3f59bSSebastian Grimberg     fprintf(stream, "CeedBasis (%s on a %s element): dim=%" CeedInt_FMT " P=%" CeedInt_FMT " Q=%" CeedInt_FMT "\n", CeedFESpaces[fe_space],
13432b730f8bSJeremy L Thompson             CeedElemTopologies[topo], basis->dim, basis->P, basis->Q);
134450c301a5SRezgar Shakeri   }
1345ea61e9acSJeremy L Thompson   // Print quadrature data, interpolation/gradient/divergence/curl of the basis
13466402da51SJeremy L Thompson   if (basis->is_tensor_basis) {  // tensor basis
13472b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_ref_1d, stream));
13482b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_weight_1d, stream));
13492b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("interp1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, basis->interp_1d, stream));
13502b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("grad1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, basis->grad_1d, stream));
135150c301a5SRezgar Shakeri   } else {  // non-tensor basis
13522b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("qref", "\t% 12.8f", 1, basis->Q * basis->dim, basis->q_ref_1d, stream));
13532b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->q_weight_1d, stream));
1354c4e3f59bSSebastian Grimberg     CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp));
1355c4e3f59bSSebastian Grimberg     CeedCall(CeedScalarView("interp", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->interp, stream));
135650c301a5SRezgar Shakeri     if (basis->grad) {
1357c4e3f59bSSebastian Grimberg       CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp));
1358c4e3f59bSSebastian Grimberg       CeedCall(CeedScalarView("grad", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->grad, stream));
13597a982d89SJeremy L. Thompson     }
136050c301a5SRezgar Shakeri     if (basis->div) {
1361c4e3f59bSSebastian Grimberg       CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_DIV, &q_comp));
1362c4e3f59bSSebastian Grimberg       CeedCall(CeedScalarView("div", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->div, stream));
1363c4e3f59bSSebastian Grimberg     }
1364c4e3f59bSSebastian Grimberg     if (basis->curl) {
1365c4e3f59bSSebastian Grimberg       CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_CURL, &q_comp));
1366c4e3f59bSSebastian Grimberg       CeedCall(CeedScalarView("curl", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->curl, stream));
136750c301a5SRezgar Shakeri     }
136850c301a5SRezgar Shakeri   }
1369e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
13707a982d89SJeremy L. Thompson }
13717a982d89SJeremy L. Thompson 
13727a982d89SJeremy L. Thompson /**
13737a982d89SJeremy L. Thompson   @brief Apply basis evaluation from nodes to quadrature points or vice versa
13747a982d89SJeremy L. Thompson 
1375ea61e9acSJeremy L Thompson   @param[in]  basis      CeedBasis to evaluate
1376ea61e9acSJeremy L Thompson   @param[in]  num_elem   The number of elements to apply the basis evaluation to;
1377ea61e9acSJeremy L Thompson                            the backend will specify the ordering in CeedElemRestrictionCreateBlocked()
1378ea61e9acSJeremy L Thompson   @param[in]  t_mode    \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature points;
1379ea61e9acSJeremy L Thompson                           \ref CEED_TRANSPOSE to apply the transpose, mapping from quadrature points to nodes
1380ea61e9acSJeremy L Thompson   @param[in]  eval_mode \ref CEED_EVAL_NONE to use values directly,
13817a982d89SJeremy L. Thompson                           \ref CEED_EVAL_INTERP to use interpolated values,
13827a982d89SJeremy L. Thompson                           \ref CEED_EVAL_GRAD to use gradients,
1383c4e3f59bSSebastian Grimberg                           \ref CEED_EVAL_DIV to use divergence,
1384c4e3f59bSSebastian Grimberg                           \ref CEED_EVAL_CURL to use curl,
13857a982d89SJeremy L. Thompson                           \ref CEED_EVAL_WEIGHT to use quadrature weights.
13867a982d89SJeremy L. Thompson   @param[in]  u        Input CeedVector
13877a982d89SJeremy L. Thompson   @param[out] v        Output CeedVector
13887a982d89SJeremy L. Thompson 
13897a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
13907a982d89SJeremy L. Thompson 
13917a982d89SJeremy L. Thompson   @ref User
13927a982d89SJeremy L. Thompson **/
13932b730f8bSJeremy L Thompson int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
13941f9221feSJeremy L Thompson   CeedSize u_length = 0, v_length;
1395c4e3f59bSSebastian Grimberg   CeedInt  dim, num_comp, q_comp, num_nodes, num_qpts;
13962b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetDimension(basis, &dim));
13972b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
1398c4e3f59bSSebastian Grimberg   CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp));
13992b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumNodes(basis, &num_nodes));
14002b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts));
14012b730f8bSJeremy L Thompson   CeedCall(CeedVectorGetLength(v, &v_length));
1402c8c3fa7dSJeremy L Thompson   if (u) CeedCall(CeedVectorGetLength(u, &u_length));
14037a982d89SJeremy L. Thompson 
14046574a04fSJeremy L Thompson   CeedCheck(basis->Apply, basis->ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisApply");
1405e15f9bd0SJeremy L Thompson 
1406e15f9bd0SJeremy L Thompson   // Check compatibility of topological and geometrical dimensions
14076574a04fSJeremy L Thompson   CeedCheck((t_mode == CEED_TRANSPOSE && v_length % num_nodes == 0 && u_length % num_qpts == 0) ||
14086574a04fSJeremy L Thompson                 (t_mode == CEED_NOTRANSPOSE && u_length % num_nodes == 0 && v_length % num_qpts == 0),
14096574a04fSJeremy L Thompson             basis->ceed, CEED_ERROR_DIMENSION, "Length of input/output vectors incompatible with basis dimensions");
14107a982d89SJeremy L. Thompson 
1411e15f9bd0SJeremy L Thompson   // Check vector lengths to prevent out of bounds issues
14126574a04fSJeremy L Thompson   bool good_dims = true;
1413d1d35e2fSjeremylt   switch (eval_mode) {
1414e15f9bd0SJeremy L Thompson     case CEED_EVAL_NONE:
14152b730f8bSJeremy L Thompson     case CEED_EVAL_INTERP:
14162b730f8bSJeremy L Thompson     case CEED_EVAL_GRAD:
1417c4e3f59bSSebastian Grimberg     case CEED_EVAL_DIV:
1418c4e3f59bSSebastian Grimberg     case CEED_EVAL_CURL:
14196574a04fSJeremy L Thompson       good_dims =
14206574a04fSJeremy L Thompson           ((t_mode == CEED_TRANSPOSE && u_length >= num_elem * num_comp * num_qpts * q_comp && v_length >= num_elem * num_comp * num_nodes) ||
14216574a04fSJeremy L Thompson            (t_mode == CEED_NOTRANSPOSE && v_length >= num_elem * num_qpts * num_comp * q_comp && u_length >= num_elem * num_comp * num_nodes));
1422e15f9bd0SJeremy L Thompson       break;
1423e15f9bd0SJeremy L Thompson     case CEED_EVAL_WEIGHT:
14246574a04fSJeremy L Thompson       good_dims = v_length >= num_elem * num_qpts;
1425e15f9bd0SJeremy L Thompson       break;
1426e15f9bd0SJeremy L Thompson   }
14276574a04fSJeremy L Thompson   CeedCheck(good_dims, basis->ceed, CEED_ERROR_DIMENSION, "Input/output vectors too short for basis and evaluation mode");
1428e15f9bd0SJeremy L Thompson 
14292b730f8bSJeremy L Thompson   CeedCall(basis->Apply(basis, num_elem, t_mode, eval_mode, u, v));
1430e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
14317a982d89SJeremy L. Thompson }
14327a982d89SJeremy L. Thompson 
14337a982d89SJeremy L. Thompson /**
1434c8c3fa7dSJeremy L Thompson   @brief Apply basis evaluation from nodes to arbitrary points
1435c8c3fa7dSJeremy L Thompson 
1436c8c3fa7dSJeremy L Thompson   @param[in]  basis      CeedBasis to evaluate
1437c8c3fa7dSJeremy L Thompson   @param[in]  num_points The number of points to apply the basis evaluation to
1438c8c3fa7dSJeremy L Thompson   @param[in]  t_mode    \ref CEED_NOTRANSPOSE to evaluate from nodes to points;
1439c8c3fa7dSJeremy L Thompson                           \ref CEED_TRANSPOSE to apply the transpose, mapping from points to nodes
1440c8c3fa7dSJeremy L Thompson   @param[in]  eval_mode \ref CEED_EVAL_INTERP to use interpolated values,
1441c8c3fa7dSJeremy L Thompson                           \ref CEED_EVAL_GRAD to use gradients
1442c8c3fa7dSJeremy L Thompson   @param[in]  x_ref    CeedVector holding reference coordinates of each point
1443c8c3fa7dSJeremy L Thompson   @param[in]  u        Input CeedVector, of length `num_nodes * num_comp` for `CEED_NOTRANSPOSE`
1444c8c3fa7dSJeremy L Thompson   @param[out] v        Output CeedVector, of length `num_points * num_q_comp` for `CEED_NOTRANSPOSE` with `CEED_EVAL_INTERP`
1445c8c3fa7dSJeremy L Thompson 
1446c8c3fa7dSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1447c8c3fa7dSJeremy L Thompson 
1448c8c3fa7dSJeremy L Thompson   @ref User
1449c8c3fa7dSJeremy L Thompson **/
1450c8c3fa7dSJeremy L Thompson int CeedBasisApplyAtPoints(CeedBasis basis, CeedInt num_points, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u,
1451c8c3fa7dSJeremy L Thompson                            CeedVector v) {
1452c8c3fa7dSJeremy L Thompson   CeedSize x_length = 0, u_length = 0, v_length;
1453c8c3fa7dSJeremy L Thompson   CeedInt  dim, num_comp, num_q_comp, num_nodes, P_1d = 1, Q_1d = 1;
1454c8c3fa7dSJeremy L Thompson 
1455c8c3fa7dSJeremy L Thompson   CeedCall(CeedBasisGetDimension(basis, &dim));
1456c8c3fa7dSJeremy L Thompson   CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d));
1457c8c3fa7dSJeremy L Thompson   CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
1458c8c3fa7dSJeremy L Thompson   CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
1459c8c3fa7dSJeremy L Thompson   CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &num_q_comp));
1460c8c3fa7dSJeremy L Thompson   CeedCall(CeedBasisGetNumNodes(basis, &num_nodes));
1461c8c3fa7dSJeremy L Thompson   CeedCall(CeedVectorGetLength(x_ref, &x_length));
1462c8c3fa7dSJeremy L Thompson   CeedCall(CeedVectorGetLength(v, &v_length));
1463c8c3fa7dSJeremy L Thompson   CeedCall(CeedVectorGetLength(u, &u_length));
1464c8c3fa7dSJeremy L Thompson 
1465c8c3fa7dSJeremy L Thompson   // Check compatibility of topological and geometrical dimensions
1466c8c3fa7dSJeremy L Thompson   CeedCheck((t_mode == CEED_TRANSPOSE && v_length % num_nodes == 0) || (t_mode == CEED_NOTRANSPOSE && u_length % num_nodes == 0), basis->ceed,
1467c8c3fa7dSJeremy L Thompson             CEED_ERROR_DIMENSION, "Length of input/output vectors incompatible with basis dimensions and number of points");
1468c8c3fa7dSJeremy L Thompson 
1469c8c3fa7dSJeremy L Thompson   // Check compatibility coordinates vector
1470c8c3fa7dSJeremy L Thompson   CeedCheck(x_length >= num_points * dim, basis->ceed, CEED_ERROR_DIMENSION,
1471c8c3fa7dSJeremy L Thompson             "Length of reference coordinate vector incompatible with basis dimension and number of points");
1472c8c3fa7dSJeremy L Thompson 
1473c8c3fa7dSJeremy L Thompson   // Check vector lengths to prevent out of bounds issues
1474c8c3fa7dSJeremy L Thompson   bool good_dims = false;
1475c8c3fa7dSJeremy L Thompson   switch (eval_mode) {
1476c8c3fa7dSJeremy L Thompson     case CEED_EVAL_INTERP:
1477c8c3fa7dSJeremy L Thompson       good_dims = ((t_mode == CEED_TRANSPOSE && (u_length >= num_points * num_q_comp || v_length >= num_nodes * num_comp)) ||
1478c8c3fa7dSJeremy L Thompson                    (t_mode == CEED_NOTRANSPOSE && (v_length >= num_points * num_q_comp || u_length >= num_nodes * num_comp)));
1479c8c3fa7dSJeremy L Thompson       break;
1480c8c3fa7dSJeremy L Thompson     case CEED_EVAL_GRAD:
1481edfbf3a6SJeremy L Thompson       good_dims = ((t_mode == CEED_TRANSPOSE && (u_length >= num_points * num_q_comp * dim || v_length >= num_nodes * num_comp)) ||
1482edfbf3a6SJeremy L Thompson                    (t_mode == CEED_NOTRANSPOSE && (v_length >= num_points * num_q_comp * dim || u_length >= num_nodes * num_comp)));
1483edfbf3a6SJeremy L Thompson       break;
1484c8c3fa7dSJeremy L Thompson     case CEED_EVAL_NONE:
1485c8c3fa7dSJeremy L Thompson     case CEED_EVAL_WEIGHT:
1486c8c3fa7dSJeremy L Thompson     case CEED_EVAL_DIV:
1487c8c3fa7dSJeremy L Thompson     case CEED_EVAL_CURL:
1488c8c3fa7dSJeremy L Thompson       // LCOV_EXCL_START
1489c8c3fa7dSJeremy L Thompson       return CeedError(basis->ceed, CEED_ERROR_UNSUPPORTED, "Evaluation at arbitrary points not supported for %s", CeedEvalModes[eval_mode]);
1490c8c3fa7dSJeremy L Thompson       // LCOV_EXCL_STOP
1491c8c3fa7dSJeremy L Thompson   }
1492c8c3fa7dSJeremy L Thompson   CeedCheck(good_dims, basis->ceed, CEED_ERROR_DIMENSION, "Input/output vectors too short for basis and evaluation mode");
1493c8c3fa7dSJeremy L Thompson 
1494c8c3fa7dSJeremy L Thompson   // Backend method
1495c8c3fa7dSJeremy L Thompson   if (basis->ApplyAtPoints) {
1496c8c3fa7dSJeremy L Thompson     CeedCall(basis->ApplyAtPoints(basis, num_points, t_mode, eval_mode, x_ref, u, v));
1497c8c3fa7dSJeremy L Thompson     return CEED_ERROR_SUCCESS;
1498c8c3fa7dSJeremy L Thompson   }
1499c8c3fa7dSJeremy L Thompson 
1500c8c3fa7dSJeremy L Thompson   // Default implementation
1501c8c3fa7dSJeremy L Thompson   CeedCheck(basis->is_tensor_basis, basis->ceed, CEED_ERROR_UNSUPPORTED, "Evaluation at arbitrary points only supported for tensor product bases");
1502edfbf3a6SJeremy L Thompson   CeedCheck(eval_mode == CEED_EVAL_INTERP || t_mode == CEED_NOTRANSPOSE, basis->ceed, CEED_ERROR_UNSUPPORTED, "%s evaluation only supported for %s",
1503edfbf3a6SJeremy L Thompson             CeedEvalModes[eval_mode], CeedTransposeModes[CEED_NOTRANSPOSE]);
1504c8c3fa7dSJeremy L Thompson   if (!basis->basis_chebyshev) {
1505c8c3fa7dSJeremy L Thompson     // Build matrix mapping from quadrature point values to Chebyshev coefficients
1506c8c3fa7dSJeremy L Thompson     CeedScalar       *tau, *C, *I, *chebyshev_coeffs_1d;
1507c8c3fa7dSJeremy L Thompson     const CeedScalar *q_ref_1d;
1508c8c3fa7dSJeremy L Thompson 
1509c8c3fa7dSJeremy L Thompson     // Build coefficient matrix
1510c8c3fa7dSJeremy L Thompson     // -- Note: Clang-tidy needs this check because it does not understand the is_tensor_basis check above
1511c8c3fa7dSJeremy L Thompson     CeedCheck(P_1d > 0 && Q_1d > 0, basis->ceed, CEED_ERROR_INCOMPATIBLE, "Basis dimensions are malformed");
1512c8c3fa7dSJeremy L Thompson     CeedCall(CeedCalloc(Q_1d * Q_1d, &C));
1513c8c3fa7dSJeremy L Thompson     CeedCall(CeedBasisGetQRef(basis, &q_ref_1d));
1514*3778dbaaSJeremy L Thompson     for (CeedInt i = 0; i < Q_1d; i++) CeedCall(CeedChebyshevPolynomialsAtPoint(q_ref_1d[i], Q_1d, &C[i * Q_1d]));
1515c8c3fa7dSJeremy L Thompson 
1516c8c3fa7dSJeremy L Thompson     // Inverse of coefficient matrix
1517c8c3fa7dSJeremy L Thompson     CeedCall(CeedCalloc(Q_1d * Q_1d, &chebyshev_coeffs_1d));
1518c8c3fa7dSJeremy L Thompson     CeedCall(CeedCalloc(Q_1d * Q_1d, &I));
1519c8c3fa7dSJeremy L Thompson     CeedCall(CeedCalloc(Q_1d, &tau));
1520c8c3fa7dSJeremy L Thompson     // -- QR Factorization, C = Q R
1521c8c3fa7dSJeremy L Thompson     CeedCall(CeedQRFactorization(basis->ceed, C, tau, Q_1d, Q_1d));
1522c8c3fa7dSJeremy L Thompson     // -- chebyshev_coeffs_1d = R_inv Q^T
1523c8c3fa7dSJeremy L Thompson     for (CeedInt i = 0; i < Q_1d; i++) I[i * Q_1d + i] = 1.0;
1524c8c3fa7dSJeremy L Thompson     // ---- Apply R_inv, chebyshev_coeffs_1d = I R_inv
1525c8c3fa7dSJeremy L Thompson     for (CeedInt i = 0; i < Q_1d; i++) {  // Row i
1526c8c3fa7dSJeremy L Thompson       chebyshev_coeffs_1d[Q_1d * i] = I[Q_1d * i] / C[0];
1527c8c3fa7dSJeremy L Thompson       for (CeedInt j = 1; j < Q_1d; j++) {  // Column j
1528c8c3fa7dSJeremy L Thompson         chebyshev_coeffs_1d[j + Q_1d * i] = I[j + Q_1d * i];
1529c8c3fa7dSJeremy L Thompson         for (CeedInt k = 0; k < j; k++) chebyshev_coeffs_1d[j + Q_1d * i] -= C[j + Q_1d * k] * chebyshev_coeffs_1d[k + Q_1d * i];
1530c8c3fa7dSJeremy L Thompson         chebyshev_coeffs_1d[j + Q_1d * i] /= C[j + Q_1d * j];
1531c8c3fa7dSJeremy L Thompson       }
1532c8c3fa7dSJeremy L Thompson     }
1533c8c3fa7dSJeremy L Thompson     // ---- Apply Q^T, chebyshev_coeffs_1d = R_inv Q^T
1534c8c3fa7dSJeremy L Thompson     CeedCall(CeedHouseholderApplyQ(chebyshev_coeffs_1d, C, tau, CEED_NOTRANSPOSE, Q_1d, Q_1d, Q_1d, 1, Q_1d));
1535c8c3fa7dSJeremy L Thompson 
1536c8c3fa7dSJeremy L Thompson     // Build basis mapping from nodes to Chebyshev coefficients
1537c8c3fa7dSJeremy L Thompson     CeedScalar       *chebyshev_interp_1d, *chebyshev_grad_1d, *chebyshev_q_weight_1d;
1538c8c3fa7dSJeremy L Thompson     const CeedScalar *interp_1d;
1539c8c3fa7dSJeremy L Thompson 
1540c8c3fa7dSJeremy L Thompson     CeedCall(CeedCalloc(Q_1d * Q_1d, &chebyshev_interp_1d));
1541c8c3fa7dSJeremy L Thompson     CeedCall(CeedCalloc(Q_1d * Q_1d, &chebyshev_grad_1d));
1542c8c3fa7dSJeremy L Thompson     CeedCall(CeedCalloc(Q_1d, &chebyshev_q_weight_1d));
1543c8c3fa7dSJeremy L Thompson     CeedCall(CeedBasisGetInterp1D(basis, &interp_1d));
1544c8c3fa7dSJeremy L Thompson     CeedCall(CeedMatrixMatrixMultiply(basis->ceed, chebyshev_coeffs_1d, interp_1d, chebyshev_interp_1d, Q_1d, P_1d, Q_1d));
1545c8c3fa7dSJeremy L Thompson 
1546c8c3fa7dSJeremy L Thompson     CeedCall(CeedVectorCreate(basis->ceed, num_comp * CeedIntPow(Q_1d, dim), &basis->vec_chebyshev));
1547c8c3fa7dSJeremy L Thompson     CeedCall(CeedBasisCreateTensorH1(basis->ceed, dim, num_comp, Q_1d, Q_1d, chebyshev_interp_1d, chebyshev_grad_1d, q_ref_1d, chebyshev_q_weight_1d,
1548c8c3fa7dSJeremy L Thompson                                      &basis->basis_chebyshev));
1549c8c3fa7dSJeremy L Thompson 
1550c8c3fa7dSJeremy L Thompson     // Cleanup
1551c8c3fa7dSJeremy L Thompson     CeedCall(CeedFree(&C));
1552c8c3fa7dSJeremy L Thompson     CeedCall(CeedFree(&chebyshev_coeffs_1d));
1553c8c3fa7dSJeremy L Thompson     CeedCall(CeedFree(&I));
1554c8c3fa7dSJeremy L Thompson     CeedCall(CeedFree(&tau));
1555c8c3fa7dSJeremy L Thompson     CeedCall(CeedFree(&chebyshev_interp_1d));
1556c8c3fa7dSJeremy L Thompson     CeedCall(CeedFree(&chebyshev_grad_1d));
1557c8c3fa7dSJeremy L Thompson     CeedCall(CeedFree(&chebyshev_q_weight_1d));
1558c8c3fa7dSJeremy L Thompson   }
1559c8c3fa7dSJeremy L Thompson 
1560c8c3fa7dSJeremy L Thompson   // Create TensorContract object if needed, such as a basis from the GPU backends
1561c8c3fa7dSJeremy L Thompson   if (!basis->contract) {
1562c8c3fa7dSJeremy L Thompson     Ceed      ceed_ref;
1563c8c3fa7dSJeremy L Thompson     CeedBasis basis_ref;
1564c8c3fa7dSJeremy L Thompson 
1565c8c3fa7dSJeremy L Thompson     CeedCall(CeedInit("/cpu/self", &ceed_ref));
1566c8c3fa7dSJeremy L Thompson     // Only need matching tensor contraction dimensions, any type of basis will work
1567c8c3fa7dSJeremy L Thompson     CeedCall(CeedBasisCreateTensorH1Lagrange(ceed_ref, dim, num_comp, Q_1d, Q_1d, CEED_GAUSS, &basis_ref));
1568c8c3fa7dSJeremy L Thompson     CeedCall(CeedTensorContractReference(basis_ref->contract));
1569c8c3fa7dSJeremy L Thompson     basis->contract = basis_ref->contract;
1570c8c3fa7dSJeremy L Thompson     CeedCall(CeedBasisDestroy(&basis_ref));
1571c8c3fa7dSJeremy L Thompson     CeedCall(CeedDestroy(&ceed_ref));
1572c8c3fa7dSJeremy L Thompson   }
1573c8c3fa7dSJeremy L Thompson 
1574c8c3fa7dSJeremy L Thompson   // Basis evaluation
1575c8c3fa7dSJeremy L Thompson   switch (t_mode) {
1576c8c3fa7dSJeremy L Thompson     case CEED_NOTRANSPOSE: {
1577c8c3fa7dSJeremy L Thompson       // Nodes to arbitrary points
1578c8c3fa7dSJeremy L Thompson       CeedScalar       *v_array;
1579c8c3fa7dSJeremy L Thompson       const CeedScalar *chebyshev_coeffs, *x_array_read;
1580c8c3fa7dSJeremy L Thompson 
1581c8c3fa7dSJeremy L Thompson       // -- Interpolate to Chebyshev coefficients
1582c8c3fa7dSJeremy L Thompson       CeedCall(CeedBasisApply(basis->basis_chebyshev, 1, CEED_NOTRANSPOSE, CEED_EVAL_INTERP, u, basis->vec_chebyshev));
1583c8c3fa7dSJeremy L Thompson 
1584c8c3fa7dSJeremy L Thompson       // -- Evaluate Chebyshev polynomials at arbitrary points
1585c8c3fa7dSJeremy L Thompson       CeedCall(CeedVectorGetArrayRead(basis->vec_chebyshev, CEED_MEM_HOST, &chebyshev_coeffs));
1586c8c3fa7dSJeremy L Thompson       CeedCall(CeedVectorGetArrayRead(x_ref, CEED_MEM_HOST, &x_array_read));
1587c8c3fa7dSJeremy L Thompson       CeedCall(CeedVectorGetArrayWrite(v, CEED_MEM_HOST, &v_array));
1588edfbf3a6SJeremy L Thompson       switch (eval_mode) {
1589edfbf3a6SJeremy L Thompson         case CEED_EVAL_INTERP: {
1590c8c3fa7dSJeremy L Thompson           CeedScalar tmp[2][num_comp * CeedIntPow(Q_1d, dim)], chebyshev_x[Q_1d];
1591c8c3fa7dSJeremy L Thompson 
1592c8c3fa7dSJeremy L Thompson           // ---- Values at point
1593c8c3fa7dSJeremy L Thompson           for (CeedInt p = 0; p < num_points; p++) {
1594c8c3fa7dSJeremy L Thompson             CeedInt pre = num_comp * CeedIntPow(Q_1d, dim - 1), post = 1;
1595c8c3fa7dSJeremy L Thompson 
1596c8c3fa7dSJeremy L Thompson             for (CeedInt d = dim - 1; d >= 0; d--) {
1597*3778dbaaSJeremy L Thompson               // ------ Tensor contract with current Chebyshev polynomial values
1598*3778dbaaSJeremy L Thompson               CeedCall(CeedChebyshevPolynomialsAtPoint(x_array_read[p * dim + d], Q_1d, chebyshev_x));
1599c8c3fa7dSJeremy L Thompson               CeedCall(CeedTensorContractApply(basis->contract, pre, Q_1d, post, 1, chebyshev_x, t_mode, false,
1600c8c3fa7dSJeremy L Thompson                                                d == (dim - 1) ? chebyshev_coeffs : tmp[d % 2], d == 0 ? &v_array[p * num_comp] : tmp[(d + 1) % 2]));
1601c8c3fa7dSJeremy L Thompson               pre /= Q_1d;
1602c8c3fa7dSJeremy L Thompson               post *= 1;
1603c8c3fa7dSJeremy L Thompson             }
1604c8c3fa7dSJeremy L Thompson           }
1605edfbf3a6SJeremy L Thompson           break;
1606edfbf3a6SJeremy L Thompson         }
1607edfbf3a6SJeremy L Thompson         case CEED_EVAL_GRAD: {
1608edfbf3a6SJeremy L Thompson           CeedScalar tmp[2][num_comp * CeedIntPow(Q_1d, dim)], chebyshev_x[Q_1d];
1609edfbf3a6SJeremy L Thompson 
1610edfbf3a6SJeremy L Thompson           // ---- Values at point
1611edfbf3a6SJeremy L Thompson           for (CeedInt p = 0; p < num_points; p++) {
1612edfbf3a6SJeremy L Thompson             // Dim**2 contractions, apply grad when pass == dim
1613edfbf3a6SJeremy L Thompson             for (CeedInt pass = dim - 1; pass >= 0; pass--) {
1614edfbf3a6SJeremy L Thompson               CeedInt pre = num_comp * CeedIntPow(Q_1d, dim - 1), post = 1;
1615edfbf3a6SJeremy L Thompson 
1616edfbf3a6SJeremy L Thompson               for (CeedInt d = dim - 1; d >= 0; d--) {
1617*3778dbaaSJeremy L Thompson                 // ------ Tensor contract with current Chebyshev polynomial values
1618*3778dbaaSJeremy L Thompson                 if (pass == d) CeedCall(CeedChebyshevDerivativeAtPoint(x_array_read[p * dim + d], Q_1d, chebyshev_x));
1619*3778dbaaSJeremy L Thompson                 else CeedCall(CeedChebyshevPolynomialsAtPoint(x_array_read[p * dim + d], Q_1d, chebyshev_x));
1620edfbf3a6SJeremy L Thompson                 CeedCall(CeedTensorContractApply(basis->contract, pre, Q_1d, post, 1, chebyshev_x, t_mode, false,
1621edfbf3a6SJeremy L Thompson                                                  d == (dim - 1) ? chebyshev_coeffs : tmp[d % 2],
1622edfbf3a6SJeremy L Thompson                                                  d == 0 ? &v_array[p * num_comp * dim + pass] : tmp[(d + 1) % 2]));
1623edfbf3a6SJeremy L Thompson                 pre /= Q_1d;
1624edfbf3a6SJeremy L Thompson                 post *= 1;
1625edfbf3a6SJeremy L Thompson               }
1626edfbf3a6SJeremy L Thompson             }
1627edfbf3a6SJeremy L Thompson           }
1628edfbf3a6SJeremy L Thompson           break;
1629edfbf3a6SJeremy L Thompson         }
1630edfbf3a6SJeremy L Thompson         default:
1631edfbf3a6SJeremy L Thompson           // Nothing to do, this won't occur
1632edfbf3a6SJeremy L Thompson           break;
1633c8c3fa7dSJeremy L Thompson       }
1634c8c3fa7dSJeremy L Thompson       CeedCall(CeedVectorRestoreArrayRead(basis->vec_chebyshev, &chebyshev_coeffs));
1635c8c3fa7dSJeremy L Thompson       CeedCall(CeedVectorRestoreArrayRead(x_ref, &x_array_read));
1636c8c3fa7dSJeremy L Thompson       CeedCall(CeedVectorRestoreArray(v, &v_array));
1637c8c3fa7dSJeremy L Thompson       break;
1638c8c3fa7dSJeremy L Thompson     }
16392a94f45fSJeremy L Thompson     case CEED_TRANSPOSE: {
1640*3778dbaaSJeremy L Thompson       // Note: No switch on e_mode here because only CEED_EVAL_INTERP is supported at this time
16412a94f45fSJeremy L Thompson       // Arbitrary points to nodes
16422a94f45fSJeremy L Thompson       CeedScalar       *chebyshev_coeffs;
16432a94f45fSJeremy L Thompson       const CeedScalar *u_array, *x_array_read;
16442a94f45fSJeremy L Thompson 
16452a94f45fSJeremy L Thompson       // -- Transpose of evaluaton of Chebyshev polynomials at arbitrary points
16462a94f45fSJeremy L Thompson       CeedCall(CeedVectorGetArrayWrite(basis->vec_chebyshev, CEED_MEM_HOST, &chebyshev_coeffs));
16472a94f45fSJeremy L Thompson       CeedCall(CeedVectorGetArrayRead(x_ref, CEED_MEM_HOST, &x_array_read));
16482a94f45fSJeremy L Thompson       CeedCall(CeedVectorGetArrayRead(u, CEED_MEM_HOST, &u_array));
16492a94f45fSJeremy L Thompson       {
16502a94f45fSJeremy L Thompson         CeedScalar tmp[2][num_comp * CeedIntPow(Q_1d, dim)], chebyshev_x[Q_1d];
16512a94f45fSJeremy L Thompson 
16522a94f45fSJeremy L Thompson         // ---- Values at point
16532a94f45fSJeremy L Thompson         for (CeedInt p = 0; p < num_points; p++) {
16542a94f45fSJeremy L Thompson           CeedInt pre = num_comp * 1, post = 1;
16552a94f45fSJeremy L Thompson 
16562a94f45fSJeremy L Thompson           for (CeedInt d = dim - 1; d >= 0; d--) {
1657*3778dbaaSJeremy L Thompson             // ------ Tensor contract with current Chebyshev polynomial values
1658*3778dbaaSJeremy L Thompson             CeedCall(CeedChebyshevPolynomialsAtPoint(x_array_read[p * dim + d], Q_1d, chebyshev_x));
16592a94f45fSJeremy L Thompson             CeedCall(CeedTensorContractApply(basis->contract, pre, 1, post, Q_1d, chebyshev_x, t_mode, p > 0 && d == 0,
16602a94f45fSJeremy L Thompson                                              d == (dim - 1) ? &u_array[p * num_comp] : tmp[d % 2], d == 0 ? chebyshev_coeffs : tmp[(d + 1) % 2]));
16612a94f45fSJeremy L Thompson             pre /= 1;
16622a94f45fSJeremy L Thompson             post *= Q_1d;
16632a94f45fSJeremy L Thompson           }
16642a94f45fSJeremy L Thompson         }
16652a94f45fSJeremy L Thompson       }
16662a94f45fSJeremy L Thompson       CeedCall(CeedVectorRestoreArray(basis->vec_chebyshev, &chebyshev_coeffs));
16672a94f45fSJeremy L Thompson       CeedCall(CeedVectorRestoreArrayRead(x_ref, &x_array_read));
16682a94f45fSJeremy L Thompson       CeedCall(CeedVectorRestoreArrayRead(u, &u_array));
16692a94f45fSJeremy L Thompson 
16702a94f45fSJeremy L Thompson       // -- Interpolate transpose from Chebyshev coefficients
16712a94f45fSJeremy L Thompson       CeedCall(CeedBasisApply(basis->basis_chebyshev, 1, CEED_TRANSPOSE, CEED_EVAL_INTERP, basis->vec_chebyshev, v));
16722a94f45fSJeremy L Thompson       break;
16732a94f45fSJeremy L Thompson     }
1674c8c3fa7dSJeremy L Thompson   }
1675c8c3fa7dSJeremy L Thompson 
1676c8c3fa7dSJeremy L Thompson   return CEED_ERROR_SUCCESS;
1677c8c3fa7dSJeremy L Thompson }
1678c8c3fa7dSJeremy L Thompson 
1679c8c3fa7dSJeremy L Thompson /**
1680b7c9bbdaSJeremy L Thompson   @brief Get Ceed associated with a CeedBasis
1681b7c9bbdaSJeremy L Thompson 
1682ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
1683b7c9bbdaSJeremy L Thompson   @param[out] ceed  Variable to store Ceed
1684b7c9bbdaSJeremy L Thompson 
1685b7c9bbdaSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1686b7c9bbdaSJeremy L Thompson 
1687b7c9bbdaSJeremy L Thompson   @ref Advanced
1688b7c9bbdaSJeremy L Thompson **/
1689b7c9bbdaSJeremy L Thompson int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) {
1690b7c9bbdaSJeremy L Thompson   *ceed = basis->ceed;
1691b7c9bbdaSJeremy L Thompson   return CEED_ERROR_SUCCESS;
1692b7c9bbdaSJeremy L Thompson }
1693b7c9bbdaSJeremy L Thompson 
1694b7c9bbdaSJeremy L Thompson /**
16959d007619Sjeremylt   @brief Get dimension for given CeedBasis
16969d007619Sjeremylt 
1697ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
16989d007619Sjeremylt   @param[out] dim   Variable to store dimension of basis
16999d007619Sjeremylt 
17009d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
17019d007619Sjeremylt 
1702b7c9bbdaSJeremy L Thompson   @ref Advanced
17039d007619Sjeremylt **/
17049d007619Sjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) {
17059d007619Sjeremylt   *dim = basis->dim;
1706e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
17079d007619Sjeremylt }
17089d007619Sjeremylt 
17099d007619Sjeremylt /**
1710d99fa3c5SJeremy L Thompson   @brief Get topology for given CeedBasis
1711d99fa3c5SJeremy L Thompson 
1712ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
1713d99fa3c5SJeremy L Thompson   @param[out] topo  Variable to store topology of basis
1714d99fa3c5SJeremy L Thompson 
1715d99fa3c5SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1716d99fa3c5SJeremy L Thompson 
1717b7c9bbdaSJeremy L Thompson   @ref Advanced
1718d99fa3c5SJeremy L Thompson **/
1719d99fa3c5SJeremy L Thompson int CeedBasisGetTopology(CeedBasis basis, CeedElemTopology *topo) {
1720d99fa3c5SJeremy L Thompson   *topo = basis->topo;
1721e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1722d99fa3c5SJeremy L Thompson }
1723d99fa3c5SJeremy L Thompson 
1724d99fa3c5SJeremy L Thompson /**
17259d007619Sjeremylt   @brief Get number of components for given CeedBasis
17269d007619Sjeremylt 
1727ea61e9acSJeremy L Thompson   @param[in]  basis    CeedBasis
1728d1d35e2fSjeremylt   @param[out] num_comp Variable to store number of components of basis
17299d007619Sjeremylt 
17309d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
17319d007619Sjeremylt 
1732b7c9bbdaSJeremy L Thompson   @ref Advanced
17339d007619Sjeremylt **/
1734d1d35e2fSjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *num_comp) {
1735d1d35e2fSjeremylt   *num_comp = basis->num_comp;
1736e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
17379d007619Sjeremylt }
17389d007619Sjeremylt 
17399d007619Sjeremylt /**
17409d007619Sjeremylt   @brief Get total number of nodes (in dim dimensions) of a CeedBasis
17419d007619Sjeremylt 
1742ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
17439d007619Sjeremylt   @param[out] P     Variable to store number of nodes
17449d007619Sjeremylt 
17459d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
17469d007619Sjeremylt 
17479d007619Sjeremylt   @ref Utility
17489d007619Sjeremylt **/
17499d007619Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) {
17509d007619Sjeremylt   *P = basis->P;
1751e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
17529d007619Sjeremylt }
17539d007619Sjeremylt 
17549d007619Sjeremylt /**
17559d007619Sjeremylt   @brief Get total number of nodes (in 1 dimension) of a CeedBasis
17569d007619Sjeremylt 
1757ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
1758d1d35e2fSjeremylt   @param[out] P_1d  Variable to store number of nodes
17599d007619Sjeremylt 
17609d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
17619d007619Sjeremylt 
1762b7c9bbdaSJeremy L Thompson   @ref Advanced
17639d007619Sjeremylt **/
1764d1d35e2fSjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P_1d) {
17656402da51SJeremy L Thompson   CeedCheck(basis->is_tensor_basis, basis->ceed, CEED_ERROR_MINOR, "Cannot supply P_1d for non-tensor basis");
1766d1d35e2fSjeremylt   *P_1d = basis->P_1d;
1767e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
17689d007619Sjeremylt }
17699d007619Sjeremylt 
17709d007619Sjeremylt /**
17719d007619Sjeremylt   @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis
17729d007619Sjeremylt 
1773ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
17749d007619Sjeremylt   @param[out] Q     Variable to store number of quadrature points
17759d007619Sjeremylt 
17769d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
17779d007619Sjeremylt 
17789d007619Sjeremylt   @ref Utility
17799d007619Sjeremylt **/
17809d007619Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) {
17819d007619Sjeremylt   *Q = basis->Q;
1782e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
17839d007619Sjeremylt }
17849d007619Sjeremylt 
17859d007619Sjeremylt /**
17869d007619Sjeremylt   @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis
17879d007619Sjeremylt 
1788ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
1789d1d35e2fSjeremylt   @param[out] Q_1d  Variable to store number of quadrature points
17909d007619Sjeremylt 
17919d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
17929d007619Sjeremylt 
1793b7c9bbdaSJeremy L Thompson   @ref Advanced
17949d007619Sjeremylt **/
1795d1d35e2fSjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q_1d) {
17966402da51SJeremy L Thompson   CeedCheck(basis->is_tensor_basis, basis->ceed, CEED_ERROR_MINOR, "Cannot supply Q_1d for non-tensor basis");
1797d1d35e2fSjeremylt   *Q_1d = basis->Q_1d;
1798e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
17999d007619Sjeremylt }
18009d007619Sjeremylt 
18019d007619Sjeremylt /**
1802ea61e9acSJeremy L Thompson   @brief Get reference coordinates of quadrature points (in dim dimensions) of a CeedBasis
18039d007619Sjeremylt 
1804ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
1805d1d35e2fSjeremylt   @param[out] q_ref Variable to store reference coordinates of quadrature points
18069d007619Sjeremylt 
18079d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
18089d007619Sjeremylt 
1809b7c9bbdaSJeremy L Thompson   @ref Advanced
18109d007619Sjeremylt **/
1811d1d35e2fSjeremylt int CeedBasisGetQRef(CeedBasis basis, const CeedScalar **q_ref) {
1812d1d35e2fSjeremylt   *q_ref = basis->q_ref_1d;
1813e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
18149d007619Sjeremylt }
18159d007619Sjeremylt 
18169d007619Sjeremylt /**
1817ea61e9acSJeremy L Thompson   @brief Get quadrature weights of quadrature points (in dim dimensions) of a CeedBasis
18189d007619Sjeremylt 
1819ea61e9acSJeremy L Thompson   @param[in]  basis    CeedBasis
1820d1d35e2fSjeremylt   @param[out] q_weight Variable to store quadrature weights
18219d007619Sjeremylt 
18229d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
18239d007619Sjeremylt 
1824b7c9bbdaSJeremy L Thompson   @ref Advanced
18259d007619Sjeremylt **/
1826d1d35e2fSjeremylt int CeedBasisGetQWeights(CeedBasis basis, const CeedScalar **q_weight) {
1827d1d35e2fSjeremylt   *q_weight = basis->q_weight_1d;
1828e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
18299d007619Sjeremylt }
18309d007619Sjeremylt 
18319d007619Sjeremylt /**
18329d007619Sjeremylt   @brief Get interpolation matrix of a CeedBasis
18339d007619Sjeremylt 
1834ea61e9acSJeremy L Thompson   @param[in]  basis  CeedBasis
18359d007619Sjeremylt   @param[out] interp Variable to store interpolation matrix
18369d007619Sjeremylt 
18379d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
18389d007619Sjeremylt 
1839b7c9bbdaSJeremy L Thompson   @ref Advanced
18409d007619Sjeremylt **/
18416c58de82SJeremy L Thompson int CeedBasisGetInterp(CeedBasis basis, const CeedScalar **interp) {
18426402da51SJeremy L Thompson   if (!basis->interp && basis->is_tensor_basis) {
18439d007619Sjeremylt     // Allocate
18442b730f8bSJeremy L Thompson     CeedCall(CeedMalloc(basis->Q * basis->P, &basis->interp));
18459d007619Sjeremylt 
18469d007619Sjeremylt     // Initialize
18472b730f8bSJeremy L Thompson     for (CeedInt i = 0; i < basis->Q * basis->P; i++) basis->interp[i] = 1.0;
18489d007619Sjeremylt 
18499d007619Sjeremylt     // Calculate
18502b730f8bSJeremy L Thompson     for (CeedInt d = 0; d < basis->dim; d++) {
18512b730f8bSJeremy L Thompson       for (CeedInt qpt = 0; qpt < basis->Q; qpt++) {
18529d007619Sjeremylt         for (CeedInt node = 0; node < basis->P; node++) {
1853d1d35e2fSjeremylt           CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d;
1854d1d35e2fSjeremylt           CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d;
1855d1d35e2fSjeremylt           basis->interp[qpt * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p];
18569d007619Sjeremylt         }
18579d007619Sjeremylt       }
18582b730f8bSJeremy L Thompson     }
18592b730f8bSJeremy L Thompson   }
18609d007619Sjeremylt   *interp = basis->interp;
1861e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
18629d007619Sjeremylt }
18639d007619Sjeremylt 
18649d007619Sjeremylt /**
18659d007619Sjeremylt   @brief Get 1D interpolation matrix of a tensor product CeedBasis
18669d007619Sjeremylt 
1867ea61e9acSJeremy L Thompson   @param[in]  basis     CeedBasis
1868d1d35e2fSjeremylt   @param[out] interp_1d Variable to store interpolation matrix
18699d007619Sjeremylt 
18709d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
18719d007619Sjeremylt 
18729d007619Sjeremylt   @ref Backend
18739d007619Sjeremylt **/
1874d1d35e2fSjeremylt int CeedBasisGetInterp1D(CeedBasis basis, const CeedScalar **interp_1d) {
18756402da51SJeremy L Thompson   CeedCheck(basis->is_tensor_basis, basis->ceed, CEED_ERROR_MINOR, "CeedBasis is not a tensor product basis.");
1876d1d35e2fSjeremylt   *interp_1d = basis->interp_1d;
1877e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
18789d007619Sjeremylt }
18799d007619Sjeremylt 
18809d007619Sjeremylt /**
18819d007619Sjeremylt   @brief Get gradient matrix of a CeedBasis
18829d007619Sjeremylt 
1883ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
18849d007619Sjeremylt   @param[out] grad  Variable to store gradient matrix
18859d007619Sjeremylt 
18869d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
18879d007619Sjeremylt 
1888b7c9bbdaSJeremy L Thompson   @ref Advanced
18899d007619Sjeremylt **/
18906c58de82SJeremy L Thompson int CeedBasisGetGrad(CeedBasis basis, const CeedScalar **grad) {
18916402da51SJeremy L Thompson   if (!basis->grad && basis->is_tensor_basis) {
18929d007619Sjeremylt     // Allocate
18932b730f8bSJeremy L Thompson     CeedCall(CeedMalloc(basis->dim * basis->Q * basis->P, &basis->grad));
18949d007619Sjeremylt 
18959d007619Sjeremylt     // Initialize
18962b730f8bSJeremy L Thompson     for (CeedInt i = 0; i < basis->dim * basis->Q * basis->P; i++) basis->grad[i] = 1.0;
18979d007619Sjeremylt 
18989d007619Sjeremylt     // Calculate
18992b730f8bSJeremy L Thompson     for (CeedInt d = 0; d < basis->dim; d++) {
19002b730f8bSJeremy L Thompson       for (CeedInt i = 0; i < basis->dim; i++) {
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;
19052b730f8bSJeremy L Thompson             if (i == d) basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->grad_1d[q * basis->P_1d + p];
19062b730f8bSJeremy L Thompson             else basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p];
19072b730f8bSJeremy L Thompson           }
19082b730f8bSJeremy L Thompson         }
19092b730f8bSJeremy L Thompson       }
19109d007619Sjeremylt     }
19119d007619Sjeremylt   }
19129d007619Sjeremylt   *grad = basis->grad;
1913e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
19149d007619Sjeremylt }
19159d007619Sjeremylt 
19169d007619Sjeremylt /**
19179d007619Sjeremylt   @brief Get 1D gradient matrix of a tensor product CeedBasis
19189d007619Sjeremylt 
1919ea61e9acSJeremy L Thompson   @param[in]  basis   CeedBasis
1920d1d35e2fSjeremylt   @param[out] grad_1d Variable to store gradient matrix
19219d007619Sjeremylt 
19229d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
19239d007619Sjeremylt 
1924b7c9bbdaSJeremy L Thompson   @ref Advanced
19259d007619Sjeremylt **/
1926d1d35e2fSjeremylt int CeedBasisGetGrad1D(CeedBasis basis, const CeedScalar **grad_1d) {
19276402da51SJeremy L Thompson   CeedCheck(basis->is_tensor_basis, basis->ceed, CEED_ERROR_MINOR, "CeedBasis is not a tensor product basis.");
1928d1d35e2fSjeremylt   *grad_1d = basis->grad_1d;
1929e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
19309d007619Sjeremylt }
19319d007619Sjeremylt 
19329d007619Sjeremylt /**
193350c301a5SRezgar Shakeri   @brief Get divergence matrix of a CeedBasis
193450c301a5SRezgar Shakeri 
1935ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
193650c301a5SRezgar Shakeri   @param[out] div   Variable to store divergence matrix
193750c301a5SRezgar Shakeri 
193850c301a5SRezgar Shakeri   @return An error code: 0 - success, otherwise - failure
193950c301a5SRezgar Shakeri 
194050c301a5SRezgar Shakeri   @ref Advanced
194150c301a5SRezgar Shakeri **/
194250c301a5SRezgar Shakeri int CeedBasisGetDiv(CeedBasis basis, const CeedScalar **div) {
19436574a04fSJeremy L Thompson   CeedCheck(basis->div, basis->ceed, CEED_ERROR_MINOR, "CeedBasis does not have divergence matrix.");
194450c301a5SRezgar Shakeri   *div = basis->div;
194550c301a5SRezgar Shakeri   return CEED_ERROR_SUCCESS;
194650c301a5SRezgar Shakeri }
194750c301a5SRezgar Shakeri 
194850c301a5SRezgar Shakeri /**
1949c4e3f59bSSebastian Grimberg   @brief Get curl matrix of a CeedBasis
1950c4e3f59bSSebastian Grimberg 
1951c4e3f59bSSebastian Grimberg   @param[in]  basis CeedBasis
1952c4e3f59bSSebastian Grimberg   @param[out] curl  Variable to store curl matrix
1953c4e3f59bSSebastian Grimberg 
1954c4e3f59bSSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
1955c4e3f59bSSebastian Grimberg 
1956c4e3f59bSSebastian Grimberg   @ref Advanced
1957c4e3f59bSSebastian Grimberg **/
1958c4e3f59bSSebastian Grimberg int CeedBasisGetCurl(CeedBasis basis, const CeedScalar **curl) {
19596574a04fSJeremy L Thompson   CeedCheck(basis->curl, basis->ceed, CEED_ERROR_MINOR, "CeedBasis does not have curl matrix.");
1960c4e3f59bSSebastian Grimberg   *curl = basis->curl;
1961c4e3f59bSSebastian Grimberg   return CEED_ERROR_SUCCESS;
1962c4e3f59bSSebastian Grimberg }
1963c4e3f59bSSebastian Grimberg 
1964c4e3f59bSSebastian Grimberg /**
19657a982d89SJeremy L. Thompson   @brief Destroy a CeedBasis
19667a982d89SJeremy L. Thompson 
1967ea61e9acSJeremy L Thompson   @param[in,out] basis CeedBasis to destroy
19687a982d89SJeremy L. Thompson 
19697a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
19707a982d89SJeremy L. Thompson 
19717a982d89SJeremy L. Thompson   @ref User
19727a982d89SJeremy L. Thompson **/
19737a982d89SJeremy L. Thompson int CeedBasisDestroy(CeedBasis *basis) {
19747425e127SJeremy L Thompson   if (!*basis || *basis == CEED_BASIS_COLLOCATED || --(*basis)->ref_count > 0) {
1975ad6481ceSJeremy L Thompson     *basis = NULL;
1976ad6481ceSJeremy L Thompson     return CEED_ERROR_SUCCESS;
1977ad6481ceSJeremy L Thompson   }
19782b730f8bSJeremy L Thompson   if ((*basis)->Destroy) CeedCall((*basis)->Destroy(*basis));
19799831d45aSJeremy L Thompson   CeedCall(CeedTensorContractDestroy(&(*basis)->contract));
1980c4e3f59bSSebastian Grimberg   CeedCall(CeedFree(&(*basis)->q_ref_1d));
1981c4e3f59bSSebastian Grimberg   CeedCall(CeedFree(&(*basis)->q_weight_1d));
19822b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->interp));
19832b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->interp_1d));
19842b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->grad));
19852b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->grad_1d));
1986c4e3f59bSSebastian Grimberg   CeedCall(CeedFree(&(*basis)->div));
1987c4e3f59bSSebastian Grimberg   CeedCall(CeedFree(&(*basis)->curl));
1988c8c3fa7dSJeremy L Thompson   CeedCall(CeedVectorDestroy(&(*basis)->vec_chebyshev));
1989c8c3fa7dSJeremy L Thompson   CeedCall(CeedBasisDestroy(&(*basis)->basis_chebyshev));
19902b730f8bSJeremy L Thompson   CeedCall(CeedDestroy(&(*basis)->ceed));
19912b730f8bSJeremy L Thompson   CeedCall(CeedFree(basis));
1992e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
19937a982d89SJeremy L. Thompson }
19947a982d89SJeremy L. Thompson 
19957a982d89SJeremy L. Thompson /**
1996b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre quadrature
1997b11c1e72Sjeremylt 
1998ea61e9acSJeremy L Thompson   @param[in]  Q           Number of quadrature points (integrates polynomials of degree 2*Q-1 exactly)
1999d1d35e2fSjeremylt   @param[out] q_ref_1d    Array of length Q to hold the abscissa on [-1, 1]
2000d1d35e2fSjeremylt   @param[out] q_weight_1d Array of length Q to hold the weights
2001b11c1e72Sjeremylt 
2002b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
2003dfdf5a53Sjeremylt 
2004dfdf5a53Sjeremylt   @ref Utility
2005b11c1e72Sjeremylt **/
20062b730f8bSJeremy L Thompson int CeedGaussQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) {
2007d7b241e6Sjeremylt   // Allocate
2008d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0 * atan(1.0);
2009d1d35e2fSjeremylt   // Build q_ref_1d, q_weight_1d
201092ae7e47SJeremy L Thompson   for (CeedInt i = 0; i <= Q / 2; i++) {
2011d7b241e6Sjeremylt     // Guess
2012d7b241e6Sjeremylt     xi = cos(PI * (CeedScalar)(2 * i + 1) / ((CeedScalar)(2 * Q)));
2013d7b241e6Sjeremylt     // Pn(xi)
2014d7b241e6Sjeremylt     P0 = 1.0;
2015d7b241e6Sjeremylt     P1 = xi;
2016d7b241e6Sjeremylt     P2 = 0.0;
201792ae7e47SJeremy L Thompson     for (CeedInt j = 2; j <= Q; j++) {
2018d7b241e6Sjeremylt       P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
2019d7b241e6Sjeremylt       P0 = P1;
2020d7b241e6Sjeremylt       P1 = P2;
2021d7b241e6Sjeremylt     }
2022d7b241e6Sjeremylt     // First Newton Step
2023d7b241e6Sjeremylt     dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
2024d7b241e6Sjeremylt     xi  = xi - P2 / dP2;
2025d7b241e6Sjeremylt     // Newton to convergence
202692ae7e47SJeremy L Thompson     for (CeedInt k = 0; k < 100 && fabs(P2) > 10 * CEED_EPSILON; k++) {
2027d7b241e6Sjeremylt       P0 = 1.0;
2028d7b241e6Sjeremylt       P1 = xi;
202992ae7e47SJeremy L Thompson       for (CeedInt j = 2; j <= Q; j++) {
2030d7b241e6Sjeremylt         P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
2031d7b241e6Sjeremylt         P0 = P1;
2032d7b241e6Sjeremylt         P1 = P2;
2033d7b241e6Sjeremylt       }
2034d7b241e6Sjeremylt       dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
2035d7b241e6Sjeremylt       xi  = xi - P2 / dP2;
2036d7b241e6Sjeremylt     }
2037d7b241e6Sjeremylt     // Save xi, wi
2038d7b241e6Sjeremylt     wi                     = 2.0 / ((1.0 - xi * xi) * dP2 * dP2);
2039d1d35e2fSjeremylt     q_weight_1d[i]         = wi;
2040d1d35e2fSjeremylt     q_weight_1d[Q - 1 - i] = wi;
2041d1d35e2fSjeremylt     q_ref_1d[i]            = -xi;
2042d1d35e2fSjeremylt     q_ref_1d[Q - 1 - i]    = xi;
2043d7b241e6Sjeremylt   }
2044e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2045d7b241e6Sjeremylt }
2046d7b241e6Sjeremylt 
2047b11c1e72Sjeremylt /**
2048b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre-Lobatto quadrature
2049b11c1e72Sjeremylt 
2050ea61e9acSJeremy L Thompson   @param[in]  Q           Number of quadrature points (integrates polynomials of degree 2*Q-3 exactly)
2051d1d35e2fSjeremylt   @param[out] q_ref_1d    Array of length Q to hold the abscissa on [-1, 1]
2052d1d35e2fSjeremylt   @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 CeedLobattoQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) {
2059d7b241e6Sjeremylt   // Allocate
2060d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0 * atan(1.0);
2061d1d35e2fSjeremylt   // Build q_ref_1d, q_weight_1d
2062d7b241e6Sjeremylt   // Set endpoints
20636574a04fSJeremy L Thompson   CeedCheck(Q > 1, NULL, CEED_ERROR_DIMENSION, "Cannot create Lobatto quadrature with Q=%" CeedInt_FMT " < 2 points", Q);
2064d7b241e6Sjeremylt   wi = 2.0 / ((CeedScalar)(Q * (Q - 1)));
2065d1d35e2fSjeremylt   if (q_weight_1d) {
2066d1d35e2fSjeremylt     q_weight_1d[0]     = wi;
2067d1d35e2fSjeremylt     q_weight_1d[Q - 1] = wi;
2068d7b241e6Sjeremylt   }
2069d1d35e2fSjeremylt   q_ref_1d[0]     = -1.0;
2070d1d35e2fSjeremylt   q_ref_1d[Q - 1] = 1.0;
2071d7b241e6Sjeremylt   // Interior
207292ae7e47SJeremy L Thompson   for (CeedInt i = 1; i <= (Q - 1) / 2; i++) {
2073d7b241e6Sjeremylt     // Guess
2074d7b241e6Sjeremylt     xi = cos(PI * (CeedScalar)(i) / (CeedScalar)(Q - 1));
2075d7b241e6Sjeremylt     // Pn(xi)
2076d7b241e6Sjeremylt     P0 = 1.0;
2077d7b241e6Sjeremylt     P1 = xi;
2078d7b241e6Sjeremylt     P2 = 0.0;
207992ae7e47SJeremy L Thompson     for (CeedInt j = 2; j < Q; j++) {
2080d7b241e6Sjeremylt       P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
2081d7b241e6Sjeremylt       P0 = P1;
2082d7b241e6Sjeremylt       P1 = P2;
2083d7b241e6Sjeremylt     }
2084d7b241e6Sjeremylt     // First Newton step
2085d7b241e6Sjeremylt     dP2  = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
2086d7b241e6Sjeremylt     d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi);
2087d7b241e6Sjeremylt     xi   = xi - dP2 / d2P2;
2088d7b241e6Sjeremylt     // Newton to convergence
208992ae7e47SJeremy L Thompson     for (CeedInt k = 0; k < 100 && fabs(dP2) > 10 * CEED_EPSILON; k++) {
2090d7b241e6Sjeremylt       P0 = 1.0;
2091d7b241e6Sjeremylt       P1 = xi;
209292ae7e47SJeremy L Thompson       for (CeedInt j = 2; j < Q; j++) {
2093d7b241e6Sjeremylt         P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
2094d7b241e6Sjeremylt         P0 = P1;
2095d7b241e6Sjeremylt         P1 = P2;
2096d7b241e6Sjeremylt       }
2097d7b241e6Sjeremylt       dP2  = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
2098d7b241e6Sjeremylt       d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi);
2099d7b241e6Sjeremylt       xi   = xi - dP2 / d2P2;
2100d7b241e6Sjeremylt     }
2101d7b241e6Sjeremylt     // Save xi, wi
2102d7b241e6Sjeremylt     wi = 2.0 / (((CeedScalar)(Q * (Q - 1))) * P2 * P2);
2103d1d35e2fSjeremylt     if (q_weight_1d) {
2104d1d35e2fSjeremylt       q_weight_1d[i]         = wi;
2105d1d35e2fSjeremylt       q_weight_1d[Q - 1 - i] = wi;
2106d7b241e6Sjeremylt     }
2107d1d35e2fSjeremylt     q_ref_1d[i]         = -xi;
2108d1d35e2fSjeremylt     q_ref_1d[Q - 1 - i] = xi;
2109d7b241e6Sjeremylt   }
2110e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2111d7b241e6Sjeremylt }
2112d7b241e6Sjeremylt 
2113d7b241e6Sjeremylt /// @}
2114