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) { 115*edf04919SJeremy L Thompson if (m > 1) { 116*edf04919SJeremy L Thompson fprintf(stream, " %s:\n", name); 117*edf04919SJeremy L Thompson } else { 118*edf04919SJeremy L Thompson char padded_name[12]; 119*edf04919SJeremy L Thompson 120*edf04919SJeremy L Thompson snprintf(padded_name, 11, "%s:", name); 121*edf04919SJeremy L Thompson fprintf(stream, " %-10s", padded_name); 122*edf04919SJeremy L Thompson } 12392ae7e47SJeremy L Thompson for (CeedInt i = 0; i < m; i++) { 124*edf04919SJeremy L Thompson if (m > 1) fprintf(stream, " [%" CeedInt_FMT "]", i); 1252b730f8bSJeremy 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); 1267a982d89SJeremy L. Thompson fputs("\n", stream); 1277a982d89SJeremy L. Thompson } 128e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1297a982d89SJeremy L. Thompson } 1307a982d89SJeremy L. Thompson 131a76a04e7SJeremy L Thompson /** 132ea61e9acSJeremy L Thompson @brief Create the interpolation and gradient matrices for projection from the nodes of `basis_from` to the nodes of `basis_to`. 133ba59ac12SSebastian Grimberg 13415ad3917SSebastian Grimberg The interpolation is given by `interp_project = interp_to^+ * interp_from`, where the pseudoinverse `interp_to^+` is given by QR factorization. 13515ad3917SSebastian 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. 13615ad3917SSebastian Grimberg 137ba59ac12SSebastian Grimberg Note: `basis_from` and `basis_to` must have compatible quadrature spaces. 138a76a04e7SJeremy L Thompson 139a76a04e7SJeremy L Thompson @param[in] basis_from CeedBasis to project from 140a76a04e7SJeremy L Thompson @param[in] basis_to CeedBasis to project to 141ea61e9acSJeremy L Thompson @param[out] interp_project Address of the variable where the newly created interpolation matrix will be stored. 142ea61e9acSJeremy L Thompson @param[out] grad_project Address of the variable where the newly created gradient matrix will be stored. 143a76a04e7SJeremy L Thompson 144a76a04e7SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 145a76a04e7SJeremy L Thompson 146a76a04e7SJeremy L Thompson @ref Developer 147a76a04e7SJeremy L Thompson **/ 1482b730f8bSJeremy L Thompson static int CeedBasisCreateProjectionMatrices(CeedBasis basis_from, CeedBasis basis_to, CeedScalar **interp_project, CeedScalar **grad_project) { 149a76a04e7SJeremy L Thompson Ceed ceed; 1502b730f8bSJeremy L Thompson CeedCall(CeedBasisGetCeed(basis_to, &ceed)); 151a76a04e7SJeremy L Thompson 152a76a04e7SJeremy L Thompson // Check for compatible quadrature spaces 153a76a04e7SJeremy L Thompson CeedInt Q_to, Q_from; 1542b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints(basis_to, &Q_to)); 1552b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints(basis_from, &Q_from)); 1566574a04fSJeremy L Thompson CeedCheck(Q_to == Q_from, ceed, CEED_ERROR_DIMENSION, "Bases must have compatible quadrature spaces"); 157a76a04e7SJeremy L Thompson 15814556e63SJeremy L Thompson // Check for matching tensor or non-tensor 159a76a04e7SJeremy L Thompson CeedInt P_to, P_from, Q = Q_to; 160a76a04e7SJeremy L Thompson bool is_tensor_to, is_tensor_from; 1612b730f8bSJeremy L Thompson CeedCall(CeedBasisIsTensor(basis_to, &is_tensor_to)); 1622b730f8bSJeremy L Thompson CeedCall(CeedBasisIsTensor(basis_from, &is_tensor_from)); 1636574a04fSJeremy L Thompson CeedCheck(is_tensor_to == is_tensor_from, ceed, CEED_ERROR_MINOR, "Bases must both be tensor or non-tensor"); 1646574a04fSJeremy L Thompson if (is_tensor_to) { 1652b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_to)); 1662b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_from)); 1672b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints1D(basis_from, &Q)); 1686574a04fSJeremy L Thompson } else { 1692b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis_to, &P_to)); 1702b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis_from, &P_from)); 171a76a04e7SJeremy L Thompson } 172a76a04e7SJeremy L Thompson 17315ad3917SSebastian Grimberg // Check for matching FE space 17415ad3917SSebastian Grimberg CeedFESpace fe_space_to, fe_space_from; 17515ad3917SSebastian Grimberg CeedCall(CeedBasisGetFESpace(basis_to, &fe_space_to)); 17615ad3917SSebastian Grimberg CeedCall(CeedBasisGetFESpace(basis_from, &fe_space_from)); 1776574a04fSJeremy L Thompson CeedCheck(fe_space_to == fe_space_from, ceed, CEED_ERROR_MINOR, "Bases must both be the same FE space type"); 17815ad3917SSebastian Grimberg 17914556e63SJeremy L Thompson // Get source matrices 18015ad3917SSebastian Grimberg CeedInt dim, q_comp = 1; 18115ad3917SSebastian Grimberg const CeedScalar *interp_to_source = NULL, *interp_from_source = NULL; 18214556e63SJeremy L Thompson CeedScalar *interp_to, *interp_from, *tau; 1832b730f8bSJeremy L Thompson CeedCall(CeedBasisGetDimension(basis_to, &dim)); 184a76a04e7SJeremy L Thompson if (is_tensor_to) { 1852b730f8bSJeremy L Thompson CeedCall(CeedBasisGetInterp1D(basis_to, &interp_to_source)); 1862b730f8bSJeremy L Thompson CeedCall(CeedBasisGetInterp1D(basis_from, &interp_from_source)); 187a76a04e7SJeremy L Thompson } else { 18815ad3917SSebastian Grimberg CeedCall(CeedBasisGetNumQuadratureComponents(basis_from, CEED_EVAL_INTERP, &q_comp)); 1892b730f8bSJeremy L Thompson CeedCall(CeedBasisGetInterp(basis_to, &interp_to_source)); 1902b730f8bSJeremy L Thompson CeedCall(CeedBasisGetInterp(basis_from, &interp_from_source)); 19115ad3917SSebastian Grimberg } 19215ad3917SSebastian Grimberg CeedCall(CeedMalloc(Q * P_from * q_comp, &interp_from)); 19315ad3917SSebastian Grimberg CeedCall(CeedMalloc(Q * P_to * q_comp, &interp_to)); 19415ad3917SSebastian Grimberg CeedCall(CeedCalloc(P_to * P_from, interp_project)); 19515ad3917SSebastian Grimberg CeedCall(CeedMalloc(Q * q_comp, &tau)); 19615ad3917SSebastian Grimberg 19715ad3917SSebastian Grimberg // `grad_project = interp_to^+ * grad_from` is computed for the H^1 space case so the 198de05fbb2SSebastian Grimberg // projection basis will have a gradient operation (allocated even if not H^1 for the 199de05fbb2SSebastian Grimberg // basis construction later on) 20015ad3917SSebastian Grimberg const CeedScalar *grad_from_source = NULL; 20115ad3917SSebastian Grimberg if (fe_space_to == CEED_FE_SPACE_H1) { 20215ad3917SSebastian Grimberg if (is_tensor_to) { 20315ad3917SSebastian Grimberg CeedCall(CeedBasisGetGrad1D(basis_from, &grad_from_source)); 20415ad3917SSebastian Grimberg } else { 2052b730f8bSJeremy L Thompson CeedCall(CeedBasisGetGrad(basis_from, &grad_from_source)); 206a76a04e7SJeremy L Thompson } 207de05fbb2SSebastian Grimberg } 20815ad3917SSebastian Grimberg CeedCall(CeedCalloc(P_to * P_from * (is_tensor_to ? 1 : dim), grad_project)); 20915ad3917SSebastian Grimberg 21015ad3917SSebastian Grimberg // QR Factorization, interp_to = Q R 21115ad3917SSebastian Grimberg memcpy(interp_to, interp_to_source, Q * P_to * q_comp * sizeof(interp_to_source[0])); 21215ad3917SSebastian Grimberg CeedCall(CeedQRFactorization(ceed, interp_to, tau, Q * q_comp, P_to)); 213a76a04e7SJeremy L Thompson 21414556e63SJeremy L Thompson // Build matrices 21515ad3917SSebastian Grimberg CeedInt num_matrices = 1 + (fe_space_to == CEED_FE_SPACE_H1) * (is_tensor_to ? 1 : dim); 21614556e63SJeremy L Thompson CeedScalar *input_from[num_matrices], *output_project[num_matrices]; 21714556e63SJeremy L Thompson input_from[0] = (CeedScalar *)interp_from_source; 21814556e63SJeremy L Thompson output_project[0] = *interp_project; 21914556e63SJeremy L Thompson for (CeedInt m = 1; m < num_matrices; m++) { 22014556e63SJeremy L Thompson input_from[m] = (CeedScalar *)&grad_from_source[(m - 1) * Q * P_from]; 22102af4036SJeremy L Thompson output_project[m] = &((*grad_project)[(m - 1) * P_to * P_from]); 22214556e63SJeremy L Thompson } 22314556e63SJeremy L Thompson for (CeedInt m = 0; m < num_matrices; m++) { 22415ad3917SSebastian Grimberg // Apply Q^T, interp_from = Q^T interp_from 22515ad3917SSebastian Grimberg memcpy(interp_from, input_from[m], Q * P_from * q_comp * sizeof(input_from[m][0])); 22615ad3917SSebastian Grimberg CeedCall(CeedHouseholderApplyQ(interp_from, interp_to, tau, CEED_TRANSPOSE, Q * q_comp, P_from, P_to, P_from, 1)); 227a76a04e7SJeremy L Thompson 22815ad3917SSebastian Grimberg // Apply Rinv, output_project = Rinv interp_from 229a76a04e7SJeremy L Thompson for (CeedInt j = 0; j < P_from; j++) { // Column j 2302b730f8bSJeremy 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]; 231a76a04e7SJeremy L Thompson for (CeedInt i = P_to - 2; i >= 0; i--) { // Row i 23214556e63SJeremy L Thompson output_project[m][j + P_from * i] = interp_from[j + P_from * i]; 233a76a04e7SJeremy L Thompson for (CeedInt k = i + 1; k < P_to; k++) { 2342b730f8bSJeremy L Thompson output_project[m][j + P_from * i] -= interp_to[k + P_to * i] * output_project[m][j + P_from * k]; 235a76a04e7SJeremy L Thompson } 23614556e63SJeremy L Thompson output_project[m][j + P_from * i] /= interp_to[i + P_to * i]; 237a76a04e7SJeremy L Thompson } 238a76a04e7SJeremy L Thompson } 23914556e63SJeremy L Thompson } 24014556e63SJeremy L Thompson 24114556e63SJeremy L Thompson // Cleanup 2422b730f8bSJeremy L Thompson CeedCall(CeedFree(&tau)); 2432b730f8bSJeremy L Thompson CeedCall(CeedFree(&interp_to)); 2442b730f8bSJeremy L Thompson CeedCall(CeedFree(&interp_from)); 245a76a04e7SJeremy L Thompson 246a76a04e7SJeremy L Thompson return CEED_ERROR_SUCCESS; 247a76a04e7SJeremy L Thompson } 248a76a04e7SJeremy L Thompson 2497a982d89SJeremy L. Thompson /// @} 2507a982d89SJeremy L. Thompson 2517a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 2527a982d89SJeremy L. Thompson /// Ceed Backend API 2537a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 2547a982d89SJeremy L. Thompson /// @addtogroup CeedBasisBackend 2557a982d89SJeremy L. Thompson /// @{ 2567a982d89SJeremy L. Thompson 2577a982d89SJeremy L. Thompson /** 2587a982d89SJeremy L. Thompson @brief Return collocated grad matrix 2597a982d89SJeremy L. Thompson 260ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 261ea61e9acSJeremy L Thompson @param[out] collo_grad_1d Row-major (Q_1d * Q_1d) matrix expressing derivatives of basis functions at quadrature points 2627a982d89SJeremy L. Thompson 2637a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 2647a982d89SJeremy L. Thompson 2657a982d89SJeremy L. Thompson @ref Backend 2667a982d89SJeremy L. Thompson **/ 267d1d35e2fSjeremylt int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *collo_grad_1d) { 2687a982d89SJeremy L. Thompson Ceed ceed; 2692b730f8bSJeremy L Thompson CeedInt P_1d = (basis)->P_1d, Q_1d = (basis)->Q_1d; 27078464608Sjeremylt CeedScalar *interp_1d, *grad_1d, *tau; 2717a982d89SJeremy L. Thompson 2722b730f8bSJeremy L Thompson CeedCall(CeedMalloc(Q_1d * P_1d, &interp_1d)); 2732b730f8bSJeremy L Thompson CeedCall(CeedMalloc(Q_1d * P_1d, &grad_1d)); 2742b730f8bSJeremy L Thompson CeedCall(CeedMalloc(Q_1d, &tau)); 275d1d35e2fSjeremylt memcpy(interp_1d, (basis)->interp_1d, Q_1d * P_1d * sizeof(basis)->interp_1d[0]); 276d1d35e2fSjeremylt memcpy(grad_1d, (basis)->grad_1d, Q_1d * P_1d * sizeof(basis)->interp_1d[0]); 2777a982d89SJeremy L. Thompson 278d1d35e2fSjeremylt // QR Factorization, interp_1d = Q R 2792b730f8bSJeremy L Thompson CeedCall(CeedBasisGetCeed(basis, &ceed)); 2802b730f8bSJeremy L Thompson CeedCall(CeedQRFactorization(ceed, interp_1d, tau, Q_1d, P_1d)); 281ea61e9acSJeremy 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. 2827a982d89SJeremy L. Thompson 283c8c3fa7dSJeremy L Thompson // Apply R_inv, collo_grad_1d = grad_1d R_inv 284c8c3fa7dSJeremy L Thompson for (CeedInt i = 0; i < Q_1d; i++) { // Row i 285d1d35e2fSjeremylt collo_grad_1d[Q_1d * i] = grad_1d[P_1d * i] / interp_1d[0]; 286c8c3fa7dSJeremy L Thompson for (CeedInt j = 1; j < P_1d; j++) { // Column j 287d1d35e2fSjeremylt collo_grad_1d[j + Q_1d * i] = grad_1d[j + P_1d * i]; 288c8c3fa7dSJeremy 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]; 289d1d35e2fSjeremylt collo_grad_1d[j + Q_1d * i] /= interp_1d[j + P_1d * j]; 2907a982d89SJeremy L. Thompson } 291c8c3fa7dSJeremy L Thompson for (CeedInt j = P_1d; j < Q_1d; j++) collo_grad_1d[j + Q_1d * i] = 0; 2927a982d89SJeremy L. Thompson } 2937a982d89SJeremy L. Thompson 29415ad3917SSebastian Grimberg // Apply Q^T, collo_grad_1d = collo_grad_1d Q^T 2952b730f8bSJeremy L Thompson CeedCall(CeedHouseholderApplyQ(collo_grad_1d, interp_1d, tau, CEED_NOTRANSPOSE, Q_1d, Q_1d, P_1d, 1, Q_1d)); 2967a982d89SJeremy L. Thompson 2972b730f8bSJeremy L Thompson CeedCall(CeedFree(&interp_1d)); 2982b730f8bSJeremy L Thompson CeedCall(CeedFree(&grad_1d)); 2992b730f8bSJeremy L Thompson CeedCall(CeedFree(&tau)); 300e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3017a982d89SJeremy L. Thompson } 3027a982d89SJeremy L. Thompson 3037a982d89SJeremy L. Thompson /** 3047a982d89SJeremy L. Thompson @brief Get tensor status for given CeedBasis 3057a982d89SJeremy L. Thompson 306ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 307d1d35e2fSjeremylt @param[out] is_tensor Variable to store tensor status 3087a982d89SJeremy L. Thompson 3097a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3107a982d89SJeremy L. Thompson 3117a982d89SJeremy L. Thompson @ref Backend 3127a982d89SJeremy L. Thompson **/ 313d1d35e2fSjeremylt int CeedBasisIsTensor(CeedBasis basis, bool *is_tensor) { 3146402da51SJeremy L Thompson *is_tensor = basis->is_tensor_basis; 315e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3167a982d89SJeremy L. Thompson } 3177a982d89SJeremy L. Thompson 3187a982d89SJeremy L. Thompson /** 3197a982d89SJeremy L. Thompson @brief Get backend data of a CeedBasis 3207a982d89SJeremy L. Thompson 321ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 3227a982d89SJeremy L. Thompson @param[out] data Variable to store data 3237a982d89SJeremy L. Thompson 3247a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3257a982d89SJeremy L. Thompson 3267a982d89SJeremy L. Thompson @ref Backend 3277a982d89SJeremy L. Thompson **/ 328777ff853SJeremy L Thompson int CeedBasisGetData(CeedBasis basis, void *data) { 329777ff853SJeremy L Thompson *(void **)data = basis->data; 330e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3317a982d89SJeremy L. Thompson } 3327a982d89SJeremy L. Thompson 3337a982d89SJeremy L. Thompson /** 3347a982d89SJeremy L. Thompson @brief Set backend data of a CeedBasis 3357a982d89SJeremy L. Thompson 336ea61e9acSJeremy L Thompson @param[in,out] basis CeedBasis 337ea61e9acSJeremy L Thompson @param[in] data Data to set 3387a982d89SJeremy L. Thompson 3397a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3407a982d89SJeremy L. Thompson 3417a982d89SJeremy L. Thompson @ref Backend 3427a982d89SJeremy L. Thompson **/ 343777ff853SJeremy L Thompson int CeedBasisSetData(CeedBasis basis, void *data) { 344777ff853SJeremy L Thompson basis->data = data; 345e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3467a982d89SJeremy L. Thompson } 3477a982d89SJeremy L. Thompson 3487a982d89SJeremy L. Thompson /** 34934359f16Sjeremylt @brief Increment the reference counter for a CeedBasis 35034359f16Sjeremylt 351ea61e9acSJeremy L Thompson @param[in,out] basis Basis to increment the reference counter 35234359f16Sjeremylt 35334359f16Sjeremylt @return An error code: 0 - success, otherwise - failure 35434359f16Sjeremylt 35534359f16Sjeremylt @ref Backend 35634359f16Sjeremylt **/ 3579560d06aSjeremylt int CeedBasisReference(CeedBasis basis) { 35834359f16Sjeremylt basis->ref_count++; 35934359f16Sjeremylt return CEED_ERROR_SUCCESS; 36034359f16Sjeremylt } 36134359f16Sjeremylt 36234359f16Sjeremylt /** 363c4e3f59bSSebastian Grimberg @brief Get number of Q-vector components for given CeedBasis 364c4e3f59bSSebastian Grimberg 365c4e3f59bSSebastian Grimberg @param[in] basis CeedBasis 366c4e3f59bSSebastian Grimberg @param[in] eval_mode \ref CEED_EVAL_INTERP to use interpolated values, 367c4e3f59bSSebastian Grimberg \ref CEED_EVAL_GRAD to use gradients, 368c4e3f59bSSebastian Grimberg \ref CEED_EVAL_DIV to use divergence, 369c4e3f59bSSebastian Grimberg \ref CEED_EVAL_CURL to use curl. 370c4e3f59bSSebastian Grimberg @param[out] q_comp Variable to store number of Q-vector components of basis 371c4e3f59bSSebastian Grimberg 372c4e3f59bSSebastian Grimberg @return An error code: 0 - success, otherwise - failure 373c4e3f59bSSebastian Grimberg 374c4e3f59bSSebastian Grimberg @ref Backend 375c4e3f59bSSebastian Grimberg **/ 376c4e3f59bSSebastian Grimberg int CeedBasisGetNumQuadratureComponents(CeedBasis basis, CeedEvalMode eval_mode, CeedInt *q_comp) { 377c4e3f59bSSebastian Grimberg switch (eval_mode) { 378c4e3f59bSSebastian Grimberg case CEED_EVAL_INTERP: 379c4e3f59bSSebastian Grimberg *q_comp = (basis->fe_space == CEED_FE_SPACE_H1) ? 1 : basis->dim; 380c4e3f59bSSebastian Grimberg break; 381c4e3f59bSSebastian Grimberg case CEED_EVAL_GRAD: 382c4e3f59bSSebastian Grimberg *q_comp = basis->dim; 383c4e3f59bSSebastian Grimberg break; 384c4e3f59bSSebastian Grimberg case CEED_EVAL_DIV: 385c4e3f59bSSebastian Grimberg *q_comp = 1; 386c4e3f59bSSebastian Grimberg break; 387c4e3f59bSSebastian Grimberg case CEED_EVAL_CURL: 388c4e3f59bSSebastian Grimberg *q_comp = (basis->dim < 3) ? 1 : basis->dim; 389c4e3f59bSSebastian Grimberg break; 390c4e3f59bSSebastian Grimberg case CEED_EVAL_NONE: 391c4e3f59bSSebastian Grimberg case CEED_EVAL_WEIGHT: 392352a5e7cSSebastian Grimberg *q_comp = 1; 393c4e3f59bSSebastian Grimberg break; 394c4e3f59bSSebastian Grimberg } 395c4e3f59bSSebastian Grimberg return CEED_ERROR_SUCCESS; 396c4e3f59bSSebastian Grimberg } 397c4e3f59bSSebastian Grimberg 398c4e3f59bSSebastian Grimberg /** 3996e15d496SJeremy L Thompson @brief Estimate number of FLOPs required to apply CeedBasis in t_mode and eval_mode 4006e15d496SJeremy L Thompson 401ea61e9acSJeremy L Thompson @param[in] basis Basis to estimate FLOPs for 402ea61e9acSJeremy L Thompson @param[in] t_mode Apply basis or transpose 403ea61e9acSJeremy L Thompson @param[in] eval_mode Basis evaluation mode 404ea61e9acSJeremy L Thompson @param[out] flops Address of variable to hold FLOPs estimate 4056e15d496SJeremy L Thompson 4066e15d496SJeremy L Thompson @ref Backend 4076e15d496SJeremy L Thompson **/ 4082b730f8bSJeremy L Thompson int CeedBasisGetFlopsEstimate(CeedBasis basis, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedSize *flops) { 4096e15d496SJeremy L Thompson bool is_tensor; 4106e15d496SJeremy L Thompson 4112b730f8bSJeremy L Thompson CeedCall(CeedBasisIsTensor(basis, &is_tensor)); 4126e15d496SJeremy L Thompson if (is_tensor) { 4136e15d496SJeremy L Thompson CeedInt dim, num_comp, P_1d, Q_1d; 4142b730f8bSJeremy L Thompson CeedCall(CeedBasisGetDimension(basis, &dim)); 4152b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumComponents(basis, &num_comp)); 4162b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d)); 4172b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 4186e15d496SJeremy L Thompson if (t_mode == CEED_TRANSPOSE) { 4192b730f8bSJeremy L Thompson P_1d = Q_1d; 4202b730f8bSJeremy L Thompson Q_1d = P_1d; 4216e15d496SJeremy L Thompson } 4226e15d496SJeremy L Thompson CeedInt tensor_flops = 0, pre = num_comp * CeedIntPow(P_1d, dim - 1), post = 1; 4236e15d496SJeremy L Thompson for (CeedInt d = 0; d < dim; d++) { 4246e15d496SJeremy L Thompson tensor_flops += 2 * pre * P_1d * post * Q_1d; 4256e15d496SJeremy L Thompson pre /= P_1d; 4266e15d496SJeremy L Thompson post *= Q_1d; 4276e15d496SJeremy L Thompson } 4286e15d496SJeremy L Thompson switch (eval_mode) { 4292b730f8bSJeremy L Thompson case CEED_EVAL_NONE: 4302b730f8bSJeremy L Thompson *flops = 0; 4312b730f8bSJeremy L Thompson break; 4322b730f8bSJeremy L Thompson case CEED_EVAL_INTERP: 4332b730f8bSJeremy L Thompson *flops = tensor_flops; 4342b730f8bSJeremy L Thompson break; 4352b730f8bSJeremy L Thompson case CEED_EVAL_GRAD: 4362b730f8bSJeremy L Thompson *flops = tensor_flops * 2; 4372b730f8bSJeremy L Thompson break; 4386e15d496SJeremy L Thompson case CEED_EVAL_DIV: 4396e15d496SJeremy L Thompson case CEED_EVAL_CURL: 4406574a04fSJeremy L Thompson // LCOV_EXCL_START 4416574a04fSJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_INCOMPATIBLE, "Tensor basis evaluation for %s not supported", CeedEvalModes[eval_mode]); 4422b730f8bSJeremy L Thompson break; 4436e15d496SJeremy L Thompson // LCOV_EXCL_STOP 4442b730f8bSJeremy L Thompson case CEED_EVAL_WEIGHT: 4452b730f8bSJeremy L Thompson *flops = dim * CeedIntPow(Q_1d, dim); 4462b730f8bSJeremy L Thompson break; 4476e15d496SJeremy L Thompson } 4486e15d496SJeremy L Thompson } else { 449c4e3f59bSSebastian Grimberg CeedInt dim, num_comp, q_comp, num_nodes, num_qpts; 4502b730f8bSJeremy L Thompson CeedCall(CeedBasisGetDimension(basis, &dim)); 4512b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumComponents(basis, &num_comp)); 452c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp)); 4532b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis, &num_nodes)); 4542b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts)); 4556e15d496SJeremy L Thompson switch (eval_mode) { 4562b730f8bSJeremy L Thompson case CEED_EVAL_NONE: 4572b730f8bSJeremy L Thompson *flops = 0; 4582b730f8bSJeremy L Thompson break; 4592b730f8bSJeremy L Thompson case CEED_EVAL_INTERP: 4602b730f8bSJeremy L Thompson case CEED_EVAL_GRAD: 4612b730f8bSJeremy L Thompson case CEED_EVAL_DIV: 4622b730f8bSJeremy L Thompson case CEED_EVAL_CURL: 463c4e3f59bSSebastian Grimberg *flops = num_nodes * num_qpts * num_comp * q_comp; 4642b730f8bSJeremy L Thompson break; 4652b730f8bSJeremy L Thompson case CEED_EVAL_WEIGHT: 4662b730f8bSJeremy L Thompson *flops = 0; 4672b730f8bSJeremy L Thompson break; 4686e15d496SJeremy L Thompson } 4696e15d496SJeremy L Thompson } 4706e15d496SJeremy L Thompson 4716e15d496SJeremy L Thompson return CEED_ERROR_SUCCESS; 4726e15d496SJeremy L Thompson } 4736e15d496SJeremy L Thompson 4746e15d496SJeremy L Thompson /** 475c4e3f59bSSebastian Grimberg @brief Get CeedFESpace for a CeedBasis 476c4e3f59bSSebastian Grimberg 477c4e3f59bSSebastian Grimberg @param[in] basis CeedBasis 478c4e3f59bSSebastian Grimberg @param[out] fe_space Variable to store CeedFESpace 479c4e3f59bSSebastian Grimberg 480c4e3f59bSSebastian Grimberg @return An error code: 0 - success, otherwise - failure 481c4e3f59bSSebastian Grimberg 482c4e3f59bSSebastian Grimberg @ref Backend 483c4e3f59bSSebastian Grimberg **/ 484c4e3f59bSSebastian Grimberg int CeedBasisGetFESpace(CeedBasis basis, CeedFESpace *fe_space) { 485c4e3f59bSSebastian Grimberg *fe_space = basis->fe_space; 486c4e3f59bSSebastian Grimberg return CEED_ERROR_SUCCESS; 487c4e3f59bSSebastian Grimberg } 488c4e3f59bSSebastian Grimberg 489c4e3f59bSSebastian Grimberg /** 4907a982d89SJeremy L. Thompson @brief Get dimension for given CeedElemTopology 4917a982d89SJeremy L. Thompson 492ea61e9acSJeremy L Thompson @param[in] topo CeedElemTopology 4937a982d89SJeremy L. Thompson @param[out] dim Variable to store dimension of topology 4947a982d89SJeremy L. Thompson 4957a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 4967a982d89SJeremy L. Thompson 4977a982d89SJeremy L. Thompson @ref Backend 4987a982d89SJeremy L. Thompson **/ 4997a982d89SJeremy L. Thompson int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) { 5007a982d89SJeremy L. Thompson *dim = (CeedInt)topo >> 16; 501e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 5027a982d89SJeremy L. Thompson } 5037a982d89SJeremy L. Thompson 5047a982d89SJeremy L. Thompson /** 5057a982d89SJeremy L. Thompson @brief Get CeedTensorContract of a CeedBasis 5067a982d89SJeremy L. Thompson 507ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 5087a982d89SJeremy L. Thompson @param[out] contract Variable to store CeedTensorContract 5097a982d89SJeremy L. Thompson 5107a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 5117a982d89SJeremy L. Thompson 5127a982d89SJeremy L. Thompson @ref Backend 5137a982d89SJeremy L. Thompson **/ 5147a982d89SJeremy L. Thompson int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) { 5157a982d89SJeremy L. Thompson *contract = basis->contract; 516e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 5177a982d89SJeremy L. Thompson } 5187a982d89SJeremy L. Thompson 5197a982d89SJeremy L. Thompson /** 5207a982d89SJeremy L. Thompson @brief Set CeedTensorContract of a CeedBasis 5217a982d89SJeremy L. Thompson 522ea61e9acSJeremy L Thompson @param[in,out] basis CeedBasis 523ea61e9acSJeremy L Thompson @param[in] contract CeedTensorContract to set 5247a982d89SJeremy L. Thompson 5257a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 5267a982d89SJeremy L. Thompson 5277a982d89SJeremy L. Thompson @ref Backend 5287a982d89SJeremy L. Thompson **/ 52934359f16Sjeremylt int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract contract) { 53034359f16Sjeremylt basis->contract = contract; 5312b730f8bSJeremy L Thompson CeedCall(CeedTensorContractReference(contract)); 532e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 5337a982d89SJeremy L. Thompson } 5347a982d89SJeremy L. Thompson 5357a982d89SJeremy L. Thompson /** 5367a982d89SJeremy L. Thompson @brief Return a reference implementation of matrix multiplication C = A B. 537ba59ac12SSebastian Grimberg 538ba59ac12SSebastian Grimberg Note: This is a reference implementation for CPU CeedScalar pointers that is not intended for high performance. 5397a982d89SJeremy L. Thompson 540ea61e9acSJeremy L Thompson @param[in] ceed Ceed context for error handling 541d1d35e2fSjeremylt @param[in] mat_A Row-major matrix A 542d1d35e2fSjeremylt @param[in] mat_B Row-major matrix B 543d1d35e2fSjeremylt @param[out] mat_C Row-major output matrix C 544ea61e9acSJeremy L Thompson @param[in] m Number of rows of C 545ea61e9acSJeremy L Thompson @param[in] n Number of columns of C 546ea61e9acSJeremy L Thompson @param[in] kk Number of columns of A/rows of B 5477a982d89SJeremy L. Thompson 5487a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 5497a982d89SJeremy L. Thompson 5507a982d89SJeremy L. Thompson @ref Utility 5517a982d89SJeremy L. Thompson **/ 5522b730f8bSJeremy L Thompson int CeedMatrixMatrixMultiply(Ceed ceed, const CeedScalar *mat_A, const CeedScalar *mat_B, CeedScalar *mat_C, CeedInt m, CeedInt n, CeedInt kk) { 5532b730f8bSJeremy L Thompson for (CeedInt i = 0; i < m; i++) { 5547a982d89SJeremy L. Thompson for (CeedInt j = 0; j < n; j++) { 5557a982d89SJeremy L. Thompson CeedScalar sum = 0; 5562b730f8bSJeremy L Thompson for (CeedInt k = 0; k < kk; k++) sum += mat_A[k + i * kk] * mat_B[j + k * n]; 557d1d35e2fSjeremylt mat_C[j + i * n] = sum; 5587a982d89SJeremy L. Thompson } 5592b730f8bSJeremy L Thompson } 560e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 5617a982d89SJeremy L. Thompson } 5627a982d89SJeremy L. Thompson 563ba59ac12SSebastian Grimberg /** 564ba59ac12SSebastian Grimberg @brief Return QR Factorization of a matrix 565ba59ac12SSebastian Grimberg 566ba59ac12SSebastian Grimberg @param[in] ceed Ceed context for error handling 567ba59ac12SSebastian Grimberg @param[in,out] mat Row-major matrix to be factorized in place 568ba59ac12SSebastian Grimberg @param[in,out] tau Vector of length m of scaling factors 569ba59ac12SSebastian Grimberg @param[in] m Number of rows 570ba59ac12SSebastian Grimberg @param[in] n Number of columns 571ba59ac12SSebastian Grimberg 572ba59ac12SSebastian Grimberg @return An error code: 0 - success, otherwise - failure 573ba59ac12SSebastian Grimberg 574ba59ac12SSebastian Grimberg @ref Utility 575ba59ac12SSebastian Grimberg **/ 576ba59ac12SSebastian Grimberg int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau, CeedInt m, CeedInt n) { 577ba59ac12SSebastian Grimberg CeedScalar v[m]; 578ba59ac12SSebastian Grimberg 579ba59ac12SSebastian Grimberg // Check matrix shape 5806574a04fSJeremy L Thompson CeedCheck(n <= m, ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute QR factorization with n > m"); 581ba59ac12SSebastian Grimberg 582ba59ac12SSebastian Grimberg for (CeedInt i = 0; i < n; i++) { 583ba59ac12SSebastian Grimberg if (i >= m - 1) { // last row of matrix, no reflection needed 584ba59ac12SSebastian Grimberg tau[i] = 0.; 585ba59ac12SSebastian Grimberg break; 586ba59ac12SSebastian Grimberg } 587ba59ac12SSebastian Grimberg // Calculate Householder vector, magnitude 588ba59ac12SSebastian Grimberg CeedScalar sigma = 0.0; 589ba59ac12SSebastian Grimberg v[i] = mat[i + n * i]; 590ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < m; j++) { 591ba59ac12SSebastian Grimberg v[j] = mat[i + n * j]; 592ba59ac12SSebastian Grimberg sigma += v[j] * v[j]; 593ba59ac12SSebastian Grimberg } 594ba59ac12SSebastian Grimberg CeedScalar norm = sqrt(v[i] * v[i] + sigma); // norm of v[i:m] 595ba59ac12SSebastian Grimberg CeedScalar R_ii = -copysign(norm, v[i]); 596ba59ac12SSebastian Grimberg v[i] -= R_ii; 597ba59ac12SSebastian Grimberg // norm of v[i:m] after modification above and scaling below 598ba59ac12SSebastian Grimberg // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 599ba59ac12SSebastian Grimberg // tau = 2 / (norm*norm) 600ba59ac12SSebastian Grimberg tau[i] = 2 * v[i] * v[i] / (v[i] * v[i] + sigma); 601ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < m; j++) v[j] /= v[i]; 602ba59ac12SSebastian Grimberg 603ba59ac12SSebastian Grimberg // Apply Householder reflector to lower right panel 604ba59ac12SSebastian Grimberg CeedHouseholderReflect(&mat[i * n + i + 1], &v[i], tau[i], m - i, n - i - 1, n, 1); 605ba59ac12SSebastian Grimberg // Save v 606ba59ac12SSebastian Grimberg mat[i + n * i] = R_ii; 607ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < m; j++) mat[i + n * j] = v[j]; 608ba59ac12SSebastian Grimberg } 609ba59ac12SSebastian Grimberg return CEED_ERROR_SUCCESS; 610ba59ac12SSebastian Grimberg } 611ba59ac12SSebastian Grimberg 612ba59ac12SSebastian Grimberg /** 613ba59ac12SSebastian Grimberg @brief Apply Householder Q matrix 614ba59ac12SSebastian Grimberg 615ba59ac12SSebastian Grimberg Compute mat_A = mat_Q mat_A, where mat_Q is mxm and mat_A is mxn. 616ba59ac12SSebastian Grimberg 617ba59ac12SSebastian Grimberg @param[in,out] mat_A Matrix to apply Householder Q to, in place 618ba59ac12SSebastian Grimberg @param[in] mat_Q Householder Q matrix 619ba59ac12SSebastian Grimberg @param[in] tau Householder scaling factors 620ba59ac12SSebastian Grimberg @param[in] t_mode Transpose mode for application 621ba59ac12SSebastian Grimberg @param[in] m Number of rows in A 622ba59ac12SSebastian Grimberg @param[in] n Number of columns in A 623ba59ac12SSebastian Grimberg @param[in] k Number of elementary reflectors in Q, k<m 624ba59ac12SSebastian Grimberg @param[in] row Row stride in A 625ba59ac12SSebastian Grimberg @param[in] col Col stride in A 626ba59ac12SSebastian Grimberg 627ba59ac12SSebastian Grimberg @return An error code: 0 - success, otherwise - failure 628ba59ac12SSebastian Grimberg 629c4e3f59bSSebastian Grimberg @ref Utility 630ba59ac12SSebastian Grimberg **/ 631ba59ac12SSebastian Grimberg int CeedHouseholderApplyQ(CeedScalar *mat_A, const CeedScalar *mat_Q, const CeedScalar *tau, CeedTransposeMode t_mode, CeedInt m, CeedInt n, 632ba59ac12SSebastian Grimberg CeedInt k, CeedInt row, CeedInt col) { 633ba59ac12SSebastian Grimberg CeedScalar *v; 634ba59ac12SSebastian Grimberg CeedCall(CeedMalloc(m, &v)); 635ba59ac12SSebastian Grimberg for (CeedInt ii = 0; ii < k; ii++) { 636ba59ac12SSebastian Grimberg CeedInt i = t_mode == CEED_TRANSPOSE ? ii : k - 1 - ii; 637ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < m; j++) v[j] = mat_Q[j * k + i]; 638ba59ac12SSebastian Grimberg // Apply Householder reflector (I - tau v v^T) collo_grad_1d^T 639ba59ac12SSebastian Grimberg CeedCall(CeedHouseholderReflect(&mat_A[i * row], &v[i], tau[i], m - i, n, row, col)); 640ba59ac12SSebastian Grimberg } 641ba59ac12SSebastian Grimberg CeedCall(CeedFree(&v)); 642ba59ac12SSebastian Grimberg return CEED_ERROR_SUCCESS; 643ba59ac12SSebastian Grimberg } 644ba59ac12SSebastian Grimberg 645ba59ac12SSebastian Grimberg /** 646ba59ac12SSebastian Grimberg @brief Return symmetric Schur decomposition of the symmetric matrix mat via symmetric QR factorization 647ba59ac12SSebastian Grimberg 648ba59ac12SSebastian Grimberg @param[in] ceed Ceed context for error handling 649ba59ac12SSebastian Grimberg @param[in,out] mat Row-major matrix to be factorized in place 650ba59ac12SSebastian Grimberg @param[out] lambda Vector of length n of eigenvalues 651ba59ac12SSebastian Grimberg @param[in] n Number of rows/columns 652ba59ac12SSebastian Grimberg 653ba59ac12SSebastian Grimberg @return An error code: 0 - success, otherwise - failure 654ba59ac12SSebastian Grimberg 655ba59ac12SSebastian Grimberg @ref Utility 656ba59ac12SSebastian Grimberg **/ 6572c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOff 6582c2ea1dbSJeremy L Thompson int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat, CeedScalar *lambda, CeedInt n) { 659ba59ac12SSebastian Grimberg // Check bounds for clang-tidy 6606574a04fSJeremy L Thompson CeedCheck(n > 1, ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute symmetric Schur decomposition of scalars"); 661ba59ac12SSebastian Grimberg 662ba59ac12SSebastian Grimberg CeedScalar v[n - 1], tau[n - 1], mat_T[n * n]; 663ba59ac12SSebastian Grimberg 664ba59ac12SSebastian Grimberg // Copy mat to mat_T and set mat to I 665ba59ac12SSebastian Grimberg memcpy(mat_T, mat, n * n * sizeof(mat[0])); 666ba59ac12SSebastian Grimberg for (CeedInt i = 0; i < n; i++) { 667ba59ac12SSebastian Grimberg for (CeedInt j = 0; j < n; j++) mat[j + n * i] = (i == j) ? 1 : 0; 668ba59ac12SSebastian Grimberg } 669ba59ac12SSebastian Grimberg 670ba59ac12SSebastian Grimberg // Reduce to tridiagonal 671ba59ac12SSebastian Grimberg for (CeedInt i = 0; i < n - 1; i++) { 672ba59ac12SSebastian Grimberg // Calculate Householder vector, magnitude 673ba59ac12SSebastian Grimberg CeedScalar sigma = 0.0; 674ba59ac12SSebastian Grimberg v[i] = mat_T[i + n * (i + 1)]; 675ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < n - 1; j++) { 676ba59ac12SSebastian Grimberg v[j] = mat_T[i + n * (j + 1)]; 677ba59ac12SSebastian Grimberg sigma += v[j] * v[j]; 678ba59ac12SSebastian Grimberg } 679ba59ac12SSebastian Grimberg CeedScalar norm = sqrt(v[i] * v[i] + sigma); // norm of v[i:n-1] 680ba59ac12SSebastian Grimberg CeedScalar R_ii = -copysign(norm, v[i]); 681ba59ac12SSebastian Grimberg v[i] -= R_ii; 682ba59ac12SSebastian Grimberg // norm of v[i:m] after modification above and scaling below 683ba59ac12SSebastian Grimberg // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 684ba59ac12SSebastian Grimberg // tau = 2 / (norm*norm) 685ba59ac12SSebastian Grimberg tau[i] = i == n - 2 ? 2 : 2 * v[i] * v[i] / (v[i] * v[i] + sigma); 686ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < n - 1; j++) v[j] /= v[i]; 687ba59ac12SSebastian Grimberg 688ba59ac12SSebastian Grimberg // Update sub and super diagonal 689ba59ac12SSebastian Grimberg for (CeedInt j = i + 2; j < n; j++) { 690ba59ac12SSebastian Grimberg mat_T[i + n * j] = 0; 691ba59ac12SSebastian Grimberg mat_T[j + n * i] = 0; 692ba59ac12SSebastian Grimberg } 693ba59ac12SSebastian Grimberg // Apply symmetric Householder reflector to lower right panel 694ba59ac12SSebastian Grimberg CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1); 695ba59ac12SSebastian Grimberg CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), 1, n); 696ba59ac12SSebastian Grimberg 697ba59ac12SSebastian Grimberg // Save v 698ba59ac12SSebastian Grimberg mat_T[i + n * (i + 1)] = R_ii; 699ba59ac12SSebastian Grimberg mat_T[(i + 1) + n * i] = R_ii; 700ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < n - 1; j++) { 701ba59ac12SSebastian Grimberg mat_T[i + n * (j + 1)] = v[j]; 702ba59ac12SSebastian Grimberg } 703ba59ac12SSebastian Grimberg } 704ba59ac12SSebastian Grimberg // Backwards accumulation of Q 705ba59ac12SSebastian Grimberg for (CeedInt i = n - 2; i >= 0; i--) { 706ba59ac12SSebastian Grimberg if (tau[i] > 0.0) { 707ba59ac12SSebastian Grimberg v[i] = 1; 708ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < n - 1; j++) { 709ba59ac12SSebastian Grimberg v[j] = mat_T[i + n * (j + 1)]; 710ba59ac12SSebastian Grimberg mat_T[i + n * (j + 1)] = 0; 711ba59ac12SSebastian Grimberg } 712ba59ac12SSebastian Grimberg CeedHouseholderReflect(&mat[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1); 713ba59ac12SSebastian Grimberg } 714ba59ac12SSebastian Grimberg } 715ba59ac12SSebastian Grimberg 716ba59ac12SSebastian Grimberg // Reduce sub and super diagonal 717ba59ac12SSebastian Grimberg CeedInt p = 0, q = 0, itr = 0, max_itr = n * n * n * n; 718ba59ac12SSebastian Grimberg CeedScalar tol = CEED_EPSILON; 719ba59ac12SSebastian Grimberg 720ba59ac12SSebastian Grimberg while (itr < max_itr) { 721ba59ac12SSebastian Grimberg // Update p, q, size of reduced portions of diagonal 722ba59ac12SSebastian Grimberg p = 0; 723ba59ac12SSebastian Grimberg q = 0; 724ba59ac12SSebastian Grimberg for (CeedInt i = n - 2; i >= 0; i--) { 725ba59ac12SSebastian Grimberg if (fabs(mat_T[i + n * (i + 1)]) < tol) q += 1; 726ba59ac12SSebastian Grimberg else break; 727ba59ac12SSebastian Grimberg } 728ba59ac12SSebastian Grimberg for (CeedInt i = 0; i < n - q - 1; i++) { 729ba59ac12SSebastian Grimberg if (fabs(mat_T[i + n * (i + 1)]) < tol) p += 1; 730ba59ac12SSebastian Grimberg else break; 731ba59ac12SSebastian Grimberg } 732ba59ac12SSebastian Grimberg if (q == n - 1) break; // Finished reducing 733ba59ac12SSebastian Grimberg 734ba59ac12SSebastian Grimberg // Reduce tridiagonal portion 735ba59ac12SSebastian Grimberg CeedScalar t_nn = mat_T[(n - 1 - q) + n * (n - 1 - q)], t_nnm1 = mat_T[(n - 2 - q) + n * (n - 1 - q)]; 736ba59ac12SSebastian Grimberg CeedScalar d = (mat_T[(n - 2 - q) + n * (n - 2 - q)] - t_nn) / 2; 737ba59ac12SSebastian Grimberg CeedScalar mu = t_nn - t_nnm1 * t_nnm1 / (d + copysign(sqrt(d * d + t_nnm1 * t_nnm1), d)); 738ba59ac12SSebastian Grimberg CeedScalar x = mat_T[p + n * p] - mu; 739ba59ac12SSebastian Grimberg CeedScalar z = mat_T[p + n * (p + 1)]; 740ba59ac12SSebastian Grimberg for (CeedInt k = p; k < n - q - 1; k++) { 741ba59ac12SSebastian Grimberg // Compute Givens rotation 742ba59ac12SSebastian Grimberg CeedScalar c = 1, s = 0; 743ba59ac12SSebastian Grimberg if (fabs(z) > tol) { 744ba59ac12SSebastian Grimberg if (fabs(z) > fabs(x)) { 745ba59ac12SSebastian Grimberg CeedScalar tau = -x / z; 746ba59ac12SSebastian Grimberg s = 1 / sqrt(1 + tau * tau), c = s * tau; 747ba59ac12SSebastian Grimberg } else { 748ba59ac12SSebastian Grimberg CeedScalar tau = -z / x; 749ba59ac12SSebastian Grimberg c = 1 / sqrt(1 + tau * tau), s = c * tau; 750ba59ac12SSebastian Grimberg } 751ba59ac12SSebastian Grimberg } 752ba59ac12SSebastian Grimberg 753ba59ac12SSebastian Grimberg // Apply Givens rotation to T 754ba59ac12SSebastian Grimberg CeedGivensRotation(mat_T, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n); 755ba59ac12SSebastian Grimberg CeedGivensRotation(mat_T, c, s, CEED_TRANSPOSE, k, k + 1, n, n); 756ba59ac12SSebastian Grimberg 757ba59ac12SSebastian Grimberg // Apply Givens rotation to Q 758ba59ac12SSebastian Grimberg CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n); 759ba59ac12SSebastian Grimberg 760ba59ac12SSebastian Grimberg // Update x, z 761ba59ac12SSebastian Grimberg if (k < n - q - 2) { 762ba59ac12SSebastian Grimberg x = mat_T[k + n * (k + 1)]; 763ba59ac12SSebastian Grimberg z = mat_T[k + n * (k + 2)]; 764ba59ac12SSebastian Grimberg } 765ba59ac12SSebastian Grimberg } 766ba59ac12SSebastian Grimberg itr++; 767ba59ac12SSebastian Grimberg } 768ba59ac12SSebastian Grimberg 769ba59ac12SSebastian Grimberg // Save eigenvalues 770ba59ac12SSebastian Grimberg for (CeedInt i = 0; i < n; i++) lambda[i] = mat_T[i + n * i]; 771ba59ac12SSebastian Grimberg 772ba59ac12SSebastian Grimberg // Check convergence 7736574a04fSJeremy L Thompson CeedCheck(itr < max_itr || q > n, ceed, CEED_ERROR_MINOR, "Symmetric QR failed to converge"); 774ba59ac12SSebastian Grimberg return CEED_ERROR_SUCCESS; 775ba59ac12SSebastian Grimberg } 7762c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOn 777ba59ac12SSebastian Grimberg 778ba59ac12SSebastian Grimberg /** 779ba59ac12SSebastian Grimberg @brief Return Simultaneous Diagonalization of two matrices. 780ba59ac12SSebastian Grimberg 781ba59ac12SSebastian Grimberg This solves the generalized eigenvalue problem A x = lambda B x, where A and B are symmetric and B is positive definite. 782ba59ac12SSebastian Grimberg We generate the matrix X and vector Lambda such that X^T A X = Lambda and X^T B X = I. 783ba59ac12SSebastian Grimberg This is equivalent to the LAPACK routine 'sygv' with TYPE = 1. 784ba59ac12SSebastian Grimberg 785ba59ac12SSebastian Grimberg @param[in] ceed Ceed context for error handling 786ba59ac12SSebastian Grimberg @param[in] mat_A Row-major matrix to be factorized with eigenvalues 787ba59ac12SSebastian Grimberg @param[in] mat_B Row-major matrix to be factorized to identity 788ba59ac12SSebastian Grimberg @param[out] mat_X Row-major orthogonal matrix 789ba59ac12SSebastian Grimberg @param[out] lambda Vector of length n of generalized eigenvalues 790ba59ac12SSebastian Grimberg @param[in] n Number of rows/columns 791ba59ac12SSebastian Grimberg 792ba59ac12SSebastian Grimberg @return An error code: 0 - success, otherwise - failure 793ba59ac12SSebastian Grimberg 794ba59ac12SSebastian Grimberg @ref Utility 795ba59ac12SSebastian Grimberg **/ 7962c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOff 7972c2ea1dbSJeremy L Thompson int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *mat_A, CeedScalar *mat_B, CeedScalar *mat_X, CeedScalar *lambda, CeedInt n) { 798ba59ac12SSebastian Grimberg CeedScalar *mat_C, *mat_G, *vec_D; 799ba59ac12SSebastian Grimberg CeedCall(CeedCalloc(n * n, &mat_C)); 800ba59ac12SSebastian Grimberg CeedCall(CeedCalloc(n * n, &mat_G)); 801ba59ac12SSebastian Grimberg CeedCall(CeedCalloc(n, &vec_D)); 802ba59ac12SSebastian Grimberg 803ba59ac12SSebastian Grimberg // Compute B = G D G^T 804ba59ac12SSebastian Grimberg memcpy(mat_G, mat_B, n * n * sizeof(mat_B[0])); 805ba59ac12SSebastian Grimberg CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_G, vec_D, n)); 806ba59ac12SSebastian Grimberg 807ba59ac12SSebastian Grimberg // Sort eigenvalues 808ba59ac12SSebastian Grimberg for (CeedInt i = n - 1; i >= 0; i--) { 809ba59ac12SSebastian Grimberg for (CeedInt j = 0; j < i; j++) { 810ba59ac12SSebastian Grimberg if (fabs(vec_D[j]) > fabs(vec_D[j + 1])) { 811ba59ac12SSebastian Grimberg CeedScalar temp; 812ba59ac12SSebastian Grimberg temp = vec_D[j]; 813ba59ac12SSebastian Grimberg vec_D[j] = vec_D[j + 1]; 814ba59ac12SSebastian Grimberg vec_D[j + 1] = temp; 815ba59ac12SSebastian Grimberg for (CeedInt k = 0; k < n; k++) { 816ba59ac12SSebastian Grimberg temp = mat_G[k * n + j]; 817ba59ac12SSebastian Grimberg mat_G[k * n + j] = mat_G[k * n + j + 1]; 818ba59ac12SSebastian Grimberg mat_G[k * n + j + 1] = temp; 819ba59ac12SSebastian Grimberg } 820ba59ac12SSebastian Grimberg } 821ba59ac12SSebastian Grimberg } 822ba59ac12SSebastian Grimberg } 823ba59ac12SSebastian Grimberg 824ba59ac12SSebastian Grimberg // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T 825ba59ac12SSebastian Grimberg // = D^-1/2 G^T A G D^-1/2 826ba59ac12SSebastian Grimberg // -- D = D^-1/2 827ba59ac12SSebastian Grimberg for (CeedInt i = 0; i < n; i++) vec_D[i] = 1. / sqrt(vec_D[i]); 828ba59ac12SSebastian Grimberg // -- G = G D^-1/2 829ba59ac12SSebastian Grimberg // -- C = D^-1/2 G^T 830ba59ac12SSebastian Grimberg for (CeedInt i = 0; i < n; i++) { 831ba59ac12SSebastian Grimberg for (CeedInt j = 0; j < n; j++) { 832ba59ac12SSebastian Grimberg mat_G[i * n + j] *= vec_D[j]; 833ba59ac12SSebastian Grimberg mat_C[j * n + i] = mat_G[i * n + j]; 834ba59ac12SSebastian Grimberg } 835ba59ac12SSebastian Grimberg } 836ba59ac12SSebastian Grimberg // -- X = (D^-1/2 G^T) A 837ba59ac12SSebastian Grimberg CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_C, (const CeedScalar *)mat_A, mat_X, n, n, n)); 838ba59ac12SSebastian Grimberg // -- C = (D^-1/2 G^T A) (G D^-1/2) 839ba59ac12SSebastian Grimberg CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_X, (const CeedScalar *)mat_G, mat_C, n, n, n)); 840ba59ac12SSebastian Grimberg 841ba59ac12SSebastian Grimberg // Compute Q^T C Q = lambda 842ba59ac12SSebastian Grimberg CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_C, lambda, n)); 843ba59ac12SSebastian Grimberg 844ba59ac12SSebastian Grimberg // Sort eigenvalues 845ba59ac12SSebastian Grimberg for (CeedInt i = n - 1; i >= 0; i--) { 846ba59ac12SSebastian Grimberg for (CeedInt j = 0; j < i; j++) { 847ba59ac12SSebastian Grimberg if (fabs(lambda[j]) > fabs(lambda[j + 1])) { 848ba59ac12SSebastian Grimberg CeedScalar temp; 849ba59ac12SSebastian Grimberg temp = lambda[j]; 850ba59ac12SSebastian Grimberg lambda[j] = lambda[j + 1]; 851ba59ac12SSebastian Grimberg lambda[j + 1] = temp; 852ba59ac12SSebastian Grimberg for (CeedInt k = 0; k < n; k++) { 853ba59ac12SSebastian Grimberg temp = mat_C[k * n + j]; 854ba59ac12SSebastian Grimberg mat_C[k * n + j] = mat_C[k * n + j + 1]; 855ba59ac12SSebastian Grimberg mat_C[k * n + j + 1] = temp; 856ba59ac12SSebastian Grimberg } 857ba59ac12SSebastian Grimberg } 858ba59ac12SSebastian Grimberg } 859ba59ac12SSebastian Grimberg } 860ba59ac12SSebastian Grimberg 861ba59ac12SSebastian Grimberg // Set X = (G D^1/2)^-T Q 862ba59ac12SSebastian Grimberg // = G D^-1/2 Q 863ba59ac12SSebastian Grimberg CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_G, (const CeedScalar *)mat_C, mat_X, n, n, n)); 864ba59ac12SSebastian Grimberg 865ba59ac12SSebastian Grimberg // Cleanup 866ba59ac12SSebastian Grimberg CeedCall(CeedFree(&mat_C)); 867ba59ac12SSebastian Grimberg CeedCall(CeedFree(&mat_G)); 868ba59ac12SSebastian Grimberg CeedCall(CeedFree(&vec_D)); 869ba59ac12SSebastian Grimberg return CEED_ERROR_SUCCESS; 870ba59ac12SSebastian Grimberg } 8712c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOn 872ba59ac12SSebastian Grimberg 8737a982d89SJeremy L. Thompson /// @} 8747a982d89SJeremy L. Thompson 8757a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 8767a982d89SJeremy L. Thompson /// CeedBasis Public API 8777a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 8787a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser 879d7b241e6Sjeremylt /// @{ 880d7b241e6Sjeremylt 881b11c1e72Sjeremylt /** 882ba59ac12SSebastian Grimberg @brief Create a tensor-product basis for H^1 discretizations 883b11c1e72Sjeremylt 884ea61e9acSJeremy L Thompson @param[in] ceed Ceed object where the CeedBasis will be created 885ea61e9acSJeremy L Thompson @param[in] dim Topological dimension 886ea61e9acSJeremy L Thompson @param[in] num_comp Number of field components (1 for scalar fields) 887ea61e9acSJeremy L Thompson @param[in] P_1d Number of nodes in one dimension 888ea61e9acSJeremy L Thompson @param[in] Q_1d Number of quadrature points in one dimension 889ea61e9acSJeremy L Thompson @param[in] interp_1d Row-major (Q_1d * P_1d) matrix expressing the values of nodal basis functions at quadrature points 890ea61e9acSJeremy L Thompson @param[in] grad_1d Row-major (Q_1d * P_1d) matrix expressing derivatives of nodal basis functions at quadrature points 891ea61e9acSJeremy 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] 892ea61e9acSJeremy L Thompson @param[in] q_weight_1d Array of length Q_1d holding the quadrature weights on the reference element 893ea61e9acSJeremy L Thompson @param[out] basis Address of the variable where the newly created CeedBasis will be stored. 894b11c1e72Sjeremylt 895b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 896dfdf5a53Sjeremylt 8977a982d89SJeremy L. Thompson @ref User 898b11c1e72Sjeremylt **/ 8992b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, 9002b730f8bSJeremy L Thompson const CeedScalar *grad_1d, const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis *basis) { 9015fe0d4faSjeremylt if (!ceed->BasisCreateTensorH1) { 9025fe0d4faSjeremylt Ceed delegate; 9036574a04fSJeremy L Thompson 9042b730f8bSJeremy L Thompson CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis")); 9056574a04fSJeremy L Thompson CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisCreateTensorH1"); 9062b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateTensorH1(delegate, dim, num_comp, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis)); 907e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 9085fe0d4faSjeremylt } 909e15f9bd0SJeremy L Thompson 9106574a04fSJeremy L Thompson CeedCheck(dim > 0, ceed, CEED_ERROR_DIMENSION, "Basis dimension must be a positive value"); 9116574a04fSJeremy L Thompson CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component"); 9126574a04fSJeremy L Thompson CeedCheck(P_1d > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node"); 9136574a04fSJeremy L Thompson CeedCheck(Q_1d > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point"); 914227444bfSJeremy L Thompson 9152b730f8bSJeremy L Thompson CeedElemTopology topo = dim == 1 ? CEED_TOPOLOGY_LINE : dim == 2 ? CEED_TOPOLOGY_QUAD : CEED_TOPOLOGY_HEX; 916e15f9bd0SJeremy L Thompson 9172b730f8bSJeremy L Thompson CeedCall(CeedCalloc(1, basis)); 918db002c03SJeremy L Thompson CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed)); 919d1d35e2fSjeremylt (*basis)->ref_count = 1; 9206402da51SJeremy L Thompson (*basis)->is_tensor_basis = true; 921d7b241e6Sjeremylt (*basis)->dim = dim; 922d99fa3c5SJeremy L Thompson (*basis)->topo = topo; 923d1d35e2fSjeremylt (*basis)->num_comp = num_comp; 924d1d35e2fSjeremylt (*basis)->P_1d = P_1d; 925d1d35e2fSjeremylt (*basis)->Q_1d = Q_1d; 926d1d35e2fSjeremylt (*basis)->P = CeedIntPow(P_1d, dim); 927d1d35e2fSjeremylt (*basis)->Q = CeedIntPow(Q_1d, dim); 928c4e3f59bSSebastian Grimberg (*basis)->fe_space = CEED_FE_SPACE_H1; 9292b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q_1d, &(*basis)->q_ref_1d)); 9302b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q_1d, &(*basis)->q_weight_1d)); 931ff3a0f91SJeremy L Thompson if (q_ref_1d) memcpy((*basis)->q_ref_1d, q_ref_1d, Q_1d * sizeof(q_ref_1d[0])); 9322b730f8bSJeremy L Thompson if (q_weight_1d) memcpy((*basis)->q_weight_1d, q_weight_1d, Q_1d * sizeof(q_weight_1d[0])); 9332b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->interp_1d)); 9342b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->grad_1d)); 9352b730f8bSJeremy L Thompson if (interp_1d) memcpy((*basis)->interp_1d, interp_1d, Q_1d * P_1d * sizeof(interp_1d[0])); 936ff3a0f91SJeremy L Thompson if (grad_1d) memcpy((*basis)->grad_1d, grad_1d, Q_1d * P_1d * sizeof(grad_1d[0])); 9372b730f8bSJeremy L Thompson CeedCall(ceed->BasisCreateTensorH1(dim, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, *basis)); 938e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 939d7b241e6Sjeremylt } 940d7b241e6Sjeremylt 941b11c1e72Sjeremylt /** 94295bb1877Svaleriabarra @brief Create a tensor-product Lagrange basis 943b11c1e72Sjeremylt 944ea61e9acSJeremy L Thompson @param[in] ceed Ceed object where the CeedBasis will be created 945ea61e9acSJeremy L Thompson @param[in] dim Topological dimension of element 946ea61e9acSJeremy L Thompson @param[in] num_comp Number of field components (1 for scalar fields) 947ea61e9acSJeremy L Thompson @param[in] P Number of Gauss-Lobatto nodes in one dimension. 948ea61e9acSJeremy L Thompson The polynomial degree of the resulting Q_k element is k=P-1. 949ea61e9acSJeremy L Thompson @param[in] Q Number of quadrature points in one dimension. 950ea61e9acSJeremy L Thompson @param[in] quad_mode Distribution of the Q quadrature points (affects order of accuracy for the quadrature) 951ea61e9acSJeremy L Thompson @param[out] basis Address of the variable where the newly created CeedBasis will be stored. 952b11c1e72Sjeremylt 953b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 954dfdf5a53Sjeremylt 9557a982d89SJeremy L. Thompson @ref User 956b11c1e72Sjeremylt **/ 9572b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P, CeedInt Q, CeedQuadMode quad_mode, CeedBasis *basis) { 958d7b241e6Sjeremylt // Allocate 959c8c3fa7dSJeremy L Thompson int ierr = CEED_ERROR_SUCCESS; 9602b730f8bSJeremy L Thompson CeedScalar c1, c2, c3, c4, dx, *nodes, *interp_1d, *grad_1d, *q_ref_1d, *q_weight_1d; 9614d537eeaSYohann 9626574a04fSJeremy L Thompson CeedCheck(dim > 0, ceed, CEED_ERROR_DIMENSION, "Basis dimension must be a positive value"); 9636574a04fSJeremy L Thompson CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component"); 9646574a04fSJeremy L Thompson CeedCheck(P > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node"); 9656574a04fSJeremy L Thompson CeedCheck(Q > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point"); 966227444bfSJeremy L Thompson 967e15f9bd0SJeremy L Thompson // Get Nodes and Weights 9682b730f8bSJeremy L Thompson CeedCall(CeedCalloc(P * Q, &interp_1d)); 9692b730f8bSJeremy L Thompson CeedCall(CeedCalloc(P * Q, &grad_1d)); 9702b730f8bSJeremy L Thompson CeedCall(CeedCalloc(P, &nodes)); 9712b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q, &q_ref_1d)); 9722b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q, &q_weight_1d)); 9732b730f8bSJeremy L Thompson if (CeedLobattoQuadrature(P, nodes, NULL) != CEED_ERROR_SUCCESS) goto cleanup; 974d1d35e2fSjeremylt switch (quad_mode) { 975d7b241e6Sjeremylt case CEED_GAUSS: 976d1d35e2fSjeremylt ierr = CeedGaussQuadrature(Q, q_ref_1d, q_weight_1d); 977d7b241e6Sjeremylt break; 978d7b241e6Sjeremylt case CEED_GAUSS_LOBATTO: 979d1d35e2fSjeremylt ierr = CeedLobattoQuadrature(Q, q_ref_1d, q_weight_1d); 980d7b241e6Sjeremylt break; 981d7b241e6Sjeremylt } 9822b730f8bSJeremy L Thompson if (ierr != CEED_ERROR_SUCCESS) goto cleanup; 983e15f9bd0SJeremy L Thompson 984d7b241e6Sjeremylt // Build B, D matrix 985d7b241e6Sjeremylt // Fornberg, 1998 986c8c3fa7dSJeremy L Thompson for (CeedInt i = 0; i < Q; i++) { 987d7b241e6Sjeremylt c1 = 1.0; 988d1d35e2fSjeremylt c3 = nodes[0] - q_ref_1d[i]; 989d1d35e2fSjeremylt interp_1d[i * P + 0] = 1.0; 990c8c3fa7dSJeremy L Thompson for (CeedInt j = 1; j < P; j++) { 991d7b241e6Sjeremylt c2 = 1.0; 992d7b241e6Sjeremylt c4 = c3; 993d1d35e2fSjeremylt c3 = nodes[j] - q_ref_1d[i]; 994c8c3fa7dSJeremy L Thompson for (CeedInt k = 0; k < j; k++) { 995d7b241e6Sjeremylt dx = nodes[j] - nodes[k]; 996d7b241e6Sjeremylt c2 *= dx; 997d7b241e6Sjeremylt if (k == j - 1) { 998d1d35e2fSjeremylt grad_1d[i * P + j] = c1 * (interp_1d[i * P + k] - c4 * grad_1d[i * P + k]) / c2; 999d1d35e2fSjeremylt interp_1d[i * P + j] = -c1 * c4 * interp_1d[i * P + k] / c2; 1000d7b241e6Sjeremylt } 1001d1d35e2fSjeremylt grad_1d[i * P + k] = (c3 * grad_1d[i * P + k] - interp_1d[i * P + k]) / dx; 1002d1d35e2fSjeremylt interp_1d[i * P + k] = c3 * interp_1d[i * P + k] / dx; 1003d7b241e6Sjeremylt } 1004d7b241e6Sjeremylt c1 = c2; 1005d7b241e6Sjeremylt } 1006d7b241e6Sjeremylt } 10079ac7b42eSJeremy L Thompson // Pass to CeedBasisCreateTensorH1 10082b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P, Q, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis)); 1009e15f9bd0SJeremy L Thompson cleanup: 10102b730f8bSJeremy L Thompson CeedCall(CeedFree(&interp_1d)); 10112b730f8bSJeremy L Thompson CeedCall(CeedFree(&grad_1d)); 10122b730f8bSJeremy L Thompson CeedCall(CeedFree(&nodes)); 10132b730f8bSJeremy L Thompson CeedCall(CeedFree(&q_ref_1d)); 10142b730f8bSJeremy L Thompson CeedCall(CeedFree(&q_weight_1d)); 1015e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1016d7b241e6Sjeremylt } 1017d7b241e6Sjeremylt 1018b11c1e72Sjeremylt /** 1019ba59ac12SSebastian Grimberg @brief Create a non tensor-product basis for H^1 discretizations 1020a8de75f0Sjeremylt 1021ea61e9acSJeremy L Thompson @param[in] ceed Ceed object where the CeedBasis will be created 1022ea61e9acSJeremy L Thompson @param[in] topo Topology of element, e.g. hypercube, simplex, ect 1023ea61e9acSJeremy L Thompson @param[in] num_comp Number of field components (1 for scalar fields) 1024ea61e9acSJeremy L Thompson @param[in] num_nodes Total number of nodes 1025ea61e9acSJeremy L Thompson @param[in] num_qpts Total number of quadrature points 1026ea61e9acSJeremy L Thompson @param[in] interp Row-major (num_qpts * num_nodes) matrix expressing the values of nodal basis functions at quadrature points 1027c4e3f59bSSebastian Grimberg @param[in] grad Row-major (dim * num_qpts * num_nodes) matrix expressing derivatives of nodal basis functions at quadrature points 10289fe083eeSJeremy L Thompson @param[in] q_ref Array of length num_qpts * dim holding the locations of quadrature points on the reference element 1029ea61e9acSJeremy L Thompson @param[in] q_weight Array of length num_qpts holding the quadrature weights on the reference element 1030ea61e9acSJeremy L Thompson @param[out] basis Address of the variable where the newly created CeedBasis will be stored. 1031a8de75f0Sjeremylt 1032a8de75f0Sjeremylt @return An error code: 0 - success, otherwise - failure 1033a8de75f0Sjeremylt 10347a982d89SJeremy L. Thompson @ref User 1035a8de75f0Sjeremylt **/ 10362b730f8bSJeremy L Thompson int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 10372b730f8bSJeremy L Thompson const CeedScalar *grad, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) { 1038d1d35e2fSjeremylt CeedInt P = num_nodes, Q = num_qpts, dim = 0; 1039a8de75f0Sjeremylt 10405fe0d4faSjeremylt if (!ceed->BasisCreateH1) { 10415fe0d4faSjeremylt Ceed delegate; 10426574a04fSJeremy L Thompson 10432b730f8bSJeremy L Thompson CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis")); 10446574a04fSJeremy L Thompson CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisCreateH1"); 10452b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateH1(delegate, topo, num_comp, num_nodes, num_qpts, interp, grad, q_ref, q_weight, basis)); 1046e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 10475fe0d4faSjeremylt } 10485fe0d4faSjeremylt 10496574a04fSJeremy L Thompson CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component"); 10506574a04fSJeremy L Thompson CeedCheck(num_nodes > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node"); 10516574a04fSJeremy L Thompson CeedCheck(num_qpts > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point"); 1052227444bfSJeremy L Thompson 10532b730f8bSJeremy L Thompson CeedCall(CeedBasisGetTopologyDimension(topo, &dim)); 1054a8de75f0Sjeremylt 1055db002c03SJeremy L Thompson CeedCall(CeedCalloc(1, basis)); 1056db002c03SJeremy L Thompson CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed)); 1057d1d35e2fSjeremylt (*basis)->ref_count = 1; 10586402da51SJeremy L Thompson (*basis)->is_tensor_basis = false; 1059a8de75f0Sjeremylt (*basis)->dim = dim; 1060d99fa3c5SJeremy L Thompson (*basis)->topo = topo; 1061d1d35e2fSjeremylt (*basis)->num_comp = num_comp; 1062a8de75f0Sjeremylt (*basis)->P = P; 1063a8de75f0Sjeremylt (*basis)->Q = Q; 1064c4e3f59bSSebastian Grimberg (*basis)->fe_space = CEED_FE_SPACE_H1; 10652b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q * dim, &(*basis)->q_ref_1d)); 10662b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q, &(*basis)->q_weight_1d)); 1067ff3a0f91SJeremy L Thompson if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0])); 1068ff3a0f91SJeremy L Thompson if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0])); 10692b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q * P, &(*basis)->interp)); 10702b730f8bSJeremy L Thompson CeedCall(CeedCalloc(dim * Q * P, &(*basis)->grad)); 1071ff3a0f91SJeremy L Thompson if (interp) memcpy((*basis)->interp, interp, Q * P * sizeof(interp[0])); 1072ff3a0f91SJeremy L Thompson if (grad) memcpy((*basis)->grad, grad, dim * Q * P * sizeof(grad[0])); 10732b730f8bSJeremy L Thompson CeedCall(ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, q_ref, q_weight, *basis)); 1074e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1075a8de75f0Sjeremylt } 1076a8de75f0Sjeremylt 1077a8de75f0Sjeremylt /** 1078859c15bbSJames Wright @brief Create a non tensor-product basis for \f$H(\mathrm{div})\f$ discretizations 107950c301a5SRezgar Shakeri 1080ea61e9acSJeremy L Thompson @param[in] ceed Ceed object where the CeedBasis will be created 1081ea61e9acSJeremy 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 1082ea61e9acSJeremy L Thompson @param[in] num_comp Number of components (usually 1 for vectors in H(div) bases) 1083ea61e9acSJeremy L Thompson @param[in] num_nodes Total number of nodes (dofs per element) 1084ea61e9acSJeremy L Thompson @param[in] num_qpts Total number of quadrature points 1085c4e3f59bSSebastian Grimberg @param[in] interp Row-major (dim * num_qpts * num_nodes) matrix expressing the values of basis functions at quadrature points 1086c4e3f59bSSebastian Grimberg @param[in] div Row-major (num_qpts * num_nodes) matrix expressing divergence of basis functions at quadrature points 10879fe083eeSJeremy L Thompson @param[in] q_ref Array of length num_qpts * dim holding the locations of quadrature points on the reference element 1088ea61e9acSJeremy L Thompson @param[in] q_weight Array of length num_qpts holding the quadrature weights on the reference element 1089ea61e9acSJeremy L Thompson @param[out] basis Address of the variable where the newly created CeedBasis will be stored. 109050c301a5SRezgar Shakeri 109150c301a5SRezgar Shakeri @return An error code: 0 - success, otherwise - failure 109250c301a5SRezgar Shakeri 109350c301a5SRezgar Shakeri @ref User 109450c301a5SRezgar Shakeri **/ 10952b730f8bSJeremy L Thompson int CeedBasisCreateHdiv(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 10962b730f8bSJeremy L Thompson const CeedScalar *div, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) { 109750c301a5SRezgar Shakeri CeedInt Q = num_qpts, P = num_nodes, dim = 0; 1098c4e3f59bSSebastian Grimberg 109950c301a5SRezgar Shakeri if (!ceed->BasisCreateHdiv) { 110050c301a5SRezgar Shakeri Ceed delegate; 11016574a04fSJeremy L Thompson 11022b730f8bSJeremy L Thompson CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis")); 11036574a04fSJeremy L Thompson CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateHdiv"); 11042b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateHdiv(delegate, topo, num_comp, num_nodes, num_qpts, interp, div, q_ref, q_weight, basis)); 110550c301a5SRezgar Shakeri return CEED_ERROR_SUCCESS; 110650c301a5SRezgar Shakeri } 110750c301a5SRezgar Shakeri 11086574a04fSJeremy L Thompson CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component"); 11096574a04fSJeremy L Thompson CeedCheck(num_nodes > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node"); 11106574a04fSJeremy L Thompson CeedCheck(num_qpts > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point"); 1111227444bfSJeremy L Thompson 1112c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetTopologyDimension(topo, &dim)); 1113c4e3f59bSSebastian Grimberg 1114db002c03SJeremy L Thompson CeedCall(CeedCalloc(1, basis)); 1115db002c03SJeremy L Thompson CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed)); 111650c301a5SRezgar Shakeri (*basis)->ref_count = 1; 11176402da51SJeremy L Thompson (*basis)->is_tensor_basis = false; 111850c301a5SRezgar Shakeri (*basis)->dim = dim; 111950c301a5SRezgar Shakeri (*basis)->topo = topo; 112050c301a5SRezgar Shakeri (*basis)->num_comp = num_comp; 112150c301a5SRezgar Shakeri (*basis)->P = P; 112250c301a5SRezgar Shakeri (*basis)->Q = Q; 1123c4e3f59bSSebastian Grimberg (*basis)->fe_space = CEED_FE_SPACE_HDIV; 11242b730f8bSJeremy L Thompson CeedCall(CeedMalloc(Q * dim, &(*basis)->q_ref_1d)); 11252b730f8bSJeremy L Thompson CeedCall(CeedMalloc(Q, &(*basis)->q_weight_1d)); 112650c301a5SRezgar Shakeri if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0])); 112750c301a5SRezgar Shakeri if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0])); 11282b730f8bSJeremy L Thompson CeedCall(CeedMalloc(dim * Q * P, &(*basis)->interp)); 11292b730f8bSJeremy L Thompson CeedCall(CeedMalloc(Q * P, &(*basis)->div)); 113050c301a5SRezgar Shakeri if (interp) memcpy((*basis)->interp, interp, dim * Q * P * sizeof(interp[0])); 113150c301a5SRezgar Shakeri if (div) memcpy((*basis)->div, div, Q * P * sizeof(div[0])); 11322b730f8bSJeremy L Thompson CeedCall(ceed->BasisCreateHdiv(topo, dim, P, Q, interp, div, q_ref, q_weight, *basis)); 113350c301a5SRezgar Shakeri return CEED_ERROR_SUCCESS; 113450c301a5SRezgar Shakeri } 113550c301a5SRezgar Shakeri 113650c301a5SRezgar Shakeri /** 11374385fb7fSSebastian Grimberg @brief Create a non tensor-product basis for \f$H(\mathrm{curl})\f$ discretizations 1138c4e3f59bSSebastian Grimberg 1139c4e3f59bSSebastian Grimberg @param[in] ceed Ceed object where the CeedBasis will be created 1140c4e3f59bSSebastian Grimberg @param[in] topo Topology of element (`CEED_TOPOLOGY_QUAD`, `CEED_TOPOLOGY_PRISM`, etc.), dimension of which is used in some array sizes below 1141c4e3f59bSSebastian Grimberg @param[in] num_comp Number of components (usually 1 for vectors in H(curl) bases) 1142c4e3f59bSSebastian Grimberg @param[in] num_nodes Total number of nodes (dofs per element) 1143c4e3f59bSSebastian Grimberg @param[in] num_qpts Total number of quadrature points 1144c4e3f59bSSebastian Grimberg @param[in] interp Row-major (dim * num_qpts * num_nodes) matrix expressing the values of basis functions at quadrature points 1145c4e3f59bSSebastian 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 1146c4e3f59bSSebastian Grimberg quadrature points 1147c4e3f59bSSebastian Grimberg @param[in] q_ref Array of length num_qpts * dim holding the locations of quadrature points on the reference element 1148c4e3f59bSSebastian Grimberg @param[in] q_weight Array of length num_qpts holding the quadrature weights on the reference element 1149c4e3f59bSSebastian Grimberg @param[out] basis Address of the variable where the newly created CeedBasis will be stored. 1150c4e3f59bSSebastian Grimberg 1151c4e3f59bSSebastian Grimberg @return An error code: 0 - success, otherwise - failure 1152c4e3f59bSSebastian Grimberg 1153c4e3f59bSSebastian Grimberg @ref User 1154c4e3f59bSSebastian Grimberg **/ 1155c4e3f59bSSebastian Grimberg int CeedBasisCreateHcurl(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 1156c4e3f59bSSebastian Grimberg const CeedScalar *curl, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) { 1157c4e3f59bSSebastian Grimberg CeedInt Q = num_qpts, P = num_nodes, dim = 0, curl_comp = 0; 1158c4e3f59bSSebastian Grimberg 1159c4e3f59bSSebastian Grimberg if (!ceed->BasisCreateHdiv) { 1160c4e3f59bSSebastian Grimberg Ceed delegate; 11616574a04fSJeremy L Thompson 1162c4e3f59bSSebastian Grimberg CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis")); 11636574a04fSJeremy L Thompson CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateHcurl"); 1164c4e3f59bSSebastian Grimberg CeedCall(CeedBasisCreateHcurl(delegate, topo, num_comp, num_nodes, num_qpts, interp, curl, q_ref, q_weight, basis)); 1165c4e3f59bSSebastian Grimberg return CEED_ERROR_SUCCESS; 1166c4e3f59bSSebastian Grimberg } 1167c4e3f59bSSebastian Grimberg 11686574a04fSJeremy L Thompson CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component"); 11696574a04fSJeremy L Thompson CeedCheck(num_nodes > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node"); 11706574a04fSJeremy L Thompson CeedCheck(num_qpts > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point"); 1171c4e3f59bSSebastian Grimberg 1172c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetTopologyDimension(topo, &dim)); 1173c4e3f59bSSebastian Grimberg curl_comp = (dim < 3) ? 1 : dim; 1174c4e3f59bSSebastian Grimberg 1175db002c03SJeremy L Thompson CeedCall(CeedCalloc(1, basis)); 1176db002c03SJeremy L Thompson CeedCall(CeedReferenceCopy(ceed, &(*basis)->ceed)); 1177c4e3f59bSSebastian Grimberg (*basis)->ref_count = 1; 11786402da51SJeremy L Thompson (*basis)->is_tensor_basis = false; 1179c4e3f59bSSebastian Grimberg (*basis)->dim = dim; 1180c4e3f59bSSebastian Grimberg (*basis)->topo = topo; 1181c4e3f59bSSebastian Grimberg (*basis)->num_comp = num_comp; 1182c4e3f59bSSebastian Grimberg (*basis)->P = P; 1183c4e3f59bSSebastian Grimberg (*basis)->Q = Q; 1184c4e3f59bSSebastian Grimberg (*basis)->fe_space = CEED_FE_SPACE_HCURL; 1185c4e3f59bSSebastian Grimberg CeedCall(CeedMalloc(Q * dim, &(*basis)->q_ref_1d)); 1186c4e3f59bSSebastian Grimberg CeedCall(CeedMalloc(Q, &(*basis)->q_weight_1d)); 1187c4e3f59bSSebastian Grimberg if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0])); 1188c4e3f59bSSebastian Grimberg if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0])); 1189c4e3f59bSSebastian Grimberg CeedCall(CeedMalloc(dim * Q * P, &(*basis)->interp)); 1190c4e3f59bSSebastian Grimberg CeedCall(CeedMalloc(curl_comp * Q * P, &(*basis)->curl)); 1191c4e3f59bSSebastian Grimberg if (interp) memcpy((*basis)->interp, interp, dim * Q * P * sizeof(interp[0])); 1192c4e3f59bSSebastian Grimberg if (curl) memcpy((*basis)->curl, curl, curl_comp * Q * P * sizeof(curl[0])); 1193c4e3f59bSSebastian Grimberg CeedCall(ceed->BasisCreateHcurl(topo, dim, P, Q, interp, curl, q_ref, q_weight, *basis)); 1194c4e3f59bSSebastian Grimberg return CEED_ERROR_SUCCESS; 1195c4e3f59bSSebastian Grimberg } 1196c4e3f59bSSebastian Grimberg 1197c4e3f59bSSebastian Grimberg /** 1198ea61e9acSJeremy L Thompson @brief Create a CeedBasis for projection from the nodes of `basis_from` to the nodes of `basis_to`. 1199ba59ac12SSebastian Grimberg 12009fd66db6SSebastian Grimberg Only `CEED_EVAL_INTERP` will be valid for the new basis, `basis_project`. 12019fd66db6SSebastian Grimberg For H^1 spaces, `CEED_EVAL_GRAD` will also be valid. 1202de05fbb2SSebastian Grimberg The interpolation is given by `interp_project = interp_to^+ * interp_from`, where the pseudoinverse `interp_to^+` is given by QR 12039fd66db6SSebastian Grimberg factorization. 12049fd66db6SSebastian Grimberg The gradient (for the H^1 case) is given by `grad_project = interp_to^+ * grad_from`. 120515ad3917SSebastian Grimberg 120615ad3917SSebastian Grimberg Note: `basis_from` and `basis_to` must have compatible quadrature spaces. 120715ad3917SSebastian Grimberg 12089fd66db6SSebastian Grimberg Note: `basis_project` will have the same number of components as `basis_from`, regardless of the number of components that `basis_to` has. 12099fd66db6SSebastian Grimberg If `basis_from` has 3 components and `basis_to` has 5 components, then `basis_project` will have 3 components. 1210f113e5dcSJeremy L Thompson 1211f113e5dcSJeremy L Thompson @param[in] basis_from CeedBasis to prolong from 1212446e7af4SJeremy L Thompson @param[in] basis_to CeedBasis to prolong to 1213ea61e9acSJeremy L Thompson @param[out] basis_project Address of the variable where the newly created CeedBasis will be stored. 1214f113e5dcSJeremy L Thompson 1215f113e5dcSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 1216f113e5dcSJeremy L Thompson 1217f113e5dcSJeremy L Thompson @ref User 1218f113e5dcSJeremy L Thompson **/ 12192b730f8bSJeremy L Thompson int CeedBasisCreateProjection(CeedBasis basis_from, CeedBasis basis_to, CeedBasis *basis_project) { 1220f113e5dcSJeremy L Thompson Ceed ceed; 12212b730f8bSJeremy L Thompson CeedCall(CeedBasisGetCeed(basis_to, &ceed)); 1222f113e5dcSJeremy L Thompson 1223ecc88aebSJeremy L Thompson // Create projection matrix 122414556e63SJeremy L Thompson CeedScalar *interp_project, *grad_project; 12252b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateProjectionMatrices(basis_from, basis_to, &interp_project, &grad_project)); 1226f113e5dcSJeremy L Thompson 1227f113e5dcSJeremy L Thompson // Build basis 1228f113e5dcSJeremy L Thompson bool is_tensor; 1229f113e5dcSJeremy L Thompson CeedInt dim, num_comp; 123014556e63SJeremy L Thompson CeedScalar *q_ref, *q_weight; 12312b730f8bSJeremy L Thompson CeedCall(CeedBasisIsTensor(basis_to, &is_tensor)); 12322b730f8bSJeremy L Thompson CeedCall(CeedBasisGetDimension(basis_to, &dim)); 12332b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumComponents(basis_from, &num_comp)); 1234f113e5dcSJeremy L Thompson if (is_tensor) { 1235f113e5dcSJeremy L Thompson CeedInt P_1d_to, P_1d_from; 12362b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_1d_from)); 12372b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_1d_to)); 12382b730f8bSJeremy L Thompson CeedCall(CeedCalloc(P_1d_to, &q_ref)); 12392b730f8bSJeremy L Thompson CeedCall(CeedCalloc(P_1d_to, &q_weight)); 12402b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P_1d_from, P_1d_to, interp_project, grad_project, q_ref, q_weight, basis_project)); 1241f113e5dcSJeremy L Thompson } else { 1242de05fbb2SSebastian Grimberg // Even if basis_to and basis_from are not H1, the resulting basis is H1 for interpolation to work 1243f113e5dcSJeremy L Thompson CeedElemTopology topo; 12442b730f8bSJeremy L Thompson CeedCall(CeedBasisGetTopology(basis_to, &topo)); 1245f113e5dcSJeremy L Thompson CeedInt num_nodes_to, num_nodes_from; 12462b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis_from, &num_nodes_from)); 12472b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis_to, &num_nodes_to)); 12482b730f8bSJeremy L Thompson CeedCall(CeedCalloc(num_nodes_to * dim, &q_ref)); 12492b730f8bSJeremy L Thompson CeedCall(CeedCalloc(num_nodes_to, &q_weight)); 12502b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateH1(ceed, topo, num_comp, num_nodes_from, num_nodes_to, interp_project, grad_project, q_ref, q_weight, basis_project)); 1251f113e5dcSJeremy L Thompson } 1252f113e5dcSJeremy L Thompson 1253f113e5dcSJeremy L Thompson // Cleanup 12542b730f8bSJeremy L Thompson CeedCall(CeedFree(&interp_project)); 12552b730f8bSJeremy L Thompson CeedCall(CeedFree(&grad_project)); 12562b730f8bSJeremy L Thompson CeedCall(CeedFree(&q_ref)); 12572b730f8bSJeremy L Thompson CeedCall(CeedFree(&q_weight)); 1258f113e5dcSJeremy L Thompson 1259f113e5dcSJeremy L Thompson return CEED_ERROR_SUCCESS; 1260f113e5dcSJeremy L Thompson } 1261f113e5dcSJeremy L Thompson 1262f113e5dcSJeremy L Thompson /** 1263ea61e9acSJeremy L Thompson @brief Copy the pointer to a CeedBasis. 12649560d06aSjeremylt 1265512bb800SJeremy 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. 1266512bb800SJeremy L Thompson This CeedBasis will be destroyed if `basis_copy` is the only reference to this CeedBasis. 1267ea61e9acSJeremy L Thompson 1268ea61e9acSJeremy L Thompson @param[in] basis CeedBasis to copy reference to 1269ea61e9acSJeremy L Thompson @param[in,out] basis_copy Variable to store copied reference 12709560d06aSjeremylt 12719560d06aSjeremylt @return An error code: 0 - success, otherwise - failure 12729560d06aSjeremylt 12739560d06aSjeremylt @ref User 12749560d06aSjeremylt **/ 12759560d06aSjeremylt int CeedBasisReferenceCopy(CeedBasis basis, CeedBasis *basis_copy) { 1276393ac2cdSJeremy L Thompson if (basis != CEED_BASIS_COLLOCATED) CeedCall(CeedBasisReference(basis)); 12772b730f8bSJeremy L Thompson CeedCall(CeedBasisDestroy(basis_copy)); 12789560d06aSjeremylt *basis_copy = basis; 12799560d06aSjeremylt return CEED_ERROR_SUCCESS; 12809560d06aSjeremylt } 12819560d06aSjeremylt 12829560d06aSjeremylt /** 12837a982d89SJeremy L. Thompson @brief View a CeedBasis 12847a982d89SJeremy L. Thompson 1285ea61e9acSJeremy L Thompson @param[in] basis CeedBasis to view 1286ea61e9acSJeremy L Thompson @param[in] stream Stream to view to, e.g., stdout 12877a982d89SJeremy L. Thompson 12887a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 12897a982d89SJeremy L. Thompson 12907a982d89SJeremy L. Thompson @ref User 12917a982d89SJeremy L. Thompson **/ 12927a982d89SJeremy L. Thompson int CeedBasisView(CeedBasis basis, FILE *stream) { 129350c301a5SRezgar Shakeri CeedElemTopology topo = basis->topo; 1294c4e3f59bSSebastian Grimberg CeedFESpace fe_space = basis->fe_space; 1295c4e3f59bSSebastian Grimberg CeedInt q_comp = 0; 12962b730f8bSJeremy L Thompson 129750c301a5SRezgar Shakeri // Print FE space and element topology of the basis 1298*edf04919SJeremy L Thompson fprintf(stream, "CeedBasis in a %s on a %s element\n", CeedFESpaces[fe_space], CeedElemTopologies[topo]); 12996402da51SJeremy L Thompson if (basis->is_tensor_basis) { 1300*edf04919SJeremy L Thompson fprintf(stream, " P: %" CeedInt_FMT "\n Q: %" CeedInt_FMT "\n", basis->P_1d, basis->Q_1d); 130150c301a5SRezgar Shakeri } else { 1302*edf04919SJeremy L Thompson fprintf(stream, " P: %" CeedInt_FMT "\n Q: %" CeedInt_FMT "\n", basis->P, basis->Q); 130350c301a5SRezgar Shakeri } 1304*edf04919SJeremy L Thompson fprintf(stream, " dimension: %" CeedInt_FMT "\n field components: %" CeedInt_FMT "\n", basis->dim, basis->num_comp); 1305ea61e9acSJeremy L Thompson // Print quadrature data, interpolation/gradient/divergence/curl of the basis 13066402da51SJeremy L Thompson if (basis->is_tensor_basis) { // tensor basis 13072b730f8bSJeremy L Thompson CeedCall(CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_ref_1d, stream)); 13082b730f8bSJeremy L Thompson CeedCall(CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_weight_1d, stream)); 13092b730f8bSJeremy L Thompson CeedCall(CeedScalarView("interp1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, basis->interp_1d, stream)); 13102b730f8bSJeremy L Thompson CeedCall(CeedScalarView("grad1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, basis->grad_1d, stream)); 131150c301a5SRezgar Shakeri } else { // non-tensor basis 13122b730f8bSJeremy L Thompson CeedCall(CeedScalarView("qref", "\t% 12.8f", 1, basis->Q * basis->dim, basis->q_ref_1d, stream)); 13132b730f8bSJeremy L Thompson CeedCall(CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->q_weight_1d, stream)); 1314c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp)); 1315c4e3f59bSSebastian Grimberg CeedCall(CeedScalarView("interp", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->interp, stream)); 131650c301a5SRezgar Shakeri if (basis->grad) { 1317c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp)); 1318c4e3f59bSSebastian Grimberg CeedCall(CeedScalarView("grad", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->grad, stream)); 13197a982d89SJeremy L. Thompson } 132050c301a5SRezgar Shakeri if (basis->div) { 1321c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_DIV, &q_comp)); 1322c4e3f59bSSebastian Grimberg CeedCall(CeedScalarView("div", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->div, stream)); 1323c4e3f59bSSebastian Grimberg } 1324c4e3f59bSSebastian Grimberg if (basis->curl) { 1325c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_CURL, &q_comp)); 1326c4e3f59bSSebastian Grimberg CeedCall(CeedScalarView("curl", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->curl, stream)); 132750c301a5SRezgar Shakeri } 132850c301a5SRezgar Shakeri } 1329e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 13307a982d89SJeremy L. Thompson } 13317a982d89SJeremy L. Thompson 13327a982d89SJeremy L. Thompson /** 13337a982d89SJeremy L. Thompson @brief Apply basis evaluation from nodes to quadrature points or vice versa 13347a982d89SJeremy L. Thompson 1335ea61e9acSJeremy L Thompson @param[in] basis CeedBasis to evaluate 1336ea61e9acSJeremy L Thompson @param[in] num_elem The number of elements to apply the basis evaluation to; 1337ea61e9acSJeremy L Thompson the backend will specify the ordering in CeedElemRestrictionCreateBlocked() 1338ea61e9acSJeremy L Thompson @param[in] t_mode \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature points; 1339ea61e9acSJeremy L Thompson \ref CEED_TRANSPOSE to apply the transpose, mapping from quadrature points to nodes 1340ea61e9acSJeremy L Thompson @param[in] eval_mode \ref CEED_EVAL_NONE to use values directly, 13417a982d89SJeremy L. Thompson \ref CEED_EVAL_INTERP to use interpolated values, 13427a982d89SJeremy L. Thompson \ref CEED_EVAL_GRAD to use gradients, 1343c4e3f59bSSebastian Grimberg \ref CEED_EVAL_DIV to use divergence, 1344c4e3f59bSSebastian Grimberg \ref CEED_EVAL_CURL to use curl, 13457a982d89SJeremy L. Thompson \ref CEED_EVAL_WEIGHT to use quadrature weights. 13467a982d89SJeremy L. Thompson @param[in] u Input CeedVector 13477a982d89SJeremy L. Thompson @param[out] v Output CeedVector 13487a982d89SJeremy L. Thompson 13497a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 13507a982d89SJeremy L. Thompson 13517a982d89SJeremy L. Thompson @ref User 13527a982d89SJeremy L. Thompson **/ 13532b730f8bSJeremy L Thompson int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, CeedVector v) { 13541f9221feSJeremy L Thompson CeedSize u_length = 0, v_length; 1355c4e3f59bSSebastian Grimberg CeedInt dim, num_comp, q_comp, num_nodes, num_qpts; 13562b730f8bSJeremy L Thompson CeedCall(CeedBasisGetDimension(basis, &dim)); 13572b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumComponents(basis, &num_comp)); 1358c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp)); 13592b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis, &num_nodes)); 13602b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts)); 13612b730f8bSJeremy L Thompson CeedCall(CeedVectorGetLength(v, &v_length)); 1362c8c3fa7dSJeremy L Thompson if (u) CeedCall(CeedVectorGetLength(u, &u_length)); 13637a982d89SJeremy L. Thompson 13646574a04fSJeremy L Thompson CeedCheck(basis->Apply, basis->ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisApply"); 1365e15f9bd0SJeremy L Thompson 1366e15f9bd0SJeremy L Thompson // Check compatibility of topological and geometrical dimensions 13676574a04fSJeremy L Thompson CeedCheck((t_mode == CEED_TRANSPOSE && v_length % num_nodes == 0 && u_length % num_qpts == 0) || 13686574a04fSJeremy L Thompson (t_mode == CEED_NOTRANSPOSE && u_length % num_nodes == 0 && v_length % num_qpts == 0), 13696574a04fSJeremy L Thompson basis->ceed, CEED_ERROR_DIMENSION, "Length of input/output vectors incompatible with basis dimensions"); 13707a982d89SJeremy L. Thompson 1371e15f9bd0SJeremy L Thompson // Check vector lengths to prevent out of bounds issues 13726574a04fSJeremy L Thompson bool good_dims = true; 1373d1d35e2fSjeremylt switch (eval_mode) { 1374e15f9bd0SJeremy L Thompson case CEED_EVAL_NONE: 13752b730f8bSJeremy L Thompson case CEED_EVAL_INTERP: 13762b730f8bSJeremy L Thompson case CEED_EVAL_GRAD: 1377c4e3f59bSSebastian Grimberg case CEED_EVAL_DIV: 1378c4e3f59bSSebastian Grimberg case CEED_EVAL_CURL: 13796574a04fSJeremy L Thompson good_dims = 13806574a04fSJeremy L Thompson ((t_mode == CEED_TRANSPOSE && u_length >= num_elem * num_comp * num_qpts * q_comp && v_length >= num_elem * num_comp * num_nodes) || 13816574a04fSJeremy L Thompson (t_mode == CEED_NOTRANSPOSE && v_length >= num_elem * num_qpts * num_comp * q_comp && u_length >= num_elem * num_comp * num_nodes)); 1382e15f9bd0SJeremy L Thompson break; 1383e15f9bd0SJeremy L Thompson case CEED_EVAL_WEIGHT: 13846574a04fSJeremy L Thompson good_dims = v_length >= num_elem * num_qpts; 1385e15f9bd0SJeremy L Thompson break; 1386e15f9bd0SJeremy L Thompson } 13876574a04fSJeremy L Thompson CeedCheck(good_dims, basis->ceed, CEED_ERROR_DIMENSION, "Input/output vectors too short for basis and evaluation mode"); 1388e15f9bd0SJeremy L Thompson 13892b730f8bSJeremy L Thompson CeedCall(basis->Apply(basis, num_elem, t_mode, eval_mode, u, v)); 1390e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 13917a982d89SJeremy L. Thompson } 13927a982d89SJeremy L. Thompson 13937a982d89SJeremy L. Thompson /** 1394c8c3fa7dSJeremy L Thompson @brief Apply basis evaluation from nodes to arbitrary points 1395c8c3fa7dSJeremy L Thompson 1396c8c3fa7dSJeremy L Thompson @param[in] basis CeedBasis to evaluate 1397c8c3fa7dSJeremy L Thompson @param[in] num_points The number of points to apply the basis evaluation to 1398c8c3fa7dSJeremy L Thompson @param[in] t_mode \ref CEED_NOTRANSPOSE to evaluate from nodes to points; 1399c8c3fa7dSJeremy L Thompson \ref CEED_TRANSPOSE to apply the transpose, mapping from points to nodes 1400c8c3fa7dSJeremy L Thompson @param[in] eval_mode \ref CEED_EVAL_INTERP to use interpolated values, 1401c8c3fa7dSJeremy L Thompson \ref CEED_EVAL_GRAD to use gradients 1402c8c3fa7dSJeremy L Thompson @param[in] x_ref CeedVector holding reference coordinates of each point 1403c8c3fa7dSJeremy L Thompson @param[in] u Input CeedVector, of length `num_nodes * num_comp` for `CEED_NOTRANSPOSE` 1404c8c3fa7dSJeremy L Thompson @param[out] v Output CeedVector, of length `num_points * num_q_comp` for `CEED_NOTRANSPOSE` with `CEED_EVAL_INTERP` 1405c8c3fa7dSJeremy L Thompson 1406c8c3fa7dSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 1407c8c3fa7dSJeremy L Thompson 1408c8c3fa7dSJeremy L Thompson @ref User 1409c8c3fa7dSJeremy L Thompson **/ 1410c8c3fa7dSJeremy L Thompson int CeedBasisApplyAtPoints(CeedBasis basis, CeedInt num_points, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, 1411c8c3fa7dSJeremy L Thompson CeedVector v) { 1412c8c3fa7dSJeremy L Thompson CeedSize x_length = 0, u_length = 0, v_length; 1413c8c3fa7dSJeremy L Thompson CeedInt dim, num_comp, num_q_comp, num_nodes, P_1d = 1, Q_1d = 1; 1414c8c3fa7dSJeremy L Thompson 1415c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisGetDimension(basis, &dim)); 1416c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d)); 1417c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 1418c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisGetNumComponents(basis, &num_comp)); 1419c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &num_q_comp)); 1420c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis, &num_nodes)); 1421c8c3fa7dSJeremy L Thompson CeedCall(CeedVectorGetLength(x_ref, &x_length)); 1422c8c3fa7dSJeremy L Thompson CeedCall(CeedVectorGetLength(v, &v_length)); 1423c8c3fa7dSJeremy L Thompson CeedCall(CeedVectorGetLength(u, &u_length)); 1424c8c3fa7dSJeremy L Thompson 1425c8c3fa7dSJeremy L Thompson // Check compatibility of topological and geometrical dimensions 1426c8c3fa7dSJeremy L Thompson CeedCheck((t_mode == CEED_TRANSPOSE && v_length % num_nodes == 0) || (t_mode == CEED_NOTRANSPOSE && u_length % num_nodes == 0), basis->ceed, 1427c8c3fa7dSJeremy L Thompson CEED_ERROR_DIMENSION, "Length of input/output vectors incompatible with basis dimensions and number of points"); 1428c8c3fa7dSJeremy L Thompson 1429c8c3fa7dSJeremy L Thompson // Check compatibility coordinates vector 1430c8c3fa7dSJeremy L Thompson CeedCheck(x_length >= num_points * dim, basis->ceed, CEED_ERROR_DIMENSION, 1431c8c3fa7dSJeremy L Thompson "Length of reference coordinate vector incompatible with basis dimension and number of points"); 1432c8c3fa7dSJeremy L Thompson 1433c8c3fa7dSJeremy L Thompson // Check vector lengths to prevent out of bounds issues 1434c8c3fa7dSJeremy L Thompson bool good_dims = false; 1435c8c3fa7dSJeremy L Thompson switch (eval_mode) { 1436c8c3fa7dSJeremy L Thompson case CEED_EVAL_INTERP: 1437c8c3fa7dSJeremy L Thompson good_dims = ((t_mode == CEED_TRANSPOSE && (u_length >= num_points * num_q_comp || v_length >= num_nodes * num_comp)) || 1438c8c3fa7dSJeremy L Thompson (t_mode == CEED_NOTRANSPOSE && (v_length >= num_points * num_q_comp || u_length >= num_nodes * num_comp))); 1439c8c3fa7dSJeremy L Thompson break; 1440c8c3fa7dSJeremy L Thompson case CEED_EVAL_GRAD: 1441c8c3fa7dSJeremy L Thompson case CEED_EVAL_NONE: 1442c8c3fa7dSJeremy L Thompson case CEED_EVAL_WEIGHT: 1443c8c3fa7dSJeremy L Thompson case CEED_EVAL_DIV: 1444c8c3fa7dSJeremy L Thompson case CEED_EVAL_CURL: 1445c8c3fa7dSJeremy L Thompson // LCOV_EXCL_START 1446c8c3fa7dSJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_UNSUPPORTED, "Evaluation at arbitrary points not supported for %s", CeedEvalModes[eval_mode]); 1447c8c3fa7dSJeremy L Thompson // LCOV_EXCL_STOP 1448c8c3fa7dSJeremy L Thompson } 1449c8c3fa7dSJeremy L Thompson CeedCheck(good_dims, basis->ceed, CEED_ERROR_DIMENSION, "Input/output vectors too short for basis and evaluation mode"); 1450c8c3fa7dSJeremy L Thompson 1451c8c3fa7dSJeremy L Thompson // Backend method 1452c8c3fa7dSJeremy L Thompson if (basis->ApplyAtPoints) { 1453c8c3fa7dSJeremy L Thompson CeedCall(basis->ApplyAtPoints(basis, num_points, t_mode, eval_mode, x_ref, u, v)); 1454c8c3fa7dSJeremy L Thompson return CEED_ERROR_SUCCESS; 1455c8c3fa7dSJeremy L Thompson } 1456c8c3fa7dSJeremy L Thompson 1457c8c3fa7dSJeremy L Thompson // Default implementation 1458c8c3fa7dSJeremy L Thompson CeedCheck(basis->is_tensor_basis, basis->ceed, CEED_ERROR_UNSUPPORTED, "Evaluation at arbitrary points only supported for tensor product bases"); 1459c8c3fa7dSJeremy L Thompson if (!basis->basis_chebyshev) { 1460c8c3fa7dSJeremy L Thompson // Build matrix mapping from quadrature point values to Chebyshev coefficients 1461c8c3fa7dSJeremy L Thompson CeedScalar *tau, *C, *I, *chebyshev_coeffs_1d; 1462c8c3fa7dSJeremy L Thompson const CeedScalar *q_ref_1d; 1463c8c3fa7dSJeremy L Thompson 1464c8c3fa7dSJeremy L Thompson // Build coefficient matrix 1465c8c3fa7dSJeremy L Thompson // -- Note: Clang-tidy needs this check because it does not understand the is_tensor_basis check above 1466c8c3fa7dSJeremy L Thompson CeedCheck(P_1d > 0 && Q_1d > 0, basis->ceed, CEED_ERROR_INCOMPATIBLE, "Basis dimensions are malformed"); 1467c8c3fa7dSJeremy L Thompson CeedCall(CeedCalloc(Q_1d * Q_1d, &C)); 1468c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisGetQRef(basis, &q_ref_1d)); 1469c8c3fa7dSJeremy L Thompson for (CeedInt i = 0; i < Q_1d; i++) { 1470c8c3fa7dSJeremy L Thompson const CeedScalar x = q_ref_1d[i]; 1471c8c3fa7dSJeremy L Thompson 1472c8c3fa7dSJeremy L Thompson C[i * Q_1d + 0] = 1.0; 1473c8c3fa7dSJeremy L Thompson C[i * Q_1d + 1] = 2 * x; 1474c8c3fa7dSJeremy 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]; 1475c8c3fa7dSJeremy L Thompson } 1476c8c3fa7dSJeremy L Thompson 1477c8c3fa7dSJeremy L Thompson // Inverse of coefficient matrix 1478c8c3fa7dSJeremy L Thompson CeedCall(CeedCalloc(Q_1d * Q_1d, &chebyshev_coeffs_1d)); 1479c8c3fa7dSJeremy L Thompson CeedCall(CeedCalloc(Q_1d * Q_1d, &I)); 1480c8c3fa7dSJeremy L Thompson CeedCall(CeedCalloc(Q_1d, &tau)); 1481c8c3fa7dSJeremy L Thompson // -- QR Factorization, C = Q R 1482c8c3fa7dSJeremy L Thompson CeedCall(CeedQRFactorization(basis->ceed, C, tau, Q_1d, Q_1d)); 1483c8c3fa7dSJeremy L Thompson // -- chebyshev_coeffs_1d = R_inv Q^T 1484c8c3fa7dSJeremy L Thompson for (CeedInt i = 0; i < Q_1d; i++) I[i * Q_1d + i] = 1.0; 1485c8c3fa7dSJeremy L Thompson // ---- Apply R_inv, chebyshev_coeffs_1d = I R_inv 1486c8c3fa7dSJeremy L Thompson for (CeedInt i = 0; i < Q_1d; i++) { // Row i 1487c8c3fa7dSJeremy L Thompson chebyshev_coeffs_1d[Q_1d * i] = I[Q_1d * i] / C[0]; 1488c8c3fa7dSJeremy L Thompson for (CeedInt j = 1; j < Q_1d; j++) { // Column j 1489c8c3fa7dSJeremy L Thompson chebyshev_coeffs_1d[j + Q_1d * i] = I[j + Q_1d * i]; 1490c8c3fa7dSJeremy 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]; 1491c8c3fa7dSJeremy L Thompson chebyshev_coeffs_1d[j + Q_1d * i] /= C[j + Q_1d * j]; 1492c8c3fa7dSJeremy L Thompson } 1493c8c3fa7dSJeremy L Thompson } 1494c8c3fa7dSJeremy L Thompson // ---- Apply Q^T, chebyshev_coeffs_1d = R_inv Q^T 1495c8c3fa7dSJeremy L Thompson CeedCall(CeedHouseholderApplyQ(chebyshev_coeffs_1d, C, tau, CEED_NOTRANSPOSE, Q_1d, Q_1d, Q_1d, 1, Q_1d)); 1496c8c3fa7dSJeremy L Thompson 1497c8c3fa7dSJeremy L Thompson // Build basis mapping from nodes to Chebyshev coefficients 1498c8c3fa7dSJeremy L Thompson CeedScalar *chebyshev_interp_1d, *chebyshev_grad_1d, *chebyshev_q_weight_1d; 1499c8c3fa7dSJeremy L Thompson const CeedScalar *interp_1d; 1500c8c3fa7dSJeremy L Thompson 1501c8c3fa7dSJeremy L Thompson CeedCall(CeedCalloc(Q_1d * Q_1d, &chebyshev_interp_1d)); 1502c8c3fa7dSJeremy L Thompson CeedCall(CeedCalloc(Q_1d * Q_1d, &chebyshev_grad_1d)); 1503c8c3fa7dSJeremy L Thompson CeedCall(CeedCalloc(Q_1d, &chebyshev_q_weight_1d)); 1504c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisGetInterp1D(basis, &interp_1d)); 1505c8c3fa7dSJeremy L Thompson CeedCall(CeedMatrixMatrixMultiply(basis->ceed, chebyshev_coeffs_1d, interp_1d, chebyshev_interp_1d, Q_1d, P_1d, Q_1d)); 1506c8c3fa7dSJeremy L Thompson 1507c8c3fa7dSJeremy L Thompson CeedCall(CeedVectorCreate(basis->ceed, num_comp * CeedIntPow(Q_1d, dim), &basis->vec_chebyshev)); 1508c8c3fa7dSJeremy 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, 1509c8c3fa7dSJeremy L Thompson &basis->basis_chebyshev)); 1510c8c3fa7dSJeremy L Thompson 1511c8c3fa7dSJeremy L Thompson // Cleanup 1512c8c3fa7dSJeremy L Thompson CeedCall(CeedFree(&C)); 1513c8c3fa7dSJeremy L Thompson CeedCall(CeedFree(&chebyshev_coeffs_1d)); 1514c8c3fa7dSJeremy L Thompson CeedCall(CeedFree(&I)); 1515c8c3fa7dSJeremy L Thompson CeedCall(CeedFree(&tau)); 1516c8c3fa7dSJeremy L Thompson CeedCall(CeedFree(&chebyshev_interp_1d)); 1517c8c3fa7dSJeremy L Thompson CeedCall(CeedFree(&chebyshev_grad_1d)); 1518c8c3fa7dSJeremy L Thompson CeedCall(CeedFree(&chebyshev_q_weight_1d)); 1519c8c3fa7dSJeremy L Thompson } 1520c8c3fa7dSJeremy L Thompson 1521c8c3fa7dSJeremy L Thompson // Create TensorContract object if needed, such as a basis from the GPU backends 1522c8c3fa7dSJeremy L Thompson if (!basis->contract) { 1523c8c3fa7dSJeremy L Thompson Ceed ceed_ref; 1524c8c3fa7dSJeremy L Thompson CeedBasis basis_ref; 1525c8c3fa7dSJeremy L Thompson 1526c8c3fa7dSJeremy L Thompson CeedCall(CeedInit("/cpu/self", &ceed_ref)); 1527c8c3fa7dSJeremy L Thompson // Only need matching tensor contraction dimensions, any type of basis will work 1528c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisCreateTensorH1Lagrange(ceed_ref, dim, num_comp, Q_1d, Q_1d, CEED_GAUSS, &basis_ref)); 1529c8c3fa7dSJeremy L Thompson CeedCall(CeedTensorContractReference(basis_ref->contract)); 1530c8c3fa7dSJeremy L Thompson basis->contract = basis_ref->contract; 1531c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisDestroy(&basis_ref)); 1532c8c3fa7dSJeremy L Thompson CeedCall(CeedDestroy(&ceed_ref)); 1533c8c3fa7dSJeremy L Thompson } 1534c8c3fa7dSJeremy L Thompson 1535c8c3fa7dSJeremy L Thompson // Basis evaluation 1536c8c3fa7dSJeremy L Thompson switch (t_mode) { 1537c8c3fa7dSJeremy L Thompson case CEED_NOTRANSPOSE: { 1538c8c3fa7dSJeremy L Thompson // Nodes to arbitrary points 1539c8c3fa7dSJeremy L Thompson CeedScalar *v_array; 1540c8c3fa7dSJeremy L Thompson const CeedScalar *chebyshev_coeffs, *x_array_read; 1541c8c3fa7dSJeremy L Thompson 1542c8c3fa7dSJeremy L Thompson // -- Interpolate to Chebyshev coefficients 1543c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisApply(basis->basis_chebyshev, 1, CEED_NOTRANSPOSE, CEED_EVAL_INTERP, u, basis->vec_chebyshev)); 1544c8c3fa7dSJeremy L Thompson 1545c8c3fa7dSJeremy L Thompson // -- Evaluate Chebyshev polynomials at arbitrary points 1546c8c3fa7dSJeremy L Thompson CeedCall(CeedVectorGetArrayRead(basis->vec_chebyshev, CEED_MEM_HOST, &chebyshev_coeffs)); 1547c8c3fa7dSJeremy L Thompson CeedCall(CeedVectorGetArrayRead(x_ref, CEED_MEM_HOST, &x_array_read)); 1548c8c3fa7dSJeremy L Thompson CeedCall(CeedVectorGetArrayWrite(v, CEED_MEM_HOST, &v_array)); 1549c8c3fa7dSJeremy L Thompson { 1550c8c3fa7dSJeremy L Thompson CeedScalar tmp[2][num_comp * CeedIntPow(Q_1d, dim)], chebyshev_x[Q_1d]; 1551c8c3fa7dSJeremy L Thompson 1552c8c3fa7dSJeremy L Thompson // ---- Values at point 1553c8c3fa7dSJeremy L Thompson for (CeedInt p = 0; p < num_points; p++) { 1554c8c3fa7dSJeremy L Thompson CeedInt pre = num_comp * CeedIntPow(Q_1d, dim - 1), post = 1; 1555c8c3fa7dSJeremy L Thompson 1556c8c3fa7dSJeremy L Thompson for (CeedInt d = dim - 1; d >= 0; d--) { 1557c8c3fa7dSJeremy L Thompson // ------ Compute Chebyshev polynomial values 1558c8c3fa7dSJeremy L Thompson { 1559c8c3fa7dSJeremy L Thompson const CeedScalar x = x_array_read[p * dim + d]; 1560c8c3fa7dSJeremy L Thompson 1561c8c3fa7dSJeremy L Thompson chebyshev_x[0] = 1.0; 1562c8c3fa7dSJeremy L Thompson chebyshev_x[1] = 2 * x; 1563c8c3fa7dSJeremy L Thompson for (CeedInt j = 2; j < Q_1d; j++) chebyshev_x[j] = 2 * x * chebyshev_x[j - 1] - chebyshev_x[j - 2]; 1564c8c3fa7dSJeremy L Thompson } 1565c8c3fa7dSJeremy L Thompson // ------ Tensor contract 1566c8c3fa7dSJeremy L Thompson CeedCall(CeedTensorContractApply(basis->contract, pre, Q_1d, post, 1, chebyshev_x, t_mode, false, 1567c8c3fa7dSJeremy L Thompson d == (dim - 1) ? chebyshev_coeffs : tmp[d % 2], d == 0 ? &v_array[p * num_comp] : tmp[(d + 1) % 2])); 1568c8c3fa7dSJeremy L Thompson pre /= Q_1d; 1569c8c3fa7dSJeremy L Thompson post *= 1; 1570c8c3fa7dSJeremy L Thompson } 1571c8c3fa7dSJeremy L Thompson } 1572c8c3fa7dSJeremy L Thompson } 1573c8c3fa7dSJeremy L Thompson CeedCall(CeedVectorRestoreArrayRead(basis->vec_chebyshev, &chebyshev_coeffs)); 1574c8c3fa7dSJeremy L Thompson CeedCall(CeedVectorRestoreArrayRead(x_ref, &x_array_read)); 1575c8c3fa7dSJeremy L Thompson CeedCall(CeedVectorRestoreArray(v, &v_array)); 1576c8c3fa7dSJeremy L Thompson break; 1577c8c3fa7dSJeremy L Thompson } 15782a94f45fSJeremy L Thompson case CEED_TRANSPOSE: { 15792a94f45fSJeremy L Thompson // Arbitrary points to nodes 15802a94f45fSJeremy L Thompson CeedScalar *chebyshev_coeffs; 15812a94f45fSJeremy L Thompson const CeedScalar *u_array, *x_array_read; 15822a94f45fSJeremy L Thompson 15832a94f45fSJeremy L Thompson // -- Transpose of evaluaton of Chebyshev polynomials at arbitrary points 15842a94f45fSJeremy L Thompson CeedCall(CeedVectorGetArrayWrite(basis->vec_chebyshev, CEED_MEM_HOST, &chebyshev_coeffs)); 15852a94f45fSJeremy L Thompson CeedCall(CeedVectorGetArrayRead(x_ref, CEED_MEM_HOST, &x_array_read)); 15862a94f45fSJeremy L Thompson CeedCall(CeedVectorGetArrayRead(u, CEED_MEM_HOST, &u_array)); 15872a94f45fSJeremy L Thompson { 15882a94f45fSJeremy L Thompson CeedScalar tmp[2][num_comp * CeedIntPow(Q_1d, dim)], chebyshev_x[Q_1d]; 15892a94f45fSJeremy L Thompson 15902a94f45fSJeremy L Thompson // ---- Values at point 15912a94f45fSJeremy L Thompson for (CeedInt p = 0; p < num_points; p++) { 15922a94f45fSJeremy L Thompson CeedInt pre = num_comp * 1, post = 1; 15932a94f45fSJeremy L Thompson 15942a94f45fSJeremy L Thompson for (CeedInt d = dim - 1; d >= 0; d--) { 15952a94f45fSJeremy L Thompson // ------ Compute Chebyshev polynomial values 15962a94f45fSJeremy L Thompson { 15972a94f45fSJeremy L Thompson const CeedScalar x = x_array_read[p * dim + d]; 15982a94f45fSJeremy L Thompson 15992a94f45fSJeremy L Thompson chebyshev_x[0] = 1.0; 16002a94f45fSJeremy L Thompson chebyshev_x[1] = 2 * x; 16012a94f45fSJeremy L Thompson for (CeedInt j = 2; j < Q_1d; j++) chebyshev_x[j] = 2 * x * chebyshev_x[j - 1] - chebyshev_x[j - 2]; 16022a94f45fSJeremy L Thompson } 16032a94f45fSJeremy L Thompson // ------ Tensor contract 16042a94f45fSJeremy L Thompson CeedCall(CeedTensorContractApply(basis->contract, pre, 1, post, Q_1d, chebyshev_x, t_mode, p > 0 && d == 0, 16052a94f45fSJeremy L Thompson d == (dim - 1) ? &u_array[p * num_comp] : tmp[d % 2], d == 0 ? chebyshev_coeffs : tmp[(d + 1) % 2])); 16062a94f45fSJeremy L Thompson pre /= 1; 16072a94f45fSJeremy L Thompson post *= Q_1d; 16082a94f45fSJeremy L Thompson } 16092a94f45fSJeremy L Thompson } 16102a94f45fSJeremy L Thompson } 16112a94f45fSJeremy L Thompson CeedCall(CeedVectorRestoreArray(basis->vec_chebyshev, &chebyshev_coeffs)); 16122a94f45fSJeremy L Thompson CeedCall(CeedVectorRestoreArrayRead(x_ref, &x_array_read)); 16132a94f45fSJeremy L Thompson CeedCall(CeedVectorRestoreArrayRead(u, &u_array)); 16142a94f45fSJeremy L Thompson 16152a94f45fSJeremy L Thompson // -- Interpolate transpose from Chebyshev coefficients 16162a94f45fSJeremy L Thompson CeedCall(CeedBasisApply(basis->basis_chebyshev, 1, CEED_TRANSPOSE, CEED_EVAL_INTERP, basis->vec_chebyshev, v)); 16172a94f45fSJeremy L Thompson break; 16182a94f45fSJeremy L Thompson } 1619c8c3fa7dSJeremy L Thompson } 1620c8c3fa7dSJeremy L Thompson 1621c8c3fa7dSJeremy L Thompson return CEED_ERROR_SUCCESS; 1622c8c3fa7dSJeremy L Thompson } 1623c8c3fa7dSJeremy L Thompson 1624c8c3fa7dSJeremy L Thompson /** 1625b7c9bbdaSJeremy L Thompson @brief Get Ceed associated with a CeedBasis 1626b7c9bbdaSJeremy L Thompson 1627ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 1628b7c9bbdaSJeremy L Thompson @param[out] ceed Variable to store Ceed 1629b7c9bbdaSJeremy L Thompson 1630b7c9bbdaSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 1631b7c9bbdaSJeremy L Thompson 1632b7c9bbdaSJeremy L Thompson @ref Advanced 1633b7c9bbdaSJeremy L Thompson **/ 1634b7c9bbdaSJeremy L Thompson int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) { 1635b7c9bbdaSJeremy L Thompson *ceed = basis->ceed; 1636b7c9bbdaSJeremy L Thompson return CEED_ERROR_SUCCESS; 1637b7c9bbdaSJeremy L Thompson } 1638b7c9bbdaSJeremy L Thompson 1639b7c9bbdaSJeremy L Thompson /** 16409d007619Sjeremylt @brief Get dimension for given CeedBasis 16419d007619Sjeremylt 1642ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 16439d007619Sjeremylt @param[out] dim Variable to store dimension of basis 16449d007619Sjeremylt 16459d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 16469d007619Sjeremylt 1647b7c9bbdaSJeremy L Thompson @ref Advanced 16489d007619Sjeremylt **/ 16499d007619Sjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) { 16509d007619Sjeremylt *dim = basis->dim; 1651e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 16529d007619Sjeremylt } 16539d007619Sjeremylt 16549d007619Sjeremylt /** 1655d99fa3c5SJeremy L Thompson @brief Get topology for given CeedBasis 1656d99fa3c5SJeremy L Thompson 1657ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 1658d99fa3c5SJeremy L Thompson @param[out] topo Variable to store topology of basis 1659d99fa3c5SJeremy L Thompson 1660d99fa3c5SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 1661d99fa3c5SJeremy L Thompson 1662b7c9bbdaSJeremy L Thompson @ref Advanced 1663d99fa3c5SJeremy L Thompson **/ 1664d99fa3c5SJeremy L Thompson int CeedBasisGetTopology(CeedBasis basis, CeedElemTopology *topo) { 1665d99fa3c5SJeremy L Thompson *topo = basis->topo; 1666e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1667d99fa3c5SJeremy L Thompson } 1668d99fa3c5SJeremy L Thompson 1669d99fa3c5SJeremy L Thompson /** 16709d007619Sjeremylt @brief Get number of components for given CeedBasis 16719d007619Sjeremylt 1672ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 1673d1d35e2fSjeremylt @param[out] num_comp Variable to store number of components of basis 16749d007619Sjeremylt 16759d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 16769d007619Sjeremylt 1677b7c9bbdaSJeremy L Thompson @ref Advanced 16789d007619Sjeremylt **/ 1679d1d35e2fSjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *num_comp) { 1680d1d35e2fSjeremylt *num_comp = basis->num_comp; 1681e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 16829d007619Sjeremylt } 16839d007619Sjeremylt 16849d007619Sjeremylt /** 16859d007619Sjeremylt @brief Get total number of nodes (in dim dimensions) of a CeedBasis 16869d007619Sjeremylt 1687ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 16889d007619Sjeremylt @param[out] P Variable to store number of nodes 16899d007619Sjeremylt 16909d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 16919d007619Sjeremylt 16929d007619Sjeremylt @ref Utility 16939d007619Sjeremylt **/ 16949d007619Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) { 16959d007619Sjeremylt *P = basis->P; 1696e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 16979d007619Sjeremylt } 16989d007619Sjeremylt 16999d007619Sjeremylt /** 17009d007619Sjeremylt @brief Get total number of nodes (in 1 dimension) of a CeedBasis 17019d007619Sjeremylt 1702ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 1703d1d35e2fSjeremylt @param[out] P_1d Variable to store number of nodes 17049d007619Sjeremylt 17059d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 17069d007619Sjeremylt 1707b7c9bbdaSJeremy L Thompson @ref Advanced 17089d007619Sjeremylt **/ 1709d1d35e2fSjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P_1d) { 17106402da51SJeremy L Thompson CeedCheck(basis->is_tensor_basis, basis->ceed, CEED_ERROR_MINOR, "Cannot supply P_1d for non-tensor basis"); 1711d1d35e2fSjeremylt *P_1d = basis->P_1d; 1712e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 17139d007619Sjeremylt } 17149d007619Sjeremylt 17159d007619Sjeremylt /** 17169d007619Sjeremylt @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis 17179d007619Sjeremylt 1718ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 17199d007619Sjeremylt @param[out] Q Variable to store number of quadrature points 17209d007619Sjeremylt 17219d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 17229d007619Sjeremylt 17239d007619Sjeremylt @ref Utility 17249d007619Sjeremylt **/ 17259d007619Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) { 17269d007619Sjeremylt *Q = basis->Q; 1727e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 17289d007619Sjeremylt } 17299d007619Sjeremylt 17309d007619Sjeremylt /** 17319d007619Sjeremylt @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis 17329d007619Sjeremylt 1733ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 1734d1d35e2fSjeremylt @param[out] Q_1d Variable to store number of quadrature points 17359d007619Sjeremylt 17369d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 17379d007619Sjeremylt 1738b7c9bbdaSJeremy L Thompson @ref Advanced 17399d007619Sjeremylt **/ 1740d1d35e2fSjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q_1d) { 17416402da51SJeremy L Thompson CeedCheck(basis->is_tensor_basis, basis->ceed, CEED_ERROR_MINOR, "Cannot supply Q_1d for non-tensor basis"); 1742d1d35e2fSjeremylt *Q_1d = basis->Q_1d; 1743e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 17449d007619Sjeremylt } 17459d007619Sjeremylt 17469d007619Sjeremylt /** 1747ea61e9acSJeremy L Thompson @brief Get reference coordinates of quadrature points (in dim dimensions) of a CeedBasis 17489d007619Sjeremylt 1749ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 1750d1d35e2fSjeremylt @param[out] q_ref Variable to store reference coordinates of quadrature points 17519d007619Sjeremylt 17529d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 17539d007619Sjeremylt 1754b7c9bbdaSJeremy L Thompson @ref Advanced 17559d007619Sjeremylt **/ 1756d1d35e2fSjeremylt int CeedBasisGetQRef(CeedBasis basis, const CeedScalar **q_ref) { 1757d1d35e2fSjeremylt *q_ref = basis->q_ref_1d; 1758e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 17599d007619Sjeremylt } 17609d007619Sjeremylt 17619d007619Sjeremylt /** 1762ea61e9acSJeremy L Thompson @brief Get quadrature weights of quadrature points (in dim dimensions) of a CeedBasis 17639d007619Sjeremylt 1764ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 1765d1d35e2fSjeremylt @param[out] q_weight Variable to store quadrature weights 17669d007619Sjeremylt 17679d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 17689d007619Sjeremylt 1769b7c9bbdaSJeremy L Thompson @ref Advanced 17709d007619Sjeremylt **/ 1771d1d35e2fSjeremylt int CeedBasisGetQWeights(CeedBasis basis, const CeedScalar **q_weight) { 1772d1d35e2fSjeremylt *q_weight = basis->q_weight_1d; 1773e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 17749d007619Sjeremylt } 17759d007619Sjeremylt 17769d007619Sjeremylt /** 17779d007619Sjeremylt @brief Get interpolation matrix of a CeedBasis 17789d007619Sjeremylt 1779ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 17809d007619Sjeremylt @param[out] interp Variable to store interpolation matrix 17819d007619Sjeremylt 17829d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 17839d007619Sjeremylt 1784b7c9bbdaSJeremy L Thompson @ref Advanced 17859d007619Sjeremylt **/ 17866c58de82SJeremy L Thompson int CeedBasisGetInterp(CeedBasis basis, const CeedScalar **interp) { 17876402da51SJeremy L Thompson if (!basis->interp && basis->is_tensor_basis) { 17889d007619Sjeremylt // Allocate 17892b730f8bSJeremy L Thompson CeedCall(CeedMalloc(basis->Q * basis->P, &basis->interp)); 17909d007619Sjeremylt 17919d007619Sjeremylt // Initialize 17922b730f8bSJeremy L Thompson for (CeedInt i = 0; i < basis->Q * basis->P; i++) basis->interp[i] = 1.0; 17939d007619Sjeremylt 17949d007619Sjeremylt // Calculate 17952b730f8bSJeremy L Thompson for (CeedInt d = 0; d < basis->dim; d++) { 17962b730f8bSJeremy L Thompson for (CeedInt qpt = 0; qpt < basis->Q; qpt++) { 17979d007619Sjeremylt for (CeedInt node = 0; node < basis->P; node++) { 1798d1d35e2fSjeremylt CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d; 1799d1d35e2fSjeremylt CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d; 1800d1d35e2fSjeremylt basis->interp[qpt * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p]; 18019d007619Sjeremylt } 18029d007619Sjeremylt } 18032b730f8bSJeremy L Thompson } 18042b730f8bSJeremy L Thompson } 18059d007619Sjeremylt *interp = basis->interp; 1806e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 18079d007619Sjeremylt } 18089d007619Sjeremylt 18099d007619Sjeremylt /** 18109d007619Sjeremylt @brief Get 1D interpolation matrix of a tensor product CeedBasis 18119d007619Sjeremylt 1812ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 1813d1d35e2fSjeremylt @param[out] interp_1d Variable to store interpolation matrix 18149d007619Sjeremylt 18159d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 18169d007619Sjeremylt 18179d007619Sjeremylt @ref Backend 18189d007619Sjeremylt **/ 1819d1d35e2fSjeremylt int CeedBasisGetInterp1D(CeedBasis basis, const CeedScalar **interp_1d) { 18206402da51SJeremy L Thompson CeedCheck(basis->is_tensor_basis, basis->ceed, CEED_ERROR_MINOR, "CeedBasis is not a tensor product basis."); 1821d1d35e2fSjeremylt *interp_1d = basis->interp_1d; 1822e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 18239d007619Sjeremylt } 18249d007619Sjeremylt 18259d007619Sjeremylt /** 18269d007619Sjeremylt @brief Get gradient matrix of a CeedBasis 18279d007619Sjeremylt 1828ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 18299d007619Sjeremylt @param[out] grad Variable to store gradient matrix 18309d007619Sjeremylt 18319d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 18329d007619Sjeremylt 1833b7c9bbdaSJeremy L Thompson @ref Advanced 18349d007619Sjeremylt **/ 18356c58de82SJeremy L Thompson int CeedBasisGetGrad(CeedBasis basis, const CeedScalar **grad) { 18366402da51SJeremy L Thompson if (!basis->grad && basis->is_tensor_basis) { 18379d007619Sjeremylt // Allocate 18382b730f8bSJeremy L Thompson CeedCall(CeedMalloc(basis->dim * basis->Q * basis->P, &basis->grad)); 18399d007619Sjeremylt 18409d007619Sjeremylt // Initialize 18412b730f8bSJeremy L Thompson for (CeedInt i = 0; i < basis->dim * basis->Q * basis->P; i++) basis->grad[i] = 1.0; 18429d007619Sjeremylt 18439d007619Sjeremylt // Calculate 18442b730f8bSJeremy L Thompson for (CeedInt d = 0; d < basis->dim; d++) { 18452b730f8bSJeremy L Thompson for (CeedInt i = 0; i < basis->dim; i++) { 18462b730f8bSJeremy L Thompson for (CeedInt qpt = 0; qpt < basis->Q; qpt++) { 18479d007619Sjeremylt for (CeedInt node = 0; node < basis->P; node++) { 1848d1d35e2fSjeremylt CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d; 1849d1d35e2fSjeremylt CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d; 18502b730f8bSJeremy L Thompson if (i == d) basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->grad_1d[q * basis->P_1d + p]; 18512b730f8bSJeremy L Thompson else basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p]; 18522b730f8bSJeremy L Thompson } 18532b730f8bSJeremy L Thompson } 18542b730f8bSJeremy L Thompson } 18559d007619Sjeremylt } 18569d007619Sjeremylt } 18579d007619Sjeremylt *grad = basis->grad; 1858e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 18599d007619Sjeremylt } 18609d007619Sjeremylt 18619d007619Sjeremylt /** 18629d007619Sjeremylt @brief Get 1D gradient matrix of a tensor product CeedBasis 18639d007619Sjeremylt 1864ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 1865d1d35e2fSjeremylt @param[out] grad_1d Variable to store gradient matrix 18669d007619Sjeremylt 18679d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 18689d007619Sjeremylt 1869b7c9bbdaSJeremy L Thompson @ref Advanced 18709d007619Sjeremylt **/ 1871d1d35e2fSjeremylt int CeedBasisGetGrad1D(CeedBasis basis, const CeedScalar **grad_1d) { 18726402da51SJeremy L Thompson CeedCheck(basis->is_tensor_basis, basis->ceed, CEED_ERROR_MINOR, "CeedBasis is not a tensor product basis."); 1873d1d35e2fSjeremylt *grad_1d = basis->grad_1d; 1874e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 18759d007619Sjeremylt } 18769d007619Sjeremylt 18779d007619Sjeremylt /** 187850c301a5SRezgar Shakeri @brief Get divergence matrix of a CeedBasis 187950c301a5SRezgar Shakeri 1880ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 188150c301a5SRezgar Shakeri @param[out] div Variable to store divergence matrix 188250c301a5SRezgar Shakeri 188350c301a5SRezgar Shakeri @return An error code: 0 - success, otherwise - failure 188450c301a5SRezgar Shakeri 188550c301a5SRezgar Shakeri @ref Advanced 188650c301a5SRezgar Shakeri **/ 188750c301a5SRezgar Shakeri int CeedBasisGetDiv(CeedBasis basis, const CeedScalar **div) { 18886574a04fSJeremy L Thompson CeedCheck(basis->div, basis->ceed, CEED_ERROR_MINOR, "CeedBasis does not have divergence matrix."); 188950c301a5SRezgar Shakeri *div = basis->div; 189050c301a5SRezgar Shakeri return CEED_ERROR_SUCCESS; 189150c301a5SRezgar Shakeri } 189250c301a5SRezgar Shakeri 189350c301a5SRezgar Shakeri /** 1894c4e3f59bSSebastian Grimberg @brief Get curl matrix of a CeedBasis 1895c4e3f59bSSebastian Grimberg 1896c4e3f59bSSebastian Grimberg @param[in] basis CeedBasis 1897c4e3f59bSSebastian Grimberg @param[out] curl Variable to store curl matrix 1898c4e3f59bSSebastian Grimberg 1899c4e3f59bSSebastian Grimberg @return An error code: 0 - success, otherwise - failure 1900c4e3f59bSSebastian Grimberg 1901c4e3f59bSSebastian Grimberg @ref Advanced 1902c4e3f59bSSebastian Grimberg **/ 1903c4e3f59bSSebastian Grimberg int CeedBasisGetCurl(CeedBasis basis, const CeedScalar **curl) { 19046574a04fSJeremy L Thompson CeedCheck(basis->curl, basis->ceed, CEED_ERROR_MINOR, "CeedBasis does not have curl matrix."); 1905c4e3f59bSSebastian Grimberg *curl = basis->curl; 1906c4e3f59bSSebastian Grimberg return CEED_ERROR_SUCCESS; 1907c4e3f59bSSebastian Grimberg } 1908c4e3f59bSSebastian Grimberg 1909c4e3f59bSSebastian Grimberg /** 19107a982d89SJeremy L. Thompson @brief Destroy a CeedBasis 19117a982d89SJeremy L. Thompson 1912ea61e9acSJeremy L Thompson @param[in,out] basis CeedBasis to destroy 19137a982d89SJeremy L. Thompson 19147a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 19157a982d89SJeremy L. Thompson 19167a982d89SJeremy L. Thompson @ref User 19177a982d89SJeremy L. Thompson **/ 19187a982d89SJeremy L. Thompson int CeedBasisDestroy(CeedBasis *basis) { 19197425e127SJeremy L Thompson if (!*basis || *basis == CEED_BASIS_COLLOCATED || --(*basis)->ref_count > 0) { 1920ad6481ceSJeremy L Thompson *basis = NULL; 1921ad6481ceSJeremy L Thompson return CEED_ERROR_SUCCESS; 1922ad6481ceSJeremy L Thompson } 19232b730f8bSJeremy L Thompson if ((*basis)->Destroy) CeedCall((*basis)->Destroy(*basis)); 19249831d45aSJeremy L Thompson CeedCall(CeedTensorContractDestroy(&(*basis)->contract)); 1925c4e3f59bSSebastian Grimberg CeedCall(CeedFree(&(*basis)->q_ref_1d)); 1926c4e3f59bSSebastian Grimberg CeedCall(CeedFree(&(*basis)->q_weight_1d)); 19272b730f8bSJeremy L Thompson CeedCall(CeedFree(&(*basis)->interp)); 19282b730f8bSJeremy L Thompson CeedCall(CeedFree(&(*basis)->interp_1d)); 19292b730f8bSJeremy L Thompson CeedCall(CeedFree(&(*basis)->grad)); 19302b730f8bSJeremy L Thompson CeedCall(CeedFree(&(*basis)->grad_1d)); 1931c4e3f59bSSebastian Grimberg CeedCall(CeedFree(&(*basis)->div)); 1932c4e3f59bSSebastian Grimberg CeedCall(CeedFree(&(*basis)->curl)); 1933c8c3fa7dSJeremy L Thompson CeedCall(CeedVectorDestroy(&(*basis)->vec_chebyshev)); 1934c8c3fa7dSJeremy L Thompson CeedCall(CeedBasisDestroy(&(*basis)->basis_chebyshev)); 19352b730f8bSJeremy L Thompson CeedCall(CeedDestroy(&(*basis)->ceed)); 19362b730f8bSJeremy L Thompson CeedCall(CeedFree(basis)); 1937e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 19387a982d89SJeremy L. Thompson } 19397a982d89SJeremy L. Thompson 19407a982d89SJeremy L. Thompson /** 1941b11c1e72Sjeremylt @brief Construct a Gauss-Legendre quadrature 1942b11c1e72Sjeremylt 1943ea61e9acSJeremy L Thompson @param[in] Q Number of quadrature points (integrates polynomials of degree 2*Q-1 exactly) 1944d1d35e2fSjeremylt @param[out] q_ref_1d Array of length Q to hold the abscissa on [-1, 1] 1945d1d35e2fSjeremylt @param[out] q_weight_1d Array of length Q to hold the weights 1946b11c1e72Sjeremylt 1947b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 1948dfdf5a53Sjeremylt 1949dfdf5a53Sjeremylt @ref Utility 1950b11c1e72Sjeremylt **/ 19512b730f8bSJeremy L Thompson int CeedGaussQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) { 1952d7b241e6Sjeremylt // Allocate 1953d7b241e6Sjeremylt CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0 * atan(1.0); 1954d1d35e2fSjeremylt // Build q_ref_1d, q_weight_1d 195592ae7e47SJeremy L Thompson for (CeedInt i = 0; i <= Q / 2; i++) { 1956d7b241e6Sjeremylt // Guess 1957d7b241e6Sjeremylt xi = cos(PI * (CeedScalar)(2 * i + 1) / ((CeedScalar)(2 * Q))); 1958d7b241e6Sjeremylt // Pn(xi) 1959d7b241e6Sjeremylt P0 = 1.0; 1960d7b241e6Sjeremylt P1 = xi; 1961d7b241e6Sjeremylt P2 = 0.0; 196292ae7e47SJeremy L Thompson for (CeedInt j = 2; j <= Q; j++) { 1963d7b241e6Sjeremylt P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j)); 1964d7b241e6Sjeremylt P0 = P1; 1965d7b241e6Sjeremylt P1 = P2; 1966d7b241e6Sjeremylt } 1967d7b241e6Sjeremylt // First Newton Step 1968d7b241e6Sjeremylt dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0); 1969d7b241e6Sjeremylt xi = xi - P2 / dP2; 1970d7b241e6Sjeremylt // Newton to convergence 197192ae7e47SJeremy L Thompson for (CeedInt k = 0; k < 100 && fabs(P2) > 10 * CEED_EPSILON; k++) { 1972d7b241e6Sjeremylt P0 = 1.0; 1973d7b241e6Sjeremylt P1 = xi; 197492ae7e47SJeremy L Thompson for (CeedInt j = 2; j <= Q; j++) { 1975d7b241e6Sjeremylt P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j)); 1976d7b241e6Sjeremylt P0 = P1; 1977d7b241e6Sjeremylt P1 = P2; 1978d7b241e6Sjeremylt } 1979d7b241e6Sjeremylt dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0); 1980d7b241e6Sjeremylt xi = xi - P2 / dP2; 1981d7b241e6Sjeremylt } 1982d7b241e6Sjeremylt // Save xi, wi 1983d7b241e6Sjeremylt wi = 2.0 / ((1.0 - xi * xi) * dP2 * dP2); 1984d1d35e2fSjeremylt q_weight_1d[i] = wi; 1985d1d35e2fSjeremylt q_weight_1d[Q - 1 - i] = wi; 1986d1d35e2fSjeremylt q_ref_1d[i] = -xi; 1987d1d35e2fSjeremylt q_ref_1d[Q - 1 - i] = xi; 1988d7b241e6Sjeremylt } 1989e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1990d7b241e6Sjeremylt } 1991d7b241e6Sjeremylt 1992b11c1e72Sjeremylt /** 1993b11c1e72Sjeremylt @brief Construct a Gauss-Legendre-Lobatto quadrature 1994b11c1e72Sjeremylt 1995ea61e9acSJeremy L Thompson @param[in] Q Number of quadrature points (integrates polynomials of degree 2*Q-3 exactly) 1996d1d35e2fSjeremylt @param[out] q_ref_1d Array of length Q to hold the abscissa on [-1, 1] 1997d1d35e2fSjeremylt @param[out] q_weight_1d Array of length Q to hold the weights 1998b11c1e72Sjeremylt 1999b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 2000dfdf5a53Sjeremylt 2001dfdf5a53Sjeremylt @ref Utility 2002b11c1e72Sjeremylt **/ 20032b730f8bSJeremy L Thompson int CeedLobattoQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) { 2004d7b241e6Sjeremylt // Allocate 2005d7b241e6Sjeremylt CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0 * atan(1.0); 2006d1d35e2fSjeremylt // Build q_ref_1d, q_weight_1d 2007d7b241e6Sjeremylt // Set endpoints 20086574a04fSJeremy L Thompson CeedCheck(Q > 1, NULL, CEED_ERROR_DIMENSION, "Cannot create Lobatto quadrature with Q=%" CeedInt_FMT " < 2 points", Q); 2009d7b241e6Sjeremylt wi = 2.0 / ((CeedScalar)(Q * (Q - 1))); 2010d1d35e2fSjeremylt if (q_weight_1d) { 2011d1d35e2fSjeremylt q_weight_1d[0] = wi; 2012d1d35e2fSjeremylt q_weight_1d[Q - 1] = wi; 2013d7b241e6Sjeremylt } 2014d1d35e2fSjeremylt q_ref_1d[0] = -1.0; 2015d1d35e2fSjeremylt q_ref_1d[Q - 1] = 1.0; 2016d7b241e6Sjeremylt // Interior 201792ae7e47SJeremy L Thompson for (CeedInt i = 1; i <= (Q - 1) / 2; i++) { 2018d7b241e6Sjeremylt // Guess 2019d7b241e6Sjeremylt xi = cos(PI * (CeedScalar)(i) / (CeedScalar)(Q - 1)); 2020d7b241e6Sjeremylt // Pn(xi) 2021d7b241e6Sjeremylt P0 = 1.0; 2022d7b241e6Sjeremylt P1 = xi; 2023d7b241e6Sjeremylt P2 = 0.0; 202492ae7e47SJeremy L Thompson for (CeedInt j = 2; j < Q; j++) { 2025d7b241e6Sjeremylt P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j)); 2026d7b241e6Sjeremylt P0 = P1; 2027d7b241e6Sjeremylt P1 = P2; 2028d7b241e6Sjeremylt } 2029d7b241e6Sjeremylt // First Newton step 2030d7b241e6Sjeremylt dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0); 2031d7b241e6Sjeremylt d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi); 2032d7b241e6Sjeremylt xi = xi - dP2 / d2P2; 2033d7b241e6Sjeremylt // Newton to convergence 203492ae7e47SJeremy L Thompson for (CeedInt k = 0; k < 100 && fabs(dP2) > 10 * CEED_EPSILON; k++) { 2035d7b241e6Sjeremylt P0 = 1.0; 2036d7b241e6Sjeremylt P1 = xi; 203792ae7e47SJeremy L Thompson for (CeedInt j = 2; j < Q; j++) { 2038d7b241e6Sjeremylt P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j)); 2039d7b241e6Sjeremylt P0 = P1; 2040d7b241e6Sjeremylt P1 = P2; 2041d7b241e6Sjeremylt } 2042d7b241e6Sjeremylt dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0); 2043d7b241e6Sjeremylt d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi); 2044d7b241e6Sjeremylt xi = xi - dP2 / d2P2; 2045d7b241e6Sjeremylt } 2046d7b241e6Sjeremylt // Save xi, wi 2047d7b241e6Sjeremylt wi = 2.0 / (((CeedScalar)(Q * (Q - 1))) * P2 * P2); 2048d1d35e2fSjeremylt if (q_weight_1d) { 2049d1d35e2fSjeremylt q_weight_1d[i] = wi; 2050d1d35e2fSjeremylt q_weight_1d[Q - 1 - i] = wi; 2051d7b241e6Sjeremylt } 2052d1d35e2fSjeremylt q_ref_1d[i] = -xi; 2053d1d35e2fSjeremylt q_ref_1d[Q - 1 - i] = xi; 2054d7b241e6Sjeremylt } 2055e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 2056d7b241e6Sjeremylt } 2057d7b241e6Sjeremylt 2058d7b241e6Sjeremylt /// @} 2059