13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3d7b241e6Sjeremylt // 43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 5d7b241e6Sjeremylt // 63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 7d7b241e6Sjeremylt 83d576824SJeremy L Thompson #include <ceed-impl.h> 949aac155SJeremy L Thompson #include <ceed.h> 102b730f8bSJeremy L Thompson #include <ceed/backend.h> 11d7b241e6Sjeremylt #include <math.h> 123d576824SJeremy L Thompson #include <stdbool.h> 13d7b241e6Sjeremylt #include <stdio.h> 14d7b241e6Sjeremylt #include <string.h> 15d7b241e6Sjeremylt 167a982d89SJeremy L. Thompson /// @file 177a982d89SJeremy L. Thompson /// Implementation of CeedBasis interfaces 187a982d89SJeremy L. Thompson 19d7b241e6Sjeremylt /// @cond DOXYGEN_SKIP 20783c99b3SValeria Barra static struct CeedBasis_private ceed_basis_collocated; 21d7b241e6Sjeremylt /// @endcond 22d7b241e6Sjeremylt 237a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser 247a982d89SJeremy L. Thompson /// @{ 257a982d89SJeremy L. Thompson 267a982d89SJeremy L. Thompson /// Indicate that the quadrature points are collocated with the nodes 277a982d89SJeremy L. Thompson const CeedBasis CEED_BASIS_COLLOCATED = &ceed_basis_collocated; 287a982d89SJeremy L. Thompson 297a982d89SJeremy L. Thompson /// @} 307a982d89SJeremy L. Thompson 317a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 327a982d89SJeremy L. Thompson /// CeedBasis Library Internal Functions 337a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 347a982d89SJeremy L. Thompson /// @addtogroup CeedBasisDeveloper 357a982d89SJeremy L. Thompson /// @{ 367a982d89SJeremy L. Thompson 377a982d89SJeremy L. Thompson /** 387a982d89SJeremy L. Thompson @brief Compute Householder reflection 397a982d89SJeremy L. Thompson 40ea61e9acSJeremy L Thompson Computes A = (I - b v v^T) A, where A is an mxn matrix indexed as A[i*row + j*col] 417a982d89SJeremy L. Thompson 427a982d89SJeremy L. Thompson @param[in,out] A Matrix to apply Householder reflection to, in place 43ea61e9acSJeremy L Thompson @param[in] v Householder vector 44ea61e9acSJeremy L Thompson @param[in] b Scaling factor 45ea61e9acSJeremy L Thompson @param[in] m Number of rows in A 46ea61e9acSJeremy L Thompson @param[in] n Number of columns in A 47ea61e9acSJeremy L Thompson @param[in] row Row stride 48ea61e9acSJeremy L Thompson @param[in] col Col stride 497a982d89SJeremy L. Thompson 507a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 517a982d89SJeremy L. Thompson 527a982d89SJeremy L. Thompson @ref Developer 537a982d89SJeremy L. Thompson **/ 542b730f8bSJeremy L Thompson static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v, CeedScalar b, CeedInt m, CeedInt n, CeedInt row, CeedInt col) { 557a982d89SJeremy L. Thompson for (CeedInt j = 0; j < n; j++) { 567a982d89SJeremy L. Thompson CeedScalar w = A[0 * row + j * col]; 572b730f8bSJeremy L Thompson for (CeedInt i = 1; i < m; i++) w += v[i] * A[i * row + j * col]; 587a982d89SJeremy L. Thompson A[0 * row + j * col] -= b * w; 592b730f8bSJeremy L Thompson for (CeedInt i = 1; i < m; i++) A[i * row + j * col] -= b * w * v[i]; 607a982d89SJeremy L. Thompson } 61e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 627a982d89SJeremy L. Thompson } 637a982d89SJeremy L. Thompson 647a982d89SJeremy L. Thompson /** 657a982d89SJeremy L. Thompson @brief Compute Givens rotation 667a982d89SJeremy L. Thompson 67ea61e9acSJeremy L Thompson Computes A = G A (or G^T A in transpose mode), where A is an mxn matrix indexed as A[i*n + j*m] 687a982d89SJeremy L. Thompson 697a982d89SJeremy L. Thompson @param[in,out] A Row major matrix to apply Givens rotation to, in place 70ea61e9acSJeremy L Thompson @param[in] c Cosine factor 71ea61e9acSJeremy L Thompson @param[in] s Sine factor 72ea61e9acSJeremy L Thompson @param[in] t_mode @ref CEED_NOTRANSPOSE to rotate the basis counter-clockwise, which has the effect of rotating columns of A clockwise; 734cc79fe7SJed Brown @ref CEED_TRANSPOSE for the opposite rotation 74ea61e9acSJeremy L Thompson @param[in] i First row/column to apply rotation 75ea61e9acSJeremy L Thompson @param[in] k Second row/column to apply rotation 76ea61e9acSJeremy L Thompson @param[in] m Number of rows in A 77ea61e9acSJeremy L Thompson @param[in] n Number of columns in A 787a982d89SJeremy L. Thompson 797a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 807a982d89SJeremy L. Thompson 817a982d89SJeremy L. Thompson @ref Developer 827a982d89SJeremy L. Thompson **/ 832b730f8bSJeremy L Thompson static int CeedGivensRotation(CeedScalar *A, CeedScalar c, CeedScalar s, CeedTransposeMode t_mode, CeedInt i, CeedInt k, CeedInt m, CeedInt n) { 84d1d35e2fSjeremylt CeedInt stride_j = 1, stride_ik = m, num_its = n; 85d1d35e2fSjeremylt if (t_mode == CEED_NOTRANSPOSE) { 862b730f8bSJeremy L Thompson stride_j = n; 872b730f8bSJeremy L Thompson stride_ik = 1; 882b730f8bSJeremy L Thompson num_its = m; 897a982d89SJeremy L. Thompson } 907a982d89SJeremy L. Thompson 917a982d89SJeremy L. Thompson // Apply rotation 92d1d35e2fSjeremylt for (CeedInt j = 0; j < num_its; j++) { 93d1d35e2fSjeremylt CeedScalar tau1 = A[i * stride_ik + j * stride_j], tau2 = A[k * stride_ik + j * stride_j]; 94d1d35e2fSjeremylt A[i * stride_ik + j * stride_j] = c * tau1 - s * tau2; 95d1d35e2fSjeremylt A[k * stride_ik + j * stride_j] = s * tau1 + c * tau2; 967a982d89SJeremy L. Thompson } 97e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 987a982d89SJeremy L. Thompson } 997a982d89SJeremy L. Thompson 1007a982d89SJeremy L. Thompson /** 1017a982d89SJeremy L. Thompson @brief View an array stored in a CeedBasis 1027a982d89SJeremy L. Thompson 1030a0da059Sjeremylt @param[in] name Name of array 104d1d35e2fSjeremylt @param[in] fp_fmt Printing format 1050a0da059Sjeremylt @param[in] m Number of rows in array 1060a0da059Sjeremylt @param[in] n Number of columns in array 1070a0da059Sjeremylt @param[in] a Array to be viewed 1080a0da059Sjeremylt @param[in] stream Stream to view to, e.g., stdout 1097a982d89SJeremy L. Thompson 1107a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 1117a982d89SJeremy L. Thompson 1127a982d89SJeremy L. Thompson @ref Developer 1137a982d89SJeremy L. Thompson **/ 1142b730f8bSJeremy L Thompson static int CeedScalarView(const char *name, const char *fp_fmt, CeedInt m, CeedInt n, const CeedScalar *a, FILE *stream) { 11592ae7e47SJeremy L Thompson for (CeedInt i = 0; i < m; i++) { 1162b730f8bSJeremy L Thompson if (m > 1) fprintf(stream, "%12s[%" CeedInt_FMT "]:", name, i); 1172b730f8bSJeremy L Thompson else fprintf(stream, "%12s:", name); 1182b730f8bSJeremy L Thompson for (CeedInt j = 0; j < n; j++) fprintf(stream, fp_fmt, fabs(a[i * n + j]) > 1E-14 ? a[i * n + j] : 0); 1197a982d89SJeremy L. Thompson fputs("\n", stream); 1207a982d89SJeremy L. Thompson } 121e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1227a982d89SJeremy L. Thompson } 1237a982d89SJeremy L. Thompson 124a76a04e7SJeremy L Thompson /** 125ea61e9acSJeremy L Thompson @brief Create the interpolation and gradient matrices for projection from the nodes of `basis_from` to the nodes of `basis_to`. 126ba59ac12SSebastian Grimberg 12715ad3917SSebastian Grimberg The interpolation is given by `interp_project = interp_to^+ * interp_from`, where the pseudoinverse `interp_to^+` is given by QR factorization. 12815ad3917SSebastian Grimberg The gradient is given by `grad_project = interp_to^+ * grad_from`, and is only computed for H^1 spaces otherwise it should not be used. 12915ad3917SSebastian Grimberg 130ba59ac12SSebastian Grimberg Note: `basis_from` and `basis_to` must have compatible quadrature spaces. 131a76a04e7SJeremy L Thompson 132a76a04e7SJeremy L Thompson @param[in] basis_from CeedBasis to project from 133a76a04e7SJeremy L Thompson @param[in] basis_to CeedBasis to project to 134ea61e9acSJeremy L Thompson @param[out] interp_project Address of the variable where the newly created interpolation matrix will be stored. 135ea61e9acSJeremy L Thompson @param[out] grad_project Address of the variable where the newly created gradient matrix will be stored. 136a76a04e7SJeremy L Thompson 137a76a04e7SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 138a76a04e7SJeremy L Thompson 139a76a04e7SJeremy L Thompson @ref Developer 140a76a04e7SJeremy L Thompson **/ 1412b730f8bSJeremy L Thompson static int CeedBasisCreateProjectionMatrices(CeedBasis basis_from, CeedBasis basis_to, CeedScalar **interp_project, CeedScalar **grad_project) { 142a76a04e7SJeremy L Thompson Ceed ceed; 1432b730f8bSJeremy L Thompson CeedCall(CeedBasisGetCeed(basis_to, &ceed)); 144a76a04e7SJeremy L Thompson 145a76a04e7SJeremy L Thompson // Check for compatible quadrature spaces 146a76a04e7SJeremy L Thompson CeedInt Q_to, Q_from; 1472b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints(basis_to, &Q_to)); 1482b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints(basis_from, &Q_from)); 149*6574a04fSJeremy L Thompson CeedCheck(Q_to == Q_from, ceed, CEED_ERROR_DIMENSION, "Bases must have compatible quadrature spaces"); 150a76a04e7SJeremy L Thompson 15114556e63SJeremy L Thompson // Check for matching tensor or non-tensor 152a76a04e7SJeremy L Thompson CeedInt P_to, P_from, Q = Q_to; 153a76a04e7SJeremy L Thompson bool is_tensor_to, is_tensor_from; 1542b730f8bSJeremy L Thompson CeedCall(CeedBasisIsTensor(basis_to, &is_tensor_to)); 1552b730f8bSJeremy L Thompson CeedCall(CeedBasisIsTensor(basis_from, &is_tensor_from)); 156*6574a04fSJeremy L Thompson CeedCheck(is_tensor_to == is_tensor_from, ceed, CEED_ERROR_MINOR, "Bases must both be tensor or non-tensor"); 157*6574a04fSJeremy L Thompson if (is_tensor_to) { 1582b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_to)); 1592b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_from)); 1602b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints1D(basis_from, &Q)); 161*6574a04fSJeremy L Thompson } else { 1622b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis_to, &P_to)); 1632b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis_from, &P_from)); 164a76a04e7SJeremy L Thompson } 165a76a04e7SJeremy L Thompson 16615ad3917SSebastian Grimberg // Check for matching FE space 16715ad3917SSebastian Grimberg CeedFESpace fe_space_to, fe_space_from; 16815ad3917SSebastian Grimberg CeedCall(CeedBasisGetFESpace(basis_to, &fe_space_to)); 16915ad3917SSebastian Grimberg CeedCall(CeedBasisGetFESpace(basis_from, &fe_space_from)); 170*6574a04fSJeremy L Thompson CeedCheck(fe_space_to == fe_space_from, ceed, CEED_ERROR_MINOR, "Bases must both be the same FE space type"); 17115ad3917SSebastian Grimberg 17214556e63SJeremy L Thompson // Get source matrices 17315ad3917SSebastian Grimberg CeedInt dim, q_comp = 1; 17415ad3917SSebastian Grimberg const CeedScalar *interp_to_source = NULL, *interp_from_source = NULL; 17514556e63SJeremy L Thompson CeedScalar *interp_to, *interp_from, *tau; 1762b730f8bSJeremy L Thompson CeedCall(CeedBasisGetDimension(basis_to, &dim)); 177a76a04e7SJeremy L Thompson if (is_tensor_to) { 1782b730f8bSJeremy L Thompson CeedCall(CeedBasisGetInterp1D(basis_to, &interp_to_source)); 1792b730f8bSJeremy L Thompson CeedCall(CeedBasisGetInterp1D(basis_from, &interp_from_source)); 180a76a04e7SJeremy L Thompson } else { 18115ad3917SSebastian Grimberg CeedCall(CeedBasisGetNumQuadratureComponents(basis_from, CEED_EVAL_INTERP, &q_comp)); 1822b730f8bSJeremy L Thompson CeedCall(CeedBasisGetInterp(basis_to, &interp_to_source)); 1832b730f8bSJeremy L Thompson CeedCall(CeedBasisGetInterp(basis_from, &interp_from_source)); 18415ad3917SSebastian Grimberg } 18515ad3917SSebastian Grimberg CeedCall(CeedMalloc(Q * P_from * q_comp, &interp_from)); 18615ad3917SSebastian Grimberg CeedCall(CeedMalloc(Q * P_to * q_comp, &interp_to)); 18715ad3917SSebastian Grimberg CeedCall(CeedCalloc(P_to * P_from, interp_project)); 18815ad3917SSebastian Grimberg CeedCall(CeedMalloc(Q * q_comp, &tau)); 18915ad3917SSebastian Grimberg 19015ad3917SSebastian Grimberg // `grad_project = interp_to^+ * grad_from` is computed for the H^1 space case so the 19115ad3917SSebastian Grimberg // projection basis will have a gradient operation 19215ad3917SSebastian Grimberg const CeedScalar *grad_from_source = NULL; 19315ad3917SSebastian Grimberg if (fe_space_to == CEED_FE_SPACE_H1) { 19415ad3917SSebastian Grimberg if (is_tensor_to) { 19515ad3917SSebastian Grimberg CeedCall(CeedBasisGetGrad1D(basis_from, &grad_from_source)); 19615ad3917SSebastian Grimberg } else { 1972b730f8bSJeremy L Thompson CeedCall(CeedBasisGetGrad(basis_from, &grad_from_source)); 198a76a04e7SJeremy L Thompson } 19915ad3917SSebastian Grimberg CeedCall(CeedCalloc(P_to * P_from * (is_tensor_to ? 1 : dim), grad_project)); 20015ad3917SSebastian Grimberg } else { 20115ad3917SSebastian Grimberg // Projected derivate operator is only calculated for H^1 but we allocate it for the 20215ad3917SSebastian Grimberg // basis construction later 20315ad3917SSebastian Grimberg CeedInt q_comp_deriv = 1; 20415ad3917SSebastian Grimberg if (fe_space_to == CEED_FE_SPACE_HDIV) { 20515ad3917SSebastian Grimberg CeedCall(CeedBasisGetNumQuadratureComponents(basis_from, CEED_EVAL_DIV, &q_comp_deriv)); 20615ad3917SSebastian Grimberg } else if (fe_space_to == CEED_FE_SPACE_HCURL) { 20715ad3917SSebastian Grimberg CeedCall(CeedBasisGetNumQuadratureComponents(basis_from, CEED_EVAL_CURL, &q_comp_deriv)); 20815ad3917SSebastian Grimberg } 20915ad3917SSebastian Grimberg CeedCall(CeedCalloc(P_to * P_from * q_comp_deriv, grad_project)); 21015ad3917SSebastian Grimberg } 21115ad3917SSebastian Grimberg 21215ad3917SSebastian Grimberg // QR Factorization, interp_to = Q R 21315ad3917SSebastian Grimberg memcpy(interp_to, interp_to_source, Q * P_to * q_comp * sizeof(interp_to_source[0])); 21415ad3917SSebastian Grimberg CeedCall(CeedQRFactorization(ceed, interp_to, tau, Q * q_comp, P_to)); 215a76a04e7SJeremy L Thompson 21614556e63SJeremy L Thompson // Build matrices 21715ad3917SSebastian Grimberg CeedInt num_matrices = 1 + (fe_space_to == CEED_FE_SPACE_H1) * (is_tensor_to ? 1 : dim); 21814556e63SJeremy L Thompson CeedScalar *input_from[num_matrices], *output_project[num_matrices]; 21914556e63SJeremy L Thompson input_from[0] = (CeedScalar *)interp_from_source; 22014556e63SJeremy L Thompson output_project[0] = *interp_project; 22114556e63SJeremy L Thompson for (CeedInt m = 1; m < num_matrices; m++) { 22214556e63SJeremy L Thompson input_from[m] = (CeedScalar *)&grad_from_source[(m - 1) * Q * P_from]; 22302af4036SJeremy L Thompson output_project[m] = &((*grad_project)[(m - 1) * P_to * P_from]); 22414556e63SJeremy L Thompson } 22514556e63SJeremy L Thompson for (CeedInt m = 0; m < num_matrices; m++) { 22615ad3917SSebastian Grimberg // Apply Q^T, interp_from = Q^T interp_from 22715ad3917SSebastian Grimberg memcpy(interp_from, input_from[m], Q * P_from * q_comp * sizeof(input_from[m][0])); 22815ad3917SSebastian Grimberg CeedCall(CeedHouseholderApplyQ(interp_from, interp_to, tau, CEED_TRANSPOSE, Q * q_comp, P_from, P_to, P_from, 1)); 229a76a04e7SJeremy L Thompson 23015ad3917SSebastian Grimberg // Apply Rinv, output_project = Rinv interp_from 231a76a04e7SJeremy L Thompson for (CeedInt j = 0; j < P_from; j++) { // Column j 2322b730f8bSJeremy 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]; 233a76a04e7SJeremy L Thompson for (CeedInt i = P_to - 2; i >= 0; i--) { // Row i 23414556e63SJeremy L Thompson output_project[m][j + P_from * i] = interp_from[j + P_from * i]; 235a76a04e7SJeremy L Thompson for (CeedInt k = i + 1; k < P_to; k++) { 2362b730f8bSJeremy L Thompson output_project[m][j + P_from * i] -= interp_to[k + P_to * i] * output_project[m][j + P_from * k]; 237a76a04e7SJeremy L Thompson } 23814556e63SJeremy L Thompson output_project[m][j + P_from * i] /= interp_to[i + P_to * i]; 239a76a04e7SJeremy L Thompson } 240a76a04e7SJeremy L Thompson } 24114556e63SJeremy L Thompson } 24214556e63SJeremy L Thompson 24314556e63SJeremy L Thompson // Cleanup 2442b730f8bSJeremy L Thompson CeedCall(CeedFree(&tau)); 2452b730f8bSJeremy L Thompson CeedCall(CeedFree(&interp_to)); 2462b730f8bSJeremy L Thompson CeedCall(CeedFree(&interp_from)); 247a76a04e7SJeremy L Thompson 248a76a04e7SJeremy L Thompson return CEED_ERROR_SUCCESS; 249a76a04e7SJeremy L Thompson } 250a76a04e7SJeremy L Thompson 2517a982d89SJeremy L. Thompson /// @} 2527a982d89SJeremy L. Thompson 2537a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 2547a982d89SJeremy L. Thompson /// Ceed Backend API 2557a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 2567a982d89SJeremy L. Thompson /// @addtogroup CeedBasisBackend 2577a982d89SJeremy L. Thompson /// @{ 2587a982d89SJeremy L. Thompson 2597a982d89SJeremy L. Thompson /** 2607a982d89SJeremy L. Thompson @brief Return collocated grad matrix 2617a982d89SJeremy L. Thompson 262ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 263ea61e9acSJeremy L Thompson @param[out] collo_grad_1d Row-major (Q_1d * Q_1d) matrix expressing derivatives of basis functions at quadrature points 2647a982d89SJeremy L. Thompson 2657a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 2667a982d89SJeremy L. Thompson 2677a982d89SJeremy L. Thompson @ref Backend 2687a982d89SJeremy L. Thompson **/ 269d1d35e2fSjeremylt int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *collo_grad_1d) { 2707a982d89SJeremy L. Thompson int i, j, k; 2717a982d89SJeremy L. Thompson Ceed ceed; 2722b730f8bSJeremy L Thompson CeedInt P_1d = (basis)->P_1d, Q_1d = (basis)->Q_1d; 27378464608Sjeremylt CeedScalar *interp_1d, *grad_1d, *tau; 2747a982d89SJeremy L. Thompson 2752b730f8bSJeremy L Thompson CeedCall(CeedMalloc(Q_1d * P_1d, &interp_1d)); 2762b730f8bSJeremy L Thompson CeedCall(CeedMalloc(Q_1d * P_1d, &grad_1d)); 2772b730f8bSJeremy L Thompson CeedCall(CeedMalloc(Q_1d, &tau)); 278d1d35e2fSjeremylt memcpy(interp_1d, (basis)->interp_1d, Q_1d * P_1d * sizeof(basis)->interp_1d[0]); 279d1d35e2fSjeremylt memcpy(grad_1d, (basis)->grad_1d, Q_1d * P_1d * sizeof(basis)->interp_1d[0]); 2807a982d89SJeremy L. Thompson 281d1d35e2fSjeremylt // QR Factorization, interp_1d = Q R 2822b730f8bSJeremy L Thompson CeedCall(CeedBasisGetCeed(basis, &ceed)); 2832b730f8bSJeremy L Thompson CeedCall(CeedQRFactorization(ceed, interp_1d, tau, Q_1d, P_1d)); 284ea61e9acSJeremy 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. 2857a982d89SJeremy L. Thompson 286d1d35e2fSjeremylt // Apply Rinv, collo_grad_1d = grad_1d Rinv 287d1d35e2fSjeremylt for (i = 0; i < Q_1d; i++) { // Row i 288d1d35e2fSjeremylt collo_grad_1d[Q_1d * i] = grad_1d[P_1d * i] / interp_1d[0]; 289d1d35e2fSjeremylt for (j = 1; j < P_1d; j++) { // Column j 290d1d35e2fSjeremylt collo_grad_1d[j + Q_1d * i] = grad_1d[j + P_1d * i]; 2912b730f8bSJeremy L Thompson for (k = 0; k < j; k++) collo_grad_1d[j + Q_1d * i] -= interp_1d[j + P_1d * k] * collo_grad_1d[k + Q_1d * i]; 292d1d35e2fSjeremylt collo_grad_1d[j + Q_1d * i] /= interp_1d[j + P_1d * j]; 2937a982d89SJeremy L. Thompson } 2942b730f8bSJeremy L Thompson for (j = P_1d; j < Q_1d; j++) collo_grad_1d[j + Q_1d * i] = 0; 2957a982d89SJeremy L. Thompson } 2967a982d89SJeremy L. Thompson 29715ad3917SSebastian Grimberg // Apply Q^T, collo_grad_1d = collo_grad_1d Q^T 2982b730f8bSJeremy L Thompson CeedCall(CeedHouseholderApplyQ(collo_grad_1d, interp_1d, tau, CEED_NOTRANSPOSE, Q_1d, Q_1d, P_1d, 1, Q_1d)); 2997a982d89SJeremy L. Thompson 3002b730f8bSJeremy L Thompson CeedCall(CeedFree(&interp_1d)); 3012b730f8bSJeremy L Thompson CeedCall(CeedFree(&grad_1d)); 3022b730f8bSJeremy L Thompson CeedCall(CeedFree(&tau)); 303e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3047a982d89SJeremy L. Thompson } 3057a982d89SJeremy L. Thompson 3067a982d89SJeremy L. Thompson /** 3077a982d89SJeremy L. Thompson @brief Get tensor status for given CeedBasis 3087a982d89SJeremy L. Thompson 309ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 310d1d35e2fSjeremylt @param[out] is_tensor Variable to store tensor status 3117a982d89SJeremy L. Thompson 3127a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3137a982d89SJeremy L. Thompson 3147a982d89SJeremy L. Thompson @ref Backend 3157a982d89SJeremy L. Thompson **/ 316d1d35e2fSjeremylt int CeedBasisIsTensor(CeedBasis basis, bool *is_tensor) { 317d1d35e2fSjeremylt *is_tensor = basis->tensor_basis; 318e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3197a982d89SJeremy L. Thompson } 3207a982d89SJeremy L. Thompson 3217a982d89SJeremy L. Thompson /** 3227a982d89SJeremy L. Thompson @brief Get backend data of a CeedBasis 3237a982d89SJeremy L. Thompson 324ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 3257a982d89SJeremy L. Thompson @param[out] data Variable to store data 3267a982d89SJeremy L. Thompson 3277a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3287a982d89SJeremy L. Thompson 3297a982d89SJeremy L. Thompson @ref Backend 3307a982d89SJeremy L. Thompson **/ 331777ff853SJeremy L Thompson int CeedBasisGetData(CeedBasis basis, void *data) { 332777ff853SJeremy L Thompson *(void **)data = basis->data; 333e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3347a982d89SJeremy L. Thompson } 3357a982d89SJeremy L. Thompson 3367a982d89SJeremy L. Thompson /** 3377a982d89SJeremy L. Thompson @brief Set backend data of a CeedBasis 3387a982d89SJeremy L. Thompson 339ea61e9acSJeremy L Thompson @param[in,out] basis CeedBasis 340ea61e9acSJeremy L Thompson @param[in] data Data to set 3417a982d89SJeremy L. Thompson 3427a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3437a982d89SJeremy L. Thompson 3447a982d89SJeremy L. Thompson @ref Backend 3457a982d89SJeremy L. Thompson **/ 346777ff853SJeremy L Thompson int CeedBasisSetData(CeedBasis basis, void *data) { 347777ff853SJeremy L Thompson basis->data = data; 348e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3497a982d89SJeremy L. Thompson } 3507a982d89SJeremy L. Thompson 3517a982d89SJeremy L. Thompson /** 35234359f16Sjeremylt @brief Increment the reference counter for a CeedBasis 35334359f16Sjeremylt 354ea61e9acSJeremy L Thompson @param[in,out] basis Basis to increment the reference counter 35534359f16Sjeremylt 35634359f16Sjeremylt @return An error code: 0 - success, otherwise - failure 35734359f16Sjeremylt 35834359f16Sjeremylt @ref Backend 35934359f16Sjeremylt **/ 3609560d06aSjeremylt int CeedBasisReference(CeedBasis basis) { 36134359f16Sjeremylt basis->ref_count++; 36234359f16Sjeremylt return CEED_ERROR_SUCCESS; 36334359f16Sjeremylt } 36434359f16Sjeremylt 36534359f16Sjeremylt /** 366c4e3f59bSSebastian Grimberg @brief Get number of Q-vector components for given CeedBasis 367c4e3f59bSSebastian Grimberg 368c4e3f59bSSebastian Grimberg @param[in] basis CeedBasis 369c4e3f59bSSebastian Grimberg @param[in] eval_mode \ref CEED_EVAL_INTERP to use interpolated values, 370c4e3f59bSSebastian Grimberg \ref CEED_EVAL_GRAD to use gradients, 371c4e3f59bSSebastian Grimberg \ref CEED_EVAL_DIV to use divergence, 372c4e3f59bSSebastian Grimberg \ref CEED_EVAL_CURL to use curl. 373c4e3f59bSSebastian Grimberg @param[out] q_comp Variable to store number of Q-vector components of basis 374c4e3f59bSSebastian Grimberg 375c4e3f59bSSebastian Grimberg @return An error code: 0 - success, otherwise - failure 376c4e3f59bSSebastian Grimberg 377c4e3f59bSSebastian Grimberg @ref Backend 378c4e3f59bSSebastian Grimberg **/ 379c4e3f59bSSebastian Grimberg int CeedBasisGetNumQuadratureComponents(CeedBasis basis, CeedEvalMode eval_mode, CeedInt *q_comp) { 380c4e3f59bSSebastian Grimberg switch (eval_mode) { 381c4e3f59bSSebastian Grimberg case CEED_EVAL_INTERP: 382c4e3f59bSSebastian Grimberg *q_comp = (basis->fe_space == CEED_FE_SPACE_H1) ? 1 : basis->dim; 383c4e3f59bSSebastian Grimberg break; 384c4e3f59bSSebastian Grimberg case CEED_EVAL_GRAD: 385c4e3f59bSSebastian Grimberg *q_comp = basis->dim; 386c4e3f59bSSebastian Grimberg break; 387c4e3f59bSSebastian Grimberg case CEED_EVAL_DIV: 388c4e3f59bSSebastian Grimberg *q_comp = 1; 389c4e3f59bSSebastian Grimberg break; 390c4e3f59bSSebastian Grimberg case CEED_EVAL_CURL: 391c4e3f59bSSebastian Grimberg *q_comp = (basis->dim < 3) ? 1 : basis->dim; 392c4e3f59bSSebastian Grimberg break; 393c4e3f59bSSebastian Grimberg case CEED_EVAL_NONE: 394c4e3f59bSSebastian Grimberg case CEED_EVAL_WEIGHT: 395352a5e7cSSebastian Grimberg *q_comp = 1; 396c4e3f59bSSebastian Grimberg break; 397c4e3f59bSSebastian Grimberg } 398c4e3f59bSSebastian Grimberg return CEED_ERROR_SUCCESS; 399c4e3f59bSSebastian Grimberg } 400c4e3f59bSSebastian Grimberg 401c4e3f59bSSebastian Grimberg /** 4026e15d496SJeremy L Thompson @brief Estimate number of FLOPs required to apply CeedBasis in t_mode and eval_mode 4036e15d496SJeremy L Thompson 404ea61e9acSJeremy L Thompson @param[in] basis Basis to estimate FLOPs for 405ea61e9acSJeremy L Thompson @param[in] t_mode Apply basis or transpose 406ea61e9acSJeremy L Thompson @param[in] eval_mode Basis evaluation mode 407ea61e9acSJeremy L Thompson @param[out] flops Address of variable to hold FLOPs estimate 4086e15d496SJeremy L Thompson 4096e15d496SJeremy L Thompson @ref Backend 4106e15d496SJeremy L Thompson **/ 4112b730f8bSJeremy L Thompson int CeedBasisGetFlopsEstimate(CeedBasis basis, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedSize *flops) { 4126e15d496SJeremy L Thompson bool is_tensor; 4136e15d496SJeremy L Thompson 4142b730f8bSJeremy L Thompson CeedCall(CeedBasisIsTensor(basis, &is_tensor)); 4156e15d496SJeremy L Thompson if (is_tensor) { 4166e15d496SJeremy L Thompson CeedInt dim, num_comp, P_1d, Q_1d; 4172b730f8bSJeremy L Thompson CeedCall(CeedBasisGetDimension(basis, &dim)); 4182b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumComponents(basis, &num_comp)); 4192b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d)); 4202b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 4216e15d496SJeremy L Thompson if (t_mode == CEED_TRANSPOSE) { 4222b730f8bSJeremy L Thompson P_1d = Q_1d; 4232b730f8bSJeremy L Thompson Q_1d = P_1d; 4246e15d496SJeremy L Thompson } 4256e15d496SJeremy L Thompson CeedInt tensor_flops = 0, pre = num_comp * CeedIntPow(P_1d, dim - 1), post = 1; 4266e15d496SJeremy L Thompson for (CeedInt d = 0; d < dim; d++) { 4276e15d496SJeremy L Thompson tensor_flops += 2 * pre * P_1d * post * Q_1d; 4286e15d496SJeremy L Thompson pre /= P_1d; 4296e15d496SJeremy L Thompson post *= Q_1d; 4306e15d496SJeremy L Thompson } 4316e15d496SJeremy L Thompson switch (eval_mode) { 4322b730f8bSJeremy L Thompson case CEED_EVAL_NONE: 4332b730f8bSJeremy L Thompson *flops = 0; 4342b730f8bSJeremy L Thompson break; 4352b730f8bSJeremy L Thompson case CEED_EVAL_INTERP: 4362b730f8bSJeremy L Thompson *flops = tensor_flops; 4372b730f8bSJeremy L Thompson break; 4382b730f8bSJeremy L Thompson case CEED_EVAL_GRAD: 4392b730f8bSJeremy L Thompson *flops = tensor_flops * 2; 4402b730f8bSJeremy L Thompson break; 4416e15d496SJeremy L Thompson case CEED_EVAL_DIV: 4426e15d496SJeremy L Thompson case CEED_EVAL_CURL: 443*6574a04fSJeremy L Thompson // LCOV_EXCL_START 444*6574a04fSJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_INCOMPATIBLE, "Tensor basis evaluation for %s not supported", CeedEvalModes[eval_mode]); 4452b730f8bSJeremy L Thompson break; 4466e15d496SJeremy L Thompson // LCOV_EXCL_STOP 4472b730f8bSJeremy L Thompson case CEED_EVAL_WEIGHT: 4482b730f8bSJeremy L Thompson *flops = dim * CeedIntPow(Q_1d, dim); 4492b730f8bSJeremy L Thompson break; 4506e15d496SJeremy L Thompson } 4516e15d496SJeremy L Thompson } else { 452c4e3f59bSSebastian Grimberg CeedInt dim, num_comp, q_comp, num_nodes, num_qpts; 4532b730f8bSJeremy L Thompson CeedCall(CeedBasisGetDimension(basis, &dim)); 4542b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumComponents(basis, &num_comp)); 455c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp)); 4562b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis, &num_nodes)); 4572b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts)); 4586e15d496SJeremy L Thompson switch (eval_mode) { 4592b730f8bSJeremy L Thompson case CEED_EVAL_NONE: 4602b730f8bSJeremy L Thompson *flops = 0; 4612b730f8bSJeremy L Thompson break; 4622b730f8bSJeremy L Thompson case CEED_EVAL_INTERP: 4632b730f8bSJeremy L Thompson case CEED_EVAL_GRAD: 4642b730f8bSJeremy L Thompson case CEED_EVAL_DIV: 4652b730f8bSJeremy L Thompson case CEED_EVAL_CURL: 466c4e3f59bSSebastian Grimberg *flops = num_nodes * num_qpts * num_comp * q_comp; 4672b730f8bSJeremy L Thompson break; 4682b730f8bSJeremy L Thompson case CEED_EVAL_WEIGHT: 4692b730f8bSJeremy L Thompson *flops = 0; 4702b730f8bSJeremy L Thompson break; 4716e15d496SJeremy L Thompson } 4726e15d496SJeremy L Thompson } 4736e15d496SJeremy L Thompson 4746e15d496SJeremy L Thompson return CEED_ERROR_SUCCESS; 4756e15d496SJeremy L Thompson } 4766e15d496SJeremy L Thompson 4776e15d496SJeremy L Thompson /** 478c4e3f59bSSebastian Grimberg @brief Get CeedFESpace for a CeedBasis 479c4e3f59bSSebastian Grimberg 480c4e3f59bSSebastian Grimberg @param[in] basis CeedBasis 481c4e3f59bSSebastian Grimberg @param[out] fe_space Variable to store CeedFESpace 482c4e3f59bSSebastian Grimberg 483c4e3f59bSSebastian Grimberg @return An error code: 0 - success, otherwise - failure 484c4e3f59bSSebastian Grimberg 485c4e3f59bSSebastian Grimberg @ref Backend 486c4e3f59bSSebastian Grimberg **/ 487c4e3f59bSSebastian Grimberg int CeedBasisGetFESpace(CeedBasis basis, CeedFESpace *fe_space) { 488c4e3f59bSSebastian Grimberg *fe_space = basis->fe_space; 489c4e3f59bSSebastian Grimberg return CEED_ERROR_SUCCESS; 490c4e3f59bSSebastian Grimberg } 491c4e3f59bSSebastian Grimberg 492c4e3f59bSSebastian Grimberg /** 4937a982d89SJeremy L. Thompson @brief Get dimension for given CeedElemTopology 4947a982d89SJeremy L. Thompson 495ea61e9acSJeremy L Thompson @param[in] topo CeedElemTopology 4967a982d89SJeremy L. Thompson @param[out] dim Variable to store dimension of topology 4977a982d89SJeremy L. Thompson 4987a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 4997a982d89SJeremy L. Thompson 5007a982d89SJeremy L. Thompson @ref Backend 5017a982d89SJeremy L. Thompson **/ 5027a982d89SJeremy L. Thompson int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) { 5037a982d89SJeremy L. Thompson *dim = (CeedInt)topo >> 16; 504e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 5057a982d89SJeremy L. Thompson } 5067a982d89SJeremy L. Thompson 5077a982d89SJeremy L. Thompson /** 5087a982d89SJeremy L. Thompson @brief Get CeedTensorContract of a CeedBasis 5097a982d89SJeremy L. Thompson 510ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 5117a982d89SJeremy L. Thompson @param[out] contract Variable to store CeedTensorContract 5127a982d89SJeremy L. Thompson 5137a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 5147a982d89SJeremy L. Thompson 5157a982d89SJeremy L. Thompson @ref Backend 5167a982d89SJeremy L. Thompson **/ 5177a982d89SJeremy L. Thompson int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) { 5187a982d89SJeremy L. Thompson *contract = basis->contract; 519e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 5207a982d89SJeremy L. Thompson } 5217a982d89SJeremy L. Thompson 5227a982d89SJeremy L. Thompson /** 5237a982d89SJeremy L. Thompson @brief Set CeedTensorContract of a CeedBasis 5247a982d89SJeremy L. Thompson 525ea61e9acSJeremy L Thompson @param[in,out] basis CeedBasis 526ea61e9acSJeremy L Thompson @param[in] contract CeedTensorContract to set 5277a982d89SJeremy L. Thompson 5287a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 5297a982d89SJeremy L. Thompson 5307a982d89SJeremy L. Thompson @ref Backend 5317a982d89SJeremy L. Thompson **/ 53234359f16Sjeremylt int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract contract) { 53334359f16Sjeremylt basis->contract = contract; 5342b730f8bSJeremy L Thompson CeedCall(CeedTensorContractReference(contract)); 535e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 5367a982d89SJeremy L. Thompson } 5377a982d89SJeremy L. Thompson 5387a982d89SJeremy L. Thompson /** 5397a982d89SJeremy L. Thompson @brief Return a reference implementation of matrix multiplication C = A B. 540ba59ac12SSebastian Grimberg 541ba59ac12SSebastian Grimberg Note: This is a reference implementation for CPU CeedScalar pointers that is not intended for high performance. 5427a982d89SJeremy L. Thompson 543ea61e9acSJeremy L Thompson @param[in] ceed Ceed context for error handling 544d1d35e2fSjeremylt @param[in] mat_A Row-major matrix A 545d1d35e2fSjeremylt @param[in] mat_B Row-major matrix B 546d1d35e2fSjeremylt @param[out] mat_C Row-major output matrix C 547ea61e9acSJeremy L Thompson @param[in] m Number of rows of C 548ea61e9acSJeremy L Thompson @param[in] n Number of columns of C 549ea61e9acSJeremy L Thompson @param[in] kk Number of columns of A/rows of B 5507a982d89SJeremy L. Thompson 5517a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 5527a982d89SJeremy L. Thompson 5537a982d89SJeremy L. Thompson @ref Utility 5547a982d89SJeremy L. Thompson **/ 5552b730f8bSJeremy L Thompson int CeedMatrixMatrixMultiply(Ceed ceed, const CeedScalar *mat_A, const CeedScalar *mat_B, CeedScalar *mat_C, CeedInt m, CeedInt n, CeedInt kk) { 5562b730f8bSJeremy L Thompson for (CeedInt i = 0; i < m; i++) { 5577a982d89SJeremy L. Thompson for (CeedInt j = 0; j < n; j++) { 5587a982d89SJeremy L. Thompson CeedScalar sum = 0; 5592b730f8bSJeremy L Thompson for (CeedInt k = 0; k < kk; k++) sum += mat_A[k + i * kk] * mat_B[j + k * n]; 560d1d35e2fSjeremylt mat_C[j + i * n] = sum; 5617a982d89SJeremy L. Thompson } 5622b730f8bSJeremy L Thompson } 563e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 5647a982d89SJeremy L. Thompson } 5657a982d89SJeremy L. Thompson 566ba59ac12SSebastian Grimberg /** 567ba59ac12SSebastian Grimberg @brief Return QR Factorization of a matrix 568ba59ac12SSebastian Grimberg 569ba59ac12SSebastian Grimberg @param[in] ceed Ceed context for error handling 570ba59ac12SSebastian Grimberg @param[in,out] mat Row-major matrix to be factorized in place 571ba59ac12SSebastian Grimberg @param[in,out] tau Vector of length m of scaling factors 572ba59ac12SSebastian Grimberg @param[in] m Number of rows 573ba59ac12SSebastian Grimberg @param[in] n Number of columns 574ba59ac12SSebastian Grimberg 575ba59ac12SSebastian Grimberg @return An error code: 0 - success, otherwise - failure 576ba59ac12SSebastian Grimberg 577ba59ac12SSebastian Grimberg @ref Utility 578ba59ac12SSebastian Grimberg **/ 579ba59ac12SSebastian Grimberg int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau, CeedInt m, CeedInt n) { 580ba59ac12SSebastian Grimberg CeedScalar v[m]; 581ba59ac12SSebastian Grimberg 582ba59ac12SSebastian Grimberg // Check matrix shape 583*6574a04fSJeremy L Thompson CeedCheck(n <= m, ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute QR factorization with n > m"); 584ba59ac12SSebastian Grimberg 585ba59ac12SSebastian Grimberg for (CeedInt i = 0; i < n; i++) { 586ba59ac12SSebastian Grimberg if (i >= m - 1) { // last row of matrix, no reflection needed 587ba59ac12SSebastian Grimberg tau[i] = 0.; 588ba59ac12SSebastian Grimberg break; 589ba59ac12SSebastian Grimberg } 590ba59ac12SSebastian Grimberg // Calculate Householder vector, magnitude 591ba59ac12SSebastian Grimberg CeedScalar sigma = 0.0; 592ba59ac12SSebastian Grimberg v[i] = mat[i + n * i]; 593ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < m; j++) { 594ba59ac12SSebastian Grimberg v[j] = mat[i + n * j]; 595ba59ac12SSebastian Grimberg sigma += v[j] * v[j]; 596ba59ac12SSebastian Grimberg } 597ba59ac12SSebastian Grimberg CeedScalar norm = sqrt(v[i] * v[i] + sigma); // norm of v[i:m] 598ba59ac12SSebastian Grimberg CeedScalar R_ii = -copysign(norm, v[i]); 599ba59ac12SSebastian Grimberg v[i] -= R_ii; 600ba59ac12SSebastian Grimberg // norm of v[i:m] after modification above and scaling below 601ba59ac12SSebastian Grimberg // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 602ba59ac12SSebastian Grimberg // tau = 2 / (norm*norm) 603ba59ac12SSebastian Grimberg tau[i] = 2 * v[i] * v[i] / (v[i] * v[i] + sigma); 604ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < m; j++) v[j] /= v[i]; 605ba59ac12SSebastian Grimberg 606ba59ac12SSebastian Grimberg // Apply Householder reflector to lower right panel 607ba59ac12SSebastian Grimberg CeedHouseholderReflect(&mat[i * n + i + 1], &v[i], tau[i], m - i, n - i - 1, n, 1); 608ba59ac12SSebastian Grimberg // Save v 609ba59ac12SSebastian Grimberg mat[i + n * i] = R_ii; 610ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < m; j++) mat[i + n * j] = v[j]; 611ba59ac12SSebastian Grimberg } 612ba59ac12SSebastian Grimberg return CEED_ERROR_SUCCESS; 613ba59ac12SSebastian Grimberg } 614ba59ac12SSebastian Grimberg 615ba59ac12SSebastian Grimberg /** 616ba59ac12SSebastian Grimberg @brief Apply Householder Q matrix 617ba59ac12SSebastian Grimberg 618ba59ac12SSebastian Grimberg Compute mat_A = mat_Q mat_A, where mat_Q is mxm and mat_A is mxn. 619ba59ac12SSebastian Grimberg 620ba59ac12SSebastian Grimberg @param[in,out] mat_A Matrix to apply Householder Q to, in place 621ba59ac12SSebastian Grimberg @param[in] mat_Q Householder Q matrix 622ba59ac12SSebastian Grimberg @param[in] tau Householder scaling factors 623ba59ac12SSebastian Grimberg @param[in] t_mode Transpose mode for application 624ba59ac12SSebastian Grimberg @param[in] m Number of rows in A 625ba59ac12SSebastian Grimberg @param[in] n Number of columns in A 626ba59ac12SSebastian Grimberg @param[in] k Number of elementary reflectors in Q, k<m 627ba59ac12SSebastian Grimberg @param[in] row Row stride in A 628ba59ac12SSebastian Grimberg @param[in] col Col stride in A 629ba59ac12SSebastian Grimberg 630ba59ac12SSebastian Grimberg @return An error code: 0 - success, otherwise - failure 631ba59ac12SSebastian Grimberg 632c4e3f59bSSebastian Grimberg @ref Utility 633ba59ac12SSebastian Grimberg **/ 634ba59ac12SSebastian Grimberg int CeedHouseholderApplyQ(CeedScalar *mat_A, const CeedScalar *mat_Q, const CeedScalar *tau, CeedTransposeMode t_mode, CeedInt m, CeedInt n, 635ba59ac12SSebastian Grimberg CeedInt k, CeedInt row, CeedInt col) { 636ba59ac12SSebastian Grimberg CeedScalar *v; 637ba59ac12SSebastian Grimberg CeedCall(CeedMalloc(m, &v)); 638ba59ac12SSebastian Grimberg for (CeedInt ii = 0; ii < k; ii++) { 639ba59ac12SSebastian Grimberg CeedInt i = t_mode == CEED_TRANSPOSE ? ii : k - 1 - ii; 640ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < m; j++) v[j] = mat_Q[j * k + i]; 641ba59ac12SSebastian Grimberg // Apply Householder reflector (I - tau v v^T) collo_grad_1d^T 642ba59ac12SSebastian Grimberg CeedCall(CeedHouseholderReflect(&mat_A[i * row], &v[i], tau[i], m - i, n, row, col)); 643ba59ac12SSebastian Grimberg } 644ba59ac12SSebastian Grimberg CeedCall(CeedFree(&v)); 645ba59ac12SSebastian Grimberg return CEED_ERROR_SUCCESS; 646ba59ac12SSebastian Grimberg } 647ba59ac12SSebastian Grimberg 648ba59ac12SSebastian Grimberg /** 649ba59ac12SSebastian Grimberg @brief Return symmetric Schur decomposition of the symmetric matrix mat via symmetric QR factorization 650ba59ac12SSebastian Grimberg 651ba59ac12SSebastian Grimberg @param[in] ceed Ceed context for error handling 652ba59ac12SSebastian Grimberg @param[in,out] mat Row-major matrix to be factorized in place 653ba59ac12SSebastian Grimberg @param[out] lambda Vector of length n of eigenvalues 654ba59ac12SSebastian Grimberg @param[in] n Number of rows/columns 655ba59ac12SSebastian Grimberg 656ba59ac12SSebastian Grimberg @return An error code: 0 - success, otherwise - failure 657ba59ac12SSebastian Grimberg 658ba59ac12SSebastian Grimberg @ref Utility 659ba59ac12SSebastian Grimberg **/ 6602c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOff 6612c2ea1dbSJeremy L Thompson int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat, CeedScalar *lambda, CeedInt n) { 662ba59ac12SSebastian Grimberg // Check bounds for clang-tidy 663*6574a04fSJeremy L Thompson CeedCheck(n > 1, ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute symmetric Schur decomposition of scalars"); 664ba59ac12SSebastian Grimberg 665ba59ac12SSebastian Grimberg CeedScalar v[n - 1], tau[n - 1], mat_T[n * n]; 666ba59ac12SSebastian Grimberg 667ba59ac12SSebastian Grimberg // Copy mat to mat_T and set mat to I 668ba59ac12SSebastian Grimberg memcpy(mat_T, mat, n * n * sizeof(mat[0])); 669ba59ac12SSebastian Grimberg for (CeedInt i = 0; i < n; i++) { 670ba59ac12SSebastian Grimberg for (CeedInt j = 0; j < n; j++) mat[j + n * i] = (i == j) ? 1 : 0; 671ba59ac12SSebastian Grimberg } 672ba59ac12SSebastian Grimberg 673ba59ac12SSebastian Grimberg // Reduce to tridiagonal 674ba59ac12SSebastian Grimberg for (CeedInt i = 0; i < n - 1; i++) { 675ba59ac12SSebastian Grimberg // Calculate Householder vector, magnitude 676ba59ac12SSebastian Grimberg CeedScalar sigma = 0.0; 677ba59ac12SSebastian Grimberg v[i] = mat_T[i + n * (i + 1)]; 678ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < n - 1; j++) { 679ba59ac12SSebastian Grimberg v[j] = mat_T[i + n * (j + 1)]; 680ba59ac12SSebastian Grimberg sigma += v[j] * v[j]; 681ba59ac12SSebastian Grimberg } 682ba59ac12SSebastian Grimberg CeedScalar norm = sqrt(v[i] * v[i] + sigma); // norm of v[i:n-1] 683ba59ac12SSebastian Grimberg CeedScalar R_ii = -copysign(norm, v[i]); 684ba59ac12SSebastian Grimberg v[i] -= R_ii; 685ba59ac12SSebastian Grimberg // norm of v[i:m] after modification above and scaling below 686ba59ac12SSebastian Grimberg // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 687ba59ac12SSebastian Grimberg // tau = 2 / (norm*norm) 688ba59ac12SSebastian Grimberg tau[i] = i == n - 2 ? 2 : 2 * v[i] * v[i] / (v[i] * v[i] + sigma); 689ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < n - 1; j++) v[j] /= v[i]; 690ba59ac12SSebastian Grimberg 691ba59ac12SSebastian Grimberg // Update sub and super diagonal 692ba59ac12SSebastian Grimberg for (CeedInt j = i + 2; j < n; j++) { 693ba59ac12SSebastian Grimberg mat_T[i + n * j] = 0; 694ba59ac12SSebastian Grimberg mat_T[j + n * i] = 0; 695ba59ac12SSebastian Grimberg } 696ba59ac12SSebastian Grimberg // Apply symmetric Householder reflector to lower right panel 697ba59ac12SSebastian Grimberg CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1); 698ba59ac12SSebastian Grimberg CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), 1, n); 699ba59ac12SSebastian Grimberg 700ba59ac12SSebastian Grimberg // Save v 701ba59ac12SSebastian Grimberg mat_T[i + n * (i + 1)] = R_ii; 702ba59ac12SSebastian Grimberg mat_T[(i + 1) + n * i] = R_ii; 703ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < n - 1; j++) { 704ba59ac12SSebastian Grimberg mat_T[i + n * (j + 1)] = v[j]; 705ba59ac12SSebastian Grimberg } 706ba59ac12SSebastian Grimberg } 707ba59ac12SSebastian Grimberg // Backwards accumulation of Q 708ba59ac12SSebastian Grimberg for (CeedInt i = n - 2; i >= 0; i--) { 709ba59ac12SSebastian Grimberg if (tau[i] > 0.0) { 710ba59ac12SSebastian Grimberg v[i] = 1; 711ba59ac12SSebastian Grimberg for (CeedInt j = i + 1; j < n - 1; j++) { 712ba59ac12SSebastian Grimberg v[j] = mat_T[i + n * (j + 1)]; 713ba59ac12SSebastian Grimberg mat_T[i + n * (j + 1)] = 0; 714ba59ac12SSebastian Grimberg } 715ba59ac12SSebastian Grimberg CeedHouseholderReflect(&mat[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1); 716ba59ac12SSebastian Grimberg } 717ba59ac12SSebastian Grimberg } 718ba59ac12SSebastian Grimberg 719ba59ac12SSebastian Grimberg // Reduce sub and super diagonal 720ba59ac12SSebastian Grimberg CeedInt p = 0, q = 0, itr = 0, max_itr = n * n * n * n; 721ba59ac12SSebastian Grimberg CeedScalar tol = CEED_EPSILON; 722ba59ac12SSebastian Grimberg 723ba59ac12SSebastian Grimberg while (itr < max_itr) { 724ba59ac12SSebastian Grimberg // Update p, q, size of reduced portions of diagonal 725ba59ac12SSebastian Grimberg p = 0; 726ba59ac12SSebastian Grimberg q = 0; 727ba59ac12SSebastian Grimberg for (CeedInt i = n - 2; i >= 0; i--) { 728ba59ac12SSebastian Grimberg if (fabs(mat_T[i + n * (i + 1)]) < tol) q += 1; 729ba59ac12SSebastian Grimberg else break; 730ba59ac12SSebastian Grimberg } 731ba59ac12SSebastian Grimberg for (CeedInt i = 0; i < n - q - 1; i++) { 732ba59ac12SSebastian Grimberg if (fabs(mat_T[i + n * (i + 1)]) < tol) p += 1; 733ba59ac12SSebastian Grimberg else break; 734ba59ac12SSebastian Grimberg } 735ba59ac12SSebastian Grimberg if (q == n - 1) break; // Finished reducing 736ba59ac12SSebastian Grimberg 737ba59ac12SSebastian Grimberg // Reduce tridiagonal portion 738ba59ac12SSebastian Grimberg CeedScalar t_nn = mat_T[(n - 1 - q) + n * (n - 1 - q)], t_nnm1 = mat_T[(n - 2 - q) + n * (n - 1 - q)]; 739ba59ac12SSebastian Grimberg CeedScalar d = (mat_T[(n - 2 - q) + n * (n - 2 - q)] - t_nn) / 2; 740ba59ac12SSebastian Grimberg CeedScalar mu = t_nn - t_nnm1 * t_nnm1 / (d + copysign(sqrt(d * d + t_nnm1 * t_nnm1), d)); 741ba59ac12SSebastian Grimberg CeedScalar x = mat_T[p + n * p] - mu; 742ba59ac12SSebastian Grimberg CeedScalar z = mat_T[p + n * (p + 1)]; 743ba59ac12SSebastian Grimberg for (CeedInt k = p; k < n - q - 1; k++) { 744ba59ac12SSebastian Grimberg // Compute Givens rotation 745ba59ac12SSebastian Grimberg CeedScalar c = 1, s = 0; 746ba59ac12SSebastian Grimberg if (fabs(z) > tol) { 747ba59ac12SSebastian Grimberg if (fabs(z) > fabs(x)) { 748ba59ac12SSebastian Grimberg CeedScalar tau = -x / z; 749ba59ac12SSebastian Grimberg s = 1 / sqrt(1 + tau * tau), c = s * tau; 750ba59ac12SSebastian Grimberg } else { 751ba59ac12SSebastian Grimberg CeedScalar tau = -z / x; 752ba59ac12SSebastian Grimberg c = 1 / sqrt(1 + tau * tau), s = c * tau; 753ba59ac12SSebastian Grimberg } 754ba59ac12SSebastian Grimberg } 755ba59ac12SSebastian Grimberg 756ba59ac12SSebastian Grimberg // Apply Givens rotation to T 757ba59ac12SSebastian Grimberg CeedGivensRotation(mat_T, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n); 758ba59ac12SSebastian Grimberg CeedGivensRotation(mat_T, c, s, CEED_TRANSPOSE, k, k + 1, n, n); 759ba59ac12SSebastian Grimberg 760ba59ac12SSebastian Grimberg // Apply Givens rotation to Q 761ba59ac12SSebastian Grimberg CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n); 762ba59ac12SSebastian Grimberg 763ba59ac12SSebastian Grimberg // Update x, z 764ba59ac12SSebastian Grimberg if (k < n - q - 2) { 765ba59ac12SSebastian Grimberg x = mat_T[k + n * (k + 1)]; 766ba59ac12SSebastian Grimberg z = mat_T[k + n * (k + 2)]; 767ba59ac12SSebastian Grimberg } 768ba59ac12SSebastian Grimberg } 769ba59ac12SSebastian Grimberg itr++; 770ba59ac12SSebastian Grimberg } 771ba59ac12SSebastian Grimberg 772ba59ac12SSebastian Grimberg // Save eigenvalues 773ba59ac12SSebastian Grimberg for (CeedInt i = 0; i < n; i++) lambda[i] = mat_T[i + n * i]; 774ba59ac12SSebastian Grimberg 775ba59ac12SSebastian Grimberg // Check convergence 776*6574a04fSJeremy L Thompson CeedCheck(itr < max_itr || q > n, ceed, CEED_ERROR_MINOR, "Symmetric QR failed to converge"); 777ba59ac12SSebastian Grimberg return CEED_ERROR_SUCCESS; 778ba59ac12SSebastian Grimberg } 7792c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOn 780ba59ac12SSebastian Grimberg 781ba59ac12SSebastian Grimberg /** 782ba59ac12SSebastian Grimberg @brief Return Simultaneous Diagonalization of two matrices. 783ba59ac12SSebastian Grimberg 784ba59ac12SSebastian Grimberg This solves the generalized eigenvalue problem A x = lambda B x, where A and B are symmetric and B is positive definite. 785ba59ac12SSebastian Grimberg We generate the matrix X and vector Lambda such that X^T A X = Lambda and X^T B X = I. 786ba59ac12SSebastian Grimberg This is equivalent to the LAPACK routine 'sygv' with TYPE = 1. 787ba59ac12SSebastian Grimberg 788ba59ac12SSebastian Grimberg @param[in] ceed Ceed context for error handling 789ba59ac12SSebastian Grimberg @param[in] mat_A Row-major matrix to be factorized with eigenvalues 790ba59ac12SSebastian Grimberg @param[in] mat_B Row-major matrix to be factorized to identity 791ba59ac12SSebastian Grimberg @param[out] mat_X Row-major orthogonal matrix 792ba59ac12SSebastian Grimberg @param[out] lambda Vector of length n of generalized eigenvalues 793ba59ac12SSebastian Grimberg @param[in] n Number of rows/columns 794ba59ac12SSebastian Grimberg 795ba59ac12SSebastian Grimberg @return An error code: 0 - success, otherwise - failure 796ba59ac12SSebastian Grimberg 797ba59ac12SSebastian Grimberg @ref Utility 798ba59ac12SSebastian Grimberg **/ 7992c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOff 8002c2ea1dbSJeremy L Thompson int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *mat_A, CeedScalar *mat_B, CeedScalar *mat_X, CeedScalar *lambda, CeedInt n) { 801ba59ac12SSebastian Grimberg CeedScalar *mat_C, *mat_G, *vec_D; 802ba59ac12SSebastian Grimberg CeedCall(CeedCalloc(n * n, &mat_C)); 803ba59ac12SSebastian Grimberg CeedCall(CeedCalloc(n * n, &mat_G)); 804ba59ac12SSebastian Grimberg CeedCall(CeedCalloc(n, &vec_D)); 805ba59ac12SSebastian Grimberg 806ba59ac12SSebastian Grimberg // Compute B = G D G^T 807ba59ac12SSebastian Grimberg memcpy(mat_G, mat_B, n * n * sizeof(mat_B[0])); 808ba59ac12SSebastian Grimberg CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_G, vec_D, n)); 809ba59ac12SSebastian Grimberg 810ba59ac12SSebastian Grimberg // Sort eigenvalues 811ba59ac12SSebastian Grimberg for (CeedInt i = n - 1; i >= 0; i--) { 812ba59ac12SSebastian Grimberg for (CeedInt j = 0; j < i; j++) { 813ba59ac12SSebastian Grimberg if (fabs(vec_D[j]) > fabs(vec_D[j + 1])) { 814ba59ac12SSebastian Grimberg CeedScalar temp; 815ba59ac12SSebastian Grimberg temp = vec_D[j]; 816ba59ac12SSebastian Grimberg vec_D[j] = vec_D[j + 1]; 817ba59ac12SSebastian Grimberg vec_D[j + 1] = temp; 818ba59ac12SSebastian Grimberg for (CeedInt k = 0; k < n; k++) { 819ba59ac12SSebastian Grimberg temp = mat_G[k * n + j]; 820ba59ac12SSebastian Grimberg mat_G[k * n + j] = mat_G[k * n + j + 1]; 821ba59ac12SSebastian Grimberg mat_G[k * n + j + 1] = temp; 822ba59ac12SSebastian Grimberg } 823ba59ac12SSebastian Grimberg } 824ba59ac12SSebastian Grimberg } 825ba59ac12SSebastian Grimberg } 826ba59ac12SSebastian Grimberg 827ba59ac12SSebastian Grimberg // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T 828ba59ac12SSebastian Grimberg // = D^-1/2 G^T A G D^-1/2 829ba59ac12SSebastian Grimberg // -- D = D^-1/2 830ba59ac12SSebastian Grimberg for (CeedInt i = 0; i < n; i++) vec_D[i] = 1. / sqrt(vec_D[i]); 831ba59ac12SSebastian Grimberg // -- G = G D^-1/2 832ba59ac12SSebastian Grimberg // -- C = D^-1/2 G^T 833ba59ac12SSebastian Grimberg for (CeedInt i = 0; i < n; i++) { 834ba59ac12SSebastian Grimberg for (CeedInt j = 0; j < n; j++) { 835ba59ac12SSebastian Grimberg mat_G[i * n + j] *= vec_D[j]; 836ba59ac12SSebastian Grimberg mat_C[j * n + i] = mat_G[i * n + j]; 837ba59ac12SSebastian Grimberg } 838ba59ac12SSebastian Grimberg } 839ba59ac12SSebastian Grimberg // -- X = (D^-1/2 G^T) A 840ba59ac12SSebastian Grimberg CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_C, (const CeedScalar *)mat_A, mat_X, n, n, n)); 841ba59ac12SSebastian Grimberg // -- C = (D^-1/2 G^T A) (G D^-1/2) 842ba59ac12SSebastian Grimberg CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_X, (const CeedScalar *)mat_G, mat_C, n, n, n)); 843ba59ac12SSebastian Grimberg 844ba59ac12SSebastian Grimberg // Compute Q^T C Q = lambda 845ba59ac12SSebastian Grimberg CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_C, lambda, n)); 846ba59ac12SSebastian Grimberg 847ba59ac12SSebastian Grimberg // Sort eigenvalues 848ba59ac12SSebastian Grimberg for (CeedInt i = n - 1; i >= 0; i--) { 849ba59ac12SSebastian Grimberg for (CeedInt j = 0; j < i; j++) { 850ba59ac12SSebastian Grimberg if (fabs(lambda[j]) > fabs(lambda[j + 1])) { 851ba59ac12SSebastian Grimberg CeedScalar temp; 852ba59ac12SSebastian Grimberg temp = lambda[j]; 853ba59ac12SSebastian Grimberg lambda[j] = lambda[j + 1]; 854ba59ac12SSebastian Grimberg lambda[j + 1] = temp; 855ba59ac12SSebastian Grimberg for (CeedInt k = 0; k < n; k++) { 856ba59ac12SSebastian Grimberg temp = mat_C[k * n + j]; 857ba59ac12SSebastian Grimberg mat_C[k * n + j] = mat_C[k * n + j + 1]; 858ba59ac12SSebastian Grimberg mat_C[k * n + j + 1] = temp; 859ba59ac12SSebastian Grimberg } 860ba59ac12SSebastian Grimberg } 861ba59ac12SSebastian Grimberg } 862ba59ac12SSebastian Grimberg } 863ba59ac12SSebastian Grimberg 864ba59ac12SSebastian Grimberg // Set X = (G D^1/2)^-T Q 865ba59ac12SSebastian Grimberg // = G D^-1/2 Q 866ba59ac12SSebastian Grimberg CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_G, (const CeedScalar *)mat_C, mat_X, n, n, n)); 867ba59ac12SSebastian Grimberg 868ba59ac12SSebastian Grimberg // Cleanup 869ba59ac12SSebastian Grimberg CeedCall(CeedFree(&mat_C)); 870ba59ac12SSebastian Grimberg CeedCall(CeedFree(&mat_G)); 871ba59ac12SSebastian Grimberg CeedCall(CeedFree(&vec_D)); 872ba59ac12SSebastian Grimberg return CEED_ERROR_SUCCESS; 873ba59ac12SSebastian Grimberg } 8742c2ea1dbSJeremy L Thompson CeedPragmaOptimizeOn 875ba59ac12SSebastian Grimberg 8767a982d89SJeremy L. Thompson /// @} 8777a982d89SJeremy L. Thompson 8787a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 8797a982d89SJeremy L. Thompson /// CeedBasis Public API 8807a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 8817a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser 882d7b241e6Sjeremylt /// @{ 883d7b241e6Sjeremylt 884b11c1e72Sjeremylt /** 885ba59ac12SSebastian Grimberg @brief Create a tensor-product basis for H^1 discretizations 886b11c1e72Sjeremylt 887ea61e9acSJeremy L Thompson @param[in] ceed Ceed object where the CeedBasis will be created 888ea61e9acSJeremy L Thompson @param[in] dim Topological dimension 889ea61e9acSJeremy L Thompson @param[in] num_comp Number of field components (1 for scalar fields) 890ea61e9acSJeremy L Thompson @param[in] P_1d Number of nodes in one dimension 891ea61e9acSJeremy L Thompson @param[in] Q_1d Number of quadrature points in one dimension 892ea61e9acSJeremy L Thompson @param[in] interp_1d Row-major (Q_1d * P_1d) matrix expressing the values of nodal basis functions at quadrature points 893ea61e9acSJeremy L Thompson @param[in] grad_1d Row-major (Q_1d * P_1d) matrix expressing derivatives of nodal basis functions at quadrature points 894ea61e9acSJeremy 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] 895ea61e9acSJeremy L Thompson @param[in] q_weight_1d Array of length Q_1d holding the quadrature weights on the reference element 896ea61e9acSJeremy L Thompson @param[out] basis Address of the variable where the newly created CeedBasis will be stored. 897b11c1e72Sjeremylt 898b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 899dfdf5a53Sjeremylt 9007a982d89SJeremy L. Thompson @ref User 901b11c1e72Sjeremylt **/ 9022b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, 9032b730f8bSJeremy L Thompson const CeedScalar *grad_1d, const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis *basis) { 9045fe0d4faSjeremylt if (!ceed->BasisCreateTensorH1) { 9055fe0d4faSjeremylt Ceed delegate; 906*6574a04fSJeremy L Thompson 9072b730f8bSJeremy L Thompson CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis")); 908*6574a04fSJeremy L Thompson CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisCreateTensorH1"); 9092b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateTensorH1(delegate, dim, num_comp, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis)); 910e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 9115fe0d4faSjeremylt } 912e15f9bd0SJeremy L Thompson 913*6574a04fSJeremy L Thompson CeedCheck(dim > 0, ceed, CEED_ERROR_DIMENSION, "Basis dimension must be a positive value"); 914*6574a04fSJeremy L Thompson CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component"); 915*6574a04fSJeremy L Thompson CeedCheck(P_1d > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node"); 916*6574a04fSJeremy L Thompson CeedCheck(Q_1d > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point"); 917227444bfSJeremy L Thompson 9182b730f8bSJeremy L Thompson CeedElemTopology topo = dim == 1 ? CEED_TOPOLOGY_LINE : dim == 2 ? CEED_TOPOLOGY_QUAD : CEED_TOPOLOGY_HEX; 919e15f9bd0SJeremy L Thompson 9202b730f8bSJeremy L Thompson CeedCall(CeedCalloc(1, basis)); 921d7b241e6Sjeremylt (*basis)->ceed = ceed; 9222b730f8bSJeremy L Thompson CeedCall(CeedReference(ceed)); 923d1d35e2fSjeremylt (*basis)->ref_count = 1; 924d1d35e2fSjeremylt (*basis)->tensor_basis = 1; 925d7b241e6Sjeremylt (*basis)->dim = dim; 926d99fa3c5SJeremy L Thompson (*basis)->topo = topo; 927d1d35e2fSjeremylt (*basis)->num_comp = num_comp; 928d1d35e2fSjeremylt (*basis)->P_1d = P_1d; 929d1d35e2fSjeremylt (*basis)->Q_1d = Q_1d; 930d1d35e2fSjeremylt (*basis)->P = CeedIntPow(P_1d, dim); 931d1d35e2fSjeremylt (*basis)->Q = CeedIntPow(Q_1d, dim); 932c4e3f59bSSebastian Grimberg (*basis)->fe_space = CEED_FE_SPACE_H1; 9332b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q_1d, &(*basis)->q_ref_1d)); 9342b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q_1d, &(*basis)->q_weight_1d)); 935ff3a0f91SJeremy L Thompson if (q_ref_1d) memcpy((*basis)->q_ref_1d, q_ref_1d, Q_1d * sizeof(q_ref_1d[0])); 9362b730f8bSJeremy L Thompson if (q_weight_1d) memcpy((*basis)->q_weight_1d, q_weight_1d, Q_1d * sizeof(q_weight_1d[0])); 9372b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->interp_1d)); 9382b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->grad_1d)); 9392b730f8bSJeremy L Thompson if (interp_1d) memcpy((*basis)->interp_1d, interp_1d, Q_1d * P_1d * sizeof(interp_1d[0])); 940ff3a0f91SJeremy L Thompson if (grad_1d) memcpy((*basis)->grad_1d, grad_1d, Q_1d * P_1d * sizeof(grad_1d[0])); 9412b730f8bSJeremy L Thompson CeedCall(ceed->BasisCreateTensorH1(dim, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, *basis)); 942e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 943d7b241e6Sjeremylt } 944d7b241e6Sjeremylt 945b11c1e72Sjeremylt /** 94695bb1877Svaleriabarra @brief Create a tensor-product Lagrange basis 947b11c1e72Sjeremylt 948ea61e9acSJeremy L Thompson @param[in] ceed Ceed object where the CeedBasis will be created 949ea61e9acSJeremy L Thompson @param[in] dim Topological dimension of element 950ea61e9acSJeremy L Thompson @param[in] num_comp Number of field components (1 for scalar fields) 951ea61e9acSJeremy L Thompson @param[in] P Number of Gauss-Lobatto nodes in one dimension. 952ea61e9acSJeremy L Thompson The polynomial degree of the resulting Q_k element is k=P-1. 953ea61e9acSJeremy L Thompson @param[in] Q Number of quadrature points in one dimension. 954ea61e9acSJeremy L Thompson @param[in] quad_mode Distribution of the Q quadrature points (affects order of accuracy for the quadrature) 955ea61e9acSJeremy L Thompson @param[out] basis Address of the variable where the newly created CeedBasis will be stored. 956b11c1e72Sjeremylt 957b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 958dfdf5a53Sjeremylt 9597a982d89SJeremy L. Thompson @ref User 960b11c1e72Sjeremylt **/ 9612b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P, CeedInt Q, CeedQuadMode quad_mode, CeedBasis *basis) { 962d7b241e6Sjeremylt // Allocate 9632b730f8bSJeremy L Thompson int ierr = CEED_ERROR_SUCCESS, i, j, k; 9642b730f8bSJeremy L Thompson CeedScalar c1, c2, c3, c4, dx, *nodes, *interp_1d, *grad_1d, *q_ref_1d, *q_weight_1d; 9654d537eeaSYohann 966*6574a04fSJeremy L Thompson CeedCheck(dim > 0, ceed, CEED_ERROR_DIMENSION, "Basis dimension must be a positive value"); 967*6574a04fSJeremy L Thompson CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component"); 968*6574a04fSJeremy L Thompson CeedCheck(P > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node"); 969*6574a04fSJeremy L Thompson CeedCheck(Q > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point"); 970227444bfSJeremy L Thompson 971e15f9bd0SJeremy L Thompson // Get Nodes and Weights 9722b730f8bSJeremy L Thompson CeedCall(CeedCalloc(P * Q, &interp_1d)); 9732b730f8bSJeremy L Thompson CeedCall(CeedCalloc(P * Q, &grad_1d)); 9742b730f8bSJeremy L Thompson CeedCall(CeedCalloc(P, &nodes)); 9752b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q, &q_ref_1d)); 9762b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q, &q_weight_1d)); 9772b730f8bSJeremy L Thompson if (CeedLobattoQuadrature(P, nodes, NULL) != CEED_ERROR_SUCCESS) goto cleanup; 978d1d35e2fSjeremylt switch (quad_mode) { 979d7b241e6Sjeremylt case CEED_GAUSS: 980d1d35e2fSjeremylt ierr = CeedGaussQuadrature(Q, q_ref_1d, q_weight_1d); 981d7b241e6Sjeremylt break; 982d7b241e6Sjeremylt case CEED_GAUSS_LOBATTO: 983d1d35e2fSjeremylt ierr = CeedLobattoQuadrature(Q, q_ref_1d, q_weight_1d); 984d7b241e6Sjeremylt break; 985d7b241e6Sjeremylt } 9862b730f8bSJeremy L Thompson if (ierr != CEED_ERROR_SUCCESS) goto cleanup; 987e15f9bd0SJeremy L Thompson 988d7b241e6Sjeremylt // Build B, D matrix 989d7b241e6Sjeremylt // Fornberg, 1998 990d7b241e6Sjeremylt for (i = 0; i < Q; i++) { 991d7b241e6Sjeremylt c1 = 1.0; 992d1d35e2fSjeremylt c3 = nodes[0] - q_ref_1d[i]; 993d1d35e2fSjeremylt interp_1d[i * P + 0] = 1.0; 994d7b241e6Sjeremylt for (j = 1; j < P; j++) { 995d7b241e6Sjeremylt c2 = 1.0; 996d7b241e6Sjeremylt c4 = c3; 997d1d35e2fSjeremylt c3 = nodes[j] - q_ref_1d[i]; 998d7b241e6Sjeremylt for (k = 0; k < j; k++) { 999d7b241e6Sjeremylt dx = nodes[j] - nodes[k]; 1000d7b241e6Sjeremylt c2 *= dx; 1001d7b241e6Sjeremylt if (k == j - 1) { 1002d1d35e2fSjeremylt grad_1d[i * P + j] = c1 * (interp_1d[i * P + k] - c4 * grad_1d[i * P + k]) / c2; 1003d1d35e2fSjeremylt interp_1d[i * P + j] = -c1 * c4 * interp_1d[i * P + k] / c2; 1004d7b241e6Sjeremylt } 1005d1d35e2fSjeremylt grad_1d[i * P + k] = (c3 * grad_1d[i * P + k] - interp_1d[i * P + k]) / dx; 1006d1d35e2fSjeremylt interp_1d[i * P + k] = c3 * interp_1d[i * P + k] / dx; 1007d7b241e6Sjeremylt } 1008d7b241e6Sjeremylt c1 = c2; 1009d7b241e6Sjeremylt } 1010d7b241e6Sjeremylt } 10119ac7b42eSJeremy L Thompson // Pass to CeedBasisCreateTensorH1 10122b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P, Q, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis)); 1013e15f9bd0SJeremy L Thompson cleanup: 10142b730f8bSJeremy L Thompson CeedCall(CeedFree(&interp_1d)); 10152b730f8bSJeremy L Thompson CeedCall(CeedFree(&grad_1d)); 10162b730f8bSJeremy L Thompson CeedCall(CeedFree(&nodes)); 10172b730f8bSJeremy L Thompson CeedCall(CeedFree(&q_ref_1d)); 10182b730f8bSJeremy L Thompson CeedCall(CeedFree(&q_weight_1d)); 1019e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1020d7b241e6Sjeremylt } 1021d7b241e6Sjeremylt 1022b11c1e72Sjeremylt /** 1023ba59ac12SSebastian Grimberg @brief Create a non tensor-product basis for H^1 discretizations 1024a8de75f0Sjeremylt 1025ea61e9acSJeremy L Thompson @param[in] ceed Ceed object where the CeedBasis will be created 1026ea61e9acSJeremy L Thompson @param[in] topo Topology of element, e.g. hypercube, simplex, ect 1027ea61e9acSJeremy L Thompson @param[in] num_comp Number of field components (1 for scalar fields) 1028ea61e9acSJeremy L Thompson @param[in] num_nodes Total number of nodes 1029ea61e9acSJeremy L Thompson @param[in] num_qpts Total number of quadrature points 1030ea61e9acSJeremy L Thompson @param[in] interp Row-major (num_qpts * num_nodes) matrix expressing the values of nodal basis functions at quadrature points 1031c4e3f59bSSebastian Grimberg @param[in] grad Row-major (dim * num_qpts * num_nodes) matrix expressing derivatives of nodal basis functions at quadrature points 10329fe083eeSJeremy L Thompson @param[in] q_ref Array of length num_qpts * dim holding the locations of quadrature points on the reference element 1033ea61e9acSJeremy L Thompson @param[in] q_weight Array of length num_qpts holding the quadrature weights on the reference element 1034ea61e9acSJeremy L Thompson @param[out] basis Address of the variable where the newly created CeedBasis will be stored. 1035a8de75f0Sjeremylt 1036a8de75f0Sjeremylt @return An error code: 0 - success, otherwise - failure 1037a8de75f0Sjeremylt 10387a982d89SJeremy L. Thompson @ref User 1039a8de75f0Sjeremylt **/ 10402b730f8bSJeremy L Thompson int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 10412b730f8bSJeremy L Thompson const CeedScalar *grad, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) { 1042d1d35e2fSjeremylt CeedInt P = num_nodes, Q = num_qpts, dim = 0; 1043a8de75f0Sjeremylt 10445fe0d4faSjeremylt if (!ceed->BasisCreateH1) { 10455fe0d4faSjeremylt Ceed delegate; 1046*6574a04fSJeremy L Thompson 10472b730f8bSJeremy L Thompson CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis")); 1048*6574a04fSJeremy L Thompson CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisCreateH1"); 10492b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateH1(delegate, topo, num_comp, num_nodes, num_qpts, interp, grad, q_ref, q_weight, basis)); 1050e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 10515fe0d4faSjeremylt } 10525fe0d4faSjeremylt 1053*6574a04fSJeremy L Thompson CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component"); 1054*6574a04fSJeremy L Thompson CeedCheck(num_nodes > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node"); 1055*6574a04fSJeremy L Thompson CeedCheck(num_qpts > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point"); 1056227444bfSJeremy L Thompson 10572b730f8bSJeremy L Thompson CeedCall(CeedCalloc(1, basis)); 1058a8de75f0Sjeremylt 10592b730f8bSJeremy L Thompson CeedCall(CeedBasisGetTopologyDimension(topo, &dim)); 1060a8de75f0Sjeremylt 1061a8de75f0Sjeremylt (*basis)->ceed = ceed; 10622b730f8bSJeremy L Thompson CeedCall(CeedReference(ceed)); 1063d1d35e2fSjeremylt (*basis)->ref_count = 1; 1064d1d35e2fSjeremylt (*basis)->tensor_basis = 0; 1065a8de75f0Sjeremylt (*basis)->dim = dim; 1066d99fa3c5SJeremy L Thompson (*basis)->topo = topo; 1067d1d35e2fSjeremylt (*basis)->num_comp = num_comp; 1068a8de75f0Sjeremylt (*basis)->P = P; 1069a8de75f0Sjeremylt (*basis)->Q = Q; 1070c4e3f59bSSebastian Grimberg (*basis)->fe_space = CEED_FE_SPACE_H1; 10712b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q * dim, &(*basis)->q_ref_1d)); 10722b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q, &(*basis)->q_weight_1d)); 1073ff3a0f91SJeremy L Thompson if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0])); 1074ff3a0f91SJeremy L Thompson if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0])); 10752b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q * P, &(*basis)->interp)); 10762b730f8bSJeremy L Thompson CeedCall(CeedCalloc(dim * Q * P, &(*basis)->grad)); 1077ff3a0f91SJeremy L Thompson if (interp) memcpy((*basis)->interp, interp, Q * P * sizeof(interp[0])); 1078ff3a0f91SJeremy L Thompson if (grad) memcpy((*basis)->grad, grad, dim * Q * P * sizeof(grad[0])); 10792b730f8bSJeremy L Thompson CeedCall(ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, q_ref, q_weight, *basis)); 1080e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1081a8de75f0Sjeremylt } 1082a8de75f0Sjeremylt 1083a8de75f0Sjeremylt /** 1084859c15bbSJames Wright @brief Create a non tensor-product basis for \f$H(\mathrm{div})\f$ discretizations 108550c301a5SRezgar Shakeri 1086ea61e9acSJeremy L Thompson @param[in] ceed Ceed object where the CeedBasis will be created 1087ea61e9acSJeremy 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 1088ea61e9acSJeremy L Thompson @param[in] num_comp Number of components (usually 1 for vectors in H(div) bases) 1089ea61e9acSJeremy L Thompson @param[in] num_nodes Total number of nodes (dofs per element) 1090ea61e9acSJeremy L Thompson @param[in] num_qpts Total number of quadrature points 1091c4e3f59bSSebastian Grimberg @param[in] interp Row-major (dim * num_qpts * num_nodes) matrix expressing the values of basis functions at quadrature points 1092c4e3f59bSSebastian Grimberg @param[in] div Row-major (num_qpts * num_nodes) matrix expressing divergence of basis functions at quadrature points 10939fe083eeSJeremy L Thompson @param[in] q_ref Array of length num_qpts * dim holding the locations of quadrature points on the reference element 1094ea61e9acSJeremy L Thompson @param[in] q_weight Array of length num_qpts holding the quadrature weights on the reference element 1095ea61e9acSJeremy L Thompson @param[out] basis Address of the variable where the newly created CeedBasis will be stored. 109650c301a5SRezgar Shakeri 109750c301a5SRezgar Shakeri @return An error code: 0 - success, otherwise - failure 109850c301a5SRezgar Shakeri 109950c301a5SRezgar Shakeri @ref User 110050c301a5SRezgar Shakeri **/ 11012b730f8bSJeremy L Thompson int CeedBasisCreateHdiv(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 11022b730f8bSJeremy L Thompson const CeedScalar *div, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) { 110350c301a5SRezgar Shakeri CeedInt Q = num_qpts, P = num_nodes, dim = 0; 1104c4e3f59bSSebastian Grimberg 110550c301a5SRezgar Shakeri if (!ceed->BasisCreateHdiv) { 110650c301a5SRezgar Shakeri Ceed delegate; 1107*6574a04fSJeremy L Thompson 11082b730f8bSJeremy L Thompson CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis")); 1109*6574a04fSJeremy L Thompson CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateHdiv"); 11102b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateHdiv(delegate, topo, num_comp, num_nodes, num_qpts, interp, div, q_ref, q_weight, basis)); 111150c301a5SRezgar Shakeri return CEED_ERROR_SUCCESS; 111250c301a5SRezgar Shakeri } 111350c301a5SRezgar Shakeri 1114*6574a04fSJeremy L Thompson CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component"); 1115*6574a04fSJeremy L Thompson CeedCheck(num_nodes > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node"); 1116*6574a04fSJeremy L Thompson CeedCheck(num_qpts > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point"); 1117227444bfSJeremy L Thompson 11182b730f8bSJeremy L Thompson CeedCall(CeedCalloc(1, basis)); 111950c301a5SRezgar Shakeri 1120c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetTopologyDimension(topo, &dim)); 1121c4e3f59bSSebastian Grimberg 112250c301a5SRezgar Shakeri (*basis)->ceed = ceed; 11232b730f8bSJeremy L Thompson CeedCall(CeedReference(ceed)); 112450c301a5SRezgar Shakeri (*basis)->ref_count = 1; 112550c301a5SRezgar Shakeri (*basis)->tensor_basis = 0; 112650c301a5SRezgar Shakeri (*basis)->dim = dim; 112750c301a5SRezgar Shakeri (*basis)->topo = topo; 112850c301a5SRezgar Shakeri (*basis)->num_comp = num_comp; 112950c301a5SRezgar Shakeri (*basis)->P = P; 113050c301a5SRezgar Shakeri (*basis)->Q = Q; 1131c4e3f59bSSebastian Grimberg (*basis)->fe_space = CEED_FE_SPACE_HDIV; 11322b730f8bSJeremy L Thompson CeedCall(CeedMalloc(Q * dim, &(*basis)->q_ref_1d)); 11332b730f8bSJeremy L Thompson CeedCall(CeedMalloc(Q, &(*basis)->q_weight_1d)); 113450c301a5SRezgar Shakeri if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0])); 113550c301a5SRezgar Shakeri if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0])); 11362b730f8bSJeremy L Thompson CeedCall(CeedMalloc(dim * Q * P, &(*basis)->interp)); 11372b730f8bSJeremy L Thompson CeedCall(CeedMalloc(Q * P, &(*basis)->div)); 113850c301a5SRezgar Shakeri if (interp) memcpy((*basis)->interp, interp, dim * Q * P * sizeof(interp[0])); 113950c301a5SRezgar Shakeri if (div) memcpy((*basis)->div, div, Q * P * sizeof(div[0])); 11402b730f8bSJeremy L Thompson CeedCall(ceed->BasisCreateHdiv(topo, dim, P, Q, interp, div, q_ref, q_weight, *basis)); 114150c301a5SRezgar Shakeri return CEED_ERROR_SUCCESS; 114250c301a5SRezgar Shakeri } 114350c301a5SRezgar Shakeri 114450c301a5SRezgar Shakeri /** 1145c4e3f59bSSebastian Grimberg @brief Create a non tensor-product basis for H(curl) discretizations 1146c4e3f59bSSebastian Grimberg 1147c4e3f59bSSebastian Grimberg @param[in] ceed Ceed object where the CeedBasis will be created 1148c4e3f59bSSebastian Grimberg @param[in] topo Topology of element (`CEED_TOPOLOGY_QUAD`, `CEED_TOPOLOGY_PRISM`, etc.), dimension of which is used in some array sizes below 1149c4e3f59bSSebastian Grimberg @param[in] num_comp Number of components (usually 1 for vectors in H(curl) bases) 1150c4e3f59bSSebastian Grimberg @param[in] num_nodes Total number of nodes (dofs per element) 1151c4e3f59bSSebastian Grimberg @param[in] num_qpts Total number of quadrature points 1152c4e3f59bSSebastian Grimberg @param[in] interp Row-major (dim * num_qpts * num_nodes) matrix expressing the values of basis functions at quadrature points 1153c4e3f59bSSebastian 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 1154c4e3f59bSSebastian Grimberg quadrature points 1155c4e3f59bSSebastian Grimberg @param[in] q_ref Array of length num_qpts * dim holding the locations of quadrature points on the reference element 1156c4e3f59bSSebastian Grimberg @param[in] q_weight Array of length num_qpts holding the quadrature weights on the reference element 1157c4e3f59bSSebastian Grimberg @param[out] basis Address of the variable where the newly created CeedBasis will be stored. 1158c4e3f59bSSebastian Grimberg 1159c4e3f59bSSebastian Grimberg @return An error code: 0 - success, otherwise - failure 1160c4e3f59bSSebastian Grimberg 1161c4e3f59bSSebastian Grimberg @ref User 1162c4e3f59bSSebastian Grimberg **/ 1163c4e3f59bSSebastian Grimberg int CeedBasisCreateHcurl(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 1164c4e3f59bSSebastian Grimberg const CeedScalar *curl, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) { 1165c4e3f59bSSebastian Grimberg CeedInt Q = num_qpts, P = num_nodes, dim = 0, curl_comp = 0; 1166c4e3f59bSSebastian Grimberg 1167c4e3f59bSSebastian Grimberg if (!ceed->BasisCreateHdiv) { 1168c4e3f59bSSebastian Grimberg Ceed delegate; 1169*6574a04fSJeremy L Thompson 1170c4e3f59bSSebastian Grimberg CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis")); 1171*6574a04fSJeremy L Thompson CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateHcurl"); 1172c4e3f59bSSebastian Grimberg CeedCall(CeedBasisCreateHcurl(delegate, topo, num_comp, num_nodes, num_qpts, interp, curl, q_ref, q_weight, basis)); 1173c4e3f59bSSebastian Grimberg return CEED_ERROR_SUCCESS; 1174c4e3f59bSSebastian Grimberg } 1175c4e3f59bSSebastian Grimberg 1176*6574a04fSJeremy L Thompson CeedCheck(num_comp > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component"); 1177*6574a04fSJeremy L Thompson CeedCheck(num_nodes > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node"); 1178*6574a04fSJeremy L Thompson CeedCheck(num_qpts > 0, ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point"); 1179c4e3f59bSSebastian Grimberg 1180c4e3f59bSSebastian Grimberg CeedCall(CeedCalloc(1, basis)); 1181c4e3f59bSSebastian Grimberg 1182c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetTopologyDimension(topo, &dim)); 1183c4e3f59bSSebastian Grimberg curl_comp = (dim < 3) ? 1 : dim; 1184c4e3f59bSSebastian Grimberg 1185c4e3f59bSSebastian Grimberg (*basis)->ceed = ceed; 1186c4e3f59bSSebastian Grimberg CeedCall(CeedReference(ceed)); 1187c4e3f59bSSebastian Grimberg (*basis)->ref_count = 1; 1188c4e3f59bSSebastian Grimberg (*basis)->tensor_basis = 0; 1189c4e3f59bSSebastian Grimberg (*basis)->dim = dim; 1190c4e3f59bSSebastian Grimberg (*basis)->topo = topo; 1191c4e3f59bSSebastian Grimberg (*basis)->num_comp = num_comp; 1192c4e3f59bSSebastian Grimberg (*basis)->P = P; 1193c4e3f59bSSebastian Grimberg (*basis)->Q = Q; 1194c4e3f59bSSebastian Grimberg (*basis)->fe_space = CEED_FE_SPACE_HCURL; 1195c4e3f59bSSebastian Grimberg CeedCall(CeedMalloc(Q * dim, &(*basis)->q_ref_1d)); 1196c4e3f59bSSebastian Grimberg CeedCall(CeedMalloc(Q, &(*basis)->q_weight_1d)); 1197c4e3f59bSSebastian Grimberg if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0])); 1198c4e3f59bSSebastian Grimberg if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0])); 1199c4e3f59bSSebastian Grimberg CeedCall(CeedMalloc(dim * Q * P, &(*basis)->interp)); 1200c4e3f59bSSebastian Grimberg CeedCall(CeedMalloc(curl_comp * Q * P, &(*basis)->curl)); 1201c4e3f59bSSebastian Grimberg if (interp) memcpy((*basis)->interp, interp, dim * Q * P * sizeof(interp[0])); 1202c4e3f59bSSebastian Grimberg if (curl) memcpy((*basis)->curl, curl, curl_comp * Q * P * sizeof(curl[0])); 1203c4e3f59bSSebastian Grimberg CeedCall(ceed->BasisCreateHcurl(topo, dim, P, Q, interp, curl, q_ref, q_weight, *basis)); 1204c4e3f59bSSebastian Grimberg return CEED_ERROR_SUCCESS; 1205c4e3f59bSSebastian Grimberg } 1206c4e3f59bSSebastian Grimberg 1207c4e3f59bSSebastian Grimberg /** 1208ea61e9acSJeremy L Thompson @brief Create a CeedBasis for projection from the nodes of `basis_from` to the nodes of `basis_to`. 1209ba59ac12SSebastian Grimberg 121015ad3917SSebastian Grimberg Only `CEED_EVAL_INTERP` will be valid for the new basis, `basis_project`. For H^1 spaces, `CEED_EVAL_GRAD` will also be valid. 1211ea61e9acSJeremy L Thompson The interpolation is given by `interp_project = interp_to^+ * interp_from`, where the pesudoinverse `interp_to^+` is given by QR 121215ad3917SSebastian Grimberg factorization. The gradient (for the H^1 case) is given by `grad_project = interp_to^+ * grad_from`. 121315ad3917SSebastian Grimberg 121415ad3917SSebastian Grimberg Note: `basis_from` and `basis_to` must have compatible quadrature spaces. 121515ad3917SSebastian Grimberg 1216ba59ac12SSebastian Grimberg Note: `basis_project` will have the same number of components as `basis_from`, regardless of the number of components that `basis_to` has. If 1217ea61e9acSJeremy L Thompson `basis_from` has 3 components and `basis_to` has 5 components, then `basis_project` will have 3 components. 1218f113e5dcSJeremy L Thompson 1219f113e5dcSJeremy L Thompson @param[in] basis_from CeedBasis to prolong from 1220446e7af4SJeremy L Thompson @param[in] basis_to CeedBasis to prolong to 1221ea61e9acSJeremy L Thompson @param[out] basis_project Address of the variable where the newly created CeedBasis will be stored. 1222f113e5dcSJeremy L Thompson 1223f113e5dcSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 1224f113e5dcSJeremy L Thompson 1225f113e5dcSJeremy L Thompson @ref User 1226f113e5dcSJeremy L Thompson **/ 12272b730f8bSJeremy L Thompson int CeedBasisCreateProjection(CeedBasis basis_from, CeedBasis basis_to, CeedBasis *basis_project) { 1228f113e5dcSJeremy L Thompson Ceed ceed; 12292b730f8bSJeremy L Thompson CeedCall(CeedBasisGetCeed(basis_to, &ceed)); 1230f113e5dcSJeremy L Thompson 1231ecc88aebSJeremy L Thompson // Create projection matrix 123214556e63SJeremy L Thompson CeedScalar *interp_project, *grad_project; 12332b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateProjectionMatrices(basis_from, basis_to, &interp_project, &grad_project)); 1234f113e5dcSJeremy L Thompson 1235f113e5dcSJeremy L Thompson // Build basis 1236f113e5dcSJeremy L Thompson bool is_tensor; 123715ad3917SSebastian Grimberg CeedFESpace fe_space; 1238f113e5dcSJeremy L Thompson CeedInt dim, num_comp; 123914556e63SJeremy L Thompson CeedScalar *q_ref, *q_weight; 12402b730f8bSJeremy L Thompson CeedCall(CeedBasisIsTensor(basis_to, &is_tensor)); 124115ad3917SSebastian Grimberg CeedCall(CeedBasisGetFESpace(basis_to, &fe_space)); 12422b730f8bSJeremy L Thompson CeedCall(CeedBasisGetDimension(basis_to, &dim)); 12432b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumComponents(basis_from, &num_comp)); 1244f113e5dcSJeremy L Thompson if (is_tensor) { 1245f113e5dcSJeremy L Thompson CeedInt P_1d_to, P_1d_from; 12462b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_1d_from)); 12472b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_1d_to)); 12482b730f8bSJeremy L Thompson CeedCall(CeedCalloc(P_1d_to, &q_ref)); 12492b730f8bSJeremy L Thompson CeedCall(CeedCalloc(P_1d_to, &q_weight)); 12502b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P_1d_from, P_1d_to, interp_project, grad_project, q_ref, q_weight, basis_project)); 1251f113e5dcSJeremy L Thompson } else { 1252f113e5dcSJeremy L Thompson CeedElemTopology topo; 12532b730f8bSJeremy L Thompson CeedCall(CeedBasisGetTopology(basis_to, &topo)); 1254f113e5dcSJeremy L Thompson CeedInt num_nodes_to, num_nodes_from; 12552b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis_from, &num_nodes_from)); 12562b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis_to, &num_nodes_to)); 12572b730f8bSJeremy L Thompson CeedCall(CeedCalloc(num_nodes_to * dim, &q_ref)); 12582b730f8bSJeremy L Thompson CeedCall(CeedCalloc(num_nodes_to, &q_weight)); 125915ad3917SSebastian Grimberg switch (fe_space) { 126015ad3917SSebastian Grimberg case CEED_FE_SPACE_H1: 12612b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateH1(ceed, topo, num_comp, num_nodes_from, num_nodes_to, interp_project, grad_project, q_ref, q_weight, basis_project)); 126215ad3917SSebastian Grimberg break; 126315ad3917SSebastian Grimberg case CEED_FE_SPACE_HDIV: 126415ad3917SSebastian Grimberg CeedCall( 126515ad3917SSebastian Grimberg CeedBasisCreateHdiv(ceed, topo, num_comp, num_nodes_from, num_nodes_to, interp_project, grad_project, q_ref, q_weight, basis_project)); 126615ad3917SSebastian Grimberg break; 126715ad3917SSebastian Grimberg case CEED_FE_SPACE_HCURL: 126815ad3917SSebastian Grimberg CeedCall( 126915ad3917SSebastian Grimberg CeedBasisCreateHcurl(ceed, topo, num_comp, num_nodes_from, num_nodes_to, interp_project, grad_project, q_ref, q_weight, basis_project)); 127015ad3917SSebastian Grimberg break; 127115ad3917SSebastian Grimberg } 1272f113e5dcSJeremy L Thompson } 1273f113e5dcSJeremy L Thompson 1274f113e5dcSJeremy L Thompson // Cleanup 12752b730f8bSJeremy L Thompson CeedCall(CeedFree(&interp_project)); 12762b730f8bSJeremy L Thompson CeedCall(CeedFree(&grad_project)); 12772b730f8bSJeremy L Thompson CeedCall(CeedFree(&q_ref)); 12782b730f8bSJeremy L Thompson CeedCall(CeedFree(&q_weight)); 1279f113e5dcSJeremy L Thompson 1280f113e5dcSJeremy L Thompson return CEED_ERROR_SUCCESS; 1281f113e5dcSJeremy L Thompson } 1282f113e5dcSJeremy L Thompson 1283f113e5dcSJeremy L Thompson /** 1284ea61e9acSJeremy L Thompson @brief Copy the pointer to a CeedBasis. 12859560d06aSjeremylt 1286512bb800SJeremy 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. 1287512bb800SJeremy L Thompson This CeedBasis will be destroyed if `basis_copy` is the only reference to this CeedBasis. 1288ea61e9acSJeremy L Thompson 1289ea61e9acSJeremy L Thompson @param[in] basis CeedBasis to copy reference to 1290ea61e9acSJeremy L Thompson @param[in,out] basis_copy Variable to store copied reference 12919560d06aSjeremylt 12929560d06aSjeremylt @return An error code: 0 - success, otherwise - failure 12939560d06aSjeremylt 12949560d06aSjeremylt @ref User 12959560d06aSjeremylt **/ 12969560d06aSjeremylt int CeedBasisReferenceCopy(CeedBasis basis, CeedBasis *basis_copy) { 1297393ac2cdSJeremy L Thompson if (basis != CEED_BASIS_COLLOCATED) CeedCall(CeedBasisReference(basis)); 12982b730f8bSJeremy L Thompson CeedCall(CeedBasisDestroy(basis_copy)); 12999560d06aSjeremylt *basis_copy = basis; 13009560d06aSjeremylt return CEED_ERROR_SUCCESS; 13019560d06aSjeremylt } 13029560d06aSjeremylt 13039560d06aSjeremylt /** 13047a982d89SJeremy L. Thompson @brief View a CeedBasis 13057a982d89SJeremy L. Thompson 1306ea61e9acSJeremy L Thompson @param[in] basis CeedBasis to view 1307ea61e9acSJeremy L Thompson @param[in] stream Stream to view to, e.g., stdout 13087a982d89SJeremy L. Thompson 13097a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 13107a982d89SJeremy L. Thompson 13117a982d89SJeremy L. Thompson @ref User 13127a982d89SJeremy L. Thompson **/ 13137a982d89SJeremy L. Thompson int CeedBasisView(CeedBasis basis, FILE *stream) { 131450c301a5SRezgar Shakeri CeedElemTopology topo = basis->topo; 1315c4e3f59bSSebastian Grimberg CeedFESpace fe_space = basis->fe_space; 1316c4e3f59bSSebastian Grimberg CeedInt q_comp = 0; 13172b730f8bSJeremy L Thompson 131850c301a5SRezgar Shakeri // Print FE space and element topology of the basis 1319d1d35e2fSjeremylt if (basis->tensor_basis) { 1320c4e3f59bSSebastian Grimberg fprintf(stream, "CeedBasis (%s on a %s element): dim=%" CeedInt_FMT " P=%" CeedInt_FMT " Q=%" CeedInt_FMT "\n", CeedFESpaces[fe_space], 13212b730f8bSJeremy L Thompson CeedElemTopologies[topo], basis->dim, basis->P_1d, basis->Q_1d); 132250c301a5SRezgar Shakeri } else { 1323c4e3f59bSSebastian Grimberg fprintf(stream, "CeedBasis (%s on a %s element): dim=%" CeedInt_FMT " P=%" CeedInt_FMT " Q=%" CeedInt_FMT "\n", CeedFESpaces[fe_space], 13242b730f8bSJeremy L Thompson CeedElemTopologies[topo], basis->dim, basis->P, basis->Q); 132550c301a5SRezgar Shakeri } 1326ea61e9acSJeremy L Thompson // Print quadrature data, interpolation/gradient/divergence/curl of the basis 132750c301a5SRezgar Shakeri if (basis->tensor_basis) { // tensor basis 13282b730f8bSJeremy L Thompson CeedCall(CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_ref_1d, stream)); 13292b730f8bSJeremy L Thompson CeedCall(CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_weight_1d, stream)); 13302b730f8bSJeremy L Thompson CeedCall(CeedScalarView("interp1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, basis->interp_1d, stream)); 13312b730f8bSJeremy L Thompson CeedCall(CeedScalarView("grad1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, basis->grad_1d, stream)); 133250c301a5SRezgar Shakeri } else { // non-tensor basis 13332b730f8bSJeremy L Thompson CeedCall(CeedScalarView("qref", "\t% 12.8f", 1, basis->Q * basis->dim, basis->q_ref_1d, stream)); 13342b730f8bSJeremy L Thompson CeedCall(CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->q_weight_1d, stream)); 1335c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp)); 1336c4e3f59bSSebastian Grimberg CeedCall(CeedScalarView("interp", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->interp, stream)); 133750c301a5SRezgar Shakeri if (basis->grad) { 1338c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp)); 1339c4e3f59bSSebastian Grimberg CeedCall(CeedScalarView("grad", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->grad, stream)); 13407a982d89SJeremy L. Thompson } 134150c301a5SRezgar Shakeri if (basis->div) { 1342c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_DIV, &q_comp)); 1343c4e3f59bSSebastian Grimberg CeedCall(CeedScalarView("div", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->div, stream)); 1344c4e3f59bSSebastian Grimberg } 1345c4e3f59bSSebastian Grimberg if (basis->curl) { 1346c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_CURL, &q_comp)); 1347c4e3f59bSSebastian Grimberg CeedCall(CeedScalarView("curl", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->curl, stream)); 134850c301a5SRezgar Shakeri } 134950c301a5SRezgar Shakeri } 1350e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 13517a982d89SJeremy L. Thompson } 13527a982d89SJeremy L. Thompson 13537a982d89SJeremy L. Thompson /** 13547a982d89SJeremy L. Thompson @brief Apply basis evaluation from nodes to quadrature points or vice versa 13557a982d89SJeremy L. Thompson 1356ea61e9acSJeremy L Thompson @param[in] basis CeedBasis to evaluate 1357ea61e9acSJeremy L Thompson @param[in] num_elem The number of elements to apply the basis evaluation to; 1358ea61e9acSJeremy L Thompson the backend will specify the ordering in CeedElemRestrictionCreateBlocked() 1359ea61e9acSJeremy L Thompson @param[in] t_mode \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature points; 1360ea61e9acSJeremy L Thompson \ref CEED_TRANSPOSE to apply the transpose, mapping from quadrature points to nodes 1361ea61e9acSJeremy L Thompson @param[in] eval_mode \ref CEED_EVAL_NONE to use values directly, 13627a982d89SJeremy L. Thompson \ref CEED_EVAL_INTERP to use interpolated values, 13637a982d89SJeremy L. Thompson \ref CEED_EVAL_GRAD to use gradients, 1364c4e3f59bSSebastian Grimberg \ref CEED_EVAL_DIV to use divergence, 1365c4e3f59bSSebastian Grimberg \ref CEED_EVAL_CURL to use curl, 13667a982d89SJeremy L. Thompson \ref CEED_EVAL_WEIGHT to use quadrature weights. 13677a982d89SJeremy L. Thompson @param[in] u Input CeedVector 13687a982d89SJeremy L. Thompson @param[out] v Output CeedVector 13697a982d89SJeremy L. Thompson 13707a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 13717a982d89SJeremy L. Thompson 13727a982d89SJeremy L. Thompson @ref User 13737a982d89SJeremy L. Thompson **/ 13742b730f8bSJeremy L Thompson int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, CeedVector v) { 13751f9221feSJeremy L Thompson CeedSize u_length = 0, v_length; 1376c4e3f59bSSebastian Grimberg CeedInt dim, num_comp, q_comp, num_nodes, num_qpts; 13772b730f8bSJeremy L Thompson CeedCall(CeedBasisGetDimension(basis, &dim)); 13782b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumComponents(basis, &num_comp)); 1379c4e3f59bSSebastian Grimberg CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp)); 13802b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis, &num_nodes)); 13812b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts)); 13822b730f8bSJeremy L Thompson CeedCall(CeedVectorGetLength(v, &v_length)); 13837a982d89SJeremy L. Thompson if (u) { 13842b730f8bSJeremy L Thompson CeedCall(CeedVectorGetLength(u, &u_length)); 13857a982d89SJeremy L. Thompson } 13867a982d89SJeremy L. Thompson 1387*6574a04fSJeremy L Thompson CeedCheck(basis->Apply, basis->ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisApply"); 1388e15f9bd0SJeremy L Thompson 1389e15f9bd0SJeremy L Thompson // Check compatibility of topological and geometrical dimensions 1390*6574a04fSJeremy L Thompson CeedCheck((t_mode == CEED_TRANSPOSE && v_length % num_nodes == 0 && u_length % num_qpts == 0) || 1391*6574a04fSJeremy L Thompson (t_mode == CEED_NOTRANSPOSE && u_length % num_nodes == 0 && v_length % num_qpts == 0), 1392*6574a04fSJeremy L Thompson basis->ceed, CEED_ERROR_DIMENSION, "Length of input/output vectors incompatible with basis dimensions"); 13937a982d89SJeremy L. Thompson 1394e15f9bd0SJeremy L Thompson // Check vector lengths to prevent out of bounds issues 1395*6574a04fSJeremy L Thompson bool good_dims = true; 1396d1d35e2fSjeremylt switch (eval_mode) { 1397e15f9bd0SJeremy L Thompson case CEED_EVAL_NONE: 13982b730f8bSJeremy L Thompson case CEED_EVAL_INTERP: 13992b730f8bSJeremy L Thompson case CEED_EVAL_GRAD: 1400c4e3f59bSSebastian Grimberg case CEED_EVAL_DIV: 1401c4e3f59bSSebastian Grimberg case CEED_EVAL_CURL: 1402*6574a04fSJeremy L Thompson good_dims = 1403*6574a04fSJeremy L Thompson ((t_mode == CEED_TRANSPOSE && u_length >= num_elem * num_comp * num_qpts * q_comp && v_length >= num_elem * num_comp * num_nodes) || 1404*6574a04fSJeremy L Thompson (t_mode == CEED_NOTRANSPOSE && v_length >= num_elem * num_qpts * num_comp * q_comp && u_length >= num_elem * num_comp * num_nodes)); 1405e15f9bd0SJeremy L Thompson break; 1406e15f9bd0SJeremy L Thompson case CEED_EVAL_WEIGHT: 1407*6574a04fSJeremy L Thompson good_dims = v_length >= num_elem * num_qpts; 1408e15f9bd0SJeremy L Thompson break; 1409e15f9bd0SJeremy L Thompson } 1410*6574a04fSJeremy L Thompson CeedCheck(good_dims, basis->ceed, CEED_ERROR_DIMENSION, "Input/output vectors too short for basis and evaluation mode"); 1411e15f9bd0SJeremy L Thompson 14122b730f8bSJeremy L Thompson CeedCall(basis->Apply(basis, num_elem, t_mode, eval_mode, u, v)); 1413e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 14147a982d89SJeremy L. Thompson } 14157a982d89SJeremy L. Thompson 14167a982d89SJeremy L. Thompson /** 1417b7c9bbdaSJeremy L Thompson @brief Get Ceed associated with a CeedBasis 1418b7c9bbdaSJeremy L Thompson 1419ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 1420b7c9bbdaSJeremy L Thompson @param[out] ceed Variable to store Ceed 1421b7c9bbdaSJeremy L Thompson 1422b7c9bbdaSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 1423b7c9bbdaSJeremy L Thompson 1424b7c9bbdaSJeremy L Thompson @ref Advanced 1425b7c9bbdaSJeremy L Thompson **/ 1426b7c9bbdaSJeremy L Thompson int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) { 1427b7c9bbdaSJeremy L Thompson *ceed = basis->ceed; 1428b7c9bbdaSJeremy L Thompson return CEED_ERROR_SUCCESS; 1429b7c9bbdaSJeremy L Thompson } 1430b7c9bbdaSJeremy L Thompson 1431b7c9bbdaSJeremy L Thompson /** 14329d007619Sjeremylt @brief Get dimension for given CeedBasis 14339d007619Sjeremylt 1434ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 14359d007619Sjeremylt @param[out] dim Variable to store dimension of basis 14369d007619Sjeremylt 14379d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 14389d007619Sjeremylt 1439b7c9bbdaSJeremy L Thompson @ref Advanced 14409d007619Sjeremylt **/ 14419d007619Sjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) { 14429d007619Sjeremylt *dim = basis->dim; 1443e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 14449d007619Sjeremylt } 14459d007619Sjeremylt 14469d007619Sjeremylt /** 1447d99fa3c5SJeremy L Thompson @brief Get topology for given CeedBasis 1448d99fa3c5SJeremy L Thompson 1449ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 1450d99fa3c5SJeremy L Thompson @param[out] topo Variable to store topology of basis 1451d99fa3c5SJeremy L Thompson 1452d99fa3c5SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 1453d99fa3c5SJeremy L Thompson 1454b7c9bbdaSJeremy L Thompson @ref Advanced 1455d99fa3c5SJeremy L Thompson **/ 1456d99fa3c5SJeremy L Thompson int CeedBasisGetTopology(CeedBasis basis, CeedElemTopology *topo) { 1457d99fa3c5SJeremy L Thompson *topo = basis->topo; 1458e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1459d99fa3c5SJeremy L Thompson } 1460d99fa3c5SJeremy L Thompson 1461d99fa3c5SJeremy L Thompson /** 14629d007619Sjeremylt @brief Get number of components for given CeedBasis 14639d007619Sjeremylt 1464ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 1465d1d35e2fSjeremylt @param[out] num_comp Variable to store number of components of basis 14669d007619Sjeremylt 14679d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 14689d007619Sjeremylt 1469b7c9bbdaSJeremy L Thompson @ref Advanced 14709d007619Sjeremylt **/ 1471d1d35e2fSjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *num_comp) { 1472d1d35e2fSjeremylt *num_comp = basis->num_comp; 1473e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 14749d007619Sjeremylt } 14759d007619Sjeremylt 14769d007619Sjeremylt /** 14779d007619Sjeremylt @brief Get total number of nodes (in dim dimensions) of a CeedBasis 14789d007619Sjeremylt 1479ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 14809d007619Sjeremylt @param[out] P Variable to store number of nodes 14819d007619Sjeremylt 14829d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 14839d007619Sjeremylt 14849d007619Sjeremylt @ref Utility 14859d007619Sjeremylt **/ 14869d007619Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) { 14879d007619Sjeremylt *P = basis->P; 1488e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 14899d007619Sjeremylt } 14909d007619Sjeremylt 14919d007619Sjeremylt /** 14929d007619Sjeremylt @brief Get total number of nodes (in 1 dimension) of a CeedBasis 14939d007619Sjeremylt 1494ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 1495d1d35e2fSjeremylt @param[out] P_1d Variable to store number of nodes 14969d007619Sjeremylt 14979d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 14989d007619Sjeremylt 1499b7c9bbdaSJeremy L Thompson @ref Advanced 15009d007619Sjeremylt **/ 1501d1d35e2fSjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P_1d) { 1502*6574a04fSJeremy L Thompson CeedCheck(basis->tensor_basis, basis->ceed, CEED_ERROR_MINOR, "Cannot supply P_1d for non-tensor basis"); 1503d1d35e2fSjeremylt *P_1d = basis->P_1d; 1504e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 15059d007619Sjeremylt } 15069d007619Sjeremylt 15079d007619Sjeremylt /** 15089d007619Sjeremylt @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis 15099d007619Sjeremylt 1510ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 15119d007619Sjeremylt @param[out] Q Variable to store number of quadrature points 15129d007619Sjeremylt 15139d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 15149d007619Sjeremylt 15159d007619Sjeremylt @ref Utility 15169d007619Sjeremylt **/ 15179d007619Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) { 15189d007619Sjeremylt *Q = basis->Q; 1519e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 15209d007619Sjeremylt } 15219d007619Sjeremylt 15229d007619Sjeremylt /** 15239d007619Sjeremylt @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis 15249d007619Sjeremylt 1525ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 1526d1d35e2fSjeremylt @param[out] Q_1d Variable to store number of quadrature points 15279d007619Sjeremylt 15289d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 15299d007619Sjeremylt 1530b7c9bbdaSJeremy L Thompson @ref Advanced 15319d007619Sjeremylt **/ 1532d1d35e2fSjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q_1d) { 1533*6574a04fSJeremy L Thompson CeedCheck(basis->tensor_basis, basis->ceed, CEED_ERROR_MINOR, "Cannot supply Q_1d for non-tensor basis"); 1534d1d35e2fSjeremylt *Q_1d = basis->Q_1d; 1535e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 15369d007619Sjeremylt } 15379d007619Sjeremylt 15389d007619Sjeremylt /** 1539ea61e9acSJeremy L Thompson @brief Get reference coordinates of quadrature points (in dim dimensions) of a CeedBasis 15409d007619Sjeremylt 1541ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 1542d1d35e2fSjeremylt @param[out] q_ref Variable to store reference coordinates of quadrature points 15439d007619Sjeremylt 15449d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 15459d007619Sjeremylt 1546b7c9bbdaSJeremy L Thompson @ref Advanced 15479d007619Sjeremylt **/ 1548d1d35e2fSjeremylt int CeedBasisGetQRef(CeedBasis basis, const CeedScalar **q_ref) { 1549d1d35e2fSjeremylt *q_ref = basis->q_ref_1d; 1550e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 15519d007619Sjeremylt } 15529d007619Sjeremylt 15539d007619Sjeremylt /** 1554ea61e9acSJeremy L Thompson @brief Get quadrature weights of quadrature points (in dim dimensions) of a CeedBasis 15559d007619Sjeremylt 1556ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 1557d1d35e2fSjeremylt @param[out] q_weight Variable to store quadrature weights 15589d007619Sjeremylt 15599d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 15609d007619Sjeremylt 1561b7c9bbdaSJeremy L Thompson @ref Advanced 15629d007619Sjeremylt **/ 1563d1d35e2fSjeremylt int CeedBasisGetQWeights(CeedBasis basis, const CeedScalar **q_weight) { 1564d1d35e2fSjeremylt *q_weight = basis->q_weight_1d; 1565e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 15669d007619Sjeremylt } 15679d007619Sjeremylt 15689d007619Sjeremylt /** 15699d007619Sjeremylt @brief Get interpolation matrix of a CeedBasis 15709d007619Sjeremylt 1571ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 15729d007619Sjeremylt @param[out] interp Variable to store interpolation matrix 15739d007619Sjeremylt 15749d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 15759d007619Sjeremylt 1576b7c9bbdaSJeremy L Thompson @ref Advanced 15779d007619Sjeremylt **/ 15786c58de82SJeremy L Thompson int CeedBasisGetInterp(CeedBasis basis, const CeedScalar **interp) { 1579d1d35e2fSjeremylt if (!basis->interp && basis->tensor_basis) { 15809d007619Sjeremylt // Allocate 15812b730f8bSJeremy L Thompson CeedCall(CeedMalloc(basis->Q * basis->P, &basis->interp)); 15829d007619Sjeremylt 15839d007619Sjeremylt // Initialize 15842b730f8bSJeremy L Thompson for (CeedInt i = 0; i < basis->Q * basis->P; i++) basis->interp[i] = 1.0; 15859d007619Sjeremylt 15869d007619Sjeremylt // Calculate 15872b730f8bSJeremy L Thompson for (CeedInt d = 0; d < basis->dim; d++) { 15882b730f8bSJeremy L Thompson for (CeedInt qpt = 0; qpt < basis->Q; qpt++) { 15899d007619Sjeremylt for (CeedInt node = 0; node < basis->P; node++) { 1590d1d35e2fSjeremylt CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d; 1591d1d35e2fSjeremylt CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d; 1592d1d35e2fSjeremylt basis->interp[qpt * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p]; 15939d007619Sjeremylt } 15949d007619Sjeremylt } 15952b730f8bSJeremy L Thompson } 15962b730f8bSJeremy L Thompson } 15979d007619Sjeremylt *interp = basis->interp; 1598e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 15999d007619Sjeremylt } 16009d007619Sjeremylt 16019d007619Sjeremylt /** 16029d007619Sjeremylt @brief Get 1D interpolation matrix of a tensor product CeedBasis 16039d007619Sjeremylt 1604ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 1605d1d35e2fSjeremylt @param[out] interp_1d Variable to store interpolation matrix 16069d007619Sjeremylt 16079d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 16089d007619Sjeremylt 16099d007619Sjeremylt @ref Backend 16109d007619Sjeremylt **/ 1611d1d35e2fSjeremylt int CeedBasisGetInterp1D(CeedBasis basis, const CeedScalar **interp_1d) { 1612*6574a04fSJeremy L Thompson CeedCheck(basis->tensor_basis, basis->ceed, CEED_ERROR_MINOR, "CeedBasis is not a tensor product basis."); 1613d1d35e2fSjeremylt *interp_1d = basis->interp_1d; 1614e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 16159d007619Sjeremylt } 16169d007619Sjeremylt 16179d007619Sjeremylt /** 16189d007619Sjeremylt @brief Get gradient matrix of a CeedBasis 16199d007619Sjeremylt 1620ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 16219d007619Sjeremylt @param[out] grad Variable to store gradient matrix 16229d007619Sjeremylt 16239d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 16249d007619Sjeremylt 1625b7c9bbdaSJeremy L Thompson @ref Advanced 16269d007619Sjeremylt **/ 16276c58de82SJeremy L Thompson int CeedBasisGetGrad(CeedBasis basis, const CeedScalar **grad) { 1628d1d35e2fSjeremylt if (!basis->grad && basis->tensor_basis) { 16299d007619Sjeremylt // Allocate 16302b730f8bSJeremy L Thompson CeedCall(CeedMalloc(basis->dim * basis->Q * basis->P, &basis->grad)); 16319d007619Sjeremylt 16329d007619Sjeremylt // Initialize 16332b730f8bSJeremy L Thompson for (CeedInt i = 0; i < basis->dim * basis->Q * basis->P; i++) basis->grad[i] = 1.0; 16349d007619Sjeremylt 16359d007619Sjeremylt // Calculate 16362b730f8bSJeremy L Thompson for (CeedInt d = 0; d < basis->dim; d++) { 16372b730f8bSJeremy L Thompson for (CeedInt i = 0; i < basis->dim; i++) { 16382b730f8bSJeremy L Thompson for (CeedInt qpt = 0; qpt < basis->Q; qpt++) { 16399d007619Sjeremylt for (CeedInt node = 0; node < basis->P; node++) { 1640d1d35e2fSjeremylt CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d; 1641d1d35e2fSjeremylt CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d; 16422b730f8bSJeremy L Thompson if (i == d) basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->grad_1d[q * basis->P_1d + p]; 16432b730f8bSJeremy L Thompson else basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p]; 16442b730f8bSJeremy L Thompson } 16452b730f8bSJeremy L Thompson } 16462b730f8bSJeremy L Thompson } 16479d007619Sjeremylt } 16489d007619Sjeremylt } 16499d007619Sjeremylt *grad = basis->grad; 1650e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 16519d007619Sjeremylt } 16529d007619Sjeremylt 16539d007619Sjeremylt /** 16549d007619Sjeremylt @brief Get 1D gradient matrix of a tensor product CeedBasis 16559d007619Sjeremylt 1656ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 1657d1d35e2fSjeremylt @param[out] grad_1d Variable to store gradient matrix 16589d007619Sjeremylt 16599d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 16609d007619Sjeremylt 1661b7c9bbdaSJeremy L Thompson @ref Advanced 16629d007619Sjeremylt **/ 1663d1d35e2fSjeremylt int CeedBasisGetGrad1D(CeedBasis basis, const CeedScalar **grad_1d) { 1664*6574a04fSJeremy L Thompson CeedCheck(basis->tensor_basis, basis->ceed, CEED_ERROR_MINOR, "CeedBasis is not a tensor product basis."); 1665d1d35e2fSjeremylt *grad_1d = basis->grad_1d; 1666e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 16679d007619Sjeremylt } 16689d007619Sjeremylt 16699d007619Sjeremylt /** 167050c301a5SRezgar Shakeri @brief Get divergence matrix of a CeedBasis 167150c301a5SRezgar Shakeri 1672ea61e9acSJeremy L Thompson @param[in] basis CeedBasis 167350c301a5SRezgar Shakeri @param[out] div Variable to store divergence matrix 167450c301a5SRezgar Shakeri 167550c301a5SRezgar Shakeri @return An error code: 0 - success, otherwise - failure 167650c301a5SRezgar Shakeri 167750c301a5SRezgar Shakeri @ref Advanced 167850c301a5SRezgar Shakeri **/ 167950c301a5SRezgar Shakeri int CeedBasisGetDiv(CeedBasis basis, const CeedScalar **div) { 1680*6574a04fSJeremy L Thompson CeedCheck(basis->div, basis->ceed, CEED_ERROR_MINOR, "CeedBasis does not have divergence matrix."); 168150c301a5SRezgar Shakeri *div = basis->div; 168250c301a5SRezgar Shakeri return CEED_ERROR_SUCCESS; 168350c301a5SRezgar Shakeri } 168450c301a5SRezgar Shakeri 168550c301a5SRezgar Shakeri /** 1686c4e3f59bSSebastian Grimberg @brief Get curl matrix of a CeedBasis 1687c4e3f59bSSebastian Grimberg 1688c4e3f59bSSebastian Grimberg @param[in] basis CeedBasis 1689c4e3f59bSSebastian Grimberg @param[out] curl Variable to store curl matrix 1690c4e3f59bSSebastian Grimberg 1691c4e3f59bSSebastian Grimberg @return An error code: 0 - success, otherwise - failure 1692c4e3f59bSSebastian Grimberg 1693c4e3f59bSSebastian Grimberg @ref Advanced 1694c4e3f59bSSebastian Grimberg **/ 1695c4e3f59bSSebastian Grimberg int CeedBasisGetCurl(CeedBasis basis, const CeedScalar **curl) { 1696*6574a04fSJeremy L Thompson CeedCheck(basis->curl, basis->ceed, CEED_ERROR_MINOR, "CeedBasis does not have curl matrix."); 1697c4e3f59bSSebastian Grimberg *curl = basis->curl; 1698c4e3f59bSSebastian Grimberg return CEED_ERROR_SUCCESS; 1699c4e3f59bSSebastian Grimberg } 1700c4e3f59bSSebastian Grimberg 1701c4e3f59bSSebastian Grimberg /** 17027a982d89SJeremy L. Thompson @brief Destroy a CeedBasis 17037a982d89SJeremy L. Thompson 1704ea61e9acSJeremy L Thompson @param[in,out] basis CeedBasis to destroy 17057a982d89SJeremy L. Thompson 17067a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 17077a982d89SJeremy L. Thompson 17087a982d89SJeremy L. Thompson @ref User 17097a982d89SJeremy L. Thompson **/ 17107a982d89SJeremy L. Thompson int CeedBasisDestroy(CeedBasis *basis) { 17117425e127SJeremy L Thompson if (!*basis || *basis == CEED_BASIS_COLLOCATED || --(*basis)->ref_count > 0) { 1712ad6481ceSJeremy L Thompson *basis = NULL; 1713ad6481ceSJeremy L Thompson return CEED_ERROR_SUCCESS; 1714ad6481ceSJeremy L Thompson } 17152b730f8bSJeremy L Thompson if ((*basis)->Destroy) CeedCall((*basis)->Destroy(*basis)); 17162b730f8bSJeremy L Thompson if ((*basis)->contract) CeedCall(CeedTensorContractDestroy(&(*basis)->contract)); 1717c4e3f59bSSebastian Grimberg CeedCall(CeedFree(&(*basis)->q_ref_1d)); 1718c4e3f59bSSebastian Grimberg CeedCall(CeedFree(&(*basis)->q_weight_1d)); 17192b730f8bSJeremy L Thompson CeedCall(CeedFree(&(*basis)->interp)); 17202b730f8bSJeremy L Thompson CeedCall(CeedFree(&(*basis)->interp_1d)); 17212b730f8bSJeremy L Thompson CeedCall(CeedFree(&(*basis)->grad)); 17222b730f8bSJeremy L Thompson CeedCall(CeedFree(&(*basis)->grad_1d)); 1723c4e3f59bSSebastian Grimberg CeedCall(CeedFree(&(*basis)->div)); 1724c4e3f59bSSebastian Grimberg CeedCall(CeedFree(&(*basis)->curl)); 17252b730f8bSJeremy L Thompson CeedCall(CeedDestroy(&(*basis)->ceed)); 17262b730f8bSJeremy L Thompson CeedCall(CeedFree(basis)); 1727e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 17287a982d89SJeremy L. Thompson } 17297a982d89SJeremy L. Thompson 17307a982d89SJeremy L. Thompson /** 1731b11c1e72Sjeremylt @brief Construct a Gauss-Legendre quadrature 1732b11c1e72Sjeremylt 1733ea61e9acSJeremy L Thompson @param[in] Q Number of quadrature points (integrates polynomials of degree 2*Q-1 exactly) 1734d1d35e2fSjeremylt @param[out] q_ref_1d Array of length Q to hold the abscissa on [-1, 1] 1735d1d35e2fSjeremylt @param[out] q_weight_1d Array of length Q to hold the weights 1736b11c1e72Sjeremylt 1737b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 1738dfdf5a53Sjeremylt 1739dfdf5a53Sjeremylt @ref Utility 1740b11c1e72Sjeremylt **/ 17412b730f8bSJeremy L Thompson int CeedGaussQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) { 1742d7b241e6Sjeremylt // Allocate 1743d7b241e6Sjeremylt CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0 * atan(1.0); 1744d1d35e2fSjeremylt // Build q_ref_1d, q_weight_1d 174592ae7e47SJeremy L Thompson for (CeedInt i = 0; i <= Q / 2; i++) { 1746d7b241e6Sjeremylt // Guess 1747d7b241e6Sjeremylt xi = cos(PI * (CeedScalar)(2 * i + 1) / ((CeedScalar)(2 * Q))); 1748d7b241e6Sjeremylt // Pn(xi) 1749d7b241e6Sjeremylt P0 = 1.0; 1750d7b241e6Sjeremylt P1 = xi; 1751d7b241e6Sjeremylt P2 = 0.0; 175292ae7e47SJeremy L Thompson for (CeedInt j = 2; j <= Q; j++) { 1753d7b241e6Sjeremylt P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j)); 1754d7b241e6Sjeremylt P0 = P1; 1755d7b241e6Sjeremylt P1 = P2; 1756d7b241e6Sjeremylt } 1757d7b241e6Sjeremylt // First Newton Step 1758d7b241e6Sjeremylt dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0); 1759d7b241e6Sjeremylt xi = xi - P2 / dP2; 1760d7b241e6Sjeremylt // Newton to convergence 176192ae7e47SJeremy L Thompson for (CeedInt k = 0; k < 100 && fabs(P2) > 10 * CEED_EPSILON; k++) { 1762d7b241e6Sjeremylt P0 = 1.0; 1763d7b241e6Sjeremylt P1 = xi; 176492ae7e47SJeremy L Thompson for (CeedInt j = 2; j <= Q; j++) { 1765d7b241e6Sjeremylt P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j)); 1766d7b241e6Sjeremylt P0 = P1; 1767d7b241e6Sjeremylt P1 = P2; 1768d7b241e6Sjeremylt } 1769d7b241e6Sjeremylt dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0); 1770d7b241e6Sjeremylt xi = xi - P2 / dP2; 1771d7b241e6Sjeremylt } 1772d7b241e6Sjeremylt // Save xi, wi 1773d7b241e6Sjeremylt wi = 2.0 / ((1.0 - xi * xi) * dP2 * dP2); 1774d1d35e2fSjeremylt q_weight_1d[i] = wi; 1775d1d35e2fSjeremylt q_weight_1d[Q - 1 - i] = wi; 1776d1d35e2fSjeremylt q_ref_1d[i] = -xi; 1777d1d35e2fSjeremylt q_ref_1d[Q - 1 - i] = xi; 1778d7b241e6Sjeremylt } 1779e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1780d7b241e6Sjeremylt } 1781d7b241e6Sjeremylt 1782b11c1e72Sjeremylt /** 1783b11c1e72Sjeremylt @brief Construct a Gauss-Legendre-Lobatto quadrature 1784b11c1e72Sjeremylt 1785ea61e9acSJeremy L Thompson @param[in] Q Number of quadrature points (integrates polynomials of degree 2*Q-3 exactly) 1786d1d35e2fSjeremylt @param[out] q_ref_1d Array of length Q to hold the abscissa on [-1, 1] 1787d1d35e2fSjeremylt @param[out] q_weight_1d Array of length Q to hold the weights 1788b11c1e72Sjeremylt 1789b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 1790dfdf5a53Sjeremylt 1791dfdf5a53Sjeremylt @ref Utility 1792b11c1e72Sjeremylt **/ 17932b730f8bSJeremy L Thompson int CeedLobattoQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) { 1794d7b241e6Sjeremylt // Allocate 1795d7b241e6Sjeremylt CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0 * atan(1.0); 1796d1d35e2fSjeremylt // Build q_ref_1d, q_weight_1d 1797d7b241e6Sjeremylt // Set endpoints 1798*6574a04fSJeremy L Thompson CeedCheck(Q > 1, NULL, CEED_ERROR_DIMENSION, "Cannot create Lobatto quadrature with Q=%" CeedInt_FMT " < 2 points", Q); 1799d7b241e6Sjeremylt wi = 2.0 / ((CeedScalar)(Q * (Q - 1))); 1800d1d35e2fSjeremylt if (q_weight_1d) { 1801d1d35e2fSjeremylt q_weight_1d[0] = wi; 1802d1d35e2fSjeremylt q_weight_1d[Q - 1] = wi; 1803d7b241e6Sjeremylt } 1804d1d35e2fSjeremylt q_ref_1d[0] = -1.0; 1805d1d35e2fSjeremylt q_ref_1d[Q - 1] = 1.0; 1806d7b241e6Sjeremylt // Interior 180792ae7e47SJeremy L Thompson for (CeedInt i = 1; i <= (Q - 1) / 2; i++) { 1808d7b241e6Sjeremylt // Guess 1809d7b241e6Sjeremylt xi = cos(PI * (CeedScalar)(i) / (CeedScalar)(Q - 1)); 1810d7b241e6Sjeremylt // Pn(xi) 1811d7b241e6Sjeremylt P0 = 1.0; 1812d7b241e6Sjeremylt P1 = xi; 1813d7b241e6Sjeremylt P2 = 0.0; 181492ae7e47SJeremy L Thompson for (CeedInt j = 2; j < Q; j++) { 1815d7b241e6Sjeremylt P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j)); 1816d7b241e6Sjeremylt P0 = P1; 1817d7b241e6Sjeremylt P1 = P2; 1818d7b241e6Sjeremylt } 1819d7b241e6Sjeremylt // First Newton step 1820d7b241e6Sjeremylt dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0); 1821d7b241e6Sjeremylt d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi); 1822d7b241e6Sjeremylt xi = xi - dP2 / d2P2; 1823d7b241e6Sjeremylt // Newton to convergence 182492ae7e47SJeremy L Thompson for (CeedInt k = 0; k < 100 && fabs(dP2) > 10 * CEED_EPSILON; k++) { 1825d7b241e6Sjeremylt P0 = 1.0; 1826d7b241e6Sjeremylt P1 = xi; 182792ae7e47SJeremy L Thompson for (CeedInt j = 2; j < Q; j++) { 1828d7b241e6Sjeremylt P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j)); 1829d7b241e6Sjeremylt P0 = P1; 1830d7b241e6Sjeremylt P1 = P2; 1831d7b241e6Sjeremylt } 1832d7b241e6Sjeremylt dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0); 1833d7b241e6Sjeremylt d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi); 1834d7b241e6Sjeremylt xi = xi - dP2 / d2P2; 1835d7b241e6Sjeremylt } 1836d7b241e6Sjeremylt // Save xi, wi 1837d7b241e6Sjeremylt wi = 2.0 / (((CeedScalar)(Q * (Q - 1))) * P2 * P2); 1838d1d35e2fSjeremylt if (q_weight_1d) { 1839d1d35e2fSjeremylt q_weight_1d[i] = wi; 1840d1d35e2fSjeremylt q_weight_1d[Q - 1 - i] = wi; 1841d7b241e6Sjeremylt } 1842d1d35e2fSjeremylt q_ref_1d[i] = -xi; 1843d1d35e2fSjeremylt q_ref_1d[Q - 1 - i] = xi; 1844d7b241e6Sjeremylt } 1845e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1846d7b241e6Sjeremylt } 1847d7b241e6Sjeremylt 1848d7b241e6Sjeremylt /// @} 1849