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 8ec3da8bcSJed Brown #include <ceed/ceed.h> 9ec3da8bcSJed Brown #include <ceed/backend.h> 103d576824SJeremy L Thompson #include <ceed-impl.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 **/ 557a982d89SJeremy L. Thompson static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v, 567a982d89SJeremy L. Thompson CeedScalar b, CeedInt m, CeedInt n, 577a982d89SJeremy L. Thompson CeedInt row, CeedInt col) { 587a982d89SJeremy L. Thompson for (CeedInt j=0; j<n; j++) { 597a982d89SJeremy L. Thompson CeedScalar w = A[0*row + j*col]; 607a982d89SJeremy L. Thompson for (CeedInt i=1; i<m; i++) 617a982d89SJeremy L. Thompson w += v[i] * A[i*row + j*col]; 627a982d89SJeremy L. Thompson A[0*row + j*col] -= b * w; 637a982d89SJeremy L. Thompson for (CeedInt i=1; i<m; i++) 647a982d89SJeremy L. Thompson A[i*row + j*col] -= b * w * v[i]; 657a982d89SJeremy L. Thompson } 66e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 677a982d89SJeremy L. Thompson } 687a982d89SJeremy L. Thompson 697a982d89SJeremy L. Thompson /** 707a982d89SJeremy L. Thompson @brief Apply Householder Q matrix 717a982d89SJeremy L. Thompson 727a982d89SJeremy L. Thompson Compute A = Q A where Q is mxm and A is mxn. 737a982d89SJeremy L. Thompson 747a982d89SJeremy L. Thompson @param[in,out] A Matrix to apply Householder Q to, in place 757a982d89SJeremy L. Thompson @param Q Householder Q matrix 767a982d89SJeremy L. Thompson @param tau Householder scaling factors 77d1d35e2fSjeremylt @param t_mode Transpose mode for application 787a982d89SJeremy L. Thompson @param m Number of rows in A 797a982d89SJeremy L. Thompson @param n Number of columns in A 807a982d89SJeremy L. Thompson @param k Number of elementary reflectors in Q, k<m 817a982d89SJeremy L. Thompson @param row Row stride in A 827a982d89SJeremy L. Thompson @param col Col stride in A 837a982d89SJeremy L. Thompson 847a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 857a982d89SJeremy L. Thompson 867a982d89SJeremy L. Thompson @ref Developer 877a982d89SJeremy L. Thompson **/ 88d99fa3c5SJeremy L Thompson int CeedHouseholderApplyQ(CeedScalar *A, const CeedScalar *Q, 89d1d35e2fSjeremylt const CeedScalar *tau, CeedTransposeMode t_mode, 907a982d89SJeremy L. Thompson CeedInt m, CeedInt n, CeedInt k, 917a982d89SJeremy L. Thompson CeedInt row, CeedInt col) { 92e15f9bd0SJeremy L Thompson int ierr; 9378464608Sjeremylt CeedScalar *v; 9478464608Sjeremylt ierr = CeedMalloc(m, &v); CeedChk(ierr); 957a982d89SJeremy L. Thompson for (CeedInt ii=0; ii<k; ii++) { 96d1d35e2fSjeremylt CeedInt i = t_mode == CEED_TRANSPOSE ? ii : k-1-ii; 977a982d89SJeremy L. Thompson for (CeedInt j=i+1; j<m; j++) 987a982d89SJeremy L. Thompson v[j] = Q[j*k+i]; 99d1d35e2fSjeremylt // Apply Householder reflector (I - tau v v^T) collo_grad_1d^T 100e15f9bd0SJeremy L Thompson ierr = CeedHouseholderReflect(&A[i*row], &v[i], tau[i], m-i, n, row, col); 101e15f9bd0SJeremy L Thompson CeedChk(ierr); 1027a982d89SJeremy L. Thompson } 10378464608Sjeremylt ierr = CeedFree(&v); CeedChk(ierr); 104e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1057a982d89SJeremy L. Thompson } 1067a982d89SJeremy L. Thompson 1077a982d89SJeremy L. Thompson /** 1087a982d89SJeremy L. Thompson @brief Compute Givens rotation 1097a982d89SJeremy L. Thompson 1107a982d89SJeremy L. Thompson Computes A = G A (or G^T A in transpose mode) 1117a982d89SJeremy L. Thompson where A is an mxn matrix indexed as A[i*n + j*m] 1127a982d89SJeremy L. Thompson 1137a982d89SJeremy L. Thompson @param[in,out] A Row major matrix to apply Givens rotation to, in place 1147a982d89SJeremy L. Thompson @param c Cosine factor 1157a982d89SJeremy L. Thompson @param s Sine factor 116d1d35e2fSjeremylt @param t_mode @ref CEED_NOTRANSPOSE to rotate the basis counter-clockwise, 1174c4400c7SValeria Barra which has the effect of rotating columns of A clockwise; 1184cc79fe7SJed Brown @ref CEED_TRANSPOSE for the opposite rotation 1197a982d89SJeremy L. Thompson @param i First row/column to apply rotation 1207a982d89SJeremy L. Thompson @param k Second row/column to apply rotation 1217a982d89SJeremy L. Thompson @param m Number of rows in A 1227a982d89SJeremy L. Thompson @param n Number of columns in A 1237a982d89SJeremy L. Thompson 1247a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 1257a982d89SJeremy L. Thompson 1267a982d89SJeremy L. Thompson @ref Developer 1277a982d89SJeremy L. Thompson **/ 1287a982d89SJeremy L. Thompson static int CeedGivensRotation(CeedScalar *A, CeedScalar c, CeedScalar s, 129d1d35e2fSjeremylt CeedTransposeMode t_mode, CeedInt i, CeedInt k, 1307a982d89SJeremy L. Thompson CeedInt m, CeedInt n) { 131d1d35e2fSjeremylt CeedInt stride_j = 1, stride_ik = m, num_its = n; 132d1d35e2fSjeremylt if (t_mode == CEED_NOTRANSPOSE) { 133d1d35e2fSjeremylt stride_j = n; stride_ik = 1; num_its = m; 1347a982d89SJeremy L. Thompson } 1357a982d89SJeremy L. Thompson 1367a982d89SJeremy L. Thompson // Apply rotation 137d1d35e2fSjeremylt for (CeedInt j=0; j<num_its; j++) { 138d1d35e2fSjeremylt CeedScalar tau1 = A[i*stride_ik+j*stride_j], tau2 = A[k*stride_ik+j*stride_j]; 139d1d35e2fSjeremylt A[i*stride_ik+j*stride_j] = c*tau1 - s*tau2; 140d1d35e2fSjeremylt A[k*stride_ik+j*stride_j] = s*tau1 + c*tau2; 1417a982d89SJeremy L. Thompson } 142e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1437a982d89SJeremy L. Thompson } 1447a982d89SJeremy L. Thompson 1457a982d89SJeremy L. Thompson /** 1467a982d89SJeremy L. Thompson @brief View an array stored in a CeedBasis 1477a982d89SJeremy L. Thompson 1480a0da059Sjeremylt @param[in] name Name of array 149d1d35e2fSjeremylt @param[in] fp_fmt Printing format 1500a0da059Sjeremylt @param[in] m Number of rows in array 1510a0da059Sjeremylt @param[in] n Number of columns in array 1520a0da059Sjeremylt @param[in] a Array to be viewed 1530a0da059Sjeremylt @param[in] stream Stream to view to, e.g., stdout 1547a982d89SJeremy L. Thompson 1557a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 1567a982d89SJeremy L. Thompson 1577a982d89SJeremy L. Thompson @ref Developer 1587a982d89SJeremy L. Thompson **/ 159d1d35e2fSjeremylt static int CeedScalarView(const char *name, const char *fp_fmt, CeedInt m, 1607a982d89SJeremy L. Thompson CeedInt n, const CeedScalar *a, FILE *stream) { 16192ae7e47SJeremy L Thompson for (CeedInt i=0; i<m; i++) { 1627a982d89SJeremy L. Thompson if (m > 1) 163990fdeb6SJeremy L Thompson fprintf(stream, "%12s[%" CeedInt_FMT "]:", name, i); 1647a982d89SJeremy L. Thompson else 1657a982d89SJeremy L. Thompson fprintf(stream, "%12s:", name); 16692ae7e47SJeremy L Thompson for (CeedInt j=0; j<n; j++) 167d1d35e2fSjeremylt fprintf(stream, fp_fmt, fabs(a[i*n+j]) > 1E-14 ? a[i*n+j] : 0); 1687a982d89SJeremy L. Thompson fputs("\n", stream); 1697a982d89SJeremy L. Thompson } 170e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1717a982d89SJeremy L. Thompson } 1727a982d89SJeremy L. Thompson 173a76a04e7SJeremy L Thompson /** 174*14556e63SJeremy L Thompson @brief Create the interpolation and gradient matrices for projection from 175*14556e63SJeremy L Thompson the nodes of `basis_from` to the nodes of `basis_to`. 176*14556e63SJeremy L Thompson The interpolation is given by `interp_project = interp_to^+ * interp_from`, 177*14556e63SJeremy L Thompson where the pesudoinverse `interp_to^+` is given by QR factorization. 178*14556e63SJeremy L Thompson The gradient is given by `grad_project = interp_to^+ * grad_from`. 179a76a04e7SJeremy L Thompson Note: `basis_from` and `basis_to` must have compatible quadrature 180a76a04e7SJeremy L Thompson spaces. 181a76a04e7SJeremy L Thompson 182a76a04e7SJeremy L Thompson @param[in] basis_from CeedBasis to project from 183a76a04e7SJeremy L Thompson @param[in] basis_to CeedBasis to project to 184a76a04e7SJeremy L Thompson @param[out] interp_project Address of the variable where the newly created 185*14556e63SJeremy L Thompson interpolation matrix will be stored. 186*14556e63SJeremy L Thompson @param[out] grad_project Address of the variable where the newly created 187*14556e63SJeremy L Thompson gradient matrix will be stored. 188a76a04e7SJeremy L Thompson 189a76a04e7SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 190a76a04e7SJeremy L Thompson 191a76a04e7SJeremy L Thompson @ref Developer 192a76a04e7SJeremy L Thompson **/ 193*14556e63SJeremy L Thompson static int CeedBasisCreateProjectionMatrices(CeedBasis basis_from, 194*14556e63SJeremy L Thompson CeedBasis basis_to, CeedScalar **interp_project, CeedScalar **grad_project) { 195a76a04e7SJeremy L Thompson int ierr; 196a76a04e7SJeremy L Thompson Ceed ceed; 197a76a04e7SJeremy L Thompson ierr = CeedBasisGetCeed(basis_to, &ceed); CeedChk(ierr); 198a76a04e7SJeremy L Thompson 199a76a04e7SJeremy L Thompson // Check for compatible quadrature spaces 200a76a04e7SJeremy L Thompson CeedInt Q_to, Q_from; 201a76a04e7SJeremy L Thompson ierr = CeedBasisGetNumQuadraturePoints(basis_to, &Q_to); CeedChk(ierr); 202a76a04e7SJeremy L Thompson ierr = CeedBasisGetNumQuadraturePoints(basis_from, &Q_from); CeedChk(ierr); 203a76a04e7SJeremy L Thompson if (Q_to != Q_from) 204a76a04e7SJeremy L Thompson // LCOV_EXCL_START 205a76a04e7SJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, 206a76a04e7SJeremy L Thompson "Bases must have compatible quadrature spaces"); 207a76a04e7SJeremy L Thompson // LCOV_EXCL_STOP 208a76a04e7SJeremy L Thompson 209*14556e63SJeremy L Thompson // Check for matching tensor or non-tensor 210a76a04e7SJeremy L Thompson CeedInt P_to, P_from, Q = Q_to; 211a76a04e7SJeremy L Thompson bool is_tensor_to, is_tensor_from; 212a76a04e7SJeremy L Thompson ierr = CeedBasisIsTensor(basis_to, &is_tensor_to); CeedChk(ierr); 213a76a04e7SJeremy L Thompson ierr = CeedBasisIsTensor(basis_from, &is_tensor_from); CeedChk(ierr); 214a76a04e7SJeremy L Thompson if (is_tensor_to && is_tensor_from) { 215a76a04e7SJeremy L Thompson ierr = CeedBasisGetNumNodes1D(basis_to, &P_to); CeedChk(ierr); 216a76a04e7SJeremy L Thompson ierr = CeedBasisGetNumNodes1D(basis_from, &P_from); CeedChk(ierr); 217a76a04e7SJeremy L Thompson ierr = CeedBasisGetNumQuadraturePoints1D(basis_from, &Q); CeedChk(ierr); 218a76a04e7SJeremy L Thompson } else if (!is_tensor_to && !is_tensor_from) { 219a76a04e7SJeremy L Thompson ierr = CeedBasisGetNumNodes(basis_to, &P_to); CeedChk(ierr); 220a76a04e7SJeremy L Thompson ierr = CeedBasisGetNumNodes(basis_from, &P_from); CeedChk(ierr); 221a76a04e7SJeremy L Thompson } else { 222a76a04e7SJeremy L Thompson // LCOV_EXCL_START 223a76a04e7SJeremy L Thompson return CeedError(ceed, CEED_ERROR_MINOR, 224a76a04e7SJeremy L Thompson "Bases must both be tensor or non-tensor"); 225a76a04e7SJeremy L Thompson // LCOV_EXCL_STOP 226a76a04e7SJeremy L Thompson } 227a76a04e7SJeremy L Thompson 228*14556e63SJeremy L Thompson // Get source matrices 229*14556e63SJeremy L Thompson CeedInt dim; 230*14556e63SJeremy L Thompson CeedScalar *interp_to, *interp_from, *tau; 231*14556e63SJeremy L Thompson ierr = CeedBasisGetDimension(basis_to, &dim); CeedChk(ierr); 232a76a04e7SJeremy L Thompson ierr = CeedMalloc(Q * P_from, &interp_from); CeedChk(ierr); 233a76a04e7SJeremy L Thompson ierr = CeedMalloc(Q * P_to, &interp_to); CeedChk(ierr); 234a76a04e7SJeremy L Thompson ierr = CeedCalloc(P_to * P_from, interp_project); CeedChk(ierr); 235*14556e63SJeremy L Thompson ierr = CeedCalloc(P_to * P_from * (is_tensor_to ? 1 : dim), grad_project); 236*14556e63SJeremy L Thompson CeedChk(ierr); 237a76a04e7SJeremy L Thompson ierr = CeedMalloc(Q, &tau); CeedChk(ierr); 238*14556e63SJeremy L Thompson const CeedScalar *interp_to_source = NULL, *interp_from_source = NULL, 239*14556e63SJeremy L Thompson *grad_from_source; 240a76a04e7SJeremy L Thompson if (is_tensor_to) { 241a76a04e7SJeremy L Thompson ierr = CeedBasisGetInterp1D(basis_to, &interp_to_source); CeedChk(ierr); 242a76a04e7SJeremy L Thompson ierr = CeedBasisGetInterp1D(basis_from, &interp_from_source); CeedChk(ierr); 243*14556e63SJeremy L Thompson ierr = CeedBasisGetGrad1D(basis_from, &grad_from_source); CeedChk(ierr); 244a76a04e7SJeremy L Thompson } else { 245a76a04e7SJeremy L Thompson ierr = CeedBasisGetInterp(basis_to, &interp_to_source); CeedChk(ierr); 246a76a04e7SJeremy L Thompson ierr = CeedBasisGetInterp(basis_from, &interp_from_source); CeedChk(ierr); 247*14556e63SJeremy L Thompson ierr = CeedBasisGetGrad(basis_from, &grad_from_source); CeedChk(ierr); 248a76a04e7SJeremy L Thompson } 249a76a04e7SJeremy L Thompson 250*14556e63SJeremy L Thompson // Build matrices 251*14556e63SJeremy L Thompson CeedInt num_matrices = 1 + (is_tensor_to ? 1 : dim); 252*14556e63SJeremy L Thompson CeedScalar *input_from[num_matrices], *output_project[num_matrices]; 253*14556e63SJeremy L Thompson input_from[0] = (CeedScalar *)interp_from_source; 254*14556e63SJeremy L Thompson output_project[0] = *interp_project; 255*14556e63SJeremy L Thompson for (CeedInt m = 1; m < num_matrices; m++) { 256*14556e63SJeremy L Thompson input_from[m] = (CeedScalar *)&grad_from_source[(m - 1) * Q * P_from]; 257*14556e63SJeremy L Thompson output_project[m] = &(*grad_project[(m - 1) * P_to * P_from]); 258*14556e63SJeremy L Thompson } 259*14556e63SJeremy L Thompson for (CeedInt m = 0; m < num_matrices; m++) { 260a76a04e7SJeremy L Thompson // -- QR Factorization, interp_to = Q R 261*14556e63SJeremy L Thompson memcpy(interp_to, interp_to_source, Q * P_to * sizeof(interp_to_source[0])); 262a76a04e7SJeremy L Thompson ierr = CeedQRFactorization(ceed, interp_to, tau, Q, P_to); CeedChk(ierr); 263a76a04e7SJeremy L Thompson 264a76a04e7SJeremy L Thompson // -- Apply Qtranspose, interp_to = Qtranspose interp_from 265*14556e63SJeremy L Thompson memcpy(interp_from, input_from[m], Q * P_from * sizeof(input_from[m][0])); 266a76a04e7SJeremy L Thompson ierr = CeedHouseholderApplyQ(interp_from, interp_to, tau, CEED_TRANSPOSE, 267a76a04e7SJeremy L Thompson Q, P_from, P_to, P_from, 1); CeedChk(ierr); 268a76a04e7SJeremy L Thompson 269a76a04e7SJeremy L Thompson // -- Apply Rinv, interp_project = Rinv interp_c 270a76a04e7SJeremy L Thompson for (CeedInt j = 0; j < P_from; j++) { // Column j 271*14556e63SJeremy L Thompson output_project[m][j + P_from * (P_to - 1)] = interp_from[j + P_from * 272a76a04e7SJeremy L Thompson (P_to - 1)] / interp_to[P_to * P_to - 1]; 273a76a04e7SJeremy L Thompson for (CeedInt i = P_to - 2; i >= 0; i--) { // Row i 274*14556e63SJeremy L Thompson output_project[m][j + P_from * i] = interp_from[j + P_from * i]; 275a76a04e7SJeremy L Thompson for (CeedInt k = i+1; k < P_to; k++) { 276*14556e63SJeremy L Thompson output_project[m][j + P_from * i] -= interp_to[k + P_to * i]* 277*14556e63SJeremy L Thompson output_project[m][j + P_from * k]; 278a76a04e7SJeremy L Thompson } 279*14556e63SJeremy L Thompson output_project[m][j + P_from * i] /= interp_to[i + P_to * i]; 280a76a04e7SJeremy L Thompson } 281a76a04e7SJeremy L Thompson } 282*14556e63SJeremy L Thompson } 283*14556e63SJeremy L Thompson 284*14556e63SJeremy L Thompson // Cleanup 285a76a04e7SJeremy L Thompson ierr = CeedFree(&tau); CeedChk(ierr); 286a76a04e7SJeremy L Thompson ierr = CeedFree(&interp_to); CeedChk(ierr); 287a76a04e7SJeremy L Thompson ierr = CeedFree(&interp_from); CeedChk(ierr); 288a76a04e7SJeremy L Thompson 289a76a04e7SJeremy L Thompson return CEED_ERROR_SUCCESS; 290a76a04e7SJeremy L Thompson } 291a76a04e7SJeremy L Thompson 2927a982d89SJeremy L. Thompson /// @} 2937a982d89SJeremy L. Thompson 2947a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 2957a982d89SJeremy L. Thompson /// Ceed Backend API 2967a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 2977a982d89SJeremy L. Thompson /// @addtogroup CeedBasisBackend 2987a982d89SJeremy L. Thompson /// @{ 2997a982d89SJeremy L. Thompson 3007a982d89SJeremy L. Thompson /** 3017a982d89SJeremy L. Thompson @brief Return collocated grad matrix 3027a982d89SJeremy L. Thompson 3037a982d89SJeremy L. Thompson @param basis CeedBasis 304d1d35e2fSjeremylt @param[out] collo_grad_1d Row-major (Q_1d * Q_1d) matrix expressing derivatives of 3057a982d89SJeremy L. Thompson basis functions at quadrature points 3067a982d89SJeremy L. Thompson 3077a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3087a982d89SJeremy L. Thompson 3097a982d89SJeremy L. Thompson @ref Backend 3107a982d89SJeremy L. Thompson **/ 311d1d35e2fSjeremylt int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *collo_grad_1d) { 3127a982d89SJeremy L. Thompson int i, j, k; 3137a982d89SJeremy L. Thompson Ceed ceed; 314d1d35e2fSjeremylt CeedInt ierr, P_1d=(basis)->P_1d, Q_1d=(basis)->Q_1d; 31578464608Sjeremylt CeedScalar *interp_1d, *grad_1d, *tau; 3167a982d89SJeremy L. Thompson 317d1d35e2fSjeremylt ierr = CeedMalloc(Q_1d*P_1d, &interp_1d); CeedChk(ierr); 318d1d35e2fSjeremylt ierr = CeedMalloc(Q_1d*P_1d, &grad_1d); CeedChk(ierr); 31978464608Sjeremylt ierr = CeedMalloc(Q_1d, &tau); CeedChk(ierr); 320d1d35e2fSjeremylt memcpy(interp_1d, (basis)->interp_1d, Q_1d*P_1d*sizeof(basis)->interp_1d[0]); 321d1d35e2fSjeremylt memcpy(grad_1d, (basis)->grad_1d, Q_1d*P_1d*sizeof(basis)->interp_1d[0]); 3227a982d89SJeremy L. Thompson 323d1d35e2fSjeremylt // QR Factorization, interp_1d = Q R 3247a982d89SJeremy L. Thompson ierr = CeedBasisGetCeed(basis, &ceed); CeedChk(ierr); 325d1d35e2fSjeremylt ierr = CeedQRFactorization(ceed, interp_1d, tau, Q_1d, P_1d); CeedChk(ierr); 326e15f9bd0SJeremy L Thompson // Note: This function is for backend use, so all errors are terminal 327e15f9bd0SJeremy L Thompson // and we do not need to clean up memory on failure. 3287a982d89SJeremy L. Thompson 329d1d35e2fSjeremylt // Apply Rinv, collo_grad_1d = grad_1d Rinv 330d1d35e2fSjeremylt for (i=0; i<Q_1d; i++) { // Row i 331d1d35e2fSjeremylt collo_grad_1d[Q_1d*i] = grad_1d[P_1d*i]/interp_1d[0]; 332d1d35e2fSjeremylt for (j=1; j<P_1d; j++) { // Column j 333d1d35e2fSjeremylt collo_grad_1d[j+Q_1d*i] = grad_1d[j+P_1d*i]; 3347a982d89SJeremy L. Thompson for (k=0; k<j; k++) 335d1d35e2fSjeremylt collo_grad_1d[j+Q_1d*i] -= interp_1d[j+P_1d*k]*collo_grad_1d[k+Q_1d*i]; 336d1d35e2fSjeremylt collo_grad_1d[j+Q_1d*i] /= interp_1d[j+P_1d*j]; 3377a982d89SJeremy L. Thompson } 338d1d35e2fSjeremylt for (j=P_1d; j<Q_1d; j++) 339d1d35e2fSjeremylt collo_grad_1d[j+Q_1d*i] = 0; 3407a982d89SJeremy L. Thompson } 3417a982d89SJeremy L. Thompson 342673160d7Sjeremylt // Apply Qtranspose, collo_grad = collo_grad Q_transpose 343d1d35e2fSjeremylt ierr = CeedHouseholderApplyQ(collo_grad_1d, interp_1d, tau, CEED_NOTRANSPOSE, 344d1d35e2fSjeremylt Q_1d, Q_1d, P_1d, 1, Q_1d); CeedChk(ierr); 3457a982d89SJeremy L. Thompson 346d1d35e2fSjeremylt ierr = CeedFree(&interp_1d); CeedChk(ierr); 347d1d35e2fSjeremylt ierr = CeedFree(&grad_1d); CeedChk(ierr); 34878464608Sjeremylt ierr = CeedFree(&tau); CeedChk(ierr); 349e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3507a982d89SJeremy L. Thompson } 3517a982d89SJeremy L. Thompson 3527a982d89SJeremy L. Thompson /** 3537a982d89SJeremy L. Thompson @brief Get tensor status for given CeedBasis 3547a982d89SJeremy L. Thompson 3557a982d89SJeremy L. Thompson @param basis CeedBasis 356d1d35e2fSjeremylt @param[out] is_tensor Variable to store tensor status 3577a982d89SJeremy L. Thompson 3587a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3597a982d89SJeremy L. Thompson 3607a982d89SJeremy L. Thompson @ref Backend 3617a982d89SJeremy L. Thompson **/ 362d1d35e2fSjeremylt int CeedBasisIsTensor(CeedBasis basis, bool *is_tensor) { 363d1d35e2fSjeremylt *is_tensor = basis->tensor_basis; 364e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3657a982d89SJeremy L. Thompson } 3667a982d89SJeremy L. Thompson 3677a982d89SJeremy L. Thompson /** 3687a982d89SJeremy L. Thompson @brief Get backend data of a CeedBasis 3697a982d89SJeremy L. Thompson 3707a982d89SJeremy L. Thompson @param basis CeedBasis 3717a982d89SJeremy L. Thompson @param[out] data Variable to store data 3727a982d89SJeremy L. Thompson 3737a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3747a982d89SJeremy L. Thompson 3757a982d89SJeremy L. Thompson @ref Backend 3767a982d89SJeremy L. Thompson **/ 377777ff853SJeremy L Thompson int CeedBasisGetData(CeedBasis basis, void *data) { 378777ff853SJeremy L Thompson *(void **)data = basis->data; 379e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3807a982d89SJeremy L. Thompson } 3817a982d89SJeremy L. Thompson 3827a982d89SJeremy L. Thompson /** 3837a982d89SJeremy L. Thompson @brief Set backend data of a CeedBasis 3847a982d89SJeremy L. Thompson 3857a982d89SJeremy L. Thompson @param[out] basis CeedBasis 3867a982d89SJeremy L. Thompson @param data Data to set 3877a982d89SJeremy L. Thompson 3887a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3897a982d89SJeremy L. Thompson 3907a982d89SJeremy L. Thompson @ref Backend 3917a982d89SJeremy L. Thompson **/ 392777ff853SJeremy L Thompson int CeedBasisSetData(CeedBasis basis, void *data) { 393777ff853SJeremy L Thompson basis->data = data; 394e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3957a982d89SJeremy L. Thompson } 3967a982d89SJeremy L. Thompson 3977a982d89SJeremy L. Thompson /** 39834359f16Sjeremylt @brief Increment the reference counter for a CeedBasis 39934359f16Sjeremylt 40034359f16Sjeremylt @param basis Basis to increment the reference counter 40134359f16Sjeremylt 40234359f16Sjeremylt @return An error code: 0 - success, otherwise - failure 40334359f16Sjeremylt 40434359f16Sjeremylt @ref Backend 40534359f16Sjeremylt **/ 4069560d06aSjeremylt int CeedBasisReference(CeedBasis basis) { 40734359f16Sjeremylt basis->ref_count++; 40834359f16Sjeremylt return CEED_ERROR_SUCCESS; 40934359f16Sjeremylt } 41034359f16Sjeremylt 41134359f16Sjeremylt /** 4126e15d496SJeremy L Thompson @brief Estimate number of FLOPs required to apply CeedBasis in t_mode and eval_mode 4136e15d496SJeremy L Thompson 4146e15d496SJeremy L Thompson @param basis Basis to estimate FLOPs for 4156e15d496SJeremy L Thompson @param t_mode Apply basis or transpose 4166e15d496SJeremy L Thompson @param eval_mode Basis evaluation mode 4176e15d496SJeremy L Thompson @param flops Address of variable to hold FLOPs estimate 4186e15d496SJeremy L Thompson 4196e15d496SJeremy L Thompson @ref Backend 4206e15d496SJeremy L Thompson **/ 4216e15d496SJeremy L Thompson int CeedBasisGetFlopsEstimate(CeedBasis basis, CeedTransposeMode t_mode, 4229d36ca50SJeremy L Thompson CeedEvalMode eval_mode, CeedSize *flops) { 4236e15d496SJeremy L Thompson int ierr; 4246e15d496SJeremy L Thompson bool is_tensor; 4256e15d496SJeremy L Thompson 4266e15d496SJeremy L Thompson ierr = CeedBasisIsTensor(basis, &is_tensor); CeedChk(ierr); 4276e15d496SJeremy L Thompson if (is_tensor) { 4286e15d496SJeremy L Thompson CeedInt dim, num_comp, P_1d, Q_1d; 4296e15d496SJeremy L Thompson ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr); 4306e15d496SJeremy L Thompson ierr = CeedBasisGetNumComponents(basis, &num_comp); CeedChk(ierr); 4316e15d496SJeremy L Thompson ierr = CeedBasisGetNumNodes1D(basis, &P_1d); CeedChk(ierr); 4326e15d496SJeremy L Thompson ierr = CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d); CeedChk(ierr); 4336e15d496SJeremy L Thompson if (t_mode == CEED_TRANSPOSE) { 4346e15d496SJeremy L Thompson P_1d = Q_1d; Q_1d = P_1d; 4356e15d496SJeremy L Thompson } 4366e15d496SJeremy L Thompson CeedInt tensor_flops = 0, pre = num_comp * CeedIntPow(P_1d, dim-1), post = 1; 4376e15d496SJeremy L Thompson for (CeedInt d = 0; d < dim; d++) { 4386e15d496SJeremy L Thompson tensor_flops += 2 * pre * P_1d * post * Q_1d; 4396e15d496SJeremy L Thompson pre /= P_1d; 4406e15d496SJeremy L Thompson post *= Q_1d; 4416e15d496SJeremy L Thompson } 4426e15d496SJeremy L Thompson switch (eval_mode) { 4436e15d496SJeremy L Thompson case CEED_EVAL_NONE: *flops = 0; break; 4446e15d496SJeremy L Thompson case CEED_EVAL_INTERP: *flops = tensor_flops; break; 4456e15d496SJeremy L Thompson case CEED_EVAL_GRAD: *flops = tensor_flops * 2; break; 4466e15d496SJeremy L Thompson case CEED_EVAL_DIV: 4476e15d496SJeremy L Thompson // LCOV_EXCL_START 4486e15d496SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_INCOMPATIBLE, 4496e15d496SJeremy L Thompson "Tensor CEED_EVAL_DIV not supported"); break; 4506e15d496SJeremy L Thompson case CEED_EVAL_CURL: 4516e15d496SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_INCOMPATIBLE, 4526e15d496SJeremy L Thompson "Tensor CEED_EVAL_CURL not supported"); break; 4536e15d496SJeremy L Thompson // LCOV_EXCL_STOP 4546e15d496SJeremy L Thompson case CEED_EVAL_WEIGHT: *flops = dim * CeedIntPow(Q_1d, dim); break; 4556e15d496SJeremy L Thompson } 4566e15d496SJeremy L Thompson } else { 4576e15d496SJeremy L Thompson CeedInt dim, num_comp, num_nodes, num_qpts, Q_comp; 4586e15d496SJeremy L Thompson ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr); 4596e15d496SJeremy L Thompson ierr = CeedBasisGetNumComponents(basis, &num_comp); CeedChk(ierr); 4606e15d496SJeremy L Thompson ierr = CeedBasisGetNumNodes(basis, &num_nodes); CeedChk(ierr); 4616e15d496SJeremy L Thompson ierr = CeedBasisGetNumQuadraturePoints(basis, &num_qpts); CeedChk(ierr); 4626e15d496SJeremy L Thompson ierr = CeedBasisGetNumQuadratureComponents(basis, &Q_comp); CeedChk(ierr); 4636e15d496SJeremy L Thompson switch (eval_mode) { 4646e15d496SJeremy L Thompson case CEED_EVAL_NONE: *flops = 0; break; 4656e15d496SJeremy L Thompson case CEED_EVAL_INTERP: *flops = num_nodes * num_qpts * num_comp; break; 4666e15d496SJeremy L Thompson case CEED_EVAL_GRAD: *flops = num_nodes * num_qpts * num_comp * dim; break; 4676e15d496SJeremy L Thompson case CEED_EVAL_DIV: *flops = num_nodes * num_qpts; break; 4686e15d496SJeremy L Thompson case CEED_EVAL_CURL: *flops = num_nodes * num_qpts * dim; break; 4696e15d496SJeremy L Thompson case CEED_EVAL_WEIGHT: *flops = 0; break; 4706e15d496SJeremy L Thompson } 4716e15d496SJeremy L Thompson } 4726e15d496SJeremy L Thompson 4736e15d496SJeremy L Thompson return CEED_ERROR_SUCCESS; 4746e15d496SJeremy L Thompson } 4756e15d496SJeremy L Thompson 4766e15d496SJeremy L Thompson /** 4777a982d89SJeremy L. Thompson @brief Get dimension for given CeedElemTopology 4787a982d89SJeremy L. Thompson 4797a982d89SJeremy L. Thompson @param topo CeedElemTopology 4807a982d89SJeremy L. Thompson @param[out] dim Variable to store dimension of topology 4817a982d89SJeremy L. Thompson 4827a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 4837a982d89SJeremy L. Thompson 4847a982d89SJeremy L. Thompson @ref Backend 4857a982d89SJeremy L. Thompson **/ 4867a982d89SJeremy L. Thompson int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) { 4877a982d89SJeremy L. Thompson *dim = (CeedInt) topo >> 16; 488e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 4897a982d89SJeremy L. Thompson } 4907a982d89SJeremy L. Thompson 4917a982d89SJeremy L. Thompson /** 4927a982d89SJeremy L. Thompson @brief Get CeedTensorContract of a CeedBasis 4937a982d89SJeremy L. Thompson 4947a982d89SJeremy L. Thompson @param basis CeedBasis 4957a982d89SJeremy L. Thompson @param[out] contract Variable to store CeedTensorContract 4967a982d89SJeremy L. Thompson 4977a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 4987a982d89SJeremy L. Thompson 4997a982d89SJeremy L. Thompson @ref Backend 5007a982d89SJeremy L. Thompson **/ 5017a982d89SJeremy L. Thompson int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) { 5027a982d89SJeremy L. Thompson *contract = basis->contract; 503e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 5047a982d89SJeremy L. Thompson } 5057a982d89SJeremy L. Thompson 5067a982d89SJeremy L. Thompson /** 5077a982d89SJeremy L. Thompson @brief Set CeedTensorContract of a CeedBasis 5087a982d89SJeremy L. Thompson 5097a982d89SJeremy L. Thompson @param[out] basis CeedBasis 5107a982d89SJeremy L. Thompson @param contract CeedTensorContract to set 5117a982d89SJeremy L. Thompson 5127a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 5137a982d89SJeremy L. Thompson 5147a982d89SJeremy L. Thompson @ref Backend 5157a982d89SJeremy L. Thompson **/ 51634359f16Sjeremylt int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract contract) { 51734359f16Sjeremylt int ierr; 51834359f16Sjeremylt basis->contract = contract; 5199560d06aSjeremylt ierr = CeedTensorContractReference(contract); CeedChk(ierr); 520e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 5217a982d89SJeremy L. Thompson } 5227a982d89SJeremy L. Thompson 5237a982d89SJeremy L. Thompson /** 5247a982d89SJeremy L. Thompson @brief Return a reference implementation of matrix multiplication C = A B. 5257a982d89SJeremy L. Thompson Note, this is a reference implementation for CPU CeedScalar pointers 5267a982d89SJeremy L. Thompson that is not intended for high performance. 5277a982d89SJeremy L. Thompson 5287a982d89SJeremy L. Thompson @param ceed A Ceed context for error handling 529d1d35e2fSjeremylt @param[in] mat_A Row-major matrix A 530d1d35e2fSjeremylt @param[in] mat_B Row-major matrix B 531d1d35e2fSjeremylt @param[out] mat_C Row-major output matrix C 5327a982d89SJeremy L. Thompson @param m Number of rows of C 5337a982d89SJeremy L. Thompson @param n Number of columns of C 5347a982d89SJeremy L. Thompson @param kk Number of columns of A/rows of B 5357a982d89SJeremy L. Thompson 5367a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 5377a982d89SJeremy L. Thompson 5387a982d89SJeremy L. Thompson @ref Utility 5397a982d89SJeremy L. Thompson **/ 540ed9e99e6SJeremy L Thompson int CeedMatrixMatrixMultiply(Ceed ceed, const CeedScalar *mat_A, 541ed9e99e6SJeremy L Thompson const CeedScalar *mat_B, CeedScalar *mat_C, 542ed9e99e6SJeremy L Thompson CeedInt m, CeedInt n, CeedInt kk) { 5437a982d89SJeremy L. Thompson for (CeedInt i=0; i<m; i++) 5447a982d89SJeremy L. Thompson for (CeedInt j=0; j<n; j++) { 5457a982d89SJeremy L. Thompson CeedScalar sum = 0; 5467a982d89SJeremy L. Thompson for (CeedInt k=0; k<kk; k++) 547d1d35e2fSjeremylt sum += mat_A[k+i*kk]*mat_B[j+k*n]; 548d1d35e2fSjeremylt mat_C[j+i*n] = sum; 5497a982d89SJeremy L. Thompson } 550e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 5517a982d89SJeremy L. Thompson } 5527a982d89SJeremy L. Thompson 5537a982d89SJeremy L. Thompson /// @} 5547a982d89SJeremy L. Thompson 5557a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 5567a982d89SJeremy L. Thompson /// CeedBasis Public API 5577a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 5587a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser 559d7b241e6Sjeremylt /// @{ 560d7b241e6Sjeremylt 561b11c1e72Sjeremylt /** 56295bb1877Svaleriabarra @brief Create a tensor-product basis for H^1 discretizations 563b11c1e72Sjeremylt 564b11c1e72Sjeremylt @param ceed A Ceed object where the CeedBasis will be created 565b11c1e72Sjeremylt @param dim Topological dimension 566d1d35e2fSjeremylt @param num_comp Number of field components (1 for scalar fields) 567d1d35e2fSjeremylt @param P_1d Number of nodes in one dimension 568d1d35e2fSjeremylt @param Q_1d Number of quadrature points in one dimension 569d1d35e2fSjeremylt @param interp_1d Row-major (Q_1d * P_1d) matrix expressing the values of nodal 570b11c1e72Sjeremylt basis functions at quadrature points 571d1d35e2fSjeremylt @param grad_1d Row-major (Q_1d * P_1d) matrix expressing derivatives of nodal 572b11c1e72Sjeremylt basis functions at quadrature points 573d1d35e2fSjeremylt @param q_ref_1d Array of length Q_1d holding the locations of quadrature points 574b11c1e72Sjeremylt on the 1D reference element [-1, 1] 575d1d35e2fSjeremylt @param q_weight_1d Array of length Q_1d holding the quadrature weights on the 576b11c1e72Sjeremylt reference element 577b11c1e72Sjeremylt @param[out] basis Address of the variable where the newly created 578b11c1e72Sjeremylt CeedBasis will be stored. 579b11c1e72Sjeremylt 580b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 581dfdf5a53Sjeremylt 5827a982d89SJeremy L. Thompson @ref User 583b11c1e72Sjeremylt **/ 584d1d35e2fSjeremylt int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt num_comp, 585d1d35e2fSjeremylt CeedInt P_1d, CeedInt Q_1d, 586d1d35e2fSjeremylt const CeedScalar *interp_1d, 587d1d35e2fSjeremylt const CeedScalar *grad_1d, const CeedScalar *q_ref_1d, 588d1d35e2fSjeremylt const CeedScalar *q_weight_1d, CeedBasis *basis) { 589d7b241e6Sjeremylt int ierr; 590d7b241e6Sjeremylt 5915fe0d4faSjeremylt if (!ceed->BasisCreateTensorH1) { 5925fe0d4faSjeremylt Ceed delegate; 593aefd8378Sjeremylt ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr); 5945fe0d4faSjeremylt 5955fe0d4faSjeremylt if (!delegate) 596c042f62fSJeremy L Thompson // LCOV_EXCL_START 597e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_UNSUPPORTED, 598e15f9bd0SJeremy L Thompson "Backend does not support BasisCreateTensorH1"); 599c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 6005fe0d4faSjeremylt 601d1d35e2fSjeremylt ierr = CeedBasisCreateTensorH1(delegate, dim, num_comp, P_1d, 602d1d35e2fSjeremylt Q_1d, interp_1d, grad_1d, q_ref_1d, 603d1d35e2fSjeremylt q_weight_1d, basis); CeedChk(ierr); 604e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 6055fe0d4faSjeremylt } 606e15f9bd0SJeremy L Thompson 607e15f9bd0SJeremy L Thompson if (dim < 1) 608e15f9bd0SJeremy L Thompson // LCOV_EXCL_START 609e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, 610e15f9bd0SJeremy L Thompson "Basis dimension must be a positive value"); 611e15f9bd0SJeremy L Thompson // LCOV_EXCL_STOP 612227444bfSJeremy L Thompson 613227444bfSJeremy L Thompson if (num_comp < 1) 614227444bfSJeremy L Thompson // LCOV_EXCL_START 615227444bfSJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, 616227444bfSJeremy L Thompson "Basis must have at least 1 component"); 617227444bfSJeremy L Thompson // LCOV_EXCL_STOP 618227444bfSJeremy L Thompson 619227444bfSJeremy L Thompson if (P_1d < 1) 620227444bfSJeremy L Thompson // LCOV_EXCL_START 621227444bfSJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, 622227444bfSJeremy L Thompson "Basis must have at least 1 node"); 623227444bfSJeremy L Thompson // LCOV_EXCL_STOP 624227444bfSJeremy L Thompson 625227444bfSJeremy L Thompson if (Q_1d < 1) 626227444bfSJeremy L Thompson // LCOV_EXCL_START 627227444bfSJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, 628227444bfSJeremy L Thompson "Basis must have at least 1 quadrature point"); 629227444bfSJeremy L Thompson // LCOV_EXCL_STOP 630227444bfSJeremy L Thompson 63150c301a5SRezgar Shakeri CeedElemTopology topo = dim == 1 ? CEED_TOPOLOGY_LINE 63250c301a5SRezgar Shakeri : dim == 2 ? CEED_TOPOLOGY_QUAD 63350c301a5SRezgar Shakeri : CEED_TOPOLOGY_HEX; 634e15f9bd0SJeremy L Thompson 635d7b241e6Sjeremylt ierr = CeedCalloc(1, basis); CeedChk(ierr); 636d7b241e6Sjeremylt (*basis)->ceed = ceed; 6379560d06aSjeremylt ierr = CeedReference(ceed); CeedChk(ierr); 638d1d35e2fSjeremylt (*basis)->ref_count = 1; 639d1d35e2fSjeremylt (*basis)->tensor_basis = 1; 640d7b241e6Sjeremylt (*basis)->dim = dim; 641d99fa3c5SJeremy L Thompson (*basis)->topo = topo; 642d1d35e2fSjeremylt (*basis)->num_comp = num_comp; 643d1d35e2fSjeremylt (*basis)->P_1d = P_1d; 644d1d35e2fSjeremylt (*basis)->Q_1d = Q_1d; 645d1d35e2fSjeremylt (*basis)->P = CeedIntPow(P_1d, dim); 646d1d35e2fSjeremylt (*basis)->Q = CeedIntPow(Q_1d, dim); 64750c301a5SRezgar Shakeri (*basis)->Q_comp = 1; 64850c301a5SRezgar Shakeri (*basis)->basis_space = 1; // 1 for H^1 space 649ff3a0f91SJeremy L Thompson ierr = CeedCalloc(Q_1d, &(*basis)->q_ref_1d); CeedChk(ierr); 650ff3a0f91SJeremy L Thompson ierr = CeedCalloc(Q_1d, &(*basis)->q_weight_1d); CeedChk(ierr); 651ff3a0f91SJeremy L Thompson if (q_ref_1d) memcpy((*basis)->q_ref_1d, q_ref_1d, Q_1d*sizeof(q_ref_1d[0])); 652ff3a0f91SJeremy L Thompson if (q_weight_1d) memcpy((*basis)->q_weight_1d, q_weight_1d, 653ff3a0f91SJeremy L Thompson Q_1d*sizeof(q_weight_1d[0])); 654ff3a0f91SJeremy L Thompson ierr = CeedCalloc(Q_1d*P_1d, &(*basis)->interp_1d); CeedChk(ierr); 655ff3a0f91SJeremy L Thompson ierr = CeedCalloc(Q_1d*P_1d, &(*basis)->grad_1d); CeedChk(ierr); 656ff3a0f91SJeremy L Thompson if (interp_1d) memcpy((*basis)->interp_1d, interp_1d, 657ff3a0f91SJeremy L Thompson Q_1d*P_1d*sizeof(interp_1d[0])); 658ff3a0f91SJeremy L Thompson if (grad_1d) memcpy((*basis)->grad_1d, grad_1d, Q_1d*P_1d*sizeof(grad_1d[0])); 659d1d35e2fSjeremylt ierr = ceed->BasisCreateTensorH1(dim, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, 660d1d35e2fSjeremylt q_weight_1d, *basis); CeedChk(ierr); 661e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 662d7b241e6Sjeremylt } 663d7b241e6Sjeremylt 664b11c1e72Sjeremylt /** 66595bb1877Svaleriabarra @brief Create a tensor-product Lagrange basis 666b11c1e72Sjeremylt 667b11c1e72Sjeremylt @param ceed A Ceed object where the CeedBasis will be created 668b11c1e72Sjeremylt @param dim Topological dimension of element 669d1d35e2fSjeremylt @param num_comp Number of field components (1 for scalar fields) 670b11c1e72Sjeremylt @param P Number of Gauss-Lobatto nodes in one dimension. The 671b11c1e72Sjeremylt polynomial degree of the resulting Q_k element is k=P-1. 672b11c1e72Sjeremylt @param Q Number of quadrature points in one dimension. 673d1d35e2fSjeremylt @param quad_mode Distribution of the Q quadrature points (affects order of 674b11c1e72Sjeremylt accuracy for the quadrature) 675b11c1e72Sjeremylt @param[out] basis Address of the variable where the newly created 676b11c1e72Sjeremylt CeedBasis will be stored. 677b11c1e72Sjeremylt 678b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 679dfdf5a53Sjeremylt 6807a982d89SJeremy L. Thompson @ref User 681b11c1e72Sjeremylt **/ 682d1d35e2fSjeremylt int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt num_comp, 683d1d35e2fSjeremylt CeedInt P, CeedInt Q, CeedQuadMode quad_mode, 684692c2638Sjeremylt CeedBasis *basis) { 685d7b241e6Sjeremylt // Allocate 686e15f9bd0SJeremy L Thompson int ierr, ierr2, i, j, k; 687d1d35e2fSjeremylt CeedScalar c1, c2, c3, c4, dx, *nodes, *interp_1d, *grad_1d, *q_ref_1d, 688d1d35e2fSjeremylt *q_weight_1d; 6894d537eeaSYohann 6904d537eeaSYohann if (dim < 1) 691c042f62fSJeremy L Thompson // LCOV_EXCL_START 692e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, 693e15f9bd0SJeremy L Thompson "Basis dimension must be a positive value"); 694c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 6954d537eeaSYohann 696227444bfSJeremy L Thompson if (num_comp < 1) 697227444bfSJeremy L Thompson // LCOV_EXCL_START 698227444bfSJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, 699227444bfSJeremy L Thompson "Basis must have at least 1 component"); 700227444bfSJeremy L Thompson // LCOV_EXCL_STOP 701227444bfSJeremy L Thompson 702227444bfSJeremy L Thompson if (P < 1) 703227444bfSJeremy L Thompson // LCOV_EXCL_START 704227444bfSJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, 705227444bfSJeremy L Thompson "Basis must have at least 1 node"); 706227444bfSJeremy L Thompson // LCOV_EXCL_STOP 707227444bfSJeremy L Thompson 708227444bfSJeremy L Thompson if (Q < 1) 709227444bfSJeremy L Thompson // LCOV_EXCL_START 710227444bfSJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, 711227444bfSJeremy L Thompson "Basis must have at least 1 quadrature point"); 712227444bfSJeremy L Thompson // LCOV_EXCL_STOP 713227444bfSJeremy L Thompson 714e15f9bd0SJeremy L Thompson // Get Nodes and Weights 715d1d35e2fSjeremylt ierr = CeedCalloc(P*Q, &interp_1d); CeedChk(ierr); 716d1d35e2fSjeremylt ierr = CeedCalloc(P*Q, &grad_1d); CeedChk(ierr); 717d7b241e6Sjeremylt ierr = CeedCalloc(P, &nodes); CeedChk(ierr); 718d1d35e2fSjeremylt ierr = CeedCalloc(Q, &q_ref_1d); CeedChk(ierr); 719d1d35e2fSjeremylt ierr = CeedCalloc(Q, &q_weight_1d); CeedChk(ierr); 720e15f9bd0SJeremy L Thompson ierr = CeedLobattoQuadrature(P, nodes, NULL); 721e15f9bd0SJeremy L Thompson if (ierr) { goto cleanup; } CeedChk(ierr); 722d1d35e2fSjeremylt switch (quad_mode) { 723d7b241e6Sjeremylt case CEED_GAUSS: 724d1d35e2fSjeremylt ierr = CeedGaussQuadrature(Q, q_ref_1d, q_weight_1d); 725d7b241e6Sjeremylt break; 726d7b241e6Sjeremylt case CEED_GAUSS_LOBATTO: 727d1d35e2fSjeremylt ierr = CeedLobattoQuadrature(Q, q_ref_1d, q_weight_1d); 728d7b241e6Sjeremylt break; 729d7b241e6Sjeremylt } 730e15f9bd0SJeremy L Thompson if (ierr) { goto cleanup; } CeedChk(ierr); 731e15f9bd0SJeremy L Thompson 732d7b241e6Sjeremylt // Build B, D matrix 733d7b241e6Sjeremylt // Fornberg, 1998 734d7b241e6Sjeremylt for (i = 0; i < Q; i++) { 735d7b241e6Sjeremylt c1 = 1.0; 736d1d35e2fSjeremylt c3 = nodes[0] - q_ref_1d[i]; 737d1d35e2fSjeremylt interp_1d[i*P+0] = 1.0; 738d7b241e6Sjeremylt for (j = 1; j < P; j++) { 739d7b241e6Sjeremylt c2 = 1.0; 740d7b241e6Sjeremylt c4 = c3; 741d1d35e2fSjeremylt c3 = nodes[j] - q_ref_1d[i]; 742d7b241e6Sjeremylt for (k = 0; k < j; k++) { 743d7b241e6Sjeremylt dx = nodes[j] - nodes[k]; 744d7b241e6Sjeremylt c2 *= dx; 745d7b241e6Sjeremylt if (k == j - 1) { 746d1d35e2fSjeremylt grad_1d[i*P + j] = c1*(interp_1d[i*P + k] - c4*grad_1d[i*P + k]) / c2; 747d1d35e2fSjeremylt interp_1d[i*P + j] = - c1*c4*interp_1d[i*P + k] / c2; 748d7b241e6Sjeremylt } 749d1d35e2fSjeremylt grad_1d[i*P + k] = (c3*grad_1d[i*P + k] - interp_1d[i*P + k]) / dx; 750d1d35e2fSjeremylt interp_1d[i*P + k] = c3*interp_1d[i*P + k] / dx; 751d7b241e6Sjeremylt } 752d7b241e6Sjeremylt c1 = c2; 753d7b241e6Sjeremylt } 754d7b241e6Sjeremylt } 7559ac7b42eSJeremy L Thompson // Pass to CeedBasisCreateTensorH1 756d1d35e2fSjeremylt ierr = CeedBasisCreateTensorH1(ceed, dim, num_comp, P, Q, interp_1d, grad_1d, 7579ac7b42eSJeremy L Thompson q_ref_1d, q_weight_1d, basis); CeedChk(ierr); 758e15f9bd0SJeremy L Thompson cleanup: 759d1d35e2fSjeremylt ierr2 = CeedFree(&interp_1d); CeedChk(ierr2); 760d1d35e2fSjeremylt ierr2 = CeedFree(&grad_1d); CeedChk(ierr2); 761e15f9bd0SJeremy L Thompson ierr2 = CeedFree(&nodes); CeedChk(ierr2); 762d1d35e2fSjeremylt ierr2 = CeedFree(&q_ref_1d); CeedChk(ierr2); 763d1d35e2fSjeremylt ierr2 = CeedFree(&q_weight_1d); CeedChk(ierr2); 764e15f9bd0SJeremy L Thompson CeedChk(ierr); 765e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 766d7b241e6Sjeremylt } 767d7b241e6Sjeremylt 768b11c1e72Sjeremylt /** 76995bb1877Svaleriabarra @brief Create a non tensor-product basis for H^1 discretizations 770a8de75f0Sjeremylt 771a8de75f0Sjeremylt @param ceed A Ceed object where the CeedBasis will be created 772a8de75f0Sjeremylt @param topo Topology of element, e.g. hypercube, simplex, ect 773d1d35e2fSjeremylt @param num_comp Number of field components (1 for scalar fields) 774d1d35e2fSjeremylt @param num_nodes Total number of nodes 775d1d35e2fSjeremylt @param num_qpts Total number of quadrature points 776d1d35e2fSjeremylt @param interp Row-major (num_qpts * num_nodes) matrix expressing the values of 7778795c945Sjeremylt nodal basis functions at quadrature points 778d1d35e2fSjeremylt @param grad Row-major (num_qpts * dim * num_nodes) matrix expressing 7798795c945Sjeremylt derivatives of nodal basis functions at quadrature points 780d1d35e2fSjeremylt @param q_ref Array of length num_qpts holding the locations of quadrature 78150c301a5SRezgar Shakeri points on the reference element 782d1d35e2fSjeremylt @param q_weight Array of length num_qpts holding the quadrature weights on the 783a8de75f0Sjeremylt reference element 784a8de75f0Sjeremylt @param[out] basis Address of the variable where the newly created 785a8de75f0Sjeremylt CeedBasis will be stored. 786a8de75f0Sjeremylt 787a8de75f0Sjeremylt @return An error code: 0 - success, otherwise - failure 788a8de75f0Sjeremylt 7897a982d89SJeremy L. Thompson @ref User 790a8de75f0Sjeremylt **/ 791d1d35e2fSjeremylt int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, 792d1d35e2fSjeremylt CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 793d1d35e2fSjeremylt const CeedScalar *grad, const CeedScalar *q_ref, 794d1d35e2fSjeremylt const CeedScalar *q_weight, CeedBasis *basis) { 795a8de75f0Sjeremylt int ierr; 796d1d35e2fSjeremylt CeedInt P = num_nodes, Q = num_qpts, dim = 0; 797a8de75f0Sjeremylt 7985fe0d4faSjeremylt if (!ceed->BasisCreateH1) { 7995fe0d4faSjeremylt Ceed delegate; 800aefd8378Sjeremylt ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr); 8015fe0d4faSjeremylt 8025fe0d4faSjeremylt if (!delegate) 803c042f62fSJeremy L Thompson // LCOV_EXCL_START 804e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_UNSUPPORTED, 805e15f9bd0SJeremy L Thompson "Backend does not support BasisCreateH1"); 806c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 8075fe0d4faSjeremylt 808d1d35e2fSjeremylt ierr = CeedBasisCreateH1(delegate, topo, num_comp, num_nodes, 809d1d35e2fSjeremylt num_qpts, interp, grad, q_ref, 810d1d35e2fSjeremylt q_weight, basis); CeedChk(ierr); 811e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 8125fe0d4faSjeremylt } 8135fe0d4faSjeremylt 814227444bfSJeremy L Thompson if (num_comp < 1) 815227444bfSJeremy L Thompson // LCOV_EXCL_START 816227444bfSJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, 817227444bfSJeremy L Thompson "Basis must have at least 1 component"); 818227444bfSJeremy L Thompson // LCOV_EXCL_STOP 819227444bfSJeremy L Thompson 820227444bfSJeremy L Thompson if (num_nodes < 1) 821227444bfSJeremy L Thompson // LCOV_EXCL_START 822227444bfSJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, 823227444bfSJeremy L Thompson "Basis must have at least 1 node"); 824227444bfSJeremy L Thompson // LCOV_EXCL_STOP 825227444bfSJeremy L Thompson 826227444bfSJeremy L Thompson if (num_qpts < 1) 827227444bfSJeremy L Thompson // LCOV_EXCL_START 828227444bfSJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, 829227444bfSJeremy L Thompson "Basis must have at least 1 quadrature point"); 830227444bfSJeremy L Thompson // LCOV_EXCL_STOP 831227444bfSJeremy L Thompson 832a8de75f0Sjeremylt ierr = CeedCalloc(1, basis); CeedChk(ierr); 833a8de75f0Sjeremylt 834a8de75f0Sjeremylt ierr = CeedBasisGetTopologyDimension(topo, &dim); CeedChk(ierr); 835a8de75f0Sjeremylt 836a8de75f0Sjeremylt (*basis)->ceed = ceed; 8379560d06aSjeremylt ierr = CeedReference(ceed); CeedChk(ierr); 838d1d35e2fSjeremylt (*basis)->ref_count = 1; 839d1d35e2fSjeremylt (*basis)->tensor_basis = 0; 840a8de75f0Sjeremylt (*basis)->dim = dim; 841d99fa3c5SJeremy L Thompson (*basis)->topo = topo; 842d1d35e2fSjeremylt (*basis)->num_comp = num_comp; 843a8de75f0Sjeremylt (*basis)->P = P; 844a8de75f0Sjeremylt (*basis)->Q = Q; 84550c301a5SRezgar Shakeri (*basis)->Q_comp = 1; 84650c301a5SRezgar Shakeri (*basis)->basis_space = 1; // 1 for H^1 space 847ff3a0f91SJeremy L Thompson ierr = CeedCalloc(Q*dim, &(*basis)->q_ref_1d); CeedChk(ierr); 848ff3a0f91SJeremy L Thompson ierr = CeedCalloc(Q, &(*basis)->q_weight_1d); CeedChk(ierr); 849ff3a0f91SJeremy L Thompson if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q*dim*sizeof(q_ref[0])); 850ff3a0f91SJeremy L Thompson if(q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q*sizeof(q_weight[0])); 851ff3a0f91SJeremy L Thompson ierr = CeedCalloc(Q*P, &(*basis)->interp); CeedChk(ierr); 852ff3a0f91SJeremy L Thompson ierr = CeedCalloc(dim*Q*P, &(*basis)->grad); CeedChk(ierr); 853ff3a0f91SJeremy L Thompson if(interp) memcpy((*basis)->interp, interp, Q*P*sizeof(interp[0])); 854ff3a0f91SJeremy L Thompson if(grad) memcpy((*basis)->grad, grad, dim*Q*P*sizeof(grad[0])); 855d1d35e2fSjeremylt ierr = ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, q_ref, 856d1d35e2fSjeremylt q_weight, *basis); CeedChk(ierr); 857e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 858a8de75f0Sjeremylt } 859a8de75f0Sjeremylt 860a8de75f0Sjeremylt /** 86150c301a5SRezgar Shakeri @brief Create a non tensor-product basis for H(div) discretizations 86250c301a5SRezgar Shakeri 86350c301a5SRezgar Shakeri @param ceed A Ceed object where the CeedBasis will be created 86450c301a5SRezgar Shakeri @param topo Topology of element (`CEED_TOPOLOGY_QUAD`, `CEED_TOPOLOGY_PRISM`, etc.), 86550c301a5SRezgar Shakeri dimension of which is used in some array sizes below 86650c301a5SRezgar Shakeri @param num_comp Number of components (usually 1 for vectors in H(div) bases) 86750c301a5SRezgar Shakeri @param num_nodes Total number of nodes (dofs per element) 86850c301a5SRezgar Shakeri @param num_qpts Total number of quadrature points 86950c301a5SRezgar Shakeri @param interp Row-major (dim*num_qpts * num_nodes) matrix expressing the values of 87050c301a5SRezgar Shakeri nodal basis functions at quadrature points 87150c301a5SRezgar Shakeri @param div Row-major (num_qpts * num_nodes) matrix expressing 87250c301a5SRezgar Shakeri divergence of nodal basis functions at quadrature points 87350c301a5SRezgar Shakeri @param q_ref Array of length num_qpts holding the locations of quadrature 87450c301a5SRezgar Shakeri points on the reference element 87550c301a5SRezgar Shakeri @param q_weight Array of length num_qpts holding the quadrature weights on the 87650c301a5SRezgar Shakeri reference element 87750c301a5SRezgar Shakeri @param[out] basis Address of the variable where the newly created 87850c301a5SRezgar Shakeri CeedBasis will be stored. 87950c301a5SRezgar Shakeri 88050c301a5SRezgar Shakeri @return An error code: 0 - success, otherwise - failure 88150c301a5SRezgar Shakeri 88250c301a5SRezgar Shakeri @ref User 88350c301a5SRezgar Shakeri **/ 88450c301a5SRezgar Shakeri int CeedBasisCreateHdiv(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, 88550c301a5SRezgar Shakeri CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 88650c301a5SRezgar Shakeri const CeedScalar *div, const CeedScalar *q_ref, 88750c301a5SRezgar Shakeri const CeedScalar *q_weight, CeedBasis *basis) { 88850c301a5SRezgar Shakeri int ierr; 88950c301a5SRezgar Shakeri CeedInt Q = num_qpts, P = num_nodes, dim = 0; 89050c301a5SRezgar Shakeri ierr = CeedBasisGetTopologyDimension(topo, &dim); CeedChk(ierr); 89150c301a5SRezgar Shakeri if (!ceed->BasisCreateHdiv) { 89250c301a5SRezgar Shakeri Ceed delegate; 89350c301a5SRezgar Shakeri ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr); 89450c301a5SRezgar Shakeri 89550c301a5SRezgar Shakeri if (!delegate) 89650c301a5SRezgar Shakeri // LCOV_EXCL_START 89750c301a5SRezgar Shakeri return CeedError(ceed, CEED_ERROR_UNSUPPORTED, 89850c301a5SRezgar Shakeri "Backend does not implement BasisCreateHdiv"); 89950c301a5SRezgar Shakeri // LCOV_EXCL_STOP 90050c301a5SRezgar Shakeri 90150c301a5SRezgar Shakeri ierr = CeedBasisCreateHdiv(delegate, topo, num_comp, num_nodes, 90250c301a5SRezgar Shakeri num_qpts, interp, div, q_ref, 90350c301a5SRezgar Shakeri q_weight, basis); CeedChk(ierr); 90450c301a5SRezgar Shakeri return CEED_ERROR_SUCCESS; 90550c301a5SRezgar Shakeri } 90650c301a5SRezgar Shakeri 907227444bfSJeremy L Thompson if (num_comp < 1) 908227444bfSJeremy L Thompson // LCOV_EXCL_START 909227444bfSJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, 910227444bfSJeremy L Thompson "Basis must have at least 1 component"); 911227444bfSJeremy L Thompson // LCOV_EXCL_STOP 912227444bfSJeremy L Thompson 913227444bfSJeremy L Thompson if (num_nodes < 1) 914227444bfSJeremy L Thompson // LCOV_EXCL_START 915227444bfSJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, 916227444bfSJeremy L Thompson "Basis must have at least 1 node"); 917227444bfSJeremy L Thompson // LCOV_EXCL_STOP 918227444bfSJeremy L Thompson 919227444bfSJeremy L Thompson if (num_qpts < 1) 920227444bfSJeremy L Thompson // LCOV_EXCL_START 921227444bfSJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, 922227444bfSJeremy L Thompson "Basis must have at least 1 quadrature point"); 923227444bfSJeremy L Thompson // LCOV_EXCL_STOP 924227444bfSJeremy L Thompson 92550c301a5SRezgar Shakeri ierr = CeedCalloc(1, basis); CeedChk(ierr); 92650c301a5SRezgar Shakeri 92750c301a5SRezgar Shakeri (*basis)->ceed = ceed; 92850c301a5SRezgar Shakeri ierr = CeedReference(ceed); CeedChk(ierr); 92950c301a5SRezgar Shakeri (*basis)->ref_count = 1; 93050c301a5SRezgar Shakeri (*basis)->tensor_basis = 0; 93150c301a5SRezgar Shakeri (*basis)->dim = dim; 93250c301a5SRezgar Shakeri (*basis)->topo = topo; 93350c301a5SRezgar Shakeri (*basis)->num_comp = num_comp; 93450c301a5SRezgar Shakeri (*basis)->P = P; 93550c301a5SRezgar Shakeri (*basis)->Q = Q; 93650c301a5SRezgar Shakeri (*basis)->Q_comp = dim; 93750c301a5SRezgar Shakeri (*basis)->basis_space = 2; // 2 for H(div) space 93850c301a5SRezgar Shakeri ierr = CeedMalloc(Q*dim, &(*basis)->q_ref_1d); CeedChk(ierr); 93950c301a5SRezgar Shakeri ierr = CeedMalloc(Q, &(*basis)->q_weight_1d); CeedChk(ierr); 94050c301a5SRezgar Shakeri if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q*dim*sizeof(q_ref[0])); 94150c301a5SRezgar Shakeri if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q*sizeof(q_weight[0])); 94250c301a5SRezgar Shakeri ierr = CeedMalloc(dim*Q*P, &(*basis)->interp); CeedChk(ierr); 94350c301a5SRezgar Shakeri ierr = CeedMalloc(Q*P, &(*basis)->div); CeedChk(ierr); 94450c301a5SRezgar Shakeri if (interp) memcpy((*basis)->interp, interp, dim*Q*P*sizeof(interp[0])); 94550c301a5SRezgar Shakeri if (div) memcpy((*basis)->div, div, Q*P*sizeof(div[0])); 94650c301a5SRezgar Shakeri ierr = ceed->BasisCreateHdiv(topo, dim, P, Q, interp, div, q_ref, 94750c301a5SRezgar Shakeri q_weight, *basis); CeedChk(ierr); 94850c301a5SRezgar Shakeri return CEED_ERROR_SUCCESS; 94950c301a5SRezgar Shakeri } 95050c301a5SRezgar Shakeri 95150c301a5SRezgar Shakeri /** 952f113e5dcSJeremy L Thompson @brief Create a CeedBasis for projection from the nodes of `basis_from` 953*14556e63SJeremy L Thompson to the nodes of `basis_to`. Only `CEED_EVAL_INTERP` and 954*14556e63SJeremy L Thompson `CEED_EVAL_GRAD` will be valid for the new basis, `basis_project`. 955*14556e63SJeremy L Thompson The interpolation is given by `interp_project = interp_to^+ * interp_from`, 956*14556e63SJeremy L Thompson where the pesudoinverse `interp_to^+` is given by QR factorization. 957*14556e63SJeremy L Thompson The gradient is given by `grad_project = interp_to^+ * grad_from`. 958f113e5dcSJeremy L Thompson Note: `basis_from` and `basis_to` must have compatible quadrature 959f113e5dcSJeremy L Thompson spaces. 960446e7af4SJeremy L Thompson Note: `basis_project` will have the same number of components as 961446e7af4SJeremy L Thompson `basis_from`, regardless of the number of components that 962446e7af4SJeremy L Thompson `basis_to` has. If `basis_from` has 3 components and `basis_to` 963446e7af4SJeremy L Thompson has 5 components, then `basis_project` will have 3 components. 964f113e5dcSJeremy L Thompson 965f113e5dcSJeremy L Thompson @param[in] basis_from CeedBasis to prolong from 966446e7af4SJeremy L Thompson @param[in] basis_to CeedBasis to prolong to 967f113e5dcSJeremy L Thompson @param[out] basis_project Address of the variable where the newly created 968f113e5dcSJeremy L Thompson CeedBasis will be stored. 969f113e5dcSJeremy L Thompson 970f113e5dcSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 971f113e5dcSJeremy L Thompson 972f113e5dcSJeremy L Thompson @ref User 973f113e5dcSJeremy L Thompson **/ 974446e7af4SJeremy L Thompson int CeedBasisCreateProjection(CeedBasis basis_from, CeedBasis basis_to, 975f113e5dcSJeremy L Thompson CeedBasis *basis_project) { 976f113e5dcSJeremy L Thompson int ierr; 977f113e5dcSJeremy L Thompson Ceed ceed; 978f113e5dcSJeremy L Thompson ierr = CeedBasisGetCeed(basis_to, &ceed); CeedChk(ierr); 979f113e5dcSJeremy L Thompson 980f113e5dcSJeremy L Thompson // Create projectior matrix 981*14556e63SJeremy L Thompson CeedScalar *interp_project, *grad_project; 982*14556e63SJeremy L Thompson ierr = CeedBasisCreateProjectionMatrices(basis_from, basis_to, 983*14556e63SJeremy L Thompson &interp_project, &grad_project); 984*14556e63SJeremy L Thompson CeedChk(ierr); 985f113e5dcSJeremy L Thompson 986f113e5dcSJeremy L Thompson // Build basis 987f113e5dcSJeremy L Thompson bool is_tensor; 988f113e5dcSJeremy L Thompson CeedInt dim, num_comp; 989*14556e63SJeremy L Thompson CeedScalar *q_ref, *q_weight; 990f113e5dcSJeremy L Thompson ierr = CeedBasisIsTensor(basis_to, &is_tensor); CeedChk(ierr); 991f113e5dcSJeremy L Thompson ierr = CeedBasisGetDimension(basis_to, &dim); CeedChk(ierr); 992446e7af4SJeremy L Thompson ierr = CeedBasisGetNumComponents(basis_from, &num_comp); CeedChk(ierr); 993f113e5dcSJeremy L Thompson if (is_tensor) { 994f113e5dcSJeremy L Thompson CeedInt P_1d_to, P_1d_from; 995f113e5dcSJeremy L Thompson ierr = CeedBasisGetNumNodes1D(basis_from, &P_1d_from); CeedChk(ierr); 996f113e5dcSJeremy L Thompson ierr = CeedBasisGetNumNodes1D(basis_to, &P_1d_to); CeedChk(ierr); 997f113e5dcSJeremy L Thompson ierr = CeedCalloc(P_1d_to, &q_ref); CeedChk(ierr); 998f113e5dcSJeremy L Thompson ierr = CeedCalloc(P_1d_to, &q_weight); CeedChk(ierr); 999f113e5dcSJeremy L Thompson ierr = CeedBasisCreateTensorH1(ceed, dim, num_comp, P_1d_from, P_1d_to, 1000*14556e63SJeremy L Thompson interp_project, grad_project, q_ref, q_weight, basis_project); 1001f113e5dcSJeremy L Thompson CeedChk(ierr); 1002f113e5dcSJeremy L Thompson } else { 1003f113e5dcSJeremy L Thompson CeedElemTopology topo; 1004f113e5dcSJeremy L Thompson ierr = CeedBasisGetTopology(basis_to, &topo); CeedChk(ierr); 1005f113e5dcSJeremy L Thompson CeedInt num_nodes_to, num_nodes_from; 1006f113e5dcSJeremy L Thompson ierr = CeedBasisGetNumNodes(basis_from, &num_nodes_from); CeedChk(ierr); 1007f113e5dcSJeremy L Thompson ierr = CeedBasisGetNumNodes(basis_to, &num_nodes_to); CeedChk(ierr); 1008f113e5dcSJeremy L Thompson ierr = CeedCalloc(num_nodes_to * dim, &q_ref); CeedChk(ierr); 1009f113e5dcSJeremy L Thompson ierr = CeedCalloc(num_nodes_to, &q_weight); CeedChk(ierr); 1010f113e5dcSJeremy L Thompson ierr = CeedBasisCreateH1(ceed, topo, num_comp, num_nodes_from, num_nodes_to, 1011*14556e63SJeremy L Thompson interp_project, grad_project, q_ref, q_weight, basis_project); 1012f113e5dcSJeremy L Thompson CeedChk(ierr); 1013f113e5dcSJeremy L Thompson } 1014f113e5dcSJeremy L Thompson 1015f113e5dcSJeremy L Thompson // Cleanup 1016f113e5dcSJeremy L Thompson ierr = CeedFree(&interp_project); CeedChk(ierr); 1017*14556e63SJeremy L Thompson ierr = CeedFree(&grad_project); CeedChk(ierr); 1018f113e5dcSJeremy L Thompson ierr = CeedFree(&q_ref); CeedChk(ierr); 1019f113e5dcSJeremy L Thompson ierr = CeedFree(&q_weight); CeedChk(ierr); 1020f113e5dcSJeremy L Thompson 1021f113e5dcSJeremy L Thompson return CEED_ERROR_SUCCESS; 1022f113e5dcSJeremy L Thompson } 1023f113e5dcSJeremy L Thompson 1024f113e5dcSJeremy L Thompson /** 10259560d06aSjeremylt @brief Copy the pointer to a CeedBasis. Both pointers should 10269560d06aSjeremylt be destroyed with `CeedBasisDestroy()`; 10279560d06aSjeremylt Note: If `*basis_copy` is non-NULL, then it is assumed that 10289560d06aSjeremylt `*basis_copy` is a pointer to a CeedBasis. This CeedBasis 10299560d06aSjeremylt will be destroyed if `*basis_copy` is the only 10309560d06aSjeremylt reference to this CeedBasis. 10319560d06aSjeremylt 10329560d06aSjeremylt @param basis CeedBasis to copy reference to 10339560d06aSjeremylt @param[out] basis_copy Variable to store copied reference 10349560d06aSjeremylt 10359560d06aSjeremylt @return An error code: 0 - success, otherwise - failure 10369560d06aSjeremylt 10379560d06aSjeremylt @ref User 10389560d06aSjeremylt **/ 10399560d06aSjeremylt int CeedBasisReferenceCopy(CeedBasis basis, CeedBasis *basis_copy) { 10409560d06aSjeremylt int ierr; 10419560d06aSjeremylt 10429560d06aSjeremylt ierr = CeedBasisReference(basis); CeedChk(ierr); 10439560d06aSjeremylt ierr = CeedBasisDestroy(basis_copy); CeedChk(ierr); 10449560d06aSjeremylt *basis_copy = basis; 10459560d06aSjeremylt return CEED_ERROR_SUCCESS; 10469560d06aSjeremylt } 10479560d06aSjeremylt 10489560d06aSjeremylt /** 10497a982d89SJeremy L. Thompson @brief View a CeedBasis 10507a982d89SJeremy L. Thompson 10517a982d89SJeremy L. Thompson @param basis CeedBasis to view 10527a982d89SJeremy L. Thompson @param stream Stream to view to, e.g., stdout 10537a982d89SJeremy L. Thompson 10547a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 10557a982d89SJeremy L. Thompson 10567a982d89SJeremy L. Thompson @ref User 10577a982d89SJeremy L. Thompson **/ 10587a982d89SJeremy L. Thompson int CeedBasisView(CeedBasis basis, FILE *stream) { 10597a982d89SJeremy L. Thompson int ierr; 106050c301a5SRezgar Shakeri CeedFESpace FE_space = basis->basis_space; 106150c301a5SRezgar Shakeri CeedElemTopology topo = basis->topo; 106250c301a5SRezgar Shakeri // Print FE space and element topology of the basis 1063d1d35e2fSjeremylt if (basis->tensor_basis) { 1064990fdeb6SJeremy L Thompson fprintf(stream, "CeedBasis (%s on a %s element): dim=%" CeedInt_FMT " P=%" 1065990fdeb6SJeremy L Thompson CeedInt_FMT " Q=%" CeedInt_FMT "\n", 106650c301a5SRezgar Shakeri CeedFESpaces[FE_space], CeedElemTopologies[topo], 106750c301a5SRezgar Shakeri basis->dim, basis->P_1d, basis->Q_1d); 106850c301a5SRezgar Shakeri } else { 1069990fdeb6SJeremy L Thompson fprintf(stream, "CeedBasis (%s on a %s element): dim=%" CeedInt_FMT " P=%" 1070990fdeb6SJeremy L Thompson CeedInt_FMT " Q=%" CeedInt_FMT "\n", 107150c301a5SRezgar Shakeri CeedFESpaces[FE_space], CeedElemTopologies[topo], 107250c301a5SRezgar Shakeri basis->dim, basis->P, basis->Q); 107350c301a5SRezgar Shakeri } 107450c301a5SRezgar Shakeri // Print quadrature data, interpolation/gradient/divergene/curl of the basis 107550c301a5SRezgar Shakeri if (basis->tensor_basis) { // tensor basis 1076d1d35e2fSjeremylt ierr = CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_ref_1d, 10777a982d89SJeremy L. Thompson stream); CeedChk(ierr); 1078d1d35e2fSjeremylt ierr = CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q_1d, 1079d1d35e2fSjeremylt basis->q_weight_1d, stream); CeedChk(ierr); 1080d1d35e2fSjeremylt ierr = CeedScalarView("interp1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, 1081d1d35e2fSjeremylt basis->interp_1d, stream); CeedChk(ierr); 1082d1d35e2fSjeremylt ierr = CeedScalarView("grad1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, 1083d1d35e2fSjeremylt basis->grad_1d, stream); CeedChk(ierr); 108450c301a5SRezgar Shakeri } else { // non-tensor basis 10857a982d89SJeremy L. Thompson ierr = CeedScalarView("qref", "\t% 12.8f", 1, basis->Q*basis->dim, 1086d1d35e2fSjeremylt basis->q_ref_1d, 10877a982d89SJeremy L. Thompson stream); CeedChk(ierr); 1088d1d35e2fSjeremylt ierr = CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->q_weight_1d, 10897a982d89SJeremy L. Thompson stream); CeedChk(ierr); 109050c301a5SRezgar Shakeri ierr = CeedScalarView("interp", "\t% 12.8f", basis->Q_comp*basis->Q, basis->P, 10917a982d89SJeremy L. Thompson basis->interp, stream); CeedChk(ierr); 109250c301a5SRezgar Shakeri if (basis->grad) { 10937a982d89SJeremy L. Thompson ierr = CeedScalarView("grad", "\t% 12.8f", basis->dim*basis->Q, basis->P, 10947a982d89SJeremy L. Thompson basis->grad, stream); CeedChk(ierr); 10957a982d89SJeremy L. Thompson } 109650c301a5SRezgar Shakeri if (basis->div) { 109750c301a5SRezgar Shakeri ierr = CeedScalarView("div", "\t% 12.8f", basis->Q, basis->P, 109850c301a5SRezgar Shakeri basis->div, stream); CeedChk(ierr); 109950c301a5SRezgar Shakeri } 110050c301a5SRezgar Shakeri } 1101e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 11027a982d89SJeremy L. Thompson } 11037a982d89SJeremy L. Thompson 11047a982d89SJeremy L. Thompson /** 11057a982d89SJeremy L. Thompson @brief Apply basis evaluation from nodes to quadrature points or vice versa 11067a982d89SJeremy L. Thompson 11077a982d89SJeremy L. Thompson @param basis CeedBasis to evaluate 1108d1d35e2fSjeremylt @param num_elem The number of elements to apply the basis evaluation to; 11097a982d89SJeremy L. Thompson the backend will specify the ordering in 11104cc79fe7SJed Brown CeedElemRestrictionCreateBlocked() 1111d1d35e2fSjeremylt @param t_mode \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature 11127a982d89SJeremy L. Thompson points, \ref CEED_TRANSPOSE to apply the transpose, mapping 11137a982d89SJeremy L. Thompson from quadrature points to nodes 1114d1d35e2fSjeremylt @param eval_mode \ref CEED_EVAL_NONE to use values directly, 11157a982d89SJeremy L. Thompson \ref CEED_EVAL_INTERP to use interpolated values, 11167a982d89SJeremy L. Thompson \ref CEED_EVAL_GRAD to use gradients, 11177a982d89SJeremy L. Thompson \ref CEED_EVAL_WEIGHT to use quadrature weights. 11187a982d89SJeremy L. Thompson @param[in] u Input CeedVector 11197a982d89SJeremy L. Thompson @param[out] v Output CeedVector 11207a982d89SJeremy L. Thompson 11217a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 11227a982d89SJeremy L. Thompson 11237a982d89SJeremy L. Thompson @ref User 11247a982d89SJeremy L. Thompson **/ 1125d1d35e2fSjeremylt int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode, 1126d1d35e2fSjeremylt CeedEvalMode eval_mode, CeedVector u, CeedVector v) { 11277a982d89SJeremy L. Thompson int ierr; 11281f9221feSJeremy L Thompson CeedSize u_length = 0, v_length; 11291f9221feSJeremy L Thompson CeedInt dim, num_comp, num_nodes, num_qpts; 1130e15f9bd0SJeremy L Thompson ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr); 1131d1d35e2fSjeremylt ierr = CeedBasisGetNumComponents(basis, &num_comp); CeedChk(ierr); 1132d1d35e2fSjeremylt ierr = CeedBasisGetNumNodes(basis, &num_nodes); CeedChk(ierr); 1133d1d35e2fSjeremylt ierr = CeedBasisGetNumQuadraturePoints(basis, &num_qpts); CeedChk(ierr); 1134d1d35e2fSjeremylt ierr = CeedVectorGetLength(v, &v_length); CeedChk(ierr); 11357a982d89SJeremy L. Thompson if (u) { 1136d1d35e2fSjeremylt ierr = CeedVectorGetLength(u, &u_length); CeedChk(ierr); 11377a982d89SJeremy L. Thompson } 11387a982d89SJeremy L. Thompson 1139e15f9bd0SJeremy L Thompson if (!basis->Apply) 1140e15f9bd0SJeremy L Thompson // LCOV_EXCL_START 1141e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_UNSUPPORTED, 1142e15f9bd0SJeremy L Thompson "Backend does not support BasisApply"); 1143e15f9bd0SJeremy L Thompson // LCOV_EXCL_STOP 1144e15f9bd0SJeremy L Thompson 1145e15f9bd0SJeremy L Thompson // Check compatibility of topological and geometrical dimensions 1146d1d35e2fSjeremylt if ((t_mode == CEED_TRANSPOSE && (v_length%num_nodes != 0 || 1147d1d35e2fSjeremylt u_length%num_qpts != 0)) || 1148d1d35e2fSjeremylt (t_mode == CEED_NOTRANSPOSE && (u_length%num_nodes != 0 || 1149d1d35e2fSjeremylt v_length%num_qpts != 0))) 11508229195eSjeremylt // LCOV_EXCL_START 1151e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_DIMENSION, 1152e15f9bd0SJeremy L Thompson "Length of input/output vectors " 11537a982d89SJeremy L. Thompson "incompatible with basis dimensions"); 11548229195eSjeremylt // LCOV_EXCL_STOP 11557a982d89SJeremy L. Thompson 1156e15f9bd0SJeremy L Thompson // Check vector lengths to prevent out of bounds issues 1157d1d35e2fSjeremylt bool bad_dims = false; 1158d1d35e2fSjeremylt switch (eval_mode) { 1159e15f9bd0SJeremy L Thompson case CEED_EVAL_NONE: 1160d1d35e2fSjeremylt case CEED_EVAL_INTERP: bad_dims = 1161d1d35e2fSjeremylt ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts || 1162d1d35e2fSjeremylt v_length < num_elem*num_comp*num_nodes)) || 1163d1d35e2fSjeremylt (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp || 1164d1d35e2fSjeremylt u_length < num_elem*num_comp*num_nodes))); 1165e15f9bd0SJeremy L Thompson break; 1166d1d35e2fSjeremylt case CEED_EVAL_GRAD: bad_dims = 1167d1d35e2fSjeremylt ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts*dim || 1168d1d35e2fSjeremylt v_length < num_elem*num_comp*num_nodes)) || 1169d1d35e2fSjeremylt (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp*dim || 1170d1d35e2fSjeremylt u_length < num_elem*num_comp*num_nodes))); 1171e15f9bd0SJeremy L Thompson break; 1172e15f9bd0SJeremy L Thompson case CEED_EVAL_WEIGHT: 1173d1d35e2fSjeremylt bad_dims = v_length < num_elem*num_qpts; 1174e15f9bd0SJeremy L Thompson break; 1175e15f9bd0SJeremy L Thompson // LCOV_EXCL_START 1176d1d35e2fSjeremylt case CEED_EVAL_DIV: bad_dims = 1177d1d35e2fSjeremylt ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts || 1178d1d35e2fSjeremylt v_length < num_elem*num_comp*num_nodes)) || 1179d1d35e2fSjeremylt (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp || 1180d1d35e2fSjeremylt u_length < num_elem*num_comp*num_nodes))); 1181e15f9bd0SJeremy L Thompson break; 1182d1d35e2fSjeremylt case CEED_EVAL_CURL: bad_dims = 1183d1d35e2fSjeremylt ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts || 1184d1d35e2fSjeremylt v_length < num_elem*num_comp*num_nodes)) || 1185d1d35e2fSjeremylt (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp || 1186d1d35e2fSjeremylt u_length < num_elem*num_comp*num_nodes))); 1187e15f9bd0SJeremy L Thompson break; 1188e15f9bd0SJeremy L Thompson // LCOV_EXCL_STOP 1189e15f9bd0SJeremy L Thompson } 1190d1d35e2fSjeremylt if (bad_dims) 1191e15f9bd0SJeremy L Thompson // LCOV_EXCL_START 1192e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_DIMENSION, 1193d1d35e2fSjeremylt "Input/output vectors too short for basis and evaluation mode"); 1194e15f9bd0SJeremy L Thompson // LCOV_EXCL_STOP 1195e15f9bd0SJeremy L Thompson 1196d1d35e2fSjeremylt ierr = basis->Apply(basis, num_elem, t_mode, eval_mode, u, v); CeedChk(ierr); 1197e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 11987a982d89SJeremy L. Thompson } 11997a982d89SJeremy L. Thompson 12007a982d89SJeremy L. Thompson /** 1201b7c9bbdaSJeremy L Thompson @brief Get Ceed associated with a CeedBasis 1202b7c9bbdaSJeremy L Thompson 1203b7c9bbdaSJeremy L Thompson @param basis CeedBasis 1204b7c9bbdaSJeremy L Thompson @param[out] ceed Variable to store Ceed 1205b7c9bbdaSJeremy L Thompson 1206b7c9bbdaSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 1207b7c9bbdaSJeremy L Thompson 1208b7c9bbdaSJeremy L Thompson @ref Advanced 1209b7c9bbdaSJeremy L Thompson **/ 1210b7c9bbdaSJeremy L Thompson int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) { 1211b7c9bbdaSJeremy L Thompson *ceed = basis->ceed; 1212b7c9bbdaSJeremy L Thompson return CEED_ERROR_SUCCESS; 1213b7c9bbdaSJeremy L Thompson } 1214b7c9bbdaSJeremy L Thompson 1215b7c9bbdaSJeremy L Thompson /** 12169d007619Sjeremylt @brief Get dimension for given CeedBasis 12179d007619Sjeremylt 12189d007619Sjeremylt @param basis CeedBasis 12199d007619Sjeremylt @param[out] dim Variable to store dimension of basis 12209d007619Sjeremylt 12219d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 12229d007619Sjeremylt 1223b7c9bbdaSJeremy L Thompson @ref Advanced 12249d007619Sjeremylt **/ 12259d007619Sjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) { 12269d007619Sjeremylt *dim = basis->dim; 1227e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 12289d007619Sjeremylt } 12299d007619Sjeremylt 12309d007619Sjeremylt /** 1231d99fa3c5SJeremy L Thompson @brief Get topology for given CeedBasis 1232d99fa3c5SJeremy L Thompson 1233d99fa3c5SJeremy L Thompson @param basis CeedBasis 1234d99fa3c5SJeremy L Thompson @param[out] topo Variable to store topology of basis 1235d99fa3c5SJeremy L Thompson 1236d99fa3c5SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 1237d99fa3c5SJeremy L Thompson 1238b7c9bbdaSJeremy L Thompson @ref Advanced 1239d99fa3c5SJeremy L Thompson **/ 1240d99fa3c5SJeremy L Thompson int CeedBasisGetTopology(CeedBasis basis, CeedElemTopology *topo) { 1241d99fa3c5SJeremy L Thompson *topo = basis->topo; 1242e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1243d99fa3c5SJeremy L Thompson } 1244d99fa3c5SJeremy L Thompson 1245d99fa3c5SJeremy L Thompson /** 124650c301a5SRezgar Shakeri @brief Get number of Q-vector components for given CeedBasis 124750c301a5SRezgar Shakeri 124850c301a5SRezgar Shakeri @param basis CeedBasis 124950c301a5SRezgar Shakeri @param[out] Q_comp Variable to store number of Q-vector components of basis 125050c301a5SRezgar Shakeri 125150c301a5SRezgar Shakeri @return An error code: 0 - success, otherwise - failure 125250c301a5SRezgar Shakeri 125350c301a5SRezgar Shakeri @ref Advanced 125450c301a5SRezgar Shakeri **/ 125550c301a5SRezgar Shakeri int CeedBasisGetNumQuadratureComponents(CeedBasis basis, CeedInt *Q_comp) { 125650c301a5SRezgar Shakeri *Q_comp = basis->Q_comp; 125750c301a5SRezgar Shakeri return CEED_ERROR_SUCCESS; 125850c301a5SRezgar Shakeri } 125950c301a5SRezgar Shakeri 126050c301a5SRezgar Shakeri /** 12619d007619Sjeremylt @brief Get number of components for given CeedBasis 12629d007619Sjeremylt 12639d007619Sjeremylt @param basis CeedBasis 1264d1d35e2fSjeremylt @param[out] num_comp Variable to store number of components of basis 12659d007619Sjeremylt 12669d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 12679d007619Sjeremylt 1268b7c9bbdaSJeremy L Thompson @ref Advanced 12699d007619Sjeremylt **/ 1270d1d35e2fSjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *num_comp) { 1271d1d35e2fSjeremylt *num_comp = basis->num_comp; 1272e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 12739d007619Sjeremylt } 12749d007619Sjeremylt 12759d007619Sjeremylt /** 12769d007619Sjeremylt @brief Get total number of nodes (in dim dimensions) of a CeedBasis 12779d007619Sjeremylt 12789d007619Sjeremylt @param basis CeedBasis 12799d007619Sjeremylt @param[out] P Variable to store number of nodes 12809d007619Sjeremylt 12819d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 12829d007619Sjeremylt 12839d007619Sjeremylt @ref Utility 12849d007619Sjeremylt **/ 12859d007619Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) { 12869d007619Sjeremylt *P = basis->P; 1287e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 12889d007619Sjeremylt } 12899d007619Sjeremylt 12909d007619Sjeremylt /** 12919d007619Sjeremylt @brief Get total number of nodes (in 1 dimension) of a CeedBasis 12929d007619Sjeremylt 12939d007619Sjeremylt @param basis CeedBasis 1294d1d35e2fSjeremylt @param[out] P_1d Variable to store number of nodes 12959d007619Sjeremylt 12969d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 12979d007619Sjeremylt 1298b7c9bbdaSJeremy L Thompson @ref Advanced 12999d007619Sjeremylt **/ 1300d1d35e2fSjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P_1d) { 1301d1d35e2fSjeremylt if (!basis->tensor_basis) 13029d007619Sjeremylt // LCOV_EXCL_START 1303e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_MINOR, 1304d1d35e2fSjeremylt "Cannot supply P_1d for non-tensor basis"); 13059d007619Sjeremylt // LCOV_EXCL_STOP 13069d007619Sjeremylt 1307d1d35e2fSjeremylt *P_1d = basis->P_1d; 1308e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 13099d007619Sjeremylt } 13109d007619Sjeremylt 13119d007619Sjeremylt /** 13129d007619Sjeremylt @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis 13139d007619Sjeremylt 13149d007619Sjeremylt @param basis CeedBasis 13159d007619Sjeremylt @param[out] Q Variable to store number of quadrature points 13169d007619Sjeremylt 13179d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 13189d007619Sjeremylt 13199d007619Sjeremylt @ref Utility 13209d007619Sjeremylt **/ 13219d007619Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) { 13229d007619Sjeremylt *Q = basis->Q; 1323e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 13249d007619Sjeremylt } 13259d007619Sjeremylt 13269d007619Sjeremylt /** 13279d007619Sjeremylt @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis 13289d007619Sjeremylt 13299d007619Sjeremylt @param basis CeedBasis 1330d1d35e2fSjeremylt @param[out] Q_1d Variable to store number of quadrature points 13319d007619Sjeremylt 13329d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 13339d007619Sjeremylt 1334b7c9bbdaSJeremy L Thompson @ref Advanced 13359d007619Sjeremylt **/ 1336d1d35e2fSjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q_1d) { 1337d1d35e2fSjeremylt if (!basis->tensor_basis) 13389d007619Sjeremylt // LCOV_EXCL_START 1339e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_MINOR, 1340d1d35e2fSjeremylt "Cannot supply Q_1d for non-tensor basis"); 13419d007619Sjeremylt // LCOV_EXCL_STOP 13429d007619Sjeremylt 1343d1d35e2fSjeremylt *Q_1d = basis->Q_1d; 1344e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 13459d007619Sjeremylt } 13469d007619Sjeremylt 13479d007619Sjeremylt /** 13489d007619Sjeremylt @brief Get reference coordinates of quadrature points (in dim dimensions) 13499d007619Sjeremylt of a CeedBasis 13509d007619Sjeremylt 13519d007619Sjeremylt @param basis CeedBasis 1352d1d35e2fSjeremylt @param[out] q_ref Variable to store reference coordinates of quadrature points 13539d007619Sjeremylt 13549d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 13559d007619Sjeremylt 1356b7c9bbdaSJeremy L Thompson @ref Advanced 13579d007619Sjeremylt **/ 1358d1d35e2fSjeremylt int CeedBasisGetQRef(CeedBasis basis, const CeedScalar **q_ref) { 1359d1d35e2fSjeremylt *q_ref = basis->q_ref_1d; 1360e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 13619d007619Sjeremylt } 13629d007619Sjeremylt 13639d007619Sjeremylt /** 13649d007619Sjeremylt @brief Get quadrature weights of quadrature points (in dim dimensions) 13659d007619Sjeremylt of a CeedBasis 13669d007619Sjeremylt 13679d007619Sjeremylt @param basis CeedBasis 1368d1d35e2fSjeremylt @param[out] q_weight Variable to store quadrature weights 13699d007619Sjeremylt 13709d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 13719d007619Sjeremylt 1372b7c9bbdaSJeremy L Thompson @ref Advanced 13739d007619Sjeremylt **/ 1374d1d35e2fSjeremylt int CeedBasisGetQWeights(CeedBasis basis, const CeedScalar **q_weight) { 1375d1d35e2fSjeremylt *q_weight = basis->q_weight_1d; 1376e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 13779d007619Sjeremylt } 13789d007619Sjeremylt 13799d007619Sjeremylt /** 13809d007619Sjeremylt @brief Get interpolation matrix of a CeedBasis 13819d007619Sjeremylt 13829d007619Sjeremylt @param basis CeedBasis 13839d007619Sjeremylt @param[out] interp Variable to store interpolation matrix 13849d007619Sjeremylt 13859d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 13869d007619Sjeremylt 1387b7c9bbdaSJeremy L Thompson @ref Advanced 13889d007619Sjeremylt **/ 13896c58de82SJeremy L Thompson int CeedBasisGetInterp(CeedBasis basis, const CeedScalar **interp) { 1390d1d35e2fSjeremylt if (!basis->interp && basis->tensor_basis) { 13919d007619Sjeremylt // Allocate 13929d007619Sjeremylt int ierr; 13939d007619Sjeremylt ierr = CeedMalloc(basis->Q*basis->P, &basis->interp); CeedChk(ierr); 13949d007619Sjeremylt 13959d007619Sjeremylt // Initialize 13969d007619Sjeremylt for (CeedInt i=0; i<basis->Q*basis->P; i++) 13979d007619Sjeremylt basis->interp[i] = 1.0; 13989d007619Sjeremylt 13999d007619Sjeremylt // Calculate 14009d007619Sjeremylt for (CeedInt d=0; d<basis->dim; d++) 14019d007619Sjeremylt for (CeedInt qpt=0; qpt<basis->Q; qpt++) 14029d007619Sjeremylt for (CeedInt node=0; node<basis->P; node++) { 1403d1d35e2fSjeremylt CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d; 1404d1d35e2fSjeremylt CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d; 1405d1d35e2fSjeremylt basis->interp[qpt*(basis->P)+node] *= basis->interp_1d[q*basis->P_1d+p]; 14069d007619Sjeremylt } 14079d007619Sjeremylt } 14089d007619Sjeremylt *interp = basis->interp; 1409e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 14109d007619Sjeremylt } 14119d007619Sjeremylt 14129d007619Sjeremylt /** 14139d007619Sjeremylt @brief Get 1D interpolation matrix of a tensor product CeedBasis 14149d007619Sjeremylt 14159d007619Sjeremylt @param basis CeedBasis 1416d1d35e2fSjeremylt @param[out] interp_1d Variable to store interpolation matrix 14179d007619Sjeremylt 14189d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 14199d007619Sjeremylt 14209d007619Sjeremylt @ref Backend 14219d007619Sjeremylt **/ 1422d1d35e2fSjeremylt int CeedBasisGetInterp1D(CeedBasis basis, const CeedScalar **interp_1d) { 1423d1d35e2fSjeremylt if (!basis->tensor_basis) 14249d007619Sjeremylt // LCOV_EXCL_START 1425e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_MINOR, 1426e15f9bd0SJeremy L Thompson "CeedBasis is not a tensor product basis."); 14279d007619Sjeremylt // LCOV_EXCL_STOP 14289d007619Sjeremylt 1429d1d35e2fSjeremylt *interp_1d = basis->interp_1d; 1430e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 14319d007619Sjeremylt } 14329d007619Sjeremylt 14339d007619Sjeremylt /** 14349d007619Sjeremylt @brief Get gradient matrix of a CeedBasis 14359d007619Sjeremylt 14369d007619Sjeremylt @param basis CeedBasis 14379d007619Sjeremylt @param[out] grad Variable to store gradient matrix 14389d007619Sjeremylt 14399d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 14409d007619Sjeremylt 1441b7c9bbdaSJeremy L Thompson @ref Advanced 14429d007619Sjeremylt **/ 14436c58de82SJeremy L Thompson int CeedBasisGetGrad(CeedBasis basis, const CeedScalar **grad) { 1444d1d35e2fSjeremylt if (!basis->grad && basis->tensor_basis) { 14459d007619Sjeremylt // Allocate 14469d007619Sjeremylt int ierr; 14479d007619Sjeremylt ierr = CeedMalloc(basis->dim*basis->Q*basis->P, &basis->grad); 14489d007619Sjeremylt CeedChk(ierr); 14499d007619Sjeremylt 14509d007619Sjeremylt // Initialize 14519d007619Sjeremylt for (CeedInt i=0; i<basis->dim*basis->Q*basis->P; i++) 14529d007619Sjeremylt basis->grad[i] = 1.0; 14539d007619Sjeremylt 14549d007619Sjeremylt // Calculate 14559d007619Sjeremylt for (CeedInt d=0; d<basis->dim; d++) 14569d007619Sjeremylt for (CeedInt i=0; i<basis->dim; i++) 14579d007619Sjeremylt for (CeedInt qpt=0; qpt<basis->Q; qpt++) 14589d007619Sjeremylt for (CeedInt node=0; node<basis->P; node++) { 1459d1d35e2fSjeremylt CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d; 1460d1d35e2fSjeremylt CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d; 14619d007619Sjeremylt if (i == d) 14629d007619Sjeremylt basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *= 1463d1d35e2fSjeremylt basis->grad_1d[q*basis->P_1d+p]; 14649d007619Sjeremylt else 14659d007619Sjeremylt basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *= 1466d1d35e2fSjeremylt basis->interp_1d[q*basis->P_1d+p]; 14679d007619Sjeremylt } 14689d007619Sjeremylt } 14699d007619Sjeremylt *grad = basis->grad; 1470e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 14719d007619Sjeremylt } 14729d007619Sjeremylt 14739d007619Sjeremylt /** 14749d007619Sjeremylt @brief Get 1D gradient matrix of a tensor product CeedBasis 14759d007619Sjeremylt 14769d007619Sjeremylt @param basis CeedBasis 1477d1d35e2fSjeremylt @param[out] grad_1d Variable to store gradient matrix 14789d007619Sjeremylt 14799d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 14809d007619Sjeremylt 1481b7c9bbdaSJeremy L Thompson @ref Advanced 14829d007619Sjeremylt **/ 1483d1d35e2fSjeremylt int CeedBasisGetGrad1D(CeedBasis basis, const CeedScalar **grad_1d) { 1484d1d35e2fSjeremylt if (!basis->tensor_basis) 14859d007619Sjeremylt // LCOV_EXCL_START 1486e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_MINOR, 1487e15f9bd0SJeremy L Thompson "CeedBasis is not a tensor product basis."); 14889d007619Sjeremylt // LCOV_EXCL_STOP 14899d007619Sjeremylt 1490d1d35e2fSjeremylt *grad_1d = basis->grad_1d; 1491e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 14929d007619Sjeremylt } 14939d007619Sjeremylt 14949d007619Sjeremylt /** 149550c301a5SRezgar Shakeri @brief Get divergence matrix of a CeedBasis 149650c301a5SRezgar Shakeri 149750c301a5SRezgar Shakeri @param basis CeedBasis 149850c301a5SRezgar Shakeri @param[out] div Variable to store divergence matrix 149950c301a5SRezgar Shakeri 150050c301a5SRezgar Shakeri @return An error code: 0 - success, otherwise - failure 150150c301a5SRezgar Shakeri 150250c301a5SRezgar Shakeri @ref Advanced 150350c301a5SRezgar Shakeri **/ 150450c301a5SRezgar Shakeri int CeedBasisGetDiv(CeedBasis basis, const CeedScalar **div) { 150550c301a5SRezgar Shakeri if (!basis->div) 150650c301a5SRezgar Shakeri // LCOV_EXCL_START 150750c301a5SRezgar Shakeri return CeedError(basis->ceed, CEED_ERROR_MINOR, 150850c301a5SRezgar Shakeri "CeedBasis does not have divergence matrix."); 150950c301a5SRezgar Shakeri // LCOV_EXCL_STOP 151050c301a5SRezgar Shakeri 151150c301a5SRezgar Shakeri *div = basis->div; 151250c301a5SRezgar Shakeri return CEED_ERROR_SUCCESS; 151350c301a5SRezgar Shakeri } 151450c301a5SRezgar Shakeri 151550c301a5SRezgar Shakeri /** 15167a982d89SJeremy L. Thompson @brief Destroy a CeedBasis 15177a982d89SJeremy L. Thompson 15187a982d89SJeremy L. Thompson @param basis CeedBasis to destroy 15197a982d89SJeremy L. Thompson 15207a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 15217a982d89SJeremy L. Thompson 15227a982d89SJeremy L. Thompson @ref User 15237a982d89SJeremy L. Thompson **/ 15247a982d89SJeremy L. Thompson int CeedBasisDestroy(CeedBasis *basis) { 15257a982d89SJeremy L. Thompson int ierr; 15267a982d89SJeremy L. Thompson 1527d1d35e2fSjeremylt if (!*basis || --(*basis)->ref_count > 0) return CEED_ERROR_SUCCESS; 15287a982d89SJeremy L. Thompson if ((*basis)->Destroy) { 15297a982d89SJeremy L. Thompson ierr = (*basis)->Destroy(*basis); CeedChk(ierr); 15307a982d89SJeremy L. Thompson } 153134359f16Sjeremylt if ((*basis)->contract) { 153234359f16Sjeremylt ierr = CeedTensorContractDestroy(&(*basis)->contract); CeedChk(ierr); 153334359f16Sjeremylt } 15347a982d89SJeremy L. Thompson ierr = CeedFree(&(*basis)->interp); CeedChk(ierr); 1535d1d35e2fSjeremylt ierr = CeedFree(&(*basis)->interp_1d); CeedChk(ierr); 15367a982d89SJeremy L. Thompson ierr = CeedFree(&(*basis)->grad); CeedChk(ierr); 153750c301a5SRezgar Shakeri ierr = CeedFree(&(*basis)->div); CeedChk(ierr); 1538d1d35e2fSjeremylt ierr = CeedFree(&(*basis)->grad_1d); CeedChk(ierr); 1539d1d35e2fSjeremylt ierr = CeedFree(&(*basis)->q_ref_1d); CeedChk(ierr); 1540d1d35e2fSjeremylt ierr = CeedFree(&(*basis)->q_weight_1d); CeedChk(ierr); 15417a982d89SJeremy L. Thompson ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr); 15427a982d89SJeremy L. Thompson ierr = CeedFree(basis); CeedChk(ierr); 1543e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 15447a982d89SJeremy L. Thompson } 15457a982d89SJeremy L. Thompson 15467a982d89SJeremy L. Thompson /** 1547b11c1e72Sjeremylt @brief Construct a Gauss-Legendre quadrature 1548b11c1e72Sjeremylt 1549b11c1e72Sjeremylt @param Q Number of quadrature points (integrates polynomials of 1550b11c1e72Sjeremylt degree 2*Q-1 exactly) 1551d1d35e2fSjeremylt @param[out] q_ref_1d Array of length Q to hold the abscissa on [-1, 1] 1552d1d35e2fSjeremylt @param[out] q_weight_1d Array of length Q to hold the weights 1553b11c1e72Sjeremylt 1554b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 1555dfdf5a53Sjeremylt 1556dfdf5a53Sjeremylt @ref Utility 1557b11c1e72Sjeremylt **/ 1558d1d35e2fSjeremylt int CeedGaussQuadrature(CeedInt Q, CeedScalar *q_ref_1d, 1559d1d35e2fSjeremylt CeedScalar *q_weight_1d) { 1560d7b241e6Sjeremylt // Allocate 1561d7b241e6Sjeremylt CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0*atan(1.0); 1562d1d35e2fSjeremylt // Build q_ref_1d, q_weight_1d 156392ae7e47SJeremy L Thompson for (CeedInt i = 0; i <= Q/2; i++) { 1564d7b241e6Sjeremylt // Guess 1565d7b241e6Sjeremylt xi = cos(PI*(CeedScalar)(2*i+1)/((CeedScalar)(2*Q))); 1566d7b241e6Sjeremylt // Pn(xi) 1567d7b241e6Sjeremylt P0 = 1.0; 1568d7b241e6Sjeremylt P1 = xi; 1569d7b241e6Sjeremylt P2 = 0.0; 157092ae7e47SJeremy L Thompson for (CeedInt j = 2; j <= Q; j++) { 1571d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 1572d7b241e6Sjeremylt P0 = P1; 1573d7b241e6Sjeremylt P1 = P2; 1574d7b241e6Sjeremylt } 1575d7b241e6Sjeremylt // First Newton Step 1576d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 1577d7b241e6Sjeremylt xi = xi-P2/dP2; 1578d7b241e6Sjeremylt // Newton to convergence 157992ae7e47SJeremy L Thompson for (CeedInt k=0; k<100 && fabs(P2)>10*CEED_EPSILON; k++) { 1580d7b241e6Sjeremylt P0 = 1.0; 1581d7b241e6Sjeremylt P1 = xi; 158292ae7e47SJeremy L Thompson for (CeedInt j = 2; j <= Q; j++) { 1583d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 1584d7b241e6Sjeremylt P0 = P1; 1585d7b241e6Sjeremylt P1 = P2; 1586d7b241e6Sjeremylt } 1587d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 1588d7b241e6Sjeremylt xi = xi-P2/dP2; 1589d7b241e6Sjeremylt } 1590d7b241e6Sjeremylt // Save xi, wi 1591d7b241e6Sjeremylt wi = 2.0/((1.0-xi*xi)*dP2*dP2); 1592d1d35e2fSjeremylt q_weight_1d[i] = wi; 1593d1d35e2fSjeremylt q_weight_1d[Q-1-i] = wi; 1594d1d35e2fSjeremylt q_ref_1d[i] = -xi; 1595d1d35e2fSjeremylt q_ref_1d[Q-1-i]= xi; 1596d7b241e6Sjeremylt } 1597e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1598d7b241e6Sjeremylt } 1599d7b241e6Sjeremylt 1600b11c1e72Sjeremylt /** 1601b11c1e72Sjeremylt @brief Construct a Gauss-Legendre-Lobatto quadrature 1602b11c1e72Sjeremylt 1603b11c1e72Sjeremylt @param Q Number of quadrature points (integrates polynomials of 1604b11c1e72Sjeremylt degree 2*Q-3 exactly) 1605d1d35e2fSjeremylt @param[out] q_ref_1d Array of length Q to hold the abscissa on [-1, 1] 1606d1d35e2fSjeremylt @param[out] q_weight_1d Array of length Q to hold the weights 1607b11c1e72Sjeremylt 1608b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 1609dfdf5a53Sjeremylt 1610dfdf5a53Sjeremylt @ref Utility 1611b11c1e72Sjeremylt **/ 1612d1d35e2fSjeremylt int CeedLobattoQuadrature(CeedInt Q, CeedScalar *q_ref_1d, 1613d1d35e2fSjeremylt CeedScalar *q_weight_1d) { 1614d7b241e6Sjeremylt // Allocate 1615d7b241e6Sjeremylt CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0*atan(1.0); 1616d1d35e2fSjeremylt // Build q_ref_1d, q_weight_1d 1617d7b241e6Sjeremylt // Set endpoints 161830a100c3SJed Brown if (Q < 2) 1619b0d62198Sjeremylt // LCOV_EXCL_START 1620e15f9bd0SJeremy L Thompson return CeedError(NULL, CEED_ERROR_DIMENSION, 1621990fdeb6SJeremy L Thompson "Cannot create Lobatto quadrature with Q=%" CeedInt_FMT " < 2 points", Q); 1622b0d62198Sjeremylt // LCOV_EXCL_STOP 1623d7b241e6Sjeremylt wi = 2.0/((CeedScalar)(Q*(Q-1))); 1624d1d35e2fSjeremylt if (q_weight_1d) { 1625d1d35e2fSjeremylt q_weight_1d[0] = wi; 1626d1d35e2fSjeremylt q_weight_1d[Q-1] = wi; 1627d7b241e6Sjeremylt } 1628d1d35e2fSjeremylt q_ref_1d[0] = -1.0; 1629d1d35e2fSjeremylt q_ref_1d[Q-1] = 1.0; 1630d7b241e6Sjeremylt // Interior 163192ae7e47SJeremy L Thompson for (CeedInt i = 1; i <= (Q-1)/2; i++) { 1632d7b241e6Sjeremylt // Guess 1633d7b241e6Sjeremylt xi = cos(PI*(CeedScalar)(i)/(CeedScalar)(Q-1)); 1634d7b241e6Sjeremylt // Pn(xi) 1635d7b241e6Sjeremylt P0 = 1.0; 1636d7b241e6Sjeremylt P1 = xi; 1637d7b241e6Sjeremylt P2 = 0.0; 163892ae7e47SJeremy L Thompson for (CeedInt j = 2; j < Q; j++) { 1639d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 1640d7b241e6Sjeremylt P0 = P1; 1641d7b241e6Sjeremylt P1 = P2; 1642d7b241e6Sjeremylt } 1643d7b241e6Sjeremylt // First Newton step 1644d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 1645d7b241e6Sjeremylt d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi); 1646d7b241e6Sjeremylt xi = xi-dP2/d2P2; 1647d7b241e6Sjeremylt // Newton to convergence 164892ae7e47SJeremy L Thompson for (CeedInt k=0; k<100 && fabs(dP2)>10*CEED_EPSILON; k++) { 1649d7b241e6Sjeremylt P0 = 1.0; 1650d7b241e6Sjeremylt P1 = xi; 165192ae7e47SJeremy L Thompson for (CeedInt j = 2; j < Q; j++) { 1652d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 1653d7b241e6Sjeremylt P0 = P1; 1654d7b241e6Sjeremylt P1 = P2; 1655d7b241e6Sjeremylt } 1656d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 1657d7b241e6Sjeremylt d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi); 1658d7b241e6Sjeremylt xi = xi-dP2/d2P2; 1659d7b241e6Sjeremylt } 1660d7b241e6Sjeremylt // Save xi, wi 1661d7b241e6Sjeremylt wi = 2.0/(((CeedScalar)(Q*(Q-1)))*P2*P2); 1662d1d35e2fSjeremylt if (q_weight_1d) { 1663d1d35e2fSjeremylt q_weight_1d[i] = wi; 1664d1d35e2fSjeremylt q_weight_1d[Q-1-i] = wi; 1665d7b241e6Sjeremylt } 1666d1d35e2fSjeremylt q_ref_1d[i] = -xi; 1667d1d35e2fSjeremylt q_ref_1d[Q-1-i]= xi; 1668d7b241e6Sjeremylt } 1669e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1670d7b241e6Sjeremylt } 1671d7b241e6Sjeremylt 1672dfdf5a53Sjeremylt /** 167395bb1877Svaleriabarra @brief Return QR Factorization of a matrix 1674b11c1e72Sjeremylt 167577645d7bSjeremylt @param ceed A Ceed context for error handling 167652bfb9bbSJeremy L Thompson @param[in,out] mat Row-major matrix to be factorized in place 167752bfb9bbSJeremy L Thompson @param[in,out] tau Vector of length m of scaling factors 1678b11c1e72Sjeremylt @param m Number of rows 1679b11c1e72Sjeremylt @param n Number of columns 1680b11c1e72Sjeremylt 1681b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 1682dfdf5a53Sjeremylt 1683dfdf5a53Sjeremylt @ref Utility 1684b11c1e72Sjeremylt **/ 1685a7bd39daSjeremylt int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau, 1686d7b241e6Sjeremylt CeedInt m, CeedInt n) { 1687d7b241e6Sjeremylt CeedScalar v[m]; 1688d7b241e6Sjeremylt 1689a7bd39daSjeremylt // Check m >= n 1690a7bd39daSjeremylt if (n > m) 1691c042f62fSJeremy L Thompson // LCOV_EXCL_START 1692e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_UNSUPPORTED, 1693e15f9bd0SJeremy L Thompson "Cannot compute QR factorization with n > m"); 1694c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 1695a7bd39daSjeremylt 169652bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) { 1697bde37e8cSJed Brown if (i >= m-1) { // last row of matrix, no reflection needed 1698bde37e8cSJed Brown tau[i] = 0.; 1699bde37e8cSJed Brown break; 1700bde37e8cSJed Brown } 1701d7b241e6Sjeremylt // Calculate Householder vector, magnitude 1702d7b241e6Sjeremylt CeedScalar sigma = 0.0; 1703d7b241e6Sjeremylt v[i] = mat[i+n*i]; 170452bfb9bbSJeremy L Thompson for (CeedInt j=i+1; j<m; j++) { 1705d7b241e6Sjeremylt v[j] = mat[i+n*j]; 1706d7b241e6Sjeremylt sigma += v[j] * v[j]; 1707d7b241e6Sjeremylt } 1708d7b241e6Sjeremylt CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m] 1709673160d7Sjeremylt CeedScalar R_ii = -copysign(norm, v[i]); 1710673160d7Sjeremylt v[i] -= R_ii; 1711d7b241e6Sjeremylt // norm of v[i:m] after modification above and scaling below 1712d7b241e6Sjeremylt // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 1713d7b241e6Sjeremylt // tau = 2 / (norm*norm) 1714d7b241e6Sjeremylt tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma); 17151d102b48SJeremy L Thompson for (CeedInt j=i+1; j<m; j++) 17161d102b48SJeremy L Thompson v[j] /= v[i]; 1717d7b241e6Sjeremylt 1718d7b241e6Sjeremylt // Apply Householder reflector to lower right panel 1719d7b241e6Sjeremylt CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1); 1720d7b241e6Sjeremylt // Save v 1721673160d7Sjeremylt mat[i+n*i] = R_ii; 17221d102b48SJeremy L Thompson for (CeedInt j=i+1; j<m; j++) 1723d7b241e6Sjeremylt mat[i+n*j] = v[j]; 1724d7b241e6Sjeremylt } 1725e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1726d7b241e6Sjeremylt } 1727d7b241e6Sjeremylt 1728b11c1e72Sjeremylt /** 172952bfb9bbSJeremy L Thompson @brief Return symmetric Schur decomposition of the symmetric matrix mat via 173052bfb9bbSJeremy L Thompson symmetric QR factorization 173152bfb9bbSJeremy L Thompson 173277645d7bSjeremylt @param ceed A Ceed context for error handling 173352bfb9bbSJeremy L Thompson @param[in,out] mat Row-major matrix to be factorized in place 1734460bf743SValeria Barra @param[out] lambda Vector of length n of eigenvalues 173552bfb9bbSJeremy L Thompson @param n Number of rows/columns 173652bfb9bbSJeremy L Thompson 173752bfb9bbSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 173852bfb9bbSJeremy L Thompson 173952bfb9bbSJeremy L Thompson @ref Utility 174052bfb9bbSJeremy L Thompson **/ 174103d18186Sjeremylt CeedPragmaOptimizeOff 174252bfb9bbSJeremy L Thompson int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat, 174352bfb9bbSJeremy L Thompson CeedScalar *lambda, CeedInt n) { 174452bfb9bbSJeremy L Thompson // Check bounds for clang-tidy 174552bfb9bbSJeremy L Thompson if (n<2) 1746c042f62fSJeremy L Thompson // LCOV_EXCL_START 1747e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_UNSUPPORTED, 1748c042f62fSJeremy L Thompson "Cannot compute symmetric Schur decomposition of scalars"); 1749c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 175052bfb9bbSJeremy L Thompson 1751673160d7Sjeremylt CeedScalar v[n-1], tau[n-1], mat_T[n*n]; 175252bfb9bbSJeremy L Thompson 1753673160d7Sjeremylt // Copy mat to mat_T and set mat to I 1754673160d7Sjeremylt memcpy(mat_T, mat, n*n*sizeof(mat[0])); 175552bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) 175652bfb9bbSJeremy L Thompson for (CeedInt j=0; j<n; j++) 175752bfb9bbSJeremy L Thompson mat[j+n*i] = (i==j) ? 1 : 0; 175852bfb9bbSJeremy L Thompson 175952bfb9bbSJeremy L Thompson // Reduce to tridiagonal 176052bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n-1; i++) { 176152bfb9bbSJeremy L Thompson // Calculate Householder vector, magnitude 176252bfb9bbSJeremy L Thompson CeedScalar sigma = 0.0; 1763673160d7Sjeremylt v[i] = mat_T[i+n*(i+1)]; 176452bfb9bbSJeremy L Thompson for (CeedInt j=i+1; j<n-1; j++) { 1765673160d7Sjeremylt v[j] = mat_T[i+n*(j+1)]; 176652bfb9bbSJeremy L Thompson sigma += v[j] * v[j]; 176752bfb9bbSJeremy L Thompson } 176852bfb9bbSJeremy L Thompson CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:n-1] 1769673160d7Sjeremylt CeedScalar R_ii = -copysign(norm, v[i]); 1770673160d7Sjeremylt v[i] -= R_ii; 177152bfb9bbSJeremy L Thompson // norm of v[i:m] after modification above and scaling below 177252bfb9bbSJeremy L Thompson // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 177352bfb9bbSJeremy L Thompson // tau = 2 / (norm*norm) 177480a9ef05SNatalie Beams tau[i] = i == n - 2 ? 2 : 2 * v[i]*v[i] / (v[i]*v[i] + sigma); 1775fb551037Sjeremylt for (CeedInt j=i+1; j<n-1; j++) 1776fb551037Sjeremylt v[j] /= v[i]; 177752bfb9bbSJeremy L Thompson 177852bfb9bbSJeremy L Thompson // Update sub and super diagonal 177952bfb9bbSJeremy L Thompson for (CeedInt j=i+2; j<n; j++) { 1780673160d7Sjeremylt mat_T[i+n*j] = 0; mat_T[j+n*i] = 0; 178152bfb9bbSJeremy L Thompson } 178252bfb9bbSJeremy L Thompson // Apply symmetric Householder reflector to lower right panel 1783673160d7Sjeremylt CeedHouseholderReflect(&mat_T[(i+1)+n*(i+1)], &v[i], tau[i], 178452bfb9bbSJeremy L Thompson n-(i+1), n-(i+1), n, 1); 1785673160d7Sjeremylt CeedHouseholderReflect(&mat_T[(i+1)+n*(i+1)], &v[i], tau[i], 178652bfb9bbSJeremy L Thompson n-(i+1), n-(i+1), 1, n); 1787673160d7Sjeremylt 178852bfb9bbSJeremy L Thompson // Save v 1789673160d7Sjeremylt mat_T[i+n*(i+1)] = R_ii; 1790673160d7Sjeremylt mat_T[(i+1)+n*i] = R_ii; 179152bfb9bbSJeremy L Thompson for (CeedInt j=i+1; j<n-1; j++) { 1792673160d7Sjeremylt mat_T[i+n*(j+1)] = v[j]; 179352bfb9bbSJeremy L Thompson } 179452bfb9bbSJeremy L Thompson } 179552bfb9bbSJeremy L Thompson // Backwards accumulation of Q 179652bfb9bbSJeremy L Thompson for (CeedInt i=n-2; i>=0; i--) { 179785cf89eaSjeremylt if (tau[i] > 0.0) { 179852bfb9bbSJeremy L Thompson v[i] = 1; 179952bfb9bbSJeremy L Thompson for (CeedInt j=i+1; j<n-1; j++) { 1800673160d7Sjeremylt v[j] = mat_T[i+n*(j+1)]; 1801673160d7Sjeremylt mat_T[i+n*(j+1)] = 0; 180252bfb9bbSJeremy L Thompson } 180352bfb9bbSJeremy L Thompson CeedHouseholderReflect(&mat[(i+1)+n*(i+1)], &v[i], tau[i], 180452bfb9bbSJeremy L Thompson n-(i+1), n-(i+1), n, 1); 180552bfb9bbSJeremy L Thompson } 180685cf89eaSjeremylt } 180752bfb9bbSJeremy L Thompson 180852bfb9bbSJeremy L Thompson // Reduce sub and super diagonal 1809673160d7Sjeremylt CeedInt p = 0, q = 0, itr = 0, max_itr = n*n*n*n; 1810673160d7Sjeremylt CeedScalar tol = CEED_EPSILON; 181152bfb9bbSJeremy L Thompson 1812673160d7Sjeremylt while (itr < max_itr) { 181352bfb9bbSJeremy L Thompson // Update p, q, size of reduced portions of diagonal 181452bfb9bbSJeremy L Thompson p = 0; q = 0; 181552bfb9bbSJeremy L Thompson for (CeedInt i=n-2; i>=0; i--) { 1816673160d7Sjeremylt if (fabs(mat_T[i+n*(i+1)]) < tol) 181752bfb9bbSJeremy L Thompson q += 1; 181852bfb9bbSJeremy L Thompson else 181952bfb9bbSJeremy L Thompson break; 182052bfb9bbSJeremy L Thompson } 1821673160d7Sjeremylt for (CeedInt i=0; i<n-q-1; i++) { 1822673160d7Sjeremylt if (fabs(mat_T[i+n*(i+1)]) < tol) 182352bfb9bbSJeremy L Thompson p += 1; 182452bfb9bbSJeremy L Thompson else 182552bfb9bbSJeremy L Thompson break; 182652bfb9bbSJeremy L Thompson } 182752bfb9bbSJeremy L Thompson if (q == n-1) break; // Finished reducing 182852bfb9bbSJeremy L Thompson 182952bfb9bbSJeremy L Thompson // Reduce tridiagonal portion 1830673160d7Sjeremylt CeedScalar t_nn = mat_T[(n-1-q)+n*(n-1-q)], 1831673160d7Sjeremylt t_nnm1 = mat_T[(n-2-q)+n*(n-1-q)]; 1832673160d7Sjeremylt CeedScalar d = (mat_T[(n-2-q)+n*(n-2-q)] - t_nn)/2; 1833673160d7Sjeremylt CeedScalar mu = t_nn - t_nnm1*t_nnm1 / 1834673160d7Sjeremylt (d + copysign(sqrt(d*d + t_nnm1*t_nnm1), d)); 1835673160d7Sjeremylt CeedScalar x = mat_T[p+n*p] - mu; 1836673160d7Sjeremylt CeedScalar z = mat_T[p+n*(p+1)]; 1837673160d7Sjeremylt for (CeedInt k=p; k<n-q-1; k++) { 183852bfb9bbSJeremy L Thompson // Compute Givens rotation 183952bfb9bbSJeremy L Thompson CeedScalar c = 1, s = 0; 184052bfb9bbSJeremy L Thompson if (fabs(z) > tol) { 184152bfb9bbSJeremy L Thompson if (fabs(z) > fabs(x)) { 184252bfb9bbSJeremy L Thompson CeedScalar tau = -x/z; 184352bfb9bbSJeremy L Thompson s = 1/sqrt(1+tau*tau), c = s*tau; 184452bfb9bbSJeremy L Thompson } else { 184552bfb9bbSJeremy L Thompson CeedScalar tau = -z/x; 184652bfb9bbSJeremy L Thompson c = 1/sqrt(1+tau*tau), s = c*tau; 184752bfb9bbSJeremy L Thompson } 184852bfb9bbSJeremy L Thompson } 184952bfb9bbSJeremy L Thompson 185052bfb9bbSJeremy L Thompson // Apply Givens rotation to T 1851673160d7Sjeremylt CeedGivensRotation(mat_T, c, s, CEED_NOTRANSPOSE, k, k+1, n, n); 1852673160d7Sjeremylt CeedGivensRotation(mat_T, c, s, CEED_TRANSPOSE, k, k+1, n, n); 185352bfb9bbSJeremy L Thompson 185452bfb9bbSJeremy L Thompson // Apply Givens rotation to Q 185552bfb9bbSJeremy L Thompson CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k+1, n, n); 185652bfb9bbSJeremy L Thompson 185752bfb9bbSJeremy L Thompson // Update x, z 185852bfb9bbSJeremy L Thompson if (k < n-q-2) { 1859673160d7Sjeremylt x = mat_T[k+n*(k+1)]; 1860673160d7Sjeremylt z = mat_T[k+n*(k+2)]; 186152bfb9bbSJeremy L Thompson } 186252bfb9bbSJeremy L Thompson } 186352bfb9bbSJeremy L Thompson itr++; 186452bfb9bbSJeremy L Thompson } 1865673160d7Sjeremylt 186652bfb9bbSJeremy L Thompson // Save eigenvalues 186752bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) 1868673160d7Sjeremylt lambda[i] = mat_T[i+n*i]; 186952bfb9bbSJeremy L Thompson 187052bfb9bbSJeremy L Thompson // Check convergence 1871673160d7Sjeremylt if (itr == max_itr && q < n-1) 1872c042f62fSJeremy L Thompson // LCOV_EXCL_START 1873e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_MINOR, 1874e15f9bd0SJeremy L Thompson "Symmetric QR failed to converge"); 1875c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 1876e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 187752bfb9bbSJeremy L Thompson } 187803d18186Sjeremylt CeedPragmaOptimizeOn 187952bfb9bbSJeremy L Thompson 188052bfb9bbSJeremy L Thompson /** 188152bfb9bbSJeremy L Thompson @brief Return Simultaneous Diagonalization of two matrices. This solves the 188252bfb9bbSJeremy L Thompson generalized eigenvalue problem A x = lambda B x, where A and B 188352bfb9bbSJeremy L Thompson are symmetric and B is positive definite. We generate the matrix X 188452bfb9bbSJeremy L Thompson and vector Lambda such that X^T A X = Lambda and X^T B X = I. This 188552bfb9bbSJeremy L Thompson is equivalent to the LAPACK routine 'sygv' with TYPE = 1. 188652bfb9bbSJeremy L Thompson 188777645d7bSjeremylt @param ceed A Ceed context for error handling 1888d1d35e2fSjeremylt @param[in] mat_A Row-major matrix to be factorized with eigenvalues 1889d1d35e2fSjeremylt @param[in] mat_B Row-major matrix to be factorized to identity 1890d3331725Sjeremylt @param[out] mat_X Row-major orthogonal matrix 1891460bf743SValeria Barra @param[out] lambda Vector of length n of generalized eigenvalues 189252bfb9bbSJeremy L Thompson @param n Number of rows/columns 189352bfb9bbSJeremy L Thompson 189452bfb9bbSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 189552bfb9bbSJeremy L Thompson 189652bfb9bbSJeremy L Thompson @ref Utility 189752bfb9bbSJeremy L Thompson **/ 189803d18186Sjeremylt CeedPragmaOptimizeOff 1899d1d35e2fSjeremylt int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *mat_A, 1900d3331725Sjeremylt CeedScalar *mat_B, CeedScalar *mat_X, 190152bfb9bbSJeremy L Thompson CeedScalar *lambda, CeedInt n) { 190252bfb9bbSJeremy L Thompson int ierr; 1903d3331725Sjeremylt CeedScalar *mat_C, *mat_G, *vec_D; 190478464608Sjeremylt ierr = CeedCalloc(n*n, &mat_C); CeedChk(ierr); 190578464608Sjeremylt ierr = CeedCalloc(n*n, &mat_G); CeedChk(ierr); 1906d3331725Sjeremylt ierr = CeedCalloc(n, &vec_D); CeedChk(ierr); 190752bfb9bbSJeremy L Thompson 190852bfb9bbSJeremy L Thompson // Compute B = G D G^T 190978464608Sjeremylt memcpy(mat_G, mat_B, n*n*sizeof(mat_B[0])); 1910d3331725Sjeremylt ierr = CeedSymmetricSchurDecomposition(ceed, mat_G, vec_D, n); CeedChk(ierr); 191152bfb9bbSJeremy L Thompson 191285cf89eaSjeremylt // Sort eigenvalues 191385cf89eaSjeremylt for (CeedInt i=n-1; i>=0; i--) 191485cf89eaSjeremylt for (CeedInt j=0; j<i; j++) { 191585cf89eaSjeremylt if (fabs(vec_D[j]) > fabs(vec_D[j+1])) { 191685cf89eaSjeremylt CeedScalar temp; 191785cf89eaSjeremylt temp = vec_D[j]; vec_D[j] = vec_D[j+1]; vec_D[j+1] = temp; 191885cf89eaSjeremylt for (CeedInt k=0; k<n; k++) { 191985cf89eaSjeremylt temp = mat_G[k*n+j]; mat_G[k*n+j] = mat_G[k*n+j+1]; mat_G[k*n+j+1] = temp; 192085cf89eaSjeremylt } 192185cf89eaSjeremylt } 192285cf89eaSjeremylt } 192385cf89eaSjeremylt 1924fb551037Sjeremylt // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T 1925fb551037Sjeremylt // = D^-1/2 G^T A G D^-1/2 1926d3331725Sjeremylt // -- D = D^-1/2 192752bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) 1928d3331725Sjeremylt vec_D[i] = 1./sqrt(vec_D[i]); 1929d3331725Sjeremylt // -- G = G D^-1/2 1930d3331725Sjeremylt // -- C = D^-1/2 G^T 1931d3331725Sjeremylt for (CeedInt i=0; i<n; i++) 1932d3331725Sjeremylt for (CeedInt j=0; j<n; j++) { 1933673160d7Sjeremylt mat_G[i*n+j] *= vec_D[j]; 1934673160d7Sjeremylt mat_C[j*n+i] = mat_G[i*n+j]; 1935d3331725Sjeremylt } 1936673160d7Sjeremylt // -- X = (D^-1/2 G^T) A 1937ed9e99e6SJeremy L Thompson ierr = CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_C, 1938d3331725Sjeremylt (const CeedScalar *)mat_A, mat_X, n, n, n); 19399289e5bfSjeremylt CeedChk(ierr); 1940673160d7Sjeremylt // -- C = (D^-1/2 G^T A) (G D^-1/2) 1941ed9e99e6SJeremy L Thompson ierr = CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_X, 194278464608Sjeremylt (const CeedScalar *)mat_G, mat_C, n, n, n); 19439289e5bfSjeremylt CeedChk(ierr); 194452bfb9bbSJeremy L Thompson 194552bfb9bbSJeremy L Thompson // Compute Q^T C Q = lambda 1946d1d35e2fSjeremylt ierr = CeedSymmetricSchurDecomposition(ceed, mat_C, lambda, n); CeedChk(ierr); 194752bfb9bbSJeremy L Thompson 194885cf89eaSjeremylt // Sort eigenvalues 194985cf89eaSjeremylt for (CeedInt i=n-1; i>=0; i--) 195085cf89eaSjeremylt for (CeedInt j=0; j<i; j++) { 195185cf89eaSjeremylt if (fabs(lambda[j]) > fabs(lambda[j+1])) { 195285cf89eaSjeremylt CeedScalar temp; 195385cf89eaSjeremylt temp = lambda[j]; lambda[j] = lambda[j+1]; lambda[j+1] = temp; 195485cf89eaSjeremylt for (CeedInt k=0; k<n; k++) { 195585cf89eaSjeremylt temp = mat_C[k*n+j]; mat_C[k*n+j] = mat_C[k*n+j+1]; mat_C[k*n+j+1] = temp; 195685cf89eaSjeremylt } 195785cf89eaSjeremylt } 195885cf89eaSjeremylt } 195985cf89eaSjeremylt 1960d3331725Sjeremylt // Set X = (G D^1/2)^-T Q 1961fb551037Sjeremylt // = G D^-1/2 Q 1962ed9e99e6SJeremy L Thompson ierr = CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_G, 1963d3331725Sjeremylt (const CeedScalar *)mat_C, mat_X, n, n, n); 19649289e5bfSjeremylt CeedChk(ierr); 196578464608Sjeremylt 196678464608Sjeremylt // Cleanup 196778464608Sjeremylt ierr = CeedFree(&mat_C); CeedChk(ierr); 196878464608Sjeremylt ierr = CeedFree(&mat_G); CeedChk(ierr); 1969d3331725Sjeremylt ierr = CeedFree(&vec_D); CeedChk(ierr); 1970e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 197152bfb9bbSJeremy L Thompson } 197203d18186Sjeremylt CeedPragmaOptimizeOn 197352bfb9bbSJeremy L Thompson 1974d7b241e6Sjeremylt /// @} 1975