xref: /libCEED/interface/ceed-basis.c (revision ba59ac12ef3579eea938de09e332a3a266c3cc08)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3d7b241e6Sjeremylt //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5d7b241e6Sjeremylt //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7d7b241e6Sjeremylt 
83d576824SJeremy L Thompson #include <ceed-impl.h>
92b730f8bSJeremy L Thompson #include <ceed/backend.h>
102b730f8bSJeremy L Thompson #include <ceed/ceed.h>
11d7b241e6Sjeremylt #include <math.h>
123d576824SJeremy L Thompson #include <stdbool.h>
13d7b241e6Sjeremylt #include <stdio.h>
14d7b241e6Sjeremylt #include <string.h>
15d7b241e6Sjeremylt 
167a982d89SJeremy L. Thompson /// @file
177a982d89SJeremy L. Thompson /// Implementation of CeedBasis interfaces
187a982d89SJeremy L. Thompson 
19d7b241e6Sjeremylt /// @cond DOXYGEN_SKIP
20783c99b3SValeria Barra static struct CeedBasis_private ceed_basis_collocated;
21d7b241e6Sjeremylt /// @endcond
22d7b241e6Sjeremylt 
237a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser
247a982d89SJeremy L. Thompson /// @{
257a982d89SJeremy L. Thompson 
267a982d89SJeremy L. Thompson /// Indicate that the quadrature points are collocated with the nodes
277a982d89SJeremy L. Thompson const CeedBasis CEED_BASIS_COLLOCATED = &ceed_basis_collocated;
287a982d89SJeremy L. Thompson 
297a982d89SJeremy L. Thompson /// @}
307a982d89SJeremy L. Thompson 
317a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
327a982d89SJeremy L. Thompson /// CeedBasis Library Internal Functions
337a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
347a982d89SJeremy L. Thompson /// @addtogroup CeedBasisDeveloper
357a982d89SJeremy L. Thompson /// @{
367a982d89SJeremy L. Thompson 
377a982d89SJeremy L. Thompson /**
387a982d89SJeremy L. Thompson   @brief Compute Householder reflection
397a982d89SJeremy L. Thompson 
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`.
126*ba59ac12SSebastian Grimberg 
127*ba59ac12SSebastian 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`.
129*ba59ac12SSebastian 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 /**
3466e15d496SJeremy L Thompson   @brief Estimate number of FLOPs required to apply CeedBasis in t_mode and eval_mode
3476e15d496SJeremy L Thompson 
348ea61e9acSJeremy L Thompson   @param[in]  basis     Basis to estimate FLOPs for
349ea61e9acSJeremy L Thompson   @param[in]  t_mode    Apply basis or transpose
350ea61e9acSJeremy L Thompson   @param[in]  eval_mode Basis evaluation mode
351ea61e9acSJeremy L Thompson   @param[out] flops     Address of variable to hold FLOPs estimate
3526e15d496SJeremy L Thompson 
3536e15d496SJeremy L Thompson   @ref Backend
3546e15d496SJeremy L Thompson **/
3552b730f8bSJeremy L Thompson int CeedBasisGetFlopsEstimate(CeedBasis basis, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedSize *flops) {
3566e15d496SJeremy L Thompson   bool is_tensor;
3576e15d496SJeremy L Thompson 
3582b730f8bSJeremy L Thompson   CeedCall(CeedBasisIsTensor(basis, &is_tensor));
3596e15d496SJeremy L Thompson   if (is_tensor) {
3606e15d496SJeremy L Thompson     CeedInt dim, num_comp, P_1d, Q_1d;
3612b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetDimension(basis, &dim));
3622b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
3632b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d));
3642b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
3656e15d496SJeremy L Thompson     if (t_mode == CEED_TRANSPOSE) {
3662b730f8bSJeremy L Thompson       P_1d = Q_1d;
3672b730f8bSJeremy L Thompson       Q_1d = P_1d;
3686e15d496SJeremy L Thompson     }
3696e15d496SJeremy L Thompson     CeedInt tensor_flops = 0, pre = num_comp * CeedIntPow(P_1d, dim - 1), post = 1;
3706e15d496SJeremy L Thompson     for (CeedInt d = 0; d < dim; d++) {
3716e15d496SJeremy L Thompson       tensor_flops += 2 * pre * P_1d * post * Q_1d;
3726e15d496SJeremy L Thompson       pre /= P_1d;
3736e15d496SJeremy L Thompson       post *= Q_1d;
3746e15d496SJeremy L Thompson     }
3756e15d496SJeremy L Thompson     switch (eval_mode) {
3762b730f8bSJeremy L Thompson       case CEED_EVAL_NONE:
3772b730f8bSJeremy L Thompson         *flops = 0;
3782b730f8bSJeremy L Thompson         break;
3792b730f8bSJeremy L Thompson       case CEED_EVAL_INTERP:
3802b730f8bSJeremy L Thompson         *flops = tensor_flops;
3812b730f8bSJeremy L Thompson         break;
3822b730f8bSJeremy L Thompson       case CEED_EVAL_GRAD:
3832b730f8bSJeremy L Thompson         *flops = tensor_flops * 2;
3842b730f8bSJeremy L Thompson         break;
3856e15d496SJeremy L Thompson       case CEED_EVAL_DIV:
3866e15d496SJeremy L Thompson         // LCOV_EXCL_START
3872b730f8bSJeremy L Thompson         return CeedError(basis->ceed, CEED_ERROR_INCOMPATIBLE, "Tensor CEED_EVAL_DIV not supported");
3882b730f8bSJeremy L Thompson         break;
3896e15d496SJeremy L Thompson       case CEED_EVAL_CURL:
3902b730f8bSJeremy L Thompson         return CeedError(basis->ceed, CEED_ERROR_INCOMPATIBLE, "Tensor CEED_EVAL_CURL not supported");
3912b730f8bSJeremy L Thompson         break;
3926e15d496SJeremy L Thompson       // LCOV_EXCL_STOP
3932b730f8bSJeremy L Thompson       case CEED_EVAL_WEIGHT:
3942b730f8bSJeremy L Thompson         *flops = dim * CeedIntPow(Q_1d, dim);
3952b730f8bSJeremy L Thompson         break;
3966e15d496SJeremy L Thompson     }
3976e15d496SJeremy L Thompson   } else {
3986e15d496SJeremy L Thompson     CeedInt dim, num_comp, num_nodes, num_qpts, Q_comp;
3992b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetDimension(basis, &dim));
4002b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
4012b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis, &num_nodes));
4022b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts));
4032b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumQuadratureComponents(basis, &Q_comp));
4046e15d496SJeremy L Thompson     switch (eval_mode) {
4052b730f8bSJeremy L Thompson       case CEED_EVAL_NONE:
4062b730f8bSJeremy L Thompson         *flops = 0;
4072b730f8bSJeremy L Thompson         break;
4082b730f8bSJeremy L Thompson       case CEED_EVAL_INTERP:
4092b730f8bSJeremy L Thompson         *flops = num_nodes * num_qpts * num_comp;
4102b730f8bSJeremy L Thompson         break;
4112b730f8bSJeremy L Thompson       case CEED_EVAL_GRAD:
4122b730f8bSJeremy L Thompson         *flops = num_nodes * num_qpts * num_comp * dim;
4132b730f8bSJeremy L Thompson         break;
4142b730f8bSJeremy L Thompson       case CEED_EVAL_DIV:
4152b730f8bSJeremy L Thompson         *flops = num_nodes * num_qpts;
4162b730f8bSJeremy L Thompson         break;
4172b730f8bSJeremy L Thompson       case CEED_EVAL_CURL:
4182b730f8bSJeremy L Thompson         *flops = num_nodes * num_qpts * dim;
4192b730f8bSJeremy L Thompson         break;
4202b730f8bSJeremy L Thompson       case CEED_EVAL_WEIGHT:
4212b730f8bSJeremy L Thompson         *flops = 0;
4222b730f8bSJeremy L Thompson         break;
4236e15d496SJeremy L Thompson     }
4246e15d496SJeremy L Thompson   }
4256e15d496SJeremy L Thompson 
4266e15d496SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4276e15d496SJeremy L Thompson }
4286e15d496SJeremy L Thompson 
4296e15d496SJeremy L Thompson /**
4307a982d89SJeremy L. Thompson   @brief Get dimension for given CeedElemTopology
4317a982d89SJeremy L. Thompson 
432ea61e9acSJeremy L Thompson   @param[in]  topo CeedElemTopology
4337a982d89SJeremy L. Thompson   @param[out] dim  Variable to store dimension of topology
4347a982d89SJeremy L. Thompson 
4357a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
4367a982d89SJeremy L. Thompson 
4377a982d89SJeremy L. Thompson   @ref Backend
4387a982d89SJeremy L. Thompson **/
4397a982d89SJeremy L. Thompson int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) {
4407a982d89SJeremy L. Thompson   *dim = (CeedInt)topo >> 16;
441e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4427a982d89SJeremy L. Thompson }
4437a982d89SJeremy L. Thompson 
4447a982d89SJeremy L. Thompson /**
4457a982d89SJeremy L. Thompson   @brief Get CeedTensorContract of a CeedBasis
4467a982d89SJeremy L. Thompson 
447ea61e9acSJeremy L Thompson   @param[in]  basis     CeedBasis
4487a982d89SJeremy L. Thompson   @param[out] contract  Variable to store CeedTensorContract
4497a982d89SJeremy L. Thompson 
4507a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
4517a982d89SJeremy L. Thompson 
4527a982d89SJeremy L. Thompson   @ref Backend
4537a982d89SJeremy L. Thompson **/
4547a982d89SJeremy L. Thompson int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) {
4557a982d89SJeremy L. Thompson   *contract = basis->contract;
456e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4577a982d89SJeremy L. Thompson }
4587a982d89SJeremy L. Thompson 
4597a982d89SJeremy L. Thompson /**
4607a982d89SJeremy L. Thompson   @brief Set CeedTensorContract of a CeedBasis
4617a982d89SJeremy L. Thompson 
462ea61e9acSJeremy L Thompson   @param[in,out] basis    CeedBasis
463ea61e9acSJeremy L Thompson   @param[in]     contract CeedTensorContract to set
4647a982d89SJeremy L. Thompson 
4657a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
4667a982d89SJeremy L. Thompson 
4677a982d89SJeremy L. Thompson   @ref Backend
4687a982d89SJeremy L. Thompson **/
46934359f16Sjeremylt int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract contract) {
47034359f16Sjeremylt   basis->contract = contract;
4712b730f8bSJeremy L Thompson   CeedCall(CeedTensorContractReference(contract));
472e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4737a982d89SJeremy L. Thompson }
4747a982d89SJeremy L. Thompson 
4757a982d89SJeremy L. Thompson /**
4767a982d89SJeremy L. Thompson   @brief Return a reference implementation of matrix multiplication C = A B.
477*ba59ac12SSebastian Grimberg 
478*ba59ac12SSebastian Grimberg   Note: This is a reference implementation for CPU CeedScalar pointers that is not intended for high performance.
4797a982d89SJeremy L. Thompson 
480ea61e9acSJeremy L Thompson   @param[in]  ceed  Ceed context for error handling
481d1d35e2fSjeremylt   @param[in]  mat_A Row-major matrix A
482d1d35e2fSjeremylt   @param[in]  mat_B Row-major matrix B
483d1d35e2fSjeremylt   @param[out] mat_C Row-major output matrix C
484ea61e9acSJeremy L Thompson   @param[in]  m     Number of rows of C
485ea61e9acSJeremy L Thompson   @param[in]  n     Number of columns of C
486ea61e9acSJeremy L Thompson   @param[in]  kk    Number of columns of A/rows of B
4877a982d89SJeremy L. Thompson 
4887a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
4897a982d89SJeremy L. Thompson 
4907a982d89SJeremy L. Thompson   @ref Utility
4917a982d89SJeremy L. Thompson **/
4922b730f8bSJeremy L Thompson int CeedMatrixMatrixMultiply(Ceed ceed, const CeedScalar *mat_A, const CeedScalar *mat_B, CeedScalar *mat_C, CeedInt m, CeedInt n, CeedInt kk) {
4932b730f8bSJeremy L Thompson   for (CeedInt i = 0; i < m; i++) {
4947a982d89SJeremy L. Thompson     for (CeedInt j = 0; j < n; j++) {
4957a982d89SJeremy L. Thompson       CeedScalar sum = 0;
4962b730f8bSJeremy L Thompson       for (CeedInt k = 0; k < kk; k++) sum += mat_A[k + i * kk] * mat_B[j + k * n];
497d1d35e2fSjeremylt       mat_C[j + i * n] = sum;
4987a982d89SJeremy L. Thompson     }
4992b730f8bSJeremy L Thompson   }
500e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5017a982d89SJeremy L. Thompson }
5027a982d89SJeremy L. Thompson 
503*ba59ac12SSebastian Grimberg /**
504*ba59ac12SSebastian Grimberg   @brief Return QR Factorization of a matrix
505*ba59ac12SSebastian Grimberg 
506*ba59ac12SSebastian Grimberg   @param[in]     ceed Ceed context for error handling
507*ba59ac12SSebastian Grimberg   @param[in,out] mat  Row-major matrix to be factorized in place
508*ba59ac12SSebastian Grimberg   @param[in,out] tau  Vector of length m of scaling factors
509*ba59ac12SSebastian Grimberg   @param[in]     m    Number of rows
510*ba59ac12SSebastian Grimberg   @param[in]     n    Number of columns
511*ba59ac12SSebastian Grimberg 
512*ba59ac12SSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
513*ba59ac12SSebastian Grimberg 
514*ba59ac12SSebastian Grimberg   @ref Utility
515*ba59ac12SSebastian Grimberg **/
516*ba59ac12SSebastian Grimberg int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau, CeedInt m, CeedInt n) {
517*ba59ac12SSebastian Grimberg   CeedScalar v[m];
518*ba59ac12SSebastian Grimberg 
519*ba59ac12SSebastian Grimberg   // Check matrix shape
520*ba59ac12SSebastian Grimberg   if (n > m) {
521*ba59ac12SSebastian Grimberg     // LCOV_EXCL_START
522*ba59ac12SSebastian Grimberg     return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute QR factorization with n > m");
523*ba59ac12SSebastian Grimberg     // LCOV_EXCL_STOP
524*ba59ac12SSebastian Grimberg   }
525*ba59ac12SSebastian Grimberg 
526*ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) {
527*ba59ac12SSebastian Grimberg     if (i >= m - 1) {  // last row of matrix, no reflection needed
528*ba59ac12SSebastian Grimberg       tau[i] = 0.;
529*ba59ac12SSebastian Grimberg       break;
530*ba59ac12SSebastian Grimberg     }
531*ba59ac12SSebastian Grimberg     // Calculate Householder vector, magnitude
532*ba59ac12SSebastian Grimberg     CeedScalar sigma = 0.0;
533*ba59ac12SSebastian Grimberg     v[i]             = mat[i + n * i];
534*ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < m; j++) {
535*ba59ac12SSebastian Grimberg       v[j] = mat[i + n * j];
536*ba59ac12SSebastian Grimberg       sigma += v[j] * v[j];
537*ba59ac12SSebastian Grimberg     }
538*ba59ac12SSebastian Grimberg     CeedScalar norm = sqrt(v[i] * v[i] + sigma);  // norm of v[i:m]
539*ba59ac12SSebastian Grimberg     CeedScalar R_ii = -copysign(norm, v[i]);
540*ba59ac12SSebastian Grimberg     v[i] -= R_ii;
541*ba59ac12SSebastian Grimberg     // norm of v[i:m] after modification above and scaling below
542*ba59ac12SSebastian Grimberg     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
543*ba59ac12SSebastian Grimberg     //   tau = 2 / (norm*norm)
544*ba59ac12SSebastian Grimberg     tau[i] = 2 * v[i] * v[i] / (v[i] * v[i] + sigma);
545*ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < m; j++) v[j] /= v[i];
546*ba59ac12SSebastian Grimberg 
547*ba59ac12SSebastian Grimberg     // Apply Householder reflector to lower right panel
548*ba59ac12SSebastian Grimberg     CeedHouseholderReflect(&mat[i * n + i + 1], &v[i], tau[i], m - i, n - i - 1, n, 1);
549*ba59ac12SSebastian Grimberg     // Save v
550*ba59ac12SSebastian Grimberg     mat[i + n * i] = R_ii;
551*ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < m; j++) mat[i + n * j] = v[j];
552*ba59ac12SSebastian Grimberg   }
553*ba59ac12SSebastian Grimberg   return CEED_ERROR_SUCCESS;
554*ba59ac12SSebastian Grimberg }
555*ba59ac12SSebastian Grimberg 
556*ba59ac12SSebastian Grimberg /**
557*ba59ac12SSebastian Grimberg   @brief Apply Householder Q matrix
558*ba59ac12SSebastian Grimberg 
559*ba59ac12SSebastian Grimberg   Compute mat_A = mat_Q mat_A, where mat_Q is mxm and mat_A is mxn.
560*ba59ac12SSebastian Grimberg 
561*ba59ac12SSebastian Grimberg   @param[in,out] mat_A  Matrix to apply Householder Q to, in place
562*ba59ac12SSebastian Grimberg   @param[in]     mat_Q  Householder Q matrix
563*ba59ac12SSebastian Grimberg   @param[in]     tau    Householder scaling factors
564*ba59ac12SSebastian Grimberg   @param[in]     t_mode Transpose mode for application
565*ba59ac12SSebastian Grimberg   @param[in]     m      Number of rows in A
566*ba59ac12SSebastian Grimberg   @param[in]     n      Number of columns in A
567*ba59ac12SSebastian Grimberg   @param[in]     k      Number of elementary reflectors in Q, k<m
568*ba59ac12SSebastian Grimberg   @param[in]     row    Row stride in A
569*ba59ac12SSebastian Grimberg   @param[in]     col    Col stride in A
570*ba59ac12SSebastian Grimberg 
571*ba59ac12SSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
572*ba59ac12SSebastian Grimberg 
573*ba59ac12SSebastian Grimberg   @ref Developer
574*ba59ac12SSebastian Grimberg **/
575*ba59ac12SSebastian Grimberg int CeedHouseholderApplyQ(CeedScalar *mat_A, const CeedScalar *mat_Q, const CeedScalar *tau, CeedTransposeMode t_mode, CeedInt m, CeedInt n,
576*ba59ac12SSebastian Grimberg                           CeedInt k, CeedInt row, CeedInt col) {
577*ba59ac12SSebastian Grimberg   CeedScalar *v;
578*ba59ac12SSebastian Grimberg   CeedCall(CeedMalloc(m, &v));
579*ba59ac12SSebastian Grimberg   for (CeedInt ii = 0; ii < k; ii++) {
580*ba59ac12SSebastian Grimberg     CeedInt i = t_mode == CEED_TRANSPOSE ? ii : k - 1 - ii;
581*ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < m; j++) v[j] = mat_Q[j * k + i];
582*ba59ac12SSebastian Grimberg     // Apply Householder reflector (I - tau v v^T) collo_grad_1d^T
583*ba59ac12SSebastian Grimberg     CeedCall(CeedHouseholderReflect(&mat_A[i * row], &v[i], tau[i], m - i, n, row, col));
584*ba59ac12SSebastian Grimberg   }
585*ba59ac12SSebastian Grimberg   CeedCall(CeedFree(&v));
586*ba59ac12SSebastian Grimberg   return CEED_ERROR_SUCCESS;
587*ba59ac12SSebastian Grimberg }
588*ba59ac12SSebastian Grimberg 
589*ba59ac12SSebastian Grimberg /**
590*ba59ac12SSebastian Grimberg   @brief Return symmetric Schur decomposition of the symmetric matrix mat via symmetric QR factorization
591*ba59ac12SSebastian Grimberg 
592*ba59ac12SSebastian Grimberg   @param[in]     ceed   Ceed context for error handling
593*ba59ac12SSebastian Grimberg   @param[in,out] mat    Row-major matrix to be factorized in place
594*ba59ac12SSebastian Grimberg   @param[out]    lambda Vector of length n of eigenvalues
595*ba59ac12SSebastian Grimberg   @param[in]     n      Number of rows/columns
596*ba59ac12SSebastian Grimberg 
597*ba59ac12SSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
598*ba59ac12SSebastian Grimberg 
599*ba59ac12SSebastian Grimberg   @ref Utility
600*ba59ac12SSebastian Grimberg **/
601*ba59ac12SSebastian Grimberg CeedPragmaOptimizeOff int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat, CeedScalar *lambda, CeedInt n) {
602*ba59ac12SSebastian Grimberg   // Check bounds for clang-tidy
603*ba59ac12SSebastian Grimberg   if (n < 2) {
604*ba59ac12SSebastian Grimberg     // LCOV_EXCL_START
605*ba59ac12SSebastian Grimberg     return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute symmetric Schur decomposition of scalars");
606*ba59ac12SSebastian Grimberg     // LCOV_EXCL_STOP
607*ba59ac12SSebastian Grimberg   }
608*ba59ac12SSebastian Grimberg 
609*ba59ac12SSebastian Grimberg   CeedScalar v[n - 1], tau[n - 1], mat_T[n * n];
610*ba59ac12SSebastian Grimberg 
611*ba59ac12SSebastian Grimberg   // Copy mat to mat_T and set mat to I
612*ba59ac12SSebastian Grimberg   memcpy(mat_T, mat, n * n * sizeof(mat[0]));
613*ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) {
614*ba59ac12SSebastian Grimberg     for (CeedInt j = 0; j < n; j++) mat[j + n * i] = (i == j) ? 1 : 0;
615*ba59ac12SSebastian Grimberg   }
616*ba59ac12SSebastian Grimberg 
617*ba59ac12SSebastian Grimberg   // Reduce to tridiagonal
618*ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n - 1; i++) {
619*ba59ac12SSebastian Grimberg     // Calculate Householder vector, magnitude
620*ba59ac12SSebastian Grimberg     CeedScalar sigma = 0.0;
621*ba59ac12SSebastian Grimberg     v[i]             = mat_T[i + n * (i + 1)];
622*ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < n - 1; j++) {
623*ba59ac12SSebastian Grimberg       v[j] = mat_T[i + n * (j + 1)];
624*ba59ac12SSebastian Grimberg       sigma += v[j] * v[j];
625*ba59ac12SSebastian Grimberg     }
626*ba59ac12SSebastian Grimberg     CeedScalar norm = sqrt(v[i] * v[i] + sigma);  // norm of v[i:n-1]
627*ba59ac12SSebastian Grimberg     CeedScalar R_ii = -copysign(norm, v[i]);
628*ba59ac12SSebastian Grimberg     v[i] -= R_ii;
629*ba59ac12SSebastian Grimberg     // norm of v[i:m] after modification above and scaling below
630*ba59ac12SSebastian Grimberg     //   norm = sqrt(v[i]*v[i] + sigma) / v[i];
631*ba59ac12SSebastian Grimberg     //   tau = 2 / (norm*norm)
632*ba59ac12SSebastian Grimberg     tau[i] = i == n - 2 ? 2 : 2 * v[i] * v[i] / (v[i] * v[i] + sigma);
633*ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < n - 1; j++) v[j] /= v[i];
634*ba59ac12SSebastian Grimberg 
635*ba59ac12SSebastian Grimberg     // Update sub and super diagonal
636*ba59ac12SSebastian Grimberg     for (CeedInt j = i + 2; j < n; j++) {
637*ba59ac12SSebastian Grimberg       mat_T[i + n * j] = 0;
638*ba59ac12SSebastian Grimberg       mat_T[j + n * i] = 0;
639*ba59ac12SSebastian Grimberg     }
640*ba59ac12SSebastian Grimberg     // Apply symmetric Householder reflector to lower right panel
641*ba59ac12SSebastian Grimberg     CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1);
642*ba59ac12SSebastian Grimberg     CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), 1, n);
643*ba59ac12SSebastian Grimberg 
644*ba59ac12SSebastian Grimberg     // Save v
645*ba59ac12SSebastian Grimberg     mat_T[i + n * (i + 1)] = R_ii;
646*ba59ac12SSebastian Grimberg     mat_T[(i + 1) + n * i] = R_ii;
647*ba59ac12SSebastian Grimberg     for (CeedInt j = i + 1; j < n - 1; j++) {
648*ba59ac12SSebastian Grimberg       mat_T[i + n * (j + 1)] = v[j];
649*ba59ac12SSebastian Grimberg     }
650*ba59ac12SSebastian Grimberg   }
651*ba59ac12SSebastian Grimberg   // Backwards accumulation of Q
652*ba59ac12SSebastian Grimberg   for (CeedInt i = n - 2; i >= 0; i--) {
653*ba59ac12SSebastian Grimberg     if (tau[i] > 0.0) {
654*ba59ac12SSebastian Grimberg       v[i] = 1;
655*ba59ac12SSebastian Grimberg       for (CeedInt j = i + 1; j < n - 1; j++) {
656*ba59ac12SSebastian Grimberg         v[j]                   = mat_T[i + n * (j + 1)];
657*ba59ac12SSebastian Grimberg         mat_T[i + n * (j + 1)] = 0;
658*ba59ac12SSebastian Grimberg       }
659*ba59ac12SSebastian Grimberg       CeedHouseholderReflect(&mat[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1);
660*ba59ac12SSebastian Grimberg     }
661*ba59ac12SSebastian Grimberg   }
662*ba59ac12SSebastian Grimberg 
663*ba59ac12SSebastian Grimberg   // Reduce sub and super diagonal
664*ba59ac12SSebastian Grimberg   CeedInt    p = 0, q = 0, itr = 0, max_itr = n * n * n * n;
665*ba59ac12SSebastian Grimberg   CeedScalar tol = CEED_EPSILON;
666*ba59ac12SSebastian Grimberg 
667*ba59ac12SSebastian Grimberg   while (itr < max_itr) {
668*ba59ac12SSebastian Grimberg     // Update p, q, size of reduced portions of diagonal
669*ba59ac12SSebastian Grimberg     p = 0;
670*ba59ac12SSebastian Grimberg     q = 0;
671*ba59ac12SSebastian Grimberg     for (CeedInt i = n - 2; i >= 0; i--) {
672*ba59ac12SSebastian Grimberg       if (fabs(mat_T[i + n * (i + 1)]) < tol) q += 1;
673*ba59ac12SSebastian Grimberg       else break;
674*ba59ac12SSebastian Grimberg     }
675*ba59ac12SSebastian Grimberg     for (CeedInt i = 0; i < n - q - 1; i++) {
676*ba59ac12SSebastian Grimberg       if (fabs(mat_T[i + n * (i + 1)]) < tol) p += 1;
677*ba59ac12SSebastian Grimberg       else break;
678*ba59ac12SSebastian Grimberg     }
679*ba59ac12SSebastian Grimberg     if (q == n - 1) break;  // Finished reducing
680*ba59ac12SSebastian Grimberg 
681*ba59ac12SSebastian Grimberg     // Reduce tridiagonal portion
682*ba59ac12SSebastian Grimberg     CeedScalar t_nn = mat_T[(n - 1 - q) + n * (n - 1 - q)], t_nnm1 = mat_T[(n - 2 - q) + n * (n - 1 - q)];
683*ba59ac12SSebastian Grimberg     CeedScalar d  = (mat_T[(n - 2 - q) + n * (n - 2 - q)] - t_nn) / 2;
684*ba59ac12SSebastian Grimberg     CeedScalar mu = t_nn - t_nnm1 * t_nnm1 / (d + copysign(sqrt(d * d + t_nnm1 * t_nnm1), d));
685*ba59ac12SSebastian Grimberg     CeedScalar x  = mat_T[p + n * p] - mu;
686*ba59ac12SSebastian Grimberg     CeedScalar z  = mat_T[p + n * (p + 1)];
687*ba59ac12SSebastian Grimberg     for (CeedInt k = p; k < n - q - 1; k++) {
688*ba59ac12SSebastian Grimberg       // Compute Givens rotation
689*ba59ac12SSebastian Grimberg       CeedScalar c = 1, s = 0;
690*ba59ac12SSebastian Grimberg       if (fabs(z) > tol) {
691*ba59ac12SSebastian Grimberg         if (fabs(z) > fabs(x)) {
692*ba59ac12SSebastian Grimberg           CeedScalar tau = -x / z;
693*ba59ac12SSebastian Grimberg           s = 1 / sqrt(1 + tau * tau), c = s * tau;
694*ba59ac12SSebastian Grimberg         } else {
695*ba59ac12SSebastian Grimberg           CeedScalar tau = -z / x;
696*ba59ac12SSebastian Grimberg           c = 1 / sqrt(1 + tau * tau), s = c * tau;
697*ba59ac12SSebastian Grimberg         }
698*ba59ac12SSebastian Grimberg       }
699*ba59ac12SSebastian Grimberg 
700*ba59ac12SSebastian Grimberg       // Apply Givens rotation to T
701*ba59ac12SSebastian Grimberg       CeedGivensRotation(mat_T, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n);
702*ba59ac12SSebastian Grimberg       CeedGivensRotation(mat_T, c, s, CEED_TRANSPOSE, k, k + 1, n, n);
703*ba59ac12SSebastian Grimberg 
704*ba59ac12SSebastian Grimberg       // Apply Givens rotation to Q
705*ba59ac12SSebastian Grimberg       CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n);
706*ba59ac12SSebastian Grimberg 
707*ba59ac12SSebastian Grimberg       // Update x, z
708*ba59ac12SSebastian Grimberg       if (k < n - q - 2) {
709*ba59ac12SSebastian Grimberg         x = mat_T[k + n * (k + 1)];
710*ba59ac12SSebastian Grimberg         z = mat_T[k + n * (k + 2)];
711*ba59ac12SSebastian Grimberg       }
712*ba59ac12SSebastian Grimberg     }
713*ba59ac12SSebastian Grimberg     itr++;
714*ba59ac12SSebastian Grimberg   }
715*ba59ac12SSebastian Grimberg 
716*ba59ac12SSebastian Grimberg   // Save eigenvalues
717*ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) lambda[i] = mat_T[i + n * i];
718*ba59ac12SSebastian Grimberg 
719*ba59ac12SSebastian Grimberg   // Check convergence
720*ba59ac12SSebastian Grimberg   if (itr == max_itr && q < n - 1) {
721*ba59ac12SSebastian Grimberg     // LCOV_EXCL_START
722*ba59ac12SSebastian Grimberg     return CeedError(ceed, CEED_ERROR_MINOR, "Symmetric QR failed to converge");
723*ba59ac12SSebastian Grimberg     // LCOV_EXCL_STOP
724*ba59ac12SSebastian Grimberg   }
725*ba59ac12SSebastian Grimberg   return CEED_ERROR_SUCCESS;
726*ba59ac12SSebastian Grimberg }
727*ba59ac12SSebastian Grimberg CeedPragmaOptimizeOn;
728*ba59ac12SSebastian Grimberg 
729*ba59ac12SSebastian Grimberg /**
730*ba59ac12SSebastian Grimberg   @brief Return Simultaneous Diagonalization of two matrices.
731*ba59ac12SSebastian Grimberg 
732*ba59ac12SSebastian Grimberg   This solves the generalized eigenvalue problem A x = lambda B x, where A and B are symmetric and B is positive definite.
733*ba59ac12SSebastian Grimberg   We generate the matrix X and vector Lambda such that X^T A X = Lambda and X^T B X = I.
734*ba59ac12SSebastian Grimberg   This is equivalent to the LAPACK routine 'sygv' with TYPE = 1.
735*ba59ac12SSebastian Grimberg 
736*ba59ac12SSebastian Grimberg   @param[in]  ceed   Ceed context for error handling
737*ba59ac12SSebastian Grimberg   @param[in]  mat_A  Row-major matrix to be factorized with eigenvalues
738*ba59ac12SSebastian Grimberg   @param[in]  mat_B  Row-major matrix to be factorized to identity
739*ba59ac12SSebastian Grimberg   @param[out] mat_X  Row-major orthogonal matrix
740*ba59ac12SSebastian Grimberg   @param[out] lambda Vector of length n of generalized eigenvalues
741*ba59ac12SSebastian Grimberg   @param[in]  n      Number of rows/columns
742*ba59ac12SSebastian Grimberg 
743*ba59ac12SSebastian Grimberg   @return An error code: 0 - success, otherwise - failure
744*ba59ac12SSebastian Grimberg 
745*ba59ac12SSebastian Grimberg   @ref Utility
746*ba59ac12SSebastian Grimberg **/
747*ba59ac12SSebastian Grimberg CeedPragmaOptimizeOff int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *mat_A, CeedScalar *mat_B, CeedScalar *mat_X, CeedScalar *lambda,
748*ba59ac12SSebastian Grimberg                                                           CeedInt n) {
749*ba59ac12SSebastian Grimberg   CeedScalar *mat_C, *mat_G, *vec_D;
750*ba59ac12SSebastian Grimberg   CeedCall(CeedCalloc(n * n, &mat_C));
751*ba59ac12SSebastian Grimberg   CeedCall(CeedCalloc(n * n, &mat_G));
752*ba59ac12SSebastian Grimberg   CeedCall(CeedCalloc(n, &vec_D));
753*ba59ac12SSebastian Grimberg 
754*ba59ac12SSebastian Grimberg   // Compute B = G D G^T
755*ba59ac12SSebastian Grimberg   memcpy(mat_G, mat_B, n * n * sizeof(mat_B[0]));
756*ba59ac12SSebastian Grimberg   CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_G, vec_D, n));
757*ba59ac12SSebastian Grimberg 
758*ba59ac12SSebastian Grimberg   // Sort eigenvalues
759*ba59ac12SSebastian Grimberg   for (CeedInt i = n - 1; i >= 0; i--) {
760*ba59ac12SSebastian Grimberg     for (CeedInt j = 0; j < i; j++) {
761*ba59ac12SSebastian Grimberg       if (fabs(vec_D[j]) > fabs(vec_D[j + 1])) {
762*ba59ac12SSebastian Grimberg         CeedScalar temp;
763*ba59ac12SSebastian Grimberg         temp         = vec_D[j];
764*ba59ac12SSebastian Grimberg         vec_D[j]     = vec_D[j + 1];
765*ba59ac12SSebastian Grimberg         vec_D[j + 1] = temp;
766*ba59ac12SSebastian Grimberg         for (CeedInt k = 0; k < n; k++) {
767*ba59ac12SSebastian Grimberg           temp                 = mat_G[k * n + j];
768*ba59ac12SSebastian Grimberg           mat_G[k * n + j]     = mat_G[k * n + j + 1];
769*ba59ac12SSebastian Grimberg           mat_G[k * n + j + 1] = temp;
770*ba59ac12SSebastian Grimberg         }
771*ba59ac12SSebastian Grimberg       }
772*ba59ac12SSebastian Grimberg     }
773*ba59ac12SSebastian Grimberg   }
774*ba59ac12SSebastian Grimberg 
775*ba59ac12SSebastian Grimberg   // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T
776*ba59ac12SSebastian Grimberg   //           = D^-1/2 G^T A G D^-1/2
777*ba59ac12SSebastian Grimberg   // -- D = D^-1/2
778*ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) vec_D[i] = 1. / sqrt(vec_D[i]);
779*ba59ac12SSebastian Grimberg   // -- G = G D^-1/2
780*ba59ac12SSebastian Grimberg   // -- C = D^-1/2 G^T
781*ba59ac12SSebastian Grimberg   for (CeedInt i = 0; i < n; i++) {
782*ba59ac12SSebastian Grimberg     for (CeedInt j = 0; j < n; j++) {
783*ba59ac12SSebastian Grimberg       mat_G[i * n + j] *= vec_D[j];
784*ba59ac12SSebastian Grimberg       mat_C[j * n + i] = mat_G[i * n + j];
785*ba59ac12SSebastian Grimberg     }
786*ba59ac12SSebastian Grimberg   }
787*ba59ac12SSebastian Grimberg   // -- X = (D^-1/2 G^T) A
788*ba59ac12SSebastian Grimberg   CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_C, (const CeedScalar *)mat_A, mat_X, n, n, n));
789*ba59ac12SSebastian Grimberg   // -- C = (D^-1/2 G^T A) (G D^-1/2)
790*ba59ac12SSebastian Grimberg   CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_X, (const CeedScalar *)mat_G, mat_C, n, n, n));
791*ba59ac12SSebastian Grimberg 
792*ba59ac12SSebastian Grimberg   // Compute Q^T C Q = lambda
793*ba59ac12SSebastian Grimberg   CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_C, lambda, n));
794*ba59ac12SSebastian Grimberg 
795*ba59ac12SSebastian Grimberg   // Sort eigenvalues
796*ba59ac12SSebastian Grimberg   for (CeedInt i = n - 1; i >= 0; i--) {
797*ba59ac12SSebastian Grimberg     for (CeedInt j = 0; j < i; j++) {
798*ba59ac12SSebastian Grimberg       if (fabs(lambda[j]) > fabs(lambda[j + 1])) {
799*ba59ac12SSebastian Grimberg         CeedScalar temp;
800*ba59ac12SSebastian Grimberg         temp          = lambda[j];
801*ba59ac12SSebastian Grimberg         lambda[j]     = lambda[j + 1];
802*ba59ac12SSebastian Grimberg         lambda[j + 1] = temp;
803*ba59ac12SSebastian Grimberg         for (CeedInt k = 0; k < n; k++) {
804*ba59ac12SSebastian Grimberg           temp                 = mat_C[k * n + j];
805*ba59ac12SSebastian Grimberg           mat_C[k * n + j]     = mat_C[k * n + j + 1];
806*ba59ac12SSebastian Grimberg           mat_C[k * n + j + 1] = temp;
807*ba59ac12SSebastian Grimberg         }
808*ba59ac12SSebastian Grimberg       }
809*ba59ac12SSebastian Grimberg     }
810*ba59ac12SSebastian Grimberg   }
811*ba59ac12SSebastian Grimberg 
812*ba59ac12SSebastian Grimberg   // Set X = (G D^1/2)^-T Q
813*ba59ac12SSebastian Grimberg   //       = G D^-1/2 Q
814*ba59ac12SSebastian Grimberg   CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_G, (const CeedScalar *)mat_C, mat_X, n, n, n));
815*ba59ac12SSebastian Grimberg 
816*ba59ac12SSebastian Grimberg   // Cleanup
817*ba59ac12SSebastian Grimberg   CeedCall(CeedFree(&mat_C));
818*ba59ac12SSebastian Grimberg   CeedCall(CeedFree(&mat_G));
819*ba59ac12SSebastian Grimberg   CeedCall(CeedFree(&vec_D));
820*ba59ac12SSebastian Grimberg   return CEED_ERROR_SUCCESS;
821*ba59ac12SSebastian Grimberg }
822*ba59ac12SSebastian Grimberg CeedPragmaOptimizeOn;
823*ba59ac12SSebastian Grimberg 
8247a982d89SJeremy L. Thompson /// @}
8257a982d89SJeremy L. Thompson 
8267a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
8277a982d89SJeremy L. Thompson /// CeedBasis Public API
8287a982d89SJeremy L. Thompson /// ----------------------------------------------------------------------------
8297a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser
830d7b241e6Sjeremylt /// @{
831d7b241e6Sjeremylt 
832b11c1e72Sjeremylt /**
833*ba59ac12SSebastian Grimberg   @brief Create a tensor-product basis for H^1 discretizations
834b11c1e72Sjeremylt 
835ea61e9acSJeremy L Thompson   @param[in]  ceed        Ceed object where the CeedBasis will be created
836ea61e9acSJeremy L Thompson   @param[in]  dim         Topological dimension
837ea61e9acSJeremy L Thompson   @param[in]  num_comp    Number of field components (1 for scalar fields)
838ea61e9acSJeremy L Thompson   @param[in]  P_1d        Number of nodes in one dimension
839ea61e9acSJeremy L Thompson   @param[in]  Q_1d        Number of quadrature points in one dimension
840ea61e9acSJeremy L Thompson   @param[in]  interp_1d   Row-major (Q_1d * P_1d) matrix expressing the values of nodal basis functions at quadrature points
841ea61e9acSJeremy L Thompson   @param[in]  grad_1d     Row-major (Q_1d * P_1d) matrix expressing derivatives of nodal basis functions at quadrature points
842ea61e9acSJeremy 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]
843ea61e9acSJeremy L Thompson   @param[in]  q_weight_1d Array of length Q_1d holding the quadrature weights on the reference element
844ea61e9acSJeremy L Thompson   @param[out] basis       Address of the variable where the newly created CeedBasis will be stored.
845b11c1e72Sjeremylt 
846b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
847dfdf5a53Sjeremylt 
8487a982d89SJeremy L. Thompson   @ref User
849b11c1e72Sjeremylt **/
8502b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d,
8512b730f8bSJeremy L Thompson                             const CeedScalar *grad_1d, const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis *basis) {
8525fe0d4faSjeremylt   if (!ceed->BasisCreateTensorH1) {
8535fe0d4faSjeremylt     Ceed delegate;
8542b730f8bSJeremy L Thompson     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
8555fe0d4faSjeremylt 
8562b730f8bSJeremy L Thompson     if (!delegate) {
857c042f62fSJeremy L Thompson       // LCOV_EXCL_START
8582b730f8bSJeremy L Thompson       return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisCreateTensorH1");
859c042f62fSJeremy L Thompson       // LCOV_EXCL_STOP
8602b730f8bSJeremy L Thompson     }
8615fe0d4faSjeremylt 
8622b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateTensorH1(delegate, dim, num_comp, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis));
863e15f9bd0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
8645fe0d4faSjeremylt   }
865e15f9bd0SJeremy L Thompson 
8662b730f8bSJeremy L Thompson   if (dim < 1) {
867e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
8682b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis dimension must be a positive value");
869e15f9bd0SJeremy L Thompson     // LCOV_EXCL_STOP
8702b730f8bSJeremy L Thompson   }
871227444bfSJeremy L Thompson 
8722b730f8bSJeremy L Thompson   if (num_comp < 1) {
873227444bfSJeremy L Thompson     // LCOV_EXCL_START
8742b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component");
875227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
8762b730f8bSJeremy L Thompson   }
877227444bfSJeremy L Thompson 
8782b730f8bSJeremy L Thompson   if (P_1d < 1) {
879227444bfSJeremy L Thompson     // LCOV_EXCL_START
8802b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node");
881227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
8822b730f8bSJeremy L Thompson   }
883227444bfSJeremy L Thompson 
8842b730f8bSJeremy L Thompson   if (Q_1d < 1) {
885227444bfSJeremy L Thompson     // LCOV_EXCL_START
8862b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point");
887227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
8882b730f8bSJeremy L Thompson   }
889227444bfSJeremy L Thompson 
8902b730f8bSJeremy L Thompson   CeedElemTopology topo = dim == 1 ? CEED_TOPOLOGY_LINE : dim == 2 ? CEED_TOPOLOGY_QUAD : CEED_TOPOLOGY_HEX;
891e15f9bd0SJeremy L Thompson 
8922b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(1, basis));
893d7b241e6Sjeremylt   (*basis)->ceed = ceed;
8942b730f8bSJeremy L Thompson   CeedCall(CeedReference(ceed));
895d1d35e2fSjeremylt   (*basis)->ref_count    = 1;
896d1d35e2fSjeremylt   (*basis)->tensor_basis = 1;
897d7b241e6Sjeremylt   (*basis)->dim          = dim;
898d99fa3c5SJeremy L Thompson   (*basis)->topo         = topo;
899d1d35e2fSjeremylt   (*basis)->num_comp     = num_comp;
900d1d35e2fSjeremylt   (*basis)->P_1d         = P_1d;
901d1d35e2fSjeremylt   (*basis)->Q_1d         = Q_1d;
902d1d35e2fSjeremylt   (*basis)->P            = CeedIntPow(P_1d, dim);
903d1d35e2fSjeremylt   (*basis)->Q            = CeedIntPow(Q_1d, dim);
90450c301a5SRezgar Shakeri   (*basis)->Q_comp       = 1;
90550c301a5SRezgar Shakeri   (*basis)->basis_space  = 1;  // 1 for H^1 space
9062b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d, &(*basis)->q_ref_1d));
9072b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d, &(*basis)->q_weight_1d));
908ff3a0f91SJeremy L Thompson   if (q_ref_1d) memcpy((*basis)->q_ref_1d, q_ref_1d, Q_1d * sizeof(q_ref_1d[0]));
9092b730f8bSJeremy L Thompson   if (q_weight_1d) memcpy((*basis)->q_weight_1d, q_weight_1d, Q_1d * sizeof(q_weight_1d[0]));
9102b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->interp_1d));
9112b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->grad_1d));
9122b730f8bSJeremy L Thompson   if (interp_1d) memcpy((*basis)->interp_1d, interp_1d, Q_1d * P_1d * sizeof(interp_1d[0]));
913ff3a0f91SJeremy L Thompson   if (grad_1d) memcpy((*basis)->grad_1d, grad_1d, Q_1d * P_1d * sizeof(grad_1d[0]));
9142b730f8bSJeremy L Thompson   CeedCall(ceed->BasisCreateTensorH1(dim, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, *basis));
915e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
916d7b241e6Sjeremylt }
917d7b241e6Sjeremylt 
918b11c1e72Sjeremylt /**
91995bb1877Svaleriabarra   @brief Create a tensor-product Lagrange basis
920b11c1e72Sjeremylt 
921ea61e9acSJeremy L Thompson   @param[in]  ceed      Ceed object where the CeedBasis will be created
922ea61e9acSJeremy L Thompson   @param[in]  dim       Topological dimension of element
923ea61e9acSJeremy L Thompson   @param[in]  num_comp  Number of field components (1 for scalar fields)
924ea61e9acSJeremy L Thompson   @param[in]  P         Number of Gauss-Lobatto nodes in one dimension.
925ea61e9acSJeremy L Thompson                           The polynomial degree of the resulting Q_k element is k=P-1.
926ea61e9acSJeremy L Thompson   @param[in]  Q         Number of quadrature points in one dimension.
927ea61e9acSJeremy L Thompson   @param[in]  quad_mode Distribution of the Q quadrature points (affects order of accuracy for the quadrature)
928ea61e9acSJeremy L Thompson   @param[out] basis     Address of the variable where the newly created CeedBasis will be stored.
929b11c1e72Sjeremylt 
930b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
931dfdf5a53Sjeremylt 
9327a982d89SJeremy L. Thompson   @ref User
933b11c1e72Sjeremylt **/
9342b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P, CeedInt Q, CeedQuadMode quad_mode, CeedBasis *basis) {
935d7b241e6Sjeremylt   // Allocate
9362b730f8bSJeremy L Thompson   int        ierr = CEED_ERROR_SUCCESS, i, j, k;
9372b730f8bSJeremy L Thompson   CeedScalar c1, c2, c3, c4, dx, *nodes, *interp_1d, *grad_1d, *q_ref_1d, *q_weight_1d;
9384d537eeaSYohann 
9392b730f8bSJeremy L Thompson   if (dim < 1) {
940c042f62fSJeremy L Thompson     // LCOV_EXCL_START
9412b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis dimension must be a positive value");
942c042f62fSJeremy L Thompson     // LCOV_EXCL_STOP
9432b730f8bSJeremy L Thompson   }
9444d537eeaSYohann 
9452b730f8bSJeremy L Thompson   if (num_comp < 1) {
946227444bfSJeremy L Thompson     // LCOV_EXCL_START
9472b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component");
948227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
9492b730f8bSJeremy L Thompson   }
950227444bfSJeremy L Thompson 
9512b730f8bSJeremy L Thompson   if (P < 1) {
952227444bfSJeremy L Thompson     // LCOV_EXCL_START
9532b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node");
954227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
9552b730f8bSJeremy L Thompson   }
956227444bfSJeremy L Thompson 
9572b730f8bSJeremy L Thompson   if (Q < 1) {
958227444bfSJeremy L Thompson     // LCOV_EXCL_START
9592b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point");
960227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
9612b730f8bSJeremy L Thompson   }
962227444bfSJeremy L Thompson 
963e15f9bd0SJeremy L Thompson   // Get Nodes and Weights
9642b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(P * Q, &interp_1d));
9652b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(P * Q, &grad_1d));
9662b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(P, &nodes));
9672b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q, &q_ref_1d));
9682b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q, &q_weight_1d));
9692b730f8bSJeremy L Thompson   if (CeedLobattoQuadrature(P, nodes, NULL) != CEED_ERROR_SUCCESS) goto cleanup;
970d1d35e2fSjeremylt   switch (quad_mode) {
971d7b241e6Sjeremylt     case CEED_GAUSS:
972d1d35e2fSjeremylt       ierr = CeedGaussQuadrature(Q, q_ref_1d, q_weight_1d);
973d7b241e6Sjeremylt       break;
974d7b241e6Sjeremylt     case CEED_GAUSS_LOBATTO:
975d1d35e2fSjeremylt       ierr = CeedLobattoQuadrature(Q, q_ref_1d, q_weight_1d);
976d7b241e6Sjeremylt       break;
977d7b241e6Sjeremylt   }
9782b730f8bSJeremy L Thompson   if (ierr != CEED_ERROR_SUCCESS) goto cleanup;
979e15f9bd0SJeremy L Thompson 
980d7b241e6Sjeremylt   // Build B, D matrix
981d7b241e6Sjeremylt   // Fornberg, 1998
982d7b241e6Sjeremylt   for (i = 0; i < Q; i++) {
983d7b241e6Sjeremylt     c1                   = 1.0;
984d1d35e2fSjeremylt     c3                   = nodes[0] - q_ref_1d[i];
985d1d35e2fSjeremylt     interp_1d[i * P + 0] = 1.0;
986d7b241e6Sjeremylt     for (j = 1; j < P; j++) {
987d7b241e6Sjeremylt       c2 = 1.0;
988d7b241e6Sjeremylt       c4 = c3;
989d1d35e2fSjeremylt       c3 = nodes[j] - q_ref_1d[i];
990d7b241e6Sjeremylt       for (k = 0; k < j; k++) {
991d7b241e6Sjeremylt         dx = nodes[j] - nodes[k];
992d7b241e6Sjeremylt         c2 *= dx;
993d7b241e6Sjeremylt         if (k == j - 1) {
994d1d35e2fSjeremylt           grad_1d[i * P + j]   = c1 * (interp_1d[i * P + k] - c4 * grad_1d[i * P + k]) / c2;
995d1d35e2fSjeremylt           interp_1d[i * P + j] = -c1 * c4 * interp_1d[i * P + k] / c2;
996d7b241e6Sjeremylt         }
997d1d35e2fSjeremylt         grad_1d[i * P + k]   = (c3 * grad_1d[i * P + k] - interp_1d[i * P + k]) / dx;
998d1d35e2fSjeremylt         interp_1d[i * P + k] = c3 * interp_1d[i * P + k] / dx;
999d7b241e6Sjeremylt       }
1000d7b241e6Sjeremylt       c1 = c2;
1001d7b241e6Sjeremylt     }
1002d7b241e6Sjeremylt   }
10039ac7b42eSJeremy L Thompson   // Pass to CeedBasisCreateTensorH1
10042b730f8bSJeremy L Thompson   CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P, Q, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis));
1005e15f9bd0SJeremy L Thompson cleanup:
10062b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_1d));
10072b730f8bSJeremy L Thompson   CeedCall(CeedFree(&grad_1d));
10082b730f8bSJeremy L Thompson   CeedCall(CeedFree(&nodes));
10092b730f8bSJeremy L Thompson   CeedCall(CeedFree(&q_ref_1d));
10102b730f8bSJeremy L Thompson   CeedCall(CeedFree(&q_weight_1d));
1011e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1012d7b241e6Sjeremylt }
1013d7b241e6Sjeremylt 
1014b11c1e72Sjeremylt /**
1015*ba59ac12SSebastian Grimberg   @brief Create a non tensor-product basis for H^1 discretizations
1016a8de75f0Sjeremylt 
1017ea61e9acSJeremy L Thompson   @param[in]  ceed      Ceed object where the CeedBasis will be created
1018ea61e9acSJeremy L Thompson   @param[in]  topo      Topology of element, e.g. hypercube, simplex, ect
1019ea61e9acSJeremy L Thompson   @param[in]  num_comp  Number of field components (1 for scalar fields)
1020ea61e9acSJeremy L Thompson   @param[in]  num_nodes Total number of nodes
1021ea61e9acSJeremy L Thompson   @param[in]  num_qpts  Total number of quadrature points
1022ea61e9acSJeremy L Thompson   @param[in]  interp    Row-major (num_qpts * num_nodes) matrix expressing the values of nodal basis functions at quadrature points
1023ea61e9acSJeremy L Thompson   @param[in]  grad      Row-major (num_qpts * dim * num_nodes) matrix expressing derivatives of nodal basis functions at quadrature points
10249fe083eeSJeremy L Thompson   @param[in]  q_ref     Array of length num_qpts * dim holding the locations of quadrature points on the reference element
1025ea61e9acSJeremy L Thompson   @param[in]  q_weight  Array of length num_qpts holding the quadrature weights on the reference element
1026ea61e9acSJeremy L Thompson   @param[out] basis     Address of the variable where the newly created CeedBasis will be stored.
1027a8de75f0Sjeremylt 
1028a8de75f0Sjeremylt   @return An error code: 0 - success, otherwise - failure
1029a8de75f0Sjeremylt 
10307a982d89SJeremy L. Thompson   @ref User
1031a8de75f0Sjeremylt **/
10322b730f8bSJeremy L Thompson int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
10332b730f8bSJeremy L Thompson                       const CeedScalar *grad, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) {
1034d1d35e2fSjeremylt   CeedInt P = num_nodes, Q = num_qpts, dim = 0;
1035a8de75f0Sjeremylt 
10365fe0d4faSjeremylt   if (!ceed->BasisCreateH1) {
10375fe0d4faSjeremylt     Ceed delegate;
10382b730f8bSJeremy L Thompson     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
10395fe0d4faSjeremylt 
10402b730f8bSJeremy L Thompson     if (!delegate) {
1041c042f62fSJeremy L Thompson       // LCOV_EXCL_START
10422b730f8bSJeremy L Thompson       return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisCreateH1");
1043c042f62fSJeremy L Thompson       // LCOV_EXCL_STOP
10442b730f8bSJeremy L Thompson     }
10455fe0d4faSjeremylt 
10462b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateH1(delegate, topo, num_comp, num_nodes, num_qpts, interp, grad, q_ref, q_weight, basis));
1047e15f9bd0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
10485fe0d4faSjeremylt   }
10495fe0d4faSjeremylt 
10502b730f8bSJeremy L Thompson   if (num_comp < 1) {
1051227444bfSJeremy L Thompson     // LCOV_EXCL_START
10522b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component");
1053227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
10542b730f8bSJeremy L Thompson   }
1055227444bfSJeremy L Thompson 
10562b730f8bSJeremy L Thompson   if (num_nodes < 1) {
1057227444bfSJeremy L Thompson     // LCOV_EXCL_START
10582b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node");
1059227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
10602b730f8bSJeremy L Thompson   }
1061227444bfSJeremy L Thompson 
10622b730f8bSJeremy L Thompson   if (num_qpts < 1) {
1063227444bfSJeremy L Thompson     // LCOV_EXCL_START
10642b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point");
1065227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
10662b730f8bSJeremy L Thompson   }
1067227444bfSJeremy L Thompson 
10682b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(1, basis));
1069a8de75f0Sjeremylt 
10702b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetTopologyDimension(topo, &dim));
1071a8de75f0Sjeremylt 
1072a8de75f0Sjeremylt   (*basis)->ceed = ceed;
10732b730f8bSJeremy L Thompson   CeedCall(CeedReference(ceed));
1074d1d35e2fSjeremylt   (*basis)->ref_count    = 1;
1075d1d35e2fSjeremylt   (*basis)->tensor_basis = 0;
1076a8de75f0Sjeremylt   (*basis)->dim          = dim;
1077d99fa3c5SJeremy L Thompson   (*basis)->topo         = topo;
1078d1d35e2fSjeremylt   (*basis)->num_comp     = num_comp;
1079a8de75f0Sjeremylt   (*basis)->P            = P;
1080a8de75f0Sjeremylt   (*basis)->Q            = Q;
108150c301a5SRezgar Shakeri   (*basis)->Q_comp       = 1;
108250c301a5SRezgar Shakeri   (*basis)->basis_space  = 1;  // 1 for H^1 space
10832b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q * dim, &(*basis)->q_ref_1d));
10842b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q, &(*basis)->q_weight_1d));
1085ff3a0f91SJeremy L Thompson   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0]));
1086ff3a0f91SJeremy L Thompson   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0]));
10872b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(Q * P, &(*basis)->interp));
10882b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(dim * Q * P, &(*basis)->grad));
1089ff3a0f91SJeremy L Thompson   if (interp) memcpy((*basis)->interp, interp, Q * P * sizeof(interp[0]));
1090ff3a0f91SJeremy L Thompson   if (grad) memcpy((*basis)->grad, grad, dim * Q * P * sizeof(grad[0]));
10912b730f8bSJeremy L Thompson   CeedCall(ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, q_ref, q_weight, *basis));
1092e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1093a8de75f0Sjeremylt }
1094a8de75f0Sjeremylt 
1095a8de75f0Sjeremylt /**
1096859c15bbSJames Wright   @brief Create a non tensor-product basis for \f$H(\mathrm{div})\f$ discretizations
109750c301a5SRezgar Shakeri 
1098ea61e9acSJeremy L Thompson   @param[in]  ceed      Ceed object where the CeedBasis will be created
1099ea61e9acSJeremy 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
1100ea61e9acSJeremy L Thompson   @param[in]  num_comp  Number of components (usually 1 for vectors in H(div) bases)
1101ea61e9acSJeremy L Thompson   @param[in]  num_nodes Total number of nodes (dofs per element)
1102ea61e9acSJeremy L Thompson   @param[in]  num_qpts  Total number of quadrature points
1103ea61e9acSJeremy L Thompson   @param[in]  interp    Row-major (dim*num_qpts * num_nodes) matrix expressing the values of nodal basis functions at quadrature points
1104ea61e9acSJeremy L Thompson   @param[in]  div       Row-major (num_qpts * num_nodes) matrix expressing divergence of nodal basis functions at quadrature points
11059fe083eeSJeremy L Thompson   @param[in]  q_ref     Array of length num_qpts * dim holding the locations of quadrature points on the reference element
1106ea61e9acSJeremy L Thompson   @param[in]  q_weight  Array of length num_qpts holding the quadrature weights on the reference element
1107ea61e9acSJeremy L Thompson   @param[out] basis     Address of the variable where the newly created CeedBasis will be stored.
110850c301a5SRezgar Shakeri 
110950c301a5SRezgar Shakeri   @return An error code: 0 - success, otherwise - failure
111050c301a5SRezgar Shakeri 
111150c301a5SRezgar Shakeri   @ref User
111250c301a5SRezgar Shakeri **/
11132b730f8bSJeremy L Thompson int CeedBasisCreateHdiv(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp,
11142b730f8bSJeremy L Thompson                         const CeedScalar *div, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) {
111550c301a5SRezgar Shakeri   CeedInt Q = num_qpts, P = num_nodes, dim = 0;
11162b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetTopologyDimension(topo, &dim));
111750c301a5SRezgar Shakeri   if (!ceed->BasisCreateHdiv) {
111850c301a5SRezgar Shakeri     Ceed delegate;
11192b730f8bSJeremy L Thompson     CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis"));
112050c301a5SRezgar Shakeri 
11212b730f8bSJeremy L Thompson     if (!delegate) {
112250c301a5SRezgar Shakeri       // LCOV_EXCL_START
11232b730f8bSJeremy L Thompson       return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateHdiv");
112450c301a5SRezgar Shakeri       // LCOV_EXCL_STOP
11252b730f8bSJeremy L Thompson     }
112650c301a5SRezgar Shakeri 
11272b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateHdiv(delegate, topo, num_comp, num_nodes, num_qpts, interp, div, q_ref, q_weight, basis));
112850c301a5SRezgar Shakeri     return CEED_ERROR_SUCCESS;
112950c301a5SRezgar Shakeri   }
113050c301a5SRezgar Shakeri 
11312b730f8bSJeremy L Thompson   if (num_comp < 1) {
1132227444bfSJeremy L Thompson     // LCOV_EXCL_START
11332b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component");
1134227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
11352b730f8bSJeremy L Thompson   }
1136227444bfSJeremy L Thompson 
11372b730f8bSJeremy L Thompson   if (num_nodes < 1) {
1138227444bfSJeremy L Thompson     // LCOV_EXCL_START
11392b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node");
1140227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
11412b730f8bSJeremy L Thompson   }
1142227444bfSJeremy L Thompson 
11432b730f8bSJeremy L Thompson   if (num_qpts < 1) {
1144227444bfSJeremy L Thompson     // LCOV_EXCL_START
11452b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point");
1146227444bfSJeremy L Thompson     // LCOV_EXCL_STOP
11472b730f8bSJeremy L Thompson   }
1148227444bfSJeremy L Thompson 
11492b730f8bSJeremy L Thompson   CeedCall(CeedCalloc(1, basis));
115050c301a5SRezgar Shakeri 
115150c301a5SRezgar Shakeri   (*basis)->ceed = ceed;
11522b730f8bSJeremy L Thompson   CeedCall(CeedReference(ceed));
115350c301a5SRezgar Shakeri   (*basis)->ref_count    = 1;
115450c301a5SRezgar Shakeri   (*basis)->tensor_basis = 0;
115550c301a5SRezgar Shakeri   (*basis)->dim          = dim;
115650c301a5SRezgar Shakeri   (*basis)->topo         = topo;
115750c301a5SRezgar Shakeri   (*basis)->num_comp     = num_comp;
115850c301a5SRezgar Shakeri   (*basis)->P            = P;
115950c301a5SRezgar Shakeri   (*basis)->Q            = Q;
116050c301a5SRezgar Shakeri   (*basis)->Q_comp       = dim;
116150c301a5SRezgar Shakeri   (*basis)->basis_space  = 2;  // 2 for H(div) space
11622b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q * dim, &(*basis)->q_ref_1d));
11632b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q, &(*basis)->q_weight_1d));
116450c301a5SRezgar Shakeri   if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0]));
116550c301a5SRezgar Shakeri   if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0]));
11662b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(dim * Q * P, &(*basis)->interp));
11672b730f8bSJeremy L Thompson   CeedCall(CeedMalloc(Q * P, &(*basis)->div));
116850c301a5SRezgar Shakeri   if (interp) memcpy((*basis)->interp, interp, dim * Q * P * sizeof(interp[0]));
116950c301a5SRezgar Shakeri   if (div) memcpy((*basis)->div, div, Q * P * sizeof(div[0]));
11702b730f8bSJeremy L Thompson   CeedCall(ceed->BasisCreateHdiv(topo, dim, P, Q, interp, div, q_ref, q_weight, *basis));
117150c301a5SRezgar Shakeri   return CEED_ERROR_SUCCESS;
117250c301a5SRezgar Shakeri }
117350c301a5SRezgar Shakeri 
117450c301a5SRezgar Shakeri /**
1175ea61e9acSJeremy L Thompson   @brief Create a CeedBasis for projection from the nodes of `basis_from` to the nodes of `basis_to`.
1176*ba59ac12SSebastian Grimberg 
1177ea61e9acSJeremy L Thompson   Only `CEED_EVAL_INTERP` and `CEED_EVAL_GRAD` will be valid for the new basis, `basis_project`.
1178ea61e9acSJeremy L Thompson   The interpolation is given by `interp_project = interp_to^+ * interp_from`, where the pesudoinverse `interp_to^+` is given by QR
1179ea61e9acSJeremy L Thompson factorization. The gradient is given by `grad_project = interp_to^+ * grad_from`. Note: `basis_from` and `basis_to` must have compatible quadrature
1180*ba59ac12SSebastian Grimberg spaces.
1181*ba59ac12SSebastian 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
1182ea61e9acSJeremy L Thompson `basis_from` has 3 components and `basis_to` has 5 components, then `basis_project` will have 3 components.
1183f113e5dcSJeremy L Thompson 
1184f113e5dcSJeremy L Thompson   @param[in]  basis_from    CeedBasis to prolong from
1185446e7af4SJeremy L Thompson   @param[in]  basis_to      CeedBasis to prolong to
1186ea61e9acSJeremy L Thompson   @param[out] basis_project Address of the variable where the newly created CeedBasis will be stored.
1187f113e5dcSJeremy L Thompson 
1188f113e5dcSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1189f113e5dcSJeremy L Thompson 
1190f113e5dcSJeremy L Thompson   @ref User
1191f113e5dcSJeremy L Thompson **/
11922b730f8bSJeremy L Thompson int CeedBasisCreateProjection(CeedBasis basis_from, CeedBasis basis_to, CeedBasis *basis_project) {
1193f113e5dcSJeremy L Thompson   Ceed ceed;
11942b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetCeed(basis_to, &ceed));
1195f113e5dcSJeremy L Thompson 
1196ecc88aebSJeremy L Thompson   // Create projection matrix
119714556e63SJeremy L Thompson   CeedScalar *interp_project, *grad_project;
11982b730f8bSJeremy L Thompson   CeedCall(CeedBasisCreateProjectionMatrices(basis_from, basis_to, &interp_project, &grad_project));
1199f113e5dcSJeremy L Thompson 
1200f113e5dcSJeremy L Thompson   // Build basis
1201f113e5dcSJeremy L Thompson   bool        is_tensor;
1202f113e5dcSJeremy L Thompson   CeedInt     dim, num_comp;
120314556e63SJeremy L Thompson   CeedScalar *q_ref, *q_weight;
12042b730f8bSJeremy L Thompson   CeedCall(CeedBasisIsTensor(basis_to, &is_tensor));
12052b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetDimension(basis_to, &dim));
12062b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumComponents(basis_from, &num_comp));
1207f113e5dcSJeremy L Thompson   if (is_tensor) {
1208f113e5dcSJeremy L Thompson     CeedInt P_1d_to, P_1d_from;
12092b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_1d_from));
12102b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_1d_to));
12112b730f8bSJeremy L Thompson     CeedCall(CeedCalloc(P_1d_to, &q_ref));
12122b730f8bSJeremy L Thompson     CeedCall(CeedCalloc(P_1d_to, &q_weight));
12132b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P_1d_from, P_1d_to, interp_project, grad_project, q_ref, q_weight, basis_project));
1214f113e5dcSJeremy L Thompson   } else {
1215f113e5dcSJeremy L Thompson     CeedElemTopology topo;
12162b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetTopology(basis_to, &topo));
1217f113e5dcSJeremy L Thompson     CeedInt num_nodes_to, num_nodes_from;
12182b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis_from, &num_nodes_from));
12192b730f8bSJeremy L Thompson     CeedCall(CeedBasisGetNumNodes(basis_to, &num_nodes_to));
12202b730f8bSJeremy L Thompson     CeedCall(CeedCalloc(num_nodes_to * dim, &q_ref));
12212b730f8bSJeremy L Thompson     CeedCall(CeedCalloc(num_nodes_to, &q_weight));
12222b730f8bSJeremy L Thompson     CeedCall(CeedBasisCreateH1(ceed, topo, num_comp, num_nodes_from, num_nodes_to, interp_project, grad_project, q_ref, q_weight, basis_project));
1223f113e5dcSJeremy L Thompson   }
1224f113e5dcSJeremy L Thompson 
1225f113e5dcSJeremy L Thompson   // Cleanup
12262b730f8bSJeremy L Thompson   CeedCall(CeedFree(&interp_project));
12272b730f8bSJeremy L Thompson   CeedCall(CeedFree(&grad_project));
12282b730f8bSJeremy L Thompson   CeedCall(CeedFree(&q_ref));
12292b730f8bSJeremy L Thompson   CeedCall(CeedFree(&q_weight));
1230f113e5dcSJeremy L Thompson 
1231f113e5dcSJeremy L Thompson   return CEED_ERROR_SUCCESS;
1232f113e5dcSJeremy L Thompson }
1233f113e5dcSJeremy L Thompson 
1234f113e5dcSJeremy L Thompson /**
1235ea61e9acSJeremy L Thompson   @brief Copy the pointer to a CeedBasis.
12369560d06aSjeremylt 
1237512bb800SJeremy 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.
1238512bb800SJeremy L Thompson   This CeedBasis will be destroyed if `basis_copy` is the only reference to this CeedBasis.
1239ea61e9acSJeremy L Thompson 
1240ea61e9acSJeremy L Thompson   @param[in]     basis      CeedBasis to copy reference to
1241ea61e9acSJeremy L Thompson   @param[in,out] basis_copy Variable to store copied reference
12429560d06aSjeremylt 
12439560d06aSjeremylt   @return An error code: 0 - success, otherwise - failure
12449560d06aSjeremylt 
12459560d06aSjeremylt   @ref User
12469560d06aSjeremylt **/
12479560d06aSjeremylt int CeedBasisReferenceCopy(CeedBasis basis, CeedBasis *basis_copy) {
12482b730f8bSJeremy L Thompson   CeedCall(CeedBasisReference(basis));
12492b730f8bSJeremy L Thompson   CeedCall(CeedBasisDestroy(basis_copy));
12509560d06aSjeremylt   *basis_copy = basis;
12519560d06aSjeremylt   return CEED_ERROR_SUCCESS;
12529560d06aSjeremylt }
12539560d06aSjeremylt 
12549560d06aSjeremylt /**
12557a982d89SJeremy L. Thompson   @brief View a CeedBasis
12567a982d89SJeremy L. Thompson 
1257ea61e9acSJeremy L Thompson   @param[in] basis  CeedBasis to view
1258ea61e9acSJeremy L Thompson   @param[in] stream Stream to view to, e.g., stdout
12597a982d89SJeremy L. Thompson 
12607a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
12617a982d89SJeremy L. Thompson 
12627a982d89SJeremy L. Thompson   @ref User
12637a982d89SJeremy L. Thompson **/
12647a982d89SJeremy L. Thompson int CeedBasisView(CeedBasis basis, FILE *stream) {
126550c301a5SRezgar Shakeri   CeedFESpace      FE_space = basis->basis_space;
126650c301a5SRezgar Shakeri   CeedElemTopology topo     = basis->topo;
12672b730f8bSJeremy L Thompson 
126850c301a5SRezgar Shakeri   // Print FE space and element topology of the basis
1269d1d35e2fSjeremylt   if (basis->tensor_basis) {
12702b730f8bSJeremy L Thompson     fprintf(stream, "CeedBasis (%s on a %s element): dim=%" CeedInt_FMT " P=%" CeedInt_FMT " Q=%" CeedInt_FMT "\n", CeedFESpaces[FE_space],
12712b730f8bSJeremy L Thompson             CeedElemTopologies[topo], basis->dim, basis->P_1d, basis->Q_1d);
127250c301a5SRezgar Shakeri   } else {
12732b730f8bSJeremy L Thompson     fprintf(stream, "CeedBasis (%s on a %s element): dim=%" CeedInt_FMT " P=%" CeedInt_FMT " Q=%" CeedInt_FMT "\n", CeedFESpaces[FE_space],
12742b730f8bSJeremy L Thompson             CeedElemTopologies[topo], basis->dim, basis->P, basis->Q);
127550c301a5SRezgar Shakeri   }
1276ea61e9acSJeremy L Thompson   // Print quadrature data, interpolation/gradient/divergence/curl of the basis
127750c301a5SRezgar Shakeri   if (basis->tensor_basis) {  // tensor basis
12782b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_ref_1d, stream));
12792b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_weight_1d, stream));
12802b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("interp1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, basis->interp_1d, stream));
12812b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("grad1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, basis->grad_1d, stream));
128250c301a5SRezgar Shakeri   } else {  // non-tensor basis
12832b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("qref", "\t% 12.8f", 1, basis->Q * basis->dim, basis->q_ref_1d, stream));
12842b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->q_weight_1d, stream));
12852b730f8bSJeremy L Thompson     CeedCall(CeedScalarView("interp", "\t% 12.8f", basis->Q_comp * basis->Q, basis->P, basis->interp, stream));
128650c301a5SRezgar Shakeri     if (basis->grad) {
12872b730f8bSJeremy L Thompson       CeedCall(CeedScalarView("grad", "\t% 12.8f", basis->dim * basis->Q, basis->P, basis->grad, stream));
12887a982d89SJeremy L. Thompson     }
128950c301a5SRezgar Shakeri     if (basis->div) {
12902b730f8bSJeremy L Thompson       CeedCall(CeedScalarView("div", "\t% 12.8f", basis->Q, basis->P, basis->div, stream));
129150c301a5SRezgar Shakeri     }
129250c301a5SRezgar Shakeri   }
1293e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
12947a982d89SJeremy L. Thompson }
12957a982d89SJeremy L. Thompson 
12967a982d89SJeremy L. Thompson /**
12977a982d89SJeremy L. Thompson   @brief Apply basis evaluation from nodes to quadrature points or vice versa
12987a982d89SJeremy L. Thompson 
1299ea61e9acSJeremy L Thompson   @param[in]  basis      CeedBasis to evaluate
1300ea61e9acSJeremy L Thompson   @param[in]  num_elem   The number of elements to apply the basis evaluation to;
1301ea61e9acSJeremy L Thompson                            the backend will specify the ordering in CeedElemRestrictionCreateBlocked()
1302ea61e9acSJeremy L Thompson   @param[in]  t_mode    \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature points;
1303ea61e9acSJeremy L Thompson                           \ref CEED_TRANSPOSE to apply the transpose, mapping from quadrature points to nodes
1304ea61e9acSJeremy L Thompson   @param[in]  eval_mode \ref CEED_EVAL_NONE to use values directly,
13057a982d89SJeremy L. Thompson                           \ref CEED_EVAL_INTERP to use interpolated values,
13067a982d89SJeremy L. Thompson                           \ref CEED_EVAL_GRAD to use gradients,
13077a982d89SJeremy L. Thompson                           \ref CEED_EVAL_WEIGHT to use quadrature weights.
13087a982d89SJeremy L. Thompson   @param[in]  u        Input CeedVector
13097a982d89SJeremy L. Thompson   @param[out] v        Output CeedVector
13107a982d89SJeremy L. Thompson 
13117a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
13127a982d89SJeremy L. Thompson 
13137a982d89SJeremy L. Thompson   @ref User
13147a982d89SJeremy L. Thompson **/
13152b730f8bSJeremy L Thompson int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, CeedVector v) {
13161f9221feSJeremy L Thompson   CeedSize u_length = 0, v_length;
13171f9221feSJeremy L Thompson   CeedInt  dim, num_comp, num_nodes, num_qpts;
13182b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetDimension(basis, &dim));
13192b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumComponents(basis, &num_comp));
13202b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumNodes(basis, &num_nodes));
13212b730f8bSJeremy L Thompson   CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts));
13222b730f8bSJeremy L Thompson   CeedCall(CeedVectorGetLength(v, &v_length));
13237a982d89SJeremy L. Thompson   if (u) {
13242b730f8bSJeremy L Thompson     CeedCall(CeedVectorGetLength(u, &u_length));
13257a982d89SJeremy L. Thompson   }
13267a982d89SJeremy L. Thompson 
13272b730f8bSJeremy L Thompson   if (!basis->Apply) {
1328e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
13292b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisApply");
1330e15f9bd0SJeremy L Thompson     // LCOV_EXCL_STOP
13312b730f8bSJeremy L Thompson   }
1332e15f9bd0SJeremy L Thompson 
1333e15f9bd0SJeremy L Thompson   // Check compatibility of topological and geometrical dimensions
13342b730f8bSJeremy L Thompson   if ((t_mode == CEED_TRANSPOSE && (v_length % num_nodes != 0 || u_length % num_qpts != 0)) ||
13352b730f8bSJeremy L Thompson       (t_mode == CEED_NOTRANSPOSE && (u_length % num_nodes != 0 || v_length % num_qpts != 0))) {
13368229195eSjeremylt     // LCOV_EXCL_START
13372b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_DIMENSION, "Length of input/output vectors incompatible with basis dimensions");
13388229195eSjeremylt     // LCOV_EXCL_STOP
13392b730f8bSJeremy L Thompson   }
13407a982d89SJeremy L. Thompson 
1341e15f9bd0SJeremy L Thompson   // Check vector lengths to prevent out of bounds issues
1342d1d35e2fSjeremylt   bool bad_dims = false;
1343d1d35e2fSjeremylt   switch (eval_mode) {
1344e15f9bd0SJeremy L Thompson     case CEED_EVAL_NONE:
13452b730f8bSJeremy L Thompson     case CEED_EVAL_INTERP:
13462b730f8bSJeremy L Thompson       bad_dims = ((t_mode == CEED_TRANSPOSE && (u_length < num_elem * num_comp * num_qpts || v_length < num_elem * num_comp * num_nodes)) ||
13472b730f8bSJeremy L Thompson                   (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem * num_qpts * num_comp || u_length < num_elem * num_comp * num_nodes)));
1348e15f9bd0SJeremy L Thompson       break;
13492b730f8bSJeremy L Thompson     case CEED_EVAL_GRAD:
13502b730f8bSJeremy L Thompson       bad_dims = ((t_mode == CEED_TRANSPOSE && (u_length < num_elem * num_comp * num_qpts * dim || v_length < num_elem * num_comp * num_nodes)) ||
13512b730f8bSJeremy L Thompson                   (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem * num_qpts * num_comp * dim || u_length < num_elem * num_comp * num_nodes)));
1352e15f9bd0SJeremy L Thompson       break;
1353e15f9bd0SJeremy L Thompson     case CEED_EVAL_WEIGHT:
1354d1d35e2fSjeremylt       bad_dims = v_length < num_elem * num_qpts;
1355e15f9bd0SJeremy L Thompson       break;
1356e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
13572b730f8bSJeremy L Thompson     case CEED_EVAL_DIV:
13582b730f8bSJeremy L Thompson       bad_dims = ((t_mode == CEED_TRANSPOSE && (u_length < num_elem * num_comp * num_qpts || v_length < num_elem * num_comp * num_nodes)) ||
13592b730f8bSJeremy L Thompson                   (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem * num_qpts * num_comp || u_length < num_elem * num_comp * num_nodes)));
1360e15f9bd0SJeremy L Thompson       break;
13612b730f8bSJeremy L Thompson     case CEED_EVAL_CURL:
13622b730f8bSJeremy L Thompson       bad_dims = ((t_mode == CEED_TRANSPOSE && (u_length < num_elem * num_comp * num_qpts || v_length < num_elem * num_comp * num_nodes)) ||
13632b730f8bSJeremy L Thompson                   (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem * num_qpts * num_comp || u_length < num_elem * num_comp * num_nodes)));
1364e15f9bd0SJeremy L Thompson       break;
1365e15f9bd0SJeremy L Thompson       // LCOV_EXCL_STOP
1366e15f9bd0SJeremy L Thompson   }
13672b730f8bSJeremy L Thompson   if (bad_dims) {
1368e15f9bd0SJeremy L Thompson     // LCOV_EXCL_START
13692b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_DIMENSION, "Input/output vectors too short for basis and evaluation mode");
1370e15f9bd0SJeremy L Thompson     // LCOV_EXCL_STOP
13712b730f8bSJeremy L Thompson   }
1372e15f9bd0SJeremy L Thompson 
13732b730f8bSJeremy L Thompson   CeedCall(basis->Apply(basis, num_elem, t_mode, eval_mode, u, v));
1374e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
13757a982d89SJeremy L. Thompson }
13767a982d89SJeremy L. Thompson 
13777a982d89SJeremy L. Thompson /**
1378b7c9bbdaSJeremy L Thompson   @brief Get Ceed associated with a CeedBasis
1379b7c9bbdaSJeremy L Thompson 
1380ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
1381b7c9bbdaSJeremy L Thompson   @param[out] ceed  Variable to store Ceed
1382b7c9bbdaSJeremy L Thompson 
1383b7c9bbdaSJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1384b7c9bbdaSJeremy L Thompson 
1385b7c9bbdaSJeremy L Thompson   @ref Advanced
1386b7c9bbdaSJeremy L Thompson **/
1387b7c9bbdaSJeremy L Thompson int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) {
1388b7c9bbdaSJeremy L Thompson   *ceed = basis->ceed;
1389b7c9bbdaSJeremy L Thompson   return CEED_ERROR_SUCCESS;
1390b7c9bbdaSJeremy L Thompson }
1391b7c9bbdaSJeremy L Thompson 
1392b7c9bbdaSJeremy L Thompson /**
13939d007619Sjeremylt   @brief Get dimension for given CeedBasis
13949d007619Sjeremylt 
1395ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
13969d007619Sjeremylt   @param[out] dim   Variable to store dimension of basis
13979d007619Sjeremylt 
13989d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
13999d007619Sjeremylt 
1400b7c9bbdaSJeremy L Thompson   @ref Advanced
14019d007619Sjeremylt **/
14029d007619Sjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) {
14039d007619Sjeremylt   *dim = basis->dim;
1404e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
14059d007619Sjeremylt }
14069d007619Sjeremylt 
14079d007619Sjeremylt /**
1408d99fa3c5SJeremy L Thompson   @brief Get topology for given CeedBasis
1409d99fa3c5SJeremy L Thompson 
1410ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
1411d99fa3c5SJeremy L Thompson   @param[out] topo  Variable to store topology of basis
1412d99fa3c5SJeremy L Thompson 
1413d99fa3c5SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
1414d99fa3c5SJeremy L Thompson 
1415b7c9bbdaSJeremy L Thompson   @ref Advanced
1416d99fa3c5SJeremy L Thompson **/
1417d99fa3c5SJeremy L Thompson int CeedBasisGetTopology(CeedBasis basis, CeedElemTopology *topo) {
1418d99fa3c5SJeremy L Thompson   *topo = basis->topo;
1419e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1420d99fa3c5SJeremy L Thompson }
1421d99fa3c5SJeremy L Thompson 
1422d99fa3c5SJeremy L Thompson /**
142350c301a5SRezgar Shakeri   @brief Get number of Q-vector components for given CeedBasis
142450c301a5SRezgar Shakeri 
1425ea61e9acSJeremy L Thompson   @param[in]  basis  CeedBasis
142650c301a5SRezgar Shakeri   @param[out] Q_comp Variable to store number of Q-vector components of basis
142750c301a5SRezgar Shakeri 
142850c301a5SRezgar Shakeri   @return An error code: 0 - success, otherwise - failure
142950c301a5SRezgar Shakeri 
143050c301a5SRezgar Shakeri   @ref Advanced
143150c301a5SRezgar Shakeri **/
143250c301a5SRezgar Shakeri int CeedBasisGetNumQuadratureComponents(CeedBasis basis, CeedInt *Q_comp) {
143350c301a5SRezgar Shakeri   *Q_comp = basis->Q_comp;
143450c301a5SRezgar Shakeri   return CEED_ERROR_SUCCESS;
143550c301a5SRezgar Shakeri }
143650c301a5SRezgar Shakeri 
143750c301a5SRezgar Shakeri /**
14389d007619Sjeremylt   @brief Get number of components for given CeedBasis
14399d007619Sjeremylt 
1440ea61e9acSJeremy L Thompson   @param[in]  basis    CeedBasis
1441d1d35e2fSjeremylt   @param[out] num_comp Variable to store number of components of basis
14429d007619Sjeremylt 
14439d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
14449d007619Sjeremylt 
1445b7c9bbdaSJeremy L Thompson   @ref Advanced
14469d007619Sjeremylt **/
1447d1d35e2fSjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *num_comp) {
1448d1d35e2fSjeremylt   *num_comp = basis->num_comp;
1449e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
14509d007619Sjeremylt }
14519d007619Sjeremylt 
14529d007619Sjeremylt /**
14539d007619Sjeremylt   @brief Get total number of nodes (in dim dimensions) of a CeedBasis
14549d007619Sjeremylt 
1455ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
14569d007619Sjeremylt   @param[out] P     Variable to store number of nodes
14579d007619Sjeremylt 
14589d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
14599d007619Sjeremylt 
14609d007619Sjeremylt   @ref Utility
14619d007619Sjeremylt **/
14629d007619Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) {
14639d007619Sjeremylt   *P = basis->P;
1464e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
14659d007619Sjeremylt }
14669d007619Sjeremylt 
14679d007619Sjeremylt /**
14689d007619Sjeremylt   @brief Get total number of nodes (in 1 dimension) of a CeedBasis
14699d007619Sjeremylt 
1470ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
1471d1d35e2fSjeremylt   @param[out] P_1d  Variable to store number of nodes
14729d007619Sjeremylt 
14739d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
14749d007619Sjeremylt 
1475b7c9bbdaSJeremy L Thompson   @ref Advanced
14769d007619Sjeremylt **/
1477d1d35e2fSjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P_1d) {
14782b730f8bSJeremy L Thompson   if (!basis->tensor_basis) {
14799d007619Sjeremylt     // LCOV_EXCL_START
14802b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR, "Cannot supply P_1d for non-tensor basis");
14819d007619Sjeremylt     // LCOV_EXCL_STOP
14822b730f8bSJeremy L Thompson   }
14839d007619Sjeremylt 
1484d1d35e2fSjeremylt   *P_1d = basis->P_1d;
1485e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
14869d007619Sjeremylt }
14879d007619Sjeremylt 
14889d007619Sjeremylt /**
14899d007619Sjeremylt   @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis
14909d007619Sjeremylt 
1491ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
14929d007619Sjeremylt   @param[out] Q     Variable to store number of quadrature points
14939d007619Sjeremylt 
14949d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
14959d007619Sjeremylt 
14969d007619Sjeremylt   @ref Utility
14979d007619Sjeremylt **/
14989d007619Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) {
14999d007619Sjeremylt   *Q = basis->Q;
1500e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
15019d007619Sjeremylt }
15029d007619Sjeremylt 
15039d007619Sjeremylt /**
15049d007619Sjeremylt   @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis
15059d007619Sjeremylt 
1506ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
1507d1d35e2fSjeremylt   @param[out] Q_1d  Variable to store number of quadrature points
15089d007619Sjeremylt 
15099d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
15109d007619Sjeremylt 
1511b7c9bbdaSJeremy L Thompson   @ref Advanced
15129d007619Sjeremylt **/
1513d1d35e2fSjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q_1d) {
15142b730f8bSJeremy L Thompson   if (!basis->tensor_basis) {
15159d007619Sjeremylt     // LCOV_EXCL_START
15162b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR, "Cannot supply Q_1d for non-tensor basis");
15179d007619Sjeremylt     // LCOV_EXCL_STOP
15182b730f8bSJeremy L Thompson   }
15199d007619Sjeremylt 
1520d1d35e2fSjeremylt   *Q_1d = basis->Q_1d;
1521e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
15229d007619Sjeremylt }
15239d007619Sjeremylt 
15249d007619Sjeremylt /**
1525ea61e9acSJeremy L Thompson   @brief Get reference coordinates of quadrature points (in dim dimensions) of a CeedBasis
15269d007619Sjeremylt 
1527ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
1528d1d35e2fSjeremylt   @param[out] q_ref Variable to store reference coordinates of quadrature points
15299d007619Sjeremylt 
15309d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
15319d007619Sjeremylt 
1532b7c9bbdaSJeremy L Thompson   @ref Advanced
15339d007619Sjeremylt **/
1534d1d35e2fSjeremylt int CeedBasisGetQRef(CeedBasis basis, const CeedScalar **q_ref) {
1535d1d35e2fSjeremylt   *q_ref = basis->q_ref_1d;
1536e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
15379d007619Sjeremylt }
15389d007619Sjeremylt 
15399d007619Sjeremylt /**
1540ea61e9acSJeremy L Thompson   @brief Get quadrature weights of quadrature points (in dim dimensions) of a CeedBasis
15419d007619Sjeremylt 
1542ea61e9acSJeremy L Thompson   @param[in]  basis    CeedBasis
1543d1d35e2fSjeremylt   @param[out] q_weight Variable to store quadrature weights
15449d007619Sjeremylt 
15459d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
15469d007619Sjeremylt 
1547b7c9bbdaSJeremy L Thompson   @ref Advanced
15489d007619Sjeremylt **/
1549d1d35e2fSjeremylt int CeedBasisGetQWeights(CeedBasis basis, const CeedScalar **q_weight) {
1550d1d35e2fSjeremylt   *q_weight = basis->q_weight_1d;
1551e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
15529d007619Sjeremylt }
15539d007619Sjeremylt 
15549d007619Sjeremylt /**
15559d007619Sjeremylt   @brief Get interpolation matrix of a CeedBasis
15569d007619Sjeremylt 
1557ea61e9acSJeremy L Thompson   @param[in]  basis  CeedBasis
15589d007619Sjeremylt   @param[out] interp Variable to store interpolation matrix
15599d007619Sjeremylt 
15609d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
15619d007619Sjeremylt 
1562b7c9bbdaSJeremy L Thompson   @ref Advanced
15639d007619Sjeremylt **/
15646c58de82SJeremy L Thompson int CeedBasisGetInterp(CeedBasis basis, const CeedScalar **interp) {
1565d1d35e2fSjeremylt   if (!basis->interp && basis->tensor_basis) {
15669d007619Sjeremylt     // Allocate
15672b730f8bSJeremy L Thompson     CeedCall(CeedMalloc(basis->Q * basis->P, &basis->interp));
15689d007619Sjeremylt 
15699d007619Sjeremylt     // Initialize
15702b730f8bSJeremy L Thompson     for (CeedInt i = 0; i < basis->Q * basis->P; i++) basis->interp[i] = 1.0;
15719d007619Sjeremylt 
15729d007619Sjeremylt     // Calculate
15732b730f8bSJeremy L Thompson     for (CeedInt d = 0; d < basis->dim; d++) {
15742b730f8bSJeremy L Thompson       for (CeedInt qpt = 0; qpt < basis->Q; qpt++) {
15759d007619Sjeremylt         for (CeedInt node = 0; node < basis->P; node++) {
1576d1d35e2fSjeremylt           CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d;
1577d1d35e2fSjeremylt           CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d;
1578d1d35e2fSjeremylt           basis->interp[qpt * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p];
15799d007619Sjeremylt         }
15809d007619Sjeremylt       }
15812b730f8bSJeremy L Thompson     }
15822b730f8bSJeremy L Thompson   }
15839d007619Sjeremylt   *interp = basis->interp;
1584e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
15859d007619Sjeremylt }
15869d007619Sjeremylt 
15879d007619Sjeremylt /**
15889d007619Sjeremylt   @brief Get 1D interpolation matrix of a tensor product CeedBasis
15899d007619Sjeremylt 
1590ea61e9acSJeremy L Thompson   @param[in]  basis     CeedBasis
1591d1d35e2fSjeremylt   @param[out] interp_1d Variable to store interpolation matrix
15929d007619Sjeremylt 
15939d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
15949d007619Sjeremylt 
15959d007619Sjeremylt   @ref Backend
15969d007619Sjeremylt **/
1597d1d35e2fSjeremylt int CeedBasisGetInterp1D(CeedBasis basis, const CeedScalar **interp_1d) {
15982b730f8bSJeremy L Thompson   if (!basis->tensor_basis) {
15999d007619Sjeremylt     // LCOV_EXCL_START
16002b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR, "CeedBasis is not a tensor product basis.");
16019d007619Sjeremylt     // LCOV_EXCL_STOP
16022b730f8bSJeremy L Thompson   }
16039d007619Sjeremylt 
1604d1d35e2fSjeremylt   *interp_1d = basis->interp_1d;
1605e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
16069d007619Sjeremylt }
16079d007619Sjeremylt 
16089d007619Sjeremylt /**
16099d007619Sjeremylt   @brief Get gradient matrix of a CeedBasis
16109d007619Sjeremylt 
1611ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
16129d007619Sjeremylt   @param[out] grad  Variable to store gradient matrix
16139d007619Sjeremylt 
16149d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
16159d007619Sjeremylt 
1616b7c9bbdaSJeremy L Thompson   @ref Advanced
16179d007619Sjeremylt **/
16186c58de82SJeremy L Thompson int CeedBasisGetGrad(CeedBasis basis, const CeedScalar **grad) {
1619d1d35e2fSjeremylt   if (!basis->grad && basis->tensor_basis) {
16209d007619Sjeremylt     // Allocate
16212b730f8bSJeremy L Thompson     CeedCall(CeedMalloc(basis->dim * basis->Q * basis->P, &basis->grad));
16229d007619Sjeremylt 
16239d007619Sjeremylt     // Initialize
16242b730f8bSJeremy L Thompson     for (CeedInt i = 0; i < basis->dim * basis->Q * basis->P; i++) basis->grad[i] = 1.0;
16259d007619Sjeremylt 
16269d007619Sjeremylt     // Calculate
16272b730f8bSJeremy L Thompson     for (CeedInt d = 0; d < basis->dim; d++) {
16282b730f8bSJeremy L Thompson       for (CeedInt i = 0; i < basis->dim; i++) {
16292b730f8bSJeremy L Thompson         for (CeedInt qpt = 0; qpt < basis->Q; qpt++) {
16309d007619Sjeremylt           for (CeedInt node = 0; node < basis->P; node++) {
1631d1d35e2fSjeremylt             CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d;
1632d1d35e2fSjeremylt             CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d;
16332b730f8bSJeremy L Thompson             if (i == d) basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->grad_1d[q * basis->P_1d + p];
16342b730f8bSJeremy L Thompson             else basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p];
16352b730f8bSJeremy L Thompson           }
16362b730f8bSJeremy L Thompson         }
16372b730f8bSJeremy L Thompson       }
16389d007619Sjeremylt     }
16399d007619Sjeremylt   }
16409d007619Sjeremylt   *grad = basis->grad;
1641e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
16429d007619Sjeremylt }
16439d007619Sjeremylt 
16449d007619Sjeremylt /**
16459d007619Sjeremylt   @brief Get 1D gradient matrix of a tensor product CeedBasis
16469d007619Sjeremylt 
1647ea61e9acSJeremy L Thompson   @param[in]  basis   CeedBasis
1648d1d35e2fSjeremylt   @param[out] grad_1d Variable to store gradient matrix
16499d007619Sjeremylt 
16509d007619Sjeremylt   @return An error code: 0 - success, otherwise - failure
16519d007619Sjeremylt 
1652b7c9bbdaSJeremy L Thompson   @ref Advanced
16539d007619Sjeremylt **/
1654d1d35e2fSjeremylt int CeedBasisGetGrad1D(CeedBasis basis, const CeedScalar **grad_1d) {
16552b730f8bSJeremy L Thompson   if (!basis->tensor_basis) {
16569d007619Sjeremylt     // LCOV_EXCL_START
16572b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR, "CeedBasis is not a tensor product basis.");
16589d007619Sjeremylt     // LCOV_EXCL_STOP
16592b730f8bSJeremy L Thompson   }
16609d007619Sjeremylt 
1661d1d35e2fSjeremylt   *grad_1d = basis->grad_1d;
1662e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
16639d007619Sjeremylt }
16649d007619Sjeremylt 
16659d007619Sjeremylt /**
166650c301a5SRezgar Shakeri   @brief Get divergence matrix of a CeedBasis
166750c301a5SRezgar Shakeri 
1668ea61e9acSJeremy L Thompson   @param[in]  basis CeedBasis
166950c301a5SRezgar Shakeri   @param[out] div   Variable to store divergence matrix
167050c301a5SRezgar Shakeri 
167150c301a5SRezgar Shakeri   @return An error code: 0 - success, otherwise - failure
167250c301a5SRezgar Shakeri 
167350c301a5SRezgar Shakeri   @ref Advanced
167450c301a5SRezgar Shakeri **/
167550c301a5SRezgar Shakeri int CeedBasisGetDiv(CeedBasis basis, const CeedScalar **div) {
16762b730f8bSJeremy L Thompson   if (!basis->div) {
167750c301a5SRezgar Shakeri     // LCOV_EXCL_START
16782b730f8bSJeremy L Thompson     return CeedError(basis->ceed, CEED_ERROR_MINOR, "CeedBasis does not have divergence matrix.");
167950c301a5SRezgar Shakeri     // LCOV_EXCL_STOP
16802b730f8bSJeremy L Thompson   }
168150c301a5SRezgar Shakeri 
168250c301a5SRezgar Shakeri   *div = basis->div;
168350c301a5SRezgar Shakeri   return CEED_ERROR_SUCCESS;
168450c301a5SRezgar Shakeri }
168550c301a5SRezgar Shakeri 
168650c301a5SRezgar Shakeri /**
16877a982d89SJeremy L. Thompson   @brief Destroy a CeedBasis
16887a982d89SJeremy L. Thompson 
1689ea61e9acSJeremy L Thompson   @param[in,out] basis CeedBasis to destroy
16907a982d89SJeremy L. Thompson 
16917a982d89SJeremy L. Thompson   @return An error code: 0 - success, otherwise - failure
16927a982d89SJeremy L. Thompson 
16937a982d89SJeremy L. Thompson   @ref User
16947a982d89SJeremy L. Thompson **/
16957a982d89SJeremy L. Thompson int CeedBasisDestroy(CeedBasis *basis) {
1696ad6481ceSJeremy L Thompson   if (!*basis || --(*basis)->ref_count > 0) {
1697ad6481ceSJeremy L Thompson     *basis = NULL;
1698ad6481ceSJeremy L Thompson     return CEED_ERROR_SUCCESS;
1699ad6481ceSJeremy L Thompson   }
17002b730f8bSJeremy L Thompson   if ((*basis)->Destroy) CeedCall((*basis)->Destroy(*basis));
17012b730f8bSJeremy L Thompson   if ((*basis)->contract) CeedCall(CeedTensorContractDestroy(&(*basis)->contract));
17022b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->interp));
17032b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->interp_1d));
17042b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->grad));
17052b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->div));
17062b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->grad_1d));
17072b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->q_ref_1d));
17082b730f8bSJeremy L Thompson   CeedCall(CeedFree(&(*basis)->q_weight_1d));
17092b730f8bSJeremy L Thompson   CeedCall(CeedDestroy(&(*basis)->ceed));
17102b730f8bSJeremy L Thompson   CeedCall(CeedFree(basis));
1711e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
17127a982d89SJeremy L. Thompson }
17137a982d89SJeremy L. Thompson 
17147a982d89SJeremy L. Thompson /**
1715b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre quadrature
1716b11c1e72Sjeremylt 
1717ea61e9acSJeremy L Thompson   @param[in]  Q           Number of quadrature points (integrates polynomials of degree 2*Q-1 exactly)
1718d1d35e2fSjeremylt   @param[out] q_ref_1d    Array of length Q to hold the abscissa on [-1, 1]
1719d1d35e2fSjeremylt   @param[out] q_weight_1d Array of length Q to hold the weights
1720b11c1e72Sjeremylt 
1721b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1722dfdf5a53Sjeremylt 
1723dfdf5a53Sjeremylt   @ref Utility
1724b11c1e72Sjeremylt **/
17252b730f8bSJeremy L Thompson int CeedGaussQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) {
1726d7b241e6Sjeremylt   // Allocate
1727d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0 * atan(1.0);
1728d1d35e2fSjeremylt   // Build q_ref_1d, q_weight_1d
172992ae7e47SJeremy L Thompson   for (CeedInt i = 0; i <= Q / 2; i++) {
1730d7b241e6Sjeremylt     // Guess
1731d7b241e6Sjeremylt     xi = cos(PI * (CeedScalar)(2 * i + 1) / ((CeedScalar)(2 * Q)));
1732d7b241e6Sjeremylt     // Pn(xi)
1733d7b241e6Sjeremylt     P0 = 1.0;
1734d7b241e6Sjeremylt     P1 = xi;
1735d7b241e6Sjeremylt     P2 = 0.0;
173692ae7e47SJeremy L Thompson     for (CeedInt j = 2; j <= Q; j++) {
1737d7b241e6Sjeremylt       P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
1738d7b241e6Sjeremylt       P0 = P1;
1739d7b241e6Sjeremylt       P1 = P2;
1740d7b241e6Sjeremylt     }
1741d7b241e6Sjeremylt     // First Newton Step
1742d7b241e6Sjeremylt     dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
1743d7b241e6Sjeremylt     xi  = xi - P2 / dP2;
1744d7b241e6Sjeremylt     // Newton to convergence
174592ae7e47SJeremy L Thompson     for (CeedInt k = 0; k < 100 && fabs(P2) > 10 * CEED_EPSILON; k++) {
1746d7b241e6Sjeremylt       P0 = 1.0;
1747d7b241e6Sjeremylt       P1 = xi;
174892ae7e47SJeremy L Thompson       for (CeedInt j = 2; j <= Q; j++) {
1749d7b241e6Sjeremylt         P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
1750d7b241e6Sjeremylt         P0 = P1;
1751d7b241e6Sjeremylt         P1 = P2;
1752d7b241e6Sjeremylt       }
1753d7b241e6Sjeremylt       dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
1754d7b241e6Sjeremylt       xi  = xi - P2 / dP2;
1755d7b241e6Sjeremylt     }
1756d7b241e6Sjeremylt     // Save xi, wi
1757d7b241e6Sjeremylt     wi                     = 2.0 / ((1.0 - xi * xi) * dP2 * dP2);
1758d1d35e2fSjeremylt     q_weight_1d[i]         = wi;
1759d1d35e2fSjeremylt     q_weight_1d[Q - 1 - i] = wi;
1760d1d35e2fSjeremylt     q_ref_1d[i]            = -xi;
1761d1d35e2fSjeremylt     q_ref_1d[Q - 1 - i]    = xi;
1762d7b241e6Sjeremylt   }
1763e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1764d7b241e6Sjeremylt }
1765d7b241e6Sjeremylt 
1766b11c1e72Sjeremylt /**
1767b11c1e72Sjeremylt   @brief Construct a Gauss-Legendre-Lobatto quadrature
1768b11c1e72Sjeremylt 
1769ea61e9acSJeremy L Thompson   @param[in]  Q           Number of quadrature points (integrates polynomials of degree 2*Q-3 exactly)
1770d1d35e2fSjeremylt   @param[out] q_ref_1d    Array of length Q to hold the abscissa on [-1, 1]
1771d1d35e2fSjeremylt   @param[out] q_weight_1d Array of length Q to hold the weights
1772b11c1e72Sjeremylt 
1773b11c1e72Sjeremylt   @return An error code: 0 - success, otherwise - failure
1774dfdf5a53Sjeremylt 
1775dfdf5a53Sjeremylt   @ref Utility
1776b11c1e72Sjeremylt **/
17772b730f8bSJeremy L Thompson int CeedLobattoQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) {
1778d7b241e6Sjeremylt   // Allocate
1779d7b241e6Sjeremylt   CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0 * atan(1.0);
1780d1d35e2fSjeremylt   // Build q_ref_1d, q_weight_1d
1781d7b241e6Sjeremylt   // Set endpoints
17822b730f8bSJeremy L Thompson   if (Q < 2) {
1783b0d62198Sjeremylt     // LCOV_EXCL_START
17842b730f8bSJeremy L Thompson     return CeedError(NULL, CEED_ERROR_DIMENSION, "Cannot create Lobatto quadrature with Q=%" CeedInt_FMT " < 2 points", Q);
1785b0d62198Sjeremylt     // LCOV_EXCL_STOP
17862b730f8bSJeremy L Thompson   }
1787d7b241e6Sjeremylt   wi = 2.0 / ((CeedScalar)(Q * (Q - 1)));
1788d1d35e2fSjeremylt   if (q_weight_1d) {
1789d1d35e2fSjeremylt     q_weight_1d[0]     = wi;
1790d1d35e2fSjeremylt     q_weight_1d[Q - 1] = wi;
1791d7b241e6Sjeremylt   }
1792d1d35e2fSjeremylt   q_ref_1d[0]     = -1.0;
1793d1d35e2fSjeremylt   q_ref_1d[Q - 1] = 1.0;
1794d7b241e6Sjeremylt   // Interior
179592ae7e47SJeremy L Thompson   for (CeedInt i = 1; i <= (Q - 1) / 2; i++) {
1796d7b241e6Sjeremylt     // Guess
1797d7b241e6Sjeremylt     xi = cos(PI * (CeedScalar)(i) / (CeedScalar)(Q - 1));
1798d7b241e6Sjeremylt     // Pn(xi)
1799d7b241e6Sjeremylt     P0 = 1.0;
1800d7b241e6Sjeremylt     P1 = xi;
1801d7b241e6Sjeremylt     P2 = 0.0;
180292ae7e47SJeremy L Thompson     for (CeedInt j = 2; j < Q; j++) {
1803d7b241e6Sjeremylt       P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
1804d7b241e6Sjeremylt       P0 = P1;
1805d7b241e6Sjeremylt       P1 = P2;
1806d7b241e6Sjeremylt     }
1807d7b241e6Sjeremylt     // First Newton step
1808d7b241e6Sjeremylt     dP2  = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
1809d7b241e6Sjeremylt     d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi);
1810d7b241e6Sjeremylt     xi   = xi - dP2 / d2P2;
1811d7b241e6Sjeremylt     // Newton to convergence
181292ae7e47SJeremy L Thompson     for (CeedInt k = 0; k < 100 && fabs(dP2) > 10 * CEED_EPSILON; k++) {
1813d7b241e6Sjeremylt       P0 = 1.0;
1814d7b241e6Sjeremylt       P1 = xi;
181592ae7e47SJeremy L Thompson       for (CeedInt j = 2; j < Q; j++) {
1816d7b241e6Sjeremylt         P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j));
1817d7b241e6Sjeremylt         P0 = P1;
1818d7b241e6Sjeremylt         P1 = P2;
1819d7b241e6Sjeremylt       }
1820d7b241e6Sjeremylt       dP2  = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0);
1821d7b241e6Sjeremylt       d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi);
1822d7b241e6Sjeremylt       xi   = xi - dP2 / d2P2;
1823d7b241e6Sjeremylt     }
1824d7b241e6Sjeremylt     // Save xi, wi
1825d7b241e6Sjeremylt     wi = 2.0 / (((CeedScalar)(Q * (Q - 1))) * P2 * P2);
1826d1d35e2fSjeremylt     if (q_weight_1d) {
1827d1d35e2fSjeremylt       q_weight_1d[i]         = wi;
1828d1d35e2fSjeremylt       q_weight_1d[Q - 1 - i] = wi;
1829d7b241e6Sjeremylt     }
1830d1d35e2fSjeremylt     q_ref_1d[i]         = -xi;
1831d1d35e2fSjeremylt     q_ref_1d[Q - 1 - i] = xi;
1832d7b241e6Sjeremylt   }
1833e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1834d7b241e6Sjeremylt }
1835d7b241e6Sjeremylt 
1836d7b241e6Sjeremylt /// @}
1837