xref: /libCEED/rust/libceed-sys/c-src/interface/ceed-basis.c (revision ea61e9ac44808524e4667c1525a05976f536c19c)
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>
92b730f8bSJeremy L Thompson #include <ceed/backend.h>
102b730f8bSJeremy L Thompson #include <ceed/ceed.h>
11d7b241e6Sjeremylt #include <math.h>
123d576824SJeremy L Thompson #include <stdbool.h>
13d7b241e6Sjeremylt #include <stdio.h>
14d7b241e6Sjeremylt #include <string.h>
15d7b241e6Sjeremylt 
167a982d89SJeremy L. Thompson /// @file
177a982d89SJeremy L. Thompson /// Implementation of CeedBasis interfaces
187a982d89SJeremy L. Thompson 
19d7b241e6Sjeremylt /// @cond DOXYGEN_SKIP
20783c99b3SValeria Barra static struct CeedBasis_private ceed_basis_collocated;
21d7b241e6Sjeremylt /// @endcond
22d7b241e6Sjeremylt 
237a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser
247a982d89SJeremy L. Thompson /// @{
257a982d89SJeremy L. Thompson 
267a982d89SJeremy L. Thompson /// Indicate that the quadrature points are collocated with the nodes
277a982d89SJeremy L. Thompson const CeedBasis CEED_BASIS_COLLOCATED = &ceed_basis_collocated;
287a982d89SJeremy L. Thompson 
297a982d89SJeremy L. Thompson /// @}
307a982d89SJeremy L. Thompson 
317a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
327a982d89SJeremy L. Thompson /// CeedBasis Library Internal Functions
337a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
347a982d89SJeremy L. Thompson /// @addtogroup CeedBasisDeveloper
357a982d89SJeremy L. Thompson /// @{
367a982d89SJeremy L. Thompson 
377a982d89SJeremy L. Thompson /**
387a982d89SJeremy L. Thompson   @brief Compute Householder reflection
397a982d89SJeremy L. Thompson 
40*ea61e9acSJeremy L Thompson     Computes A = (I - b v v^T) A, where A is an mxn matrix indexed as A[i*row + j*col]
417a982d89SJeremy L. Thompson 
427a982d89SJeremy L. Thompson   @param[in,out] A   Matrix to apply Householder reflection to, in place
43*ea61e9acSJeremy L Thompson   @param[in]     v   Householder vector
44*ea61e9acSJeremy L Thompson   @param[in]     b   Scaling factor
45*ea61e9acSJeremy L Thompson   @param[in]     m   Number of rows in A
46*ea61e9acSJeremy L Thompson   @param[in]     n   Number of columns in A
47*ea61e9acSJeremy L Thompson   @param[in]     row Row stride
48*ea61e9acSJeremy L Thompson   @param[in]     col Col stride
497a982d89SJeremy L. Thompson 
507a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
517a982d89SJeremy L. Thompson 
527a982d89SJeremy L. Thompson   @ref Developer
537a982d89SJeremy L. Thompson **/
542b730f8bSJeremy L Thompson static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v, CeedScalar b, CeedInt m, CeedInt n, CeedInt row, CeedInt col) {
557a982d89SJeremy L. Thompson   for (CeedInt j = 0; j < n; j++) {
567a982d89SJeremy L. Thompson     CeedScalar w = A[0 * row + j * col];
572b730f8bSJeremy L Thompson     for (CeedInt i = 1; i < m; i++) w += v[i] * A[i * row + j * col];
587a982d89SJeremy L. Thompson     A[0 * row + j * col] -= b * w;
592b730f8bSJeremy L Thompson     for (CeedInt i = 1; i < m; i++) A[i * row + j * col] -= b * w * v[i];
607a982d89SJeremy L. Thompson   }
61e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
627a982d89SJeremy L. Thompson }
637a982d89SJeremy L. Thompson 
647a982d89SJeremy L. Thompson /**
657a982d89SJeremy L. Thompson   @brief Apply Householder Q matrix
667a982d89SJeremy L. Thompson 
67*ea61e9acSJeremy L Thompson     Compute A = Q A, where Q is mxm and A is mxn.
687a982d89SJeremy L. Thompson 
697a982d89SJeremy L. Thompson   @param[in,out] A      Matrix to apply Householder Q to, in place
70*ea61e9acSJeremy L Thompson   @param[in]     Q      Householder Q matrix
71*ea61e9acSJeremy L Thompson   @param[in]     tau    Householder scaling factors
72*ea61e9acSJeremy L Thompson   @param[in]     t_mode Transpose mode for application
73*ea61e9acSJeremy L Thompson   @param[in]     m      Number of rows in A
74*ea61e9acSJeremy L Thompson   @param[in]     n      Number of columns in A
75*ea61e9acSJeremy L Thompson   @param[in]     k      Number of elementary reflectors in Q, k<m
76*ea61e9acSJeremy L Thompson   @param[in]     row    Row stride in A
77*ea61e9acSJeremy L Thompson   @param[in]     col    Col stride in A
787a982d89SJeremy L. Thompson 
797a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
807a982d89SJeremy L. Thompson 
817a982d89SJeremy L. Thompson   @ref Developer
827a982d89SJeremy L. Thompson **/
832b730f8bSJeremy L Thompson int CeedHouseholderApplyQ(CeedScalar *A, const CeedScalar *Q, const CeedScalar *tau, CeedTransposeMode t_mode, CeedInt m, CeedInt n, CeedInt k,
847a982d89SJeremy L. Thompson                           CeedInt row, CeedInt col) {
8578464608Sjeremylt   CeedScalar *v;
862b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(m, &v));
877a982d89SJeremy L. Thompson   for (CeedInt ii = 0; ii < k; ii++) {
88d1d35e2fSjeremylt     CeedInt i = t_mode == CEED_TRANSPOSE ? ii : k - 1 - ii;
892b730f8bSJeremy L Thompson     for (CeedInt j = i + 1; j < m; j++) v[j] = Q[j * k + i];
90d1d35e2fSjeremylt     // Apply Householder reflector (I - tau v v^T) collo_grad_1d^T
912b730f8bSJeremy L Thompson     CeedCall(CeedHouseholderReflect(&A[i * row], &v[i], tau[i], m - i, n, row, col));
927a982d89SJeremy L. Thompson   }
932b730f8bSJeremy L Thompson   CeedCall(CeedFree(&v));
94e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
957a982d89SJeremy L. Thompson }
967a982d89SJeremy L. Thompson 
977a982d89SJeremy L. Thompson /**
987a982d89SJeremy L. Thompson   @brief Compute Givens rotation
997a982d89SJeremy L. Thompson 
100*ea61e9acSJeremy 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]
1017a982d89SJeremy L. Thompson 
1027a982d89SJeremy L. Thompson   @param[in,out] A      Row major matrix to apply Givens rotation to, in place
103*ea61e9acSJeremy L Thompson   @param[in]     c      Cosine factor
104*ea61e9acSJeremy L Thompson   @param[in]     s      Sine factor
105*ea61e9acSJeremy 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;
1064cc79fe7SJed Brown                           @ref CEED_TRANSPOSE for the opposite rotation
107*ea61e9acSJeremy L Thompson   @param[in]     i      First row/column to apply rotation
108*ea61e9acSJeremy L Thompson   @param[in]     k      Second row/column to apply rotation
109*ea61e9acSJeremy L Thompson   @param[in]     m      Number of rows in A
110*ea61e9acSJeremy L Thompson   @param[in]     n      Number of columns in A
1117a982d89SJeremy L. Thompson 
1127a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
1137a982d89SJeremy L. Thompson 
1147a982d89SJeremy L. Thompson   @ref Developer
1157a982d89SJeremy L. Thompson **/
1162b730f8bSJeremy L Thompson static int CeedGivensRotation(CeedScalar *A, CeedScalar c, CeedScalar s, CeedTransposeMode t_mode, CeedInt i, CeedInt k, CeedInt m, CeedInt n) {
117d1d35e2fSjeremylt   CeedInt stride_j = 1, stride_ik = m, num_its = n;
118d1d35e2fSjeremylt   if (t_mode == CEED_NOTRANSPOSE) {
1192b730f8bSJeremy L Thompson     stride_j  = n;
1202b730f8bSJeremy L Thompson     stride_ik = 1;
1212b730f8bSJeremy L Thompson     num_its   = m;
1227a982d89SJeremy L. Thompson   }
1237a982d89SJeremy L. Thompson 
1247a982d89SJeremy L. Thompson   // Apply rotation
125d1d35e2fSjeremylt   for (CeedInt j = 0; j < num_its; j++) {
126d1d35e2fSjeremylt     CeedScalar tau1 = A[i * stride_ik + j * stride_j], tau2 = A[k * stride_ik + j * stride_j];
127d1d35e2fSjeremylt     A[i * stride_ik + j * stride_j] = c * tau1 - s * tau2;
128d1d35e2fSjeremylt     A[k * stride_ik + j * stride_j] = s * tau1 + c * tau2;
1297a982d89SJeremy L. Thompson   }
130e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1317a982d89SJeremy L. Thompson }
1327a982d89SJeremy L. Thompson 
1337a982d89SJeremy L. Thompson /**
1347a982d89SJeremy L. Thompson   @brief View an array stored in a CeedBasis
1357a982d89SJeremy L. Thompson 
1360a0da059Sjeremylt   @param[in] name   Name of array
137d1d35e2fSjeremylt   @param[in] fp_fmt Printing format
1380a0da059Sjeremylt   @param[in] m      Number of rows in array
1390a0da059Sjeremylt   @param[in] n      Number of columns in array
1400a0da059Sjeremylt   @param[in] a      Array to be viewed
1410a0da059Sjeremylt   @param[in] stream Stream to view to, e.g., stdout
1427a982d89SJeremy L. Thompson 
1437a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
1447a982d89SJeremy L. Thompson 
1457a982d89SJeremy L. Thompson   @ref Developer
1467a982d89SJeremy L. Thompson **/
1472b730f8bSJeremy L Thompson static int CeedScalarView(const char *name, const char *fp_fmt, CeedInt m, CeedInt n, const CeedScalar *a, FILE *stream) {
14892ae7e47SJeremy L Thompson   for (CeedInt i = 0; i < m; i++) {
1492b730f8bSJeremy L Thompson     if (m > 1) fprintf(stream, "%12s[%" CeedInt_FMT "]:", name, i);
1502b730f8bSJeremy L Thompson     else fprintf(stream, "%12s:", name);
1512b730f8bSJeremy 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);
1527a982d89SJeremy L. Thompson     fputs("\n", stream);
1537a982d89SJeremy L. Thompson   }
154e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1557a982d89SJeremy L. Thompson }
1567a982d89SJeremy L. Thompson 
157a76a04e7SJeremy L Thompson /**
158*ea61e9acSJeremy L Thompson   @brief Create the interpolation and gradient matrices for projection from the nodes of `basis_from` to the nodes of `basis_to`.
159*ea61e9acSJeremy L Thompson            The interpolation is given by `interp_project = interp_to^+ * interp_from`, where the pesudoinverse `interp_to^+` is given by QR
160*ea61e9acSJeremy L Thompson factorization.
16114556e63SJeremy L Thompson            The gradient is given by `grad_project = interp_to^+ * grad_from`.
162a76a04e7SJeremy L Thompson            Note: `basis_from` and `basis_to` must have compatible quadrature
163a76a04e7SJeremy L Thompson spaces.
164a76a04e7SJeremy L Thompson 
165a76a04e7SJeremy L Thompson   @param[in]  basis_from     CeedBasis to project from
166a76a04e7SJeremy L Thompson   @param[in]  basis_to       CeedBasis to project to
167*ea61e9acSJeremy L Thompson   @param[out] interp_project Address of the variable where the newly created interpolation matrix will be stored.
168*ea61e9acSJeremy L Thompson   @param[out] grad_project   Address of the variable where the newly created gradient matrix will be stored.
169a76a04e7SJeremy L Thompson 
170a76a04e7SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
171a76a04e7SJeremy L Thompson 
172a76a04e7SJeremy L Thompson   @ref Developer
173a76a04e7SJeremy L Thompson **/
1742b730f8bSJeremy L Thompson static int CeedBasisCreateProjectionMatrices(CeedBasis basis_from, CeedBasis basis_to, CeedScalar **interp_project, CeedScalar **grad_project) {
175a76a04e7SJeremy L Thompson   Ceed ceed;
1762b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetCeed(basis_to, &ceed));
177a76a04e7SJeremy L Thompson 
178a76a04e7SJeremy L Thompson   // Check for compatible quadrature spaces
179a76a04e7SJeremy L Thompson   CeedInt Q_to, Q_from;
1802b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumQuadraturePoints(basis_to, &Q_to));
1812b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumQuadraturePoints(basis_from, &Q_from));
1822b730f8bSJeremy L Thompson   if (Q_to != Q_from) {
183a76a04e7SJeremy L Thompson     // LCOV_EXCL_START
1842b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Bases must have compatible quadrature spaces");
185a76a04e7SJeremy L Thompson     // LCOV_EXCL_STOP
1862b730f8bSJeremy L Thompson   }
187a76a04e7SJeremy L Thompson 
18814556e63SJeremy L Thompson   // Check for matching tensor or non-tensor
189a76a04e7SJeremy L Thompson   CeedInt P_to, P_from, Q = Q_to;
190a76a04e7SJeremy L Thompson   bool    is_tensor_to, is_tensor_from;
1912b730f8bSJeremy L Thompson   CeedCall(CeedBasisIsTensor(basis_to, &is_tensor_to));
1922b730f8bSJeremy L Thompson   CeedCall(CeedBasisIsTensor(basis_from, &is_tensor_from));
193a76a04e7SJeremy L Thompson   if (is_tensor_to && is_tensor_from) {
1942b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_to));
1952b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_from));
1962b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumQuadraturePoints1D(basis_from, &Q));
197a76a04e7SJeremy L Thompson   } else if (!is_tensor_to && !is_tensor_from) {
1982b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis_to, &P_to));
1992b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis_from, &P_from));
200a76a04e7SJeremy L Thompson   } else {
201a76a04e7SJeremy L Thompson     // LCOV_EXCL_START
2022b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_MINOR, "Bases must both be tensor or non-tensor");
203a76a04e7SJeremy L Thompson     // LCOV_EXCL_STOP
204a76a04e7SJeremy L Thompson   }
205a76a04e7SJeremy L Thompson 
20614556e63SJeremy L Thompson   // Get source matrices
20714556e63SJeremy L Thompson   CeedInt     dim;
20814556e63SJeremy L Thompson   CeedScalar *interp_to, *interp_from, *tau;
2092b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetDimension(basis_to, &dim));
2102b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q * P_from, &interp_from));
2112b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q * P_to, &interp_to));
2122b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(P_to * P_from, interp_project));
2132b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(P_to * P_from * (is_tensor_to ? 1 : dim), grad_project));
2142b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q, &tau));
2152b730f8bSJeremy L Thompson   const CeedScalar *interp_to_source = NULL, *interp_from_source = NULL, *grad_from_source;
216a76a04e7SJeremy L Thompson   if (is_tensor_to) {
2172b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetInterp1D(basis_to, &interp_to_source));
2182b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetInterp1D(basis_from, &interp_from_source));
2192b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetGrad1D(basis_from, &grad_from_source));
220a76a04e7SJeremy L Thompson   } else {
2212b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetInterp(basis_to, &interp_to_source));
2222b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetInterp(basis_from, &interp_from_source));
2232b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetGrad(basis_from, &grad_from_source));
224a76a04e7SJeremy L Thompson   }
225a76a04e7SJeremy L Thompson 
22614556e63SJeremy L Thompson   // Build matrices
22714556e63SJeremy L Thompson   CeedInt     num_matrices = 1 + (is_tensor_to ? 1 : dim);
22814556e63SJeremy L Thompson   CeedScalar *input_from[num_matrices], *output_project[num_matrices];
22914556e63SJeremy L Thompson   input_from[0]     = (CeedScalar *)interp_from_source;
23014556e63SJeremy L Thompson   output_project[0] = *interp_project;
23114556e63SJeremy L Thompson   for (CeedInt m = 1; m < num_matrices; m++) {
23214556e63SJeremy L Thompson     input_from[m]     = (CeedScalar *)&grad_from_source[(m - 1) * Q * P_from];
23302af4036SJeremy L Thompson     output_project[m] = &((*grad_project)[(m - 1) * P_to * P_from]);
23414556e63SJeremy L Thompson   }
23514556e63SJeremy L Thompson   for (CeedInt m = 0; m < num_matrices; m++) {
236a76a04e7SJeremy L Thompson     // -- QR Factorization, interp_to = Q R
23714556e63SJeremy L Thompson     memcpy(interp_to, interp_to_source, Q * P_to * sizeof(interp_to_source[0]));
2382b730f8bSJeremy L Thompson     CeedCall(CeedQRFactorization(ceed, interp_to, tau, Q, P_to));
239a76a04e7SJeremy L Thompson 
240a76a04e7SJeremy L Thompson     // -- Apply Qtranspose, interp_to = Qtranspose interp_from
24114556e63SJeremy L Thompson     memcpy(interp_from, input_from[m], Q * P_from * sizeof(input_from[m][0]));
2422b730f8bSJeremy L Thompson     CeedCall(CeedHouseholderApplyQ(interp_from, interp_to, tau, CEED_TRANSPOSE, Q, P_from, P_to, P_from, 1));
243a76a04e7SJeremy L Thompson 
244a76a04e7SJeremy L Thompson     // -- Apply Rinv, interp_project = Rinv interp_c
245a76a04e7SJeremy L Thompson     for (CeedInt j = 0; j < P_from; j++) {  // Column j
2462b730f8bSJeremy 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];
247a76a04e7SJeremy L Thompson       for (CeedInt i = P_to - 2; i >= 0; i--) {  // Row i
24814556e63SJeremy L Thompson         output_project[m][j + P_from * i] = interp_from[j + P_from * i];
249a76a04e7SJeremy L Thompson         for (CeedInt k = i + 1; k < P_to; k++) {
2502b730f8bSJeremy L Thompson           output_project[m][j + P_from * i] -= interp_to[k + P_to * i] * output_project[m][j + P_from * k];
251a76a04e7SJeremy L Thompson         }
25214556e63SJeremy L Thompson         output_project[m][j + P_from * i] /= interp_to[i + P_to * i];
253a76a04e7SJeremy L Thompson       }
254a76a04e7SJeremy L Thompson     }
25514556e63SJeremy L Thompson   }
25614556e63SJeremy L Thompson 
25714556e63SJeremy L Thompson   // Cleanup
2582b730f8bSJeremy L Thompson   CeedCall(CeedFree(&tau));
2592b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_to));
2602b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_from));
261a76a04e7SJeremy L Thompson 
262a76a04e7SJeremy L Thompson   return CEED_ERROR_SUCCESS;
263a76a04e7SJeremy L Thompson }
264a76a04e7SJeremy L Thompson 
2657a982d89SJeremy L. Thompson /// @}
2667a982d89SJeremy L. Thompson 
2677a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
2687a982d89SJeremy L. Thompson /// Ceed Backend API
2697a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
2707a982d89SJeremy L. Thompson /// @addtogroup CeedBasisBackend
2717a982d89SJeremy L. Thompson /// @{
2727a982d89SJeremy L. Thompson 
2737a982d89SJeremy L. Thompson /**
2747a982d89SJeremy L. Thompson   @brief Return collocated grad matrix
2757a982d89SJeremy L. Thompson 
276*ea61e9acSJeremy L Thompson   @param[in]  basis         CeedBasis
277*ea61e9acSJeremy L Thompson   @param[out] collo_grad_1d Row-major (Q_1d * Q_1d) matrix expressing derivatives of basis functions at quadrature points
2787a982d89SJeremy L. Thompson 
2797a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
2807a982d89SJeremy L. Thompson 
2817a982d89SJeremy L. Thompson   @ref Backend
2827a982d89SJeremy L. Thompson **/
283d1d35e2fSjeremylt int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *collo_grad_1d) {
2847a982d89SJeremy L. Thompson   int         i, j, k;
2857a982d89SJeremy L. Thompson   Ceed        ceed;
2862b730f8bSJeremy L Thompson   CeedInt     P_1d = (basis)->P_1d, Q_1d = (basis)->Q_1d;
28778464608Sjeremylt   CeedScalar *interp_1d, *grad_1d, *tau;
2887a982d89SJeremy L. Thompson 
2892b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q_1d * P_1d, &interp_1d));
2902b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q_1d * P_1d, &grad_1d));
2912b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q_1d, &tau));
292d1d35e2fSjeremylt   memcpy(interp_1d, (basis)->interp_1d, Q_1d * P_1d * sizeof(basis)->interp_1d[0]);
293d1d35e2fSjeremylt   memcpy(grad_1d, (basis)->grad_1d, Q_1d * P_1d * sizeof(basis)->interp_1d[0]);
2947a982d89SJeremy L. Thompson 
295d1d35e2fSjeremylt   // QR Factorization, interp_1d = Q R
2962b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetCeed(basis, &ceed));
2972b730f8bSJeremy L Thompson   CeedCall(CeedQRFactorization(ceed, interp_1d, tau, Q_1d, P_1d));
298*ea61e9acSJeremy 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.
2997a982d89SJeremy L. Thompson 
300d1d35e2fSjeremylt   // Apply Rinv, collo_grad_1d = grad_1d Rinv
301d1d35e2fSjeremylt   for (i = 0; i < Q_1d; i++) {  // Row i
302d1d35e2fSjeremylt     collo_grad_1d[Q_1d * i] = grad_1d[P_1d * i] / interp_1d[0];
303d1d35e2fSjeremylt     for (j = 1; j < P_1d; j++) {  // Column j
304d1d35e2fSjeremylt       collo_grad_1d[j + Q_1d * i] = grad_1d[j + P_1d * i];
3052b730f8bSJeremy L Thompson       for (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];
306d1d35e2fSjeremylt       collo_grad_1d[j + Q_1d * i] /= interp_1d[j + P_1d * j];
3077a982d89SJeremy L. Thompson     }
3082b730f8bSJeremy L Thompson     for (j = P_1d; j < Q_1d; j++) collo_grad_1d[j + Q_1d * i] = 0;
3097a982d89SJeremy L. Thompson   }
3107a982d89SJeremy L. Thompson 
311673160d7Sjeremylt   // Apply Qtranspose, collo_grad = collo_grad Q_transpose
3122b730f8bSJeremy L Thompson   CeedCall(CeedHouseholderApplyQ(collo_grad_1d, interp_1d, tau, CEED_NOTRANSPOSE, Q_1d, Q_1d, P_1d, 1, Q_1d));
3137a982d89SJeremy L. Thompson 
3142b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_1d));
3152b730f8bSJeremy L Thompson   CeedCall(CeedFree(&grad_1d));
3162b730f8bSJeremy L Thompson   CeedCall(CeedFree(&tau));
317e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3187a982d89SJeremy L. Thompson }
3197a982d89SJeremy L. Thompson 
3207a982d89SJeremy L. Thompson /**
3217a982d89SJeremy L. Thompson   @brief Get tensor status for given CeedBasis
3227a982d89SJeremy L. Thompson 
323*ea61e9acSJeremy L Thompson   @param[in]  basis     CeedBasis
324d1d35e2fSjeremylt   @param[out] is_tensor Variable to store tensor status
3257a982d89SJeremy L. Thompson 
3267a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3277a982d89SJeremy L. Thompson 
3287a982d89SJeremy L. Thompson   @ref Backend
3297a982d89SJeremy L. Thompson **/
330d1d35e2fSjeremylt int CeedBasisIsTensor(CeedBasis basis, bool *is_tensor) {
331d1d35e2fSjeremylt   *is_tensor = basis->tensor_basis;
332e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3337a982d89SJeremy L. Thompson }
3347a982d89SJeremy L. Thompson 
3357a982d89SJeremy L. Thompson /**
3367a982d89SJeremy L. Thompson   @brief Get backend data of a CeedBasis
3377a982d89SJeremy L. Thompson 
338*ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
3397a982d89SJeremy L. Thompson   @param[out] data  Variable to store data
3407a982d89SJeremy L. Thompson 
3417a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3427a982d89SJeremy L. Thompson 
3437a982d89SJeremy L. Thompson   @ref Backend
3447a982d89SJeremy L. Thompson **/
345777ff853SJeremy L Thompson int CeedBasisGetData(CeedBasis basis, void *data) {
346777ff853SJeremy L Thompson   *(void **)data = basis->data;
347e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3487a982d89SJeremy L. Thompson }
3497a982d89SJeremy L. Thompson 
3507a982d89SJeremy L. Thompson /**
3517a982d89SJeremy L. Thompson   @brief Set backend data of a CeedBasis
3527a982d89SJeremy L. Thompson 
353*ea61e9acSJeremy L Thompson   @param[in,out] basis  CeedBasis
354*ea61e9acSJeremy L Thompson   @param[in]     data   Data to set
3557a982d89SJeremy L. Thompson 
3567a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3577a982d89SJeremy L. Thompson 
3587a982d89SJeremy L. Thompson   @ref Backend
3597a982d89SJeremy L. Thompson **/
360777ff853SJeremy L Thompson int CeedBasisSetData(CeedBasis basis, void *data) {
361777ff853SJeremy L Thompson   basis->data = data;
362e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3637a982d89SJeremy L. Thompson }
3647a982d89SJeremy L. Thompson 
3657a982d89SJeremy L. Thompson /**
36634359f16Sjeremylt   @brief Increment the reference counter for a CeedBasis
36734359f16Sjeremylt 
368*ea61e9acSJeremy L Thompson   @param[in,out] basis Basis to increment the reference counter
36934359f16Sjeremylt 
37034359f16Sjeremylt   @return An error code: 0 - success, otherwise - failure
37134359f16Sjeremylt 
37234359f16Sjeremylt   @ref Backend
37334359f16Sjeremylt **/
3749560d06aSjeremylt int CeedBasisReference(CeedBasis basis) {
37534359f16Sjeremylt   basis->ref_count++;
37634359f16Sjeremylt   return CEED_ERROR_SUCCESS;
37734359f16Sjeremylt }
37834359f16Sjeremylt 
37934359f16Sjeremylt /**
3806e15d496SJeremy L Thompson   @brief Estimate number of FLOPs required to apply CeedBasis in t_mode and eval_mode
3816e15d496SJeremy L Thompson 
382*ea61e9acSJeremy L Thompson   @param[in]  basis     Basis to estimate FLOPs for
383*ea61e9acSJeremy L Thompson   @param[in]  t_mode    Apply basis or transpose
384*ea61e9acSJeremy L Thompson   @param[in]  eval_mode Basis evaluation mode
385*ea61e9acSJeremy L Thompson   @param[out] flops     Address of variable to hold FLOPs estimate
3866e15d496SJeremy L Thompson 
3876e15d496SJeremy L Thompson   @ref Backend
3886e15d496SJeremy L Thompson **/
3892b730f8bSJeremy L Thompson int CeedBasisGetFlopsEstimate(CeedBasis basis, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedSize *flops) {
3906e15d496SJeremy L Thompson   bool is_tensor;
3916e15d496SJeremy L Thompson 
3922b730f8bSJeremy L Thompson   CeedCall(CeedBasisIsTensor(basis, &is_tensor));
3936e15d496SJeremy L Thompson   if (is_tensor) {
3946e15d496SJeremy L Thompson     CeedInt dim, num_comp, P_1d, Q_1d;
3952b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetDimension(basis, &dim));
3962b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
3972b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d));
3982b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
3996e15d496SJeremy L Thompson     if (t_mode == CEED_TRANSPOSE) {
4002b730f8bSJeremy L Thompson       P_1d = Q_1d;
4012b730f8bSJeremy L Thompson       Q_1d = P_1d;
4026e15d496SJeremy L Thompson     }
4036e15d496SJeremy L Thompson     CeedInt tensor_flops = 0, pre = num_comp * CeedIntPow(P_1d, dim - 1), post = 1;
4046e15d496SJeremy L Thompson     for (CeedInt d = 0; d < dim; d++) {
4056e15d496SJeremy L Thompson       tensor_flops += 2 * pre * P_1d * post * Q_1d;
4066e15d496SJeremy L Thompson       pre /= P_1d;
4076e15d496SJeremy L Thompson       post *= Q_1d;
4086e15d496SJeremy L Thompson     }
4096e15d496SJeremy L Thompson     switch (eval_mode) {
4102b730f8bSJeremy L Thompson       case CEED_EVAL_NONE:
4112b730f8bSJeremy L Thompson         *flops = 0;
4122b730f8bSJeremy L Thompson         break;
4132b730f8bSJeremy L Thompson       case CEED_EVAL_INTERP:
4142b730f8bSJeremy L Thompson         *flops = tensor_flops;
4152b730f8bSJeremy L Thompson         break;
4162b730f8bSJeremy L Thompson       case CEED_EVAL_GRAD:
4172b730f8bSJeremy L Thompson         *flops = tensor_flops * 2;
4182b730f8bSJeremy L Thompson         break;
4196e15d496SJeremy L Thompson       case CEED_EVAL_DIV:
4206e15d496SJeremy L Thompson         // LCOV_EXCL_START
4212b730f8bSJeremy L Thompson         return CeedError(basis->ceed, CEED_ERROR_INCOMPATIBLE, "Tensor CEED_EVAL_DIV not supported");
4222b730f8bSJeremy L Thompson         break;
4236e15d496SJeremy L Thompson       case CEED_EVAL_CURL:
4242b730f8bSJeremy L Thompson         return CeedError(basis->ceed, CEED_ERROR_INCOMPATIBLE, "Tensor CEED_EVAL_CURL not supported");
4252b730f8bSJeremy L Thompson         break;
4266e15d496SJeremy L Thompson       // LCOV_EXCL_STOP
4272b730f8bSJeremy L Thompson       case CEED_EVAL_WEIGHT:
4282b730f8bSJeremy L Thompson         *flops = dim * CeedIntPow(Q_1d, dim);
4292b730f8bSJeremy L Thompson         break;
4306e15d496SJeremy L Thompson     }
4316e15d496SJeremy L Thompson   } else {
4326e15d496SJeremy L Thompson     CeedInt dim, num_comp, num_nodes, num_qpts, Q_comp;
4332b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetDimension(basis, &dim));
4342b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
4352b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis, &num_nodes));
4362b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts));
4372b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumQuadratureComponents(basis, &Q_comp));
4386e15d496SJeremy L Thompson     switch (eval_mode) {
4392b730f8bSJeremy L Thompson       case CEED_EVAL_NONE:
4402b730f8bSJeremy L Thompson         *flops = 0;
4412b730f8bSJeremy L Thompson         break;
4422b730f8bSJeremy L Thompson       case CEED_EVAL_INTERP:
4432b730f8bSJeremy L Thompson         *flops = num_nodes * num_qpts * num_comp;
4442b730f8bSJeremy L Thompson         break;
4452b730f8bSJeremy L Thompson       case CEED_EVAL_GRAD:
4462b730f8bSJeremy L Thompson         *flops = num_nodes * num_qpts * num_comp * dim;
4472b730f8bSJeremy L Thompson         break;
4482b730f8bSJeremy L Thompson       case CEED_EVAL_DIV:
4492b730f8bSJeremy L Thompson         *flops = num_nodes * num_qpts;
4502b730f8bSJeremy L Thompson         break;
4512b730f8bSJeremy L Thompson       case CEED_EVAL_CURL:
4522b730f8bSJeremy L Thompson         *flops = num_nodes * num_qpts * dim;
4532b730f8bSJeremy L Thompson         break;
4542b730f8bSJeremy L Thompson       case CEED_EVAL_WEIGHT:
4552b730f8bSJeremy L Thompson         *flops = 0;
4562b730f8bSJeremy L Thompson         break;
4576e15d496SJeremy L Thompson     }
4586e15d496SJeremy L Thompson   }
4596e15d496SJeremy L Thompson 
4606e15d496SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4616e15d496SJeremy L Thompson }
4626e15d496SJeremy L Thompson 
4636e15d496SJeremy L Thompson /**
4647a982d89SJeremy L. Thompson   @brief Get dimension for given CeedElemTopology
4657a982d89SJeremy L. Thompson 
466*ea61e9acSJeremy L Thompson   @param[in]  topo CeedElemTopology
4677a982d89SJeremy L. Thompson   @param[out] dim  Variable to store dimension of topology
4687a982d89SJeremy L. Thompson 
4697a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
4707a982d89SJeremy L. Thompson 
4717a982d89SJeremy L. Thompson   @ref Backend
4727a982d89SJeremy L. Thompson **/
4737a982d89SJeremy L. Thompson int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) {
4747a982d89SJeremy L. Thompson   *dim = (CeedInt)topo >> 16;
475e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4767a982d89SJeremy L. Thompson }
4777a982d89SJeremy L. Thompson 
4787a982d89SJeremy L. Thompson /**
4797a982d89SJeremy L. Thompson   @brief Get CeedTensorContract of a CeedBasis
4807a982d89SJeremy L. Thompson 
481*ea61e9acSJeremy L Thompson   @param[in]  basis     CeedBasis
4827a982d89SJeremy L. Thompson   @param[out] contract  Variable to store CeedTensorContract
4837a982d89SJeremy L. Thompson 
4847a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
4857a982d89SJeremy L. Thompson 
4867a982d89SJeremy L. Thompson   @ref Backend
4877a982d89SJeremy L. Thompson **/
4887a982d89SJeremy L. Thompson int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) {
4897a982d89SJeremy L. Thompson   *contract = basis->contract;
490e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4917a982d89SJeremy L. Thompson }
4927a982d89SJeremy L. Thompson 
4937a982d89SJeremy L. Thompson /**
4947a982d89SJeremy L. Thompson   @brief Set CeedTensorContract of a CeedBasis
4957a982d89SJeremy L. Thompson 
496*ea61e9acSJeremy L Thompson   @param[in,out] basis    CeedBasis
497*ea61e9acSJeremy L Thompson   @param[in]     contract CeedTensorContract to set
4987a982d89SJeremy L. Thompson 
4997a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
5007a982d89SJeremy L. Thompson 
5017a982d89SJeremy L. Thompson   @ref Backend
5027a982d89SJeremy L. Thompson **/
50334359f16Sjeremylt int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract contract) {
50434359f16Sjeremylt   basis->contract = contract;
5052b730f8bSJeremy L Thompson   CeedCall(CeedTensorContractReference(contract));
506e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5077a982d89SJeremy L. Thompson }
5087a982d89SJeremy L. Thompson 
5097a982d89SJeremy L. Thompson /**
5107a982d89SJeremy L. Thompson   @brief Return a reference implementation of matrix multiplication C = A B.
511*ea61e9acSJeremy L Thompson            Note, this is a reference implementation for CPU CeedScalar pointers that is not intended for high performance.
5127a982d89SJeremy L. Thompson 
513*ea61e9acSJeremy L Thompson   @param[in]  ceed  Ceed context for error handling
514d1d35e2fSjeremylt   @param[in]  mat_A Row-major matrix A
515d1d35e2fSjeremylt   @param[in]  mat_B Row-major matrix B
516d1d35e2fSjeremylt   @param[out] mat_C Row-major output matrix C
517*ea61e9acSJeremy L Thompson   @param[in]  m     Number of rows of C
518*ea61e9acSJeremy L Thompson   @param[in]  n     Number of columns of C
519*ea61e9acSJeremy L Thompson   @param[in]  kk    Number of columns of A/rows of B
5207a982d89SJeremy L. Thompson 
5217a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
5227a982d89SJeremy L. Thompson 
5237a982d89SJeremy L. Thompson   @ref Utility
5247a982d89SJeremy L. Thompson **/
5252b730f8bSJeremy L Thompson int CeedMatrixMatrixMultiply(Ceed ceed, const CeedScalar *mat_A, const CeedScalar *mat_B, CeedScalar *mat_C, CeedInt m, CeedInt n, CeedInt kk) {
5262b730f8bSJeremy L Thompson   for (CeedInt i = 0; i < m; i++) {
5277a982d89SJeremy L. Thompson     for (CeedInt j = 0; j < n; j++) {
5287a982d89SJeremy L. Thompson       CeedScalar sum = 0;
5292b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < kk; k++) sum += mat_A[k + i * kk] * mat_B[j + k * n];
530d1d35e2fSjeremylt       mat_C[j + i * n] = sum;
5317a982d89SJeremy L. Thompson     }
5322b730f8bSJeremy L Thompson   }
533e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5347a982d89SJeremy L. Thompson }
5357a982d89SJeremy L. Thompson 
5367a982d89SJeremy L. Thompson /// @}
5377a982d89SJeremy L. Thompson 
5387a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
5397a982d89SJeremy L. Thompson /// CeedBasis Public API
5407a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
5417a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser
542d7b241e6Sjeremylt /// @{
543d7b241e6Sjeremylt 
544b11c1e72Sjeremylt /**
54595bb1877Svaleriabarra   @brief Create a tensor-product basis for H^1 discretizations
546b11c1e72Sjeremylt 
547*ea61e9acSJeremy L Thompson   @param[in]  ceed        Ceed object where the CeedBasis will be created
548*ea61e9acSJeremy L Thompson   @param[in]  dim         Topological dimension
549*ea61e9acSJeremy L Thompson   @param[in]  num_comp    Number of field components (1 for scalar fields)
550*ea61e9acSJeremy L Thompson   @param[in]  P_1d        Number of nodes in one dimension
551*ea61e9acSJeremy L Thompson   @param[in]  Q_1d        Number of quadrature points in one dimension
552*ea61e9acSJeremy L Thompson   @param[in]  interp_1d   Row-major (Q_1d * P_1d) matrix expressing the values of nodal basis functions at quadrature points
553*ea61e9acSJeremy L Thompson   @param[in]  grad_1d     Row-major (Q_1d * P_1d) matrix expressing derivatives of nodal basis functions at quadrature points
554*ea61e9acSJeremy 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]
555*ea61e9acSJeremy L Thompson   @param[in]  q_weight_1d Array of length Q_1d holding the quadrature weights on the reference element
556*ea61e9acSJeremy L Thompson   @param[out] basis       Address of the variable where the newly created CeedBasis will be stored.
557b11c1e72Sjeremylt 
558b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
559dfdf5a53Sjeremylt 
5607a982d89SJeremy L. Thompson   @ref User
561b11c1e72Sjeremylt **/
5622b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d,
5632b730f8bSJeremy L Thompson                             const CeedScalar *grad_1d, const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis *basis) {
5645fe0d4faSjeremylt   if (!ceed->BasisCreateTensorH1) {
5655fe0d4faSjeremylt     Ceed delegate;
5662b730f8bSJeremy L Thompson     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
5675fe0d4faSjeremylt 
5682b730f8bSJeremy L Thompson     if (!delegate) {
569c042f62fSJeremy L Thompson       // LCOV_EXCL_START
5702b730f8bSJeremy L Thompson       return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisCreateTensorH1");
571c042f62fSJeremy L Thompson       // LCOV_EXCL_STOP
5722b730f8bSJeremy L Thompson     }
5735fe0d4faSjeremylt 
5742b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateTensorH1(delegate, dim, num_comp, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis));
575e15f9bd0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
5765fe0d4faSjeremylt   }
577e15f9bd0SJeremy L Thompson 
5782b730f8bSJeremy L Thompson   if (dim < 1) {
579e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
5802b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis dimension must be a positive value");
581e15f9bd0SJeremy L Thompson     // LCOV_EXCL_STOP
5822b730f8bSJeremy L Thompson   }
583227444bfSJeremy L Thompson 
5842b730f8bSJeremy L Thompson   if (num_comp < 1) {
585227444bfSJeremy L Thompson     // LCOV_EXCL_START
5862b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component");
587227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
5882b730f8bSJeremy L Thompson   }
589227444bfSJeremy L Thompson 
5902b730f8bSJeremy L Thompson   if (P_1d < 1) {
591227444bfSJeremy L Thompson     // LCOV_EXCL_START
5922b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node");
593227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
5942b730f8bSJeremy L Thompson   }
595227444bfSJeremy L Thompson 
5962b730f8bSJeremy L Thompson   if (Q_1d < 1) {
597227444bfSJeremy L Thompson     // LCOV_EXCL_START
5982b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point");
599227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
6002b730f8bSJeremy L Thompson   }
601227444bfSJeremy L Thompson 
6022b730f8bSJeremy L Thompson   CeedElemTopology topo = dim == 1 ? CEED_TOPOLOGY_LINE : dim == 2 ? CEED_TOPOLOGY_QUAD : CEED_TOPOLOGY_HEX;
603e15f9bd0SJeremy L Thompson 
6042b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(1, basis));
605d7b241e6Sjeremylt   (*basis)->ceed = ceed;
6062b730f8bSJeremy L Thompson   CeedCall(CeedReference(ceed));
607d1d35e2fSjeremylt   (*basis)->ref_count    = 1;
608d1d35e2fSjeremylt   (*basis)->tensor_basis = 1;
609d7b241e6Sjeremylt   (*basis)->dim          = dim;
610d99fa3c5SJeremy L Thompson   (*basis)->topo         = topo;
611d1d35e2fSjeremylt   (*basis)->num_comp     = num_comp;
612d1d35e2fSjeremylt   (*basis)->P_1d         = P_1d;
613d1d35e2fSjeremylt   (*basis)->Q_1d         = Q_1d;
614d1d35e2fSjeremylt   (*basis)->P            = CeedIntPow(P_1d, dim);
615d1d35e2fSjeremylt   (*basis)->Q            = CeedIntPow(Q_1d, dim);
61650c301a5SRezgar Shakeri   (*basis)->Q_comp       = 1;
61750c301a5SRezgar Shakeri   (*basis)->basis_space  = 1;  // 1 for H^1 space
6182b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d, &(*basis)->q_ref_1d));
6192b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d, &(*basis)->q_weight_1d));
620ff3a0f91SJeremy L Thompson   if (q_ref_1d) memcpy((*basis)->q_ref_1d, q_ref_1d, Q_1d * sizeof(q_ref_1d[0]));
6212b730f8bSJeremy L Thompson   if (q_weight_1d) memcpy((*basis)->q_weight_1d, q_weight_1d, Q_1d * sizeof(q_weight_1d[0]));
6222b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->interp_1d));
6232b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->grad_1d));
6242b730f8bSJeremy L Thompson   if (interp_1d) memcpy((*basis)->interp_1d, interp_1d, Q_1d * P_1d * sizeof(interp_1d[0]));
625ff3a0f91SJeremy L Thompson   if (grad_1d) memcpy((*basis)->grad_1d, grad_1d, Q_1d * P_1d * sizeof(grad_1d[0]));
6262b730f8bSJeremy L Thompson   CeedCall(ceed->BasisCreateTensorH1(dim, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, *basis));
627e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
628d7b241e6Sjeremylt }
629d7b241e6Sjeremylt 
630b11c1e72Sjeremylt /**
63195bb1877Svaleriabarra   @brief Create a tensor-product Lagrange basis
632b11c1e72Sjeremylt 
633*ea61e9acSJeremy L Thompson   @param[in]  ceed      Ceed object where the CeedBasis will be created
634*ea61e9acSJeremy L Thompson   @param[in]  dim       Topological dimension of element
635*ea61e9acSJeremy L Thompson   @param[in]  num_comp  Number of field components (1 for scalar fields)
636*ea61e9acSJeremy L Thompson   @param[in]  P         Number of Gauss-Lobatto nodes in one dimension.
637*ea61e9acSJeremy L Thompson                           The polynomial degree of the resulting Q_k element is k=P-1.
638*ea61e9acSJeremy L Thompson   @param[in]  Q         Number of quadrature points in one dimension.
639*ea61e9acSJeremy L Thompson   @param[in]  quad_mode Distribution of the Q quadrature points (affects order of accuracy for the quadrature)
640*ea61e9acSJeremy L Thompson   @param[out] basis     Address of the variable where the newly created CeedBasis will be stored.
641b11c1e72Sjeremylt 
642b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
643dfdf5a53Sjeremylt 
6447a982d89SJeremy L. Thompson   @ref User
645b11c1e72Sjeremylt **/
6462b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P, CeedInt Q, CeedQuadMode quad_mode, CeedBasis *basis) {
647d7b241e6Sjeremylt   // Allocate
6482b730f8bSJeremy L Thompson   int        ierr = CEED_ERROR_SUCCESS, i, j, k;
6492b730f8bSJeremy L Thompson   CeedScalar c1, c2, c3, c4, dx, *nodes, *interp_1d, *grad_1d, *q_ref_1d, *q_weight_1d;
6504d537eeaSYohann 
6512b730f8bSJeremy L Thompson   if (dim < 1) {
652c042f62fSJeremy L Thompson     // LCOV_EXCL_START
6532b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis dimension must be a positive value");
654c042f62fSJeremy L Thompson     // LCOV_EXCL_STOP
6552b730f8bSJeremy L Thompson   }
6564d537eeaSYohann 
6572b730f8bSJeremy L Thompson   if (num_comp < 1) {
658227444bfSJeremy L Thompson     // LCOV_EXCL_START
6592b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component");
660227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
6612b730f8bSJeremy L Thompson   }
662227444bfSJeremy L Thompson 
6632b730f8bSJeremy L Thompson   if (P < 1) {
664227444bfSJeremy L Thompson     // LCOV_EXCL_START
6652b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node");
666227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
6672b730f8bSJeremy L Thompson   }
668227444bfSJeremy L Thompson 
6692b730f8bSJeremy L Thompson   if (Q < 1) {
670227444bfSJeremy L Thompson     // LCOV_EXCL_START
6712b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point");
672227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
6732b730f8bSJeremy L Thompson   }
674227444bfSJeremy L Thompson 
675e15f9bd0SJeremy L Thompson   // Get Nodes and Weights
6762b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(P * Q, &interp_1d));
6772b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(P * Q, &grad_1d));
6782b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(P, &nodes));
6792b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q, &q_ref_1d));
6802b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q, &q_weight_1d));
6812b730f8bSJeremy L Thompson   if (CeedLobattoQuadrature(P, nodes, NULL) != CEED_ERROR_SUCCESS) goto cleanup;
682d1d35e2fSjeremylt   switch (quad_mode) {
683d7b241e6Sjeremylt     case CEED_GAUSS:
684d1d35e2fSjeremylt       ierr = CeedGaussQuadrature(Q, q_ref_1d, q_weight_1d);
685d7b241e6Sjeremylt       break;
686d7b241e6Sjeremylt     case CEED_GAUSS_LOBATTO:
687d1d35e2fSjeremylt       ierr = CeedLobattoQuadrature(Q, q_ref_1d, q_weight_1d);
688d7b241e6Sjeremylt       break;
689d7b241e6Sjeremylt   }
6902b730f8bSJeremy L Thompson   if (ierr != CEED_ERROR_SUCCESS) goto cleanup;
691e15f9bd0SJeremy L Thompson 
692d7b241e6Sjeremylt   // Build B, D matrix
693d7b241e6Sjeremylt   // Fornberg, 1998
694d7b241e6Sjeremylt   for (i = 0; i < Q; i++) {
695d7b241e6Sjeremylt     c1                   = 1.0;
696d1d35e2fSjeremylt     c3                   = nodes[0] - q_ref_1d[i];
697d1d35e2fSjeremylt     interp_1d[i * P + 0] = 1.0;
698d7b241e6Sjeremylt     for (j = 1; j < P; j++) {
699d7b241e6Sjeremylt       c2 = 1.0;
700d7b241e6Sjeremylt       c4 = c3;
701d1d35e2fSjeremylt       c3 = nodes[j] - q_ref_1d[i];
702d7b241e6Sjeremylt       for (k = 0; k < j; k++) {
703d7b241e6Sjeremylt         dx = nodes[j] - nodes[k];
704d7b241e6Sjeremylt         c2 *= dx;
705d7b241e6Sjeremylt         if (k == j - 1) {
706d1d35e2fSjeremylt           grad_1d[i * P + j]   = c1 * (interp_1d[i * P + k] - c4 * grad_1d[i * P + k]) / c2;
707d1d35e2fSjeremylt           interp_1d[i * P + j] = -c1 * c4 * interp_1d[i * P + k] / c2;
708d7b241e6Sjeremylt         }
709d1d35e2fSjeremylt         grad_1d[i * P + k]   = (c3 * grad_1d[i * P + k] - interp_1d[i * P + k]) / dx;
710d1d35e2fSjeremylt         interp_1d[i * P + k] = c3 * interp_1d[i * P + k] / dx;
711d7b241e6Sjeremylt       }
712d7b241e6Sjeremylt       c1 = c2;
713d7b241e6Sjeremylt     }
714d7b241e6Sjeremylt   }
7159ac7b42eSJeremy L Thompson   // Pass to CeedBasisCreateTensorH1
7162b730f8bSJeremy L Thompson   CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P, Q, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis));
717e15f9bd0SJeremy L Thompson cleanup:
7182b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_1d));
7192b730f8bSJeremy L Thompson   CeedCall(CeedFree(&grad_1d));
7202b730f8bSJeremy L Thompson   CeedCall(CeedFree(&nodes));
7212b730f8bSJeremy L Thompson   CeedCall(CeedFree(&q_ref_1d));
7222b730f8bSJeremy L Thompson   CeedCall(CeedFree(&q_weight_1d));
723e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
724d7b241e6Sjeremylt }
725d7b241e6Sjeremylt 
726b11c1e72Sjeremylt /**
72795bb1877Svaleriabarra   @brief Create a non tensor-product basis for H^1 discretizations
728a8de75f0Sjeremylt 
729*ea61e9acSJeremy L Thompson   @param[in]  ceed      Ceed object where the CeedBasis will be created
730*ea61e9acSJeremy L Thompson   @param[in]  topo      Topology of element, e.g. hypercube, simplex, ect
731*ea61e9acSJeremy L Thompson   @param[in]  num_comp  Number of field components (1 for scalar fields)
732*ea61e9acSJeremy L Thompson   @param[in]  num_nodes Total number of nodes
733*ea61e9acSJeremy L Thompson   @param[in]  num_qpts  Total number of quadrature points
734*ea61e9acSJeremy L Thompson   @param[in]  interp    Row-major (num_qpts * num_nodes) matrix expressing the values of nodal basis functions at quadrature points
735*ea61e9acSJeremy L Thompson   @param[in]  grad      Row-major (num_qpts * dim * num_nodes) matrix expressing derivatives of nodal basis functions at quadrature points
736*ea61e9acSJeremy L Thompson   @param[in]  q_ref     Array of length num_qpts holding the locations of quadrature points on the reference element
737*ea61e9acSJeremy L Thompson   @param[in]  q_weight  Array of length num_qpts holding the quadrature weights on the reference element
738*ea61e9acSJeremy L Thompson   @param[out] basis     Address of the variable where the newly created CeedBasis will be stored.
739a8de75f0Sjeremylt 
740a8de75f0Sjeremylt   @return An error code: 0 - success, otherwise - failure
741a8de75f0Sjeremylt 
7427a982d89SJeremy L. Thompson   @ref User
743a8de75f0Sjeremylt **/
7442b730f8bSJeremy L Thompson int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
7452b730f8bSJeremy L Thompson                       const CeedScalar *grad, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) {
746d1d35e2fSjeremylt   CeedInt P = num_nodes, Q = num_qpts, dim = 0;
747a8de75f0Sjeremylt 
7485fe0d4faSjeremylt   if (!ceed->BasisCreateH1) {
7495fe0d4faSjeremylt     Ceed delegate;
7502b730f8bSJeremy L Thompson     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
7515fe0d4faSjeremylt 
7522b730f8bSJeremy L Thompson     if (!delegate) {
753c042f62fSJeremy L Thompson       // LCOV_EXCL_START
7542b730f8bSJeremy L Thompson       return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisCreateH1");
755c042f62fSJeremy L Thompson       // LCOV_EXCL_STOP
7562b730f8bSJeremy L Thompson     }
7575fe0d4faSjeremylt 
7582b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateH1(delegate, topo, num_comp, num_nodes, num_qpts, interp, grad, q_ref, q_weight, basis));
759e15f9bd0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
7605fe0d4faSjeremylt   }
7615fe0d4faSjeremylt 
7622b730f8bSJeremy L Thompson   if (num_comp < 1) {
763227444bfSJeremy L Thompson     // LCOV_EXCL_START
7642b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component");
765227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
7662b730f8bSJeremy L Thompson   }
767227444bfSJeremy L Thompson 
7682b730f8bSJeremy L Thompson   if (num_nodes < 1) {
769227444bfSJeremy L Thompson     // LCOV_EXCL_START
7702b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node");
771227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
7722b730f8bSJeremy L Thompson   }
773227444bfSJeremy L Thompson 
7742b730f8bSJeremy L Thompson   if (num_qpts < 1) {
775227444bfSJeremy L Thompson     // LCOV_EXCL_START
7762b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point");
777227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
7782b730f8bSJeremy L Thompson   }
779227444bfSJeremy L Thompson 
7802b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(1, basis));
781a8de75f0Sjeremylt 
7822b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetTopologyDimension(topo, &dim));
783a8de75f0Sjeremylt 
784a8de75f0Sjeremylt   (*basis)->ceed = ceed;
7852b730f8bSJeremy L Thompson   CeedCall(CeedReference(ceed));
786d1d35e2fSjeremylt   (*basis)->ref_count    = 1;
787d1d35e2fSjeremylt   (*basis)->tensor_basis = 0;
788a8de75f0Sjeremylt   (*basis)->dim          = dim;
789d99fa3c5SJeremy L Thompson   (*basis)->topo         = topo;
790d1d35e2fSjeremylt   (*basis)->num_comp     = num_comp;
791a8de75f0Sjeremylt   (*basis)->P            = P;
792a8de75f0Sjeremylt   (*basis)->Q            = Q;
79350c301a5SRezgar Shakeri   (*basis)->Q_comp       = 1;
79450c301a5SRezgar Shakeri   (*basis)->basis_space  = 1;  // 1 for H^1 space
7952b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q * dim, &(*basis)->q_ref_1d));
7962b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q, &(*basis)->q_weight_1d));
797ff3a0f91SJeremy L Thompson   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0]));
798ff3a0f91SJeremy L Thompson   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0]));
7992b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q * P, &(*basis)->interp));
8002b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(dim * Q * P, &(*basis)->grad));
801ff3a0f91SJeremy L Thompson   if (interp) memcpy((*basis)->interp, interp, Q * P * sizeof(interp[0]));
802ff3a0f91SJeremy L Thompson   if (grad) memcpy((*basis)->grad, grad, dim * Q * P * sizeof(grad[0]));
8032b730f8bSJeremy L Thompson   CeedCall(ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, q_ref, q_weight, *basis));
804e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
805a8de75f0Sjeremylt }
806a8de75f0Sjeremylt 
807a8de75f0Sjeremylt /**
80850c301a5SRezgar Shakeri   @brief Create a non tensor-product basis for H(div) discretizations
80950c301a5SRezgar Shakeri 
810*ea61e9acSJeremy L Thompson   @param[in]  ceed      Ceed object where the CeedBasis will be created
811*ea61e9acSJeremy 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
812*ea61e9acSJeremy L Thompson   @param[in]  num_comp  Number of components (usually 1 for vectors in H(div) bases)
813*ea61e9acSJeremy L Thompson   @param[in]  num_nodes Total number of nodes (dofs per element)
814*ea61e9acSJeremy L Thompson   @param[in]  num_qpts  Total number of quadrature points
815*ea61e9acSJeremy L Thompson   @param[in]  interp    Row-major (dim*num_qpts * num_nodes) matrix expressing the values of nodal basis functions at quadrature points
816*ea61e9acSJeremy L Thompson   @param[in]  div       Row-major (num_qpts * num_nodes) matrix expressing divergence of nodal basis functions at quadrature points
817*ea61e9acSJeremy L Thompson   @param[in]  q_ref     Array of length num_qpts holding the locations of quadrature points on the reference element
818*ea61e9acSJeremy L Thompson   @param[in]  q_weight  Array of length num_qpts holding the quadrature weights on the reference element
819*ea61e9acSJeremy L Thompson   @param[out] basis     Address of the variable where the newly created CeedBasis will be stored.
82050c301a5SRezgar Shakeri 
82150c301a5SRezgar Shakeri   @return An error code: 0 - success, otherwise - failure
82250c301a5SRezgar Shakeri 
82350c301a5SRezgar Shakeri   @ref User
82450c301a5SRezgar Shakeri **/
8252b730f8bSJeremy L Thompson int CeedBasisCreateHdiv(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
8262b730f8bSJeremy L Thompson                         const CeedScalar *div, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) {
82750c301a5SRezgar Shakeri   CeedInt Q = num_qpts, P = num_nodes, dim = 0;
8282b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetTopologyDimension(topo, &dim));
82950c301a5SRezgar Shakeri   if (!ceed->BasisCreateHdiv) {
83050c301a5SRezgar Shakeri     Ceed delegate;
8312b730f8bSJeremy L Thompson     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
83250c301a5SRezgar Shakeri 
8332b730f8bSJeremy L Thompson     if (!delegate) {
83450c301a5SRezgar Shakeri       // LCOV_EXCL_START
8352b730f8bSJeremy L Thompson       return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateHdiv");
83650c301a5SRezgar Shakeri       // LCOV_EXCL_STOP
8372b730f8bSJeremy L Thompson     }
83850c301a5SRezgar Shakeri 
8392b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateHdiv(delegate, topo, num_comp, num_nodes, num_qpts, interp, div, q_ref, q_weight, basis));
84050c301a5SRezgar Shakeri     return CEED_ERROR_SUCCESS;
84150c301a5SRezgar Shakeri   }
84250c301a5SRezgar Shakeri 
8432b730f8bSJeremy L Thompson   if (num_comp < 1) {
844227444bfSJeremy L Thompson     // LCOV_EXCL_START
8452b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component");
846227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
8472b730f8bSJeremy L Thompson   }
848227444bfSJeremy L Thompson 
8492b730f8bSJeremy L Thompson   if (num_nodes < 1) {
850227444bfSJeremy L Thompson     // LCOV_EXCL_START
8512b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node");
852227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
8532b730f8bSJeremy L Thompson   }
854227444bfSJeremy L Thompson 
8552b730f8bSJeremy L Thompson   if (num_qpts < 1) {
856227444bfSJeremy L Thompson     // LCOV_EXCL_START
8572b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point");
858227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
8592b730f8bSJeremy L Thompson   }
860227444bfSJeremy L Thompson 
8612b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(1, basis));
86250c301a5SRezgar Shakeri 
86350c301a5SRezgar Shakeri   (*basis)->ceed = ceed;
8642b730f8bSJeremy L Thompson   CeedCall(CeedReference(ceed));
86550c301a5SRezgar Shakeri   (*basis)->ref_count    = 1;
86650c301a5SRezgar Shakeri   (*basis)->tensor_basis = 0;
86750c301a5SRezgar Shakeri   (*basis)->dim          = dim;
86850c301a5SRezgar Shakeri   (*basis)->topo         = topo;
86950c301a5SRezgar Shakeri   (*basis)->num_comp     = num_comp;
87050c301a5SRezgar Shakeri   (*basis)->P            = P;
87150c301a5SRezgar Shakeri   (*basis)->Q            = Q;
87250c301a5SRezgar Shakeri   (*basis)->Q_comp       = dim;
87350c301a5SRezgar Shakeri   (*basis)->basis_space  = 2;  // 2 for H(div) space
8742b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q * dim, &(*basis)->q_ref_1d));
8752b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q, &(*basis)->q_weight_1d));
87650c301a5SRezgar Shakeri   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0]));
87750c301a5SRezgar Shakeri   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0]));
8782b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(dim * Q * P, &(*basis)->interp));
8792b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q * P, &(*basis)->div));
88050c301a5SRezgar Shakeri   if (interp) memcpy((*basis)->interp, interp, dim * Q * P * sizeof(interp[0]));
88150c301a5SRezgar Shakeri   if (div) memcpy((*basis)->div, div, Q * P * sizeof(div[0]));
8822b730f8bSJeremy L Thompson   CeedCall(ceed->BasisCreateHdiv(topo, dim, P, Q, interp, div, q_ref, q_weight, *basis));
88350c301a5SRezgar Shakeri   return CEED_ERROR_SUCCESS;
88450c301a5SRezgar Shakeri }
88550c301a5SRezgar Shakeri 
88650c301a5SRezgar Shakeri /**
887*ea61e9acSJeremy L Thompson   @brief Create a CeedBasis for projection from the nodes of `basis_from` to the nodes of `basis_to`.
888*ea61e9acSJeremy L Thompson            Only `CEED_EVAL_INTERP` and `CEED_EVAL_GRAD` will be valid for the new basis, `basis_project`.
889*ea61e9acSJeremy L Thompson            The interpolation is given by `interp_project = interp_to^+ * interp_from`, where the pesudoinverse `interp_to^+` is given by QR
890*ea61e9acSJeremy L Thompson factorization. The gradient is given by `grad_project = interp_to^+ * grad_from`. Note: `basis_from` and `basis_to` must have compatible quadrature
891*ea61e9acSJeremy L Thompson spaces. Note: `basis_project` will have the same number of components as `basis_from`, regardless of the number of components that `basis_to` has. If
892*ea61e9acSJeremy L Thompson `basis_from` has 3 components and `basis_to` has 5 components, then `basis_project` will have 3 components.
893f113e5dcSJeremy L Thompson 
894f113e5dcSJeremy L Thompson   @param[in]  basis_from    CeedBasis to prolong from
895446e7af4SJeremy L Thompson   @param[in]  basis_to      CeedBasis to prolong to
896*ea61e9acSJeremy L Thompson   @param[out] basis_project Address of the variable where the newly created CeedBasis will be stored.
897f113e5dcSJeremy L Thompson 
898f113e5dcSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
899f113e5dcSJeremy L Thompson 
900f113e5dcSJeremy L Thompson   @ref User
901f113e5dcSJeremy L Thompson **/
9022b730f8bSJeremy L Thompson int CeedBasisCreateProjection(CeedBasis basis_from, CeedBasis basis_to, CeedBasis *basis_project) {
903f113e5dcSJeremy L Thompson   Ceed ceed;
9042b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetCeed(basis_to, &ceed));
905f113e5dcSJeremy L Thompson 
906f113e5dcSJeremy L Thompson   // Create projectior matrix
90714556e63SJeremy L Thompson   CeedScalar *interp_project, *grad_project;
9082b730f8bSJeremy L Thompson   CeedCall(CeedBasisCreateProjectionMatrices(basis_from, basis_to, &interp_project, &grad_project));
909f113e5dcSJeremy L Thompson 
910f113e5dcSJeremy L Thompson   // Build basis
911f113e5dcSJeremy L Thompson   bool        is_tensor;
912f113e5dcSJeremy L Thompson   CeedInt     dim, num_comp;
91314556e63SJeremy L Thompson   CeedScalar *q_ref, *q_weight;
9142b730f8bSJeremy L Thompson   CeedCall(CeedBasisIsTensor(basis_to, &is_tensor));
9152b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetDimension(basis_to, &dim));
9162b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumComponents(basis_from, &num_comp));
917f113e5dcSJeremy L Thompson   if (is_tensor) {
918f113e5dcSJeremy L Thompson     CeedInt P_1d_to, P_1d_from;
9192b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_1d_from));
9202b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_1d_to));
9212b730f8bSJeremy L Thompson     CeedCall(CeedCalloc(P_1d_to, &q_ref));
9222b730f8bSJeremy L Thompson     CeedCall(CeedCalloc(P_1d_to, &q_weight));
9232b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P_1d_from, P_1d_to, interp_project, grad_project, q_ref, q_weight, basis_project));
924f113e5dcSJeremy L Thompson   } else {
925f113e5dcSJeremy L Thompson     CeedElemTopology topo;
9262b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetTopology(basis_to, &topo));
927f113e5dcSJeremy L Thompson     CeedInt num_nodes_to, num_nodes_from;
9282b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis_from, &num_nodes_from));
9292b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis_to, &num_nodes_to));
9302b730f8bSJeremy L Thompson     CeedCall(CeedCalloc(num_nodes_to * dim, &q_ref));
9312b730f8bSJeremy L Thompson     CeedCall(CeedCalloc(num_nodes_to, &q_weight));
9322b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateH1(ceed, topo, num_comp, num_nodes_from, num_nodes_to, interp_project, grad_project, q_ref, q_weight, basis_project));
933f113e5dcSJeremy L Thompson   }
934f113e5dcSJeremy L Thompson 
935f113e5dcSJeremy L Thompson   // Cleanup
9362b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_project));
9372b730f8bSJeremy L Thompson   CeedCall(CeedFree(&grad_project));
9382b730f8bSJeremy L Thompson   CeedCall(CeedFree(&q_ref));
9392b730f8bSJeremy L Thompson   CeedCall(CeedFree(&q_weight));
940f113e5dcSJeremy L Thompson 
941f113e5dcSJeremy L Thompson   return CEED_ERROR_SUCCESS;
942f113e5dcSJeremy L Thompson }
943f113e5dcSJeremy L Thompson 
944f113e5dcSJeremy L Thompson /**
945*ea61e9acSJeremy L Thompson   @brief Copy the pointer to a CeedBasis.
946*ea61e9acSJeremy L Thompson            Both pointers should be destroyed with `CeedBasisDestroy()`.
9479560d06aSjeremylt 
948*ea61e9acSJeremy L Thompson            Note: If `*basis_copy` is non-NULL, then it is assumed that `*basis_copy` is a pointer to a CeedBasis.
949*ea61e9acSJeremy L Thompson              This CeedBasis will be destroyed if `*basis_copy` is the only reference to this CeedBasis.
950*ea61e9acSJeremy L Thompson 
951*ea61e9acSJeremy L Thompson   @param[in]     basis      CeedBasis to copy reference to
952*ea61e9acSJeremy L Thompson   @param[in,out] basis_copy Variable to store copied reference
9539560d06aSjeremylt 
9549560d06aSjeremylt   @return An error code: 0 - success, otherwise - failure
9559560d06aSjeremylt 
9569560d06aSjeremylt   @ref User
9579560d06aSjeremylt **/
9589560d06aSjeremylt int CeedBasisReferenceCopy(CeedBasis basis, CeedBasis *basis_copy) {
9592b730f8bSJeremy L Thompson   CeedCall(CeedBasisReference(basis));
9602b730f8bSJeremy L Thompson   CeedCall(CeedBasisDestroy(basis_copy));
9619560d06aSjeremylt   *basis_copy = basis;
9629560d06aSjeremylt   return CEED_ERROR_SUCCESS;
9639560d06aSjeremylt }
9649560d06aSjeremylt 
9659560d06aSjeremylt /**
9667a982d89SJeremy L. Thompson   @brief View a CeedBasis
9677a982d89SJeremy L. Thompson 
968*ea61e9acSJeremy L Thompson   @param[in] basis  CeedBasis to view
969*ea61e9acSJeremy L Thompson   @param[in] stream Stream to view to, e.g., stdout
9707a982d89SJeremy L. Thompson 
9717a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
9727a982d89SJeremy L. Thompson 
9737a982d89SJeremy L. Thompson   @ref User
9747a982d89SJeremy L. Thompson **/
9757a982d89SJeremy L. Thompson int CeedBasisView(CeedBasis basis, FILE *stream) {
97650c301a5SRezgar Shakeri   CeedFESpace      FE_space = basis->basis_space;
97750c301a5SRezgar Shakeri   CeedElemTopology topo     = basis->topo;
9782b730f8bSJeremy L Thompson 
97950c301a5SRezgar Shakeri   // Print FE space and element topology of the basis
980d1d35e2fSjeremylt   if (basis->tensor_basis) {
9812b730f8bSJeremy L Thompson     fprintf(stream, "CeedBasis (%s on a %s element): dim=%" CeedInt_FMT " P=%" CeedInt_FMT " Q=%" CeedInt_FMT "\n", CeedFESpaces[FE_space],
9822b730f8bSJeremy L Thompson             CeedElemTopologies[topo], basis->dim, basis->P_1d, basis->Q_1d);
98350c301a5SRezgar Shakeri   } else {
9842b730f8bSJeremy L Thompson     fprintf(stream, "CeedBasis (%s on a %s element): dim=%" CeedInt_FMT " P=%" CeedInt_FMT " Q=%" CeedInt_FMT "\n", CeedFESpaces[FE_space],
9852b730f8bSJeremy L Thompson             CeedElemTopologies[topo], basis->dim, basis->P, basis->Q);
98650c301a5SRezgar Shakeri   }
987*ea61e9acSJeremy L Thompson   // Print quadrature data, interpolation/gradient/divergence/curl of the basis
98850c301a5SRezgar Shakeri   if (basis->tensor_basis) {  // tensor basis
9892b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_ref_1d, stream));
9902b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_weight_1d, stream));
9912b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("interp1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, basis->interp_1d, stream));
9922b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("grad1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, basis->grad_1d, stream));
99350c301a5SRezgar Shakeri   } else {  // non-tensor basis
9942b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("qref", "\t% 12.8f", 1, basis->Q * basis->dim, basis->q_ref_1d, stream));
9952b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->q_weight_1d, stream));
9962b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("interp", "\t% 12.8f", basis->Q_comp * basis->Q, basis->P, basis->interp, stream));
99750c301a5SRezgar Shakeri     if (basis->grad) {
9982b730f8bSJeremy L Thompson       CeedCall(CeedScalarView("grad", "\t% 12.8f", basis->dim * basis->Q, basis->P, basis->grad, stream));
9997a982d89SJeremy L. Thompson     }
100050c301a5SRezgar Shakeri     if (basis->div) {
10012b730f8bSJeremy L Thompson       CeedCall(CeedScalarView("div", "\t% 12.8f", basis->Q, basis->P, basis->div, stream));
100250c301a5SRezgar Shakeri     }
100350c301a5SRezgar Shakeri   }
1004e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
10057a982d89SJeremy L. Thompson }
10067a982d89SJeremy L. Thompson 
10077a982d89SJeremy L. Thompson /**
10087a982d89SJeremy L. Thompson   @brief Apply basis evaluation from nodes to quadrature points or vice versa
10097a982d89SJeremy L. Thompson 
1010*ea61e9acSJeremy L Thompson   @param[in]  basis      CeedBasis to evaluate
1011*ea61e9acSJeremy L Thompson   @param[in]  num_elem   The number of elements to apply the basis evaluation to;
1012*ea61e9acSJeremy L Thompson                            the backend will specify the ordering in CeedElemRestrictionCreateBlocked()
1013*ea61e9acSJeremy L Thompson   @param[in]  t_mode    \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature points;
1014*ea61e9acSJeremy L Thompson                           \ref CEED_TRANSPOSE to apply the transpose, mapping from quadrature points to nodes
1015*ea61e9acSJeremy L Thompson   @param[in]  eval_mode \ref CEED_EVAL_NONE to use values directly,
10167a982d89SJeremy L. Thompson                           \ref CEED_EVAL_INTERP to use interpolated values,
10177a982d89SJeremy L. Thompson                           \ref CEED_EVAL_GRAD to use gradients,
10187a982d89SJeremy L. Thompson                           \ref CEED_EVAL_WEIGHT to use quadrature weights.
10197a982d89SJeremy L. Thompson   @param[in]  u        Input CeedVector
10207a982d89SJeremy L. Thompson   @param[out] v        Output CeedVector
10217a982d89SJeremy L. Thompson 
10227a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
10237a982d89SJeremy L. Thompson 
10247a982d89SJeremy L. Thompson   @ref User
10257a982d89SJeremy L. Thompson **/
10262b730f8bSJeremy L Thompson int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
10271f9221feSJeremy L Thompson   CeedSize u_length = 0, v_length;
10281f9221feSJeremy L Thompson   CeedInt  dim, num_comp, num_nodes, num_qpts;
10292b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetDimension(basis, &dim));
10302b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
10312b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumNodes(basis, &num_nodes));
10322b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts));
10332b730f8bSJeremy L Thompson   CeedCall(CeedVectorGetLength(v, &v_length));
10347a982d89SJeremy L. Thompson   if (u) {
10352b730f8bSJeremy L Thompson     CeedCall(CeedVectorGetLength(u, &u_length));
10367a982d89SJeremy L. Thompson   }
10377a982d89SJeremy L. Thompson 
10382b730f8bSJeremy L Thompson   if (!basis->Apply) {
1039e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
10402b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisApply");
1041e15f9bd0SJeremy L Thompson     // LCOV_EXCL_STOP
10422b730f8bSJeremy L Thompson   }
1043e15f9bd0SJeremy L Thompson 
1044e15f9bd0SJeremy L Thompson   // Check compatibility of topological and geometrical dimensions
10452b730f8bSJeremy L Thompson   if ((t_mode == CEED_TRANSPOSE && (v_length % num_nodes != 0 || u_length % num_qpts != 0)) ||
10462b730f8bSJeremy L Thompson       (t_mode == CEED_NOTRANSPOSE && (u_length % num_nodes != 0 || v_length % num_qpts != 0))) {
10478229195eSjeremylt     // LCOV_EXCL_START
10482b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_DIMENSION, "Length of input/output vectors incompatible with basis dimensions");
10498229195eSjeremylt     // LCOV_EXCL_STOP
10502b730f8bSJeremy L Thompson   }
10517a982d89SJeremy L. Thompson 
1052e15f9bd0SJeremy L Thompson   // Check vector lengths to prevent out of bounds issues
1053d1d35e2fSjeremylt   bool bad_dims = false;
1054d1d35e2fSjeremylt   switch (eval_mode) {
1055e15f9bd0SJeremy L Thompson     case CEED_EVAL_NONE:
10562b730f8bSJeremy L Thompson     case CEED_EVAL_INTERP:
10572b730f8bSJeremy L Thompson       bad_dims = ((t_mode == CEED_TRANSPOSE && (u_length < num_elem * num_comp * num_qpts || v_length < num_elem * num_comp * num_nodes)) ||
10582b730f8bSJeremy L Thompson                   (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem * num_qpts * num_comp || u_length < num_elem * num_comp * num_nodes)));
1059e15f9bd0SJeremy L Thompson       break;
10602b730f8bSJeremy L Thompson     case CEED_EVAL_GRAD:
10612b730f8bSJeremy L Thompson       bad_dims = ((t_mode == CEED_TRANSPOSE && (u_length < num_elem * num_comp * num_qpts * dim || v_length < num_elem * num_comp * num_nodes)) ||
10622b730f8bSJeremy L Thompson                   (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem * num_qpts * num_comp * dim || u_length < num_elem * num_comp * num_nodes)));
1063e15f9bd0SJeremy L Thompson       break;
1064e15f9bd0SJeremy L Thompson     case CEED_EVAL_WEIGHT:
1065d1d35e2fSjeremylt       bad_dims = v_length < num_elem * num_qpts;
1066e15f9bd0SJeremy L Thompson       break;
1067e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
10682b730f8bSJeremy L Thompson     case CEED_EVAL_DIV:
10692b730f8bSJeremy L Thompson       bad_dims = ((t_mode == CEED_TRANSPOSE && (u_length < num_elem * num_comp * num_qpts || v_length < num_elem * num_comp * num_nodes)) ||
10702b730f8bSJeremy L Thompson                   (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem * num_qpts * num_comp || u_length < num_elem * num_comp * num_nodes)));
1071e15f9bd0SJeremy L Thompson       break;
10722b730f8bSJeremy L Thompson     case CEED_EVAL_CURL:
10732b730f8bSJeremy L Thompson       bad_dims = ((t_mode == CEED_TRANSPOSE && (u_length < num_elem * num_comp * num_qpts || v_length < num_elem * num_comp * num_nodes)) ||
10742b730f8bSJeremy L Thompson                   (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem * num_qpts * num_comp || u_length < num_elem * num_comp * num_nodes)));
1075e15f9bd0SJeremy L Thompson       break;
1076e15f9bd0SJeremy L Thompson       // LCOV_EXCL_STOP
1077e15f9bd0SJeremy L Thompson   }
10782b730f8bSJeremy L Thompson   if (bad_dims) {
1079e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
10802b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_DIMENSION, "Input/output vectors too short for basis and evaluation mode");
1081e15f9bd0SJeremy L Thompson     // LCOV_EXCL_STOP
10822b730f8bSJeremy L Thompson   }
1083e15f9bd0SJeremy L Thompson 
10842b730f8bSJeremy L Thompson   CeedCall(basis->Apply(basis, num_elem, t_mode, eval_mode, u, v));
1085e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
10867a982d89SJeremy L. Thompson }
10877a982d89SJeremy L. Thompson 
10887a982d89SJeremy L. Thompson /**
1089b7c9bbdaSJeremy L Thompson   @brief Get Ceed associated with a CeedBasis
1090b7c9bbdaSJeremy L Thompson 
1091*ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
1092b7c9bbdaSJeremy L Thompson   @param[out] ceed  Variable to store Ceed
1093b7c9bbdaSJeremy L Thompson 
1094b7c9bbdaSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1095b7c9bbdaSJeremy L Thompson 
1096b7c9bbdaSJeremy L Thompson   @ref Advanced
1097b7c9bbdaSJeremy L Thompson **/
1098b7c9bbdaSJeremy L Thompson int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) {
1099b7c9bbdaSJeremy L Thompson   *ceed = basis->ceed;
1100b7c9bbdaSJeremy L Thompson   return CEED_ERROR_SUCCESS;
1101b7c9bbdaSJeremy L Thompson }
1102b7c9bbdaSJeremy L Thompson 
1103b7c9bbdaSJeremy L Thompson /**
11049d007619Sjeremylt   @brief Get dimension for given CeedBasis
11059d007619Sjeremylt 
1106*ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
11079d007619Sjeremylt   @param[out] dim   Variable to store dimension of basis
11089d007619Sjeremylt 
11099d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
11109d007619Sjeremylt 
1111b7c9bbdaSJeremy L Thompson   @ref Advanced
11129d007619Sjeremylt **/
11139d007619Sjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) {
11149d007619Sjeremylt   *dim = basis->dim;
1115e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
11169d007619Sjeremylt }
11179d007619Sjeremylt 
11189d007619Sjeremylt /**
1119d99fa3c5SJeremy L Thompson   @brief Get topology for given CeedBasis
1120d99fa3c5SJeremy L Thompson 
1121*ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
1122d99fa3c5SJeremy L Thompson   @param[out] topo  Variable to store topology of basis
1123d99fa3c5SJeremy L Thompson 
1124d99fa3c5SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1125d99fa3c5SJeremy L Thompson 
1126b7c9bbdaSJeremy L Thompson   @ref Advanced
1127d99fa3c5SJeremy L Thompson **/
1128d99fa3c5SJeremy L Thompson int CeedBasisGetTopology(CeedBasis basis, CeedElemTopology *topo) {
1129d99fa3c5SJeremy L Thompson   *topo = basis->topo;
1130e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1131d99fa3c5SJeremy L Thompson }
1132d99fa3c5SJeremy L Thompson 
1133d99fa3c5SJeremy L Thompson /**
113450c301a5SRezgar Shakeri   @brief Get number of Q-vector components for given CeedBasis
113550c301a5SRezgar Shakeri 
1136*ea61e9acSJeremy L Thompson   @param[in]  basis  CeedBasis
113750c301a5SRezgar Shakeri   @param[out] Q_comp Variable to store number of Q-vector components of basis
113850c301a5SRezgar Shakeri 
113950c301a5SRezgar Shakeri   @return An error code: 0 - success, otherwise - failure
114050c301a5SRezgar Shakeri 
114150c301a5SRezgar Shakeri   @ref Advanced
114250c301a5SRezgar Shakeri **/
114350c301a5SRezgar Shakeri int CeedBasisGetNumQuadratureComponents(CeedBasis basis, CeedInt *Q_comp) {
114450c301a5SRezgar Shakeri   *Q_comp = basis->Q_comp;
114550c301a5SRezgar Shakeri   return CEED_ERROR_SUCCESS;
114650c301a5SRezgar Shakeri }
114750c301a5SRezgar Shakeri 
114850c301a5SRezgar Shakeri /**
11499d007619Sjeremylt   @brief Get number of components for given CeedBasis
11509d007619Sjeremylt 
1151*ea61e9acSJeremy L Thompson   @param[in]  basis    CeedBasis
1152d1d35e2fSjeremylt   @param[out] num_comp Variable to store number of components of basis
11539d007619Sjeremylt 
11549d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
11559d007619Sjeremylt 
1156b7c9bbdaSJeremy L Thompson   @ref Advanced
11579d007619Sjeremylt **/
1158d1d35e2fSjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *num_comp) {
1159d1d35e2fSjeremylt   *num_comp = basis->num_comp;
1160e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
11619d007619Sjeremylt }
11629d007619Sjeremylt 
11639d007619Sjeremylt /**
11649d007619Sjeremylt   @brief Get total number of nodes (in dim dimensions) of a CeedBasis
11659d007619Sjeremylt 
1166*ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
11679d007619Sjeremylt   @param[out] P     Variable to store number of nodes
11689d007619Sjeremylt 
11699d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
11709d007619Sjeremylt 
11719d007619Sjeremylt   @ref Utility
11729d007619Sjeremylt **/
11739d007619Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) {
11749d007619Sjeremylt   *P = basis->P;
1175e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
11769d007619Sjeremylt }
11779d007619Sjeremylt 
11789d007619Sjeremylt /**
11799d007619Sjeremylt   @brief Get total number of nodes (in 1 dimension) of a CeedBasis
11809d007619Sjeremylt 
1181*ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
1182d1d35e2fSjeremylt   @param[out] P_1d  Variable to store number of nodes
11839d007619Sjeremylt 
11849d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
11859d007619Sjeremylt 
1186b7c9bbdaSJeremy L Thompson   @ref Advanced
11879d007619Sjeremylt **/
1188d1d35e2fSjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P_1d) {
11892b730f8bSJeremy L Thompson   if (!basis->tensor_basis) {
11909d007619Sjeremylt     // LCOV_EXCL_START
11912b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR, "Cannot supply P_1d for non-tensor basis");
11929d007619Sjeremylt     // LCOV_EXCL_STOP
11932b730f8bSJeremy L Thompson   }
11949d007619Sjeremylt 
1195d1d35e2fSjeremylt   *P_1d = basis->P_1d;
1196e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
11979d007619Sjeremylt }
11989d007619Sjeremylt 
11999d007619Sjeremylt /**
12009d007619Sjeremylt   @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis
12019d007619Sjeremylt 
1202*ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
12039d007619Sjeremylt   @param[out] Q     Variable to store number of quadrature points
12049d007619Sjeremylt 
12059d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
12069d007619Sjeremylt 
12079d007619Sjeremylt   @ref Utility
12089d007619Sjeremylt **/
12099d007619Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) {
12109d007619Sjeremylt   *Q = basis->Q;
1211e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
12129d007619Sjeremylt }
12139d007619Sjeremylt 
12149d007619Sjeremylt /**
12159d007619Sjeremylt   @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis
12169d007619Sjeremylt 
1217*ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
1218d1d35e2fSjeremylt   @param[out] Q_1d  Variable to store number of quadrature points
12199d007619Sjeremylt 
12209d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
12219d007619Sjeremylt 
1222b7c9bbdaSJeremy L Thompson   @ref Advanced
12239d007619Sjeremylt **/
1224d1d35e2fSjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q_1d) {
12252b730f8bSJeremy L Thompson   if (!basis->tensor_basis) {
12269d007619Sjeremylt     // LCOV_EXCL_START
12272b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR, "Cannot supply Q_1d for non-tensor basis");
12289d007619Sjeremylt     // LCOV_EXCL_STOP
12292b730f8bSJeremy L Thompson   }
12309d007619Sjeremylt 
1231d1d35e2fSjeremylt   *Q_1d = basis->Q_1d;
1232e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
12339d007619Sjeremylt }
12349d007619Sjeremylt 
12359d007619Sjeremylt /**
1236*ea61e9acSJeremy L Thompson   @brief Get reference coordinates of quadrature points (in dim dimensions) of a CeedBasis
12379d007619Sjeremylt 
1238*ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
1239d1d35e2fSjeremylt   @param[out] q_ref Variable to store reference coordinates of quadrature points
12409d007619Sjeremylt 
12419d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
12429d007619Sjeremylt 
1243b7c9bbdaSJeremy L Thompson   @ref Advanced
12449d007619Sjeremylt **/
1245d1d35e2fSjeremylt int CeedBasisGetQRef(CeedBasis basis, const CeedScalar **q_ref) {
1246d1d35e2fSjeremylt   *q_ref = basis->q_ref_1d;
1247e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
12489d007619Sjeremylt }
12499d007619Sjeremylt 
12509d007619Sjeremylt /**
1251*ea61e9acSJeremy L Thompson   @brief Get quadrature weights of quadrature points (in dim dimensions) of a CeedBasis
12529d007619Sjeremylt 
1253*ea61e9acSJeremy L Thompson   @param[in]  basis    CeedBasis
1254d1d35e2fSjeremylt   @param[out] q_weight Variable to store quadrature weights
12559d007619Sjeremylt 
12569d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
12579d007619Sjeremylt 
1258b7c9bbdaSJeremy L Thompson   @ref Advanced
12599d007619Sjeremylt **/
1260d1d35e2fSjeremylt int CeedBasisGetQWeights(CeedBasis basis, const CeedScalar **q_weight) {
1261d1d35e2fSjeremylt   *q_weight = basis->q_weight_1d;
1262e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
12639d007619Sjeremylt }
12649d007619Sjeremylt 
12659d007619Sjeremylt /**
12669d007619Sjeremylt   @brief Get interpolation matrix of a CeedBasis
12679d007619Sjeremylt 
1268*ea61e9acSJeremy L Thompson   @param[in]  basis  CeedBasis
12699d007619Sjeremylt   @param[out] interp Variable to store interpolation matrix
12709d007619Sjeremylt 
12719d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
12729d007619Sjeremylt 
1273b7c9bbdaSJeremy L Thompson   @ref Advanced
12749d007619Sjeremylt **/
12756c58de82SJeremy L Thompson int CeedBasisGetInterp(CeedBasis basis, const CeedScalar **interp) {
1276d1d35e2fSjeremylt   if (!basis->interp && basis->tensor_basis) {
12779d007619Sjeremylt     // Allocate
12782b730f8bSJeremy L Thompson     CeedCall(CeedMalloc(basis->Q * basis->P, &basis->interp));
12799d007619Sjeremylt 
12809d007619Sjeremylt     // Initialize
12812b730f8bSJeremy L Thompson     for (CeedInt i = 0; i < basis->Q * basis->P; i++) basis->interp[i] = 1.0;
12829d007619Sjeremylt 
12839d007619Sjeremylt     // Calculate
12842b730f8bSJeremy L Thompson     for (CeedInt d = 0; d < basis->dim; d++) {
12852b730f8bSJeremy L Thompson       for (CeedInt qpt = 0; qpt < basis->Q; qpt++) {
12869d007619Sjeremylt         for (CeedInt node = 0; node < basis->P; node++) {
1287d1d35e2fSjeremylt           CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d;
1288d1d35e2fSjeremylt           CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d;
1289d1d35e2fSjeremylt           basis->interp[qpt * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p];
12909d007619Sjeremylt         }
12919d007619Sjeremylt       }
12922b730f8bSJeremy L Thompson     }
12932b730f8bSJeremy L Thompson   }
12949d007619Sjeremylt   *interp = basis->interp;
1295e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
12969d007619Sjeremylt }
12979d007619Sjeremylt 
12989d007619Sjeremylt /**
12999d007619Sjeremylt   @brief Get 1D interpolation matrix of a tensor product CeedBasis
13009d007619Sjeremylt 
1301*ea61e9acSJeremy L Thompson   @param[in]  basis     CeedBasis
1302d1d35e2fSjeremylt   @param[out] interp_1d Variable to store interpolation matrix
13039d007619Sjeremylt 
13049d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
13059d007619Sjeremylt 
13069d007619Sjeremylt   @ref Backend
13079d007619Sjeremylt **/
1308d1d35e2fSjeremylt int CeedBasisGetInterp1D(CeedBasis basis, const CeedScalar **interp_1d) {
13092b730f8bSJeremy L Thompson   if (!basis->tensor_basis) {
13109d007619Sjeremylt     // LCOV_EXCL_START
13112b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR, "CeedBasis is not a tensor product basis.");
13129d007619Sjeremylt     // LCOV_EXCL_STOP
13132b730f8bSJeremy L Thompson   }
13149d007619Sjeremylt 
1315d1d35e2fSjeremylt   *interp_1d = basis->interp_1d;
1316e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
13179d007619Sjeremylt }
13189d007619Sjeremylt 
13199d007619Sjeremylt /**
13209d007619Sjeremylt   @brief Get gradient matrix of a CeedBasis
13219d007619Sjeremylt 
1322*ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
13239d007619Sjeremylt   @param[out] grad  Variable to store gradient matrix
13249d007619Sjeremylt 
13259d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
13269d007619Sjeremylt 
1327b7c9bbdaSJeremy L Thompson   @ref Advanced
13289d007619Sjeremylt **/
13296c58de82SJeremy L Thompson int CeedBasisGetGrad(CeedBasis basis, const CeedScalar **grad) {
1330d1d35e2fSjeremylt   if (!basis->grad && basis->tensor_basis) {
13319d007619Sjeremylt     // Allocate
13322b730f8bSJeremy L Thompson     CeedCall(CeedMalloc(basis->dim * basis->Q * basis->P, &basis->grad));
13339d007619Sjeremylt 
13349d007619Sjeremylt     // Initialize
13352b730f8bSJeremy L Thompson     for (CeedInt i = 0; i < basis->dim * basis->Q * basis->P; i++) basis->grad[i] = 1.0;
13369d007619Sjeremylt 
13379d007619Sjeremylt     // Calculate
13382b730f8bSJeremy L Thompson     for (CeedInt d = 0; d < basis->dim; d++) {
13392b730f8bSJeremy L Thompson       for (CeedInt i = 0; i < basis->dim; i++) {
13402b730f8bSJeremy L Thompson         for (CeedInt qpt = 0; qpt < basis->Q; qpt++) {
13419d007619Sjeremylt           for (CeedInt node = 0; node < basis->P; node++) {
1342d1d35e2fSjeremylt             CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d;
1343d1d35e2fSjeremylt             CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d;
13442b730f8bSJeremy L Thompson             if (i == d) basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->grad_1d[q * basis->P_1d + p];
13452b730f8bSJeremy L Thompson             else basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p];
13462b730f8bSJeremy L Thompson           }
13472b730f8bSJeremy L Thompson         }
13482b730f8bSJeremy L Thompson       }
13499d007619Sjeremylt     }
13509d007619Sjeremylt   }
13519d007619Sjeremylt   *grad = basis->grad;
1352e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
13539d007619Sjeremylt }
13549d007619Sjeremylt 
13559d007619Sjeremylt /**
13569d007619Sjeremylt   @brief Get 1D gradient matrix of a tensor product CeedBasis
13579d007619Sjeremylt 
1358*ea61e9acSJeremy L Thompson   @param[in]  basis   CeedBasis
1359d1d35e2fSjeremylt   @param[out] grad_1d Variable to store gradient matrix
13609d007619Sjeremylt 
13619d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
13629d007619Sjeremylt 
1363b7c9bbdaSJeremy L Thompson   @ref Advanced
13649d007619Sjeremylt **/
1365d1d35e2fSjeremylt int CeedBasisGetGrad1D(CeedBasis basis, const CeedScalar **grad_1d) {
13662b730f8bSJeremy L Thompson   if (!basis->tensor_basis) {
13679d007619Sjeremylt     // LCOV_EXCL_START
13682b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR, "CeedBasis is not a tensor product basis.");
13699d007619Sjeremylt     // LCOV_EXCL_STOP
13702b730f8bSJeremy L Thompson   }
13719d007619Sjeremylt 
1372d1d35e2fSjeremylt   *grad_1d = basis->grad_1d;
1373e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
13749d007619Sjeremylt }
13759d007619Sjeremylt 
13769d007619Sjeremylt /**
137750c301a5SRezgar Shakeri   @brief Get divergence matrix of a CeedBasis
137850c301a5SRezgar Shakeri 
1379*ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
138050c301a5SRezgar Shakeri   @param[out] div   Variable to store divergence matrix
138150c301a5SRezgar Shakeri 
138250c301a5SRezgar Shakeri   @return An error code: 0 - success, otherwise - failure
138350c301a5SRezgar Shakeri 
138450c301a5SRezgar Shakeri   @ref Advanced
138550c301a5SRezgar Shakeri **/
138650c301a5SRezgar Shakeri int CeedBasisGetDiv(CeedBasis basis, const CeedScalar **div) {
13872b730f8bSJeremy L Thompson   if (!basis->div) {
138850c301a5SRezgar Shakeri     // LCOV_EXCL_START
13892b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR, "CeedBasis does not have divergence matrix.");
139050c301a5SRezgar Shakeri     // LCOV_EXCL_STOP
13912b730f8bSJeremy L Thompson   }
139250c301a5SRezgar Shakeri 
139350c301a5SRezgar Shakeri   *div = basis->div;
139450c301a5SRezgar Shakeri   return CEED_ERROR_SUCCESS;
139550c301a5SRezgar Shakeri }
139650c301a5SRezgar Shakeri 
139750c301a5SRezgar Shakeri /**
13987a982d89SJeremy L. Thompson   @brief Destroy a CeedBasis
13997a982d89SJeremy L. Thompson 
1400*ea61e9acSJeremy L Thompson   @param[in,out] basis CeedBasis to destroy
14017a982d89SJeremy L. Thompson 
14027a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
14037a982d89SJeremy L. Thompson 
14047a982d89SJeremy L. Thompson   @ref User
14057a982d89SJeremy L. Thompson **/
14067a982d89SJeremy L. Thompson int CeedBasisDestroy(CeedBasis *basis) {
1407d1d35e2fSjeremylt   if (!*basis || --(*basis)->ref_count > 0) return CEED_ERROR_SUCCESS;
14082b730f8bSJeremy L Thompson   if ((*basis)->Destroy) CeedCall((*basis)->Destroy(*basis));
14092b730f8bSJeremy L Thompson   if ((*basis)->contract) CeedCall(CeedTensorContractDestroy(&(*basis)->contract));
14102b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->interp));
14112b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->interp_1d));
14122b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->grad));
14132b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->div));
14142b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->grad_1d));
14152b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->q_ref_1d));
14162b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->q_weight_1d));
14172b730f8bSJeremy L Thompson   CeedCall(CeedDestroy(&(*basis)->ceed));
14182b730f8bSJeremy L Thompson   CeedCall(CeedFree(basis));
1419e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
14207a982d89SJeremy L. Thompson }
14217a982d89SJeremy L. Thompson 
14227a982d89SJeremy L. Thompson /**
1423b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre quadrature
1424b11c1e72Sjeremylt 
1425*ea61e9acSJeremy L Thompson   @param[in]  Q           Number of quadrature points (integrates polynomials of degree 2*Q-1 exactly)
1426d1d35e2fSjeremylt   @param[out] q_ref_1d    Array of length Q to hold the abscissa on [-1, 1]
1427d1d35e2fSjeremylt   @param[out] q_weight_1d Array of length Q to hold the weights
1428b11c1e72Sjeremylt 
1429b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1430dfdf5a53Sjeremylt 
1431dfdf5a53Sjeremylt   @ref Utility
1432b11c1e72Sjeremylt **/
14332b730f8bSJeremy L Thompson int CeedGaussQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) {
1434d7b241e6Sjeremylt   // Allocate
1435d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0 * atan(1.0);
1436d1d35e2fSjeremylt   // Build q_ref_1d, q_weight_1d
143792ae7e47SJeremy L Thompson   for (CeedInt i = 0; i <= Q / 2; i++) {
1438d7b241e6Sjeremylt     // Guess
1439d7b241e6Sjeremylt     xi = cos(PI * (CeedScalar)(2 * i + 1) / ((CeedScalar)(2 * Q)));
1440d7b241e6Sjeremylt     // Pn(xi)
1441d7b241e6Sjeremylt     P0 = 1.0;
1442d7b241e6Sjeremylt     P1 = xi;
1443d7b241e6Sjeremylt     P2 = 0.0;
144492ae7e47SJeremy L Thompson     for (CeedInt j = 2; j <= Q; j++) {
1445d7b241e6Sjeremylt       P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
1446d7b241e6Sjeremylt       P0 = P1;
1447d7b241e6Sjeremylt       P1 = P2;
1448d7b241e6Sjeremylt     }
1449d7b241e6Sjeremylt     // First Newton Step
1450d7b241e6Sjeremylt     dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
1451d7b241e6Sjeremylt     xi  = xi - P2 / dP2;
1452d7b241e6Sjeremylt     // Newton to convergence
145392ae7e47SJeremy L Thompson     for (CeedInt k = 0; k < 100 && fabs(P2) > 10 * CEED_EPSILON; k++) {
1454d7b241e6Sjeremylt       P0 = 1.0;
1455d7b241e6Sjeremylt       P1 = xi;
145692ae7e47SJeremy L Thompson       for (CeedInt j = 2; j <= Q; j++) {
1457d7b241e6Sjeremylt         P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
1458d7b241e6Sjeremylt         P0 = P1;
1459d7b241e6Sjeremylt         P1 = P2;
1460d7b241e6Sjeremylt       }
1461d7b241e6Sjeremylt       dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
1462d7b241e6Sjeremylt       xi  = xi - P2 / dP2;
1463d7b241e6Sjeremylt     }
1464d7b241e6Sjeremylt     // Save xi, wi
1465d7b241e6Sjeremylt     wi                     = 2.0 / ((1.0 - xi * xi) * dP2 * dP2);
1466d1d35e2fSjeremylt     q_weight_1d[i]         = wi;
1467d1d35e2fSjeremylt     q_weight_1d[Q - 1 - i] = wi;
1468d1d35e2fSjeremylt     q_ref_1d[i]            = -xi;
1469d1d35e2fSjeremylt     q_ref_1d[Q - 1 - i]    = xi;
1470d7b241e6Sjeremylt   }
1471e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1472d7b241e6Sjeremylt }
1473d7b241e6Sjeremylt 
1474b11c1e72Sjeremylt /**
1475b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre-Lobatto quadrature
1476b11c1e72Sjeremylt 
1477*ea61e9acSJeremy L Thompson   @param[in]  Q           Number of quadrature points (integrates polynomials of degree 2*Q-3 exactly)
1478d1d35e2fSjeremylt   @param[out] q_ref_1d    Array of length Q to hold the abscissa on [-1, 1]
1479d1d35e2fSjeremylt   @param[out] q_weight_1d Array of length Q to hold the weights
1480b11c1e72Sjeremylt 
1481b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1482dfdf5a53Sjeremylt 
1483dfdf5a53Sjeremylt   @ref Utility
1484b11c1e72Sjeremylt **/
14852b730f8bSJeremy L Thompson int CeedLobattoQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) {
1486d7b241e6Sjeremylt   // Allocate
1487d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0 * atan(1.0);
1488d1d35e2fSjeremylt   // Build q_ref_1d, q_weight_1d
1489d7b241e6Sjeremylt   // Set endpoints
14902b730f8bSJeremy L Thompson   if (Q < 2) {
1491b0d62198Sjeremylt     // LCOV_EXCL_START
14922b730f8bSJeremy L Thompson     return CeedError(NULL, CEED_ERROR_DIMENSION, "Cannot create Lobatto quadrature with Q=%" CeedInt_FMT " < 2 points", Q);
1493b0d62198Sjeremylt     // LCOV_EXCL_STOP
14942b730f8bSJeremy L Thompson   }
1495d7b241e6Sjeremylt   wi = 2.0 / ((CeedScalar)(Q * (Q - 1)));
1496d1d35e2fSjeremylt   if (q_weight_1d) {
1497d1d35e2fSjeremylt     q_weight_1d[0]     = wi;
1498d1d35e2fSjeremylt     q_weight_1d[Q - 1] = wi;
1499d7b241e6Sjeremylt   }
1500d1d35e2fSjeremylt   q_ref_1d[0]     = -1.0;
1501d1d35e2fSjeremylt   q_ref_1d[Q - 1] = 1.0;
1502d7b241e6Sjeremylt   // Interior
150392ae7e47SJeremy L Thompson   for (CeedInt i = 1; i <= (Q - 1) / 2; i++) {
1504d7b241e6Sjeremylt     // Guess
1505d7b241e6Sjeremylt     xi = cos(PI * (CeedScalar)(i) / (CeedScalar)(Q - 1));
1506d7b241e6Sjeremylt     // Pn(xi)
1507d7b241e6Sjeremylt     P0 = 1.0;
1508d7b241e6Sjeremylt     P1 = xi;
1509d7b241e6Sjeremylt     P2 = 0.0;
151092ae7e47SJeremy L Thompson     for (CeedInt j = 2; j < Q; j++) {
1511d7b241e6Sjeremylt       P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
1512d7b241e6Sjeremylt       P0 = P1;
1513d7b241e6Sjeremylt       P1 = P2;
1514d7b241e6Sjeremylt     }
1515d7b241e6Sjeremylt     // First Newton step
1516d7b241e6Sjeremylt     dP2  = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
1517d7b241e6Sjeremylt     d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi);
1518d7b241e6Sjeremylt     xi   = xi - dP2 / d2P2;
1519d7b241e6Sjeremylt     // Newton to convergence
152092ae7e47SJeremy L Thompson     for (CeedInt k = 0; k < 100 && fabs(dP2) > 10 * CEED_EPSILON; k++) {
1521d7b241e6Sjeremylt       P0 = 1.0;
1522d7b241e6Sjeremylt       P1 = xi;
152392ae7e47SJeremy L Thompson       for (CeedInt j = 2; j < Q; j++) {
1524d7b241e6Sjeremylt         P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
1525d7b241e6Sjeremylt         P0 = P1;
1526d7b241e6Sjeremylt         P1 = P2;
1527d7b241e6Sjeremylt       }
1528d7b241e6Sjeremylt       dP2  = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
1529d7b241e6Sjeremylt       d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi);
1530d7b241e6Sjeremylt       xi   = xi - dP2 / d2P2;
1531d7b241e6Sjeremylt     }
1532d7b241e6Sjeremylt     // Save xi, wi
1533d7b241e6Sjeremylt     wi = 2.0 / (((CeedScalar)(Q * (Q - 1))) * P2 * P2);
1534d1d35e2fSjeremylt     if (q_weight_1d) {
1535d1d35e2fSjeremylt       q_weight_1d[i]         = wi;
1536d1d35e2fSjeremylt       q_weight_1d[Q - 1 - i] = wi;
1537d7b241e6Sjeremylt     }
1538d1d35e2fSjeremylt     q_ref_1d[i]         = -xi;
1539d1d35e2fSjeremylt     q_ref_1d[Q - 1 - i] = xi;
1540d7b241e6Sjeremylt   }
1541e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1542d7b241e6Sjeremylt }
1543d7b241e6Sjeremylt 
1544dfdf5a53Sjeremylt /**
154595bb1877Svaleriabarra   @brief Return QR Factorization of a matrix
1546b11c1e72Sjeremylt 
1547*ea61e9acSJeremy L Thompson   @param[in]     ceed Ceed context for error handling
154852bfb9bbSJeremy L Thompson   @param[in,out] mat  Row-major matrix to be factorized in place
154952bfb9bbSJeremy L Thompson   @param[in,out] tau  Vector of length m of scaling factors
1550*ea61e9acSJeremy L Thompson   @param[in]     m    Number of rows
1551*ea61e9acSJeremy L Thompson   @param[in]     n    Number of columns
1552b11c1e72Sjeremylt 
1553b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1554dfdf5a53Sjeremylt 
1555dfdf5a53Sjeremylt   @ref Utility
1556b11c1e72Sjeremylt **/
15572b730f8bSJeremy L Thompson int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau, CeedInt m, CeedInt n) {
1558d7b241e6Sjeremylt   CeedScalar v[m];
1559d7b241e6Sjeremylt 
15602b730f8bSJeremy L Thompson   // Check matrix shape
15612b730f8bSJeremy L Thompson   if (n > m) {
1562c042f62fSJeremy L Thompson     // LCOV_EXCL_START
15632b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute QR factorization with n > m");
1564c042f62fSJeremy L Thompson     // LCOV_EXCL_STOP
15652b730f8bSJeremy L Thompson   }
1566a7bd39daSjeremylt 
156752bfb9bbSJeremy L Thompson   for (CeedInt i = 0; i < n; i++) {
1568bde37e8cSJed Brown     if (i >= m - 1) {  // last row of matrix, no reflection needed
1569bde37e8cSJed Brown       tau[i] = 0.;
1570bde37e8cSJed Brown       break;
1571bde37e8cSJed Brown     }
1572d7b241e6Sjeremylt     // Calculate Householder vector, magnitude
1573d7b241e6Sjeremylt     CeedScalar sigma = 0.0;
1574d7b241e6Sjeremylt     v[i]             = mat[i + n * i];
157552bfb9bbSJeremy L Thompson     for (CeedInt j = i + 1; j < m; j++) {
1576d7b241e6Sjeremylt       v[j] = mat[i + n * j];
1577d7b241e6Sjeremylt       sigma += v[j] * v[j];
1578d7b241e6Sjeremylt     }
1579d7b241e6Sjeremylt     CeedScalar norm = sqrt(v[i] * v[i] + sigma);  // norm of v[i:m]
1580673160d7Sjeremylt     CeedScalar R_ii = -copysign(norm, v[i]);
1581673160d7Sjeremylt     v[i] -= R_ii;
1582d7b241e6Sjeremylt     // norm of v[i:m] after modification above and scaling below
1583d7b241e6Sjeremylt     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
1584d7b241e6Sjeremylt     //   tau = 2 / (norm*norm)
1585d7b241e6Sjeremylt     tau[i] = 2 * v[i] * v[i] / (v[i] * v[i] + sigma);
15862b730f8bSJeremy L Thompson     for (CeedInt j = i + 1; j < m; j++) v[j] /= v[i];
1587d7b241e6Sjeremylt 
1588d7b241e6Sjeremylt     // Apply Householder reflector to lower right panel
1589d7b241e6Sjeremylt     CeedHouseholderReflect(&mat[i * n + i + 1], &v[i], tau[i], m - i, n - i - 1, n, 1);
1590d7b241e6Sjeremylt     // Save v
1591673160d7Sjeremylt     mat[i + n * i] = R_ii;
15922b730f8bSJeremy L Thompson     for (CeedInt j = i + 1; j < m; j++) mat[i + n * j] = v[j];
1593d7b241e6Sjeremylt   }
1594e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1595d7b241e6Sjeremylt }
1596d7b241e6Sjeremylt 
1597b11c1e72Sjeremylt /**
1598*ea61e9acSJeremy L Thompson   @brief Return symmetric Schur decomposition of the symmetric matrix mat via symmetric QR factorization
159952bfb9bbSJeremy L Thompson 
1600*ea61e9acSJeremy L Thompson   @param[in]     ceed   Ceed context for error handling
160152bfb9bbSJeremy L Thompson   @param[in,out] mat    Row-major matrix to be factorized in place
1602460bf743SValeria Barra   @param[out]    lambda Vector of length n of eigenvalues
1603*ea61e9acSJeremy L Thompson   @param[in]     n      Number of rows/columns
160452bfb9bbSJeremy L Thompson 
160552bfb9bbSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
160652bfb9bbSJeremy L Thompson 
160752bfb9bbSJeremy L Thompson   @ref Utility
160852bfb9bbSJeremy L Thompson **/
16092b730f8bSJeremy L Thompson CeedPragmaOptimizeOff int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat, CeedScalar *lambda, CeedInt n) {
161052bfb9bbSJeremy L Thompson   // Check bounds for clang-tidy
16112b730f8bSJeremy L Thompson   if (n < 2) {
1612c042f62fSJeremy L Thompson     // LCOV_EXCL_START
16132b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute symmetric Schur decomposition of scalars");
1614c042f62fSJeremy L Thompson     // LCOV_EXCL_STOP
16152b730f8bSJeremy L Thompson   }
161652bfb9bbSJeremy L Thompson 
1617673160d7Sjeremylt   CeedScalar v[n - 1], tau[n - 1], mat_T[n * n];
161852bfb9bbSJeremy L Thompson 
1619673160d7Sjeremylt   // Copy mat to mat_T and set mat to I
1620673160d7Sjeremylt   memcpy(mat_T, mat, n * n * sizeof(mat[0]));
16212b730f8bSJeremy L Thompson   for (CeedInt i = 0; i < n; i++) {
16222b730f8bSJeremy L Thompson     for (CeedInt j = 0; j < n; j++) mat[j + n * i] = (i == j) ? 1 : 0;
16232b730f8bSJeremy L Thompson   }
162452bfb9bbSJeremy L Thompson 
162552bfb9bbSJeremy L Thompson   // Reduce to tridiagonal
162652bfb9bbSJeremy L Thompson   for (CeedInt i = 0; i < n - 1; i++) {
162752bfb9bbSJeremy L Thompson     // Calculate Householder vector, magnitude
162852bfb9bbSJeremy L Thompson     CeedScalar sigma = 0.0;
1629673160d7Sjeremylt     v[i]             = mat_T[i + n * (i + 1)];
163052bfb9bbSJeremy L Thompson     for (CeedInt j = i + 1; j < n - 1; j++) {
1631673160d7Sjeremylt       v[j] = mat_T[i + n * (j + 1)];
163252bfb9bbSJeremy L Thompson       sigma += v[j] * v[j];
163352bfb9bbSJeremy L Thompson     }
163452bfb9bbSJeremy L Thompson     CeedScalar norm = sqrt(v[i] * v[i] + sigma);  // norm of v[i:n-1]
1635673160d7Sjeremylt     CeedScalar R_ii = -copysign(norm, v[i]);
1636673160d7Sjeremylt     v[i] -= R_ii;
163752bfb9bbSJeremy L Thompson     // norm of v[i:m] after modification above and scaling below
163852bfb9bbSJeremy L Thompson     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
163952bfb9bbSJeremy L Thompson     //   tau = 2 / (norm*norm)
164080a9ef05SNatalie Beams     tau[i] = i == n - 2 ? 2 : 2 * v[i] * v[i] / (v[i] * v[i] + sigma);
16412b730f8bSJeremy L Thompson     for (CeedInt j = i + 1; j < n - 1; j++) v[j] /= v[i];
164252bfb9bbSJeremy L Thompson 
164352bfb9bbSJeremy L Thompson     // Update sub and super diagonal
164452bfb9bbSJeremy L Thompson     for (CeedInt j = i + 2; j < n; j++) {
16452b730f8bSJeremy L Thompson       mat_T[i + n * j] = 0;
16462b730f8bSJeremy L Thompson       mat_T[j + n * i] = 0;
164752bfb9bbSJeremy L Thompson     }
164852bfb9bbSJeremy L Thompson     // Apply symmetric Householder reflector to lower right panel
16492b730f8bSJeremy L Thompson     CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1);
16502b730f8bSJeremy L Thompson     CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), 1, n);
1651673160d7Sjeremylt 
165252bfb9bbSJeremy L Thompson     // Save v
1653673160d7Sjeremylt     mat_T[i + n * (i + 1)] = R_ii;
1654673160d7Sjeremylt     mat_T[(i + 1) + n * i] = R_ii;
165552bfb9bbSJeremy L Thompson     for (CeedInt j = i + 1; j < n - 1; j++) {
1656673160d7Sjeremylt       mat_T[i + n * (j + 1)] = v[j];
165752bfb9bbSJeremy L Thompson     }
165852bfb9bbSJeremy L Thompson   }
165952bfb9bbSJeremy L Thompson   // Backwards accumulation of Q
166052bfb9bbSJeremy L Thompson   for (CeedInt i = n - 2; i >= 0; i--) {
166185cf89eaSjeremylt     if (tau[i] > 0.0) {
166252bfb9bbSJeremy L Thompson       v[i] = 1;
166352bfb9bbSJeremy L Thompson       for (CeedInt j = i + 1; j < n - 1; j++) {
1664673160d7Sjeremylt         v[j]                   = mat_T[i + n * (j + 1)];
1665673160d7Sjeremylt         mat_T[i + n * (j + 1)] = 0;
166652bfb9bbSJeremy L Thompson       }
16672b730f8bSJeremy L Thompson       CeedHouseholderReflect(&mat[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1);
166852bfb9bbSJeremy L Thompson     }
166985cf89eaSjeremylt   }
167052bfb9bbSJeremy L Thompson 
167152bfb9bbSJeremy L Thompson   // Reduce sub and super diagonal
1672673160d7Sjeremylt   CeedInt    p = 0, q = 0, itr = 0, max_itr = n * n * n * n;
1673673160d7Sjeremylt   CeedScalar tol = CEED_EPSILON;
167452bfb9bbSJeremy L Thompson 
1675673160d7Sjeremylt   while (itr < max_itr) {
167652bfb9bbSJeremy L Thompson     // Update p, q, size of reduced portions of diagonal
16772b730f8bSJeremy L Thompson     p = 0;
16782b730f8bSJeremy L Thompson     q = 0;
167952bfb9bbSJeremy L Thompson     for (CeedInt i = n - 2; i >= 0; i--) {
16802b730f8bSJeremy L Thompson       if (fabs(mat_T[i + n * (i + 1)]) < tol) q += 1;
16812b730f8bSJeremy L Thompson       else break;
168252bfb9bbSJeremy L Thompson     }
1683673160d7Sjeremylt     for (CeedInt i = 0; i < n - q - 1; i++) {
16842b730f8bSJeremy L Thompson       if (fabs(mat_T[i + n * (i + 1)]) < tol) p += 1;
16852b730f8bSJeremy L Thompson       else break;
168652bfb9bbSJeremy L Thompson     }
168752bfb9bbSJeremy L Thompson     if (q == n - 1) break;  // Finished reducing
168852bfb9bbSJeremy L Thompson 
168952bfb9bbSJeremy L Thompson     // Reduce tridiagonal portion
16902b730f8bSJeremy L Thompson     CeedScalar t_nn = mat_T[(n - 1 - q) + n * (n - 1 - q)], t_nnm1 = mat_T[(n - 2 - q) + n * (n - 1 - q)];
1691673160d7Sjeremylt     CeedScalar d  = (mat_T[(n - 2 - q) + n * (n - 2 - q)] - t_nn) / 2;
16922b730f8bSJeremy L Thompson     CeedScalar mu = t_nn - t_nnm1 * t_nnm1 / (d + copysign(sqrt(d * d + t_nnm1 * t_nnm1), d));
1693673160d7Sjeremylt     CeedScalar x  = mat_T[p + n * p] - mu;
1694673160d7Sjeremylt     CeedScalar z  = mat_T[p + n * (p + 1)];
1695673160d7Sjeremylt     for (CeedInt k = p; k < n - q - 1; k++) {
169652bfb9bbSJeremy L Thompson       // Compute Givens rotation
169752bfb9bbSJeremy L Thompson       CeedScalar c = 1, s = 0;
169852bfb9bbSJeremy L Thompson       if (fabs(z) > tol) {
169952bfb9bbSJeremy L Thompson         if (fabs(z) > fabs(x)) {
170052bfb9bbSJeremy L Thompson           CeedScalar tau = -x / z;
170152bfb9bbSJeremy L Thompson           s = 1 / sqrt(1 + tau * tau), c = s * tau;
170252bfb9bbSJeremy L Thompson         } else {
170352bfb9bbSJeremy L Thompson           CeedScalar tau = -z / x;
170452bfb9bbSJeremy L Thompson           c = 1 / sqrt(1 + tau * tau), s = c * tau;
170552bfb9bbSJeremy L Thompson         }
170652bfb9bbSJeremy L Thompson       }
170752bfb9bbSJeremy L Thompson 
170852bfb9bbSJeremy L Thompson       // Apply Givens rotation to T
1709673160d7Sjeremylt       CeedGivensRotation(mat_T, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n);
1710673160d7Sjeremylt       CeedGivensRotation(mat_T, c, s, CEED_TRANSPOSE, k, k + 1, n, n);
171152bfb9bbSJeremy L Thompson 
171252bfb9bbSJeremy L Thompson       // Apply Givens rotation to Q
171352bfb9bbSJeremy L Thompson       CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n);
171452bfb9bbSJeremy L Thompson 
171552bfb9bbSJeremy L Thompson       // Update x, z
171652bfb9bbSJeremy L Thompson       if (k < n - q - 2) {
1717673160d7Sjeremylt         x = mat_T[k + n * (k + 1)];
1718673160d7Sjeremylt         z = mat_T[k + n * (k + 2)];
171952bfb9bbSJeremy L Thompson       }
172052bfb9bbSJeremy L Thompson     }
172152bfb9bbSJeremy L Thompson     itr++;
172252bfb9bbSJeremy L Thompson   }
1723673160d7Sjeremylt 
172452bfb9bbSJeremy L Thompson   // Save eigenvalues
17252b730f8bSJeremy L Thompson   for (CeedInt i = 0; i < n; i++) lambda[i] = mat_T[i + n * i];
172652bfb9bbSJeremy L Thompson 
172752bfb9bbSJeremy L Thompson   // Check convergence
17282b730f8bSJeremy L Thompson   if (itr == max_itr && q < n - 1) {
1729c042f62fSJeremy L Thompson     // LCOV_EXCL_START
17302b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_MINOR, "Symmetric QR failed to converge");
1731c042f62fSJeremy L Thompson     // LCOV_EXCL_STOP
17322b730f8bSJeremy L Thompson   }
1733e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
173452bfb9bbSJeremy L Thompson }
173503d18186Sjeremylt CeedPragmaOptimizeOn
173652bfb9bbSJeremy L Thompson 
173752bfb9bbSJeremy L Thompson     /**
1738*ea61e9acSJeremy L Thompson       @brief Return Simultaneous Diagonalization of two matrices.
1739*ea61e9acSJeremy L Thompson                This solves the generalized eigenvalue problem A x = lambda B x, where A and B are symmetric and B is positive definite.
1740*ea61e9acSJeremy L Thompson                We generate the matrix X and vector Lambda such that X^T A X = Lambda and X^T B X = I.
1741*ea61e9acSJeremy L Thompson                This is equivalent to the LAPACK routine 'sygv' with TYPE = 1.
174252bfb9bbSJeremy L Thompson 
1743*ea61e9acSJeremy L Thompson       @param[in]  ceed   Ceed context for error handling
1744d1d35e2fSjeremylt       @param[in]  mat_A  Row-major matrix to be factorized with eigenvalues
1745d1d35e2fSjeremylt       @param[in]  mat_B  Row-major matrix to be factorized to identity
1746d3331725Sjeremylt       @param[out] mat_X  Row-major orthogonal matrix
1747460bf743SValeria Barra       @param[out] lambda Vector of length n of generalized eigenvalues
1748*ea61e9acSJeremy L Thompson       @param[in]  n      Number of rows/columns
174952bfb9bbSJeremy L Thompson 
175052bfb9bbSJeremy L Thompson       @return An error code: 0 - success, otherwise - failure
175152bfb9bbSJeremy L Thompson 
175252bfb9bbSJeremy L Thompson       @ref Utility
175352bfb9bbSJeremy L Thompson     **/
17542b730f8bSJeremy L Thompson     CeedPragmaOptimizeOff int
17552b730f8bSJeremy L Thompson     CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *mat_A, CeedScalar *mat_B, CeedScalar *mat_X, CeedScalar *lambda, CeedInt n) {
1756d3331725Sjeremylt   CeedScalar *mat_C, *mat_G, *vec_D;
17572b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(n * n, &mat_C));
17582b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(n * n, &mat_G));
17592b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(n, &vec_D));
176052bfb9bbSJeremy L Thompson 
176152bfb9bbSJeremy L Thompson   // Compute B = G D G^T
176278464608Sjeremylt   memcpy(mat_G, mat_B, n * n * sizeof(mat_B[0]));
17632b730f8bSJeremy L Thompson   CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_G, vec_D, n));
176452bfb9bbSJeremy L Thompson 
176585cf89eaSjeremylt   // Sort eigenvalues
17662b730f8bSJeremy L Thompson   for (CeedInt i = n - 1; i >= 0; i--) {
176785cf89eaSjeremylt     for (CeedInt j = 0; j < i; j++) {
176885cf89eaSjeremylt       if (fabs(vec_D[j]) > fabs(vec_D[j + 1])) {
176985cf89eaSjeremylt         CeedScalar temp;
17702b730f8bSJeremy L Thompson         temp         = vec_D[j];
17712b730f8bSJeremy L Thompson         vec_D[j]     = vec_D[j + 1];
17722b730f8bSJeremy L Thompson         vec_D[j + 1] = temp;
177385cf89eaSjeremylt         for (CeedInt k = 0; k < n; k++) {
17742b730f8bSJeremy L Thompson           temp                 = mat_G[k * n + j];
17752b730f8bSJeremy L Thompson           mat_G[k * n + j]     = mat_G[k * n + j + 1];
17762b730f8bSJeremy L Thompson           mat_G[k * n + j + 1] = temp;
17772b730f8bSJeremy L Thompson         }
177885cf89eaSjeremylt       }
177985cf89eaSjeremylt     }
178085cf89eaSjeremylt   }
178185cf89eaSjeremylt 
1782fb551037Sjeremylt   // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T
1783fb551037Sjeremylt   //           = D^-1/2 G^T A G D^-1/2
1784d3331725Sjeremylt   // -- D = D^-1/2
17852b730f8bSJeremy L Thompson   for (CeedInt i = 0; i < n; i++) vec_D[i] = 1. / sqrt(vec_D[i]);
1786d3331725Sjeremylt   // -- G = G D^-1/2
1787d3331725Sjeremylt   // -- C = D^-1/2 G^T
17882b730f8bSJeremy L Thompson   for (CeedInt i = 0; i < n; i++) {
1789d3331725Sjeremylt     for (CeedInt j = 0; j < n; j++) {
1790673160d7Sjeremylt       mat_G[i * n + j] *= vec_D[j];
1791673160d7Sjeremylt       mat_C[j * n + i] = mat_G[i * n + j];
1792d3331725Sjeremylt     }
17932b730f8bSJeremy L Thompson   }
1794673160d7Sjeremylt   // -- X = (D^-1/2 G^T) A
17952b730f8bSJeremy L Thompson   CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_C, (const CeedScalar *)mat_A, mat_X, n, n, n));
1796673160d7Sjeremylt   // -- C = (D^-1/2 G^T A) (G D^-1/2)
17972b730f8bSJeremy L Thompson   CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_X, (const CeedScalar *)mat_G, mat_C, n, n, n));
179852bfb9bbSJeremy L Thompson 
179952bfb9bbSJeremy L Thompson   // Compute Q^T C Q = lambda
18002b730f8bSJeremy L Thompson   CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_C, lambda, n));
180152bfb9bbSJeremy L Thompson 
180285cf89eaSjeremylt   // Sort eigenvalues
18032b730f8bSJeremy L Thompson   for (CeedInt i = n - 1; i >= 0; i--) {
180485cf89eaSjeremylt     for (CeedInt j = 0; j < i; j++) {
180585cf89eaSjeremylt       if (fabs(lambda[j]) > fabs(lambda[j + 1])) {
180685cf89eaSjeremylt         CeedScalar temp;
18072b730f8bSJeremy L Thompson         temp          = lambda[j];
18082b730f8bSJeremy L Thompson         lambda[j]     = lambda[j + 1];
18092b730f8bSJeremy L Thompson         lambda[j + 1] = temp;
181085cf89eaSjeremylt         for (CeedInt k = 0; k < n; k++) {
18112b730f8bSJeremy L Thompson           temp                 = mat_C[k * n + j];
18122b730f8bSJeremy L Thompson           mat_C[k * n + j]     = mat_C[k * n + j + 1];
18132b730f8bSJeremy L Thompson           mat_C[k * n + j + 1] = temp;
18142b730f8bSJeremy L Thompson         }
181585cf89eaSjeremylt       }
181685cf89eaSjeremylt     }
181785cf89eaSjeremylt   }
181885cf89eaSjeremylt 
1819d3331725Sjeremylt   // Set X = (G D^1/2)^-T Q
1820fb551037Sjeremylt   //       = G D^-1/2 Q
18212b730f8bSJeremy L Thompson   CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_G, (const CeedScalar *)mat_C, mat_X, n, n, n));
182278464608Sjeremylt 
182378464608Sjeremylt   // Cleanup
18242b730f8bSJeremy L Thompson   CeedCall(CeedFree(&mat_C));
18252b730f8bSJeremy L Thompson   CeedCall(CeedFree(&mat_G));
18262b730f8bSJeremy L Thompson   CeedCall(CeedFree(&vec_D));
1827e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
182852bfb9bbSJeremy L Thompson }
182903d18186Sjeremylt CeedPragmaOptimizeOn
183052bfb9bbSJeremy L Thompson 
1831d7b241e6Sjeremylt     /// @}
1832