13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3d7b241e6Sjeremylt // 43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 5d7b241e6Sjeremylt // 63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 7d7b241e6Sjeremylt 83d576824SJeremy L Thompson #include <ceed-impl.h> 949aac155SJeremy L Thompson #include <ceed.h> 102b730f8bSJeremy L Thompson #include <ceed/backend.h> 11d7b241e6Sjeremylt #include <math.h> 123d576824SJeremy L Thompson #include <stdbool.h> 13d7b241e6Sjeremylt #include <stdio.h> 14d7b241e6Sjeremylt #include <string.h> 15d7b241e6Sjeremylt 167a982d89SJeremy L. Thompson /// @file 177a982d89SJeremy L. Thompson /// Implementation of CeedBasis interfaces 187a982d89SJeremy L. Thompson 19d7b241e6Sjeremylt /// @cond DOXYGEN_SKIP 20783c99b3SValeria Barra static struct CeedBasis_private ceed_basis_collocated; 21d7b241e6Sjeremylt /// @endcond 22d7b241e6Sjeremylt 237a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser 247a982d89SJeremy L. Thompson /// @{ 257a982d89SJeremy L. Thompson 267a982d89SJeremy L. Thompson /// Indicate that the quadrature points are collocated with the nodes 277a982d89SJeremy L. Thompson const CeedBasis CEED_BASIS_COLLOCATED = &ceed_basis_collocated; 287a982d89SJeremy L. Thompson 297a982d89SJeremy L. Thompson /// @} 307a982d89SJeremy L. Thompson 317a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 327a982d89SJeremy L. Thompson /// CeedBasis Library Internal Functions 337a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 347a982d89SJeremy L. Thompson /// @addtogroup CeedBasisDeveloper 357a982d89SJeremy L. Thompson /// @{ 367a982d89SJeremy L. Thompson 377a982d89SJeremy L. Thompson /** 387a982d89SJeremy L. Thompson @brief Compute Householder reflection 397a982d89SJeremy L. Thompson 40ea61e9acSJeremy L Thompson Computes A = (I - b v v^T) A, where A is an mxn matrix indexed as A[i*row + j*col] 417a982d89SJeremy L. Thompson 427a982d89SJeremy L. Thompson @param[in,out] A Matrix to apply Householder reflection to, in place 43ea61e9acSJeremy L Thompson @param[in] v Householder vector 44ea61e9acSJeremy L Thompson @param[in] b Scaling factor 45ea61e9acSJeremy L Thompson @param[in] m Number of rows in A 46ea61e9acSJeremy L Thompson @param[in] n Number of columns in A 47ea61e9acSJeremy L Thompson @param[in] row Row stride 48ea61e9acSJeremy L Thompson @param[in] col Col stride 497a982d89SJeremy L. Thompson 507a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 517a982d89SJeremy L. Thompson 527a982d89SJeremy L. Thompson @ref Developer 537a982d89SJeremy L. Thompson **/ 542b730f8bSJeremy L Thompson static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v, CeedScalar b, CeedInt m, CeedInt n, CeedInt row, CeedInt col) { 557a982d89SJeremy L. Thompson for (CeedInt j = 0; j < n; j++) { 567a982d89SJeremy L. Thompson CeedScalar w = A[0 * row + j * col]; 572b730f8bSJeremy L Thompson for (CeedInt i = 1; i < m; i++) w += v[i] * A[i * row + j * col]; 587a982d89SJeremy L. Thompson A[0 * row + j * col] -= b * w; 592b730f8bSJeremy L Thompson for (CeedInt i = 1; i < m; i++) A[i * row + j * col] -= b * w * v[i]; 607a982d89SJeremy L. Thompson } 61e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 627a982d89SJeremy L. Thompson } 637a982d89SJeremy L. Thompson 647a982d89SJeremy L. Thompson /** 657a982d89SJeremy L. Thompson @brief Compute Givens rotation 667a982d89SJeremy L. Thompson 67ea61e9acSJeremy L Thompson Computes A = G A (or G^T A in transpose mode), where A is an mxn matrix indexed as A[i*n + j*m] 687a982d89SJeremy L. Thompson 697a982d89SJeremy L. Thompson @param[in,out] A Row major matrix to apply Givens rotation to, in place 70ea61e9acSJeremy L Thompson @param[in] c Cosine factor 71ea61e9acSJeremy L Thompson @param[in] s Sine factor 72ea61e9acSJeremy L Thompson @param[in] t_mode @ref CEED_NOTRANSPOSE to rotate the basis counter-clockwise, which has the effect of rotating columns of A clockwise; 734cc79fe7SJed Brown @ref CEED_TRANSPOSE for the opposite rotation 74ea61e9acSJeremy L Thompson @param[in] i First row/column to apply rotation 75ea61e9acSJeremy L Thompson @param[in] k Second row/column to apply rotation 76ea61e9acSJeremy L Thompson @param[in] m Number of rows in A 77ea61e9acSJeremy L Thompson @param[in] n Number of columns in A 787a982d89SJeremy L. Thompson 797a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 807a982d89SJeremy L. Thompson 817a982d89SJeremy L. Thompson @ref Developer 827a982d89SJeremy L. Thompson **/ 832b730f8bSJeremy L Thompson static int CeedGivensRotation(CeedScalar *A, CeedScalar c, CeedScalar s, CeedTransposeMode t_mode, CeedInt i, CeedInt k, CeedInt m, CeedInt n) { 84d1d35e2fSjeremylt CeedInt stride_j = 1, stride_ik = m, num_its = n; 85d1d35e2fSjeremylt if (t_mode == CEED_NOTRANSPOSE) { 862b730f8bSJeremy L Thompson stride_j = n; 872b730f8bSJeremy L Thompson stride_ik = 1; 882b730f8bSJeremy L Thompson num_its = m; 897a982d89SJeremy L. Thompson } 907a982d89SJeremy L. Thompson 917a982d89SJeremy L. Thompson // Apply rotation 92d1d35e2fSjeremylt for (CeedInt j = 0; j < num_its; j++) { 93d1d35e2fSjeremylt CeedScalar tau1 = A[i * stride_ik + j * stride_j], tau2 = A[k * stride_ik + j * stride_j]; 94d1d35e2fSjeremylt A[i * stride_ik + j * stride_j] = c * tau1 - s * tau2; 95d1d35e2fSjeremylt A[k * stride_ik + j * stride_j] = s * tau1 + c * tau2; 967a982d89SJeremy L. Thompson } 97e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 987a982d89SJeremy L. Thompson } 997a982d89SJeremy L. Thompson 1007a982d89SJeremy L. Thompson /** 1017a982d89SJeremy L. Thompson @brief View an array stored in a CeedBasis 1027a982d89SJeremy L. Thompson 1030a0da059Sjeremylt @param[in] name Name of array 104d1d35e2fSjeremylt @param[in] fp_fmt Printing format 1050a0da059Sjeremylt @param[in] m Number of rows in array 1060a0da059Sjeremylt @param[in] n Number of columns in array 1070a0da059Sjeremylt @param[in] a Array to be viewed 1080a0da059Sjeremylt @param[in] stream Stream to view to, e.g., stdout 1097a982d89SJeremy L. Thompson 1107a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 1117a982d89SJeremy L. Thompson 1127a982d89SJeremy L. Thompson @ref Developer 1137a982d89SJeremy L. Thompson **/ 1142b730f8bSJeremy L Thompson static int CeedScalarView(const char *name, const char *fp_fmt, CeedInt m, CeedInt n, const CeedScalar *a, FILE *stream) { 11592ae7e47SJeremy L Thompson for (CeedInt i = 0; i < m; i++) { 1162b730f8bSJeremy L Thompson if (m > 1) fprintf(stream, "%12s[%" CeedInt_FMT "]:", name, i); 1172b730f8bSJeremy L Thompson else fprintf(stream, "%12s:", name); 1182b730f8bSJeremy L Thompson for (CeedInt j = 0; j < n; j++) fprintf(stream, fp_fmt, fabs(a[i * n + j]) > 1E-14 ? a[i * n + j] : 0); 1197a982d89SJeremy L. Thompson fputs("\n", stream); 1207a982d89SJeremy L. Thompson } 121e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1227a982d89SJeremy L. Thompson } 1237a982d89SJeremy L. Thompson 124a76a04e7SJeremy L Thompson /** 125ea61e9acSJeremy L Thompson @brief Create the interpolation and gradient matrices for projection from the nodes of `basis_from` to the nodes of `basis_to`. 126ba59ac12SSebastian Grimberg 12715ad3917SSebastian Grimberg The interpolation is given by `interp_project = interp_to^+ * interp_from`, where the pseudoinverse `interp_to^+` is given by QR factorization. 12815ad3917SSebastian Grimberg The gradient is given by `grad_project = interp_to^+ * grad_from`, and is only computed for H^1 spaces otherwise it should not be used. 12915ad3917SSebastian Grimberg 130ba59ac12SSebastian Grimberg Note: `basis_from` and `basis_to` must have compatible quadrature spaces. 131a76a04e7SJeremy L Thompson 132a76a04e7SJeremy L Thompson @param[in] basis_from CeedBasis to project from 133a76a04e7SJeremy L Thompson @param[in] basis_to CeedBasis to project to 134ea61e9acSJeremy L Thompson @param[out] interp_project Address of the variable where the newly created interpolation matrix will be stored. 135ea61e9acSJeremy L Thompson @param[out] grad_project Address of the variable where the newly created gradient matrix will be stored. 136a76a04e7SJeremy L Thompson 137a76a04e7SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 138a76a04e7SJeremy L Thompson 139a76a04e7SJeremy L Thompson @ref Developer 140a76a04e7SJeremy L Thompson **/ 1412b730f8bSJeremy L Thompson static int CeedBasisCreateProjectionMatrices(CeedBasis basis_from, CeedBasis basis_to, CeedScalar **interp_project, CeedScalar **grad_project) { 142a76a04e7SJeremy L Thompson Ceed ceed; 1432b730f8bSJeremy L Thompson CeedCall(CeedBasisGetCeed(basis_to, &ceed)); 144a76a04e7SJeremy L Thompson 145a76a04e7SJeremy L Thompson // Check for compatible quadrature spaces 146a76a04e7SJeremy L Thompson CeedInt Q_to, Q_from; 1472b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints(basis_to, &Q_to)); 1482b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints(basis_from, &Q_from)); 1496574a04fSJeremy L Thompson CeedCheck(Q_to == Q_from, ceed, CEED_ERROR_DIMENSION, "Bases must have compatible quadrature spaces"); 150a76a04e7SJeremy L Thompson 15114556e63SJeremy L Thompson // Check for matching tensor or non-tensor 152a76a04e7SJeremy L Thompson CeedInt P_to, P_from, Q = Q_to; 153a76a04e7SJeremy L Thompson bool is_tensor_to, is_tensor_from; 1542b730f8bSJeremy L Thompson CeedCall(CeedBasisIsTensor(basis_to, &is_tensor_to)); 1552b730f8bSJeremy L Thompson CeedCall(CeedBasisIsTensor(basis_from, &is_tensor_from)); 1566574a04fSJeremy L Thompson CeedCheck(is_tensor_to == is_tensor_from, ceed, CEED_ERROR_MINOR, "Bases must both be tensor or non-tensor"); 1576574a04fSJeremy L Thompson if (is_tensor_to) { 1582b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_to)); 1592b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_from)); 1602b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints1D(basis_from, &Q)); 1616574a04fSJeremy L Thompson } else { 1622b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis_to, &P_to)); 1632b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis_from, &P_from)); 164a76a04e7SJeremy L Thompson } 165a76a04e7SJeremy L Thompson 16615ad3917SSebastian Grimberg // Check for matching FE space 16715ad3917SSebastian Grimberg CeedFESpace fe_space_to, fe_space_from; 16815ad3917SSebastian Grimberg CeedCall(CeedBasisGetFESpace(basis_to, &fe_space_to)); 16915ad3917SSebastian Grimberg CeedCall(CeedBasisGetFESpace(basis_from, &fe_space_from)); 1706574a04fSJeremy L Thompson CeedCheck(fe_space_to == fe_space_from, ceed, CEED_ERROR_MINOR, "Bases must both be the same FE space type"); 17115ad3917SSebastian Grimberg 17214556e63SJeremy L Thompson // Get source matrices 17315ad3917SSebastian Grimberg CeedInt dim, q_comp = 1; 17415ad3917SSebastian Grimberg const CeedScalar *interp_to_source = NULL, *interp_from_source = NULL; 17514556e63SJeremy L Thompson CeedScalar *interp_to, *interp_from, *tau; 1762b730f8bSJeremy L Thompson CeedCall(CeedBasisGetDimension(basis_to, &dim)); 177a76a04e7SJeremy L Thompson if (is_tensor_to) { 1782b730f8bSJeremy L Thompson CeedCall(CeedBasisGetInterp1D(basis_to, &interp_to_source)); 1792b730f8bSJeremy L Thompson CeedCall(CeedBasisGetInterp1D(basis_from, &interp_from_source)); 180a76a04e7SJeremy L Thompson } else { 18115ad3917SSebastian Grimberg CeedCall(CeedBasisGetNumQuadratureComponents(basis_from, CEED_EVAL_INTERP, &q_comp)); 1822b730f8bSJeremy L Thompson CeedCall(CeedBasisGetInterp(basis_to, &interp_to_source)); 1832b730f8bSJeremy L Thompson CeedCall(CeedBasisGetInterp(basis_from, &interp_from_source)); 18415ad3917SSebastian Grimberg } 18515ad3917SSebastian Grimberg CeedCall(CeedMalloc(Q * P_from * q_comp, &interp_from)); 18615ad3917SSebastian Grimberg CeedCall(CeedMalloc(Q * P_to * q_comp, &interp_to)); 18715ad3917SSebastian Grimberg CeedCall(CeedCalloc(P_to * P_from, interp_project)); 18815ad3917SSebastian Grimberg CeedCall(CeedMalloc(Q * q_comp, &tau)); 18915ad3917SSebastian Grimberg 19015ad3917SSebastian Grimberg // `grad_project = interp_to^+ * grad_from` is computed for the H^1 space case so the 191de05fbb2SSebastian Grimberg // projection basis will have a gradient operation (allocated even if not H^1 for the 192de05fbb2SSebastian Grimberg // basis construction later on) 19315ad3917SSebastian Grimberg const CeedScalar *grad_from_source = NULL; 19415ad3917SSebastian Grimberg if (fe_space_to == CEED_FE_SPACE_H1) { 19515ad3917SSebastian Grimberg if (is_tensor_to) { 19615ad3917SSebastian Grimberg CeedCall(CeedBasisGetGrad1D(basis_from, &grad_from_source)); 19715ad3917SSebastian Grimberg } else { 1982b730f8bSJeremy L Thompson CeedCall(CeedBasisGetGrad(basis_from, &grad_from_source)); 199a76a04e7SJeremy L Thompson } 200de05fbb2SSebastian Grimberg } 20115ad3917SSebastian Grimberg CeedCall(CeedCalloc(P_to * P_from * (is_tensor_to ? 1 : dim), grad_project)); 20215ad3917SSebastian Grimberg 20315ad3917SSebastian Grimberg // QR Factorization, interp_to = Q R 20415ad3917SSebastian Grimberg memcpy(interp_to, interp_to_source, Q * P_to * q_comp * sizeof(interp_to_source[0])); 20515ad3917SSebastian Grimberg CeedCall(CeedQRFactorization(ceed, interp_to, tau, Q * q_comp, P_to)); 206a76a04e7SJeremy L Thompson 20714556e63SJeremy L Thompson // Build matrices 20815ad3917SSebastian Grimberg CeedInt num_matrices = 1 + (fe_space_to == CEED_FE_SPACE_H1) * (is_tensor_to ? 1 : dim); 20914556e63SJeremy L Thompson CeedScalar *input_from[num_matrices], *output_project[num_matrices]; 21014556e63SJeremy L Thompson input_from[0] = (CeedScalar *)interp_from_source; 21114556e63SJeremy L Thompson output_project[0] = *interp_project; 21214556e63SJeremy L Thompson for (CeedInt m = 1; m < num_matrices; m++) { 21314556e63SJeremy L Thompson input_from[m] = (CeedScalar *)&grad_from_source[(m - 1) * Q * P_from]; 21402af4036SJeremy L Thompson output_project[m] = &((*grad_project)[(m - 1) * P_to * P_from]); 21514556e63SJeremy L Thompson } 21614556e63SJeremy L Thompson for (CeedInt m = 0; m < num_matrices; m++) { 21715ad3917SSebastian Grimberg // Apply Q^T, interp_from = Q^T interp_from 21815ad3917SSebastian Grimberg memcpy(interp_from, input_from[m], Q * P_from * q_comp * sizeof(input_from[m][0])); 21915ad3917SSebastian Grimberg CeedCall(CeedHouseholderApplyQ(interp_from, interp_to, tau, CEED_TRANSPOSE, Q * q_comp, P_from, P_to, P_from, 1)); 220a76a04e7SJeremy L Thompson 22115ad3917SSebastian Grimberg // Apply Rinv, output_project = Rinv interp_from 222a76a04e7SJeremy L Thompson for (CeedInt j = 0; j < P_from; j++) { // Column j 2232b730f8bSJeremy L Thompson output_project[m][j + P_from * (P_to - 1)] = interp_from[j + P_from * (P_to - 1)] / interp_to[P_to * P_to - 1]; 224a76a04e7SJeremy L Thompson for (CeedInt i = P_to - 2; i >= 0; i--) { // Row i 22514556e63SJeremy L Thompson output_project[m][j + P_from * i] = interp_from[j + P_from * i]; 226a76a04e7SJeremy L Thompson for (CeedInt k = i + 1; k < P_to; k++) { 2272b730f8bSJeremy L Thompson output_project[m][j + P_from * i] -= interp_to[k + P_to * i] * output_project[m][j + P_from * k]; 228a76a04e7SJeremy L Thompson } 22914556e63SJeremy L Thompson output_project[m][j + P_from * i] /= interp_to[i + P_to * i]; 230a76a04e7SJeremy L Thompson } 231a76a04e7SJeremy L Thompson } 23214556e63SJeremy L Thompson } 23314556e63SJeremy L Thompson 23414556e63SJeremy L Thompson // Cleanup 2352b730f8bSJeremy L Thompson CeedCall(CeedFree(&tau)); 2362b730f8bSJeremy L Thompson CeedCall(CeedFree(&interp_to)); 2372b730f8bSJeremy L Thompson CeedCall(CeedFree(&interp_from)); 238a76a04e7SJeremy L Thompson 239a76a04e7SJeremy L Thompson return CEED_ERROR_SUCCESS; 240a76a04e7SJeremy L Thompson } 241a76a04e7SJeremy L Thompson 2427a982d89SJeremy L. Thompson /// @} 2437a982d89SJeremy L. Thompson 2447a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 2457a982d89SJeremy L. Thompson /// Ceed Backend API 2467a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 2477a982d89SJeremy L. Thompson /// @addtogroup CeedBasisBackend 2487a982d89SJeremy L. Thompson /// @{ 2497a982d89SJeremy L. Thompson 2507a982d89SJeremy L. Thompson /** 2517a982d89SJeremy L. Thompson @brief Return collocated grad matrix 2527a982d89SJeremy L. Thompson 253ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 254ea61e9acSJeremy L Thompson @param[out] collo_grad_1d Row-major (Q_1d * Q_1d) matrix expressing derivatives of basis functions at quadrature points 2557a982d89SJeremy L. Thompson 2567a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 2577a982d89SJeremy L. Thompson 2587a982d89SJeremy L. Thompson @ref Backend 2597a982d89SJeremy L. Thompson **/ 260d1d35e2fSjeremylt int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *collo_grad_1d) { 2617a982d89SJeremy L. Thompson Ceed ceed; 2622b730f8bSJeremy L Thompson CeedInt P_1d = (basis)->P_1d, Q_1d = (basis)->Q_1d; 26378464608Sjeremylt CeedScalar *interp_1d, *grad_1d, *tau; 2647a982d89SJeremy L. Thompson 2652b730f8bSJeremy L Thompson CeedCall(CeedMalloc(Q_1d * P_1d, &interp_1d)); 2662b730f8bSJeremy L Thompson CeedCall(CeedMalloc(Q_1d * P_1d, &grad_1d)); 2672b730f8bSJeremy L Thompson CeedCall(CeedMalloc(Q_1d, &tau)); 268d1d35e2fSjeremylt memcpy(interp_1d, (basis)->interp_1d, Q_1d * P_1d * sizeof(basis)->interp_1d[0]); 269d1d35e2fSjeremylt memcpy(grad_1d, (basis)->grad_1d, Q_1d * P_1d * sizeof(basis)->interp_1d[0]); 2707a982d89SJeremy L. Thompson 271d1d35e2fSjeremylt // QR Factorization, interp_1d = Q R 2722b730f8bSJeremy L Thompson CeedCall(CeedBasisGetCeed(basis, &ceed)); 2732b730f8bSJeremy L Thompson CeedCall(CeedQRFactorization(ceed, interp_1d, tau, Q_1d, P_1d)); 274ea61e9acSJeremy 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. 2757a982d89SJeremy L. Thompson 276c8c3fa7dSJeremy L Thompson // Apply R_inv, collo_grad_1d = grad_1d R_inv 277c8c3fa7dSJeremy L Thompson for (CeedInt i = 0; i < Q_1d; i++) { // Row i 278d1d35e2fSjeremylt collo_grad_1d[Q_1d * i] = grad_1d[P_1d * i] / interp_1d[0]; 279c8c3fa7dSJeremy L Thompson for (CeedInt j = 1; j < P_1d; j++) { // Column j 280d1d35e2fSjeremylt collo_grad_1d[j + Q_1d * i] = grad_1d[j + P_1d * i]; 281c8c3fa7dSJeremy L Thompson for (CeedInt k = 0; k < j; k++) collo_grad_1d[j + Q_1d * i] -= interp_1d[j + P_1d * k] * collo_grad_1d[k + Q_1d * i]; 282d1d35e2fSjeremylt collo_grad_1d[j + Q_1d * i] /= interp_1d[j + P_1d * j]; 2837a982d89SJeremy L. Thompson } 284c8c3fa7dSJeremy L Thompson for (CeedInt j = P_1d; j < Q_1d; j++) collo_grad_1d[j + Q_1d * i] = 0; 2857a982d89SJeremy L. Thompson } 2867a982d89SJeremy L. Thompson 28715ad3917SSebastian Grimberg // Apply Q^T, collo_grad_1d = collo_grad_1d Q^T 2882b730f8bSJeremy L Thompson CeedCall(CeedHouseholderApplyQ(collo_grad_1d, interp_1d, tau, CEED_NOTRANSPOSE, Q_1d, Q_1d, P_1d, 1, Q_1d)); 2897a982d89SJeremy L. Thompson 2902b730f8bSJeremy L Thompson CeedCall(CeedFree(&interp_1d)); 2912b730f8bSJeremy L Thompson CeedCall(CeedFree(&grad_1d)); 2922b730f8bSJeremy L Thompson CeedCall(CeedFree(&tau)); 293e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 2947a982d89SJeremy L. Thompson } 2957a982d89SJeremy L. Thompson 2967a982d89SJeremy L. Thompson /** 2977a982d89SJeremy L. Thompson @brief Get tensor status for given CeedBasis 2987a982d89SJeremy L. Thompson 299ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 300d1d35e2fSjeremylt @param[out] is_tensor Variable to store tensor status 3017a982d89SJeremy L. Thompson 3027a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3037a982d89SJeremy L. Thompson 3047a982d89SJeremy L. Thompson @ref Backend 3057a982d89SJeremy L. Thompson **/ 306d1d35e2fSjeremylt int CeedBasisIsTensor(CeedBasis basis, bool *is_tensor) { 3076402da51SJeremy L Thompson *is_tensor = basis->is_tensor_basis; 308e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3097a982d89SJeremy L. Thompson } 3107a982d89SJeremy L. Thompson 3117a982d89SJeremy L. Thompson /** 3127a982d89SJeremy L. Thompson @brief Get backend data of a CeedBasis 3137a982d89SJeremy L. Thompson 314ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 3157a982d89SJeremy L. Thompson @param[out] data Variable to store data 3167a982d89SJeremy L. Thompson 3177a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3187a982d89SJeremy L. Thompson 3197a982d89SJeremy L. Thompson @ref Backend 3207a982d89SJeremy L. Thompson **/ 321777ff853SJeremy L Thompson int CeedBasisGetData(CeedBasis basis, void *data) { 322777ff853SJeremy L Thompson *(void **)data = basis->data; 323e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3247a982d89SJeremy L. Thompson } 3257a982d89SJeremy L. Thompson 3267a982d89SJeremy L. Thompson /** 3277a982d89SJeremy L. Thompson @brief Set backend data of a CeedBasis 3287a982d89SJeremy L. Thompson 329ea61e9acSJeremy L Thompson @param[in,out] basis CeedBasis 330ea61e9acSJeremy L Thompson @param[in] data Data to set 3317a982d89SJeremy L. Thompson 3327a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3337a982d89SJeremy L. Thompson 3347a982d89SJeremy L. Thompson @ref Backend 3357a982d89SJeremy L. Thompson **/ 336777ff853SJeremy L Thompson int CeedBasisSetData(CeedBasis basis, void *data) { 337777ff853SJeremy L Thompson basis->data = data; 338e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3397a982d89SJeremy L. Thompson } 3407a982d89SJeremy L. Thompson 3417a982d89SJeremy L. Thompson /** 34234359f16Sjeremylt @brief Increment the reference counter for a CeedBasis 34334359f16Sjeremylt 344ea61e9acSJeremy L Thompson @param[in,out] basis Basis to increment the reference counter 34534359f16Sjeremylt 34634359f16Sjeremylt @return An error code: 0 - success, otherwise - failure 34734359f16Sjeremylt 34834359f16Sjeremylt @ref Backend 34934359f16Sjeremylt **/ 3509560d06aSjeremylt int CeedBasisReference(CeedBasis basis) { 35134359f16Sjeremylt basis->ref_count++; 35234359f16Sjeremylt return CEED_ERROR_SUCCESS; 35334359f16Sjeremylt } 35434359f16Sjeremylt 35534359f16Sjeremylt /** 356c4e3f59bSSebastian Grimberg @brief Get number of Q-vector components for given CeedBasis 357c4e3f59bSSebastian Grimberg 358c4e3f59bSSebastian Grimberg @param[in] basis CeedBasis 359c4e3f59bSSebastian Grimberg @param[in] eval_mode \ref CEED_EVAL_INTERP to use interpolated values, 360c4e3f59bSSebastian Grimberg \ref CEED_EVAL_GRAD to use gradients, 361c4e3f59bSSebastian Grimberg \ref CEED_EVAL_DIV to use divergence, 362c4e3f59bSSebastian Grimberg \ref CEED_EVAL_CURL to use curl. 363c4e3f59bSSebastian Grimberg @param[out] q_comp Variable to store number of Q-vector components of basis 364c4e3f59bSSebastian Grimberg 365c4e3f59bSSebastian Grimberg @return An error code: 0 - success, otherwise - failure 366c4e3f59bSSebastian Grimberg 367c4e3f59bSSebastian Grimberg @ref Backend 368c4e3f59bSSebastian Grimberg **/ 369c4e3f59bSSebastian Grimberg int CeedBasisGetNumQuadratureComponents(CeedBasis basis, CeedEvalMode eval_mode, CeedInt *q_comp) { 370c4e3f59bSSebastian Grimberg switch (eval_mode) { 371c4e3f59bSSebastian Grimberg case CEED_EVAL_INTERP: 372c4e3f59bSSebastian Grimberg *q_comp = (basis->fe_space == CEED_FE_SPACE_H1) ? 1 : basis->dim; 373c4e3f59bSSebastian Grimberg break; 374c4e3f59bSSebastian Grimberg case CEED_EVAL_GRAD: 375c4e3f59bSSebastian Grimberg *q_comp = basis->dim; 376c4e3f59bSSebastian Grimberg break; 377c4e3f59bSSebastian Grimberg case CEED_EVAL_DIV: 378c4e3f59bSSebastian Grimberg *q_comp = 1; 379c4e3f59bSSebastian Grimberg break; 380c4e3f59bSSebastian Grimberg case CEED_EVAL_CURL: 381c4e3f59bSSebastian Grimberg *q_comp = (basis->dim < 3) ? 1 : basis->dim; 382c4e3f59bSSebastian Grimberg break; 383c4e3f59bSSebastian Grimberg case CEED_EVAL_NONE: 384c4e3f59bSSebastian Grimberg case CEED_EVAL_WEIGHT: 385352a5e7cSSebastian Grimberg *q_comp = 1; 386c4e3f59bSSebastian Grimberg break; 387c4e3f59bSSebastian Grimberg } 388c4e3f59bSSebastian Grimberg return CEED_ERROR_SUCCESS; 389c4e3f59bSSebastian Grimberg } 390c4e3f59bSSebastian Grimberg 391c4e3f59bSSebastian Grimberg /** 3926e15d496SJeremy L Thompson @brief Estimate number of FLOPs required to apply CeedBasis in t_mode and eval_mode 3936e15d496SJeremy L Thompson 394ea61e9acSJeremy L Thompson @param[in] basis Basis to estimate FLOPs for 395ea61e9acSJeremy L Thompson @param[in] t_mode Apply basis or transpose 396ea61e9acSJeremy L Thompson @param[in] eval_mode Basis evaluation mode 397ea61e9acSJeremy L Thompson @param[out] flops Address of variable to hold FLOPs estimate 3986e15d496SJeremy L Thompson 3996e15d496SJeremy L Thompson @ref Backend 4006e15d496SJeremy L Thompson **/ 4012b730f8bSJeremy L Thompson int CeedBasisGetFlopsEstimate(CeedBasis basis, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedSize *flops) { 4026e15d496SJeremy L Thompson bool is_tensor; 4036e15d496SJeremy L Thompson 4042b730f8bSJeremy L Thompson CeedCall(CeedBasisIsTensor(basis, &is_tensor)); 4056e15d496SJeremy L Thompson if (is_tensor) { 4066e15d496SJeremy L Thompson CeedInt dim, num_comp, P_1d, Q_1d; 4072b730f8bSJeremy L Thompson CeedCall(CeedBasisGetDimension(basis, &dim)); 4082b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumComponents(basis, &num_comp)); 4092b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d)); 4102b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 4116e15d496SJeremy L Thompson if (t_mode == CEED_TRANSPOSE) { 4122b730f8bSJeremy L Thompson P_1d = Q_1d; 4132b730f8bSJeremy L Thompson Q_1d = P_1d; 4146e15d496SJeremy L Thompson } 4156e15d496SJeremy L Thompson CeedInt tensor_flops = 0, pre = num_comp * CeedIntPow(P_1d, dim - 1), post = 1; 4166e15d496SJeremy L Thompson for (CeedInt d = 0; d < dim; d++) { 4176e15d496SJeremy L Thompson tensor_flops += 2 * pre * P_1d * post * Q_1d; 4186e15d496SJeremy L Thompson pre /= P_1d; 4196e15d496SJeremy L Thompson post *= Q_1d; 4206e15d496SJeremy L Thompson } 4216e15d496SJeremy L Thompson switch (eval_mode) { 4222b730f8bSJeremy L Thompson case CEED_EVAL_NONE: 4232b730f8bSJeremy L Thompson *flops = 0; 4242b730f8bSJeremy L Thompson break; 4252b730f8bSJeremy L Thompson case CEED_EVAL_INTERP: 4262b730f8bSJeremy L Thompson *flops = tensor_flops; 4272b730f8bSJeremy L Thompson break; 4282b730f8bSJeremy L Thompson case CEED_EVAL_GRAD: 4292b730f8bSJeremy L Thompson *flops = tensor_flops * 2; 4302b730f8bSJeremy L Thompson break; 4316e15d496SJeremy L Thompson case CEED_EVAL_DIV: 4326e15d496SJeremy L Thompson case CEED_EVAL_CURL: 4336574a04fSJeremy L Thompson // LCOV_EXCL_START 4346574a04fSJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_INCOMPATIBLE, "Tensor basis evaluation for %s not supported", CeedEvalModes[eval_mode]); 4352b730f8bSJeremy L Thompson break; 4366e15d496SJeremy L Thompson // LCOV_EXCL_STOP 4372b730f8bSJeremy L Thompson case CEED_EVAL_WEIGHT: 4382b730f8bSJeremy L Thompson *flops = dim * CeedIntPow(Q_1d, dim); 4392b730f8bSJeremy L Thompson break; 4406e15d496SJeremy L Thompson } 4416e15d496SJeremy L Thompson } else { 442c4e3f59bSSebastian Grimberg CeedInt dim, num_comp, q_comp, num_nodes, num_qpts; 4432b730f8bSJeremy L Thompson CeedCall(CeedBasisGetDimension(basis, &dim)); 4442b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumComponents(basis, &num_comp)); 445c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp)); 4462b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis, &num_nodes)); 4472b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts)); 4486e15d496SJeremy L Thompson switch (eval_mode) { 4492b730f8bSJeremy L Thompson case CEED_EVAL_NONE: 4502b730f8bSJeremy L Thompson *flops = 0; 4512b730f8bSJeremy L Thompson break; 4522b730f8bSJeremy L Thompson case CEED_EVAL_INTERP: 4532b730f8bSJeremy L Thompson case CEED_EVAL_GRAD: 4542b730f8bSJeremy L Thompson case CEED_EVAL_DIV: 4552b730f8bSJeremy L Thompson case CEED_EVAL_CURL: 456c4e3f59bSSebastian Grimberg *flops = num_nodes * num_qpts * num_comp * q_comp; 4572b730f8bSJeremy L Thompson break; 4582b730f8bSJeremy L Thompson case CEED_EVAL_WEIGHT: 4592b730f8bSJeremy L Thompson *flops = 0; 4602b730f8bSJeremy L Thompson break; 4616e15d496SJeremy L Thompson } 4626e15d496SJeremy L Thompson } 4636e15d496SJeremy L Thompson 4646e15d496SJeremy L Thompson return CEED_ERROR_SUCCESS; 4656e15d496SJeremy L Thompson } 4666e15d496SJeremy L Thompson 4676e15d496SJeremy L Thompson /** 468c4e3f59bSSebastian Grimberg @brief Get CeedFESpace for a CeedBasis 469c4e3f59bSSebastian Grimberg 470c4e3f59bSSebastian Grimberg @param[in] basis CeedBasis 471c4e3f59bSSebastian Grimberg @param[out] fe_space Variable to store CeedFESpace 472c4e3f59bSSebastian Grimberg 473c4e3f59bSSebastian Grimberg @return An error code: 0 - success, otherwise - failure 474c4e3f59bSSebastian Grimberg 475c4e3f59bSSebastian Grimberg @ref Backend 476c4e3f59bSSebastian Grimberg **/ 477c4e3f59bSSebastian Grimberg int CeedBasisGetFESpace(CeedBasis basis, CeedFESpace *fe_space) { 478c4e3f59bSSebastian Grimberg *fe_space = basis->fe_space; 479c4e3f59bSSebastian Grimberg return CEED_ERROR_SUCCESS; 480c4e3f59bSSebastian Grimberg } 481c4e3f59bSSebastian Grimberg 482c4e3f59bSSebastian Grimberg /** 4837a982d89SJeremy L. Thompson @brief Get dimension for given CeedElemTopology 4847a982d89SJeremy L. Thompson 485ea61e9acSJeremy L Thompson @param[in] topo CeedElemTopology 4867a982d89SJeremy L. Thompson @param[out] dim Variable to store dimension of topology 4877a982d89SJeremy L. Thompson 4887a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 4897a982d89SJeremy L. Thompson 4907a982d89SJeremy L. Thompson @ref Backend 4917a982d89SJeremy L. Thompson **/ 4927a982d89SJeremy L. Thompson int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) { 4937a982d89SJeremy L. Thompson *dim = (CeedInt)topo >> 16; 494e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 4957a982d89SJeremy L. Thompson } 4967a982d89SJeremy L. Thompson 4977a982d89SJeremy L. Thompson /** 4987a982d89SJeremy L. Thompson @brief Get CeedTensorContract of a CeedBasis 4997a982d89SJeremy L. Thompson 500ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 5017a982d89SJeremy L. Thompson @param[out] contract Variable to store CeedTensorContract 5027a982d89SJeremy L. Thompson 5037a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 5047a982d89SJeremy L. Thompson 5057a982d89SJeremy L. Thompson @ref Backend 5067a982d89SJeremy L. Thompson **/ 5077a982d89SJeremy L. Thompson int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) { 5087a982d89SJeremy L. Thompson *contract = basis->contract; 509e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 5107a982d89SJeremy L. Thompson } 5117a982d89SJeremy L. Thompson 5127a982d89SJeremy L. Thompson /** 5137a982d89SJeremy L. Thompson @brief Set CeedTensorContract of a CeedBasis 5147a982d89SJeremy L. Thompson 515ea61e9acSJeremy L Thompson @param[in,out] basis CeedBasis 516ea61e9acSJeremy L Thompson @param[in] contract CeedTensorContract to set 5177a982d89SJeremy L. Thompson 5187a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 5197a982d89SJeremy L. Thompson 5207a982d89SJeremy L. Thompson @ref Backend 5217a982d89SJeremy L. Thompson **/ 52234359f16Sjeremylt int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract contract) { 52334359f16Sjeremylt basis->contract = contract; 5242b730f8bSJeremy L Thompson CeedCall(CeedTensorContractReference(contract)); 525e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 5267a982d89SJeremy L. Thompson } 5277a982d89SJeremy L. Thompson 5287a982d89SJeremy L. Thompson /** 5297a982d89SJeremy L. Thompson @brief Return a reference implementation of matrix multiplication C = A B. 530ba59ac12SSebastian Grimberg 531ba59ac12SSebastian Grimberg Note: This is a reference implementation for CPU CeedScalar pointers that is not intended for high performance. 5327a982d89SJeremy L. Thompson 533ea61e9acSJeremy L Thompson @param[in] ceed Ceed context for error handling 534d1d35e2fSjeremylt @param[in] mat_A Row-major matrix A 535d1d35e2fSjeremylt @param[in] mat_B Row-major matrix B 536d1d35e2fSjeremylt @param[out] mat_C Row-major output matrix C 537ea61e9acSJeremy L Thompson @param[in] m Number of rows of C 538ea61e9acSJeremy L Thompson @param[in] n Number of columns of C 539ea61e9acSJeremy L Thompson @param[in] kk Number of columns of A/rows of B 5407a982d89SJeremy L. Thompson 5417a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 5427a982d89SJeremy L. Thompson 5437a982d89SJeremy L. Thompson @ref Utility 5447a982d89SJeremy L. Thompson **/ 5452b730f8bSJeremy L Thompson int CeedMatrixMatrixMultiply(Ceed ceed, const CeedScalar *mat_A, const CeedScalar *mat_B, CeedScalar *mat_C, CeedInt m, CeedInt n, CeedInt kk) { 5462b730f8bSJeremy L Thompson for (CeedInt i = 0; i < m; i++) { 5477a982d89SJeremy L. Thompson for (CeedInt j = 0; j < n; j++) { 5487a982d89SJeremy L. Thompson CeedScalar sum = 0; 5492b730f8bSJeremy L Thompson for (CeedInt k = 0; k < kk; k++) sum += mat_A[k + i * kk] * mat_B[j + k * n]; 550d1d35e2fSjeremylt mat_C[j + i * n] = sum; 5517a982d89SJeremy L. Thompson } 5522b730f8bSJeremy L Thompson } 553e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 5547a982d89SJeremy L. Thompson } 5557a982d89SJeremy L. Thompson 556ba59ac12SSebastian Grimberg /** 557ba59ac12SSebastian Grimberg @brief Return QR Factorization of a matrix 558ba59ac12SSebastian Grimberg 559ba59ac12SSebastian Grimberg @param[in] ceed Ceed context for error handling 560ba59ac12SSebastian Grimberg @param[in,out] mat Row-major matrix to be factorized in place 561ba59ac12SSebastian Grimberg @param[in,out] tau Vector of length m of scaling factors 562ba59ac12SSebastian Grimberg @param[in] m Number of rows 563ba59ac12SSebastian Grimberg @param[in] n Number of columns 564ba59ac12SSebastian Grimberg 565ba59ac12SSebastian Grimberg @return An error code: 0 - success, otherwise - failure 566ba59ac12SSebastian Grimberg 567ba59ac12SSebastian Grimberg @ref Utility 568ba59ac12SSebastian Grimberg **/ 569ba59ac12SSebastian Grimberg int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau, CeedInt m, CeedInt n) { 570ba59ac12SSebastian Grimberg CeedScalar v[m]; 571ba59ac12SSebastian Grimberg 572ba59ac12SSebastian Grimberg // Check matrix shape 5736574a04fSJeremy L Thompson CeedCheck(n <= m, ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute QR factorization with n > m"); 574ba59ac12SSebastian Grimberg 575ba59ac12SSebastian Grimberg for (CeedInt i = 0; i < n; i++) { 576ba59ac12SSebastian Grimberg if (i >= m - 1) { // last row of matrix, no reflection needed 577ba59ac12SSebastian Grimberg tau[i] = 0.; 578ba59ac12SSebastian Grimberg break; 579ba59ac12SSebastian Grimberg } 580ba59ac12SSebastian Grimberg // Calculate Householder vector, magnitude 581ba59ac12SSebastian Grimberg CeedScalar sigma = 0.0; 582ba59ac12SSebastian Grimberg v[i] = mat[i + n * i]; 583ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < m; j++) { 584ba59ac12SSebastian Grimberg v[j] = mat[i + n * j]; 585ba59ac12SSebastian Grimberg sigma += v[j] * v[j]; 586ba59ac12SSebastian Grimberg } 587ba59ac12SSebastian Grimberg CeedScalar norm = sqrt(v[i] * v[i] + sigma); // norm of v[i:m] 588ba59ac12SSebastian Grimberg CeedScalar R_ii = -copysign(norm, v[i]); 589ba59ac12SSebastian Grimberg v[i] -= R_ii; 590ba59ac12SSebastian Grimberg // norm of v[i:m] after modification above and scaling below 591ba59ac12SSebastian Grimberg // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 592ba59ac12SSebastian Grimberg // tau = 2 / (norm*norm) 593ba59ac12SSebastian Grimberg tau[i] = 2 * v[i] * v[i] / (v[i] * v[i] + sigma); 594ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < m; j++) v[j] /= v[i]; 595ba59ac12SSebastian Grimberg 596ba59ac12SSebastian Grimberg // Apply Householder reflector to lower right panel 597ba59ac12SSebastian Grimberg CeedHouseholderReflect(&mat[i * n + i + 1], &v[i], tau[i], m - i, n - i - 1, n, 1); 598ba59ac12SSebastian Grimberg // Save v 599ba59ac12SSebastian Grimberg mat[i + n * i] = R_ii; 600ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < m; j++) mat[i + n * j] = v[j]; 601ba59ac12SSebastian Grimberg } 602ba59ac12SSebastian Grimberg return CEED_ERROR_SUCCESS; 603ba59ac12SSebastian Grimberg } 604ba59ac12SSebastian Grimberg 605ba59ac12SSebastian Grimberg /** 606ba59ac12SSebastian Grimberg @brief Apply Householder Q matrix 607ba59ac12SSebastian Grimberg 608ba59ac12SSebastian Grimberg Compute mat_A = mat_Q mat_A, where mat_Q is mxm and mat_A is mxn. 609ba59ac12SSebastian Grimberg 610ba59ac12SSebastian Grimberg @param[in,out] mat_A Matrix to apply Householder Q to, in place 611ba59ac12SSebastian Grimberg @param[in] mat_Q Householder Q matrix 612ba59ac12SSebastian Grimberg @param[in] tau Householder scaling factors 613ba59ac12SSebastian Grimberg @param[in] t_mode Transpose mode for application 614ba59ac12SSebastian Grimberg @param[in] m Number of rows in A 615ba59ac12SSebastian Grimberg @param[in] n Number of columns in A 616ba59ac12SSebastian Grimberg @param[in] k Number of elementary reflectors in Q, k<m 617ba59ac12SSebastian Grimberg @param[in] row Row stride in A 618ba59ac12SSebastian Grimberg @param[in] col Col stride in A 619ba59ac12SSebastian Grimberg 620ba59ac12SSebastian Grimberg @return An error code: 0 - success, otherwise - failure 621ba59ac12SSebastian Grimberg 622c4e3f59bSSebastian Grimberg @ref Utility 623ba59ac12SSebastian Grimberg **/ 624ba59ac12SSebastian Grimberg int CeedHouseholderApplyQ(CeedScalar *mat_A, const CeedScalar *mat_Q, const CeedScalar *tau, CeedTransposeMode t_mode, CeedInt m, CeedInt n, 625ba59ac12SSebastian Grimberg CeedInt k, CeedInt row, CeedInt col) { 626ba59ac12SSebastian Grimberg CeedScalar *v; 627ba59ac12SSebastian Grimberg CeedCall(CeedMalloc(m, &v)); 628ba59ac12SSebastian Grimberg for (CeedInt ii = 0; ii < k; ii++) { 629ba59ac12SSebastian Grimberg CeedInt i = t_mode == CEED_TRANSPOSE ? ii : k - 1 - ii; 630ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < m; j++) v[j] = mat_Q[j * k + i]; 631ba59ac12SSebastian Grimberg // Apply Householder reflector (I - tau v v^T) collo_grad_1d^T 632ba59ac12SSebastian Grimberg CeedCall(CeedHouseholderReflect(&mat_A[i * row], &v[i], tau[i], m - i, n, row, col)); 633ba59ac12SSebastian Grimberg } 634ba59ac12SSebastian Grimberg CeedCall(CeedFree(&v)); 635ba59ac12SSebastian Grimberg return CEED_ERROR_SUCCESS; 636ba59ac12SSebastian Grimberg } 637ba59ac12SSebastian Grimberg 638ba59ac12SSebastian Grimberg /** 639ba59ac12SSebastian Grimberg @brief Return symmetric Schur decomposition of the symmetric matrix mat via symmetric QR factorization 640ba59ac12SSebastian Grimberg 641ba59ac12SSebastian Grimberg @param[in] ceed Ceed context for error handling 642ba59ac12SSebastian Grimberg @param[in,out] mat Row-major matrix to be factorized in place 643ba59ac12SSebastian Grimberg @param[out] lambda Vector of length n of eigenvalues 644ba59ac12SSebastian Grimberg @param[in] n Number of rows/columns 645ba59ac12SSebastian Grimberg 646ba59ac12SSebastian Grimberg @return An error code: 0 - success, otherwise - failure 647ba59ac12SSebastian Grimberg 648ba59ac12SSebastian Grimberg @ref Utility 649ba59ac12SSebastian Grimberg **/ 6502c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOff 6512c2ea1dbSJeremy L Thompson int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat, CeedScalar *lambda, CeedInt n) { 652ba59ac12SSebastian Grimberg // Check bounds for clang-tidy 6536574a04fSJeremy L Thompson CeedCheck(n > 1, ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute symmetric Schur decomposition of scalars"); 654ba59ac12SSebastian Grimberg 655ba59ac12SSebastian Grimberg CeedScalar v[n - 1], tau[n - 1], mat_T[n * n]; 656ba59ac12SSebastian Grimberg 657ba59ac12SSebastian Grimberg // Copy mat to mat_T and set mat to I 658ba59ac12SSebastian Grimberg memcpy(mat_T, mat, n * n * sizeof(mat[0])); 659ba59ac12SSebastian Grimberg for (CeedInt i = 0; i < n; i++) { 660ba59ac12SSebastian Grimberg for (CeedInt j = 0; j < n; j++) mat[j + n * i] = (i == j) ? 1 : 0; 661ba59ac12SSebastian Grimberg } 662ba59ac12SSebastian Grimberg 663ba59ac12SSebastian Grimberg // Reduce to tridiagonal 664ba59ac12SSebastian Grimberg for (CeedInt i = 0; i < n - 1; i++) { 665ba59ac12SSebastian Grimberg // Calculate Householder vector, magnitude 666ba59ac12SSebastian Grimberg CeedScalar sigma = 0.0; 667ba59ac12SSebastian Grimberg v[i] = mat_T[i + n * (i + 1)]; 668ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < n - 1; j++) { 669ba59ac12SSebastian Grimberg v[j] = mat_T[i + n * (j + 1)]; 670ba59ac12SSebastian Grimberg sigma += v[j] * v[j]; 671ba59ac12SSebastian Grimberg } 672ba59ac12SSebastian Grimberg CeedScalar norm = sqrt(v[i] * v[i] + sigma); // norm of v[i:n-1] 673ba59ac12SSebastian Grimberg CeedScalar R_ii = -copysign(norm, v[i]); 674ba59ac12SSebastian Grimberg v[i] -= R_ii; 675ba59ac12SSebastian Grimberg // norm of v[i:m] after modification above and scaling below 676ba59ac12SSebastian Grimberg // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 677ba59ac12SSebastian Grimberg // tau = 2 / (norm*norm) 678ba59ac12SSebastian Grimberg tau[i] = i == n - 2 ? 2 : 2 * v[i] * v[i] / (v[i] * v[i] + sigma); 679ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < n - 1; j++) v[j] /= v[i]; 680ba59ac12SSebastian Grimberg 681ba59ac12SSebastian Grimberg // Update sub and super diagonal 682ba59ac12SSebastian Grimberg for (CeedInt j = i + 2; j < n; j++) { 683ba59ac12SSebastian Grimberg mat_T[i + n * j] = 0; 684ba59ac12SSebastian Grimberg mat_T[j + n * i] = 0; 685ba59ac12SSebastian Grimberg } 686ba59ac12SSebastian Grimberg // Apply symmetric Householder reflector to lower right panel 687ba59ac12SSebastian Grimberg CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1); 688ba59ac12SSebastian Grimberg CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), 1, n); 689ba59ac12SSebastian Grimberg 690ba59ac12SSebastian Grimberg // Save v 691ba59ac12SSebastian Grimberg mat_T[i + n * (i + 1)] = R_ii; 692ba59ac12SSebastian Grimberg mat_T[(i + 1) + n * i] = R_ii; 693ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < n - 1; j++) { 694ba59ac12SSebastian Grimberg mat_T[i + n * (j + 1)] = v[j]; 695ba59ac12SSebastian Grimberg } 696ba59ac12SSebastian Grimberg } 697ba59ac12SSebastian Grimberg // Backwards accumulation of Q 698ba59ac12SSebastian Grimberg for (CeedInt i = n - 2; i >= 0; i--) { 699ba59ac12SSebastian Grimberg if (tau[i] > 0.0) { 700ba59ac12SSebastian Grimberg v[i] = 1; 701ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < n - 1; j++) { 702ba59ac12SSebastian Grimberg v[j] = mat_T[i + n * (j + 1)]; 703ba59ac12SSebastian Grimberg mat_T[i + n * (j + 1)] = 0; 704ba59ac12SSebastian Grimberg } 705ba59ac12SSebastian Grimberg CeedHouseholderReflect(&mat[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1); 706ba59ac12SSebastian Grimberg } 707ba59ac12SSebastian Grimberg } 708ba59ac12SSebastian Grimberg 709ba59ac12SSebastian Grimberg // Reduce sub and super diagonal 710ba59ac12SSebastian Grimberg CeedInt p = 0, q = 0, itr = 0, max_itr = n * n * n * n; 711ba59ac12SSebastian Grimberg CeedScalar tol = CEED_EPSILON; 712ba59ac12SSebastian Grimberg 713ba59ac12SSebastian Grimberg while (itr < max_itr) { 714ba59ac12SSebastian Grimberg // Update p, q, size of reduced portions of diagonal 715ba59ac12SSebastian Grimberg p = 0; 716ba59ac12SSebastian Grimberg q = 0; 717ba59ac12SSebastian Grimberg for (CeedInt i = n - 2; i >= 0; i--) { 718ba59ac12SSebastian Grimberg if (fabs(mat_T[i + n * (i + 1)]) < tol) q += 1; 719ba59ac12SSebastian Grimberg else break; 720ba59ac12SSebastian Grimberg } 721ba59ac12SSebastian Grimberg for (CeedInt i = 0; i < n - q - 1; i++) { 722ba59ac12SSebastian Grimberg if (fabs(mat_T[i + n * (i + 1)]) < tol) p += 1; 723ba59ac12SSebastian Grimberg else break; 724ba59ac12SSebastian Grimberg } 725ba59ac12SSebastian Grimberg if (q == n - 1) break; // Finished reducing 726ba59ac12SSebastian Grimberg 727ba59ac12SSebastian Grimberg // Reduce tridiagonal portion 728ba59ac12SSebastian Grimberg CeedScalar t_nn = mat_T[(n - 1 - q) + n * (n - 1 - q)], t_nnm1 = mat_T[(n - 2 - q) + n * (n - 1 - q)]; 729ba59ac12SSebastian Grimberg CeedScalar d = (mat_T[(n - 2 - q) + n * (n - 2 - q)] - t_nn) / 2; 730ba59ac12SSebastian Grimberg CeedScalar mu = t_nn - t_nnm1 * t_nnm1 / (d + copysign(sqrt(d * d + t_nnm1 * t_nnm1), d)); 731ba59ac12SSebastian Grimberg CeedScalar x = mat_T[p + n * p] - mu; 732ba59ac12SSebastian Grimberg CeedScalar z = mat_T[p + n * (p + 1)]; 733ba59ac12SSebastian Grimberg for (CeedInt k = p; k < n - q - 1; k++) { 734ba59ac12SSebastian Grimberg // Compute Givens rotation 735ba59ac12SSebastian Grimberg CeedScalar c = 1, s = 0; 736ba59ac12SSebastian Grimberg if (fabs(z) > tol) { 737ba59ac12SSebastian Grimberg if (fabs(z) > fabs(x)) { 738ba59ac12SSebastian Grimberg CeedScalar tau = -x / z; 739ba59ac12SSebastian Grimberg s = 1 / sqrt(1 + tau * tau), c = s * tau; 740ba59ac12SSebastian Grimberg } else { 741ba59ac12SSebastian Grimberg CeedScalar tau = -z / x; 742ba59ac12SSebastian Grimberg c = 1 / sqrt(1 + tau * tau), s = c * tau; 743ba59ac12SSebastian Grimberg } 744ba59ac12SSebastian Grimberg } 745ba59ac12SSebastian Grimberg 746ba59ac12SSebastian Grimberg // Apply Givens rotation to T 747ba59ac12SSebastian Grimberg CeedGivensRotation(mat_T, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n); 748ba59ac12SSebastian Grimberg CeedGivensRotation(mat_T, c, s, CEED_TRANSPOSE, k, k + 1, n, n); 749ba59ac12SSebastian Grimberg 750ba59ac12SSebastian Grimberg // Apply Givens rotation to Q 751ba59ac12SSebastian Grimberg CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n); 752ba59ac12SSebastian Grimberg 753ba59ac12SSebastian Grimberg // Update x, z 754ba59ac12SSebastian Grimberg if (k < n - q - 2) { 755ba59ac12SSebastian Grimberg x = mat_T[k + n * (k + 1)]; 756ba59ac12SSebastian Grimberg z = mat_T[k + n * (k + 2)]; 757ba59ac12SSebastian Grimberg } 758ba59ac12SSebastian Grimberg } 759ba59ac12SSebastian Grimberg itr++; 760ba59ac12SSebastian Grimberg } 761ba59ac12SSebastian Grimberg 762ba59ac12SSebastian Grimberg // Save eigenvalues 763ba59ac12SSebastian Grimberg for (CeedInt i = 0; i < n; i++) lambda[i] = mat_T[i + n * i]; 764ba59ac12SSebastian Grimberg 765ba59ac12SSebastian Grimberg // Check convergence 7666574a04fSJeremy L Thompson CeedCheck(itr < max_itr || q > n, ceed, CEED_ERROR_MINOR, "Symmetric QR failed to converge"); 767ba59ac12SSebastian Grimberg return CEED_ERROR_SUCCESS; 768ba59ac12SSebastian Grimberg } 7692c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOn 770ba59ac12SSebastian Grimberg 771ba59ac12SSebastian Grimberg /** 772ba59ac12SSebastian Grimberg @brief Return Simultaneous Diagonalization of two matrices. 773ba59ac12SSebastian Grimberg 774ba59ac12SSebastian Grimberg This solves the generalized eigenvalue problem A x = lambda B x, where A and B are symmetric and B is positive definite. 775ba59ac12SSebastian Grimberg We generate the matrix X and vector Lambda such that X^T A X = Lambda and X^T B X = I. 776ba59ac12SSebastian Grimberg This is equivalent to the LAPACK routine 'sygv' with TYPE = 1. 777ba59ac12SSebastian Grimberg 778ba59ac12SSebastian Grimberg @param[in] ceed Ceed context for error handling 779ba59ac12SSebastian Grimberg @param[in] mat_A Row-major matrix to be factorized with eigenvalues 780ba59ac12SSebastian Grimberg @param[in] mat_B Row-major matrix to be factorized to identity 781ba59ac12SSebastian Grimberg @param[out] mat_X Row-major orthogonal matrix 782ba59ac12SSebastian Grimberg @param[out] lambda Vector of length n of generalized eigenvalues 783ba59ac12SSebastian Grimberg @param[in] n Number of rows/columns 784ba59ac12SSebastian Grimberg 785ba59ac12SSebastian Grimberg @return An error code: 0 - success, otherwise - failure 786ba59ac12SSebastian Grimberg 787ba59ac12SSebastian Grimberg @ref Utility 788ba59ac12SSebastian Grimberg **/ 7892c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOff 7902c2ea1dbSJeremy L Thompson int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *mat_A, CeedScalar *mat_B, CeedScalar *mat_X, CeedScalar *lambda, CeedInt n) { 791ba59ac12SSebastian Grimberg CeedScalar *mat_C, *mat_G, *vec_D; 792ba59ac12SSebastian Grimberg CeedCall(CeedCalloc(n * n, &mat_C)); 793ba59ac12SSebastian Grimberg CeedCall(CeedCalloc(n * n, &mat_G)); 794ba59ac12SSebastian Grimberg CeedCall(CeedCalloc(n, &vec_D)); 795ba59ac12SSebastian Grimberg 796ba59ac12SSebastian Grimberg // Compute B = G D G^T 797ba59ac12SSebastian Grimberg memcpy(mat_G, mat_B, n * n * sizeof(mat_B[0])); 798ba59ac12SSebastian Grimberg CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_G, vec_D, n)); 799ba59ac12SSebastian Grimberg 800ba59ac12SSebastian Grimberg // Sort eigenvalues 801ba59ac12SSebastian Grimberg for (CeedInt i = n - 1; i >= 0; i--) { 802ba59ac12SSebastian Grimberg for (CeedInt j = 0; j < i; j++) { 803ba59ac12SSebastian Grimberg if (fabs(vec_D[j]) > fabs(vec_D[j + 1])) { 804ba59ac12SSebastian Grimberg CeedScalar temp; 805ba59ac12SSebastian Grimberg temp = vec_D[j]; 806ba59ac12SSebastian Grimberg vec_D[j] = vec_D[j + 1]; 807ba59ac12SSebastian Grimberg vec_D[j + 1] = temp; 808ba59ac12SSebastian Grimberg for (CeedInt k = 0; k < n; k++) { 809ba59ac12SSebastian Grimberg temp = mat_G[k * n + j]; 810ba59ac12SSebastian Grimberg mat_G[k * n + j] = mat_G[k * n + j + 1]; 811ba59ac12SSebastian Grimberg mat_G[k * n + j + 1] = temp; 812ba59ac12SSebastian Grimberg } 813ba59ac12SSebastian Grimberg } 814ba59ac12SSebastian Grimberg } 815ba59ac12SSebastian Grimberg } 816ba59ac12SSebastian Grimberg 817ba59ac12SSebastian Grimberg // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T 818ba59ac12SSebastian Grimberg // = D^-1/2 G^T A G D^-1/2 819ba59ac12SSebastian Grimberg // -- D = D^-1/2 820ba59ac12SSebastian Grimberg for (CeedInt i = 0; i < n; i++) vec_D[i] = 1. / sqrt(vec_D[i]); 821ba59ac12SSebastian Grimberg // -- G = G D^-1/2 822ba59ac12SSebastian Grimberg // -- C = D^-1/2 G^T 823ba59ac12SSebastian Grimberg for (CeedInt i = 0; i < n; i++) { 824ba59ac12SSebastian Grimberg for (CeedInt j = 0; j < n; j++) { 825ba59ac12SSebastian Grimberg mat_G[i * n + j] *= vec_D[j]; 826ba59ac12SSebastian Grimberg mat_C[j * n + i] = mat_G[i * n + j]; 827ba59ac12SSebastian Grimberg } 828ba59ac12SSebastian Grimberg } 829ba59ac12SSebastian Grimberg // -- X = (D^-1/2 G^T) A 830ba59ac12SSebastian Grimberg CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_C, (const CeedScalar *)mat_A, mat_X, n, n, n)); 831ba59ac12SSebastian Grimberg // -- C = (D^-1/2 G^T A) (G D^-1/2) 832ba59ac12SSebastian Grimberg CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_X, (const CeedScalar *)mat_G, mat_C, n, n, n)); 833ba59ac12SSebastian Grimberg 834ba59ac12SSebastian Grimberg // Compute Q^T C Q = lambda 835ba59ac12SSebastian Grimberg CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_C, lambda, n)); 836ba59ac12SSebastian Grimberg 837ba59ac12SSebastian Grimberg // Sort eigenvalues 838ba59ac12SSebastian Grimberg for (CeedInt i = n - 1; i >= 0; i--) { 839ba59ac12SSebastian Grimberg for (CeedInt j = 0; j < i; j++) { 840ba59ac12SSebastian Grimberg if (fabs(lambda[j]) > fabs(lambda[j + 1])) { 841ba59ac12SSebastian Grimberg CeedScalar temp; 842ba59ac12SSebastian Grimberg temp = lambda[j]; 843ba59ac12SSebastian Grimberg lambda[j] = lambda[j + 1]; 844ba59ac12SSebastian Grimberg lambda[j + 1] = temp; 845ba59ac12SSebastian Grimberg for (CeedInt k = 0; k < n; k++) { 846ba59ac12SSebastian Grimberg temp = mat_C[k * n + j]; 847ba59ac12SSebastian Grimberg mat_C[k * n + j] = mat_C[k * n + j + 1]; 848ba59ac12SSebastian Grimberg mat_C[k * n + j + 1] = temp; 849ba59ac12SSebastian Grimberg } 850ba59ac12SSebastian Grimberg } 851ba59ac12SSebastian Grimberg } 852ba59ac12SSebastian Grimberg } 853ba59ac12SSebastian Grimberg 854ba59ac12SSebastian Grimberg // Set X = (G D^1/2)^-T Q 855ba59ac12SSebastian Grimberg // = G D^-1/2 Q 856ba59ac12SSebastian Grimberg CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_G, (const CeedScalar *)mat_C, mat_X, n, n, n)); 857ba59ac12SSebastian Grimberg 858ba59ac12SSebastian Grimberg // Cleanup 859ba59ac12SSebastian Grimberg CeedCall(CeedFree(&mat_C)); 860ba59ac12SSebastian Grimberg CeedCall(CeedFree(&mat_G)); 861ba59ac12SSebastian Grimberg CeedCall(CeedFree(&vec_D)); 862ba59ac12SSebastian Grimberg return CEED_ERROR_SUCCESS; 863ba59ac12SSebastian Grimberg } 8642c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOn 865ba59ac12SSebastian Grimberg 8667a982d89SJeremy L. Thompson /// @} 8677a982d89SJeremy L. Thompson 8687a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 8697a982d89SJeremy L. Thompson /// CeedBasis Public API 8707a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 8717a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser 872d7b241e6Sjeremylt /// @{ 873d7b241e6Sjeremylt 874b11c1e72Sjeremylt /** 875ba59ac12SSebastian Grimberg @brief Create a tensor-product basis for H^1 discretizations 876b11c1e72Sjeremylt 877ea61e9acSJeremy L Thompson @param[in] ceed Ceed object where the CeedBasis will be created 878ea61e9acSJeremy L Thompson @param[in] dim Topological dimension 879ea61e9acSJeremy L Thompson @param[in] num_comp Number of field components (1 for scalar fields) 880ea61e9acSJeremy L Thompson @param[in] P_1d Number of nodes in one dimension 881ea61e9acSJeremy L Thompson @param[in] Q_1d Number of quadrature points in one dimension 882ea61e9acSJeremy L Thompson @param[in] interp_1d Row-major (Q_1d * P_1d) matrix expressing the values of nodal basis functions at quadrature points 883ea61e9acSJeremy L Thompson @param[in] grad_1d Row-major (Q_1d * P_1d) matrix expressing derivatives of nodal basis functions at quadrature points 884ea61e9acSJeremy 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] 885ea61e9acSJeremy L Thompson @param[in] q_weight_1d Array of length Q_1d holding the quadrature weights on the reference element 886ea61e9acSJeremy L Thompson @param[out] basis Address of the variable where the newly created CeedBasis will be stored. 887b11c1e72Sjeremylt 888b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 889dfdf5a53Sjeremylt 8907a982d89SJeremy L. Thompson @ref User 891b11c1e72Sjeremylt **/ 8922b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, 8932b730f8bSJeremy L Thompson const CeedScalar *grad_1d, const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis *basis) { 8945fe0d4faSjeremylt if (!ceed->BasisCreateTensorH1) { 8955fe0d4faSjeremylt Ceed delegate; 8966574a04fSJeremy L Thompson 8972b730f8bSJeremy L Thompson CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis")); 8986574a04fSJeremy L Thompson CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisCreateTensorH1"); 8992b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateTensorH1(delegate, dim, num_comp, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis)); 900e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 9015fe0d4faSjeremylt } 902e15f9bd0SJeremy L Thompson 9036574a04fSJeremy L Thompson CeedCheck(dim > 0, ceed, CEED_ERROR_DIMENSION, "Basis dimension must be a positive value"); 9046574a04fSJeremy L Thompson CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component"); 9056574a04fSJeremy L Thompson CeedCheck(P_1d > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node"); 9066574a04fSJeremy L Thompson CeedCheck(Q_1d > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point"); 907227444bfSJeremy L Thompson 9082b730f8bSJeremy L Thompson CeedElemTopology topo = dim == 1 ? CEED_TOPOLOGY_LINE : dim == 2 ? CEED_TOPOLOGY_QUAD : CEED_TOPOLOGY_HEX; 909e15f9bd0SJeremy L Thompson 9102b730f8bSJeremy L Thompson CeedCall(CeedCalloc(1, basis)); 911db002c03SJeremy L Thompson CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed)); 912d1d35e2fSjeremylt (*basis)->ref_count = 1; 9136402da51SJeremy L Thompson (*basis)->is_tensor_basis = true; 914d7b241e6Sjeremylt (*basis)->dim = dim; 915d99fa3c5SJeremy L Thompson (*basis)->topo = topo; 916d1d35e2fSjeremylt (*basis)->num_comp = num_comp; 917d1d35e2fSjeremylt (*basis)->P_1d = P_1d; 918d1d35e2fSjeremylt (*basis)->Q_1d = Q_1d; 919d1d35e2fSjeremylt (*basis)->P = CeedIntPow(P_1d, dim); 920d1d35e2fSjeremylt (*basis)->Q = CeedIntPow(Q_1d, dim); 921c4e3f59bSSebastian Grimberg (*basis)->fe_space = CEED_FE_SPACE_H1; 9222b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q_1d, &(*basis)->q_ref_1d)); 9232b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q_1d, &(*basis)->q_weight_1d)); 924ff3a0f91SJeremy L Thompson if (q_ref_1d) memcpy((*basis)->q_ref_1d, q_ref_1d, Q_1d * sizeof(q_ref_1d[0])); 9252b730f8bSJeremy L Thompson if (q_weight_1d) memcpy((*basis)->q_weight_1d, q_weight_1d, Q_1d * sizeof(q_weight_1d[0])); 9262b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->interp_1d)); 9272b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->grad_1d)); 9282b730f8bSJeremy L Thompson if (interp_1d) memcpy((*basis)->interp_1d, interp_1d, Q_1d * P_1d * sizeof(interp_1d[0])); 929ff3a0f91SJeremy L Thompson if (grad_1d) memcpy((*basis)->grad_1d, grad_1d, Q_1d * P_1d * sizeof(grad_1d[0])); 9302b730f8bSJeremy L Thompson CeedCall(ceed->BasisCreateTensorH1(dim, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, *basis)); 931e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 932d7b241e6Sjeremylt } 933d7b241e6Sjeremylt 934b11c1e72Sjeremylt /** 93595bb1877Svaleriabarra @brief Create a tensor-product Lagrange basis 936b11c1e72Sjeremylt 937ea61e9acSJeremy L Thompson @param[in] ceed Ceed object where the CeedBasis will be created 938ea61e9acSJeremy L Thompson @param[in] dim Topological dimension of element 939ea61e9acSJeremy L Thompson @param[in] num_comp Number of field components (1 for scalar fields) 940ea61e9acSJeremy L Thompson @param[in] P Number of Gauss-Lobatto nodes in one dimension. 941ea61e9acSJeremy L Thompson The polynomial degree of the resulting Q_k element is k=P-1. 942ea61e9acSJeremy L Thompson @param[in] Q Number of quadrature points in one dimension. 943ea61e9acSJeremy L Thompson @param[in] quad_mode Distribution of the Q quadrature points (affects order of accuracy for the quadrature) 944ea61e9acSJeremy L Thompson @param[out] basis Address of the variable where the newly created CeedBasis will be stored. 945b11c1e72Sjeremylt 946b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 947dfdf5a53Sjeremylt 9487a982d89SJeremy L. Thompson @ref User 949b11c1e72Sjeremylt **/ 9502b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P, CeedInt Q, CeedQuadMode quad_mode, CeedBasis *basis) { 951d7b241e6Sjeremylt // Allocate 952c8c3fa7dSJeremy L Thompson int ierr = CEED_ERROR_SUCCESS; 9532b730f8bSJeremy L Thompson CeedScalar c1, c2, c3, c4, dx, *nodes, *interp_1d, *grad_1d, *q_ref_1d, *q_weight_1d; 9544d537eeaSYohann 9556574a04fSJeremy L Thompson CeedCheck(dim > 0, ceed, CEED_ERROR_DIMENSION, "Basis dimension must be a positive value"); 9566574a04fSJeremy L Thompson CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component"); 9576574a04fSJeremy L Thompson CeedCheck(P > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node"); 9586574a04fSJeremy L Thompson CeedCheck(Q > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point"); 959227444bfSJeremy L Thompson 960e15f9bd0SJeremy L Thompson // Get Nodes and Weights 9612b730f8bSJeremy L Thompson CeedCall(CeedCalloc(P * Q, &interp_1d)); 9622b730f8bSJeremy L Thompson CeedCall(CeedCalloc(P * Q, &grad_1d)); 9632b730f8bSJeremy L Thompson CeedCall(CeedCalloc(P, &nodes)); 9642b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q, &q_ref_1d)); 9652b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q, &q_weight_1d)); 9662b730f8bSJeremy L Thompson if (CeedLobattoQuadrature(P, nodes, NULL) != CEED_ERROR_SUCCESS) goto cleanup; 967d1d35e2fSjeremylt switch (quad_mode) { 968d7b241e6Sjeremylt case CEED_GAUSS: 969d1d35e2fSjeremylt ierr = CeedGaussQuadrature(Q, q_ref_1d, q_weight_1d); 970d7b241e6Sjeremylt break; 971d7b241e6Sjeremylt case CEED_GAUSS_LOBATTO: 972d1d35e2fSjeremylt ierr = CeedLobattoQuadrature(Q, q_ref_1d, q_weight_1d); 973d7b241e6Sjeremylt break; 974d7b241e6Sjeremylt } 9752b730f8bSJeremy L Thompson if (ierr != CEED_ERROR_SUCCESS) goto cleanup; 976e15f9bd0SJeremy L Thompson 977d7b241e6Sjeremylt // Build B, D matrix 978d7b241e6Sjeremylt // Fornberg, 1998 979c8c3fa7dSJeremy L Thompson for (CeedInt i = 0; i < Q; i++) { 980d7b241e6Sjeremylt c1 = 1.0; 981d1d35e2fSjeremylt c3 = nodes[0] - q_ref_1d[i]; 982d1d35e2fSjeremylt interp_1d[i * P + 0] = 1.0; 983c8c3fa7dSJeremy L Thompson for (CeedInt j = 1; j < P; j++) { 984d7b241e6Sjeremylt c2 = 1.0; 985d7b241e6Sjeremylt c4 = c3; 986d1d35e2fSjeremylt c3 = nodes[j] - q_ref_1d[i]; 987c8c3fa7dSJeremy L Thompson for (CeedInt k = 0; k < j; k++) { 988d7b241e6Sjeremylt dx = nodes[j] - nodes[k]; 989d7b241e6Sjeremylt c2 *= dx; 990d7b241e6Sjeremylt if (k == j - 1) { 991d1d35e2fSjeremylt grad_1d[i * P + j] = c1 * (interp_1d[i * P + k] - c4 * grad_1d[i * P + k]) / c2; 992d1d35e2fSjeremylt interp_1d[i * P + j] = -c1 * c4 * interp_1d[i * P + k] / c2; 993d7b241e6Sjeremylt } 994d1d35e2fSjeremylt grad_1d[i * P + k] = (c3 * grad_1d[i * P + k] - interp_1d[i * P + k]) / dx; 995d1d35e2fSjeremylt interp_1d[i * P + k] = c3 * interp_1d[i * P + k] / dx; 996d7b241e6Sjeremylt } 997d7b241e6Sjeremylt c1 = c2; 998d7b241e6Sjeremylt } 999d7b241e6Sjeremylt } 10009ac7b42eSJeremy L Thompson // Pass to CeedBasisCreateTensorH1 10012b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P, Q, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis)); 1002e15f9bd0SJeremy L Thompson cleanup: 10032b730f8bSJeremy L Thompson CeedCall(CeedFree(&interp_1d)); 10042b730f8bSJeremy L Thompson CeedCall(CeedFree(&grad_1d)); 10052b730f8bSJeremy L Thompson CeedCall(CeedFree(&nodes)); 10062b730f8bSJeremy L Thompson CeedCall(CeedFree(&q_ref_1d)); 10072b730f8bSJeremy L Thompson CeedCall(CeedFree(&q_weight_1d)); 1008e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1009d7b241e6Sjeremylt } 1010d7b241e6Sjeremylt 1011b11c1e72Sjeremylt /** 1012ba59ac12SSebastian Grimberg @brief Create a non tensor-product basis for H^1 discretizations 1013a8de75f0Sjeremylt 1014ea61e9acSJeremy L Thompson @param[in] ceed Ceed object where the CeedBasis will be created 1015ea61e9acSJeremy L Thompson @param[in] topo Topology of element, e.g. hypercube, simplex, ect 1016ea61e9acSJeremy L Thompson @param[in] num_comp Number of field components (1 for scalar fields) 1017ea61e9acSJeremy L Thompson @param[in] num_nodes Total number of nodes 1018ea61e9acSJeremy L Thompson @param[in] num_qpts Total number of quadrature points 1019ea61e9acSJeremy L Thompson @param[in] interp Row-major (num_qpts * num_nodes) matrix expressing the values of nodal basis functions at quadrature points 1020c4e3f59bSSebastian Grimberg @param[in] grad Row-major (dim * num_qpts * num_nodes) matrix expressing derivatives of nodal basis functions at quadrature points 10219fe083eeSJeremy L Thompson @param[in] q_ref Array of length num_qpts * dim holding the locations of quadrature points on the reference element 1022ea61e9acSJeremy L Thompson @param[in] q_weight Array of length num_qpts holding the quadrature weights on the reference element 1023ea61e9acSJeremy L Thompson @param[out] basis Address of the variable where the newly created CeedBasis will be stored. 1024a8de75f0Sjeremylt 1025a8de75f0Sjeremylt @return An error code: 0 - success, otherwise - failure 1026a8de75f0Sjeremylt 10277a982d89SJeremy L. Thompson @ref User 1028a8de75f0Sjeremylt **/ 10292b730f8bSJeremy L Thompson int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 10302b730f8bSJeremy L Thompson const CeedScalar *grad, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) { 1031d1d35e2fSjeremylt CeedInt P = num_nodes, Q = num_qpts, dim = 0; 1032a8de75f0Sjeremylt 10335fe0d4faSjeremylt if (!ceed->BasisCreateH1) { 10345fe0d4faSjeremylt Ceed delegate; 10356574a04fSJeremy L Thompson 10362b730f8bSJeremy L Thompson CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis")); 10376574a04fSJeremy L Thompson CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisCreateH1"); 10382b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateH1(delegate, topo, num_comp, num_nodes, num_qpts, interp, grad, q_ref, q_weight, basis)); 1039e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 10405fe0d4faSjeremylt } 10415fe0d4faSjeremylt 10426574a04fSJeremy L Thompson CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component"); 10436574a04fSJeremy L Thompson CeedCheck(num_nodes > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node"); 10446574a04fSJeremy L Thompson CeedCheck(num_qpts > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point"); 1045227444bfSJeremy L Thompson 10462b730f8bSJeremy L Thompson CeedCall(CeedBasisGetTopologyDimension(topo, &dim)); 1047a8de75f0Sjeremylt 1048db002c03SJeremy L Thompson CeedCall(CeedCalloc(1, basis)); 1049db002c03SJeremy L Thompson CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed)); 1050d1d35e2fSjeremylt (*basis)->ref_count = 1; 10516402da51SJeremy L Thompson (*basis)->is_tensor_basis = false; 1052a8de75f0Sjeremylt (*basis)->dim = dim; 1053d99fa3c5SJeremy L Thompson (*basis)->topo = topo; 1054d1d35e2fSjeremylt (*basis)->num_comp = num_comp; 1055a8de75f0Sjeremylt (*basis)->P = P; 1056a8de75f0Sjeremylt (*basis)->Q = Q; 1057c4e3f59bSSebastian Grimberg (*basis)->fe_space = CEED_FE_SPACE_H1; 10582b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q * dim, &(*basis)->q_ref_1d)); 10592b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q, &(*basis)->q_weight_1d)); 1060ff3a0f91SJeremy L Thompson if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0])); 1061ff3a0f91SJeremy L Thompson if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0])); 10622b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q * P, &(*basis)->interp)); 10632b730f8bSJeremy L Thompson CeedCall(CeedCalloc(dim * Q * P, &(*basis)->grad)); 1064ff3a0f91SJeremy L Thompson if (interp) memcpy((*basis)->interp, interp, Q * P * sizeof(interp[0])); 1065ff3a0f91SJeremy L Thompson if (grad) memcpy((*basis)->grad, grad, dim * Q * P * sizeof(grad[0])); 10662b730f8bSJeremy L Thompson CeedCall(ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, q_ref, q_weight, *basis)); 1067e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1068a8de75f0Sjeremylt } 1069a8de75f0Sjeremylt 1070a8de75f0Sjeremylt /** 1071859c15bbSJames Wright @brief Create a non tensor-product basis for \f$H(\mathrm{div})\f$ discretizations 107250c301a5SRezgar Shakeri 1073ea61e9acSJeremy L Thompson @param[in] ceed Ceed object where the CeedBasis will be created 1074ea61e9acSJeremy 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 1075ea61e9acSJeremy L Thompson @param[in] num_comp Number of components (usually 1 for vectors in H(div) bases) 1076ea61e9acSJeremy L Thompson @param[in] num_nodes Total number of nodes (dofs per element) 1077ea61e9acSJeremy L Thompson @param[in] num_qpts Total number of quadrature points 1078c4e3f59bSSebastian Grimberg @param[in] interp Row-major (dim * num_qpts * num_nodes) matrix expressing the values of basis functions at quadrature points 1079c4e3f59bSSebastian Grimberg @param[in] div Row-major (num_qpts * num_nodes) matrix expressing divergence of basis functions at quadrature points 10809fe083eeSJeremy L Thompson @param[in] q_ref Array of length num_qpts * dim holding the locations of quadrature points on the reference element 1081ea61e9acSJeremy L Thompson @param[in] q_weight Array of length num_qpts holding the quadrature weights on the reference element 1082ea61e9acSJeremy L Thompson @param[out] basis Address of the variable where the newly created CeedBasis will be stored. 108350c301a5SRezgar Shakeri 108450c301a5SRezgar Shakeri @return An error code: 0 - success, otherwise - failure 108550c301a5SRezgar Shakeri 108650c301a5SRezgar Shakeri @ref User 108750c301a5SRezgar Shakeri **/ 10882b730f8bSJeremy L Thompson int CeedBasisCreateHdiv(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 10892b730f8bSJeremy L Thompson const CeedScalar *div, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) { 109050c301a5SRezgar Shakeri CeedInt Q = num_qpts, P = num_nodes, dim = 0; 1091c4e3f59bSSebastian Grimberg 109250c301a5SRezgar Shakeri if (!ceed->BasisCreateHdiv) { 109350c301a5SRezgar Shakeri Ceed delegate; 10946574a04fSJeremy L Thompson 10952b730f8bSJeremy L Thompson CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis")); 10966574a04fSJeremy L Thompson CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateHdiv"); 10972b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateHdiv(delegate, topo, num_comp, num_nodes, num_qpts, interp, div, q_ref, q_weight, basis)); 109850c301a5SRezgar Shakeri return CEED_ERROR_SUCCESS; 109950c301a5SRezgar Shakeri } 110050c301a5SRezgar Shakeri 11016574a04fSJeremy L Thompson CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component"); 11026574a04fSJeremy L Thompson CeedCheck(num_nodes > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node"); 11036574a04fSJeremy L Thompson CeedCheck(num_qpts > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point"); 1104227444bfSJeremy L Thompson 1105c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetTopologyDimension(topo, &dim)); 1106c4e3f59bSSebastian Grimberg 1107db002c03SJeremy L Thompson CeedCall(CeedCalloc(1, basis)); 1108db002c03SJeremy L Thompson CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed)); 110950c301a5SRezgar Shakeri (*basis)->ref_count = 1; 11106402da51SJeremy L Thompson (*basis)->is_tensor_basis = false; 111150c301a5SRezgar Shakeri (*basis)->dim = dim; 111250c301a5SRezgar Shakeri (*basis)->topo = topo; 111350c301a5SRezgar Shakeri (*basis)->num_comp = num_comp; 111450c301a5SRezgar Shakeri (*basis)->P = P; 111550c301a5SRezgar Shakeri (*basis)->Q = Q; 1116c4e3f59bSSebastian Grimberg (*basis)->fe_space = CEED_FE_SPACE_HDIV; 11172b730f8bSJeremy L Thompson CeedCall(CeedMalloc(Q * dim, &(*basis)->q_ref_1d)); 11182b730f8bSJeremy L Thompson CeedCall(CeedMalloc(Q, &(*basis)->q_weight_1d)); 111950c301a5SRezgar Shakeri if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0])); 112050c301a5SRezgar Shakeri if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0])); 11212b730f8bSJeremy L Thompson CeedCall(CeedMalloc(dim * Q * P, &(*basis)->interp)); 11222b730f8bSJeremy L Thompson CeedCall(CeedMalloc(Q * P, &(*basis)->div)); 112350c301a5SRezgar Shakeri if (interp) memcpy((*basis)->interp, interp, dim * Q * P * sizeof(interp[0])); 112450c301a5SRezgar Shakeri if (div) memcpy((*basis)->div, div, Q * P * sizeof(div[0])); 11252b730f8bSJeremy L Thompson CeedCall(ceed->BasisCreateHdiv(topo, dim, P, Q, interp, div, q_ref, q_weight, *basis)); 112650c301a5SRezgar Shakeri return CEED_ERROR_SUCCESS; 112750c301a5SRezgar Shakeri } 112850c301a5SRezgar Shakeri 112950c301a5SRezgar Shakeri /** 11304385fb7fSSebastian Grimberg @brief Create a non tensor-product basis for \f$H(\mathrm{curl})\f$ discretizations 1131c4e3f59bSSebastian Grimberg 1132c4e3f59bSSebastian Grimberg @param[in] ceed Ceed object where the CeedBasis will be created 1133c4e3f59bSSebastian Grimberg @param[in] topo Topology of element (`CEED_TOPOLOGY_QUAD`, `CEED_TOPOLOGY_PRISM`, etc.), dimension of which is used in some array sizes below 1134c4e3f59bSSebastian Grimberg @param[in] num_comp Number of components (usually 1 for vectors in H(curl) bases) 1135c4e3f59bSSebastian Grimberg @param[in] num_nodes Total number of nodes (dofs per element) 1136c4e3f59bSSebastian Grimberg @param[in] num_qpts Total number of quadrature points 1137c4e3f59bSSebastian Grimberg @param[in] interp Row-major (dim * num_qpts * num_nodes) matrix expressing the values of basis functions at quadrature points 1138c4e3f59bSSebastian Grimberg @param[in] curl Row-major (curl_comp * num_qpts * num_nodes, curl_comp = 1 if dim < 3 else dim) matrix expressing curl of basis functions at 1139c4e3f59bSSebastian Grimberg quadrature points 1140c4e3f59bSSebastian Grimberg @param[in] q_ref Array of length num_qpts * dim holding the locations of quadrature points on the reference element 1141c4e3f59bSSebastian Grimberg @param[in] q_weight Array of length num_qpts holding the quadrature weights on the reference element 1142c4e3f59bSSebastian Grimberg @param[out] basis Address of the variable where the newly created CeedBasis will be stored. 1143c4e3f59bSSebastian Grimberg 1144c4e3f59bSSebastian Grimberg @return An error code: 0 - success, otherwise - failure 1145c4e3f59bSSebastian Grimberg 1146c4e3f59bSSebastian Grimberg @ref User 1147c4e3f59bSSebastian Grimberg **/ 1148c4e3f59bSSebastian Grimberg int CeedBasisCreateHcurl(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 1149c4e3f59bSSebastian Grimberg const CeedScalar *curl, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) { 1150c4e3f59bSSebastian Grimberg CeedInt Q = num_qpts, P = num_nodes, dim = 0, curl_comp = 0; 1151c4e3f59bSSebastian Grimberg 1152c4e3f59bSSebastian Grimberg if (!ceed->BasisCreateHdiv) { 1153c4e3f59bSSebastian Grimberg Ceed delegate; 11546574a04fSJeremy L Thompson 1155c4e3f59bSSebastian Grimberg CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis")); 11566574a04fSJeremy L Thompson CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateHcurl"); 1157c4e3f59bSSebastian Grimberg CeedCall(CeedBasisCreateHcurl(delegate, topo, num_comp, num_nodes, num_qpts, interp, curl, q_ref, q_weight, basis)); 1158c4e3f59bSSebastian Grimberg return CEED_ERROR_SUCCESS; 1159c4e3f59bSSebastian Grimberg } 1160c4e3f59bSSebastian Grimberg 11616574a04fSJeremy L Thompson CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component"); 11626574a04fSJeremy L Thompson CeedCheck(num_nodes > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node"); 11636574a04fSJeremy L Thompson CeedCheck(num_qpts > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point"); 1164c4e3f59bSSebastian Grimberg 1165c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetTopologyDimension(topo, &dim)); 1166c4e3f59bSSebastian Grimberg curl_comp = (dim < 3) ? 1 : dim; 1167c4e3f59bSSebastian Grimberg 1168db002c03SJeremy L Thompson CeedCall(CeedCalloc(1, basis)); 1169db002c03SJeremy L Thompson CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed)); 1170c4e3f59bSSebastian Grimberg (*basis)->ref_count = 1; 11716402da51SJeremy L Thompson (*basis)->is_tensor_basis = false; 1172c4e3f59bSSebastian Grimberg (*basis)->dim = dim; 1173c4e3f59bSSebastian Grimberg (*basis)->topo = topo; 1174c4e3f59bSSebastian Grimberg (*basis)->num_comp = num_comp; 1175c4e3f59bSSebastian Grimberg (*basis)->P = P; 1176c4e3f59bSSebastian Grimberg (*basis)->Q = Q; 1177c4e3f59bSSebastian Grimberg (*basis)->fe_space = CEED_FE_SPACE_HCURL; 1178c4e3f59bSSebastian Grimberg CeedCall(CeedMalloc(Q * dim, &(*basis)->q_ref_1d)); 1179c4e3f59bSSebastian Grimberg CeedCall(CeedMalloc(Q, &(*basis)->q_weight_1d)); 1180c4e3f59bSSebastian Grimberg if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0])); 1181c4e3f59bSSebastian Grimberg if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0])); 1182c4e3f59bSSebastian Grimberg CeedCall(CeedMalloc(dim * Q * P, &(*basis)->interp)); 1183c4e3f59bSSebastian Grimberg CeedCall(CeedMalloc(curl_comp * Q * P, &(*basis)->curl)); 1184c4e3f59bSSebastian Grimberg if (interp) memcpy((*basis)->interp, interp, dim * Q * P * sizeof(interp[0])); 1185c4e3f59bSSebastian Grimberg if (curl) memcpy((*basis)->curl, curl, curl_comp * Q * P * sizeof(curl[0])); 1186c4e3f59bSSebastian Grimberg CeedCall(ceed->BasisCreateHcurl(topo, dim, P, Q, interp, curl, q_ref, q_weight, *basis)); 1187c4e3f59bSSebastian Grimberg return CEED_ERROR_SUCCESS; 1188c4e3f59bSSebastian Grimberg } 1189c4e3f59bSSebastian Grimberg 1190c4e3f59bSSebastian Grimberg /** 1191ea61e9acSJeremy L Thompson @brief Create a CeedBasis for projection from the nodes of `basis_from` to the nodes of `basis_to`. 1192ba59ac12SSebastian Grimberg 11939fd66db6SSebastian Grimberg Only `CEED_EVAL_INTERP` will be valid for the new basis, `basis_project`. 11949fd66db6SSebastian Grimberg For H^1 spaces, `CEED_EVAL_GRAD` will also be valid. 1195de05fbb2SSebastian Grimberg The interpolation is given by `interp_project = interp_to^+ * interp_from`, where the pseudoinverse `interp_to^+` is given by QR 11969fd66db6SSebastian Grimberg factorization. 11979fd66db6SSebastian Grimberg The gradient (for the H^1 case) is given by `grad_project = interp_to^+ * grad_from`. 119815ad3917SSebastian Grimberg 119915ad3917SSebastian Grimberg Note: `basis_from` and `basis_to` must have compatible quadrature spaces. 120015ad3917SSebastian Grimberg 12019fd66db6SSebastian Grimberg Note: `basis_project` will have the same number of components as `basis_from`, regardless of the number of components that `basis_to` has. 12029fd66db6SSebastian Grimberg If `basis_from` has 3 components and `basis_to` has 5 components, then `basis_project` will have 3 components. 1203f113e5dcSJeremy L Thompson 1204f113e5dcSJeremy L Thompson @param[in] basis_from CeedBasis to prolong from 1205446e7af4SJeremy L Thompson @param[in] basis_to CeedBasis to prolong to 1206ea61e9acSJeremy L Thompson @param[out] basis_project Address of the variable where the newly created CeedBasis will be stored. 1207f113e5dcSJeremy L Thompson 1208f113e5dcSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 1209f113e5dcSJeremy L Thompson 1210f113e5dcSJeremy L Thompson @ref User 1211f113e5dcSJeremy L Thompson **/ 12122b730f8bSJeremy L Thompson int CeedBasisCreateProjection(CeedBasis basis_from, CeedBasis basis_to, CeedBasis *basis_project) { 1213f113e5dcSJeremy L Thompson Ceed ceed; 12142b730f8bSJeremy L Thompson CeedCall(CeedBasisGetCeed(basis_to, &ceed)); 1215f113e5dcSJeremy L Thompson 1216ecc88aebSJeremy L Thompson // Create projection matrix 121714556e63SJeremy L Thompson CeedScalar *interp_project, *grad_project; 12182b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateProjectionMatrices(basis_from, basis_to, &interp_project, &grad_project)); 1219f113e5dcSJeremy L Thompson 1220f113e5dcSJeremy L Thompson // Build basis 1221f113e5dcSJeremy L Thompson bool is_tensor; 1222f113e5dcSJeremy L Thompson CeedInt dim, num_comp; 122314556e63SJeremy L Thompson CeedScalar *q_ref, *q_weight; 12242b730f8bSJeremy L Thompson CeedCall(CeedBasisIsTensor(basis_to, &is_tensor)); 12252b730f8bSJeremy L Thompson CeedCall(CeedBasisGetDimension(basis_to, &dim)); 12262b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumComponents(basis_from, &num_comp)); 1227f113e5dcSJeremy L Thompson if (is_tensor) { 1228f113e5dcSJeremy L Thompson CeedInt P_1d_to, P_1d_from; 12292b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_1d_from)); 12302b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_1d_to)); 12312b730f8bSJeremy L Thompson CeedCall(CeedCalloc(P_1d_to, &q_ref)); 12322b730f8bSJeremy L Thompson CeedCall(CeedCalloc(P_1d_to, &q_weight)); 12332b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P_1d_from, P_1d_to, interp_project, grad_project, q_ref, q_weight, basis_project)); 1234f113e5dcSJeremy L Thompson } else { 1235de05fbb2SSebastian Grimberg // Even if basis_to and basis_from are not H1, the resulting basis is H1 for interpolation to work 1236f113e5dcSJeremy L Thompson CeedElemTopology topo; 12372b730f8bSJeremy L Thompson CeedCall(CeedBasisGetTopology(basis_to, &topo)); 1238f113e5dcSJeremy L Thompson CeedInt num_nodes_to, num_nodes_from; 12392b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis_from, &num_nodes_from)); 12402b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis_to, &num_nodes_to)); 12412b730f8bSJeremy L Thompson CeedCall(CeedCalloc(num_nodes_to * dim, &q_ref)); 12422b730f8bSJeremy L Thompson CeedCall(CeedCalloc(num_nodes_to, &q_weight)); 12432b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateH1(ceed, topo, num_comp, num_nodes_from, num_nodes_to, interp_project, grad_project, q_ref, q_weight, basis_project)); 1244f113e5dcSJeremy L Thompson } 1245f113e5dcSJeremy L Thompson 1246f113e5dcSJeremy L Thompson // Cleanup 12472b730f8bSJeremy L Thompson CeedCall(CeedFree(&interp_project)); 12482b730f8bSJeremy L Thompson CeedCall(CeedFree(&grad_project)); 12492b730f8bSJeremy L Thompson CeedCall(CeedFree(&q_ref)); 12502b730f8bSJeremy L Thompson CeedCall(CeedFree(&q_weight)); 1251f113e5dcSJeremy L Thompson 1252f113e5dcSJeremy L Thompson return CEED_ERROR_SUCCESS; 1253f113e5dcSJeremy L Thompson } 1254f113e5dcSJeremy L Thompson 1255f113e5dcSJeremy L Thompson /** 1256ea61e9acSJeremy L Thompson @brief Copy the pointer to a CeedBasis. 12579560d06aSjeremylt 1258512bb800SJeremy 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. 1259512bb800SJeremy L Thompson This CeedBasis will be destroyed if `basis_copy` is the only reference to this CeedBasis. 1260ea61e9acSJeremy L Thompson 1261ea61e9acSJeremy L Thompson @param[in] basis CeedBasis to copy reference to 1262ea61e9acSJeremy L Thompson @param[in,out] basis_copy Variable to store copied reference 12639560d06aSjeremylt 12649560d06aSjeremylt @return An error code: 0 - success, otherwise - failure 12659560d06aSjeremylt 12669560d06aSjeremylt @ref User 12679560d06aSjeremylt **/ 12689560d06aSjeremylt int CeedBasisReferenceCopy(CeedBasis basis, CeedBasis *basis_copy) { 1269393ac2cdSJeremy L Thompson if (basis != CEED_BASIS_COLLOCATED) CeedCall(CeedBasisReference(basis)); 12702b730f8bSJeremy L Thompson CeedCall(CeedBasisDestroy(basis_copy)); 12719560d06aSjeremylt *basis_copy = basis; 12729560d06aSjeremylt return CEED_ERROR_SUCCESS; 12739560d06aSjeremylt } 12749560d06aSjeremylt 12759560d06aSjeremylt /** 12767a982d89SJeremy L. Thompson @brief View a CeedBasis 12777a982d89SJeremy L. Thompson 1278ea61e9acSJeremy L Thompson @param[in] basis CeedBasis to view 1279ea61e9acSJeremy L Thompson @param[in] stream Stream to view to, e.g., stdout 12807a982d89SJeremy L. Thompson 12817a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 12827a982d89SJeremy L. Thompson 12837a982d89SJeremy L. Thompson @ref User 12847a982d89SJeremy L. Thompson **/ 12857a982d89SJeremy L. Thompson int CeedBasisView(CeedBasis basis, FILE *stream) { 128650c301a5SRezgar Shakeri CeedElemTopology topo = basis->topo; 1287c4e3f59bSSebastian Grimberg CeedFESpace fe_space = basis->fe_space; 1288c4e3f59bSSebastian Grimberg CeedInt q_comp = 0; 12892b730f8bSJeremy L Thompson 129050c301a5SRezgar Shakeri // Print FE space and element topology of the basis 12916402da51SJeremy L Thompson if (basis->is_tensor_basis) { 1292c4e3f59bSSebastian Grimberg fprintf(stream, "CeedBasis (%s on a %s element): dim=%" CeedInt_FMT " P=%" CeedInt_FMT " Q=%" CeedInt_FMT "\n", CeedFESpaces[fe_space], 12932b730f8bSJeremy L Thompson CeedElemTopologies[topo], basis->dim, basis->P_1d, basis->Q_1d); 129450c301a5SRezgar Shakeri } else { 1295c4e3f59bSSebastian Grimberg fprintf(stream, "CeedBasis (%s on a %s element): dim=%" CeedInt_FMT " P=%" CeedInt_FMT " Q=%" CeedInt_FMT "\n", CeedFESpaces[fe_space], 12962b730f8bSJeremy L Thompson CeedElemTopologies[topo], basis->dim, basis->P, basis->Q); 129750c301a5SRezgar Shakeri } 1298ea61e9acSJeremy L Thompson // Print quadrature data, interpolation/gradient/divergence/curl of the basis 12996402da51SJeremy L Thompson if (basis->is_tensor_basis) { // tensor basis 13002b730f8bSJeremy L Thompson CeedCall(CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_ref_1d, stream)); 13012b730f8bSJeremy L Thompson CeedCall(CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_weight_1d, stream)); 13022b730f8bSJeremy L Thompson CeedCall(CeedScalarView("interp1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, basis->interp_1d, stream)); 13032b730f8bSJeremy L Thompson CeedCall(CeedScalarView("grad1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, basis->grad_1d, stream)); 130450c301a5SRezgar Shakeri } else { // non-tensor basis 13052b730f8bSJeremy L Thompson CeedCall(CeedScalarView("qref", "\t% 12.8f", 1, basis->Q * basis->dim, basis->q_ref_1d, stream)); 13062b730f8bSJeremy L Thompson CeedCall(CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->q_weight_1d, stream)); 1307c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp)); 1308c4e3f59bSSebastian Grimberg CeedCall(CeedScalarView("interp", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->interp, stream)); 130950c301a5SRezgar Shakeri if (basis->grad) { 1310c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp)); 1311c4e3f59bSSebastian Grimberg CeedCall(CeedScalarView("grad", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->grad, stream)); 13127a982d89SJeremy L. Thompson } 131350c301a5SRezgar Shakeri if (basis->div) { 1314c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_DIV, &q_comp)); 1315c4e3f59bSSebastian Grimberg CeedCall(CeedScalarView("div", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->div, stream)); 1316c4e3f59bSSebastian Grimberg } 1317c4e3f59bSSebastian Grimberg if (basis->curl) { 1318c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_CURL, &q_comp)); 1319c4e3f59bSSebastian Grimberg CeedCall(CeedScalarView("curl", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->curl, stream)); 132050c301a5SRezgar Shakeri } 132150c301a5SRezgar Shakeri } 1322e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 13237a982d89SJeremy L. Thompson } 13247a982d89SJeremy L. Thompson 13257a982d89SJeremy L. Thompson /** 13267a982d89SJeremy L. Thompson @brief Apply basis evaluation from nodes to quadrature points or vice versa 13277a982d89SJeremy L. Thompson 1328ea61e9acSJeremy L Thompson @param[in] basis CeedBasis to evaluate 1329ea61e9acSJeremy L Thompson @param[in] num_elem The number of elements to apply the basis evaluation to; 1330ea61e9acSJeremy L Thompson the backend will specify the ordering in CeedElemRestrictionCreateBlocked() 1331ea61e9acSJeremy L Thompson @param[in] t_mode \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature points; 1332ea61e9acSJeremy L Thompson \ref CEED_TRANSPOSE to apply the transpose, mapping from quadrature points to nodes 1333ea61e9acSJeremy L Thompson @param[in] eval_mode \ref CEED_EVAL_NONE to use values directly, 13347a982d89SJeremy L. Thompson \ref CEED_EVAL_INTERP to use interpolated values, 13357a982d89SJeremy L. Thompson \ref CEED_EVAL_GRAD to use gradients, 1336c4e3f59bSSebastian Grimberg \ref CEED_EVAL_DIV to use divergence, 1337c4e3f59bSSebastian Grimberg \ref CEED_EVAL_CURL to use curl, 13387a982d89SJeremy L. Thompson \ref CEED_EVAL_WEIGHT to use quadrature weights. 13397a982d89SJeremy L. Thompson @param[in] u Input CeedVector 13407a982d89SJeremy L. Thompson @param[out] v Output CeedVector 13417a982d89SJeremy L. Thompson 13427a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 13437a982d89SJeremy L. Thompson 13447a982d89SJeremy L. Thompson @ref User 13457a982d89SJeremy L. Thompson **/ 13462b730f8bSJeremy L Thompson int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, CeedVector v) { 13471f9221feSJeremy L Thompson CeedSize u_length = 0, v_length; 1348c4e3f59bSSebastian Grimberg CeedInt dim, num_comp, q_comp, num_nodes, num_qpts; 13492b730f8bSJeremy L Thompson CeedCall(CeedBasisGetDimension(basis, &dim)); 13502b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumComponents(basis, &num_comp)); 1351c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp)); 13522b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis, &num_nodes)); 13532b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts)); 13542b730f8bSJeremy L Thompson CeedCall(CeedVectorGetLength(v, &v_length)); 1355c8c3fa7dSJeremy L Thompson if (u) CeedCall(CeedVectorGetLength(u, &u_length)); 13567a982d89SJeremy L. Thompson 13576574a04fSJeremy L Thompson CeedCheck(basis->Apply, basis->ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisApply"); 1358e15f9bd0SJeremy L Thompson 1359e15f9bd0SJeremy L Thompson // Check compatibility of topological and geometrical dimensions 13606574a04fSJeremy L Thompson CeedCheck((t_mode == CEED_TRANSPOSE && v_length % num_nodes == 0 && u_length % num_qpts == 0) || 13616574a04fSJeremy L Thompson (t_mode == CEED_NOTRANSPOSE && u_length % num_nodes == 0 && v_length % num_qpts == 0), 13626574a04fSJeremy L Thompson basis->ceed, CEED_ERROR_DIMENSION, "Length of input/output vectors incompatible with basis dimensions"); 13637a982d89SJeremy L. Thompson 1364e15f9bd0SJeremy L Thompson // Check vector lengths to prevent out of bounds issues 13656574a04fSJeremy L Thompson bool good_dims = true; 1366d1d35e2fSjeremylt switch (eval_mode) { 1367e15f9bd0SJeremy L Thompson case CEED_EVAL_NONE: 13682b730f8bSJeremy L Thompson case CEED_EVAL_INTERP: 13692b730f8bSJeremy L Thompson case CEED_EVAL_GRAD: 1370c4e3f59bSSebastian Grimberg case CEED_EVAL_DIV: 1371c4e3f59bSSebastian Grimberg case CEED_EVAL_CURL: 13726574a04fSJeremy L Thompson good_dims = 13736574a04fSJeremy L Thompson ((t_mode == CEED_TRANSPOSE && u_length >= num_elem * num_comp * num_qpts * q_comp && v_length >= num_elem * num_comp * num_nodes) || 13746574a04fSJeremy L Thompson (t_mode == CEED_NOTRANSPOSE && v_length >= num_elem * num_qpts * num_comp * q_comp && u_length >= num_elem * num_comp * num_nodes)); 1375e15f9bd0SJeremy L Thompson break; 1376e15f9bd0SJeremy L Thompson case CEED_EVAL_WEIGHT: 13776574a04fSJeremy L Thompson good_dims = v_length >= num_elem * num_qpts; 1378e15f9bd0SJeremy L Thompson break; 1379e15f9bd0SJeremy L Thompson } 13806574a04fSJeremy L Thompson CeedCheck(good_dims, basis->ceed, CEED_ERROR_DIMENSION, "Input/output vectors too short for basis and evaluation mode"); 1381e15f9bd0SJeremy L Thompson 13822b730f8bSJeremy L Thompson CeedCall(basis->Apply(basis, num_elem, t_mode, eval_mode, u, v)); 1383e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 13847a982d89SJeremy L. Thompson } 13857a982d89SJeremy L. Thompson 13867a982d89SJeremy L. Thompson /** 1387c8c3fa7dSJeremy L Thompson @brief Apply basis evaluation from nodes to arbitrary points 1388c8c3fa7dSJeremy L Thompson 1389c8c3fa7dSJeremy L Thompson @param[in] basis CeedBasis to evaluate 1390c8c3fa7dSJeremy L Thompson @param[in] num_points The number of points to apply the basis evaluation to 1391c8c3fa7dSJeremy L Thompson @param[in] t_mode \ref CEED_NOTRANSPOSE to evaluate from nodes to points; 1392c8c3fa7dSJeremy L Thompson \ref CEED_TRANSPOSE to apply the transpose, mapping from points to nodes 1393c8c3fa7dSJeremy L Thompson @param[in] eval_mode \ref CEED_EVAL_INTERP to use interpolated values, 1394c8c3fa7dSJeremy L Thompson \ref CEED_EVAL_GRAD to use gradients 1395c8c3fa7dSJeremy L Thompson @param[in] x_ref CeedVector holding reference coordinates of each point 1396c8c3fa7dSJeremy L Thompson @param[in] u Input CeedVector, of length `num_nodes * num_comp` for `CEED_NOTRANSPOSE` 1397c8c3fa7dSJeremy L Thompson @param[out] v Output CeedVector, of length `num_points * num_q_comp` for `CEED_NOTRANSPOSE` with `CEED_EVAL_INTERP` 1398c8c3fa7dSJeremy L Thompson 1399c8c3fa7dSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 1400c8c3fa7dSJeremy L Thompson 1401c8c3fa7dSJeremy L Thompson @ref User 1402c8c3fa7dSJeremy L Thompson **/ 1403c8c3fa7dSJeremy L Thompson int CeedBasisApplyAtPoints(CeedBasis basis, CeedInt num_points, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, 1404c8c3fa7dSJeremy L Thompson CeedVector v) { 1405c8c3fa7dSJeremy L Thompson CeedSize x_length = 0, u_length = 0, v_length; 1406c8c3fa7dSJeremy L Thompson CeedInt dim, num_comp, num_q_comp, num_nodes, P_1d = 1, Q_1d = 1; 1407c8c3fa7dSJeremy L Thompson 1408c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisGetDimension(basis, &dim)); 1409c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d)); 1410c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 1411c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisGetNumComponents(basis, &num_comp)); 1412c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &num_q_comp)); 1413c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis, &num_nodes)); 1414c8c3fa7dSJeremy L Thompson CeedCall(CeedVectorGetLength(x_ref, &x_length)); 1415c8c3fa7dSJeremy L Thompson CeedCall(CeedVectorGetLength(v, &v_length)); 1416c8c3fa7dSJeremy L Thompson CeedCall(CeedVectorGetLength(u, &u_length)); 1417c8c3fa7dSJeremy L Thompson 1418c8c3fa7dSJeremy L Thompson // Check compatibility of topological and geometrical dimensions 1419c8c3fa7dSJeremy L Thompson CeedCheck((t_mode == CEED_TRANSPOSE && v_length % num_nodes == 0) || (t_mode == CEED_NOTRANSPOSE && u_length % num_nodes == 0), basis->ceed, 1420c8c3fa7dSJeremy L Thompson CEED_ERROR_DIMENSION, "Length of input/output vectors incompatible with basis dimensions and number of points"); 1421c8c3fa7dSJeremy L Thompson 1422c8c3fa7dSJeremy L Thompson // Check compatibility coordinates vector 1423c8c3fa7dSJeremy L Thompson CeedCheck(x_length >= num_points * dim, basis->ceed, CEED_ERROR_DIMENSION, 1424c8c3fa7dSJeremy L Thompson "Length of reference coordinate vector incompatible with basis dimension and number of points"); 1425c8c3fa7dSJeremy L Thompson 1426c8c3fa7dSJeremy L Thompson // Check vector lengths to prevent out of bounds issues 1427c8c3fa7dSJeremy L Thompson bool good_dims = false; 1428c8c3fa7dSJeremy L Thompson switch (eval_mode) { 1429c8c3fa7dSJeremy L Thompson case CEED_EVAL_INTERP: 1430c8c3fa7dSJeremy L Thompson good_dims = ((t_mode == CEED_TRANSPOSE && (u_length >= num_points * num_q_comp || v_length >= num_nodes * num_comp)) || 1431c8c3fa7dSJeremy L Thompson (t_mode == CEED_NOTRANSPOSE && (v_length >= num_points * num_q_comp || u_length >= num_nodes * num_comp))); 1432c8c3fa7dSJeremy L Thompson break; 1433c8c3fa7dSJeremy L Thompson case CEED_EVAL_GRAD: 1434*edfbf3a6SJeremy L Thompson good_dims = ((t_mode == CEED_TRANSPOSE && (u_length >= num_points * num_q_comp * dim || v_length >= num_nodes * num_comp)) || 1435*edfbf3a6SJeremy L Thompson (t_mode == CEED_NOTRANSPOSE && (v_length >= num_points * num_q_comp * dim || u_length >= num_nodes * num_comp))); 1436*edfbf3a6SJeremy L Thompson break; 1437c8c3fa7dSJeremy L Thompson case CEED_EVAL_NONE: 1438c8c3fa7dSJeremy L Thompson case CEED_EVAL_WEIGHT: 1439c8c3fa7dSJeremy L Thompson case CEED_EVAL_DIV: 1440c8c3fa7dSJeremy L Thompson case CEED_EVAL_CURL: 1441c8c3fa7dSJeremy L Thompson // LCOV_EXCL_START 1442c8c3fa7dSJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_UNSUPPORTED, "Evaluation at arbitrary points not supported for %s", CeedEvalModes[eval_mode]); 1443c8c3fa7dSJeremy L Thompson // LCOV_EXCL_STOP 1444c8c3fa7dSJeremy L Thompson } 1445c8c3fa7dSJeremy L Thompson CeedCheck(good_dims, basis->ceed, CEED_ERROR_DIMENSION, "Input/output vectors too short for basis and evaluation mode"); 1446c8c3fa7dSJeremy L Thompson 1447c8c3fa7dSJeremy L Thompson // Backend method 1448c8c3fa7dSJeremy L Thompson if (basis->ApplyAtPoints) { 1449c8c3fa7dSJeremy L Thompson CeedCall(basis->ApplyAtPoints(basis, num_points, t_mode, eval_mode, x_ref, u, v)); 1450c8c3fa7dSJeremy L Thompson return CEED_ERROR_SUCCESS; 1451c8c3fa7dSJeremy L Thompson } 1452c8c3fa7dSJeremy L Thompson 1453c8c3fa7dSJeremy L Thompson // Default implementation 1454c8c3fa7dSJeremy L Thompson CeedCheck(basis->is_tensor_basis, basis->ceed, CEED_ERROR_UNSUPPORTED, "Evaluation at arbitrary points only supported for tensor product bases"); 1455*edfbf3a6SJeremy L Thompson CeedCheck(eval_mode == CEED_EVAL_INTERP || t_mode == CEED_NOTRANSPOSE, basis->ceed, CEED_ERROR_UNSUPPORTED, "%s evaluation only supported for %s", 1456*edfbf3a6SJeremy L Thompson CeedEvalModes[eval_mode], CeedTransposeModes[CEED_NOTRANSPOSE]); 1457c8c3fa7dSJeremy L Thompson if (!basis->basis_chebyshev) { 1458c8c3fa7dSJeremy L Thompson // Build matrix mapping from quadrature point values to Chebyshev coefficients 1459c8c3fa7dSJeremy L Thompson CeedScalar *tau, *C, *I, *chebyshev_coeffs_1d; 1460c8c3fa7dSJeremy L Thompson const CeedScalar *q_ref_1d; 1461c8c3fa7dSJeremy L Thompson 1462c8c3fa7dSJeremy L Thompson // Build coefficient matrix 1463c8c3fa7dSJeremy L Thompson // -- Note: Clang-tidy needs this check because it does not understand the is_tensor_basis check above 1464c8c3fa7dSJeremy L Thompson CeedCheck(P_1d > 0 && Q_1d > 0, basis->ceed, CEED_ERROR_INCOMPATIBLE, "Basis dimensions are malformed"); 1465c8c3fa7dSJeremy L Thompson CeedCall(CeedCalloc(Q_1d * Q_1d, &C)); 1466c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisGetQRef(basis, &q_ref_1d)); 1467c8c3fa7dSJeremy L Thompson for (CeedInt i = 0; i < Q_1d; i++) { 1468c8c3fa7dSJeremy L Thompson const CeedScalar x = q_ref_1d[i]; 1469c8c3fa7dSJeremy L Thompson 1470c8c3fa7dSJeremy L Thompson C[i * Q_1d + 0] = 1.0; 1471c8c3fa7dSJeremy L Thompson C[i * Q_1d + 1] = 2 * x; 1472c8c3fa7dSJeremy L Thompson for (CeedInt j = 2; j < Q_1d; j++) C[i * Q_1d + j] = 2 * x * C[i * Q_1d + j - 1] - C[i * Q_1d + j - 2]; 1473c8c3fa7dSJeremy L Thompson } 1474c8c3fa7dSJeremy L Thompson 1475c8c3fa7dSJeremy L Thompson // Inverse of coefficient matrix 1476c8c3fa7dSJeremy L Thompson CeedCall(CeedCalloc(Q_1d * Q_1d, &chebyshev_coeffs_1d)); 1477c8c3fa7dSJeremy L Thompson CeedCall(CeedCalloc(Q_1d * Q_1d, &I)); 1478c8c3fa7dSJeremy L Thompson CeedCall(CeedCalloc(Q_1d, &tau)); 1479c8c3fa7dSJeremy L Thompson // -- QR Factorization, C = Q R 1480c8c3fa7dSJeremy L Thompson CeedCall(CeedQRFactorization(basis->ceed, C, tau, Q_1d, Q_1d)); 1481c8c3fa7dSJeremy L Thompson // -- chebyshev_coeffs_1d = R_inv Q^T 1482c8c3fa7dSJeremy L Thompson for (CeedInt i = 0; i < Q_1d; i++) I[i * Q_1d + i] = 1.0; 1483c8c3fa7dSJeremy L Thompson // ---- Apply R_inv, chebyshev_coeffs_1d = I R_inv 1484c8c3fa7dSJeremy L Thompson for (CeedInt i = 0; i < Q_1d; i++) { // Row i 1485c8c3fa7dSJeremy L Thompson chebyshev_coeffs_1d[Q_1d * i] = I[Q_1d * i] / C[0]; 1486c8c3fa7dSJeremy L Thompson for (CeedInt j = 1; j < Q_1d; j++) { // Column j 1487c8c3fa7dSJeremy L Thompson chebyshev_coeffs_1d[j + Q_1d * i] = I[j + Q_1d * i]; 1488c8c3fa7dSJeremy L Thompson for (CeedInt k = 0; k < j; k++) chebyshev_coeffs_1d[j + Q_1d * i] -= C[j + Q_1d * k] * chebyshev_coeffs_1d[k + Q_1d * i]; 1489c8c3fa7dSJeremy L Thompson chebyshev_coeffs_1d[j + Q_1d * i] /= C[j + Q_1d * j]; 1490c8c3fa7dSJeremy L Thompson } 1491c8c3fa7dSJeremy L Thompson } 1492c8c3fa7dSJeremy L Thompson // ---- Apply Q^T, chebyshev_coeffs_1d = R_inv Q^T 1493c8c3fa7dSJeremy L Thompson CeedCall(CeedHouseholderApplyQ(chebyshev_coeffs_1d, C, tau, CEED_NOTRANSPOSE, Q_1d, Q_1d, Q_1d, 1, Q_1d)); 1494c8c3fa7dSJeremy L Thompson 1495c8c3fa7dSJeremy L Thompson // Build basis mapping from nodes to Chebyshev coefficients 1496c8c3fa7dSJeremy L Thompson CeedScalar *chebyshev_interp_1d, *chebyshev_grad_1d, *chebyshev_q_weight_1d; 1497c8c3fa7dSJeremy L Thompson const CeedScalar *interp_1d; 1498c8c3fa7dSJeremy L Thompson 1499c8c3fa7dSJeremy L Thompson CeedCall(CeedCalloc(Q_1d * Q_1d, &chebyshev_interp_1d)); 1500c8c3fa7dSJeremy L Thompson CeedCall(CeedCalloc(Q_1d * Q_1d, &chebyshev_grad_1d)); 1501c8c3fa7dSJeremy L Thompson CeedCall(CeedCalloc(Q_1d, &chebyshev_q_weight_1d)); 1502c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisGetInterp1D(basis, &interp_1d)); 1503c8c3fa7dSJeremy L Thompson CeedCall(CeedMatrixMatrixMultiply(basis->ceed, chebyshev_coeffs_1d, interp_1d, chebyshev_interp_1d, Q_1d, P_1d, Q_1d)); 1504c8c3fa7dSJeremy L Thompson 1505c8c3fa7dSJeremy L Thompson CeedCall(CeedVectorCreate(basis->ceed, num_comp * CeedIntPow(Q_1d, dim), &basis->vec_chebyshev)); 1506c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisCreateTensorH1(basis->ceed, dim, num_comp, Q_1d, Q_1d, chebyshev_interp_1d, chebyshev_grad_1d, q_ref_1d, chebyshev_q_weight_1d, 1507c8c3fa7dSJeremy L Thompson &basis->basis_chebyshev)); 1508c8c3fa7dSJeremy L Thompson 1509c8c3fa7dSJeremy L Thompson // Cleanup 1510c8c3fa7dSJeremy L Thompson CeedCall(CeedFree(&C)); 1511c8c3fa7dSJeremy L Thompson CeedCall(CeedFree(&chebyshev_coeffs_1d)); 1512c8c3fa7dSJeremy L Thompson CeedCall(CeedFree(&I)); 1513c8c3fa7dSJeremy L Thompson CeedCall(CeedFree(&tau)); 1514c8c3fa7dSJeremy L Thompson CeedCall(CeedFree(&chebyshev_interp_1d)); 1515c8c3fa7dSJeremy L Thompson CeedCall(CeedFree(&chebyshev_grad_1d)); 1516c8c3fa7dSJeremy L Thompson CeedCall(CeedFree(&chebyshev_q_weight_1d)); 1517c8c3fa7dSJeremy L Thompson } 1518c8c3fa7dSJeremy L Thompson 1519c8c3fa7dSJeremy L Thompson // Create TensorContract object if needed, such as a basis from the GPU backends 1520c8c3fa7dSJeremy L Thompson if (!basis->contract) { 1521c8c3fa7dSJeremy L Thompson Ceed ceed_ref; 1522c8c3fa7dSJeremy L Thompson CeedBasis basis_ref; 1523c8c3fa7dSJeremy L Thompson 1524c8c3fa7dSJeremy L Thompson CeedCall(CeedInit("/cpu/self", &ceed_ref)); 1525c8c3fa7dSJeremy L Thompson // Only need matching tensor contraction dimensions, any type of basis will work 1526c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisCreateTensorH1Lagrange(ceed_ref, dim, num_comp, Q_1d, Q_1d, CEED_GAUSS, &basis_ref)); 1527c8c3fa7dSJeremy L Thompson CeedCall(CeedTensorContractReference(basis_ref->contract)); 1528c8c3fa7dSJeremy L Thompson basis->contract = basis_ref->contract; 1529c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisDestroy(&basis_ref)); 1530c8c3fa7dSJeremy L Thompson CeedCall(CeedDestroy(&ceed_ref)); 1531c8c3fa7dSJeremy L Thompson } 1532c8c3fa7dSJeremy L Thompson 1533c8c3fa7dSJeremy L Thompson // Basis evaluation 1534c8c3fa7dSJeremy L Thompson switch (t_mode) { 1535c8c3fa7dSJeremy L Thompson case CEED_NOTRANSPOSE: { 1536c8c3fa7dSJeremy L Thompson // Nodes to arbitrary points 1537c8c3fa7dSJeremy L Thompson CeedScalar *v_array; 1538c8c3fa7dSJeremy L Thompson const CeedScalar *chebyshev_coeffs, *x_array_read; 1539c8c3fa7dSJeremy L Thompson 1540c8c3fa7dSJeremy L Thompson // -- Interpolate to Chebyshev coefficients 1541c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisApply(basis->basis_chebyshev, 1, CEED_NOTRANSPOSE, CEED_EVAL_INTERP, u, basis->vec_chebyshev)); 1542c8c3fa7dSJeremy L Thompson 1543c8c3fa7dSJeremy L Thompson // -- Evaluate Chebyshev polynomials at arbitrary points 1544c8c3fa7dSJeremy L Thompson CeedCall(CeedVectorGetArrayRead(basis->vec_chebyshev, CEED_MEM_HOST, &chebyshev_coeffs)); 1545c8c3fa7dSJeremy L Thompson CeedCall(CeedVectorGetArrayRead(x_ref, CEED_MEM_HOST, &x_array_read)); 1546c8c3fa7dSJeremy L Thompson CeedCall(CeedVectorGetArrayWrite(v, CEED_MEM_HOST, &v_array)); 1547*edfbf3a6SJeremy L Thompson switch (eval_mode) { 1548*edfbf3a6SJeremy L Thompson case CEED_EVAL_INTERP: { 1549c8c3fa7dSJeremy L Thompson CeedScalar tmp[2][num_comp * CeedIntPow(Q_1d, dim)], chebyshev_x[Q_1d]; 1550c8c3fa7dSJeremy L Thompson 1551c8c3fa7dSJeremy L Thompson // ---- Values at point 1552c8c3fa7dSJeremy L Thompson for (CeedInt p = 0; p < num_points; p++) { 1553c8c3fa7dSJeremy L Thompson CeedInt pre = num_comp * CeedIntPow(Q_1d, dim - 1), post = 1; 1554c8c3fa7dSJeremy L Thompson 1555c8c3fa7dSJeremy L Thompson for (CeedInt d = dim - 1; d >= 0; d--) { 1556c8c3fa7dSJeremy L Thompson // ------ Compute Chebyshev polynomial values 1557c8c3fa7dSJeremy L Thompson { 1558c8c3fa7dSJeremy L Thompson const CeedScalar x = x_array_read[p * dim + d]; 1559c8c3fa7dSJeremy L Thompson 1560c8c3fa7dSJeremy L Thompson chebyshev_x[0] = 1.0; 1561c8c3fa7dSJeremy L Thompson chebyshev_x[1] = 2 * x; 1562c8c3fa7dSJeremy L Thompson for (CeedInt j = 2; j < Q_1d; j++) chebyshev_x[j] = 2 * x * chebyshev_x[j - 1] - chebyshev_x[j - 2]; 1563c8c3fa7dSJeremy L Thompson } 1564c8c3fa7dSJeremy L Thompson // ------ Tensor contract 1565c8c3fa7dSJeremy L Thompson CeedCall(CeedTensorContractApply(basis->contract, pre, Q_1d, post, 1, chebyshev_x, t_mode, false, 1566c8c3fa7dSJeremy L Thompson d == (dim - 1) ? chebyshev_coeffs : tmp[d % 2], d == 0 ? &v_array[p * num_comp] : tmp[(d + 1) % 2])); 1567c8c3fa7dSJeremy L Thompson pre /= Q_1d; 1568c8c3fa7dSJeremy L Thompson post *= 1; 1569c8c3fa7dSJeremy L Thompson } 1570c8c3fa7dSJeremy L Thompson } 1571*edfbf3a6SJeremy L Thompson break; 1572*edfbf3a6SJeremy L Thompson } 1573*edfbf3a6SJeremy L Thompson case CEED_EVAL_GRAD: { 1574*edfbf3a6SJeremy L Thompson CeedScalar tmp[2][num_comp * CeedIntPow(Q_1d, dim)], chebyshev_x[Q_1d]; 1575*edfbf3a6SJeremy L Thompson 1576*edfbf3a6SJeremy L Thompson // ---- Values at point 1577*edfbf3a6SJeremy L Thompson for (CeedInt p = 0; p < num_points; p++) { 1578*edfbf3a6SJeremy L Thompson // Dim**2 contractions, apply grad when pass == dim 1579*edfbf3a6SJeremy L Thompson for (CeedInt pass = dim - 1; pass >= 0; pass--) { 1580*edfbf3a6SJeremy L Thompson CeedInt pre = num_comp * CeedIntPow(Q_1d, dim - 1), post = 1; 1581*edfbf3a6SJeremy L Thompson 1582*edfbf3a6SJeremy L Thompson for (CeedInt d = dim - 1; d >= 0; d--) { 1583*edfbf3a6SJeremy L Thompson const CeedScalar x = x_array_read[p * dim + d]; 1584*edfbf3a6SJeremy L Thompson 1585*edfbf3a6SJeremy L Thompson // ------ Compute Chebyshev polynomial values 1586*edfbf3a6SJeremy L Thompson if (pass == d) { 1587*edfbf3a6SJeremy L Thompson // chebyshev_temp holds the original series values 1588*edfbf3a6SJeremy L Thompson // chebyshev_x holds the derivative of the original series 1589*edfbf3a6SJeremy L Thompson CeedScalar chebyshev_temp[3]; 1590*edfbf3a6SJeremy L Thompson 1591*edfbf3a6SJeremy L Thompson chebyshev_temp[1] = 1.0; 1592*edfbf3a6SJeremy L Thompson chebyshev_temp[2] = 2 * x; 1593*edfbf3a6SJeremy L Thompson chebyshev_x[0] = 0.0; 1594*edfbf3a6SJeremy L Thompson chebyshev_x[1] = 2.0; 1595*edfbf3a6SJeremy L Thompson for (CeedInt j = 2; j < Q_1d; j++) { 1596*edfbf3a6SJeremy L Thompson chebyshev_temp[0] = chebyshev_temp[1]; 1597*edfbf3a6SJeremy L Thompson chebyshev_temp[1] = chebyshev_temp[2]; 1598*edfbf3a6SJeremy L Thompson chebyshev_temp[2] = 2 * x * chebyshev_temp[1] - chebyshev_temp[0]; 1599*edfbf3a6SJeremy L Thompson chebyshev_x[j] = 2 * x * chebyshev_x[j - 1] + 2 * chebyshev_temp[1] - chebyshev_x[j - 2]; 1600*edfbf3a6SJeremy L Thompson } 1601*edfbf3a6SJeremy L Thompson } else { 1602*edfbf3a6SJeremy L Thompson chebyshev_x[0] = 1.0; 1603*edfbf3a6SJeremy L Thompson chebyshev_x[1] = 2 * x; 1604*edfbf3a6SJeremy L Thompson for (CeedInt j = 2; j < Q_1d; j++) chebyshev_x[j] = 2 * x * chebyshev_x[j - 1] - chebyshev_x[j - 2]; 1605*edfbf3a6SJeremy L Thompson } 1606*edfbf3a6SJeremy L Thompson // ------ Tensor contract 1607*edfbf3a6SJeremy L Thompson CeedCall(CeedTensorContractApply(basis->contract, pre, Q_1d, post, 1, chebyshev_x, t_mode, false, 1608*edfbf3a6SJeremy L Thompson d == (dim - 1) ? chebyshev_coeffs : tmp[d % 2], 1609*edfbf3a6SJeremy L Thompson d == 0 ? &v_array[p * num_comp * dim + pass] : tmp[(d + 1) % 2])); 1610*edfbf3a6SJeremy L Thompson pre /= Q_1d; 1611*edfbf3a6SJeremy L Thompson post *= 1; 1612*edfbf3a6SJeremy L Thompson } 1613*edfbf3a6SJeremy L Thompson } 1614*edfbf3a6SJeremy L Thompson } 1615*edfbf3a6SJeremy L Thompson break; 1616*edfbf3a6SJeremy L Thompson } 1617*edfbf3a6SJeremy L Thompson default: 1618*edfbf3a6SJeremy L Thompson // Nothing to do, this won't occur 1619*edfbf3a6SJeremy L Thompson break; 1620c8c3fa7dSJeremy L Thompson } 1621c8c3fa7dSJeremy L Thompson CeedCall(CeedVectorRestoreArrayRead(basis->vec_chebyshev, &chebyshev_coeffs)); 1622c8c3fa7dSJeremy L Thompson CeedCall(CeedVectorRestoreArrayRead(x_ref, &x_array_read)); 1623c8c3fa7dSJeremy L Thompson CeedCall(CeedVectorRestoreArray(v, &v_array)); 1624c8c3fa7dSJeremy L Thompson break; 1625c8c3fa7dSJeremy L Thompson } 16262a94f45fSJeremy L Thompson case CEED_TRANSPOSE: { 16272a94f45fSJeremy L Thompson // Arbitrary points to nodes 16282a94f45fSJeremy L Thompson CeedScalar *chebyshev_coeffs; 16292a94f45fSJeremy L Thompson const CeedScalar *u_array, *x_array_read; 16302a94f45fSJeremy L Thompson 16312a94f45fSJeremy L Thompson // -- Transpose of evaluaton of Chebyshev polynomials at arbitrary points 16322a94f45fSJeremy L Thompson CeedCall(CeedVectorGetArrayWrite(basis->vec_chebyshev, CEED_MEM_HOST, &chebyshev_coeffs)); 16332a94f45fSJeremy L Thompson CeedCall(CeedVectorGetArrayRead(x_ref, CEED_MEM_HOST, &x_array_read)); 16342a94f45fSJeremy L Thompson CeedCall(CeedVectorGetArrayRead(u, CEED_MEM_HOST, &u_array)); 16352a94f45fSJeremy L Thompson { 16362a94f45fSJeremy L Thompson CeedScalar tmp[2][num_comp * CeedIntPow(Q_1d, dim)], chebyshev_x[Q_1d]; 16372a94f45fSJeremy L Thompson 16382a94f45fSJeremy L Thompson // ---- Values at point 16392a94f45fSJeremy L Thompson for (CeedInt p = 0; p < num_points; p++) { 16402a94f45fSJeremy L Thompson CeedInt pre = num_comp * 1, post = 1; 16412a94f45fSJeremy L Thompson 16422a94f45fSJeremy L Thompson for (CeedInt d = dim - 1; d >= 0; d--) { 16432a94f45fSJeremy L Thompson // ------ Compute Chebyshev polynomial values 16442a94f45fSJeremy L Thompson { 16452a94f45fSJeremy L Thompson const CeedScalar x = x_array_read[p * dim + d]; 16462a94f45fSJeremy L Thompson 16472a94f45fSJeremy L Thompson chebyshev_x[0] = 1.0; 16482a94f45fSJeremy L Thompson chebyshev_x[1] = 2 * x; 16492a94f45fSJeremy L Thompson for (CeedInt j = 2; j < Q_1d; j++) chebyshev_x[j] = 2 * x * chebyshev_x[j - 1] - chebyshev_x[j - 2]; 16502a94f45fSJeremy L Thompson } 16512a94f45fSJeremy L Thompson // ------ Tensor contract 16522a94f45fSJeremy L Thompson CeedCall(CeedTensorContractApply(basis->contract, pre, 1, post, Q_1d, chebyshev_x, t_mode, p > 0 && d == 0, 16532a94f45fSJeremy L Thompson d == (dim - 1) ? &u_array[p * num_comp] : tmp[d % 2], d == 0 ? chebyshev_coeffs : tmp[(d + 1) % 2])); 16542a94f45fSJeremy L Thompson pre /= 1; 16552a94f45fSJeremy L Thompson post *= Q_1d; 16562a94f45fSJeremy L Thompson } 16572a94f45fSJeremy L Thompson } 16582a94f45fSJeremy L Thompson } 16592a94f45fSJeremy L Thompson CeedCall(CeedVectorRestoreArray(basis->vec_chebyshev, &chebyshev_coeffs)); 16602a94f45fSJeremy L Thompson CeedCall(CeedVectorRestoreArrayRead(x_ref, &x_array_read)); 16612a94f45fSJeremy L Thompson CeedCall(CeedVectorRestoreArrayRead(u, &u_array)); 16622a94f45fSJeremy L Thompson 16632a94f45fSJeremy L Thompson // -- Interpolate transpose from Chebyshev coefficients 16642a94f45fSJeremy L Thompson CeedCall(CeedBasisApply(basis->basis_chebyshev, 1, CEED_TRANSPOSE, CEED_EVAL_INTERP, basis->vec_chebyshev, v)); 16652a94f45fSJeremy L Thompson break; 16662a94f45fSJeremy L Thompson } 1667c8c3fa7dSJeremy L Thompson } 1668c8c3fa7dSJeremy L Thompson 1669c8c3fa7dSJeremy L Thompson return CEED_ERROR_SUCCESS; 1670c8c3fa7dSJeremy L Thompson } 1671c8c3fa7dSJeremy L Thompson 1672c8c3fa7dSJeremy L Thompson /** 1673b7c9bbdaSJeremy L Thompson @brief Get Ceed associated with a CeedBasis 1674b7c9bbdaSJeremy L Thompson 1675ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 1676b7c9bbdaSJeremy L Thompson @param[out] ceed Variable to store Ceed 1677b7c9bbdaSJeremy L Thompson 1678b7c9bbdaSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 1679b7c9bbdaSJeremy L Thompson 1680b7c9bbdaSJeremy L Thompson @ref Advanced 1681b7c9bbdaSJeremy L Thompson **/ 1682b7c9bbdaSJeremy L Thompson int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) { 1683b7c9bbdaSJeremy L Thompson *ceed = basis->ceed; 1684b7c9bbdaSJeremy L Thompson return CEED_ERROR_SUCCESS; 1685b7c9bbdaSJeremy L Thompson } 1686b7c9bbdaSJeremy L Thompson 1687b7c9bbdaSJeremy L Thompson /** 16889d007619Sjeremylt @brief Get dimension for given CeedBasis 16899d007619Sjeremylt 1690ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 16919d007619Sjeremylt @param[out] dim Variable to store dimension of basis 16929d007619Sjeremylt 16939d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 16949d007619Sjeremylt 1695b7c9bbdaSJeremy L Thompson @ref Advanced 16969d007619Sjeremylt **/ 16979d007619Sjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) { 16989d007619Sjeremylt *dim = basis->dim; 1699e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 17009d007619Sjeremylt } 17019d007619Sjeremylt 17029d007619Sjeremylt /** 1703d99fa3c5SJeremy L Thompson @brief Get topology for given CeedBasis 1704d99fa3c5SJeremy L Thompson 1705ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 1706d99fa3c5SJeremy L Thompson @param[out] topo Variable to store topology of basis 1707d99fa3c5SJeremy L Thompson 1708d99fa3c5SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 1709d99fa3c5SJeremy L Thompson 1710b7c9bbdaSJeremy L Thompson @ref Advanced 1711d99fa3c5SJeremy L Thompson **/ 1712d99fa3c5SJeremy L Thompson int CeedBasisGetTopology(CeedBasis basis, CeedElemTopology *topo) { 1713d99fa3c5SJeremy L Thompson *topo = basis->topo; 1714e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1715d99fa3c5SJeremy L Thompson } 1716d99fa3c5SJeremy L Thompson 1717d99fa3c5SJeremy L Thompson /** 17189d007619Sjeremylt @brief Get number of components for given CeedBasis 17199d007619Sjeremylt 1720ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 1721d1d35e2fSjeremylt @param[out] num_comp Variable to store number of components of basis 17229d007619Sjeremylt 17239d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 17249d007619Sjeremylt 1725b7c9bbdaSJeremy L Thompson @ref Advanced 17269d007619Sjeremylt **/ 1727d1d35e2fSjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *num_comp) { 1728d1d35e2fSjeremylt *num_comp = basis->num_comp; 1729e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 17309d007619Sjeremylt } 17319d007619Sjeremylt 17329d007619Sjeremylt /** 17339d007619Sjeremylt @brief Get total number of nodes (in dim dimensions) of a CeedBasis 17349d007619Sjeremylt 1735ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 17369d007619Sjeremylt @param[out] P Variable to store number of nodes 17379d007619Sjeremylt 17389d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 17399d007619Sjeremylt 17409d007619Sjeremylt @ref Utility 17419d007619Sjeremylt **/ 17429d007619Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) { 17439d007619Sjeremylt *P = basis->P; 1744e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 17459d007619Sjeremylt } 17469d007619Sjeremylt 17479d007619Sjeremylt /** 17489d007619Sjeremylt @brief Get total number of nodes (in 1 dimension) of a CeedBasis 17499d007619Sjeremylt 1750ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 1751d1d35e2fSjeremylt @param[out] P_1d Variable to store number of nodes 17529d007619Sjeremylt 17539d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 17549d007619Sjeremylt 1755b7c9bbdaSJeremy L Thompson @ref Advanced 17569d007619Sjeremylt **/ 1757d1d35e2fSjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P_1d) { 17586402da51SJeremy L Thompson CeedCheck(basis->is_tensor_basis, basis->ceed, CEED_ERROR_MINOR, "Cannot supply P_1d for non-tensor basis"); 1759d1d35e2fSjeremylt *P_1d = basis->P_1d; 1760e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 17619d007619Sjeremylt } 17629d007619Sjeremylt 17639d007619Sjeremylt /** 17649d007619Sjeremylt @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis 17659d007619Sjeremylt 1766ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 17679d007619Sjeremylt @param[out] Q Variable to store number of quadrature points 17689d007619Sjeremylt 17699d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 17709d007619Sjeremylt 17719d007619Sjeremylt @ref Utility 17729d007619Sjeremylt **/ 17739d007619Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) { 17749d007619Sjeremylt *Q = basis->Q; 1775e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 17769d007619Sjeremylt } 17779d007619Sjeremylt 17789d007619Sjeremylt /** 17799d007619Sjeremylt @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis 17809d007619Sjeremylt 1781ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 1782d1d35e2fSjeremylt @param[out] Q_1d Variable to store number of quadrature points 17839d007619Sjeremylt 17849d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 17859d007619Sjeremylt 1786b7c9bbdaSJeremy L Thompson @ref Advanced 17879d007619Sjeremylt **/ 1788d1d35e2fSjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q_1d) { 17896402da51SJeremy L Thompson CeedCheck(basis->is_tensor_basis, basis->ceed, CEED_ERROR_MINOR, "Cannot supply Q_1d for non-tensor basis"); 1790d1d35e2fSjeremylt *Q_1d = basis->Q_1d; 1791e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 17929d007619Sjeremylt } 17939d007619Sjeremylt 17949d007619Sjeremylt /** 1795ea61e9acSJeremy L Thompson @brief Get reference coordinates of quadrature points (in dim dimensions) of a CeedBasis 17969d007619Sjeremylt 1797ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 1798d1d35e2fSjeremylt @param[out] q_ref Variable to store reference coordinates of quadrature points 17999d007619Sjeremylt 18009d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 18019d007619Sjeremylt 1802b7c9bbdaSJeremy L Thompson @ref Advanced 18039d007619Sjeremylt **/ 1804d1d35e2fSjeremylt int CeedBasisGetQRef(CeedBasis basis, const CeedScalar **q_ref) { 1805d1d35e2fSjeremylt *q_ref = basis->q_ref_1d; 1806e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 18079d007619Sjeremylt } 18089d007619Sjeremylt 18099d007619Sjeremylt /** 1810ea61e9acSJeremy L Thompson @brief Get quadrature weights of quadrature points (in dim dimensions) of a CeedBasis 18119d007619Sjeremylt 1812ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 1813d1d35e2fSjeremylt @param[out] q_weight Variable to store quadrature weights 18149d007619Sjeremylt 18159d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 18169d007619Sjeremylt 1817b7c9bbdaSJeremy L Thompson @ref Advanced 18189d007619Sjeremylt **/ 1819d1d35e2fSjeremylt int CeedBasisGetQWeights(CeedBasis basis, const CeedScalar **q_weight) { 1820d1d35e2fSjeremylt *q_weight = basis->q_weight_1d; 1821e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 18229d007619Sjeremylt } 18239d007619Sjeremylt 18249d007619Sjeremylt /** 18259d007619Sjeremylt @brief Get interpolation matrix of a CeedBasis 18269d007619Sjeremylt 1827ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 18289d007619Sjeremylt @param[out] interp Variable to store interpolation matrix 18299d007619Sjeremylt 18309d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 18319d007619Sjeremylt 1832b7c9bbdaSJeremy L Thompson @ref Advanced 18339d007619Sjeremylt **/ 18346c58de82SJeremy L Thompson int CeedBasisGetInterp(CeedBasis basis, const CeedScalar **interp) { 18356402da51SJeremy L Thompson if (!basis->interp && basis->is_tensor_basis) { 18369d007619Sjeremylt // Allocate 18372b730f8bSJeremy L Thompson CeedCall(CeedMalloc(basis->Q * basis->P, &basis->interp)); 18389d007619Sjeremylt 18399d007619Sjeremylt // Initialize 18402b730f8bSJeremy L Thompson for (CeedInt i = 0; i < basis->Q * basis->P; i++) basis->interp[i] = 1.0; 18419d007619Sjeremylt 18429d007619Sjeremylt // Calculate 18432b730f8bSJeremy L Thompson for (CeedInt d = 0; d < basis->dim; d++) { 18442b730f8bSJeremy L Thompson for (CeedInt qpt = 0; qpt < basis->Q; qpt++) { 18459d007619Sjeremylt for (CeedInt node = 0; node < basis->P; node++) { 1846d1d35e2fSjeremylt CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d; 1847d1d35e2fSjeremylt CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d; 1848d1d35e2fSjeremylt basis->interp[qpt * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p]; 18499d007619Sjeremylt } 18509d007619Sjeremylt } 18512b730f8bSJeremy L Thompson } 18522b730f8bSJeremy L Thompson } 18539d007619Sjeremylt *interp = basis->interp; 1854e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 18559d007619Sjeremylt } 18569d007619Sjeremylt 18579d007619Sjeremylt /** 18589d007619Sjeremylt @brief Get 1D interpolation matrix of a tensor product CeedBasis 18599d007619Sjeremylt 1860ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 1861d1d35e2fSjeremylt @param[out] interp_1d Variable to store interpolation matrix 18629d007619Sjeremylt 18639d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 18649d007619Sjeremylt 18659d007619Sjeremylt @ref Backend 18669d007619Sjeremylt **/ 1867d1d35e2fSjeremylt int CeedBasisGetInterp1D(CeedBasis basis, const CeedScalar **interp_1d) { 18686402da51SJeremy L Thompson CeedCheck(basis->is_tensor_basis, basis->ceed, CEED_ERROR_MINOR, "CeedBasis is not a tensor product basis."); 1869d1d35e2fSjeremylt *interp_1d = basis->interp_1d; 1870e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 18719d007619Sjeremylt } 18729d007619Sjeremylt 18739d007619Sjeremylt /** 18749d007619Sjeremylt @brief Get gradient matrix of a CeedBasis 18759d007619Sjeremylt 1876ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 18779d007619Sjeremylt @param[out] grad Variable to store gradient matrix 18789d007619Sjeremylt 18799d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 18809d007619Sjeremylt 1881b7c9bbdaSJeremy L Thompson @ref Advanced 18829d007619Sjeremylt **/ 18836c58de82SJeremy L Thompson int CeedBasisGetGrad(CeedBasis basis, const CeedScalar **grad) { 18846402da51SJeremy L Thompson if (!basis->grad && basis->is_tensor_basis) { 18859d007619Sjeremylt // Allocate 18862b730f8bSJeremy L Thompson CeedCall(CeedMalloc(basis->dim * basis->Q * basis->P, &basis->grad)); 18879d007619Sjeremylt 18889d007619Sjeremylt // Initialize 18892b730f8bSJeremy L Thompson for (CeedInt i = 0; i < basis->dim * basis->Q * basis->P; i++) basis->grad[i] = 1.0; 18909d007619Sjeremylt 18919d007619Sjeremylt // Calculate 18922b730f8bSJeremy L Thompson for (CeedInt d = 0; d < basis->dim; d++) { 18932b730f8bSJeremy L Thompson for (CeedInt i = 0; i < basis->dim; i++) { 18942b730f8bSJeremy L Thompson for (CeedInt qpt = 0; qpt < basis->Q; qpt++) { 18959d007619Sjeremylt for (CeedInt node = 0; node < basis->P; node++) { 1896d1d35e2fSjeremylt CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d; 1897d1d35e2fSjeremylt CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d; 18982b730f8bSJeremy L Thompson if (i == d) basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->grad_1d[q * basis->P_1d + p]; 18992b730f8bSJeremy L Thompson else basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p]; 19002b730f8bSJeremy L Thompson } 19012b730f8bSJeremy L Thompson } 19022b730f8bSJeremy L Thompson } 19039d007619Sjeremylt } 19049d007619Sjeremylt } 19059d007619Sjeremylt *grad = basis->grad; 1906e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 19079d007619Sjeremylt } 19089d007619Sjeremylt 19099d007619Sjeremylt /** 19109d007619Sjeremylt @brief Get 1D gradient matrix of a tensor product CeedBasis 19119d007619Sjeremylt 1912ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 1913d1d35e2fSjeremylt @param[out] grad_1d Variable to store gradient matrix 19149d007619Sjeremylt 19159d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 19169d007619Sjeremylt 1917b7c9bbdaSJeremy L Thompson @ref Advanced 19189d007619Sjeremylt **/ 1919d1d35e2fSjeremylt int CeedBasisGetGrad1D(CeedBasis basis, const CeedScalar **grad_1d) { 19206402da51SJeremy L Thompson CeedCheck(basis->is_tensor_basis, basis->ceed, CEED_ERROR_MINOR, "CeedBasis is not a tensor product basis."); 1921d1d35e2fSjeremylt *grad_1d = basis->grad_1d; 1922e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 19239d007619Sjeremylt } 19249d007619Sjeremylt 19259d007619Sjeremylt /** 192650c301a5SRezgar Shakeri @brief Get divergence matrix of a CeedBasis 192750c301a5SRezgar Shakeri 1928ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 192950c301a5SRezgar Shakeri @param[out] div Variable to store divergence matrix 193050c301a5SRezgar Shakeri 193150c301a5SRezgar Shakeri @return An error code: 0 - success, otherwise - failure 193250c301a5SRezgar Shakeri 193350c301a5SRezgar Shakeri @ref Advanced 193450c301a5SRezgar Shakeri **/ 193550c301a5SRezgar Shakeri int CeedBasisGetDiv(CeedBasis basis, const CeedScalar **div) { 19366574a04fSJeremy L Thompson CeedCheck(basis->div, basis->ceed, CEED_ERROR_MINOR, "CeedBasis does not have divergence matrix."); 193750c301a5SRezgar Shakeri *div = basis->div; 193850c301a5SRezgar Shakeri return CEED_ERROR_SUCCESS; 193950c301a5SRezgar Shakeri } 194050c301a5SRezgar Shakeri 194150c301a5SRezgar Shakeri /** 1942c4e3f59bSSebastian Grimberg @brief Get curl matrix of a CeedBasis 1943c4e3f59bSSebastian Grimberg 1944c4e3f59bSSebastian Grimberg @param[in] basis CeedBasis 1945c4e3f59bSSebastian Grimberg @param[out] curl Variable to store curl matrix 1946c4e3f59bSSebastian Grimberg 1947c4e3f59bSSebastian Grimberg @return An error code: 0 - success, otherwise - failure 1948c4e3f59bSSebastian Grimberg 1949c4e3f59bSSebastian Grimberg @ref Advanced 1950c4e3f59bSSebastian Grimberg **/ 1951c4e3f59bSSebastian Grimberg int CeedBasisGetCurl(CeedBasis basis, const CeedScalar **curl) { 19526574a04fSJeremy L Thompson CeedCheck(basis->curl, basis->ceed, CEED_ERROR_MINOR, "CeedBasis does not have curl matrix."); 1953c4e3f59bSSebastian Grimberg *curl = basis->curl; 1954c4e3f59bSSebastian Grimberg return CEED_ERROR_SUCCESS; 1955c4e3f59bSSebastian Grimberg } 1956c4e3f59bSSebastian Grimberg 1957c4e3f59bSSebastian Grimberg /** 19587a982d89SJeremy L. Thompson @brief Destroy a CeedBasis 19597a982d89SJeremy L. Thompson 1960ea61e9acSJeremy L Thompson @param[in,out] basis CeedBasis to destroy 19617a982d89SJeremy L. Thompson 19627a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 19637a982d89SJeremy L. Thompson 19647a982d89SJeremy L. Thompson @ref User 19657a982d89SJeremy L. Thompson **/ 19667a982d89SJeremy L. Thompson int CeedBasisDestroy(CeedBasis *basis) { 19677425e127SJeremy L Thompson if (!*basis || *basis == CEED_BASIS_COLLOCATED || --(*basis)->ref_count > 0) { 1968ad6481ceSJeremy L Thompson *basis = NULL; 1969ad6481ceSJeremy L Thompson return CEED_ERROR_SUCCESS; 1970ad6481ceSJeremy L Thompson } 19712b730f8bSJeremy L Thompson if ((*basis)->Destroy) CeedCall((*basis)->Destroy(*basis)); 19729831d45aSJeremy L Thompson CeedCall(CeedTensorContractDestroy(&(*basis)->contract)); 1973c4e3f59bSSebastian Grimberg CeedCall(CeedFree(&(*basis)->q_ref_1d)); 1974c4e3f59bSSebastian Grimberg CeedCall(CeedFree(&(*basis)->q_weight_1d)); 19752b730f8bSJeremy L Thompson CeedCall(CeedFree(&(*basis)->interp)); 19762b730f8bSJeremy L Thompson CeedCall(CeedFree(&(*basis)->interp_1d)); 19772b730f8bSJeremy L Thompson CeedCall(CeedFree(&(*basis)->grad)); 19782b730f8bSJeremy L Thompson CeedCall(CeedFree(&(*basis)->grad_1d)); 1979c4e3f59bSSebastian Grimberg CeedCall(CeedFree(&(*basis)->div)); 1980c4e3f59bSSebastian Grimberg CeedCall(CeedFree(&(*basis)->curl)); 1981c8c3fa7dSJeremy L Thompson CeedCall(CeedVectorDestroy(&(*basis)->vec_chebyshev)); 1982c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisDestroy(&(*basis)->basis_chebyshev)); 19832b730f8bSJeremy L Thompson CeedCall(CeedDestroy(&(*basis)->ceed)); 19842b730f8bSJeremy L Thompson CeedCall(CeedFree(basis)); 1985e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 19867a982d89SJeremy L. Thompson } 19877a982d89SJeremy L. Thompson 19887a982d89SJeremy L. Thompson /** 1989b11c1e72Sjeremylt @brief Construct a Gauss-Legendre quadrature 1990b11c1e72Sjeremylt 1991ea61e9acSJeremy L Thompson @param[in] Q Number of quadrature points (integrates polynomials of degree 2*Q-1 exactly) 1992d1d35e2fSjeremylt @param[out] q_ref_1d Array of length Q to hold the abscissa on [-1, 1] 1993d1d35e2fSjeremylt @param[out] q_weight_1d Array of length Q to hold the weights 1994b11c1e72Sjeremylt 1995b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 1996dfdf5a53Sjeremylt 1997dfdf5a53Sjeremylt @ref Utility 1998b11c1e72Sjeremylt **/ 19992b730f8bSJeremy L Thompson int CeedGaussQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) { 2000d7b241e6Sjeremylt // Allocate 2001d7b241e6Sjeremylt CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0 * atan(1.0); 2002d1d35e2fSjeremylt // Build q_ref_1d, q_weight_1d 200392ae7e47SJeremy L Thompson for (CeedInt i = 0; i <= Q / 2; i++) { 2004d7b241e6Sjeremylt // Guess 2005d7b241e6Sjeremylt xi = cos(PI * (CeedScalar)(2 * i + 1) / ((CeedScalar)(2 * Q))); 2006d7b241e6Sjeremylt // Pn(xi) 2007d7b241e6Sjeremylt P0 = 1.0; 2008d7b241e6Sjeremylt P1 = xi; 2009d7b241e6Sjeremylt P2 = 0.0; 201092ae7e47SJeremy L Thompson for (CeedInt j = 2; j <= Q; j++) { 2011d7b241e6Sjeremylt P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j)); 2012d7b241e6Sjeremylt P0 = P1; 2013d7b241e6Sjeremylt P1 = P2; 2014d7b241e6Sjeremylt } 2015d7b241e6Sjeremylt // First Newton Step 2016d7b241e6Sjeremylt dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0); 2017d7b241e6Sjeremylt xi = xi - P2 / dP2; 2018d7b241e6Sjeremylt // Newton to convergence 201992ae7e47SJeremy L Thompson for (CeedInt k = 0; k < 100 && fabs(P2) > 10 * CEED_EPSILON; k++) { 2020d7b241e6Sjeremylt P0 = 1.0; 2021d7b241e6Sjeremylt P1 = xi; 202292ae7e47SJeremy L Thompson for (CeedInt j = 2; j <= Q; j++) { 2023d7b241e6Sjeremylt P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j)); 2024d7b241e6Sjeremylt P0 = P1; 2025d7b241e6Sjeremylt P1 = P2; 2026d7b241e6Sjeremylt } 2027d7b241e6Sjeremylt dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0); 2028d7b241e6Sjeremylt xi = xi - P2 / dP2; 2029d7b241e6Sjeremylt } 2030d7b241e6Sjeremylt // Save xi, wi 2031d7b241e6Sjeremylt wi = 2.0 / ((1.0 - xi * xi) * dP2 * dP2); 2032d1d35e2fSjeremylt q_weight_1d[i] = wi; 2033d1d35e2fSjeremylt q_weight_1d[Q - 1 - i] = wi; 2034d1d35e2fSjeremylt q_ref_1d[i] = -xi; 2035d1d35e2fSjeremylt q_ref_1d[Q - 1 - i] = xi; 2036d7b241e6Sjeremylt } 2037e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 2038d7b241e6Sjeremylt } 2039d7b241e6Sjeremylt 2040b11c1e72Sjeremylt /** 2041b11c1e72Sjeremylt @brief Construct a Gauss-Legendre-Lobatto quadrature 2042b11c1e72Sjeremylt 2043ea61e9acSJeremy L Thompson @param[in] Q Number of quadrature points (integrates polynomials of degree 2*Q-3 exactly) 2044d1d35e2fSjeremylt @param[out] q_ref_1d Array of length Q to hold the abscissa on [-1, 1] 2045d1d35e2fSjeremylt @param[out] q_weight_1d Array of length Q to hold the weights 2046b11c1e72Sjeremylt 2047b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 2048dfdf5a53Sjeremylt 2049dfdf5a53Sjeremylt @ref Utility 2050b11c1e72Sjeremylt **/ 20512b730f8bSJeremy L Thompson int CeedLobattoQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) { 2052d7b241e6Sjeremylt // Allocate 2053d7b241e6Sjeremylt CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0 * atan(1.0); 2054d1d35e2fSjeremylt // Build q_ref_1d, q_weight_1d 2055d7b241e6Sjeremylt // Set endpoints 20566574a04fSJeremy L Thompson CeedCheck(Q > 1, NULL, CEED_ERROR_DIMENSION, "Cannot create Lobatto quadrature with Q=%" CeedInt_FMT " < 2 points", Q); 2057d7b241e6Sjeremylt wi = 2.0 / ((CeedScalar)(Q * (Q - 1))); 2058d1d35e2fSjeremylt if (q_weight_1d) { 2059d1d35e2fSjeremylt q_weight_1d[0] = wi; 2060d1d35e2fSjeremylt q_weight_1d[Q - 1] = wi; 2061d7b241e6Sjeremylt } 2062d1d35e2fSjeremylt q_ref_1d[0] = -1.0; 2063d1d35e2fSjeremylt q_ref_1d[Q - 1] = 1.0; 2064d7b241e6Sjeremylt // Interior 206592ae7e47SJeremy L Thompson for (CeedInt i = 1; i <= (Q - 1) / 2; i++) { 2066d7b241e6Sjeremylt // Guess 2067d7b241e6Sjeremylt xi = cos(PI * (CeedScalar)(i) / (CeedScalar)(Q - 1)); 2068d7b241e6Sjeremylt // Pn(xi) 2069d7b241e6Sjeremylt P0 = 1.0; 2070d7b241e6Sjeremylt P1 = xi; 2071d7b241e6Sjeremylt P2 = 0.0; 207292ae7e47SJeremy L Thompson for (CeedInt j = 2; j < Q; j++) { 2073d7b241e6Sjeremylt P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j)); 2074d7b241e6Sjeremylt P0 = P1; 2075d7b241e6Sjeremylt P1 = P2; 2076d7b241e6Sjeremylt } 2077d7b241e6Sjeremylt // First Newton step 2078d7b241e6Sjeremylt dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0); 2079d7b241e6Sjeremylt d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi); 2080d7b241e6Sjeremylt xi = xi - dP2 / d2P2; 2081d7b241e6Sjeremylt // Newton to convergence 208292ae7e47SJeremy L Thompson for (CeedInt k = 0; k < 100 && fabs(dP2) > 10 * CEED_EPSILON; k++) { 2083d7b241e6Sjeremylt P0 = 1.0; 2084d7b241e6Sjeremylt P1 = xi; 208592ae7e47SJeremy L Thompson for (CeedInt j = 2; j < Q; j++) { 2086d7b241e6Sjeremylt P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j)); 2087d7b241e6Sjeremylt P0 = P1; 2088d7b241e6Sjeremylt P1 = P2; 2089d7b241e6Sjeremylt } 2090d7b241e6Sjeremylt dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0); 2091d7b241e6Sjeremylt d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi); 2092d7b241e6Sjeremylt xi = xi - dP2 / d2P2; 2093d7b241e6Sjeremylt } 2094d7b241e6Sjeremylt // Save xi, wi 2095d7b241e6Sjeremylt wi = 2.0 / (((CeedScalar)(Q * (Q - 1))) * P2 * P2); 2096d1d35e2fSjeremylt if (q_weight_1d) { 2097d1d35e2fSjeremylt q_weight_1d[i] = wi; 2098d1d35e2fSjeremylt q_weight_1d[Q - 1 - i] = wi; 2099d7b241e6Sjeremylt } 2100d1d35e2fSjeremylt q_ref_1d[i] = -xi; 2101d1d35e2fSjeremylt q_ref_1d[Q - 1 - i] = xi; 2102d7b241e6Sjeremylt } 2103e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 2104d7b241e6Sjeremylt } 2105d7b241e6Sjeremylt 2106d7b241e6Sjeremylt /// @} 2107