xref: /libCEED/interface/ceed-basis.c (revision c4e3f59b2ea5a0c95cc0118aa5026c447cce3092)
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 
127ba59ac12SSebastian Grimberg   The interpolation is given by `interp_project = interp_to^+ * interp_from`, where the pesudoinverse `interp_to^+` is given by QR factorization.
12814556e63SJeremy L Thompson   The gradient is given by `grad_project = interp_to^+ * grad_from`.
129ba59ac12SSebastian Grimberg   Note: `basis_from` and `basis_to` must have compatible quadrature spaces.
130a76a04e7SJeremy L Thompson 
131a76a04e7SJeremy L Thompson   @param[in]  basis_from     CeedBasis to project from
132a76a04e7SJeremy L Thompson   @param[in]  basis_to       CeedBasis to project to
133ea61e9acSJeremy L Thompson   @param[out] interp_project Address of the variable where the newly created interpolation matrix will be stored.
134ea61e9acSJeremy L Thompson   @param[out] grad_project   Address of the variable where the newly created gradient matrix will be stored.
135a76a04e7SJeremy L Thompson 
136a76a04e7SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
137a76a04e7SJeremy L Thompson 
138a76a04e7SJeremy L Thompson   @ref Developer
139a76a04e7SJeremy L Thompson **/
1402b730f8bSJeremy L Thompson static int CeedBasisCreateProjectionMatrices(CeedBasis basis_from, CeedBasis basis_to, CeedScalar **interp_project, CeedScalar **grad_project) {
141a76a04e7SJeremy L Thompson   Ceed ceed;
1422b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetCeed(basis_to, &ceed));
143a76a04e7SJeremy L Thompson 
144a76a04e7SJeremy L Thompson   // Check for compatible quadrature spaces
145a76a04e7SJeremy L Thompson   CeedInt Q_to, Q_from;
1462b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumQuadraturePoints(basis_to, &Q_to));
1472b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumQuadraturePoints(basis_from, &Q_from));
1482b730f8bSJeremy L Thompson   if (Q_to != Q_from) {
149a76a04e7SJeremy L Thompson     // LCOV_EXCL_START
1502b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Bases must have compatible quadrature spaces");
151a76a04e7SJeremy L Thompson     // LCOV_EXCL_STOP
1522b730f8bSJeremy L Thompson   }
153a76a04e7SJeremy L Thompson 
15414556e63SJeremy L Thompson   // Check for matching tensor or non-tensor
155a76a04e7SJeremy L Thompson   CeedInt P_to, P_from, Q = Q_to;
156a76a04e7SJeremy L Thompson   bool    is_tensor_to, is_tensor_from;
1572b730f8bSJeremy L Thompson   CeedCall(CeedBasisIsTensor(basis_to, &is_tensor_to));
1582b730f8bSJeremy L Thompson   CeedCall(CeedBasisIsTensor(basis_from, &is_tensor_from));
159a76a04e7SJeremy L Thompson   if (is_tensor_to && is_tensor_from) {
1602b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_to));
1612b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_from));
1622b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumQuadraturePoints1D(basis_from, &Q));
163a76a04e7SJeremy L Thompson   } else if (!is_tensor_to && !is_tensor_from) {
1642b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis_to, &P_to));
1652b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis_from, &P_from));
166a76a04e7SJeremy L Thompson   } else {
167a76a04e7SJeremy L Thompson     // LCOV_EXCL_START
1682b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_MINOR, "Bases must both be tensor or non-tensor");
169a76a04e7SJeremy L Thompson     // LCOV_EXCL_STOP
170a76a04e7SJeremy L Thompson   }
171a76a04e7SJeremy L Thompson 
17214556e63SJeremy L Thompson   // Get source matrices
17314556e63SJeremy L Thompson   CeedInt     dim;
17414556e63SJeremy L Thompson   CeedScalar *interp_to, *interp_from, *tau;
1752b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetDimension(basis_to, &dim));
1762b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q * P_from, &interp_from));
1772b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q * P_to, &interp_to));
1782b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(P_to * P_from, interp_project));
1792b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(P_to * P_from * (is_tensor_to ? 1 : dim), grad_project));
1802b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q, &tau));
1812b730f8bSJeremy L Thompson   const CeedScalar *interp_to_source = NULL, *interp_from_source = NULL, *grad_from_source;
182a76a04e7SJeremy L Thompson   if (is_tensor_to) {
1832b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetInterp1D(basis_to, &interp_to_source));
1842b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetInterp1D(basis_from, &interp_from_source));
1852b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetGrad1D(basis_from, &grad_from_source));
186a76a04e7SJeremy L Thompson   } else {
1872b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetInterp(basis_to, &interp_to_source));
1882b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetInterp(basis_from, &interp_from_source));
1892b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetGrad(basis_from, &grad_from_source));
190a76a04e7SJeremy L Thompson   }
191a76a04e7SJeremy L Thompson 
19214556e63SJeremy L Thompson   // Build matrices
19314556e63SJeremy L Thompson   CeedInt     num_matrices = 1 + (is_tensor_to ? 1 : dim);
19414556e63SJeremy L Thompson   CeedScalar *input_from[num_matrices], *output_project[num_matrices];
19514556e63SJeremy L Thompson   input_from[0]     = (CeedScalar *)interp_from_source;
19614556e63SJeremy L Thompson   output_project[0] = *interp_project;
19714556e63SJeremy L Thompson   for (CeedInt m = 1; m < num_matrices; m++) {
19814556e63SJeremy L Thompson     input_from[m]     = (CeedScalar *)&grad_from_source[(m - 1) * Q * P_from];
19902af4036SJeremy L Thompson     output_project[m] = &((*grad_project)[(m - 1) * P_to * P_from]);
20014556e63SJeremy L Thompson   }
20114556e63SJeremy L Thompson   for (CeedInt m = 0; m < num_matrices; m++) {
202a76a04e7SJeremy L Thompson     // -- QR Factorization, interp_to = Q R
20314556e63SJeremy L Thompson     memcpy(interp_to, interp_to_source, Q * P_to * sizeof(interp_to_source[0]));
2042b730f8bSJeremy L Thompson     CeedCall(CeedQRFactorization(ceed, interp_to, tau, Q, P_to));
205a76a04e7SJeremy L Thompson 
206a76a04e7SJeremy L Thompson     // -- Apply Qtranspose, interp_to = Qtranspose interp_from
20714556e63SJeremy L Thompson     memcpy(interp_from, input_from[m], Q * P_from * sizeof(input_from[m][0]));
2082b730f8bSJeremy L Thompson     CeedCall(CeedHouseholderApplyQ(interp_from, interp_to, tau, CEED_TRANSPOSE, Q, P_from, P_to, P_from, 1));
209a76a04e7SJeremy L Thompson 
210a76a04e7SJeremy L Thompson     // -- Apply Rinv, interp_project = Rinv interp_c
211a76a04e7SJeremy L Thompson     for (CeedInt j = 0; j < P_from; j++) {  // Column j
2122b730f8bSJeremy 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];
213a76a04e7SJeremy L Thompson       for (CeedInt i = P_to - 2; i >= 0; i--) {  // Row i
21414556e63SJeremy L Thompson         output_project[m][j + P_from * i] = interp_from[j + P_from * i];
215a76a04e7SJeremy L Thompson         for (CeedInt k = i + 1; k < P_to; k++) {
2162b730f8bSJeremy L Thompson           output_project[m][j + P_from * i] -= interp_to[k + P_to * i] * output_project[m][j + P_from * k];
217a76a04e7SJeremy L Thompson         }
21814556e63SJeremy L Thompson         output_project[m][j + P_from * i] /= interp_to[i + P_to * i];
219a76a04e7SJeremy L Thompson       }
220a76a04e7SJeremy L Thompson     }
22114556e63SJeremy L Thompson   }
22214556e63SJeremy L Thompson 
22314556e63SJeremy L Thompson   // Cleanup
2242b730f8bSJeremy L Thompson   CeedCall(CeedFree(&tau));
2252b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_to));
2262b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_from));
227a76a04e7SJeremy L Thompson 
228a76a04e7SJeremy L Thompson   return CEED_ERROR_SUCCESS;
229a76a04e7SJeremy L Thompson }
230a76a04e7SJeremy L Thompson 
2317a982d89SJeremy L. Thompson /// @}
2327a982d89SJeremy L. Thompson 
2337a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
2347a982d89SJeremy L. Thompson /// Ceed Backend API
2357a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
2367a982d89SJeremy L. Thompson /// @addtogroup CeedBasisBackend
2377a982d89SJeremy L. Thompson /// @{
2387a982d89SJeremy L. Thompson 
2397a982d89SJeremy L. Thompson /**
2407a982d89SJeremy L. Thompson   @brief Return collocated grad matrix
2417a982d89SJeremy L. Thompson 
242ea61e9acSJeremy L Thompson   @param[in]  basis         CeedBasis
243ea61e9acSJeremy L Thompson   @param[out] collo_grad_1d Row-major (Q_1d * Q_1d) matrix expressing derivatives of basis functions at quadrature points
2447a982d89SJeremy L. Thompson 
2457a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
2467a982d89SJeremy L. Thompson 
2477a982d89SJeremy L. Thompson   @ref Backend
2487a982d89SJeremy L. Thompson **/
249d1d35e2fSjeremylt int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *collo_grad_1d) {
2507a982d89SJeremy L. Thompson   int         i, j, k;
2517a982d89SJeremy L. Thompson   Ceed        ceed;
2522b730f8bSJeremy L Thompson   CeedInt     P_1d = (basis)->P_1d, Q_1d = (basis)->Q_1d;
25378464608Sjeremylt   CeedScalar *interp_1d, *grad_1d, *tau;
2547a982d89SJeremy L. Thompson 
2552b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q_1d * P_1d, &interp_1d));
2562b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q_1d * P_1d, &grad_1d));
2572b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q_1d, &tau));
258d1d35e2fSjeremylt   memcpy(interp_1d, (basis)->interp_1d, Q_1d * P_1d * sizeof(basis)->interp_1d[0]);
259d1d35e2fSjeremylt   memcpy(grad_1d, (basis)->grad_1d, Q_1d * P_1d * sizeof(basis)->interp_1d[0]);
2607a982d89SJeremy L. Thompson 
261d1d35e2fSjeremylt   // QR Factorization, interp_1d = Q R
2622b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetCeed(basis, &ceed));
2632b730f8bSJeremy L Thompson   CeedCall(CeedQRFactorization(ceed, interp_1d, tau, Q_1d, P_1d));
264ea61e9acSJeremy 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.
2657a982d89SJeremy L. Thompson 
266d1d35e2fSjeremylt   // Apply Rinv, collo_grad_1d = grad_1d Rinv
267d1d35e2fSjeremylt   for (i = 0; i < Q_1d; i++) {  // Row i
268d1d35e2fSjeremylt     collo_grad_1d[Q_1d * i] = grad_1d[P_1d * i] / interp_1d[0];
269d1d35e2fSjeremylt     for (j = 1; j < P_1d; j++) {  // Column j
270d1d35e2fSjeremylt       collo_grad_1d[j + Q_1d * i] = grad_1d[j + P_1d * i];
2712b730f8bSJeremy 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];
272d1d35e2fSjeremylt       collo_grad_1d[j + Q_1d * i] /= interp_1d[j + P_1d * j];
2737a982d89SJeremy L. Thompson     }
2742b730f8bSJeremy L Thompson     for (j = P_1d; j < Q_1d; j++) collo_grad_1d[j + Q_1d * i] = 0;
2757a982d89SJeremy L. Thompson   }
2767a982d89SJeremy L. Thompson 
277673160d7Sjeremylt   // Apply Qtranspose, collo_grad = collo_grad Q_transpose
2782b730f8bSJeremy L Thompson   CeedCall(CeedHouseholderApplyQ(collo_grad_1d, interp_1d, tau, CEED_NOTRANSPOSE, Q_1d, Q_1d, P_1d, 1, Q_1d));
2797a982d89SJeremy L. Thompson 
2802b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_1d));
2812b730f8bSJeremy L Thompson   CeedCall(CeedFree(&grad_1d));
2822b730f8bSJeremy L Thompson   CeedCall(CeedFree(&tau));
283e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2847a982d89SJeremy L. Thompson }
2857a982d89SJeremy L. Thompson 
2867a982d89SJeremy L. Thompson /**
2877a982d89SJeremy L. Thompson   @brief Get tensor status for given CeedBasis
2887a982d89SJeremy L. Thompson 
289ea61e9acSJeremy L Thompson   @param[in]  basis     CeedBasis
290d1d35e2fSjeremylt   @param[out] is_tensor Variable to store tensor status
2917a982d89SJeremy L. Thompson 
2927a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
2937a982d89SJeremy L. Thompson 
2947a982d89SJeremy L. Thompson   @ref Backend
2957a982d89SJeremy L. Thompson **/
296d1d35e2fSjeremylt int CeedBasisIsTensor(CeedBasis basis, bool *is_tensor) {
297d1d35e2fSjeremylt   *is_tensor = basis->tensor_basis;
298e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2997a982d89SJeremy L. Thompson }
3007a982d89SJeremy L. Thompson 
3017a982d89SJeremy L. Thompson /**
3027a982d89SJeremy L. Thompson   @brief Get backend data of a CeedBasis
3037a982d89SJeremy L. Thompson 
304ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
3057a982d89SJeremy L. Thompson   @param[out] data  Variable to store data
3067a982d89SJeremy L. Thompson 
3077a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3087a982d89SJeremy L. Thompson 
3097a982d89SJeremy L. Thompson   @ref Backend
3107a982d89SJeremy L. Thompson **/
311777ff853SJeremy L Thompson int CeedBasisGetData(CeedBasis basis, void *data) {
312777ff853SJeremy L Thompson   *(void **)data = basis->data;
313e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3147a982d89SJeremy L. Thompson }
3157a982d89SJeremy L. Thompson 
3167a982d89SJeremy L. Thompson /**
3177a982d89SJeremy L. Thompson   @brief Set backend data of a CeedBasis
3187a982d89SJeremy L. Thompson 
319ea61e9acSJeremy L Thompson   @param[in,out] basis  CeedBasis
320ea61e9acSJeremy L Thompson   @param[in]     data   Data to set
3217a982d89SJeremy L. Thompson 
3227a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
3237a982d89SJeremy L. Thompson 
3247a982d89SJeremy L. Thompson   @ref Backend
3257a982d89SJeremy L. Thompson **/
326777ff853SJeremy L Thompson int CeedBasisSetData(CeedBasis basis, void *data) {
327777ff853SJeremy L Thompson   basis->data = data;
328e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3297a982d89SJeremy L. Thompson }
3307a982d89SJeremy L. Thompson 
3317a982d89SJeremy L. Thompson /**
33234359f16Sjeremylt   @brief Increment the reference counter for a CeedBasis
33334359f16Sjeremylt 
334ea61e9acSJeremy L Thompson   @param[in,out] basis Basis to increment the reference counter
33534359f16Sjeremylt 
33634359f16Sjeremylt   @return An error code: 0 - success, otherwise - failure
33734359f16Sjeremylt 
33834359f16Sjeremylt   @ref Backend
33934359f16Sjeremylt **/
3409560d06aSjeremylt int CeedBasisReference(CeedBasis basis) {
34134359f16Sjeremylt   basis->ref_count++;
34234359f16Sjeremylt   return CEED_ERROR_SUCCESS;
34334359f16Sjeremylt }
34434359f16Sjeremylt 
34534359f16Sjeremylt /**
346*c4e3f59bSSebastian Grimberg   @brief Get number of Q-vector components for given CeedBasis
347*c4e3f59bSSebastian Grimberg 
348*c4e3f59bSSebastian Grimberg   @param[in]  basis  CeedBasis
349*c4e3f59bSSebastian Grimberg   @param[in]  eval_mode \ref CEED_EVAL_INTERP to use interpolated values,
350*c4e3f59bSSebastian Grimberg                           \ref CEED_EVAL_GRAD to use gradients,
351*c4e3f59bSSebastian Grimberg                           \ref CEED_EVAL_DIV to use divergence,
352*c4e3f59bSSebastian Grimberg                           \ref CEED_EVAL_CURL to use curl.
353*c4e3f59bSSebastian Grimberg   @param[out] q_comp Variable to store number of Q-vector components of basis
354*c4e3f59bSSebastian Grimberg 
355*c4e3f59bSSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
356*c4e3f59bSSebastian Grimberg 
357*c4e3f59bSSebastian Grimberg   @ref Backend
358*c4e3f59bSSebastian Grimberg **/
359*c4e3f59bSSebastian Grimberg int CeedBasisGetNumQuadratureComponents(CeedBasis basis, CeedEvalMode eval_mode, CeedInt *q_comp) {
360*c4e3f59bSSebastian Grimberg   switch (eval_mode) {
361*c4e3f59bSSebastian Grimberg     case CEED_EVAL_INTERP:
362*c4e3f59bSSebastian Grimberg       *q_comp = (basis->fe_space == CEED_FE_SPACE_H1) ? 1 : basis->dim;
363*c4e3f59bSSebastian Grimberg       break;
364*c4e3f59bSSebastian Grimberg     case CEED_EVAL_GRAD:
365*c4e3f59bSSebastian Grimberg       *q_comp = basis->dim;
366*c4e3f59bSSebastian Grimberg       break;
367*c4e3f59bSSebastian Grimberg     case CEED_EVAL_DIV:
368*c4e3f59bSSebastian Grimberg       *q_comp = 1;
369*c4e3f59bSSebastian Grimberg       break;
370*c4e3f59bSSebastian Grimberg     case CEED_EVAL_CURL:
371*c4e3f59bSSebastian Grimberg       *q_comp = (basis->dim < 3) ? 1 : basis->dim;
372*c4e3f59bSSebastian Grimberg       break;
373*c4e3f59bSSebastian Grimberg     case CEED_EVAL_NONE:
374*c4e3f59bSSebastian Grimberg     case CEED_EVAL_WEIGHT:
375*c4e3f59bSSebastian Grimberg       *q_comp = 0;
376*c4e3f59bSSebastian Grimberg       break;
377*c4e3f59bSSebastian Grimberg   }
378*c4e3f59bSSebastian Grimberg   return CEED_ERROR_SUCCESS;
379*c4e3f59bSSebastian Grimberg }
380*c4e3f59bSSebastian Grimberg 
381*c4e3f59bSSebastian Grimberg /**
3826e15d496SJeremy L Thompson   @brief Estimate number of FLOPs required to apply CeedBasis in t_mode and eval_mode
3836e15d496SJeremy L Thompson 
384ea61e9acSJeremy L Thompson   @param[in]  basis     Basis to estimate FLOPs for
385ea61e9acSJeremy L Thompson   @param[in]  t_mode    Apply basis or transpose
386ea61e9acSJeremy L Thompson   @param[in]  eval_mode Basis evaluation mode
387ea61e9acSJeremy L Thompson   @param[out] flops     Address of variable to hold FLOPs estimate
3886e15d496SJeremy L Thompson 
3896e15d496SJeremy L Thompson   @ref Backend
3906e15d496SJeremy L Thompson **/
3912b730f8bSJeremy L Thompson int CeedBasisGetFlopsEstimate(CeedBasis basis, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedSize *flops) {
3926e15d496SJeremy L Thompson   bool is_tensor;
3936e15d496SJeremy L Thompson 
3942b730f8bSJeremy L Thompson   CeedCall(CeedBasisIsTensor(basis, &is_tensor));
3956e15d496SJeremy L Thompson   if (is_tensor) {
3966e15d496SJeremy L Thompson     CeedInt dim, num_comp, P_1d, Q_1d;
3972b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetDimension(basis, &dim));
3982b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
3992b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d));
4002b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
4016e15d496SJeremy L Thompson     if (t_mode == CEED_TRANSPOSE) {
4022b730f8bSJeremy L Thompson       P_1d = Q_1d;
4032b730f8bSJeremy L Thompson       Q_1d = P_1d;
4046e15d496SJeremy L Thompson     }
4056e15d496SJeremy L Thompson     CeedInt tensor_flops = 0, pre = num_comp * CeedIntPow(P_1d, dim - 1), post = 1;
4066e15d496SJeremy L Thompson     for (CeedInt d = 0; d < dim; d++) {
4076e15d496SJeremy L Thompson       tensor_flops += 2 * pre * P_1d * post * Q_1d;
4086e15d496SJeremy L Thompson       pre /= P_1d;
4096e15d496SJeremy L Thompson       post *= Q_1d;
4106e15d496SJeremy L Thompson     }
4116e15d496SJeremy L Thompson     switch (eval_mode) {
4122b730f8bSJeremy L Thompson       case CEED_EVAL_NONE:
4132b730f8bSJeremy L Thompson         *flops = 0;
4142b730f8bSJeremy L Thompson         break;
4152b730f8bSJeremy L Thompson       case CEED_EVAL_INTERP:
4162b730f8bSJeremy L Thompson         *flops = tensor_flops;
4172b730f8bSJeremy L Thompson         break;
4182b730f8bSJeremy L Thompson       case CEED_EVAL_GRAD:
4192b730f8bSJeremy L Thompson         *flops = tensor_flops * 2;
4202b730f8bSJeremy L Thompson         break;
4216e15d496SJeremy L Thompson       case CEED_EVAL_DIV:
4226e15d496SJeremy L Thompson         // LCOV_EXCL_START
4232b730f8bSJeremy L Thompson         return CeedError(basis->ceed, CEED_ERROR_INCOMPATIBLE, "Tensor CEED_EVAL_DIV not supported");
4242b730f8bSJeremy L Thompson         break;
4256e15d496SJeremy L Thompson       case CEED_EVAL_CURL:
4262b730f8bSJeremy L Thompson         return CeedError(basis->ceed, CEED_ERROR_INCOMPATIBLE, "Tensor CEED_EVAL_CURL not supported");
4272b730f8bSJeremy L Thompson         break;
4286e15d496SJeremy L Thompson       // LCOV_EXCL_STOP
4292b730f8bSJeremy L Thompson       case CEED_EVAL_WEIGHT:
4302b730f8bSJeremy L Thompson         *flops = dim * CeedIntPow(Q_1d, dim);
4312b730f8bSJeremy L Thompson         break;
4326e15d496SJeremy L Thompson     }
4336e15d496SJeremy L Thompson   } else {
434*c4e3f59bSSebastian Grimberg     CeedInt dim, num_comp, q_comp, num_nodes, num_qpts;
4352b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetDimension(basis, &dim));
4362b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
437*c4e3f59bSSebastian Grimberg     CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp));
4382b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis, &num_nodes));
4392b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts));
4406e15d496SJeremy L Thompson     switch (eval_mode) {
4412b730f8bSJeremy L Thompson       case CEED_EVAL_NONE:
4422b730f8bSJeremy L Thompson         *flops = 0;
4432b730f8bSJeremy L Thompson         break;
4442b730f8bSJeremy L Thompson       case CEED_EVAL_INTERP:
4452b730f8bSJeremy L Thompson       case CEED_EVAL_GRAD:
4462b730f8bSJeremy L Thompson       case CEED_EVAL_DIV:
4472b730f8bSJeremy L Thompson       case CEED_EVAL_CURL:
448*c4e3f59bSSebastian Grimberg         *flops = num_nodes * num_qpts * num_comp * q_comp;
4492b730f8bSJeremy L Thompson         break;
4502b730f8bSJeremy L Thompson       case CEED_EVAL_WEIGHT:
4512b730f8bSJeremy L Thompson         *flops = 0;
4522b730f8bSJeremy L Thompson         break;
4536e15d496SJeremy L Thompson     }
4546e15d496SJeremy L Thompson   }
4556e15d496SJeremy L Thompson 
4566e15d496SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4576e15d496SJeremy L Thompson }
4586e15d496SJeremy L Thompson 
4596e15d496SJeremy L Thompson /**
460*c4e3f59bSSebastian Grimberg   @brief Get CeedFESpace for a CeedBasis
461*c4e3f59bSSebastian Grimberg 
462*c4e3f59bSSebastian Grimberg   @param[in]  basis     CeedBasis
463*c4e3f59bSSebastian Grimberg   @param[out] fe_space  Variable to store CeedFESpace
464*c4e3f59bSSebastian Grimberg 
465*c4e3f59bSSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
466*c4e3f59bSSebastian Grimberg 
467*c4e3f59bSSebastian Grimberg   @ref Backend
468*c4e3f59bSSebastian Grimberg **/
469*c4e3f59bSSebastian Grimberg int CeedBasisGetFESpace(CeedBasis basis, CeedFESpace *fe_space) {
470*c4e3f59bSSebastian Grimberg   *fe_space = basis->fe_space;
471*c4e3f59bSSebastian Grimberg   return CEED_ERROR_SUCCESS;
472*c4e3f59bSSebastian Grimberg }
473*c4e3f59bSSebastian Grimberg 
474*c4e3f59bSSebastian Grimberg /**
4757a982d89SJeremy L. Thompson   @brief Get dimension for given CeedElemTopology
4767a982d89SJeremy L. Thompson 
477ea61e9acSJeremy L Thompson   @param[in]  topo CeedElemTopology
4787a982d89SJeremy L. Thompson   @param[out] dim  Variable to store dimension of topology
4797a982d89SJeremy L. Thompson 
4807a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
4817a982d89SJeremy L. Thompson 
4827a982d89SJeremy L. Thompson   @ref Backend
4837a982d89SJeremy L. Thompson **/
4847a982d89SJeremy L. Thompson int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) {
4857a982d89SJeremy L. Thompson   *dim = (CeedInt)topo >> 16;
486e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4877a982d89SJeremy L. Thompson }
4887a982d89SJeremy L. Thompson 
4897a982d89SJeremy L. Thompson /**
4907a982d89SJeremy L. Thompson   @brief Get CeedTensorContract of a CeedBasis
4917a982d89SJeremy L. Thompson 
492ea61e9acSJeremy L Thompson   @param[in]  basis     CeedBasis
4937a982d89SJeremy L. Thompson   @param[out] contract  Variable to store CeedTensorContract
4947a982d89SJeremy L. Thompson 
4957a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
4967a982d89SJeremy L. Thompson 
4977a982d89SJeremy L. Thompson   @ref Backend
4987a982d89SJeremy L. Thompson **/
4997a982d89SJeremy L. Thompson int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) {
5007a982d89SJeremy L. Thompson   *contract = basis->contract;
501e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5027a982d89SJeremy L. Thompson }
5037a982d89SJeremy L. Thompson 
5047a982d89SJeremy L. Thompson /**
5057a982d89SJeremy L. Thompson   @brief Set CeedTensorContract of a CeedBasis
5067a982d89SJeremy L. Thompson 
507ea61e9acSJeremy L Thompson   @param[in,out] basis    CeedBasis
508ea61e9acSJeremy L Thompson   @param[in]     contract CeedTensorContract to set
5097a982d89SJeremy L. Thompson 
5107a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
5117a982d89SJeremy L. Thompson 
5127a982d89SJeremy L. Thompson   @ref Backend
5137a982d89SJeremy L. Thompson **/
51434359f16Sjeremylt int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract contract) {
51534359f16Sjeremylt   basis->contract = contract;
5162b730f8bSJeremy L Thompson   CeedCall(CeedTensorContractReference(contract));
517e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5187a982d89SJeremy L. Thompson }
5197a982d89SJeremy L. Thompson 
5207a982d89SJeremy L. Thompson /**
5217a982d89SJeremy L. Thompson   @brief Return a reference implementation of matrix multiplication C = A B.
522ba59ac12SSebastian Grimberg 
523ba59ac12SSebastian Grimberg   Note: This is a reference implementation for CPU CeedScalar pointers that is not intended for high performance.
5247a982d89SJeremy L. Thompson 
525ea61e9acSJeremy L Thompson   @param[in]  ceed  Ceed context for error handling
526d1d35e2fSjeremylt   @param[in]  mat_A Row-major matrix A
527d1d35e2fSjeremylt   @param[in]  mat_B Row-major matrix B
528d1d35e2fSjeremylt   @param[out] mat_C Row-major output matrix C
529ea61e9acSJeremy L Thompson   @param[in]  m     Number of rows of C
530ea61e9acSJeremy L Thompson   @param[in]  n     Number of columns of C
531ea61e9acSJeremy L Thompson   @param[in]  kk    Number of columns of A/rows of B
5327a982d89SJeremy L. Thompson 
5337a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
5347a982d89SJeremy L. Thompson 
5357a982d89SJeremy L. Thompson   @ref Utility
5367a982d89SJeremy L. Thompson **/
5372b730f8bSJeremy L Thompson int CeedMatrixMatrixMultiply(Ceed ceed, const CeedScalar *mat_A, const CeedScalar *mat_B, CeedScalar *mat_C, CeedInt m, CeedInt n, CeedInt kk) {
5382b730f8bSJeremy L Thompson   for (CeedInt i = 0; i < m; i++) {
5397a982d89SJeremy L. Thompson     for (CeedInt j = 0; j < n; j++) {
5407a982d89SJeremy L. Thompson       CeedScalar sum = 0;
5412b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < kk; k++) sum += mat_A[k + i * kk] * mat_B[j + k * n];
542d1d35e2fSjeremylt       mat_C[j + i * n] = sum;
5437a982d89SJeremy L. Thompson     }
5442b730f8bSJeremy L Thompson   }
545e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5467a982d89SJeremy L. Thompson }
5477a982d89SJeremy L. Thompson 
548ba59ac12SSebastian Grimberg /**
549ba59ac12SSebastian Grimberg   @brief Return QR Factorization of a matrix
550ba59ac12SSebastian Grimberg 
551ba59ac12SSebastian Grimberg   @param[in]     ceed Ceed context for error handling
552ba59ac12SSebastian Grimberg   @param[in,out] mat  Row-major matrix to be factorized in place
553ba59ac12SSebastian Grimberg   @param[in,out] tau  Vector of length m of scaling factors
554ba59ac12SSebastian Grimberg   @param[in]     m    Number of rows
555ba59ac12SSebastian Grimberg   @param[in]     n    Number of columns
556ba59ac12SSebastian Grimberg 
557ba59ac12SSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
558ba59ac12SSebastian Grimberg 
559ba59ac12SSebastian Grimberg   @ref Utility
560ba59ac12SSebastian Grimberg **/
561ba59ac12SSebastian Grimberg int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau, CeedInt m, CeedInt n) {
562ba59ac12SSebastian Grimberg   CeedScalar v[m];
563ba59ac12SSebastian Grimberg 
564ba59ac12SSebastian Grimberg   // Check matrix shape
565ba59ac12SSebastian Grimberg   if (n > m) {
566ba59ac12SSebastian Grimberg     // LCOV_EXCL_START
567ba59ac12SSebastian Grimberg     return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute QR factorization with n > m");
568ba59ac12SSebastian Grimberg     // LCOV_EXCL_STOP
569ba59ac12SSebastian Grimberg   }
570ba59ac12SSebastian Grimberg 
571ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) {
572ba59ac12SSebastian Grimberg     if (i >= m - 1) {  // last row of matrix, no reflection needed
573ba59ac12SSebastian Grimberg       tau[i] = 0.;
574ba59ac12SSebastian Grimberg       break;
575ba59ac12SSebastian Grimberg     }
576ba59ac12SSebastian Grimberg     // Calculate Householder vector, magnitude
577ba59ac12SSebastian Grimberg     CeedScalar sigma = 0.0;
578ba59ac12SSebastian Grimberg     v[i]             = mat[i + n * i];
579ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < m; j++) {
580ba59ac12SSebastian Grimberg       v[j] = mat[i + n * j];
581ba59ac12SSebastian Grimberg       sigma += v[j] * v[j];
582ba59ac12SSebastian Grimberg     }
583ba59ac12SSebastian Grimberg     CeedScalar norm = sqrt(v[i] * v[i] + sigma);  // norm of v[i:m]
584ba59ac12SSebastian Grimberg     CeedScalar R_ii = -copysign(norm, v[i]);
585ba59ac12SSebastian Grimberg     v[i] -= R_ii;
586ba59ac12SSebastian Grimberg     // norm of v[i:m] after modification above and scaling below
587ba59ac12SSebastian Grimberg     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
588ba59ac12SSebastian Grimberg     //   tau = 2 / (norm*norm)
589ba59ac12SSebastian Grimberg     tau[i] = 2 * v[i] * v[i] / (v[i] * v[i] + sigma);
590ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < m; j++) v[j] /= v[i];
591ba59ac12SSebastian Grimberg 
592ba59ac12SSebastian Grimberg     // Apply Householder reflector to lower right panel
593ba59ac12SSebastian Grimberg     CeedHouseholderReflect(&mat[i * n + i + 1], &v[i], tau[i], m - i, n - i - 1, n, 1);
594ba59ac12SSebastian Grimberg     // Save v
595ba59ac12SSebastian Grimberg     mat[i + n * i] = R_ii;
596ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < m; j++) mat[i + n * j] = v[j];
597ba59ac12SSebastian Grimberg   }
598ba59ac12SSebastian Grimberg   return CEED_ERROR_SUCCESS;
599ba59ac12SSebastian Grimberg }
600ba59ac12SSebastian Grimberg 
601ba59ac12SSebastian Grimberg /**
602ba59ac12SSebastian Grimberg   @brief Apply Householder Q matrix
603ba59ac12SSebastian Grimberg 
604ba59ac12SSebastian Grimberg   Compute mat_A = mat_Q mat_A, where mat_Q is mxm and mat_A is mxn.
605ba59ac12SSebastian Grimberg 
606ba59ac12SSebastian Grimberg   @param[in,out] mat_A  Matrix to apply Householder Q to, in place
607ba59ac12SSebastian Grimberg   @param[in]     mat_Q  Householder Q matrix
608ba59ac12SSebastian Grimberg   @param[in]     tau    Householder scaling factors
609ba59ac12SSebastian Grimberg   @param[in]     t_mode Transpose mode for application
610ba59ac12SSebastian Grimberg   @param[in]     m      Number of rows in A
611ba59ac12SSebastian Grimberg   @param[in]     n      Number of columns in A
612ba59ac12SSebastian Grimberg   @param[in]     k      Number of elementary reflectors in Q, k<m
613ba59ac12SSebastian Grimberg   @param[in]     row    Row stride in A
614ba59ac12SSebastian Grimberg   @param[in]     col    Col stride in A
615ba59ac12SSebastian Grimberg 
616ba59ac12SSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
617ba59ac12SSebastian Grimberg 
618*c4e3f59bSSebastian Grimberg   @ref Utility
619ba59ac12SSebastian Grimberg **/
620ba59ac12SSebastian Grimberg int CeedHouseholderApplyQ(CeedScalar *mat_A, const CeedScalar *mat_Q, const CeedScalar *tau, CeedTransposeMode t_mode, CeedInt m, CeedInt n,
621ba59ac12SSebastian Grimberg                           CeedInt k, CeedInt row, CeedInt col) {
622ba59ac12SSebastian Grimberg   CeedScalar *v;
623ba59ac12SSebastian Grimberg   CeedCall(CeedMalloc(m, &v));
624ba59ac12SSebastian Grimberg   for (CeedInt ii = 0; ii < k; ii++) {
625ba59ac12SSebastian Grimberg     CeedInt i = t_mode == CEED_TRANSPOSE ? ii : k - 1 - ii;
626ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < m; j++) v[j] = mat_Q[j * k + i];
627ba59ac12SSebastian Grimberg     // Apply Householder reflector (I - tau v v^T) collo_grad_1d^T
628ba59ac12SSebastian Grimberg     CeedCall(CeedHouseholderReflect(&mat_A[i * row], &v[i], tau[i], m - i, n, row, col));
629ba59ac12SSebastian Grimberg   }
630ba59ac12SSebastian Grimberg   CeedCall(CeedFree(&v));
631ba59ac12SSebastian Grimberg   return CEED_ERROR_SUCCESS;
632ba59ac12SSebastian Grimberg }
633ba59ac12SSebastian Grimberg 
634ba59ac12SSebastian Grimberg /**
635ba59ac12SSebastian Grimberg   @brief Return symmetric Schur decomposition of the symmetric matrix mat via symmetric QR factorization
636ba59ac12SSebastian Grimberg 
637ba59ac12SSebastian Grimberg   @param[in]     ceed   Ceed context for error handling
638ba59ac12SSebastian Grimberg   @param[in,out] mat    Row-major matrix to be factorized in place
639ba59ac12SSebastian Grimberg   @param[out]    lambda Vector of length n of eigenvalues
640ba59ac12SSebastian Grimberg   @param[in]     n      Number of rows/columns
641ba59ac12SSebastian Grimberg 
642ba59ac12SSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
643ba59ac12SSebastian Grimberg 
644ba59ac12SSebastian Grimberg   @ref Utility
645ba59ac12SSebastian Grimberg **/
646ba59ac12SSebastian Grimberg CeedPragmaOptimizeOff int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat, CeedScalar *lambda, CeedInt n) {
647ba59ac12SSebastian Grimberg   // Check bounds for clang-tidy
648ba59ac12SSebastian Grimberg   if (n < 2) {
649ba59ac12SSebastian Grimberg     // LCOV_EXCL_START
650ba59ac12SSebastian Grimberg     return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute symmetric Schur decomposition of scalars");
651ba59ac12SSebastian Grimberg     // LCOV_EXCL_STOP
652ba59ac12SSebastian Grimberg   }
653ba59ac12SSebastian Grimberg 
654ba59ac12SSebastian Grimberg   CeedScalar v[n - 1], tau[n - 1], mat_T[n * n];
655ba59ac12SSebastian Grimberg 
656ba59ac12SSebastian Grimberg   // Copy mat to mat_T and set mat to I
657ba59ac12SSebastian Grimberg   memcpy(mat_T, mat, n * n * sizeof(mat[0]));
658ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) {
659ba59ac12SSebastian Grimberg     for (CeedInt j = 0; j < n; j++) mat[j + n * i] = (i == j) ? 1 : 0;
660ba59ac12SSebastian Grimberg   }
661ba59ac12SSebastian Grimberg 
662ba59ac12SSebastian Grimberg   // Reduce to tridiagonal
663ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n - 1; i++) {
664ba59ac12SSebastian Grimberg     // Calculate Householder vector, magnitude
665ba59ac12SSebastian Grimberg     CeedScalar sigma = 0.0;
666ba59ac12SSebastian Grimberg     v[i]             = mat_T[i + n * (i + 1)];
667ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < n - 1; j++) {
668ba59ac12SSebastian Grimberg       v[j] = mat_T[i + n * (j + 1)];
669ba59ac12SSebastian Grimberg       sigma += v[j] * v[j];
670ba59ac12SSebastian Grimberg     }
671ba59ac12SSebastian Grimberg     CeedScalar norm = sqrt(v[i] * v[i] + sigma);  // norm of v[i:n-1]
672ba59ac12SSebastian Grimberg     CeedScalar R_ii = -copysign(norm, v[i]);
673ba59ac12SSebastian Grimberg     v[i] -= R_ii;
674ba59ac12SSebastian Grimberg     // norm of v[i:m] after modification above and scaling below
675ba59ac12SSebastian Grimberg     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
676ba59ac12SSebastian Grimberg     //   tau = 2 / (norm*norm)
677ba59ac12SSebastian Grimberg     tau[i] = i == n - 2 ? 2 : 2 * v[i] * v[i] / (v[i] * v[i] + sigma);
678ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < n - 1; j++) v[j] /= v[i];
679ba59ac12SSebastian Grimberg 
680ba59ac12SSebastian Grimberg     // Update sub and super diagonal
681ba59ac12SSebastian Grimberg     for (CeedInt j = i + 2; j < n; j++) {
682ba59ac12SSebastian Grimberg       mat_T[i + n * j] = 0;
683ba59ac12SSebastian Grimberg       mat_T[j + n * i] = 0;
684ba59ac12SSebastian Grimberg     }
685ba59ac12SSebastian Grimberg     // Apply symmetric Householder reflector to lower right panel
686ba59ac12SSebastian Grimberg     CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1);
687ba59ac12SSebastian Grimberg     CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), 1, n);
688ba59ac12SSebastian Grimberg 
689ba59ac12SSebastian Grimberg     // Save v
690ba59ac12SSebastian Grimberg     mat_T[i + n * (i + 1)] = R_ii;
691ba59ac12SSebastian Grimberg     mat_T[(i + 1) + n * i] = R_ii;
692ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < n - 1; j++) {
693ba59ac12SSebastian Grimberg       mat_T[i + n * (j + 1)] = v[j];
694ba59ac12SSebastian Grimberg     }
695ba59ac12SSebastian Grimberg   }
696ba59ac12SSebastian Grimberg   // Backwards accumulation of Q
697ba59ac12SSebastian Grimberg   for (CeedInt i = n - 2; i >= 0; i--) {
698ba59ac12SSebastian Grimberg     if (tau[i] > 0.0) {
699ba59ac12SSebastian Grimberg       v[i] = 1;
700ba59ac12SSebastian Grimberg       for (CeedInt j = i + 1; j < n - 1; j++) {
701ba59ac12SSebastian Grimberg         v[j]                   = mat_T[i + n * (j + 1)];
702ba59ac12SSebastian Grimberg         mat_T[i + n * (j + 1)] = 0;
703ba59ac12SSebastian Grimberg       }
704ba59ac12SSebastian Grimberg       CeedHouseholderReflect(&mat[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1);
705ba59ac12SSebastian Grimberg     }
706ba59ac12SSebastian Grimberg   }
707ba59ac12SSebastian Grimberg 
708ba59ac12SSebastian Grimberg   // Reduce sub and super diagonal
709ba59ac12SSebastian Grimberg   CeedInt    p = 0, q = 0, itr = 0, max_itr = n * n * n * n;
710ba59ac12SSebastian Grimberg   CeedScalar tol = CEED_EPSILON;
711ba59ac12SSebastian Grimberg 
712ba59ac12SSebastian Grimberg   while (itr < max_itr) {
713ba59ac12SSebastian Grimberg     // Update p, q, size of reduced portions of diagonal
714ba59ac12SSebastian Grimberg     p = 0;
715ba59ac12SSebastian Grimberg     q = 0;
716ba59ac12SSebastian Grimberg     for (CeedInt i = n - 2; i >= 0; i--) {
717ba59ac12SSebastian Grimberg       if (fabs(mat_T[i + n * (i + 1)]) < tol) q += 1;
718ba59ac12SSebastian Grimberg       else break;
719ba59ac12SSebastian Grimberg     }
720ba59ac12SSebastian Grimberg     for (CeedInt i = 0; i < n - q - 1; i++) {
721ba59ac12SSebastian Grimberg       if (fabs(mat_T[i + n * (i + 1)]) < tol) p += 1;
722ba59ac12SSebastian Grimberg       else break;
723ba59ac12SSebastian Grimberg     }
724ba59ac12SSebastian Grimberg     if (q == n - 1) break;  // Finished reducing
725ba59ac12SSebastian Grimberg 
726ba59ac12SSebastian Grimberg     // Reduce tridiagonal portion
727ba59ac12SSebastian Grimberg     CeedScalar t_nn = mat_T[(n - 1 - q) + n * (n - 1 - q)], t_nnm1 = mat_T[(n - 2 - q) + n * (n - 1 - q)];
728ba59ac12SSebastian Grimberg     CeedScalar d  = (mat_T[(n - 2 - q) + n * (n - 2 - q)] - t_nn) / 2;
729ba59ac12SSebastian Grimberg     CeedScalar mu = t_nn - t_nnm1 * t_nnm1 / (d + copysign(sqrt(d * d + t_nnm1 * t_nnm1), d));
730ba59ac12SSebastian Grimberg     CeedScalar x  = mat_T[p + n * p] - mu;
731ba59ac12SSebastian Grimberg     CeedScalar z  = mat_T[p + n * (p + 1)];
732ba59ac12SSebastian Grimberg     for (CeedInt k = p; k < n - q - 1; k++) {
733ba59ac12SSebastian Grimberg       // Compute Givens rotation
734ba59ac12SSebastian Grimberg       CeedScalar c = 1, s = 0;
735ba59ac12SSebastian Grimberg       if (fabs(z) > tol) {
736ba59ac12SSebastian Grimberg         if (fabs(z) > fabs(x)) {
737ba59ac12SSebastian Grimberg           CeedScalar tau = -x / z;
738ba59ac12SSebastian Grimberg           s = 1 / sqrt(1 + tau * tau), c = s * tau;
739ba59ac12SSebastian Grimberg         } else {
740ba59ac12SSebastian Grimberg           CeedScalar tau = -z / x;
741ba59ac12SSebastian Grimberg           c = 1 / sqrt(1 + tau * tau), s = c * tau;
742ba59ac12SSebastian Grimberg         }
743ba59ac12SSebastian Grimberg       }
744ba59ac12SSebastian Grimberg 
745ba59ac12SSebastian Grimberg       // Apply Givens rotation to T
746ba59ac12SSebastian Grimberg       CeedGivensRotation(mat_T, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n);
747ba59ac12SSebastian Grimberg       CeedGivensRotation(mat_T, c, s, CEED_TRANSPOSE, k, k + 1, n, n);
748ba59ac12SSebastian Grimberg 
749ba59ac12SSebastian Grimberg       // Apply Givens rotation to Q
750ba59ac12SSebastian Grimberg       CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n);
751ba59ac12SSebastian Grimberg 
752ba59ac12SSebastian Grimberg       // Update x, z
753ba59ac12SSebastian Grimberg       if (k < n - q - 2) {
754ba59ac12SSebastian Grimberg         x = mat_T[k + n * (k + 1)];
755ba59ac12SSebastian Grimberg         z = mat_T[k + n * (k + 2)];
756ba59ac12SSebastian Grimberg       }
757ba59ac12SSebastian Grimberg     }
758ba59ac12SSebastian Grimberg     itr++;
759ba59ac12SSebastian Grimberg   }
760ba59ac12SSebastian Grimberg 
761ba59ac12SSebastian Grimberg   // Save eigenvalues
762ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) lambda[i] = mat_T[i + n * i];
763ba59ac12SSebastian Grimberg 
764ba59ac12SSebastian Grimberg   // Check convergence
765ba59ac12SSebastian Grimberg   if (itr == max_itr && q < n - 1) {
766ba59ac12SSebastian Grimberg     // LCOV_EXCL_START
767ba59ac12SSebastian Grimberg     return CeedError(ceed, CEED_ERROR_MINOR, "Symmetric QR failed to converge");
768ba59ac12SSebastian Grimberg     // LCOV_EXCL_STOP
769ba59ac12SSebastian Grimberg   }
770ba59ac12SSebastian Grimberg   return CEED_ERROR_SUCCESS;
771ba59ac12SSebastian Grimberg }
772ba59ac12SSebastian Grimberg CeedPragmaOptimizeOn;
773ba59ac12SSebastian Grimberg 
774ba59ac12SSebastian Grimberg /**
775ba59ac12SSebastian Grimberg   @brief Return Simultaneous Diagonalization of two matrices.
776ba59ac12SSebastian Grimberg 
777ba59ac12SSebastian Grimberg   This solves the generalized eigenvalue problem A x = lambda B x, where A and B are symmetric and B is positive definite.
778ba59ac12SSebastian Grimberg   We generate the matrix X and vector Lambda such that X^T A X = Lambda and X^T B X = I.
779ba59ac12SSebastian Grimberg   This is equivalent to the LAPACK routine 'sygv' with TYPE = 1.
780ba59ac12SSebastian Grimberg 
781ba59ac12SSebastian Grimberg   @param[in]  ceed   Ceed context for error handling
782ba59ac12SSebastian Grimberg   @param[in]  mat_A  Row-major matrix to be factorized with eigenvalues
783ba59ac12SSebastian Grimberg   @param[in]  mat_B  Row-major matrix to be factorized to identity
784ba59ac12SSebastian Grimberg   @param[out] mat_X  Row-major orthogonal matrix
785ba59ac12SSebastian Grimberg   @param[out] lambda Vector of length n of generalized eigenvalues
786ba59ac12SSebastian Grimberg   @param[in]  n      Number of rows/columns
787ba59ac12SSebastian Grimberg 
788ba59ac12SSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
789ba59ac12SSebastian Grimberg 
790ba59ac12SSebastian Grimberg   @ref Utility
791ba59ac12SSebastian Grimberg **/
792ba59ac12SSebastian Grimberg CeedPragmaOptimizeOff int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *mat_A, CeedScalar *mat_B, CeedScalar *mat_X, CeedScalar *lambda,
793ba59ac12SSebastian Grimberg                                                           CeedInt n) {
794ba59ac12SSebastian Grimberg   CeedScalar *mat_C, *mat_G, *vec_D;
795ba59ac12SSebastian Grimberg   CeedCall(CeedCalloc(n * n, &mat_C));
796ba59ac12SSebastian Grimberg   CeedCall(CeedCalloc(n * n, &mat_G));
797ba59ac12SSebastian Grimberg   CeedCall(CeedCalloc(n, &vec_D));
798ba59ac12SSebastian Grimberg 
799ba59ac12SSebastian Grimberg   // Compute B = G D G^T
800ba59ac12SSebastian Grimberg   memcpy(mat_G, mat_B, n * n * sizeof(mat_B[0]));
801ba59ac12SSebastian Grimberg   CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_G, vec_D, n));
802ba59ac12SSebastian Grimberg 
803ba59ac12SSebastian Grimberg   // Sort eigenvalues
804ba59ac12SSebastian Grimberg   for (CeedInt i = n - 1; i >= 0; i--) {
805ba59ac12SSebastian Grimberg     for (CeedInt j = 0; j < i; j++) {
806ba59ac12SSebastian Grimberg       if (fabs(vec_D[j]) > fabs(vec_D[j + 1])) {
807ba59ac12SSebastian Grimberg         CeedScalar temp;
808ba59ac12SSebastian Grimberg         temp         = vec_D[j];
809ba59ac12SSebastian Grimberg         vec_D[j]     = vec_D[j + 1];
810ba59ac12SSebastian Grimberg         vec_D[j + 1] = temp;
811ba59ac12SSebastian Grimberg         for (CeedInt k = 0; k < n; k++) {
812ba59ac12SSebastian Grimberg           temp                 = mat_G[k * n + j];
813ba59ac12SSebastian Grimberg           mat_G[k * n + j]     = mat_G[k * n + j + 1];
814ba59ac12SSebastian Grimberg           mat_G[k * n + j + 1] = temp;
815ba59ac12SSebastian Grimberg         }
816ba59ac12SSebastian Grimberg       }
817ba59ac12SSebastian Grimberg     }
818ba59ac12SSebastian Grimberg   }
819ba59ac12SSebastian Grimberg 
820ba59ac12SSebastian Grimberg   // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T
821ba59ac12SSebastian Grimberg   //           = D^-1/2 G^T A G D^-1/2
822ba59ac12SSebastian Grimberg   // -- D = D^-1/2
823ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) vec_D[i] = 1. / sqrt(vec_D[i]);
824ba59ac12SSebastian Grimberg   // -- G = G D^-1/2
825ba59ac12SSebastian Grimberg   // -- C = D^-1/2 G^T
826ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) {
827ba59ac12SSebastian Grimberg     for (CeedInt j = 0; j < n; j++) {
828ba59ac12SSebastian Grimberg       mat_G[i * n + j] *= vec_D[j];
829ba59ac12SSebastian Grimberg       mat_C[j * n + i] = mat_G[i * n + j];
830ba59ac12SSebastian Grimberg     }
831ba59ac12SSebastian Grimberg   }
832ba59ac12SSebastian Grimberg   // -- X = (D^-1/2 G^T) A
833ba59ac12SSebastian Grimberg   CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_C, (const CeedScalar *)mat_A, mat_X, n, n, n));
834ba59ac12SSebastian Grimberg   // -- C = (D^-1/2 G^T A) (G D^-1/2)
835ba59ac12SSebastian Grimberg   CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_X, (const CeedScalar *)mat_G, mat_C, n, n, n));
836ba59ac12SSebastian Grimberg 
837ba59ac12SSebastian Grimberg   // Compute Q^T C Q = lambda
838ba59ac12SSebastian Grimberg   CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_C, lambda, n));
839ba59ac12SSebastian Grimberg 
840ba59ac12SSebastian Grimberg   // Sort eigenvalues
841ba59ac12SSebastian Grimberg   for (CeedInt i = n - 1; i >= 0; i--) {
842ba59ac12SSebastian Grimberg     for (CeedInt j = 0; j < i; j++) {
843ba59ac12SSebastian Grimberg       if (fabs(lambda[j]) > fabs(lambda[j + 1])) {
844ba59ac12SSebastian Grimberg         CeedScalar temp;
845ba59ac12SSebastian Grimberg         temp          = lambda[j];
846ba59ac12SSebastian Grimberg         lambda[j]     = lambda[j + 1];
847ba59ac12SSebastian Grimberg         lambda[j + 1] = temp;
848ba59ac12SSebastian Grimberg         for (CeedInt k = 0; k < n; k++) {
849ba59ac12SSebastian Grimberg           temp                 = mat_C[k * n + j];
850ba59ac12SSebastian Grimberg           mat_C[k * n + j]     = mat_C[k * n + j + 1];
851ba59ac12SSebastian Grimberg           mat_C[k * n + j + 1] = temp;
852ba59ac12SSebastian Grimberg         }
853ba59ac12SSebastian Grimberg       }
854ba59ac12SSebastian Grimberg     }
855ba59ac12SSebastian Grimberg   }
856ba59ac12SSebastian Grimberg 
857ba59ac12SSebastian Grimberg   // Set X = (G D^1/2)^-T Q
858ba59ac12SSebastian Grimberg   //       = G D^-1/2 Q
859ba59ac12SSebastian Grimberg   CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_G, (const CeedScalar *)mat_C, mat_X, n, n, n));
860ba59ac12SSebastian Grimberg 
861ba59ac12SSebastian Grimberg   // Cleanup
862ba59ac12SSebastian Grimberg   CeedCall(CeedFree(&mat_C));
863ba59ac12SSebastian Grimberg   CeedCall(CeedFree(&mat_G));
864ba59ac12SSebastian Grimberg   CeedCall(CeedFree(&vec_D));
865ba59ac12SSebastian Grimberg   return CEED_ERROR_SUCCESS;
866ba59ac12SSebastian Grimberg }
867ba59ac12SSebastian Grimberg CeedPragmaOptimizeOn;
868ba59ac12SSebastian Grimberg 
8697a982d89SJeremy L. Thompson /// @}
8707a982d89SJeremy L. Thompson 
8717a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
8727a982d89SJeremy L. Thompson /// CeedBasis Public API
8737a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
8747a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser
875d7b241e6Sjeremylt /// @{
876d7b241e6Sjeremylt 
877b11c1e72Sjeremylt /**
878ba59ac12SSebastian Grimberg   @brief Create a tensor-product basis for H^1 discretizations
879b11c1e72Sjeremylt 
880ea61e9acSJeremy L Thompson   @param[in]  ceed        Ceed object where the CeedBasis will be created
881ea61e9acSJeremy L Thompson   @param[in]  dim         Topological dimension
882ea61e9acSJeremy L Thompson   @param[in]  num_comp    Number of field components (1 for scalar fields)
883ea61e9acSJeremy L Thompson   @param[in]  P_1d        Number of nodes in one dimension
884ea61e9acSJeremy L Thompson   @param[in]  Q_1d        Number of quadrature points in one dimension
885ea61e9acSJeremy L Thompson   @param[in]  interp_1d   Row-major (Q_1d * P_1d) matrix expressing the values of nodal basis functions at quadrature points
886ea61e9acSJeremy L Thompson   @param[in]  grad_1d     Row-major (Q_1d * P_1d) matrix expressing derivatives of nodal basis functions at quadrature points
887ea61e9acSJeremy 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]
888ea61e9acSJeremy L Thompson   @param[in]  q_weight_1d Array of length Q_1d holding the quadrature weights on the reference element
889ea61e9acSJeremy L Thompson   @param[out] basis       Address of the variable where the newly created CeedBasis will be stored.
890b11c1e72Sjeremylt 
891b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
892dfdf5a53Sjeremylt 
8937a982d89SJeremy L. Thompson   @ref User
894b11c1e72Sjeremylt **/
8952b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d,
8962b730f8bSJeremy L Thompson                             const CeedScalar *grad_1d, const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis *basis) {
8975fe0d4faSjeremylt   if (!ceed->BasisCreateTensorH1) {
8985fe0d4faSjeremylt     Ceed delegate;
8992b730f8bSJeremy L Thompson     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
9005fe0d4faSjeremylt 
9012b730f8bSJeremy L Thompson     if (!delegate) {
902c042f62fSJeremy L Thompson       // LCOV_EXCL_START
9032b730f8bSJeremy L Thompson       return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisCreateTensorH1");
904c042f62fSJeremy L Thompson       // LCOV_EXCL_STOP
9052b730f8bSJeremy L Thompson     }
9065fe0d4faSjeremylt 
9072b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateTensorH1(delegate, dim, num_comp, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis));
908e15f9bd0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
9095fe0d4faSjeremylt   }
910e15f9bd0SJeremy L Thompson 
9112b730f8bSJeremy L Thompson   if (dim < 1) {
912e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
9132b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis dimension must be a positive value");
914e15f9bd0SJeremy L Thompson     // LCOV_EXCL_STOP
9152b730f8bSJeremy L Thompson   }
916227444bfSJeremy L Thompson 
9172b730f8bSJeremy L Thompson   if (num_comp < 1) {
918227444bfSJeremy L Thompson     // LCOV_EXCL_START
9192b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component");
920227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
9212b730f8bSJeremy L Thompson   }
922227444bfSJeremy L Thompson 
9232b730f8bSJeremy L Thompson   if (P_1d < 1) {
924227444bfSJeremy L Thompson     // LCOV_EXCL_START
9252b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node");
926227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
9272b730f8bSJeremy L Thompson   }
928227444bfSJeremy L Thompson 
9292b730f8bSJeremy L Thompson   if (Q_1d < 1) {
930227444bfSJeremy L Thompson     // LCOV_EXCL_START
9312b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point");
932227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
9332b730f8bSJeremy L Thompson   }
934227444bfSJeremy L Thompson 
9352b730f8bSJeremy L Thompson   CeedElemTopology topo = dim == 1 ? CEED_TOPOLOGY_LINE : dim == 2 ? CEED_TOPOLOGY_QUAD : CEED_TOPOLOGY_HEX;
936e15f9bd0SJeremy L Thompson 
9372b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(1, basis));
938d7b241e6Sjeremylt   (*basis)->ceed = ceed;
9392b730f8bSJeremy L Thompson   CeedCall(CeedReference(ceed));
940d1d35e2fSjeremylt   (*basis)->ref_count    = 1;
941d1d35e2fSjeremylt   (*basis)->tensor_basis = 1;
942d7b241e6Sjeremylt   (*basis)->dim          = dim;
943d99fa3c5SJeremy L Thompson   (*basis)->topo         = topo;
944d1d35e2fSjeremylt   (*basis)->num_comp     = num_comp;
945d1d35e2fSjeremylt   (*basis)->P_1d         = P_1d;
946d1d35e2fSjeremylt   (*basis)->Q_1d         = Q_1d;
947d1d35e2fSjeremylt   (*basis)->P            = CeedIntPow(P_1d, dim);
948d1d35e2fSjeremylt   (*basis)->Q            = CeedIntPow(Q_1d, dim);
949*c4e3f59bSSebastian Grimberg   (*basis)->fe_space     = CEED_FE_SPACE_H1;
9502b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d, &(*basis)->q_ref_1d));
9512b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d, &(*basis)->q_weight_1d));
952ff3a0f91SJeremy L Thompson   if (q_ref_1d) memcpy((*basis)->q_ref_1d, q_ref_1d, Q_1d * sizeof(q_ref_1d[0]));
9532b730f8bSJeremy L Thompson   if (q_weight_1d) memcpy((*basis)->q_weight_1d, q_weight_1d, Q_1d * sizeof(q_weight_1d[0]));
9542b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->interp_1d));
9552b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->grad_1d));
9562b730f8bSJeremy L Thompson   if (interp_1d) memcpy((*basis)->interp_1d, interp_1d, Q_1d * P_1d * sizeof(interp_1d[0]));
957ff3a0f91SJeremy L Thompson   if (grad_1d) memcpy((*basis)->grad_1d, grad_1d, Q_1d * P_1d * sizeof(grad_1d[0]));
9582b730f8bSJeremy L Thompson   CeedCall(ceed->BasisCreateTensorH1(dim, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, *basis));
959e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
960d7b241e6Sjeremylt }
961d7b241e6Sjeremylt 
962b11c1e72Sjeremylt /**
96395bb1877Svaleriabarra   @brief Create a tensor-product Lagrange basis
964b11c1e72Sjeremylt 
965ea61e9acSJeremy L Thompson   @param[in]  ceed      Ceed object where the CeedBasis will be created
966ea61e9acSJeremy L Thompson   @param[in]  dim       Topological dimension of element
967ea61e9acSJeremy L Thompson   @param[in]  num_comp  Number of field components (1 for scalar fields)
968ea61e9acSJeremy L Thompson   @param[in]  P         Number of Gauss-Lobatto nodes in one dimension.
969ea61e9acSJeremy L Thompson                           The polynomial degree of the resulting Q_k element is k=P-1.
970ea61e9acSJeremy L Thompson   @param[in]  Q         Number of quadrature points in one dimension.
971ea61e9acSJeremy L Thompson   @param[in]  quad_mode Distribution of the Q quadrature points (affects order of accuracy for the quadrature)
972ea61e9acSJeremy L Thompson   @param[out] basis     Address of the variable where the newly created CeedBasis will be stored.
973b11c1e72Sjeremylt 
974b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
975dfdf5a53Sjeremylt 
9767a982d89SJeremy L. Thompson   @ref User
977b11c1e72Sjeremylt **/
9782b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P, CeedInt Q, CeedQuadMode quad_mode, CeedBasis *basis) {
979d7b241e6Sjeremylt   // Allocate
9802b730f8bSJeremy L Thompson   int        ierr = CEED_ERROR_SUCCESS, i, j, k;
9812b730f8bSJeremy L Thompson   CeedScalar c1, c2, c3, c4, dx, *nodes, *interp_1d, *grad_1d, *q_ref_1d, *q_weight_1d;
9824d537eeaSYohann 
9832b730f8bSJeremy L Thompson   if (dim < 1) {
984c042f62fSJeremy L Thompson     // LCOV_EXCL_START
9852b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis dimension must be a positive value");
986c042f62fSJeremy L Thompson     // LCOV_EXCL_STOP
9872b730f8bSJeremy L Thompson   }
9884d537eeaSYohann 
9892b730f8bSJeremy L Thompson   if (num_comp < 1) {
990227444bfSJeremy L Thompson     // LCOV_EXCL_START
9912b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component");
992227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
9932b730f8bSJeremy L Thompson   }
994227444bfSJeremy L Thompson 
9952b730f8bSJeremy L Thompson   if (P < 1) {
996227444bfSJeremy L Thompson     // LCOV_EXCL_START
9972b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node");
998227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
9992b730f8bSJeremy L Thompson   }
1000227444bfSJeremy L Thompson 
10012b730f8bSJeremy L Thompson   if (Q < 1) {
1002227444bfSJeremy L Thompson     // LCOV_EXCL_START
10032b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point");
1004227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
10052b730f8bSJeremy L Thompson   }
1006227444bfSJeremy L Thompson 
1007e15f9bd0SJeremy L Thompson   // Get Nodes and Weights
10082b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(P * Q, &interp_1d));
10092b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(P * Q, &grad_1d));
10102b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(P, &nodes));
10112b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q, &q_ref_1d));
10122b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q, &q_weight_1d));
10132b730f8bSJeremy L Thompson   if (CeedLobattoQuadrature(P, nodes, NULL) != CEED_ERROR_SUCCESS) goto cleanup;
1014d1d35e2fSjeremylt   switch (quad_mode) {
1015d7b241e6Sjeremylt     case CEED_GAUSS:
1016d1d35e2fSjeremylt       ierr = CeedGaussQuadrature(Q, q_ref_1d, q_weight_1d);
1017d7b241e6Sjeremylt       break;
1018d7b241e6Sjeremylt     case CEED_GAUSS_LOBATTO:
1019d1d35e2fSjeremylt       ierr = CeedLobattoQuadrature(Q, q_ref_1d, q_weight_1d);
1020d7b241e6Sjeremylt       break;
1021d7b241e6Sjeremylt   }
10222b730f8bSJeremy L Thompson   if (ierr != CEED_ERROR_SUCCESS) goto cleanup;
1023e15f9bd0SJeremy L Thompson 
1024d7b241e6Sjeremylt   // Build B, D matrix
1025d7b241e6Sjeremylt   // Fornberg, 1998
1026d7b241e6Sjeremylt   for (i = 0; i < Q; i++) {
1027d7b241e6Sjeremylt     c1                   = 1.0;
1028d1d35e2fSjeremylt     c3                   = nodes[0] - q_ref_1d[i];
1029d1d35e2fSjeremylt     interp_1d[i * P + 0] = 1.0;
1030d7b241e6Sjeremylt     for (j = 1; j < P; j++) {
1031d7b241e6Sjeremylt       c2 = 1.0;
1032d7b241e6Sjeremylt       c4 = c3;
1033d1d35e2fSjeremylt       c3 = nodes[j] - q_ref_1d[i];
1034d7b241e6Sjeremylt       for (k = 0; k < j; k++) {
1035d7b241e6Sjeremylt         dx = nodes[j] - nodes[k];
1036d7b241e6Sjeremylt         c2 *= dx;
1037d7b241e6Sjeremylt         if (k == j - 1) {
1038d1d35e2fSjeremylt           grad_1d[i * P + j]   = c1 * (interp_1d[i * P + k] - c4 * grad_1d[i * P + k]) / c2;
1039d1d35e2fSjeremylt           interp_1d[i * P + j] = -c1 * c4 * interp_1d[i * P + k] / c2;
1040d7b241e6Sjeremylt         }
1041d1d35e2fSjeremylt         grad_1d[i * P + k]   = (c3 * grad_1d[i * P + k] - interp_1d[i * P + k]) / dx;
1042d1d35e2fSjeremylt         interp_1d[i * P + k] = c3 * interp_1d[i * P + k] / dx;
1043d7b241e6Sjeremylt       }
1044d7b241e6Sjeremylt       c1 = c2;
1045d7b241e6Sjeremylt     }
1046d7b241e6Sjeremylt   }
10479ac7b42eSJeremy L Thompson   // Pass to CeedBasisCreateTensorH1
10482b730f8bSJeremy L Thompson   CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P, Q, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis));
1049e15f9bd0SJeremy L Thompson cleanup:
10502b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_1d));
10512b730f8bSJeremy L Thompson   CeedCall(CeedFree(&grad_1d));
10522b730f8bSJeremy L Thompson   CeedCall(CeedFree(&nodes));
10532b730f8bSJeremy L Thompson   CeedCall(CeedFree(&q_ref_1d));
10542b730f8bSJeremy L Thompson   CeedCall(CeedFree(&q_weight_1d));
1055e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1056d7b241e6Sjeremylt }
1057d7b241e6Sjeremylt 
1058b11c1e72Sjeremylt /**
1059ba59ac12SSebastian Grimberg   @brief Create a non tensor-product basis for H^1 discretizations
1060a8de75f0Sjeremylt 
1061ea61e9acSJeremy L Thompson   @param[in]  ceed      Ceed object where the CeedBasis will be created
1062ea61e9acSJeremy L Thompson   @param[in]  topo      Topology of element, e.g. hypercube, simplex, ect
1063ea61e9acSJeremy L Thompson   @param[in]  num_comp  Number of field components (1 for scalar fields)
1064ea61e9acSJeremy L Thompson   @param[in]  num_nodes Total number of nodes
1065ea61e9acSJeremy L Thompson   @param[in]  num_qpts  Total number of quadrature points
1066ea61e9acSJeremy L Thompson   @param[in]  interp    Row-major (num_qpts * num_nodes) matrix expressing the values of nodal basis functions at quadrature points
1067*c4e3f59bSSebastian Grimberg   @param[in]  grad      Row-major (dim * num_qpts * num_nodes) matrix expressing derivatives of nodal basis functions at quadrature points
10689fe083eeSJeremy L Thompson   @param[in]  q_ref     Array of length num_qpts * dim holding the locations of quadrature points on the reference element
1069ea61e9acSJeremy L Thompson   @param[in]  q_weight  Array of length num_qpts holding the quadrature weights on the reference element
1070ea61e9acSJeremy L Thompson   @param[out] basis     Address of the variable where the newly created CeedBasis will be stored.
1071a8de75f0Sjeremylt 
1072a8de75f0Sjeremylt   @return An error code: 0 - success, otherwise - failure
1073a8de75f0Sjeremylt 
10747a982d89SJeremy L. Thompson   @ref User
1075a8de75f0Sjeremylt **/
10762b730f8bSJeremy L Thompson int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
10772b730f8bSJeremy L Thompson                       const CeedScalar *grad, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) {
1078d1d35e2fSjeremylt   CeedInt P = num_nodes, Q = num_qpts, dim = 0;
1079a8de75f0Sjeremylt 
10805fe0d4faSjeremylt   if (!ceed->BasisCreateH1) {
10815fe0d4faSjeremylt     Ceed delegate;
10822b730f8bSJeremy L Thompson     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
10835fe0d4faSjeremylt 
10842b730f8bSJeremy L Thompson     if (!delegate) {
1085c042f62fSJeremy L Thompson       // LCOV_EXCL_START
10862b730f8bSJeremy L Thompson       return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisCreateH1");
1087c042f62fSJeremy L Thompson       // LCOV_EXCL_STOP
10882b730f8bSJeremy L Thompson     }
10895fe0d4faSjeremylt 
10902b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateH1(delegate, topo, num_comp, num_nodes, num_qpts, interp, grad, q_ref, q_weight, basis));
1091e15f9bd0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
10925fe0d4faSjeremylt   }
10935fe0d4faSjeremylt 
10942b730f8bSJeremy L Thompson   if (num_comp < 1) {
1095227444bfSJeremy L Thompson     // LCOV_EXCL_START
10962b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component");
1097227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
10982b730f8bSJeremy L Thompson   }
1099227444bfSJeremy L Thompson 
11002b730f8bSJeremy L Thompson   if (num_nodes < 1) {
1101227444bfSJeremy L Thompson     // LCOV_EXCL_START
11022b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node");
1103227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
11042b730f8bSJeremy L Thompson   }
1105227444bfSJeremy L Thompson 
11062b730f8bSJeremy L Thompson   if (num_qpts < 1) {
1107227444bfSJeremy L Thompson     // LCOV_EXCL_START
11082b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point");
1109227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
11102b730f8bSJeremy L Thompson   }
1111227444bfSJeremy L Thompson 
11122b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(1, basis));
1113a8de75f0Sjeremylt 
11142b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetTopologyDimension(topo, &dim));
1115a8de75f0Sjeremylt 
1116a8de75f0Sjeremylt   (*basis)->ceed = ceed;
11172b730f8bSJeremy L Thompson   CeedCall(CeedReference(ceed));
1118d1d35e2fSjeremylt   (*basis)->ref_count    = 1;
1119d1d35e2fSjeremylt   (*basis)->tensor_basis = 0;
1120a8de75f0Sjeremylt   (*basis)->dim          = dim;
1121d99fa3c5SJeremy L Thompson   (*basis)->topo         = topo;
1122d1d35e2fSjeremylt   (*basis)->num_comp     = num_comp;
1123a8de75f0Sjeremylt   (*basis)->P            = P;
1124a8de75f0Sjeremylt   (*basis)->Q            = Q;
1125*c4e3f59bSSebastian Grimberg   (*basis)->fe_space     = CEED_FE_SPACE_H1;
11262b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q * dim, &(*basis)->q_ref_1d));
11272b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q, &(*basis)->q_weight_1d));
1128ff3a0f91SJeremy L Thompson   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0]));
1129ff3a0f91SJeremy L Thompson   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0]));
11302b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q * P, &(*basis)->interp));
11312b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(dim * Q * P, &(*basis)->grad));
1132ff3a0f91SJeremy L Thompson   if (interp) memcpy((*basis)->interp, interp, Q * P * sizeof(interp[0]));
1133ff3a0f91SJeremy L Thompson   if (grad) memcpy((*basis)->grad, grad, dim * Q * P * sizeof(grad[0]));
11342b730f8bSJeremy L Thompson   CeedCall(ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, q_ref, q_weight, *basis));
1135e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1136a8de75f0Sjeremylt }
1137a8de75f0Sjeremylt 
1138a8de75f0Sjeremylt /**
1139859c15bbSJames Wright   @brief Create a non tensor-product basis for \f$H(\mathrm{div})\f$ discretizations
114050c301a5SRezgar Shakeri 
1141ea61e9acSJeremy L Thompson   @param[in]  ceed      Ceed object where the CeedBasis will be created
1142ea61e9acSJeremy 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
1143ea61e9acSJeremy L Thompson   @param[in]  num_comp  Number of components (usually 1 for vectors in H(div) bases)
1144ea61e9acSJeremy L Thompson   @param[in]  num_nodes Total number of nodes (dofs per element)
1145ea61e9acSJeremy L Thompson   @param[in]  num_qpts  Total number of quadrature points
1146*c4e3f59bSSebastian Grimberg   @param[in]  interp    Row-major (dim * num_qpts * num_nodes) matrix expressing the values of basis functions at quadrature points
1147*c4e3f59bSSebastian Grimberg   @param[in]  div       Row-major (num_qpts * num_nodes) matrix expressing divergence of basis functions at quadrature points
11489fe083eeSJeremy L Thompson   @param[in]  q_ref     Array of length num_qpts * dim holding the locations of quadrature points on the reference element
1149ea61e9acSJeremy L Thompson   @param[in]  q_weight  Array of length num_qpts holding the quadrature weights on the reference element
1150ea61e9acSJeremy L Thompson   @param[out] basis     Address of the variable where the newly created CeedBasis will be stored.
115150c301a5SRezgar Shakeri 
115250c301a5SRezgar Shakeri   @return An error code: 0 - success, otherwise - failure
115350c301a5SRezgar Shakeri 
115450c301a5SRezgar Shakeri   @ref User
115550c301a5SRezgar Shakeri **/
11562b730f8bSJeremy L Thompson int CeedBasisCreateHdiv(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
11572b730f8bSJeremy L Thompson                         const CeedScalar *div, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) {
115850c301a5SRezgar Shakeri   CeedInt Q = num_qpts, P = num_nodes, dim = 0;
1159*c4e3f59bSSebastian Grimberg 
116050c301a5SRezgar Shakeri   if (!ceed->BasisCreateHdiv) {
116150c301a5SRezgar Shakeri     Ceed delegate;
11622b730f8bSJeremy L Thompson     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
116350c301a5SRezgar Shakeri 
11642b730f8bSJeremy L Thompson     if (!delegate) {
116550c301a5SRezgar Shakeri       // LCOV_EXCL_START
11662b730f8bSJeremy L Thompson       return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateHdiv");
116750c301a5SRezgar Shakeri       // LCOV_EXCL_STOP
11682b730f8bSJeremy L Thompson     }
116950c301a5SRezgar Shakeri 
11702b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateHdiv(delegate, topo, num_comp, num_nodes, num_qpts, interp, div, q_ref, q_weight, basis));
117150c301a5SRezgar Shakeri     return CEED_ERROR_SUCCESS;
117250c301a5SRezgar Shakeri   }
117350c301a5SRezgar Shakeri 
11742b730f8bSJeremy L Thompson   if (num_comp < 1) {
1175227444bfSJeremy L Thompson     // LCOV_EXCL_START
11762b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component");
1177227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
11782b730f8bSJeremy L Thompson   }
1179227444bfSJeremy L Thompson 
11802b730f8bSJeremy L Thompson   if (num_nodes < 1) {
1181227444bfSJeremy L Thompson     // LCOV_EXCL_START
11822b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node");
1183227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
11842b730f8bSJeremy L Thompson   }
1185227444bfSJeremy L Thompson 
11862b730f8bSJeremy L Thompson   if (num_qpts < 1) {
1187227444bfSJeremy L Thompson     // LCOV_EXCL_START
11882b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point");
1189227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
11902b730f8bSJeremy L Thompson   }
1191227444bfSJeremy L Thompson 
11922b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(1, basis));
119350c301a5SRezgar Shakeri 
1194*c4e3f59bSSebastian Grimberg   CeedCall(CeedBasisGetTopologyDimension(topo, &dim));
1195*c4e3f59bSSebastian Grimberg 
119650c301a5SRezgar Shakeri   (*basis)->ceed = ceed;
11972b730f8bSJeremy L Thompson   CeedCall(CeedReference(ceed));
119850c301a5SRezgar Shakeri   (*basis)->ref_count    = 1;
119950c301a5SRezgar Shakeri   (*basis)->tensor_basis = 0;
120050c301a5SRezgar Shakeri   (*basis)->dim          = dim;
120150c301a5SRezgar Shakeri   (*basis)->topo         = topo;
120250c301a5SRezgar Shakeri   (*basis)->num_comp     = num_comp;
120350c301a5SRezgar Shakeri   (*basis)->P            = P;
120450c301a5SRezgar Shakeri   (*basis)->Q            = Q;
1205*c4e3f59bSSebastian Grimberg   (*basis)->fe_space     = CEED_FE_SPACE_HDIV;
12062b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q * dim, &(*basis)->q_ref_1d));
12072b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q, &(*basis)->q_weight_1d));
120850c301a5SRezgar Shakeri   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0]));
120950c301a5SRezgar Shakeri   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0]));
12102b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(dim * Q * P, &(*basis)->interp));
12112b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q * P, &(*basis)->div));
121250c301a5SRezgar Shakeri   if (interp) memcpy((*basis)->interp, interp, dim * Q * P * sizeof(interp[0]));
121350c301a5SRezgar Shakeri   if (div) memcpy((*basis)->div, div, Q * P * sizeof(div[0]));
12142b730f8bSJeremy L Thompson   CeedCall(ceed->BasisCreateHdiv(topo, dim, P, Q, interp, div, q_ref, q_weight, *basis));
121550c301a5SRezgar Shakeri   return CEED_ERROR_SUCCESS;
121650c301a5SRezgar Shakeri }
121750c301a5SRezgar Shakeri 
121850c301a5SRezgar Shakeri /**
1219*c4e3f59bSSebastian Grimberg   @brief Create a non tensor-product basis for H(curl) discretizations
1220*c4e3f59bSSebastian Grimberg 
1221*c4e3f59bSSebastian Grimberg   @param[in]  ceed      Ceed object where the CeedBasis will be created
1222*c4e3f59bSSebastian Grimberg   @param[in]  topo      Topology of element (`CEED_TOPOLOGY_QUAD`, `CEED_TOPOLOGY_PRISM`, etc.), dimension of which is used in some array sizes below
1223*c4e3f59bSSebastian Grimberg   @param[in]  num_comp  Number of components (usually 1 for vectors in H(curl) bases)
1224*c4e3f59bSSebastian Grimberg   @param[in]  num_nodes Total number of nodes (dofs per element)
1225*c4e3f59bSSebastian Grimberg   @param[in]  num_qpts  Total number of quadrature points
1226*c4e3f59bSSebastian Grimberg   @param[in]  interp    Row-major (dim * num_qpts * num_nodes) matrix expressing the values of basis functions at quadrature points
1227*c4e3f59bSSebastian 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
1228*c4e3f59bSSebastian Grimberg quadrature points
1229*c4e3f59bSSebastian Grimberg   @param[in]  q_ref     Array of length num_qpts * dim holding the locations of quadrature points on the reference element
1230*c4e3f59bSSebastian Grimberg   @param[in]  q_weight  Array of length num_qpts holding the quadrature weights on the reference element
1231*c4e3f59bSSebastian Grimberg   @param[out] basis     Address of the variable where the newly created CeedBasis will be stored.
1232*c4e3f59bSSebastian Grimberg 
1233*c4e3f59bSSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
1234*c4e3f59bSSebastian Grimberg 
1235*c4e3f59bSSebastian Grimberg   @ref User
1236*c4e3f59bSSebastian Grimberg **/
1237*c4e3f59bSSebastian Grimberg int CeedBasisCreateHcurl(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
1238*c4e3f59bSSebastian Grimberg                          const CeedScalar *curl, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) {
1239*c4e3f59bSSebastian Grimberg   CeedInt Q = num_qpts, P = num_nodes, dim = 0, curl_comp = 0;
1240*c4e3f59bSSebastian Grimberg 
1241*c4e3f59bSSebastian Grimberg   if (!ceed->BasisCreateHdiv) {
1242*c4e3f59bSSebastian Grimberg     Ceed delegate;
1243*c4e3f59bSSebastian Grimberg     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
1244*c4e3f59bSSebastian Grimberg 
1245*c4e3f59bSSebastian Grimberg     if (!delegate) {
1246*c4e3f59bSSebastian Grimberg       // LCOV_EXCL_START
1247*c4e3f59bSSebastian Grimberg       return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateHcurl");
1248*c4e3f59bSSebastian Grimberg       // LCOV_EXCL_STOP
1249*c4e3f59bSSebastian Grimberg     }
1250*c4e3f59bSSebastian Grimberg 
1251*c4e3f59bSSebastian Grimberg     CeedCall(CeedBasisCreateHcurl(delegate, topo, num_comp, num_nodes, num_qpts, interp, curl, q_ref, q_weight, basis));
1252*c4e3f59bSSebastian Grimberg     return CEED_ERROR_SUCCESS;
1253*c4e3f59bSSebastian Grimberg   }
1254*c4e3f59bSSebastian Grimberg 
1255*c4e3f59bSSebastian Grimberg   if (num_comp < 1) {
1256*c4e3f59bSSebastian Grimberg     // LCOV_EXCL_START
1257*c4e3f59bSSebastian Grimberg     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component");
1258*c4e3f59bSSebastian Grimberg     // LCOV_EXCL_STOP
1259*c4e3f59bSSebastian Grimberg   }
1260*c4e3f59bSSebastian Grimberg 
1261*c4e3f59bSSebastian Grimberg   if (num_nodes < 1) {
1262*c4e3f59bSSebastian Grimberg     // LCOV_EXCL_START
1263*c4e3f59bSSebastian Grimberg     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node");
1264*c4e3f59bSSebastian Grimberg     // LCOV_EXCL_STOP
1265*c4e3f59bSSebastian Grimberg   }
1266*c4e3f59bSSebastian Grimberg 
1267*c4e3f59bSSebastian Grimberg   if (num_qpts < 1) {
1268*c4e3f59bSSebastian Grimberg     // LCOV_EXCL_START
1269*c4e3f59bSSebastian Grimberg     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point");
1270*c4e3f59bSSebastian Grimberg     // LCOV_EXCL_STOP
1271*c4e3f59bSSebastian Grimberg   }
1272*c4e3f59bSSebastian Grimberg 
1273*c4e3f59bSSebastian Grimberg   CeedCall(CeedCalloc(1, basis));
1274*c4e3f59bSSebastian Grimberg 
1275*c4e3f59bSSebastian Grimberg   CeedCall(CeedBasisGetTopologyDimension(topo, &dim));
1276*c4e3f59bSSebastian Grimberg   curl_comp = (dim < 3) ? 1 : dim;
1277*c4e3f59bSSebastian Grimberg 
1278*c4e3f59bSSebastian Grimberg   (*basis)->ceed = ceed;
1279*c4e3f59bSSebastian Grimberg   CeedCall(CeedReference(ceed));
1280*c4e3f59bSSebastian Grimberg   (*basis)->ref_count    = 1;
1281*c4e3f59bSSebastian Grimberg   (*basis)->tensor_basis = 0;
1282*c4e3f59bSSebastian Grimberg   (*basis)->dim          = dim;
1283*c4e3f59bSSebastian Grimberg   (*basis)->topo         = topo;
1284*c4e3f59bSSebastian Grimberg   (*basis)->num_comp     = num_comp;
1285*c4e3f59bSSebastian Grimberg   (*basis)->P            = P;
1286*c4e3f59bSSebastian Grimberg   (*basis)->Q            = Q;
1287*c4e3f59bSSebastian Grimberg   (*basis)->fe_space     = CEED_FE_SPACE_HCURL;
1288*c4e3f59bSSebastian Grimberg   CeedCall(CeedMalloc(Q * dim, &(*basis)->q_ref_1d));
1289*c4e3f59bSSebastian Grimberg   CeedCall(CeedMalloc(Q, &(*basis)->q_weight_1d));
1290*c4e3f59bSSebastian Grimberg   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0]));
1291*c4e3f59bSSebastian Grimberg   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0]));
1292*c4e3f59bSSebastian Grimberg   CeedCall(CeedMalloc(dim * Q * P, &(*basis)->interp));
1293*c4e3f59bSSebastian Grimberg   CeedCall(CeedMalloc(curl_comp * Q * P, &(*basis)->curl));
1294*c4e3f59bSSebastian Grimberg   if (interp) memcpy((*basis)->interp, interp, dim * Q * P * sizeof(interp[0]));
1295*c4e3f59bSSebastian Grimberg   if (curl) memcpy((*basis)->curl, curl, curl_comp * Q * P * sizeof(curl[0]));
1296*c4e3f59bSSebastian Grimberg   CeedCall(ceed->BasisCreateHcurl(topo, dim, P, Q, interp, curl, q_ref, q_weight, *basis));
1297*c4e3f59bSSebastian Grimberg   return CEED_ERROR_SUCCESS;
1298*c4e3f59bSSebastian Grimberg }
1299*c4e3f59bSSebastian Grimberg 
1300*c4e3f59bSSebastian Grimberg /**
1301ea61e9acSJeremy L Thompson   @brief Create a CeedBasis for projection from the nodes of `basis_from` to the nodes of `basis_to`.
1302ba59ac12SSebastian Grimberg 
1303ea61e9acSJeremy L Thompson   Only `CEED_EVAL_INTERP` and `CEED_EVAL_GRAD` will be valid for the new basis, `basis_project`.
1304ea61e9acSJeremy L Thompson   The interpolation is given by `interp_project = interp_to^+ * interp_from`, where the pesudoinverse `interp_to^+` is given by QR
1305ea61e9acSJeremy L Thompson factorization. The gradient is given by `grad_project = interp_to^+ * grad_from`. Note: `basis_from` and `basis_to` must have compatible quadrature
1306ba59ac12SSebastian Grimberg spaces.
1307ba59ac12SSebastian 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
1308ea61e9acSJeremy L Thompson `basis_from` has 3 components and `basis_to` has 5 components, then `basis_project` will have 3 components.
1309f113e5dcSJeremy L Thompson 
1310f113e5dcSJeremy L Thompson   @param[in]  basis_from    CeedBasis to prolong from
1311446e7af4SJeremy L Thompson   @param[in]  basis_to      CeedBasis to prolong to
1312ea61e9acSJeremy L Thompson   @param[out] basis_project Address of the variable where the newly created CeedBasis will be stored.
1313f113e5dcSJeremy L Thompson 
1314f113e5dcSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1315f113e5dcSJeremy L Thompson 
1316f113e5dcSJeremy L Thompson   @ref User
1317f113e5dcSJeremy L Thompson **/
13182b730f8bSJeremy L Thompson int CeedBasisCreateProjection(CeedBasis basis_from, CeedBasis basis_to, CeedBasis *basis_project) {
1319f113e5dcSJeremy L Thompson   Ceed ceed;
13202b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetCeed(basis_to, &ceed));
1321f113e5dcSJeremy L Thompson 
1322ecc88aebSJeremy L Thompson   // Create projection matrix
132314556e63SJeremy L Thompson   CeedScalar *interp_project, *grad_project;
13242b730f8bSJeremy L Thompson   CeedCall(CeedBasisCreateProjectionMatrices(basis_from, basis_to, &interp_project, &grad_project));
1325f113e5dcSJeremy L Thompson 
1326f113e5dcSJeremy L Thompson   // Build basis
1327f113e5dcSJeremy L Thompson   bool        is_tensor;
1328f113e5dcSJeremy L Thompson   CeedInt     dim, num_comp;
132914556e63SJeremy L Thompson   CeedScalar *q_ref, *q_weight;
13302b730f8bSJeremy L Thompson   CeedCall(CeedBasisIsTensor(basis_to, &is_tensor));
13312b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetDimension(basis_to, &dim));
13322b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumComponents(basis_from, &num_comp));
1333f113e5dcSJeremy L Thompson   if (is_tensor) {
1334f113e5dcSJeremy L Thompson     CeedInt P_1d_to, P_1d_from;
13352b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_1d_from));
13362b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_1d_to));
13372b730f8bSJeremy L Thompson     CeedCall(CeedCalloc(P_1d_to, &q_ref));
13382b730f8bSJeremy L Thompson     CeedCall(CeedCalloc(P_1d_to, &q_weight));
13392b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P_1d_from, P_1d_to, interp_project, grad_project, q_ref, q_weight, basis_project));
1340f113e5dcSJeremy L Thompson   } else {
1341f113e5dcSJeremy L Thompson     CeedElemTopology topo;
13422b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetTopology(basis_to, &topo));
1343f113e5dcSJeremy L Thompson     CeedInt num_nodes_to, num_nodes_from;
13442b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis_from, &num_nodes_from));
13452b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis_to, &num_nodes_to));
13462b730f8bSJeremy L Thompson     CeedCall(CeedCalloc(num_nodes_to * dim, &q_ref));
13472b730f8bSJeremy L Thompson     CeedCall(CeedCalloc(num_nodes_to, &q_weight));
13482b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateH1(ceed, topo, num_comp, num_nodes_from, num_nodes_to, interp_project, grad_project, q_ref, q_weight, basis_project));
1349f113e5dcSJeremy L Thompson   }
1350f113e5dcSJeremy L Thompson 
1351f113e5dcSJeremy L Thompson   // Cleanup
13522b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_project));
13532b730f8bSJeremy L Thompson   CeedCall(CeedFree(&grad_project));
13542b730f8bSJeremy L Thompson   CeedCall(CeedFree(&q_ref));
13552b730f8bSJeremy L Thompson   CeedCall(CeedFree(&q_weight));
1356f113e5dcSJeremy L Thompson 
1357f113e5dcSJeremy L Thompson   return CEED_ERROR_SUCCESS;
1358f113e5dcSJeremy L Thompson }
1359f113e5dcSJeremy L Thompson 
1360f113e5dcSJeremy L Thompson /**
1361ea61e9acSJeremy L Thompson   @brief Copy the pointer to a CeedBasis.
13629560d06aSjeremylt 
1363512bb800SJeremy 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.
1364512bb800SJeremy L Thompson   This CeedBasis will be destroyed if `basis_copy` is the only reference to this CeedBasis.
1365ea61e9acSJeremy L Thompson 
1366ea61e9acSJeremy L Thompson   @param[in]     basis      CeedBasis to copy reference to
1367ea61e9acSJeremy L Thompson   @param[in,out] basis_copy Variable to store copied reference
13689560d06aSjeremylt 
13699560d06aSjeremylt   @return An error code: 0 - success, otherwise - failure
13709560d06aSjeremylt 
13719560d06aSjeremylt   @ref User
13729560d06aSjeremylt **/
13739560d06aSjeremylt int CeedBasisReferenceCopy(CeedBasis basis, CeedBasis *basis_copy) {
1374393ac2cdSJeremy L Thompson   if (basis != CEED_BASIS_COLLOCATED) CeedCall(CeedBasisReference(basis));
13752b730f8bSJeremy L Thompson   CeedCall(CeedBasisDestroy(basis_copy));
13769560d06aSjeremylt   *basis_copy = basis;
13779560d06aSjeremylt   return CEED_ERROR_SUCCESS;
13789560d06aSjeremylt }
13799560d06aSjeremylt 
13809560d06aSjeremylt /**
13817a982d89SJeremy L. Thompson   @brief View a CeedBasis
13827a982d89SJeremy L. Thompson 
1383ea61e9acSJeremy L Thompson   @param[in] basis  CeedBasis to view
1384ea61e9acSJeremy L Thompson   @param[in] stream Stream to view to, e.g., stdout
13857a982d89SJeremy L. Thompson 
13867a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
13877a982d89SJeremy L. Thompson 
13887a982d89SJeremy L. Thompson   @ref User
13897a982d89SJeremy L. Thompson **/
13907a982d89SJeremy L. Thompson int CeedBasisView(CeedBasis basis, FILE *stream) {
139150c301a5SRezgar Shakeri   CeedElemTopology topo     = basis->topo;
1392*c4e3f59bSSebastian Grimberg   CeedFESpace      fe_space = basis->fe_space;
1393*c4e3f59bSSebastian Grimberg   CeedInt          q_comp   = 0;
13942b730f8bSJeremy L Thompson 
139550c301a5SRezgar Shakeri   // Print FE space and element topology of the basis
1396d1d35e2fSjeremylt   if (basis->tensor_basis) {
1397*c4e3f59bSSebastian Grimberg     fprintf(stream, "CeedBasis (%s on a %s element): dim=%" CeedInt_FMT " P=%" CeedInt_FMT " Q=%" CeedInt_FMT "\n", CeedFESpaces[fe_space],
13982b730f8bSJeremy L Thompson             CeedElemTopologies[topo], basis->dim, basis->P_1d, basis->Q_1d);
139950c301a5SRezgar Shakeri   } else {
1400*c4e3f59bSSebastian Grimberg     fprintf(stream, "CeedBasis (%s on a %s element): dim=%" CeedInt_FMT " P=%" CeedInt_FMT " Q=%" CeedInt_FMT "\n", CeedFESpaces[fe_space],
14012b730f8bSJeremy L Thompson             CeedElemTopologies[topo], basis->dim, basis->P, basis->Q);
140250c301a5SRezgar Shakeri   }
1403ea61e9acSJeremy L Thompson   // Print quadrature data, interpolation/gradient/divergence/curl of the basis
140450c301a5SRezgar Shakeri   if (basis->tensor_basis) {  // tensor basis
14052b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_ref_1d, stream));
14062b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_weight_1d, stream));
14072b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("interp1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, basis->interp_1d, stream));
14082b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("grad1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, basis->grad_1d, stream));
140950c301a5SRezgar Shakeri   } else {  // non-tensor basis
14102b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("qref", "\t% 12.8f", 1, basis->Q * basis->dim, basis->q_ref_1d, stream));
14112b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->q_weight_1d, stream));
1412*c4e3f59bSSebastian Grimberg     CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp));
1413*c4e3f59bSSebastian Grimberg     CeedCall(CeedScalarView("interp", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->interp, stream));
141450c301a5SRezgar Shakeri     if (basis->grad) {
1415*c4e3f59bSSebastian Grimberg       CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp));
1416*c4e3f59bSSebastian Grimberg       CeedCall(CeedScalarView("grad", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->grad, stream));
14177a982d89SJeremy L. Thompson     }
141850c301a5SRezgar Shakeri     if (basis->div) {
1419*c4e3f59bSSebastian Grimberg       CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_DIV, &q_comp));
1420*c4e3f59bSSebastian Grimberg       CeedCall(CeedScalarView("div", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->div, stream));
1421*c4e3f59bSSebastian Grimberg     }
1422*c4e3f59bSSebastian Grimberg     if (basis->curl) {
1423*c4e3f59bSSebastian Grimberg       CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_CURL, &q_comp));
1424*c4e3f59bSSebastian Grimberg       CeedCall(CeedScalarView("curl", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->curl, stream));
142550c301a5SRezgar Shakeri     }
142650c301a5SRezgar Shakeri   }
1427e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
14287a982d89SJeremy L. Thompson }
14297a982d89SJeremy L. Thompson 
14307a982d89SJeremy L. Thompson /**
14317a982d89SJeremy L. Thompson   @brief Apply basis evaluation from nodes to quadrature points or vice versa
14327a982d89SJeremy L. Thompson 
1433ea61e9acSJeremy L Thompson   @param[in]  basis      CeedBasis to evaluate
1434ea61e9acSJeremy L Thompson   @param[in]  num_elem   The number of elements to apply the basis evaluation to;
1435ea61e9acSJeremy L Thompson                            the backend will specify the ordering in CeedElemRestrictionCreateBlocked()
1436ea61e9acSJeremy L Thompson   @param[in]  t_mode    \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature points;
1437ea61e9acSJeremy L Thompson                           \ref CEED_TRANSPOSE to apply the transpose, mapping from quadrature points to nodes
1438ea61e9acSJeremy L Thompson   @param[in]  eval_mode \ref CEED_EVAL_NONE to use values directly,
14397a982d89SJeremy L. Thompson                           \ref CEED_EVAL_INTERP to use interpolated values,
14407a982d89SJeremy L. Thompson                           \ref CEED_EVAL_GRAD to use gradients,
1441*c4e3f59bSSebastian Grimberg                           \ref CEED_EVAL_DIV to use divergence,
1442*c4e3f59bSSebastian Grimberg                           \ref CEED_EVAL_CURL to use curl,
14437a982d89SJeremy L. Thompson                           \ref CEED_EVAL_WEIGHT to use quadrature weights.
14447a982d89SJeremy L. Thompson   @param[in]  u        Input CeedVector
14457a982d89SJeremy L. Thompson   @param[out] v        Output CeedVector
14467a982d89SJeremy L. Thompson 
14477a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
14487a982d89SJeremy L. Thompson 
14497a982d89SJeremy L. Thompson   @ref User
14507a982d89SJeremy L. Thompson **/
14512b730f8bSJeremy L Thompson int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
14521f9221feSJeremy L Thompson   CeedSize u_length = 0, v_length;
1453*c4e3f59bSSebastian Grimberg   CeedInt  dim, num_comp, q_comp, num_nodes, num_qpts;
14542b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetDimension(basis, &dim));
14552b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
1456*c4e3f59bSSebastian Grimberg   CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp));
14572b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumNodes(basis, &num_nodes));
14582b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts));
14592b730f8bSJeremy L Thompson   CeedCall(CeedVectorGetLength(v, &v_length));
14607a982d89SJeremy L. Thompson   if (u) {
14612b730f8bSJeremy L Thompson     CeedCall(CeedVectorGetLength(u, &u_length));
14627a982d89SJeremy L. Thompson   }
14637a982d89SJeremy L. Thompson 
14642b730f8bSJeremy L Thompson   if (!basis->Apply) {
1465e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
14662b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisApply");
1467e15f9bd0SJeremy L Thompson     // LCOV_EXCL_STOP
14682b730f8bSJeremy L Thompson   }
1469e15f9bd0SJeremy L Thompson 
1470e15f9bd0SJeremy L Thompson   // Check compatibility of topological and geometrical dimensions
14712b730f8bSJeremy L Thompson   if ((t_mode == CEED_TRANSPOSE && (v_length % num_nodes != 0 || u_length % num_qpts != 0)) ||
14722b730f8bSJeremy L Thompson       (t_mode == CEED_NOTRANSPOSE && (u_length % num_nodes != 0 || v_length % num_qpts != 0))) {
14738229195eSjeremylt     // LCOV_EXCL_START
14742b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_DIMENSION, "Length of input/output vectors incompatible with basis dimensions");
14758229195eSjeremylt     // LCOV_EXCL_STOP
14762b730f8bSJeremy L Thompson   }
14777a982d89SJeremy L. Thompson 
1478e15f9bd0SJeremy L Thompson   // Check vector lengths to prevent out of bounds issues
1479d1d35e2fSjeremylt   bool bad_dims = false;
1480d1d35e2fSjeremylt   switch (eval_mode) {
1481e15f9bd0SJeremy L Thompson     case CEED_EVAL_NONE:
14822b730f8bSJeremy L Thompson     case CEED_EVAL_INTERP:
14832b730f8bSJeremy L Thompson     case CEED_EVAL_GRAD:
1484*c4e3f59bSSebastian Grimberg     case CEED_EVAL_DIV:
1485*c4e3f59bSSebastian Grimberg     case CEED_EVAL_CURL:
1486*c4e3f59bSSebastian 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)) ||
1487*c4e3f59bSSebastian Grimberg                   (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem * num_qpts * num_comp * q_comp || u_length < num_elem * num_comp * num_nodes)));
1488e15f9bd0SJeremy L Thompson       break;
1489e15f9bd0SJeremy L Thompson     case CEED_EVAL_WEIGHT:
1490d1d35e2fSjeremylt       bad_dims = v_length < num_elem * num_qpts;
1491e15f9bd0SJeremy L Thompson       break;
1492e15f9bd0SJeremy L Thompson   }
14932b730f8bSJeremy L Thompson   if (bad_dims) {
1494e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
14952b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_DIMENSION, "Input/output vectors too short for basis and evaluation mode");
1496e15f9bd0SJeremy L Thompson     // LCOV_EXCL_STOP
14972b730f8bSJeremy L Thompson   }
1498e15f9bd0SJeremy L Thompson 
14992b730f8bSJeremy L Thompson   CeedCall(basis->Apply(basis, num_elem, t_mode, eval_mode, u, v));
1500e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
15017a982d89SJeremy L. Thompson }
15027a982d89SJeremy L. Thompson 
15037a982d89SJeremy L. Thompson /**
1504b7c9bbdaSJeremy L Thompson   @brief Get Ceed associated with a CeedBasis
1505b7c9bbdaSJeremy L Thompson 
1506ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
1507b7c9bbdaSJeremy L Thompson   @param[out] ceed  Variable to store Ceed
1508b7c9bbdaSJeremy L Thompson 
1509b7c9bbdaSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1510b7c9bbdaSJeremy L Thompson 
1511b7c9bbdaSJeremy L Thompson   @ref Advanced
1512b7c9bbdaSJeremy L Thompson **/
1513b7c9bbdaSJeremy L Thompson int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) {
1514b7c9bbdaSJeremy L Thompson   *ceed = basis->ceed;
1515b7c9bbdaSJeremy L Thompson   return CEED_ERROR_SUCCESS;
1516b7c9bbdaSJeremy L Thompson }
1517b7c9bbdaSJeremy L Thompson 
1518b7c9bbdaSJeremy L Thompson /**
15199d007619Sjeremylt   @brief Get dimension for given CeedBasis
15209d007619Sjeremylt 
1521ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
15229d007619Sjeremylt   @param[out] dim   Variable to store dimension of basis
15239d007619Sjeremylt 
15249d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
15259d007619Sjeremylt 
1526b7c9bbdaSJeremy L Thompson   @ref Advanced
15279d007619Sjeremylt **/
15289d007619Sjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) {
15299d007619Sjeremylt   *dim = basis->dim;
1530e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
15319d007619Sjeremylt }
15329d007619Sjeremylt 
15339d007619Sjeremylt /**
1534d99fa3c5SJeremy L Thompson   @brief Get topology for given CeedBasis
1535d99fa3c5SJeremy L Thompson 
1536ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
1537d99fa3c5SJeremy L Thompson   @param[out] topo  Variable to store topology of basis
1538d99fa3c5SJeremy L Thompson 
1539d99fa3c5SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1540d99fa3c5SJeremy L Thompson 
1541b7c9bbdaSJeremy L Thompson   @ref Advanced
1542d99fa3c5SJeremy L Thompson **/
1543d99fa3c5SJeremy L Thompson int CeedBasisGetTopology(CeedBasis basis, CeedElemTopology *topo) {
1544d99fa3c5SJeremy L Thompson   *topo = basis->topo;
1545e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1546d99fa3c5SJeremy L Thompson }
1547d99fa3c5SJeremy L Thompson 
1548d99fa3c5SJeremy L Thompson /**
15499d007619Sjeremylt   @brief Get number of components for given CeedBasis
15509d007619Sjeremylt 
1551ea61e9acSJeremy L Thompson   @param[in]  basis    CeedBasis
1552d1d35e2fSjeremylt   @param[out] num_comp Variable to store number of components of basis
15539d007619Sjeremylt 
15549d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
15559d007619Sjeremylt 
1556b7c9bbdaSJeremy L Thompson   @ref Advanced
15579d007619Sjeremylt **/
1558d1d35e2fSjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *num_comp) {
1559d1d35e2fSjeremylt   *num_comp = basis->num_comp;
1560e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
15619d007619Sjeremylt }
15629d007619Sjeremylt 
15639d007619Sjeremylt /**
15649d007619Sjeremylt   @brief Get total number of nodes (in dim dimensions) of a CeedBasis
15659d007619Sjeremylt 
1566ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
15679d007619Sjeremylt   @param[out] P     Variable to store number of nodes
15689d007619Sjeremylt 
15699d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
15709d007619Sjeremylt 
15719d007619Sjeremylt   @ref Utility
15729d007619Sjeremylt **/
15739d007619Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) {
15749d007619Sjeremylt   *P = basis->P;
1575e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
15769d007619Sjeremylt }
15779d007619Sjeremylt 
15789d007619Sjeremylt /**
15799d007619Sjeremylt   @brief Get total number of nodes (in 1 dimension) of a CeedBasis
15809d007619Sjeremylt 
1581ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
1582d1d35e2fSjeremylt   @param[out] P_1d  Variable to store number of nodes
15839d007619Sjeremylt 
15849d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
15859d007619Sjeremylt 
1586b7c9bbdaSJeremy L Thompson   @ref Advanced
15879d007619Sjeremylt **/
1588d1d35e2fSjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P_1d) {
15892b730f8bSJeremy L Thompson   if (!basis->tensor_basis) {
15909d007619Sjeremylt     // LCOV_EXCL_START
15912b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR, "Cannot supply P_1d for non-tensor basis");
15929d007619Sjeremylt     // LCOV_EXCL_STOP
15932b730f8bSJeremy L Thompson   }
15949d007619Sjeremylt 
1595d1d35e2fSjeremylt   *P_1d = basis->P_1d;
1596e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
15979d007619Sjeremylt }
15989d007619Sjeremylt 
15999d007619Sjeremylt /**
16009d007619Sjeremylt   @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis
16019d007619Sjeremylt 
1602ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
16039d007619Sjeremylt   @param[out] Q     Variable to store number of quadrature points
16049d007619Sjeremylt 
16059d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
16069d007619Sjeremylt 
16079d007619Sjeremylt   @ref Utility
16089d007619Sjeremylt **/
16099d007619Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) {
16109d007619Sjeremylt   *Q = basis->Q;
1611e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
16129d007619Sjeremylt }
16139d007619Sjeremylt 
16149d007619Sjeremylt /**
16159d007619Sjeremylt   @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis
16169d007619Sjeremylt 
1617ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
1618d1d35e2fSjeremylt   @param[out] Q_1d  Variable to store number of quadrature points
16199d007619Sjeremylt 
16209d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
16219d007619Sjeremylt 
1622b7c9bbdaSJeremy L Thompson   @ref Advanced
16239d007619Sjeremylt **/
1624d1d35e2fSjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q_1d) {
16252b730f8bSJeremy L Thompson   if (!basis->tensor_basis) {
16269d007619Sjeremylt     // LCOV_EXCL_START
16272b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR, "Cannot supply Q_1d for non-tensor basis");
16289d007619Sjeremylt     // LCOV_EXCL_STOP
16292b730f8bSJeremy L Thompson   }
16309d007619Sjeremylt 
1631d1d35e2fSjeremylt   *Q_1d = basis->Q_1d;
1632e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
16339d007619Sjeremylt }
16349d007619Sjeremylt 
16359d007619Sjeremylt /**
1636ea61e9acSJeremy L Thompson   @brief Get reference coordinates of quadrature points (in dim dimensions) of a CeedBasis
16379d007619Sjeremylt 
1638ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
1639d1d35e2fSjeremylt   @param[out] q_ref Variable to store reference coordinates of quadrature points
16409d007619Sjeremylt 
16419d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
16429d007619Sjeremylt 
1643b7c9bbdaSJeremy L Thompson   @ref Advanced
16449d007619Sjeremylt **/
1645d1d35e2fSjeremylt int CeedBasisGetQRef(CeedBasis basis, const CeedScalar **q_ref) {
1646d1d35e2fSjeremylt   *q_ref = basis->q_ref_1d;
1647e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
16489d007619Sjeremylt }
16499d007619Sjeremylt 
16509d007619Sjeremylt /**
1651ea61e9acSJeremy L Thompson   @brief Get quadrature weights of quadrature points (in dim dimensions) of a CeedBasis
16529d007619Sjeremylt 
1653ea61e9acSJeremy L Thompson   @param[in]  basis    CeedBasis
1654d1d35e2fSjeremylt   @param[out] q_weight Variable to store quadrature weights
16559d007619Sjeremylt 
16569d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
16579d007619Sjeremylt 
1658b7c9bbdaSJeremy L Thompson   @ref Advanced
16599d007619Sjeremylt **/
1660d1d35e2fSjeremylt int CeedBasisGetQWeights(CeedBasis basis, const CeedScalar **q_weight) {
1661d1d35e2fSjeremylt   *q_weight = basis->q_weight_1d;
1662e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
16639d007619Sjeremylt }
16649d007619Sjeremylt 
16659d007619Sjeremylt /**
16669d007619Sjeremylt   @brief Get interpolation matrix of a CeedBasis
16679d007619Sjeremylt 
1668ea61e9acSJeremy L Thompson   @param[in]  basis  CeedBasis
16699d007619Sjeremylt   @param[out] interp Variable to store interpolation matrix
16709d007619Sjeremylt 
16719d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
16729d007619Sjeremylt 
1673b7c9bbdaSJeremy L Thompson   @ref Advanced
16749d007619Sjeremylt **/
16756c58de82SJeremy L Thompson int CeedBasisGetInterp(CeedBasis basis, const CeedScalar **interp) {
1676d1d35e2fSjeremylt   if (!basis->interp && basis->tensor_basis) {
16779d007619Sjeremylt     // Allocate
16782b730f8bSJeremy L Thompson     CeedCall(CeedMalloc(basis->Q * basis->P, &basis->interp));
16799d007619Sjeremylt 
16809d007619Sjeremylt     // Initialize
16812b730f8bSJeremy L Thompson     for (CeedInt i = 0; i < basis->Q * basis->P; i++) basis->interp[i] = 1.0;
16829d007619Sjeremylt 
16839d007619Sjeremylt     // Calculate
16842b730f8bSJeremy L Thompson     for (CeedInt d = 0; d < basis->dim; d++) {
16852b730f8bSJeremy L Thompson       for (CeedInt qpt = 0; qpt < basis->Q; qpt++) {
16869d007619Sjeremylt         for (CeedInt node = 0; node < basis->P; node++) {
1687d1d35e2fSjeremylt           CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d;
1688d1d35e2fSjeremylt           CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d;
1689d1d35e2fSjeremylt           basis->interp[qpt * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p];
16909d007619Sjeremylt         }
16919d007619Sjeremylt       }
16922b730f8bSJeremy L Thompson     }
16932b730f8bSJeremy L Thompson   }
16949d007619Sjeremylt   *interp = basis->interp;
1695e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
16969d007619Sjeremylt }
16979d007619Sjeremylt 
16989d007619Sjeremylt /**
16999d007619Sjeremylt   @brief Get 1D interpolation matrix of a tensor product CeedBasis
17009d007619Sjeremylt 
1701ea61e9acSJeremy L Thompson   @param[in]  basis     CeedBasis
1702d1d35e2fSjeremylt   @param[out] interp_1d Variable to store interpolation matrix
17039d007619Sjeremylt 
17049d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
17059d007619Sjeremylt 
17069d007619Sjeremylt   @ref Backend
17079d007619Sjeremylt **/
1708d1d35e2fSjeremylt int CeedBasisGetInterp1D(CeedBasis basis, const CeedScalar **interp_1d) {
17092b730f8bSJeremy L Thompson   if (!basis->tensor_basis) {
17109d007619Sjeremylt     // LCOV_EXCL_START
17112b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR, "CeedBasis is not a tensor product basis.");
17129d007619Sjeremylt     // LCOV_EXCL_STOP
17132b730f8bSJeremy L Thompson   }
17149d007619Sjeremylt 
1715d1d35e2fSjeremylt   *interp_1d = basis->interp_1d;
1716e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
17179d007619Sjeremylt }
17189d007619Sjeremylt 
17199d007619Sjeremylt /**
17209d007619Sjeremylt   @brief Get gradient matrix of a CeedBasis
17219d007619Sjeremylt 
1722ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
17239d007619Sjeremylt   @param[out] grad  Variable to store gradient matrix
17249d007619Sjeremylt 
17259d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
17269d007619Sjeremylt 
1727b7c9bbdaSJeremy L Thompson   @ref Advanced
17289d007619Sjeremylt **/
17296c58de82SJeremy L Thompson int CeedBasisGetGrad(CeedBasis basis, const CeedScalar **grad) {
1730d1d35e2fSjeremylt   if (!basis->grad && basis->tensor_basis) {
17319d007619Sjeremylt     // Allocate
17322b730f8bSJeremy L Thompson     CeedCall(CeedMalloc(basis->dim * basis->Q * basis->P, &basis->grad));
17339d007619Sjeremylt 
17349d007619Sjeremylt     // Initialize
17352b730f8bSJeremy L Thompson     for (CeedInt i = 0; i < basis->dim * basis->Q * basis->P; i++) basis->grad[i] = 1.0;
17369d007619Sjeremylt 
17379d007619Sjeremylt     // Calculate
17382b730f8bSJeremy L Thompson     for (CeedInt d = 0; d < basis->dim; d++) {
17392b730f8bSJeremy L Thompson       for (CeedInt i = 0; i < basis->dim; i++) {
17402b730f8bSJeremy L Thompson         for (CeedInt qpt = 0; qpt < basis->Q; qpt++) {
17419d007619Sjeremylt           for (CeedInt node = 0; node < basis->P; node++) {
1742d1d35e2fSjeremylt             CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d;
1743d1d35e2fSjeremylt             CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d;
17442b730f8bSJeremy L Thompson             if (i == d) basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->grad_1d[q * basis->P_1d + p];
17452b730f8bSJeremy L Thompson             else basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p];
17462b730f8bSJeremy L Thompson           }
17472b730f8bSJeremy L Thompson         }
17482b730f8bSJeremy L Thompson       }
17499d007619Sjeremylt     }
17509d007619Sjeremylt   }
17519d007619Sjeremylt   *grad = basis->grad;
1752e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
17539d007619Sjeremylt }
17549d007619Sjeremylt 
17559d007619Sjeremylt /**
17569d007619Sjeremylt   @brief Get 1D gradient matrix of a tensor product CeedBasis
17579d007619Sjeremylt 
1758ea61e9acSJeremy L Thompson   @param[in]  basis   CeedBasis
1759d1d35e2fSjeremylt   @param[out] grad_1d Variable to store gradient matrix
17609d007619Sjeremylt 
17619d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
17629d007619Sjeremylt 
1763b7c9bbdaSJeremy L Thompson   @ref Advanced
17649d007619Sjeremylt **/
1765d1d35e2fSjeremylt int CeedBasisGetGrad1D(CeedBasis basis, const CeedScalar **grad_1d) {
17662b730f8bSJeremy L Thompson   if (!basis->tensor_basis) {
17679d007619Sjeremylt     // LCOV_EXCL_START
17682b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR, "CeedBasis is not a tensor product basis.");
17699d007619Sjeremylt     // LCOV_EXCL_STOP
17702b730f8bSJeremy L Thompson   }
17719d007619Sjeremylt 
1772d1d35e2fSjeremylt   *grad_1d = basis->grad_1d;
1773e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
17749d007619Sjeremylt }
17759d007619Sjeremylt 
17769d007619Sjeremylt /**
177750c301a5SRezgar Shakeri   @brief Get divergence matrix of a CeedBasis
177850c301a5SRezgar Shakeri 
1779ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
178050c301a5SRezgar Shakeri   @param[out] div   Variable to store divergence matrix
178150c301a5SRezgar Shakeri 
178250c301a5SRezgar Shakeri   @return An error code: 0 - success, otherwise - failure
178350c301a5SRezgar Shakeri 
178450c301a5SRezgar Shakeri   @ref Advanced
178550c301a5SRezgar Shakeri **/
178650c301a5SRezgar Shakeri int CeedBasisGetDiv(CeedBasis basis, const CeedScalar **div) {
17872b730f8bSJeremy L Thompson   if (!basis->div) {
178850c301a5SRezgar Shakeri     // LCOV_EXCL_START
17892b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR, "CeedBasis does not have divergence matrix.");
179050c301a5SRezgar Shakeri     // LCOV_EXCL_STOP
17912b730f8bSJeremy L Thompson   }
179250c301a5SRezgar Shakeri 
179350c301a5SRezgar Shakeri   *div = basis->div;
179450c301a5SRezgar Shakeri   return CEED_ERROR_SUCCESS;
179550c301a5SRezgar Shakeri }
179650c301a5SRezgar Shakeri 
179750c301a5SRezgar Shakeri /**
1798*c4e3f59bSSebastian Grimberg   @brief Get curl matrix of a CeedBasis
1799*c4e3f59bSSebastian Grimberg 
1800*c4e3f59bSSebastian Grimberg   @param[in]  basis CeedBasis
1801*c4e3f59bSSebastian Grimberg   @param[out] curl  Variable to store curl matrix
1802*c4e3f59bSSebastian Grimberg 
1803*c4e3f59bSSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
1804*c4e3f59bSSebastian Grimberg 
1805*c4e3f59bSSebastian Grimberg   @ref Advanced
1806*c4e3f59bSSebastian Grimberg **/
1807*c4e3f59bSSebastian Grimberg int CeedBasisGetCurl(CeedBasis basis, const CeedScalar **curl) {
1808*c4e3f59bSSebastian Grimberg   if (!basis->curl) {
1809*c4e3f59bSSebastian Grimberg     // LCOV_EXCL_START
1810*c4e3f59bSSebastian Grimberg     return CeedError(basis->ceed, CEED_ERROR_MINOR, "CeedBasis does not have curl matrix.");
1811*c4e3f59bSSebastian Grimberg     // LCOV_EXCL_STOP
1812*c4e3f59bSSebastian Grimberg   }
1813*c4e3f59bSSebastian Grimberg 
1814*c4e3f59bSSebastian Grimberg   *curl = basis->curl;
1815*c4e3f59bSSebastian Grimberg   return CEED_ERROR_SUCCESS;
1816*c4e3f59bSSebastian Grimberg }
1817*c4e3f59bSSebastian Grimberg 
1818*c4e3f59bSSebastian Grimberg /**
18197a982d89SJeremy L. Thompson   @brief Destroy a CeedBasis
18207a982d89SJeremy L. Thompson 
1821ea61e9acSJeremy L Thompson   @param[in,out] basis CeedBasis to destroy
18227a982d89SJeremy L. Thompson 
18237a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
18247a982d89SJeremy L. Thompson 
18257a982d89SJeremy L. Thompson   @ref User
18267a982d89SJeremy L. Thompson **/
18277a982d89SJeremy L. Thompson int CeedBasisDestroy(CeedBasis *basis) {
18287425e127SJeremy L Thompson   if (!*basis || *basis == CEED_BASIS_COLLOCATED || --(*basis)->ref_count > 0) {
1829ad6481ceSJeremy L Thompson     *basis = NULL;
1830ad6481ceSJeremy L Thompson     return CEED_ERROR_SUCCESS;
1831ad6481ceSJeremy L Thompson   }
18322b730f8bSJeremy L Thompson   if ((*basis)->Destroy) CeedCall((*basis)->Destroy(*basis));
18332b730f8bSJeremy L Thompson   if ((*basis)->contract) CeedCall(CeedTensorContractDestroy(&(*basis)->contract));
1834*c4e3f59bSSebastian Grimberg   CeedCall(CeedFree(&(*basis)->q_ref_1d));
1835*c4e3f59bSSebastian Grimberg   CeedCall(CeedFree(&(*basis)->q_weight_1d));
18362b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->interp));
18372b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->interp_1d));
18382b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->grad));
18392b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->grad_1d));
1840*c4e3f59bSSebastian Grimberg   CeedCall(CeedFree(&(*basis)->div));
1841*c4e3f59bSSebastian Grimberg   CeedCall(CeedFree(&(*basis)->curl));
18422b730f8bSJeremy L Thompson   CeedCall(CeedDestroy(&(*basis)->ceed));
18432b730f8bSJeremy L Thompson   CeedCall(CeedFree(basis));
1844e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
18457a982d89SJeremy L. Thompson }
18467a982d89SJeremy L. Thompson 
18477a982d89SJeremy L. Thompson /**
1848b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre quadrature
1849b11c1e72Sjeremylt 
1850ea61e9acSJeremy L Thompson   @param[in]  Q           Number of quadrature points (integrates polynomials of degree 2*Q-1 exactly)
1851d1d35e2fSjeremylt   @param[out] q_ref_1d    Array of length Q to hold the abscissa on [-1, 1]
1852d1d35e2fSjeremylt   @param[out] q_weight_1d Array of length Q to hold the weights
1853b11c1e72Sjeremylt 
1854b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1855dfdf5a53Sjeremylt 
1856dfdf5a53Sjeremylt   @ref Utility
1857b11c1e72Sjeremylt **/
18582b730f8bSJeremy L Thompson int CeedGaussQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) {
1859d7b241e6Sjeremylt   // Allocate
1860d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0 * atan(1.0);
1861d1d35e2fSjeremylt   // Build q_ref_1d, q_weight_1d
186292ae7e47SJeremy L Thompson   for (CeedInt i = 0; i <= Q / 2; i++) {
1863d7b241e6Sjeremylt     // Guess
1864d7b241e6Sjeremylt     xi = cos(PI * (CeedScalar)(2 * i + 1) / ((CeedScalar)(2 * Q)));
1865d7b241e6Sjeremylt     // Pn(xi)
1866d7b241e6Sjeremylt     P0 = 1.0;
1867d7b241e6Sjeremylt     P1 = xi;
1868d7b241e6Sjeremylt     P2 = 0.0;
186992ae7e47SJeremy L Thompson     for (CeedInt j = 2; j <= Q; j++) {
1870d7b241e6Sjeremylt       P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
1871d7b241e6Sjeremylt       P0 = P1;
1872d7b241e6Sjeremylt       P1 = P2;
1873d7b241e6Sjeremylt     }
1874d7b241e6Sjeremylt     // First Newton Step
1875d7b241e6Sjeremylt     dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
1876d7b241e6Sjeremylt     xi  = xi - P2 / dP2;
1877d7b241e6Sjeremylt     // Newton to convergence
187892ae7e47SJeremy L Thompson     for (CeedInt k = 0; k < 100 && fabs(P2) > 10 * CEED_EPSILON; k++) {
1879d7b241e6Sjeremylt       P0 = 1.0;
1880d7b241e6Sjeremylt       P1 = xi;
188192ae7e47SJeremy L Thompson       for (CeedInt j = 2; j <= Q; j++) {
1882d7b241e6Sjeremylt         P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
1883d7b241e6Sjeremylt         P0 = P1;
1884d7b241e6Sjeremylt         P1 = P2;
1885d7b241e6Sjeremylt       }
1886d7b241e6Sjeremylt       dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
1887d7b241e6Sjeremylt       xi  = xi - P2 / dP2;
1888d7b241e6Sjeremylt     }
1889d7b241e6Sjeremylt     // Save xi, wi
1890d7b241e6Sjeremylt     wi                     = 2.0 / ((1.0 - xi * xi) * dP2 * dP2);
1891d1d35e2fSjeremylt     q_weight_1d[i]         = wi;
1892d1d35e2fSjeremylt     q_weight_1d[Q - 1 - i] = wi;
1893d1d35e2fSjeremylt     q_ref_1d[i]            = -xi;
1894d1d35e2fSjeremylt     q_ref_1d[Q - 1 - i]    = xi;
1895d7b241e6Sjeremylt   }
1896e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1897d7b241e6Sjeremylt }
1898d7b241e6Sjeremylt 
1899b11c1e72Sjeremylt /**
1900b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre-Lobatto quadrature
1901b11c1e72Sjeremylt 
1902ea61e9acSJeremy L Thompson   @param[in]  Q           Number of quadrature points (integrates polynomials of degree 2*Q-3 exactly)
1903d1d35e2fSjeremylt   @param[out] q_ref_1d    Array of length Q to hold the abscissa on [-1, 1]
1904d1d35e2fSjeremylt   @param[out] q_weight_1d Array of length Q to hold the weights
1905b11c1e72Sjeremylt 
1906b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1907dfdf5a53Sjeremylt 
1908dfdf5a53Sjeremylt   @ref Utility
1909b11c1e72Sjeremylt **/
19102b730f8bSJeremy L Thompson int CeedLobattoQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) {
1911d7b241e6Sjeremylt   // Allocate
1912d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0 * atan(1.0);
1913d1d35e2fSjeremylt   // Build q_ref_1d, q_weight_1d
1914d7b241e6Sjeremylt   // Set endpoints
19152b730f8bSJeremy L Thompson   if (Q < 2) {
1916b0d62198Sjeremylt     // LCOV_EXCL_START
19172b730f8bSJeremy L Thompson     return CeedError(NULL, CEED_ERROR_DIMENSION, "Cannot create Lobatto quadrature with Q=%" CeedInt_FMT " < 2 points", Q);
1918b0d62198Sjeremylt     // LCOV_EXCL_STOP
19192b730f8bSJeremy L Thompson   }
1920d7b241e6Sjeremylt   wi = 2.0 / ((CeedScalar)(Q * (Q - 1)));
1921d1d35e2fSjeremylt   if (q_weight_1d) {
1922d1d35e2fSjeremylt     q_weight_1d[0]     = wi;
1923d1d35e2fSjeremylt     q_weight_1d[Q - 1] = wi;
1924d7b241e6Sjeremylt   }
1925d1d35e2fSjeremylt   q_ref_1d[0]     = -1.0;
1926d1d35e2fSjeremylt   q_ref_1d[Q - 1] = 1.0;
1927d7b241e6Sjeremylt   // Interior
192892ae7e47SJeremy L Thompson   for (CeedInt i = 1; i <= (Q - 1) / 2; i++) {
1929d7b241e6Sjeremylt     // Guess
1930d7b241e6Sjeremylt     xi = cos(PI * (CeedScalar)(i) / (CeedScalar)(Q - 1));
1931d7b241e6Sjeremylt     // Pn(xi)
1932d7b241e6Sjeremylt     P0 = 1.0;
1933d7b241e6Sjeremylt     P1 = xi;
1934d7b241e6Sjeremylt     P2 = 0.0;
193592ae7e47SJeremy L Thompson     for (CeedInt j = 2; j < Q; j++) {
1936d7b241e6Sjeremylt       P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
1937d7b241e6Sjeremylt       P0 = P1;
1938d7b241e6Sjeremylt       P1 = P2;
1939d7b241e6Sjeremylt     }
1940d7b241e6Sjeremylt     // First Newton step
1941d7b241e6Sjeremylt     dP2  = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
1942d7b241e6Sjeremylt     d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi);
1943d7b241e6Sjeremylt     xi   = xi - dP2 / d2P2;
1944d7b241e6Sjeremylt     // Newton to convergence
194592ae7e47SJeremy L Thompson     for (CeedInt k = 0; k < 100 && fabs(dP2) > 10 * CEED_EPSILON; k++) {
1946d7b241e6Sjeremylt       P0 = 1.0;
1947d7b241e6Sjeremylt       P1 = xi;
194892ae7e47SJeremy L Thompson       for (CeedInt j = 2; j < Q; j++) {
1949d7b241e6Sjeremylt         P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
1950d7b241e6Sjeremylt         P0 = P1;
1951d7b241e6Sjeremylt         P1 = P2;
1952d7b241e6Sjeremylt       }
1953d7b241e6Sjeremylt       dP2  = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
1954d7b241e6Sjeremylt       d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi);
1955d7b241e6Sjeremylt       xi   = xi - dP2 / d2P2;
1956d7b241e6Sjeremylt     }
1957d7b241e6Sjeremylt     // Save xi, wi
1958d7b241e6Sjeremylt     wi = 2.0 / (((CeedScalar)(Q * (Q - 1))) * P2 * P2);
1959d1d35e2fSjeremylt     if (q_weight_1d) {
1960d1d35e2fSjeremylt       q_weight_1d[i]         = wi;
1961d1d35e2fSjeremylt       q_weight_1d[Q - 1 - i] = wi;
1962d7b241e6Sjeremylt     }
1963d1d35e2fSjeremylt     q_ref_1d[i]         = -xi;
1964d1d35e2fSjeremylt     q_ref_1d[Q - 1 - i] = xi;
1965d7b241e6Sjeremylt   }
1966e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1967d7b241e6Sjeremylt }
1968d7b241e6Sjeremylt 
1969d7b241e6Sjeremylt /// @}
1970