15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors. 23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3d7b241e6Sjeremylt // 43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 5d7b241e6Sjeremylt // 63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 7d7b241e6Sjeremylt 83d576824SJeremy L Thompson #include <ceed-impl.h> 949aac155SJeremy L Thompson #include <ceed.h> 102b730f8bSJeremy L Thompson #include <ceed/backend.h> 11d7b241e6Sjeremylt #include <math.h> 123d576824SJeremy L Thompson #include <stdbool.h> 13d7b241e6Sjeremylt #include <stdio.h> 14d7b241e6Sjeremylt #include <string.h> 15d7b241e6Sjeremylt 167a982d89SJeremy L. Thompson /// @file 177a982d89SJeremy L. Thompson /// Implementation of CeedBasis interfaces 187a982d89SJeremy L. Thompson 19d7b241e6Sjeremylt /// @cond DOXYGEN_SKIP 20356036faSJeremy L Thompson static struct CeedBasis_private ceed_basis_none; 21d7b241e6Sjeremylt /// @endcond 22d7b241e6Sjeremylt 237a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser 247a982d89SJeremy L. Thompson /// @{ 257a982d89SJeremy L. Thompson 26ca94c3ddSJeremy L Thompson /// Argument for @ref CeedOperatorSetField() indicating that the field does not require a `CeedBasis` 27356036faSJeremy L Thompson const CeedBasis CEED_BASIS_NONE = &ceed_basis_none; 28356036faSJeremy 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 /** 383778dbaaSJeremy L Thompson @brief Compute Chebyshev polynomial values at a point 393778dbaaSJeremy L Thompson 403778dbaaSJeremy L Thompson @param[in] x Coordinate to evaluate Chebyshev polynomials at 41ca94c3ddSJeremy L Thompson @param[in] n Number of Chebyshev polynomials to evaluate, `n >= 2` 423778dbaaSJeremy L Thompson @param[out] chebyshev_x Array of Chebyshev polynomial values 433778dbaaSJeremy L Thompson 443778dbaaSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 453778dbaaSJeremy L Thompson 463778dbaaSJeremy L Thompson @ref Developer 473778dbaaSJeremy L Thompson **/ 483778dbaaSJeremy L Thompson static int CeedChebyshevPolynomialsAtPoint(CeedScalar x, CeedInt n, CeedScalar *chebyshev_x) { 493778dbaaSJeremy L Thompson chebyshev_x[0] = 1.0; 503778dbaaSJeremy L Thompson chebyshev_x[1] = 2 * x; 513778dbaaSJeremy L Thompson for (CeedInt i = 2; i < n; i++) chebyshev_x[i] = 2 * x * chebyshev_x[i - 1] - chebyshev_x[i - 2]; 523778dbaaSJeremy L Thompson return CEED_ERROR_SUCCESS; 533778dbaaSJeremy L Thompson } 543778dbaaSJeremy L Thompson 553778dbaaSJeremy L Thompson /** 563778dbaaSJeremy L Thompson @brief Compute values of the derivative of Chebyshev polynomials at a point 573778dbaaSJeremy L Thompson 583778dbaaSJeremy L Thompson @param[in] x Coordinate to evaluate derivative of Chebyshev polynomials at 59ca94c3ddSJeremy L Thompson @param[in] n Number of Chebyshev polynomials to evaluate, `n >= 2` 606cec60aaSJed Brown @param[out] chebyshev_dx Array of Chebyshev polynomial derivative values 613778dbaaSJeremy L Thompson 623778dbaaSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 633778dbaaSJeremy L Thompson 643778dbaaSJeremy L Thompson @ref Developer 653778dbaaSJeremy L Thompson **/ 663778dbaaSJeremy L Thompson static int CeedChebyshevDerivativeAtPoint(CeedScalar x, CeedInt n, CeedScalar *chebyshev_dx) { 673778dbaaSJeremy L Thompson CeedScalar chebyshev_x[3]; 683778dbaaSJeremy L Thompson 693778dbaaSJeremy L Thompson chebyshev_x[1] = 1.0; 703778dbaaSJeremy L Thompson chebyshev_x[2] = 2 * x; 713778dbaaSJeremy L Thompson chebyshev_dx[0] = 0.0; 723778dbaaSJeremy L Thompson chebyshev_dx[1] = 2.0; 733778dbaaSJeremy L Thompson for (CeedInt i = 2; i < n; i++) { 743778dbaaSJeremy L Thompson chebyshev_x[0] = chebyshev_x[1]; 753778dbaaSJeremy L Thompson chebyshev_x[1] = chebyshev_x[2]; 763778dbaaSJeremy L Thompson chebyshev_x[2] = 2 * x * chebyshev_x[1] - chebyshev_x[0]; 773778dbaaSJeremy L Thompson chebyshev_dx[i] = 2 * x * chebyshev_dx[i - 1] + 2 * chebyshev_x[1] - chebyshev_dx[i - 2]; 783778dbaaSJeremy L Thompson } 793778dbaaSJeremy L Thompson return CEED_ERROR_SUCCESS; 803778dbaaSJeremy L Thompson } 813778dbaaSJeremy L Thompson 823778dbaaSJeremy L Thompson /** 83ca94c3ddSJeremy L Thompson @brief Compute Householder reflection. 847a982d89SJeremy L. Thompson 85ca94c3ddSJeremy L Thompson Computes \f$A = (I - b v v^T) A\f$, where \f$A\f$ is an \f$m \times n\f$ matrix indexed as `A[i*row + j*col]`. 867a982d89SJeremy L. Thompson 877a982d89SJeremy L. Thompson @param[in,out] A Matrix to apply Householder reflection to, in place 88ea61e9acSJeremy L Thompson @param[in] v Householder vector 89ea61e9acSJeremy L Thompson @param[in] b Scaling factor 90ca94c3ddSJeremy L Thompson @param[in] m Number of rows in `A` 91ca94c3ddSJeremy L Thompson @param[in] n Number of columns in `A` 92ea61e9acSJeremy L Thompson @param[in] row Row stride 93ea61e9acSJeremy L Thompson @param[in] col Col stride 947a982d89SJeremy L. Thompson 957a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 967a982d89SJeremy L. Thompson 977a982d89SJeremy L. Thompson @ref Developer 987a982d89SJeremy L. Thompson **/ 992b730f8bSJeremy L Thompson static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v, CeedScalar b, CeedInt m, CeedInt n, CeedInt row, CeedInt col) { 1007a982d89SJeremy L. Thompson for (CeedInt j = 0; j < n; j++) { 1017a982d89SJeremy L. Thompson CeedScalar w = A[0 * row + j * col]; 1021c66c397SJeremy L Thompson 1032b730f8bSJeremy L Thompson for (CeedInt i = 1; i < m; i++) w += v[i] * A[i * row + j * col]; 1047a982d89SJeremy L. Thompson A[0 * row + j * col] -= b * w; 1052b730f8bSJeremy L Thompson for (CeedInt i = 1; i < m; i++) A[i * row + j * col] -= b * w * v[i]; 1067a982d89SJeremy L. Thompson } 107e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1087a982d89SJeremy L. Thompson } 1097a982d89SJeremy L. Thompson 1107a982d89SJeremy L. Thompson /** 1117a982d89SJeremy L. Thompson @brief Compute Givens rotation 1127a982d89SJeremy L. Thompson 113ca94c3ddSJeremy L Thompson Computes \f$A = G A\f$ (or \f$G^T A\f$ in transpose mode), where \f$A\f$ is an \f$m \times n\f$ matrix indexed as `A[i*n + j*m]`. 1147a982d89SJeremy L. Thompson 1157a982d89SJeremy L. Thompson @param[in,out] A Row major matrix to apply Givens rotation to, in place 116ea61e9acSJeremy L Thompson @param[in] c Cosine factor 117ea61e9acSJeremy L Thompson @param[in] s Sine factor 118ca94c3ddSJeremy 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; 1194cc79fe7SJed Brown @ref CEED_TRANSPOSE for the opposite rotation 120ea61e9acSJeremy L Thompson @param[in] i First row/column to apply rotation 121ea61e9acSJeremy L Thompson @param[in] k Second row/column to apply rotation 122ca94c3ddSJeremy L Thompson @param[in] m Number of rows in `A` 123ca94c3ddSJeremy L Thompson @param[in] n Number of columns in `A` 1247a982d89SJeremy L. Thompson 1257a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 1267a982d89SJeremy L. Thompson 1277a982d89SJeremy L. Thompson @ref Developer 1287a982d89SJeremy L. Thompson **/ 1292b730f8bSJeremy L Thompson static int CeedGivensRotation(CeedScalar *A, CeedScalar c, CeedScalar s, CeedTransposeMode t_mode, CeedInt i, CeedInt k, CeedInt m, CeedInt n) { 130d1d35e2fSjeremylt CeedInt stride_j = 1, stride_ik = m, num_its = n; 1311c66c397SJeremy L Thompson 132d1d35e2fSjeremylt if (t_mode == CEED_NOTRANSPOSE) { 1332b730f8bSJeremy L Thompson stride_j = n; 1342b730f8bSJeremy L Thompson stride_ik = 1; 1352b730f8bSJeremy L Thompson num_its = m; 1367a982d89SJeremy L. Thompson } 1377a982d89SJeremy L. Thompson 1387a982d89SJeremy L. Thompson // Apply rotation 139d1d35e2fSjeremylt for (CeedInt j = 0; j < num_its; j++) { 140d1d35e2fSjeremylt CeedScalar tau1 = A[i * stride_ik + j * stride_j], tau2 = A[k * stride_ik + j * stride_j]; 1411c66c397SJeremy L Thompson 142d1d35e2fSjeremylt A[i * stride_ik + j * stride_j] = c * tau1 - s * tau2; 143d1d35e2fSjeremylt A[k * stride_ik + j * stride_j] = s * tau1 + c * tau2; 1447a982d89SJeremy L. Thompson } 145e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1467a982d89SJeremy L. Thompson } 1477a982d89SJeremy L. Thompson 1487a982d89SJeremy L. Thompson /** 149ca94c3ddSJeremy L Thompson @brief View an array stored in a `CeedBasis` 1507a982d89SJeremy L. Thompson 1510a0da059Sjeremylt @param[in] name Name of array 152d1d35e2fSjeremylt @param[in] fp_fmt Printing format 1530a0da059Sjeremylt @param[in] m Number of rows in array 1540a0da059Sjeremylt @param[in] n Number of columns in array 1550a0da059Sjeremylt @param[in] a Array to be viewed 156ca94c3ddSJeremy L Thompson @param[in] stream Stream to view to, e.g., `stdout` 1577a982d89SJeremy L. Thompson 1587a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 1597a982d89SJeremy L. Thompson 1607a982d89SJeremy L. Thompson @ref Developer 1617a982d89SJeremy L. Thompson **/ 1622b730f8bSJeremy L Thompson static int CeedScalarView(const char *name, const char *fp_fmt, CeedInt m, CeedInt n, const CeedScalar *a, FILE *stream) { 163edf04919SJeremy L Thompson if (m > 1) { 164edf04919SJeremy L Thompson fprintf(stream, " %s:\n", name); 165edf04919SJeremy L Thompson } else { 166edf04919SJeremy L Thompson char padded_name[12]; 167edf04919SJeremy L Thompson 168edf04919SJeremy L Thompson snprintf(padded_name, 11, "%s:", name); 169edf04919SJeremy L Thompson fprintf(stream, " %-10s", padded_name); 170edf04919SJeremy L Thompson } 17192ae7e47SJeremy L Thompson for (CeedInt i = 0; i < m; i++) { 172edf04919SJeremy L Thompson if (m > 1) fprintf(stream, " [%" CeedInt_FMT "]", i); 1732b730f8bSJeremy 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); 1747a982d89SJeremy L. Thompson fputs("\n", stream); 1757a982d89SJeremy L. Thompson } 176e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1777a982d89SJeremy L. Thompson } 1787a982d89SJeremy L. Thompson 179a76a04e7SJeremy L Thompson /** 180ea61e9acSJeremy L Thompson @brief Create the interpolation and gradient matrices for projection from the nodes of `basis_from` to the nodes of `basis_to`. 181ba59ac12SSebastian Grimberg 18215ad3917SSebastian Grimberg The interpolation is given by `interp_project = interp_to^+ * interp_from`, where the pseudoinverse `interp_to^+` is given by QR factorization. 183ca94c3ddSJeremy L Thompson The gradient is given by `grad_project = interp_to^+ * grad_from`, and is only computed for \f$H^1\f$ spaces otherwise it should not be used. 18415ad3917SSebastian Grimberg 185ba59ac12SSebastian Grimberg Note: `basis_from` and `basis_to` must have compatible quadrature spaces. 186a76a04e7SJeremy L Thompson 187ca94c3ddSJeremy L Thompson @param[in] basis_from `CeedBasis` to project from 188ca94c3ddSJeremy L Thompson @param[in] basis_to `CeedBasis` to project to 189ca94c3ddSJeremy L Thompson @param[out] interp_project Address of the variable where the newly created interpolation matrix will be stored 190ca94c3ddSJeremy L Thompson @param[out] grad_project Address of the variable where the newly created gradient matrix will be stored 191a76a04e7SJeremy L Thompson 192a76a04e7SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 193a76a04e7SJeremy L Thompson 194a76a04e7SJeremy L Thompson @ref Developer 195a76a04e7SJeremy L Thompson **/ 1962b730f8bSJeremy L Thompson static int CeedBasisCreateProjectionMatrices(CeedBasis basis_from, CeedBasis basis_to, CeedScalar **interp_project, CeedScalar **grad_project) { 197a76a04e7SJeremy L Thompson Ceed ceed; 198e104ad11SJames Wright bool are_both_tensor; 1991c66c397SJeremy L Thompson CeedInt Q, Q_to, Q_from, P_to, P_from; 2001c66c397SJeremy L Thompson 2012b730f8bSJeremy L Thompson CeedCall(CeedBasisGetCeed(basis_to, &ceed)); 202a76a04e7SJeremy L Thompson 203a76a04e7SJeremy L Thompson // Check for compatible quadrature spaces 2042b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints(basis_to, &Q_to)); 2052b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints(basis_from, &Q_from)); 2066574a04fSJeremy L Thompson CeedCheck(Q_to == Q_from, ceed, CEED_ERROR_DIMENSION, "Bases must have compatible quadrature spaces"); 2071c66c397SJeremy L Thompson Q = Q_to; 208a76a04e7SJeremy L Thompson 20914556e63SJeremy L Thompson // Check for matching tensor or non-tensor 210e104ad11SJames Wright { 211e104ad11SJames Wright bool is_tensor_to, is_tensor_from; 212e104ad11SJames Wright 2132b730f8bSJeremy L Thompson CeedCall(CeedBasisIsTensor(basis_to, &is_tensor_to)); 2142b730f8bSJeremy L Thompson CeedCall(CeedBasisIsTensor(basis_from, &is_tensor_from)); 215e104ad11SJames Wright are_both_tensor = is_tensor_to && is_tensor_from; 216e104ad11SJames Wright } 217e104ad11SJames Wright if (are_both_tensor) { 2182b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_to)); 2192b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_from)); 2202b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints1D(basis_from, &Q)); 2216574a04fSJeremy L Thompson } else { 2222b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis_to, &P_to)); 2232b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis_from, &P_from)); 224a76a04e7SJeremy L Thompson } 225a76a04e7SJeremy L Thompson 22615ad3917SSebastian Grimberg // Check for matching FE space 22715ad3917SSebastian Grimberg CeedFESpace fe_space_to, fe_space_from; 22815ad3917SSebastian Grimberg CeedCall(CeedBasisGetFESpace(basis_to, &fe_space_to)); 22915ad3917SSebastian Grimberg CeedCall(CeedBasisGetFESpace(basis_from, &fe_space_from)); 2306574a04fSJeremy L Thompson CeedCheck(fe_space_to == fe_space_from, ceed, CEED_ERROR_MINOR, "Bases must both be the same FE space type"); 23115ad3917SSebastian Grimberg 23214556e63SJeremy L Thompson // Get source matrices 23315ad3917SSebastian Grimberg CeedInt dim, q_comp = 1; 2342247a93fSRezgar Shakeri CeedScalar *interp_to_inv, *interp_from; 2351c66c397SJeremy L Thompson const CeedScalar *interp_to_source = NULL, *interp_from_source = NULL, *grad_from_source = NULL; 2361c66c397SJeremy L Thompson 2372b730f8bSJeremy L Thompson CeedCall(CeedBasisGetDimension(basis_to, &dim)); 238e104ad11SJames Wright if (are_both_tensor) { 2392b730f8bSJeremy L Thompson CeedCall(CeedBasisGetInterp1D(basis_to, &interp_to_source)); 2402b730f8bSJeremy L Thompson CeedCall(CeedBasisGetInterp1D(basis_from, &interp_from_source)); 241a76a04e7SJeremy L Thompson } else { 24215ad3917SSebastian Grimberg CeedCall(CeedBasisGetNumQuadratureComponents(basis_from, CEED_EVAL_INTERP, &q_comp)); 2432b730f8bSJeremy L Thompson CeedCall(CeedBasisGetInterp(basis_to, &interp_to_source)); 2442b730f8bSJeremy L Thompson CeedCall(CeedBasisGetInterp(basis_from, &interp_from_source)); 24515ad3917SSebastian Grimberg } 24615ad3917SSebastian Grimberg CeedCall(CeedMalloc(Q * P_from * q_comp, &interp_from)); 24715ad3917SSebastian Grimberg CeedCall(CeedCalloc(P_to * P_from, interp_project)); 24815ad3917SSebastian Grimberg 24915ad3917SSebastian Grimberg // `grad_project = interp_to^+ * grad_from` is computed for the H^1 space case so the 250de05fbb2SSebastian Grimberg // projection basis will have a gradient operation (allocated even if not H^1 for the 251de05fbb2SSebastian Grimberg // basis construction later on) 25215ad3917SSebastian Grimberg if (fe_space_to == CEED_FE_SPACE_H1) { 253e104ad11SJames Wright if (are_both_tensor) { 25415ad3917SSebastian Grimberg CeedCall(CeedBasisGetGrad1D(basis_from, &grad_from_source)); 25515ad3917SSebastian Grimberg } else { 2562b730f8bSJeremy L Thompson CeedCall(CeedBasisGetGrad(basis_from, &grad_from_source)); 257a76a04e7SJeremy L Thompson } 258de05fbb2SSebastian Grimberg } 259e104ad11SJames Wright CeedCall(CeedCalloc(P_to * P_from * (are_both_tensor ? 1 : dim), grad_project)); 26015ad3917SSebastian Grimberg 2612247a93fSRezgar Shakeri // Compute interp_to^+, pseudoinverse of interp_to 2622247a93fSRezgar Shakeri CeedCall(CeedCalloc(Q * q_comp * P_to, &interp_to_inv)); 2631203703bSJeremy L Thompson CeedCall(CeedMatrixPseudoinverse(ceed, interp_to_source, Q * q_comp, P_to, interp_to_inv)); 26414556e63SJeremy L Thompson // Build matrices 265e104ad11SJames Wright CeedInt num_matrices = 1 + (fe_space_to == CEED_FE_SPACE_H1) * (are_both_tensor ? 1 : dim); 26614556e63SJeremy L Thompson CeedScalar *input_from[num_matrices], *output_project[num_matrices]; 2671c66c397SJeremy L Thompson 26814556e63SJeremy L Thompson input_from[0] = (CeedScalar *)interp_from_source; 26914556e63SJeremy L Thompson output_project[0] = *interp_project; 27014556e63SJeremy L Thompson for (CeedInt m = 1; m < num_matrices; m++) { 27114556e63SJeremy L Thompson input_from[m] = (CeedScalar *)&grad_from_source[(m - 1) * Q * P_from]; 27202af4036SJeremy L Thompson output_project[m] = &((*grad_project)[(m - 1) * P_to * P_from]); 27314556e63SJeremy L Thompson } 27414556e63SJeremy L Thompson for (CeedInt m = 0; m < num_matrices; m++) { 2752247a93fSRezgar Shakeri // output_project = interp_to^+ * interp_from 27615ad3917SSebastian Grimberg memcpy(interp_from, input_from[m], Q * P_from * q_comp * sizeof(input_from[m][0])); 2772247a93fSRezgar Shakeri CeedCall(CeedMatrixMatrixMultiply(ceed, interp_to_inv, input_from[m], output_project[m], P_to, P_from, Q * q_comp)); 2782247a93fSRezgar Shakeri // Round zero to machine precision 2792247a93fSRezgar Shakeri for (CeedInt i = 0; i < P_to * P_from; i++) { 2802247a93fSRezgar Shakeri if (fabs(output_project[m][i]) < 10 * CEED_EPSILON) output_project[m][i] = 0.0; 281a76a04e7SJeremy L Thompson } 28214556e63SJeremy L Thompson } 28314556e63SJeremy L Thompson 28414556e63SJeremy L Thompson // Cleanup 2852247a93fSRezgar Shakeri CeedCall(CeedFree(&interp_to_inv)); 2862b730f8bSJeremy L Thompson CeedCall(CeedFree(&interp_from)); 287a76a04e7SJeremy L Thompson return CEED_ERROR_SUCCESS; 288a76a04e7SJeremy L Thompson } 289a76a04e7SJeremy L Thompson 2907a982d89SJeremy L. Thompson /// @} 2917a982d89SJeremy L. Thompson 2927a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 2937a982d89SJeremy L. Thompson /// Ceed Backend API 2947a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 2957a982d89SJeremy L. Thompson /// @addtogroup CeedBasisBackend 2967a982d89SJeremy L. Thompson /// @{ 2977a982d89SJeremy L. Thompson 2987a982d89SJeremy L. Thompson /** 299ca94c3ddSJeremy L Thompson @brief Return collocated gradient matrix 3007a982d89SJeremy L. Thompson 301ca94c3ddSJeremy L Thompson @param[in] basis `CeedBasis` 302ca94c3ddSJeremy L Thompson @param[out] collo_grad_1d Row-major (`Q_1d * Q_1d`) matrix expressing derivatives of basis functions at quadrature points 3037a982d89SJeremy L. Thompson 3047a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3057a982d89SJeremy L. Thompson 3067a982d89SJeremy L. Thompson @ref Backend 3077a982d89SJeremy L. Thompson **/ 308d1d35e2fSjeremylt int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *collo_grad_1d) { 3097a982d89SJeremy L. Thompson Ceed ceed; 3102247a93fSRezgar Shakeri CeedInt P_1d, Q_1d; 3112247a93fSRezgar Shakeri CeedScalar *interp_1d_pinv; 3121203703bSJeremy L Thompson const CeedScalar *grad_1d, *interp_1d; 3131203703bSJeremy L Thompson 314ea61e9acSJeremy 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. 3152247a93fSRezgar Shakeri CeedCall(CeedBasisGetCeed(basis, &ceed)); 3162247a93fSRezgar Shakeri CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d)); 3172247a93fSRezgar Shakeri CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 3187a982d89SJeremy L. Thompson 3192247a93fSRezgar Shakeri // Compute interp_1d^+, pseudoinverse of interp_1d 3202247a93fSRezgar Shakeri CeedCall(CeedCalloc(P_1d * Q_1d, &interp_1d_pinv)); 3211203703bSJeremy L Thompson CeedCall(CeedBasisGetInterp1D(basis, &interp_1d)); 3221203703bSJeremy L Thompson CeedCall(CeedMatrixPseudoinverse(ceed, interp_1d, Q_1d, P_1d, interp_1d_pinv)); 3231203703bSJeremy L Thompson CeedCall(CeedBasisGetGrad1D(basis, &grad_1d)); 3241203703bSJeremy L Thompson CeedCall(CeedMatrixMatrixMultiply(ceed, grad_1d, (const CeedScalar *)interp_1d_pinv, collo_grad_1d, Q_1d, Q_1d, P_1d)); 3257a982d89SJeremy L. Thompson 3262247a93fSRezgar Shakeri CeedCall(CeedFree(&interp_1d_pinv)); 327e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3287a982d89SJeremy L. Thompson } 3297a982d89SJeremy L. Thompson 3307a982d89SJeremy L. Thompson /** 331b0cc4569SJeremy L Thompson @brief Return 1D interpolation matrix to Chebyshev polynomial coefficients on quadrature space 332b0cc4569SJeremy L Thompson 333b0cc4569SJeremy L Thompson @param[in] basis `CeedBasis` 334b0cc4569SJeremy L Thompson @param[out] chebyshev_interp_1d Row-major (`P_1d * Q_1d`) matrix interpolating from basis nodes to Chebyshev polynomial coefficients 335b0cc4569SJeremy L Thompson 336b0cc4569SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 337b0cc4569SJeremy L Thompson 338b0cc4569SJeremy L Thompson @ref Backend 339b0cc4569SJeremy L Thompson **/ 340b0cc4569SJeremy L Thompson int CeedBasisGetChebyshevInterp1D(CeedBasis basis, CeedScalar *chebyshev_interp_1d) { 341b0cc4569SJeremy L Thompson CeedInt P_1d, Q_1d; 342b0cc4569SJeremy L Thompson CeedScalar *C, *chebyshev_coeffs_1d_inv; 343b0cc4569SJeremy L Thompson const CeedScalar *interp_1d, *q_ref_1d; 344b0cc4569SJeremy L Thompson Ceed ceed; 345b0cc4569SJeremy L Thompson 346b0cc4569SJeremy L Thompson CeedCall(CeedBasisGetCeed(basis, &ceed)); 347b0cc4569SJeremy L Thompson CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d)); 348b0cc4569SJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 349b0cc4569SJeremy L Thompson 350b0cc4569SJeremy L Thompson // Build coefficient matrix 351bd83cbc5SJeremy L Thompson // -- Note: Clang-tidy needs this check 352bd83cbc5SJeremy L Thompson CeedCheck((P_1d > 0) && (Q_1d > 0), ceed, CEED_ERROR_INCOMPATIBLE, "CeedBasis dimensions are malformed"); 353b0cc4569SJeremy L Thompson CeedCall(CeedCalloc(Q_1d * Q_1d, &C)); 354b0cc4569SJeremy L Thompson CeedCall(CeedBasisGetQRef(basis, &q_ref_1d)); 355b0cc4569SJeremy L Thompson for (CeedInt i = 0; i < Q_1d; i++) CeedCall(CeedChebyshevPolynomialsAtPoint(q_ref_1d[i], Q_1d, &C[i * Q_1d])); 356b0cc4569SJeremy L Thompson 357b0cc4569SJeremy L Thompson // Compute C^+, pseudoinverse of coefficient matrix 358b0cc4569SJeremy L Thompson CeedCall(CeedCalloc(Q_1d * Q_1d, &chebyshev_coeffs_1d_inv)); 359b0cc4569SJeremy L Thompson CeedCall(CeedMatrixPseudoinverse(ceed, C, Q_1d, Q_1d, chebyshev_coeffs_1d_inv)); 360b0cc4569SJeremy L Thompson 361b0cc4569SJeremy L Thompson // Build mapping from nodes to Chebyshev coefficients 362b0cc4569SJeremy L Thompson CeedCall(CeedBasisGetInterp1D(basis, &interp_1d)); 363b0cc4569SJeremy L Thompson CeedCall(CeedMatrixMatrixMultiply(ceed, chebyshev_coeffs_1d_inv, interp_1d, chebyshev_interp_1d, Q_1d, P_1d, Q_1d)); 364b0cc4569SJeremy L Thompson 365b0cc4569SJeremy L Thompson // Cleanup 366b0cc4569SJeremy L Thompson CeedCall(CeedFree(&C)); 367b0cc4569SJeremy L Thompson CeedCall(CeedFree(&chebyshev_coeffs_1d_inv)); 368b0cc4569SJeremy L Thompson return CEED_ERROR_SUCCESS; 369b0cc4569SJeremy L Thompson } 370b0cc4569SJeremy L Thompson 371b0cc4569SJeremy L Thompson /** 372ca94c3ddSJeremy L Thompson @brief Get tensor status for given `CeedBasis` 3737a982d89SJeremy L. Thompson 374ca94c3ddSJeremy L Thompson @param[in] basis `CeedBasis` 375d1d35e2fSjeremylt @param[out] is_tensor Variable to store tensor status 3767a982d89SJeremy L. Thompson 3777a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3787a982d89SJeremy L. Thompson 3797a982d89SJeremy L. Thompson @ref Backend 3807a982d89SJeremy L. Thompson **/ 381d1d35e2fSjeremylt int CeedBasisIsTensor(CeedBasis basis, bool *is_tensor) { 3826402da51SJeremy L Thompson *is_tensor = basis->is_tensor_basis; 383e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3847a982d89SJeremy L. Thompson } 3857a982d89SJeremy L. Thompson 3867a982d89SJeremy L. Thompson /** 387ca94c3ddSJeremy L Thompson @brief Get backend data of a `CeedBasis` 3887a982d89SJeremy L. Thompson 389ca94c3ddSJeremy L Thompson @param[in] basis `CeedBasis` 3907a982d89SJeremy L. Thompson @param[out] data Variable to store data 3917a982d89SJeremy L. Thompson 3927a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3937a982d89SJeremy L. Thompson 3947a982d89SJeremy L. Thompson @ref Backend 3957a982d89SJeremy L. Thompson **/ 396777ff853SJeremy L Thompson int CeedBasisGetData(CeedBasis basis, void *data) { 397777ff853SJeremy L Thompson *(void **)data = basis->data; 398e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3997a982d89SJeremy L. Thompson } 4007a982d89SJeremy L. Thompson 4017a982d89SJeremy L. Thompson /** 402ca94c3ddSJeremy L Thompson @brief Set backend data of a `CeedBasis` 4037a982d89SJeremy L. Thompson 404ca94c3ddSJeremy L Thompson @param[in,out] basis `CeedBasis` 405ea61e9acSJeremy L Thompson @param[in] data Data to set 4067a982d89SJeremy L. Thompson 4077a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 4087a982d89SJeremy L. Thompson 4097a982d89SJeremy L. Thompson @ref Backend 4107a982d89SJeremy L. Thompson **/ 411777ff853SJeremy L Thompson int CeedBasisSetData(CeedBasis basis, void *data) { 412777ff853SJeremy L Thompson basis->data = data; 413e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 4147a982d89SJeremy L. Thompson } 4157a982d89SJeremy L. Thompson 4167a982d89SJeremy L. Thompson /** 417ca94c3ddSJeremy L Thompson @brief Increment the reference counter for a `CeedBasis` 41834359f16Sjeremylt 419ca94c3ddSJeremy L Thompson @param[in,out] basis `CeedBasis` to increment the reference counter 42034359f16Sjeremylt 42134359f16Sjeremylt @return An error code: 0 - success, otherwise - failure 42234359f16Sjeremylt 42334359f16Sjeremylt @ref Backend 42434359f16Sjeremylt **/ 4259560d06aSjeremylt int CeedBasisReference(CeedBasis basis) { 42634359f16Sjeremylt basis->ref_count++; 42734359f16Sjeremylt return CEED_ERROR_SUCCESS; 42834359f16Sjeremylt } 42934359f16Sjeremylt 43034359f16Sjeremylt /** 431ca94c3ddSJeremy L Thompson @brief Get number of Q-vector components for given `CeedBasis` 432c4e3f59bSSebastian Grimberg 433ca94c3ddSJeremy L Thompson @param[in] basis `CeedBasis` 434ca94c3ddSJeremy L Thompson @param[in] eval_mode @ref CEED_EVAL_INTERP to use interpolated values, 435ca94c3ddSJeremy L Thompson @ref CEED_EVAL_GRAD to use gradients, 436ca94c3ddSJeremy L Thompson @ref CEED_EVAL_DIV to use divergence, 437ca94c3ddSJeremy L Thompson @ref CEED_EVAL_CURL to use curl 438c4e3f59bSSebastian Grimberg @param[out] q_comp Variable to store number of Q-vector components of basis 439c4e3f59bSSebastian Grimberg 440c4e3f59bSSebastian Grimberg @return An error code: 0 - success, otherwise - failure 441c4e3f59bSSebastian Grimberg 442c4e3f59bSSebastian Grimberg @ref Backend 443c4e3f59bSSebastian Grimberg **/ 444c4e3f59bSSebastian Grimberg int CeedBasisGetNumQuadratureComponents(CeedBasis basis, CeedEvalMode eval_mode, CeedInt *q_comp) { 4451203703bSJeremy L Thompson CeedInt dim; 4461203703bSJeremy L Thompson 4471203703bSJeremy L Thompson CeedCall(CeedBasisGetDimension(basis, &dim)); 448c4e3f59bSSebastian Grimberg switch (eval_mode) { 4491203703bSJeremy L Thompson case CEED_EVAL_INTERP: { 4501203703bSJeremy L Thompson CeedFESpace fe_space; 4511203703bSJeremy L Thompson 4521203703bSJeremy L Thompson CeedCall(CeedBasisGetFESpace(basis, &fe_space)); 4531203703bSJeremy L Thompson *q_comp = (fe_space == CEED_FE_SPACE_H1) ? 1 : dim; 4541203703bSJeremy L Thompson } break; 455c4e3f59bSSebastian Grimberg case CEED_EVAL_GRAD: 4561203703bSJeremy L Thompson *q_comp = dim; 457c4e3f59bSSebastian Grimberg break; 458c4e3f59bSSebastian Grimberg case CEED_EVAL_DIV: 459c4e3f59bSSebastian Grimberg *q_comp = 1; 460c4e3f59bSSebastian Grimberg break; 461c4e3f59bSSebastian Grimberg case CEED_EVAL_CURL: 4621203703bSJeremy L Thompson *q_comp = (dim < 3) ? 1 : dim; 463c4e3f59bSSebastian Grimberg break; 464c4e3f59bSSebastian Grimberg case CEED_EVAL_NONE: 465c4e3f59bSSebastian Grimberg case CEED_EVAL_WEIGHT: 466352a5e7cSSebastian Grimberg *q_comp = 1; 467c4e3f59bSSebastian Grimberg break; 468c4e3f59bSSebastian Grimberg } 469c4e3f59bSSebastian Grimberg return CEED_ERROR_SUCCESS; 470c4e3f59bSSebastian Grimberg } 471c4e3f59bSSebastian Grimberg 472c4e3f59bSSebastian Grimberg /** 473ca94c3ddSJeremy L Thompson @brief Estimate number of FLOPs required to apply `CeedBasis` in `t_mode` and `eval_mode` 4746e15d496SJeremy L Thompson 475ca94c3ddSJeremy L Thompson @param[in] basis `CeedBasis` to estimate FLOPs for 476ea61e9acSJeremy L Thompson @param[in] t_mode Apply basis or transpose 477ca94c3ddSJeremy L Thompson @param[in] eval_mode @ref CeedEvalMode 478ea61e9acSJeremy L Thompson @param[out] flops Address of variable to hold FLOPs estimate 4796e15d496SJeremy L Thompson 4806e15d496SJeremy L Thompson @ref Backend 4816e15d496SJeremy L Thompson **/ 4822b730f8bSJeremy L Thompson int CeedBasisGetFlopsEstimate(CeedBasis basis, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedSize *flops) { 4836e15d496SJeremy L Thompson bool is_tensor; 4846e15d496SJeremy L Thompson 4852b730f8bSJeremy L Thompson CeedCall(CeedBasisIsTensor(basis, &is_tensor)); 4866e15d496SJeremy L Thompson if (is_tensor) { 4876e15d496SJeremy L Thompson CeedInt dim, num_comp, P_1d, Q_1d; 4881c66c397SJeremy L Thompson 4892b730f8bSJeremy L Thompson CeedCall(CeedBasisGetDimension(basis, &dim)); 4902b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumComponents(basis, &num_comp)); 4912b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d)); 4922b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 4936e15d496SJeremy L Thompson if (t_mode == CEED_TRANSPOSE) { 4942b730f8bSJeremy L Thompson P_1d = Q_1d; 4952b730f8bSJeremy L Thompson Q_1d = P_1d; 4966e15d496SJeremy L Thompson } 4976e15d496SJeremy L Thompson CeedInt tensor_flops = 0, pre = num_comp * CeedIntPow(P_1d, dim - 1), post = 1; 4986e15d496SJeremy L Thompson for (CeedInt d = 0; d < dim; d++) { 4996e15d496SJeremy L Thompson tensor_flops += 2 * pre * P_1d * post * Q_1d; 5006e15d496SJeremy L Thompson pre /= P_1d; 5016e15d496SJeremy L Thompson post *= Q_1d; 5026e15d496SJeremy L Thompson } 5036e15d496SJeremy L Thompson switch (eval_mode) { 5042b730f8bSJeremy L Thompson case CEED_EVAL_NONE: 5052b730f8bSJeremy L Thompson *flops = 0; 5062b730f8bSJeremy L Thompson break; 5072b730f8bSJeremy L Thompson case CEED_EVAL_INTERP: 5082b730f8bSJeremy L Thompson *flops = tensor_flops; 5092b730f8bSJeremy L Thompson break; 5102b730f8bSJeremy L Thompson case CEED_EVAL_GRAD: 5112b730f8bSJeremy L Thompson *flops = tensor_flops * 2; 5122b730f8bSJeremy L Thompson break; 5136e15d496SJeremy L Thompson case CEED_EVAL_DIV: 5141203703bSJeremy L Thompson case CEED_EVAL_CURL: { 5156574a04fSJeremy L Thompson // LCOV_EXCL_START 5166e536b99SJeremy L Thompson return CeedError(CeedBasisReturnCeed(basis), CEED_ERROR_INCOMPATIBLE, "Tensor basis evaluation for %s not supported", 5176e536b99SJeremy L Thompson CeedEvalModes[eval_mode]); 5182b730f8bSJeremy L Thompson break; 5196e15d496SJeremy L Thompson // LCOV_EXCL_STOP 5201203703bSJeremy L Thompson } 5212b730f8bSJeremy L Thompson case CEED_EVAL_WEIGHT: 5222b730f8bSJeremy L Thompson *flops = dim * CeedIntPow(Q_1d, dim); 5232b730f8bSJeremy L Thompson break; 5246e15d496SJeremy L Thompson } 5256e15d496SJeremy L Thompson } else { 526c4e3f59bSSebastian Grimberg CeedInt dim, num_comp, q_comp, num_nodes, num_qpts; 5271c66c397SJeremy L Thompson 5282b730f8bSJeremy L Thompson CeedCall(CeedBasisGetDimension(basis, &dim)); 5292b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumComponents(basis, &num_comp)); 530c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp)); 5312b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis, &num_nodes)); 5322b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts)); 5336e15d496SJeremy L Thompson switch (eval_mode) { 5342b730f8bSJeremy L Thompson case CEED_EVAL_NONE: 5352b730f8bSJeremy L Thompson *flops = 0; 5362b730f8bSJeremy L Thompson break; 5372b730f8bSJeremy L Thompson case CEED_EVAL_INTERP: 5382b730f8bSJeremy L Thompson case CEED_EVAL_GRAD: 5392b730f8bSJeremy L Thompson case CEED_EVAL_DIV: 5402b730f8bSJeremy L Thompson case CEED_EVAL_CURL: 541c4e3f59bSSebastian Grimberg *flops = num_nodes * num_qpts * num_comp * q_comp; 5422b730f8bSJeremy L Thompson break; 5432b730f8bSJeremy L Thompson case CEED_EVAL_WEIGHT: 5442b730f8bSJeremy L Thompson *flops = 0; 5452b730f8bSJeremy L Thompson break; 5466e15d496SJeremy L Thompson } 5476e15d496SJeremy L Thompson } 5486e15d496SJeremy L Thompson return CEED_ERROR_SUCCESS; 5496e15d496SJeremy L Thompson } 5506e15d496SJeremy L Thompson 5516e15d496SJeremy L Thompson /** 552ca94c3ddSJeremy L Thompson @brief Get `CeedFESpace` for a `CeedBasis` 553c4e3f59bSSebastian Grimberg 554ca94c3ddSJeremy L Thompson @param[in] basis `CeedBasis` 555ca94c3ddSJeremy L Thompson @param[out] fe_space Variable to store `CeedFESpace` 556c4e3f59bSSebastian Grimberg 557c4e3f59bSSebastian Grimberg @return An error code: 0 - success, otherwise - failure 558c4e3f59bSSebastian Grimberg 559c4e3f59bSSebastian Grimberg @ref Backend 560c4e3f59bSSebastian Grimberg **/ 561c4e3f59bSSebastian Grimberg int CeedBasisGetFESpace(CeedBasis basis, CeedFESpace *fe_space) { 562c4e3f59bSSebastian Grimberg *fe_space = basis->fe_space; 563c4e3f59bSSebastian Grimberg return CEED_ERROR_SUCCESS; 564c4e3f59bSSebastian Grimberg } 565c4e3f59bSSebastian Grimberg 566c4e3f59bSSebastian Grimberg /** 567ca94c3ddSJeremy L Thompson @brief Get dimension for given `CeedElemTopology` 5687a982d89SJeremy L. Thompson 569ca94c3ddSJeremy L Thompson @param[in] topo `CeedElemTopology` 5707a982d89SJeremy L. Thompson @param[out] dim Variable to store dimension of topology 5717a982d89SJeremy L. Thompson 5727a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 5737a982d89SJeremy L. Thompson 5747a982d89SJeremy L. Thompson @ref Backend 5757a982d89SJeremy L. Thompson **/ 5767a982d89SJeremy L. Thompson int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) { 5777a982d89SJeremy L. Thompson *dim = (CeedInt)topo >> 16; 578e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 5797a982d89SJeremy L. Thompson } 5807a982d89SJeremy L. Thompson 5817a982d89SJeremy L. Thompson /** 582ca94c3ddSJeremy L Thompson @brief Get `CeedTensorContract` of a `CeedBasis` 5837a982d89SJeremy L. Thompson 584ca94c3ddSJeremy L Thompson @param[in] basis `CeedBasis` 585ca94c3ddSJeremy L Thompson @param[out] contract Variable to store `CeedTensorContract` 5867a982d89SJeremy L. Thompson 5877a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 5887a982d89SJeremy L. Thompson 5897a982d89SJeremy L. Thompson @ref Backend 5907a982d89SJeremy L. Thompson **/ 5917a982d89SJeremy L. Thompson int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) { 5927a982d89SJeremy L. Thompson *contract = basis->contract; 593e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 5947a982d89SJeremy L. Thompson } 5957a982d89SJeremy L. Thompson 5967a982d89SJeremy L. Thompson /** 597ca94c3ddSJeremy L Thompson @brief Set `CeedTensorContract` of a `CeedBasis` 5987a982d89SJeremy L. Thompson 599ca94c3ddSJeremy L Thompson @param[in,out] basis `CeedBasis` 600ca94c3ddSJeremy L Thompson @param[in] contract `CeedTensorContract` to set 6017a982d89SJeremy L. Thompson 6027a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 6037a982d89SJeremy L. Thompson 6047a982d89SJeremy L. Thompson @ref Backend 6057a982d89SJeremy L. Thompson **/ 60634359f16Sjeremylt int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract contract) { 60734359f16Sjeremylt basis->contract = contract; 6082b730f8bSJeremy L Thompson CeedCall(CeedTensorContractReference(contract)); 609e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 6107a982d89SJeremy L. Thompson } 6117a982d89SJeremy L. Thompson 6127a982d89SJeremy L. Thompson /** 613ca94c3ddSJeremy L Thompson @brief Return a reference implementation of matrix multiplication \f$C = A B\f$. 614ba59ac12SSebastian Grimberg 615ca94c3ddSJeremy L Thompson Note: This is a reference implementation for CPU `CeedScalar` pointers that is not intended for high performance. 6167a982d89SJeremy L. Thompson 617ca94c3ddSJeremy L Thompson @param[in] ceed `Ceed` context for error handling 618ca94c3ddSJeremy L Thompson @param[in] mat_A Row-major matrix `A` 619ca94c3ddSJeremy L Thompson @param[in] mat_B Row-major matrix `B` 620ca94c3ddSJeremy L Thompson @param[out] mat_C Row-major output matrix `C` 621ca94c3ddSJeremy L Thompson @param[in] m Number of rows of `C` 622ca94c3ddSJeremy L Thompson @param[in] n Number of columns of `C` 623ca94c3ddSJeremy L Thompson @param[in] kk Number of columns of `A`/rows of `B` 6247a982d89SJeremy L. Thompson 6257a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 6267a982d89SJeremy L. Thompson 6277a982d89SJeremy L. Thompson @ref Utility 6287a982d89SJeremy L. Thompson **/ 6292b730f8bSJeremy L Thompson int CeedMatrixMatrixMultiply(Ceed ceed, const CeedScalar *mat_A, const CeedScalar *mat_B, CeedScalar *mat_C, CeedInt m, CeedInt n, CeedInt kk) { 6302b730f8bSJeremy L Thompson for (CeedInt i = 0; i < m; i++) { 6317a982d89SJeremy L. Thompson for (CeedInt j = 0; j < n; j++) { 6327a982d89SJeremy L. Thompson CeedScalar sum = 0; 6331c66c397SJeremy L Thompson 6342b730f8bSJeremy L Thompson for (CeedInt k = 0; k < kk; k++) sum += mat_A[k + i * kk] * mat_B[j + k * n]; 635d1d35e2fSjeremylt mat_C[j + i * n] = sum; 6367a982d89SJeremy L. Thompson } 6372b730f8bSJeremy L Thompson } 638e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 6397a982d89SJeremy L. Thompson } 6407a982d89SJeremy L. Thompson 641ba59ac12SSebastian Grimberg /** 642ba59ac12SSebastian Grimberg @brief Return QR Factorization of a matrix 643ba59ac12SSebastian Grimberg 644ca94c3ddSJeremy L Thompson @param[in] ceed `Ceed` context for error handling 645ba59ac12SSebastian Grimberg @param[in,out] mat Row-major matrix to be factorized in place 646ca94c3ddSJeremy L Thompson @param[in,out] tau Vector of length `m` of scaling factors 647ba59ac12SSebastian Grimberg @param[in] m Number of rows 648ba59ac12SSebastian Grimberg @param[in] n Number of columns 649ba59ac12SSebastian Grimberg 650ba59ac12SSebastian Grimberg @return An error code: 0 - success, otherwise - failure 651ba59ac12SSebastian Grimberg 652ba59ac12SSebastian Grimberg @ref Utility 653ba59ac12SSebastian Grimberg **/ 654ba59ac12SSebastian Grimberg int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau, CeedInt m, CeedInt n) { 655ba59ac12SSebastian Grimberg CeedScalar v[m]; 656ba59ac12SSebastian Grimberg 657ba59ac12SSebastian Grimberg // Check matrix shape 6586574a04fSJeremy L Thompson CeedCheck(n <= m, ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute QR factorization with n > m"); 659ba59ac12SSebastian Grimberg 660ba59ac12SSebastian Grimberg for (CeedInt i = 0; i < n; i++) { 6611c66c397SJeremy L Thompson CeedScalar sigma = 0.0; 6621c66c397SJeremy L Thompson 663ba59ac12SSebastian Grimberg if (i >= m - 1) { // last row of matrix, no reflection needed 664ba59ac12SSebastian Grimberg tau[i] = 0.; 665ba59ac12SSebastian Grimberg break; 666ba59ac12SSebastian Grimberg } 667ba59ac12SSebastian Grimberg // Calculate Householder vector, magnitude 668ba59ac12SSebastian Grimberg v[i] = mat[i + n * i]; 669ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < m; j++) { 670ba59ac12SSebastian Grimberg v[j] = mat[i + n * j]; 671ba59ac12SSebastian Grimberg sigma += v[j] * v[j]; 672ba59ac12SSebastian Grimberg } 6731c66c397SJeremy L Thompson const CeedScalar norm = sqrt(v[i] * v[i] + sigma); // norm of v[i:m] 6741c66c397SJeremy L Thompson const CeedScalar R_ii = -copysign(norm, v[i]); 6751c66c397SJeremy L Thompson 676ba59ac12SSebastian Grimberg v[i] -= R_ii; 677ba59ac12SSebastian Grimberg // norm of v[i:m] after modification above and scaling below 678ba59ac12SSebastian Grimberg // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 679ba59ac12SSebastian Grimberg // tau = 2 / (norm*norm) 680ba59ac12SSebastian Grimberg tau[i] = 2 * v[i] * v[i] / (v[i] * v[i] + sigma); 681ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < m; j++) v[j] /= v[i]; 682ba59ac12SSebastian Grimberg 683ba59ac12SSebastian Grimberg // Apply Householder reflector to lower right panel 684ba59ac12SSebastian Grimberg CeedHouseholderReflect(&mat[i * n + i + 1], &v[i], tau[i], m - i, n - i - 1, n, 1); 685ba59ac12SSebastian Grimberg // Save v 686ba59ac12SSebastian Grimberg mat[i + n * i] = R_ii; 687ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < m; j++) mat[i + n * j] = v[j]; 688ba59ac12SSebastian Grimberg } 689ba59ac12SSebastian Grimberg return CEED_ERROR_SUCCESS; 690ba59ac12SSebastian Grimberg } 691ba59ac12SSebastian Grimberg 692ba59ac12SSebastian Grimberg /** 693ba59ac12SSebastian Grimberg @brief Apply Householder Q matrix 694ba59ac12SSebastian Grimberg 695ca94c3ddSJeremy L Thompson Compute `mat_A = mat_Q mat_A`, where `mat_Q` is \f$m \times m\f$ and `mat_A` is \f$m \times n\f$. 696ba59ac12SSebastian Grimberg 697ba59ac12SSebastian Grimberg @param[in,out] mat_A Matrix to apply Householder Q to, in place 698ba59ac12SSebastian Grimberg @param[in] mat_Q Householder Q matrix 699ba59ac12SSebastian Grimberg @param[in] tau Householder scaling factors 700ba59ac12SSebastian Grimberg @param[in] t_mode Transpose mode for application 701ca94c3ddSJeremy L Thompson @param[in] m Number of rows in `A` 702ca94c3ddSJeremy L Thompson @param[in] n Number of columns in `A` 703ca94c3ddSJeremy L Thompson @param[in] k Number of elementary reflectors in Q, `k < m` 704ca94c3ddSJeremy L Thompson @param[in] row Row stride in `A` 705ca94c3ddSJeremy L Thompson @param[in] col Col stride in `A` 706ba59ac12SSebastian Grimberg 707ba59ac12SSebastian Grimberg @return An error code: 0 - success, otherwise - failure 708ba59ac12SSebastian Grimberg 709c4e3f59bSSebastian Grimberg @ref Utility 710ba59ac12SSebastian Grimberg **/ 711ba59ac12SSebastian Grimberg int CeedHouseholderApplyQ(CeedScalar *mat_A, const CeedScalar *mat_Q, const CeedScalar *tau, CeedTransposeMode t_mode, CeedInt m, CeedInt n, 712ba59ac12SSebastian Grimberg CeedInt k, CeedInt row, CeedInt col) { 713ba59ac12SSebastian Grimberg CeedScalar *v; 7141c66c397SJeremy L Thompson 715ba59ac12SSebastian Grimberg CeedCall(CeedMalloc(m, &v)); 716ba59ac12SSebastian Grimberg for (CeedInt ii = 0; ii < k; ii++) { 717ba59ac12SSebastian Grimberg CeedInt i = t_mode == CEED_TRANSPOSE ? ii : k - 1 - ii; 718ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < m; j++) v[j] = mat_Q[j * k + i]; 719ba59ac12SSebastian Grimberg // Apply Householder reflector (I - tau v v^T) collo_grad_1d^T 720ba59ac12SSebastian Grimberg CeedCall(CeedHouseholderReflect(&mat_A[i * row], &v[i], tau[i], m - i, n, row, col)); 721ba59ac12SSebastian Grimberg } 722ba59ac12SSebastian Grimberg CeedCall(CeedFree(&v)); 723ba59ac12SSebastian Grimberg return CEED_ERROR_SUCCESS; 724ba59ac12SSebastian Grimberg } 725ba59ac12SSebastian Grimberg 726ba59ac12SSebastian Grimberg /** 7272247a93fSRezgar Shakeri @brief Return pseudoinverse of a matrix 7282247a93fSRezgar Shakeri 7292247a93fSRezgar Shakeri @param[in] ceed Ceed context for error handling 7302247a93fSRezgar Shakeri @param[in] mat Row-major matrix to compute pseudoinverse of 7312247a93fSRezgar Shakeri @param[in] m Number of rows 7322247a93fSRezgar Shakeri @param[in] n Number of columns 7332247a93fSRezgar Shakeri @param[out] mat_pinv Row-major pseudoinverse matrix 7342247a93fSRezgar Shakeri 7352247a93fSRezgar Shakeri @return An error code: 0 - success, otherwise - failure 7362247a93fSRezgar Shakeri 7372247a93fSRezgar Shakeri @ref Utility 7382247a93fSRezgar Shakeri **/ 7391203703bSJeremy L Thompson int CeedMatrixPseudoinverse(Ceed ceed, const CeedScalar *mat, CeedInt m, CeedInt n, CeedScalar *mat_pinv) { 7402247a93fSRezgar Shakeri CeedScalar *tau, *I, *mat_copy; 7412247a93fSRezgar Shakeri 7422247a93fSRezgar Shakeri CeedCall(CeedCalloc(m, &tau)); 7432247a93fSRezgar Shakeri CeedCall(CeedCalloc(m * m, &I)); 7442247a93fSRezgar Shakeri CeedCall(CeedCalloc(m * n, &mat_copy)); 7452247a93fSRezgar Shakeri memcpy(mat_copy, mat, m * n * sizeof mat[0]); 7462247a93fSRezgar Shakeri 7472247a93fSRezgar Shakeri // QR Factorization, mat = Q R 7482247a93fSRezgar Shakeri CeedCall(CeedQRFactorization(ceed, mat_copy, tau, m, n)); 7492247a93fSRezgar Shakeri 7502247a93fSRezgar Shakeri // -- Apply Q^T, I = Q^T * I 7512247a93fSRezgar Shakeri for (CeedInt i = 0; i < m; i++) I[i * m + i] = 1.0; 7522247a93fSRezgar Shakeri CeedCall(CeedHouseholderApplyQ(I, mat_copy, tau, CEED_TRANSPOSE, m, m, n, m, 1)); 7532247a93fSRezgar Shakeri // -- Apply R_inv, mat_pinv = R_inv * Q^T 7542247a93fSRezgar Shakeri for (CeedInt j = 0; j < m; j++) { // Column j 7552247a93fSRezgar Shakeri mat_pinv[j + m * (n - 1)] = I[j + m * (n - 1)] / mat_copy[n * n - 1]; 7562247a93fSRezgar Shakeri for (CeedInt i = n - 2; i >= 0; i--) { // Row i 7572247a93fSRezgar Shakeri mat_pinv[j + m * i] = I[j + m * i]; 7582247a93fSRezgar Shakeri for (CeedInt k = i + 1; k < n; k++) mat_pinv[j + m * i] -= mat_copy[k + n * i] * mat_pinv[j + m * k]; 7592247a93fSRezgar Shakeri mat_pinv[j + m * i] /= mat_copy[i + n * i]; 7602247a93fSRezgar Shakeri } 7612247a93fSRezgar Shakeri } 7622247a93fSRezgar Shakeri 7632247a93fSRezgar Shakeri // Cleanup 7642247a93fSRezgar Shakeri CeedCall(CeedFree(&I)); 7652247a93fSRezgar Shakeri CeedCall(CeedFree(&tau)); 7662247a93fSRezgar Shakeri CeedCall(CeedFree(&mat_copy)); 7672247a93fSRezgar Shakeri return CEED_ERROR_SUCCESS; 7682247a93fSRezgar Shakeri } 7692247a93fSRezgar Shakeri 7702247a93fSRezgar Shakeri /** 771ba59ac12SSebastian Grimberg @brief Return symmetric Schur decomposition of the symmetric matrix mat via symmetric QR factorization 772ba59ac12SSebastian Grimberg 773ca94c3ddSJeremy L Thompson @param[in] ceed `Ceed` context for error handling 774ba59ac12SSebastian Grimberg @param[in,out] mat Row-major matrix to be factorized in place 775ba59ac12SSebastian Grimberg @param[out] lambda Vector of length n of eigenvalues 776ba59ac12SSebastian Grimberg @param[in] n Number of rows/columns 777ba59ac12SSebastian Grimberg 778ba59ac12SSebastian Grimberg @return An error code: 0 - success, otherwise - failure 779ba59ac12SSebastian Grimberg 780ba59ac12SSebastian Grimberg @ref Utility 781ba59ac12SSebastian Grimberg **/ 7822c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOff 7832c2ea1dbSJeremy L Thompson int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat, CeedScalar *lambda, CeedInt n) { 784ba59ac12SSebastian Grimberg // Check bounds for clang-tidy 7856574a04fSJeremy L Thompson CeedCheck(n > 1, ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute symmetric Schur decomposition of scalars"); 786ba59ac12SSebastian Grimberg 787ba59ac12SSebastian Grimberg CeedScalar v[n - 1], tau[n - 1], mat_T[n * n]; 788ba59ac12SSebastian Grimberg 789ba59ac12SSebastian Grimberg // Copy mat to mat_T and set mat to I 790ba59ac12SSebastian Grimberg memcpy(mat_T, mat, n * n * sizeof(mat[0])); 791ba59ac12SSebastian Grimberg for (CeedInt i = 0; i < n; i++) { 792ba59ac12SSebastian Grimberg for (CeedInt j = 0; j < n; j++) mat[j + n * i] = (i == j) ? 1 : 0; 793ba59ac12SSebastian Grimberg } 794ba59ac12SSebastian Grimberg 795ba59ac12SSebastian Grimberg // Reduce to tridiagonal 796ba59ac12SSebastian Grimberg for (CeedInt i = 0; i < n - 1; i++) { 797ba59ac12SSebastian Grimberg // Calculate Householder vector, magnitude 798ba59ac12SSebastian Grimberg CeedScalar sigma = 0.0; 7991c66c397SJeremy L Thompson 800ba59ac12SSebastian Grimberg v[i] = mat_T[i + n * (i + 1)]; 801ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < n - 1; j++) { 802ba59ac12SSebastian Grimberg v[j] = mat_T[i + n * (j + 1)]; 803ba59ac12SSebastian Grimberg sigma += v[j] * v[j]; 804ba59ac12SSebastian Grimberg } 8051c66c397SJeremy L Thompson const CeedScalar norm = sqrt(v[i] * v[i] + sigma); // norm of v[i:n-1] 8061c66c397SJeremy L Thompson const CeedScalar R_ii = -copysign(norm, v[i]); 8071c66c397SJeremy L Thompson 808ba59ac12SSebastian Grimberg v[i] -= R_ii; 809ba59ac12SSebastian Grimberg // norm of v[i:m] after modification above and scaling below 810ba59ac12SSebastian Grimberg // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 811ba59ac12SSebastian Grimberg // tau = 2 / (norm*norm) 812ba59ac12SSebastian Grimberg tau[i] = i == n - 2 ? 2 : 2 * v[i] * v[i] / (v[i] * v[i] + sigma); 813ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < n - 1; j++) v[j] /= v[i]; 814ba59ac12SSebastian Grimberg 815ba59ac12SSebastian Grimberg // Update sub and super diagonal 816ba59ac12SSebastian Grimberg for (CeedInt j = i + 2; j < n; j++) { 817ba59ac12SSebastian Grimberg mat_T[i + n * j] = 0; 818ba59ac12SSebastian Grimberg mat_T[j + n * i] = 0; 819ba59ac12SSebastian Grimberg } 820ba59ac12SSebastian Grimberg // Apply symmetric Householder reflector to lower right panel 821ba59ac12SSebastian Grimberg CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1); 822ba59ac12SSebastian Grimberg CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), 1, n); 823ba59ac12SSebastian Grimberg 824ba59ac12SSebastian Grimberg // Save v 825ba59ac12SSebastian Grimberg mat_T[i + n * (i + 1)] = R_ii; 826ba59ac12SSebastian Grimberg mat_T[(i + 1) + n * i] = R_ii; 827ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < n - 1; j++) { 828ba59ac12SSebastian Grimberg mat_T[i + n * (j + 1)] = v[j]; 829ba59ac12SSebastian Grimberg } 830ba59ac12SSebastian Grimberg } 831ba59ac12SSebastian Grimberg // Backwards accumulation of Q 832ba59ac12SSebastian Grimberg for (CeedInt i = n - 2; i >= 0; i--) { 833ba59ac12SSebastian Grimberg if (tau[i] > 0.0) { 834ba59ac12SSebastian Grimberg v[i] = 1; 835ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < n - 1; j++) { 836ba59ac12SSebastian Grimberg v[j] = mat_T[i + n * (j + 1)]; 837ba59ac12SSebastian Grimberg mat_T[i + n * (j + 1)] = 0; 838ba59ac12SSebastian Grimberg } 839ba59ac12SSebastian Grimberg CeedHouseholderReflect(&mat[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1); 840ba59ac12SSebastian Grimberg } 841ba59ac12SSebastian Grimberg } 842ba59ac12SSebastian Grimberg 843ba59ac12SSebastian Grimberg // Reduce sub and super diagonal 844ba59ac12SSebastian Grimberg CeedInt p = 0, q = 0, itr = 0, max_itr = n * n * n * n; 845ba59ac12SSebastian Grimberg CeedScalar tol = CEED_EPSILON; 846ba59ac12SSebastian Grimberg 847ba59ac12SSebastian Grimberg while (itr < max_itr) { 848ba59ac12SSebastian Grimberg // Update p, q, size of reduced portions of diagonal 849ba59ac12SSebastian Grimberg p = 0; 850ba59ac12SSebastian Grimberg q = 0; 851ba59ac12SSebastian Grimberg for (CeedInt i = n - 2; i >= 0; i--) { 852ba59ac12SSebastian Grimberg if (fabs(mat_T[i + n * (i + 1)]) < tol) q += 1; 853ba59ac12SSebastian Grimberg else break; 854ba59ac12SSebastian Grimberg } 855ba59ac12SSebastian Grimberg for (CeedInt i = 0; i < n - q - 1; i++) { 856ba59ac12SSebastian Grimberg if (fabs(mat_T[i + n * (i + 1)]) < tol) p += 1; 857ba59ac12SSebastian Grimberg else break; 858ba59ac12SSebastian Grimberg } 859ba59ac12SSebastian Grimberg if (q == n - 1) break; // Finished reducing 860ba59ac12SSebastian Grimberg 861ba59ac12SSebastian Grimberg // Reduce tridiagonal portion 862ba59ac12SSebastian Grimberg CeedScalar t_nn = mat_T[(n - 1 - q) + n * (n - 1 - q)], t_nnm1 = mat_T[(n - 2 - q) + n * (n - 1 - q)]; 863ba59ac12SSebastian Grimberg CeedScalar d = (mat_T[(n - 2 - q) + n * (n - 2 - q)] - t_nn) / 2; 864ba59ac12SSebastian Grimberg CeedScalar mu = t_nn - t_nnm1 * t_nnm1 / (d + copysign(sqrt(d * d + t_nnm1 * t_nnm1), d)); 865ba59ac12SSebastian Grimberg CeedScalar x = mat_T[p + n * p] - mu; 866ba59ac12SSebastian Grimberg CeedScalar z = mat_T[p + n * (p + 1)]; 8671c66c397SJeremy L Thompson 868ba59ac12SSebastian Grimberg for (CeedInt k = p; k < n - q - 1; k++) { 869ba59ac12SSebastian Grimberg // Compute Givens rotation 870ba59ac12SSebastian Grimberg CeedScalar c = 1, s = 0; 8711c66c397SJeremy L Thompson 872ba59ac12SSebastian Grimberg if (fabs(z) > tol) { 873ba59ac12SSebastian Grimberg if (fabs(z) > fabs(x)) { 8741c66c397SJeremy L Thompson const CeedScalar tau = -x / z; 8751c66c397SJeremy L Thompson 8761c66c397SJeremy L Thompson s = 1 / sqrt(1 + tau * tau); 8771c66c397SJeremy L Thompson c = s * tau; 878ba59ac12SSebastian Grimberg } else { 8791c66c397SJeremy L Thompson const CeedScalar tau = -z / x; 8801c66c397SJeremy L Thompson 8811c66c397SJeremy L Thompson c = 1 / sqrt(1 + tau * tau); 8821c66c397SJeremy L Thompson s = c * tau; 883ba59ac12SSebastian Grimberg } 884ba59ac12SSebastian Grimberg } 885ba59ac12SSebastian Grimberg 886ba59ac12SSebastian Grimberg // Apply Givens rotation to T 887ba59ac12SSebastian Grimberg CeedGivensRotation(mat_T, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n); 888ba59ac12SSebastian Grimberg CeedGivensRotation(mat_T, c, s, CEED_TRANSPOSE, k, k + 1, n, n); 889ba59ac12SSebastian Grimberg 890ba59ac12SSebastian Grimberg // Apply Givens rotation to Q 891ba59ac12SSebastian Grimberg CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n); 892ba59ac12SSebastian Grimberg 893ba59ac12SSebastian Grimberg // Update x, z 894ba59ac12SSebastian Grimberg if (k < n - q - 2) { 895ba59ac12SSebastian Grimberg x = mat_T[k + n * (k + 1)]; 896ba59ac12SSebastian Grimberg z = mat_T[k + n * (k + 2)]; 897ba59ac12SSebastian Grimberg } 898ba59ac12SSebastian Grimberg } 899ba59ac12SSebastian Grimberg itr++; 900ba59ac12SSebastian Grimberg } 901ba59ac12SSebastian Grimberg 902ba59ac12SSebastian Grimberg // Save eigenvalues 903ba59ac12SSebastian Grimberg for (CeedInt i = 0; i < n; i++) lambda[i] = mat_T[i + n * i]; 904ba59ac12SSebastian Grimberg 905ba59ac12SSebastian Grimberg // Check convergence 9066574a04fSJeremy L Thompson CeedCheck(itr < max_itr || q > n, ceed, CEED_ERROR_MINOR, "Symmetric QR failed to converge"); 907ba59ac12SSebastian Grimberg return CEED_ERROR_SUCCESS; 908ba59ac12SSebastian Grimberg } 9092c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOn 910ba59ac12SSebastian Grimberg 911ba59ac12SSebastian Grimberg /** 912ba59ac12SSebastian Grimberg @brief Return Simultaneous Diagonalization of two matrices. 913ba59ac12SSebastian Grimberg 914ca94c3ddSJeremy L Thompson This solves the generalized eigenvalue problem `A x = lambda B x`, where `A` and `B` are symmetric and `B` is positive definite. 915ca94c3ddSJeremy L Thompson We generate the matrix `X` and vector `Lambda` such that `X^T A X = Lambda` and `X^T B X = I`. 916ca94c3ddSJeremy L Thompson This is equivalent to the LAPACK routine 'sygv' with `TYPE = 1`. 917ba59ac12SSebastian Grimberg 918ca94c3ddSJeremy L Thompson @param[in] ceed `Ceed` context for error handling 919ba59ac12SSebastian Grimberg @param[in] mat_A Row-major matrix to be factorized with eigenvalues 920ba59ac12SSebastian Grimberg @param[in] mat_B Row-major matrix to be factorized to identity 921ba59ac12SSebastian Grimberg @param[out] mat_X Row-major orthogonal matrix 922ca94c3ddSJeremy L Thompson @param[out] lambda Vector of length `n` of generalized eigenvalues 923ba59ac12SSebastian Grimberg @param[in] n Number of rows/columns 924ba59ac12SSebastian Grimberg 925ba59ac12SSebastian Grimberg @return An error code: 0 - success, otherwise - failure 926ba59ac12SSebastian Grimberg 927ba59ac12SSebastian Grimberg @ref Utility 928ba59ac12SSebastian Grimberg **/ 9292c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOff 9302c2ea1dbSJeremy L Thompson int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *mat_A, CeedScalar *mat_B, CeedScalar *mat_X, CeedScalar *lambda, CeedInt n) { 931ba59ac12SSebastian Grimberg CeedScalar *mat_C, *mat_G, *vec_D; 9321c66c397SJeremy L Thompson 933ba59ac12SSebastian Grimberg CeedCall(CeedCalloc(n * n, &mat_C)); 934ba59ac12SSebastian Grimberg CeedCall(CeedCalloc(n * n, &mat_G)); 935ba59ac12SSebastian Grimberg CeedCall(CeedCalloc(n, &vec_D)); 936ba59ac12SSebastian Grimberg 937ba59ac12SSebastian Grimberg // Compute B = G D G^T 938ba59ac12SSebastian Grimberg memcpy(mat_G, mat_B, n * n * sizeof(mat_B[0])); 939ba59ac12SSebastian Grimberg CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_G, vec_D, n)); 940ba59ac12SSebastian Grimberg 941ba59ac12SSebastian Grimberg // Sort eigenvalues 942ba59ac12SSebastian Grimberg for (CeedInt i = n - 1; i >= 0; i--) { 943ba59ac12SSebastian Grimberg for (CeedInt j = 0; j < i; j++) { 944ba59ac12SSebastian Grimberg if (fabs(vec_D[j]) > fabs(vec_D[j + 1])) { 9451c66c397SJeremy L Thompson CeedScalarSwap(vec_D[j], vec_D[j + 1]); 9461c66c397SJeremy L Thompson for (CeedInt k = 0; k < n; k++) CeedScalarSwap(mat_G[k * n + j], mat_G[k * n + j + 1]); 947ba59ac12SSebastian Grimberg } 948ba59ac12SSebastian Grimberg } 949ba59ac12SSebastian Grimberg } 950ba59ac12SSebastian Grimberg 951ba59ac12SSebastian Grimberg // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T 952ba59ac12SSebastian Grimberg // = D^-1/2 G^T A G D^-1/2 953ba59ac12SSebastian Grimberg // -- D = D^-1/2 954ba59ac12SSebastian Grimberg for (CeedInt i = 0; i < n; i++) vec_D[i] = 1. / sqrt(vec_D[i]); 955ba59ac12SSebastian Grimberg // -- G = G D^-1/2 956ba59ac12SSebastian Grimberg // -- C = D^-1/2 G^T 957ba59ac12SSebastian Grimberg for (CeedInt i = 0; i < n; i++) { 958ba59ac12SSebastian Grimberg for (CeedInt j = 0; j < n; j++) { 959ba59ac12SSebastian Grimberg mat_G[i * n + j] *= vec_D[j]; 960ba59ac12SSebastian Grimberg mat_C[j * n + i] = mat_G[i * n + j]; 961ba59ac12SSebastian Grimberg } 962ba59ac12SSebastian Grimberg } 963ba59ac12SSebastian Grimberg // -- X = (D^-1/2 G^T) A 964ba59ac12SSebastian Grimberg CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_C, (const CeedScalar *)mat_A, mat_X, n, n, n)); 965ba59ac12SSebastian Grimberg // -- C = (D^-1/2 G^T A) (G D^-1/2) 966ba59ac12SSebastian Grimberg CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_X, (const CeedScalar *)mat_G, mat_C, n, n, n)); 967ba59ac12SSebastian Grimberg 968ba59ac12SSebastian Grimberg // Compute Q^T C Q = lambda 969ba59ac12SSebastian Grimberg CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_C, lambda, n)); 970ba59ac12SSebastian Grimberg 971ba59ac12SSebastian Grimberg // Sort eigenvalues 972ba59ac12SSebastian Grimberg for (CeedInt i = n - 1; i >= 0; i--) { 973ba59ac12SSebastian Grimberg for (CeedInt j = 0; j < i; j++) { 974ba59ac12SSebastian Grimberg if (fabs(lambda[j]) > fabs(lambda[j + 1])) { 9751c66c397SJeremy L Thompson CeedScalarSwap(lambda[j], lambda[j + 1]); 9761c66c397SJeremy L Thompson for (CeedInt k = 0; k < n; k++) CeedScalarSwap(mat_C[k * n + j], mat_C[k * n + j + 1]); 977ba59ac12SSebastian Grimberg } 978ba59ac12SSebastian Grimberg } 979ba59ac12SSebastian Grimberg } 980ba59ac12SSebastian Grimberg 981ba59ac12SSebastian Grimberg // Set X = (G D^1/2)^-T Q 982ba59ac12SSebastian Grimberg // = G D^-1/2 Q 983ba59ac12SSebastian Grimberg CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_G, (const CeedScalar *)mat_C, mat_X, n, n, n)); 984ba59ac12SSebastian Grimberg 985ba59ac12SSebastian Grimberg // Cleanup 986ba59ac12SSebastian Grimberg CeedCall(CeedFree(&mat_C)); 987ba59ac12SSebastian Grimberg CeedCall(CeedFree(&mat_G)); 988ba59ac12SSebastian Grimberg CeedCall(CeedFree(&vec_D)); 989ba59ac12SSebastian Grimberg return CEED_ERROR_SUCCESS; 990ba59ac12SSebastian Grimberg } 9912c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOn 992ba59ac12SSebastian Grimberg 9937a982d89SJeremy L. Thompson /// @} 9947a982d89SJeremy L. Thompson 9957a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 9967a982d89SJeremy L. Thompson /// CeedBasis Public API 9977a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 9987a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser 999d7b241e6Sjeremylt /// @{ 1000d7b241e6Sjeremylt 1001b11c1e72Sjeremylt /** 1002ca94c3ddSJeremy L Thompson @brief Create a tensor-product basis for \f$H^1\f$ discretizations 1003b11c1e72Sjeremylt 1004ca94c3ddSJeremy L Thompson @param[in] ceed `Ceed` object used to create the `CeedBasis` 1005ea61e9acSJeremy L Thompson @param[in] dim Topological dimension 1006ea61e9acSJeremy L Thompson @param[in] num_comp Number of field components (1 for scalar fields) 1007ea61e9acSJeremy L Thompson @param[in] P_1d Number of nodes in one dimension 1008ea61e9acSJeremy L Thompson @param[in] Q_1d Number of quadrature points in one dimension 1009ca94c3ddSJeremy L Thompson @param[in] interp_1d Row-major (`Q_1d * P_1d`) matrix expressing the values of nodal basis functions at quadrature points 1010ca94c3ddSJeremy L Thompson @param[in] grad_1d Row-major (`Q_1d * P_1d`) matrix expressing derivatives of nodal basis functions at quadrature points 1011ca94c3ddSJeremy 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]` 1012ca94c3ddSJeremy L Thompson @param[in] q_weight_1d Array of length `Q_1d` holding the quadrature weights on the reference element 1013ca94c3ddSJeremy L Thompson @param[out] basis Address of the variable where the newly created `CeedBasis` will be stored 1014b11c1e72Sjeremylt 1015b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 1016dfdf5a53Sjeremylt 10177a982d89SJeremy L. Thompson @ref User 1018b11c1e72Sjeremylt **/ 10192b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, 10202b730f8bSJeremy L Thompson const CeedScalar *grad_1d, const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis *basis) { 10215fe0d4faSjeremylt if (!ceed->BasisCreateTensorH1) { 10225fe0d4faSjeremylt Ceed delegate; 10236574a04fSJeremy L Thompson 10242b730f8bSJeremy L Thompson CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis")); 10251ef3a2a9SJeremy L Thompson CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateTensorH1"); 10262b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateTensorH1(delegate, dim, num_comp, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis)); 1027e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 10285fe0d4faSjeremylt } 1029e15f9bd0SJeremy L Thompson 1030ca94c3ddSJeremy L Thompson CeedCheck(dim > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis dimension must be a positive value"); 1031ca94c3ddSJeremy L Thompson CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 component"); 1032ca94c3ddSJeremy L Thompson CeedCheck(P_1d > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 node"); 1033ca94c3ddSJeremy L Thompson CeedCheck(Q_1d > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 quadrature point"); 1034227444bfSJeremy L Thompson 10352b730f8bSJeremy L Thompson CeedElemTopology topo = dim == 1 ? CEED_TOPOLOGY_LINE : dim == 2 ? CEED_TOPOLOGY_QUAD : CEED_TOPOLOGY_HEX; 1036e15f9bd0SJeremy L Thompson 10372b730f8bSJeremy L Thompson CeedCall(CeedCalloc(1, basis)); 1038db002c03SJeremy L Thompson CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed)); 1039d1d35e2fSjeremylt (*basis)->ref_count = 1; 10406402da51SJeremy L Thompson (*basis)->is_tensor_basis = true; 1041d7b241e6Sjeremylt (*basis)->dim = dim; 1042d99fa3c5SJeremy L Thompson (*basis)->topo = topo; 1043d1d35e2fSjeremylt (*basis)->num_comp = num_comp; 1044d1d35e2fSjeremylt (*basis)->P_1d = P_1d; 1045d1d35e2fSjeremylt (*basis)->Q_1d = Q_1d; 1046d1d35e2fSjeremylt (*basis)->P = CeedIntPow(P_1d, dim); 1047d1d35e2fSjeremylt (*basis)->Q = CeedIntPow(Q_1d, dim); 1048c4e3f59bSSebastian Grimberg (*basis)->fe_space = CEED_FE_SPACE_H1; 10492b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q_1d, &(*basis)->q_ref_1d)); 10502b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q_1d, &(*basis)->q_weight_1d)); 1051ff3a0f91SJeremy L Thompson if (q_ref_1d) memcpy((*basis)->q_ref_1d, q_ref_1d, Q_1d * sizeof(q_ref_1d[0])); 10522b730f8bSJeremy L Thompson if (q_weight_1d) memcpy((*basis)->q_weight_1d, q_weight_1d, Q_1d * sizeof(q_weight_1d[0])); 10532b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->interp_1d)); 10542b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->grad_1d)); 10552b730f8bSJeremy L Thompson if (interp_1d) memcpy((*basis)->interp_1d, interp_1d, Q_1d * P_1d * sizeof(interp_1d[0])); 1056ff3a0f91SJeremy L Thompson if (grad_1d) memcpy((*basis)->grad_1d, grad_1d, Q_1d * P_1d * sizeof(grad_1d[0])); 10572b730f8bSJeremy L Thompson CeedCall(ceed->BasisCreateTensorH1(dim, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, *basis)); 1058e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1059d7b241e6Sjeremylt } 1060d7b241e6Sjeremylt 1061b11c1e72Sjeremylt /** 1062ca94c3ddSJeremy L Thompson @brief Create a tensor-product \f$H^1\f$ Lagrange basis 1063b11c1e72Sjeremylt 1064ca94c3ddSJeremy L Thompson @param[in] ceed `Ceed` object used to create the `CeedBasis` 1065ea61e9acSJeremy L Thompson @param[in] dim Topological dimension of element 1066ea61e9acSJeremy L Thompson @param[in] num_comp Number of field components (1 for scalar fields) 1067ea61e9acSJeremy L Thompson @param[in] P Number of Gauss-Lobatto nodes in one dimension. 1068ca94c3ddSJeremy L Thompson The polynomial degree of the resulting `Q_k` element is `k = P - 1`. 1069ea61e9acSJeremy L Thompson @param[in] Q Number of quadrature points in one dimension. 1070ca94c3ddSJeremy L Thompson @param[in] quad_mode Distribution of the `Q` quadrature points (affects order of accuracy for the quadrature) 1071ca94c3ddSJeremy L Thompson @param[out] basis Address of the variable where the newly created `CeedBasis` will be stored 1072b11c1e72Sjeremylt 1073b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 1074dfdf5a53Sjeremylt 10757a982d89SJeremy L. Thompson @ref User 1076b11c1e72Sjeremylt **/ 10772b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P, CeedInt Q, CeedQuadMode quad_mode, CeedBasis *basis) { 1078d7b241e6Sjeremylt // Allocate 1079c8c3fa7dSJeremy L Thompson int ierr = CEED_ERROR_SUCCESS; 10802b730f8bSJeremy L Thompson CeedScalar c1, c2, c3, c4, dx, *nodes, *interp_1d, *grad_1d, *q_ref_1d, *q_weight_1d; 10814d537eeaSYohann 1082ca94c3ddSJeremy L Thompson CeedCheck(dim > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis dimension must be a positive value"); 1083ca94c3ddSJeremy L Thompson CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 component"); 1084ca94c3ddSJeremy L Thompson CeedCheck(P > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 node"); 1085ca94c3ddSJeremy L Thompson CeedCheck(Q > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 quadrature point"); 1086227444bfSJeremy L Thompson 1087e15f9bd0SJeremy L Thompson // Get Nodes and Weights 10882b730f8bSJeremy L Thompson CeedCall(CeedCalloc(P * Q, &interp_1d)); 10892b730f8bSJeremy L Thompson CeedCall(CeedCalloc(P * Q, &grad_1d)); 10902b730f8bSJeremy L Thompson CeedCall(CeedCalloc(P, &nodes)); 10912b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q, &q_ref_1d)); 10922b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q, &q_weight_1d)); 10932b730f8bSJeremy L Thompson if (CeedLobattoQuadrature(P, nodes, NULL) != CEED_ERROR_SUCCESS) goto cleanup; 1094d1d35e2fSjeremylt switch (quad_mode) { 1095d7b241e6Sjeremylt case CEED_GAUSS: 1096d1d35e2fSjeremylt ierr = CeedGaussQuadrature(Q, q_ref_1d, q_weight_1d); 1097d7b241e6Sjeremylt break; 1098d7b241e6Sjeremylt case CEED_GAUSS_LOBATTO: 1099d1d35e2fSjeremylt ierr = CeedLobattoQuadrature(Q, q_ref_1d, q_weight_1d); 1100d7b241e6Sjeremylt break; 1101d7b241e6Sjeremylt } 11022b730f8bSJeremy L Thompson if (ierr != CEED_ERROR_SUCCESS) goto cleanup; 1103e15f9bd0SJeremy L Thompson 1104d7b241e6Sjeremylt // Build B, D matrix 1105d7b241e6Sjeremylt // Fornberg, 1998 1106c8c3fa7dSJeremy L Thompson for (CeedInt i = 0; i < Q; i++) { 1107d7b241e6Sjeremylt c1 = 1.0; 1108d1d35e2fSjeremylt c3 = nodes[0] - q_ref_1d[i]; 1109d1d35e2fSjeremylt interp_1d[i * P + 0] = 1.0; 1110c8c3fa7dSJeremy L Thompson for (CeedInt j = 1; j < P; j++) { 1111d7b241e6Sjeremylt c2 = 1.0; 1112d7b241e6Sjeremylt c4 = c3; 1113d1d35e2fSjeremylt c3 = nodes[j] - q_ref_1d[i]; 1114c8c3fa7dSJeremy L Thompson for (CeedInt k = 0; k < j; k++) { 1115d7b241e6Sjeremylt dx = nodes[j] - nodes[k]; 1116d7b241e6Sjeremylt c2 *= dx; 1117d7b241e6Sjeremylt if (k == j - 1) { 1118d1d35e2fSjeremylt grad_1d[i * P + j] = c1 * (interp_1d[i * P + k] - c4 * grad_1d[i * P + k]) / c2; 1119d1d35e2fSjeremylt interp_1d[i * P + j] = -c1 * c4 * interp_1d[i * P + k] / c2; 1120d7b241e6Sjeremylt } 1121d1d35e2fSjeremylt grad_1d[i * P + k] = (c3 * grad_1d[i * P + k] - interp_1d[i * P + k]) / dx; 1122d1d35e2fSjeremylt interp_1d[i * P + k] = c3 * interp_1d[i * P + k] / dx; 1123d7b241e6Sjeremylt } 1124d7b241e6Sjeremylt c1 = c2; 1125d7b241e6Sjeremylt } 1126d7b241e6Sjeremylt } 11279ac7b42eSJeremy L Thompson // Pass to CeedBasisCreateTensorH1 11282b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P, Q, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis)); 1129e15f9bd0SJeremy L Thompson cleanup: 11302b730f8bSJeremy L Thompson CeedCall(CeedFree(&interp_1d)); 11312b730f8bSJeremy L Thompson CeedCall(CeedFree(&grad_1d)); 11322b730f8bSJeremy L Thompson CeedCall(CeedFree(&nodes)); 11332b730f8bSJeremy L Thompson CeedCall(CeedFree(&q_ref_1d)); 11342b730f8bSJeremy L Thompson CeedCall(CeedFree(&q_weight_1d)); 1135e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1136d7b241e6Sjeremylt } 1137d7b241e6Sjeremylt 1138b11c1e72Sjeremylt /** 1139ca94c3ddSJeremy L Thompson @brief Create a non tensor-product basis for \f$H^1\f$ discretizations 1140a8de75f0Sjeremylt 1141ca94c3ddSJeremy L Thompson @param[in] ceed `Ceed` object used to create the `CeedBasis` 1142*e00f3be8SJames Wright @param[in] topo Topology of element, e.g. hypercube, simplex, etc 1143ea61e9acSJeremy L Thompson @param[in] num_comp Number of field components (1 for scalar fields) 1144ea61e9acSJeremy L Thompson @param[in] num_nodes Total number of nodes 1145ea61e9acSJeremy L Thompson @param[in] num_qpts Total number of quadrature points 1146ca94c3ddSJeremy L Thompson @param[in] interp Row-major (`num_qpts * num_nodes`) matrix expressing the values of nodal basis functions at quadrature points 1147ca94c3ddSJeremy L Thompson @param[in] grad Row-major (`dim * num_qpts * num_nodes`) matrix expressing derivatives of nodal basis functions at quadrature points 1148ca94c3ddSJeremy L Thompson @param[in] q_ref Array of length `num_qpts` * dim holding the locations of quadrature points on the reference element 1149ca94c3ddSJeremy L Thompson @param[in] q_weight Array of length `num_qpts` holding the quadrature weights on the reference element 1150ca94c3ddSJeremy L Thompson @param[out] basis Address of the variable where the newly created `CeedBasis` will be stored 1151a8de75f0Sjeremylt 1152a8de75f0Sjeremylt @return An error code: 0 - success, otherwise - failure 1153a8de75f0Sjeremylt 11547a982d89SJeremy L. Thompson @ref User 1155a8de75f0Sjeremylt **/ 11562b730f8bSJeremy L Thompson int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 11572b730f8bSJeremy L Thompson const CeedScalar *grad, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) { 1158d1d35e2fSjeremylt CeedInt P = num_nodes, Q = num_qpts, dim = 0; 1159a8de75f0Sjeremylt 11605fe0d4faSjeremylt if (!ceed->BasisCreateH1) { 11615fe0d4faSjeremylt Ceed delegate; 11626574a04fSJeremy L Thompson 11632b730f8bSJeremy L Thompson CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis")); 11641ef3a2a9SJeremy L Thompson CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateH1"); 11652b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateH1(delegate, topo, num_comp, num_nodes, num_qpts, interp, grad, q_ref, q_weight, basis)); 1166e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 11675fe0d4faSjeremylt } 11685fe0d4faSjeremylt 1169ca94c3ddSJeremy L Thompson CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 component"); 1170ca94c3ddSJeremy L Thompson CeedCheck(num_nodes > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 node"); 1171ca94c3ddSJeremy L Thompson CeedCheck(num_qpts > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 quadrature point"); 1172227444bfSJeremy L Thompson 11732b730f8bSJeremy L Thompson CeedCall(CeedBasisGetTopologyDimension(topo, &dim)); 1174a8de75f0Sjeremylt 1175db002c03SJeremy L Thompson CeedCall(CeedCalloc(1, basis)); 1176db002c03SJeremy L Thompson CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed)); 1177d1d35e2fSjeremylt (*basis)->ref_count = 1; 11786402da51SJeremy L Thompson (*basis)->is_tensor_basis = false; 1179a8de75f0Sjeremylt (*basis)->dim = dim; 1180d99fa3c5SJeremy L Thompson (*basis)->topo = topo; 1181d1d35e2fSjeremylt (*basis)->num_comp = num_comp; 1182a8de75f0Sjeremylt (*basis)->P = P; 1183a8de75f0Sjeremylt (*basis)->Q = Q; 1184c4e3f59bSSebastian Grimberg (*basis)->fe_space = CEED_FE_SPACE_H1; 11852b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q * dim, &(*basis)->q_ref_1d)); 11862b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q, &(*basis)->q_weight_1d)); 1187ff3a0f91SJeremy L Thompson if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0])); 1188ff3a0f91SJeremy L Thompson if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0])); 11892b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q * P, &(*basis)->interp)); 11902b730f8bSJeremy L Thompson CeedCall(CeedCalloc(dim * Q * P, &(*basis)->grad)); 1191ff3a0f91SJeremy L Thompson if (interp) memcpy((*basis)->interp, interp, Q * P * sizeof(interp[0])); 1192ff3a0f91SJeremy L Thompson if (grad) memcpy((*basis)->grad, grad, dim * Q * P * sizeof(grad[0])); 11932b730f8bSJeremy L Thompson CeedCall(ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, q_ref, q_weight, *basis)); 1194e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1195a8de75f0Sjeremylt } 1196a8de75f0Sjeremylt 1197a8de75f0Sjeremylt /** 1198859c15bbSJames Wright @brief Create a non tensor-product basis for \f$H(\mathrm{div})\f$ discretizations 119950c301a5SRezgar Shakeri 1200ca94c3ddSJeremy L Thompson @param[in] ceed `Ceed` object used to create the `CeedBasis` 1201ea61e9acSJeremy 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 1202ea61e9acSJeremy L Thompson @param[in] num_comp Number of components (usually 1 for vectors in H(div) bases) 1203ca94c3ddSJeremy L Thompson @param[in] num_nodes Total number of nodes (DoFs per element) 1204ea61e9acSJeremy L Thompson @param[in] num_qpts Total number of quadrature points 1205ca94c3ddSJeremy L Thompson @param[in] interp Row-major (`dim * num_qpts * num_nodes`) matrix expressing the values of basis functions at quadrature points 1206ca94c3ddSJeremy L Thompson @param[in] div Row-major (`num_qpts * num_nodes`) matrix expressing divergence of basis functions at quadrature points 1207ca94c3ddSJeremy L Thompson @param[in] q_ref Array of length `num_qpts` * dim holding the locations of quadrature points on the reference element 1208ca94c3ddSJeremy L Thompson @param[in] q_weight Array of length `num_qpts` holding the quadrature weights on the reference element 1209ca94c3ddSJeremy L Thompson @param[out] basis Address of the variable where the newly created `CeedBasis` will be stored 121050c301a5SRezgar Shakeri 121150c301a5SRezgar Shakeri @return An error code: 0 - success, otherwise - failure 121250c301a5SRezgar Shakeri 121350c301a5SRezgar Shakeri @ref User 121450c301a5SRezgar Shakeri **/ 12152b730f8bSJeremy L Thompson int CeedBasisCreateHdiv(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 12162b730f8bSJeremy L Thompson const CeedScalar *div, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) { 121750c301a5SRezgar Shakeri CeedInt Q = num_qpts, P = num_nodes, dim = 0; 1218c4e3f59bSSebastian Grimberg 121950c301a5SRezgar Shakeri if (!ceed->BasisCreateHdiv) { 122050c301a5SRezgar Shakeri Ceed delegate; 12216574a04fSJeremy L Thompson 12222b730f8bSJeremy L Thompson CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis")); 12236574a04fSJeremy L Thompson CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateHdiv"); 12242b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateHdiv(delegate, topo, num_comp, num_nodes, num_qpts, interp, div, q_ref, q_weight, basis)); 122550c301a5SRezgar Shakeri return CEED_ERROR_SUCCESS; 122650c301a5SRezgar Shakeri } 122750c301a5SRezgar Shakeri 1228ca94c3ddSJeremy L Thompson CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 component"); 1229ca94c3ddSJeremy L Thompson CeedCheck(num_nodes > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 node"); 1230ca94c3ddSJeremy L Thompson CeedCheck(num_qpts > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 quadrature point"); 1231227444bfSJeremy L Thompson 1232c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetTopologyDimension(topo, &dim)); 1233c4e3f59bSSebastian Grimberg 1234db002c03SJeremy L Thompson CeedCall(CeedCalloc(1, basis)); 1235db002c03SJeremy L Thompson CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed)); 123650c301a5SRezgar Shakeri (*basis)->ref_count = 1; 12376402da51SJeremy L Thompson (*basis)->is_tensor_basis = false; 123850c301a5SRezgar Shakeri (*basis)->dim = dim; 123950c301a5SRezgar Shakeri (*basis)->topo = topo; 124050c301a5SRezgar Shakeri (*basis)->num_comp = num_comp; 124150c301a5SRezgar Shakeri (*basis)->P = P; 124250c301a5SRezgar Shakeri (*basis)->Q = Q; 1243c4e3f59bSSebastian Grimberg (*basis)->fe_space = CEED_FE_SPACE_HDIV; 12442b730f8bSJeremy L Thompson CeedCall(CeedMalloc(Q * dim, &(*basis)->q_ref_1d)); 12452b730f8bSJeremy L Thompson CeedCall(CeedMalloc(Q, &(*basis)->q_weight_1d)); 124650c301a5SRezgar Shakeri if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0])); 124750c301a5SRezgar Shakeri if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0])); 12482b730f8bSJeremy L Thompson CeedCall(CeedMalloc(dim * Q * P, &(*basis)->interp)); 12492b730f8bSJeremy L Thompson CeedCall(CeedMalloc(Q * P, &(*basis)->div)); 125050c301a5SRezgar Shakeri if (interp) memcpy((*basis)->interp, interp, dim * Q * P * sizeof(interp[0])); 125150c301a5SRezgar Shakeri if (div) memcpy((*basis)->div, div, Q * P * sizeof(div[0])); 12522b730f8bSJeremy L Thompson CeedCall(ceed->BasisCreateHdiv(topo, dim, P, Q, interp, div, q_ref, q_weight, *basis)); 125350c301a5SRezgar Shakeri return CEED_ERROR_SUCCESS; 125450c301a5SRezgar Shakeri } 125550c301a5SRezgar Shakeri 125650c301a5SRezgar Shakeri /** 12574385fb7fSSebastian Grimberg @brief Create a non tensor-product basis for \f$H(\mathrm{curl})\f$ discretizations 1258c4e3f59bSSebastian Grimberg 1259ca94c3ddSJeremy L Thompson @param[in] ceed `Ceed` object used to create the `CeedBasis` 1260c4e3f59bSSebastian Grimberg @param[in] topo Topology of element (`CEED_TOPOLOGY_QUAD`, `CEED_TOPOLOGY_PRISM`, etc.), dimension of which is used in some array sizes below 1261ca94c3ddSJeremy L Thompson @param[in] num_comp Number of components (usually 1 for vectors in \f$H(\mathrm{curl})\f$ bases) 1262ca94c3ddSJeremy L Thompson @param[in] num_nodes Total number of nodes (DoFs per element) 1263c4e3f59bSSebastian Grimberg @param[in] num_qpts Total number of quadrature points 1264ca94c3ddSJeremy L Thompson @param[in] interp Row-major (`dim * num_qpts * num_nodes`) matrix expressing the values of basis functions at quadrature points 1265ca94c3ddSJeremy L Thompson @param[in] curl Row-major (`curl_comp * num_qpts * num_nodes`, `curl_comp = 1` if `dim < 3` otherwise `curl_comp = dim`) matrix expressing curl of basis functions at quadrature points 1266ca94c3ddSJeremy L Thompson @param[in] q_ref Array of length `num_qpts * dim` holding the locations of quadrature points on the reference element 1267ca94c3ddSJeremy L Thompson @param[in] q_weight Array of length `num_qpts` holding the quadrature weights on the reference element 1268ca94c3ddSJeremy L Thompson @param[out] basis Address of the variable where the newly created `CeedBasis` will be stored 1269c4e3f59bSSebastian Grimberg 1270c4e3f59bSSebastian Grimberg @return An error code: 0 - success, otherwise - failure 1271c4e3f59bSSebastian Grimberg 1272c4e3f59bSSebastian Grimberg @ref User 1273c4e3f59bSSebastian Grimberg **/ 1274c4e3f59bSSebastian Grimberg int CeedBasisCreateHcurl(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 1275c4e3f59bSSebastian Grimberg const CeedScalar *curl, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) { 1276c4e3f59bSSebastian Grimberg CeedInt Q = num_qpts, P = num_nodes, dim = 0, curl_comp = 0; 1277c4e3f59bSSebastian Grimberg 1278d075f50bSSebastian Grimberg if (!ceed->BasisCreateHcurl) { 1279c4e3f59bSSebastian Grimberg Ceed delegate; 12806574a04fSJeremy L Thompson 1281c4e3f59bSSebastian Grimberg CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis")); 12826574a04fSJeremy L Thompson CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateHcurl"); 1283c4e3f59bSSebastian Grimberg CeedCall(CeedBasisCreateHcurl(delegate, topo, num_comp, num_nodes, num_qpts, interp, curl, q_ref, q_weight, basis)); 1284c4e3f59bSSebastian Grimberg return CEED_ERROR_SUCCESS; 1285c4e3f59bSSebastian Grimberg } 1286c4e3f59bSSebastian Grimberg 1287ca94c3ddSJeremy L Thompson CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 component"); 1288ca94c3ddSJeremy L Thompson CeedCheck(num_nodes > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 node"); 1289ca94c3ddSJeremy L Thompson CeedCheck(num_qpts > 0, ceed, CEED_ERROR_DIMENSION, "CeedBasis must have at least 1 quadrature point"); 1290c4e3f59bSSebastian Grimberg 1291c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetTopologyDimension(topo, &dim)); 1292c4e3f59bSSebastian Grimberg curl_comp = (dim < 3) ? 1 : dim; 1293c4e3f59bSSebastian Grimberg 1294db002c03SJeremy L Thompson CeedCall(CeedCalloc(1, basis)); 1295db002c03SJeremy L Thompson CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed)); 1296c4e3f59bSSebastian Grimberg (*basis)->ref_count = 1; 12976402da51SJeremy L Thompson (*basis)->is_tensor_basis = false; 1298c4e3f59bSSebastian Grimberg (*basis)->dim = dim; 1299c4e3f59bSSebastian Grimberg (*basis)->topo = topo; 1300c4e3f59bSSebastian Grimberg (*basis)->num_comp = num_comp; 1301c4e3f59bSSebastian Grimberg (*basis)->P = P; 1302c4e3f59bSSebastian Grimberg (*basis)->Q = Q; 1303c4e3f59bSSebastian Grimberg (*basis)->fe_space = CEED_FE_SPACE_HCURL; 1304c4e3f59bSSebastian Grimberg CeedCall(CeedMalloc(Q * dim, &(*basis)->q_ref_1d)); 1305c4e3f59bSSebastian Grimberg CeedCall(CeedMalloc(Q, &(*basis)->q_weight_1d)); 1306c4e3f59bSSebastian Grimberg if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0])); 1307c4e3f59bSSebastian Grimberg if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0])); 1308c4e3f59bSSebastian Grimberg CeedCall(CeedMalloc(dim * Q * P, &(*basis)->interp)); 1309c4e3f59bSSebastian Grimberg CeedCall(CeedMalloc(curl_comp * Q * P, &(*basis)->curl)); 1310c4e3f59bSSebastian Grimberg if (interp) memcpy((*basis)->interp, interp, dim * Q * P * sizeof(interp[0])); 1311c4e3f59bSSebastian Grimberg if (curl) memcpy((*basis)->curl, curl, curl_comp * Q * P * sizeof(curl[0])); 1312c4e3f59bSSebastian Grimberg CeedCall(ceed->BasisCreateHcurl(topo, dim, P, Q, interp, curl, q_ref, q_weight, *basis)); 1313c4e3f59bSSebastian Grimberg return CEED_ERROR_SUCCESS; 1314c4e3f59bSSebastian Grimberg } 1315c4e3f59bSSebastian Grimberg 1316c4e3f59bSSebastian Grimberg /** 1317ca94c3ddSJeremy L Thompson @brief Create a `CeedBasis` for projection from the nodes of `basis_from` to the nodes of `basis_to`. 1318ba59ac12SSebastian Grimberg 1319ca94c3ddSJeremy L Thompson Only @ref CEED_EVAL_INTERP will be valid for the new basis, `basis_project`. 1320ca94c3ddSJeremy L Thompson For \f$H^1\f$ spaces, @ref CEED_EVAL_GRAD will also be valid. 1321ca94c3ddSJeremy L Thompson The interpolation is given by `interp_project = interp_to^+ * interp_from`, where the pseudoinverse `interp_to^+` is given by QR factorization. 1322ca94c3ddSJeremy L Thompson The gradient (for the \f$H^1\f$ case) is given by `grad_project = interp_to^+ * grad_from`. 132315ad3917SSebastian Grimberg 132415ad3917SSebastian Grimberg Note: `basis_from` and `basis_to` must have compatible quadrature spaces. 132515ad3917SSebastian Grimberg 13269fd66db6SSebastian Grimberg Note: `basis_project` will have the same number of components as `basis_from`, regardless of the number of components that `basis_to` has. 13279fd66db6SSebastian Grimberg If `basis_from` has 3 components and `basis_to` has 5 components, then `basis_project` will have 3 components. 1328f113e5dcSJeremy L Thompson 1329e104ad11SJames Wright Note: If either `basis_from` or `basis_to` are non-tensor, then `basis_project` will also be non-tensor 1330e104ad11SJames Wright 1331ca94c3ddSJeremy L Thompson @param[in] basis_from `CeedBasis` to prolong from 1332ca94c3ddSJeremy L Thompson @param[in] basis_to `CeedBasis` to prolong to 1333ca94c3ddSJeremy L Thompson @param[out] basis_project Address of the variable where the newly created `CeedBasis` will be stored 1334f113e5dcSJeremy L Thompson 1335f113e5dcSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 1336f113e5dcSJeremy L Thompson 1337f113e5dcSJeremy L Thompson @ref User 1338f113e5dcSJeremy L Thompson **/ 13392b730f8bSJeremy L Thompson int CeedBasisCreateProjection(CeedBasis basis_from, CeedBasis basis_to, CeedBasis *basis_project) { 1340f113e5dcSJeremy L Thompson Ceed ceed; 1341e104ad11SJames Wright bool create_tensor; 13421c66c397SJeremy L Thompson CeedInt dim, num_comp; 1343097cc795SJames Wright CeedScalar *interp_project, *grad_project; 13441c66c397SJeremy L Thompson 13452b730f8bSJeremy L Thompson CeedCall(CeedBasisGetCeed(basis_to, &ceed)); 1346f113e5dcSJeremy L Thompson 1347ecc88aebSJeremy L Thompson // Create projection matrix 13482b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateProjectionMatrices(basis_from, basis_to, &interp_project, &grad_project)); 1349f113e5dcSJeremy L Thompson 1350f113e5dcSJeremy L Thompson // Build basis 1351e104ad11SJames Wright { 1352e104ad11SJames Wright bool is_tensor_to, is_tensor_from; 1353e104ad11SJames Wright 1354e104ad11SJames Wright CeedCall(CeedBasisIsTensor(basis_to, &is_tensor_to)); 1355e104ad11SJames Wright CeedCall(CeedBasisIsTensor(basis_from, &is_tensor_from)); 1356e104ad11SJames Wright create_tensor = is_tensor_from && is_tensor_to; 1357e104ad11SJames Wright } 13582b730f8bSJeremy L Thompson CeedCall(CeedBasisGetDimension(basis_to, &dim)); 13592b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumComponents(basis_from, &num_comp)); 1360e104ad11SJames Wright if (create_tensor) { 1361f113e5dcSJeremy L Thompson CeedInt P_1d_to, P_1d_from; 13621c66c397SJeremy L Thompson 13632b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_1d_from)); 13642b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_1d_to)); 1365097cc795SJames Wright CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P_1d_from, P_1d_to, interp_project, grad_project, NULL, NULL, basis_project)); 1366f113e5dcSJeremy L Thompson } else { 1367de05fbb2SSebastian Grimberg // Even if basis_to and basis_from are not H1, the resulting basis is H1 for interpolation to work 1368f113e5dcSJeremy L Thompson CeedInt num_nodes_to, num_nodes_from; 13691c66c397SJeremy L Thompson CeedElemTopology topo; 13701c66c397SJeremy L Thompson 1371*e00f3be8SJames Wright CeedCall(CeedBasisGetTopology(basis_from, &topo)); 13722b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis_from, &num_nodes_from)); 13732b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis_to, &num_nodes_to)); 1374097cc795SJames Wright CeedCall(CeedBasisCreateH1(ceed, topo, num_comp, num_nodes_from, num_nodes_to, interp_project, grad_project, NULL, NULL, basis_project)); 1375f113e5dcSJeremy L Thompson } 1376f113e5dcSJeremy L Thompson 1377f113e5dcSJeremy L Thompson // Cleanup 13782b730f8bSJeremy L Thompson CeedCall(CeedFree(&interp_project)); 13792b730f8bSJeremy L Thompson CeedCall(CeedFree(&grad_project)); 1380f113e5dcSJeremy L Thompson return CEED_ERROR_SUCCESS; 1381f113e5dcSJeremy L Thompson } 1382f113e5dcSJeremy L Thompson 1383f113e5dcSJeremy L Thompson /** 1384ca94c3ddSJeremy L Thompson @brief Copy the pointer to a `CeedBasis`. 13859560d06aSjeremylt 1386ca94c3ddSJeremy 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`. 1387ca94c3ddSJeremy L Thompson This `CeedBasis` will be destroyed if `*basis_copy` is the only reference to this `CeedBasis`. 1388ea61e9acSJeremy L Thompson 1389ca94c3ddSJeremy L Thompson @param[in] basis `CeedBasis` to copy reference to 1390ea61e9acSJeremy L Thompson @param[in,out] basis_copy Variable to store copied reference 13919560d06aSjeremylt 13929560d06aSjeremylt @return An error code: 0 - success, otherwise - failure 13939560d06aSjeremylt 13949560d06aSjeremylt @ref User 13959560d06aSjeremylt **/ 13969560d06aSjeremylt int CeedBasisReferenceCopy(CeedBasis basis, CeedBasis *basis_copy) { 1397356036faSJeremy L Thompson if (basis != CEED_BASIS_NONE) CeedCall(CeedBasisReference(basis)); 13982b730f8bSJeremy L Thompson CeedCall(CeedBasisDestroy(basis_copy)); 13999560d06aSjeremylt *basis_copy = basis; 14009560d06aSjeremylt return CEED_ERROR_SUCCESS; 14019560d06aSjeremylt } 14029560d06aSjeremylt 14039560d06aSjeremylt /** 1404ca94c3ddSJeremy L Thompson @brief View a `CeedBasis` 14057a982d89SJeremy L. Thompson 1406ca94c3ddSJeremy L Thompson @param[in] basis `CeedBasis` to view 1407ca94c3ddSJeremy L Thompson @param[in] stream Stream to view to, e.g., `stdout` 14087a982d89SJeremy L. Thompson 14097a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 14107a982d89SJeremy L. Thompson 14117a982d89SJeremy L. Thompson @ref User 14127a982d89SJeremy L. Thompson **/ 14137a982d89SJeremy L. Thompson int CeedBasisView(CeedBasis basis, FILE *stream) { 14141203703bSJeremy L Thompson bool is_tensor_basis; 14151203703bSJeremy L Thompson CeedElemTopology topo; 14161203703bSJeremy L Thompson CeedFESpace fe_space; 14171203703bSJeremy L Thompson 14181203703bSJeremy L Thompson // Basis data 14191203703bSJeremy L Thompson CeedCall(CeedBasisIsTensor(basis, &is_tensor_basis)); 14201203703bSJeremy L Thompson CeedCall(CeedBasisGetTopology(basis, &topo)); 14211203703bSJeremy L Thompson CeedCall(CeedBasisGetFESpace(basis, &fe_space)); 14222b730f8bSJeremy L Thompson 142350c301a5SRezgar Shakeri // Print FE space and element topology of the basis 1424edf04919SJeremy L Thompson fprintf(stream, "CeedBasis in a %s on a %s element\n", CeedFESpaces[fe_space], CeedElemTopologies[topo]); 14251203703bSJeremy L Thompson if (is_tensor_basis) { 1426edf04919SJeremy L Thompson fprintf(stream, " P: %" CeedInt_FMT "\n Q: %" CeedInt_FMT "\n", basis->P_1d, basis->Q_1d); 142750c301a5SRezgar Shakeri } else { 1428edf04919SJeremy L Thompson fprintf(stream, " P: %" CeedInt_FMT "\n Q: %" CeedInt_FMT "\n", basis->P, basis->Q); 142950c301a5SRezgar Shakeri } 1430edf04919SJeremy L Thompson fprintf(stream, " dimension: %" CeedInt_FMT "\n field components: %" CeedInt_FMT "\n", basis->dim, basis->num_comp); 1431ea61e9acSJeremy L Thompson // Print quadrature data, interpolation/gradient/divergence/curl of the basis 14321203703bSJeremy L Thompson if (is_tensor_basis) { // tensor basis 14331203703bSJeremy L Thompson CeedInt P_1d, Q_1d; 14341203703bSJeremy L Thompson const CeedScalar *q_ref_1d, *q_weight_1d, *interp_1d, *grad_1d; 14351203703bSJeremy L Thompson 14361203703bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d)); 14371203703bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 14381203703bSJeremy L Thompson CeedCall(CeedBasisGetQRef(basis, &q_ref_1d)); 14391203703bSJeremy L Thompson CeedCall(CeedBasisGetQWeights(basis, &q_weight_1d)); 14401203703bSJeremy L Thompson CeedCall(CeedBasisGetInterp1D(basis, &interp_1d)); 14411203703bSJeremy L Thompson CeedCall(CeedBasisGetGrad1D(basis, &grad_1d)); 14421203703bSJeremy L Thompson 14431203703bSJeremy L Thompson CeedCall(CeedScalarView("qref1d", "\t% 12.8f", 1, Q_1d, q_ref_1d, stream)); 14441203703bSJeremy L Thompson CeedCall(CeedScalarView("qweight1d", "\t% 12.8f", 1, Q_1d, q_weight_1d, stream)); 14451203703bSJeremy L Thompson CeedCall(CeedScalarView("interp1d", "\t% 12.8f", Q_1d, P_1d, interp_1d, stream)); 14461203703bSJeremy L Thompson CeedCall(CeedScalarView("grad1d", "\t% 12.8f", Q_1d, P_1d, grad_1d, stream)); 144750c301a5SRezgar Shakeri } else { // non-tensor basis 14481203703bSJeremy L Thompson CeedInt P, Q, dim, q_comp; 14491203703bSJeremy L Thompson const CeedScalar *q_ref, *q_weight, *interp, *grad, *div, *curl; 14501203703bSJeremy L Thompson 14511203703bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis, &P)); 14521203703bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints(basis, &Q)); 14531203703bSJeremy L Thompson CeedCall(CeedBasisGetDimension(basis, &dim)); 14541203703bSJeremy L Thompson CeedCall(CeedBasisGetQRef(basis, &q_ref)); 14551203703bSJeremy L Thompson CeedCall(CeedBasisGetQWeights(basis, &q_weight)); 14561203703bSJeremy L Thompson CeedCall(CeedBasisGetInterp(basis, &interp)); 14571203703bSJeremy L Thompson CeedCall(CeedBasisGetGrad(basis, &grad)); 14581203703bSJeremy L Thompson CeedCall(CeedBasisGetDiv(basis, &div)); 14591203703bSJeremy L Thompson CeedCall(CeedBasisGetCurl(basis, &curl)); 14601203703bSJeremy L Thompson 14611203703bSJeremy L Thompson CeedCall(CeedScalarView("qref", "\t% 12.8f", 1, Q * dim, q_ref, stream)); 14621203703bSJeremy L Thompson CeedCall(CeedScalarView("qweight", "\t% 12.8f", 1, Q, q_weight, stream)); 1463c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp)); 14641203703bSJeremy L Thompson CeedCall(CeedScalarView("interp", "\t% 12.8f", q_comp * Q, P, interp, stream)); 14651203703bSJeremy L Thompson if (grad) { 1466c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp)); 14671203703bSJeremy L Thompson CeedCall(CeedScalarView("grad", "\t% 12.8f", q_comp * Q, P, grad, stream)); 14687a982d89SJeremy L. Thompson } 14691203703bSJeremy L Thompson if (div) { 1470c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_DIV, &q_comp)); 14711203703bSJeremy L Thompson CeedCall(CeedScalarView("div", "\t% 12.8f", q_comp * Q, P, div, stream)); 1472c4e3f59bSSebastian Grimberg } 14731203703bSJeremy L Thompson if (curl) { 1474c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_CURL, &q_comp)); 14751203703bSJeremy L Thompson CeedCall(CeedScalarView("curl", "\t% 12.8f", q_comp * Q, P, curl, stream)); 147650c301a5SRezgar Shakeri } 147750c301a5SRezgar Shakeri } 1478e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 14797a982d89SJeremy L. Thompson } 14807a982d89SJeremy L. Thompson 14817a982d89SJeremy L. Thompson /** 14827a982d89SJeremy L. Thompson @brief Apply basis evaluation from nodes to quadrature points or vice versa 14837a982d89SJeremy L. Thompson 1484ca94c3ddSJeremy L Thompson @param[in] basis `CeedBasis` to evaluate 1485ea61e9acSJeremy L Thompson @param[in] num_elem The number of elements to apply the basis evaluation to; 1486ca94c3ddSJeremy L Thompson the backend will specify the ordering in @ref CeedElemRestrictionCreate() 1487ca94c3ddSJeremy L Thompson @param[in] t_mode @ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature points; 1488ca94c3ddSJeremy L Thompson @ref CEED_TRANSPOSE to apply the transpose, mapping from quadrature points to nodes 1489ca94c3ddSJeremy L Thompson @param[in] eval_mode @ref CEED_EVAL_NONE to use values directly, 1490ca94c3ddSJeremy L Thompson @ref CEED_EVAL_INTERP to use interpolated values, 1491ca94c3ddSJeremy L Thompson @ref CEED_EVAL_GRAD to use gradients, 1492ca94c3ddSJeremy L Thompson @ref CEED_EVAL_DIV to use divergence, 1493ca94c3ddSJeremy L Thompson @ref CEED_EVAL_CURL to use curl, 1494ca94c3ddSJeremy L Thompson @ref CEED_EVAL_WEIGHT to use quadrature weights 1495ca94c3ddSJeremy L Thompson @param[in] u Input `CeedVector` 1496ca94c3ddSJeremy L Thompson @param[out] v Output `CeedVector` 14977a982d89SJeremy L. Thompson 14987a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 14997a982d89SJeremy L. Thompson 15007a982d89SJeremy L. Thompson @ref User 15017a982d89SJeremy L. Thompson **/ 15022b730f8bSJeremy L Thompson int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, CeedVector v) { 1503c4e3f59bSSebastian Grimberg CeedInt dim, num_comp, q_comp, num_nodes, num_qpts; 15041c66c397SJeremy L Thompson CeedSize u_length = 0, v_length; 15051203703bSJeremy L Thompson Ceed ceed; 15061c66c397SJeremy L Thompson 15071203703bSJeremy L Thompson CeedCall(CeedBasisGetCeed(basis, &ceed)); 15082b730f8bSJeremy L Thompson CeedCall(CeedBasisGetDimension(basis, &dim)); 15092b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumComponents(basis, &num_comp)); 1510c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp)); 15112b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis, &num_nodes)); 15122b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts)); 15132b730f8bSJeremy L Thompson CeedCall(CeedVectorGetLength(v, &v_length)); 1514c8c3fa7dSJeremy L Thompson if (u) CeedCall(CeedVectorGetLength(u, &u_length)); 15157a982d89SJeremy L. Thompson 15161203703bSJeremy L Thompson CeedCheck(basis->Apply, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support CeedBasisApply"); 1517e15f9bd0SJeremy L Thompson 1518e15f9bd0SJeremy L Thompson // Check compatibility of topological and geometrical dimensions 15196574a04fSJeremy L Thompson CeedCheck((t_mode == CEED_TRANSPOSE && v_length % num_nodes == 0 && u_length % num_qpts == 0) || 15206574a04fSJeremy L Thompson (t_mode == CEED_NOTRANSPOSE && u_length % num_nodes == 0 && v_length % num_qpts == 0), 15211203703bSJeremy L Thompson ceed, CEED_ERROR_DIMENSION, "Length of input/output vectors incompatible with basis dimensions"); 15227a982d89SJeremy L. Thompson 1523e15f9bd0SJeremy L Thompson // Check vector lengths to prevent out of bounds issues 152499e754f0SJeremy L Thompson bool has_good_dims = true; 1525d1d35e2fSjeremylt switch (eval_mode) { 1526e15f9bd0SJeremy L Thompson case CEED_EVAL_NONE: 15272b730f8bSJeremy L Thompson case CEED_EVAL_INTERP: 15282b730f8bSJeremy L Thompson case CEED_EVAL_GRAD: 1529c4e3f59bSSebastian Grimberg case CEED_EVAL_DIV: 1530c4e3f59bSSebastian Grimberg case CEED_EVAL_CURL: 153199e754f0SJeremy L Thompson has_good_dims = 15326574a04fSJeremy L Thompson ((t_mode == CEED_TRANSPOSE && u_length >= num_elem * num_comp * num_qpts * q_comp && v_length >= num_elem * num_comp * num_nodes) || 15336574a04fSJeremy L Thompson (t_mode == CEED_NOTRANSPOSE && v_length >= num_elem * num_qpts * num_comp * q_comp && u_length >= num_elem * num_comp * num_nodes)); 1534e15f9bd0SJeremy L Thompson break; 1535e15f9bd0SJeremy L Thompson case CEED_EVAL_WEIGHT: 153699e754f0SJeremy L Thompson has_good_dims = v_length >= num_elem * num_qpts; 1537e15f9bd0SJeremy L Thompson break; 1538e15f9bd0SJeremy L Thompson } 153999e754f0SJeremy L Thompson CeedCheck(has_good_dims, ceed, CEED_ERROR_DIMENSION, "Input/output vectors too short for basis and evaluation mode"); 1540e15f9bd0SJeremy L Thompson 15412b730f8bSJeremy L Thompson CeedCall(basis->Apply(basis, num_elem, t_mode, eval_mode, u, v)); 1542e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 15437a982d89SJeremy L. Thompson } 15447a982d89SJeremy L. Thompson 15457a982d89SJeremy L. Thompson /** 1546c8c3fa7dSJeremy L Thompson @brief Apply basis evaluation from nodes to arbitrary points 1547c8c3fa7dSJeremy L Thompson 1548ca94c3ddSJeremy L Thompson @param[in] basis `CeedBasis` to evaluate 1549fc0f7cc6SJeremy L Thompson @param[in] num_elem The number of elements to apply the basis evaluation to; 1550fc0f7cc6SJeremy L Thompson the backend will specify the ordering in @ref CeedElemRestrictionCreate() 1551faed4840SJeremy L Thompson @param[in] num_points Array of the number of points to apply the basis evaluation to in each element, size `num_elem` 1552ca94c3ddSJeremy L Thompson @param[in] t_mode @ref CEED_NOTRANSPOSE to evaluate from nodes to points; 1553ca94c3ddSJeremy L Thompson @ref CEED_TRANSPOSE to apply the transpose, mapping from points to nodes 1554ca94c3ddSJeremy L Thompson @param[in] eval_mode @ref CEED_EVAL_INTERP to use interpolated values, 1555ca94c3ddSJeremy L Thompson @ref CEED_EVAL_GRAD to use gradients, 1556ca94c3ddSJeremy L Thompson @ref CEED_EVAL_WEIGHT to use quadrature weights 1557ca94c3ddSJeremy L Thompson @param[in] x_ref `CeedVector` holding reference coordinates of each point 1558ca94c3ddSJeremy L Thompson @param[in] u Input `CeedVector`, of length `num_nodes * num_comp` for @ref CEED_NOTRANSPOSE 1559ca94c3ddSJeremy L Thompson @param[out] v Output `CeedVector`, of length `num_points * num_q_comp` for @ref CEED_NOTRANSPOSE with @ref CEED_EVAL_INTERP 1560c8c3fa7dSJeremy L Thompson 1561c8c3fa7dSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 1562c8c3fa7dSJeremy L Thompson 1563c8c3fa7dSJeremy L Thompson @ref User 1564c8c3fa7dSJeremy L Thompson **/ 1565fc0f7cc6SJeremy L Thompson int CeedBasisApplyAtPoints(CeedBasis basis, CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode, CeedEvalMode eval_mode, 1566fc0f7cc6SJeremy L Thompson CeedVector x_ref, CeedVector u, CeedVector v) { 15671203703bSJeremy L Thompson bool is_tensor_basis; 1568fc0f7cc6SJeremy L Thompson CeedInt dim, num_comp, num_q_comp, num_nodes, P_1d = 1, Q_1d = 1, total_num_points = 0; 15691c66c397SJeremy L Thompson CeedSize x_length = 0, u_length = 0, v_length; 15701203703bSJeremy L Thompson Ceed ceed; 1571c8c3fa7dSJeremy L Thompson 15721203703bSJeremy L Thompson CeedCall(CeedBasisGetCeed(basis, &ceed)); 1573c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisGetDimension(basis, &dim)); 1574c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d)); 1575c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 1576c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisGetNumComponents(basis, &num_comp)); 1577c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &num_q_comp)); 1578c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis, &num_nodes)); 1579c8c3fa7dSJeremy L Thompson CeedCall(CeedVectorGetLength(v, &v_length)); 1580953190f4SJeremy L Thompson if (x_ref != CEED_VECTOR_NONE) CeedCall(CeedVectorGetLength(x_ref, &x_length)); 1581953190f4SJeremy L Thompson if (u != CEED_VECTOR_NONE) CeedCall(CeedVectorGetLength(u, &u_length)); 1582c8c3fa7dSJeremy L Thompson 1583c8c3fa7dSJeremy L Thompson // Check compatibility of topological and geometrical dimensions 1584fc0f7cc6SJeremy L Thompson for (CeedInt i = 0; i < num_elem; i++) total_num_points += num_points[i]; 1585953190f4SJeremy L Thompson CeedCheck((t_mode == CEED_TRANSPOSE && v_length % num_nodes == 0) || (t_mode == CEED_NOTRANSPOSE && u_length % num_nodes == 0) || 1586953190f4SJeremy L Thompson (eval_mode == CEED_EVAL_WEIGHT), 15871203703bSJeremy L Thompson ceed, CEED_ERROR_DIMENSION, "Length of input/output vectors incompatible with basis dimensions and number of points"); 1588c8c3fa7dSJeremy L Thompson 1589c8c3fa7dSJeremy L Thompson // Check compatibility coordinates vector 1590fc0f7cc6SJeremy L Thompson CeedCheck((x_length >= total_num_points * dim) || (eval_mode == CEED_EVAL_WEIGHT), ceed, CEED_ERROR_DIMENSION, 1591c8c3fa7dSJeremy L Thompson "Length of reference coordinate vector incompatible with basis dimension and number of points"); 1592c8c3fa7dSJeremy L Thompson 1593953190f4SJeremy L Thompson // Check CEED_EVAL_WEIGHT only on CEED_NOTRANSPOSE 15941203703bSJeremy L Thompson CeedCheck(eval_mode != CEED_EVAL_WEIGHT || t_mode == CEED_NOTRANSPOSE, ceed, CEED_ERROR_UNSUPPORTED, 1595953190f4SJeremy L Thompson "CEED_EVAL_WEIGHT only supported with CEED_NOTRANSPOSE"); 1596953190f4SJeremy L Thompson 1597c8c3fa7dSJeremy L Thompson // Check vector lengths to prevent out of bounds issues 159899e754f0SJeremy L Thompson bool has_good_dims = true; 1599c8c3fa7dSJeremy L Thompson switch (eval_mode) { 1600c8c3fa7dSJeremy L Thompson case CEED_EVAL_INTERP: 1601fc0f7cc6SJeremy L Thompson has_good_dims = ((t_mode == CEED_TRANSPOSE && (u_length >= total_num_points * num_q_comp || v_length >= num_elem * num_nodes * num_comp)) || 1602fc0f7cc6SJeremy L Thompson (t_mode == CEED_NOTRANSPOSE && (v_length >= total_num_points * num_q_comp || u_length >= num_elem * num_nodes * num_comp))); 1603c8c3fa7dSJeremy L Thompson break; 1604c8c3fa7dSJeremy L Thompson case CEED_EVAL_GRAD: 1605fc0f7cc6SJeremy L Thompson has_good_dims = 1606fc0f7cc6SJeremy L Thompson ((t_mode == CEED_TRANSPOSE && (u_length >= total_num_points * num_q_comp * dim || v_length >= num_elem * num_nodes * num_comp)) || 1607fc0f7cc6SJeremy L Thompson (t_mode == CEED_NOTRANSPOSE && (v_length >= total_num_points * num_q_comp * dim || u_length >= num_elem * num_nodes * num_comp))); 1608edfbf3a6SJeremy L Thompson break; 1609c8c3fa7dSJeremy L Thompson case CEED_EVAL_WEIGHT: 1610fc0f7cc6SJeremy L Thompson has_good_dims = t_mode == CEED_NOTRANSPOSE && (v_length >= total_num_points); 1611953190f4SJeremy L Thompson break; 161299e754f0SJeremy L Thompson // LCOV_EXCL_START 1613953190f4SJeremy L Thompson case CEED_EVAL_NONE: 1614c8c3fa7dSJeremy L Thompson case CEED_EVAL_DIV: 1615c8c3fa7dSJeremy L Thompson case CEED_EVAL_CURL: 16161203703bSJeremy L Thompson return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Evaluation at arbitrary points not supported for %s", CeedEvalModes[eval_mode]); 1617c8c3fa7dSJeremy L Thompson // LCOV_EXCL_STOP 1618c8c3fa7dSJeremy L Thompson } 161999e754f0SJeremy L Thompson CeedCheck(has_good_dims, ceed, CEED_ERROR_DIMENSION, "Input/output vectors too short for basis and evaluation mode"); 1620c8c3fa7dSJeremy L Thompson 1621c8c3fa7dSJeremy L Thompson // Backend method 1622c8c3fa7dSJeremy L Thompson if (basis->ApplyAtPoints) { 1623fc0f7cc6SJeremy L Thompson CeedCall(basis->ApplyAtPoints(basis, num_elem, num_points, t_mode, eval_mode, x_ref, u, v)); 1624c8c3fa7dSJeremy L Thompson return CEED_ERROR_SUCCESS; 1625c8c3fa7dSJeremy L Thompson } 1626c8c3fa7dSJeremy L Thompson 1627c8c3fa7dSJeremy L Thompson // Default implementation 16281203703bSJeremy L Thompson CeedCall(CeedBasisIsTensor(basis, &is_tensor_basis)); 16291203703bSJeremy L Thompson CeedCheck(is_tensor_basis, ceed, CEED_ERROR_UNSUPPORTED, "Evaluation at arbitrary points only supported for tensor product bases"); 1630fc0f7cc6SJeremy L Thompson CeedCheck(num_elem == 1, ceed, CEED_ERROR_UNSUPPORTED, "Evaluation at arbitrary points only supported for a single element at a time"); 1631953190f4SJeremy L Thompson if (eval_mode == CEED_EVAL_WEIGHT) { 1632953190f4SJeremy L Thompson CeedCall(CeedVectorSetValue(v, 1.0)); 1633953190f4SJeremy L Thompson return CEED_ERROR_SUCCESS; 1634953190f4SJeremy L Thompson } 1635c8c3fa7dSJeremy L Thompson if (!basis->basis_chebyshev) { 1636c8c3fa7dSJeremy L Thompson // Build basis mapping from nodes to Chebyshev coefficients 1637c8c3fa7dSJeremy L Thompson CeedScalar *chebyshev_interp_1d, *chebyshev_grad_1d, *chebyshev_q_weight_1d; 1638b0cc4569SJeremy L Thompson const CeedScalar *q_ref_1d; 1639c8c3fa7dSJeremy L Thompson 164071a83b88SJeremy L Thompson CeedCall(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d)); 164171a83b88SJeremy L Thompson CeedCall(CeedCalloc(P_1d * Q_1d, &chebyshev_grad_1d)); 1642c8c3fa7dSJeremy L Thompson CeedCall(CeedCalloc(Q_1d, &chebyshev_q_weight_1d)); 1643b0cc4569SJeremy L Thompson CeedCall(CeedBasisGetQRef(basis, &q_ref_1d)); 1644b0cc4569SJeremy L Thompson CeedCall(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d)); 1645c8c3fa7dSJeremy L Thompson 16461203703bSJeremy L Thompson CeedCall(CeedVectorCreate(ceed, num_comp * CeedIntPow(Q_1d, dim), &basis->vec_chebyshev)); 16471203703bSJeremy L Thompson CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P_1d, Q_1d, chebyshev_interp_1d, chebyshev_grad_1d, q_ref_1d, chebyshev_q_weight_1d, 1648c8c3fa7dSJeremy L Thompson &basis->basis_chebyshev)); 1649c8c3fa7dSJeremy L Thompson 1650c8c3fa7dSJeremy L Thompson // Cleanup 1651c8c3fa7dSJeremy L Thompson CeedCall(CeedFree(&chebyshev_interp_1d)); 1652c8c3fa7dSJeremy L Thompson CeedCall(CeedFree(&chebyshev_grad_1d)); 1653c8c3fa7dSJeremy L Thompson CeedCall(CeedFree(&chebyshev_q_weight_1d)); 1654c8c3fa7dSJeremy L Thompson } 1655c8c3fa7dSJeremy L Thompson 1656c8c3fa7dSJeremy L Thompson // Create TensorContract object if needed, such as a basis from the GPU backends 1657c8c3fa7dSJeremy L Thompson if (!basis->contract) { 1658c8c3fa7dSJeremy L Thompson Ceed ceed_ref; 1659585a562dSJeremy L Thompson CeedBasis basis_ref = NULL; 1660c8c3fa7dSJeremy L Thompson 1661c8c3fa7dSJeremy L Thompson CeedCall(CeedInit("/cpu/self", &ceed_ref)); 1662c8c3fa7dSJeremy L Thompson // Only need matching tensor contraction dimensions, any type of basis will work 166371a83b88SJeremy L Thompson CeedCall(CeedBasisCreateTensorH1Lagrange(ceed_ref, dim, num_comp, P_1d, Q_1d, CEED_GAUSS, &basis_ref)); 1664585a562dSJeremy L Thompson // Note - clang-tidy doesn't know basis_ref->contract must be valid here 16651203703bSJeremy L Thompson CeedCheck(basis_ref && basis_ref->contract, ceed, CEED_ERROR_UNSUPPORTED, "Reference CPU ceed failed to create a tensor contraction object"); 1666585a562dSJeremy L Thompson CeedCall(CeedTensorContractReferenceCopy(basis_ref->contract, &basis->contract)); 1667c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisDestroy(&basis_ref)); 1668c8c3fa7dSJeremy L Thompson CeedCall(CeedDestroy(&ceed_ref)); 1669c8c3fa7dSJeremy L Thompson } 1670c8c3fa7dSJeremy L Thompson 1671c8c3fa7dSJeremy L Thompson // Basis evaluation 1672c8c3fa7dSJeremy L Thompson switch (t_mode) { 1673c8c3fa7dSJeremy L Thompson case CEED_NOTRANSPOSE: { 1674c8c3fa7dSJeremy L Thompson // Nodes to arbitrary points 1675c8c3fa7dSJeremy L Thompson CeedScalar *v_array; 1676c8c3fa7dSJeremy L Thompson const CeedScalar *chebyshev_coeffs, *x_array_read; 1677c8c3fa7dSJeremy L Thompson 1678c8c3fa7dSJeremy L Thompson // -- Interpolate to Chebyshev coefficients 1679c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisApply(basis->basis_chebyshev, 1, CEED_NOTRANSPOSE, CEED_EVAL_INTERP, u, basis->vec_chebyshev)); 1680c8c3fa7dSJeremy L Thompson 1681c8c3fa7dSJeremy L Thompson // -- Evaluate Chebyshev polynomials at arbitrary points 1682c8c3fa7dSJeremy L Thompson CeedCall(CeedVectorGetArrayRead(basis->vec_chebyshev, CEED_MEM_HOST, &chebyshev_coeffs)); 1683c8c3fa7dSJeremy L Thompson CeedCall(CeedVectorGetArrayRead(x_ref, CEED_MEM_HOST, &x_array_read)); 1684c8c3fa7dSJeremy L Thompson CeedCall(CeedVectorGetArrayWrite(v, CEED_MEM_HOST, &v_array)); 1685edfbf3a6SJeremy L Thompson switch (eval_mode) { 1686edfbf3a6SJeremy L Thompson case CEED_EVAL_INTERP: { 1687c8c3fa7dSJeremy L Thompson CeedScalar tmp[2][num_comp * CeedIntPow(Q_1d, dim)], chebyshev_x[Q_1d]; 1688c8c3fa7dSJeremy L Thompson 1689c8c3fa7dSJeremy L Thompson // ---- Values at point 1690fc0f7cc6SJeremy L Thompson for (CeedInt p = 0; p < total_num_points; p++) { 1691c8c3fa7dSJeremy L Thompson CeedInt pre = num_comp * CeedIntPow(Q_1d, dim - 1), post = 1; 1692c8c3fa7dSJeremy L Thompson 169353ef2869SZach Atkins for (CeedInt d = 0; d < dim; d++) { 16943778dbaaSJeremy L Thompson // ------ Tensor contract with current Chebyshev polynomial values 1695fc0f7cc6SJeremy L Thompson CeedCall(CeedChebyshevPolynomialsAtPoint(x_array_read[d * total_num_points + p], Q_1d, chebyshev_x)); 1696c8c3fa7dSJeremy L Thompson CeedCall(CeedTensorContractApply(basis->contract, pre, Q_1d, post, 1, chebyshev_x, t_mode, false, 16974608bdaaSJeremy L Thompson d == 0 ? chebyshev_coeffs : tmp[d % 2], tmp[(d + 1) % 2])); 1698c8c3fa7dSJeremy L Thompson pre /= Q_1d; 1699c8c3fa7dSJeremy L Thompson post *= 1; 1700c8c3fa7dSJeremy L Thompson } 1701fc0f7cc6SJeremy L Thompson for (CeedInt c = 0; c < num_comp; c++) v_array[c * total_num_points + p] = tmp[dim % 2][c]; 1702c8c3fa7dSJeremy L Thompson } 1703edfbf3a6SJeremy L Thompson break; 1704edfbf3a6SJeremy L Thompson } 1705edfbf3a6SJeremy L Thompson case CEED_EVAL_GRAD: { 1706edfbf3a6SJeremy L Thompson CeedScalar tmp[2][num_comp * CeedIntPow(Q_1d, dim)], chebyshev_x[Q_1d]; 1707edfbf3a6SJeremy L Thompson 1708edfbf3a6SJeremy L Thompson // ---- Values at point 1709fc0f7cc6SJeremy L Thompson for (CeedInt p = 0; p < total_num_points; p++) { 1710edfbf3a6SJeremy L Thompson // Dim**2 contractions, apply grad when pass == dim 171153ef2869SZach Atkins for (CeedInt pass = 0; pass < dim; pass++) { 1712edfbf3a6SJeremy L Thompson CeedInt pre = num_comp * CeedIntPow(Q_1d, dim - 1), post = 1; 1713edfbf3a6SJeremy L Thompson 171453ef2869SZach Atkins for (CeedInt d = 0; d < dim; d++) { 17153778dbaaSJeremy L Thompson // ------ Tensor contract with current Chebyshev polynomial values 1716fc0f7cc6SJeremy L Thompson if (pass == d) CeedCall(CeedChebyshevDerivativeAtPoint(x_array_read[d * total_num_points + p], Q_1d, chebyshev_x)); 1717fc0f7cc6SJeremy L Thompson else CeedCall(CeedChebyshevPolynomialsAtPoint(x_array_read[d * total_num_points + p], Q_1d, chebyshev_x)); 1718edfbf3a6SJeremy L Thompson CeedCall(CeedTensorContractApply(basis->contract, pre, Q_1d, post, 1, chebyshev_x, t_mode, false, 17194608bdaaSJeremy L Thompson d == 0 ? chebyshev_coeffs : tmp[d % 2], tmp[(d + 1) % 2])); 1720edfbf3a6SJeremy L Thompson pre /= Q_1d; 1721edfbf3a6SJeremy L Thompson post *= 1; 1722edfbf3a6SJeremy L Thompson } 1723fc0f7cc6SJeremy L Thompson for (CeedInt c = 0; c < num_comp; c++) v_array[(pass * num_comp + c) * total_num_points + p] = tmp[dim % 2][c]; 1724edfbf3a6SJeremy L Thompson } 1725edfbf3a6SJeremy L Thompson } 1726edfbf3a6SJeremy L Thompson break; 1727edfbf3a6SJeremy L Thompson } 1728edfbf3a6SJeremy L Thompson default: 1729953190f4SJeremy L Thompson // Nothing to do, excluded above 1730edfbf3a6SJeremy L Thompson break; 1731c8c3fa7dSJeremy L Thompson } 1732c8c3fa7dSJeremy L Thompson CeedCall(CeedVectorRestoreArrayRead(basis->vec_chebyshev, &chebyshev_coeffs)); 1733c8c3fa7dSJeremy L Thompson CeedCall(CeedVectorRestoreArrayRead(x_ref, &x_array_read)); 1734c8c3fa7dSJeremy L Thompson CeedCall(CeedVectorRestoreArray(v, &v_array)); 1735c8c3fa7dSJeremy L Thompson break; 1736c8c3fa7dSJeremy L Thompson } 17372a94f45fSJeremy L Thompson case CEED_TRANSPOSE: { 17383778dbaaSJeremy L Thompson // Note: No switch on e_mode here because only CEED_EVAL_INTERP is supported at this time 17392a94f45fSJeremy L Thompson // Arbitrary points to nodes 17402a94f45fSJeremy L Thompson CeedScalar *chebyshev_coeffs; 17412a94f45fSJeremy L Thompson const CeedScalar *u_array, *x_array_read; 17422a94f45fSJeremy L Thompson 17431c66c397SJeremy L Thompson // -- Transpose of evaluation of Chebyshev polynomials at arbitrary points 17442a94f45fSJeremy L Thompson CeedCall(CeedVectorGetArrayWrite(basis->vec_chebyshev, CEED_MEM_HOST, &chebyshev_coeffs)); 17452a94f45fSJeremy L Thompson CeedCall(CeedVectorGetArrayRead(x_ref, CEED_MEM_HOST, &x_array_read)); 17462a94f45fSJeremy L Thompson CeedCall(CeedVectorGetArrayRead(u, CEED_MEM_HOST, &u_array)); 1747038a8942SZach Atkins 1748038a8942SZach Atkins switch (eval_mode) { 1749038a8942SZach Atkins case CEED_EVAL_INTERP: { 17502a94f45fSJeremy L Thompson CeedScalar tmp[2][num_comp * CeedIntPow(Q_1d, dim)], chebyshev_x[Q_1d]; 17512a94f45fSJeremy L Thompson 17522a94f45fSJeremy L Thompson // ---- Values at point 1753fc0f7cc6SJeremy L Thompson for (CeedInt p = 0; p < total_num_points; p++) { 17542a94f45fSJeremy L Thompson CeedInt pre = num_comp * 1, post = 1; 17552a94f45fSJeremy L Thompson 1756fc0f7cc6SJeremy L Thompson for (CeedInt c = 0; c < num_comp; c++) tmp[0][c] = u_array[c * total_num_points + p]; 175753ef2869SZach Atkins for (CeedInt d = 0; d < dim; d++) { 17583778dbaaSJeremy L Thompson // ------ Tensor contract with current Chebyshev polynomial values 1759fc0f7cc6SJeremy L Thompson CeedCall(CeedChebyshevPolynomialsAtPoint(x_array_read[d * total_num_points + p], Q_1d, chebyshev_x)); 17604608bdaaSJeremy L Thompson CeedCall(CeedTensorContractApply(basis->contract, pre, 1, post, Q_1d, chebyshev_x, t_mode, p > 0 && d == (dim - 1), tmp[d % 2], 17614608bdaaSJeremy L Thompson d == (dim - 1) ? chebyshev_coeffs : tmp[(d + 1) % 2])); 17622a94f45fSJeremy L Thompson pre /= 1; 17632a94f45fSJeremy L Thompson post *= Q_1d; 17642a94f45fSJeremy L Thompson } 17652a94f45fSJeremy L Thompson } 1766038a8942SZach Atkins break; 1767038a8942SZach Atkins } 1768038a8942SZach Atkins case CEED_EVAL_GRAD: { 1769038a8942SZach Atkins CeedScalar tmp[2][num_comp * CeedIntPow(Q_1d, dim)], chebyshev_x[Q_1d]; 1770038a8942SZach Atkins 1771038a8942SZach Atkins // ---- Values at point 1772fc0f7cc6SJeremy L Thompson for (CeedInt p = 0; p < total_num_points; p++) { 1773038a8942SZach Atkins // Dim**2 contractions, apply grad when pass == dim 1774038a8942SZach Atkins for (CeedInt pass = 0; pass < dim; pass++) { 1775038a8942SZach Atkins CeedInt pre = num_comp * 1, post = 1; 1776038a8942SZach Atkins 1777fc0f7cc6SJeremy L Thompson for (CeedInt c = 0; c < num_comp; c++) tmp[0][c] = u_array[(pass * num_comp + c) * total_num_points + p]; 1778038a8942SZach Atkins for (CeedInt d = 0; d < dim; d++) { 1779038a8942SZach Atkins // ------ Tensor contract with current Chebyshev polynomial values 1780fc0f7cc6SJeremy L Thompson if (pass == d) CeedCall(CeedChebyshevDerivativeAtPoint(x_array_read[d * total_num_points + p], Q_1d, chebyshev_x)); 1781fc0f7cc6SJeremy L Thompson else CeedCall(CeedChebyshevPolynomialsAtPoint(x_array_read[d * total_num_points + p], Q_1d, chebyshev_x)); 17824608bdaaSJeremy L Thompson CeedCall(CeedTensorContractApply(basis->contract, pre, 1, post, Q_1d, chebyshev_x, t_mode, 17834608bdaaSJeremy L Thompson (p > 0 || (p == 0 && pass > 0)) && d == (dim - 1), tmp[d % 2], 17844608bdaaSJeremy L Thompson d == (dim - 1) ? chebyshev_coeffs : tmp[(d + 1) % 2])); 1785038a8942SZach Atkins pre /= 1; 1786038a8942SZach Atkins post *= Q_1d; 1787038a8942SZach Atkins } 1788038a8942SZach Atkins } 1789038a8942SZach Atkins } 1790038a8942SZach Atkins break; 1791038a8942SZach Atkins } 1792038a8942SZach Atkins default: 1793038a8942SZach Atkins // Nothing to do, excluded above 1794038a8942SZach Atkins break; 17952a94f45fSJeremy L Thompson } 17962a94f45fSJeremy L Thompson CeedCall(CeedVectorRestoreArray(basis->vec_chebyshev, &chebyshev_coeffs)); 17972a94f45fSJeremy L Thompson CeedCall(CeedVectorRestoreArrayRead(x_ref, &x_array_read)); 17982a94f45fSJeremy L Thompson CeedCall(CeedVectorRestoreArrayRead(u, &u_array)); 17992a94f45fSJeremy L Thompson 18002a94f45fSJeremy L Thompson // -- Interpolate transpose from Chebyshev coefficients 18012a94f45fSJeremy L Thompson CeedCall(CeedBasisApply(basis->basis_chebyshev, 1, CEED_TRANSPOSE, CEED_EVAL_INTERP, basis->vec_chebyshev, v)); 18022a94f45fSJeremy L Thompson break; 18032a94f45fSJeremy L Thompson } 1804c8c3fa7dSJeremy L Thompson } 1805c8c3fa7dSJeremy L Thompson return CEED_ERROR_SUCCESS; 1806c8c3fa7dSJeremy L Thompson } 1807c8c3fa7dSJeremy L Thompson 1808c8c3fa7dSJeremy L Thompson /** 18096e536b99SJeremy L Thompson @brief Get the `Ceed` associated with a `CeedBasis` 1810b7c9bbdaSJeremy L Thompson 1811ca94c3ddSJeremy L Thompson @param[in] basis `CeedBasis` 1812ca94c3ddSJeremy L Thompson @param[out] ceed Variable to store `Ceed` 1813b7c9bbdaSJeremy L Thompson 1814b7c9bbdaSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 1815b7c9bbdaSJeremy L Thompson 1816b7c9bbdaSJeremy L Thompson @ref Advanced 1817b7c9bbdaSJeremy L Thompson **/ 1818b7c9bbdaSJeremy L Thompson int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) { 18196e536b99SJeremy L Thompson *ceed = CeedBasisReturnCeed(basis); 1820b7c9bbdaSJeremy L Thompson return CEED_ERROR_SUCCESS; 1821b7c9bbdaSJeremy L Thompson } 1822b7c9bbdaSJeremy L Thompson 1823b7c9bbdaSJeremy L Thompson /** 18246e536b99SJeremy L Thompson @brief Return the `Ceed` associated with a `CeedBasis` 18256e536b99SJeremy L Thompson 18266e536b99SJeremy L Thompson @param[in] basis `CeedBasis` 18276e536b99SJeremy L Thompson 18286e536b99SJeremy L Thompson @return `Ceed` associated with the `basis` 18296e536b99SJeremy L Thompson 18306e536b99SJeremy L Thompson @ref Advanced 18316e536b99SJeremy L Thompson **/ 18326e536b99SJeremy L Thompson Ceed CeedBasisReturnCeed(CeedBasis basis) { return basis->ceed; } 18336e536b99SJeremy L Thompson 18346e536b99SJeremy L Thompson /** 1835ca94c3ddSJeremy L Thompson @brief Get dimension for given `CeedBasis` 18369d007619Sjeremylt 1837ca94c3ddSJeremy L Thompson @param[in] basis `CeedBasis` 18389d007619Sjeremylt @param[out] dim Variable to store dimension of basis 18399d007619Sjeremylt 18409d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 18419d007619Sjeremylt 1842b7c9bbdaSJeremy L Thompson @ref Advanced 18439d007619Sjeremylt **/ 18449d007619Sjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) { 18459d007619Sjeremylt *dim = basis->dim; 1846e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 18479d007619Sjeremylt } 18489d007619Sjeremylt 18499d007619Sjeremylt /** 1850ca94c3ddSJeremy L Thompson @brief Get topology for given `CeedBasis` 1851d99fa3c5SJeremy L Thompson 1852ca94c3ddSJeremy L Thompson @param[in] basis `CeedBasis` 1853d99fa3c5SJeremy L Thompson @param[out] topo Variable to store topology of basis 1854d99fa3c5SJeremy L Thompson 1855d99fa3c5SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 1856d99fa3c5SJeremy L Thompson 1857b7c9bbdaSJeremy L Thompson @ref Advanced 1858d99fa3c5SJeremy L Thompson **/ 1859d99fa3c5SJeremy L Thompson int CeedBasisGetTopology(CeedBasis basis, CeedElemTopology *topo) { 1860d99fa3c5SJeremy L Thompson *topo = basis->topo; 1861e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1862d99fa3c5SJeremy L Thompson } 1863d99fa3c5SJeremy L Thompson 1864d99fa3c5SJeremy L Thompson /** 1865ca94c3ddSJeremy L Thompson @brief Get number of components for given `CeedBasis` 18669d007619Sjeremylt 1867ca94c3ddSJeremy L Thompson @param[in] basis `CeedBasis` 1868ca94c3ddSJeremy L Thompson @param[out] num_comp Variable to store number of components 18699d007619Sjeremylt 18709d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 18719d007619Sjeremylt 1872b7c9bbdaSJeremy L Thompson @ref Advanced 18739d007619Sjeremylt **/ 1874d1d35e2fSjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *num_comp) { 1875d1d35e2fSjeremylt *num_comp = basis->num_comp; 1876e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 18779d007619Sjeremylt } 18789d007619Sjeremylt 18799d007619Sjeremylt /** 1880ca94c3ddSJeremy L Thompson @brief Get total number of nodes (in `dim` dimensions) of a `CeedBasis` 18819d007619Sjeremylt 1882ca94c3ddSJeremy L Thompson @param[in] basis `CeedBasis` 18839d007619Sjeremylt @param[out] P Variable to store number of nodes 18849d007619Sjeremylt 18859d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 18869d007619Sjeremylt 18879d007619Sjeremylt @ref Utility 18889d007619Sjeremylt **/ 18899d007619Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) { 18909d007619Sjeremylt *P = basis->P; 1891e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 18929d007619Sjeremylt } 18939d007619Sjeremylt 18949d007619Sjeremylt /** 1895ca94c3ddSJeremy L Thompson @brief Get total number of nodes (in 1 dimension) of a `CeedBasis` 18969d007619Sjeremylt 1897ca94c3ddSJeremy L Thompson @param[in] basis `CeedBasis` 1898d1d35e2fSjeremylt @param[out] P_1d Variable to store number of nodes 18999d007619Sjeremylt 19009d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 19019d007619Sjeremylt 1902b7c9bbdaSJeremy L Thompson @ref Advanced 19039d007619Sjeremylt **/ 1904d1d35e2fSjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P_1d) { 19056e536b99SJeremy L Thompson CeedCheck(basis->is_tensor_basis, CeedBasisReturnCeed(basis), CEED_ERROR_MINOR, "Cannot supply P_1d for non-tensor CeedBasis"); 1906d1d35e2fSjeremylt *P_1d = basis->P_1d; 1907e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 19089d007619Sjeremylt } 19099d007619Sjeremylt 19109d007619Sjeremylt /** 1911ca94c3ddSJeremy L Thompson @brief Get total number of quadrature points (in `dim` dimensions) of a `CeedBasis` 19129d007619Sjeremylt 1913ca94c3ddSJeremy L Thompson @param[in] basis `CeedBasis` 19149d007619Sjeremylt @param[out] Q Variable to store number of quadrature points 19159d007619Sjeremylt 19169d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 19179d007619Sjeremylt 19189d007619Sjeremylt @ref Utility 19199d007619Sjeremylt **/ 19209d007619Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) { 19219d007619Sjeremylt *Q = basis->Q; 1922e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 19239d007619Sjeremylt } 19249d007619Sjeremylt 19259d007619Sjeremylt /** 1926ca94c3ddSJeremy L Thompson @brief Get total number of quadrature points (in 1 dimension) of a `CeedBasis` 19279d007619Sjeremylt 1928ca94c3ddSJeremy L Thompson @param[in] basis `CeedBasis` 1929d1d35e2fSjeremylt @param[out] Q_1d Variable to store number of quadrature points 19309d007619Sjeremylt 19319d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 19329d007619Sjeremylt 1933b7c9bbdaSJeremy L Thompson @ref Advanced 19349d007619Sjeremylt **/ 1935d1d35e2fSjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q_1d) { 19366e536b99SJeremy L Thompson CeedCheck(basis->is_tensor_basis, CeedBasisReturnCeed(basis), CEED_ERROR_MINOR, "Cannot supply Q_1d for non-tensor CeedBasis"); 1937d1d35e2fSjeremylt *Q_1d = basis->Q_1d; 1938e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 19399d007619Sjeremylt } 19409d007619Sjeremylt 19419d007619Sjeremylt /** 1942ca94c3ddSJeremy L Thompson @brief Get reference coordinates of quadrature points (in `dim` dimensions) of a `CeedBasis` 19439d007619Sjeremylt 1944ca94c3ddSJeremy L Thompson @param[in] basis `CeedBasis` 1945d1d35e2fSjeremylt @param[out] q_ref Variable to store reference coordinates of quadrature points 19469d007619Sjeremylt 19479d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 19489d007619Sjeremylt 1949b7c9bbdaSJeremy L Thompson @ref Advanced 19509d007619Sjeremylt **/ 1951d1d35e2fSjeremylt int CeedBasisGetQRef(CeedBasis basis, const CeedScalar **q_ref) { 1952d1d35e2fSjeremylt *q_ref = basis->q_ref_1d; 1953e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 19549d007619Sjeremylt } 19559d007619Sjeremylt 19569d007619Sjeremylt /** 1957ca94c3ddSJeremy L Thompson @brief Get quadrature weights of quadrature points (in `dim` dimensions) of a `CeedBasis` 19589d007619Sjeremylt 1959ca94c3ddSJeremy L Thompson @param[in] basis `CeedBasis` 1960d1d35e2fSjeremylt @param[out] q_weight Variable to store quadrature weights 19619d007619Sjeremylt 19629d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 19639d007619Sjeremylt 1964b7c9bbdaSJeremy L Thompson @ref Advanced 19659d007619Sjeremylt **/ 1966d1d35e2fSjeremylt int CeedBasisGetQWeights(CeedBasis basis, const CeedScalar **q_weight) { 1967d1d35e2fSjeremylt *q_weight = basis->q_weight_1d; 1968e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 19699d007619Sjeremylt } 19709d007619Sjeremylt 19719d007619Sjeremylt /** 1972ca94c3ddSJeremy L Thompson @brief Get interpolation matrix of a `CeedBasis` 19739d007619Sjeremylt 1974ca94c3ddSJeremy L Thompson @param[in] basis `CeedBasis` 19759d007619Sjeremylt @param[out] interp Variable to store interpolation matrix 19769d007619Sjeremylt 19779d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 19789d007619Sjeremylt 1979b7c9bbdaSJeremy L Thompson @ref Advanced 19809d007619Sjeremylt **/ 19816c58de82SJeremy L Thompson int CeedBasisGetInterp(CeedBasis basis, const CeedScalar **interp) { 19826402da51SJeremy L Thompson if (!basis->interp && basis->is_tensor_basis) { 19839d007619Sjeremylt // Allocate 19842b730f8bSJeremy L Thompson CeedCall(CeedMalloc(basis->Q * basis->P, &basis->interp)); 19859d007619Sjeremylt 19869d007619Sjeremylt // Initialize 19872b730f8bSJeremy L Thompson for (CeedInt i = 0; i < basis->Q * basis->P; i++) basis->interp[i] = 1.0; 19889d007619Sjeremylt 19899d007619Sjeremylt // Calculate 19902b730f8bSJeremy L Thompson for (CeedInt d = 0; d < basis->dim; d++) { 19912b730f8bSJeremy L Thompson for (CeedInt qpt = 0; qpt < basis->Q; qpt++) { 19929d007619Sjeremylt for (CeedInt node = 0; node < basis->P; node++) { 1993d1d35e2fSjeremylt CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d; 1994d1d35e2fSjeremylt CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d; 19951c66c397SJeremy L Thompson 1996d1d35e2fSjeremylt basis->interp[qpt * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p]; 19979d007619Sjeremylt } 19989d007619Sjeremylt } 19992b730f8bSJeremy L Thompson } 20002b730f8bSJeremy L Thompson } 20019d007619Sjeremylt *interp = basis->interp; 2002e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 20039d007619Sjeremylt } 20049d007619Sjeremylt 20059d007619Sjeremylt /** 2006ca94c3ddSJeremy L Thompson @brief Get 1D interpolation matrix of a tensor product `CeedBasis` 20079d007619Sjeremylt 2008ca94c3ddSJeremy L Thompson @param[in] basis `CeedBasis` 2009d1d35e2fSjeremylt @param[out] interp_1d Variable to store interpolation matrix 20109d007619Sjeremylt 20119d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 20129d007619Sjeremylt 20139d007619Sjeremylt @ref Backend 20149d007619Sjeremylt **/ 2015d1d35e2fSjeremylt int CeedBasisGetInterp1D(CeedBasis basis, const CeedScalar **interp_1d) { 20161203703bSJeremy L Thompson bool is_tensor_basis; 20171203703bSJeremy L Thompson 20181203703bSJeremy L Thompson CeedCall(CeedBasisIsTensor(basis, &is_tensor_basis)); 20196e536b99SJeremy L Thompson CeedCheck(is_tensor_basis, CeedBasisReturnCeed(basis), CEED_ERROR_MINOR, "CeedBasis is not a tensor product CeedBasis"); 2020d1d35e2fSjeremylt *interp_1d = basis->interp_1d; 2021e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 20229d007619Sjeremylt } 20239d007619Sjeremylt 20249d007619Sjeremylt /** 2025ca94c3ddSJeremy L Thompson @brief Get gradient matrix of a `CeedBasis` 20269d007619Sjeremylt 2027ca94c3ddSJeremy L Thompson @param[in] basis `CeedBasis` 20289d007619Sjeremylt @param[out] grad Variable to store gradient matrix 20299d007619Sjeremylt 20309d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 20319d007619Sjeremylt 2032b7c9bbdaSJeremy L Thompson @ref Advanced 20339d007619Sjeremylt **/ 20346c58de82SJeremy L Thompson int CeedBasisGetGrad(CeedBasis basis, const CeedScalar **grad) { 20356402da51SJeremy L Thompson if (!basis->grad && basis->is_tensor_basis) { 20369d007619Sjeremylt // Allocate 20372b730f8bSJeremy L Thompson CeedCall(CeedMalloc(basis->dim * basis->Q * basis->P, &basis->grad)); 20389d007619Sjeremylt 20399d007619Sjeremylt // Initialize 20402b730f8bSJeremy L Thompson for (CeedInt i = 0; i < basis->dim * basis->Q * basis->P; i++) basis->grad[i] = 1.0; 20419d007619Sjeremylt 20429d007619Sjeremylt // Calculate 20432b730f8bSJeremy L Thompson for (CeedInt d = 0; d < basis->dim; d++) { 20442b730f8bSJeremy L Thompson for (CeedInt i = 0; i < basis->dim; i++) { 20452b730f8bSJeremy L Thompson for (CeedInt qpt = 0; qpt < basis->Q; qpt++) { 20469d007619Sjeremylt for (CeedInt node = 0; node < basis->P; node++) { 2047d1d35e2fSjeremylt CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d; 2048d1d35e2fSjeremylt CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d; 20491c66c397SJeremy L Thompson 20502b730f8bSJeremy L Thompson if (i == d) basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->grad_1d[q * basis->P_1d + p]; 20512b730f8bSJeremy L Thompson else basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p]; 20522b730f8bSJeremy L Thompson } 20532b730f8bSJeremy L Thompson } 20542b730f8bSJeremy L Thompson } 20559d007619Sjeremylt } 20569d007619Sjeremylt } 20579d007619Sjeremylt *grad = basis->grad; 2058e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 20599d007619Sjeremylt } 20609d007619Sjeremylt 20619d007619Sjeremylt /** 2062ca94c3ddSJeremy L Thompson @brief Get 1D gradient matrix of a tensor product `CeedBasis` 20639d007619Sjeremylt 2064ca94c3ddSJeremy L Thompson @param[in] basis `CeedBasis` 2065d1d35e2fSjeremylt @param[out] grad_1d Variable to store gradient matrix 20669d007619Sjeremylt 20679d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 20689d007619Sjeremylt 2069b7c9bbdaSJeremy L Thompson @ref Advanced 20709d007619Sjeremylt **/ 2071d1d35e2fSjeremylt int CeedBasisGetGrad1D(CeedBasis basis, const CeedScalar **grad_1d) { 20721203703bSJeremy L Thompson bool is_tensor_basis; 20731203703bSJeremy L Thompson 20741203703bSJeremy L Thompson CeedCall(CeedBasisIsTensor(basis, &is_tensor_basis)); 20756e536b99SJeremy L Thompson CeedCheck(is_tensor_basis, CeedBasisReturnCeed(basis), CEED_ERROR_MINOR, "CeedBasis is not a tensor product CeedBasis"); 2076d1d35e2fSjeremylt *grad_1d = basis->grad_1d; 2077e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 20789d007619Sjeremylt } 20799d007619Sjeremylt 20809d007619Sjeremylt /** 2081ca94c3ddSJeremy L Thompson @brief Get divergence matrix of a `CeedBasis` 208250c301a5SRezgar Shakeri 2083ca94c3ddSJeremy L Thompson @param[in] basis `CeedBasis` 208450c301a5SRezgar Shakeri @param[out] div Variable to store divergence matrix 208550c301a5SRezgar Shakeri 208650c301a5SRezgar Shakeri @return An error code: 0 - success, otherwise - failure 208750c301a5SRezgar Shakeri 208850c301a5SRezgar Shakeri @ref Advanced 208950c301a5SRezgar Shakeri **/ 209050c301a5SRezgar Shakeri int CeedBasisGetDiv(CeedBasis basis, const CeedScalar **div) { 209150c301a5SRezgar Shakeri *div = basis->div; 209250c301a5SRezgar Shakeri return CEED_ERROR_SUCCESS; 209350c301a5SRezgar Shakeri } 209450c301a5SRezgar Shakeri 209550c301a5SRezgar Shakeri /** 2096ca94c3ddSJeremy L Thompson @brief Get curl matrix of a `CeedBasis` 2097c4e3f59bSSebastian Grimberg 2098ca94c3ddSJeremy L Thompson @param[in] basis `CeedBasis` 2099c4e3f59bSSebastian Grimberg @param[out] curl Variable to store curl matrix 2100c4e3f59bSSebastian Grimberg 2101c4e3f59bSSebastian Grimberg @return An error code: 0 - success, otherwise - failure 2102c4e3f59bSSebastian Grimberg 2103c4e3f59bSSebastian Grimberg @ref Advanced 2104c4e3f59bSSebastian Grimberg **/ 2105c4e3f59bSSebastian Grimberg int CeedBasisGetCurl(CeedBasis basis, const CeedScalar **curl) { 2106c4e3f59bSSebastian Grimberg *curl = basis->curl; 2107c4e3f59bSSebastian Grimberg return CEED_ERROR_SUCCESS; 2108c4e3f59bSSebastian Grimberg } 2109c4e3f59bSSebastian Grimberg 2110c4e3f59bSSebastian Grimberg /** 2111ca94c3ddSJeremy L Thompson @brief Destroy a @ref CeedBasis 21127a982d89SJeremy L. Thompson 2113ca94c3ddSJeremy L Thompson @param[in,out] basis `CeedBasis` to destroy 21147a982d89SJeremy L. Thompson 21157a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 21167a982d89SJeremy L. Thompson 21177a982d89SJeremy L. Thompson @ref User 21187a982d89SJeremy L. Thompson **/ 21197a982d89SJeremy L. Thompson int CeedBasisDestroy(CeedBasis *basis) { 2120356036faSJeremy L Thompson if (!*basis || *basis == CEED_BASIS_NONE || --(*basis)->ref_count > 0) { 2121ad6481ceSJeremy L Thompson *basis = NULL; 2122ad6481ceSJeremy L Thompson return CEED_ERROR_SUCCESS; 2123ad6481ceSJeremy L Thompson } 21242b730f8bSJeremy L Thompson if ((*basis)->Destroy) CeedCall((*basis)->Destroy(*basis)); 21259831d45aSJeremy L Thompson CeedCall(CeedTensorContractDestroy(&(*basis)->contract)); 2126c4e3f59bSSebastian Grimberg CeedCall(CeedFree(&(*basis)->q_ref_1d)); 2127c4e3f59bSSebastian Grimberg CeedCall(CeedFree(&(*basis)->q_weight_1d)); 21282b730f8bSJeremy L Thompson CeedCall(CeedFree(&(*basis)->interp)); 21292b730f8bSJeremy L Thompson CeedCall(CeedFree(&(*basis)->interp_1d)); 21302b730f8bSJeremy L Thompson CeedCall(CeedFree(&(*basis)->grad)); 21312b730f8bSJeremy L Thompson CeedCall(CeedFree(&(*basis)->grad_1d)); 2132c4e3f59bSSebastian Grimberg CeedCall(CeedFree(&(*basis)->div)); 2133c4e3f59bSSebastian Grimberg CeedCall(CeedFree(&(*basis)->curl)); 2134c8c3fa7dSJeremy L Thompson CeedCall(CeedVectorDestroy(&(*basis)->vec_chebyshev)); 2135c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisDestroy(&(*basis)->basis_chebyshev)); 21362b730f8bSJeremy L Thompson CeedCall(CeedDestroy(&(*basis)->ceed)); 21372b730f8bSJeremy L Thompson CeedCall(CeedFree(basis)); 2138e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 21397a982d89SJeremy L. Thompson } 21407a982d89SJeremy L. Thompson 21417a982d89SJeremy L. Thompson /** 2142b11c1e72Sjeremylt @brief Construct a Gauss-Legendre quadrature 2143b11c1e72Sjeremylt 2144ca94c3ddSJeremy L Thompson @param[in] Q Number of quadrature points (integrates polynomials of degree `2*Q-1` exactly) 2145ca94c3ddSJeremy L Thompson @param[out] q_ref_1d Array of length `Q` to hold the abscissa on `[-1, 1]` 2146ca94c3ddSJeremy L Thompson @param[out] q_weight_1d Array of length `Q` to hold the weights 2147b11c1e72Sjeremylt 2148b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 2149dfdf5a53Sjeremylt 2150dfdf5a53Sjeremylt @ref Utility 2151b11c1e72Sjeremylt **/ 21522b730f8bSJeremy L Thompson int CeedGaussQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) { 2153d7b241e6Sjeremylt CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0 * atan(1.0); 21541c66c397SJeremy L Thompson 2155d1d35e2fSjeremylt // Build q_ref_1d, q_weight_1d 215692ae7e47SJeremy L Thompson for (CeedInt i = 0; i <= Q / 2; i++) { 2157d7b241e6Sjeremylt // Guess 2158d7b241e6Sjeremylt xi = cos(PI * (CeedScalar)(2 * i + 1) / ((CeedScalar)(2 * Q))); 2159d7b241e6Sjeremylt // Pn(xi) 2160d7b241e6Sjeremylt P0 = 1.0; 2161d7b241e6Sjeremylt P1 = xi; 2162d7b241e6Sjeremylt P2 = 0.0; 216392ae7e47SJeremy L Thompson for (CeedInt j = 2; j <= Q; j++) { 2164d7b241e6Sjeremylt P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j)); 2165d7b241e6Sjeremylt P0 = P1; 2166d7b241e6Sjeremylt P1 = P2; 2167d7b241e6Sjeremylt } 2168d7b241e6Sjeremylt // First Newton Step 2169d7b241e6Sjeremylt dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0); 2170d7b241e6Sjeremylt xi = xi - P2 / dP2; 2171d7b241e6Sjeremylt // Newton to convergence 217292ae7e47SJeremy L Thompson for (CeedInt k = 0; k < 100 && fabs(P2) > 10 * CEED_EPSILON; k++) { 2173d7b241e6Sjeremylt P0 = 1.0; 2174d7b241e6Sjeremylt P1 = xi; 217592ae7e47SJeremy L Thompson for (CeedInt j = 2; j <= Q; j++) { 2176d7b241e6Sjeremylt P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j)); 2177d7b241e6Sjeremylt P0 = P1; 2178d7b241e6Sjeremylt P1 = P2; 2179d7b241e6Sjeremylt } 2180d7b241e6Sjeremylt dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0); 2181d7b241e6Sjeremylt xi = xi - P2 / dP2; 2182d7b241e6Sjeremylt } 2183d7b241e6Sjeremylt // Save xi, wi 2184d7b241e6Sjeremylt wi = 2.0 / ((1.0 - xi * xi) * dP2 * dP2); 2185d1d35e2fSjeremylt q_weight_1d[i] = wi; 2186d1d35e2fSjeremylt q_weight_1d[Q - 1 - i] = wi; 2187d1d35e2fSjeremylt q_ref_1d[i] = -xi; 2188d1d35e2fSjeremylt q_ref_1d[Q - 1 - i] = xi; 2189d7b241e6Sjeremylt } 2190e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 2191d7b241e6Sjeremylt } 2192d7b241e6Sjeremylt 2193b11c1e72Sjeremylt /** 2194b11c1e72Sjeremylt @brief Construct a Gauss-Legendre-Lobatto quadrature 2195b11c1e72Sjeremylt 2196ca94c3ddSJeremy L Thompson @param[in] Q Number of quadrature points (integrates polynomials of degree `2*Q-3` exactly) 2197ca94c3ddSJeremy L Thompson @param[out] q_ref_1d Array of length `Q` to hold the abscissa on `[-1, 1]` 2198ca94c3ddSJeremy L Thompson @param[out] q_weight_1d Array of length `Q` to hold the weights 2199b11c1e72Sjeremylt 2200b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 2201dfdf5a53Sjeremylt 2202dfdf5a53Sjeremylt @ref Utility 2203b11c1e72Sjeremylt **/ 22042b730f8bSJeremy L Thompson int CeedLobattoQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) { 2205d7b241e6Sjeremylt CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0 * atan(1.0); 22061c66c397SJeremy L Thompson 2207d1d35e2fSjeremylt // Build q_ref_1d, q_weight_1d 2208d7b241e6Sjeremylt // Set endpoints 22096574a04fSJeremy L Thompson CeedCheck(Q > 1, NULL, CEED_ERROR_DIMENSION, "Cannot create Lobatto quadrature with Q=%" CeedInt_FMT " < 2 points", Q); 2210d7b241e6Sjeremylt wi = 2.0 / ((CeedScalar)(Q * (Q - 1))); 2211d1d35e2fSjeremylt if (q_weight_1d) { 2212d1d35e2fSjeremylt q_weight_1d[0] = wi; 2213d1d35e2fSjeremylt q_weight_1d[Q - 1] = wi; 2214d7b241e6Sjeremylt } 2215d1d35e2fSjeremylt q_ref_1d[0] = -1.0; 2216d1d35e2fSjeremylt q_ref_1d[Q - 1] = 1.0; 2217d7b241e6Sjeremylt // Interior 221892ae7e47SJeremy L Thompson for (CeedInt i = 1; i <= (Q - 1) / 2; i++) { 2219d7b241e6Sjeremylt // Guess 2220d7b241e6Sjeremylt xi = cos(PI * (CeedScalar)(i) / (CeedScalar)(Q - 1)); 2221d7b241e6Sjeremylt // Pn(xi) 2222d7b241e6Sjeremylt P0 = 1.0; 2223d7b241e6Sjeremylt P1 = xi; 2224d7b241e6Sjeremylt P2 = 0.0; 222592ae7e47SJeremy L Thompson for (CeedInt j = 2; j < Q; j++) { 2226d7b241e6Sjeremylt P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j)); 2227d7b241e6Sjeremylt P0 = P1; 2228d7b241e6Sjeremylt P1 = P2; 2229d7b241e6Sjeremylt } 2230d7b241e6Sjeremylt // First Newton step 2231d7b241e6Sjeremylt dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0); 2232d7b241e6Sjeremylt d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi); 2233d7b241e6Sjeremylt xi = xi - dP2 / d2P2; 2234d7b241e6Sjeremylt // Newton to convergence 223592ae7e47SJeremy L Thompson for (CeedInt k = 0; k < 100 && fabs(dP2) > 10 * CEED_EPSILON; k++) { 2236d7b241e6Sjeremylt P0 = 1.0; 2237d7b241e6Sjeremylt P1 = xi; 223892ae7e47SJeremy L Thompson for (CeedInt j = 2; j < Q; j++) { 2239d7b241e6Sjeremylt P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j)); 2240d7b241e6Sjeremylt P0 = P1; 2241d7b241e6Sjeremylt P1 = P2; 2242d7b241e6Sjeremylt } 2243d7b241e6Sjeremylt dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0); 2244d7b241e6Sjeremylt d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi); 2245d7b241e6Sjeremylt xi = xi - dP2 / d2P2; 2246d7b241e6Sjeremylt } 2247d7b241e6Sjeremylt // Save xi, wi 2248d7b241e6Sjeremylt wi = 2.0 / (((CeedScalar)(Q * (Q - 1))) * P2 * P2); 2249d1d35e2fSjeremylt if (q_weight_1d) { 2250d1d35e2fSjeremylt q_weight_1d[i] = wi; 2251d1d35e2fSjeremylt q_weight_1d[Q - 1 - i] = wi; 2252d7b241e6Sjeremylt } 2253d1d35e2fSjeremylt q_ref_1d[i] = -xi; 2254d1d35e2fSjeremylt q_ref_1d[Q - 1 - i] = xi; 2255d7b241e6Sjeremylt } 2256e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 2257d7b241e6Sjeremylt } 2258d7b241e6Sjeremylt 2259d7b241e6Sjeremylt /// @} 2260