1d7b241e6Sjeremylt // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 2d7b241e6Sjeremylt // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 3d7b241e6Sjeremylt // reserved. See files LICENSE and NOTICE for details. 4d7b241e6Sjeremylt // 5d7b241e6Sjeremylt // This file is part of CEED, a collection of benchmarks, miniapps, software 6d7b241e6Sjeremylt // libraries and APIs for efficient high-order finite element and spectral 7d7b241e6Sjeremylt // element discretizations for exascale applications. For more information and 8d7b241e6Sjeremylt // source code availability see http://github.com/ceed. 9d7b241e6Sjeremylt // 10d7b241e6Sjeremylt // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 11d7b241e6Sjeremylt // a collaborative effort of two U.S. Department of Energy organizations (Office 12d7b241e6Sjeremylt // of Science and the National Nuclear Security Administration) responsible for 13d7b241e6Sjeremylt // the planning and preparation of a capable exascale ecosystem, including 14d7b241e6Sjeremylt // software, applications, hardware, advanced system engineering and early 15d7b241e6Sjeremylt // testbed platforms, in support of the nation's exascale computing imperative. 16d7b241e6Sjeremylt 17ec3da8bcSJed Brown #include <ceed/ceed.h> 18ec3da8bcSJed Brown #include <ceed/backend.h> 193d576824SJeremy L Thompson #include <ceed-impl.h> 20d7b241e6Sjeremylt #include <math.h> 213d576824SJeremy L Thompson #include <stdbool.h> 22d7b241e6Sjeremylt #include <stdio.h> 23d7b241e6Sjeremylt #include <string.h> 24d7b241e6Sjeremylt 257a982d89SJeremy L. Thompson /// @file 267a982d89SJeremy L. Thompson /// Implementation of CeedBasis interfaces 277a982d89SJeremy L. Thompson 28d7b241e6Sjeremylt /// @cond DOXYGEN_SKIP 29783c99b3SValeria Barra static struct CeedBasis_private ceed_basis_collocated; 30d7b241e6Sjeremylt /// @endcond 31d7b241e6Sjeremylt 327a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser 337a982d89SJeremy L. Thompson /// @{ 347a982d89SJeremy L. Thompson 357a982d89SJeremy L. Thompson /// Indicate that the quadrature points are collocated with the nodes 367a982d89SJeremy L. Thompson const CeedBasis CEED_BASIS_COLLOCATED = &ceed_basis_collocated; 377a982d89SJeremy L. Thompson 387a982d89SJeremy L. Thompson /// @} 397a982d89SJeremy L. Thompson 407a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 417a982d89SJeremy L. Thompson /// CeedBasis Library Internal Functions 427a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 437a982d89SJeremy L. Thompson /// @addtogroup CeedBasisDeveloper 447a982d89SJeremy L. Thompson /// @{ 457a982d89SJeremy L. Thompson 467a982d89SJeremy L. Thompson /** 477a982d89SJeremy L. Thompson @brief Compute Householder reflection 487a982d89SJeremy L. Thompson 497a982d89SJeremy L. Thompson Computes A = (I - b v v^T) A 507a982d89SJeremy L. Thompson where A is an mxn matrix indexed as A[i*row + j*col] 517a982d89SJeremy L. Thompson 527a982d89SJeremy L. Thompson @param[in,out] A Matrix to apply Householder reflection to, in place 537a982d89SJeremy L. Thompson @param v Householder vector 547a982d89SJeremy L. Thompson @param b Scaling factor 557a982d89SJeremy L. Thompson @param m Number of rows in A 567a982d89SJeremy L. Thompson @param n Number of columns in A 577a982d89SJeremy L. Thompson @param row Row stride 587a982d89SJeremy L. Thompson @param col Col stride 597a982d89SJeremy L. Thompson 607a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 617a982d89SJeremy L. Thompson 627a982d89SJeremy L. Thompson @ref Developer 637a982d89SJeremy L. Thompson **/ 647a982d89SJeremy L. Thompson static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v, 657a982d89SJeremy L. Thompson CeedScalar b, CeedInt m, CeedInt n, 667a982d89SJeremy L. Thompson CeedInt row, CeedInt col) { 677a982d89SJeremy L. Thompson for (CeedInt j=0; j<n; j++) { 687a982d89SJeremy L. Thompson CeedScalar w = A[0*row + j*col]; 697a982d89SJeremy L. Thompson for (CeedInt i=1; i<m; i++) 707a982d89SJeremy L. Thompson w += v[i] * A[i*row + j*col]; 717a982d89SJeremy L. Thompson A[0*row + j*col] -= b * w; 727a982d89SJeremy L. Thompson for (CeedInt i=1; i<m; i++) 737a982d89SJeremy L. Thompson A[i*row + j*col] -= b * w * v[i]; 747a982d89SJeremy L. Thompson } 75e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 767a982d89SJeremy L. Thompson } 777a982d89SJeremy L. Thompson 787a982d89SJeremy L. Thompson /** 797a982d89SJeremy L. Thompson @brief Apply Householder Q matrix 807a982d89SJeremy L. Thompson 817a982d89SJeremy L. Thompson Compute A = Q A where Q is mxm and A is mxn. 827a982d89SJeremy L. Thompson 837a982d89SJeremy L. Thompson @param[in,out] A Matrix to apply Householder Q to, in place 847a982d89SJeremy L. Thompson @param Q Householder Q matrix 857a982d89SJeremy L. Thompson @param tau Householder scaling factors 86d1d35e2fSjeremylt @param t_mode Transpose mode for application 877a982d89SJeremy L. Thompson @param m Number of rows in A 887a982d89SJeremy L. Thompson @param n Number of columns in A 897a982d89SJeremy L. Thompson @param k Number of elementary reflectors in Q, k<m 907a982d89SJeremy L. Thompson @param row Row stride in A 917a982d89SJeremy L. Thompson @param col Col stride in A 927a982d89SJeremy L. Thompson 937a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 947a982d89SJeremy L. Thompson 957a982d89SJeremy L. Thompson @ref Developer 967a982d89SJeremy L. Thompson **/ 97d99fa3c5SJeremy L Thompson int CeedHouseholderApplyQ(CeedScalar *A, const CeedScalar *Q, 98d1d35e2fSjeremylt const CeedScalar *tau, CeedTransposeMode t_mode, 997a982d89SJeremy L. Thompson CeedInt m, CeedInt n, CeedInt k, 1007a982d89SJeremy L. Thompson CeedInt row, CeedInt col) { 101e15f9bd0SJeremy L Thompson int ierr; 10278464608Sjeremylt CeedScalar *v; 10378464608Sjeremylt ierr = CeedMalloc(m, &v); CeedChk(ierr); 1047a982d89SJeremy L. Thompson for (CeedInt ii=0; ii<k; ii++) { 105d1d35e2fSjeremylt CeedInt i = t_mode == CEED_TRANSPOSE ? ii : k-1-ii; 1067a982d89SJeremy L. Thompson for (CeedInt j=i+1; j<m; j++) 1077a982d89SJeremy L. Thompson v[j] = Q[j*k+i]; 108d1d35e2fSjeremylt // Apply Householder reflector (I - tau v v^T) collo_grad_1d^T 109e15f9bd0SJeremy L Thompson ierr = CeedHouseholderReflect(&A[i*row], &v[i], tau[i], m-i, n, row, col); 110e15f9bd0SJeremy L Thompson CeedChk(ierr); 1117a982d89SJeremy L. Thompson } 11278464608Sjeremylt ierr = CeedFree(&v); CeedChk(ierr); 113e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1147a982d89SJeremy L. Thompson } 1157a982d89SJeremy L. Thompson 1167a982d89SJeremy L. Thompson /** 1177a982d89SJeremy L. Thompson @brief Compute Givens rotation 1187a982d89SJeremy L. Thompson 1197a982d89SJeremy L. Thompson Computes A = G A (or G^T A in transpose mode) 1207a982d89SJeremy L. Thompson where A is an mxn matrix indexed as A[i*n + j*m] 1217a982d89SJeremy L. Thompson 1227a982d89SJeremy L. Thompson @param[in,out] A Row major matrix to apply Givens rotation to, in place 1237a982d89SJeremy L. Thompson @param c Cosine factor 1247a982d89SJeremy L. Thompson @param s Sine factor 125d1d35e2fSjeremylt @param t_mode @ref CEED_NOTRANSPOSE to rotate the basis counter-clockwise, 1264c4400c7SValeria Barra which has the effect of rotating columns of A clockwise; 1274cc79fe7SJed Brown @ref CEED_TRANSPOSE for the opposite rotation 1287a982d89SJeremy L. Thompson @param i First row/column to apply rotation 1297a982d89SJeremy L. Thompson @param k Second row/column to apply rotation 1307a982d89SJeremy L. Thompson @param m Number of rows in A 1317a982d89SJeremy L. Thompson @param n Number of columns in A 1327a982d89SJeremy L. Thompson 1337a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 1347a982d89SJeremy L. Thompson 1357a982d89SJeremy L. Thompson @ref Developer 1367a982d89SJeremy L. Thompson **/ 1377a982d89SJeremy L. Thompson static int CeedGivensRotation(CeedScalar *A, CeedScalar c, CeedScalar s, 138d1d35e2fSjeremylt CeedTransposeMode t_mode, CeedInt i, CeedInt k, 1397a982d89SJeremy L. Thompson CeedInt m, CeedInt n) { 140d1d35e2fSjeremylt CeedInt stride_j = 1, stride_ik = m, num_its = n; 141d1d35e2fSjeremylt if (t_mode == CEED_NOTRANSPOSE) { 142d1d35e2fSjeremylt stride_j = n; stride_ik = 1; num_its = m; 1437a982d89SJeremy L. Thompson } 1447a982d89SJeremy L. Thompson 1457a982d89SJeremy L. Thompson // Apply rotation 146d1d35e2fSjeremylt for (CeedInt j=0; j<num_its; j++) { 147d1d35e2fSjeremylt CeedScalar tau1 = A[i*stride_ik+j*stride_j], tau2 = A[k*stride_ik+j*stride_j]; 148d1d35e2fSjeremylt A[i*stride_ik+j*stride_j] = c*tau1 - s*tau2; 149d1d35e2fSjeremylt A[k*stride_ik+j*stride_j] = s*tau1 + c*tau2; 1507a982d89SJeremy L. Thompson } 151e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1527a982d89SJeremy L. Thompson } 1537a982d89SJeremy L. Thompson 1547a982d89SJeremy L. Thompson /** 1557a982d89SJeremy L. Thompson @brief View an array stored in a CeedBasis 1567a982d89SJeremy L. Thompson 1570a0da059Sjeremylt @param[in] name Name of array 158d1d35e2fSjeremylt @param[in] fp_fmt Printing format 1590a0da059Sjeremylt @param[in] m Number of rows in array 1600a0da059Sjeremylt @param[in] n Number of columns in array 1610a0da059Sjeremylt @param[in] a Array to be viewed 1620a0da059Sjeremylt @param[in] stream Stream to view to, e.g., stdout 1637a982d89SJeremy L. Thompson 1647a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 1657a982d89SJeremy L. Thompson 1667a982d89SJeremy L. Thompson @ref Developer 1677a982d89SJeremy L. Thompson **/ 168d1d35e2fSjeremylt static int CeedScalarView(const char *name, const char *fp_fmt, CeedInt m, 1697a982d89SJeremy L. Thompson CeedInt n, const CeedScalar *a, FILE *stream) { 1707a982d89SJeremy L. Thompson for (int i=0; i<m; i++) { 1717a982d89SJeremy L. Thompson if (m > 1) 1727a982d89SJeremy L. Thompson fprintf(stream, "%12s[%d]:", name, i); 1737a982d89SJeremy L. Thompson else 1747a982d89SJeremy L. Thompson fprintf(stream, "%12s:", name); 1757a982d89SJeremy L. Thompson for (int j=0; j<n; j++) 176d1d35e2fSjeremylt fprintf(stream, fp_fmt, fabs(a[i*n+j]) > 1E-14 ? a[i*n+j] : 0); 1777a982d89SJeremy L. Thompson fputs("\n", stream); 1787a982d89SJeremy L. Thompson } 179e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1807a982d89SJeremy L. Thompson } 1817a982d89SJeremy L. Thompson 1827a982d89SJeremy L. Thompson /// @} 1837a982d89SJeremy L. Thompson 1847a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 1857a982d89SJeremy L. Thompson /// Ceed Backend API 1867a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 1877a982d89SJeremy L. Thompson /// @addtogroup CeedBasisBackend 1887a982d89SJeremy L. Thompson /// @{ 1897a982d89SJeremy L. Thompson 1907a982d89SJeremy L. Thompson /** 1917a982d89SJeremy L. Thompson @brief Return collocated grad matrix 1927a982d89SJeremy L. Thompson 1937a982d89SJeremy L. Thompson @param basis CeedBasis 194d1d35e2fSjeremylt @param[out] collo_grad_1d Row-major (Q_1d * Q_1d) matrix expressing derivatives of 1957a982d89SJeremy L. Thompson basis functions at quadrature points 1967a982d89SJeremy L. Thompson 1977a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 1987a982d89SJeremy L. Thompson 1997a982d89SJeremy L. Thompson @ref Backend 2007a982d89SJeremy L. Thompson **/ 201d1d35e2fSjeremylt int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *collo_grad_1d) { 2027a982d89SJeremy L. Thompson int i, j, k; 2037a982d89SJeremy L. Thompson Ceed ceed; 204d1d35e2fSjeremylt CeedInt ierr, P_1d=(basis)->P_1d, Q_1d=(basis)->Q_1d; 20578464608Sjeremylt CeedScalar *interp_1d, *grad_1d, *tau; 2067a982d89SJeremy L. Thompson 207d1d35e2fSjeremylt ierr = CeedMalloc(Q_1d*P_1d, &interp_1d); CeedChk(ierr); 208d1d35e2fSjeremylt ierr = CeedMalloc(Q_1d*P_1d, &grad_1d); CeedChk(ierr); 20978464608Sjeremylt ierr = CeedMalloc(Q_1d, &tau); CeedChk(ierr); 210d1d35e2fSjeremylt memcpy(interp_1d, (basis)->interp_1d, Q_1d*P_1d*sizeof(basis)->interp_1d[0]); 211d1d35e2fSjeremylt memcpy(grad_1d, (basis)->grad_1d, Q_1d*P_1d*sizeof(basis)->interp_1d[0]); 2127a982d89SJeremy L. Thompson 213d1d35e2fSjeremylt // QR Factorization, interp_1d = Q R 2147a982d89SJeremy L. Thompson ierr = CeedBasisGetCeed(basis, &ceed); CeedChk(ierr); 215d1d35e2fSjeremylt ierr = CeedQRFactorization(ceed, interp_1d, tau, Q_1d, P_1d); CeedChk(ierr); 216e15f9bd0SJeremy L Thompson // Note: This function is for backend use, so all errors are terminal 217e15f9bd0SJeremy L Thompson // and we do not need to clean up memory on failure. 2187a982d89SJeremy L. Thompson 219d1d35e2fSjeremylt // Apply Rinv, collo_grad_1d = grad_1d Rinv 220d1d35e2fSjeremylt for (i=0; i<Q_1d; i++) { // Row i 221d1d35e2fSjeremylt collo_grad_1d[Q_1d*i] = grad_1d[P_1d*i]/interp_1d[0]; 222d1d35e2fSjeremylt for (j=1; j<P_1d; j++) { // Column j 223d1d35e2fSjeremylt collo_grad_1d[j+Q_1d*i] = grad_1d[j+P_1d*i]; 2247a982d89SJeremy L. Thompson for (k=0; k<j; k++) 225d1d35e2fSjeremylt collo_grad_1d[j+Q_1d*i] -= interp_1d[j+P_1d*k]*collo_grad_1d[k+Q_1d*i]; 226d1d35e2fSjeremylt collo_grad_1d[j+Q_1d*i] /= interp_1d[j+P_1d*j]; 2277a982d89SJeremy L. Thompson } 228d1d35e2fSjeremylt for (j=P_1d; j<Q_1d; j++) 229d1d35e2fSjeremylt collo_grad_1d[j+Q_1d*i] = 0; 2307a982d89SJeremy L. Thompson } 2317a982d89SJeremy L. Thompson 232673160d7Sjeremylt // Apply Qtranspose, collo_grad = collo_grad Q_transpose 233d1d35e2fSjeremylt ierr = CeedHouseholderApplyQ(collo_grad_1d, interp_1d, tau, CEED_NOTRANSPOSE, 234d1d35e2fSjeremylt Q_1d, Q_1d, P_1d, 1, Q_1d); CeedChk(ierr); 2357a982d89SJeremy L. Thompson 236d1d35e2fSjeremylt ierr = CeedFree(&interp_1d); CeedChk(ierr); 237d1d35e2fSjeremylt ierr = CeedFree(&grad_1d); CeedChk(ierr); 23878464608Sjeremylt ierr = CeedFree(&tau); CeedChk(ierr); 239e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 2407a982d89SJeremy L. Thompson } 2417a982d89SJeremy L. Thompson 2427a982d89SJeremy L. Thompson /** 2437a982d89SJeremy L. Thompson @brief Get tensor status for given CeedBasis 2447a982d89SJeremy L. Thompson 2457a982d89SJeremy L. Thompson @param basis CeedBasis 246d1d35e2fSjeremylt @param[out] is_tensor Variable to store tensor status 2477a982d89SJeremy L. Thompson 2487a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 2497a982d89SJeremy L. Thompson 2507a982d89SJeremy L. Thompson @ref Backend 2517a982d89SJeremy L. Thompson **/ 252d1d35e2fSjeremylt int CeedBasisIsTensor(CeedBasis basis, bool *is_tensor) { 253d1d35e2fSjeremylt *is_tensor = basis->tensor_basis; 254e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 2557a982d89SJeremy L. Thompson } 2567a982d89SJeremy L. Thompson 2577a982d89SJeremy L. Thompson /** 2587a982d89SJeremy L. Thompson @brief Get backend data of a CeedBasis 2597a982d89SJeremy L. Thompson 2607a982d89SJeremy L. Thompson @param basis CeedBasis 2617a982d89SJeremy L. Thompson @param[out] data Variable to store data 2627a982d89SJeremy L. Thompson 2637a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 2647a982d89SJeremy L. Thompson 2657a982d89SJeremy L. Thompson @ref Backend 2667a982d89SJeremy L. Thompson **/ 267777ff853SJeremy L Thompson int CeedBasisGetData(CeedBasis basis, void *data) { 268777ff853SJeremy L Thompson *(void **)data = basis->data; 269e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 2707a982d89SJeremy L. Thompson } 2717a982d89SJeremy L. Thompson 2727a982d89SJeremy L. Thompson /** 2737a982d89SJeremy L. Thompson @brief Set backend data of a CeedBasis 2747a982d89SJeremy L. Thompson 2757a982d89SJeremy L. Thompson @param[out] basis CeedBasis 2767a982d89SJeremy L. Thompson @param data Data to set 2777a982d89SJeremy L. Thompson 2787a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 2797a982d89SJeremy L. Thompson 2807a982d89SJeremy L. Thompson @ref Backend 2817a982d89SJeremy L. Thompson **/ 282777ff853SJeremy L Thompson int CeedBasisSetData(CeedBasis basis, void *data) { 283777ff853SJeremy L Thompson basis->data = data; 284e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 2857a982d89SJeremy L. Thompson } 2867a982d89SJeremy L. Thompson 2877a982d89SJeremy L. Thompson /** 28834359f16Sjeremylt @brief Increment the reference counter for a CeedBasis 28934359f16Sjeremylt 29034359f16Sjeremylt @param basis Basis to increment the reference counter 29134359f16Sjeremylt 29234359f16Sjeremylt @return An error code: 0 - success, otherwise - failure 29334359f16Sjeremylt 29434359f16Sjeremylt @ref Backend 29534359f16Sjeremylt **/ 2969560d06aSjeremylt int CeedBasisReference(CeedBasis basis) { 29734359f16Sjeremylt basis->ref_count++; 29834359f16Sjeremylt return CEED_ERROR_SUCCESS; 29934359f16Sjeremylt } 30034359f16Sjeremylt 30134359f16Sjeremylt /** 3027a982d89SJeremy L. Thompson @brief Get dimension for given CeedElemTopology 3037a982d89SJeremy L. Thompson 3047a982d89SJeremy L. Thompson @param topo CeedElemTopology 3057a982d89SJeremy L. Thompson @param[out] dim Variable to store dimension of topology 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 **/ 3117a982d89SJeremy L. Thompson int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) { 3127a982d89SJeremy L. Thompson *dim = (CeedInt) topo >> 16; 313e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3147a982d89SJeremy L. Thompson } 3157a982d89SJeremy L. Thompson 3167a982d89SJeremy L. Thompson /** 3177a982d89SJeremy L. Thompson @brief Get CeedTensorContract of a CeedBasis 3187a982d89SJeremy L. Thompson 3197a982d89SJeremy L. Thompson @param basis CeedBasis 3207a982d89SJeremy L. Thompson @param[out] contract Variable to store CeedTensorContract 3217a982d89SJeremy L. Thompson 3227a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3237a982d89SJeremy L. Thompson 3247a982d89SJeremy L. Thompson @ref Backend 3257a982d89SJeremy L. Thompson **/ 3267a982d89SJeremy L. Thompson int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) { 3277a982d89SJeremy L. Thompson *contract = basis->contract; 328e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3297a982d89SJeremy L. Thompson } 3307a982d89SJeremy L. Thompson 3317a982d89SJeremy L. Thompson /** 3327a982d89SJeremy L. Thompson @brief Set CeedTensorContract of a CeedBasis 3337a982d89SJeremy L. Thompson 3347a982d89SJeremy L. Thompson @param[out] basis CeedBasis 3357a982d89SJeremy L. Thompson @param contract CeedTensorContract to set 3367a982d89SJeremy L. Thompson 3377a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3387a982d89SJeremy L. Thompson 3397a982d89SJeremy L. Thompson @ref Backend 3407a982d89SJeremy L. Thompson **/ 34134359f16Sjeremylt int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract contract) { 34234359f16Sjeremylt int ierr; 34334359f16Sjeremylt basis->contract = contract; 3449560d06aSjeremylt ierr = CeedTensorContractReference(contract); CeedChk(ierr); 345e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3467a982d89SJeremy L. Thompson } 3477a982d89SJeremy L. Thompson 3487a982d89SJeremy L. Thompson /** 3497a982d89SJeremy L. Thompson @brief Return a reference implementation of matrix multiplication C = A B. 3507a982d89SJeremy L. Thompson Note, this is a reference implementation for CPU CeedScalar pointers 3517a982d89SJeremy L. Thompson that is not intended for high performance. 3527a982d89SJeremy L. Thompson 3537a982d89SJeremy L. Thompson @param ceed A Ceed context for error handling 354d1d35e2fSjeremylt @param[in] mat_A Row-major matrix A 355d1d35e2fSjeremylt @param[in] mat_B Row-major matrix B 356d1d35e2fSjeremylt @param[out] mat_C Row-major output matrix C 3577a982d89SJeremy L. Thompson @param m Number of rows of C 3587a982d89SJeremy L. Thompson @param n Number of columns of C 3597a982d89SJeremy L. Thompson @param kk Number of columns of A/rows of B 3607a982d89SJeremy L. Thompson 3617a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3627a982d89SJeremy L. Thompson 3637a982d89SJeremy L. Thompson @ref Utility 3647a982d89SJeremy L. Thompson **/ 365d1d35e2fSjeremylt int CeedMatrixMultiply(Ceed ceed, const CeedScalar *mat_A, 366d1d35e2fSjeremylt const CeedScalar *mat_B, CeedScalar *mat_C, CeedInt m, 3677a982d89SJeremy L. Thompson CeedInt n, CeedInt kk) { 3687a982d89SJeremy L. Thompson for (CeedInt i=0; i<m; i++) 3697a982d89SJeremy L. Thompson for (CeedInt j=0; j<n; j++) { 3707a982d89SJeremy L. Thompson CeedScalar sum = 0; 3717a982d89SJeremy L. Thompson for (CeedInt k=0; k<kk; k++) 372d1d35e2fSjeremylt sum += mat_A[k+i*kk]*mat_B[j+k*n]; 373d1d35e2fSjeremylt mat_C[j+i*n] = sum; 3747a982d89SJeremy L. Thompson } 375e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3767a982d89SJeremy L. Thompson } 3777a982d89SJeremy L. Thompson 3787a982d89SJeremy L. Thompson /// @} 3797a982d89SJeremy L. Thompson 3807a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 3817a982d89SJeremy L. Thompson /// CeedBasis Public API 3827a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 3837a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser 384d7b241e6Sjeremylt /// @{ 385d7b241e6Sjeremylt 386b11c1e72Sjeremylt /** 38795bb1877Svaleriabarra @brief Create a tensor-product basis for H^1 discretizations 388b11c1e72Sjeremylt 389b11c1e72Sjeremylt @param ceed A Ceed object where the CeedBasis will be created 390b11c1e72Sjeremylt @param dim Topological dimension 391d1d35e2fSjeremylt @param num_comp Number of field components (1 for scalar fields) 392d1d35e2fSjeremylt @param P_1d Number of nodes in one dimension 393d1d35e2fSjeremylt @param Q_1d Number of quadrature points in one dimension 394d1d35e2fSjeremylt @param interp_1d Row-major (Q_1d * P_1d) matrix expressing the values of nodal 395b11c1e72Sjeremylt basis functions at quadrature points 396d1d35e2fSjeremylt @param grad_1d Row-major (Q_1d * P_1d) matrix expressing derivatives of nodal 397b11c1e72Sjeremylt basis functions at quadrature points 398d1d35e2fSjeremylt @param q_ref_1d Array of length Q_1d holding the locations of quadrature points 399b11c1e72Sjeremylt on the 1D reference element [-1, 1] 400d1d35e2fSjeremylt @param q_weight_1d Array of length Q_1d holding the quadrature weights on the 401b11c1e72Sjeremylt reference element 402b11c1e72Sjeremylt @param[out] basis Address of the variable where the newly created 403b11c1e72Sjeremylt CeedBasis will be stored. 404b11c1e72Sjeremylt 405b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 406dfdf5a53Sjeremylt 4077a982d89SJeremy L. Thompson @ref User 408b11c1e72Sjeremylt **/ 409d1d35e2fSjeremylt int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt num_comp, 410d1d35e2fSjeremylt CeedInt P_1d, CeedInt Q_1d, 411d1d35e2fSjeremylt const CeedScalar *interp_1d, 412d1d35e2fSjeremylt const CeedScalar *grad_1d, const CeedScalar *q_ref_1d, 413d1d35e2fSjeremylt const CeedScalar *q_weight_1d, CeedBasis *basis) { 414d7b241e6Sjeremylt int ierr; 415d7b241e6Sjeremylt 4165fe0d4faSjeremylt if (!ceed->BasisCreateTensorH1) { 4175fe0d4faSjeremylt Ceed delegate; 418aefd8378Sjeremylt ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr); 4195fe0d4faSjeremylt 4205fe0d4faSjeremylt if (!delegate) 421c042f62fSJeremy L Thompson // LCOV_EXCL_START 422e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_UNSUPPORTED, 423e15f9bd0SJeremy L Thompson "Backend does not support BasisCreateTensorH1"); 424c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 4255fe0d4faSjeremylt 426d1d35e2fSjeremylt ierr = CeedBasisCreateTensorH1(delegate, dim, num_comp, P_1d, 427d1d35e2fSjeremylt Q_1d, interp_1d, grad_1d, q_ref_1d, 428d1d35e2fSjeremylt q_weight_1d, basis); CeedChk(ierr); 429e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 4305fe0d4faSjeremylt } 431e15f9bd0SJeremy L Thompson 432e15f9bd0SJeremy L Thompson if (dim<1) 433e15f9bd0SJeremy L Thompson // LCOV_EXCL_START 434e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, 435e15f9bd0SJeremy L Thompson "Basis dimension must be a positive value"); 436e15f9bd0SJeremy L Thompson // LCOV_EXCL_STOP 43750c301a5SRezgar Shakeri CeedElemTopology topo = dim == 1 ? CEED_TOPOLOGY_LINE 43850c301a5SRezgar Shakeri : dim == 2 ? CEED_TOPOLOGY_QUAD 43950c301a5SRezgar Shakeri : CEED_TOPOLOGY_HEX; 440e15f9bd0SJeremy L Thompson 441d7b241e6Sjeremylt ierr = CeedCalloc(1, basis); CeedChk(ierr); 442d7b241e6Sjeremylt (*basis)->ceed = ceed; 4439560d06aSjeremylt ierr = CeedReference(ceed); CeedChk(ierr); 444d1d35e2fSjeremylt (*basis)->ref_count = 1; 445d1d35e2fSjeremylt (*basis)->tensor_basis = 1; 446d7b241e6Sjeremylt (*basis)->dim = dim; 447d99fa3c5SJeremy L Thompson (*basis)->topo = topo; 448d1d35e2fSjeremylt (*basis)->num_comp = num_comp; 449d1d35e2fSjeremylt (*basis)->P_1d = P_1d; 450d1d35e2fSjeremylt (*basis)->Q_1d = Q_1d; 451d1d35e2fSjeremylt (*basis)->P = CeedIntPow(P_1d, dim); 452d1d35e2fSjeremylt (*basis)->Q = CeedIntPow(Q_1d, dim); 45350c301a5SRezgar Shakeri (*basis)->Q_comp = 1; 45450c301a5SRezgar Shakeri (*basis)->basis_space = 1; // 1 for H^1 space 455ff3a0f91SJeremy L Thompson ierr = CeedCalloc(Q_1d, &(*basis)->q_ref_1d); CeedChk(ierr); 456ff3a0f91SJeremy L Thompson ierr = CeedCalloc(Q_1d, &(*basis)->q_weight_1d); CeedChk(ierr); 457ff3a0f91SJeremy L Thompson if (q_ref_1d) memcpy((*basis)->q_ref_1d, q_ref_1d, Q_1d*sizeof(q_ref_1d[0])); 458ff3a0f91SJeremy L Thompson if (q_weight_1d) memcpy((*basis)->q_weight_1d, q_weight_1d, 459ff3a0f91SJeremy L Thompson Q_1d*sizeof(q_weight_1d[0])); 460ff3a0f91SJeremy L Thompson ierr = CeedCalloc(Q_1d*P_1d, &(*basis)->interp_1d); CeedChk(ierr); 461ff3a0f91SJeremy L Thompson ierr = CeedCalloc(Q_1d*P_1d, &(*basis)->grad_1d); CeedChk(ierr); 462ff3a0f91SJeremy L Thompson if (interp_1d) memcpy((*basis)->interp_1d, interp_1d, 463ff3a0f91SJeremy L Thompson Q_1d*P_1d*sizeof(interp_1d[0])); 464ff3a0f91SJeremy L Thompson if (grad_1d) memcpy((*basis)->grad_1d, grad_1d, Q_1d*P_1d*sizeof(grad_1d[0])); 465d1d35e2fSjeremylt ierr = ceed->BasisCreateTensorH1(dim, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, 466d1d35e2fSjeremylt q_weight_1d, *basis); CeedChk(ierr); 467e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 468d7b241e6Sjeremylt } 469d7b241e6Sjeremylt 470b11c1e72Sjeremylt /** 47195bb1877Svaleriabarra @brief Create a tensor-product Lagrange basis 472b11c1e72Sjeremylt 473b11c1e72Sjeremylt @param ceed A Ceed object where the CeedBasis will be created 474b11c1e72Sjeremylt @param dim Topological dimension of element 475d1d35e2fSjeremylt @param num_comp Number of field components (1 for scalar fields) 476b11c1e72Sjeremylt @param P Number of Gauss-Lobatto nodes in one dimension. The 477b11c1e72Sjeremylt polynomial degree of the resulting Q_k element is k=P-1. 478b11c1e72Sjeremylt @param Q Number of quadrature points in one dimension. 479d1d35e2fSjeremylt @param quad_mode Distribution of the Q quadrature points (affects order of 480b11c1e72Sjeremylt accuracy for the quadrature) 481b11c1e72Sjeremylt @param[out] basis Address of the variable where the newly created 482b11c1e72Sjeremylt CeedBasis will be stored. 483b11c1e72Sjeremylt 484b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 485dfdf5a53Sjeremylt 4867a982d89SJeremy L. Thompson @ref User 487b11c1e72Sjeremylt **/ 488d1d35e2fSjeremylt int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt num_comp, 489d1d35e2fSjeremylt CeedInt P, CeedInt Q, CeedQuadMode quad_mode, 490692c2638Sjeremylt CeedBasis *basis) { 491d7b241e6Sjeremylt // Allocate 492e15f9bd0SJeremy L Thompson int ierr, ierr2, i, j, k; 493d1d35e2fSjeremylt CeedScalar c1, c2, c3, c4, dx, *nodes, *interp_1d, *grad_1d, *q_ref_1d, 494d1d35e2fSjeremylt *q_weight_1d; 4954d537eeaSYohann 4964d537eeaSYohann if (dim<1) 497c042f62fSJeremy L Thompson // LCOV_EXCL_START 498e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, 499e15f9bd0SJeremy L Thompson "Basis dimension must be a positive value"); 500c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 5014d537eeaSYohann 502e15f9bd0SJeremy L Thompson // Get Nodes and Weights 503d1d35e2fSjeremylt ierr = CeedCalloc(P*Q, &interp_1d); CeedChk(ierr); 504d1d35e2fSjeremylt ierr = CeedCalloc(P*Q, &grad_1d); CeedChk(ierr); 505d7b241e6Sjeremylt ierr = CeedCalloc(P, &nodes); CeedChk(ierr); 506d1d35e2fSjeremylt ierr = CeedCalloc(Q, &q_ref_1d); CeedChk(ierr); 507d1d35e2fSjeremylt ierr = CeedCalloc(Q, &q_weight_1d); CeedChk(ierr); 508e15f9bd0SJeremy L Thompson ierr = CeedLobattoQuadrature(P, nodes, NULL); 509e15f9bd0SJeremy L Thompson if (ierr) { goto cleanup; } CeedChk(ierr); 510d1d35e2fSjeremylt switch (quad_mode) { 511d7b241e6Sjeremylt case CEED_GAUSS: 512d1d35e2fSjeremylt ierr = CeedGaussQuadrature(Q, q_ref_1d, q_weight_1d); 513d7b241e6Sjeremylt break; 514d7b241e6Sjeremylt case CEED_GAUSS_LOBATTO: 515d1d35e2fSjeremylt ierr = CeedLobattoQuadrature(Q, q_ref_1d, q_weight_1d); 516d7b241e6Sjeremylt break; 517d7b241e6Sjeremylt } 518e15f9bd0SJeremy L Thompson if (ierr) { goto cleanup; } CeedChk(ierr); 519e15f9bd0SJeremy L Thompson 520d7b241e6Sjeremylt // Build B, D matrix 521d7b241e6Sjeremylt // Fornberg, 1998 522d7b241e6Sjeremylt for (i = 0; i < Q; i++) { 523d7b241e6Sjeremylt c1 = 1.0; 524d1d35e2fSjeremylt c3 = nodes[0] - q_ref_1d[i]; 525d1d35e2fSjeremylt interp_1d[i*P+0] = 1.0; 526d7b241e6Sjeremylt for (j = 1; j < P; j++) { 527d7b241e6Sjeremylt c2 = 1.0; 528d7b241e6Sjeremylt c4 = c3; 529d1d35e2fSjeremylt c3 = nodes[j] - q_ref_1d[i]; 530d7b241e6Sjeremylt for (k = 0; k < j; k++) { 531d7b241e6Sjeremylt dx = nodes[j] - nodes[k]; 532d7b241e6Sjeremylt c2 *= dx; 533d7b241e6Sjeremylt if (k == j - 1) { 534d1d35e2fSjeremylt grad_1d[i*P + j] = c1*(interp_1d[i*P + k] - c4*grad_1d[i*P + k]) / c2; 535d1d35e2fSjeremylt interp_1d[i*P + j] = - c1*c4*interp_1d[i*P + k] / c2; 536d7b241e6Sjeremylt } 537d1d35e2fSjeremylt grad_1d[i*P + k] = (c3*grad_1d[i*P + k] - interp_1d[i*P + k]) / dx; 538d1d35e2fSjeremylt interp_1d[i*P + k] = c3*interp_1d[i*P + k] / dx; 539d7b241e6Sjeremylt } 540d7b241e6Sjeremylt c1 = c2; 541d7b241e6Sjeremylt } 542d7b241e6Sjeremylt } 5439ac7b42eSJeremy L Thompson // Pass to CeedBasisCreateTensorH1 544d1d35e2fSjeremylt ierr = CeedBasisCreateTensorH1(ceed, dim, num_comp, P, Q, interp_1d, grad_1d, 5459ac7b42eSJeremy L Thompson q_ref_1d, q_weight_1d, basis); CeedChk(ierr); 546e15f9bd0SJeremy L Thompson cleanup: 547d1d35e2fSjeremylt ierr2 = CeedFree(&interp_1d); CeedChk(ierr2); 548d1d35e2fSjeremylt ierr2 = CeedFree(&grad_1d); CeedChk(ierr2); 549e15f9bd0SJeremy L Thompson ierr2 = CeedFree(&nodes); CeedChk(ierr2); 550d1d35e2fSjeremylt ierr2 = CeedFree(&q_ref_1d); CeedChk(ierr2); 551d1d35e2fSjeremylt ierr2 = CeedFree(&q_weight_1d); CeedChk(ierr2); 552e15f9bd0SJeremy L Thompson CeedChk(ierr); 553e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 554d7b241e6Sjeremylt } 555d7b241e6Sjeremylt 556b11c1e72Sjeremylt /** 55795bb1877Svaleriabarra @brief Create a non tensor-product basis for H^1 discretizations 558a8de75f0Sjeremylt 559a8de75f0Sjeremylt @param ceed A Ceed object where the CeedBasis will be created 560a8de75f0Sjeremylt @param topo Topology of element, e.g. hypercube, simplex, ect 561d1d35e2fSjeremylt @param num_comp Number of field components (1 for scalar fields) 562d1d35e2fSjeremylt @param num_nodes Total number of nodes 563d1d35e2fSjeremylt @param num_qpts Total number of quadrature points 564d1d35e2fSjeremylt @param interp Row-major (num_qpts * num_nodes) matrix expressing the values of 5658795c945Sjeremylt nodal basis functions at quadrature points 566d1d35e2fSjeremylt @param grad Row-major (num_qpts * dim * num_nodes) matrix expressing 5678795c945Sjeremylt derivatives of nodal basis functions at quadrature points 568d1d35e2fSjeremylt @param q_ref Array of length num_qpts holding the locations of quadrature 56950c301a5SRezgar Shakeri points on the reference element 570d1d35e2fSjeremylt @param q_weight Array of length num_qpts holding the quadrature weights on the 571a8de75f0Sjeremylt reference element 572a8de75f0Sjeremylt @param[out] basis Address of the variable where the newly created 573a8de75f0Sjeremylt CeedBasis will be stored. 574a8de75f0Sjeremylt 575a8de75f0Sjeremylt @return An error code: 0 - success, otherwise - failure 576a8de75f0Sjeremylt 5777a982d89SJeremy L. Thompson @ref User 578a8de75f0Sjeremylt **/ 579d1d35e2fSjeremylt int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, 580d1d35e2fSjeremylt CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 581d1d35e2fSjeremylt const CeedScalar *grad, const CeedScalar *q_ref, 582d1d35e2fSjeremylt const CeedScalar *q_weight, CeedBasis *basis) { 583a8de75f0Sjeremylt int ierr; 584d1d35e2fSjeremylt CeedInt P = num_nodes, Q = num_qpts, dim = 0; 585a8de75f0Sjeremylt 5865fe0d4faSjeremylt if (!ceed->BasisCreateH1) { 5875fe0d4faSjeremylt Ceed delegate; 588aefd8378Sjeremylt ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr); 5895fe0d4faSjeremylt 5905fe0d4faSjeremylt if (!delegate) 591c042f62fSJeremy L Thompson // LCOV_EXCL_START 592e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_UNSUPPORTED, 593e15f9bd0SJeremy L Thompson "Backend does not support BasisCreateH1"); 594c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 5955fe0d4faSjeremylt 596d1d35e2fSjeremylt ierr = CeedBasisCreateH1(delegate, topo, num_comp, num_nodes, 597d1d35e2fSjeremylt num_qpts, interp, grad, q_ref, 598d1d35e2fSjeremylt q_weight, basis); CeedChk(ierr); 599e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 6005fe0d4faSjeremylt } 6015fe0d4faSjeremylt 602a8de75f0Sjeremylt ierr = CeedCalloc(1, basis); CeedChk(ierr); 603a8de75f0Sjeremylt 604a8de75f0Sjeremylt ierr = CeedBasisGetTopologyDimension(topo, &dim); CeedChk(ierr); 605a8de75f0Sjeremylt 606a8de75f0Sjeremylt (*basis)->ceed = ceed; 6079560d06aSjeremylt ierr = CeedReference(ceed); CeedChk(ierr); 608d1d35e2fSjeremylt (*basis)->ref_count = 1; 609d1d35e2fSjeremylt (*basis)->tensor_basis = 0; 610a8de75f0Sjeremylt (*basis)->dim = dim; 611d99fa3c5SJeremy L Thompson (*basis)->topo = topo; 612d1d35e2fSjeremylt (*basis)->num_comp = num_comp; 613a8de75f0Sjeremylt (*basis)->P = P; 614a8de75f0Sjeremylt (*basis)->Q = Q; 61550c301a5SRezgar Shakeri (*basis)->Q_comp = 1; 61650c301a5SRezgar Shakeri (*basis)->basis_space = 1; // 1 for H^1 space 617ff3a0f91SJeremy L Thompson ierr = CeedCalloc(Q*dim, &(*basis)->q_ref_1d); CeedChk(ierr); 618ff3a0f91SJeremy L Thompson ierr = CeedCalloc(Q, &(*basis)->q_weight_1d); CeedChk(ierr); 619ff3a0f91SJeremy L Thompson if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q*dim*sizeof(q_ref[0])); 620ff3a0f91SJeremy L Thompson if(q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q*sizeof(q_weight[0])); 621ff3a0f91SJeremy L Thompson ierr = CeedCalloc(Q*P, &(*basis)->interp); CeedChk(ierr); 622ff3a0f91SJeremy L Thompson ierr = CeedCalloc(dim*Q*P, &(*basis)->grad); CeedChk(ierr); 623ff3a0f91SJeremy L Thompson if(interp) memcpy((*basis)->interp, interp, Q*P*sizeof(interp[0])); 624ff3a0f91SJeremy L Thompson if(grad) memcpy((*basis)->grad, grad, dim*Q*P*sizeof(grad[0])); 625d1d35e2fSjeremylt ierr = ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, q_ref, 626d1d35e2fSjeremylt q_weight, *basis); CeedChk(ierr); 627e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 628a8de75f0Sjeremylt } 629a8de75f0Sjeremylt 630a8de75f0Sjeremylt /** 63150c301a5SRezgar Shakeri @brief Create a non tensor-product basis for H(div) discretizations 63250c301a5SRezgar Shakeri 63350c301a5SRezgar Shakeri @param ceed A Ceed object where the CeedBasis will be created 63450c301a5SRezgar Shakeri @param topo Topology of element (`CEED_TOPOLOGY_QUAD`, `CEED_TOPOLOGY_PRISM`, etc.), 63550c301a5SRezgar Shakeri dimension of which is used in some array sizes below 63650c301a5SRezgar Shakeri @param num_comp Number of components (usually 1 for vectors in H(div) bases) 63750c301a5SRezgar Shakeri @param num_nodes Total number of nodes (dofs per element) 63850c301a5SRezgar Shakeri @param num_qpts Total number of quadrature points 63950c301a5SRezgar Shakeri @param interp Row-major (dim*num_qpts * num_nodes) matrix expressing the values of 64050c301a5SRezgar Shakeri nodal basis functions at quadrature points 64150c301a5SRezgar Shakeri @param div Row-major (num_qpts * num_nodes) matrix expressing 64250c301a5SRezgar Shakeri divergence of nodal basis functions at quadrature points 64350c301a5SRezgar Shakeri @param q_ref Array of length num_qpts holding the locations of quadrature 64450c301a5SRezgar Shakeri points on the reference element 64550c301a5SRezgar Shakeri @param q_weight Array of length num_qpts holding the quadrature weights on the 64650c301a5SRezgar Shakeri reference element 64750c301a5SRezgar Shakeri @param[out] basis Address of the variable where the newly created 64850c301a5SRezgar Shakeri CeedBasis will be stored. 64950c301a5SRezgar Shakeri 65050c301a5SRezgar Shakeri @return An error code: 0 - success, otherwise - failure 65150c301a5SRezgar Shakeri 65250c301a5SRezgar Shakeri @ref User 65350c301a5SRezgar Shakeri **/ 65450c301a5SRezgar Shakeri int CeedBasisCreateHdiv(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, 65550c301a5SRezgar Shakeri CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 65650c301a5SRezgar Shakeri const CeedScalar *div, const CeedScalar *q_ref, 65750c301a5SRezgar Shakeri const CeedScalar *q_weight, CeedBasis *basis) { 65850c301a5SRezgar Shakeri int ierr; 65950c301a5SRezgar Shakeri CeedInt Q = num_qpts, P = num_nodes, dim = 0; 66050c301a5SRezgar Shakeri ierr = CeedBasisGetTopologyDimension(topo, &dim); CeedChk(ierr); 66150c301a5SRezgar Shakeri if (!ceed->BasisCreateHdiv) { 66250c301a5SRezgar Shakeri Ceed delegate; 66350c301a5SRezgar Shakeri ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr); 66450c301a5SRezgar Shakeri 66550c301a5SRezgar Shakeri if (!delegate) 66650c301a5SRezgar Shakeri // LCOV_EXCL_START 66750c301a5SRezgar Shakeri return CeedError(ceed, CEED_ERROR_UNSUPPORTED, 66850c301a5SRezgar Shakeri "Backend does not implement BasisCreateHdiv"); 66950c301a5SRezgar Shakeri // LCOV_EXCL_STOP 67050c301a5SRezgar Shakeri 67150c301a5SRezgar Shakeri ierr = CeedBasisCreateHdiv(delegate, topo, num_comp, num_nodes, 67250c301a5SRezgar Shakeri num_qpts, interp, div, q_ref, 67350c301a5SRezgar Shakeri q_weight, basis); CeedChk(ierr); 67450c301a5SRezgar Shakeri return CEED_ERROR_SUCCESS; 67550c301a5SRezgar Shakeri } 67650c301a5SRezgar Shakeri 67750c301a5SRezgar Shakeri ierr = CeedCalloc(1, basis); CeedChk(ierr); 67850c301a5SRezgar Shakeri 67950c301a5SRezgar Shakeri (*basis)->ceed = ceed; 68050c301a5SRezgar Shakeri ierr = CeedReference(ceed); CeedChk(ierr); 68150c301a5SRezgar Shakeri (*basis)->ref_count = 1; 68250c301a5SRezgar Shakeri (*basis)->tensor_basis = 0; 68350c301a5SRezgar Shakeri (*basis)->dim = dim; 68450c301a5SRezgar Shakeri (*basis)->topo = topo; 68550c301a5SRezgar Shakeri (*basis)->num_comp = num_comp; 68650c301a5SRezgar Shakeri (*basis)->P = P; 68750c301a5SRezgar Shakeri (*basis)->Q = Q; 68850c301a5SRezgar Shakeri (*basis)->Q_comp = dim; 68950c301a5SRezgar Shakeri (*basis)->basis_space = 2; // 2 for H(div) space 69050c301a5SRezgar Shakeri ierr = CeedMalloc(Q*dim, &(*basis)->q_ref_1d); CeedChk(ierr); 69150c301a5SRezgar Shakeri ierr = CeedMalloc(Q, &(*basis)->q_weight_1d); CeedChk(ierr); 69250c301a5SRezgar Shakeri if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q*dim*sizeof(q_ref[0])); 69350c301a5SRezgar Shakeri if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q*sizeof(q_weight[0])); 69450c301a5SRezgar Shakeri ierr = CeedMalloc(dim*Q*P, &(*basis)->interp); CeedChk(ierr); 69550c301a5SRezgar Shakeri ierr = CeedMalloc(Q*P, &(*basis)->div); CeedChk(ierr); 69650c301a5SRezgar Shakeri if (interp) memcpy((*basis)->interp, interp, dim*Q*P*sizeof(interp[0])); 69750c301a5SRezgar Shakeri if (div) memcpy((*basis)->div, div, Q*P*sizeof(div[0])); 69850c301a5SRezgar Shakeri ierr = ceed->BasisCreateHdiv(topo, dim, P, Q, interp, div, q_ref, 69950c301a5SRezgar Shakeri q_weight, *basis); CeedChk(ierr); 70050c301a5SRezgar Shakeri return CEED_ERROR_SUCCESS; 70150c301a5SRezgar Shakeri } 70250c301a5SRezgar Shakeri 70350c301a5SRezgar Shakeri /** 7049560d06aSjeremylt @brief Copy the pointer to a CeedBasis. Both pointers should 7059560d06aSjeremylt be destroyed with `CeedBasisDestroy()`; 7069560d06aSjeremylt Note: If `*basis_copy` is non-NULL, then it is assumed that 7079560d06aSjeremylt `*basis_copy` is a pointer to a CeedBasis. This CeedBasis 7089560d06aSjeremylt will be destroyed if `*basis_copy` is the only 7099560d06aSjeremylt reference to this CeedBasis. 7109560d06aSjeremylt 7119560d06aSjeremylt @param basis CeedBasis to copy reference to 7129560d06aSjeremylt @param[out] basis_copy Variable to store copied reference 7139560d06aSjeremylt 7149560d06aSjeremylt @return An error code: 0 - success, otherwise - failure 7159560d06aSjeremylt 7169560d06aSjeremylt @ref User 7179560d06aSjeremylt **/ 7189560d06aSjeremylt int CeedBasisReferenceCopy(CeedBasis basis, CeedBasis *basis_copy) { 7199560d06aSjeremylt int ierr; 7209560d06aSjeremylt 7219560d06aSjeremylt ierr = CeedBasisReference(basis); CeedChk(ierr); 7229560d06aSjeremylt ierr = CeedBasisDestroy(basis_copy); CeedChk(ierr); 7239560d06aSjeremylt *basis_copy = basis; 7249560d06aSjeremylt return CEED_ERROR_SUCCESS; 7259560d06aSjeremylt } 7269560d06aSjeremylt 7279560d06aSjeremylt /** 7287a982d89SJeremy L. Thompson @brief View a CeedBasis 7297a982d89SJeremy L. Thompson 7307a982d89SJeremy L. Thompson @param basis CeedBasis to view 7317a982d89SJeremy L. Thompson @param stream Stream to view to, e.g., stdout 7327a982d89SJeremy L. Thompson 7337a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 7347a982d89SJeremy L. Thompson 7357a982d89SJeremy L. Thompson @ref User 7367a982d89SJeremy L. Thompson **/ 7377a982d89SJeremy L. Thompson int CeedBasisView(CeedBasis basis, FILE *stream) { 7387a982d89SJeremy L. Thompson int ierr; 73950c301a5SRezgar Shakeri CeedFESpace FE_space = basis->basis_space; 74050c301a5SRezgar Shakeri CeedElemTopology topo = basis->topo; 74150c301a5SRezgar Shakeri // Print FE space and element topology of the basis 742d1d35e2fSjeremylt if (basis->tensor_basis) { 74350c301a5SRezgar Shakeri fprintf(stream, "CeedBasis (%s on a %s element): dim=%d P=%d Q=%d\n", 74450c301a5SRezgar Shakeri CeedFESpaces[FE_space], CeedElemTopologies[topo], 74550c301a5SRezgar Shakeri basis->dim, basis->P_1d, basis->Q_1d); 74650c301a5SRezgar Shakeri } else { 74750c301a5SRezgar Shakeri fprintf(stream, "CeedBasis (%s on a %s element): dim=%d P=%d Q=%d\n", 74850c301a5SRezgar Shakeri CeedFESpaces[FE_space], CeedElemTopologies[topo], 74950c301a5SRezgar Shakeri basis->dim, basis->P, basis->Q); 75050c301a5SRezgar Shakeri } 75150c301a5SRezgar Shakeri // Print quadrature data, interpolation/gradient/divergene/curl of the basis 75250c301a5SRezgar Shakeri if (basis->tensor_basis) { // tensor basis 753d1d35e2fSjeremylt ierr = CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_ref_1d, 7547a982d89SJeremy L. Thompson stream); CeedChk(ierr); 755d1d35e2fSjeremylt ierr = CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q_1d, 756d1d35e2fSjeremylt basis->q_weight_1d, stream); CeedChk(ierr); 757d1d35e2fSjeremylt ierr = CeedScalarView("interp1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, 758d1d35e2fSjeremylt basis->interp_1d, stream); CeedChk(ierr); 759d1d35e2fSjeremylt ierr = CeedScalarView("grad1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, 760d1d35e2fSjeremylt basis->grad_1d, stream); CeedChk(ierr); 76150c301a5SRezgar Shakeri } else { // non-tensor basis 7627a982d89SJeremy L. Thompson ierr = CeedScalarView("qref", "\t% 12.8f", 1, basis->Q*basis->dim, 763d1d35e2fSjeremylt basis->q_ref_1d, 7647a982d89SJeremy L. Thompson stream); CeedChk(ierr); 765d1d35e2fSjeremylt ierr = CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->q_weight_1d, 7667a982d89SJeremy L. Thompson stream); CeedChk(ierr); 76750c301a5SRezgar Shakeri ierr = CeedScalarView("interp", "\t% 12.8f", basis->Q_comp*basis->Q, basis->P, 7687a982d89SJeremy L. Thompson basis->interp, stream); CeedChk(ierr); 76950c301a5SRezgar Shakeri if (basis->grad) { 7707a982d89SJeremy L. Thompson ierr = CeedScalarView("grad", "\t% 12.8f", basis->dim*basis->Q, basis->P, 7717a982d89SJeremy L. Thompson basis->grad, stream); CeedChk(ierr); 7727a982d89SJeremy L. Thompson } 77350c301a5SRezgar Shakeri if (basis->div) { 77450c301a5SRezgar Shakeri ierr = CeedScalarView("div", "\t% 12.8f", basis->Q, basis->P, 77550c301a5SRezgar Shakeri basis->div, stream); CeedChk(ierr); 77650c301a5SRezgar Shakeri } 77750c301a5SRezgar Shakeri } 778e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 7797a982d89SJeremy L. Thompson } 7807a982d89SJeremy L. Thompson 7817a982d89SJeremy L. Thompson /** 7827a982d89SJeremy L. Thompson @brief Apply basis evaluation from nodes to quadrature points or vice versa 7837a982d89SJeremy L. Thompson 7847a982d89SJeremy L. Thompson @param basis CeedBasis to evaluate 785d1d35e2fSjeremylt @param num_elem The number of elements to apply the basis evaluation to; 7867a982d89SJeremy L. Thompson the backend will specify the ordering in 7874cc79fe7SJed Brown CeedElemRestrictionCreateBlocked() 788d1d35e2fSjeremylt @param t_mode \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature 7897a982d89SJeremy L. Thompson points, \ref CEED_TRANSPOSE to apply the transpose, mapping 7907a982d89SJeremy L. Thompson from quadrature points to nodes 791d1d35e2fSjeremylt @param eval_mode \ref CEED_EVAL_NONE to use values directly, 7927a982d89SJeremy L. Thompson \ref CEED_EVAL_INTERP to use interpolated values, 7937a982d89SJeremy L. Thompson \ref CEED_EVAL_GRAD to use gradients, 7947a982d89SJeremy L. Thompson \ref CEED_EVAL_WEIGHT to use quadrature weights. 7957a982d89SJeremy L. Thompson @param[in] u Input CeedVector 7967a982d89SJeremy L. Thompson @param[out] v Output CeedVector 7977a982d89SJeremy L. Thompson 7987a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 7997a982d89SJeremy L. Thompson 8007a982d89SJeremy L. Thompson @ref User 8017a982d89SJeremy L. Thompson **/ 802d1d35e2fSjeremylt int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode, 803d1d35e2fSjeremylt CeedEvalMode eval_mode, CeedVector u, CeedVector v) { 8047a982d89SJeremy L. Thompson int ierr; 805*1f9221feSJeremy L Thompson CeedSize u_length = 0, v_length; 806*1f9221feSJeremy L Thompson CeedInt dim, num_comp, num_nodes, num_qpts; 807e15f9bd0SJeremy L Thompson ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr); 808d1d35e2fSjeremylt ierr = CeedBasisGetNumComponents(basis, &num_comp); CeedChk(ierr); 809d1d35e2fSjeremylt ierr = CeedBasisGetNumNodes(basis, &num_nodes); CeedChk(ierr); 810d1d35e2fSjeremylt ierr = CeedBasisGetNumQuadraturePoints(basis, &num_qpts); CeedChk(ierr); 811d1d35e2fSjeremylt ierr = CeedVectorGetLength(v, &v_length); CeedChk(ierr); 8127a982d89SJeremy L. Thompson if (u) { 813d1d35e2fSjeremylt ierr = CeedVectorGetLength(u, &u_length); CeedChk(ierr); 8147a982d89SJeremy L. Thompson } 8157a982d89SJeremy L. Thompson 816e15f9bd0SJeremy L Thompson if (!basis->Apply) 817e15f9bd0SJeremy L Thompson // LCOV_EXCL_START 818e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_UNSUPPORTED, 819e15f9bd0SJeremy L Thompson "Backend does not support BasisApply"); 820e15f9bd0SJeremy L Thompson // LCOV_EXCL_STOP 821e15f9bd0SJeremy L Thompson 822e15f9bd0SJeremy L Thompson // Check compatibility of topological and geometrical dimensions 823d1d35e2fSjeremylt if ((t_mode == CEED_TRANSPOSE && (v_length%num_nodes != 0 || 824d1d35e2fSjeremylt u_length%num_qpts != 0)) || 825d1d35e2fSjeremylt (t_mode == CEED_NOTRANSPOSE && (u_length%num_nodes != 0 || 826d1d35e2fSjeremylt v_length%num_qpts != 0))) 8278229195eSjeremylt // LCOV_EXCL_START 828e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_DIMENSION, 829e15f9bd0SJeremy L Thompson "Length of input/output vectors " 8307a982d89SJeremy L. Thompson "incompatible with basis dimensions"); 8318229195eSjeremylt // LCOV_EXCL_STOP 8327a982d89SJeremy L. Thompson 833e15f9bd0SJeremy L Thompson // Check vector lengths to prevent out of bounds issues 834d1d35e2fSjeremylt bool bad_dims = false; 835d1d35e2fSjeremylt switch (eval_mode) { 836e15f9bd0SJeremy L Thompson case CEED_EVAL_NONE: 837d1d35e2fSjeremylt case CEED_EVAL_INTERP: bad_dims = 838d1d35e2fSjeremylt ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts || 839d1d35e2fSjeremylt v_length < num_elem*num_comp*num_nodes)) || 840d1d35e2fSjeremylt (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp || 841d1d35e2fSjeremylt u_length < num_elem*num_comp*num_nodes))); 842e15f9bd0SJeremy L Thompson break; 843d1d35e2fSjeremylt case CEED_EVAL_GRAD: bad_dims = 844d1d35e2fSjeremylt ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts*dim || 845d1d35e2fSjeremylt v_length < num_elem*num_comp*num_nodes)) || 846d1d35e2fSjeremylt (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp*dim || 847d1d35e2fSjeremylt u_length < num_elem*num_comp*num_nodes))); 848e15f9bd0SJeremy L Thompson break; 849e15f9bd0SJeremy L Thompson case CEED_EVAL_WEIGHT: 850d1d35e2fSjeremylt bad_dims = v_length < num_elem*num_qpts; 851e15f9bd0SJeremy L Thompson break; 852e15f9bd0SJeremy L Thompson // LCOV_EXCL_START 853d1d35e2fSjeremylt case CEED_EVAL_DIV: bad_dims = 854d1d35e2fSjeremylt ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts || 855d1d35e2fSjeremylt v_length < num_elem*num_comp*num_nodes)) || 856d1d35e2fSjeremylt (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp || 857d1d35e2fSjeremylt u_length < num_elem*num_comp*num_nodes))); 858e15f9bd0SJeremy L Thompson break; 859d1d35e2fSjeremylt case CEED_EVAL_CURL: bad_dims = 860d1d35e2fSjeremylt ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts || 861d1d35e2fSjeremylt v_length < num_elem*num_comp*num_nodes)) || 862d1d35e2fSjeremylt (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp || 863d1d35e2fSjeremylt u_length < num_elem*num_comp*num_nodes))); 864e15f9bd0SJeremy L Thompson break; 865e15f9bd0SJeremy L Thompson // LCOV_EXCL_STOP 866e15f9bd0SJeremy L Thompson } 867d1d35e2fSjeremylt if (bad_dims) 868e15f9bd0SJeremy L Thompson // LCOV_EXCL_START 869e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_DIMENSION, 870d1d35e2fSjeremylt "Input/output vectors too short for basis and evaluation mode"); 871e15f9bd0SJeremy L Thompson // LCOV_EXCL_STOP 872e15f9bd0SJeremy L Thompson 873d1d35e2fSjeremylt ierr = basis->Apply(basis, num_elem, t_mode, eval_mode, u, v); CeedChk(ierr); 874e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 8757a982d89SJeremy L. Thompson } 8767a982d89SJeremy L. Thompson 8777a982d89SJeremy L. Thompson /** 878b7c9bbdaSJeremy L Thompson @brief Get Ceed associated with a CeedBasis 879b7c9bbdaSJeremy L Thompson 880b7c9bbdaSJeremy L Thompson @param basis CeedBasis 881b7c9bbdaSJeremy L Thompson @param[out] ceed Variable to store Ceed 882b7c9bbdaSJeremy L Thompson 883b7c9bbdaSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 884b7c9bbdaSJeremy L Thompson 885b7c9bbdaSJeremy L Thompson @ref Advanced 886b7c9bbdaSJeremy L Thompson **/ 887b7c9bbdaSJeremy L Thompson int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) { 888b7c9bbdaSJeremy L Thompson *ceed = basis->ceed; 889b7c9bbdaSJeremy L Thompson return CEED_ERROR_SUCCESS; 890b7c9bbdaSJeremy L Thompson } 891b7c9bbdaSJeremy L Thompson 892b7c9bbdaSJeremy L Thompson /** 8939d007619Sjeremylt @brief Get dimension for given CeedBasis 8949d007619Sjeremylt 8959d007619Sjeremylt @param basis CeedBasis 8969d007619Sjeremylt @param[out] dim Variable to store dimension of basis 8979d007619Sjeremylt 8989d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 8999d007619Sjeremylt 900b7c9bbdaSJeremy L Thompson @ref Advanced 9019d007619Sjeremylt **/ 9029d007619Sjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) { 9039d007619Sjeremylt *dim = basis->dim; 904e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 9059d007619Sjeremylt } 9069d007619Sjeremylt 9079d007619Sjeremylt /** 908d99fa3c5SJeremy L Thompson @brief Get topology for given CeedBasis 909d99fa3c5SJeremy L Thompson 910d99fa3c5SJeremy L Thompson @param basis CeedBasis 911d99fa3c5SJeremy L Thompson @param[out] topo Variable to store topology of basis 912d99fa3c5SJeremy L Thompson 913d99fa3c5SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 914d99fa3c5SJeremy L Thompson 915b7c9bbdaSJeremy L Thompson @ref Advanced 916d99fa3c5SJeremy L Thompson **/ 917d99fa3c5SJeremy L Thompson int CeedBasisGetTopology(CeedBasis basis, CeedElemTopology *topo) { 918d99fa3c5SJeremy L Thompson *topo = basis->topo; 919e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 920d99fa3c5SJeremy L Thompson } 921d99fa3c5SJeremy L Thompson 922d99fa3c5SJeremy L Thompson /** 92350c301a5SRezgar Shakeri @brief Get number of Q-vector components for given CeedBasis 92450c301a5SRezgar Shakeri 92550c301a5SRezgar Shakeri @param basis CeedBasis 92650c301a5SRezgar Shakeri @param[out] Q_comp Variable to store number of Q-vector components of basis 92750c301a5SRezgar Shakeri 92850c301a5SRezgar Shakeri @return An error code: 0 - success, otherwise - failure 92950c301a5SRezgar Shakeri 93050c301a5SRezgar Shakeri @ref Advanced 93150c301a5SRezgar Shakeri **/ 93250c301a5SRezgar Shakeri int CeedBasisGetNumQuadratureComponents(CeedBasis basis, CeedInt *Q_comp) { 93350c301a5SRezgar Shakeri *Q_comp = basis->Q_comp; 93450c301a5SRezgar Shakeri return CEED_ERROR_SUCCESS; 93550c301a5SRezgar Shakeri } 93650c301a5SRezgar Shakeri 93750c301a5SRezgar Shakeri /** 9389d007619Sjeremylt @brief Get number of components for given CeedBasis 9399d007619Sjeremylt 9409d007619Sjeremylt @param basis CeedBasis 941d1d35e2fSjeremylt @param[out] num_comp Variable to store number of components of basis 9429d007619Sjeremylt 9439d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 9449d007619Sjeremylt 945b7c9bbdaSJeremy L Thompson @ref Advanced 9469d007619Sjeremylt **/ 947d1d35e2fSjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *num_comp) { 948d1d35e2fSjeremylt *num_comp = basis->num_comp; 949e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 9509d007619Sjeremylt } 9519d007619Sjeremylt 9529d007619Sjeremylt /** 9539d007619Sjeremylt @brief Get total number of nodes (in dim dimensions) of a CeedBasis 9549d007619Sjeremylt 9559d007619Sjeremylt @param basis CeedBasis 9569d007619Sjeremylt @param[out] P Variable to store number of nodes 9579d007619Sjeremylt 9589d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 9599d007619Sjeremylt 9609d007619Sjeremylt @ref Utility 9619d007619Sjeremylt **/ 9629d007619Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) { 9639d007619Sjeremylt *P = basis->P; 964e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 9659d007619Sjeremylt } 9669d007619Sjeremylt 9679d007619Sjeremylt /** 9689d007619Sjeremylt @brief Get total number of nodes (in 1 dimension) of a CeedBasis 9699d007619Sjeremylt 9709d007619Sjeremylt @param basis CeedBasis 971d1d35e2fSjeremylt @param[out] P_1d Variable to store number of nodes 9729d007619Sjeremylt 9739d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 9749d007619Sjeremylt 975b7c9bbdaSJeremy L Thompson @ref Advanced 9769d007619Sjeremylt **/ 977d1d35e2fSjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P_1d) { 978d1d35e2fSjeremylt if (!basis->tensor_basis) 9799d007619Sjeremylt // LCOV_EXCL_START 980e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_MINOR, 981d1d35e2fSjeremylt "Cannot supply P_1d for non-tensor basis"); 9829d007619Sjeremylt // LCOV_EXCL_STOP 9839d007619Sjeremylt 984d1d35e2fSjeremylt *P_1d = basis->P_1d; 985e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 9869d007619Sjeremylt } 9879d007619Sjeremylt 9889d007619Sjeremylt /** 9899d007619Sjeremylt @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis 9909d007619Sjeremylt 9919d007619Sjeremylt @param basis CeedBasis 9929d007619Sjeremylt @param[out] Q Variable to store number of quadrature points 9939d007619Sjeremylt 9949d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 9959d007619Sjeremylt 9969d007619Sjeremylt @ref Utility 9979d007619Sjeremylt **/ 9989d007619Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) { 9999d007619Sjeremylt *Q = basis->Q; 1000e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 10019d007619Sjeremylt } 10029d007619Sjeremylt 10039d007619Sjeremylt /** 10049d007619Sjeremylt @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis 10059d007619Sjeremylt 10069d007619Sjeremylt @param basis CeedBasis 1007d1d35e2fSjeremylt @param[out] Q_1d Variable to store number of quadrature points 10089d007619Sjeremylt 10099d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 10109d007619Sjeremylt 1011b7c9bbdaSJeremy L Thompson @ref Advanced 10129d007619Sjeremylt **/ 1013d1d35e2fSjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q_1d) { 1014d1d35e2fSjeremylt if (!basis->tensor_basis) 10159d007619Sjeremylt // LCOV_EXCL_START 1016e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_MINOR, 1017d1d35e2fSjeremylt "Cannot supply Q_1d for non-tensor basis"); 10189d007619Sjeremylt // LCOV_EXCL_STOP 10199d007619Sjeremylt 1020d1d35e2fSjeremylt *Q_1d = basis->Q_1d; 1021e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 10229d007619Sjeremylt } 10239d007619Sjeremylt 10249d007619Sjeremylt /** 10259d007619Sjeremylt @brief Get reference coordinates of quadrature points (in dim dimensions) 10269d007619Sjeremylt of a CeedBasis 10279d007619Sjeremylt 10289d007619Sjeremylt @param basis CeedBasis 1029d1d35e2fSjeremylt @param[out] q_ref Variable to store reference coordinates of quadrature points 10309d007619Sjeremylt 10319d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 10329d007619Sjeremylt 1033b7c9bbdaSJeremy L Thompson @ref Advanced 10349d007619Sjeremylt **/ 1035d1d35e2fSjeremylt int CeedBasisGetQRef(CeedBasis basis, const CeedScalar **q_ref) { 1036d1d35e2fSjeremylt *q_ref = basis->q_ref_1d; 1037e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 10389d007619Sjeremylt } 10399d007619Sjeremylt 10409d007619Sjeremylt /** 10419d007619Sjeremylt @brief Get quadrature weights of quadrature points (in dim dimensions) 10429d007619Sjeremylt of a CeedBasis 10439d007619Sjeremylt 10449d007619Sjeremylt @param basis CeedBasis 1045d1d35e2fSjeremylt @param[out] q_weight Variable to store quadrature weights 10469d007619Sjeremylt 10479d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 10489d007619Sjeremylt 1049b7c9bbdaSJeremy L Thompson @ref Advanced 10509d007619Sjeremylt **/ 1051d1d35e2fSjeremylt int CeedBasisGetQWeights(CeedBasis basis, const CeedScalar **q_weight) { 1052d1d35e2fSjeremylt *q_weight = basis->q_weight_1d; 1053e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 10549d007619Sjeremylt } 10559d007619Sjeremylt 10569d007619Sjeremylt /** 10579d007619Sjeremylt @brief Get interpolation matrix of a CeedBasis 10589d007619Sjeremylt 10599d007619Sjeremylt @param basis CeedBasis 10609d007619Sjeremylt @param[out] interp Variable to store interpolation matrix 10619d007619Sjeremylt 10629d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 10639d007619Sjeremylt 1064b7c9bbdaSJeremy L Thompson @ref Advanced 10659d007619Sjeremylt **/ 10666c58de82SJeremy L Thompson int CeedBasisGetInterp(CeedBasis basis, const CeedScalar **interp) { 1067d1d35e2fSjeremylt if (!basis->interp && basis->tensor_basis) { 10689d007619Sjeremylt // Allocate 10699d007619Sjeremylt int ierr; 10709d007619Sjeremylt ierr = CeedMalloc(basis->Q*basis->P, &basis->interp); CeedChk(ierr); 10719d007619Sjeremylt 10729d007619Sjeremylt // Initialize 10739d007619Sjeremylt for (CeedInt i=0; i<basis->Q*basis->P; i++) 10749d007619Sjeremylt basis->interp[i] = 1.0; 10759d007619Sjeremylt 10769d007619Sjeremylt // Calculate 10779d007619Sjeremylt for (CeedInt d=0; d<basis->dim; d++) 10789d007619Sjeremylt for (CeedInt qpt=0; qpt<basis->Q; qpt++) 10799d007619Sjeremylt for (CeedInt node=0; node<basis->P; node++) { 1080d1d35e2fSjeremylt CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d; 1081d1d35e2fSjeremylt CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d; 1082d1d35e2fSjeremylt basis->interp[qpt*(basis->P)+node] *= basis->interp_1d[q*basis->P_1d+p]; 10839d007619Sjeremylt } 10849d007619Sjeremylt } 10859d007619Sjeremylt *interp = basis->interp; 1086e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 10879d007619Sjeremylt } 10889d007619Sjeremylt 10899d007619Sjeremylt /** 10909d007619Sjeremylt @brief Get 1D interpolation matrix of a tensor product CeedBasis 10919d007619Sjeremylt 10929d007619Sjeremylt @param basis CeedBasis 1093d1d35e2fSjeremylt @param[out] interp_1d Variable to store interpolation matrix 10949d007619Sjeremylt 10959d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 10969d007619Sjeremylt 10979d007619Sjeremylt @ref Backend 10989d007619Sjeremylt **/ 1099d1d35e2fSjeremylt int CeedBasisGetInterp1D(CeedBasis basis, const CeedScalar **interp_1d) { 1100d1d35e2fSjeremylt if (!basis->tensor_basis) 11019d007619Sjeremylt // LCOV_EXCL_START 1102e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_MINOR, 1103e15f9bd0SJeremy L Thompson "CeedBasis is not a tensor product basis."); 11049d007619Sjeremylt // LCOV_EXCL_STOP 11059d007619Sjeremylt 1106d1d35e2fSjeremylt *interp_1d = basis->interp_1d; 1107e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 11089d007619Sjeremylt } 11099d007619Sjeremylt 11109d007619Sjeremylt /** 11119d007619Sjeremylt @brief Get gradient matrix of a CeedBasis 11129d007619Sjeremylt 11139d007619Sjeremylt @param basis CeedBasis 11149d007619Sjeremylt @param[out] grad Variable to store gradient matrix 11159d007619Sjeremylt 11169d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 11179d007619Sjeremylt 1118b7c9bbdaSJeremy L Thompson @ref Advanced 11199d007619Sjeremylt **/ 11206c58de82SJeremy L Thompson int CeedBasisGetGrad(CeedBasis basis, const CeedScalar **grad) { 1121d1d35e2fSjeremylt if (!basis->grad && basis->tensor_basis) { 11229d007619Sjeremylt // Allocate 11239d007619Sjeremylt int ierr; 11249d007619Sjeremylt ierr = CeedMalloc(basis->dim*basis->Q*basis->P, &basis->grad); 11259d007619Sjeremylt CeedChk(ierr); 11269d007619Sjeremylt 11279d007619Sjeremylt // Initialize 11289d007619Sjeremylt for (CeedInt i=0; i<basis->dim*basis->Q*basis->P; i++) 11299d007619Sjeremylt basis->grad[i] = 1.0; 11309d007619Sjeremylt 11319d007619Sjeremylt // Calculate 11329d007619Sjeremylt for (CeedInt d=0; d<basis->dim; d++) 11339d007619Sjeremylt for (CeedInt i=0; i<basis->dim; i++) 11349d007619Sjeremylt for (CeedInt qpt=0; qpt<basis->Q; qpt++) 11359d007619Sjeremylt for (CeedInt node=0; node<basis->P; node++) { 1136d1d35e2fSjeremylt CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d; 1137d1d35e2fSjeremylt CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d; 11389d007619Sjeremylt if (i == d) 11399d007619Sjeremylt basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *= 1140d1d35e2fSjeremylt basis->grad_1d[q*basis->P_1d+p]; 11419d007619Sjeremylt else 11429d007619Sjeremylt basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *= 1143d1d35e2fSjeremylt basis->interp_1d[q*basis->P_1d+p]; 11449d007619Sjeremylt } 11459d007619Sjeremylt } 11469d007619Sjeremylt *grad = basis->grad; 1147e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 11489d007619Sjeremylt } 11499d007619Sjeremylt 11509d007619Sjeremylt /** 11519d007619Sjeremylt @brief Get 1D gradient matrix of a tensor product CeedBasis 11529d007619Sjeremylt 11539d007619Sjeremylt @param basis CeedBasis 1154d1d35e2fSjeremylt @param[out] grad_1d Variable to store gradient matrix 11559d007619Sjeremylt 11569d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 11579d007619Sjeremylt 1158b7c9bbdaSJeremy L Thompson @ref Advanced 11599d007619Sjeremylt **/ 1160d1d35e2fSjeremylt int CeedBasisGetGrad1D(CeedBasis basis, const CeedScalar **grad_1d) { 1161d1d35e2fSjeremylt if (!basis->tensor_basis) 11629d007619Sjeremylt // LCOV_EXCL_START 1163e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_MINOR, 1164e15f9bd0SJeremy L Thompson "CeedBasis is not a tensor product basis."); 11659d007619Sjeremylt // LCOV_EXCL_STOP 11669d007619Sjeremylt 1167d1d35e2fSjeremylt *grad_1d = basis->grad_1d; 1168e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 11699d007619Sjeremylt } 11709d007619Sjeremylt 11719d007619Sjeremylt /** 117250c301a5SRezgar Shakeri @brief Get divergence matrix of a CeedBasis 117350c301a5SRezgar Shakeri 117450c301a5SRezgar Shakeri @param basis CeedBasis 117550c301a5SRezgar Shakeri @param[out] div Variable to store divergence matrix 117650c301a5SRezgar Shakeri 117750c301a5SRezgar Shakeri @return An error code: 0 - success, otherwise - failure 117850c301a5SRezgar Shakeri 117950c301a5SRezgar Shakeri @ref Advanced 118050c301a5SRezgar Shakeri **/ 118150c301a5SRezgar Shakeri int CeedBasisGetDiv(CeedBasis basis, const CeedScalar **div) { 118250c301a5SRezgar Shakeri if (!basis->div) 118350c301a5SRezgar Shakeri // LCOV_EXCL_START 118450c301a5SRezgar Shakeri return CeedError(basis->ceed, CEED_ERROR_MINOR, 118550c301a5SRezgar Shakeri "CeedBasis does not have divergence matrix."); 118650c301a5SRezgar Shakeri // LCOV_EXCL_STOP 118750c301a5SRezgar Shakeri 118850c301a5SRezgar Shakeri *div = basis->div; 118950c301a5SRezgar Shakeri return CEED_ERROR_SUCCESS; 119050c301a5SRezgar Shakeri } 119150c301a5SRezgar Shakeri 119250c301a5SRezgar Shakeri /** 11937a982d89SJeremy L. Thompson @brief Destroy a CeedBasis 11947a982d89SJeremy L. Thompson 11957a982d89SJeremy L. Thompson @param basis CeedBasis to destroy 11967a982d89SJeremy L. Thompson 11977a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 11987a982d89SJeremy L. Thompson 11997a982d89SJeremy L. Thompson @ref User 12007a982d89SJeremy L. Thompson **/ 12017a982d89SJeremy L. Thompson int CeedBasisDestroy(CeedBasis *basis) { 12027a982d89SJeremy L. Thompson int ierr; 12037a982d89SJeremy L. Thompson 1204d1d35e2fSjeremylt if (!*basis || --(*basis)->ref_count > 0) return CEED_ERROR_SUCCESS; 12057a982d89SJeremy L. Thompson if ((*basis)->Destroy) { 12067a982d89SJeremy L. Thompson ierr = (*basis)->Destroy(*basis); CeedChk(ierr); 12077a982d89SJeremy L. Thompson } 120834359f16Sjeremylt if ((*basis)->contract) { 120934359f16Sjeremylt ierr = CeedTensorContractDestroy(&(*basis)->contract); CeedChk(ierr); 121034359f16Sjeremylt } 12117a982d89SJeremy L. Thompson ierr = CeedFree(&(*basis)->interp); CeedChk(ierr); 1212d1d35e2fSjeremylt ierr = CeedFree(&(*basis)->interp_1d); CeedChk(ierr); 12137a982d89SJeremy L. Thompson ierr = CeedFree(&(*basis)->grad); CeedChk(ierr); 121450c301a5SRezgar Shakeri ierr = CeedFree(&(*basis)->div); CeedChk(ierr); 1215d1d35e2fSjeremylt ierr = CeedFree(&(*basis)->grad_1d); CeedChk(ierr); 1216d1d35e2fSjeremylt ierr = CeedFree(&(*basis)->q_ref_1d); CeedChk(ierr); 1217d1d35e2fSjeremylt ierr = CeedFree(&(*basis)->q_weight_1d); CeedChk(ierr); 12187a982d89SJeremy L. Thompson ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr); 12197a982d89SJeremy L. Thompson ierr = CeedFree(basis); CeedChk(ierr); 1220e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 12217a982d89SJeremy L. Thompson } 12227a982d89SJeremy L. Thompson 12237a982d89SJeremy L. Thompson /** 1224b11c1e72Sjeremylt @brief Construct a Gauss-Legendre quadrature 1225b11c1e72Sjeremylt 1226b11c1e72Sjeremylt @param Q Number of quadrature points (integrates polynomials of 1227b11c1e72Sjeremylt degree 2*Q-1 exactly) 1228d1d35e2fSjeremylt @param[out] q_ref_1d Array of length Q to hold the abscissa on [-1, 1] 1229d1d35e2fSjeremylt @param[out] q_weight_1d Array of length Q to hold the weights 1230b11c1e72Sjeremylt 1231b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 1232dfdf5a53Sjeremylt 1233dfdf5a53Sjeremylt @ref Utility 1234b11c1e72Sjeremylt **/ 1235d1d35e2fSjeremylt int CeedGaussQuadrature(CeedInt Q, CeedScalar *q_ref_1d, 1236d1d35e2fSjeremylt CeedScalar *q_weight_1d) { 1237d7b241e6Sjeremylt // Allocate 1238d7b241e6Sjeremylt CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0*atan(1.0); 1239d1d35e2fSjeremylt // Build q_ref_1d, q_weight_1d 1240d7b241e6Sjeremylt for (int i = 0; i <= Q/2; i++) { 1241d7b241e6Sjeremylt // Guess 1242d7b241e6Sjeremylt xi = cos(PI*(CeedScalar)(2*i+1)/((CeedScalar)(2*Q))); 1243d7b241e6Sjeremylt // Pn(xi) 1244d7b241e6Sjeremylt P0 = 1.0; 1245d7b241e6Sjeremylt P1 = xi; 1246d7b241e6Sjeremylt P2 = 0.0; 1247d7b241e6Sjeremylt for (int j = 2; j <= Q; j++) { 1248d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 1249d7b241e6Sjeremylt P0 = P1; 1250d7b241e6Sjeremylt P1 = P2; 1251d7b241e6Sjeremylt } 1252d7b241e6Sjeremylt // First Newton Step 1253d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 1254d7b241e6Sjeremylt xi = xi-P2/dP2; 1255d7b241e6Sjeremylt // Newton to convergence 12560e4d4210Sjeremylt for (int k=0; k<100 && fabs(P2)>10*CEED_EPSILON; k++) { 1257d7b241e6Sjeremylt P0 = 1.0; 1258d7b241e6Sjeremylt P1 = xi; 1259d7b241e6Sjeremylt for (int j = 2; j <= Q; j++) { 1260d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 1261d7b241e6Sjeremylt P0 = P1; 1262d7b241e6Sjeremylt P1 = P2; 1263d7b241e6Sjeremylt } 1264d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 1265d7b241e6Sjeremylt xi = xi-P2/dP2; 1266d7b241e6Sjeremylt } 1267d7b241e6Sjeremylt // Save xi, wi 1268d7b241e6Sjeremylt wi = 2.0/((1.0-xi*xi)*dP2*dP2); 1269d1d35e2fSjeremylt q_weight_1d[i] = wi; 1270d1d35e2fSjeremylt q_weight_1d[Q-1-i] = wi; 1271d1d35e2fSjeremylt q_ref_1d[i] = -xi; 1272d1d35e2fSjeremylt q_ref_1d[Q-1-i]= xi; 1273d7b241e6Sjeremylt } 1274e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1275d7b241e6Sjeremylt } 1276d7b241e6Sjeremylt 1277b11c1e72Sjeremylt /** 1278b11c1e72Sjeremylt @brief Construct a Gauss-Legendre-Lobatto quadrature 1279b11c1e72Sjeremylt 1280b11c1e72Sjeremylt @param Q Number of quadrature points (integrates polynomials of 1281b11c1e72Sjeremylt degree 2*Q-3 exactly) 1282d1d35e2fSjeremylt @param[out] q_ref_1d Array of length Q to hold the abscissa on [-1, 1] 1283d1d35e2fSjeremylt @param[out] q_weight_1d Array of length Q to hold the weights 1284b11c1e72Sjeremylt 1285b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 1286dfdf5a53Sjeremylt 1287dfdf5a53Sjeremylt @ref Utility 1288b11c1e72Sjeremylt **/ 1289d1d35e2fSjeremylt int CeedLobattoQuadrature(CeedInt Q, CeedScalar *q_ref_1d, 1290d1d35e2fSjeremylt CeedScalar *q_weight_1d) { 1291d7b241e6Sjeremylt // Allocate 1292d7b241e6Sjeremylt CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0*atan(1.0); 1293d1d35e2fSjeremylt // Build q_ref_1d, q_weight_1d 1294d7b241e6Sjeremylt // Set endpoints 129530a100c3SJed Brown if (Q < 2) 1296b0d62198Sjeremylt // LCOV_EXCL_START 1297e15f9bd0SJeremy L Thompson return CeedError(NULL, CEED_ERROR_DIMENSION, 12987ed52d01Sjeremylt "Cannot create Lobatto quadrature with Q=%d < 2 points", Q); 1299b0d62198Sjeremylt // LCOV_EXCL_STOP 1300d7b241e6Sjeremylt wi = 2.0/((CeedScalar)(Q*(Q-1))); 1301d1d35e2fSjeremylt if (q_weight_1d) { 1302d1d35e2fSjeremylt q_weight_1d[0] = wi; 1303d1d35e2fSjeremylt q_weight_1d[Q-1] = wi; 1304d7b241e6Sjeremylt } 1305d1d35e2fSjeremylt q_ref_1d[0] = -1.0; 1306d1d35e2fSjeremylt q_ref_1d[Q-1] = 1.0; 1307d7b241e6Sjeremylt // Interior 1308d7b241e6Sjeremylt for (int i = 1; i <= (Q-1)/2; i++) { 1309d7b241e6Sjeremylt // Guess 1310d7b241e6Sjeremylt xi = cos(PI*(CeedScalar)(i)/(CeedScalar)(Q-1)); 1311d7b241e6Sjeremylt // Pn(xi) 1312d7b241e6Sjeremylt P0 = 1.0; 1313d7b241e6Sjeremylt P1 = xi; 1314d7b241e6Sjeremylt P2 = 0.0; 1315d7b241e6Sjeremylt for (int j = 2; j < Q; j++) { 1316d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 1317d7b241e6Sjeremylt P0 = P1; 1318d7b241e6Sjeremylt P1 = P2; 1319d7b241e6Sjeremylt } 1320d7b241e6Sjeremylt // First Newton step 1321d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 1322d7b241e6Sjeremylt d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi); 1323d7b241e6Sjeremylt xi = xi-dP2/d2P2; 1324d7b241e6Sjeremylt // Newton to convergence 13250e4d4210Sjeremylt for (int k=0; k<100 && fabs(dP2)>10*CEED_EPSILON; k++) { 1326d7b241e6Sjeremylt P0 = 1.0; 1327d7b241e6Sjeremylt P1 = xi; 1328d7b241e6Sjeremylt for (int j = 2; j < Q; j++) { 1329d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 1330d7b241e6Sjeremylt P0 = P1; 1331d7b241e6Sjeremylt P1 = P2; 1332d7b241e6Sjeremylt } 1333d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 1334d7b241e6Sjeremylt d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi); 1335d7b241e6Sjeremylt xi = xi-dP2/d2P2; 1336d7b241e6Sjeremylt } 1337d7b241e6Sjeremylt // Save xi, wi 1338d7b241e6Sjeremylt wi = 2.0/(((CeedScalar)(Q*(Q-1)))*P2*P2); 1339d1d35e2fSjeremylt if (q_weight_1d) { 1340d1d35e2fSjeremylt q_weight_1d[i] = wi; 1341d1d35e2fSjeremylt q_weight_1d[Q-1-i] = wi; 1342d7b241e6Sjeremylt } 1343d1d35e2fSjeremylt q_ref_1d[i] = -xi; 1344d1d35e2fSjeremylt q_ref_1d[Q-1-i]= xi; 1345d7b241e6Sjeremylt } 1346e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1347d7b241e6Sjeremylt } 1348d7b241e6Sjeremylt 1349dfdf5a53Sjeremylt /** 135095bb1877Svaleriabarra @brief Return QR Factorization of a matrix 1351b11c1e72Sjeremylt 135277645d7bSjeremylt @param ceed A Ceed context for error handling 135352bfb9bbSJeremy L Thompson @param[in,out] mat Row-major matrix to be factorized in place 135452bfb9bbSJeremy L Thompson @param[in,out] tau Vector of length m of scaling factors 1355b11c1e72Sjeremylt @param m Number of rows 1356b11c1e72Sjeremylt @param n Number of columns 1357b11c1e72Sjeremylt 1358b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 1359dfdf5a53Sjeremylt 1360dfdf5a53Sjeremylt @ref Utility 1361b11c1e72Sjeremylt **/ 1362a7bd39daSjeremylt int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau, 1363d7b241e6Sjeremylt CeedInt m, CeedInt n) { 1364d7b241e6Sjeremylt CeedScalar v[m]; 1365d7b241e6Sjeremylt 1366a7bd39daSjeremylt // Check m >= n 1367a7bd39daSjeremylt if (n > m) 1368c042f62fSJeremy L Thompson // LCOV_EXCL_START 1369e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_UNSUPPORTED, 1370e15f9bd0SJeremy L Thompson "Cannot compute QR factorization with n > m"); 1371c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 1372a7bd39daSjeremylt 137352bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) { 1374bde37e8cSJed Brown if (i >= m-1) { // last row of matrix, no reflection needed 1375bde37e8cSJed Brown tau[i] = 0.; 1376bde37e8cSJed Brown break; 1377bde37e8cSJed Brown } 1378d7b241e6Sjeremylt // Calculate Householder vector, magnitude 1379d7b241e6Sjeremylt CeedScalar sigma = 0.0; 1380d7b241e6Sjeremylt v[i] = mat[i+n*i]; 138152bfb9bbSJeremy L Thompson for (CeedInt j=i+1; j<m; j++) { 1382d7b241e6Sjeremylt v[j] = mat[i+n*j]; 1383d7b241e6Sjeremylt sigma += v[j] * v[j]; 1384d7b241e6Sjeremylt } 1385d7b241e6Sjeremylt CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m] 1386673160d7Sjeremylt CeedScalar R_ii = -copysign(norm, v[i]); 1387673160d7Sjeremylt v[i] -= R_ii; 1388d7b241e6Sjeremylt // norm of v[i:m] after modification above and scaling below 1389d7b241e6Sjeremylt // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 1390d7b241e6Sjeremylt // tau = 2 / (norm*norm) 1391d7b241e6Sjeremylt tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma); 13921d102b48SJeremy L Thompson for (CeedInt j=i+1; j<m; j++) 13931d102b48SJeremy L Thompson v[j] /= v[i]; 1394d7b241e6Sjeremylt 1395d7b241e6Sjeremylt // Apply Householder reflector to lower right panel 1396d7b241e6Sjeremylt CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1); 1397d7b241e6Sjeremylt // Save v 1398673160d7Sjeremylt mat[i+n*i] = R_ii; 13991d102b48SJeremy L Thompson for (CeedInt j=i+1; j<m; j++) 1400d7b241e6Sjeremylt mat[i+n*j] = v[j]; 1401d7b241e6Sjeremylt } 1402e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1403d7b241e6Sjeremylt } 1404d7b241e6Sjeremylt 1405b11c1e72Sjeremylt /** 140652bfb9bbSJeremy L Thompson @brief Return symmetric Schur decomposition of the symmetric matrix mat via 140752bfb9bbSJeremy L Thompson symmetric QR factorization 140852bfb9bbSJeremy L Thompson 140977645d7bSjeremylt @param ceed A Ceed context for error handling 141052bfb9bbSJeremy L Thompson @param[in,out] mat Row-major matrix to be factorized in place 1411460bf743SValeria Barra @param[out] lambda Vector of length n of eigenvalues 141252bfb9bbSJeremy L Thompson @param n Number of rows/columns 141352bfb9bbSJeremy L Thompson 141452bfb9bbSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 141552bfb9bbSJeremy L Thompson 141652bfb9bbSJeremy L Thompson @ref Utility 141752bfb9bbSJeremy L Thompson **/ 141803d18186Sjeremylt CeedPragmaOptimizeOff 141952bfb9bbSJeremy L Thompson int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat, 142052bfb9bbSJeremy L Thompson CeedScalar *lambda, CeedInt n) { 142152bfb9bbSJeremy L Thompson // Check bounds for clang-tidy 142252bfb9bbSJeremy L Thompson if (n<2) 1423c042f62fSJeremy L Thompson // LCOV_EXCL_START 1424e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_UNSUPPORTED, 1425c042f62fSJeremy L Thompson "Cannot compute symmetric Schur decomposition of scalars"); 1426c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 142752bfb9bbSJeremy L Thompson 1428673160d7Sjeremylt CeedScalar v[n-1], tau[n-1], mat_T[n*n]; 142952bfb9bbSJeremy L Thompson 1430673160d7Sjeremylt // Copy mat to mat_T and set mat to I 1431673160d7Sjeremylt memcpy(mat_T, mat, n*n*sizeof(mat[0])); 143252bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) 143352bfb9bbSJeremy L Thompson for (CeedInt j=0; j<n; j++) 143452bfb9bbSJeremy L Thompson mat[j+n*i] = (i==j) ? 1 : 0; 143552bfb9bbSJeremy L Thompson 143652bfb9bbSJeremy L Thompson // Reduce to tridiagonal 143752bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n-1; i++) { 143852bfb9bbSJeremy L Thompson // Calculate Householder vector, magnitude 143952bfb9bbSJeremy L Thompson CeedScalar sigma = 0.0; 1440673160d7Sjeremylt v[i] = mat_T[i+n*(i+1)]; 144152bfb9bbSJeremy L Thompson for (CeedInt j=i+1; j<n-1; j++) { 1442673160d7Sjeremylt v[j] = mat_T[i+n*(j+1)]; 144352bfb9bbSJeremy L Thompson sigma += v[j] * v[j]; 144452bfb9bbSJeremy L Thompson } 144552bfb9bbSJeremy L Thompson CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:n-1] 1446673160d7Sjeremylt CeedScalar R_ii = -copysign(norm, v[i]); 1447673160d7Sjeremylt v[i] -= R_ii; 144852bfb9bbSJeremy L Thompson // norm of v[i:m] after modification above and scaling below 144952bfb9bbSJeremy L Thompson // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 145052bfb9bbSJeremy L Thompson // tau = 2 / (norm*norm) 145180a9ef05SNatalie Beams tau[i] = i == n - 2 ? 2 : 2 * v[i]*v[i] / (v[i]*v[i] + sigma); 1452fb551037Sjeremylt for (CeedInt j=i+1; j<n-1; j++) 1453fb551037Sjeremylt v[j] /= v[i]; 145452bfb9bbSJeremy L Thompson 145552bfb9bbSJeremy L Thompson // Update sub and super diagonal 145652bfb9bbSJeremy L Thompson for (CeedInt j=i+2; j<n; j++) { 1457673160d7Sjeremylt mat_T[i+n*j] = 0; mat_T[j+n*i] = 0; 145852bfb9bbSJeremy L Thompson } 145952bfb9bbSJeremy L Thompson // Apply symmetric Householder reflector to lower right panel 1460673160d7Sjeremylt CeedHouseholderReflect(&mat_T[(i+1)+n*(i+1)], &v[i], tau[i], 146152bfb9bbSJeremy L Thompson n-(i+1), n-(i+1), n, 1); 1462673160d7Sjeremylt CeedHouseholderReflect(&mat_T[(i+1)+n*(i+1)], &v[i], tau[i], 146352bfb9bbSJeremy L Thompson n-(i+1), n-(i+1), 1, n); 1464673160d7Sjeremylt 146552bfb9bbSJeremy L Thompson // Save v 1466673160d7Sjeremylt mat_T[i+n*(i+1)] = R_ii; 1467673160d7Sjeremylt mat_T[(i+1)+n*i] = R_ii; 146852bfb9bbSJeremy L Thompson for (CeedInt j=i+1; j<n-1; j++) { 1469673160d7Sjeremylt mat_T[i+n*(j+1)] = v[j]; 147052bfb9bbSJeremy L Thompson } 147152bfb9bbSJeremy L Thompson } 147252bfb9bbSJeremy L Thompson // Backwards accumulation of Q 147352bfb9bbSJeremy L Thompson for (CeedInt i=n-2; i>=0; i--) { 147485cf89eaSjeremylt if (tau[i] > 0.0) { 147552bfb9bbSJeremy L Thompson v[i] = 1; 147652bfb9bbSJeremy L Thompson for (CeedInt j=i+1; j<n-1; j++) { 1477673160d7Sjeremylt v[j] = mat_T[i+n*(j+1)]; 1478673160d7Sjeremylt mat_T[i+n*(j+1)] = 0; 147952bfb9bbSJeremy L Thompson } 148052bfb9bbSJeremy L Thompson CeedHouseholderReflect(&mat[(i+1)+n*(i+1)], &v[i], tau[i], 148152bfb9bbSJeremy L Thompson n-(i+1), n-(i+1), n, 1); 148252bfb9bbSJeremy L Thompson } 148385cf89eaSjeremylt } 148452bfb9bbSJeremy L Thompson 148552bfb9bbSJeremy L Thompson // Reduce sub and super diagonal 1486673160d7Sjeremylt CeedInt p = 0, q = 0, itr = 0, max_itr = n*n*n*n; 1487673160d7Sjeremylt CeedScalar tol = CEED_EPSILON; 148852bfb9bbSJeremy L Thompson 1489673160d7Sjeremylt while (itr < max_itr) { 149052bfb9bbSJeremy L Thompson // Update p, q, size of reduced portions of diagonal 149152bfb9bbSJeremy L Thompson p = 0; q = 0; 149252bfb9bbSJeremy L Thompson for (CeedInt i=n-2; i>=0; i--) { 1493673160d7Sjeremylt if (fabs(mat_T[i+n*(i+1)]) < tol) 149452bfb9bbSJeremy L Thompson q += 1; 149552bfb9bbSJeremy L Thompson else 149652bfb9bbSJeremy L Thompson break; 149752bfb9bbSJeremy L Thompson } 1498673160d7Sjeremylt for (CeedInt i=0; i<n-q-1; i++) { 1499673160d7Sjeremylt if (fabs(mat_T[i+n*(i+1)]) < tol) 150052bfb9bbSJeremy L Thompson p += 1; 150152bfb9bbSJeremy L Thompson else 150252bfb9bbSJeremy L Thompson break; 150352bfb9bbSJeremy L Thompson } 150452bfb9bbSJeremy L Thompson if (q == n-1) break; // Finished reducing 150552bfb9bbSJeremy L Thompson 150652bfb9bbSJeremy L Thompson // Reduce tridiagonal portion 1507673160d7Sjeremylt CeedScalar t_nn = mat_T[(n-1-q)+n*(n-1-q)], 1508673160d7Sjeremylt t_nnm1 = mat_T[(n-2-q)+n*(n-1-q)]; 1509673160d7Sjeremylt CeedScalar d = (mat_T[(n-2-q)+n*(n-2-q)] - t_nn)/2; 1510673160d7Sjeremylt CeedScalar mu = t_nn - t_nnm1*t_nnm1 / 1511673160d7Sjeremylt (d + copysign(sqrt(d*d + t_nnm1*t_nnm1), d)); 1512673160d7Sjeremylt CeedScalar x = mat_T[p+n*p] - mu; 1513673160d7Sjeremylt CeedScalar z = mat_T[p+n*(p+1)]; 1514673160d7Sjeremylt for (CeedInt k=p; k<n-q-1; k++) { 151552bfb9bbSJeremy L Thompson // Compute Givens rotation 151652bfb9bbSJeremy L Thompson CeedScalar c = 1, s = 0; 151752bfb9bbSJeremy L Thompson if (fabs(z) > tol) { 151852bfb9bbSJeremy L Thompson if (fabs(z) > fabs(x)) { 151952bfb9bbSJeremy L Thompson CeedScalar tau = -x/z; 152052bfb9bbSJeremy L Thompson s = 1/sqrt(1+tau*tau), c = s*tau; 152152bfb9bbSJeremy L Thompson } else { 152252bfb9bbSJeremy L Thompson CeedScalar tau = -z/x; 152352bfb9bbSJeremy L Thompson c = 1/sqrt(1+tau*tau), s = c*tau; 152452bfb9bbSJeremy L Thompson } 152552bfb9bbSJeremy L Thompson } 152652bfb9bbSJeremy L Thompson 152752bfb9bbSJeremy L Thompson // Apply Givens rotation to T 1528673160d7Sjeremylt CeedGivensRotation(mat_T, c, s, CEED_NOTRANSPOSE, k, k+1, n, n); 1529673160d7Sjeremylt CeedGivensRotation(mat_T, c, s, CEED_TRANSPOSE, k, k+1, n, n); 153052bfb9bbSJeremy L Thompson 153152bfb9bbSJeremy L Thompson // Apply Givens rotation to Q 153252bfb9bbSJeremy L Thompson CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k+1, n, n); 153352bfb9bbSJeremy L Thompson 153452bfb9bbSJeremy L Thompson // Update x, z 153552bfb9bbSJeremy L Thompson if (k < n-q-2) { 1536673160d7Sjeremylt x = mat_T[k+n*(k+1)]; 1537673160d7Sjeremylt z = mat_T[k+n*(k+2)]; 153852bfb9bbSJeremy L Thompson } 153952bfb9bbSJeremy L Thompson } 154052bfb9bbSJeremy L Thompson itr++; 154152bfb9bbSJeremy L Thompson } 1542673160d7Sjeremylt 154352bfb9bbSJeremy L Thompson // Save eigenvalues 154452bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) 1545673160d7Sjeremylt lambda[i] = mat_T[i+n*i]; 154652bfb9bbSJeremy L Thompson 154752bfb9bbSJeremy L Thompson // Check convergence 1548673160d7Sjeremylt if (itr == max_itr && q < n-1) 1549c042f62fSJeremy L Thompson // LCOV_EXCL_START 1550e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_MINOR, 1551e15f9bd0SJeremy L Thompson "Symmetric QR failed to converge"); 1552c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 1553e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 155452bfb9bbSJeremy L Thompson } 155503d18186Sjeremylt CeedPragmaOptimizeOn 155652bfb9bbSJeremy L Thompson 155752bfb9bbSJeremy L Thompson /** 155852bfb9bbSJeremy L Thompson @brief Return Simultaneous Diagonalization of two matrices. This solves the 155952bfb9bbSJeremy L Thompson generalized eigenvalue problem A x = lambda B x, where A and B 156052bfb9bbSJeremy L Thompson are symmetric and B is positive definite. We generate the matrix X 156152bfb9bbSJeremy L Thompson and vector Lambda such that X^T A X = Lambda and X^T B X = I. This 156252bfb9bbSJeremy L Thompson is equivalent to the LAPACK routine 'sygv' with TYPE = 1. 156352bfb9bbSJeremy L Thompson 156477645d7bSjeremylt @param ceed A Ceed context for error handling 1565d1d35e2fSjeremylt @param[in] mat_A Row-major matrix to be factorized with eigenvalues 1566d1d35e2fSjeremylt @param[in] mat_B Row-major matrix to be factorized to identity 1567d3331725Sjeremylt @param[out] mat_X Row-major orthogonal matrix 1568460bf743SValeria Barra @param[out] lambda Vector of length n of generalized eigenvalues 156952bfb9bbSJeremy L Thompson @param n Number of rows/columns 157052bfb9bbSJeremy L Thompson 157152bfb9bbSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 157252bfb9bbSJeremy L Thompson 157352bfb9bbSJeremy L Thompson @ref Utility 157452bfb9bbSJeremy L Thompson **/ 157503d18186Sjeremylt CeedPragmaOptimizeOff 1576d1d35e2fSjeremylt int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *mat_A, 1577d3331725Sjeremylt CeedScalar *mat_B, CeedScalar *mat_X, 157852bfb9bbSJeremy L Thompson CeedScalar *lambda, CeedInt n) { 157952bfb9bbSJeremy L Thompson int ierr; 1580d3331725Sjeremylt CeedScalar *mat_C, *mat_G, *vec_D; 158178464608Sjeremylt ierr = CeedCalloc(n*n, &mat_C); CeedChk(ierr); 158278464608Sjeremylt ierr = CeedCalloc(n*n, &mat_G); CeedChk(ierr); 1583d3331725Sjeremylt ierr = CeedCalloc(n, &vec_D); CeedChk(ierr); 158452bfb9bbSJeremy L Thompson 158552bfb9bbSJeremy L Thompson // Compute B = G D G^T 158678464608Sjeremylt memcpy(mat_G, mat_B, n*n*sizeof(mat_B[0])); 1587d3331725Sjeremylt ierr = CeedSymmetricSchurDecomposition(ceed, mat_G, vec_D, n); CeedChk(ierr); 158852bfb9bbSJeremy L Thompson 158985cf89eaSjeremylt // Sort eigenvalues 159085cf89eaSjeremylt for (CeedInt i=n-1; i>=0; i--) 159185cf89eaSjeremylt for (CeedInt j=0; j<i; j++) { 159285cf89eaSjeremylt if (fabs(vec_D[j]) > fabs(vec_D[j+1])) { 159385cf89eaSjeremylt CeedScalar temp; 159485cf89eaSjeremylt temp = vec_D[j]; vec_D[j] = vec_D[j+1]; vec_D[j+1] = temp; 159585cf89eaSjeremylt for (CeedInt k=0; k<n; k++) { 159685cf89eaSjeremylt temp = mat_G[k*n+j]; mat_G[k*n+j] = mat_G[k*n+j+1]; mat_G[k*n+j+1] = temp; 159785cf89eaSjeremylt } 159885cf89eaSjeremylt } 159985cf89eaSjeremylt } 160085cf89eaSjeremylt 1601fb551037Sjeremylt // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T 1602fb551037Sjeremylt // = D^-1/2 G^T A G D^-1/2 1603d3331725Sjeremylt // -- D = D^-1/2 160452bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) 1605d3331725Sjeremylt vec_D[i] = 1./sqrt(vec_D[i]); 1606d3331725Sjeremylt // -- G = G D^-1/2 1607d3331725Sjeremylt // -- C = D^-1/2 G^T 1608d3331725Sjeremylt for (CeedInt i=0; i<n; i++) 1609d3331725Sjeremylt for (CeedInt j=0; j<n; j++) { 1610673160d7Sjeremylt mat_G[i*n+j] *= vec_D[j]; 1611673160d7Sjeremylt mat_C[j*n+i] = mat_G[i*n+j]; 1612d3331725Sjeremylt } 1613673160d7Sjeremylt // -- X = (D^-1/2 G^T) A 1614d1d35e2fSjeremylt ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)mat_C, 1615d3331725Sjeremylt (const CeedScalar *)mat_A, mat_X, n, n, n); 16169289e5bfSjeremylt CeedChk(ierr); 1617673160d7Sjeremylt // -- C = (D^-1/2 G^T A) (G D^-1/2) 1618d3331725Sjeremylt ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)mat_X, 161978464608Sjeremylt (const CeedScalar *)mat_G, mat_C, n, n, n); 16209289e5bfSjeremylt CeedChk(ierr); 162152bfb9bbSJeremy L Thompson 162252bfb9bbSJeremy L Thompson // Compute Q^T C Q = lambda 1623d1d35e2fSjeremylt ierr = CeedSymmetricSchurDecomposition(ceed, mat_C, lambda, n); CeedChk(ierr); 162452bfb9bbSJeremy L Thompson 162585cf89eaSjeremylt // Sort eigenvalues 162685cf89eaSjeremylt for (CeedInt i=n-1; i>=0; i--) 162785cf89eaSjeremylt for (CeedInt j=0; j<i; j++) { 162885cf89eaSjeremylt if (fabs(lambda[j]) > fabs(lambda[j+1])) { 162985cf89eaSjeremylt CeedScalar temp; 163085cf89eaSjeremylt temp = lambda[j]; lambda[j] = lambda[j+1]; lambda[j+1] = temp; 163185cf89eaSjeremylt for (CeedInt k=0; k<n; k++) { 163285cf89eaSjeremylt temp = mat_C[k*n+j]; mat_C[k*n+j] = mat_C[k*n+j+1]; mat_C[k*n+j+1] = temp; 163385cf89eaSjeremylt } 163485cf89eaSjeremylt } 163585cf89eaSjeremylt } 163685cf89eaSjeremylt 1637d3331725Sjeremylt // Set X = (G D^1/2)^-T Q 1638fb551037Sjeremylt // = G D^-1/2 Q 163978464608Sjeremylt ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)mat_G, 1640d3331725Sjeremylt (const CeedScalar *)mat_C, mat_X, n, n, n); 16419289e5bfSjeremylt CeedChk(ierr); 164278464608Sjeremylt 164378464608Sjeremylt // Cleanup 164478464608Sjeremylt ierr = CeedFree(&mat_C); CeedChk(ierr); 164578464608Sjeremylt ierr = CeedFree(&mat_G); CeedChk(ierr); 1646d3331725Sjeremylt ierr = CeedFree(&vec_D); CeedChk(ierr); 1647e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 164852bfb9bbSJeremy L Thompson } 164903d18186Sjeremylt CeedPragmaOptimizeOn 165052bfb9bbSJeremy L Thompson 1651d7b241e6Sjeremylt /// @} 1652