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> 9*2b730f8bSJeremy L Thompson #include <ceed/backend.h> 10*2b730f8bSJeremy L Thompson #include <ceed/ceed.h> 11d7b241e6Sjeremylt #include <math.h> 123d576824SJeremy L Thompson #include <stdbool.h> 13d7b241e6Sjeremylt #include <stdio.h> 14d7b241e6Sjeremylt #include <string.h> 15d7b241e6Sjeremylt 167a982d89SJeremy L. Thompson /// @file 177a982d89SJeremy L. Thompson /// Implementation of CeedBasis interfaces 187a982d89SJeremy L. Thompson 19d7b241e6Sjeremylt /// @cond DOXYGEN_SKIP 20783c99b3SValeria Barra static struct CeedBasis_private ceed_basis_collocated; 21d7b241e6Sjeremylt /// @endcond 22d7b241e6Sjeremylt 237a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser 247a982d89SJeremy L. Thompson /// @{ 257a982d89SJeremy L. Thompson 267a982d89SJeremy L. Thompson /// Indicate that the quadrature points are collocated with the nodes 277a982d89SJeremy L. Thompson const CeedBasis CEED_BASIS_COLLOCATED = &ceed_basis_collocated; 287a982d89SJeremy L. Thompson 297a982d89SJeremy L. Thompson /// @} 307a982d89SJeremy L. Thompson 317a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 327a982d89SJeremy L. Thompson /// CeedBasis Library Internal Functions 337a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 347a982d89SJeremy L. Thompson /// @addtogroup CeedBasisDeveloper 357a982d89SJeremy L. Thompson /// @{ 367a982d89SJeremy L. Thompson 377a982d89SJeremy L. Thompson /** 387a982d89SJeremy L. Thompson @brief Compute Householder reflection 397a982d89SJeremy L. Thompson 407a982d89SJeremy L. Thompson Computes A = (I - b v v^T) A 417a982d89SJeremy L. Thompson where A is an mxn matrix indexed as A[i*row + j*col] 427a982d89SJeremy L. Thompson 437a982d89SJeremy L. Thompson @param[in,out] A Matrix to apply Householder reflection to, in place 447a982d89SJeremy L. Thompson @param v Householder vector 457a982d89SJeremy L. Thompson @param b Scaling factor 467a982d89SJeremy L. Thompson @param m Number of rows in A 477a982d89SJeremy L. Thompson @param n Number of columns in A 487a982d89SJeremy L. Thompson @param row Row stride 497a982d89SJeremy L. Thompson @param col Col stride 507a982d89SJeremy L. Thompson 517a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 527a982d89SJeremy L. Thompson 537a982d89SJeremy L. Thompson @ref Developer 547a982d89SJeremy L. Thompson **/ 55*2b730f8bSJeremy L Thompson static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v, CeedScalar b, CeedInt m, CeedInt n, CeedInt row, CeedInt col) { 567a982d89SJeremy L. Thompson for (CeedInt j = 0; j < n; j++) { 577a982d89SJeremy L. Thompson CeedScalar w = A[0 * row + j * col]; 58*2b730f8bSJeremy L Thompson for (CeedInt i = 1; i < m; i++) w += v[i] * A[i * row + j * col]; 597a982d89SJeremy L. Thompson A[0 * row + j * col] -= b * w; 60*2b730f8bSJeremy L Thompson for (CeedInt i = 1; i < m; i++) A[i * row + j * col] -= b * w * v[i]; 617a982d89SJeremy L. Thompson } 62e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 637a982d89SJeremy L. Thompson } 647a982d89SJeremy L. Thompson 657a982d89SJeremy L. Thompson /** 667a982d89SJeremy L. Thompson @brief Apply Householder Q matrix 677a982d89SJeremy L. Thompson 687a982d89SJeremy L. Thompson Compute A = Q A where Q is mxm and A is mxn. 697a982d89SJeremy L. Thompson 707a982d89SJeremy L. Thompson @param[in,out] A Matrix to apply Householder Q to, in place 717a982d89SJeremy L. Thompson @param Q Householder Q matrix 727a982d89SJeremy L. Thompson @param tau Householder scaling factors 73d1d35e2fSjeremylt @param t_mode Transpose mode for application 747a982d89SJeremy L. Thompson @param m Number of rows in A 757a982d89SJeremy L. Thompson @param n Number of columns in A 767a982d89SJeremy L. Thompson @param k Number of elementary reflectors in Q, k<m 777a982d89SJeremy L. Thompson @param row Row stride in A 787a982d89SJeremy L. Thompson @param col Col stride in A 797a982d89SJeremy L. Thompson 807a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 817a982d89SJeremy L. Thompson 827a982d89SJeremy L. Thompson @ref Developer 837a982d89SJeremy L. Thompson **/ 84*2b730f8bSJeremy L Thompson int CeedHouseholderApplyQ(CeedScalar *A, const CeedScalar *Q, const CeedScalar *tau, CeedTransposeMode t_mode, CeedInt m, CeedInt n, CeedInt k, 857a982d89SJeremy L. Thompson CeedInt row, CeedInt col) { 8678464608Sjeremylt CeedScalar *v; 87*2b730f8bSJeremy L Thompson CeedCall(CeedMalloc(m, &v)); 887a982d89SJeremy L. Thompson for (CeedInt ii = 0; ii < k; ii++) { 89d1d35e2fSjeremylt CeedInt i = t_mode == CEED_TRANSPOSE ? ii : k - 1 - ii; 90*2b730f8bSJeremy L Thompson for (CeedInt j = i + 1; j < m; j++) v[j] = Q[j * k + i]; 91d1d35e2fSjeremylt // Apply Householder reflector (I - tau v v^T) collo_grad_1d^T 92*2b730f8bSJeremy L Thompson CeedCall(CeedHouseholderReflect(&A[i * row], &v[i], tau[i], m - i, n, row, col)); 937a982d89SJeremy L. Thompson } 94*2b730f8bSJeremy L Thompson CeedCall(CeedFree(&v)); 95e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 967a982d89SJeremy L. Thompson } 977a982d89SJeremy L. Thompson 987a982d89SJeremy L. Thompson /** 997a982d89SJeremy L. Thompson @brief Compute Givens rotation 1007a982d89SJeremy L. Thompson 1017a982d89SJeremy L. Thompson Computes A = G A (or G^T A in transpose mode) 1027a982d89SJeremy L. Thompson where A is an mxn matrix indexed as A[i*n + j*m] 1037a982d89SJeremy L. Thompson 1047a982d89SJeremy L. Thompson @param[in,out] A Row major matrix to apply Givens rotation to, in place 1057a982d89SJeremy L. Thompson @param c Cosine factor 1067a982d89SJeremy L. Thompson @param s Sine factor 107d1d35e2fSjeremylt @param t_mode @ref CEED_NOTRANSPOSE to rotate the basis counter-clockwise, 1084c4400c7SValeria Barra which has the effect of rotating columns of A clockwise; 1094cc79fe7SJed Brown @ref CEED_TRANSPOSE for the opposite rotation 1107a982d89SJeremy L. Thompson @param i First row/column to apply rotation 1117a982d89SJeremy L. Thompson @param k Second row/column to apply rotation 1127a982d89SJeremy L. Thompson @param m Number of rows in A 1137a982d89SJeremy L. Thompson @param n Number of columns in A 1147a982d89SJeremy L. Thompson 1157a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 1167a982d89SJeremy L. Thompson 1177a982d89SJeremy L. Thompson @ref Developer 1187a982d89SJeremy L. Thompson **/ 119*2b730f8bSJeremy L Thompson static int CeedGivensRotation(CeedScalar *A, CeedScalar c, CeedScalar s, CeedTransposeMode t_mode, CeedInt i, CeedInt k, CeedInt m, CeedInt n) { 120d1d35e2fSjeremylt CeedInt stride_j = 1, stride_ik = m, num_its = n; 121d1d35e2fSjeremylt if (t_mode == CEED_NOTRANSPOSE) { 122*2b730f8bSJeremy L Thompson stride_j = n; 123*2b730f8bSJeremy L Thompson stride_ik = 1; 124*2b730f8bSJeremy L Thompson num_its = m; 1257a982d89SJeremy L. Thompson } 1267a982d89SJeremy L. Thompson 1277a982d89SJeremy L. Thompson // Apply rotation 128d1d35e2fSjeremylt for (CeedInt j = 0; j < num_its; j++) { 129d1d35e2fSjeremylt CeedScalar tau1 = A[i * stride_ik + j * stride_j], tau2 = A[k * stride_ik + j * stride_j]; 130d1d35e2fSjeremylt A[i * stride_ik + j * stride_j] = c * tau1 - s * tau2; 131d1d35e2fSjeremylt A[k * stride_ik + j * stride_j] = s * tau1 + c * tau2; 1327a982d89SJeremy L. Thompson } 133e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1347a982d89SJeremy L. Thompson } 1357a982d89SJeremy L. Thompson 1367a982d89SJeremy L. Thompson /** 1377a982d89SJeremy L. Thompson @brief View an array stored in a CeedBasis 1387a982d89SJeremy L. Thompson 1390a0da059Sjeremylt @param[in] name Name of array 140d1d35e2fSjeremylt @param[in] fp_fmt Printing format 1410a0da059Sjeremylt @param[in] m Number of rows in array 1420a0da059Sjeremylt @param[in] n Number of columns in array 1430a0da059Sjeremylt @param[in] a Array to be viewed 1440a0da059Sjeremylt @param[in] stream Stream to view to, e.g., stdout 1457a982d89SJeremy L. Thompson 1467a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 1477a982d89SJeremy L. Thompson 1487a982d89SJeremy L. Thompson @ref Developer 1497a982d89SJeremy L. Thompson **/ 150*2b730f8bSJeremy L Thompson static int CeedScalarView(const char *name, const char *fp_fmt, CeedInt m, CeedInt n, const CeedScalar *a, FILE *stream) { 15192ae7e47SJeremy L Thompson for (CeedInt i = 0; i < m; i++) { 152*2b730f8bSJeremy L Thompson if (m > 1) fprintf(stream, "%12s[%" CeedInt_FMT "]:", name, i); 153*2b730f8bSJeremy L Thompson else fprintf(stream, "%12s:", name); 154*2b730f8bSJeremy 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); 1557a982d89SJeremy L. Thompson fputs("\n", stream); 1567a982d89SJeremy L. Thompson } 157e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1587a982d89SJeremy L. Thompson } 1597a982d89SJeremy L. Thompson 160a76a04e7SJeremy L Thompson /** 16114556e63SJeremy L Thompson @brief Create the interpolation and gradient matrices for projection from 16214556e63SJeremy L Thompson the nodes of `basis_from` to the nodes of `basis_to`. 16314556e63SJeremy L Thompson The interpolation is given by `interp_project = interp_to^+ * interp_from`, 16414556e63SJeremy L Thompson where the pesudoinverse `interp_to^+` is given by QR factorization. 16514556e63SJeremy L Thompson The gradient is given by `grad_project = interp_to^+ * grad_from`. 166a76a04e7SJeremy L Thompson Note: `basis_from` and `basis_to` must have compatible quadrature 167a76a04e7SJeremy L Thompson spaces. 168a76a04e7SJeremy L Thompson 169a76a04e7SJeremy L Thompson @param[in] basis_from CeedBasis to project from 170a76a04e7SJeremy L Thompson @param[in] basis_to CeedBasis to project to 171a76a04e7SJeremy L Thompson @param[out] interp_project Address of the variable where the newly created 17214556e63SJeremy L Thompson interpolation matrix will be stored. 17314556e63SJeremy L Thompson @param[out] grad_project Address of the variable where the newly created 17414556e63SJeremy L Thompson gradient matrix will be stored. 175a76a04e7SJeremy L Thompson 176a76a04e7SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 177a76a04e7SJeremy L Thompson 178a76a04e7SJeremy L Thompson @ref Developer 179a76a04e7SJeremy L Thompson **/ 180*2b730f8bSJeremy L Thompson static int CeedBasisCreateProjectionMatrices(CeedBasis basis_from, CeedBasis basis_to, CeedScalar **interp_project, CeedScalar **grad_project) { 181a76a04e7SJeremy L Thompson Ceed ceed; 182*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetCeed(basis_to, &ceed)); 183a76a04e7SJeremy L Thompson 184a76a04e7SJeremy L Thompson // Check for compatible quadrature spaces 185a76a04e7SJeremy L Thompson CeedInt Q_to, Q_from; 186*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints(basis_to, &Q_to)); 187*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints(basis_from, &Q_from)); 188*2b730f8bSJeremy L Thompson if (Q_to != Q_from) { 189a76a04e7SJeremy L Thompson // LCOV_EXCL_START 190*2b730f8bSJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, "Bases must have compatible quadrature spaces"); 191a76a04e7SJeremy L Thompson // LCOV_EXCL_STOP 192*2b730f8bSJeremy L Thompson } 193a76a04e7SJeremy L Thompson 19414556e63SJeremy L Thompson // Check for matching tensor or non-tensor 195a76a04e7SJeremy L Thompson CeedInt P_to, P_from, Q = Q_to; 196a76a04e7SJeremy L Thompson bool is_tensor_to, is_tensor_from; 197*2b730f8bSJeremy L Thompson CeedCall(CeedBasisIsTensor(basis_to, &is_tensor_to)); 198*2b730f8bSJeremy L Thompson CeedCall(CeedBasisIsTensor(basis_from, &is_tensor_from)); 199a76a04e7SJeremy L Thompson if (is_tensor_to && is_tensor_from) { 200*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_to)); 201*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_from)); 202*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints1D(basis_from, &Q)); 203a76a04e7SJeremy L Thompson } else if (!is_tensor_to && !is_tensor_from) { 204*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis_to, &P_to)); 205*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis_from, &P_from)); 206a76a04e7SJeremy L Thompson } else { 207a76a04e7SJeremy L Thompson // LCOV_EXCL_START 208*2b730f8bSJeremy L Thompson return CeedError(ceed, CEED_ERROR_MINOR, "Bases must both be tensor or non-tensor"); 209a76a04e7SJeremy L Thompson // LCOV_EXCL_STOP 210a76a04e7SJeremy L Thompson } 211a76a04e7SJeremy L Thompson 21214556e63SJeremy L Thompson // Get source matrices 21314556e63SJeremy L Thompson CeedInt dim; 21414556e63SJeremy L Thompson CeedScalar *interp_to, *interp_from, *tau; 215*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetDimension(basis_to, &dim)); 216*2b730f8bSJeremy L Thompson CeedCall(CeedMalloc(Q * P_from, &interp_from)); 217*2b730f8bSJeremy L Thompson CeedCall(CeedMalloc(Q * P_to, &interp_to)); 218*2b730f8bSJeremy L Thompson CeedCall(CeedCalloc(P_to * P_from, interp_project)); 219*2b730f8bSJeremy L Thompson CeedCall(CeedCalloc(P_to * P_from * (is_tensor_to ? 1 : dim), grad_project)); 220*2b730f8bSJeremy L Thompson CeedCall(CeedMalloc(Q, &tau)); 221*2b730f8bSJeremy L Thompson const CeedScalar *interp_to_source = NULL, *interp_from_source = NULL, *grad_from_source; 222a76a04e7SJeremy L Thompson if (is_tensor_to) { 223*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetInterp1D(basis_to, &interp_to_source)); 224*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetInterp1D(basis_from, &interp_from_source)); 225*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetGrad1D(basis_from, &grad_from_source)); 226a76a04e7SJeremy L Thompson } else { 227*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetInterp(basis_to, &interp_to_source)); 228*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetInterp(basis_from, &interp_from_source)); 229*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetGrad(basis_from, &grad_from_source)); 230a76a04e7SJeremy L Thompson } 231a76a04e7SJeremy L Thompson 23214556e63SJeremy L Thompson // Build matrices 23314556e63SJeremy L Thompson CeedInt num_matrices = 1 + (is_tensor_to ? 1 : dim); 23414556e63SJeremy L Thompson CeedScalar *input_from[num_matrices], *output_project[num_matrices]; 23514556e63SJeremy L Thompson input_from[0] = (CeedScalar *)interp_from_source; 23614556e63SJeremy L Thompson output_project[0] = *interp_project; 23714556e63SJeremy L Thompson for (CeedInt m = 1; m < num_matrices; m++) { 23814556e63SJeremy L Thompson input_from[m] = (CeedScalar *)&grad_from_source[(m - 1) * Q * P_from]; 23902af4036SJeremy L Thompson output_project[m] = &((*grad_project)[(m - 1) * P_to * P_from]); 24014556e63SJeremy L Thompson } 24114556e63SJeremy L Thompson for (CeedInt m = 0; m < num_matrices; m++) { 242a76a04e7SJeremy L Thompson // -- QR Factorization, interp_to = Q R 24314556e63SJeremy L Thompson memcpy(interp_to, interp_to_source, Q * P_to * sizeof(interp_to_source[0])); 244*2b730f8bSJeremy L Thompson CeedCall(CeedQRFactorization(ceed, interp_to, tau, Q, P_to)); 245a76a04e7SJeremy L Thompson 246a76a04e7SJeremy L Thompson // -- Apply Qtranspose, interp_to = Qtranspose interp_from 24714556e63SJeremy L Thompson memcpy(interp_from, input_from[m], Q * P_from * sizeof(input_from[m][0])); 248*2b730f8bSJeremy L Thompson CeedCall(CeedHouseholderApplyQ(interp_from, interp_to, tau, CEED_TRANSPOSE, Q, P_from, P_to, P_from, 1)); 249a76a04e7SJeremy L Thompson 250a76a04e7SJeremy L Thompson // -- Apply Rinv, interp_project = Rinv interp_c 251a76a04e7SJeremy L Thompson for (CeedInt j = 0; j < P_from; j++) { // Column j 252*2b730f8bSJeremy 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]; 253a76a04e7SJeremy L Thompson for (CeedInt i = P_to - 2; i >= 0; i--) { // Row i 25414556e63SJeremy L Thompson output_project[m][j + P_from * i] = interp_from[j + P_from * i]; 255a76a04e7SJeremy L Thompson for (CeedInt k = i + 1; k < P_to; k++) { 256*2b730f8bSJeremy L Thompson output_project[m][j + P_from * i] -= interp_to[k + P_to * i] * output_project[m][j + P_from * k]; 257a76a04e7SJeremy L Thompson } 25814556e63SJeremy L Thompson output_project[m][j + P_from * i] /= interp_to[i + P_to * i]; 259a76a04e7SJeremy L Thompson } 260a76a04e7SJeremy L Thompson } 26114556e63SJeremy L Thompson } 26214556e63SJeremy L Thompson 26314556e63SJeremy L Thompson // Cleanup 264*2b730f8bSJeremy L Thompson CeedCall(CeedFree(&tau)); 265*2b730f8bSJeremy L Thompson CeedCall(CeedFree(&interp_to)); 266*2b730f8bSJeremy L Thompson CeedCall(CeedFree(&interp_from)); 267a76a04e7SJeremy L Thompson 268a76a04e7SJeremy L Thompson return CEED_ERROR_SUCCESS; 269a76a04e7SJeremy L Thompson } 270a76a04e7SJeremy L Thompson 2717a982d89SJeremy L. Thompson /// @} 2727a982d89SJeremy L. Thompson 2737a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 2747a982d89SJeremy L. Thompson /// Ceed Backend API 2757a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 2767a982d89SJeremy L. Thompson /// @addtogroup CeedBasisBackend 2777a982d89SJeremy L. Thompson /// @{ 2787a982d89SJeremy L. Thompson 2797a982d89SJeremy L. Thompson /** 2807a982d89SJeremy L. Thompson @brief Return collocated grad matrix 2817a982d89SJeremy L. Thompson 2827a982d89SJeremy L. Thompson @param basis CeedBasis 283d1d35e2fSjeremylt @param[out] collo_grad_1d Row-major (Q_1d * Q_1d) matrix expressing derivatives of 2847a982d89SJeremy L. Thompson basis functions at quadrature points 2857a982d89SJeremy L. Thompson 2867a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 2877a982d89SJeremy L. Thompson 2887a982d89SJeremy L. Thompson @ref Backend 2897a982d89SJeremy L. Thompson **/ 290d1d35e2fSjeremylt int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *collo_grad_1d) { 2917a982d89SJeremy L. Thompson int i, j, k; 2927a982d89SJeremy L. Thompson Ceed ceed; 293*2b730f8bSJeremy L Thompson CeedInt P_1d = (basis)->P_1d, Q_1d = (basis)->Q_1d; 29478464608Sjeremylt CeedScalar *interp_1d, *grad_1d, *tau; 2957a982d89SJeremy L. Thompson 296*2b730f8bSJeremy L Thompson CeedCall(CeedMalloc(Q_1d * P_1d, &interp_1d)); 297*2b730f8bSJeremy L Thompson CeedCall(CeedMalloc(Q_1d * P_1d, &grad_1d)); 298*2b730f8bSJeremy L Thompson CeedCall(CeedMalloc(Q_1d, &tau)); 299d1d35e2fSjeremylt memcpy(interp_1d, (basis)->interp_1d, Q_1d * P_1d * sizeof(basis)->interp_1d[0]); 300d1d35e2fSjeremylt memcpy(grad_1d, (basis)->grad_1d, Q_1d * P_1d * sizeof(basis)->interp_1d[0]); 3017a982d89SJeremy L. Thompson 302d1d35e2fSjeremylt // QR Factorization, interp_1d = Q R 303*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetCeed(basis, &ceed)); 304*2b730f8bSJeremy L Thompson CeedCall(CeedQRFactorization(ceed, interp_1d, tau, Q_1d, P_1d)); 305e15f9bd0SJeremy L Thompson // Note: This function is for backend use, so all errors are terminal 306e15f9bd0SJeremy L Thompson // and we do not need to clean up memory on failure. 3077a982d89SJeremy L. Thompson 308d1d35e2fSjeremylt // Apply Rinv, collo_grad_1d = grad_1d Rinv 309d1d35e2fSjeremylt for (i = 0; i < Q_1d; i++) { // Row i 310d1d35e2fSjeremylt collo_grad_1d[Q_1d * i] = grad_1d[P_1d * i] / interp_1d[0]; 311d1d35e2fSjeremylt for (j = 1; j < P_1d; j++) { // Column j 312d1d35e2fSjeremylt collo_grad_1d[j + Q_1d * i] = grad_1d[j + P_1d * i]; 313*2b730f8bSJeremy 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]; 314d1d35e2fSjeremylt collo_grad_1d[j + Q_1d * i] /= interp_1d[j + P_1d * j]; 3157a982d89SJeremy L. Thompson } 316*2b730f8bSJeremy L Thompson for (j = P_1d; j < Q_1d; j++) collo_grad_1d[j + Q_1d * i] = 0; 3177a982d89SJeremy L. Thompson } 3187a982d89SJeremy L. Thompson 319673160d7Sjeremylt // Apply Qtranspose, collo_grad = collo_grad Q_transpose 320*2b730f8bSJeremy L Thompson CeedCall(CeedHouseholderApplyQ(collo_grad_1d, interp_1d, tau, CEED_NOTRANSPOSE, Q_1d, Q_1d, P_1d, 1, Q_1d)); 3217a982d89SJeremy L. Thompson 322*2b730f8bSJeremy L Thompson CeedCall(CeedFree(&interp_1d)); 323*2b730f8bSJeremy L Thompson CeedCall(CeedFree(&grad_1d)); 324*2b730f8bSJeremy L Thompson CeedCall(CeedFree(&tau)); 325e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3267a982d89SJeremy L. Thompson } 3277a982d89SJeremy L. Thompson 3287a982d89SJeremy L. Thompson /** 3297a982d89SJeremy L. Thompson @brief Get tensor status for given CeedBasis 3307a982d89SJeremy L. Thompson 3317a982d89SJeremy L. Thompson @param basis CeedBasis 332d1d35e2fSjeremylt @param[out] is_tensor Variable to store tensor status 3337a982d89SJeremy L. Thompson 3347a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3357a982d89SJeremy L. Thompson 3367a982d89SJeremy L. Thompson @ref Backend 3377a982d89SJeremy L. Thompson **/ 338d1d35e2fSjeremylt int CeedBasisIsTensor(CeedBasis basis, bool *is_tensor) { 339d1d35e2fSjeremylt *is_tensor = basis->tensor_basis; 340e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3417a982d89SJeremy L. Thompson } 3427a982d89SJeremy L. Thompson 3437a982d89SJeremy L. Thompson /** 3447a982d89SJeremy L. Thompson @brief Get backend data of a CeedBasis 3457a982d89SJeremy L. Thompson 3467a982d89SJeremy L. Thompson @param basis CeedBasis 3477a982d89SJeremy L. Thompson @param[out] data Variable to store data 3487a982d89SJeremy L. Thompson 3497a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3507a982d89SJeremy L. Thompson 3517a982d89SJeremy L. Thompson @ref Backend 3527a982d89SJeremy L. Thompson **/ 353777ff853SJeremy L Thompson int CeedBasisGetData(CeedBasis basis, void *data) { 354777ff853SJeremy L Thompson *(void **)data = basis->data; 355e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3567a982d89SJeremy L. Thompson } 3577a982d89SJeremy L. Thompson 3587a982d89SJeremy L. Thompson /** 3597a982d89SJeremy L. Thompson @brief Set backend data of a CeedBasis 3607a982d89SJeremy L. Thompson 3617a982d89SJeremy L. Thompson @param[out] basis CeedBasis 3627a982d89SJeremy L. Thompson @param data Data to set 3637a982d89SJeremy L. Thompson 3647a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3657a982d89SJeremy L. Thompson 3667a982d89SJeremy L. Thompson @ref Backend 3677a982d89SJeremy L. Thompson **/ 368777ff853SJeremy L Thompson int CeedBasisSetData(CeedBasis basis, void *data) { 369777ff853SJeremy L Thompson basis->data = data; 370e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3717a982d89SJeremy L. Thompson } 3727a982d89SJeremy L. Thompson 3737a982d89SJeremy L. Thompson /** 37434359f16Sjeremylt @brief Increment the reference counter for a CeedBasis 37534359f16Sjeremylt 37634359f16Sjeremylt @param basis Basis to increment the reference counter 37734359f16Sjeremylt 37834359f16Sjeremylt @return An error code: 0 - success, otherwise - failure 37934359f16Sjeremylt 38034359f16Sjeremylt @ref Backend 38134359f16Sjeremylt **/ 3829560d06aSjeremylt int CeedBasisReference(CeedBasis basis) { 38334359f16Sjeremylt basis->ref_count++; 38434359f16Sjeremylt return CEED_ERROR_SUCCESS; 38534359f16Sjeremylt } 38634359f16Sjeremylt 38734359f16Sjeremylt /** 3886e15d496SJeremy L Thompson @brief Estimate number of FLOPs required to apply CeedBasis in t_mode and eval_mode 3896e15d496SJeremy L Thompson 3906e15d496SJeremy L Thompson @param basis Basis to estimate FLOPs for 3916e15d496SJeremy L Thompson @param t_mode Apply basis or transpose 3926e15d496SJeremy L Thompson @param eval_mode Basis evaluation mode 3936e15d496SJeremy L Thompson @param flops Address of variable to hold FLOPs estimate 3946e15d496SJeremy L Thompson 3956e15d496SJeremy L Thompson @ref Backend 3966e15d496SJeremy L Thompson **/ 397*2b730f8bSJeremy L Thompson int CeedBasisGetFlopsEstimate(CeedBasis basis, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedSize *flops) { 3986e15d496SJeremy L Thompson bool is_tensor; 3996e15d496SJeremy L Thompson 400*2b730f8bSJeremy L Thompson CeedCall(CeedBasisIsTensor(basis, &is_tensor)); 4016e15d496SJeremy L Thompson if (is_tensor) { 4026e15d496SJeremy L Thompson CeedInt dim, num_comp, P_1d, Q_1d; 403*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetDimension(basis, &dim)); 404*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumComponents(basis, &num_comp)); 405*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d)); 406*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 4076e15d496SJeremy L Thompson if (t_mode == CEED_TRANSPOSE) { 408*2b730f8bSJeremy L Thompson P_1d = Q_1d; 409*2b730f8bSJeremy L Thompson Q_1d = P_1d; 4106e15d496SJeremy L Thompson } 4116e15d496SJeremy L Thompson CeedInt tensor_flops = 0, pre = num_comp * CeedIntPow(P_1d, dim - 1), post = 1; 4126e15d496SJeremy L Thompson for (CeedInt d = 0; d < dim; d++) { 4136e15d496SJeremy L Thompson tensor_flops += 2 * pre * P_1d * post * Q_1d; 4146e15d496SJeremy L Thompson pre /= P_1d; 4156e15d496SJeremy L Thompson post *= Q_1d; 4166e15d496SJeremy L Thompson } 4176e15d496SJeremy L Thompson switch (eval_mode) { 418*2b730f8bSJeremy L Thompson case CEED_EVAL_NONE: 419*2b730f8bSJeremy L Thompson *flops = 0; 420*2b730f8bSJeremy L Thompson break; 421*2b730f8bSJeremy L Thompson case CEED_EVAL_INTERP: 422*2b730f8bSJeremy L Thompson *flops = tensor_flops; 423*2b730f8bSJeremy L Thompson break; 424*2b730f8bSJeremy L Thompson case CEED_EVAL_GRAD: 425*2b730f8bSJeremy L Thompson *flops = tensor_flops * 2; 426*2b730f8bSJeremy L Thompson break; 4276e15d496SJeremy L Thompson case CEED_EVAL_DIV: 4286e15d496SJeremy L Thompson // LCOV_EXCL_START 429*2b730f8bSJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_INCOMPATIBLE, "Tensor CEED_EVAL_DIV not supported"); 430*2b730f8bSJeremy L Thompson break; 4316e15d496SJeremy L Thompson case CEED_EVAL_CURL: 432*2b730f8bSJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_INCOMPATIBLE, "Tensor CEED_EVAL_CURL not supported"); 433*2b730f8bSJeremy L Thompson break; 4346e15d496SJeremy L Thompson // LCOV_EXCL_STOP 435*2b730f8bSJeremy L Thompson case CEED_EVAL_WEIGHT: 436*2b730f8bSJeremy L Thompson *flops = dim * CeedIntPow(Q_1d, dim); 437*2b730f8bSJeremy L Thompson break; 4386e15d496SJeremy L Thompson } 4396e15d496SJeremy L Thompson } else { 4406e15d496SJeremy L Thompson CeedInt dim, num_comp, num_nodes, num_qpts, Q_comp; 441*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetDimension(basis, &dim)); 442*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumComponents(basis, &num_comp)); 443*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis, &num_nodes)); 444*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts)); 445*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadratureComponents(basis, &Q_comp)); 4466e15d496SJeremy L Thompson switch (eval_mode) { 447*2b730f8bSJeremy L Thompson case CEED_EVAL_NONE: 448*2b730f8bSJeremy L Thompson *flops = 0; 449*2b730f8bSJeremy L Thompson break; 450*2b730f8bSJeremy L Thompson case CEED_EVAL_INTERP: 451*2b730f8bSJeremy L Thompson *flops = num_nodes * num_qpts * num_comp; 452*2b730f8bSJeremy L Thompson break; 453*2b730f8bSJeremy L Thompson case CEED_EVAL_GRAD: 454*2b730f8bSJeremy L Thompson *flops = num_nodes * num_qpts * num_comp * dim; 455*2b730f8bSJeremy L Thompson break; 456*2b730f8bSJeremy L Thompson case CEED_EVAL_DIV: 457*2b730f8bSJeremy L Thompson *flops = num_nodes * num_qpts; 458*2b730f8bSJeremy L Thompson break; 459*2b730f8bSJeremy L Thompson case CEED_EVAL_CURL: 460*2b730f8bSJeremy L Thompson *flops = num_nodes * num_qpts * dim; 461*2b730f8bSJeremy L Thompson break; 462*2b730f8bSJeremy L Thompson case CEED_EVAL_WEIGHT: 463*2b730f8bSJeremy L Thompson *flops = 0; 464*2b730f8bSJeremy L Thompson break; 4656e15d496SJeremy L Thompson } 4666e15d496SJeremy L Thompson } 4676e15d496SJeremy L Thompson 4686e15d496SJeremy L Thompson return CEED_ERROR_SUCCESS; 4696e15d496SJeremy L Thompson } 4706e15d496SJeremy L Thompson 4716e15d496SJeremy L Thompson /** 4727a982d89SJeremy L. Thompson @brief Get dimension for given CeedElemTopology 4737a982d89SJeremy L. Thompson 4747a982d89SJeremy L. Thompson @param topo CeedElemTopology 4757a982d89SJeremy L. Thompson @param[out] dim Variable to store dimension of topology 4767a982d89SJeremy L. Thompson 4777a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 4787a982d89SJeremy L. Thompson 4797a982d89SJeremy L. Thompson @ref Backend 4807a982d89SJeremy L. Thompson **/ 4817a982d89SJeremy L. Thompson int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) { 4827a982d89SJeremy L. Thompson *dim = (CeedInt)topo >> 16; 483e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 4847a982d89SJeremy L. Thompson } 4857a982d89SJeremy L. Thompson 4867a982d89SJeremy L. Thompson /** 4877a982d89SJeremy L. Thompson @brief Get CeedTensorContract of a CeedBasis 4887a982d89SJeremy L. Thompson 4897a982d89SJeremy L. Thompson @param basis CeedBasis 4907a982d89SJeremy L. Thompson @param[out] contract Variable to store CeedTensorContract 4917a982d89SJeremy L. Thompson 4927a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 4937a982d89SJeremy L. Thompson 4947a982d89SJeremy L. Thompson @ref Backend 4957a982d89SJeremy L. Thompson **/ 4967a982d89SJeremy L. Thompson int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) { 4977a982d89SJeremy L. Thompson *contract = basis->contract; 498e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 4997a982d89SJeremy L. Thompson } 5007a982d89SJeremy L. Thompson 5017a982d89SJeremy L. Thompson /** 5027a982d89SJeremy L. Thompson @brief Set CeedTensorContract of a CeedBasis 5037a982d89SJeremy L. Thompson 5047a982d89SJeremy L. Thompson @param[out] basis CeedBasis 5057a982d89SJeremy L. Thompson @param contract CeedTensorContract to set 5067a982d89SJeremy L. Thompson 5077a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 5087a982d89SJeremy L. Thompson 5097a982d89SJeremy L. Thompson @ref Backend 5107a982d89SJeremy L. Thompson **/ 51134359f16Sjeremylt int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract contract) { 51234359f16Sjeremylt basis->contract = contract; 513*2b730f8bSJeremy L Thompson CeedCall(CeedTensorContractReference(contract)); 514e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 5157a982d89SJeremy L. Thompson } 5167a982d89SJeremy L. Thompson 5177a982d89SJeremy L. Thompson /** 5187a982d89SJeremy L. Thompson @brief Return a reference implementation of matrix multiplication C = A B. 5197a982d89SJeremy L. Thompson Note, this is a reference implementation for CPU CeedScalar pointers 5207a982d89SJeremy L. Thompson that is not intended for high performance. 5217a982d89SJeremy L. Thompson 5227a982d89SJeremy L. Thompson @param ceed A Ceed context for error handling 523d1d35e2fSjeremylt @param[in] mat_A Row-major matrix A 524d1d35e2fSjeremylt @param[in] mat_B Row-major matrix B 525d1d35e2fSjeremylt @param[out] mat_C Row-major output matrix C 5267a982d89SJeremy L. Thompson @param m Number of rows of C 5277a982d89SJeremy L. Thompson @param n Number of columns of C 5287a982d89SJeremy L. Thompson @param kk Number of columns of A/rows of B 5297a982d89SJeremy L. Thompson 5307a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 5317a982d89SJeremy L. Thompson 5327a982d89SJeremy L. Thompson @ref Utility 5337a982d89SJeremy L. Thompson **/ 534*2b730f8bSJeremy L Thompson int CeedMatrixMatrixMultiply(Ceed ceed, const CeedScalar *mat_A, const CeedScalar *mat_B, CeedScalar *mat_C, CeedInt m, CeedInt n, CeedInt kk) { 535*2b730f8bSJeremy L Thompson for (CeedInt i = 0; i < m; i++) { 5367a982d89SJeremy L. Thompson for (CeedInt j = 0; j < n; j++) { 5377a982d89SJeremy L. Thompson CeedScalar sum = 0; 538*2b730f8bSJeremy L Thompson for (CeedInt k = 0; k < kk; k++) sum += mat_A[k + i * kk] * mat_B[j + k * n]; 539d1d35e2fSjeremylt mat_C[j + i * n] = sum; 5407a982d89SJeremy L. Thompson } 541*2b730f8bSJeremy L Thompson } 542e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 5437a982d89SJeremy L. Thompson } 5447a982d89SJeremy L. Thompson 5457a982d89SJeremy L. Thompson /// @} 5467a982d89SJeremy L. Thompson 5477a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 5487a982d89SJeremy L. Thompson /// CeedBasis Public API 5497a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 5507a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser 551d7b241e6Sjeremylt /// @{ 552d7b241e6Sjeremylt 553b11c1e72Sjeremylt /** 55495bb1877Svaleriabarra @brief Create a tensor-product basis for H^1 discretizations 555b11c1e72Sjeremylt 556b11c1e72Sjeremylt @param ceed A Ceed object where the CeedBasis will be created 557b11c1e72Sjeremylt @param dim Topological dimension 558d1d35e2fSjeremylt @param num_comp Number of field components (1 for scalar fields) 559d1d35e2fSjeremylt @param P_1d Number of nodes in one dimension 560d1d35e2fSjeremylt @param Q_1d Number of quadrature points in one dimension 561d1d35e2fSjeremylt @param interp_1d Row-major (Q_1d * P_1d) matrix expressing the values of nodal 562b11c1e72Sjeremylt basis functions at quadrature points 563d1d35e2fSjeremylt @param grad_1d Row-major (Q_1d * P_1d) matrix expressing derivatives of nodal 564b11c1e72Sjeremylt basis functions at quadrature points 565d1d35e2fSjeremylt @param q_ref_1d Array of length Q_1d holding the locations of quadrature points 566b11c1e72Sjeremylt on the 1D reference element [-1, 1] 567d1d35e2fSjeremylt @param q_weight_1d Array of length Q_1d holding the quadrature weights on the 568b11c1e72Sjeremylt reference element 569b11c1e72Sjeremylt @param[out] basis Address of the variable where the newly created 570b11c1e72Sjeremylt CeedBasis will be stored. 571b11c1e72Sjeremylt 572b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 573dfdf5a53Sjeremylt 5747a982d89SJeremy L. Thompson @ref User 575b11c1e72Sjeremylt **/ 576*2b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, 577*2b730f8bSJeremy L Thompson const CeedScalar *grad_1d, const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis *basis) { 5785fe0d4faSjeremylt if (!ceed->BasisCreateTensorH1) { 5795fe0d4faSjeremylt Ceed delegate; 580*2b730f8bSJeremy L Thompson CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis")); 5815fe0d4faSjeremylt 582*2b730f8bSJeremy L Thompson if (!delegate) { 583c042f62fSJeremy L Thompson // LCOV_EXCL_START 584*2b730f8bSJeremy L Thompson return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisCreateTensorH1"); 585c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 586*2b730f8bSJeremy L Thompson } 5875fe0d4faSjeremylt 588*2b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateTensorH1(delegate, dim, num_comp, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis)); 589e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 5905fe0d4faSjeremylt } 591e15f9bd0SJeremy L Thompson 592*2b730f8bSJeremy L Thompson if (dim < 1) { 593e15f9bd0SJeremy L Thompson // LCOV_EXCL_START 594*2b730f8bSJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis dimension must be a positive value"); 595e15f9bd0SJeremy L Thompson // LCOV_EXCL_STOP 596*2b730f8bSJeremy L Thompson } 597227444bfSJeremy L Thompson 598*2b730f8bSJeremy L Thompson if (num_comp < 1) { 599227444bfSJeremy L Thompson // LCOV_EXCL_START 600*2b730f8bSJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component"); 601227444bfSJeremy L Thompson // LCOV_EXCL_STOP 602*2b730f8bSJeremy L Thompson } 603227444bfSJeremy L Thompson 604*2b730f8bSJeremy L Thompson if (P_1d < 1) { 605227444bfSJeremy L Thompson // LCOV_EXCL_START 606*2b730f8bSJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node"); 607227444bfSJeremy L Thompson // LCOV_EXCL_STOP 608*2b730f8bSJeremy L Thompson } 609227444bfSJeremy L Thompson 610*2b730f8bSJeremy L Thompson if (Q_1d < 1) { 611227444bfSJeremy L Thompson // LCOV_EXCL_START 612*2b730f8bSJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point"); 613227444bfSJeremy L Thompson // LCOV_EXCL_STOP 614*2b730f8bSJeremy L Thompson } 615227444bfSJeremy L Thompson 616*2b730f8bSJeremy L Thompson CeedElemTopology topo = dim == 1 ? CEED_TOPOLOGY_LINE : dim == 2 ? CEED_TOPOLOGY_QUAD : CEED_TOPOLOGY_HEX; 617e15f9bd0SJeremy L Thompson 618*2b730f8bSJeremy L Thompson CeedCall(CeedCalloc(1, basis)); 619d7b241e6Sjeremylt (*basis)->ceed = ceed; 620*2b730f8bSJeremy L Thompson CeedCall(CeedReference(ceed)); 621d1d35e2fSjeremylt (*basis)->ref_count = 1; 622d1d35e2fSjeremylt (*basis)->tensor_basis = 1; 623d7b241e6Sjeremylt (*basis)->dim = dim; 624d99fa3c5SJeremy L Thompson (*basis)->topo = topo; 625d1d35e2fSjeremylt (*basis)->num_comp = num_comp; 626d1d35e2fSjeremylt (*basis)->P_1d = P_1d; 627d1d35e2fSjeremylt (*basis)->Q_1d = Q_1d; 628d1d35e2fSjeremylt (*basis)->P = CeedIntPow(P_1d, dim); 629d1d35e2fSjeremylt (*basis)->Q = CeedIntPow(Q_1d, dim); 63050c301a5SRezgar Shakeri (*basis)->Q_comp = 1; 63150c301a5SRezgar Shakeri (*basis)->basis_space = 1; // 1 for H^1 space 632*2b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q_1d, &(*basis)->q_ref_1d)); 633*2b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q_1d, &(*basis)->q_weight_1d)); 634ff3a0f91SJeremy L Thompson if (q_ref_1d) memcpy((*basis)->q_ref_1d, q_ref_1d, Q_1d * sizeof(q_ref_1d[0])); 635*2b730f8bSJeremy L Thompson if (q_weight_1d) memcpy((*basis)->q_weight_1d, q_weight_1d, Q_1d * sizeof(q_weight_1d[0])); 636*2b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->interp_1d)); 637*2b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->grad_1d)); 638*2b730f8bSJeremy L Thompson if (interp_1d) memcpy((*basis)->interp_1d, interp_1d, Q_1d * P_1d * sizeof(interp_1d[0])); 639ff3a0f91SJeremy L Thompson if (grad_1d) memcpy((*basis)->grad_1d, grad_1d, Q_1d * P_1d * sizeof(grad_1d[0])); 640*2b730f8bSJeremy L Thompson CeedCall(ceed->BasisCreateTensorH1(dim, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, *basis)); 641e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 642d7b241e6Sjeremylt } 643d7b241e6Sjeremylt 644b11c1e72Sjeremylt /** 64595bb1877Svaleriabarra @brief Create a tensor-product Lagrange basis 646b11c1e72Sjeremylt 647b11c1e72Sjeremylt @param ceed A Ceed object where the CeedBasis will be created 648b11c1e72Sjeremylt @param dim Topological dimension of element 649d1d35e2fSjeremylt @param num_comp Number of field components (1 for scalar fields) 650b11c1e72Sjeremylt @param P Number of Gauss-Lobatto nodes in one dimension. The 651b11c1e72Sjeremylt polynomial degree of the resulting Q_k element is k=P-1. 652b11c1e72Sjeremylt @param Q Number of quadrature points in one dimension. 653d1d35e2fSjeremylt @param quad_mode Distribution of the Q quadrature points (affects order of 654b11c1e72Sjeremylt accuracy for the quadrature) 655b11c1e72Sjeremylt @param[out] basis Address of the variable where the newly created 656b11c1e72Sjeremylt CeedBasis will be stored. 657b11c1e72Sjeremylt 658b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 659dfdf5a53Sjeremylt 6607a982d89SJeremy L. Thompson @ref User 661b11c1e72Sjeremylt **/ 662*2b730f8bSJeremy L Thompson int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P, CeedInt Q, CeedQuadMode quad_mode, CeedBasis *basis) { 663d7b241e6Sjeremylt // Allocate 664*2b730f8bSJeremy L Thompson int ierr = CEED_ERROR_SUCCESS, i, j, k; 665*2b730f8bSJeremy L Thompson CeedScalar c1, c2, c3, c4, dx, *nodes, *interp_1d, *grad_1d, *q_ref_1d, *q_weight_1d; 6664d537eeaSYohann 667*2b730f8bSJeremy L Thompson if (dim < 1) { 668c042f62fSJeremy L Thompson // LCOV_EXCL_START 669*2b730f8bSJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis dimension must be a positive value"); 670c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 671*2b730f8bSJeremy L Thompson } 6724d537eeaSYohann 673*2b730f8bSJeremy L Thompson if (num_comp < 1) { 674227444bfSJeremy L Thompson // LCOV_EXCL_START 675*2b730f8bSJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component"); 676227444bfSJeremy L Thompson // LCOV_EXCL_STOP 677*2b730f8bSJeremy L Thompson } 678227444bfSJeremy L Thompson 679*2b730f8bSJeremy L Thompson if (P < 1) { 680227444bfSJeremy L Thompson // LCOV_EXCL_START 681*2b730f8bSJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node"); 682227444bfSJeremy L Thompson // LCOV_EXCL_STOP 683*2b730f8bSJeremy L Thompson } 684227444bfSJeremy L Thompson 685*2b730f8bSJeremy L Thompson if (Q < 1) { 686227444bfSJeremy L Thompson // LCOV_EXCL_START 687*2b730f8bSJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point"); 688227444bfSJeremy L Thompson // LCOV_EXCL_STOP 689*2b730f8bSJeremy L Thompson } 690227444bfSJeremy L Thompson 691e15f9bd0SJeremy L Thompson // Get Nodes and Weights 692*2b730f8bSJeremy L Thompson CeedCall(CeedCalloc(P * Q, &interp_1d)); 693*2b730f8bSJeremy L Thompson CeedCall(CeedCalloc(P * Q, &grad_1d)); 694*2b730f8bSJeremy L Thompson CeedCall(CeedCalloc(P, &nodes)); 695*2b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q, &q_ref_1d)); 696*2b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q, &q_weight_1d)); 697*2b730f8bSJeremy L Thompson if (CeedLobattoQuadrature(P, nodes, NULL) != CEED_ERROR_SUCCESS) goto cleanup; 698d1d35e2fSjeremylt switch (quad_mode) { 699d7b241e6Sjeremylt case CEED_GAUSS: 700d1d35e2fSjeremylt ierr = CeedGaussQuadrature(Q, q_ref_1d, q_weight_1d); 701d7b241e6Sjeremylt break; 702d7b241e6Sjeremylt case CEED_GAUSS_LOBATTO: 703d1d35e2fSjeremylt ierr = CeedLobattoQuadrature(Q, q_ref_1d, q_weight_1d); 704d7b241e6Sjeremylt break; 705d7b241e6Sjeremylt } 706*2b730f8bSJeremy L Thompson if (ierr != CEED_ERROR_SUCCESS) goto cleanup; 707e15f9bd0SJeremy L Thompson 708d7b241e6Sjeremylt // Build B, D matrix 709d7b241e6Sjeremylt // Fornberg, 1998 710d7b241e6Sjeremylt for (i = 0; i < Q; i++) { 711d7b241e6Sjeremylt c1 = 1.0; 712d1d35e2fSjeremylt c3 = nodes[0] - q_ref_1d[i]; 713d1d35e2fSjeremylt interp_1d[i * P + 0] = 1.0; 714d7b241e6Sjeremylt for (j = 1; j < P; j++) { 715d7b241e6Sjeremylt c2 = 1.0; 716d7b241e6Sjeremylt c4 = c3; 717d1d35e2fSjeremylt c3 = nodes[j] - q_ref_1d[i]; 718d7b241e6Sjeremylt for (k = 0; k < j; k++) { 719d7b241e6Sjeremylt dx = nodes[j] - nodes[k]; 720d7b241e6Sjeremylt c2 *= dx; 721d7b241e6Sjeremylt if (k == j - 1) { 722d1d35e2fSjeremylt grad_1d[i * P + j] = c1 * (interp_1d[i * P + k] - c4 * grad_1d[i * P + k]) / c2; 723d1d35e2fSjeremylt interp_1d[i * P + j] = -c1 * c4 * interp_1d[i * P + k] / c2; 724d7b241e6Sjeremylt } 725d1d35e2fSjeremylt grad_1d[i * P + k] = (c3 * grad_1d[i * P + k] - interp_1d[i * P + k]) / dx; 726d1d35e2fSjeremylt interp_1d[i * P + k] = c3 * interp_1d[i * P + k] / dx; 727d7b241e6Sjeremylt } 728d7b241e6Sjeremylt c1 = c2; 729d7b241e6Sjeremylt } 730d7b241e6Sjeremylt } 7319ac7b42eSJeremy L Thompson // Pass to CeedBasisCreateTensorH1 732*2b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P, Q, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis)); 733e15f9bd0SJeremy L Thompson cleanup: 734*2b730f8bSJeremy L Thompson CeedCall(CeedFree(&interp_1d)); 735*2b730f8bSJeremy L Thompson CeedCall(CeedFree(&grad_1d)); 736*2b730f8bSJeremy L Thompson CeedCall(CeedFree(&nodes)); 737*2b730f8bSJeremy L Thompson CeedCall(CeedFree(&q_ref_1d)); 738*2b730f8bSJeremy L Thompson CeedCall(CeedFree(&q_weight_1d)); 739e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 740d7b241e6Sjeremylt } 741d7b241e6Sjeremylt 742b11c1e72Sjeremylt /** 74395bb1877Svaleriabarra @brief Create a non tensor-product basis for H^1 discretizations 744a8de75f0Sjeremylt 745a8de75f0Sjeremylt @param ceed A Ceed object where the CeedBasis will be created 746a8de75f0Sjeremylt @param topo Topology of element, e.g. hypercube, simplex, ect 747d1d35e2fSjeremylt @param num_comp Number of field components (1 for scalar fields) 748d1d35e2fSjeremylt @param num_nodes Total number of nodes 749d1d35e2fSjeremylt @param num_qpts Total number of quadrature points 750d1d35e2fSjeremylt @param interp Row-major (num_qpts * num_nodes) matrix expressing the values of 7518795c945Sjeremylt nodal basis functions at quadrature points 752d1d35e2fSjeremylt @param grad Row-major (num_qpts * dim * num_nodes) matrix expressing 7538795c945Sjeremylt derivatives of nodal basis functions at quadrature points 754d1d35e2fSjeremylt @param q_ref Array of length num_qpts holding the locations of quadrature 75550c301a5SRezgar Shakeri points on the reference element 756d1d35e2fSjeremylt @param q_weight Array of length num_qpts holding the quadrature weights on the 757a8de75f0Sjeremylt reference element 758a8de75f0Sjeremylt @param[out] basis Address of the variable where the newly created 759a8de75f0Sjeremylt CeedBasis will be stored. 760a8de75f0Sjeremylt 761a8de75f0Sjeremylt @return An error code: 0 - success, otherwise - failure 762a8de75f0Sjeremylt 7637a982d89SJeremy L. Thompson @ref User 764a8de75f0Sjeremylt **/ 765*2b730f8bSJeremy L Thompson int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 766*2b730f8bSJeremy L Thompson const CeedScalar *grad, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) { 767d1d35e2fSjeremylt CeedInt P = num_nodes, Q = num_qpts, dim = 0; 768a8de75f0Sjeremylt 7695fe0d4faSjeremylt if (!ceed->BasisCreateH1) { 7705fe0d4faSjeremylt Ceed delegate; 771*2b730f8bSJeremy L Thompson CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis")); 7725fe0d4faSjeremylt 773*2b730f8bSJeremy L Thompson if (!delegate) { 774c042f62fSJeremy L Thompson // LCOV_EXCL_START 775*2b730f8bSJeremy L Thompson return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisCreateH1"); 776c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 777*2b730f8bSJeremy L Thompson } 7785fe0d4faSjeremylt 779*2b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateH1(delegate, topo, num_comp, num_nodes, num_qpts, interp, grad, q_ref, q_weight, basis)); 780e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 7815fe0d4faSjeremylt } 7825fe0d4faSjeremylt 783*2b730f8bSJeremy L Thompson if (num_comp < 1) { 784227444bfSJeremy L Thompson // LCOV_EXCL_START 785*2b730f8bSJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component"); 786227444bfSJeremy L Thompson // LCOV_EXCL_STOP 787*2b730f8bSJeremy L Thompson } 788227444bfSJeremy L Thompson 789*2b730f8bSJeremy L Thompson if (num_nodes < 1) { 790227444bfSJeremy L Thompson // LCOV_EXCL_START 791*2b730f8bSJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node"); 792227444bfSJeremy L Thompson // LCOV_EXCL_STOP 793*2b730f8bSJeremy L Thompson } 794227444bfSJeremy L Thompson 795*2b730f8bSJeremy L Thompson if (num_qpts < 1) { 796227444bfSJeremy L Thompson // LCOV_EXCL_START 797*2b730f8bSJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point"); 798227444bfSJeremy L Thompson // LCOV_EXCL_STOP 799*2b730f8bSJeremy L Thompson } 800227444bfSJeremy L Thompson 801*2b730f8bSJeremy L Thompson CeedCall(CeedCalloc(1, basis)); 802a8de75f0Sjeremylt 803*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetTopologyDimension(topo, &dim)); 804a8de75f0Sjeremylt 805a8de75f0Sjeremylt (*basis)->ceed = ceed; 806*2b730f8bSJeremy L Thompson CeedCall(CeedReference(ceed)); 807d1d35e2fSjeremylt (*basis)->ref_count = 1; 808d1d35e2fSjeremylt (*basis)->tensor_basis = 0; 809a8de75f0Sjeremylt (*basis)->dim = dim; 810d99fa3c5SJeremy L Thompson (*basis)->topo = topo; 811d1d35e2fSjeremylt (*basis)->num_comp = num_comp; 812a8de75f0Sjeremylt (*basis)->P = P; 813a8de75f0Sjeremylt (*basis)->Q = Q; 81450c301a5SRezgar Shakeri (*basis)->Q_comp = 1; 81550c301a5SRezgar Shakeri (*basis)->basis_space = 1; // 1 for H^1 space 816*2b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q * dim, &(*basis)->q_ref_1d)); 817*2b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q, &(*basis)->q_weight_1d)); 818ff3a0f91SJeremy L Thompson if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0])); 819ff3a0f91SJeremy L Thompson if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0])); 820*2b730f8bSJeremy L Thompson CeedCall(CeedCalloc(Q * P, &(*basis)->interp)); 821*2b730f8bSJeremy L Thompson CeedCall(CeedCalloc(dim * Q * P, &(*basis)->grad)); 822ff3a0f91SJeremy L Thompson if (interp) memcpy((*basis)->interp, interp, Q * P * sizeof(interp[0])); 823ff3a0f91SJeremy L Thompson if (grad) memcpy((*basis)->grad, grad, dim * Q * P * sizeof(grad[0])); 824*2b730f8bSJeremy L Thompson CeedCall(ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, q_ref, q_weight, *basis)); 825e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 826a8de75f0Sjeremylt } 827a8de75f0Sjeremylt 828a8de75f0Sjeremylt /** 82950c301a5SRezgar Shakeri @brief Create a non tensor-product basis for H(div) discretizations 83050c301a5SRezgar Shakeri 83150c301a5SRezgar Shakeri @param ceed A Ceed object where the CeedBasis will be created 83250c301a5SRezgar Shakeri @param topo Topology of element (`CEED_TOPOLOGY_QUAD`, `CEED_TOPOLOGY_PRISM`, etc.), 83350c301a5SRezgar Shakeri dimension of which is used in some array sizes below 83450c301a5SRezgar Shakeri @param num_comp Number of components (usually 1 for vectors in H(div) bases) 83550c301a5SRezgar Shakeri @param num_nodes Total number of nodes (dofs per element) 83650c301a5SRezgar Shakeri @param num_qpts Total number of quadrature points 83750c301a5SRezgar Shakeri @param interp Row-major (dim*num_qpts * num_nodes) matrix expressing the values of 83850c301a5SRezgar Shakeri nodal basis functions at quadrature points 83950c301a5SRezgar Shakeri @param div Row-major (num_qpts * num_nodes) matrix expressing 84050c301a5SRezgar Shakeri divergence of nodal basis functions at quadrature points 84150c301a5SRezgar Shakeri @param q_ref Array of length num_qpts holding the locations of quadrature 84250c301a5SRezgar Shakeri points on the reference element 84350c301a5SRezgar Shakeri @param q_weight Array of length num_qpts holding the quadrature weights on the 84450c301a5SRezgar Shakeri reference element 84550c301a5SRezgar Shakeri @param[out] basis Address of the variable where the newly created 84650c301a5SRezgar Shakeri CeedBasis will be stored. 84750c301a5SRezgar Shakeri 84850c301a5SRezgar Shakeri @return An error code: 0 - success, otherwise - failure 84950c301a5SRezgar Shakeri 85050c301a5SRezgar Shakeri @ref User 85150c301a5SRezgar Shakeri **/ 852*2b730f8bSJeremy L Thompson int CeedBasisCreateHdiv(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 853*2b730f8bSJeremy L Thompson const CeedScalar *div, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) { 85450c301a5SRezgar Shakeri CeedInt Q = num_qpts, P = num_nodes, dim = 0; 855*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetTopologyDimension(topo, &dim)); 85650c301a5SRezgar Shakeri if (!ceed->BasisCreateHdiv) { 85750c301a5SRezgar Shakeri Ceed delegate; 858*2b730f8bSJeremy L Thompson CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis")); 85950c301a5SRezgar Shakeri 860*2b730f8bSJeremy L Thompson if (!delegate) { 86150c301a5SRezgar Shakeri // LCOV_EXCL_START 862*2b730f8bSJeremy L Thompson return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateHdiv"); 86350c301a5SRezgar Shakeri // LCOV_EXCL_STOP 864*2b730f8bSJeremy L Thompson } 86550c301a5SRezgar Shakeri 866*2b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateHdiv(delegate, topo, num_comp, num_nodes, num_qpts, interp, div, q_ref, q_weight, basis)); 86750c301a5SRezgar Shakeri return CEED_ERROR_SUCCESS; 86850c301a5SRezgar Shakeri } 86950c301a5SRezgar Shakeri 870*2b730f8bSJeremy L Thompson if (num_comp < 1) { 871227444bfSJeremy L Thompson // LCOV_EXCL_START 872*2b730f8bSJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component"); 873227444bfSJeremy L Thompson // LCOV_EXCL_STOP 874*2b730f8bSJeremy L Thompson } 875227444bfSJeremy L Thompson 876*2b730f8bSJeremy L Thompson if (num_nodes < 1) { 877227444bfSJeremy L Thompson // LCOV_EXCL_START 878*2b730f8bSJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node"); 879227444bfSJeremy L Thompson // LCOV_EXCL_STOP 880*2b730f8bSJeremy L Thompson } 881227444bfSJeremy L Thompson 882*2b730f8bSJeremy L Thompson if (num_qpts < 1) { 883227444bfSJeremy L Thompson // LCOV_EXCL_START 884*2b730f8bSJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point"); 885227444bfSJeremy L Thompson // LCOV_EXCL_STOP 886*2b730f8bSJeremy L Thompson } 887227444bfSJeremy L Thompson 888*2b730f8bSJeremy L Thompson CeedCall(CeedCalloc(1, basis)); 88950c301a5SRezgar Shakeri 89050c301a5SRezgar Shakeri (*basis)->ceed = ceed; 891*2b730f8bSJeremy L Thompson CeedCall(CeedReference(ceed)); 89250c301a5SRezgar Shakeri (*basis)->ref_count = 1; 89350c301a5SRezgar Shakeri (*basis)->tensor_basis = 0; 89450c301a5SRezgar Shakeri (*basis)->dim = dim; 89550c301a5SRezgar Shakeri (*basis)->topo = topo; 89650c301a5SRezgar Shakeri (*basis)->num_comp = num_comp; 89750c301a5SRezgar Shakeri (*basis)->P = P; 89850c301a5SRezgar Shakeri (*basis)->Q = Q; 89950c301a5SRezgar Shakeri (*basis)->Q_comp = dim; 90050c301a5SRezgar Shakeri (*basis)->basis_space = 2; // 2 for H(div) space 901*2b730f8bSJeremy L Thompson CeedCall(CeedMalloc(Q * dim, &(*basis)->q_ref_1d)); 902*2b730f8bSJeremy L Thompson CeedCall(CeedMalloc(Q, &(*basis)->q_weight_1d)); 90350c301a5SRezgar Shakeri if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0])); 90450c301a5SRezgar Shakeri if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0])); 905*2b730f8bSJeremy L Thompson CeedCall(CeedMalloc(dim * Q * P, &(*basis)->interp)); 906*2b730f8bSJeremy L Thompson CeedCall(CeedMalloc(Q * P, &(*basis)->div)); 90750c301a5SRezgar Shakeri if (interp) memcpy((*basis)->interp, interp, dim * Q * P * sizeof(interp[0])); 90850c301a5SRezgar Shakeri if (div) memcpy((*basis)->div, div, Q * P * sizeof(div[0])); 909*2b730f8bSJeremy L Thompson CeedCall(ceed->BasisCreateHdiv(topo, dim, P, Q, interp, div, q_ref, q_weight, *basis)); 91050c301a5SRezgar Shakeri return CEED_ERROR_SUCCESS; 91150c301a5SRezgar Shakeri } 91250c301a5SRezgar Shakeri 91350c301a5SRezgar Shakeri /** 914f113e5dcSJeremy L Thompson @brief Create a CeedBasis for projection from the nodes of `basis_from` 91514556e63SJeremy L Thompson to the nodes of `basis_to`. Only `CEED_EVAL_INTERP` and 91614556e63SJeremy L Thompson `CEED_EVAL_GRAD` will be valid for the new basis, `basis_project`. 91714556e63SJeremy L Thompson The interpolation is given by `interp_project = interp_to^+ * interp_from`, 91814556e63SJeremy L Thompson where the pesudoinverse `interp_to^+` is given by QR factorization. 91914556e63SJeremy L Thompson The gradient is given by `grad_project = interp_to^+ * grad_from`. 920f113e5dcSJeremy L Thompson Note: `basis_from` and `basis_to` must have compatible quadrature 921f113e5dcSJeremy L Thompson spaces. 922446e7af4SJeremy L Thompson Note: `basis_project` will have the same number of components as 923446e7af4SJeremy L Thompson `basis_from`, regardless of the number of components that 924446e7af4SJeremy L Thompson `basis_to` has. If `basis_from` has 3 components and `basis_to` 925446e7af4SJeremy L Thompson has 5 components, then `basis_project` will have 3 components. 926f113e5dcSJeremy L Thompson 927f113e5dcSJeremy L Thompson @param[in] basis_from CeedBasis to prolong from 928446e7af4SJeremy L Thompson @param[in] basis_to CeedBasis to prolong to 929f113e5dcSJeremy L Thompson @param[out] basis_project Address of the variable where the newly created 930f113e5dcSJeremy L Thompson CeedBasis will be stored. 931f113e5dcSJeremy L Thompson 932f113e5dcSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 933f113e5dcSJeremy L Thompson 934f113e5dcSJeremy L Thompson @ref User 935f113e5dcSJeremy L Thompson **/ 936*2b730f8bSJeremy L Thompson int CeedBasisCreateProjection(CeedBasis basis_from, CeedBasis basis_to, CeedBasis *basis_project) { 937f113e5dcSJeremy L Thompson Ceed ceed; 938*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetCeed(basis_to, &ceed)); 939f113e5dcSJeremy L Thompson 940f113e5dcSJeremy L Thompson // Create projectior matrix 94114556e63SJeremy L Thompson CeedScalar *interp_project, *grad_project; 942*2b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateProjectionMatrices(basis_from, basis_to, &interp_project, &grad_project)); 943f113e5dcSJeremy L Thompson 944f113e5dcSJeremy L Thompson // Build basis 945f113e5dcSJeremy L Thompson bool is_tensor; 946f113e5dcSJeremy L Thompson CeedInt dim, num_comp; 94714556e63SJeremy L Thompson CeedScalar *q_ref, *q_weight; 948*2b730f8bSJeremy L Thompson CeedCall(CeedBasisIsTensor(basis_to, &is_tensor)); 949*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetDimension(basis_to, &dim)); 950*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumComponents(basis_from, &num_comp)); 951f113e5dcSJeremy L Thompson if (is_tensor) { 952f113e5dcSJeremy L Thompson CeedInt P_1d_to, P_1d_from; 953*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_1d_from)); 954*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_1d_to)); 955*2b730f8bSJeremy L Thompson CeedCall(CeedCalloc(P_1d_to, &q_ref)); 956*2b730f8bSJeremy L Thompson CeedCall(CeedCalloc(P_1d_to, &q_weight)); 957*2b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P_1d_from, P_1d_to, interp_project, grad_project, q_ref, q_weight, basis_project)); 958f113e5dcSJeremy L Thompson } else { 959f113e5dcSJeremy L Thompson CeedElemTopology topo; 960*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetTopology(basis_to, &topo)); 961f113e5dcSJeremy L Thompson CeedInt num_nodes_to, num_nodes_from; 962*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis_from, &num_nodes_from)); 963*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis_to, &num_nodes_to)); 964*2b730f8bSJeremy L Thompson CeedCall(CeedCalloc(num_nodes_to * dim, &q_ref)); 965*2b730f8bSJeremy L Thompson CeedCall(CeedCalloc(num_nodes_to, &q_weight)); 966*2b730f8bSJeremy L Thompson CeedCall(CeedBasisCreateH1(ceed, topo, num_comp, num_nodes_from, num_nodes_to, interp_project, grad_project, q_ref, q_weight, basis_project)); 967f113e5dcSJeremy L Thompson } 968f113e5dcSJeremy L Thompson 969f113e5dcSJeremy L Thompson // Cleanup 970*2b730f8bSJeremy L Thompson CeedCall(CeedFree(&interp_project)); 971*2b730f8bSJeremy L Thompson CeedCall(CeedFree(&grad_project)); 972*2b730f8bSJeremy L Thompson CeedCall(CeedFree(&q_ref)); 973*2b730f8bSJeremy L Thompson CeedCall(CeedFree(&q_weight)); 974f113e5dcSJeremy L Thompson 975f113e5dcSJeremy L Thompson return CEED_ERROR_SUCCESS; 976f113e5dcSJeremy L Thompson } 977f113e5dcSJeremy L Thompson 978f113e5dcSJeremy L Thompson /** 9799560d06aSjeremylt @brief Copy the pointer to a CeedBasis. Both pointers should 9809560d06aSjeremylt be destroyed with `CeedBasisDestroy()`; 9819560d06aSjeremylt Note: If `*basis_copy` is non-NULL, then it is assumed that 9829560d06aSjeremylt `*basis_copy` is a pointer to a CeedBasis. This CeedBasis 9839560d06aSjeremylt will be destroyed if `*basis_copy` is the only 9849560d06aSjeremylt reference to this CeedBasis. 9859560d06aSjeremylt 9869560d06aSjeremylt @param basis CeedBasis to copy reference to 9879560d06aSjeremylt @param[out] basis_copy Variable to store copied reference 9889560d06aSjeremylt 9899560d06aSjeremylt @return An error code: 0 - success, otherwise - failure 9909560d06aSjeremylt 9919560d06aSjeremylt @ref User 9929560d06aSjeremylt **/ 9939560d06aSjeremylt int CeedBasisReferenceCopy(CeedBasis basis, CeedBasis *basis_copy) { 994*2b730f8bSJeremy L Thompson CeedCall(CeedBasisReference(basis)); 995*2b730f8bSJeremy L Thompson CeedCall(CeedBasisDestroy(basis_copy)); 9969560d06aSjeremylt *basis_copy = basis; 9979560d06aSjeremylt return CEED_ERROR_SUCCESS; 9989560d06aSjeremylt } 9999560d06aSjeremylt 10009560d06aSjeremylt /** 10017a982d89SJeremy L. Thompson @brief View a CeedBasis 10027a982d89SJeremy L. Thompson 10037a982d89SJeremy L. Thompson @param basis CeedBasis to view 10047a982d89SJeremy L. Thompson @param stream Stream to view to, e.g., stdout 10057a982d89SJeremy L. Thompson 10067a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 10077a982d89SJeremy L. Thompson 10087a982d89SJeremy L. Thompson @ref User 10097a982d89SJeremy L. Thompson **/ 10107a982d89SJeremy L. Thompson int CeedBasisView(CeedBasis basis, FILE *stream) { 101150c301a5SRezgar Shakeri CeedFESpace FE_space = basis->basis_space; 101250c301a5SRezgar Shakeri CeedElemTopology topo = basis->topo; 1013*2b730f8bSJeremy L Thompson 101450c301a5SRezgar Shakeri // Print FE space and element topology of the basis 1015d1d35e2fSjeremylt if (basis->tensor_basis) { 1016*2b730f8bSJeremy L Thompson fprintf(stream, "CeedBasis (%s on a %s element): dim=%" CeedInt_FMT " P=%" CeedInt_FMT " Q=%" CeedInt_FMT "\n", CeedFESpaces[FE_space], 1017*2b730f8bSJeremy L Thompson CeedElemTopologies[topo], basis->dim, basis->P_1d, basis->Q_1d); 101850c301a5SRezgar Shakeri } else { 1019*2b730f8bSJeremy L Thompson fprintf(stream, "CeedBasis (%s on a %s element): dim=%" CeedInt_FMT " P=%" CeedInt_FMT " Q=%" CeedInt_FMT "\n", CeedFESpaces[FE_space], 1020*2b730f8bSJeremy L Thompson CeedElemTopologies[topo], basis->dim, basis->P, basis->Q); 102150c301a5SRezgar Shakeri } 102250c301a5SRezgar Shakeri // Print quadrature data, interpolation/gradient/divergene/curl of the basis 102350c301a5SRezgar Shakeri if (basis->tensor_basis) { // tensor basis 1024*2b730f8bSJeremy L Thompson CeedCall(CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_ref_1d, stream)); 1025*2b730f8bSJeremy L Thompson CeedCall(CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_weight_1d, stream)); 1026*2b730f8bSJeremy L Thompson CeedCall(CeedScalarView("interp1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, basis->interp_1d, stream)); 1027*2b730f8bSJeremy L Thompson CeedCall(CeedScalarView("grad1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, basis->grad_1d, stream)); 102850c301a5SRezgar Shakeri } else { // non-tensor basis 1029*2b730f8bSJeremy L Thompson CeedCall(CeedScalarView("qref", "\t% 12.8f", 1, basis->Q * basis->dim, basis->q_ref_1d, stream)); 1030*2b730f8bSJeremy L Thompson CeedCall(CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->q_weight_1d, stream)); 1031*2b730f8bSJeremy L Thompson CeedCall(CeedScalarView("interp", "\t% 12.8f", basis->Q_comp * basis->Q, basis->P, basis->interp, stream)); 103250c301a5SRezgar Shakeri if (basis->grad) { 1033*2b730f8bSJeremy L Thompson CeedCall(CeedScalarView("grad", "\t% 12.8f", basis->dim * basis->Q, basis->P, basis->grad, stream)); 10347a982d89SJeremy L. Thompson } 103550c301a5SRezgar Shakeri if (basis->div) { 1036*2b730f8bSJeremy L Thompson CeedCall(CeedScalarView("div", "\t% 12.8f", basis->Q, basis->P, basis->div, stream)); 103750c301a5SRezgar Shakeri } 103850c301a5SRezgar Shakeri } 1039e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 10407a982d89SJeremy L. Thompson } 10417a982d89SJeremy L. Thompson 10427a982d89SJeremy L. Thompson /** 10437a982d89SJeremy L. Thompson @brief Apply basis evaluation from nodes to quadrature points or vice versa 10447a982d89SJeremy L. Thompson 10457a982d89SJeremy L. Thompson @param basis CeedBasis to evaluate 1046d1d35e2fSjeremylt @param num_elem The number of elements to apply the basis evaluation to; 10477a982d89SJeremy L. Thompson the backend will specify the ordering in 10484cc79fe7SJed Brown CeedElemRestrictionCreateBlocked() 1049d1d35e2fSjeremylt @param t_mode \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature 10507a982d89SJeremy L. Thompson points, \ref CEED_TRANSPOSE to apply the transpose, mapping 10517a982d89SJeremy L. Thompson from quadrature points to nodes 1052d1d35e2fSjeremylt @param eval_mode \ref CEED_EVAL_NONE to use values directly, 10537a982d89SJeremy L. Thompson \ref CEED_EVAL_INTERP to use interpolated values, 10547a982d89SJeremy L. Thompson \ref CEED_EVAL_GRAD to use gradients, 10557a982d89SJeremy L. Thompson \ref CEED_EVAL_WEIGHT to use quadrature weights. 10567a982d89SJeremy L. Thompson @param[in] u Input CeedVector 10577a982d89SJeremy L. Thompson @param[out] v Output CeedVector 10587a982d89SJeremy L. Thompson 10597a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 10607a982d89SJeremy L. Thompson 10617a982d89SJeremy L. Thompson @ref User 10627a982d89SJeremy L. Thompson **/ 1063*2b730f8bSJeremy L Thompson int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, CeedVector v) { 10641f9221feSJeremy L Thompson CeedSize u_length = 0, v_length; 10651f9221feSJeremy L Thompson CeedInt dim, num_comp, num_nodes, num_qpts; 1066*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetDimension(basis, &dim)); 1067*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumComponents(basis, &num_comp)); 1068*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumNodes(basis, &num_nodes)); 1069*2b730f8bSJeremy L Thompson CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts)); 1070*2b730f8bSJeremy L Thompson CeedCall(CeedVectorGetLength(v, &v_length)); 10717a982d89SJeremy L. Thompson if (u) { 1072*2b730f8bSJeremy L Thompson CeedCall(CeedVectorGetLength(u, &u_length)); 10737a982d89SJeremy L. Thompson } 10747a982d89SJeremy L. Thompson 1075*2b730f8bSJeremy L Thompson if (!basis->Apply) { 1076e15f9bd0SJeremy L Thompson // LCOV_EXCL_START 1077*2b730f8bSJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisApply"); 1078e15f9bd0SJeremy L Thompson // LCOV_EXCL_STOP 1079*2b730f8bSJeremy L Thompson } 1080e15f9bd0SJeremy L Thompson 1081e15f9bd0SJeremy L Thompson // Check compatibility of topological and geometrical dimensions 1082*2b730f8bSJeremy L Thompson if ((t_mode == CEED_TRANSPOSE && (v_length % num_nodes != 0 || u_length % num_qpts != 0)) || 1083*2b730f8bSJeremy L Thompson (t_mode == CEED_NOTRANSPOSE && (u_length % num_nodes != 0 || v_length % num_qpts != 0))) { 10848229195eSjeremylt // LCOV_EXCL_START 1085*2b730f8bSJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_DIMENSION, "Length of input/output vectors incompatible with basis dimensions"); 10868229195eSjeremylt // LCOV_EXCL_STOP 1087*2b730f8bSJeremy L Thompson } 10887a982d89SJeremy L. Thompson 1089e15f9bd0SJeremy L Thompson // Check vector lengths to prevent out of bounds issues 1090d1d35e2fSjeremylt bool bad_dims = false; 1091d1d35e2fSjeremylt switch (eval_mode) { 1092e15f9bd0SJeremy L Thompson case CEED_EVAL_NONE: 1093*2b730f8bSJeremy L Thompson case CEED_EVAL_INTERP: 1094*2b730f8bSJeremy L Thompson bad_dims = ((t_mode == CEED_TRANSPOSE && (u_length < num_elem * num_comp * num_qpts || v_length < num_elem * num_comp * num_nodes)) || 1095*2b730f8bSJeremy L Thompson (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem * num_qpts * num_comp || u_length < num_elem * num_comp * num_nodes))); 1096e15f9bd0SJeremy L Thompson break; 1097*2b730f8bSJeremy L Thompson case CEED_EVAL_GRAD: 1098*2b730f8bSJeremy L Thompson bad_dims = ((t_mode == CEED_TRANSPOSE && (u_length < num_elem * num_comp * num_qpts * dim || v_length < num_elem * num_comp * num_nodes)) || 1099*2b730f8bSJeremy L Thompson (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem * num_qpts * num_comp * dim || u_length < num_elem * num_comp * num_nodes))); 1100e15f9bd0SJeremy L Thompson break; 1101e15f9bd0SJeremy L Thompson case CEED_EVAL_WEIGHT: 1102d1d35e2fSjeremylt bad_dims = v_length < num_elem * num_qpts; 1103e15f9bd0SJeremy L Thompson break; 1104e15f9bd0SJeremy L Thompson // LCOV_EXCL_START 1105*2b730f8bSJeremy L Thompson case CEED_EVAL_DIV: 1106*2b730f8bSJeremy L Thompson bad_dims = ((t_mode == CEED_TRANSPOSE && (u_length < num_elem * num_comp * num_qpts || v_length < num_elem * num_comp * num_nodes)) || 1107*2b730f8bSJeremy L Thompson (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem * num_qpts * num_comp || u_length < num_elem * num_comp * num_nodes))); 1108e15f9bd0SJeremy L Thompson break; 1109*2b730f8bSJeremy L Thompson case CEED_EVAL_CURL: 1110*2b730f8bSJeremy L Thompson bad_dims = ((t_mode == CEED_TRANSPOSE && (u_length < num_elem * num_comp * num_qpts || v_length < num_elem * num_comp * num_nodes)) || 1111*2b730f8bSJeremy L Thompson (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem * num_qpts * num_comp || u_length < num_elem * num_comp * num_nodes))); 1112e15f9bd0SJeremy L Thompson break; 1113e15f9bd0SJeremy L Thompson // LCOV_EXCL_STOP 1114e15f9bd0SJeremy L Thompson } 1115*2b730f8bSJeremy L Thompson if (bad_dims) { 1116e15f9bd0SJeremy L Thompson // LCOV_EXCL_START 1117*2b730f8bSJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_DIMENSION, "Input/output vectors too short for basis and evaluation mode"); 1118e15f9bd0SJeremy L Thompson // LCOV_EXCL_STOP 1119*2b730f8bSJeremy L Thompson } 1120e15f9bd0SJeremy L Thompson 1121*2b730f8bSJeremy L Thompson CeedCall(basis->Apply(basis, num_elem, t_mode, eval_mode, u, v)); 1122e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 11237a982d89SJeremy L. Thompson } 11247a982d89SJeremy L. Thompson 11257a982d89SJeremy L. Thompson /** 1126b7c9bbdaSJeremy L Thompson @brief Get Ceed associated with a CeedBasis 1127b7c9bbdaSJeremy L Thompson 1128b7c9bbdaSJeremy L Thompson @param basis CeedBasis 1129b7c9bbdaSJeremy L Thompson @param[out] ceed Variable to store Ceed 1130b7c9bbdaSJeremy L Thompson 1131b7c9bbdaSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 1132b7c9bbdaSJeremy L Thompson 1133b7c9bbdaSJeremy L Thompson @ref Advanced 1134b7c9bbdaSJeremy L Thompson **/ 1135b7c9bbdaSJeremy L Thompson int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) { 1136b7c9bbdaSJeremy L Thompson *ceed = basis->ceed; 1137b7c9bbdaSJeremy L Thompson return CEED_ERROR_SUCCESS; 1138b7c9bbdaSJeremy L Thompson } 1139b7c9bbdaSJeremy L Thompson 1140b7c9bbdaSJeremy L Thompson /** 11419d007619Sjeremylt @brief Get dimension for given CeedBasis 11429d007619Sjeremylt 11439d007619Sjeremylt @param basis CeedBasis 11449d007619Sjeremylt @param[out] dim Variable to store dimension of basis 11459d007619Sjeremylt 11469d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 11479d007619Sjeremylt 1148b7c9bbdaSJeremy L Thompson @ref Advanced 11499d007619Sjeremylt **/ 11509d007619Sjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) { 11519d007619Sjeremylt *dim = basis->dim; 1152e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 11539d007619Sjeremylt } 11549d007619Sjeremylt 11559d007619Sjeremylt /** 1156d99fa3c5SJeremy L Thompson @brief Get topology for given CeedBasis 1157d99fa3c5SJeremy L Thompson 1158d99fa3c5SJeremy L Thompson @param basis CeedBasis 1159d99fa3c5SJeremy L Thompson @param[out] topo Variable to store topology of basis 1160d99fa3c5SJeremy L Thompson 1161d99fa3c5SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 1162d99fa3c5SJeremy L Thompson 1163b7c9bbdaSJeremy L Thompson @ref Advanced 1164d99fa3c5SJeremy L Thompson **/ 1165d99fa3c5SJeremy L Thompson int CeedBasisGetTopology(CeedBasis basis, CeedElemTopology *topo) { 1166d99fa3c5SJeremy L Thompson *topo = basis->topo; 1167e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1168d99fa3c5SJeremy L Thompson } 1169d99fa3c5SJeremy L Thompson 1170d99fa3c5SJeremy L Thompson /** 117150c301a5SRezgar Shakeri @brief Get number of Q-vector components for given CeedBasis 117250c301a5SRezgar Shakeri 117350c301a5SRezgar Shakeri @param basis CeedBasis 117450c301a5SRezgar Shakeri @param[out] Q_comp Variable to store number of Q-vector components of basis 117550c301a5SRezgar Shakeri 117650c301a5SRezgar Shakeri @return An error code: 0 - success, otherwise - failure 117750c301a5SRezgar Shakeri 117850c301a5SRezgar Shakeri @ref Advanced 117950c301a5SRezgar Shakeri **/ 118050c301a5SRezgar Shakeri int CeedBasisGetNumQuadratureComponents(CeedBasis basis, CeedInt *Q_comp) { 118150c301a5SRezgar Shakeri *Q_comp = basis->Q_comp; 118250c301a5SRezgar Shakeri return CEED_ERROR_SUCCESS; 118350c301a5SRezgar Shakeri } 118450c301a5SRezgar Shakeri 118550c301a5SRezgar Shakeri /** 11869d007619Sjeremylt @brief Get number of components for given CeedBasis 11879d007619Sjeremylt 11889d007619Sjeremylt @param basis CeedBasis 1189d1d35e2fSjeremylt @param[out] num_comp Variable to store number of components of basis 11909d007619Sjeremylt 11919d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 11929d007619Sjeremylt 1193b7c9bbdaSJeremy L Thompson @ref Advanced 11949d007619Sjeremylt **/ 1195d1d35e2fSjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *num_comp) { 1196d1d35e2fSjeremylt *num_comp = basis->num_comp; 1197e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 11989d007619Sjeremylt } 11999d007619Sjeremylt 12009d007619Sjeremylt /** 12019d007619Sjeremylt @brief Get total number of nodes (in dim dimensions) of a CeedBasis 12029d007619Sjeremylt 12039d007619Sjeremylt @param basis CeedBasis 12049d007619Sjeremylt @param[out] P Variable to store number of nodes 12059d007619Sjeremylt 12069d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 12079d007619Sjeremylt 12089d007619Sjeremylt @ref Utility 12099d007619Sjeremylt **/ 12109d007619Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) { 12119d007619Sjeremylt *P = basis->P; 1212e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 12139d007619Sjeremylt } 12149d007619Sjeremylt 12159d007619Sjeremylt /** 12169d007619Sjeremylt @brief Get total number of nodes (in 1 dimension) of a CeedBasis 12179d007619Sjeremylt 12189d007619Sjeremylt @param basis CeedBasis 1219d1d35e2fSjeremylt @param[out] P_1d Variable to store number of nodes 12209d007619Sjeremylt 12219d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 12229d007619Sjeremylt 1223b7c9bbdaSJeremy L Thompson @ref Advanced 12249d007619Sjeremylt **/ 1225d1d35e2fSjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P_1d) { 1226*2b730f8bSJeremy L Thompson if (!basis->tensor_basis) { 12279d007619Sjeremylt // LCOV_EXCL_START 1228*2b730f8bSJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_MINOR, "Cannot supply P_1d for non-tensor basis"); 12299d007619Sjeremylt // LCOV_EXCL_STOP 1230*2b730f8bSJeremy L Thompson } 12319d007619Sjeremylt 1232d1d35e2fSjeremylt *P_1d = basis->P_1d; 1233e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 12349d007619Sjeremylt } 12359d007619Sjeremylt 12369d007619Sjeremylt /** 12379d007619Sjeremylt @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis 12389d007619Sjeremylt 12399d007619Sjeremylt @param basis CeedBasis 12409d007619Sjeremylt @param[out] Q Variable to store number of quadrature points 12419d007619Sjeremylt 12429d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 12439d007619Sjeremylt 12449d007619Sjeremylt @ref Utility 12459d007619Sjeremylt **/ 12469d007619Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) { 12479d007619Sjeremylt *Q = basis->Q; 1248e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 12499d007619Sjeremylt } 12509d007619Sjeremylt 12519d007619Sjeremylt /** 12529d007619Sjeremylt @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis 12539d007619Sjeremylt 12549d007619Sjeremylt @param basis CeedBasis 1255d1d35e2fSjeremylt @param[out] Q_1d Variable to store number of quadrature points 12569d007619Sjeremylt 12579d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 12589d007619Sjeremylt 1259b7c9bbdaSJeremy L Thompson @ref Advanced 12609d007619Sjeremylt **/ 1261d1d35e2fSjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q_1d) { 1262*2b730f8bSJeremy L Thompson if (!basis->tensor_basis) { 12639d007619Sjeremylt // LCOV_EXCL_START 1264*2b730f8bSJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_MINOR, "Cannot supply Q_1d for non-tensor basis"); 12659d007619Sjeremylt // LCOV_EXCL_STOP 1266*2b730f8bSJeremy L Thompson } 12679d007619Sjeremylt 1268d1d35e2fSjeremylt *Q_1d = basis->Q_1d; 1269e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 12709d007619Sjeremylt } 12719d007619Sjeremylt 12729d007619Sjeremylt /** 12739d007619Sjeremylt @brief Get reference coordinates of quadrature points (in dim dimensions) 12749d007619Sjeremylt of a CeedBasis 12759d007619Sjeremylt 12769d007619Sjeremylt @param basis CeedBasis 1277d1d35e2fSjeremylt @param[out] q_ref Variable to store reference coordinates of quadrature points 12789d007619Sjeremylt 12799d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 12809d007619Sjeremylt 1281b7c9bbdaSJeremy L Thompson @ref Advanced 12829d007619Sjeremylt **/ 1283d1d35e2fSjeremylt int CeedBasisGetQRef(CeedBasis basis, const CeedScalar **q_ref) { 1284d1d35e2fSjeremylt *q_ref = basis->q_ref_1d; 1285e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 12869d007619Sjeremylt } 12879d007619Sjeremylt 12889d007619Sjeremylt /** 12899d007619Sjeremylt @brief Get quadrature weights of quadrature points (in dim dimensions) 12909d007619Sjeremylt of a CeedBasis 12919d007619Sjeremylt 12929d007619Sjeremylt @param basis CeedBasis 1293d1d35e2fSjeremylt @param[out] q_weight Variable to store quadrature weights 12949d007619Sjeremylt 12959d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 12969d007619Sjeremylt 1297b7c9bbdaSJeremy L Thompson @ref Advanced 12989d007619Sjeremylt **/ 1299d1d35e2fSjeremylt int CeedBasisGetQWeights(CeedBasis basis, const CeedScalar **q_weight) { 1300d1d35e2fSjeremylt *q_weight = basis->q_weight_1d; 1301e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 13029d007619Sjeremylt } 13039d007619Sjeremylt 13049d007619Sjeremylt /** 13059d007619Sjeremylt @brief Get interpolation matrix of a CeedBasis 13069d007619Sjeremylt 13079d007619Sjeremylt @param basis CeedBasis 13089d007619Sjeremylt @param[out] interp Variable to store interpolation matrix 13099d007619Sjeremylt 13109d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 13119d007619Sjeremylt 1312b7c9bbdaSJeremy L Thompson @ref Advanced 13139d007619Sjeremylt **/ 13146c58de82SJeremy L Thompson int CeedBasisGetInterp(CeedBasis basis, const CeedScalar **interp) { 1315d1d35e2fSjeremylt if (!basis->interp && basis->tensor_basis) { 13169d007619Sjeremylt // Allocate 1317*2b730f8bSJeremy L Thompson CeedCall(CeedMalloc(basis->Q * basis->P, &basis->interp)); 13189d007619Sjeremylt 13199d007619Sjeremylt // Initialize 1320*2b730f8bSJeremy L Thompson for (CeedInt i = 0; i < basis->Q * basis->P; i++) basis->interp[i] = 1.0; 13219d007619Sjeremylt 13229d007619Sjeremylt // Calculate 1323*2b730f8bSJeremy L Thompson for (CeedInt d = 0; d < basis->dim; d++) { 1324*2b730f8bSJeremy L Thompson for (CeedInt qpt = 0; qpt < basis->Q; qpt++) { 13259d007619Sjeremylt for (CeedInt node = 0; node < basis->P; node++) { 1326d1d35e2fSjeremylt CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d; 1327d1d35e2fSjeremylt CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d; 1328d1d35e2fSjeremylt basis->interp[qpt * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p]; 13299d007619Sjeremylt } 13309d007619Sjeremylt } 1331*2b730f8bSJeremy L Thompson } 1332*2b730f8bSJeremy L Thompson } 13339d007619Sjeremylt *interp = basis->interp; 1334e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 13359d007619Sjeremylt } 13369d007619Sjeremylt 13379d007619Sjeremylt /** 13389d007619Sjeremylt @brief Get 1D interpolation matrix of a tensor product CeedBasis 13399d007619Sjeremylt 13409d007619Sjeremylt @param basis CeedBasis 1341d1d35e2fSjeremylt @param[out] interp_1d Variable to store interpolation matrix 13429d007619Sjeremylt 13439d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 13449d007619Sjeremylt 13459d007619Sjeremylt @ref Backend 13469d007619Sjeremylt **/ 1347d1d35e2fSjeremylt int CeedBasisGetInterp1D(CeedBasis basis, const CeedScalar **interp_1d) { 1348*2b730f8bSJeremy L Thompson if (!basis->tensor_basis) { 13499d007619Sjeremylt // LCOV_EXCL_START 1350*2b730f8bSJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_MINOR, "CeedBasis is not a tensor product basis."); 13519d007619Sjeremylt // LCOV_EXCL_STOP 1352*2b730f8bSJeremy L Thompson } 13539d007619Sjeremylt 1354d1d35e2fSjeremylt *interp_1d = basis->interp_1d; 1355e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 13569d007619Sjeremylt } 13579d007619Sjeremylt 13589d007619Sjeremylt /** 13599d007619Sjeremylt @brief Get gradient matrix of a CeedBasis 13609d007619Sjeremylt 13619d007619Sjeremylt @param basis CeedBasis 13629d007619Sjeremylt @param[out] grad Variable to store gradient matrix 13639d007619Sjeremylt 13649d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 13659d007619Sjeremylt 1366b7c9bbdaSJeremy L Thompson @ref Advanced 13679d007619Sjeremylt **/ 13686c58de82SJeremy L Thompson int CeedBasisGetGrad(CeedBasis basis, const CeedScalar **grad) { 1369d1d35e2fSjeremylt if (!basis->grad && basis->tensor_basis) { 13709d007619Sjeremylt // Allocate 1371*2b730f8bSJeremy L Thompson CeedCall(CeedMalloc(basis->dim * basis->Q * basis->P, &basis->grad)); 13729d007619Sjeremylt 13739d007619Sjeremylt // Initialize 1374*2b730f8bSJeremy L Thompson for (CeedInt i = 0; i < basis->dim * basis->Q * basis->P; i++) basis->grad[i] = 1.0; 13759d007619Sjeremylt 13769d007619Sjeremylt // Calculate 1377*2b730f8bSJeremy L Thompson for (CeedInt d = 0; d < basis->dim; d++) { 1378*2b730f8bSJeremy L Thompson for (CeedInt i = 0; i < basis->dim; i++) { 1379*2b730f8bSJeremy L Thompson for (CeedInt qpt = 0; qpt < basis->Q; qpt++) { 13809d007619Sjeremylt for (CeedInt node = 0; node < basis->P; node++) { 1381d1d35e2fSjeremylt CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d; 1382d1d35e2fSjeremylt CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d; 1383*2b730f8bSJeremy L Thompson if (i == d) basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->grad_1d[q * basis->P_1d + p]; 1384*2b730f8bSJeremy L Thompson else basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p]; 1385*2b730f8bSJeremy L Thompson } 1386*2b730f8bSJeremy L Thompson } 1387*2b730f8bSJeremy L Thompson } 13889d007619Sjeremylt } 13899d007619Sjeremylt } 13909d007619Sjeremylt *grad = basis->grad; 1391e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 13929d007619Sjeremylt } 13939d007619Sjeremylt 13949d007619Sjeremylt /** 13959d007619Sjeremylt @brief Get 1D gradient matrix of a tensor product CeedBasis 13969d007619Sjeremylt 13979d007619Sjeremylt @param basis CeedBasis 1398d1d35e2fSjeremylt @param[out] grad_1d Variable to store gradient matrix 13999d007619Sjeremylt 14009d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 14019d007619Sjeremylt 1402b7c9bbdaSJeremy L Thompson @ref Advanced 14039d007619Sjeremylt **/ 1404d1d35e2fSjeremylt int CeedBasisGetGrad1D(CeedBasis basis, const CeedScalar **grad_1d) { 1405*2b730f8bSJeremy L Thompson if (!basis->tensor_basis) { 14069d007619Sjeremylt // LCOV_EXCL_START 1407*2b730f8bSJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_MINOR, "CeedBasis is not a tensor product basis."); 14089d007619Sjeremylt // LCOV_EXCL_STOP 1409*2b730f8bSJeremy L Thompson } 14109d007619Sjeremylt 1411d1d35e2fSjeremylt *grad_1d = basis->grad_1d; 1412e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 14139d007619Sjeremylt } 14149d007619Sjeremylt 14159d007619Sjeremylt /** 141650c301a5SRezgar Shakeri @brief Get divergence matrix of a CeedBasis 141750c301a5SRezgar Shakeri 141850c301a5SRezgar Shakeri @param basis CeedBasis 141950c301a5SRezgar Shakeri @param[out] div Variable to store divergence matrix 142050c301a5SRezgar Shakeri 142150c301a5SRezgar Shakeri @return An error code: 0 - success, otherwise - failure 142250c301a5SRezgar Shakeri 142350c301a5SRezgar Shakeri @ref Advanced 142450c301a5SRezgar Shakeri **/ 142550c301a5SRezgar Shakeri int CeedBasisGetDiv(CeedBasis basis, const CeedScalar **div) { 1426*2b730f8bSJeremy L Thompson if (!basis->div) { 142750c301a5SRezgar Shakeri // LCOV_EXCL_START 1428*2b730f8bSJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_MINOR, "CeedBasis does not have divergence matrix."); 142950c301a5SRezgar Shakeri // LCOV_EXCL_STOP 1430*2b730f8bSJeremy L Thompson } 143150c301a5SRezgar Shakeri 143250c301a5SRezgar Shakeri *div = basis->div; 143350c301a5SRezgar Shakeri return CEED_ERROR_SUCCESS; 143450c301a5SRezgar Shakeri } 143550c301a5SRezgar Shakeri 143650c301a5SRezgar Shakeri /** 14377a982d89SJeremy L. Thompson @brief Destroy a CeedBasis 14387a982d89SJeremy L. Thompson 14397a982d89SJeremy L. Thompson @param basis CeedBasis to destroy 14407a982d89SJeremy L. Thompson 14417a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 14427a982d89SJeremy L. Thompson 14437a982d89SJeremy L. Thompson @ref User 14447a982d89SJeremy L. Thompson **/ 14457a982d89SJeremy L. Thompson int CeedBasisDestroy(CeedBasis *basis) { 1446d1d35e2fSjeremylt if (!*basis || --(*basis)->ref_count > 0) return CEED_ERROR_SUCCESS; 1447*2b730f8bSJeremy L Thompson if ((*basis)->Destroy) CeedCall((*basis)->Destroy(*basis)); 1448*2b730f8bSJeremy L Thompson if ((*basis)->contract) CeedCall(CeedTensorContractDestroy(&(*basis)->contract)); 1449*2b730f8bSJeremy L Thompson CeedCall(CeedFree(&(*basis)->interp)); 1450*2b730f8bSJeremy L Thompson CeedCall(CeedFree(&(*basis)->interp_1d)); 1451*2b730f8bSJeremy L Thompson CeedCall(CeedFree(&(*basis)->grad)); 1452*2b730f8bSJeremy L Thompson CeedCall(CeedFree(&(*basis)->div)); 1453*2b730f8bSJeremy L Thompson CeedCall(CeedFree(&(*basis)->grad_1d)); 1454*2b730f8bSJeremy L Thompson CeedCall(CeedFree(&(*basis)->q_ref_1d)); 1455*2b730f8bSJeremy L Thompson CeedCall(CeedFree(&(*basis)->q_weight_1d)); 1456*2b730f8bSJeremy L Thompson CeedCall(CeedDestroy(&(*basis)->ceed)); 1457*2b730f8bSJeremy L Thompson CeedCall(CeedFree(basis)); 1458e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 14597a982d89SJeremy L. Thompson } 14607a982d89SJeremy L. Thompson 14617a982d89SJeremy L. Thompson /** 1462b11c1e72Sjeremylt @brief Construct a Gauss-Legendre quadrature 1463b11c1e72Sjeremylt 1464b11c1e72Sjeremylt @param Q Number of quadrature points (integrates polynomials of 1465b11c1e72Sjeremylt degree 2*Q-1 exactly) 1466d1d35e2fSjeremylt @param[out] q_ref_1d Array of length Q to hold the abscissa on [-1, 1] 1467d1d35e2fSjeremylt @param[out] q_weight_1d Array of length Q to hold the weights 1468b11c1e72Sjeremylt 1469b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 1470dfdf5a53Sjeremylt 1471dfdf5a53Sjeremylt @ref Utility 1472b11c1e72Sjeremylt **/ 1473*2b730f8bSJeremy L Thompson int CeedGaussQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) { 1474d7b241e6Sjeremylt // Allocate 1475d7b241e6Sjeremylt CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0 * atan(1.0); 1476d1d35e2fSjeremylt // Build q_ref_1d, q_weight_1d 147792ae7e47SJeremy L Thompson for (CeedInt i = 0; i <= Q / 2; i++) { 1478d7b241e6Sjeremylt // Guess 1479d7b241e6Sjeremylt xi = cos(PI * (CeedScalar)(2 * i + 1) / ((CeedScalar)(2 * Q))); 1480d7b241e6Sjeremylt // Pn(xi) 1481d7b241e6Sjeremylt P0 = 1.0; 1482d7b241e6Sjeremylt P1 = xi; 1483d7b241e6Sjeremylt P2 = 0.0; 148492ae7e47SJeremy L Thompson for (CeedInt j = 2; j <= Q; j++) { 1485d7b241e6Sjeremylt P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j)); 1486d7b241e6Sjeremylt P0 = P1; 1487d7b241e6Sjeremylt P1 = P2; 1488d7b241e6Sjeremylt } 1489d7b241e6Sjeremylt // First Newton Step 1490d7b241e6Sjeremylt dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0); 1491d7b241e6Sjeremylt xi = xi - P2 / dP2; 1492d7b241e6Sjeremylt // Newton to convergence 149392ae7e47SJeremy L Thompson for (CeedInt k = 0; k < 100 && fabs(P2) > 10 * CEED_EPSILON; k++) { 1494d7b241e6Sjeremylt P0 = 1.0; 1495d7b241e6Sjeremylt P1 = xi; 149692ae7e47SJeremy L Thompson for (CeedInt j = 2; j <= Q; j++) { 1497d7b241e6Sjeremylt P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j)); 1498d7b241e6Sjeremylt P0 = P1; 1499d7b241e6Sjeremylt P1 = P2; 1500d7b241e6Sjeremylt } 1501d7b241e6Sjeremylt dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0); 1502d7b241e6Sjeremylt xi = xi - P2 / dP2; 1503d7b241e6Sjeremylt } 1504d7b241e6Sjeremylt // Save xi, wi 1505d7b241e6Sjeremylt wi = 2.0 / ((1.0 - xi * xi) * dP2 * dP2); 1506d1d35e2fSjeremylt q_weight_1d[i] = wi; 1507d1d35e2fSjeremylt q_weight_1d[Q - 1 - i] = wi; 1508d1d35e2fSjeremylt q_ref_1d[i] = -xi; 1509d1d35e2fSjeremylt q_ref_1d[Q - 1 - i] = xi; 1510d7b241e6Sjeremylt } 1511e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1512d7b241e6Sjeremylt } 1513d7b241e6Sjeremylt 1514b11c1e72Sjeremylt /** 1515b11c1e72Sjeremylt @brief Construct a Gauss-Legendre-Lobatto quadrature 1516b11c1e72Sjeremylt 1517b11c1e72Sjeremylt @param Q Number of quadrature points (integrates polynomials of 1518b11c1e72Sjeremylt degree 2*Q-3 exactly) 1519d1d35e2fSjeremylt @param[out] q_ref_1d Array of length Q to hold the abscissa on [-1, 1] 1520d1d35e2fSjeremylt @param[out] q_weight_1d Array of length Q to hold the weights 1521b11c1e72Sjeremylt 1522b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 1523dfdf5a53Sjeremylt 1524dfdf5a53Sjeremylt @ref Utility 1525b11c1e72Sjeremylt **/ 1526*2b730f8bSJeremy L Thompson int CeedLobattoQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) { 1527d7b241e6Sjeremylt // Allocate 1528d7b241e6Sjeremylt CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0 * atan(1.0); 1529d1d35e2fSjeremylt // Build q_ref_1d, q_weight_1d 1530d7b241e6Sjeremylt // Set endpoints 1531*2b730f8bSJeremy L Thompson if (Q < 2) { 1532b0d62198Sjeremylt // LCOV_EXCL_START 1533*2b730f8bSJeremy L Thompson return CeedError(NULL, CEED_ERROR_DIMENSION, "Cannot create Lobatto quadrature with Q=%" CeedInt_FMT " < 2 points", Q); 1534b0d62198Sjeremylt // LCOV_EXCL_STOP 1535*2b730f8bSJeremy L Thompson } 1536d7b241e6Sjeremylt wi = 2.0 / ((CeedScalar)(Q * (Q - 1))); 1537d1d35e2fSjeremylt if (q_weight_1d) { 1538d1d35e2fSjeremylt q_weight_1d[0] = wi; 1539d1d35e2fSjeremylt q_weight_1d[Q - 1] = wi; 1540d7b241e6Sjeremylt } 1541d1d35e2fSjeremylt q_ref_1d[0] = -1.0; 1542d1d35e2fSjeremylt q_ref_1d[Q - 1] = 1.0; 1543d7b241e6Sjeremylt // Interior 154492ae7e47SJeremy L Thompson for (CeedInt i = 1; i <= (Q - 1) / 2; i++) { 1545d7b241e6Sjeremylt // Guess 1546d7b241e6Sjeremylt xi = cos(PI * (CeedScalar)(i) / (CeedScalar)(Q - 1)); 1547d7b241e6Sjeremylt // Pn(xi) 1548d7b241e6Sjeremylt P0 = 1.0; 1549d7b241e6Sjeremylt P1 = xi; 1550d7b241e6Sjeremylt P2 = 0.0; 155192ae7e47SJeremy L Thompson for (CeedInt j = 2; j < Q; j++) { 1552d7b241e6Sjeremylt P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j)); 1553d7b241e6Sjeremylt P0 = P1; 1554d7b241e6Sjeremylt P1 = P2; 1555d7b241e6Sjeremylt } 1556d7b241e6Sjeremylt // First Newton step 1557d7b241e6Sjeremylt dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0); 1558d7b241e6Sjeremylt d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi); 1559d7b241e6Sjeremylt xi = xi - dP2 / d2P2; 1560d7b241e6Sjeremylt // Newton to convergence 156192ae7e47SJeremy L Thompson for (CeedInt k = 0; k < 100 && fabs(dP2) > 10 * CEED_EPSILON; k++) { 1562d7b241e6Sjeremylt P0 = 1.0; 1563d7b241e6Sjeremylt P1 = xi; 156492ae7e47SJeremy L Thompson for (CeedInt j = 2; j < Q; j++) { 1565d7b241e6Sjeremylt P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j)); 1566d7b241e6Sjeremylt P0 = P1; 1567d7b241e6Sjeremylt P1 = P2; 1568d7b241e6Sjeremylt } 1569d7b241e6Sjeremylt dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0); 1570d7b241e6Sjeremylt d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi); 1571d7b241e6Sjeremylt xi = xi - dP2 / d2P2; 1572d7b241e6Sjeremylt } 1573d7b241e6Sjeremylt // Save xi, wi 1574d7b241e6Sjeremylt wi = 2.0 / (((CeedScalar)(Q * (Q - 1))) * P2 * P2); 1575d1d35e2fSjeremylt if (q_weight_1d) { 1576d1d35e2fSjeremylt q_weight_1d[i] = wi; 1577d1d35e2fSjeremylt q_weight_1d[Q - 1 - i] = wi; 1578d7b241e6Sjeremylt } 1579d1d35e2fSjeremylt q_ref_1d[i] = -xi; 1580d1d35e2fSjeremylt q_ref_1d[Q - 1 - i] = xi; 1581d7b241e6Sjeremylt } 1582e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1583d7b241e6Sjeremylt } 1584d7b241e6Sjeremylt 1585dfdf5a53Sjeremylt /** 158695bb1877Svaleriabarra @brief Return QR Factorization of a matrix 1587b11c1e72Sjeremylt 158877645d7bSjeremylt @param ceed A Ceed context for error handling 158952bfb9bbSJeremy L Thompson @param[in,out] mat Row-major matrix to be factorized in place 159052bfb9bbSJeremy L Thompson @param[in,out] tau Vector of length m of scaling factors 1591b11c1e72Sjeremylt @param m Number of rows 1592b11c1e72Sjeremylt @param n Number of columns 1593b11c1e72Sjeremylt 1594b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 1595dfdf5a53Sjeremylt 1596dfdf5a53Sjeremylt @ref Utility 1597b11c1e72Sjeremylt **/ 1598*2b730f8bSJeremy L Thompson int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau, CeedInt m, CeedInt n) { 1599d7b241e6Sjeremylt CeedScalar v[m]; 1600d7b241e6Sjeremylt 1601*2b730f8bSJeremy L Thompson // Check matrix shape 1602*2b730f8bSJeremy L Thompson if (n > m) { 1603c042f62fSJeremy L Thompson // LCOV_EXCL_START 1604*2b730f8bSJeremy L Thompson return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute QR factorization with n > m"); 1605c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 1606*2b730f8bSJeremy L Thompson } 1607a7bd39daSjeremylt 160852bfb9bbSJeremy L Thompson for (CeedInt i = 0; i < n; i++) { 1609bde37e8cSJed Brown if (i >= m - 1) { // last row of matrix, no reflection needed 1610bde37e8cSJed Brown tau[i] = 0.; 1611bde37e8cSJed Brown break; 1612bde37e8cSJed Brown } 1613d7b241e6Sjeremylt // Calculate Householder vector, magnitude 1614d7b241e6Sjeremylt CeedScalar sigma = 0.0; 1615d7b241e6Sjeremylt v[i] = mat[i + n * i]; 161652bfb9bbSJeremy L Thompson for (CeedInt j = i + 1; j < m; j++) { 1617d7b241e6Sjeremylt v[j] = mat[i + n * j]; 1618d7b241e6Sjeremylt sigma += v[j] * v[j]; 1619d7b241e6Sjeremylt } 1620d7b241e6Sjeremylt CeedScalar norm = sqrt(v[i] * v[i] + sigma); // norm of v[i:m] 1621673160d7Sjeremylt CeedScalar R_ii = -copysign(norm, v[i]); 1622673160d7Sjeremylt v[i] -= R_ii; 1623d7b241e6Sjeremylt // norm of v[i:m] after modification above and scaling below 1624d7b241e6Sjeremylt // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 1625d7b241e6Sjeremylt // tau = 2 / (norm*norm) 1626d7b241e6Sjeremylt tau[i] = 2 * v[i] * v[i] / (v[i] * v[i] + sigma); 1627*2b730f8bSJeremy L Thompson for (CeedInt j = i + 1; j < m; j++) v[j] /= v[i]; 1628d7b241e6Sjeremylt 1629d7b241e6Sjeremylt // Apply Householder reflector to lower right panel 1630d7b241e6Sjeremylt CeedHouseholderReflect(&mat[i * n + i + 1], &v[i], tau[i], m - i, n - i - 1, n, 1); 1631d7b241e6Sjeremylt // Save v 1632673160d7Sjeremylt mat[i + n * i] = R_ii; 1633*2b730f8bSJeremy L Thompson for (CeedInt j = i + 1; j < m; j++) mat[i + n * j] = v[j]; 1634d7b241e6Sjeremylt } 1635e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1636d7b241e6Sjeremylt } 1637d7b241e6Sjeremylt 1638b11c1e72Sjeremylt /** 163952bfb9bbSJeremy L Thompson @brief Return symmetric Schur decomposition of the symmetric matrix mat via 164052bfb9bbSJeremy L Thompson symmetric QR factorization 164152bfb9bbSJeremy L Thompson 164277645d7bSjeremylt @param ceed A Ceed context for error handling 164352bfb9bbSJeremy L Thompson @param[in,out] mat Row-major matrix to be factorized in place 1644460bf743SValeria Barra @param[out] lambda Vector of length n of eigenvalues 164552bfb9bbSJeremy L Thompson @param n Number of rows/columns 164652bfb9bbSJeremy L Thompson 164752bfb9bbSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 164852bfb9bbSJeremy L Thompson 164952bfb9bbSJeremy L Thompson @ref Utility 165052bfb9bbSJeremy L Thompson **/ 1651*2b730f8bSJeremy L Thompson CeedPragmaOptimizeOff int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat, CeedScalar *lambda, CeedInt n) { 165252bfb9bbSJeremy L Thompson // Check bounds for clang-tidy 1653*2b730f8bSJeremy L Thompson if (n < 2) { 1654c042f62fSJeremy L Thompson // LCOV_EXCL_START 1655*2b730f8bSJeremy L Thompson return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute symmetric Schur decomposition of scalars"); 1656c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 1657*2b730f8bSJeremy L Thompson } 165852bfb9bbSJeremy L Thompson 1659673160d7Sjeremylt CeedScalar v[n - 1], tau[n - 1], mat_T[n * n]; 166052bfb9bbSJeremy L Thompson 1661673160d7Sjeremylt // Copy mat to mat_T and set mat to I 1662673160d7Sjeremylt memcpy(mat_T, mat, n * n * sizeof(mat[0])); 1663*2b730f8bSJeremy L Thompson for (CeedInt i = 0; i < n; i++) { 1664*2b730f8bSJeremy L Thompson for (CeedInt j = 0; j < n; j++) mat[j + n * i] = (i == j) ? 1 : 0; 1665*2b730f8bSJeremy L Thompson } 166652bfb9bbSJeremy L Thompson 166752bfb9bbSJeremy L Thompson // Reduce to tridiagonal 166852bfb9bbSJeremy L Thompson for (CeedInt i = 0; i < n - 1; i++) { 166952bfb9bbSJeremy L Thompson // Calculate Householder vector, magnitude 167052bfb9bbSJeremy L Thompson CeedScalar sigma = 0.0; 1671673160d7Sjeremylt v[i] = mat_T[i + n * (i + 1)]; 167252bfb9bbSJeremy L Thompson for (CeedInt j = i + 1; j < n - 1; j++) { 1673673160d7Sjeremylt v[j] = mat_T[i + n * (j + 1)]; 167452bfb9bbSJeremy L Thompson sigma += v[j] * v[j]; 167552bfb9bbSJeremy L Thompson } 167652bfb9bbSJeremy L Thompson CeedScalar norm = sqrt(v[i] * v[i] + sigma); // norm of v[i:n-1] 1677673160d7Sjeremylt CeedScalar R_ii = -copysign(norm, v[i]); 1678673160d7Sjeremylt v[i] -= R_ii; 167952bfb9bbSJeremy L Thompson // norm of v[i:m] after modification above and scaling below 168052bfb9bbSJeremy L Thompson // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 168152bfb9bbSJeremy L Thompson // tau = 2 / (norm*norm) 168280a9ef05SNatalie Beams tau[i] = i == n - 2 ? 2 : 2 * v[i] * v[i] / (v[i] * v[i] + sigma); 1683*2b730f8bSJeremy L Thompson for (CeedInt j = i + 1; j < n - 1; j++) v[j] /= v[i]; 168452bfb9bbSJeremy L Thompson 168552bfb9bbSJeremy L Thompson // Update sub and super diagonal 168652bfb9bbSJeremy L Thompson for (CeedInt j = i + 2; j < n; j++) { 1687*2b730f8bSJeremy L Thompson mat_T[i + n * j] = 0; 1688*2b730f8bSJeremy L Thompson mat_T[j + n * i] = 0; 168952bfb9bbSJeremy L Thompson } 169052bfb9bbSJeremy L Thompson // Apply symmetric Householder reflector to lower right panel 1691*2b730f8bSJeremy L Thompson CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1); 1692*2b730f8bSJeremy L Thompson CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), 1, n); 1693673160d7Sjeremylt 169452bfb9bbSJeremy L Thompson // Save v 1695673160d7Sjeremylt mat_T[i + n * (i + 1)] = R_ii; 1696673160d7Sjeremylt mat_T[(i + 1) + n * i] = R_ii; 169752bfb9bbSJeremy L Thompson for (CeedInt j = i + 1; j < n - 1; j++) { 1698673160d7Sjeremylt mat_T[i + n * (j + 1)] = v[j]; 169952bfb9bbSJeremy L Thompson } 170052bfb9bbSJeremy L Thompson } 170152bfb9bbSJeremy L Thompson // Backwards accumulation of Q 170252bfb9bbSJeremy L Thompson for (CeedInt i = n - 2; i >= 0; i--) { 170385cf89eaSjeremylt if (tau[i] > 0.0) { 170452bfb9bbSJeremy L Thompson v[i] = 1; 170552bfb9bbSJeremy L Thompson for (CeedInt j = i + 1; j < n - 1; j++) { 1706673160d7Sjeremylt v[j] = mat_T[i + n * (j + 1)]; 1707673160d7Sjeremylt mat_T[i + n * (j + 1)] = 0; 170852bfb9bbSJeremy L Thompson } 1709*2b730f8bSJeremy L Thompson CeedHouseholderReflect(&mat[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1); 171052bfb9bbSJeremy L Thompson } 171185cf89eaSjeremylt } 171252bfb9bbSJeremy L Thompson 171352bfb9bbSJeremy L Thompson // Reduce sub and super diagonal 1714673160d7Sjeremylt CeedInt p = 0, q = 0, itr = 0, max_itr = n * n * n * n; 1715673160d7Sjeremylt CeedScalar tol = CEED_EPSILON; 171652bfb9bbSJeremy L Thompson 1717673160d7Sjeremylt while (itr < max_itr) { 171852bfb9bbSJeremy L Thompson // Update p, q, size of reduced portions of diagonal 1719*2b730f8bSJeremy L Thompson p = 0; 1720*2b730f8bSJeremy L Thompson q = 0; 172152bfb9bbSJeremy L Thompson for (CeedInt i = n - 2; i >= 0; i--) { 1722*2b730f8bSJeremy L Thompson if (fabs(mat_T[i + n * (i + 1)]) < tol) q += 1; 1723*2b730f8bSJeremy L Thompson else break; 172452bfb9bbSJeremy L Thompson } 1725673160d7Sjeremylt for (CeedInt i = 0; i < n - q - 1; i++) { 1726*2b730f8bSJeremy L Thompson if (fabs(mat_T[i + n * (i + 1)]) < tol) p += 1; 1727*2b730f8bSJeremy L Thompson else break; 172852bfb9bbSJeremy L Thompson } 172952bfb9bbSJeremy L Thompson if (q == n - 1) break; // Finished reducing 173052bfb9bbSJeremy L Thompson 173152bfb9bbSJeremy L Thompson // Reduce tridiagonal portion 1732*2b730f8bSJeremy L Thompson CeedScalar t_nn = mat_T[(n - 1 - q) + n * (n - 1 - q)], t_nnm1 = mat_T[(n - 2 - q) + n * (n - 1 - q)]; 1733673160d7Sjeremylt CeedScalar d = (mat_T[(n - 2 - q) + n * (n - 2 - q)] - t_nn) / 2; 1734*2b730f8bSJeremy L Thompson CeedScalar mu = t_nn - t_nnm1 * t_nnm1 / (d + copysign(sqrt(d * d + t_nnm1 * t_nnm1), d)); 1735673160d7Sjeremylt CeedScalar x = mat_T[p + n * p] - mu; 1736673160d7Sjeremylt CeedScalar z = mat_T[p + n * (p + 1)]; 1737673160d7Sjeremylt for (CeedInt k = p; k < n - q - 1; k++) { 173852bfb9bbSJeremy L Thompson // Compute Givens rotation 173952bfb9bbSJeremy L Thompson CeedScalar c = 1, s = 0; 174052bfb9bbSJeremy L Thompson if (fabs(z) > tol) { 174152bfb9bbSJeremy L Thompson if (fabs(z) > fabs(x)) { 174252bfb9bbSJeremy L Thompson CeedScalar tau = -x / z; 174352bfb9bbSJeremy L Thompson s = 1 / sqrt(1 + tau * tau), c = s * tau; 174452bfb9bbSJeremy L Thompson } else { 174552bfb9bbSJeremy L Thompson CeedScalar tau = -z / x; 174652bfb9bbSJeremy L Thompson c = 1 / sqrt(1 + tau * tau), s = c * tau; 174752bfb9bbSJeremy L Thompson } 174852bfb9bbSJeremy L Thompson } 174952bfb9bbSJeremy L Thompson 175052bfb9bbSJeremy L Thompson // Apply Givens rotation to T 1751673160d7Sjeremylt CeedGivensRotation(mat_T, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n); 1752673160d7Sjeremylt CeedGivensRotation(mat_T, c, s, CEED_TRANSPOSE, k, k + 1, n, n); 175352bfb9bbSJeremy L Thompson 175452bfb9bbSJeremy L Thompson // Apply Givens rotation to Q 175552bfb9bbSJeremy L Thompson CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n); 175652bfb9bbSJeremy L Thompson 175752bfb9bbSJeremy L Thompson // Update x, z 175852bfb9bbSJeremy L Thompson if (k < n - q - 2) { 1759673160d7Sjeremylt x = mat_T[k + n * (k + 1)]; 1760673160d7Sjeremylt z = mat_T[k + n * (k + 2)]; 176152bfb9bbSJeremy L Thompson } 176252bfb9bbSJeremy L Thompson } 176352bfb9bbSJeremy L Thompson itr++; 176452bfb9bbSJeremy L Thompson } 1765673160d7Sjeremylt 176652bfb9bbSJeremy L Thompson // Save eigenvalues 1767*2b730f8bSJeremy L Thompson for (CeedInt i = 0; i < n; i++) lambda[i] = mat_T[i + n * i]; 176852bfb9bbSJeremy L Thompson 176952bfb9bbSJeremy L Thompson // Check convergence 1770*2b730f8bSJeremy L Thompson if (itr == max_itr && q < n - 1) { 1771c042f62fSJeremy L Thompson // LCOV_EXCL_START 1772*2b730f8bSJeremy L Thompson return CeedError(ceed, CEED_ERROR_MINOR, "Symmetric QR failed to converge"); 1773c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 1774*2b730f8bSJeremy L Thompson } 1775e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 177652bfb9bbSJeremy L Thompson } 177703d18186Sjeremylt CeedPragmaOptimizeOn 177852bfb9bbSJeremy L Thompson 177952bfb9bbSJeremy L Thompson /** 178052bfb9bbSJeremy L Thompson @brief Return Simultaneous Diagonalization of two matrices. This solves the 178152bfb9bbSJeremy L Thompson generalized eigenvalue problem A x = lambda B x, where A and B 178252bfb9bbSJeremy L Thompson are symmetric and B is positive definite. We generate the matrix X 178352bfb9bbSJeremy L Thompson and vector Lambda such that X^T A X = Lambda and X^T B X = I. This 178452bfb9bbSJeremy L Thompson is equivalent to the LAPACK routine 'sygv' with TYPE = 1. 178552bfb9bbSJeremy L Thompson 178677645d7bSjeremylt @param ceed A Ceed context for error handling 1787d1d35e2fSjeremylt @param[in] mat_A Row-major matrix to be factorized with eigenvalues 1788d1d35e2fSjeremylt @param[in] mat_B Row-major matrix to be factorized to identity 1789d3331725Sjeremylt @param[out] mat_X Row-major orthogonal matrix 1790460bf743SValeria Barra @param[out] lambda Vector of length n of generalized eigenvalues 179152bfb9bbSJeremy L Thompson @param n Number of rows/columns 179252bfb9bbSJeremy L Thompson 179352bfb9bbSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 179452bfb9bbSJeremy L Thompson 179552bfb9bbSJeremy L Thompson @ref Utility 179652bfb9bbSJeremy L Thompson **/ 1797*2b730f8bSJeremy L Thompson CeedPragmaOptimizeOff int 1798*2b730f8bSJeremy L Thompson CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *mat_A, CeedScalar *mat_B, CeedScalar *mat_X, CeedScalar *lambda, CeedInt n) { 1799d3331725Sjeremylt CeedScalar *mat_C, *mat_G, *vec_D; 1800*2b730f8bSJeremy L Thompson CeedCall(CeedCalloc(n * n, &mat_C)); 1801*2b730f8bSJeremy L Thompson CeedCall(CeedCalloc(n * n, &mat_G)); 1802*2b730f8bSJeremy L Thompson CeedCall(CeedCalloc(n, &vec_D)); 180352bfb9bbSJeremy L Thompson 180452bfb9bbSJeremy L Thompson // Compute B = G D G^T 180578464608Sjeremylt memcpy(mat_G, mat_B, n * n * sizeof(mat_B[0])); 1806*2b730f8bSJeremy L Thompson CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_G, vec_D, n)); 180752bfb9bbSJeremy L Thompson 180885cf89eaSjeremylt // Sort eigenvalues 1809*2b730f8bSJeremy L Thompson for (CeedInt i = n - 1; i >= 0; i--) { 181085cf89eaSjeremylt for (CeedInt j = 0; j < i; j++) { 181185cf89eaSjeremylt if (fabs(vec_D[j]) > fabs(vec_D[j + 1])) { 181285cf89eaSjeremylt CeedScalar temp; 1813*2b730f8bSJeremy L Thompson temp = vec_D[j]; 1814*2b730f8bSJeremy L Thompson vec_D[j] = vec_D[j + 1]; 1815*2b730f8bSJeremy L Thompson vec_D[j + 1] = temp; 181685cf89eaSjeremylt for (CeedInt k = 0; k < n; k++) { 1817*2b730f8bSJeremy L Thompson temp = mat_G[k * n + j]; 1818*2b730f8bSJeremy L Thompson mat_G[k * n + j] = mat_G[k * n + j + 1]; 1819*2b730f8bSJeremy L Thompson mat_G[k * n + j + 1] = temp; 1820*2b730f8bSJeremy L Thompson } 182185cf89eaSjeremylt } 182285cf89eaSjeremylt } 182385cf89eaSjeremylt } 182485cf89eaSjeremylt 1825fb551037Sjeremylt // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T 1826fb551037Sjeremylt // = D^-1/2 G^T A G D^-1/2 1827d3331725Sjeremylt // -- D = D^-1/2 1828*2b730f8bSJeremy L Thompson for (CeedInt i = 0; i < n; i++) vec_D[i] = 1. / sqrt(vec_D[i]); 1829d3331725Sjeremylt // -- G = G D^-1/2 1830d3331725Sjeremylt // -- C = D^-1/2 G^T 1831*2b730f8bSJeremy L Thompson for (CeedInt i = 0; i < n; i++) { 1832d3331725Sjeremylt for (CeedInt j = 0; j < n; j++) { 1833673160d7Sjeremylt mat_G[i * n + j] *= vec_D[j]; 1834673160d7Sjeremylt mat_C[j * n + i] = mat_G[i * n + j]; 1835d3331725Sjeremylt } 1836*2b730f8bSJeremy L Thompson } 1837673160d7Sjeremylt // -- X = (D^-1/2 G^T) A 1838*2b730f8bSJeremy L Thompson CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_C, (const CeedScalar *)mat_A, mat_X, n, n, n)); 1839673160d7Sjeremylt // -- C = (D^-1/2 G^T A) (G D^-1/2) 1840*2b730f8bSJeremy L Thompson CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_X, (const CeedScalar *)mat_G, mat_C, n, n, n)); 184152bfb9bbSJeremy L Thompson 184252bfb9bbSJeremy L Thompson // Compute Q^T C Q = lambda 1843*2b730f8bSJeremy L Thompson CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_C, lambda, n)); 184452bfb9bbSJeremy L Thompson 184585cf89eaSjeremylt // Sort eigenvalues 1846*2b730f8bSJeremy L Thompson for (CeedInt i = n - 1; i >= 0; i--) { 184785cf89eaSjeremylt for (CeedInt j = 0; j < i; j++) { 184885cf89eaSjeremylt if (fabs(lambda[j]) > fabs(lambda[j + 1])) { 184985cf89eaSjeremylt CeedScalar temp; 1850*2b730f8bSJeremy L Thompson temp = lambda[j]; 1851*2b730f8bSJeremy L Thompson lambda[j] = lambda[j + 1]; 1852*2b730f8bSJeremy L Thompson lambda[j + 1] = temp; 185385cf89eaSjeremylt for (CeedInt k = 0; k < n; k++) { 1854*2b730f8bSJeremy L Thompson temp = mat_C[k * n + j]; 1855*2b730f8bSJeremy L Thompson mat_C[k * n + j] = mat_C[k * n + j + 1]; 1856*2b730f8bSJeremy L Thompson mat_C[k * n + j + 1] = temp; 1857*2b730f8bSJeremy L Thompson } 185885cf89eaSjeremylt } 185985cf89eaSjeremylt } 186085cf89eaSjeremylt } 186185cf89eaSjeremylt 1862d3331725Sjeremylt // Set X = (G D^1/2)^-T Q 1863fb551037Sjeremylt // = G D^-1/2 Q 1864*2b730f8bSJeremy L Thompson CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_G, (const CeedScalar *)mat_C, mat_X, n, n, n)); 186578464608Sjeremylt 186678464608Sjeremylt // Cleanup 1867*2b730f8bSJeremy L Thompson CeedCall(CeedFree(&mat_C)); 1868*2b730f8bSJeremy L Thompson CeedCall(CeedFree(&mat_G)); 1869*2b730f8bSJeremy L Thompson CeedCall(CeedFree(&vec_D)); 1870e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 187152bfb9bbSJeremy L Thompson } 187203d18186Sjeremylt CeedPragmaOptimizeOn 187352bfb9bbSJeremy L Thompson 1874d7b241e6Sjeremylt /// @} 1875