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 Ceed associated with a CeedBasis 2447a982d89SJeremy L. Thompson 2457a982d89SJeremy L. Thompson @param basis CeedBasis 2467a982d89SJeremy L. Thompson @param[out] ceed Variable to store Ceed 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 **/ 2527a982d89SJeremy L. Thompson int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) { 2537a982d89SJeremy L. Thompson *ceed = basis->ceed; 254e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 2557a982d89SJeremy L. Thompson } 2567a982d89SJeremy L. Thompson 2577a982d89SJeremy L. Thompson /** 2587a982d89SJeremy L. Thompson @brief Get tensor status for given CeedBasis 2597a982d89SJeremy L. Thompson 2607a982d89SJeremy L. Thompson @param basis CeedBasis 261d1d35e2fSjeremylt @param[out] is_tensor Variable to store tensor status 2627a982d89SJeremy L. Thompson 2637a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 2647a982d89SJeremy L. Thompson 2657a982d89SJeremy L. Thompson @ref Backend 2667a982d89SJeremy L. Thompson **/ 267d1d35e2fSjeremylt int CeedBasisIsTensor(CeedBasis basis, bool *is_tensor) { 268d1d35e2fSjeremylt *is_tensor = basis->tensor_basis; 269e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 2707a982d89SJeremy L. Thompson } 2717a982d89SJeremy L. Thompson 2727a982d89SJeremy L. Thompson /** 2737a982d89SJeremy L. Thompson @brief Get backend data of a CeedBasis 2747a982d89SJeremy L. Thompson 2757a982d89SJeremy L. Thompson @param basis CeedBasis 2767a982d89SJeremy L. Thompson @param[out] data Variable to store data 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 CeedBasisGetData(CeedBasis basis, void *data) { 283777ff853SJeremy L Thompson *(void **)data = basis->data; 284e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 2857a982d89SJeremy L. Thompson } 2867a982d89SJeremy L. Thompson 2877a982d89SJeremy L. Thompson /** 2887a982d89SJeremy L. Thompson @brief Set backend data of a CeedBasis 2897a982d89SJeremy L. Thompson 2907a982d89SJeremy L. Thompson @param[out] basis CeedBasis 2917a982d89SJeremy L. Thompson @param data Data to set 2927a982d89SJeremy L. Thompson 2937a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 2947a982d89SJeremy L. Thompson 2957a982d89SJeremy L. Thompson @ref Backend 2967a982d89SJeremy L. Thompson **/ 297777ff853SJeremy L Thompson int CeedBasisSetData(CeedBasis basis, void *data) { 298777ff853SJeremy L Thompson basis->data = data; 299e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3007a982d89SJeremy L. Thompson } 3017a982d89SJeremy L. Thompson 3027a982d89SJeremy L. Thompson /** 30334359f16Sjeremylt @brief Increment the reference counter for a CeedBasis 30434359f16Sjeremylt 30534359f16Sjeremylt @param basis Basis to increment the reference counter 30634359f16Sjeremylt 30734359f16Sjeremylt @return An error code: 0 - success, otherwise - failure 30834359f16Sjeremylt 30934359f16Sjeremylt @ref Backend 31034359f16Sjeremylt **/ 3119560d06aSjeremylt int CeedBasisReference(CeedBasis basis) { 31234359f16Sjeremylt basis->ref_count++; 31334359f16Sjeremylt return CEED_ERROR_SUCCESS; 31434359f16Sjeremylt } 31534359f16Sjeremylt 31634359f16Sjeremylt /** 3177a982d89SJeremy L. Thompson @brief Get dimension for given CeedElemTopology 3187a982d89SJeremy L. Thompson 3197a982d89SJeremy L. Thompson @param topo CeedElemTopology 3207a982d89SJeremy L. Thompson @param[out] dim Variable to store dimension of topology 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 CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) { 3277a982d89SJeremy L. Thompson *dim = (CeedInt) topo >> 16; 328e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3297a982d89SJeremy L. Thompson } 3307a982d89SJeremy L. Thompson 3317a982d89SJeremy L. Thompson /** 3327a982d89SJeremy L. Thompson @brief Get CeedTensorContract of a CeedBasis 3337a982d89SJeremy L. Thompson 3347a982d89SJeremy L. Thompson @param basis CeedBasis 3357a982d89SJeremy L. Thompson @param[out] contract Variable to store CeedTensorContract 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 **/ 3417a982d89SJeremy L. Thompson int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) { 3427a982d89SJeremy L. Thompson *contract = basis->contract; 343e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3447a982d89SJeremy L. Thompson } 3457a982d89SJeremy L. Thompson 3467a982d89SJeremy L. Thompson /** 3477a982d89SJeremy L. Thompson @brief Set CeedTensorContract of a CeedBasis 3487a982d89SJeremy L. Thompson 3497a982d89SJeremy L. Thompson @param[out] basis CeedBasis 3507a982d89SJeremy L. Thompson @param contract CeedTensorContract to set 3517a982d89SJeremy L. Thompson 3527a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3537a982d89SJeremy L. Thompson 3547a982d89SJeremy L. Thompson @ref Backend 3557a982d89SJeremy L. Thompson **/ 35634359f16Sjeremylt int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract contract) { 35734359f16Sjeremylt int ierr; 35834359f16Sjeremylt basis->contract = contract; 3599560d06aSjeremylt ierr = CeedTensorContractReference(contract); CeedChk(ierr); 360e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3617a982d89SJeremy L. Thompson } 3627a982d89SJeremy L. Thompson 3637a982d89SJeremy L. Thompson /** 3647a982d89SJeremy L. Thompson @brief Return a reference implementation of matrix multiplication C = A B. 3657a982d89SJeremy L. Thompson Note, this is a reference implementation for CPU CeedScalar pointers 3667a982d89SJeremy L. Thompson that is not intended for high performance. 3677a982d89SJeremy L. Thompson 3687a982d89SJeremy L. Thompson @param ceed A Ceed context for error handling 369d1d35e2fSjeremylt @param[in] mat_A Row-major matrix A 370d1d35e2fSjeremylt @param[in] mat_B Row-major matrix B 371d1d35e2fSjeremylt @param[out] mat_C Row-major output matrix C 3727a982d89SJeremy L. Thompson @param m Number of rows of C 3737a982d89SJeremy L. Thompson @param n Number of columns of C 3747a982d89SJeremy L. Thompson @param kk Number of columns of A/rows of B 3757a982d89SJeremy L. Thompson 3767a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3777a982d89SJeremy L. Thompson 3787a982d89SJeremy L. Thompson @ref Utility 3797a982d89SJeremy L. Thompson **/ 380d1d35e2fSjeremylt int CeedMatrixMultiply(Ceed ceed, const CeedScalar *mat_A, 381d1d35e2fSjeremylt const CeedScalar *mat_B, CeedScalar *mat_C, CeedInt m, 3827a982d89SJeremy L. Thompson CeedInt n, CeedInt kk) { 3837a982d89SJeremy L. Thompson for (CeedInt i=0; i<m; i++) 3847a982d89SJeremy L. Thompson for (CeedInt j=0; j<n; j++) { 3857a982d89SJeremy L. Thompson CeedScalar sum = 0; 3867a982d89SJeremy L. Thompson for (CeedInt k=0; k<kk; k++) 387d1d35e2fSjeremylt sum += mat_A[k+i*kk]*mat_B[j+k*n]; 388d1d35e2fSjeremylt mat_C[j+i*n] = sum; 3897a982d89SJeremy L. Thompson } 390e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3917a982d89SJeremy L. Thompson } 3927a982d89SJeremy L. Thompson 3937a982d89SJeremy L. Thompson /// @} 3947a982d89SJeremy L. Thompson 3957a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 3967a982d89SJeremy L. Thompson /// CeedBasis Public API 3977a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 3987a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser 399d7b241e6Sjeremylt /// @{ 400d7b241e6Sjeremylt 401b11c1e72Sjeremylt /** 40295bb1877Svaleriabarra @brief Create a tensor-product basis for H^1 discretizations 403b11c1e72Sjeremylt 404b11c1e72Sjeremylt @param ceed A Ceed object where the CeedBasis will be created 405b11c1e72Sjeremylt @param dim Topological dimension 406d1d35e2fSjeremylt @param num_comp Number of field components (1 for scalar fields) 407d1d35e2fSjeremylt @param P_1d Number of nodes in one dimension 408d1d35e2fSjeremylt @param Q_1d Number of quadrature points in one dimension 409d1d35e2fSjeremylt @param interp_1d Row-major (Q_1d * P_1d) matrix expressing the values of nodal 410b11c1e72Sjeremylt basis functions at quadrature points 411d1d35e2fSjeremylt @param grad_1d Row-major (Q_1d * P_1d) matrix expressing derivatives of nodal 412b11c1e72Sjeremylt basis functions at quadrature points 413d1d35e2fSjeremylt @param q_ref_1d Array of length Q_1d holding the locations of quadrature points 414b11c1e72Sjeremylt on the 1D reference element [-1, 1] 415d1d35e2fSjeremylt @param q_weight_1d Array of length Q_1d holding the quadrature weights on the 416b11c1e72Sjeremylt reference element 417b11c1e72Sjeremylt @param[out] basis Address of the variable where the newly created 418b11c1e72Sjeremylt CeedBasis will be stored. 419b11c1e72Sjeremylt 420b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 421dfdf5a53Sjeremylt 4227a982d89SJeremy L. Thompson @ref User 423b11c1e72Sjeremylt **/ 424d1d35e2fSjeremylt int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt num_comp, 425d1d35e2fSjeremylt CeedInt P_1d, CeedInt Q_1d, 426d1d35e2fSjeremylt const CeedScalar *interp_1d, 427d1d35e2fSjeremylt const CeedScalar *grad_1d, const CeedScalar *q_ref_1d, 428d1d35e2fSjeremylt const CeedScalar *q_weight_1d, CeedBasis *basis) { 429d7b241e6Sjeremylt int ierr; 430d7b241e6Sjeremylt 4315fe0d4faSjeremylt if (!ceed->BasisCreateTensorH1) { 4325fe0d4faSjeremylt Ceed delegate; 433aefd8378Sjeremylt ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr); 4345fe0d4faSjeremylt 4355fe0d4faSjeremylt if (!delegate) 436c042f62fSJeremy L Thompson // LCOV_EXCL_START 437e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_UNSUPPORTED, 438e15f9bd0SJeremy L Thompson "Backend does not support BasisCreateTensorH1"); 439c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 4405fe0d4faSjeremylt 441d1d35e2fSjeremylt ierr = CeedBasisCreateTensorH1(delegate, dim, num_comp, P_1d, 442d1d35e2fSjeremylt Q_1d, interp_1d, grad_1d, q_ref_1d, 443d1d35e2fSjeremylt q_weight_1d, basis); CeedChk(ierr); 444e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 4455fe0d4faSjeremylt } 446e15f9bd0SJeremy L Thompson 447e15f9bd0SJeremy L Thompson if (dim<1) 448e15f9bd0SJeremy L Thompson // LCOV_EXCL_START 449e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, 450e15f9bd0SJeremy L Thompson "Basis dimension must be a positive value"); 451e15f9bd0SJeremy L Thompson // LCOV_EXCL_STOP 452d1d35e2fSjeremylt CeedElemTopology topo = dim == 1 ? CEED_LINE 453d1d35e2fSjeremylt : dim == 2 ? CEED_QUAD 454d1d35e2fSjeremylt : CEED_HEX; 455e15f9bd0SJeremy L Thompson 456d7b241e6Sjeremylt ierr = CeedCalloc(1, basis); CeedChk(ierr); 457d7b241e6Sjeremylt (*basis)->ceed = ceed; 4589560d06aSjeremylt ierr = CeedReference(ceed); CeedChk(ierr); 459d1d35e2fSjeremylt (*basis)->ref_count = 1; 460d1d35e2fSjeremylt (*basis)->tensor_basis = 1; 461d7b241e6Sjeremylt (*basis)->dim = dim; 462d99fa3c5SJeremy L Thompson (*basis)->topo = topo; 463d1d35e2fSjeremylt (*basis)->num_comp = num_comp; 464d1d35e2fSjeremylt (*basis)->P_1d = P_1d; 465d1d35e2fSjeremylt (*basis)->Q_1d = Q_1d; 466d1d35e2fSjeremylt (*basis)->P = CeedIntPow(P_1d, dim); 467d1d35e2fSjeremylt (*basis)->Q = CeedIntPow(Q_1d, dim); 468d1d35e2fSjeremylt ierr = CeedMalloc(Q_1d,&(*basis)->q_ref_1d); CeedChk(ierr); 469d1d35e2fSjeremylt ierr = CeedMalloc(Q_1d,&(*basis)->q_weight_1d); CeedChk(ierr); 470d1d35e2fSjeremylt memcpy((*basis)->q_ref_1d, q_ref_1d, Q_1d*sizeof(q_ref_1d[0])); 471d1d35e2fSjeremylt memcpy((*basis)->q_weight_1d, q_weight_1d, Q_1d*sizeof(q_weight_1d[0])); 472d1d35e2fSjeremylt ierr = CeedMalloc(Q_1d*P_1d,&(*basis)->interp_1d); CeedChk(ierr); 473d1d35e2fSjeremylt ierr = CeedMalloc(Q_1d*P_1d,&(*basis)->grad_1d); CeedChk(ierr); 474d1d35e2fSjeremylt memcpy((*basis)->interp_1d, interp_1d, Q_1d*P_1d*sizeof(interp_1d[0])); 475d1d35e2fSjeremylt memcpy((*basis)->grad_1d, grad_1d, Q_1d*P_1d*sizeof(grad_1d[0])); 476d1d35e2fSjeremylt ierr = ceed->BasisCreateTensorH1(dim, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, 477d1d35e2fSjeremylt q_weight_1d, *basis); CeedChk(ierr); 478e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 479d7b241e6Sjeremylt } 480d7b241e6Sjeremylt 481b11c1e72Sjeremylt /** 48295bb1877Svaleriabarra @brief Create a tensor-product Lagrange basis 483b11c1e72Sjeremylt 484b11c1e72Sjeremylt @param ceed A Ceed object where the CeedBasis will be created 485b11c1e72Sjeremylt @param dim Topological dimension of element 486d1d35e2fSjeremylt @param num_comp Number of field components (1 for scalar fields) 487b11c1e72Sjeremylt @param P Number of Gauss-Lobatto nodes in one dimension. The 488b11c1e72Sjeremylt polynomial degree of the resulting Q_k element is k=P-1. 489b11c1e72Sjeremylt @param Q Number of quadrature points in one dimension. 490d1d35e2fSjeremylt @param quad_mode Distribution of the Q quadrature points (affects order of 491b11c1e72Sjeremylt accuracy for the quadrature) 492b11c1e72Sjeremylt @param[out] basis Address of the variable where the newly created 493b11c1e72Sjeremylt CeedBasis will be stored. 494b11c1e72Sjeremylt 495b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 496dfdf5a53Sjeremylt 4977a982d89SJeremy L. Thompson @ref User 498b11c1e72Sjeremylt **/ 499d1d35e2fSjeremylt int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt num_comp, 500d1d35e2fSjeremylt CeedInt P, CeedInt Q, CeedQuadMode quad_mode, 501692c2638Sjeremylt CeedBasis *basis) { 502d7b241e6Sjeremylt // Allocate 503e15f9bd0SJeremy L Thompson int ierr, ierr2, i, j, k; 504d1d35e2fSjeremylt CeedScalar c1, c2, c3, c4, dx, *nodes, *interp_1d, *grad_1d, *q_ref_1d, 505d1d35e2fSjeremylt *q_weight_1d; 5064d537eeaSYohann 5074d537eeaSYohann if (dim<1) 508c042f62fSJeremy L Thompson // LCOV_EXCL_START 509e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, 510e15f9bd0SJeremy L Thompson "Basis dimension must be a positive value"); 511c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 5124d537eeaSYohann 513e15f9bd0SJeremy L Thompson // Get Nodes and Weights 514d1d35e2fSjeremylt ierr = CeedCalloc(P*Q, &interp_1d); CeedChk(ierr); 515d1d35e2fSjeremylt ierr = CeedCalloc(P*Q, &grad_1d); CeedChk(ierr); 516d7b241e6Sjeremylt ierr = CeedCalloc(P, &nodes); CeedChk(ierr); 517d1d35e2fSjeremylt ierr = CeedCalloc(Q, &q_ref_1d); CeedChk(ierr); 518d1d35e2fSjeremylt ierr = CeedCalloc(Q, &q_weight_1d); CeedChk(ierr); 519e15f9bd0SJeremy L Thompson ierr = CeedLobattoQuadrature(P, nodes, NULL); 520e15f9bd0SJeremy L Thompson if (ierr) { goto cleanup; } CeedChk(ierr); 521d1d35e2fSjeremylt switch (quad_mode) { 522d7b241e6Sjeremylt case CEED_GAUSS: 523d1d35e2fSjeremylt ierr = CeedGaussQuadrature(Q, q_ref_1d, q_weight_1d); 524d7b241e6Sjeremylt break; 525d7b241e6Sjeremylt case CEED_GAUSS_LOBATTO: 526d1d35e2fSjeremylt ierr = CeedLobattoQuadrature(Q, q_ref_1d, q_weight_1d); 527d7b241e6Sjeremylt break; 528d7b241e6Sjeremylt } 529e15f9bd0SJeremy L Thompson if (ierr) { goto cleanup; } CeedChk(ierr); 530e15f9bd0SJeremy L Thompson 531d7b241e6Sjeremylt // Build B, D matrix 532d7b241e6Sjeremylt // Fornberg, 1998 533d7b241e6Sjeremylt for (i = 0; i < Q; i++) { 534d7b241e6Sjeremylt c1 = 1.0; 535d1d35e2fSjeremylt c3 = nodes[0] - q_ref_1d[i]; 536d1d35e2fSjeremylt interp_1d[i*P+0] = 1.0; 537d7b241e6Sjeremylt for (j = 1; j < P; j++) { 538d7b241e6Sjeremylt c2 = 1.0; 539d7b241e6Sjeremylt c4 = c3; 540d1d35e2fSjeremylt c3 = nodes[j] - q_ref_1d[i]; 541d7b241e6Sjeremylt for (k = 0; k < j; k++) { 542d7b241e6Sjeremylt dx = nodes[j] - nodes[k]; 543d7b241e6Sjeremylt c2 *= dx; 544d7b241e6Sjeremylt if (k == j - 1) { 545d1d35e2fSjeremylt grad_1d[i*P + j] = c1*(interp_1d[i*P + k] - c4*grad_1d[i*P + k]) / c2; 546d1d35e2fSjeremylt interp_1d[i*P + j] = - c1*c4*interp_1d[i*P + k] / c2; 547d7b241e6Sjeremylt } 548d1d35e2fSjeremylt grad_1d[i*P + k] = (c3*grad_1d[i*P + k] - interp_1d[i*P + k]) / dx; 549d1d35e2fSjeremylt interp_1d[i*P + k] = c3*interp_1d[i*P + k] / dx; 550d7b241e6Sjeremylt } 551d7b241e6Sjeremylt c1 = c2; 552d7b241e6Sjeremylt } 553d7b241e6Sjeremylt } 5549ac7b42eSJeremy L Thompson // Pass to CeedBasisCreateTensorH1 555d1d35e2fSjeremylt ierr = CeedBasisCreateTensorH1(ceed, dim, num_comp, P, Q, interp_1d, grad_1d, 5569ac7b42eSJeremy L Thompson q_ref_1d, q_weight_1d, basis); CeedChk(ierr); 557e15f9bd0SJeremy L Thompson cleanup: 558d1d35e2fSjeremylt ierr2 = CeedFree(&interp_1d); CeedChk(ierr2); 559d1d35e2fSjeremylt ierr2 = CeedFree(&grad_1d); CeedChk(ierr2); 560e15f9bd0SJeremy L Thompson ierr2 = CeedFree(&nodes); CeedChk(ierr2); 561d1d35e2fSjeremylt ierr2 = CeedFree(&q_ref_1d); CeedChk(ierr2); 562d1d35e2fSjeremylt ierr2 = CeedFree(&q_weight_1d); CeedChk(ierr2); 563e15f9bd0SJeremy L Thompson CeedChk(ierr); 564e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 565d7b241e6Sjeremylt } 566d7b241e6Sjeremylt 567b11c1e72Sjeremylt /** 56895bb1877Svaleriabarra @brief Create a non tensor-product basis for H^1 discretizations 569a8de75f0Sjeremylt 570a8de75f0Sjeremylt @param ceed A Ceed object where the CeedBasis will be created 571a8de75f0Sjeremylt @param topo Topology of element, e.g. hypercube, simplex, ect 572d1d35e2fSjeremylt @param num_comp Number of field components (1 for scalar fields) 573d1d35e2fSjeremylt @param num_nodes Total number of nodes 574d1d35e2fSjeremylt @param num_qpts Total number of quadrature points 575d1d35e2fSjeremylt @param interp Row-major (num_qpts * num_nodes) matrix expressing the values of 5768795c945Sjeremylt nodal basis functions at quadrature points 577d1d35e2fSjeremylt @param grad Row-major (num_qpts * dim * num_nodes) matrix expressing 5788795c945Sjeremylt derivatives of nodal basis functions at quadrature points 579d1d35e2fSjeremylt @param q_ref Array of length num_qpts holding the locations of quadrature 5808795c945Sjeremylt points on the reference element [-1, 1] 581d1d35e2fSjeremylt @param q_weight Array of length num_qpts holding the quadrature weights on the 582a8de75f0Sjeremylt reference element 583a8de75f0Sjeremylt @param[out] basis Address of the variable where the newly created 584a8de75f0Sjeremylt CeedBasis will be stored. 585a8de75f0Sjeremylt 586a8de75f0Sjeremylt @return An error code: 0 - success, otherwise - failure 587a8de75f0Sjeremylt 5887a982d89SJeremy L. Thompson @ref User 589a8de75f0Sjeremylt **/ 590d1d35e2fSjeremylt int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, 591d1d35e2fSjeremylt CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 592d1d35e2fSjeremylt const CeedScalar *grad, const CeedScalar *q_ref, 593d1d35e2fSjeremylt const CeedScalar *q_weight, CeedBasis *basis) { 594a8de75f0Sjeremylt int ierr; 595d1d35e2fSjeremylt CeedInt P = num_nodes, Q = num_qpts, dim = 0; 596a8de75f0Sjeremylt 5975fe0d4faSjeremylt if (!ceed->BasisCreateH1) { 5985fe0d4faSjeremylt Ceed delegate; 599aefd8378Sjeremylt ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr); 6005fe0d4faSjeremylt 6015fe0d4faSjeremylt if (!delegate) 602c042f62fSJeremy L Thompson // LCOV_EXCL_START 603e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_UNSUPPORTED, 604e15f9bd0SJeremy L Thompson "Backend does not support BasisCreateH1"); 605c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 6065fe0d4faSjeremylt 607d1d35e2fSjeremylt ierr = CeedBasisCreateH1(delegate, topo, num_comp, num_nodes, 608d1d35e2fSjeremylt num_qpts, interp, grad, q_ref, 609d1d35e2fSjeremylt q_weight, basis); CeedChk(ierr); 610e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 6115fe0d4faSjeremylt } 6125fe0d4faSjeremylt 613a8de75f0Sjeremylt ierr = CeedCalloc(1,basis); CeedChk(ierr); 614a8de75f0Sjeremylt 615a8de75f0Sjeremylt ierr = CeedBasisGetTopologyDimension(topo, &dim); CeedChk(ierr); 616a8de75f0Sjeremylt 617a8de75f0Sjeremylt (*basis)->ceed = ceed; 6189560d06aSjeremylt ierr = CeedReference(ceed); CeedChk(ierr); 619d1d35e2fSjeremylt (*basis)->ref_count = 1; 620d1d35e2fSjeremylt (*basis)->tensor_basis = 0; 621a8de75f0Sjeremylt (*basis)->dim = dim; 622d99fa3c5SJeremy L Thompson (*basis)->topo = topo; 623d1d35e2fSjeremylt (*basis)->num_comp = num_comp; 624a8de75f0Sjeremylt (*basis)->P = P; 625a8de75f0Sjeremylt (*basis)->Q = Q; 626d1d35e2fSjeremylt ierr = CeedMalloc(Q*dim,&(*basis)->q_ref_1d); CeedChk(ierr); 627d1d35e2fSjeremylt ierr = CeedMalloc(Q,&(*basis)->q_weight_1d); CeedChk(ierr); 628d1d35e2fSjeremylt memcpy((*basis)->q_ref_1d, q_ref, Q*dim*sizeof(q_ref[0])); 629d1d35e2fSjeremylt memcpy((*basis)->q_weight_1d, q_weight, Q*sizeof(q_weight[0])); 63000f91b2bSjeremylt ierr = CeedMalloc(Q*P, &(*basis)->interp); CeedChk(ierr); 63100f91b2bSjeremylt ierr = CeedMalloc(dim*Q*P, &(*basis)->grad); CeedChk(ierr); 63200f91b2bSjeremylt memcpy((*basis)->interp, interp, Q*P*sizeof(interp[0])); 63300f91b2bSjeremylt memcpy((*basis)->grad, grad, dim*Q*P*sizeof(grad[0])); 634d1d35e2fSjeremylt ierr = ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, q_ref, 635d1d35e2fSjeremylt q_weight, *basis); CeedChk(ierr); 636e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 637a8de75f0Sjeremylt } 638a8de75f0Sjeremylt 639a8de75f0Sjeremylt /** 6409560d06aSjeremylt @brief Copy the pointer to a CeedBasis. Both pointers should 6419560d06aSjeremylt be destroyed with `CeedBasisDestroy()`; 6429560d06aSjeremylt Note: If `*basis_copy` is non-NULL, then it is assumed that 6439560d06aSjeremylt `*basis_copy` is a pointer to a CeedBasis. This CeedBasis 6449560d06aSjeremylt will be destroyed if `*basis_copy` is the only 6459560d06aSjeremylt reference to this CeedBasis. 6469560d06aSjeremylt 6479560d06aSjeremylt @param basis CeedBasis to copy reference to 6489560d06aSjeremylt @param[out] basis_copy Variable to store copied reference 6499560d06aSjeremylt 6509560d06aSjeremylt @return An error code: 0 - success, otherwise - failure 6519560d06aSjeremylt 6529560d06aSjeremylt @ref User 6539560d06aSjeremylt **/ 6549560d06aSjeremylt int CeedBasisReferenceCopy(CeedBasis basis, CeedBasis *basis_copy) { 6559560d06aSjeremylt int ierr; 6569560d06aSjeremylt 6579560d06aSjeremylt ierr = CeedBasisReference(basis); CeedChk(ierr); 6589560d06aSjeremylt ierr = CeedBasisDestroy(basis_copy); CeedChk(ierr); 6599560d06aSjeremylt *basis_copy = basis; 6609560d06aSjeremylt return CEED_ERROR_SUCCESS; 6619560d06aSjeremylt } 6629560d06aSjeremylt 6639560d06aSjeremylt /** 6647a982d89SJeremy L. Thompson @brief View a CeedBasis 6657a982d89SJeremy L. Thompson 6667a982d89SJeremy L. Thompson @param basis CeedBasis to view 6677a982d89SJeremy L. Thompson @param stream Stream to view to, e.g., stdout 6687a982d89SJeremy L. Thompson 6697a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 6707a982d89SJeremy L. Thompson 6717a982d89SJeremy L. Thompson @ref User 6727a982d89SJeremy L. Thompson **/ 6737a982d89SJeremy L. Thompson int CeedBasisView(CeedBasis basis, FILE *stream) { 6747a982d89SJeremy L. Thompson int ierr; 6757a982d89SJeremy L. Thompson 676d1d35e2fSjeremylt if (basis->tensor_basis) { 677d1d35e2fSjeremylt fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P_1d, 678d1d35e2fSjeremylt basis->Q_1d); 679d1d35e2fSjeremylt ierr = CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_ref_1d, 6807a982d89SJeremy L. Thompson stream); CeedChk(ierr); 681d1d35e2fSjeremylt ierr = CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q_1d, 682d1d35e2fSjeremylt basis->q_weight_1d, stream); CeedChk(ierr); 683d1d35e2fSjeremylt ierr = CeedScalarView("interp1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, 684d1d35e2fSjeremylt basis->interp_1d, stream); CeedChk(ierr); 685d1d35e2fSjeremylt ierr = CeedScalarView("grad1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, 686d1d35e2fSjeremylt basis->grad_1d, stream); CeedChk(ierr); 6877a982d89SJeremy L. Thompson } else { 6887a982d89SJeremy L. Thompson fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P, 6897a982d89SJeremy L. Thompson basis->Q); 6907a982d89SJeremy L. Thompson ierr = CeedScalarView("qref", "\t% 12.8f", 1, basis->Q*basis->dim, 691d1d35e2fSjeremylt basis->q_ref_1d, 6927a982d89SJeremy L. Thompson stream); CeedChk(ierr); 693d1d35e2fSjeremylt ierr = CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->q_weight_1d, 6947a982d89SJeremy L. Thompson stream); CeedChk(ierr); 6957a982d89SJeremy L. Thompson ierr = CeedScalarView("interp", "\t% 12.8f", basis->Q, basis->P, 6967a982d89SJeremy L. Thompson basis->interp, stream); CeedChk(ierr); 6977a982d89SJeremy L. Thompson ierr = CeedScalarView("grad", "\t% 12.8f", basis->dim*basis->Q, basis->P, 6987a982d89SJeremy L. Thompson basis->grad, stream); CeedChk(ierr); 6997a982d89SJeremy L. Thompson } 700e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 7017a982d89SJeremy L. Thompson } 7027a982d89SJeremy L. Thompson 7037a982d89SJeremy L. Thompson /** 7047a982d89SJeremy L. Thompson @brief Apply basis evaluation from nodes to quadrature points or vice versa 7057a982d89SJeremy L. Thompson 7067a982d89SJeremy L. Thompson @param basis CeedBasis to evaluate 707d1d35e2fSjeremylt @param num_elem The number of elements to apply the basis evaluation to; 7087a982d89SJeremy L. Thompson the backend will specify the ordering in 7094cc79fe7SJed Brown CeedElemRestrictionCreateBlocked() 710d1d35e2fSjeremylt @param t_mode \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature 7117a982d89SJeremy L. Thompson points, \ref CEED_TRANSPOSE to apply the transpose, mapping 7127a982d89SJeremy L. Thompson from quadrature points to nodes 713d1d35e2fSjeremylt @param eval_mode \ref CEED_EVAL_NONE to use values directly, 7147a982d89SJeremy L. Thompson \ref CEED_EVAL_INTERP to use interpolated values, 7157a982d89SJeremy L. Thompson \ref CEED_EVAL_GRAD to use gradients, 7167a982d89SJeremy L. Thompson \ref CEED_EVAL_WEIGHT to use quadrature weights. 7177a982d89SJeremy L. Thompson @param[in] u Input CeedVector 7187a982d89SJeremy L. Thompson @param[out] v Output CeedVector 7197a982d89SJeremy L. Thompson 7207a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 7217a982d89SJeremy L. Thompson 7227a982d89SJeremy L. Thompson @ref User 7237a982d89SJeremy L. Thompson **/ 724d1d35e2fSjeremylt int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode, 725d1d35e2fSjeremylt CeedEvalMode eval_mode, CeedVector u, CeedVector v) { 7267a982d89SJeremy L. Thompson int ierr; 727d1d35e2fSjeremylt CeedInt u_length = 0, v_length, dim, num_comp, num_nodes, num_qpts; 728e15f9bd0SJeremy L Thompson ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr); 729d1d35e2fSjeremylt ierr = CeedBasisGetNumComponents(basis, &num_comp); CeedChk(ierr); 730d1d35e2fSjeremylt ierr = CeedBasisGetNumNodes(basis, &num_nodes); CeedChk(ierr); 731d1d35e2fSjeremylt ierr = CeedBasisGetNumQuadraturePoints(basis, &num_qpts); CeedChk(ierr); 732d1d35e2fSjeremylt ierr = CeedVectorGetLength(v, &v_length); CeedChk(ierr); 7337a982d89SJeremy L. Thompson if (u) { 734d1d35e2fSjeremylt ierr = CeedVectorGetLength(u, &u_length); CeedChk(ierr); 7357a982d89SJeremy L. Thompson } 7367a982d89SJeremy L. Thompson 737e15f9bd0SJeremy L Thompson if (!basis->Apply) 738e15f9bd0SJeremy L Thompson // LCOV_EXCL_START 739e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_UNSUPPORTED, 740e15f9bd0SJeremy L Thompson "Backend does not support BasisApply"); 741e15f9bd0SJeremy L Thompson // LCOV_EXCL_STOP 742e15f9bd0SJeremy L Thompson 743e15f9bd0SJeremy L Thompson // Check compatibility of topological and geometrical dimensions 744d1d35e2fSjeremylt if ((t_mode == CEED_TRANSPOSE && (v_length%num_nodes != 0 || 745d1d35e2fSjeremylt u_length%num_qpts != 0)) || 746d1d35e2fSjeremylt (t_mode == CEED_NOTRANSPOSE && (u_length%num_nodes != 0 || 747d1d35e2fSjeremylt v_length%num_qpts != 0))) 7488229195eSjeremylt // LCOV_EXCL_START 749e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_DIMENSION, 750e15f9bd0SJeremy L Thompson "Length of input/output vectors " 7517a982d89SJeremy L. Thompson "incompatible with basis dimensions"); 7528229195eSjeremylt // LCOV_EXCL_STOP 7537a982d89SJeremy L. Thompson 754e15f9bd0SJeremy L Thompson // Check vector lengths to prevent out of bounds issues 755d1d35e2fSjeremylt bool bad_dims = false; 756d1d35e2fSjeremylt switch (eval_mode) { 757e15f9bd0SJeremy L Thompson case CEED_EVAL_NONE: 758d1d35e2fSjeremylt case CEED_EVAL_INTERP: bad_dims = 759d1d35e2fSjeremylt ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts || 760d1d35e2fSjeremylt v_length < num_elem*num_comp*num_nodes)) || 761d1d35e2fSjeremylt (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp || 762d1d35e2fSjeremylt u_length < num_elem*num_comp*num_nodes))); 763e15f9bd0SJeremy L Thompson break; 764d1d35e2fSjeremylt case CEED_EVAL_GRAD: bad_dims = 765d1d35e2fSjeremylt ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts*dim || 766d1d35e2fSjeremylt v_length < num_elem*num_comp*num_nodes)) || 767d1d35e2fSjeremylt (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp*dim || 768d1d35e2fSjeremylt u_length < num_elem*num_comp*num_nodes))); 769e15f9bd0SJeremy L Thompson break; 770e15f9bd0SJeremy L Thompson case CEED_EVAL_WEIGHT: 771d1d35e2fSjeremylt bad_dims = v_length < num_elem*num_qpts; 772e15f9bd0SJeremy L Thompson break; 773e15f9bd0SJeremy L Thompson // LCOV_EXCL_START 774d1d35e2fSjeremylt case CEED_EVAL_DIV: bad_dims = 775d1d35e2fSjeremylt ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts || 776d1d35e2fSjeremylt v_length < num_elem*num_comp*num_nodes)) || 777d1d35e2fSjeremylt (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp || 778d1d35e2fSjeremylt u_length < num_elem*num_comp*num_nodes))); 779e15f9bd0SJeremy L Thompson break; 780d1d35e2fSjeremylt case CEED_EVAL_CURL: bad_dims = 781d1d35e2fSjeremylt ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts || 782d1d35e2fSjeremylt v_length < num_elem*num_comp*num_nodes)) || 783d1d35e2fSjeremylt (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp || 784d1d35e2fSjeremylt u_length < num_elem*num_comp*num_nodes))); 785e15f9bd0SJeremy L Thompson break; 786e15f9bd0SJeremy L Thompson // LCOV_EXCL_STOP 787e15f9bd0SJeremy L Thompson } 788d1d35e2fSjeremylt if (bad_dims) 789e15f9bd0SJeremy L Thompson // LCOV_EXCL_START 790e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_DIMENSION, 791d1d35e2fSjeremylt "Input/output vectors too short for basis and evaluation mode"); 792e15f9bd0SJeremy L Thompson // LCOV_EXCL_STOP 793e15f9bd0SJeremy L Thompson 794d1d35e2fSjeremylt ierr = basis->Apply(basis, num_elem, t_mode, eval_mode, u, v); CeedChk(ierr); 795e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 7967a982d89SJeremy L. Thompson } 7977a982d89SJeremy L. Thompson 7987a982d89SJeremy L. Thompson /** 7999d007619Sjeremylt @brief Get dimension for given CeedBasis 8009d007619Sjeremylt 8019d007619Sjeremylt @param basis CeedBasis 8029d007619Sjeremylt @param[out] dim Variable to store dimension of basis 8039d007619Sjeremylt 8049d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 8059d007619Sjeremylt 8069d007619Sjeremylt @ref Backend 8079d007619Sjeremylt **/ 8089d007619Sjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) { 8099d007619Sjeremylt *dim = basis->dim; 810e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 8119d007619Sjeremylt } 8129d007619Sjeremylt 8139d007619Sjeremylt /** 814d99fa3c5SJeremy L Thompson @brief Get topology for given CeedBasis 815d99fa3c5SJeremy L Thompson 816d99fa3c5SJeremy L Thompson @param basis CeedBasis 817d99fa3c5SJeremy L Thompson @param[out] topo Variable to store topology of basis 818d99fa3c5SJeremy L Thompson 819d99fa3c5SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 820d99fa3c5SJeremy L Thompson 821d99fa3c5SJeremy L Thompson @ref Backend 822d99fa3c5SJeremy L Thompson **/ 823d99fa3c5SJeremy L Thompson int CeedBasisGetTopology(CeedBasis basis, CeedElemTopology *topo) { 824d99fa3c5SJeremy L Thompson *topo = basis->topo; 825e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 826d99fa3c5SJeremy L Thompson } 827d99fa3c5SJeremy L Thompson 828d99fa3c5SJeremy L Thompson /** 8299d007619Sjeremylt @brief Get number of components for given CeedBasis 8309d007619Sjeremylt 8319d007619Sjeremylt @param basis CeedBasis 832d1d35e2fSjeremylt @param[out] num_comp Variable to store number of components of basis 8339d007619Sjeremylt 8349d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 8359d007619Sjeremylt 8369d007619Sjeremylt @ref Backend 8379d007619Sjeremylt **/ 838d1d35e2fSjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *num_comp) { 839d1d35e2fSjeremylt *num_comp = basis->num_comp; 840e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 8419d007619Sjeremylt } 8429d007619Sjeremylt 8439d007619Sjeremylt /** 8449d007619Sjeremylt @brief Get total number of nodes (in dim dimensions) of a CeedBasis 8459d007619Sjeremylt 8469d007619Sjeremylt @param basis CeedBasis 8479d007619Sjeremylt @param[out] P Variable to store number of nodes 8489d007619Sjeremylt 8499d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 8509d007619Sjeremylt 8519d007619Sjeremylt @ref Utility 8529d007619Sjeremylt **/ 8539d007619Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) { 8549d007619Sjeremylt *P = basis->P; 855e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 8569d007619Sjeremylt } 8579d007619Sjeremylt 8589d007619Sjeremylt /** 8599d007619Sjeremylt @brief Get total number of nodes (in 1 dimension) of a CeedBasis 8609d007619Sjeremylt 8619d007619Sjeremylt @param basis CeedBasis 862d1d35e2fSjeremylt @param[out] P_1d Variable to store number of nodes 8639d007619Sjeremylt 8649d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 8659d007619Sjeremylt 8669d007619Sjeremylt @ref Backend 8679d007619Sjeremylt **/ 868d1d35e2fSjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P_1d) { 869d1d35e2fSjeremylt if (!basis->tensor_basis) 8709d007619Sjeremylt // LCOV_EXCL_START 871e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_MINOR, 872d1d35e2fSjeremylt "Cannot supply P_1d for non-tensor basis"); 8739d007619Sjeremylt // LCOV_EXCL_STOP 8749d007619Sjeremylt 875d1d35e2fSjeremylt *P_1d = basis->P_1d; 876e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 8779d007619Sjeremylt } 8789d007619Sjeremylt 8799d007619Sjeremylt /** 8809d007619Sjeremylt @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis 8819d007619Sjeremylt 8829d007619Sjeremylt @param basis CeedBasis 8839d007619Sjeremylt @param[out] Q Variable to store number of quadrature points 8849d007619Sjeremylt 8859d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 8869d007619Sjeremylt 8879d007619Sjeremylt @ref Utility 8889d007619Sjeremylt **/ 8899d007619Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) { 8909d007619Sjeremylt *Q = basis->Q; 891e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 8929d007619Sjeremylt } 8939d007619Sjeremylt 8949d007619Sjeremylt /** 8959d007619Sjeremylt @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis 8969d007619Sjeremylt 8979d007619Sjeremylt @param basis CeedBasis 898d1d35e2fSjeremylt @param[out] Q_1d Variable to store number of quadrature points 8999d007619Sjeremylt 9009d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 9019d007619Sjeremylt 9029d007619Sjeremylt @ref Backend 9039d007619Sjeremylt **/ 904d1d35e2fSjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q_1d) { 905d1d35e2fSjeremylt if (!basis->tensor_basis) 9069d007619Sjeremylt // LCOV_EXCL_START 907e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_MINOR, 908d1d35e2fSjeremylt "Cannot supply Q_1d for non-tensor basis"); 9099d007619Sjeremylt // LCOV_EXCL_STOP 9109d007619Sjeremylt 911d1d35e2fSjeremylt *Q_1d = basis->Q_1d; 912e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 9139d007619Sjeremylt } 9149d007619Sjeremylt 9159d007619Sjeremylt /** 9169d007619Sjeremylt @brief Get reference coordinates of quadrature points (in dim dimensions) 9179d007619Sjeremylt of a CeedBasis 9189d007619Sjeremylt 9199d007619Sjeremylt @param basis CeedBasis 920d1d35e2fSjeremylt @param[out] q_ref Variable to store reference coordinates of quadrature points 9219d007619Sjeremylt 9229d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 9239d007619Sjeremylt 9249d007619Sjeremylt @ref Backend 9259d007619Sjeremylt **/ 926d1d35e2fSjeremylt int CeedBasisGetQRef(CeedBasis basis, const CeedScalar **q_ref) { 927d1d35e2fSjeremylt *q_ref = basis->q_ref_1d; 928e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 9299d007619Sjeremylt } 9309d007619Sjeremylt 9319d007619Sjeremylt /** 9329d007619Sjeremylt @brief Get quadrature weights of quadrature points (in dim dimensions) 9339d007619Sjeremylt of a CeedBasis 9349d007619Sjeremylt 9359d007619Sjeremylt @param basis CeedBasis 936d1d35e2fSjeremylt @param[out] q_weight Variable to store quadrature weights 9379d007619Sjeremylt 9389d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 9399d007619Sjeremylt 9409d007619Sjeremylt @ref Backend 9419d007619Sjeremylt **/ 942d1d35e2fSjeremylt int CeedBasisGetQWeights(CeedBasis basis, const CeedScalar **q_weight) { 943d1d35e2fSjeremylt *q_weight = basis->q_weight_1d; 944e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 9459d007619Sjeremylt } 9469d007619Sjeremylt 9479d007619Sjeremylt /** 9489d007619Sjeremylt @brief Get interpolation matrix of a CeedBasis 9499d007619Sjeremylt 9509d007619Sjeremylt @param basis CeedBasis 9519d007619Sjeremylt @param[out] interp Variable to store interpolation matrix 9529d007619Sjeremylt 9539d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 9549d007619Sjeremylt 9559d007619Sjeremylt @ref Backend 9569d007619Sjeremylt **/ 9576c58de82SJeremy L Thompson int CeedBasisGetInterp(CeedBasis basis, const CeedScalar **interp) { 958d1d35e2fSjeremylt if (!basis->interp && basis->tensor_basis) { 9599d007619Sjeremylt // Allocate 9609d007619Sjeremylt int ierr; 9619d007619Sjeremylt ierr = CeedMalloc(basis->Q*basis->P, &basis->interp); CeedChk(ierr); 9629d007619Sjeremylt 9639d007619Sjeremylt // Initialize 9649d007619Sjeremylt for (CeedInt i=0; i<basis->Q*basis->P; i++) 9659d007619Sjeremylt basis->interp[i] = 1.0; 9669d007619Sjeremylt 9679d007619Sjeremylt // Calculate 9689d007619Sjeremylt for (CeedInt d=0; d<basis->dim; d++) 9699d007619Sjeremylt for (CeedInt qpt=0; qpt<basis->Q; qpt++) 9709d007619Sjeremylt for (CeedInt node=0; node<basis->P; node++) { 971d1d35e2fSjeremylt CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d; 972d1d35e2fSjeremylt CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d; 973d1d35e2fSjeremylt basis->interp[qpt*(basis->P)+node] *= basis->interp_1d[q*basis->P_1d+p]; 9749d007619Sjeremylt } 9759d007619Sjeremylt } 9769d007619Sjeremylt *interp = basis->interp; 977e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 9789d007619Sjeremylt } 9799d007619Sjeremylt 9809d007619Sjeremylt /** 9819d007619Sjeremylt @brief Get 1D interpolation matrix of a tensor product CeedBasis 9829d007619Sjeremylt 9839d007619Sjeremylt @param basis CeedBasis 984d1d35e2fSjeremylt @param[out] interp_1d Variable to store interpolation matrix 9859d007619Sjeremylt 9869d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 9879d007619Sjeremylt 9889d007619Sjeremylt @ref Backend 9899d007619Sjeremylt **/ 990d1d35e2fSjeremylt int CeedBasisGetInterp1D(CeedBasis basis, const CeedScalar **interp_1d) { 991d1d35e2fSjeremylt if (!basis->tensor_basis) 9929d007619Sjeremylt // LCOV_EXCL_START 993e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_MINOR, 994e15f9bd0SJeremy L Thompson "CeedBasis is not a tensor product basis."); 9959d007619Sjeremylt // LCOV_EXCL_STOP 9969d007619Sjeremylt 997d1d35e2fSjeremylt *interp_1d = basis->interp_1d; 998e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 9999d007619Sjeremylt } 10009d007619Sjeremylt 10019d007619Sjeremylt /** 10029d007619Sjeremylt @brief Get gradient matrix of a CeedBasis 10039d007619Sjeremylt 10049d007619Sjeremylt @param basis CeedBasis 10059d007619Sjeremylt @param[out] grad Variable to store gradient matrix 10069d007619Sjeremylt 10079d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 10089d007619Sjeremylt 10099d007619Sjeremylt @ref Backend 10109d007619Sjeremylt **/ 10116c58de82SJeremy L Thompson int CeedBasisGetGrad(CeedBasis basis, const CeedScalar **grad) { 1012d1d35e2fSjeremylt if (!basis->grad && basis->tensor_basis) { 10139d007619Sjeremylt // Allocate 10149d007619Sjeremylt int ierr; 10159d007619Sjeremylt ierr = CeedMalloc(basis->dim*basis->Q*basis->P, &basis->grad); 10169d007619Sjeremylt CeedChk(ierr); 10179d007619Sjeremylt 10189d007619Sjeremylt // Initialize 10199d007619Sjeremylt for (CeedInt i=0; i<basis->dim*basis->Q*basis->P; i++) 10209d007619Sjeremylt basis->grad[i] = 1.0; 10219d007619Sjeremylt 10229d007619Sjeremylt // Calculate 10239d007619Sjeremylt for (CeedInt d=0; d<basis->dim; d++) 10249d007619Sjeremylt for (CeedInt i=0; i<basis->dim; i++) 10259d007619Sjeremylt for (CeedInt qpt=0; qpt<basis->Q; qpt++) 10269d007619Sjeremylt for (CeedInt node=0; node<basis->P; node++) { 1027d1d35e2fSjeremylt CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d; 1028d1d35e2fSjeremylt CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d; 10299d007619Sjeremylt if (i == d) 10309d007619Sjeremylt basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *= 1031d1d35e2fSjeremylt basis->grad_1d[q*basis->P_1d+p]; 10329d007619Sjeremylt else 10339d007619Sjeremylt basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *= 1034d1d35e2fSjeremylt basis->interp_1d[q*basis->P_1d+p]; 10359d007619Sjeremylt } 10369d007619Sjeremylt } 10379d007619Sjeremylt *grad = basis->grad; 1038e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 10399d007619Sjeremylt } 10409d007619Sjeremylt 10419d007619Sjeremylt /** 10429d007619Sjeremylt @brief Get 1D gradient matrix of a tensor product CeedBasis 10439d007619Sjeremylt 10449d007619Sjeremylt @param basis CeedBasis 1045d1d35e2fSjeremylt @param[out] grad_1d Variable to store gradient matrix 10469d007619Sjeremylt 10479d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 10489d007619Sjeremylt 10499d007619Sjeremylt @ref Backend 10509d007619Sjeremylt **/ 1051d1d35e2fSjeremylt int CeedBasisGetGrad1D(CeedBasis basis, const CeedScalar **grad_1d) { 1052d1d35e2fSjeremylt if (!basis->tensor_basis) 10539d007619Sjeremylt // LCOV_EXCL_START 1054e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_MINOR, 1055e15f9bd0SJeremy L Thompson "CeedBasis is not a tensor product basis."); 10569d007619Sjeremylt // LCOV_EXCL_STOP 10579d007619Sjeremylt 1058d1d35e2fSjeremylt *grad_1d = basis->grad_1d; 1059e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 10609d007619Sjeremylt } 10619d007619Sjeremylt 10629d007619Sjeremylt /** 10637a982d89SJeremy L. Thompson @brief Destroy a CeedBasis 10647a982d89SJeremy L. Thompson 10657a982d89SJeremy L. Thompson @param basis CeedBasis to destroy 10667a982d89SJeremy L. Thompson 10677a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 10687a982d89SJeremy L. Thompson 10697a982d89SJeremy L. Thompson @ref User 10707a982d89SJeremy L. Thompson **/ 10717a982d89SJeremy L. Thompson int CeedBasisDestroy(CeedBasis *basis) { 10727a982d89SJeremy L. Thompson int ierr; 10737a982d89SJeremy L. Thompson 1074d1d35e2fSjeremylt if (!*basis || --(*basis)->ref_count > 0) return CEED_ERROR_SUCCESS; 10757a982d89SJeremy L. Thompson if ((*basis)->Destroy) { 10767a982d89SJeremy L. Thompson ierr = (*basis)->Destroy(*basis); CeedChk(ierr); 10777a982d89SJeremy L. Thompson } 107834359f16Sjeremylt if ((*basis)->contract) { 107934359f16Sjeremylt ierr = CeedTensorContractDestroy(&(*basis)->contract); CeedChk(ierr); 108034359f16Sjeremylt } 10817a982d89SJeremy L. Thompson ierr = CeedFree(&(*basis)->interp); CeedChk(ierr); 1082d1d35e2fSjeremylt ierr = CeedFree(&(*basis)->interp_1d); CeedChk(ierr); 10837a982d89SJeremy L. Thompson ierr = CeedFree(&(*basis)->grad); CeedChk(ierr); 1084d1d35e2fSjeremylt ierr = CeedFree(&(*basis)->grad_1d); CeedChk(ierr); 1085d1d35e2fSjeremylt ierr = CeedFree(&(*basis)->q_ref_1d); CeedChk(ierr); 1086d1d35e2fSjeremylt ierr = CeedFree(&(*basis)->q_weight_1d); CeedChk(ierr); 10877a982d89SJeremy L. Thompson ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr); 10887a982d89SJeremy L. Thompson ierr = CeedFree(basis); CeedChk(ierr); 1089e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 10907a982d89SJeremy L. Thompson } 10917a982d89SJeremy L. Thompson 10927a982d89SJeremy L. Thompson /** 1093b11c1e72Sjeremylt @brief Construct a Gauss-Legendre quadrature 1094b11c1e72Sjeremylt 1095b11c1e72Sjeremylt @param Q Number of quadrature points (integrates polynomials of 1096b11c1e72Sjeremylt degree 2*Q-1 exactly) 1097d1d35e2fSjeremylt @param[out] q_ref_1d Array of length Q to hold the abscissa on [-1, 1] 1098d1d35e2fSjeremylt @param[out] q_weight_1d Array of length Q to hold the weights 1099b11c1e72Sjeremylt 1100b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 1101dfdf5a53Sjeremylt 1102dfdf5a53Sjeremylt @ref Utility 1103b11c1e72Sjeremylt **/ 1104d1d35e2fSjeremylt int CeedGaussQuadrature(CeedInt Q, CeedScalar *q_ref_1d, 1105d1d35e2fSjeremylt CeedScalar *q_weight_1d) { 1106d7b241e6Sjeremylt // Allocate 1107d7b241e6Sjeremylt CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0*atan(1.0); 1108d1d35e2fSjeremylt // Build q_ref_1d, q_weight_1d 1109d7b241e6Sjeremylt for (int i = 0; i <= Q/2; i++) { 1110d7b241e6Sjeremylt // Guess 1111d7b241e6Sjeremylt xi = cos(PI*(CeedScalar)(2*i+1)/((CeedScalar)(2*Q))); 1112d7b241e6Sjeremylt // Pn(xi) 1113d7b241e6Sjeremylt P0 = 1.0; 1114d7b241e6Sjeremylt P1 = xi; 1115d7b241e6Sjeremylt P2 = 0.0; 1116d7b241e6Sjeremylt for (int j = 2; j <= Q; j++) { 1117d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 1118d7b241e6Sjeremylt P0 = P1; 1119d7b241e6Sjeremylt P1 = P2; 1120d7b241e6Sjeremylt } 1121d7b241e6Sjeremylt // First Newton Step 1122d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 1123d7b241e6Sjeremylt xi = xi-P2/dP2; 1124d7b241e6Sjeremylt // Newton to convergence 11250e4d4210Sjeremylt for (int k=0; k<100 && fabs(P2)>10*CEED_EPSILON; k++) { 1126d7b241e6Sjeremylt P0 = 1.0; 1127d7b241e6Sjeremylt P1 = xi; 1128d7b241e6Sjeremylt for (int j = 2; j <= Q; j++) { 1129d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 1130d7b241e6Sjeremylt P0 = P1; 1131d7b241e6Sjeremylt P1 = P2; 1132d7b241e6Sjeremylt } 1133d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 1134d7b241e6Sjeremylt xi = xi-P2/dP2; 1135d7b241e6Sjeremylt } 1136d7b241e6Sjeremylt // Save xi, wi 1137d7b241e6Sjeremylt wi = 2.0/((1.0-xi*xi)*dP2*dP2); 1138d1d35e2fSjeremylt q_weight_1d[i] = wi; 1139d1d35e2fSjeremylt q_weight_1d[Q-1-i] = wi; 1140d1d35e2fSjeremylt q_ref_1d[i] = -xi; 1141d1d35e2fSjeremylt q_ref_1d[Q-1-i]= xi; 1142d7b241e6Sjeremylt } 1143e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1144d7b241e6Sjeremylt } 1145d7b241e6Sjeremylt 1146b11c1e72Sjeremylt /** 1147b11c1e72Sjeremylt @brief Construct a Gauss-Legendre-Lobatto quadrature 1148b11c1e72Sjeremylt 1149b11c1e72Sjeremylt @param Q Number of quadrature points (integrates polynomials of 1150b11c1e72Sjeremylt degree 2*Q-3 exactly) 1151d1d35e2fSjeremylt @param[out] q_ref_1d Array of length Q to hold the abscissa on [-1, 1] 1152d1d35e2fSjeremylt @param[out] q_weight_1d Array of length Q to hold the weights 1153b11c1e72Sjeremylt 1154b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 1155dfdf5a53Sjeremylt 1156dfdf5a53Sjeremylt @ref Utility 1157b11c1e72Sjeremylt **/ 1158d1d35e2fSjeremylt int CeedLobattoQuadrature(CeedInt Q, CeedScalar *q_ref_1d, 1159d1d35e2fSjeremylt CeedScalar *q_weight_1d) { 1160d7b241e6Sjeremylt // Allocate 1161d7b241e6Sjeremylt CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0*atan(1.0); 1162d1d35e2fSjeremylt // Build q_ref_1d, q_weight_1d 1163d7b241e6Sjeremylt // Set endpoints 116430a100c3SJed Brown if (Q < 2) 1165b0d62198Sjeremylt // LCOV_EXCL_START 1166e15f9bd0SJeremy L Thompson return CeedError(NULL, CEED_ERROR_DIMENSION, 11677ed52d01Sjeremylt "Cannot create Lobatto quadrature with Q=%d < 2 points", Q); 1168b0d62198Sjeremylt // LCOV_EXCL_STOP 1169d7b241e6Sjeremylt wi = 2.0/((CeedScalar)(Q*(Q-1))); 1170d1d35e2fSjeremylt if (q_weight_1d) { 1171d1d35e2fSjeremylt q_weight_1d[0] = wi; 1172d1d35e2fSjeremylt q_weight_1d[Q-1] = wi; 1173d7b241e6Sjeremylt } 1174d1d35e2fSjeremylt q_ref_1d[0] = -1.0; 1175d1d35e2fSjeremylt q_ref_1d[Q-1] = 1.0; 1176d7b241e6Sjeremylt // Interior 1177d7b241e6Sjeremylt for (int i = 1; i <= (Q-1)/2; i++) { 1178d7b241e6Sjeremylt // Guess 1179d7b241e6Sjeremylt xi = cos(PI*(CeedScalar)(i)/(CeedScalar)(Q-1)); 1180d7b241e6Sjeremylt // Pn(xi) 1181d7b241e6Sjeremylt P0 = 1.0; 1182d7b241e6Sjeremylt P1 = xi; 1183d7b241e6Sjeremylt P2 = 0.0; 1184d7b241e6Sjeremylt for (int j = 2; j < Q; j++) { 1185d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 1186d7b241e6Sjeremylt P0 = P1; 1187d7b241e6Sjeremylt P1 = P2; 1188d7b241e6Sjeremylt } 1189d7b241e6Sjeremylt // First Newton step 1190d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 1191d7b241e6Sjeremylt d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi); 1192d7b241e6Sjeremylt xi = xi-dP2/d2P2; 1193d7b241e6Sjeremylt // Newton to convergence 11940e4d4210Sjeremylt for (int k=0; k<100 && fabs(dP2)>10*CEED_EPSILON; k++) { 1195d7b241e6Sjeremylt P0 = 1.0; 1196d7b241e6Sjeremylt P1 = xi; 1197d7b241e6Sjeremylt for (int j = 2; j < Q; j++) { 1198d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 1199d7b241e6Sjeremylt P0 = P1; 1200d7b241e6Sjeremylt P1 = P2; 1201d7b241e6Sjeremylt } 1202d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 1203d7b241e6Sjeremylt d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi); 1204d7b241e6Sjeremylt xi = xi-dP2/d2P2; 1205d7b241e6Sjeremylt } 1206d7b241e6Sjeremylt // Save xi, wi 1207d7b241e6Sjeremylt wi = 2.0/(((CeedScalar)(Q*(Q-1)))*P2*P2); 1208d1d35e2fSjeremylt if (q_weight_1d) { 1209d1d35e2fSjeremylt q_weight_1d[i] = wi; 1210d1d35e2fSjeremylt q_weight_1d[Q-1-i] = wi; 1211d7b241e6Sjeremylt } 1212d1d35e2fSjeremylt q_ref_1d[i] = -xi; 1213d1d35e2fSjeremylt q_ref_1d[Q-1-i]= xi; 1214d7b241e6Sjeremylt } 1215e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1216d7b241e6Sjeremylt } 1217d7b241e6Sjeremylt 1218dfdf5a53Sjeremylt /** 121995bb1877Svaleriabarra @brief Return QR Factorization of a matrix 1220b11c1e72Sjeremylt 122177645d7bSjeremylt @param ceed A Ceed context for error handling 122252bfb9bbSJeremy L Thompson @param[in,out] mat Row-major matrix to be factorized in place 122352bfb9bbSJeremy L Thompson @param[in,out] tau Vector of length m of scaling factors 1224b11c1e72Sjeremylt @param m Number of rows 1225b11c1e72Sjeremylt @param n Number of columns 1226b11c1e72Sjeremylt 1227b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 1228dfdf5a53Sjeremylt 1229dfdf5a53Sjeremylt @ref Utility 1230b11c1e72Sjeremylt **/ 1231a7bd39daSjeremylt int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau, 1232d7b241e6Sjeremylt CeedInt m, CeedInt n) { 1233d7b241e6Sjeremylt CeedScalar v[m]; 1234d7b241e6Sjeremylt 1235a7bd39daSjeremylt // Check m >= n 1236a7bd39daSjeremylt if (n > m) 1237c042f62fSJeremy L Thompson // LCOV_EXCL_START 1238e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_UNSUPPORTED, 1239e15f9bd0SJeremy L Thompson "Cannot compute QR factorization with n > m"); 1240c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 1241a7bd39daSjeremylt 124252bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) { 1243bde37e8cSJed Brown if (i >= m-1) { // last row of matrix, no reflection needed 1244bde37e8cSJed Brown tau[i] = 0.; 1245bde37e8cSJed Brown break; 1246bde37e8cSJed Brown } 1247d7b241e6Sjeremylt // Calculate Householder vector, magnitude 1248d7b241e6Sjeremylt CeedScalar sigma = 0.0; 1249d7b241e6Sjeremylt v[i] = mat[i+n*i]; 125052bfb9bbSJeremy L Thompson for (CeedInt j=i+1; j<m; j++) { 1251d7b241e6Sjeremylt v[j] = mat[i+n*j]; 1252d7b241e6Sjeremylt sigma += v[j] * v[j]; 1253d7b241e6Sjeremylt } 1254d7b241e6Sjeremylt CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m] 1255673160d7Sjeremylt CeedScalar R_ii = -copysign(norm, v[i]); 1256673160d7Sjeremylt v[i] -= R_ii; 1257d7b241e6Sjeremylt // norm of v[i:m] after modification above and scaling below 1258d7b241e6Sjeremylt // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 1259d7b241e6Sjeremylt // tau = 2 / (norm*norm) 1260d7b241e6Sjeremylt tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma); 12611d102b48SJeremy L Thompson for (CeedInt j=i+1; j<m; j++) 12621d102b48SJeremy L Thompson v[j] /= v[i]; 1263d7b241e6Sjeremylt 1264d7b241e6Sjeremylt // Apply Householder reflector to lower right panel 1265d7b241e6Sjeremylt CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1); 1266d7b241e6Sjeremylt // Save v 1267673160d7Sjeremylt mat[i+n*i] = R_ii; 12681d102b48SJeremy L Thompson for (CeedInt j=i+1; j<m; j++) 1269d7b241e6Sjeremylt mat[i+n*j] = v[j]; 1270d7b241e6Sjeremylt } 1271e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1272d7b241e6Sjeremylt } 1273d7b241e6Sjeremylt 1274b11c1e72Sjeremylt /** 127552bfb9bbSJeremy L Thompson @brief Return symmetric Schur decomposition of the symmetric matrix mat via 127652bfb9bbSJeremy L Thompson symmetric QR factorization 127752bfb9bbSJeremy L Thompson 127877645d7bSjeremylt @param ceed A Ceed context for error handling 127952bfb9bbSJeremy L Thompson @param[in,out] mat Row-major matrix to be factorized in place 1280460bf743SValeria Barra @param[out] lambda Vector of length n of eigenvalues 128152bfb9bbSJeremy L Thompson @param n Number of rows/columns 128252bfb9bbSJeremy L Thompson 128352bfb9bbSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 128452bfb9bbSJeremy L Thompson 128552bfb9bbSJeremy L Thompson @ref Utility 128652bfb9bbSJeremy L Thompson **/ 128703d18186Sjeremylt CeedPragmaOptimizeOff 128852bfb9bbSJeremy L Thompson int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat, 128952bfb9bbSJeremy L Thompson CeedScalar *lambda, CeedInt n) { 129052bfb9bbSJeremy L Thompson // Check bounds for clang-tidy 129152bfb9bbSJeremy L Thompson if (n<2) 1292c042f62fSJeremy L Thompson // LCOV_EXCL_START 1293e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_UNSUPPORTED, 1294c042f62fSJeremy L Thompson "Cannot compute symmetric Schur decomposition of scalars"); 1295c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 129652bfb9bbSJeremy L Thompson 1297673160d7Sjeremylt CeedScalar v[n-1], tau[n-1], mat_T[n*n]; 129852bfb9bbSJeremy L Thompson 1299673160d7Sjeremylt // Copy mat to mat_T and set mat to I 1300673160d7Sjeremylt memcpy(mat_T, mat, n*n*sizeof(mat[0])); 130152bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) 130252bfb9bbSJeremy L Thompson for (CeedInt j=0; j<n; j++) 130352bfb9bbSJeremy L Thompson mat[j+n*i] = (i==j) ? 1 : 0; 130452bfb9bbSJeremy L Thompson 130552bfb9bbSJeremy L Thompson // Reduce to tridiagonal 130652bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n-1; i++) { 130752bfb9bbSJeremy L Thompson // Calculate Householder vector, magnitude 130852bfb9bbSJeremy L Thompson CeedScalar sigma = 0.0; 1309673160d7Sjeremylt v[i] = mat_T[i+n*(i+1)]; 131052bfb9bbSJeremy L Thompson for (CeedInt j=i+1; j<n-1; j++) { 1311673160d7Sjeremylt v[j] = mat_T[i+n*(j+1)]; 131252bfb9bbSJeremy L Thompson sigma += v[j] * v[j]; 131352bfb9bbSJeremy L Thompson } 131452bfb9bbSJeremy L Thompson CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:n-1] 1315673160d7Sjeremylt CeedScalar R_ii = -copysign(norm, v[i]); 1316673160d7Sjeremylt v[i] -= R_ii; 131752bfb9bbSJeremy L Thompson // norm of v[i:m] after modification above and scaling below 131852bfb9bbSJeremy L Thompson // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 131952bfb9bbSJeremy L Thompson // tau = 2 / (norm*norm) 1320*80a9ef05SNatalie Beams tau[i] = i == n - 2 ? 2 : 2 * v[i]*v[i] / (v[i]*v[i] + sigma); 1321fb551037Sjeremylt for (CeedInt j=i+1; j<n-1; j++) 1322fb551037Sjeremylt v[j] /= v[i]; 132352bfb9bbSJeremy L Thompson 132452bfb9bbSJeremy L Thompson // Update sub and super diagonal 132552bfb9bbSJeremy L Thompson for (CeedInt j=i+2; j<n; j++) { 1326673160d7Sjeremylt mat_T[i+n*j] = 0; mat_T[j+n*i] = 0; 132752bfb9bbSJeremy L Thompson } 132852bfb9bbSJeremy L Thompson // Apply symmetric Householder reflector to lower right panel 1329673160d7Sjeremylt CeedHouseholderReflect(&mat_T[(i+1)+n*(i+1)], &v[i], tau[i], 133052bfb9bbSJeremy L Thompson n-(i+1), n-(i+1), n, 1); 1331673160d7Sjeremylt CeedHouseholderReflect(&mat_T[(i+1)+n*(i+1)], &v[i], tau[i], 133252bfb9bbSJeremy L Thompson n-(i+1), n-(i+1), 1, n); 1333673160d7Sjeremylt 133452bfb9bbSJeremy L Thompson // Save v 1335673160d7Sjeremylt mat_T[i+n*(i+1)] = R_ii; 1336673160d7Sjeremylt mat_T[(i+1)+n*i] = R_ii; 133752bfb9bbSJeremy L Thompson for (CeedInt j=i+1; j<n-1; j++) { 1338673160d7Sjeremylt mat_T[i+n*(j+1)] = v[j]; 133952bfb9bbSJeremy L Thompson } 134052bfb9bbSJeremy L Thompson } 134152bfb9bbSJeremy L Thompson // Backwards accumulation of Q 134252bfb9bbSJeremy L Thompson for (CeedInt i=n-2; i>=0; i--) { 134385cf89eaSjeremylt if (tau[i] > 0.0) { 134452bfb9bbSJeremy L Thompson v[i] = 1; 134552bfb9bbSJeremy L Thompson for (CeedInt j=i+1; j<n-1; j++) { 1346673160d7Sjeremylt v[j] = mat_T[i+n*(j+1)]; 1347673160d7Sjeremylt mat_T[i+n*(j+1)] = 0; 134852bfb9bbSJeremy L Thompson } 134952bfb9bbSJeremy L Thompson CeedHouseholderReflect(&mat[(i+1)+n*(i+1)], &v[i], tau[i], 135052bfb9bbSJeremy L Thompson n-(i+1), n-(i+1), n, 1); 135152bfb9bbSJeremy L Thompson } 135285cf89eaSjeremylt } 135352bfb9bbSJeremy L Thompson 135452bfb9bbSJeremy L Thompson // Reduce sub and super diagonal 1355673160d7Sjeremylt CeedInt p = 0, q = 0, itr = 0, max_itr = n*n*n*n; 1356673160d7Sjeremylt CeedScalar tol = CEED_EPSILON; 135752bfb9bbSJeremy L Thompson 1358673160d7Sjeremylt while (itr < max_itr) { 135952bfb9bbSJeremy L Thompson // Update p, q, size of reduced portions of diagonal 136052bfb9bbSJeremy L Thompson p = 0; q = 0; 136152bfb9bbSJeremy L Thompson for (CeedInt i=n-2; i>=0; i--) { 1362673160d7Sjeremylt if (fabs(mat_T[i+n*(i+1)]) < tol) 136352bfb9bbSJeremy L Thompson q += 1; 136452bfb9bbSJeremy L Thompson else 136552bfb9bbSJeremy L Thompson break; 136652bfb9bbSJeremy L Thompson } 1367673160d7Sjeremylt for (CeedInt i=0; i<n-q-1; i++) { 1368673160d7Sjeremylt if (fabs(mat_T[i+n*(i+1)]) < tol) 136952bfb9bbSJeremy L Thompson p += 1; 137052bfb9bbSJeremy L Thompson else 137152bfb9bbSJeremy L Thompson break; 137252bfb9bbSJeremy L Thompson } 137352bfb9bbSJeremy L Thompson if (q == n-1) break; // Finished reducing 137452bfb9bbSJeremy L Thompson 137552bfb9bbSJeremy L Thompson // Reduce tridiagonal portion 1376673160d7Sjeremylt CeedScalar t_nn = mat_T[(n-1-q)+n*(n-1-q)], 1377673160d7Sjeremylt t_nnm1 = mat_T[(n-2-q)+n*(n-1-q)]; 1378673160d7Sjeremylt CeedScalar d = (mat_T[(n-2-q)+n*(n-2-q)] - t_nn)/2; 1379673160d7Sjeremylt CeedScalar mu = t_nn - t_nnm1*t_nnm1 / 1380673160d7Sjeremylt (d + copysign(sqrt(d*d + t_nnm1*t_nnm1), d)); 1381673160d7Sjeremylt CeedScalar x = mat_T[p+n*p] - mu; 1382673160d7Sjeremylt CeedScalar z = mat_T[p+n*(p+1)]; 1383673160d7Sjeremylt for (CeedInt k=p; k<n-q-1; k++) { 138452bfb9bbSJeremy L Thompson // Compute Givens rotation 138552bfb9bbSJeremy L Thompson CeedScalar c = 1, s = 0; 138652bfb9bbSJeremy L Thompson if (fabs(z) > tol) { 138752bfb9bbSJeremy L Thompson if (fabs(z) > fabs(x)) { 138852bfb9bbSJeremy L Thompson CeedScalar tau = -x/z; 138952bfb9bbSJeremy L Thompson s = 1/sqrt(1+tau*tau), c = s*tau; 139052bfb9bbSJeremy L Thompson } else { 139152bfb9bbSJeremy L Thompson CeedScalar tau = -z/x; 139252bfb9bbSJeremy L Thompson c = 1/sqrt(1+tau*tau), s = c*tau; 139352bfb9bbSJeremy L Thompson } 139452bfb9bbSJeremy L Thompson } 139552bfb9bbSJeremy L Thompson 139652bfb9bbSJeremy L Thompson // Apply Givens rotation to T 1397673160d7Sjeremylt CeedGivensRotation(mat_T, c, s, CEED_NOTRANSPOSE, k, k+1, n, n); 1398673160d7Sjeremylt CeedGivensRotation(mat_T, c, s, CEED_TRANSPOSE, k, k+1, n, n); 139952bfb9bbSJeremy L Thompson 140052bfb9bbSJeremy L Thompson // Apply Givens rotation to Q 140152bfb9bbSJeremy L Thompson CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k+1, n, n); 140252bfb9bbSJeremy L Thompson 140352bfb9bbSJeremy L Thompson // Update x, z 140452bfb9bbSJeremy L Thompson if (k < n-q-2) { 1405673160d7Sjeremylt x = mat_T[k+n*(k+1)]; 1406673160d7Sjeremylt z = mat_T[k+n*(k+2)]; 140752bfb9bbSJeremy L Thompson } 140852bfb9bbSJeremy L Thompson } 140952bfb9bbSJeremy L Thompson itr++; 141052bfb9bbSJeremy L Thompson } 1411673160d7Sjeremylt 141252bfb9bbSJeremy L Thompson // Save eigenvalues 141352bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) 1414673160d7Sjeremylt lambda[i] = mat_T[i+n*i]; 141552bfb9bbSJeremy L Thompson 141652bfb9bbSJeremy L Thompson // Check convergence 1417673160d7Sjeremylt if (itr == max_itr && q < n-1) 1418c042f62fSJeremy L Thompson // LCOV_EXCL_START 1419e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_MINOR, 1420e15f9bd0SJeremy L Thompson "Symmetric QR failed to converge"); 1421c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 1422e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 142352bfb9bbSJeremy L Thompson } 142403d18186Sjeremylt CeedPragmaOptimizeOn 142552bfb9bbSJeremy L Thompson 142652bfb9bbSJeremy L Thompson /** 142752bfb9bbSJeremy L Thompson @brief Return Simultaneous Diagonalization of two matrices. This solves the 142852bfb9bbSJeremy L Thompson generalized eigenvalue problem A x = lambda B x, where A and B 142952bfb9bbSJeremy L Thompson are symmetric and B is positive definite. We generate the matrix X 143052bfb9bbSJeremy L Thompson and vector Lambda such that X^T A X = Lambda and X^T B X = I. This 143152bfb9bbSJeremy L Thompson is equivalent to the LAPACK routine 'sygv' with TYPE = 1. 143252bfb9bbSJeremy L Thompson 143377645d7bSjeremylt @param ceed A Ceed context for error handling 1434d1d35e2fSjeremylt @param[in] mat_A Row-major matrix to be factorized with eigenvalues 1435d1d35e2fSjeremylt @param[in] mat_B Row-major matrix to be factorized to identity 1436d3331725Sjeremylt @param[out] mat_X Row-major orthogonal matrix 1437460bf743SValeria Barra @param[out] lambda Vector of length n of generalized eigenvalues 143852bfb9bbSJeremy L Thompson @param n Number of rows/columns 143952bfb9bbSJeremy L Thompson 144052bfb9bbSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 144152bfb9bbSJeremy L Thompson 144252bfb9bbSJeremy L Thompson @ref Utility 144352bfb9bbSJeremy L Thompson **/ 144403d18186Sjeremylt CeedPragmaOptimizeOff 1445d1d35e2fSjeremylt int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *mat_A, 1446d3331725Sjeremylt CeedScalar *mat_B, CeedScalar *mat_X, 144752bfb9bbSJeremy L Thompson CeedScalar *lambda, CeedInt n) { 144852bfb9bbSJeremy L Thompson int ierr; 1449d3331725Sjeremylt CeedScalar *mat_C, *mat_G, *vec_D; 145078464608Sjeremylt ierr = CeedCalloc(n*n, &mat_C); CeedChk(ierr); 145178464608Sjeremylt ierr = CeedCalloc(n*n, &mat_G); CeedChk(ierr); 1452d3331725Sjeremylt ierr = CeedCalloc(n, &vec_D); CeedChk(ierr); 145352bfb9bbSJeremy L Thompson 145452bfb9bbSJeremy L Thompson // Compute B = G D G^T 145578464608Sjeremylt memcpy(mat_G, mat_B, n*n*sizeof(mat_B[0])); 1456d3331725Sjeremylt ierr = CeedSymmetricSchurDecomposition(ceed, mat_G, vec_D, n); CeedChk(ierr); 145752bfb9bbSJeremy L Thompson 145885cf89eaSjeremylt // Sort eigenvalues 145985cf89eaSjeremylt for (CeedInt i=n-1; i>=0; i--) 146085cf89eaSjeremylt for (CeedInt j=0; j<i; j++) { 146185cf89eaSjeremylt if (fabs(vec_D[j]) > fabs(vec_D[j+1])) { 146285cf89eaSjeremylt CeedScalar temp; 146385cf89eaSjeremylt temp = vec_D[j]; vec_D[j] = vec_D[j+1]; vec_D[j+1] = temp; 146485cf89eaSjeremylt for (CeedInt k=0; k<n; k++) { 146585cf89eaSjeremylt temp = mat_G[k*n+j]; mat_G[k*n+j] = mat_G[k*n+j+1]; mat_G[k*n+j+1] = temp; 146685cf89eaSjeremylt } 146785cf89eaSjeremylt } 146885cf89eaSjeremylt } 146985cf89eaSjeremylt 1470fb551037Sjeremylt // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T 1471fb551037Sjeremylt // = D^-1/2 G^T A G D^-1/2 1472d3331725Sjeremylt // -- D = D^-1/2 147352bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) 1474d3331725Sjeremylt vec_D[i] = 1./sqrt(vec_D[i]); 1475d3331725Sjeremylt // -- G = G D^-1/2 1476d3331725Sjeremylt // -- C = D^-1/2 G^T 1477d3331725Sjeremylt for (CeedInt i=0; i<n; i++) 1478d3331725Sjeremylt for (CeedInt j=0; j<n; j++) { 1479673160d7Sjeremylt mat_G[i*n+j] *= vec_D[j]; 1480673160d7Sjeremylt mat_C[j*n+i] = mat_G[i*n+j]; 1481d3331725Sjeremylt } 1482673160d7Sjeremylt // -- X = (D^-1/2 G^T) A 1483d1d35e2fSjeremylt ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)mat_C, 1484d3331725Sjeremylt (const CeedScalar *)mat_A, mat_X, n, n, n); 14859289e5bfSjeremylt CeedChk(ierr); 1486673160d7Sjeremylt // -- C = (D^-1/2 G^T A) (G D^-1/2) 1487d3331725Sjeremylt ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)mat_X, 148878464608Sjeremylt (const CeedScalar *)mat_G, mat_C, n, n, n); 14899289e5bfSjeremylt CeedChk(ierr); 149052bfb9bbSJeremy L Thompson 149152bfb9bbSJeremy L Thompson // Compute Q^T C Q = lambda 1492d1d35e2fSjeremylt ierr = CeedSymmetricSchurDecomposition(ceed, mat_C, lambda, n); CeedChk(ierr); 149352bfb9bbSJeremy L Thompson 149485cf89eaSjeremylt // Sort eigenvalues 149585cf89eaSjeremylt for (CeedInt i=n-1; i>=0; i--) 149685cf89eaSjeremylt for (CeedInt j=0; j<i; j++) { 149785cf89eaSjeremylt if (fabs(lambda[j]) > fabs(lambda[j+1])) { 149885cf89eaSjeremylt CeedScalar temp; 149985cf89eaSjeremylt temp = lambda[j]; lambda[j] = lambda[j+1]; lambda[j+1] = temp; 150085cf89eaSjeremylt for (CeedInt k=0; k<n; k++) { 150185cf89eaSjeremylt temp = mat_C[k*n+j]; mat_C[k*n+j] = mat_C[k*n+j+1]; mat_C[k*n+j+1] = temp; 150285cf89eaSjeremylt } 150385cf89eaSjeremylt } 150485cf89eaSjeremylt } 150585cf89eaSjeremylt 1506d3331725Sjeremylt // Set X = (G D^1/2)^-T Q 1507fb551037Sjeremylt // = G D^-1/2 Q 150878464608Sjeremylt ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)mat_G, 1509d3331725Sjeremylt (const CeedScalar *)mat_C, mat_X, n, n, n); 15109289e5bfSjeremylt CeedChk(ierr); 151178464608Sjeremylt 151278464608Sjeremylt // Cleanup 151378464608Sjeremylt ierr = CeedFree(&mat_C); CeedChk(ierr); 151478464608Sjeremylt ierr = CeedFree(&mat_G); CeedChk(ierr); 1515d3331725Sjeremylt ierr = CeedFree(&vec_D); CeedChk(ierr); 1516e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 151752bfb9bbSJeremy L Thompson } 151803d18186Sjeremylt CeedPragmaOptimizeOn 151952bfb9bbSJeremy L Thompson 1520d7b241e6Sjeremylt /// @} 1521