xref: /libCEED/interface/ceed-basis.c (revision 15ad391736ade9c8ee9580dd2a112afd9b87264b)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3d7b241e6Sjeremylt //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5d7b241e6Sjeremylt //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7d7b241e6Sjeremylt 
83d576824SJeremy L Thompson #include <ceed-impl.h>
949aac155SJeremy L Thompson #include <ceed.h>
102b730f8bSJeremy L Thompson #include <ceed/backend.h>
11d7b241e6Sjeremylt #include <math.h>
123d576824SJeremy L Thompson #include <stdbool.h>
13d7b241e6Sjeremylt #include <stdio.h>
14d7b241e6Sjeremylt #include <string.h>
15d7b241e6Sjeremylt 
167a982d89SJeremy L. Thompson /// @file
177a982d89SJeremy L. Thompson /// Implementation of CeedBasis interfaces
187a982d89SJeremy L. Thompson 
19d7b241e6Sjeremylt /// @cond DOXYGEN_SKIP
20783c99b3SValeria Barra static struct CeedBasis_private ceed_basis_collocated;
21d7b241e6Sjeremylt /// @endcond
22d7b241e6Sjeremylt 
237a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser
247a982d89SJeremy L. Thompson /// @{
257a982d89SJeremy L. Thompson 
267a982d89SJeremy L. Thompson /// Indicate that the quadrature points are collocated with the nodes
277a982d89SJeremy L. Thompson const CeedBasis CEED_BASIS_COLLOCATED = &ceed_basis_collocated;
287a982d89SJeremy L. Thompson 
297a982d89SJeremy L. Thompson /// @}
307a982d89SJeremy L. Thompson 
317a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
327a982d89SJeremy L. Thompson /// CeedBasis Library Internal Functions
337a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
347a982d89SJeremy L. Thompson /// @addtogroup CeedBasisDeveloper
357a982d89SJeremy L. Thompson /// @{
367a982d89SJeremy L. Thompson 
377a982d89SJeremy L. Thompson /**
387a982d89SJeremy L. Thompson   @brief Compute Householder reflection
397a982d89SJeremy L. Thompson 
40ea61e9acSJeremy 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
43ea61e9acSJeremy L Thompson   @param[in]     v   Householder vector
44ea61e9acSJeremy L Thompson   @param[in]     b   Scaling factor
45ea61e9acSJeremy L Thompson   @param[in]     m   Number of rows in A
46ea61e9acSJeremy L Thompson   @param[in]     n   Number of columns in A
47ea61e9acSJeremy L Thompson   @param[in]     row Row stride
48ea61e9acSJeremy 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 Compute Givens rotation
667a982d89SJeremy L. Thompson 
67ea61e9acSJeremy 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]
687a982d89SJeremy L. Thompson 
697a982d89SJeremy L. Thompson   @param[in,out] A      Row major matrix to apply Givens rotation to, in place
70ea61e9acSJeremy L Thompson   @param[in]     c      Cosine factor
71ea61e9acSJeremy L Thompson   @param[in]     s      Sine factor
72ea61e9acSJeremy 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;
734cc79fe7SJed Brown                           @ref CEED_TRANSPOSE for the opposite rotation
74ea61e9acSJeremy L Thompson   @param[in]     i      First row/column to apply rotation
75ea61e9acSJeremy L Thompson   @param[in]     k      Second row/column to apply rotation
76ea61e9acSJeremy L Thompson   @param[in]     m      Number of rows in A
77ea61e9acSJeremy L Thompson   @param[in]     n      Number of columns 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 static int CeedGivensRotation(CeedScalar *A, CeedScalar c, CeedScalar s, CeedTransposeMode t_mode, CeedInt i, CeedInt k, CeedInt m, CeedInt n) {
84d1d35e2fSjeremylt   CeedInt stride_j = 1, stride_ik = m, num_its = n;
85d1d35e2fSjeremylt   if (t_mode == CEED_NOTRANSPOSE) {
862b730f8bSJeremy L Thompson     stride_j  = n;
872b730f8bSJeremy L Thompson     stride_ik = 1;
882b730f8bSJeremy L Thompson     num_its   = m;
897a982d89SJeremy L. Thompson   }
907a982d89SJeremy L. Thompson 
917a982d89SJeremy L. Thompson   // Apply rotation
92d1d35e2fSjeremylt   for (CeedInt j = 0; j < num_its; j++) {
93d1d35e2fSjeremylt     CeedScalar tau1 = A[i * stride_ik + j * stride_j], tau2 = A[k * stride_ik + j * stride_j];
94d1d35e2fSjeremylt     A[i * stride_ik + j * stride_j] = c * tau1 - s * tau2;
95d1d35e2fSjeremylt     A[k * stride_ik + j * stride_j] = s * tau1 + c * tau2;
967a982d89SJeremy L. Thompson   }
97e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
987a982d89SJeremy L. Thompson }
997a982d89SJeremy L. Thompson 
1007a982d89SJeremy L. Thompson /**
1017a982d89SJeremy L. Thompson   @brief View an array stored in a CeedBasis
1027a982d89SJeremy L. Thompson 
1030a0da059Sjeremylt   @param[in] name   Name of array
104d1d35e2fSjeremylt   @param[in] fp_fmt Printing format
1050a0da059Sjeremylt   @param[in] m      Number of rows in array
1060a0da059Sjeremylt   @param[in] n      Number of columns in array
1070a0da059Sjeremylt   @param[in] a      Array to be viewed
1080a0da059Sjeremylt   @param[in] stream Stream to view to, e.g., stdout
1097a982d89SJeremy L. Thompson 
1107a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
1117a982d89SJeremy L. Thompson 
1127a982d89SJeremy L. Thompson   @ref Developer
1137a982d89SJeremy L. Thompson **/
1142b730f8bSJeremy L Thompson static int CeedScalarView(const char *name, const char *fp_fmt, CeedInt m, CeedInt n, const CeedScalar *a, FILE *stream) {
11592ae7e47SJeremy L Thompson   for (CeedInt i = 0; i < m; i++) {
1162b730f8bSJeremy L Thompson     if (m > 1) fprintf(stream, "%12s[%" CeedInt_FMT "]:", name, i);
1172b730f8bSJeremy L Thompson     else fprintf(stream, "%12s:", name);
1182b730f8bSJeremy 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);
1197a982d89SJeremy L. Thompson     fputs("\n", stream);
1207a982d89SJeremy L. Thompson   }
121e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1227a982d89SJeremy L. Thompson }
1237a982d89SJeremy L. Thompson 
124a76a04e7SJeremy L Thompson /**
125ea61e9acSJeremy L Thompson   @brief Create the interpolation and gradient matrices for projection from the nodes of `basis_from` to the nodes of `basis_to`.
126ba59ac12SSebastian Grimberg 
127*15ad3917SSebastian Grimberg   The interpolation is given by `interp_project = interp_to^+ * interp_from`, where the pseudoinverse `interp_to^+` is given by QR factorization.
128*15ad3917SSebastian Grimberg   The gradient is given by `grad_project = interp_to^+ * grad_from`, and is only computed for H^1 spaces otherwise it should not be used.
129*15ad3917SSebastian Grimberg 
130ba59ac12SSebastian Grimberg   Note: `basis_from` and `basis_to` must have compatible quadrature spaces.
131a76a04e7SJeremy L Thompson 
132a76a04e7SJeremy L Thompson   @param[in]  basis_from     CeedBasis to project from
133a76a04e7SJeremy L Thompson   @param[in]  basis_to       CeedBasis to project to
134ea61e9acSJeremy L Thompson   @param[out] interp_project Address of the variable where the newly created interpolation matrix will be stored.
135ea61e9acSJeremy L Thompson   @param[out] grad_project   Address of the variable where the newly created gradient matrix will be stored.
136a76a04e7SJeremy L Thompson 
137a76a04e7SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
138a76a04e7SJeremy L Thompson 
139a76a04e7SJeremy L Thompson   @ref Developer
140a76a04e7SJeremy L Thompson **/
1412b730f8bSJeremy L Thompson static int CeedBasisCreateProjectionMatrices(CeedBasis basis_from, CeedBasis basis_to, CeedScalar **interp_project, CeedScalar **grad_project) {
142a76a04e7SJeremy L Thompson   Ceed ceed;
1432b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetCeed(basis_to, &ceed));
144a76a04e7SJeremy L Thompson 
145a76a04e7SJeremy L Thompson   // Check for compatible quadrature spaces
146a76a04e7SJeremy L Thompson   CeedInt Q_to, Q_from;
1472b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumQuadraturePoints(basis_to, &Q_to));
1482b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumQuadraturePoints(basis_from, &Q_from));
1492b730f8bSJeremy L Thompson   if (Q_to != Q_from) {
150a76a04e7SJeremy L Thompson     // LCOV_EXCL_START
1512b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Bases must have compatible quadrature spaces");
152a76a04e7SJeremy L Thompson     // LCOV_EXCL_STOP
1532b730f8bSJeremy L Thompson   }
154a76a04e7SJeremy L Thompson 
15514556e63SJeremy L Thompson   // Check for matching tensor or non-tensor
156a76a04e7SJeremy L Thompson   CeedInt P_to, P_from, Q = Q_to;
157a76a04e7SJeremy L Thompson   bool    is_tensor_to, is_tensor_from;
1582b730f8bSJeremy L Thompson   CeedCall(CeedBasisIsTensor(basis_to, &is_tensor_to));
1592b730f8bSJeremy L Thompson   CeedCall(CeedBasisIsTensor(basis_from, &is_tensor_from));
160a76a04e7SJeremy L Thompson   if (is_tensor_to && is_tensor_from) {
1612b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_to));
1622b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_from));
1632b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumQuadraturePoints1D(basis_from, &Q));
164a76a04e7SJeremy L Thompson   } else if (!is_tensor_to && !is_tensor_from) {
1652b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis_to, &P_to));
1662b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis_from, &P_from));
167a76a04e7SJeremy L Thompson   } else {
168a76a04e7SJeremy L Thompson     // LCOV_EXCL_START
1692b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_MINOR, "Bases must both be tensor or non-tensor");
170a76a04e7SJeremy L Thompson     // LCOV_EXCL_STOP
171a76a04e7SJeremy L Thompson   }
172a76a04e7SJeremy L Thompson 
173*15ad3917SSebastian Grimberg   // Check for matching FE space
174*15ad3917SSebastian Grimberg   CeedFESpace fe_space_to, fe_space_from;
175*15ad3917SSebastian Grimberg   CeedCall(CeedBasisGetFESpace(basis_to, &fe_space_to));
176*15ad3917SSebastian Grimberg   CeedCall(CeedBasisGetFESpace(basis_from, &fe_space_from));
177*15ad3917SSebastian Grimberg   if (fe_space_to != fe_space_from) {
178*15ad3917SSebastian Grimberg     // LCOV_EXCL_START
179*15ad3917SSebastian Grimberg     return CeedError(ceed, CEED_ERROR_MINOR, "Bases must both be the same FE space type");
180*15ad3917SSebastian Grimberg     // LCOV_EXCL_STOP
181*15ad3917SSebastian Grimberg   }
182*15ad3917SSebastian Grimberg 
18314556e63SJeremy L Thompson   // Get source matrices
184*15ad3917SSebastian Grimberg   CeedInt           dim, q_comp = 1;
185*15ad3917SSebastian Grimberg   const CeedScalar *interp_to_source = NULL, *interp_from_source = NULL;
18614556e63SJeremy L Thompson   CeedScalar       *interp_to, *interp_from, *tau;
1872b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetDimension(basis_to, &dim));
188a76a04e7SJeremy L Thompson   if (is_tensor_to) {
1892b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetInterp1D(basis_to, &interp_to_source));
1902b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetInterp1D(basis_from, &interp_from_source));
191a76a04e7SJeremy L Thompson   } else {
192*15ad3917SSebastian Grimberg     CeedCall(CeedBasisGetNumQuadratureComponents(basis_from, CEED_EVAL_INTERP, &q_comp));
1932b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetInterp(basis_to, &interp_to_source));
1942b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetInterp(basis_from, &interp_from_source));
195*15ad3917SSebastian Grimberg   }
196*15ad3917SSebastian Grimberg   CeedCall(CeedMalloc(Q * P_from * q_comp, &interp_from));
197*15ad3917SSebastian Grimberg   CeedCall(CeedMalloc(Q * P_to * q_comp, &interp_to));
198*15ad3917SSebastian Grimberg   CeedCall(CeedCalloc(P_to * P_from, interp_project));
199*15ad3917SSebastian Grimberg   CeedCall(CeedMalloc(Q * q_comp, &tau));
200*15ad3917SSebastian Grimberg 
201*15ad3917SSebastian Grimberg   // `grad_project = interp_to^+ * grad_from` is computed for the H^1 space case so the
202*15ad3917SSebastian Grimberg   // projection basis will have a gradient operation
203*15ad3917SSebastian Grimberg   const CeedScalar *grad_from_source = NULL;
204*15ad3917SSebastian Grimberg   if (fe_space_to == CEED_FE_SPACE_H1) {
205*15ad3917SSebastian Grimberg     if (is_tensor_to) {
206*15ad3917SSebastian Grimberg       CeedCall(CeedBasisGetGrad1D(basis_from, &grad_from_source));
207*15ad3917SSebastian Grimberg     } else {
2082b730f8bSJeremy L Thompson       CeedCall(CeedBasisGetGrad(basis_from, &grad_from_source));
209a76a04e7SJeremy L Thompson     }
210*15ad3917SSebastian Grimberg     CeedCall(CeedCalloc(P_to * P_from * (is_tensor_to ? 1 : dim), grad_project));
211*15ad3917SSebastian Grimberg   } else {
212*15ad3917SSebastian Grimberg     // Projected derivate operator is only calculated for H^1 but we allocate it for the
213*15ad3917SSebastian Grimberg     // basis construction later
214*15ad3917SSebastian Grimberg     CeedInt q_comp_deriv = 1;
215*15ad3917SSebastian Grimberg     if (fe_space_to == CEED_FE_SPACE_HDIV) {
216*15ad3917SSebastian Grimberg       CeedCall(CeedBasisGetNumQuadratureComponents(basis_from, CEED_EVAL_DIV, &q_comp_deriv));
217*15ad3917SSebastian Grimberg     } else if (fe_space_to == CEED_FE_SPACE_HCURL) {
218*15ad3917SSebastian Grimberg       CeedCall(CeedBasisGetNumQuadratureComponents(basis_from, CEED_EVAL_CURL, &q_comp_deriv));
219*15ad3917SSebastian Grimberg     }
220*15ad3917SSebastian Grimberg     CeedCall(CeedCalloc(P_to * P_from * q_comp_deriv, grad_project));
221*15ad3917SSebastian Grimberg   }
222*15ad3917SSebastian Grimberg 
223*15ad3917SSebastian Grimberg   // QR Factorization, interp_to = Q R
224*15ad3917SSebastian Grimberg   memcpy(interp_to, interp_to_source, Q * P_to * q_comp * sizeof(interp_to_source[0]));
225*15ad3917SSebastian Grimberg   CeedCall(CeedQRFactorization(ceed, interp_to, tau, Q * q_comp, P_to));
226a76a04e7SJeremy L Thompson 
22714556e63SJeremy L Thompson   // Build matrices
228*15ad3917SSebastian Grimberg   CeedInt     num_matrices = 1 + (fe_space_to == CEED_FE_SPACE_H1) * (is_tensor_to ? 1 : dim);
22914556e63SJeremy L Thompson   CeedScalar *input_from[num_matrices], *output_project[num_matrices];
23014556e63SJeremy L Thompson   input_from[0]     = (CeedScalar *)interp_from_source;
23114556e63SJeremy L Thompson   output_project[0] = *interp_project;
23214556e63SJeremy L Thompson   for (CeedInt m = 1; m < num_matrices; m++) {
23314556e63SJeremy L Thompson     input_from[m]     = (CeedScalar *)&grad_from_source[(m - 1) * Q * P_from];
23402af4036SJeremy L Thompson     output_project[m] = &((*grad_project)[(m - 1) * P_to * P_from]);
23514556e63SJeremy L Thompson   }
23614556e63SJeremy L Thompson   for (CeedInt m = 0; m < num_matrices; m++) {
237*15ad3917SSebastian Grimberg     // Apply Q^T, interp_from = Q^T interp_from
238*15ad3917SSebastian Grimberg     memcpy(interp_from, input_from[m], Q * P_from * q_comp * sizeof(input_from[m][0]));
239*15ad3917SSebastian Grimberg     CeedCall(CeedHouseholderApplyQ(interp_from, interp_to, tau, CEED_TRANSPOSE, Q * q_comp, P_from, P_to, P_from, 1));
240a76a04e7SJeremy L Thompson 
241*15ad3917SSebastian Grimberg     // Apply Rinv, output_project = Rinv interp_from
242a76a04e7SJeremy L Thompson     for (CeedInt j = 0; j < P_from; j++) {  // Column j
2432b730f8bSJeremy 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];
244a76a04e7SJeremy L Thompson       for (CeedInt i = P_to - 2; i >= 0; i--) {  // Row i
24514556e63SJeremy L Thompson         output_project[m][j + P_from * i] = interp_from[j + P_from * i];
246a76a04e7SJeremy L Thompson         for (CeedInt k = i + 1; k < P_to; k++) {
2472b730f8bSJeremy L Thompson           output_project[m][j + P_from * i] -= interp_to[k + P_to * i] * output_project[m][j + P_from * k];
248a76a04e7SJeremy L Thompson         }
24914556e63SJeremy L Thompson         output_project[m][j + P_from * i] /= interp_to[i + P_to * i];
250a76a04e7SJeremy L Thompson       }
251a76a04e7SJeremy L Thompson     }
25214556e63SJeremy L Thompson   }
25314556e63SJeremy L Thompson 
25414556e63SJeremy L Thompson   // Cleanup
2552b730f8bSJeremy L Thompson   CeedCall(CeedFree(&tau));
2562b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_to));
2572b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_from));
258a76a04e7SJeremy L Thompson 
259a76a04e7SJeremy L Thompson   return CEED_ERROR_SUCCESS;
260a76a04e7SJeremy L Thompson }
261a76a04e7SJeremy L Thompson 
2627a982d89SJeremy L. Thompson /// @}
2637a982d89SJeremy L. Thompson 
2647a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
2657a982d89SJeremy L. Thompson /// Ceed Backend API
2667a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
2677a982d89SJeremy L. Thompson /// @addtogroup CeedBasisBackend
2687a982d89SJeremy L. Thompson /// @{
2697a982d89SJeremy L. Thompson 
2707a982d89SJeremy L. Thompson /**
2717a982d89SJeremy L. Thompson   @brief Return collocated grad matrix
2727a982d89SJeremy L. Thompson 
273ea61e9acSJeremy L Thompson   @param[in]  basis         CeedBasis
274ea61e9acSJeremy L Thompson   @param[out] collo_grad_1d Row-major (Q_1d * Q_1d) matrix expressing derivatives of basis functions at quadrature points
2757a982d89SJeremy L. Thompson 
2767a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
2777a982d89SJeremy L. Thompson 
2787a982d89SJeremy L. Thompson   @ref Backend
2797a982d89SJeremy L. Thompson **/
280d1d35e2fSjeremylt int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *collo_grad_1d) {
2817a982d89SJeremy L. Thompson   int         i, j, k;
2827a982d89SJeremy L. Thompson   Ceed        ceed;
2832b730f8bSJeremy L Thompson   CeedInt     P_1d = (basis)->P_1d, Q_1d = (basis)->Q_1d;
28478464608Sjeremylt   CeedScalar *interp_1d, *grad_1d, *tau;
2857a982d89SJeremy L. Thompson 
2862b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q_1d * P_1d, &interp_1d));
2872b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q_1d * P_1d, &grad_1d));
2882b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q_1d, &tau));
289d1d35e2fSjeremylt   memcpy(interp_1d, (basis)->interp_1d, Q_1d * P_1d * sizeof(basis)->interp_1d[0]);
290d1d35e2fSjeremylt   memcpy(grad_1d, (basis)->grad_1d, Q_1d * P_1d * sizeof(basis)->interp_1d[0]);
2917a982d89SJeremy L. Thompson 
292d1d35e2fSjeremylt   // QR Factorization, interp_1d = Q R
2932b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetCeed(basis, &ceed));
2942b730f8bSJeremy L Thompson   CeedCall(CeedQRFactorization(ceed, interp_1d, tau, Q_1d, P_1d));
295ea61e9acSJeremy 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.
2967a982d89SJeremy L. Thompson 
297d1d35e2fSjeremylt   // Apply Rinv, collo_grad_1d = grad_1d Rinv
298d1d35e2fSjeremylt   for (i = 0; i < Q_1d; i++) {  // Row i
299d1d35e2fSjeremylt     collo_grad_1d[Q_1d * i] = grad_1d[P_1d * i] / interp_1d[0];
300d1d35e2fSjeremylt     for (j = 1; j < P_1d; j++) {  // Column j
301d1d35e2fSjeremylt       collo_grad_1d[j + Q_1d * i] = grad_1d[j + P_1d * i];
3022b730f8bSJeremy 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];
303d1d35e2fSjeremylt       collo_grad_1d[j + Q_1d * i] /= interp_1d[j + P_1d * j];
3047a982d89SJeremy L. Thompson     }
3052b730f8bSJeremy L Thompson     for (j = P_1d; j < Q_1d; j++) collo_grad_1d[j + Q_1d * i] = 0;
3067a982d89SJeremy L. Thompson   }
3077a982d89SJeremy L. Thompson 
308*15ad3917SSebastian Grimberg   // Apply Q^T, collo_grad_1d = collo_grad_1d Q^T
3092b730f8bSJeremy L Thompson   CeedCall(CeedHouseholderApplyQ(collo_grad_1d, interp_1d, tau, CEED_NOTRANSPOSE, Q_1d, Q_1d, P_1d, 1, Q_1d));
3107a982d89SJeremy L. Thompson 
3112b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_1d));
3122b730f8bSJeremy L Thompson   CeedCall(CeedFree(&grad_1d));
3132b730f8bSJeremy L Thompson   CeedCall(CeedFree(&tau));
314e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3157a982d89SJeremy L. Thompson }
3167a982d89SJeremy L. Thompson 
3177a982d89SJeremy L. Thompson /**
3187a982d89SJeremy L. Thompson   @brief Get tensor status for given CeedBasis
3197a982d89SJeremy L. Thompson 
320ea61e9acSJeremy L Thompson   @param[in]  basis     CeedBasis
321d1d35e2fSjeremylt   @param[out] is_tensor Variable to store tensor status
3227a982d89SJeremy L. Thompson 
3237a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3247a982d89SJeremy L. Thompson 
3257a982d89SJeremy L. Thompson   @ref Backend
3267a982d89SJeremy L. Thompson **/
327d1d35e2fSjeremylt int CeedBasisIsTensor(CeedBasis basis, bool *is_tensor) {
328d1d35e2fSjeremylt   *is_tensor = basis->tensor_basis;
329e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3307a982d89SJeremy L. Thompson }
3317a982d89SJeremy L. Thompson 
3327a982d89SJeremy L. Thompson /**
3337a982d89SJeremy L. Thompson   @brief Get backend data of a CeedBasis
3347a982d89SJeremy L. Thompson 
335ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
3367a982d89SJeremy L. Thompson   @param[out] data  Variable to store data
3377a982d89SJeremy L. Thompson 
3387a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3397a982d89SJeremy L. Thompson 
3407a982d89SJeremy L. Thompson   @ref Backend
3417a982d89SJeremy L. Thompson **/
342777ff853SJeremy L Thompson int CeedBasisGetData(CeedBasis basis, void *data) {
343777ff853SJeremy L Thompson   *(void **)data = basis->data;
344e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3457a982d89SJeremy L. Thompson }
3467a982d89SJeremy L. Thompson 
3477a982d89SJeremy L. Thompson /**
3487a982d89SJeremy L. Thompson   @brief Set backend data of a CeedBasis
3497a982d89SJeremy L. Thompson 
350ea61e9acSJeremy L Thompson   @param[in,out] basis  CeedBasis
351ea61e9acSJeremy L Thompson   @param[in]     data   Data to set
3527a982d89SJeremy L. Thompson 
3537a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3547a982d89SJeremy L. Thompson 
3557a982d89SJeremy L. Thompson   @ref Backend
3567a982d89SJeremy L. Thompson **/
357777ff853SJeremy L Thompson int CeedBasisSetData(CeedBasis basis, void *data) {
358777ff853SJeremy L Thompson   basis->data = data;
359e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3607a982d89SJeremy L. Thompson }
3617a982d89SJeremy L. Thompson 
3627a982d89SJeremy L. Thompson /**
36334359f16Sjeremylt   @brief Increment the reference counter for a CeedBasis
36434359f16Sjeremylt 
365ea61e9acSJeremy L Thompson   @param[in,out] basis Basis to increment the reference counter
36634359f16Sjeremylt 
36734359f16Sjeremylt   @return An error code: 0 - success, otherwise - failure
36834359f16Sjeremylt 
36934359f16Sjeremylt   @ref Backend
37034359f16Sjeremylt **/
3719560d06aSjeremylt int CeedBasisReference(CeedBasis basis) {
37234359f16Sjeremylt   basis->ref_count++;
37334359f16Sjeremylt   return CEED_ERROR_SUCCESS;
37434359f16Sjeremylt }
37534359f16Sjeremylt 
37634359f16Sjeremylt /**
377c4e3f59bSSebastian Grimberg   @brief Get number of Q-vector components for given CeedBasis
378c4e3f59bSSebastian Grimberg 
379c4e3f59bSSebastian Grimberg   @param[in]  basis  CeedBasis
380c4e3f59bSSebastian Grimberg   @param[in]  eval_mode \ref CEED_EVAL_INTERP to use interpolated values,
381c4e3f59bSSebastian Grimberg                           \ref CEED_EVAL_GRAD to use gradients,
382c4e3f59bSSebastian Grimberg                           \ref CEED_EVAL_DIV to use divergence,
383c4e3f59bSSebastian Grimberg                           \ref CEED_EVAL_CURL to use curl.
384c4e3f59bSSebastian Grimberg   @param[out] q_comp Variable to store number of Q-vector components of basis
385c4e3f59bSSebastian Grimberg 
386c4e3f59bSSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
387c4e3f59bSSebastian Grimberg 
388c4e3f59bSSebastian Grimberg   @ref Backend
389c4e3f59bSSebastian Grimberg **/
390c4e3f59bSSebastian Grimberg int CeedBasisGetNumQuadratureComponents(CeedBasis basis, CeedEvalMode eval_mode, CeedInt *q_comp) {
391c4e3f59bSSebastian Grimberg   switch (eval_mode) {
392c4e3f59bSSebastian Grimberg     case CEED_EVAL_INTERP:
393c4e3f59bSSebastian Grimberg       *q_comp = (basis->fe_space == CEED_FE_SPACE_H1) ? 1 : basis->dim;
394c4e3f59bSSebastian Grimberg       break;
395c4e3f59bSSebastian Grimberg     case CEED_EVAL_GRAD:
396c4e3f59bSSebastian Grimberg       *q_comp = basis->dim;
397c4e3f59bSSebastian Grimberg       break;
398c4e3f59bSSebastian Grimberg     case CEED_EVAL_DIV:
399c4e3f59bSSebastian Grimberg       *q_comp = 1;
400c4e3f59bSSebastian Grimberg       break;
401c4e3f59bSSebastian Grimberg     case CEED_EVAL_CURL:
402c4e3f59bSSebastian Grimberg       *q_comp = (basis->dim < 3) ? 1 : basis->dim;
403c4e3f59bSSebastian Grimberg       break;
404c4e3f59bSSebastian Grimberg     case CEED_EVAL_NONE:
405c4e3f59bSSebastian Grimberg     case CEED_EVAL_WEIGHT:
406c4e3f59bSSebastian Grimberg       *q_comp = 0;
407c4e3f59bSSebastian Grimberg       break;
408c4e3f59bSSebastian Grimberg   }
409c4e3f59bSSebastian Grimberg   return CEED_ERROR_SUCCESS;
410c4e3f59bSSebastian Grimberg }
411c4e3f59bSSebastian Grimberg 
412c4e3f59bSSebastian Grimberg /**
4136e15d496SJeremy L Thompson   @brief Estimate number of FLOPs required to apply CeedBasis in t_mode and eval_mode
4146e15d496SJeremy L Thompson 
415ea61e9acSJeremy L Thompson   @param[in]  basis     Basis to estimate FLOPs for
416ea61e9acSJeremy L Thompson   @param[in]  t_mode    Apply basis or transpose
417ea61e9acSJeremy L Thompson   @param[in]  eval_mode Basis evaluation mode
418ea61e9acSJeremy L Thompson   @param[out] flops     Address of variable to hold FLOPs estimate
4196e15d496SJeremy L Thompson 
4206e15d496SJeremy L Thompson   @ref Backend
4216e15d496SJeremy L Thompson **/
4222b730f8bSJeremy L Thompson int CeedBasisGetFlopsEstimate(CeedBasis basis, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedSize *flops) {
4236e15d496SJeremy L Thompson   bool is_tensor;
4246e15d496SJeremy L Thompson 
4252b730f8bSJeremy L Thompson   CeedCall(CeedBasisIsTensor(basis, &is_tensor));
4266e15d496SJeremy L Thompson   if (is_tensor) {
4276e15d496SJeremy L Thompson     CeedInt dim, num_comp, P_1d, Q_1d;
4282b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetDimension(basis, &dim));
4292b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
4302b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d));
4312b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
4326e15d496SJeremy L Thompson     if (t_mode == CEED_TRANSPOSE) {
4332b730f8bSJeremy L Thompson       P_1d = Q_1d;
4342b730f8bSJeremy L Thompson       Q_1d = P_1d;
4356e15d496SJeremy L Thompson     }
4366e15d496SJeremy L Thompson     CeedInt tensor_flops = 0, pre = num_comp * CeedIntPow(P_1d, dim - 1), post = 1;
4376e15d496SJeremy L Thompson     for (CeedInt d = 0; d < dim; d++) {
4386e15d496SJeremy L Thompson       tensor_flops += 2 * pre * P_1d * post * Q_1d;
4396e15d496SJeremy L Thompson       pre /= P_1d;
4406e15d496SJeremy L Thompson       post *= Q_1d;
4416e15d496SJeremy L Thompson     }
4426e15d496SJeremy L Thompson     switch (eval_mode) {
4432b730f8bSJeremy L Thompson       case CEED_EVAL_NONE:
4442b730f8bSJeremy L Thompson         *flops = 0;
4452b730f8bSJeremy L Thompson         break;
4462b730f8bSJeremy L Thompson       case CEED_EVAL_INTERP:
4472b730f8bSJeremy L Thompson         *flops = tensor_flops;
4482b730f8bSJeremy L Thompson         break;
4492b730f8bSJeremy L Thompson       case CEED_EVAL_GRAD:
4502b730f8bSJeremy L Thompson         *flops = tensor_flops * 2;
4512b730f8bSJeremy L Thompson         break;
4526e15d496SJeremy L Thompson       case CEED_EVAL_DIV:
4536e15d496SJeremy L Thompson         // LCOV_EXCL_START
4542b730f8bSJeremy L Thompson         return CeedError(basis->ceed, CEED_ERROR_INCOMPATIBLE, "Tensor CEED_EVAL_DIV not supported");
4552b730f8bSJeremy L Thompson         break;
4566e15d496SJeremy L Thompson       case CEED_EVAL_CURL:
4572b730f8bSJeremy L Thompson         return CeedError(basis->ceed, CEED_ERROR_INCOMPATIBLE, "Tensor CEED_EVAL_CURL not supported");
4582b730f8bSJeremy L Thompson         break;
4596e15d496SJeremy L Thompson       // LCOV_EXCL_STOP
4602b730f8bSJeremy L Thompson       case CEED_EVAL_WEIGHT:
4612b730f8bSJeremy L Thompson         *flops = dim * CeedIntPow(Q_1d, dim);
4622b730f8bSJeremy L Thompson         break;
4636e15d496SJeremy L Thompson     }
4646e15d496SJeremy L Thompson   } else {
465c4e3f59bSSebastian Grimberg     CeedInt dim, num_comp, q_comp, num_nodes, num_qpts;
4662b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetDimension(basis, &dim));
4672b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
468c4e3f59bSSebastian Grimberg     CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp));
4692b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis, &num_nodes));
4702b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts));
4716e15d496SJeremy L Thompson     switch (eval_mode) {
4722b730f8bSJeremy L Thompson       case CEED_EVAL_NONE:
4732b730f8bSJeremy L Thompson         *flops = 0;
4742b730f8bSJeremy L Thompson         break;
4752b730f8bSJeremy L Thompson       case CEED_EVAL_INTERP:
4762b730f8bSJeremy L Thompson       case CEED_EVAL_GRAD:
4772b730f8bSJeremy L Thompson       case CEED_EVAL_DIV:
4782b730f8bSJeremy L Thompson       case CEED_EVAL_CURL:
479c4e3f59bSSebastian Grimberg         *flops = num_nodes * num_qpts * num_comp * q_comp;
4802b730f8bSJeremy L Thompson         break;
4812b730f8bSJeremy L Thompson       case CEED_EVAL_WEIGHT:
4822b730f8bSJeremy L Thompson         *flops = 0;
4832b730f8bSJeremy L Thompson         break;
4846e15d496SJeremy L Thompson     }
4856e15d496SJeremy L Thompson   }
4866e15d496SJeremy L Thompson 
4876e15d496SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4886e15d496SJeremy L Thompson }
4896e15d496SJeremy L Thompson 
4906e15d496SJeremy L Thompson /**
491c4e3f59bSSebastian Grimberg   @brief Get CeedFESpace for a CeedBasis
492c4e3f59bSSebastian Grimberg 
493c4e3f59bSSebastian Grimberg   @param[in]  basis     CeedBasis
494c4e3f59bSSebastian Grimberg   @param[out] fe_space  Variable to store CeedFESpace
495c4e3f59bSSebastian Grimberg 
496c4e3f59bSSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
497c4e3f59bSSebastian Grimberg 
498c4e3f59bSSebastian Grimberg   @ref Backend
499c4e3f59bSSebastian Grimberg **/
500c4e3f59bSSebastian Grimberg int CeedBasisGetFESpace(CeedBasis basis, CeedFESpace *fe_space) {
501c4e3f59bSSebastian Grimberg   *fe_space = basis->fe_space;
502c4e3f59bSSebastian Grimberg   return CEED_ERROR_SUCCESS;
503c4e3f59bSSebastian Grimberg }
504c4e3f59bSSebastian Grimberg 
505c4e3f59bSSebastian Grimberg /**
5067a982d89SJeremy L. Thompson   @brief Get dimension for given CeedElemTopology
5077a982d89SJeremy L. Thompson 
508ea61e9acSJeremy L Thompson   @param[in]  topo CeedElemTopology
5097a982d89SJeremy L. Thompson   @param[out] dim  Variable to store dimension of topology
5107a982d89SJeremy L. Thompson 
5117a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
5127a982d89SJeremy L. Thompson 
5137a982d89SJeremy L. Thompson   @ref Backend
5147a982d89SJeremy L. Thompson **/
5157a982d89SJeremy L. Thompson int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) {
5167a982d89SJeremy L. Thompson   *dim = (CeedInt)topo >> 16;
517e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5187a982d89SJeremy L. Thompson }
5197a982d89SJeremy L. Thompson 
5207a982d89SJeremy L. Thompson /**
5217a982d89SJeremy L. Thompson   @brief Get CeedTensorContract of a CeedBasis
5227a982d89SJeremy L. Thompson 
523ea61e9acSJeremy L Thompson   @param[in]  basis     CeedBasis
5247a982d89SJeremy L. Thompson   @param[out] contract  Variable to store CeedTensorContract
5257a982d89SJeremy L. Thompson 
5267a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
5277a982d89SJeremy L. Thompson 
5287a982d89SJeremy L. Thompson   @ref Backend
5297a982d89SJeremy L. Thompson **/
5307a982d89SJeremy L. Thompson int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) {
5317a982d89SJeremy L. Thompson   *contract = basis->contract;
532e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5337a982d89SJeremy L. Thompson }
5347a982d89SJeremy L. Thompson 
5357a982d89SJeremy L. Thompson /**
5367a982d89SJeremy L. Thompson   @brief Set CeedTensorContract of a CeedBasis
5377a982d89SJeremy L. Thompson 
538ea61e9acSJeremy L Thompson   @param[in,out] basis    CeedBasis
539ea61e9acSJeremy L Thompson   @param[in]     contract CeedTensorContract to set
5407a982d89SJeremy L. Thompson 
5417a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
5427a982d89SJeremy L. Thompson 
5437a982d89SJeremy L. Thompson   @ref Backend
5447a982d89SJeremy L. Thompson **/
54534359f16Sjeremylt int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract contract) {
54634359f16Sjeremylt   basis->contract = contract;
5472b730f8bSJeremy L Thompson   CeedCall(CeedTensorContractReference(contract));
548e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5497a982d89SJeremy L. Thompson }
5507a982d89SJeremy L. Thompson 
5517a982d89SJeremy L. Thompson /**
5527a982d89SJeremy L. Thompson   @brief Return a reference implementation of matrix multiplication C = A B.
553ba59ac12SSebastian Grimberg 
554ba59ac12SSebastian Grimberg   Note: This is a reference implementation for CPU CeedScalar pointers that is not intended for high performance.
5557a982d89SJeremy L. Thompson 
556ea61e9acSJeremy L Thompson   @param[in]  ceed  Ceed context for error handling
557d1d35e2fSjeremylt   @param[in]  mat_A Row-major matrix A
558d1d35e2fSjeremylt   @param[in]  mat_B Row-major matrix B
559d1d35e2fSjeremylt   @param[out] mat_C Row-major output matrix C
560ea61e9acSJeremy L Thompson   @param[in]  m     Number of rows of C
561ea61e9acSJeremy L Thompson   @param[in]  n     Number of columns of C
562ea61e9acSJeremy L Thompson   @param[in]  kk    Number of columns of A/rows of B
5637a982d89SJeremy L. Thompson 
5647a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
5657a982d89SJeremy L. Thompson 
5667a982d89SJeremy L. Thompson   @ref Utility
5677a982d89SJeremy L. Thompson **/
5682b730f8bSJeremy L Thompson int CeedMatrixMatrixMultiply(Ceed ceed, const CeedScalar *mat_A, const CeedScalar *mat_B, CeedScalar *mat_C, CeedInt m, CeedInt n, CeedInt kk) {
5692b730f8bSJeremy L Thompson   for (CeedInt i = 0; i < m; i++) {
5707a982d89SJeremy L. Thompson     for (CeedInt j = 0; j < n; j++) {
5717a982d89SJeremy L. Thompson       CeedScalar sum = 0;
5722b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < kk; k++) sum += mat_A[k + i * kk] * mat_B[j + k * n];
573d1d35e2fSjeremylt       mat_C[j + i * n] = sum;
5747a982d89SJeremy L. Thompson     }
5752b730f8bSJeremy L Thompson   }
576e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5777a982d89SJeremy L. Thompson }
5787a982d89SJeremy L. Thompson 
579ba59ac12SSebastian Grimberg /**
580ba59ac12SSebastian Grimberg   @brief Return QR Factorization of a matrix
581ba59ac12SSebastian Grimberg 
582ba59ac12SSebastian Grimberg   @param[in]     ceed Ceed context for error handling
583ba59ac12SSebastian Grimberg   @param[in,out] mat  Row-major matrix to be factorized in place
584ba59ac12SSebastian Grimberg   @param[in,out] tau  Vector of length m of scaling factors
585ba59ac12SSebastian Grimberg   @param[in]     m    Number of rows
586ba59ac12SSebastian Grimberg   @param[in]     n    Number of columns
587ba59ac12SSebastian Grimberg 
588ba59ac12SSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
589ba59ac12SSebastian Grimberg 
590ba59ac12SSebastian Grimberg   @ref Utility
591ba59ac12SSebastian Grimberg **/
592ba59ac12SSebastian Grimberg int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau, CeedInt m, CeedInt n) {
593ba59ac12SSebastian Grimberg   CeedScalar v[m];
594ba59ac12SSebastian Grimberg 
595ba59ac12SSebastian Grimberg   // Check matrix shape
596ba59ac12SSebastian Grimberg   if (n > m) {
597ba59ac12SSebastian Grimberg     // LCOV_EXCL_START
598ba59ac12SSebastian Grimberg     return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute QR factorization with n > m");
599ba59ac12SSebastian Grimberg     // LCOV_EXCL_STOP
600ba59ac12SSebastian Grimberg   }
601ba59ac12SSebastian Grimberg 
602ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) {
603ba59ac12SSebastian Grimberg     if (i >= m - 1) {  // last row of matrix, no reflection needed
604ba59ac12SSebastian Grimberg       tau[i] = 0.;
605ba59ac12SSebastian Grimberg       break;
606ba59ac12SSebastian Grimberg     }
607ba59ac12SSebastian Grimberg     // Calculate Householder vector, magnitude
608ba59ac12SSebastian Grimberg     CeedScalar sigma = 0.0;
609ba59ac12SSebastian Grimberg     v[i]             = mat[i + n * i];
610ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < m; j++) {
611ba59ac12SSebastian Grimberg       v[j] = mat[i + n * j];
612ba59ac12SSebastian Grimberg       sigma += v[j] * v[j];
613ba59ac12SSebastian Grimberg     }
614ba59ac12SSebastian Grimberg     CeedScalar norm = sqrt(v[i] * v[i] + sigma);  // norm of v[i:m]
615ba59ac12SSebastian Grimberg     CeedScalar R_ii = -copysign(norm, v[i]);
616ba59ac12SSebastian Grimberg     v[i] -= R_ii;
617ba59ac12SSebastian Grimberg     // norm of v[i:m] after modification above and scaling below
618ba59ac12SSebastian Grimberg     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
619ba59ac12SSebastian Grimberg     //   tau = 2 / (norm*norm)
620ba59ac12SSebastian Grimberg     tau[i] = 2 * v[i] * v[i] / (v[i] * v[i] + sigma);
621ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < m; j++) v[j] /= v[i];
622ba59ac12SSebastian Grimberg 
623ba59ac12SSebastian Grimberg     // Apply Householder reflector to lower right panel
624ba59ac12SSebastian Grimberg     CeedHouseholderReflect(&mat[i * n + i + 1], &v[i], tau[i], m - i, n - i - 1, n, 1);
625ba59ac12SSebastian Grimberg     // Save v
626ba59ac12SSebastian Grimberg     mat[i + n * i] = R_ii;
627ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < m; j++) mat[i + n * j] = v[j];
628ba59ac12SSebastian Grimberg   }
629ba59ac12SSebastian Grimberg   return CEED_ERROR_SUCCESS;
630ba59ac12SSebastian Grimberg }
631ba59ac12SSebastian Grimberg 
632ba59ac12SSebastian Grimberg /**
633ba59ac12SSebastian Grimberg   @brief Apply Householder Q matrix
634ba59ac12SSebastian Grimberg 
635ba59ac12SSebastian Grimberg   Compute mat_A = mat_Q mat_A, where mat_Q is mxm and mat_A is mxn.
636ba59ac12SSebastian Grimberg 
637ba59ac12SSebastian Grimberg   @param[in,out] mat_A  Matrix to apply Householder Q to, in place
638ba59ac12SSebastian Grimberg   @param[in]     mat_Q  Householder Q matrix
639ba59ac12SSebastian Grimberg   @param[in]     tau    Householder scaling factors
640ba59ac12SSebastian Grimberg   @param[in]     t_mode Transpose mode for application
641ba59ac12SSebastian Grimberg   @param[in]     m      Number of rows in A
642ba59ac12SSebastian Grimberg   @param[in]     n      Number of columns in A
643ba59ac12SSebastian Grimberg   @param[in]     k      Number of elementary reflectors in Q, k<m
644ba59ac12SSebastian Grimberg   @param[in]     row    Row stride in A
645ba59ac12SSebastian Grimberg   @param[in]     col    Col stride in A
646ba59ac12SSebastian Grimberg 
647ba59ac12SSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
648ba59ac12SSebastian Grimberg 
649c4e3f59bSSebastian Grimberg   @ref Utility
650ba59ac12SSebastian Grimberg **/
651ba59ac12SSebastian Grimberg int CeedHouseholderApplyQ(CeedScalar *mat_A, const CeedScalar *mat_Q, const CeedScalar *tau, CeedTransposeMode t_mode, CeedInt m, CeedInt n,
652ba59ac12SSebastian Grimberg                           CeedInt k, CeedInt row, CeedInt col) {
653ba59ac12SSebastian Grimberg   CeedScalar *v;
654ba59ac12SSebastian Grimberg   CeedCall(CeedMalloc(m, &v));
655ba59ac12SSebastian Grimberg   for (CeedInt ii = 0; ii < k; ii++) {
656ba59ac12SSebastian Grimberg     CeedInt i = t_mode == CEED_TRANSPOSE ? ii : k - 1 - ii;
657ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < m; j++) v[j] = mat_Q[j * k + i];
658ba59ac12SSebastian Grimberg     // Apply Householder reflector (I - tau v v^T) collo_grad_1d^T
659ba59ac12SSebastian Grimberg     CeedCall(CeedHouseholderReflect(&mat_A[i * row], &v[i], tau[i], m - i, n, row, col));
660ba59ac12SSebastian Grimberg   }
661ba59ac12SSebastian Grimberg   CeedCall(CeedFree(&v));
662ba59ac12SSebastian Grimberg   return CEED_ERROR_SUCCESS;
663ba59ac12SSebastian Grimberg }
664ba59ac12SSebastian Grimberg 
665ba59ac12SSebastian Grimberg /**
666ba59ac12SSebastian Grimberg   @brief Return symmetric Schur decomposition of the symmetric matrix mat via symmetric QR factorization
667ba59ac12SSebastian Grimberg 
668ba59ac12SSebastian Grimberg   @param[in]     ceed   Ceed context for error handling
669ba59ac12SSebastian Grimberg   @param[in,out] mat    Row-major matrix to be factorized in place
670ba59ac12SSebastian Grimberg   @param[out]    lambda Vector of length n of eigenvalues
671ba59ac12SSebastian Grimberg   @param[in]     n      Number of rows/columns
672ba59ac12SSebastian Grimberg 
673ba59ac12SSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
674ba59ac12SSebastian Grimberg 
675ba59ac12SSebastian Grimberg   @ref Utility
676ba59ac12SSebastian Grimberg **/
677ba59ac12SSebastian Grimberg CeedPragmaOptimizeOff int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat, CeedScalar *lambda, CeedInt n) {
678ba59ac12SSebastian Grimberg   // Check bounds for clang-tidy
679ba59ac12SSebastian Grimberg   if (n < 2) {
680ba59ac12SSebastian Grimberg     // LCOV_EXCL_START
681ba59ac12SSebastian Grimberg     return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute symmetric Schur decomposition of scalars");
682ba59ac12SSebastian Grimberg     // LCOV_EXCL_STOP
683ba59ac12SSebastian Grimberg   }
684ba59ac12SSebastian Grimberg 
685ba59ac12SSebastian Grimberg   CeedScalar v[n - 1], tau[n - 1], mat_T[n * n];
686ba59ac12SSebastian Grimberg 
687ba59ac12SSebastian Grimberg   // Copy mat to mat_T and set mat to I
688ba59ac12SSebastian Grimberg   memcpy(mat_T, mat, n * n * sizeof(mat[0]));
689ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) {
690ba59ac12SSebastian Grimberg     for (CeedInt j = 0; j < n; j++) mat[j + n * i] = (i == j) ? 1 : 0;
691ba59ac12SSebastian Grimberg   }
692ba59ac12SSebastian Grimberg 
693ba59ac12SSebastian Grimberg   // Reduce to tridiagonal
694ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n - 1; i++) {
695ba59ac12SSebastian Grimberg     // Calculate Householder vector, magnitude
696ba59ac12SSebastian Grimberg     CeedScalar sigma = 0.0;
697ba59ac12SSebastian Grimberg     v[i]             = mat_T[i + n * (i + 1)];
698ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < n - 1; j++) {
699ba59ac12SSebastian Grimberg       v[j] = mat_T[i + n * (j + 1)];
700ba59ac12SSebastian Grimberg       sigma += v[j] * v[j];
701ba59ac12SSebastian Grimberg     }
702ba59ac12SSebastian Grimberg     CeedScalar norm = sqrt(v[i] * v[i] + sigma);  // norm of v[i:n-1]
703ba59ac12SSebastian Grimberg     CeedScalar R_ii = -copysign(norm, v[i]);
704ba59ac12SSebastian Grimberg     v[i] -= R_ii;
705ba59ac12SSebastian Grimberg     // norm of v[i:m] after modification above and scaling below
706ba59ac12SSebastian Grimberg     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
707ba59ac12SSebastian Grimberg     //   tau = 2 / (norm*norm)
708ba59ac12SSebastian Grimberg     tau[i] = i == n - 2 ? 2 : 2 * v[i] * v[i] / (v[i] * v[i] + sigma);
709ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < n - 1; j++) v[j] /= v[i];
710ba59ac12SSebastian Grimberg 
711ba59ac12SSebastian Grimberg     // Update sub and super diagonal
712ba59ac12SSebastian Grimberg     for (CeedInt j = i + 2; j < n; j++) {
713ba59ac12SSebastian Grimberg       mat_T[i + n * j] = 0;
714ba59ac12SSebastian Grimberg       mat_T[j + n * i] = 0;
715ba59ac12SSebastian Grimberg     }
716ba59ac12SSebastian Grimberg     // Apply symmetric Householder reflector to lower right panel
717ba59ac12SSebastian Grimberg     CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1);
718ba59ac12SSebastian Grimberg     CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), 1, n);
719ba59ac12SSebastian Grimberg 
720ba59ac12SSebastian Grimberg     // Save v
721ba59ac12SSebastian Grimberg     mat_T[i + n * (i + 1)] = R_ii;
722ba59ac12SSebastian Grimberg     mat_T[(i + 1) + n * i] = R_ii;
723ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < n - 1; j++) {
724ba59ac12SSebastian Grimberg       mat_T[i + n * (j + 1)] = v[j];
725ba59ac12SSebastian Grimberg     }
726ba59ac12SSebastian Grimberg   }
727ba59ac12SSebastian Grimberg   // Backwards accumulation of Q
728ba59ac12SSebastian Grimberg   for (CeedInt i = n - 2; i >= 0; i--) {
729ba59ac12SSebastian Grimberg     if (tau[i] > 0.0) {
730ba59ac12SSebastian Grimberg       v[i] = 1;
731ba59ac12SSebastian Grimberg       for (CeedInt j = i + 1; j < n - 1; j++) {
732ba59ac12SSebastian Grimberg         v[j]                   = mat_T[i + n * (j + 1)];
733ba59ac12SSebastian Grimberg         mat_T[i + n * (j + 1)] = 0;
734ba59ac12SSebastian Grimberg       }
735ba59ac12SSebastian Grimberg       CeedHouseholderReflect(&mat[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1);
736ba59ac12SSebastian Grimberg     }
737ba59ac12SSebastian Grimberg   }
738ba59ac12SSebastian Grimberg 
739ba59ac12SSebastian Grimberg   // Reduce sub and super diagonal
740ba59ac12SSebastian Grimberg   CeedInt    p = 0, q = 0, itr = 0, max_itr = n * n * n * n;
741ba59ac12SSebastian Grimberg   CeedScalar tol = CEED_EPSILON;
742ba59ac12SSebastian Grimberg 
743ba59ac12SSebastian Grimberg   while (itr < max_itr) {
744ba59ac12SSebastian Grimberg     // Update p, q, size of reduced portions of diagonal
745ba59ac12SSebastian Grimberg     p = 0;
746ba59ac12SSebastian Grimberg     q = 0;
747ba59ac12SSebastian Grimberg     for (CeedInt i = n - 2; i >= 0; i--) {
748ba59ac12SSebastian Grimberg       if (fabs(mat_T[i + n * (i + 1)]) < tol) q += 1;
749ba59ac12SSebastian Grimberg       else break;
750ba59ac12SSebastian Grimberg     }
751ba59ac12SSebastian Grimberg     for (CeedInt i = 0; i < n - q - 1; i++) {
752ba59ac12SSebastian Grimberg       if (fabs(mat_T[i + n * (i + 1)]) < tol) p += 1;
753ba59ac12SSebastian Grimberg       else break;
754ba59ac12SSebastian Grimberg     }
755ba59ac12SSebastian Grimberg     if (q == n - 1) break;  // Finished reducing
756ba59ac12SSebastian Grimberg 
757ba59ac12SSebastian Grimberg     // Reduce tridiagonal portion
758ba59ac12SSebastian Grimberg     CeedScalar t_nn = mat_T[(n - 1 - q) + n * (n - 1 - q)], t_nnm1 = mat_T[(n - 2 - q) + n * (n - 1 - q)];
759ba59ac12SSebastian Grimberg     CeedScalar d  = (mat_T[(n - 2 - q) + n * (n - 2 - q)] - t_nn) / 2;
760ba59ac12SSebastian Grimberg     CeedScalar mu = t_nn - t_nnm1 * t_nnm1 / (d + copysign(sqrt(d * d + t_nnm1 * t_nnm1), d));
761ba59ac12SSebastian Grimberg     CeedScalar x  = mat_T[p + n * p] - mu;
762ba59ac12SSebastian Grimberg     CeedScalar z  = mat_T[p + n * (p + 1)];
763ba59ac12SSebastian Grimberg     for (CeedInt k = p; k < n - q - 1; k++) {
764ba59ac12SSebastian Grimberg       // Compute Givens rotation
765ba59ac12SSebastian Grimberg       CeedScalar c = 1, s = 0;
766ba59ac12SSebastian Grimberg       if (fabs(z) > tol) {
767ba59ac12SSebastian Grimberg         if (fabs(z) > fabs(x)) {
768ba59ac12SSebastian Grimberg           CeedScalar tau = -x / z;
769ba59ac12SSebastian Grimberg           s = 1 / sqrt(1 + tau * tau), c = s * tau;
770ba59ac12SSebastian Grimberg         } else {
771ba59ac12SSebastian Grimberg           CeedScalar tau = -z / x;
772ba59ac12SSebastian Grimberg           c = 1 / sqrt(1 + tau * tau), s = c * tau;
773ba59ac12SSebastian Grimberg         }
774ba59ac12SSebastian Grimberg       }
775ba59ac12SSebastian Grimberg 
776ba59ac12SSebastian Grimberg       // Apply Givens rotation to T
777ba59ac12SSebastian Grimberg       CeedGivensRotation(mat_T, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n);
778ba59ac12SSebastian Grimberg       CeedGivensRotation(mat_T, c, s, CEED_TRANSPOSE, k, k + 1, n, n);
779ba59ac12SSebastian Grimberg 
780ba59ac12SSebastian Grimberg       // Apply Givens rotation to Q
781ba59ac12SSebastian Grimberg       CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n);
782ba59ac12SSebastian Grimberg 
783ba59ac12SSebastian Grimberg       // Update x, z
784ba59ac12SSebastian Grimberg       if (k < n - q - 2) {
785ba59ac12SSebastian Grimberg         x = mat_T[k + n * (k + 1)];
786ba59ac12SSebastian Grimberg         z = mat_T[k + n * (k + 2)];
787ba59ac12SSebastian Grimberg       }
788ba59ac12SSebastian Grimberg     }
789ba59ac12SSebastian Grimberg     itr++;
790ba59ac12SSebastian Grimberg   }
791ba59ac12SSebastian Grimberg 
792ba59ac12SSebastian Grimberg   // Save eigenvalues
793ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) lambda[i] = mat_T[i + n * i];
794ba59ac12SSebastian Grimberg 
795ba59ac12SSebastian Grimberg   // Check convergence
796ba59ac12SSebastian Grimberg   if (itr == max_itr && q < n - 1) {
797ba59ac12SSebastian Grimberg     // LCOV_EXCL_START
798ba59ac12SSebastian Grimberg     return CeedError(ceed, CEED_ERROR_MINOR, "Symmetric QR failed to converge");
799ba59ac12SSebastian Grimberg     // LCOV_EXCL_STOP
800ba59ac12SSebastian Grimberg   }
801ba59ac12SSebastian Grimberg   return CEED_ERROR_SUCCESS;
802ba59ac12SSebastian Grimberg }
803ba59ac12SSebastian Grimberg CeedPragmaOptimizeOn;
804ba59ac12SSebastian Grimberg 
805ba59ac12SSebastian Grimberg /**
806ba59ac12SSebastian Grimberg   @brief Return Simultaneous Diagonalization of two matrices.
807ba59ac12SSebastian Grimberg 
808ba59ac12SSebastian Grimberg   This solves the generalized eigenvalue problem A x = lambda B x, where A and B are symmetric and B is positive definite.
809ba59ac12SSebastian Grimberg   We generate the matrix X and vector Lambda such that X^T A X = Lambda and X^T B X = I.
810ba59ac12SSebastian Grimberg   This is equivalent to the LAPACK routine 'sygv' with TYPE = 1.
811ba59ac12SSebastian Grimberg 
812ba59ac12SSebastian Grimberg   @param[in]  ceed   Ceed context for error handling
813ba59ac12SSebastian Grimberg   @param[in]  mat_A  Row-major matrix to be factorized with eigenvalues
814ba59ac12SSebastian Grimberg   @param[in]  mat_B  Row-major matrix to be factorized to identity
815ba59ac12SSebastian Grimberg   @param[out] mat_X  Row-major orthogonal matrix
816ba59ac12SSebastian Grimberg   @param[out] lambda Vector of length n of generalized eigenvalues
817ba59ac12SSebastian Grimberg   @param[in]  n      Number of rows/columns
818ba59ac12SSebastian Grimberg 
819ba59ac12SSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
820ba59ac12SSebastian Grimberg 
821ba59ac12SSebastian Grimberg   @ref Utility
822ba59ac12SSebastian Grimberg **/
823ba59ac12SSebastian Grimberg CeedPragmaOptimizeOff int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *mat_A, CeedScalar *mat_B, CeedScalar *mat_X, CeedScalar *lambda,
824ba59ac12SSebastian Grimberg                                                           CeedInt n) {
825ba59ac12SSebastian Grimberg   CeedScalar *mat_C, *mat_G, *vec_D;
826ba59ac12SSebastian Grimberg   CeedCall(CeedCalloc(n * n, &mat_C));
827ba59ac12SSebastian Grimberg   CeedCall(CeedCalloc(n * n, &mat_G));
828ba59ac12SSebastian Grimberg   CeedCall(CeedCalloc(n, &vec_D));
829ba59ac12SSebastian Grimberg 
830ba59ac12SSebastian Grimberg   // Compute B = G D G^T
831ba59ac12SSebastian Grimberg   memcpy(mat_G, mat_B, n * n * sizeof(mat_B[0]));
832ba59ac12SSebastian Grimberg   CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_G, vec_D, n));
833ba59ac12SSebastian Grimberg 
834ba59ac12SSebastian Grimberg   // Sort eigenvalues
835ba59ac12SSebastian Grimberg   for (CeedInt i = n - 1; i >= 0; i--) {
836ba59ac12SSebastian Grimberg     for (CeedInt j = 0; j < i; j++) {
837ba59ac12SSebastian Grimberg       if (fabs(vec_D[j]) > fabs(vec_D[j + 1])) {
838ba59ac12SSebastian Grimberg         CeedScalar temp;
839ba59ac12SSebastian Grimberg         temp         = vec_D[j];
840ba59ac12SSebastian Grimberg         vec_D[j]     = vec_D[j + 1];
841ba59ac12SSebastian Grimberg         vec_D[j + 1] = temp;
842ba59ac12SSebastian Grimberg         for (CeedInt k = 0; k < n; k++) {
843ba59ac12SSebastian Grimberg           temp                 = mat_G[k * n + j];
844ba59ac12SSebastian Grimberg           mat_G[k * n + j]     = mat_G[k * n + j + 1];
845ba59ac12SSebastian Grimberg           mat_G[k * n + j + 1] = temp;
846ba59ac12SSebastian Grimberg         }
847ba59ac12SSebastian Grimberg       }
848ba59ac12SSebastian Grimberg     }
849ba59ac12SSebastian Grimberg   }
850ba59ac12SSebastian Grimberg 
851ba59ac12SSebastian Grimberg   // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T
852ba59ac12SSebastian Grimberg   //           = D^-1/2 G^T A G D^-1/2
853ba59ac12SSebastian Grimberg   // -- D = D^-1/2
854ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) vec_D[i] = 1. / sqrt(vec_D[i]);
855ba59ac12SSebastian Grimberg   // -- G = G D^-1/2
856ba59ac12SSebastian Grimberg   // -- C = D^-1/2 G^T
857ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) {
858ba59ac12SSebastian Grimberg     for (CeedInt j = 0; j < n; j++) {
859ba59ac12SSebastian Grimberg       mat_G[i * n + j] *= vec_D[j];
860ba59ac12SSebastian Grimberg       mat_C[j * n + i] = mat_G[i * n + j];
861ba59ac12SSebastian Grimberg     }
862ba59ac12SSebastian Grimberg   }
863ba59ac12SSebastian Grimberg   // -- X = (D^-1/2 G^T) A
864ba59ac12SSebastian Grimberg   CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_C, (const CeedScalar *)mat_A, mat_X, n, n, n));
865ba59ac12SSebastian Grimberg   // -- C = (D^-1/2 G^T A) (G D^-1/2)
866ba59ac12SSebastian Grimberg   CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_X, (const CeedScalar *)mat_G, mat_C, n, n, n));
867ba59ac12SSebastian Grimberg 
868ba59ac12SSebastian Grimberg   // Compute Q^T C Q = lambda
869ba59ac12SSebastian Grimberg   CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_C, lambda, n));
870ba59ac12SSebastian Grimberg 
871ba59ac12SSebastian Grimberg   // Sort eigenvalues
872ba59ac12SSebastian Grimberg   for (CeedInt i = n - 1; i >= 0; i--) {
873ba59ac12SSebastian Grimberg     for (CeedInt j = 0; j < i; j++) {
874ba59ac12SSebastian Grimberg       if (fabs(lambda[j]) > fabs(lambda[j + 1])) {
875ba59ac12SSebastian Grimberg         CeedScalar temp;
876ba59ac12SSebastian Grimberg         temp          = lambda[j];
877ba59ac12SSebastian Grimberg         lambda[j]     = lambda[j + 1];
878ba59ac12SSebastian Grimberg         lambda[j + 1] = temp;
879ba59ac12SSebastian Grimberg         for (CeedInt k = 0; k < n; k++) {
880ba59ac12SSebastian Grimberg           temp                 = mat_C[k * n + j];
881ba59ac12SSebastian Grimberg           mat_C[k * n + j]     = mat_C[k * n + j + 1];
882ba59ac12SSebastian Grimberg           mat_C[k * n + j + 1] = temp;
883ba59ac12SSebastian Grimberg         }
884ba59ac12SSebastian Grimberg       }
885ba59ac12SSebastian Grimberg     }
886ba59ac12SSebastian Grimberg   }
887ba59ac12SSebastian Grimberg 
888ba59ac12SSebastian Grimberg   // Set X = (G D^1/2)^-T Q
889ba59ac12SSebastian Grimberg   //       = G D^-1/2 Q
890ba59ac12SSebastian Grimberg   CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_G, (const CeedScalar *)mat_C, mat_X, n, n, n));
891ba59ac12SSebastian Grimberg 
892ba59ac12SSebastian Grimberg   // Cleanup
893ba59ac12SSebastian Grimberg   CeedCall(CeedFree(&mat_C));
894ba59ac12SSebastian Grimberg   CeedCall(CeedFree(&mat_G));
895ba59ac12SSebastian Grimberg   CeedCall(CeedFree(&vec_D));
896ba59ac12SSebastian Grimberg   return CEED_ERROR_SUCCESS;
897ba59ac12SSebastian Grimberg }
898ba59ac12SSebastian Grimberg CeedPragmaOptimizeOn;
899ba59ac12SSebastian Grimberg 
9007a982d89SJeremy L. Thompson /// @}
9017a982d89SJeremy L. Thompson 
9027a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
9037a982d89SJeremy L. Thompson /// CeedBasis Public API
9047a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
9057a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser
906d7b241e6Sjeremylt /// @{
907d7b241e6Sjeremylt 
908b11c1e72Sjeremylt /**
909ba59ac12SSebastian Grimberg   @brief Create a tensor-product basis for H^1 discretizations
910b11c1e72Sjeremylt 
911ea61e9acSJeremy L Thompson   @param[in]  ceed        Ceed object where the CeedBasis will be created
912ea61e9acSJeremy L Thompson   @param[in]  dim         Topological dimension
913ea61e9acSJeremy L Thompson   @param[in]  num_comp    Number of field components (1 for scalar fields)
914ea61e9acSJeremy L Thompson   @param[in]  P_1d        Number of nodes in one dimension
915ea61e9acSJeremy L Thompson   @param[in]  Q_1d        Number of quadrature points in one dimension
916ea61e9acSJeremy L Thompson   @param[in]  interp_1d   Row-major (Q_1d * P_1d) matrix expressing the values of nodal basis functions at quadrature points
917ea61e9acSJeremy L Thompson   @param[in]  grad_1d     Row-major (Q_1d * P_1d) matrix expressing derivatives of nodal basis functions at quadrature points
918ea61e9acSJeremy 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]
919ea61e9acSJeremy L Thompson   @param[in]  q_weight_1d Array of length Q_1d holding the quadrature weights on the reference element
920ea61e9acSJeremy L Thompson   @param[out] basis       Address of the variable where the newly created CeedBasis will be stored.
921b11c1e72Sjeremylt 
922b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
923dfdf5a53Sjeremylt 
9247a982d89SJeremy L. Thompson   @ref User
925b11c1e72Sjeremylt **/
9262b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d,
9272b730f8bSJeremy L Thompson                             const CeedScalar *grad_1d, const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis *basis) {
9285fe0d4faSjeremylt   if (!ceed->BasisCreateTensorH1) {
9295fe0d4faSjeremylt     Ceed delegate;
9302b730f8bSJeremy L Thompson     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
9315fe0d4faSjeremylt 
9322b730f8bSJeremy L Thompson     if (!delegate) {
933c042f62fSJeremy L Thompson       // LCOV_EXCL_START
9342b730f8bSJeremy L Thompson       return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisCreateTensorH1");
935c042f62fSJeremy L Thompson       // LCOV_EXCL_STOP
9362b730f8bSJeremy L Thompson     }
9375fe0d4faSjeremylt 
9382b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateTensorH1(delegate, dim, num_comp, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis));
939e15f9bd0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
9405fe0d4faSjeremylt   }
941e15f9bd0SJeremy L Thompson 
9422b730f8bSJeremy L Thompson   if (dim < 1) {
943e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
9442b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis dimension must be a positive value");
945e15f9bd0SJeremy L Thompson     // LCOV_EXCL_STOP
9462b730f8bSJeremy L Thompson   }
947227444bfSJeremy L Thompson 
9482b730f8bSJeremy L Thompson   if (num_comp < 1) {
949227444bfSJeremy L Thompson     // LCOV_EXCL_START
9502b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component");
951227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
9522b730f8bSJeremy L Thompson   }
953227444bfSJeremy L Thompson 
9542b730f8bSJeremy L Thompson   if (P_1d < 1) {
955227444bfSJeremy L Thompson     // LCOV_EXCL_START
9562b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node");
957227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
9582b730f8bSJeremy L Thompson   }
959227444bfSJeremy L Thompson 
9602b730f8bSJeremy L Thompson   if (Q_1d < 1) {
961227444bfSJeremy L Thompson     // LCOV_EXCL_START
9622b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point");
963227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
9642b730f8bSJeremy L Thompson   }
965227444bfSJeremy L Thompson 
9662b730f8bSJeremy L Thompson   CeedElemTopology topo = dim == 1 ? CEED_TOPOLOGY_LINE : dim == 2 ? CEED_TOPOLOGY_QUAD : CEED_TOPOLOGY_HEX;
967e15f9bd0SJeremy L Thompson 
9682b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(1, basis));
969d7b241e6Sjeremylt   (*basis)->ceed = ceed;
9702b730f8bSJeremy L Thompson   CeedCall(CeedReference(ceed));
971d1d35e2fSjeremylt   (*basis)->ref_count    = 1;
972d1d35e2fSjeremylt   (*basis)->tensor_basis = 1;
973d7b241e6Sjeremylt   (*basis)->dim          = dim;
974d99fa3c5SJeremy L Thompson   (*basis)->topo         = topo;
975d1d35e2fSjeremylt   (*basis)->num_comp     = num_comp;
976d1d35e2fSjeremylt   (*basis)->P_1d         = P_1d;
977d1d35e2fSjeremylt   (*basis)->Q_1d         = Q_1d;
978d1d35e2fSjeremylt   (*basis)->P            = CeedIntPow(P_1d, dim);
979d1d35e2fSjeremylt   (*basis)->Q            = CeedIntPow(Q_1d, dim);
980c4e3f59bSSebastian Grimberg   (*basis)->fe_space     = CEED_FE_SPACE_H1;
9812b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d, &(*basis)->q_ref_1d));
9822b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d, &(*basis)->q_weight_1d));
983ff3a0f91SJeremy L Thompson   if (q_ref_1d) memcpy((*basis)->q_ref_1d, q_ref_1d, Q_1d * sizeof(q_ref_1d[0]));
9842b730f8bSJeremy L Thompson   if (q_weight_1d) memcpy((*basis)->q_weight_1d, q_weight_1d, Q_1d * sizeof(q_weight_1d[0]));
9852b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->interp_1d));
9862b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->grad_1d));
9872b730f8bSJeremy L Thompson   if (interp_1d) memcpy((*basis)->interp_1d, interp_1d, Q_1d * P_1d * sizeof(interp_1d[0]));
988ff3a0f91SJeremy L Thompson   if (grad_1d) memcpy((*basis)->grad_1d, grad_1d, Q_1d * P_1d * sizeof(grad_1d[0]));
9892b730f8bSJeremy L Thompson   CeedCall(ceed->BasisCreateTensorH1(dim, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, *basis));
990e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
991d7b241e6Sjeremylt }
992d7b241e6Sjeremylt 
993b11c1e72Sjeremylt /**
99495bb1877Svaleriabarra   @brief Create a tensor-product Lagrange basis
995b11c1e72Sjeremylt 
996ea61e9acSJeremy L Thompson   @param[in]  ceed      Ceed object where the CeedBasis will be created
997ea61e9acSJeremy L Thompson   @param[in]  dim       Topological dimension of element
998ea61e9acSJeremy L Thompson   @param[in]  num_comp  Number of field components (1 for scalar fields)
999ea61e9acSJeremy L Thompson   @param[in]  P         Number of Gauss-Lobatto nodes in one dimension.
1000ea61e9acSJeremy L Thompson                           The polynomial degree of the resulting Q_k element is k=P-1.
1001ea61e9acSJeremy L Thompson   @param[in]  Q         Number of quadrature points in one dimension.
1002ea61e9acSJeremy L Thompson   @param[in]  quad_mode Distribution of the Q quadrature points (affects order of accuracy for the quadrature)
1003ea61e9acSJeremy L Thompson   @param[out] basis     Address of the variable where the newly created CeedBasis will be stored.
1004b11c1e72Sjeremylt 
1005b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1006dfdf5a53Sjeremylt 
10077a982d89SJeremy L. Thompson   @ref User
1008b11c1e72Sjeremylt **/
10092b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P, CeedInt Q, CeedQuadMode quad_mode, CeedBasis *basis) {
1010d7b241e6Sjeremylt   // Allocate
10112b730f8bSJeremy L Thompson   int        ierr = CEED_ERROR_SUCCESS, i, j, k;
10122b730f8bSJeremy L Thompson   CeedScalar c1, c2, c3, c4, dx, *nodes, *interp_1d, *grad_1d, *q_ref_1d, *q_weight_1d;
10134d537eeaSYohann 
10142b730f8bSJeremy L Thompson   if (dim < 1) {
1015c042f62fSJeremy L Thompson     // LCOV_EXCL_START
10162b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis dimension must be a positive value");
1017c042f62fSJeremy L Thompson     // LCOV_EXCL_STOP
10182b730f8bSJeremy L Thompson   }
10194d537eeaSYohann 
10202b730f8bSJeremy L Thompson   if (num_comp < 1) {
1021227444bfSJeremy L Thompson     // LCOV_EXCL_START
10222b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component");
1023227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
10242b730f8bSJeremy L Thompson   }
1025227444bfSJeremy L Thompson 
10262b730f8bSJeremy L Thompson   if (P < 1) {
1027227444bfSJeremy L Thompson     // LCOV_EXCL_START
10282b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node");
1029227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
10302b730f8bSJeremy L Thompson   }
1031227444bfSJeremy L Thompson 
10322b730f8bSJeremy L Thompson   if (Q < 1) {
1033227444bfSJeremy L Thompson     // LCOV_EXCL_START
10342b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point");
1035227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
10362b730f8bSJeremy L Thompson   }
1037227444bfSJeremy L Thompson 
1038e15f9bd0SJeremy L Thompson   // Get Nodes and Weights
10392b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(P * Q, &interp_1d));
10402b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(P * Q, &grad_1d));
10412b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(P, &nodes));
10422b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q, &q_ref_1d));
10432b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q, &q_weight_1d));
10442b730f8bSJeremy L Thompson   if (CeedLobattoQuadrature(P, nodes, NULL) != CEED_ERROR_SUCCESS) goto cleanup;
1045d1d35e2fSjeremylt   switch (quad_mode) {
1046d7b241e6Sjeremylt     case CEED_GAUSS:
1047d1d35e2fSjeremylt       ierr = CeedGaussQuadrature(Q, q_ref_1d, q_weight_1d);
1048d7b241e6Sjeremylt       break;
1049d7b241e6Sjeremylt     case CEED_GAUSS_LOBATTO:
1050d1d35e2fSjeremylt       ierr = CeedLobattoQuadrature(Q, q_ref_1d, q_weight_1d);
1051d7b241e6Sjeremylt       break;
1052d7b241e6Sjeremylt   }
10532b730f8bSJeremy L Thompson   if (ierr != CEED_ERROR_SUCCESS) goto cleanup;
1054e15f9bd0SJeremy L Thompson 
1055d7b241e6Sjeremylt   // Build B, D matrix
1056d7b241e6Sjeremylt   // Fornberg, 1998
1057d7b241e6Sjeremylt   for (i = 0; i < Q; i++) {
1058d7b241e6Sjeremylt     c1                   = 1.0;
1059d1d35e2fSjeremylt     c3                   = nodes[0] - q_ref_1d[i];
1060d1d35e2fSjeremylt     interp_1d[i * P + 0] = 1.0;
1061d7b241e6Sjeremylt     for (j = 1; j < P; j++) {
1062d7b241e6Sjeremylt       c2 = 1.0;
1063d7b241e6Sjeremylt       c4 = c3;
1064d1d35e2fSjeremylt       c3 = nodes[j] - q_ref_1d[i];
1065d7b241e6Sjeremylt       for (k = 0; k < j; k++) {
1066d7b241e6Sjeremylt         dx = nodes[j] - nodes[k];
1067d7b241e6Sjeremylt         c2 *= dx;
1068d7b241e6Sjeremylt         if (k == j - 1) {
1069d1d35e2fSjeremylt           grad_1d[i * P + j]   = c1 * (interp_1d[i * P + k] - c4 * grad_1d[i * P + k]) / c2;
1070d1d35e2fSjeremylt           interp_1d[i * P + j] = -c1 * c4 * interp_1d[i * P + k] / c2;
1071d7b241e6Sjeremylt         }
1072d1d35e2fSjeremylt         grad_1d[i * P + k]   = (c3 * grad_1d[i * P + k] - interp_1d[i * P + k]) / dx;
1073d1d35e2fSjeremylt         interp_1d[i * P + k] = c3 * interp_1d[i * P + k] / dx;
1074d7b241e6Sjeremylt       }
1075d7b241e6Sjeremylt       c1 = c2;
1076d7b241e6Sjeremylt     }
1077d7b241e6Sjeremylt   }
10789ac7b42eSJeremy L Thompson   // Pass to CeedBasisCreateTensorH1
10792b730f8bSJeremy L Thompson   CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P, Q, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis));
1080e15f9bd0SJeremy L Thompson cleanup:
10812b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_1d));
10822b730f8bSJeremy L Thompson   CeedCall(CeedFree(&grad_1d));
10832b730f8bSJeremy L Thompson   CeedCall(CeedFree(&nodes));
10842b730f8bSJeremy L Thompson   CeedCall(CeedFree(&q_ref_1d));
10852b730f8bSJeremy L Thompson   CeedCall(CeedFree(&q_weight_1d));
1086e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1087d7b241e6Sjeremylt }
1088d7b241e6Sjeremylt 
1089b11c1e72Sjeremylt /**
1090ba59ac12SSebastian Grimberg   @brief Create a non tensor-product basis for H^1 discretizations
1091a8de75f0Sjeremylt 
1092ea61e9acSJeremy L Thompson   @param[in]  ceed      Ceed object where the CeedBasis will be created
1093ea61e9acSJeremy L Thompson   @param[in]  topo      Topology of element, e.g. hypercube, simplex, ect
1094ea61e9acSJeremy L Thompson   @param[in]  num_comp  Number of field components (1 for scalar fields)
1095ea61e9acSJeremy L Thompson   @param[in]  num_nodes Total number of nodes
1096ea61e9acSJeremy L Thompson   @param[in]  num_qpts  Total number of quadrature points
1097ea61e9acSJeremy L Thompson   @param[in]  interp    Row-major (num_qpts * num_nodes) matrix expressing the values of nodal basis functions at quadrature points
1098c4e3f59bSSebastian Grimberg   @param[in]  grad      Row-major (dim * num_qpts * num_nodes) matrix expressing derivatives of nodal basis functions at quadrature points
10999fe083eeSJeremy L Thompson   @param[in]  q_ref     Array of length num_qpts * dim holding the locations of quadrature points on the reference element
1100ea61e9acSJeremy L Thompson   @param[in]  q_weight  Array of length num_qpts holding the quadrature weights on the reference element
1101ea61e9acSJeremy L Thompson   @param[out] basis     Address of the variable where the newly created CeedBasis will be stored.
1102a8de75f0Sjeremylt 
1103a8de75f0Sjeremylt   @return An error code: 0 - success, otherwise - failure
1104a8de75f0Sjeremylt 
11057a982d89SJeremy L. Thompson   @ref User
1106a8de75f0Sjeremylt **/
11072b730f8bSJeremy L Thompson int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
11082b730f8bSJeremy L Thompson                       const CeedScalar *grad, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) {
1109d1d35e2fSjeremylt   CeedInt P = num_nodes, Q = num_qpts, dim = 0;
1110a8de75f0Sjeremylt 
11115fe0d4faSjeremylt   if (!ceed->BasisCreateH1) {
11125fe0d4faSjeremylt     Ceed delegate;
11132b730f8bSJeremy L Thompson     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
11145fe0d4faSjeremylt 
11152b730f8bSJeremy L Thompson     if (!delegate) {
1116c042f62fSJeremy L Thompson       // LCOV_EXCL_START
11172b730f8bSJeremy L Thompson       return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisCreateH1");
1118c042f62fSJeremy L Thompson       // LCOV_EXCL_STOP
11192b730f8bSJeremy L Thompson     }
11205fe0d4faSjeremylt 
11212b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateH1(delegate, topo, num_comp, num_nodes, num_qpts, interp, grad, q_ref, q_weight, basis));
1122e15f9bd0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
11235fe0d4faSjeremylt   }
11245fe0d4faSjeremylt 
11252b730f8bSJeremy L Thompson   if (num_comp < 1) {
1126227444bfSJeremy L Thompson     // LCOV_EXCL_START
11272b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component");
1128227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
11292b730f8bSJeremy L Thompson   }
1130227444bfSJeremy L Thompson 
11312b730f8bSJeremy L Thompson   if (num_nodes < 1) {
1132227444bfSJeremy L Thompson     // LCOV_EXCL_START
11332b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node");
1134227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
11352b730f8bSJeremy L Thompson   }
1136227444bfSJeremy L Thompson 
11372b730f8bSJeremy L Thompson   if (num_qpts < 1) {
1138227444bfSJeremy L Thompson     // LCOV_EXCL_START
11392b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point");
1140227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
11412b730f8bSJeremy L Thompson   }
1142227444bfSJeremy L Thompson 
11432b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(1, basis));
1144a8de75f0Sjeremylt 
11452b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetTopologyDimension(topo, &dim));
1146a8de75f0Sjeremylt 
1147a8de75f0Sjeremylt   (*basis)->ceed = ceed;
11482b730f8bSJeremy L Thompson   CeedCall(CeedReference(ceed));
1149d1d35e2fSjeremylt   (*basis)->ref_count    = 1;
1150d1d35e2fSjeremylt   (*basis)->tensor_basis = 0;
1151a8de75f0Sjeremylt   (*basis)->dim          = dim;
1152d99fa3c5SJeremy L Thompson   (*basis)->topo         = topo;
1153d1d35e2fSjeremylt   (*basis)->num_comp     = num_comp;
1154a8de75f0Sjeremylt   (*basis)->P            = P;
1155a8de75f0Sjeremylt   (*basis)->Q            = Q;
1156c4e3f59bSSebastian Grimberg   (*basis)->fe_space     = CEED_FE_SPACE_H1;
11572b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q * dim, &(*basis)->q_ref_1d));
11582b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q, &(*basis)->q_weight_1d));
1159ff3a0f91SJeremy L Thompson   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0]));
1160ff3a0f91SJeremy L Thompson   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0]));
11612b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q * P, &(*basis)->interp));
11622b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(dim * Q * P, &(*basis)->grad));
1163ff3a0f91SJeremy L Thompson   if (interp) memcpy((*basis)->interp, interp, Q * P * sizeof(interp[0]));
1164ff3a0f91SJeremy L Thompson   if (grad) memcpy((*basis)->grad, grad, dim * Q * P * sizeof(grad[0]));
11652b730f8bSJeremy L Thompson   CeedCall(ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, q_ref, q_weight, *basis));
1166e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1167a8de75f0Sjeremylt }
1168a8de75f0Sjeremylt 
1169a8de75f0Sjeremylt /**
1170859c15bbSJames Wright   @brief Create a non tensor-product basis for \f$H(\mathrm{div})\f$ discretizations
117150c301a5SRezgar Shakeri 
1172ea61e9acSJeremy L Thompson   @param[in]  ceed      Ceed object where the CeedBasis will be created
1173ea61e9acSJeremy 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
1174ea61e9acSJeremy L Thompson   @param[in]  num_comp  Number of components (usually 1 for vectors in H(div) bases)
1175ea61e9acSJeremy L Thompson   @param[in]  num_nodes Total number of nodes (dofs per element)
1176ea61e9acSJeremy L Thompson   @param[in]  num_qpts  Total number of quadrature points
1177c4e3f59bSSebastian Grimberg   @param[in]  interp    Row-major (dim * num_qpts * num_nodes) matrix expressing the values of basis functions at quadrature points
1178c4e3f59bSSebastian Grimberg   @param[in]  div       Row-major (num_qpts * num_nodes) matrix expressing divergence of basis functions at quadrature points
11799fe083eeSJeremy L Thompson   @param[in]  q_ref     Array of length num_qpts * dim holding the locations of quadrature points on the reference element
1180ea61e9acSJeremy L Thompson   @param[in]  q_weight  Array of length num_qpts holding the quadrature weights on the reference element
1181ea61e9acSJeremy L Thompson   @param[out] basis     Address of the variable where the newly created CeedBasis will be stored.
118250c301a5SRezgar Shakeri 
118350c301a5SRezgar Shakeri   @return An error code: 0 - success, otherwise - failure
118450c301a5SRezgar Shakeri 
118550c301a5SRezgar Shakeri   @ref User
118650c301a5SRezgar Shakeri **/
11872b730f8bSJeremy L Thompson int CeedBasisCreateHdiv(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
11882b730f8bSJeremy L Thompson                         const CeedScalar *div, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) {
118950c301a5SRezgar Shakeri   CeedInt Q = num_qpts, P = num_nodes, dim = 0;
1190c4e3f59bSSebastian Grimberg 
119150c301a5SRezgar Shakeri   if (!ceed->BasisCreateHdiv) {
119250c301a5SRezgar Shakeri     Ceed delegate;
11932b730f8bSJeremy L Thompson     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
119450c301a5SRezgar Shakeri 
11952b730f8bSJeremy L Thompson     if (!delegate) {
119650c301a5SRezgar Shakeri       // LCOV_EXCL_START
11972b730f8bSJeremy L Thompson       return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateHdiv");
119850c301a5SRezgar Shakeri       // LCOV_EXCL_STOP
11992b730f8bSJeremy L Thompson     }
120050c301a5SRezgar Shakeri 
12012b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateHdiv(delegate, topo, num_comp, num_nodes, num_qpts, interp, div, q_ref, q_weight, basis));
120250c301a5SRezgar Shakeri     return CEED_ERROR_SUCCESS;
120350c301a5SRezgar Shakeri   }
120450c301a5SRezgar Shakeri 
12052b730f8bSJeremy L Thompson   if (num_comp < 1) {
1206227444bfSJeremy L Thompson     // LCOV_EXCL_START
12072b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component");
1208227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
12092b730f8bSJeremy L Thompson   }
1210227444bfSJeremy L Thompson 
12112b730f8bSJeremy L Thompson   if (num_nodes < 1) {
1212227444bfSJeremy L Thompson     // LCOV_EXCL_START
12132b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node");
1214227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
12152b730f8bSJeremy L Thompson   }
1216227444bfSJeremy L Thompson 
12172b730f8bSJeremy L Thompson   if (num_qpts < 1) {
1218227444bfSJeremy L Thompson     // LCOV_EXCL_START
12192b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point");
1220227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
12212b730f8bSJeremy L Thompson   }
1222227444bfSJeremy L Thompson 
12232b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(1, basis));
122450c301a5SRezgar Shakeri 
1225c4e3f59bSSebastian Grimberg   CeedCall(CeedBasisGetTopologyDimension(topo, &dim));
1226c4e3f59bSSebastian Grimberg 
122750c301a5SRezgar Shakeri   (*basis)->ceed = ceed;
12282b730f8bSJeremy L Thompson   CeedCall(CeedReference(ceed));
122950c301a5SRezgar Shakeri   (*basis)->ref_count    = 1;
123050c301a5SRezgar Shakeri   (*basis)->tensor_basis = 0;
123150c301a5SRezgar Shakeri   (*basis)->dim          = dim;
123250c301a5SRezgar Shakeri   (*basis)->topo         = topo;
123350c301a5SRezgar Shakeri   (*basis)->num_comp     = num_comp;
123450c301a5SRezgar Shakeri   (*basis)->P            = P;
123550c301a5SRezgar Shakeri   (*basis)->Q            = Q;
1236c4e3f59bSSebastian Grimberg   (*basis)->fe_space     = CEED_FE_SPACE_HDIV;
12372b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q * dim, &(*basis)->q_ref_1d));
12382b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q, &(*basis)->q_weight_1d));
123950c301a5SRezgar Shakeri   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0]));
124050c301a5SRezgar Shakeri   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0]));
12412b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(dim * Q * P, &(*basis)->interp));
12422b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q * P, &(*basis)->div));
124350c301a5SRezgar Shakeri   if (interp) memcpy((*basis)->interp, interp, dim * Q * P * sizeof(interp[0]));
124450c301a5SRezgar Shakeri   if (div) memcpy((*basis)->div, div, Q * P * sizeof(div[0]));
12452b730f8bSJeremy L Thompson   CeedCall(ceed->BasisCreateHdiv(topo, dim, P, Q, interp, div, q_ref, q_weight, *basis));
124650c301a5SRezgar Shakeri   return CEED_ERROR_SUCCESS;
124750c301a5SRezgar Shakeri }
124850c301a5SRezgar Shakeri 
124950c301a5SRezgar Shakeri /**
1250c4e3f59bSSebastian Grimberg   @brief Create a non tensor-product basis for H(curl) discretizations
1251c4e3f59bSSebastian Grimberg 
1252c4e3f59bSSebastian Grimberg   @param[in]  ceed      Ceed object where the CeedBasis will be created
1253c4e3f59bSSebastian Grimberg   @param[in]  topo      Topology of element (`CEED_TOPOLOGY_QUAD`, `CEED_TOPOLOGY_PRISM`, etc.), dimension of which is used in some array sizes below
1254c4e3f59bSSebastian Grimberg   @param[in]  num_comp  Number of components (usually 1 for vectors in H(curl) bases)
1255c4e3f59bSSebastian Grimberg   @param[in]  num_nodes Total number of nodes (dofs per element)
1256c4e3f59bSSebastian Grimberg   @param[in]  num_qpts  Total number of quadrature points
1257c4e3f59bSSebastian Grimberg   @param[in]  interp    Row-major (dim * num_qpts * num_nodes) matrix expressing the values of basis functions at quadrature points
1258c4e3f59bSSebastian Grimberg   @param[in]  curl      Row-major (curl_comp * num_qpts * num_nodes, curl_comp = 1 if dim < 3 else dim) matrix expressing curl of basis functions at
1259c4e3f59bSSebastian Grimberg quadrature points
1260c4e3f59bSSebastian Grimberg   @param[in]  q_ref     Array of length num_qpts * dim holding the locations of quadrature points on the reference element
1261c4e3f59bSSebastian Grimberg   @param[in]  q_weight  Array of length num_qpts holding the quadrature weights on the reference element
1262c4e3f59bSSebastian Grimberg   @param[out] basis     Address of the variable where the newly created CeedBasis will be stored.
1263c4e3f59bSSebastian Grimberg 
1264c4e3f59bSSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
1265c4e3f59bSSebastian Grimberg 
1266c4e3f59bSSebastian Grimberg   @ref User
1267c4e3f59bSSebastian Grimberg **/
1268c4e3f59bSSebastian Grimberg int CeedBasisCreateHcurl(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
1269c4e3f59bSSebastian Grimberg                          const CeedScalar *curl, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) {
1270c4e3f59bSSebastian Grimberg   CeedInt Q = num_qpts, P = num_nodes, dim = 0, curl_comp = 0;
1271c4e3f59bSSebastian Grimberg 
1272c4e3f59bSSebastian Grimberg   if (!ceed->BasisCreateHdiv) {
1273c4e3f59bSSebastian Grimberg     Ceed delegate;
1274c4e3f59bSSebastian Grimberg     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
1275c4e3f59bSSebastian Grimberg 
1276c4e3f59bSSebastian Grimberg     if (!delegate) {
1277c4e3f59bSSebastian Grimberg       // LCOV_EXCL_START
1278c4e3f59bSSebastian Grimberg       return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateHcurl");
1279c4e3f59bSSebastian Grimberg       // LCOV_EXCL_STOP
1280c4e3f59bSSebastian Grimberg     }
1281c4e3f59bSSebastian Grimberg 
1282c4e3f59bSSebastian Grimberg     CeedCall(CeedBasisCreateHcurl(delegate, topo, num_comp, num_nodes, num_qpts, interp, curl, q_ref, q_weight, basis));
1283c4e3f59bSSebastian Grimberg     return CEED_ERROR_SUCCESS;
1284c4e3f59bSSebastian Grimberg   }
1285c4e3f59bSSebastian Grimberg 
1286c4e3f59bSSebastian Grimberg   if (num_comp < 1) {
1287c4e3f59bSSebastian Grimberg     // LCOV_EXCL_START
1288c4e3f59bSSebastian Grimberg     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component");
1289c4e3f59bSSebastian Grimberg     // LCOV_EXCL_STOP
1290c4e3f59bSSebastian Grimberg   }
1291c4e3f59bSSebastian Grimberg 
1292c4e3f59bSSebastian Grimberg   if (num_nodes < 1) {
1293c4e3f59bSSebastian Grimberg     // LCOV_EXCL_START
1294c4e3f59bSSebastian Grimberg     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node");
1295c4e3f59bSSebastian Grimberg     // LCOV_EXCL_STOP
1296c4e3f59bSSebastian Grimberg   }
1297c4e3f59bSSebastian Grimberg 
1298c4e3f59bSSebastian Grimberg   if (num_qpts < 1) {
1299c4e3f59bSSebastian Grimberg     // LCOV_EXCL_START
1300c4e3f59bSSebastian Grimberg     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point");
1301c4e3f59bSSebastian Grimberg     // LCOV_EXCL_STOP
1302c4e3f59bSSebastian Grimberg   }
1303c4e3f59bSSebastian Grimberg 
1304c4e3f59bSSebastian Grimberg   CeedCall(CeedCalloc(1, basis));
1305c4e3f59bSSebastian Grimberg 
1306c4e3f59bSSebastian Grimberg   CeedCall(CeedBasisGetTopologyDimension(topo, &dim));
1307c4e3f59bSSebastian Grimberg   curl_comp = (dim < 3) ? 1 : dim;
1308c4e3f59bSSebastian Grimberg 
1309c4e3f59bSSebastian Grimberg   (*basis)->ceed = ceed;
1310c4e3f59bSSebastian Grimberg   CeedCall(CeedReference(ceed));
1311c4e3f59bSSebastian Grimberg   (*basis)->ref_count    = 1;
1312c4e3f59bSSebastian Grimberg   (*basis)->tensor_basis = 0;
1313c4e3f59bSSebastian Grimberg   (*basis)->dim          = dim;
1314c4e3f59bSSebastian Grimberg   (*basis)->topo         = topo;
1315c4e3f59bSSebastian Grimberg   (*basis)->num_comp     = num_comp;
1316c4e3f59bSSebastian Grimberg   (*basis)->P            = P;
1317c4e3f59bSSebastian Grimberg   (*basis)->Q            = Q;
1318c4e3f59bSSebastian Grimberg   (*basis)->fe_space     = CEED_FE_SPACE_HCURL;
1319c4e3f59bSSebastian Grimberg   CeedCall(CeedMalloc(Q * dim, &(*basis)->q_ref_1d));
1320c4e3f59bSSebastian Grimberg   CeedCall(CeedMalloc(Q, &(*basis)->q_weight_1d));
1321c4e3f59bSSebastian Grimberg   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0]));
1322c4e3f59bSSebastian Grimberg   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0]));
1323c4e3f59bSSebastian Grimberg   CeedCall(CeedMalloc(dim * Q * P, &(*basis)->interp));
1324c4e3f59bSSebastian Grimberg   CeedCall(CeedMalloc(curl_comp * Q * P, &(*basis)->curl));
1325c4e3f59bSSebastian Grimberg   if (interp) memcpy((*basis)->interp, interp, dim * Q * P * sizeof(interp[0]));
1326c4e3f59bSSebastian Grimberg   if (curl) memcpy((*basis)->curl, curl, curl_comp * Q * P * sizeof(curl[0]));
1327c4e3f59bSSebastian Grimberg   CeedCall(ceed->BasisCreateHcurl(topo, dim, P, Q, interp, curl, q_ref, q_weight, *basis));
1328c4e3f59bSSebastian Grimberg   return CEED_ERROR_SUCCESS;
1329c4e3f59bSSebastian Grimberg }
1330c4e3f59bSSebastian Grimberg 
1331c4e3f59bSSebastian Grimberg /**
1332ea61e9acSJeremy L Thompson   @brief Create a CeedBasis for projection from the nodes of `basis_from` to the nodes of `basis_to`.
1333ba59ac12SSebastian Grimberg 
1334*15ad3917SSebastian Grimberg   Only `CEED_EVAL_INTERP` will be valid for the new basis, `basis_project`. For H^1 spaces, `CEED_EVAL_GRAD` will also be valid.
1335ea61e9acSJeremy L Thompson   The interpolation is given by `interp_project = interp_to^+ * interp_from`, where the pesudoinverse `interp_to^+` is given by QR
1336*15ad3917SSebastian Grimberg factorization.  The gradient (for the H^1 case) is given by `grad_project = interp_to^+ * grad_from`.
1337*15ad3917SSebastian Grimberg 
1338*15ad3917SSebastian Grimberg   Note: `basis_from` and `basis_to` must have compatible quadrature spaces.
1339*15ad3917SSebastian Grimberg 
1340ba59ac12SSebastian Grimberg   Note: `basis_project` will have the same number of components as `basis_from`, regardless of the number of components that `basis_to` has. If
1341ea61e9acSJeremy L Thompson `basis_from` has 3 components and `basis_to` has 5 components, then `basis_project` will have 3 components.
1342f113e5dcSJeremy L Thompson 
1343f113e5dcSJeremy L Thompson   @param[in]  basis_from    CeedBasis to prolong from
1344446e7af4SJeremy L Thompson   @param[in]  basis_to      CeedBasis to prolong to
1345ea61e9acSJeremy L Thompson   @param[out] basis_project Address of the variable where the newly created CeedBasis will be stored.
1346f113e5dcSJeremy L Thompson 
1347f113e5dcSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1348f113e5dcSJeremy L Thompson 
1349f113e5dcSJeremy L Thompson   @ref User
1350f113e5dcSJeremy L Thompson **/
13512b730f8bSJeremy L Thompson int CeedBasisCreateProjection(CeedBasis basis_from, CeedBasis basis_to, CeedBasis *basis_project) {
1352f113e5dcSJeremy L Thompson   Ceed ceed;
13532b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetCeed(basis_to, &ceed));
1354f113e5dcSJeremy L Thompson 
1355ecc88aebSJeremy L Thompson   // Create projection matrix
135614556e63SJeremy L Thompson   CeedScalar *interp_project, *grad_project;
13572b730f8bSJeremy L Thompson   CeedCall(CeedBasisCreateProjectionMatrices(basis_from, basis_to, &interp_project, &grad_project));
1358f113e5dcSJeremy L Thompson 
1359f113e5dcSJeremy L Thompson   // Build basis
1360f113e5dcSJeremy L Thompson   bool        is_tensor;
1361*15ad3917SSebastian Grimberg   CeedFESpace fe_space;
1362f113e5dcSJeremy L Thompson   CeedInt     dim, num_comp;
136314556e63SJeremy L Thompson   CeedScalar *q_ref, *q_weight;
13642b730f8bSJeremy L Thompson   CeedCall(CeedBasisIsTensor(basis_to, &is_tensor));
1365*15ad3917SSebastian Grimberg   CeedCall(CeedBasisGetFESpace(basis_to, &fe_space));
13662b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetDimension(basis_to, &dim));
13672b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumComponents(basis_from, &num_comp));
1368f113e5dcSJeremy L Thompson   if (is_tensor) {
1369f113e5dcSJeremy L Thompson     CeedInt P_1d_to, P_1d_from;
13702b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_1d_from));
13712b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_1d_to));
13722b730f8bSJeremy L Thompson     CeedCall(CeedCalloc(P_1d_to, &q_ref));
13732b730f8bSJeremy L Thompson     CeedCall(CeedCalloc(P_1d_to, &q_weight));
13742b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P_1d_from, P_1d_to, interp_project, grad_project, q_ref, q_weight, basis_project));
1375f113e5dcSJeremy L Thompson   } else {
1376f113e5dcSJeremy L Thompson     CeedElemTopology topo;
13772b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetTopology(basis_to, &topo));
1378f113e5dcSJeremy L Thompson     CeedInt num_nodes_to, num_nodes_from;
13792b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis_from, &num_nodes_from));
13802b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis_to, &num_nodes_to));
13812b730f8bSJeremy L Thompson     CeedCall(CeedCalloc(num_nodes_to * dim, &q_ref));
13822b730f8bSJeremy L Thompson     CeedCall(CeedCalloc(num_nodes_to, &q_weight));
1383*15ad3917SSebastian Grimberg     switch (fe_space) {
1384*15ad3917SSebastian Grimberg       case CEED_FE_SPACE_H1:
13852b730f8bSJeremy L Thompson         CeedCall(CeedBasisCreateH1(ceed, topo, num_comp, num_nodes_from, num_nodes_to, interp_project, grad_project, q_ref, q_weight, basis_project));
1386*15ad3917SSebastian Grimberg         break;
1387*15ad3917SSebastian Grimberg       case CEED_FE_SPACE_HDIV:
1388*15ad3917SSebastian Grimberg         CeedCall(
1389*15ad3917SSebastian Grimberg             CeedBasisCreateHdiv(ceed, topo, num_comp, num_nodes_from, num_nodes_to, interp_project, grad_project, q_ref, q_weight, basis_project));
1390*15ad3917SSebastian Grimberg         break;
1391*15ad3917SSebastian Grimberg       case CEED_FE_SPACE_HCURL:
1392*15ad3917SSebastian Grimberg         CeedCall(
1393*15ad3917SSebastian Grimberg             CeedBasisCreateHcurl(ceed, topo, num_comp, num_nodes_from, num_nodes_to, interp_project, grad_project, q_ref, q_weight, basis_project));
1394*15ad3917SSebastian Grimberg         break;
1395*15ad3917SSebastian Grimberg     }
1396f113e5dcSJeremy L Thompson   }
1397f113e5dcSJeremy L Thompson 
1398f113e5dcSJeremy L Thompson   // Cleanup
13992b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_project));
14002b730f8bSJeremy L Thompson   CeedCall(CeedFree(&grad_project));
14012b730f8bSJeremy L Thompson   CeedCall(CeedFree(&q_ref));
14022b730f8bSJeremy L Thompson   CeedCall(CeedFree(&q_weight));
1403f113e5dcSJeremy L Thompson 
1404f113e5dcSJeremy L Thompson   return CEED_ERROR_SUCCESS;
1405f113e5dcSJeremy L Thompson }
1406f113e5dcSJeremy L Thompson 
1407f113e5dcSJeremy L Thompson /**
1408ea61e9acSJeremy L Thompson   @brief Copy the pointer to a CeedBasis.
14099560d06aSjeremylt 
1410512bb800SJeremy L Thompson   Note: If the value of `basis_copy` passed into this function is non-NULL, then it is assumed that `basis_copy` is a pointer to a CeedBasis.
1411512bb800SJeremy L Thompson   This CeedBasis will be destroyed if `basis_copy` is the only reference to this CeedBasis.
1412ea61e9acSJeremy L Thompson 
1413ea61e9acSJeremy L Thompson   @param[in]     basis      CeedBasis to copy reference to
1414ea61e9acSJeremy L Thompson   @param[in,out] basis_copy Variable to store copied reference
14159560d06aSjeremylt 
14169560d06aSjeremylt   @return An error code: 0 - success, otherwise - failure
14179560d06aSjeremylt 
14189560d06aSjeremylt   @ref User
14199560d06aSjeremylt **/
14209560d06aSjeremylt int CeedBasisReferenceCopy(CeedBasis basis, CeedBasis *basis_copy) {
1421393ac2cdSJeremy L Thompson   if (basis != CEED_BASIS_COLLOCATED) CeedCall(CeedBasisReference(basis));
14222b730f8bSJeremy L Thompson   CeedCall(CeedBasisDestroy(basis_copy));
14239560d06aSjeremylt   *basis_copy = basis;
14249560d06aSjeremylt   return CEED_ERROR_SUCCESS;
14259560d06aSjeremylt }
14269560d06aSjeremylt 
14279560d06aSjeremylt /**
14287a982d89SJeremy L. Thompson   @brief View a CeedBasis
14297a982d89SJeremy L. Thompson 
1430ea61e9acSJeremy L Thompson   @param[in] basis  CeedBasis to view
1431ea61e9acSJeremy L Thompson   @param[in] stream Stream to view to, e.g., stdout
14327a982d89SJeremy L. Thompson 
14337a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
14347a982d89SJeremy L. Thompson 
14357a982d89SJeremy L. Thompson   @ref User
14367a982d89SJeremy L. Thompson **/
14377a982d89SJeremy L. Thompson int CeedBasisView(CeedBasis basis, FILE *stream) {
143850c301a5SRezgar Shakeri   CeedElemTopology topo     = basis->topo;
1439c4e3f59bSSebastian Grimberg   CeedFESpace      fe_space = basis->fe_space;
1440c4e3f59bSSebastian Grimberg   CeedInt          q_comp   = 0;
14412b730f8bSJeremy L Thompson 
144250c301a5SRezgar Shakeri   // Print FE space and element topology of the basis
1443d1d35e2fSjeremylt   if (basis->tensor_basis) {
1444c4e3f59bSSebastian Grimberg     fprintf(stream, "CeedBasis (%s on a %s element): dim=%" CeedInt_FMT " P=%" CeedInt_FMT " Q=%" CeedInt_FMT "\n", CeedFESpaces[fe_space],
14452b730f8bSJeremy L Thompson             CeedElemTopologies[topo], basis->dim, basis->P_1d, basis->Q_1d);
144650c301a5SRezgar Shakeri   } else {
1447c4e3f59bSSebastian Grimberg     fprintf(stream, "CeedBasis (%s on a %s element): dim=%" CeedInt_FMT " P=%" CeedInt_FMT " Q=%" CeedInt_FMT "\n", CeedFESpaces[fe_space],
14482b730f8bSJeremy L Thompson             CeedElemTopologies[topo], basis->dim, basis->P, basis->Q);
144950c301a5SRezgar Shakeri   }
1450ea61e9acSJeremy L Thompson   // Print quadrature data, interpolation/gradient/divergence/curl of the basis
145150c301a5SRezgar Shakeri   if (basis->tensor_basis) {  // tensor basis
14522b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_ref_1d, stream));
14532b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_weight_1d, stream));
14542b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("interp1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, basis->interp_1d, stream));
14552b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("grad1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, basis->grad_1d, stream));
145650c301a5SRezgar Shakeri   } else {  // non-tensor basis
14572b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("qref", "\t% 12.8f", 1, basis->Q * basis->dim, basis->q_ref_1d, stream));
14582b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->q_weight_1d, stream));
1459c4e3f59bSSebastian Grimberg     CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp));
1460c4e3f59bSSebastian Grimberg     CeedCall(CeedScalarView("interp", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->interp, stream));
146150c301a5SRezgar Shakeri     if (basis->grad) {
1462c4e3f59bSSebastian Grimberg       CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp));
1463c4e3f59bSSebastian Grimberg       CeedCall(CeedScalarView("grad", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->grad, stream));
14647a982d89SJeremy L. Thompson     }
146550c301a5SRezgar Shakeri     if (basis->div) {
1466c4e3f59bSSebastian Grimberg       CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_DIV, &q_comp));
1467c4e3f59bSSebastian Grimberg       CeedCall(CeedScalarView("div", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->div, stream));
1468c4e3f59bSSebastian Grimberg     }
1469c4e3f59bSSebastian Grimberg     if (basis->curl) {
1470c4e3f59bSSebastian Grimberg       CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_CURL, &q_comp));
1471c4e3f59bSSebastian Grimberg       CeedCall(CeedScalarView("curl", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->curl, stream));
147250c301a5SRezgar Shakeri     }
147350c301a5SRezgar Shakeri   }
1474e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
14757a982d89SJeremy L. Thompson }
14767a982d89SJeremy L. Thompson 
14777a982d89SJeremy L. Thompson /**
14787a982d89SJeremy L. Thompson   @brief Apply basis evaluation from nodes to quadrature points or vice versa
14797a982d89SJeremy L. Thompson 
1480ea61e9acSJeremy L Thompson   @param[in]  basis      CeedBasis to evaluate
1481ea61e9acSJeremy L Thompson   @param[in]  num_elem   The number of elements to apply the basis evaluation to;
1482ea61e9acSJeremy L Thompson                            the backend will specify the ordering in CeedElemRestrictionCreateBlocked()
1483ea61e9acSJeremy L Thompson   @param[in]  t_mode    \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature points;
1484ea61e9acSJeremy L Thompson                           \ref CEED_TRANSPOSE to apply the transpose, mapping from quadrature points to nodes
1485ea61e9acSJeremy L Thompson   @param[in]  eval_mode \ref CEED_EVAL_NONE to use values directly,
14867a982d89SJeremy L. Thompson                           \ref CEED_EVAL_INTERP to use interpolated values,
14877a982d89SJeremy L. Thompson                           \ref CEED_EVAL_GRAD to use gradients,
1488c4e3f59bSSebastian Grimberg                           \ref CEED_EVAL_DIV to use divergence,
1489c4e3f59bSSebastian Grimberg                           \ref CEED_EVAL_CURL to use curl,
14907a982d89SJeremy L. Thompson                           \ref CEED_EVAL_WEIGHT to use quadrature weights.
14917a982d89SJeremy L. Thompson   @param[in]  u        Input CeedVector
14927a982d89SJeremy L. Thompson   @param[out] v        Output CeedVector
14937a982d89SJeremy L. Thompson 
14947a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
14957a982d89SJeremy L. Thompson 
14967a982d89SJeremy L. Thompson   @ref User
14977a982d89SJeremy L. Thompson **/
14982b730f8bSJeremy L Thompson int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
14991f9221feSJeremy L Thompson   CeedSize u_length = 0, v_length;
1500c4e3f59bSSebastian Grimberg   CeedInt  dim, num_comp, q_comp, num_nodes, num_qpts;
15012b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetDimension(basis, &dim));
15022b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
1503c4e3f59bSSebastian Grimberg   CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp));
15042b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumNodes(basis, &num_nodes));
15052b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts));
15062b730f8bSJeremy L Thompson   CeedCall(CeedVectorGetLength(v, &v_length));
15077a982d89SJeremy L. Thompson   if (u) {
15082b730f8bSJeremy L Thompson     CeedCall(CeedVectorGetLength(u, &u_length));
15097a982d89SJeremy L. Thompson   }
15107a982d89SJeremy L. Thompson 
15112b730f8bSJeremy L Thompson   if (!basis->Apply) {
1512e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
15132b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisApply");
1514e15f9bd0SJeremy L Thompson     // LCOV_EXCL_STOP
15152b730f8bSJeremy L Thompson   }
1516e15f9bd0SJeremy L Thompson 
1517e15f9bd0SJeremy L Thompson   // Check compatibility of topological and geometrical dimensions
15182b730f8bSJeremy L Thompson   if ((t_mode == CEED_TRANSPOSE && (v_length % num_nodes != 0 || u_length % num_qpts != 0)) ||
15192b730f8bSJeremy L Thompson       (t_mode == CEED_NOTRANSPOSE && (u_length % num_nodes != 0 || v_length % num_qpts != 0))) {
15208229195eSjeremylt     // LCOV_EXCL_START
15212b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_DIMENSION, "Length of input/output vectors incompatible with basis dimensions");
15228229195eSjeremylt     // LCOV_EXCL_STOP
15232b730f8bSJeremy L Thompson   }
15247a982d89SJeremy L. Thompson 
1525e15f9bd0SJeremy L Thompson   // Check vector lengths to prevent out of bounds issues
1526d1d35e2fSjeremylt   bool bad_dims = false;
1527d1d35e2fSjeremylt   switch (eval_mode) {
1528e15f9bd0SJeremy L Thompson     case CEED_EVAL_NONE:
15292b730f8bSJeremy L Thompson     case CEED_EVAL_INTERP:
15302b730f8bSJeremy L Thompson     case CEED_EVAL_GRAD:
1531c4e3f59bSSebastian Grimberg     case CEED_EVAL_DIV:
1532c4e3f59bSSebastian Grimberg     case CEED_EVAL_CURL:
1533c4e3f59bSSebastian Grimberg       bad_dims = ((t_mode == CEED_TRANSPOSE && (u_length < num_elem * num_comp * num_qpts * q_comp || v_length < num_elem * num_comp * num_nodes)) ||
1534c4e3f59bSSebastian Grimberg                   (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem * num_qpts * num_comp * q_comp || u_length < num_elem * num_comp * num_nodes)));
1535e15f9bd0SJeremy L Thompson       break;
1536e15f9bd0SJeremy L Thompson     case CEED_EVAL_WEIGHT:
1537d1d35e2fSjeremylt       bad_dims = v_length < num_elem * num_qpts;
1538e15f9bd0SJeremy L Thompson       break;
1539e15f9bd0SJeremy L Thompson   }
15402b730f8bSJeremy L Thompson   if (bad_dims) {
1541e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
15422b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_DIMENSION, "Input/output vectors too short for basis and evaluation mode");
1543e15f9bd0SJeremy L Thompson     // LCOV_EXCL_STOP
15442b730f8bSJeremy L Thompson   }
1545e15f9bd0SJeremy L Thompson 
15462b730f8bSJeremy L Thompson   CeedCall(basis->Apply(basis, num_elem, t_mode, eval_mode, u, v));
1547e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
15487a982d89SJeremy L. Thompson }
15497a982d89SJeremy L. Thompson 
15507a982d89SJeremy L. Thompson /**
1551b7c9bbdaSJeremy L Thompson   @brief Get Ceed associated with a CeedBasis
1552b7c9bbdaSJeremy L Thompson 
1553ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
1554b7c9bbdaSJeremy L Thompson   @param[out] ceed  Variable to store Ceed
1555b7c9bbdaSJeremy L Thompson 
1556b7c9bbdaSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1557b7c9bbdaSJeremy L Thompson 
1558b7c9bbdaSJeremy L Thompson   @ref Advanced
1559b7c9bbdaSJeremy L Thompson **/
1560b7c9bbdaSJeremy L Thompson int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) {
1561b7c9bbdaSJeremy L Thompson   *ceed = basis->ceed;
1562b7c9bbdaSJeremy L Thompson   return CEED_ERROR_SUCCESS;
1563b7c9bbdaSJeremy L Thompson }
1564b7c9bbdaSJeremy L Thompson 
1565b7c9bbdaSJeremy L Thompson /**
15669d007619Sjeremylt   @brief Get dimension for given CeedBasis
15679d007619Sjeremylt 
1568ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
15699d007619Sjeremylt   @param[out] dim   Variable to store dimension of basis
15709d007619Sjeremylt 
15719d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
15729d007619Sjeremylt 
1573b7c9bbdaSJeremy L Thompson   @ref Advanced
15749d007619Sjeremylt **/
15759d007619Sjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) {
15769d007619Sjeremylt   *dim = basis->dim;
1577e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
15789d007619Sjeremylt }
15799d007619Sjeremylt 
15809d007619Sjeremylt /**
1581d99fa3c5SJeremy L Thompson   @brief Get topology for given CeedBasis
1582d99fa3c5SJeremy L Thompson 
1583ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
1584d99fa3c5SJeremy L Thompson   @param[out] topo  Variable to store topology of basis
1585d99fa3c5SJeremy L Thompson 
1586d99fa3c5SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1587d99fa3c5SJeremy L Thompson 
1588b7c9bbdaSJeremy L Thompson   @ref Advanced
1589d99fa3c5SJeremy L Thompson **/
1590d99fa3c5SJeremy L Thompson int CeedBasisGetTopology(CeedBasis basis, CeedElemTopology *topo) {
1591d99fa3c5SJeremy L Thompson   *topo = basis->topo;
1592e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1593d99fa3c5SJeremy L Thompson }
1594d99fa3c5SJeremy L Thompson 
1595d99fa3c5SJeremy L Thompson /**
15969d007619Sjeremylt   @brief Get number of components for given CeedBasis
15979d007619Sjeremylt 
1598ea61e9acSJeremy L Thompson   @param[in]  basis    CeedBasis
1599d1d35e2fSjeremylt   @param[out] num_comp Variable to store number of components of basis
16009d007619Sjeremylt 
16019d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
16029d007619Sjeremylt 
1603b7c9bbdaSJeremy L Thompson   @ref Advanced
16049d007619Sjeremylt **/
1605d1d35e2fSjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *num_comp) {
1606d1d35e2fSjeremylt   *num_comp = basis->num_comp;
1607e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
16089d007619Sjeremylt }
16099d007619Sjeremylt 
16109d007619Sjeremylt /**
16119d007619Sjeremylt   @brief Get total number of nodes (in dim dimensions) of a CeedBasis
16129d007619Sjeremylt 
1613ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
16149d007619Sjeremylt   @param[out] P     Variable to store number of nodes
16159d007619Sjeremylt 
16169d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
16179d007619Sjeremylt 
16189d007619Sjeremylt   @ref Utility
16199d007619Sjeremylt **/
16209d007619Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) {
16219d007619Sjeremylt   *P = basis->P;
1622e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
16239d007619Sjeremylt }
16249d007619Sjeremylt 
16259d007619Sjeremylt /**
16269d007619Sjeremylt   @brief Get total number of nodes (in 1 dimension) of a CeedBasis
16279d007619Sjeremylt 
1628ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
1629d1d35e2fSjeremylt   @param[out] P_1d  Variable to store number of nodes
16309d007619Sjeremylt 
16319d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
16329d007619Sjeremylt 
1633b7c9bbdaSJeremy L Thompson   @ref Advanced
16349d007619Sjeremylt **/
1635d1d35e2fSjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P_1d) {
16362b730f8bSJeremy L Thompson   if (!basis->tensor_basis) {
16379d007619Sjeremylt     // LCOV_EXCL_START
16382b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR, "Cannot supply P_1d for non-tensor basis");
16399d007619Sjeremylt     // LCOV_EXCL_STOP
16402b730f8bSJeremy L Thompson   }
16419d007619Sjeremylt 
1642d1d35e2fSjeremylt   *P_1d = basis->P_1d;
1643e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
16449d007619Sjeremylt }
16459d007619Sjeremylt 
16469d007619Sjeremylt /**
16479d007619Sjeremylt   @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis
16489d007619Sjeremylt 
1649ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
16509d007619Sjeremylt   @param[out] Q     Variable to store number of quadrature points
16519d007619Sjeremylt 
16529d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
16539d007619Sjeremylt 
16549d007619Sjeremylt   @ref Utility
16559d007619Sjeremylt **/
16569d007619Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) {
16579d007619Sjeremylt   *Q = basis->Q;
1658e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
16599d007619Sjeremylt }
16609d007619Sjeremylt 
16619d007619Sjeremylt /**
16629d007619Sjeremylt   @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis
16639d007619Sjeremylt 
1664ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
1665d1d35e2fSjeremylt   @param[out] Q_1d  Variable to store number of quadrature points
16669d007619Sjeremylt 
16679d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
16689d007619Sjeremylt 
1669b7c9bbdaSJeremy L Thompson   @ref Advanced
16709d007619Sjeremylt **/
1671d1d35e2fSjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q_1d) {
16722b730f8bSJeremy L Thompson   if (!basis->tensor_basis) {
16739d007619Sjeremylt     // LCOV_EXCL_START
16742b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR, "Cannot supply Q_1d for non-tensor basis");
16759d007619Sjeremylt     // LCOV_EXCL_STOP
16762b730f8bSJeremy L Thompson   }
16779d007619Sjeremylt 
1678d1d35e2fSjeremylt   *Q_1d = basis->Q_1d;
1679e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
16809d007619Sjeremylt }
16819d007619Sjeremylt 
16829d007619Sjeremylt /**
1683ea61e9acSJeremy L Thompson   @brief Get reference coordinates of quadrature points (in dim dimensions) of a CeedBasis
16849d007619Sjeremylt 
1685ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
1686d1d35e2fSjeremylt   @param[out] q_ref Variable to store reference coordinates of quadrature points
16879d007619Sjeremylt 
16889d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
16899d007619Sjeremylt 
1690b7c9bbdaSJeremy L Thompson   @ref Advanced
16919d007619Sjeremylt **/
1692d1d35e2fSjeremylt int CeedBasisGetQRef(CeedBasis basis, const CeedScalar **q_ref) {
1693d1d35e2fSjeremylt   *q_ref = basis->q_ref_1d;
1694e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
16959d007619Sjeremylt }
16969d007619Sjeremylt 
16979d007619Sjeremylt /**
1698ea61e9acSJeremy L Thompson   @brief Get quadrature weights of quadrature points (in dim dimensions) of a CeedBasis
16999d007619Sjeremylt 
1700ea61e9acSJeremy L Thompson   @param[in]  basis    CeedBasis
1701d1d35e2fSjeremylt   @param[out] q_weight Variable to store quadrature weights
17029d007619Sjeremylt 
17039d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
17049d007619Sjeremylt 
1705b7c9bbdaSJeremy L Thompson   @ref Advanced
17069d007619Sjeremylt **/
1707d1d35e2fSjeremylt int CeedBasisGetQWeights(CeedBasis basis, const CeedScalar **q_weight) {
1708d1d35e2fSjeremylt   *q_weight = basis->q_weight_1d;
1709e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
17109d007619Sjeremylt }
17119d007619Sjeremylt 
17129d007619Sjeremylt /**
17139d007619Sjeremylt   @brief Get interpolation matrix of a CeedBasis
17149d007619Sjeremylt 
1715ea61e9acSJeremy L Thompson   @param[in]  basis  CeedBasis
17169d007619Sjeremylt   @param[out] interp Variable to store interpolation matrix
17179d007619Sjeremylt 
17189d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
17199d007619Sjeremylt 
1720b7c9bbdaSJeremy L Thompson   @ref Advanced
17219d007619Sjeremylt **/
17226c58de82SJeremy L Thompson int CeedBasisGetInterp(CeedBasis basis, const CeedScalar **interp) {
1723d1d35e2fSjeremylt   if (!basis->interp && basis->tensor_basis) {
17249d007619Sjeremylt     // Allocate
17252b730f8bSJeremy L Thompson     CeedCall(CeedMalloc(basis->Q * basis->P, &basis->interp));
17269d007619Sjeremylt 
17279d007619Sjeremylt     // Initialize
17282b730f8bSJeremy L Thompson     for (CeedInt i = 0; i < basis->Q * basis->P; i++) basis->interp[i] = 1.0;
17299d007619Sjeremylt 
17309d007619Sjeremylt     // Calculate
17312b730f8bSJeremy L Thompson     for (CeedInt d = 0; d < basis->dim; d++) {
17322b730f8bSJeremy L Thompson       for (CeedInt qpt = 0; qpt < basis->Q; qpt++) {
17339d007619Sjeremylt         for (CeedInt node = 0; node < basis->P; node++) {
1734d1d35e2fSjeremylt           CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d;
1735d1d35e2fSjeremylt           CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d;
1736d1d35e2fSjeremylt           basis->interp[qpt * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p];
17379d007619Sjeremylt         }
17389d007619Sjeremylt       }
17392b730f8bSJeremy L Thompson     }
17402b730f8bSJeremy L Thompson   }
17419d007619Sjeremylt   *interp = basis->interp;
1742e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
17439d007619Sjeremylt }
17449d007619Sjeremylt 
17459d007619Sjeremylt /**
17469d007619Sjeremylt   @brief Get 1D interpolation matrix of a tensor product CeedBasis
17479d007619Sjeremylt 
1748ea61e9acSJeremy L Thompson   @param[in]  basis     CeedBasis
1749d1d35e2fSjeremylt   @param[out] interp_1d Variable to store interpolation matrix
17509d007619Sjeremylt 
17519d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
17529d007619Sjeremylt 
17539d007619Sjeremylt   @ref Backend
17549d007619Sjeremylt **/
1755d1d35e2fSjeremylt int CeedBasisGetInterp1D(CeedBasis basis, const CeedScalar **interp_1d) {
17562b730f8bSJeremy L Thompson   if (!basis->tensor_basis) {
17579d007619Sjeremylt     // LCOV_EXCL_START
17582b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR, "CeedBasis is not a tensor product basis.");
17599d007619Sjeremylt     // LCOV_EXCL_STOP
17602b730f8bSJeremy L Thompson   }
17619d007619Sjeremylt 
1762d1d35e2fSjeremylt   *interp_1d = basis->interp_1d;
1763e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
17649d007619Sjeremylt }
17659d007619Sjeremylt 
17669d007619Sjeremylt /**
17679d007619Sjeremylt   @brief Get gradient matrix of a CeedBasis
17689d007619Sjeremylt 
1769ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
17709d007619Sjeremylt   @param[out] grad  Variable to store gradient matrix
17719d007619Sjeremylt 
17729d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
17739d007619Sjeremylt 
1774b7c9bbdaSJeremy L Thompson   @ref Advanced
17759d007619Sjeremylt **/
17766c58de82SJeremy L Thompson int CeedBasisGetGrad(CeedBasis basis, const CeedScalar **grad) {
1777d1d35e2fSjeremylt   if (!basis->grad && basis->tensor_basis) {
17789d007619Sjeremylt     // Allocate
17792b730f8bSJeremy L Thompson     CeedCall(CeedMalloc(basis->dim * basis->Q * basis->P, &basis->grad));
17809d007619Sjeremylt 
17819d007619Sjeremylt     // Initialize
17822b730f8bSJeremy L Thompson     for (CeedInt i = 0; i < basis->dim * basis->Q * basis->P; i++) basis->grad[i] = 1.0;
17839d007619Sjeremylt 
17849d007619Sjeremylt     // Calculate
17852b730f8bSJeremy L Thompson     for (CeedInt d = 0; d < basis->dim; d++) {
17862b730f8bSJeremy L Thompson       for (CeedInt i = 0; i < basis->dim; i++) {
17872b730f8bSJeremy L Thompson         for (CeedInt qpt = 0; qpt < basis->Q; qpt++) {
17889d007619Sjeremylt           for (CeedInt node = 0; node < basis->P; node++) {
1789d1d35e2fSjeremylt             CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d;
1790d1d35e2fSjeremylt             CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d;
17912b730f8bSJeremy L Thompson             if (i == d) basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->grad_1d[q * basis->P_1d + p];
17922b730f8bSJeremy L Thompson             else basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p];
17932b730f8bSJeremy L Thompson           }
17942b730f8bSJeremy L Thompson         }
17952b730f8bSJeremy L Thompson       }
17969d007619Sjeremylt     }
17979d007619Sjeremylt   }
17989d007619Sjeremylt   *grad = basis->grad;
1799e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
18009d007619Sjeremylt }
18019d007619Sjeremylt 
18029d007619Sjeremylt /**
18039d007619Sjeremylt   @brief Get 1D gradient matrix of a tensor product CeedBasis
18049d007619Sjeremylt 
1805ea61e9acSJeremy L Thompson   @param[in]  basis   CeedBasis
1806d1d35e2fSjeremylt   @param[out] grad_1d Variable to store gradient matrix
18079d007619Sjeremylt 
18089d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
18099d007619Sjeremylt 
1810b7c9bbdaSJeremy L Thompson   @ref Advanced
18119d007619Sjeremylt **/
1812d1d35e2fSjeremylt int CeedBasisGetGrad1D(CeedBasis basis, const CeedScalar **grad_1d) {
18132b730f8bSJeremy L Thompson   if (!basis->tensor_basis) {
18149d007619Sjeremylt     // LCOV_EXCL_START
18152b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR, "CeedBasis is not a tensor product basis.");
18169d007619Sjeremylt     // LCOV_EXCL_STOP
18172b730f8bSJeremy L Thompson   }
18189d007619Sjeremylt 
1819d1d35e2fSjeremylt   *grad_1d = basis->grad_1d;
1820e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
18219d007619Sjeremylt }
18229d007619Sjeremylt 
18239d007619Sjeremylt /**
182450c301a5SRezgar Shakeri   @brief Get divergence matrix of a CeedBasis
182550c301a5SRezgar Shakeri 
1826ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
182750c301a5SRezgar Shakeri   @param[out] div   Variable to store divergence matrix
182850c301a5SRezgar Shakeri 
182950c301a5SRezgar Shakeri   @return An error code: 0 - success, otherwise - failure
183050c301a5SRezgar Shakeri 
183150c301a5SRezgar Shakeri   @ref Advanced
183250c301a5SRezgar Shakeri **/
183350c301a5SRezgar Shakeri int CeedBasisGetDiv(CeedBasis basis, const CeedScalar **div) {
18342b730f8bSJeremy L Thompson   if (!basis->div) {
183550c301a5SRezgar Shakeri     // LCOV_EXCL_START
18362b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR, "CeedBasis does not have divergence matrix.");
183750c301a5SRezgar Shakeri     // LCOV_EXCL_STOP
18382b730f8bSJeremy L Thompson   }
183950c301a5SRezgar Shakeri 
184050c301a5SRezgar Shakeri   *div = basis->div;
184150c301a5SRezgar Shakeri   return CEED_ERROR_SUCCESS;
184250c301a5SRezgar Shakeri }
184350c301a5SRezgar Shakeri 
184450c301a5SRezgar Shakeri /**
1845c4e3f59bSSebastian Grimberg   @brief Get curl matrix of a CeedBasis
1846c4e3f59bSSebastian Grimberg 
1847c4e3f59bSSebastian Grimberg   @param[in]  basis CeedBasis
1848c4e3f59bSSebastian Grimberg   @param[out] curl  Variable to store curl matrix
1849c4e3f59bSSebastian Grimberg 
1850c4e3f59bSSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
1851c4e3f59bSSebastian Grimberg 
1852c4e3f59bSSebastian Grimberg   @ref Advanced
1853c4e3f59bSSebastian Grimberg **/
1854c4e3f59bSSebastian Grimberg int CeedBasisGetCurl(CeedBasis basis, const CeedScalar **curl) {
1855c4e3f59bSSebastian Grimberg   if (!basis->curl) {
1856c4e3f59bSSebastian Grimberg     // LCOV_EXCL_START
1857c4e3f59bSSebastian Grimberg     return CeedError(basis->ceed, CEED_ERROR_MINOR, "CeedBasis does not have curl matrix.");
1858c4e3f59bSSebastian Grimberg     // LCOV_EXCL_STOP
1859c4e3f59bSSebastian Grimberg   }
1860c4e3f59bSSebastian Grimberg 
1861c4e3f59bSSebastian Grimberg   *curl = basis->curl;
1862c4e3f59bSSebastian Grimberg   return CEED_ERROR_SUCCESS;
1863c4e3f59bSSebastian Grimberg }
1864c4e3f59bSSebastian Grimberg 
1865c4e3f59bSSebastian Grimberg /**
18667a982d89SJeremy L. Thompson   @brief Destroy a CeedBasis
18677a982d89SJeremy L. Thompson 
1868ea61e9acSJeremy L Thompson   @param[in,out] basis CeedBasis to destroy
18697a982d89SJeremy L. Thompson 
18707a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
18717a982d89SJeremy L. Thompson 
18727a982d89SJeremy L. Thompson   @ref User
18737a982d89SJeremy L. Thompson **/
18747a982d89SJeremy L. Thompson int CeedBasisDestroy(CeedBasis *basis) {
18757425e127SJeremy L Thompson   if (!*basis || *basis == CEED_BASIS_COLLOCATED || --(*basis)->ref_count > 0) {
1876ad6481ceSJeremy L Thompson     *basis = NULL;
1877ad6481ceSJeremy L Thompson     return CEED_ERROR_SUCCESS;
1878ad6481ceSJeremy L Thompson   }
18792b730f8bSJeremy L Thompson   if ((*basis)->Destroy) CeedCall((*basis)->Destroy(*basis));
18802b730f8bSJeremy L Thompson   if ((*basis)->contract) CeedCall(CeedTensorContractDestroy(&(*basis)->contract));
1881c4e3f59bSSebastian Grimberg   CeedCall(CeedFree(&(*basis)->q_ref_1d));
1882c4e3f59bSSebastian Grimberg   CeedCall(CeedFree(&(*basis)->q_weight_1d));
18832b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->interp));
18842b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->interp_1d));
18852b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->grad));
18862b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->grad_1d));
1887c4e3f59bSSebastian Grimberg   CeedCall(CeedFree(&(*basis)->div));
1888c4e3f59bSSebastian Grimberg   CeedCall(CeedFree(&(*basis)->curl));
18892b730f8bSJeremy L Thompson   CeedCall(CeedDestroy(&(*basis)->ceed));
18902b730f8bSJeremy L Thompson   CeedCall(CeedFree(basis));
1891e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
18927a982d89SJeremy L. Thompson }
18937a982d89SJeremy L. Thompson 
18947a982d89SJeremy L. Thompson /**
1895b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre quadrature
1896b11c1e72Sjeremylt 
1897ea61e9acSJeremy L Thompson   @param[in]  Q           Number of quadrature points (integrates polynomials of degree 2*Q-1 exactly)
1898d1d35e2fSjeremylt   @param[out] q_ref_1d    Array of length Q to hold the abscissa on [-1, 1]
1899d1d35e2fSjeremylt   @param[out] q_weight_1d Array of length Q to hold the weights
1900b11c1e72Sjeremylt 
1901b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1902dfdf5a53Sjeremylt 
1903dfdf5a53Sjeremylt   @ref Utility
1904b11c1e72Sjeremylt **/
19052b730f8bSJeremy L Thompson int CeedGaussQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) {
1906d7b241e6Sjeremylt   // Allocate
1907d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0 * atan(1.0);
1908d1d35e2fSjeremylt   // Build q_ref_1d, q_weight_1d
190992ae7e47SJeremy L Thompson   for (CeedInt i = 0; i <= Q / 2; i++) {
1910d7b241e6Sjeremylt     // Guess
1911d7b241e6Sjeremylt     xi = cos(PI * (CeedScalar)(2 * i + 1) / ((CeedScalar)(2 * Q)));
1912d7b241e6Sjeremylt     // Pn(xi)
1913d7b241e6Sjeremylt     P0 = 1.0;
1914d7b241e6Sjeremylt     P1 = xi;
1915d7b241e6Sjeremylt     P2 = 0.0;
191692ae7e47SJeremy L Thompson     for (CeedInt j = 2; j <= Q; j++) {
1917d7b241e6Sjeremylt       P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
1918d7b241e6Sjeremylt       P0 = P1;
1919d7b241e6Sjeremylt       P1 = P2;
1920d7b241e6Sjeremylt     }
1921d7b241e6Sjeremylt     // First Newton Step
1922d7b241e6Sjeremylt     dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
1923d7b241e6Sjeremylt     xi  = xi - P2 / dP2;
1924d7b241e6Sjeremylt     // Newton to convergence
192592ae7e47SJeremy L Thompson     for (CeedInt k = 0; k < 100 && fabs(P2) > 10 * CEED_EPSILON; k++) {
1926d7b241e6Sjeremylt       P0 = 1.0;
1927d7b241e6Sjeremylt       P1 = xi;
192892ae7e47SJeremy L Thompson       for (CeedInt j = 2; j <= Q; j++) {
1929d7b241e6Sjeremylt         P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
1930d7b241e6Sjeremylt         P0 = P1;
1931d7b241e6Sjeremylt         P1 = P2;
1932d7b241e6Sjeremylt       }
1933d7b241e6Sjeremylt       dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
1934d7b241e6Sjeremylt       xi  = xi - P2 / dP2;
1935d7b241e6Sjeremylt     }
1936d7b241e6Sjeremylt     // Save xi, wi
1937d7b241e6Sjeremylt     wi                     = 2.0 / ((1.0 - xi * xi) * dP2 * dP2);
1938d1d35e2fSjeremylt     q_weight_1d[i]         = wi;
1939d1d35e2fSjeremylt     q_weight_1d[Q - 1 - i] = wi;
1940d1d35e2fSjeremylt     q_ref_1d[i]            = -xi;
1941d1d35e2fSjeremylt     q_ref_1d[Q - 1 - i]    = xi;
1942d7b241e6Sjeremylt   }
1943e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1944d7b241e6Sjeremylt }
1945d7b241e6Sjeremylt 
1946b11c1e72Sjeremylt /**
1947b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre-Lobatto quadrature
1948b11c1e72Sjeremylt 
1949ea61e9acSJeremy L Thompson   @param[in]  Q           Number of quadrature points (integrates polynomials of degree 2*Q-3 exactly)
1950d1d35e2fSjeremylt   @param[out] q_ref_1d    Array of length Q to hold the abscissa on [-1, 1]
1951d1d35e2fSjeremylt   @param[out] q_weight_1d Array of length Q to hold the weights
1952b11c1e72Sjeremylt 
1953b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1954dfdf5a53Sjeremylt 
1955dfdf5a53Sjeremylt   @ref Utility
1956b11c1e72Sjeremylt **/
19572b730f8bSJeremy L Thompson int CeedLobattoQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) {
1958d7b241e6Sjeremylt   // Allocate
1959d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0 * atan(1.0);
1960d1d35e2fSjeremylt   // Build q_ref_1d, q_weight_1d
1961d7b241e6Sjeremylt   // Set endpoints
19622b730f8bSJeremy L Thompson   if (Q < 2) {
1963b0d62198Sjeremylt     // LCOV_EXCL_START
19642b730f8bSJeremy L Thompson     return CeedError(NULL, CEED_ERROR_DIMENSION, "Cannot create Lobatto quadrature with Q=%" CeedInt_FMT " < 2 points", Q);
1965b0d62198Sjeremylt     // LCOV_EXCL_STOP
19662b730f8bSJeremy L Thompson   }
1967d7b241e6Sjeremylt   wi = 2.0 / ((CeedScalar)(Q * (Q - 1)));
1968d1d35e2fSjeremylt   if (q_weight_1d) {
1969d1d35e2fSjeremylt     q_weight_1d[0]     = wi;
1970d1d35e2fSjeremylt     q_weight_1d[Q - 1] = wi;
1971d7b241e6Sjeremylt   }
1972d1d35e2fSjeremylt   q_ref_1d[0]     = -1.0;
1973d1d35e2fSjeremylt   q_ref_1d[Q - 1] = 1.0;
1974d7b241e6Sjeremylt   // Interior
197592ae7e47SJeremy L Thompson   for (CeedInt i = 1; i <= (Q - 1) / 2; i++) {
1976d7b241e6Sjeremylt     // Guess
1977d7b241e6Sjeremylt     xi = cos(PI * (CeedScalar)(i) / (CeedScalar)(Q - 1));
1978d7b241e6Sjeremylt     // Pn(xi)
1979d7b241e6Sjeremylt     P0 = 1.0;
1980d7b241e6Sjeremylt     P1 = xi;
1981d7b241e6Sjeremylt     P2 = 0.0;
198292ae7e47SJeremy L Thompson     for (CeedInt j = 2; j < Q; j++) {
1983d7b241e6Sjeremylt       P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
1984d7b241e6Sjeremylt       P0 = P1;
1985d7b241e6Sjeremylt       P1 = P2;
1986d7b241e6Sjeremylt     }
1987d7b241e6Sjeremylt     // First Newton step
1988d7b241e6Sjeremylt     dP2  = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
1989d7b241e6Sjeremylt     d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi);
1990d7b241e6Sjeremylt     xi   = xi - dP2 / d2P2;
1991d7b241e6Sjeremylt     // Newton to convergence
199292ae7e47SJeremy L Thompson     for (CeedInt k = 0; k < 100 && fabs(dP2) > 10 * CEED_EPSILON; k++) {
1993d7b241e6Sjeremylt       P0 = 1.0;
1994d7b241e6Sjeremylt       P1 = xi;
199592ae7e47SJeremy L Thompson       for (CeedInt j = 2; j < Q; j++) {
1996d7b241e6Sjeremylt         P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
1997d7b241e6Sjeremylt         P0 = P1;
1998d7b241e6Sjeremylt         P1 = P2;
1999d7b241e6Sjeremylt       }
2000d7b241e6Sjeremylt       dP2  = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
2001d7b241e6Sjeremylt       d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi);
2002d7b241e6Sjeremylt       xi   = xi - dP2 / d2P2;
2003d7b241e6Sjeremylt     }
2004d7b241e6Sjeremylt     // Save xi, wi
2005d7b241e6Sjeremylt     wi = 2.0 / (((CeedScalar)(Q * (Q - 1))) * P2 * P2);
2006d1d35e2fSjeremylt     if (q_weight_1d) {
2007d1d35e2fSjeremylt       q_weight_1d[i]         = wi;
2008d1d35e2fSjeremylt       q_weight_1d[Q - 1 - i] = wi;
2009d7b241e6Sjeremylt     }
2010d1d35e2fSjeremylt     q_ref_1d[i]         = -xi;
2011d1d35e2fSjeremylt     q_ref_1d[Q - 1 - i] = xi;
2012d7b241e6Sjeremylt   }
2013e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2014d7b241e6Sjeremylt }
2015d7b241e6Sjeremylt 
2016d7b241e6Sjeremylt /// @}
2017