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; 1027a982d89SJeremy L. Thompson CeedScalar v[m]; 1037a982d89SJeremy L. Thompson for (CeedInt ii=0; ii<k; ii++) { 104d1d35e2fSjeremylt CeedInt i = t_mode == CEED_TRANSPOSE ? ii : k-1-ii; 1057a982d89SJeremy L. Thompson for (CeedInt j=i+1; j<m; j++) 1067a982d89SJeremy L. Thompson v[j] = Q[j*k+i]; 107d1d35e2fSjeremylt // Apply Householder reflector (I - tau v v^T) collo_grad_1d^T 108e15f9bd0SJeremy L Thompson ierr = CeedHouseholderReflect(&A[i*row], &v[i], tau[i], m-i, n, row, col); 109e15f9bd0SJeremy L Thompson CeedChk(ierr); 1107a982d89SJeremy L. Thompson } 111e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1127a982d89SJeremy L. Thompson } 1137a982d89SJeremy L. Thompson 1147a982d89SJeremy L. Thompson /** 1157a982d89SJeremy L. Thompson @brief Compute Givens rotation 1167a982d89SJeremy L. Thompson 1177a982d89SJeremy L. Thompson Computes A = G A (or G^T A in transpose mode) 1187a982d89SJeremy L. Thompson where A is an mxn matrix indexed as A[i*n + j*m] 1197a982d89SJeremy L. Thompson 1207a982d89SJeremy L. Thompson @param[in,out] A Row major matrix to apply Givens rotation to, in place 1217a982d89SJeremy L. Thompson @param c Cosine factor 1227a982d89SJeremy L. Thompson @param s Sine factor 123d1d35e2fSjeremylt @param t_mode @ref CEED_NOTRANSPOSE to rotate the basis counter-clockwise, 1244c4400c7SValeria Barra which has the effect of rotating columns of A clockwise; 1254cc79fe7SJed Brown @ref CEED_TRANSPOSE for the opposite rotation 1267a982d89SJeremy L. Thompson @param i First row/column to apply rotation 1277a982d89SJeremy L. Thompson @param k Second row/column to apply rotation 1287a982d89SJeremy L. Thompson @param m Number of rows in A 1297a982d89SJeremy L. Thompson @param n Number of columns in A 1307a982d89SJeremy L. Thompson 1317a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 1327a982d89SJeremy L. Thompson 1337a982d89SJeremy L. Thompson @ref Developer 1347a982d89SJeremy L. Thompson **/ 1357a982d89SJeremy L. Thompson static int CeedGivensRotation(CeedScalar *A, CeedScalar c, CeedScalar s, 136d1d35e2fSjeremylt CeedTransposeMode t_mode, CeedInt i, CeedInt k, 1377a982d89SJeremy L. Thompson CeedInt m, CeedInt n) { 138d1d35e2fSjeremylt CeedInt stride_j = 1, stride_ik = m, num_its = n; 139d1d35e2fSjeremylt if (t_mode == CEED_NOTRANSPOSE) { 140d1d35e2fSjeremylt stride_j = n; stride_ik = 1; num_its = m; 1417a982d89SJeremy L. Thompson } 1427a982d89SJeremy L. Thompson 1437a982d89SJeremy L. Thompson // Apply rotation 144d1d35e2fSjeremylt for (CeedInt j=0; j<num_its; j++) { 145d1d35e2fSjeremylt CeedScalar tau1 = A[i*stride_ik+j*stride_j], tau2 = A[k*stride_ik+j*stride_j]; 146d1d35e2fSjeremylt A[i*stride_ik+j*stride_j] = c*tau1 - s*tau2; 147d1d35e2fSjeremylt A[k*stride_ik+j*stride_j] = s*tau1 + c*tau2; 1487a982d89SJeremy L. Thompson } 149e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1507a982d89SJeremy L. Thompson } 1517a982d89SJeremy L. Thompson 1527a982d89SJeremy L. Thompson /** 1537a982d89SJeremy L. Thompson @brief View an array stored in a CeedBasis 1547a982d89SJeremy L. Thompson 1550a0da059Sjeremylt @param[in] name Name of array 156d1d35e2fSjeremylt @param[in] fp_fmt Printing format 1570a0da059Sjeremylt @param[in] m Number of rows in array 1580a0da059Sjeremylt @param[in] n Number of columns in array 1590a0da059Sjeremylt @param[in] a Array to be viewed 1600a0da059Sjeremylt @param[in] stream Stream to view to, e.g., stdout 1617a982d89SJeremy L. Thompson 1627a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 1637a982d89SJeremy L. Thompson 1647a982d89SJeremy L. Thompson @ref Developer 1657a982d89SJeremy L. Thompson **/ 166d1d35e2fSjeremylt static int CeedScalarView(const char *name, const char *fp_fmt, CeedInt m, 1677a982d89SJeremy L. Thompson CeedInt n, const CeedScalar *a, FILE *stream) { 1687a982d89SJeremy L. Thompson for (int i=0; i<m; i++) { 1697a982d89SJeremy L. Thompson if (m > 1) 1707a982d89SJeremy L. Thompson fprintf(stream, "%12s[%d]:", name, i); 1717a982d89SJeremy L. Thompson else 1727a982d89SJeremy L. Thompson fprintf(stream, "%12s:", name); 1737a982d89SJeremy L. Thompson for (int j=0; j<n; j++) 174d1d35e2fSjeremylt fprintf(stream, fp_fmt, fabs(a[i*n+j]) > 1E-14 ? a[i*n+j] : 0); 1757a982d89SJeremy L. Thompson fputs("\n", stream); 1767a982d89SJeremy L. Thompson } 177e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1787a982d89SJeremy L. Thompson } 1797a982d89SJeremy L. Thompson 1807a982d89SJeremy L. Thompson /// @} 1817a982d89SJeremy L. Thompson 1827a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 1837a982d89SJeremy L. Thompson /// Ceed Backend API 1847a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 1857a982d89SJeremy L. Thompson /// @addtogroup CeedBasisBackend 1867a982d89SJeremy L. Thompson /// @{ 1877a982d89SJeremy L. Thompson 1887a982d89SJeremy L. Thompson /** 1897a982d89SJeremy L. Thompson @brief Return collocated grad matrix 1907a982d89SJeremy L. Thompson 1917a982d89SJeremy L. Thompson @param basis CeedBasis 192d1d35e2fSjeremylt @param[out] collo_grad_1d Row-major (Q_1d * Q_1d) matrix expressing derivatives of 1937a982d89SJeremy L. Thompson basis functions at quadrature points 1947a982d89SJeremy L. Thompson 1957a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 1967a982d89SJeremy L. Thompson 1977a982d89SJeremy L. Thompson @ref Backend 1987a982d89SJeremy L. Thompson **/ 199d1d35e2fSjeremylt int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *collo_grad_1d) { 2007a982d89SJeremy L. Thompson int i, j, k; 2017a982d89SJeremy L. Thompson Ceed ceed; 202d1d35e2fSjeremylt CeedInt ierr, P_1d=(basis)->P_1d, Q_1d=(basis)->Q_1d; 203d1d35e2fSjeremylt CeedScalar *interp_1d, *grad_1d, tau[Q_1d]; 2047a982d89SJeremy L. Thompson 205d1d35e2fSjeremylt ierr = CeedMalloc(Q_1d*P_1d, &interp_1d); CeedChk(ierr); 206d1d35e2fSjeremylt ierr = CeedMalloc(Q_1d*P_1d, &grad_1d); CeedChk(ierr); 207d1d35e2fSjeremylt memcpy(interp_1d, (basis)->interp_1d, Q_1d*P_1d*sizeof(basis)->interp_1d[0]); 208d1d35e2fSjeremylt memcpy(grad_1d, (basis)->grad_1d, Q_1d*P_1d*sizeof(basis)->interp_1d[0]); 2097a982d89SJeremy L. Thompson 210d1d35e2fSjeremylt // QR Factorization, interp_1d = Q R 2117a982d89SJeremy L. Thompson ierr = CeedBasisGetCeed(basis, &ceed); CeedChk(ierr); 212d1d35e2fSjeremylt ierr = CeedQRFactorization(ceed, interp_1d, tau, Q_1d, P_1d); CeedChk(ierr); 213e15f9bd0SJeremy L Thompson // Note: This function is for backend use, so all errors are terminal 214e15f9bd0SJeremy L Thompson // and we do not need to clean up memory on failure. 2157a982d89SJeremy L. Thompson 216d1d35e2fSjeremylt // Apply Rinv, collo_grad_1d = grad_1d Rinv 217d1d35e2fSjeremylt for (i=0; i<Q_1d; i++) { // Row i 218d1d35e2fSjeremylt collo_grad_1d[Q_1d*i] = grad_1d[P_1d*i]/interp_1d[0]; 219d1d35e2fSjeremylt for (j=1; j<P_1d; j++) { // Column j 220d1d35e2fSjeremylt collo_grad_1d[j+Q_1d*i] = grad_1d[j+P_1d*i]; 2217a982d89SJeremy L. Thompson for (k=0; k<j; k++) 222d1d35e2fSjeremylt collo_grad_1d[j+Q_1d*i] -= interp_1d[j+P_1d*k]*collo_grad_1d[k+Q_1d*i]; 223d1d35e2fSjeremylt collo_grad_1d[j+Q_1d*i] /= interp_1d[j+P_1d*j]; 2247a982d89SJeremy L. Thompson } 225d1d35e2fSjeremylt for (j=P_1d; j<Q_1d; j++) 226d1d35e2fSjeremylt collo_grad_1d[j+Q_1d*i] = 0; 2277a982d89SJeremy L. Thompson } 2287a982d89SJeremy L. Thompson 229d1d35e2fSjeremylt // Apply Qtranspose, collograd = collo_grad Q_transpose 230d1d35e2fSjeremylt ierr = CeedHouseholderApplyQ(collo_grad_1d, interp_1d, tau, CEED_NOTRANSPOSE, 231d1d35e2fSjeremylt Q_1d, Q_1d, P_1d, 1, Q_1d); CeedChk(ierr); 2327a982d89SJeremy L. Thompson 233d1d35e2fSjeremylt ierr = CeedFree(&interp_1d); CeedChk(ierr); 234d1d35e2fSjeremylt ierr = CeedFree(&grad_1d); CeedChk(ierr); 235e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 2367a982d89SJeremy L. Thompson } 2377a982d89SJeremy L. Thompson 2387a982d89SJeremy L. Thompson /** 2397a982d89SJeremy L. Thompson @brief Get Ceed associated with a CeedBasis 2407a982d89SJeremy L. Thompson 2417a982d89SJeremy L. Thompson @param basis CeedBasis 2427a982d89SJeremy L. Thompson @param[out] ceed Variable to store Ceed 2437a982d89SJeremy L. Thompson 2447a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 2457a982d89SJeremy L. Thompson 2467a982d89SJeremy L. Thompson @ref Backend 2477a982d89SJeremy L. Thompson **/ 2487a982d89SJeremy L. Thompson int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) { 2497a982d89SJeremy L. Thompson *ceed = basis->ceed; 250e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 2517a982d89SJeremy L. Thompson } 2527a982d89SJeremy L. Thompson 2537a982d89SJeremy L. Thompson /** 2547a982d89SJeremy L. Thompson @brief Get tensor status for given CeedBasis 2557a982d89SJeremy L. Thompson 2567a982d89SJeremy L. Thompson @param basis CeedBasis 257d1d35e2fSjeremylt @param[out] is_tensor Variable to store tensor status 2587a982d89SJeremy L. Thompson 2597a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 2607a982d89SJeremy L. Thompson 2617a982d89SJeremy L. Thompson @ref Backend 2627a982d89SJeremy L. Thompson **/ 263d1d35e2fSjeremylt int CeedBasisIsTensor(CeedBasis basis, bool *is_tensor) { 264d1d35e2fSjeremylt *is_tensor = basis->tensor_basis; 265e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 2667a982d89SJeremy L. Thompson } 2677a982d89SJeremy L. Thompson 2687a982d89SJeremy L. Thompson /** 2697a982d89SJeremy L. Thompson @brief Get backend data of a CeedBasis 2707a982d89SJeremy L. Thompson 2717a982d89SJeremy L. Thompson @param basis CeedBasis 2727a982d89SJeremy L. Thompson @param[out] data Variable to store data 2737a982d89SJeremy L. Thompson 2747a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 2757a982d89SJeremy L. Thompson 2767a982d89SJeremy L. Thompson @ref Backend 2777a982d89SJeremy L. Thompson **/ 278777ff853SJeremy L Thompson int CeedBasisGetData(CeedBasis basis, void *data) { 279777ff853SJeremy L Thompson *(void **)data = basis->data; 280e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 2817a982d89SJeremy L. Thompson } 2827a982d89SJeremy L. Thompson 2837a982d89SJeremy L. Thompson /** 2847a982d89SJeremy L. Thompson @brief Set backend data of a CeedBasis 2857a982d89SJeremy L. Thompson 2867a982d89SJeremy L. Thompson @param[out] basis CeedBasis 2877a982d89SJeremy L. Thompson @param data Data to set 2887a982d89SJeremy L. Thompson 2897a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 2907a982d89SJeremy L. Thompson 2917a982d89SJeremy L. Thompson @ref Backend 2927a982d89SJeremy L. Thompson **/ 293777ff853SJeremy L Thompson int CeedBasisSetData(CeedBasis basis, void *data) { 294777ff853SJeremy L Thompson basis->data = data; 295e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 2967a982d89SJeremy L. Thompson } 2977a982d89SJeremy L. Thompson 2987a982d89SJeremy L. Thompson /** 29934359f16Sjeremylt @brief Increment the reference counter for a CeedBasis 30034359f16Sjeremylt 30134359f16Sjeremylt @param basis Basis to increment the reference counter 30234359f16Sjeremylt 30334359f16Sjeremylt @return An error code: 0 - success, otherwise - failure 30434359f16Sjeremylt 30534359f16Sjeremylt @ref Backend 30634359f16Sjeremylt **/ 307*9560d06aSjeremylt int CeedBasisReference(CeedBasis basis) { 30834359f16Sjeremylt basis->ref_count++; 30934359f16Sjeremylt return CEED_ERROR_SUCCESS; 31034359f16Sjeremylt } 31134359f16Sjeremylt 31234359f16Sjeremylt /** 3137a982d89SJeremy L. Thompson @brief Get dimension for given CeedElemTopology 3147a982d89SJeremy L. Thompson 3157a982d89SJeremy L. Thompson @param topo CeedElemTopology 3167a982d89SJeremy L. Thompson @param[out] dim Variable to store dimension of topology 3177a982d89SJeremy L. Thompson 3187a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3197a982d89SJeremy L. Thompson 3207a982d89SJeremy L. Thompson @ref Backend 3217a982d89SJeremy L. Thompson **/ 3227a982d89SJeremy L. Thompson int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) { 3237a982d89SJeremy L. Thompson *dim = (CeedInt) topo >> 16; 324e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3257a982d89SJeremy L. Thompson } 3267a982d89SJeremy L. Thompson 3277a982d89SJeremy L. Thompson /** 3287a982d89SJeremy L. Thompson @brief Get CeedTensorContract of a CeedBasis 3297a982d89SJeremy L. Thompson 3307a982d89SJeremy L. Thompson @param basis CeedBasis 3317a982d89SJeremy L. Thompson @param[out] contract Variable to store CeedTensorContract 3327a982d89SJeremy L. Thompson 3337a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3347a982d89SJeremy L. Thompson 3357a982d89SJeremy L. Thompson @ref Backend 3367a982d89SJeremy L. Thompson **/ 3377a982d89SJeremy L. Thompson int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) { 3387a982d89SJeremy L. Thompson *contract = basis->contract; 339e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3407a982d89SJeremy L. Thompson } 3417a982d89SJeremy L. Thompson 3427a982d89SJeremy L. Thompson /** 3437a982d89SJeremy L. Thompson @brief Set CeedTensorContract of a CeedBasis 3447a982d89SJeremy L. Thompson 3457a982d89SJeremy L. Thompson @param[out] basis CeedBasis 3467a982d89SJeremy L. Thompson @param contract CeedTensorContract to set 3477a982d89SJeremy L. Thompson 3487a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3497a982d89SJeremy L. Thompson 3507a982d89SJeremy L. Thompson @ref Backend 3517a982d89SJeremy L. Thompson **/ 35234359f16Sjeremylt int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract contract) { 35334359f16Sjeremylt int ierr; 35434359f16Sjeremylt basis->contract = contract; 355*9560d06aSjeremylt ierr = CeedTensorContractReference(contract); CeedChk(ierr); 356e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3577a982d89SJeremy L. Thompson } 3587a982d89SJeremy L. Thompson 3597a982d89SJeremy L. Thompson /** 3607a982d89SJeremy L. Thompson @brief Return a reference implementation of matrix multiplication C = A B. 3617a982d89SJeremy L. Thompson Note, this is a reference implementation for CPU CeedScalar pointers 3627a982d89SJeremy L. Thompson that is not intended for high performance. 3637a982d89SJeremy L. Thompson 3647a982d89SJeremy L. Thompson @param ceed A Ceed context for error handling 365d1d35e2fSjeremylt @param[in] mat_A Row-major matrix A 366d1d35e2fSjeremylt @param[in] mat_B Row-major matrix B 367d1d35e2fSjeremylt @param[out] mat_C Row-major output matrix C 3687a982d89SJeremy L. Thompson @param m Number of rows of C 3697a982d89SJeremy L. Thompson @param n Number of columns of C 3707a982d89SJeremy L. Thompson @param kk Number of columns of A/rows of B 3717a982d89SJeremy L. Thompson 3727a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3737a982d89SJeremy L. Thompson 3747a982d89SJeremy L. Thompson @ref Utility 3757a982d89SJeremy L. Thompson **/ 376d1d35e2fSjeremylt int CeedMatrixMultiply(Ceed ceed, const CeedScalar *mat_A, 377d1d35e2fSjeremylt const CeedScalar *mat_B, CeedScalar *mat_C, CeedInt m, 3787a982d89SJeremy L. Thompson CeedInt n, CeedInt kk) { 3797a982d89SJeremy L. Thompson for (CeedInt i=0; i<m; i++) 3807a982d89SJeremy L. Thompson for (CeedInt j=0; j<n; j++) { 3817a982d89SJeremy L. Thompson CeedScalar sum = 0; 3827a982d89SJeremy L. Thompson for (CeedInt k=0; k<kk; k++) 383d1d35e2fSjeremylt sum += mat_A[k+i*kk]*mat_B[j+k*n]; 384d1d35e2fSjeremylt mat_C[j+i*n] = sum; 3857a982d89SJeremy L. Thompson } 386e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3877a982d89SJeremy L. Thompson } 3887a982d89SJeremy L. Thompson 3897a982d89SJeremy L. Thompson /// @} 3907a982d89SJeremy L. Thompson 3917a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 3927a982d89SJeremy L. Thompson /// CeedBasis Public API 3937a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 3947a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser 395d7b241e6Sjeremylt /// @{ 396d7b241e6Sjeremylt 397b11c1e72Sjeremylt /** 39895bb1877Svaleriabarra @brief Create a tensor-product basis for H^1 discretizations 399b11c1e72Sjeremylt 400b11c1e72Sjeremylt @param ceed A Ceed object where the CeedBasis will be created 401b11c1e72Sjeremylt @param dim Topological dimension 402d1d35e2fSjeremylt @param num_comp Number of field components (1 for scalar fields) 403d1d35e2fSjeremylt @param P_1d Number of nodes in one dimension 404d1d35e2fSjeremylt @param Q_1d Number of quadrature points in one dimension 405d1d35e2fSjeremylt @param interp_1d Row-major (Q_1d * P_1d) matrix expressing the values of nodal 406b11c1e72Sjeremylt basis functions at quadrature points 407d1d35e2fSjeremylt @param grad_1d Row-major (Q_1d * P_1d) matrix expressing derivatives of nodal 408b11c1e72Sjeremylt basis functions at quadrature points 409d1d35e2fSjeremylt @param q_ref_1d Array of length Q_1d holding the locations of quadrature points 410b11c1e72Sjeremylt on the 1D reference element [-1, 1] 411d1d35e2fSjeremylt @param q_weight_1d Array of length Q_1d holding the quadrature weights on the 412b11c1e72Sjeremylt reference element 413b11c1e72Sjeremylt @param[out] basis Address of the variable where the newly created 414b11c1e72Sjeremylt CeedBasis will be stored. 415b11c1e72Sjeremylt 416b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 417dfdf5a53Sjeremylt 4187a982d89SJeremy L. Thompson @ref User 419b11c1e72Sjeremylt **/ 420d1d35e2fSjeremylt int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt num_comp, 421d1d35e2fSjeremylt CeedInt P_1d, CeedInt Q_1d, 422d1d35e2fSjeremylt const CeedScalar *interp_1d, 423d1d35e2fSjeremylt const CeedScalar *grad_1d, const CeedScalar *q_ref_1d, 424d1d35e2fSjeremylt const CeedScalar *q_weight_1d, CeedBasis *basis) { 425d7b241e6Sjeremylt int ierr; 426d7b241e6Sjeremylt 4275fe0d4faSjeremylt if (!ceed->BasisCreateTensorH1) { 4285fe0d4faSjeremylt Ceed delegate; 429aefd8378Sjeremylt ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr); 4305fe0d4faSjeremylt 4315fe0d4faSjeremylt if (!delegate) 432c042f62fSJeremy L Thompson // LCOV_EXCL_START 433e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_UNSUPPORTED, 434e15f9bd0SJeremy L Thompson "Backend does not support BasisCreateTensorH1"); 435c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 4365fe0d4faSjeremylt 437d1d35e2fSjeremylt ierr = CeedBasisCreateTensorH1(delegate, dim, num_comp, P_1d, 438d1d35e2fSjeremylt Q_1d, interp_1d, grad_1d, q_ref_1d, 439d1d35e2fSjeremylt q_weight_1d, basis); CeedChk(ierr); 440e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 4415fe0d4faSjeremylt } 442e15f9bd0SJeremy L Thompson 443e15f9bd0SJeremy L Thompson if (dim<1) 444e15f9bd0SJeremy L Thompson // LCOV_EXCL_START 445e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, 446e15f9bd0SJeremy L Thompson "Basis dimension must be a positive value"); 447e15f9bd0SJeremy L Thompson // LCOV_EXCL_STOP 448d1d35e2fSjeremylt CeedElemTopology topo = dim == 1 ? CEED_LINE 449d1d35e2fSjeremylt : dim == 2 ? CEED_QUAD 450d1d35e2fSjeremylt : CEED_HEX; 451e15f9bd0SJeremy L Thompson 452d7b241e6Sjeremylt ierr = CeedCalloc(1, basis); CeedChk(ierr); 453d7b241e6Sjeremylt (*basis)->ceed = ceed; 454*9560d06aSjeremylt ierr = CeedReference(ceed); CeedChk(ierr); 455d1d35e2fSjeremylt (*basis)->ref_count = 1; 456d1d35e2fSjeremylt (*basis)->tensor_basis = 1; 457d7b241e6Sjeremylt (*basis)->dim = dim; 458d99fa3c5SJeremy L Thompson (*basis)->topo = topo; 459d1d35e2fSjeremylt (*basis)->num_comp = num_comp; 460d1d35e2fSjeremylt (*basis)->P_1d = P_1d; 461d1d35e2fSjeremylt (*basis)->Q_1d = Q_1d; 462d1d35e2fSjeremylt (*basis)->P = CeedIntPow(P_1d, dim); 463d1d35e2fSjeremylt (*basis)->Q = CeedIntPow(Q_1d, dim); 464d1d35e2fSjeremylt ierr = CeedMalloc(Q_1d,&(*basis)->q_ref_1d); CeedChk(ierr); 465d1d35e2fSjeremylt ierr = CeedMalloc(Q_1d,&(*basis)->q_weight_1d); CeedChk(ierr); 466d1d35e2fSjeremylt memcpy((*basis)->q_ref_1d, q_ref_1d, Q_1d*sizeof(q_ref_1d[0])); 467d1d35e2fSjeremylt memcpy((*basis)->q_weight_1d, q_weight_1d, Q_1d*sizeof(q_weight_1d[0])); 468d1d35e2fSjeremylt ierr = CeedMalloc(Q_1d*P_1d,&(*basis)->interp_1d); CeedChk(ierr); 469d1d35e2fSjeremylt ierr = CeedMalloc(Q_1d*P_1d,&(*basis)->grad_1d); CeedChk(ierr); 470d1d35e2fSjeremylt memcpy((*basis)->interp_1d, interp_1d, Q_1d*P_1d*sizeof(interp_1d[0])); 471d1d35e2fSjeremylt memcpy((*basis)->grad_1d, grad_1d, Q_1d*P_1d*sizeof(grad_1d[0])); 472d1d35e2fSjeremylt ierr = ceed->BasisCreateTensorH1(dim, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, 473d1d35e2fSjeremylt q_weight_1d, *basis); CeedChk(ierr); 474e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 475d7b241e6Sjeremylt } 476d7b241e6Sjeremylt 477b11c1e72Sjeremylt /** 47895bb1877Svaleriabarra @brief Create a tensor-product Lagrange basis 479b11c1e72Sjeremylt 480b11c1e72Sjeremylt @param ceed A Ceed object where the CeedBasis will be created 481b11c1e72Sjeremylt @param dim Topological dimension of element 482d1d35e2fSjeremylt @param num_comp Number of field components (1 for scalar fields) 483b11c1e72Sjeremylt @param P Number of Gauss-Lobatto nodes in one dimension. The 484b11c1e72Sjeremylt polynomial degree of the resulting Q_k element is k=P-1. 485b11c1e72Sjeremylt @param Q Number of quadrature points in one dimension. 486d1d35e2fSjeremylt @param quad_mode Distribution of the Q quadrature points (affects order of 487b11c1e72Sjeremylt accuracy for the quadrature) 488b11c1e72Sjeremylt @param[out] basis Address of the variable where the newly created 489b11c1e72Sjeremylt CeedBasis will be stored. 490b11c1e72Sjeremylt 491b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 492dfdf5a53Sjeremylt 4937a982d89SJeremy L. Thompson @ref User 494b11c1e72Sjeremylt **/ 495d1d35e2fSjeremylt int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt num_comp, 496d1d35e2fSjeremylt CeedInt P, CeedInt Q, CeedQuadMode quad_mode, 497692c2638Sjeremylt CeedBasis *basis) { 498d7b241e6Sjeremylt // Allocate 499e15f9bd0SJeremy L Thompson int ierr, ierr2, i, j, k; 500d1d35e2fSjeremylt CeedScalar c1, c2, c3, c4, dx, *nodes, *interp_1d, *grad_1d, *q_ref_1d, 501d1d35e2fSjeremylt *q_weight_1d; 5024d537eeaSYohann 5034d537eeaSYohann if (dim<1) 504c042f62fSJeremy L Thompson // LCOV_EXCL_START 505e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, 506e15f9bd0SJeremy L Thompson "Basis dimension must be a positive value"); 507c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 5084d537eeaSYohann 509e15f9bd0SJeremy L Thompson // Get Nodes and Weights 510d1d35e2fSjeremylt ierr = CeedCalloc(P*Q, &interp_1d); CeedChk(ierr); 511d1d35e2fSjeremylt ierr = CeedCalloc(P*Q, &grad_1d); CeedChk(ierr); 512d7b241e6Sjeremylt ierr = CeedCalloc(P, &nodes); CeedChk(ierr); 513d1d35e2fSjeremylt ierr = CeedCalloc(Q, &q_ref_1d); CeedChk(ierr); 514d1d35e2fSjeremylt ierr = CeedCalloc(Q, &q_weight_1d); CeedChk(ierr); 515e15f9bd0SJeremy L Thompson ierr = CeedLobattoQuadrature(P, nodes, NULL); 516e15f9bd0SJeremy L Thompson if (ierr) { goto cleanup; } CeedChk(ierr); 517d1d35e2fSjeremylt switch (quad_mode) { 518d7b241e6Sjeremylt case CEED_GAUSS: 519d1d35e2fSjeremylt ierr = CeedGaussQuadrature(Q, q_ref_1d, q_weight_1d); 520d7b241e6Sjeremylt break; 521d7b241e6Sjeremylt case CEED_GAUSS_LOBATTO: 522d1d35e2fSjeremylt ierr = CeedLobattoQuadrature(Q, q_ref_1d, q_weight_1d); 523d7b241e6Sjeremylt break; 524d7b241e6Sjeremylt } 525e15f9bd0SJeremy L Thompson if (ierr) { goto cleanup; } CeedChk(ierr); 526e15f9bd0SJeremy L Thompson 527d7b241e6Sjeremylt // Build B, D matrix 528d7b241e6Sjeremylt // Fornberg, 1998 529d7b241e6Sjeremylt for (i = 0; i < Q; i++) { 530d7b241e6Sjeremylt c1 = 1.0; 531d1d35e2fSjeremylt c3 = nodes[0] - q_ref_1d[i]; 532d1d35e2fSjeremylt interp_1d[i*P+0] = 1.0; 533d7b241e6Sjeremylt for (j = 1; j < P; j++) { 534d7b241e6Sjeremylt c2 = 1.0; 535d7b241e6Sjeremylt c4 = c3; 536d1d35e2fSjeremylt c3 = nodes[j] - q_ref_1d[i]; 537d7b241e6Sjeremylt for (k = 0; k < j; k++) { 538d7b241e6Sjeremylt dx = nodes[j] - nodes[k]; 539d7b241e6Sjeremylt c2 *= dx; 540d7b241e6Sjeremylt if (k == j - 1) { 541d1d35e2fSjeremylt grad_1d[i*P + j] = c1*(interp_1d[i*P + k] - c4*grad_1d[i*P + k]) / c2; 542d1d35e2fSjeremylt interp_1d[i*P + j] = - c1*c4*interp_1d[i*P + k] / c2; 543d7b241e6Sjeremylt } 544d1d35e2fSjeremylt grad_1d[i*P + k] = (c3*grad_1d[i*P + k] - interp_1d[i*P + k]) / dx; 545d1d35e2fSjeremylt interp_1d[i*P + k] = c3*interp_1d[i*P + k] / dx; 546d7b241e6Sjeremylt } 547d7b241e6Sjeremylt c1 = c2; 548d7b241e6Sjeremylt } 549d7b241e6Sjeremylt } 550d7b241e6Sjeremylt // // Pass to CeedBasisCreateTensorH1 551d1d35e2fSjeremylt ierr = CeedBasisCreateTensorH1(ceed, dim, num_comp, P, Q, interp_1d, grad_1d, 552d1d35e2fSjeremylt q_ref_1d, 553d1d35e2fSjeremylt q_weight_1d, basis); CeedChk(ierr); 554e15f9bd0SJeremy L Thompson cleanup: 555d1d35e2fSjeremylt ierr2 = CeedFree(&interp_1d); CeedChk(ierr2); 556d1d35e2fSjeremylt ierr2 = CeedFree(&grad_1d); CeedChk(ierr2); 557e15f9bd0SJeremy L Thompson ierr2 = CeedFree(&nodes); CeedChk(ierr2); 558d1d35e2fSjeremylt ierr2 = CeedFree(&q_ref_1d); CeedChk(ierr2); 559d1d35e2fSjeremylt ierr2 = CeedFree(&q_weight_1d); CeedChk(ierr2); 560e15f9bd0SJeremy L Thompson CeedChk(ierr); 561e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 562d7b241e6Sjeremylt } 563d7b241e6Sjeremylt 564b11c1e72Sjeremylt /** 56595bb1877Svaleriabarra @brief Create a non tensor-product basis for H^1 discretizations 566a8de75f0Sjeremylt 567a8de75f0Sjeremylt @param ceed A Ceed object where the CeedBasis will be created 568a8de75f0Sjeremylt @param topo Topology of element, e.g. hypercube, simplex, ect 569d1d35e2fSjeremylt @param num_comp Number of field components (1 for scalar fields) 570d1d35e2fSjeremylt @param num_nodes Total number of nodes 571d1d35e2fSjeremylt @param num_qpts Total number of quadrature points 572d1d35e2fSjeremylt @param interp Row-major (num_qpts * num_nodes) matrix expressing the values of 5738795c945Sjeremylt nodal basis functions at quadrature points 574d1d35e2fSjeremylt @param grad Row-major (num_qpts * dim * num_nodes) matrix expressing 5758795c945Sjeremylt derivatives of nodal basis functions at quadrature points 576d1d35e2fSjeremylt @param q_ref Array of length num_qpts holding the locations of quadrature 5778795c945Sjeremylt points on the reference element [-1, 1] 578d1d35e2fSjeremylt @param q_weight Array of length num_qpts holding the quadrature weights on the 579a8de75f0Sjeremylt reference element 580a8de75f0Sjeremylt @param[out] basis Address of the variable where the newly created 581a8de75f0Sjeremylt CeedBasis will be stored. 582a8de75f0Sjeremylt 583a8de75f0Sjeremylt @return An error code: 0 - success, otherwise - failure 584a8de75f0Sjeremylt 5857a982d89SJeremy L. Thompson @ref User 586a8de75f0Sjeremylt **/ 587d1d35e2fSjeremylt int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, 588d1d35e2fSjeremylt CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 589d1d35e2fSjeremylt const CeedScalar *grad, const CeedScalar *q_ref, 590d1d35e2fSjeremylt const CeedScalar *q_weight, CeedBasis *basis) { 591a8de75f0Sjeremylt int ierr; 592d1d35e2fSjeremylt CeedInt P = num_nodes, Q = num_qpts, dim = 0; 593a8de75f0Sjeremylt 5945fe0d4faSjeremylt if (!ceed->BasisCreateH1) { 5955fe0d4faSjeremylt Ceed delegate; 596aefd8378Sjeremylt ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr); 5975fe0d4faSjeremylt 5985fe0d4faSjeremylt if (!delegate) 599c042f62fSJeremy L Thompson // LCOV_EXCL_START 600e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_UNSUPPORTED, 601e15f9bd0SJeremy L Thompson "Backend does not support BasisCreateH1"); 602c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 6035fe0d4faSjeremylt 604d1d35e2fSjeremylt ierr = CeedBasisCreateH1(delegate, topo, num_comp, num_nodes, 605d1d35e2fSjeremylt num_qpts, interp, grad, q_ref, 606d1d35e2fSjeremylt q_weight, basis); CeedChk(ierr); 607e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 6085fe0d4faSjeremylt } 6095fe0d4faSjeremylt 610a8de75f0Sjeremylt ierr = CeedCalloc(1,basis); CeedChk(ierr); 611a8de75f0Sjeremylt 612a8de75f0Sjeremylt ierr = CeedBasisGetTopologyDimension(topo, &dim); CeedChk(ierr); 613a8de75f0Sjeremylt 614a8de75f0Sjeremylt (*basis)->ceed = ceed; 615*9560d06aSjeremylt ierr = CeedReference(ceed); CeedChk(ierr); 616d1d35e2fSjeremylt (*basis)->ref_count = 1; 617d1d35e2fSjeremylt (*basis)->tensor_basis = 0; 618a8de75f0Sjeremylt (*basis)->dim = dim; 619d99fa3c5SJeremy L Thompson (*basis)->topo = topo; 620d1d35e2fSjeremylt (*basis)->num_comp = num_comp; 621a8de75f0Sjeremylt (*basis)->P = P; 622a8de75f0Sjeremylt (*basis)->Q = Q; 623d1d35e2fSjeremylt ierr = CeedMalloc(Q*dim,&(*basis)->q_ref_1d); CeedChk(ierr); 624d1d35e2fSjeremylt ierr = CeedMalloc(Q,&(*basis)->q_weight_1d); CeedChk(ierr); 625d1d35e2fSjeremylt memcpy((*basis)->q_ref_1d, q_ref, Q*dim*sizeof(q_ref[0])); 626d1d35e2fSjeremylt memcpy((*basis)->q_weight_1d, q_weight, Q*sizeof(q_weight[0])); 62700f91b2bSjeremylt ierr = CeedMalloc(Q*P, &(*basis)->interp); CeedChk(ierr); 62800f91b2bSjeremylt ierr = CeedMalloc(dim*Q*P, &(*basis)->grad); CeedChk(ierr); 62900f91b2bSjeremylt memcpy((*basis)->interp, interp, Q*P*sizeof(interp[0])); 63000f91b2bSjeremylt memcpy((*basis)->grad, grad, dim*Q*P*sizeof(grad[0])); 631d1d35e2fSjeremylt ierr = ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, q_ref, 632d1d35e2fSjeremylt q_weight, *basis); CeedChk(ierr); 633e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 634a8de75f0Sjeremylt } 635a8de75f0Sjeremylt 636a8de75f0Sjeremylt /** 637*9560d06aSjeremylt @brief Copy the pointer to a CeedBasis. Both pointers should 638*9560d06aSjeremylt be destroyed with `CeedBasisDestroy()`; 639*9560d06aSjeremylt Note: If `*basis_copy` is non-NULL, then it is assumed that 640*9560d06aSjeremylt `*basis_copy` is a pointer to a CeedBasis. This CeedBasis 641*9560d06aSjeremylt will be destroyed if `*basis_copy` is the only 642*9560d06aSjeremylt reference to this CeedBasis. 643*9560d06aSjeremylt 644*9560d06aSjeremylt @param basis CeedBasis to copy reference to 645*9560d06aSjeremylt @param[out] basis_copy Variable to store copied reference 646*9560d06aSjeremylt 647*9560d06aSjeremylt @return An error code: 0 - success, otherwise - failure 648*9560d06aSjeremylt 649*9560d06aSjeremylt @ref User 650*9560d06aSjeremylt **/ 651*9560d06aSjeremylt int CeedBasisReferenceCopy(CeedBasis basis, CeedBasis *basis_copy) { 652*9560d06aSjeremylt int ierr; 653*9560d06aSjeremylt 654*9560d06aSjeremylt ierr = CeedBasisReference(basis); CeedChk(ierr); 655*9560d06aSjeremylt ierr = CeedBasisDestroy(basis_copy); CeedChk(ierr); 656*9560d06aSjeremylt *basis_copy = basis; 657*9560d06aSjeremylt return CEED_ERROR_SUCCESS; 658*9560d06aSjeremylt } 659*9560d06aSjeremylt 660*9560d06aSjeremylt /** 6617a982d89SJeremy L. Thompson @brief View a CeedBasis 6627a982d89SJeremy L. Thompson 6637a982d89SJeremy L. Thompson @param basis CeedBasis to view 6647a982d89SJeremy L. Thompson @param stream Stream to view to, e.g., stdout 6657a982d89SJeremy L. Thompson 6667a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 6677a982d89SJeremy L. Thompson 6687a982d89SJeremy L. Thompson @ref User 6697a982d89SJeremy L. Thompson **/ 6707a982d89SJeremy L. Thompson int CeedBasisView(CeedBasis basis, FILE *stream) { 6717a982d89SJeremy L. Thompson int ierr; 6727a982d89SJeremy L. Thompson 673d1d35e2fSjeremylt if (basis->tensor_basis) { 674d1d35e2fSjeremylt fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P_1d, 675d1d35e2fSjeremylt basis->Q_1d); 676d1d35e2fSjeremylt ierr = CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_ref_1d, 6777a982d89SJeremy L. Thompson stream); CeedChk(ierr); 678d1d35e2fSjeremylt ierr = CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q_1d, 679d1d35e2fSjeremylt basis->q_weight_1d, stream); CeedChk(ierr); 680d1d35e2fSjeremylt ierr = CeedScalarView("interp1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, 681d1d35e2fSjeremylt basis->interp_1d, stream); CeedChk(ierr); 682d1d35e2fSjeremylt ierr = CeedScalarView("grad1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, 683d1d35e2fSjeremylt basis->grad_1d, stream); CeedChk(ierr); 6847a982d89SJeremy L. Thompson } else { 6857a982d89SJeremy L. Thompson fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P, 6867a982d89SJeremy L. Thompson basis->Q); 6877a982d89SJeremy L. Thompson ierr = CeedScalarView("qref", "\t% 12.8f", 1, basis->Q*basis->dim, 688d1d35e2fSjeremylt basis->q_ref_1d, 6897a982d89SJeremy L. Thompson stream); CeedChk(ierr); 690d1d35e2fSjeremylt ierr = CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->q_weight_1d, 6917a982d89SJeremy L. Thompson stream); CeedChk(ierr); 6927a982d89SJeremy L. Thompson ierr = CeedScalarView("interp", "\t% 12.8f", basis->Q, basis->P, 6937a982d89SJeremy L. Thompson basis->interp, stream); CeedChk(ierr); 6947a982d89SJeremy L. Thompson ierr = CeedScalarView("grad", "\t% 12.8f", basis->dim*basis->Q, basis->P, 6957a982d89SJeremy L. Thompson basis->grad, stream); CeedChk(ierr); 6967a982d89SJeremy L. Thompson } 697e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 6987a982d89SJeremy L. Thompson } 6997a982d89SJeremy L. Thompson 7007a982d89SJeremy L. Thompson /** 7017a982d89SJeremy L. Thompson @brief Apply basis evaluation from nodes to quadrature points or vice versa 7027a982d89SJeremy L. Thompson 7037a982d89SJeremy L. Thompson @param basis CeedBasis to evaluate 704d1d35e2fSjeremylt @param num_elem The number of elements to apply the basis evaluation to; 7057a982d89SJeremy L. Thompson the backend will specify the ordering in 7064cc79fe7SJed Brown CeedElemRestrictionCreateBlocked() 707d1d35e2fSjeremylt @param t_mode \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature 7087a982d89SJeremy L. Thompson points, \ref CEED_TRANSPOSE to apply the transpose, mapping 7097a982d89SJeremy L. Thompson from quadrature points to nodes 710d1d35e2fSjeremylt @param eval_mode \ref CEED_EVAL_NONE to use values directly, 7117a982d89SJeremy L. Thompson \ref CEED_EVAL_INTERP to use interpolated values, 7127a982d89SJeremy L. Thompson \ref CEED_EVAL_GRAD to use gradients, 7137a982d89SJeremy L. Thompson \ref CEED_EVAL_WEIGHT to use quadrature weights. 7147a982d89SJeremy L. Thompson @param[in] u Input CeedVector 7157a982d89SJeremy L. Thompson @param[out] v Output CeedVector 7167a982d89SJeremy L. Thompson 7177a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 7187a982d89SJeremy L. Thompson 7197a982d89SJeremy L. Thompson @ref User 7207a982d89SJeremy L. Thompson **/ 721d1d35e2fSjeremylt int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode, 722d1d35e2fSjeremylt CeedEvalMode eval_mode, CeedVector u, CeedVector v) { 7237a982d89SJeremy L. Thompson int ierr; 724d1d35e2fSjeremylt CeedInt u_length = 0, v_length, dim, num_comp, num_nodes, num_qpts; 725e15f9bd0SJeremy L Thompson ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr); 726d1d35e2fSjeremylt ierr = CeedBasisGetNumComponents(basis, &num_comp); CeedChk(ierr); 727d1d35e2fSjeremylt ierr = CeedBasisGetNumNodes(basis, &num_nodes); CeedChk(ierr); 728d1d35e2fSjeremylt ierr = CeedBasisGetNumQuadraturePoints(basis, &num_qpts); CeedChk(ierr); 729d1d35e2fSjeremylt ierr = CeedVectorGetLength(v, &v_length); CeedChk(ierr); 7307a982d89SJeremy L. Thompson if (u) { 731d1d35e2fSjeremylt ierr = CeedVectorGetLength(u, &u_length); CeedChk(ierr); 7327a982d89SJeremy L. Thompson } 7337a982d89SJeremy L. Thompson 734e15f9bd0SJeremy L Thompson if (!basis->Apply) 735e15f9bd0SJeremy L Thompson // LCOV_EXCL_START 736e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_UNSUPPORTED, 737e15f9bd0SJeremy L Thompson "Backend does not support BasisApply"); 738e15f9bd0SJeremy L Thompson // LCOV_EXCL_STOP 739e15f9bd0SJeremy L Thompson 740e15f9bd0SJeremy L Thompson // Check compatibility of topological and geometrical dimensions 741d1d35e2fSjeremylt if ((t_mode == CEED_TRANSPOSE && (v_length%num_nodes != 0 || 742d1d35e2fSjeremylt u_length%num_qpts != 0)) || 743d1d35e2fSjeremylt (t_mode == CEED_NOTRANSPOSE && (u_length%num_nodes != 0 || 744d1d35e2fSjeremylt v_length%num_qpts != 0))) 745e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_DIMENSION, 746e15f9bd0SJeremy L Thompson "Length of input/output vectors " 7477a982d89SJeremy L. Thompson "incompatible with basis dimensions"); 7487a982d89SJeremy L. Thompson 749e15f9bd0SJeremy L Thompson // Check vector lengths to prevent out of bounds issues 750d1d35e2fSjeremylt bool bad_dims = false; 751d1d35e2fSjeremylt switch (eval_mode) { 752e15f9bd0SJeremy L Thompson case CEED_EVAL_NONE: 753d1d35e2fSjeremylt case CEED_EVAL_INTERP: bad_dims = 754d1d35e2fSjeremylt ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts || 755d1d35e2fSjeremylt v_length < num_elem*num_comp*num_nodes)) || 756d1d35e2fSjeremylt (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp || 757d1d35e2fSjeremylt u_length < num_elem*num_comp*num_nodes))); 758e15f9bd0SJeremy L Thompson break; 759d1d35e2fSjeremylt case CEED_EVAL_GRAD: bad_dims = 760d1d35e2fSjeremylt ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts*dim || 761d1d35e2fSjeremylt v_length < num_elem*num_comp*num_nodes)) || 762d1d35e2fSjeremylt (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp*dim || 763d1d35e2fSjeremylt u_length < num_elem*num_comp*num_nodes))); 764e15f9bd0SJeremy L Thompson break; 765e15f9bd0SJeremy L Thompson case CEED_EVAL_WEIGHT: 766d1d35e2fSjeremylt bad_dims = v_length < num_elem*num_qpts; 767e15f9bd0SJeremy L Thompson break; 768e15f9bd0SJeremy L Thompson // LCOV_EXCL_START 769d1d35e2fSjeremylt case CEED_EVAL_DIV: bad_dims = 770d1d35e2fSjeremylt ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts || 771d1d35e2fSjeremylt v_length < num_elem*num_comp*num_nodes)) || 772d1d35e2fSjeremylt (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp || 773d1d35e2fSjeremylt u_length < num_elem*num_comp*num_nodes))); 774e15f9bd0SJeremy L Thompson break; 775d1d35e2fSjeremylt case CEED_EVAL_CURL: bad_dims = 776d1d35e2fSjeremylt ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts || 777d1d35e2fSjeremylt v_length < num_elem*num_comp*num_nodes)) || 778d1d35e2fSjeremylt (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp || 779d1d35e2fSjeremylt u_length < num_elem*num_comp*num_nodes))); 780e15f9bd0SJeremy L Thompson break; 781e15f9bd0SJeremy L Thompson // LCOV_EXCL_STOP 782e15f9bd0SJeremy L Thompson } 783d1d35e2fSjeremylt if (bad_dims) 784e15f9bd0SJeremy L Thompson // LCOV_EXCL_START 785e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_DIMENSION, 786d1d35e2fSjeremylt "Input/output vectors too short for basis and evaluation mode"); 787e15f9bd0SJeremy L Thompson // LCOV_EXCL_STOP 788e15f9bd0SJeremy L Thompson 789d1d35e2fSjeremylt ierr = basis->Apply(basis, num_elem, t_mode, eval_mode, u, v); CeedChk(ierr); 790e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 7917a982d89SJeremy L. Thompson } 7927a982d89SJeremy L. Thompson 7937a982d89SJeremy L. Thompson /** 7949d007619Sjeremylt @brief Get dimension for given CeedBasis 7959d007619Sjeremylt 7969d007619Sjeremylt @param basis CeedBasis 7979d007619Sjeremylt @param[out] dim Variable to store dimension of basis 7989d007619Sjeremylt 7999d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 8009d007619Sjeremylt 8019d007619Sjeremylt @ref Backend 8029d007619Sjeremylt **/ 8039d007619Sjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) { 8049d007619Sjeremylt *dim = basis->dim; 805e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 8069d007619Sjeremylt } 8079d007619Sjeremylt 8089d007619Sjeremylt /** 809d99fa3c5SJeremy L Thompson @brief Get topology for given CeedBasis 810d99fa3c5SJeremy L Thompson 811d99fa3c5SJeremy L Thompson @param basis CeedBasis 812d99fa3c5SJeremy L Thompson @param[out] topo Variable to store topology of basis 813d99fa3c5SJeremy L Thompson 814d99fa3c5SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 815d99fa3c5SJeremy L Thompson 816d99fa3c5SJeremy L Thompson @ref Backend 817d99fa3c5SJeremy L Thompson **/ 818d99fa3c5SJeremy L Thompson int CeedBasisGetTopology(CeedBasis basis, CeedElemTopology *topo) { 819d99fa3c5SJeremy L Thompson *topo = basis->topo; 820e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 821d99fa3c5SJeremy L Thompson } 822d99fa3c5SJeremy L Thompson 823d99fa3c5SJeremy L Thompson /** 8249d007619Sjeremylt @brief Get number of components for given CeedBasis 8259d007619Sjeremylt 8269d007619Sjeremylt @param basis CeedBasis 827d1d35e2fSjeremylt @param[out] num_comp Variable to store number of components of basis 8289d007619Sjeremylt 8299d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 8309d007619Sjeremylt 8319d007619Sjeremylt @ref Backend 8329d007619Sjeremylt **/ 833d1d35e2fSjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *num_comp) { 834d1d35e2fSjeremylt *num_comp = basis->num_comp; 835e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 8369d007619Sjeremylt } 8379d007619Sjeremylt 8389d007619Sjeremylt /** 8399d007619Sjeremylt @brief Get total number of nodes (in dim dimensions) of a CeedBasis 8409d007619Sjeremylt 8419d007619Sjeremylt @param basis CeedBasis 8429d007619Sjeremylt @param[out] P Variable to store number of nodes 8439d007619Sjeremylt 8449d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 8459d007619Sjeremylt 8469d007619Sjeremylt @ref Utility 8479d007619Sjeremylt **/ 8489d007619Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) { 8499d007619Sjeremylt *P = basis->P; 850e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 8519d007619Sjeremylt } 8529d007619Sjeremylt 8539d007619Sjeremylt /** 8549d007619Sjeremylt @brief Get total number of nodes (in 1 dimension) of a CeedBasis 8559d007619Sjeremylt 8569d007619Sjeremylt @param basis CeedBasis 857d1d35e2fSjeremylt @param[out] P_1d Variable to store number of nodes 8589d007619Sjeremylt 8599d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 8609d007619Sjeremylt 8619d007619Sjeremylt @ref Backend 8629d007619Sjeremylt **/ 863d1d35e2fSjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P_1d) { 864d1d35e2fSjeremylt if (!basis->tensor_basis) 8659d007619Sjeremylt // LCOV_EXCL_START 866e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_MINOR, 867d1d35e2fSjeremylt "Cannot supply P_1d for non-tensor basis"); 8689d007619Sjeremylt // LCOV_EXCL_STOP 8699d007619Sjeremylt 870d1d35e2fSjeremylt *P_1d = basis->P_1d; 871e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 8729d007619Sjeremylt } 8739d007619Sjeremylt 8749d007619Sjeremylt /** 8759d007619Sjeremylt @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis 8769d007619Sjeremylt 8779d007619Sjeremylt @param basis CeedBasis 8789d007619Sjeremylt @param[out] Q Variable to store number of quadrature points 8799d007619Sjeremylt 8809d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 8819d007619Sjeremylt 8829d007619Sjeremylt @ref Utility 8839d007619Sjeremylt **/ 8849d007619Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) { 8859d007619Sjeremylt *Q = basis->Q; 886e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 8879d007619Sjeremylt } 8889d007619Sjeremylt 8899d007619Sjeremylt /** 8909d007619Sjeremylt @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis 8919d007619Sjeremylt 8929d007619Sjeremylt @param basis CeedBasis 893d1d35e2fSjeremylt @param[out] Q_1d Variable to store number of quadrature points 8949d007619Sjeremylt 8959d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 8969d007619Sjeremylt 8979d007619Sjeremylt @ref Backend 8989d007619Sjeremylt **/ 899d1d35e2fSjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q_1d) { 900d1d35e2fSjeremylt if (!basis->tensor_basis) 9019d007619Sjeremylt // LCOV_EXCL_START 902e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_MINOR, 903d1d35e2fSjeremylt "Cannot supply Q_1d for non-tensor basis"); 9049d007619Sjeremylt // LCOV_EXCL_STOP 9059d007619Sjeremylt 906d1d35e2fSjeremylt *Q_1d = basis->Q_1d; 907e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 9089d007619Sjeremylt } 9099d007619Sjeremylt 9109d007619Sjeremylt /** 9119d007619Sjeremylt @brief Get reference coordinates of quadrature points (in dim dimensions) 9129d007619Sjeremylt of a CeedBasis 9139d007619Sjeremylt 9149d007619Sjeremylt @param basis CeedBasis 915d1d35e2fSjeremylt @param[out] q_ref Variable to store reference coordinates of quadrature points 9169d007619Sjeremylt 9179d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 9189d007619Sjeremylt 9199d007619Sjeremylt @ref Backend 9209d007619Sjeremylt **/ 921d1d35e2fSjeremylt int CeedBasisGetQRef(CeedBasis basis, const CeedScalar **q_ref) { 922d1d35e2fSjeremylt *q_ref = basis->q_ref_1d; 923e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 9249d007619Sjeremylt } 9259d007619Sjeremylt 9269d007619Sjeremylt /** 9279d007619Sjeremylt @brief Get quadrature weights of quadrature points (in dim dimensions) 9289d007619Sjeremylt of a CeedBasis 9299d007619Sjeremylt 9309d007619Sjeremylt @param basis CeedBasis 931d1d35e2fSjeremylt @param[out] q_weight Variable to store quadrature weights 9329d007619Sjeremylt 9339d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 9349d007619Sjeremylt 9359d007619Sjeremylt @ref Backend 9369d007619Sjeremylt **/ 937d1d35e2fSjeremylt int CeedBasisGetQWeights(CeedBasis basis, const CeedScalar **q_weight) { 938d1d35e2fSjeremylt *q_weight = basis->q_weight_1d; 939e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 9409d007619Sjeremylt } 9419d007619Sjeremylt 9429d007619Sjeremylt /** 9439d007619Sjeremylt @brief Get interpolation matrix of a CeedBasis 9449d007619Sjeremylt 9459d007619Sjeremylt @param basis CeedBasis 9469d007619Sjeremylt @param[out] interp Variable to store interpolation matrix 9479d007619Sjeremylt 9489d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 9499d007619Sjeremylt 9509d007619Sjeremylt @ref Backend 9519d007619Sjeremylt **/ 9526c58de82SJeremy L Thompson int CeedBasisGetInterp(CeedBasis basis, const CeedScalar **interp) { 953d1d35e2fSjeremylt if (!basis->interp && basis->tensor_basis) { 9549d007619Sjeremylt // Allocate 9559d007619Sjeremylt int ierr; 9569d007619Sjeremylt ierr = CeedMalloc(basis->Q*basis->P, &basis->interp); CeedChk(ierr); 9579d007619Sjeremylt 9589d007619Sjeremylt // Initialize 9599d007619Sjeremylt for (CeedInt i=0; i<basis->Q*basis->P; i++) 9609d007619Sjeremylt basis->interp[i] = 1.0; 9619d007619Sjeremylt 9629d007619Sjeremylt // Calculate 9639d007619Sjeremylt for (CeedInt d=0; d<basis->dim; d++) 9649d007619Sjeremylt for (CeedInt qpt=0; qpt<basis->Q; qpt++) 9659d007619Sjeremylt for (CeedInt node=0; node<basis->P; node++) { 966d1d35e2fSjeremylt CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d; 967d1d35e2fSjeremylt CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d; 968d1d35e2fSjeremylt basis->interp[qpt*(basis->P)+node] *= basis->interp_1d[q*basis->P_1d+p]; 9699d007619Sjeremylt } 9709d007619Sjeremylt } 9719d007619Sjeremylt *interp = basis->interp; 972e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 9739d007619Sjeremylt } 9749d007619Sjeremylt 9759d007619Sjeremylt /** 9769d007619Sjeremylt @brief Get 1D interpolation matrix of a tensor product CeedBasis 9779d007619Sjeremylt 9789d007619Sjeremylt @param basis CeedBasis 979d1d35e2fSjeremylt @param[out] interp_1d Variable to store interpolation matrix 9809d007619Sjeremylt 9819d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 9829d007619Sjeremylt 9839d007619Sjeremylt @ref Backend 9849d007619Sjeremylt **/ 985d1d35e2fSjeremylt int CeedBasisGetInterp1D(CeedBasis basis, const CeedScalar **interp_1d) { 986d1d35e2fSjeremylt if (!basis->tensor_basis) 9879d007619Sjeremylt // LCOV_EXCL_START 988e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_MINOR, 989e15f9bd0SJeremy L Thompson "CeedBasis is not a tensor product basis."); 9909d007619Sjeremylt // LCOV_EXCL_STOP 9919d007619Sjeremylt 992d1d35e2fSjeremylt *interp_1d = basis->interp_1d; 993e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 9949d007619Sjeremylt } 9959d007619Sjeremylt 9969d007619Sjeremylt /** 9979d007619Sjeremylt @brief Get gradient matrix of a CeedBasis 9989d007619Sjeremylt 9999d007619Sjeremylt @param basis CeedBasis 10009d007619Sjeremylt @param[out] grad Variable to store gradient matrix 10019d007619Sjeremylt 10029d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 10039d007619Sjeremylt 10049d007619Sjeremylt @ref Backend 10059d007619Sjeremylt **/ 10066c58de82SJeremy L Thompson int CeedBasisGetGrad(CeedBasis basis, const CeedScalar **grad) { 1007d1d35e2fSjeremylt if (!basis->grad && basis->tensor_basis) { 10089d007619Sjeremylt // Allocate 10099d007619Sjeremylt int ierr; 10109d007619Sjeremylt ierr = CeedMalloc(basis->dim*basis->Q*basis->P, &basis->grad); 10119d007619Sjeremylt CeedChk(ierr); 10129d007619Sjeremylt 10139d007619Sjeremylt // Initialize 10149d007619Sjeremylt for (CeedInt i=0; i<basis->dim*basis->Q*basis->P; i++) 10159d007619Sjeremylt basis->grad[i] = 1.0; 10169d007619Sjeremylt 10179d007619Sjeremylt // Calculate 10189d007619Sjeremylt for (CeedInt d=0; d<basis->dim; d++) 10199d007619Sjeremylt for (CeedInt i=0; i<basis->dim; i++) 10209d007619Sjeremylt for (CeedInt qpt=0; qpt<basis->Q; qpt++) 10219d007619Sjeremylt for (CeedInt node=0; node<basis->P; node++) { 1022d1d35e2fSjeremylt CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d; 1023d1d35e2fSjeremylt CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d; 10249d007619Sjeremylt if (i == d) 10259d007619Sjeremylt basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *= 1026d1d35e2fSjeremylt basis->grad_1d[q*basis->P_1d+p]; 10279d007619Sjeremylt else 10289d007619Sjeremylt basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *= 1029d1d35e2fSjeremylt basis->interp_1d[q*basis->P_1d+p]; 10309d007619Sjeremylt } 10319d007619Sjeremylt } 10329d007619Sjeremylt *grad = basis->grad; 1033e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 10349d007619Sjeremylt } 10359d007619Sjeremylt 10369d007619Sjeremylt /** 10379d007619Sjeremylt @brief Get 1D gradient matrix of a tensor product CeedBasis 10389d007619Sjeremylt 10399d007619Sjeremylt @param basis CeedBasis 1040d1d35e2fSjeremylt @param[out] grad_1d Variable to store gradient matrix 10419d007619Sjeremylt 10429d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 10439d007619Sjeremylt 10449d007619Sjeremylt @ref Backend 10459d007619Sjeremylt **/ 1046d1d35e2fSjeremylt int CeedBasisGetGrad1D(CeedBasis basis, const CeedScalar **grad_1d) { 1047d1d35e2fSjeremylt if (!basis->tensor_basis) 10489d007619Sjeremylt // LCOV_EXCL_START 1049e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_MINOR, 1050e15f9bd0SJeremy L Thompson "CeedBasis is not a tensor product basis."); 10519d007619Sjeremylt // LCOV_EXCL_STOP 10529d007619Sjeremylt 1053d1d35e2fSjeremylt *grad_1d = basis->grad_1d; 1054e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 10559d007619Sjeremylt } 10569d007619Sjeremylt 10579d007619Sjeremylt /** 10587a982d89SJeremy L. Thompson @brief Destroy a CeedBasis 10597a982d89SJeremy L. Thompson 10607a982d89SJeremy L. Thompson @param basis CeedBasis to destroy 10617a982d89SJeremy L. Thompson 10627a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 10637a982d89SJeremy L. Thompson 10647a982d89SJeremy L. Thompson @ref User 10657a982d89SJeremy L. Thompson **/ 10667a982d89SJeremy L. Thompson int CeedBasisDestroy(CeedBasis *basis) { 10677a982d89SJeremy L. Thompson int ierr; 10687a982d89SJeremy L. Thompson 1069d1d35e2fSjeremylt if (!*basis || --(*basis)->ref_count > 0) return CEED_ERROR_SUCCESS; 10707a982d89SJeremy L. Thompson if ((*basis)->Destroy) { 10717a982d89SJeremy L. Thompson ierr = (*basis)->Destroy(*basis); CeedChk(ierr); 10727a982d89SJeremy L. Thompson } 107334359f16Sjeremylt if ((*basis)->contract) { 107434359f16Sjeremylt ierr = CeedTensorContractDestroy(&(*basis)->contract); CeedChk(ierr); 107534359f16Sjeremylt } 10767a982d89SJeremy L. Thompson ierr = CeedFree(&(*basis)->interp); CeedChk(ierr); 1077d1d35e2fSjeremylt ierr = CeedFree(&(*basis)->interp_1d); CeedChk(ierr); 10787a982d89SJeremy L. Thompson ierr = CeedFree(&(*basis)->grad); CeedChk(ierr); 1079d1d35e2fSjeremylt ierr = CeedFree(&(*basis)->grad_1d); CeedChk(ierr); 1080d1d35e2fSjeremylt ierr = CeedFree(&(*basis)->q_ref_1d); CeedChk(ierr); 1081d1d35e2fSjeremylt ierr = CeedFree(&(*basis)->q_weight_1d); CeedChk(ierr); 10827a982d89SJeremy L. Thompson ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr); 10837a982d89SJeremy L. Thompson ierr = CeedFree(basis); CeedChk(ierr); 1084e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 10857a982d89SJeremy L. Thompson } 10867a982d89SJeremy L. Thompson 10877a982d89SJeremy L. Thompson /** 1088b11c1e72Sjeremylt @brief Construct a Gauss-Legendre quadrature 1089b11c1e72Sjeremylt 1090b11c1e72Sjeremylt @param Q Number of quadrature points (integrates polynomials of 1091b11c1e72Sjeremylt degree 2*Q-1 exactly) 1092d1d35e2fSjeremylt @param[out] q_ref_1d Array of length Q to hold the abscissa on [-1, 1] 1093d1d35e2fSjeremylt @param[out] q_weight_1d Array of length Q to hold the weights 1094b11c1e72Sjeremylt 1095b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 1096dfdf5a53Sjeremylt 1097dfdf5a53Sjeremylt @ref Utility 1098b11c1e72Sjeremylt **/ 1099d1d35e2fSjeremylt int CeedGaussQuadrature(CeedInt Q, CeedScalar *q_ref_1d, 1100d1d35e2fSjeremylt CeedScalar *q_weight_1d) { 1101d7b241e6Sjeremylt // Allocate 1102d7b241e6Sjeremylt CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0*atan(1.0); 1103d1d35e2fSjeremylt // Build q_ref_1d, q_weight_1d 1104d7b241e6Sjeremylt for (int i = 0; i <= Q/2; i++) { 1105d7b241e6Sjeremylt // Guess 1106d7b241e6Sjeremylt xi = cos(PI*(CeedScalar)(2*i+1)/((CeedScalar)(2*Q))); 1107d7b241e6Sjeremylt // Pn(xi) 1108d7b241e6Sjeremylt P0 = 1.0; 1109d7b241e6Sjeremylt P1 = xi; 1110d7b241e6Sjeremylt P2 = 0.0; 1111d7b241e6Sjeremylt for (int j = 2; j <= Q; j++) { 1112d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 1113d7b241e6Sjeremylt P0 = P1; 1114d7b241e6Sjeremylt P1 = P2; 1115d7b241e6Sjeremylt } 1116d7b241e6Sjeremylt // First Newton Step 1117d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 1118d7b241e6Sjeremylt xi = xi-P2/dP2; 1119d7b241e6Sjeremylt // Newton to convergence 11200e4d4210Sjeremylt for (int k=0; k<100 && fabs(P2)>10*CEED_EPSILON; k++) { 1121d7b241e6Sjeremylt P0 = 1.0; 1122d7b241e6Sjeremylt P1 = xi; 1123d7b241e6Sjeremylt for (int j = 2; j <= Q; j++) { 1124d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 1125d7b241e6Sjeremylt P0 = P1; 1126d7b241e6Sjeremylt P1 = P2; 1127d7b241e6Sjeremylt } 1128d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 1129d7b241e6Sjeremylt xi = xi-P2/dP2; 1130d7b241e6Sjeremylt } 1131d7b241e6Sjeremylt // Save xi, wi 1132d7b241e6Sjeremylt wi = 2.0/((1.0-xi*xi)*dP2*dP2); 1133d1d35e2fSjeremylt q_weight_1d[i] = wi; 1134d1d35e2fSjeremylt q_weight_1d[Q-1-i] = wi; 1135d1d35e2fSjeremylt q_ref_1d[i] = -xi; 1136d1d35e2fSjeremylt q_ref_1d[Q-1-i]= xi; 1137d7b241e6Sjeremylt } 1138e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1139d7b241e6Sjeremylt } 1140d7b241e6Sjeremylt 1141b11c1e72Sjeremylt /** 1142b11c1e72Sjeremylt @brief Construct a Gauss-Legendre-Lobatto quadrature 1143b11c1e72Sjeremylt 1144b11c1e72Sjeremylt @param Q Number of quadrature points (integrates polynomials of 1145b11c1e72Sjeremylt degree 2*Q-3 exactly) 1146d1d35e2fSjeremylt @param[out] q_ref_1d Array of length Q to hold the abscissa on [-1, 1] 1147d1d35e2fSjeremylt @param[out] q_weight_1d Array of length Q to hold the weights 1148b11c1e72Sjeremylt 1149b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 1150dfdf5a53Sjeremylt 1151dfdf5a53Sjeremylt @ref Utility 1152b11c1e72Sjeremylt **/ 1153d1d35e2fSjeremylt int CeedLobattoQuadrature(CeedInt Q, CeedScalar *q_ref_1d, 1154d1d35e2fSjeremylt CeedScalar *q_weight_1d) { 1155d7b241e6Sjeremylt // Allocate 1156d7b241e6Sjeremylt CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0*atan(1.0); 1157d1d35e2fSjeremylt // Build q_ref_1d, q_weight_1d 1158d7b241e6Sjeremylt // Set endpoints 115930a100c3SJed Brown if (Q < 2) 1160b0d62198Sjeremylt // LCOV_EXCL_START 1161e15f9bd0SJeremy L Thompson return CeedError(NULL, CEED_ERROR_DIMENSION, 11627ed52d01Sjeremylt "Cannot create Lobatto quadrature with Q=%d < 2 points", Q); 1163b0d62198Sjeremylt // LCOV_EXCL_STOP 1164d7b241e6Sjeremylt wi = 2.0/((CeedScalar)(Q*(Q-1))); 1165d1d35e2fSjeremylt if (q_weight_1d) { 1166d1d35e2fSjeremylt q_weight_1d[0] = wi; 1167d1d35e2fSjeremylt q_weight_1d[Q-1] = wi; 1168d7b241e6Sjeremylt } 1169d1d35e2fSjeremylt q_ref_1d[0] = -1.0; 1170d1d35e2fSjeremylt q_ref_1d[Q-1] = 1.0; 1171d7b241e6Sjeremylt // Interior 1172d7b241e6Sjeremylt for (int i = 1; i <= (Q-1)/2; i++) { 1173d7b241e6Sjeremylt // Guess 1174d7b241e6Sjeremylt xi = cos(PI*(CeedScalar)(i)/(CeedScalar)(Q-1)); 1175d7b241e6Sjeremylt // Pn(xi) 1176d7b241e6Sjeremylt P0 = 1.0; 1177d7b241e6Sjeremylt P1 = xi; 1178d7b241e6Sjeremylt P2 = 0.0; 1179d7b241e6Sjeremylt for (int j = 2; j < Q; j++) { 1180d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 1181d7b241e6Sjeremylt P0 = P1; 1182d7b241e6Sjeremylt P1 = P2; 1183d7b241e6Sjeremylt } 1184d7b241e6Sjeremylt // First Newton step 1185d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 1186d7b241e6Sjeremylt d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi); 1187d7b241e6Sjeremylt xi = xi-dP2/d2P2; 1188d7b241e6Sjeremylt // Newton to convergence 11890e4d4210Sjeremylt for (int k=0; k<100 && fabs(dP2)>10*CEED_EPSILON; k++) { 1190d7b241e6Sjeremylt P0 = 1.0; 1191d7b241e6Sjeremylt P1 = xi; 1192d7b241e6Sjeremylt for (int j = 2; j < Q; j++) { 1193d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 1194d7b241e6Sjeremylt P0 = P1; 1195d7b241e6Sjeremylt P1 = P2; 1196d7b241e6Sjeremylt } 1197d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 1198d7b241e6Sjeremylt d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi); 1199d7b241e6Sjeremylt xi = xi-dP2/d2P2; 1200d7b241e6Sjeremylt } 1201d7b241e6Sjeremylt // Save xi, wi 1202d7b241e6Sjeremylt wi = 2.0/(((CeedScalar)(Q*(Q-1)))*P2*P2); 1203d1d35e2fSjeremylt if (q_weight_1d) { 1204d1d35e2fSjeremylt q_weight_1d[i] = wi; 1205d1d35e2fSjeremylt q_weight_1d[Q-1-i] = wi; 1206d7b241e6Sjeremylt } 1207d1d35e2fSjeremylt q_ref_1d[i] = -xi; 1208d1d35e2fSjeremylt q_ref_1d[Q-1-i]= xi; 1209d7b241e6Sjeremylt } 1210e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1211d7b241e6Sjeremylt } 1212d7b241e6Sjeremylt 1213dfdf5a53Sjeremylt /** 121495bb1877Svaleriabarra @brief Return QR Factorization of a matrix 1215b11c1e72Sjeremylt 121677645d7bSjeremylt @param ceed A Ceed context for error handling 121752bfb9bbSJeremy L Thompson @param[in,out] mat Row-major matrix to be factorized in place 121852bfb9bbSJeremy L Thompson @param[in,out] tau Vector of length m of scaling factors 1219b11c1e72Sjeremylt @param m Number of rows 1220b11c1e72Sjeremylt @param n Number of columns 1221b11c1e72Sjeremylt 1222b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 1223dfdf5a53Sjeremylt 1224dfdf5a53Sjeremylt @ref Utility 1225b11c1e72Sjeremylt **/ 1226a7bd39daSjeremylt int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau, 1227d7b241e6Sjeremylt CeedInt m, CeedInt n) { 1228d7b241e6Sjeremylt CeedScalar v[m]; 1229d7b241e6Sjeremylt 1230a7bd39daSjeremylt // Check m >= n 1231a7bd39daSjeremylt if (n > m) 1232c042f62fSJeremy L Thompson // LCOV_EXCL_START 1233e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_UNSUPPORTED, 1234e15f9bd0SJeremy L Thompson "Cannot compute QR factorization with n > m"); 1235c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 1236a7bd39daSjeremylt 123752bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) { 1238d7b241e6Sjeremylt // Calculate Householder vector, magnitude 1239d7b241e6Sjeremylt CeedScalar sigma = 0.0; 1240d7b241e6Sjeremylt v[i] = mat[i+n*i]; 124152bfb9bbSJeremy L Thompson for (CeedInt j=i+1; j<m; j++) { 1242d7b241e6Sjeremylt v[j] = mat[i+n*j]; 1243d7b241e6Sjeremylt sigma += v[j] * v[j]; 1244d7b241e6Sjeremylt } 1245d7b241e6Sjeremylt CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m] 1246d7b241e6Sjeremylt CeedScalar Rii = -copysign(norm, v[i]); 1247d7b241e6Sjeremylt v[i] -= Rii; 1248d7b241e6Sjeremylt // norm of v[i:m] after modification above and scaling below 1249d7b241e6Sjeremylt // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 1250d7b241e6Sjeremylt // tau = 2 / (norm*norm) 1251d7b241e6Sjeremylt tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma); 1252fb551037Sjeremylt 12531d102b48SJeremy L Thompson for (CeedInt j=i+1; j<m; j++) 12541d102b48SJeremy L Thompson v[j] /= v[i]; 1255d7b241e6Sjeremylt 1256d7b241e6Sjeremylt // Apply Householder reflector to lower right panel 1257d7b241e6Sjeremylt CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1); 1258d7b241e6Sjeremylt // Save v 1259d7b241e6Sjeremylt mat[i+n*i] = Rii; 12601d102b48SJeremy L Thompson for (CeedInt j=i+1; j<m; j++) 1261d7b241e6Sjeremylt mat[i+n*j] = v[j]; 1262d7b241e6Sjeremylt } 1263e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1264d7b241e6Sjeremylt } 1265d7b241e6Sjeremylt 1266b11c1e72Sjeremylt /** 126752bfb9bbSJeremy L Thompson @brief Return symmetric Schur decomposition of the symmetric matrix mat via 126852bfb9bbSJeremy L Thompson symmetric QR factorization 126952bfb9bbSJeremy L Thompson 127077645d7bSjeremylt @param ceed A Ceed context for error handling 127152bfb9bbSJeremy L Thompson @param[in,out] mat Row-major matrix to be factorized in place 1272460bf743SValeria Barra @param[out] lambda Vector of length n of eigenvalues 127352bfb9bbSJeremy L Thompson @param n Number of rows/columns 127452bfb9bbSJeremy L Thompson 127552bfb9bbSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 127652bfb9bbSJeremy L Thompson 127752bfb9bbSJeremy L Thompson @ref Utility 127852bfb9bbSJeremy L Thompson **/ 127952bfb9bbSJeremy L Thompson int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat, 128052bfb9bbSJeremy L Thompson CeedScalar *lambda, CeedInt n) { 128152bfb9bbSJeremy L Thompson // Check bounds for clang-tidy 128252bfb9bbSJeremy L Thompson if (n<2) 1283c042f62fSJeremy L Thompson // LCOV_EXCL_START 1284e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_UNSUPPORTED, 1285c042f62fSJeremy L Thompson "Cannot compute symmetric Schur decomposition of scalars"); 1286c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 128752bfb9bbSJeremy L Thompson 128852bfb9bbSJeremy L Thompson CeedScalar v[n-1], tau[n-1], matT[n*n]; 128952bfb9bbSJeremy L Thompson 129052bfb9bbSJeremy L Thompson // Copy mat to matT and set mat to I 129152bfb9bbSJeremy L Thompson memcpy(matT, mat, n*n*sizeof(mat[0])); 129252bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) 129352bfb9bbSJeremy L Thompson for (CeedInt j=0; j<n; j++) 129452bfb9bbSJeremy L Thompson mat[j+n*i] = (i==j) ? 1 : 0; 129552bfb9bbSJeremy L Thompson 129652bfb9bbSJeremy L Thompson // Reduce to tridiagonal 129752bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n-1; i++) { 129852bfb9bbSJeremy L Thompson // Calculate Householder vector, magnitude 129952bfb9bbSJeremy L Thompson CeedScalar sigma = 0.0; 130052bfb9bbSJeremy L Thompson v[i] = matT[i+n*(i+1)]; 130152bfb9bbSJeremy L Thompson for (CeedInt j=i+1; j<n-1; j++) { 130252bfb9bbSJeremy L Thompson v[j] = matT[i+n*(j+1)]; 130352bfb9bbSJeremy L Thompson sigma += v[j] * v[j]; 130452bfb9bbSJeremy L Thompson } 130552bfb9bbSJeremy L Thompson CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:n-1] 130652bfb9bbSJeremy L Thompson CeedScalar Rii = -copysign(norm, v[i]); 130752bfb9bbSJeremy L Thompson v[i] -= Rii; 130852bfb9bbSJeremy L Thompson // norm of v[i:m] after modification above and scaling below 130952bfb9bbSJeremy L Thompson // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 131052bfb9bbSJeremy L Thompson // tau = 2 / (norm*norm) 13110e4d4210Sjeremylt if (sigma > 10*CEED_EPSILON) 131252bfb9bbSJeremy L Thompson tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma); 1313fb551037Sjeremylt else 1314fb551037Sjeremylt tau[i] = 0; 1315fb551037Sjeremylt 1316fb551037Sjeremylt for (CeedInt j=i+1; j<n-1; j++) 1317fb551037Sjeremylt v[j] /= v[i]; 131852bfb9bbSJeremy L Thompson 131952bfb9bbSJeremy L Thompson // Update sub and super diagonal 132052bfb9bbSJeremy L Thompson matT[i+n*(i+1)] = Rii; 132152bfb9bbSJeremy L Thompson matT[(i+1)+n*i] = Rii; 132252bfb9bbSJeremy L Thompson for (CeedInt j=i+2; j<n; j++) { 132352bfb9bbSJeremy L Thompson matT[i+n*j] = 0; matT[j+n*i] = 0; 132452bfb9bbSJeremy L Thompson } 132552bfb9bbSJeremy L Thompson // Apply symmetric Householder reflector to lower right panel 132652bfb9bbSJeremy L Thompson CeedHouseholderReflect(&matT[(i+1)+n*(i+1)], &v[i], tau[i], 132752bfb9bbSJeremy L Thompson n-(i+1), n-(i+1), n, 1); 132852bfb9bbSJeremy L Thompson CeedHouseholderReflect(&matT[(i+1)+n*(i+1)], &v[i], tau[i], 132952bfb9bbSJeremy L Thompson n-(i+1), n-(i+1), 1, n); 133052bfb9bbSJeremy L Thompson // Save v 133152bfb9bbSJeremy L Thompson for (CeedInt j=i+1; j<n-1; j++) { 133252bfb9bbSJeremy L Thompson matT[i+n*(j+1)] = v[j]; 133352bfb9bbSJeremy L Thompson } 133452bfb9bbSJeremy L Thompson } 133552bfb9bbSJeremy L Thompson // Backwards accumulation of Q 133652bfb9bbSJeremy L Thompson for (CeedInt i=n-2; i>=0; i--) { 133752bfb9bbSJeremy L Thompson v[i] = 1; 133852bfb9bbSJeremy L Thompson for (CeedInt j=i+1; j<n-1; j++) { 133952bfb9bbSJeremy L Thompson v[j] = matT[i+n*(j+1)]; 134052bfb9bbSJeremy L Thompson matT[i+n*(j+1)] = 0; 134152bfb9bbSJeremy L Thompson } 134252bfb9bbSJeremy L Thompson CeedHouseholderReflect(&mat[(i+1)+n*(i+1)], &v[i], tau[i], 134352bfb9bbSJeremy L Thompson n-(i+1), n-(i+1), n, 1); 134452bfb9bbSJeremy L Thompson } 134552bfb9bbSJeremy L Thompson 134652bfb9bbSJeremy L Thompson // Reduce sub and super diagonal 134752bfb9bbSJeremy L Thompson CeedInt p = 0, q = 0, itr = 0, maxitr = n*n*n; 13480e4d4210Sjeremylt CeedScalar tol = 10*CEED_EPSILON; 134952bfb9bbSJeremy L Thompson 135052bfb9bbSJeremy L Thompson while (q < n && itr < maxitr) { 135152bfb9bbSJeremy L Thompson // Update p, q, size of reduced portions of diagonal 135252bfb9bbSJeremy L Thompson p = 0; q = 0; 135352bfb9bbSJeremy L Thompson for (CeedInt i=n-2; i>=0; i--) { 135452bfb9bbSJeremy L Thompson if (fabs(matT[i+n*(i+1)]) < tol) 135552bfb9bbSJeremy L Thompson q += 1; 135652bfb9bbSJeremy L Thompson else 135752bfb9bbSJeremy L Thompson break; 135852bfb9bbSJeremy L Thompson } 135952bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n-1-q; i++) { 136052bfb9bbSJeremy L Thompson if (fabs(matT[i+n*(i+1)]) < tol) 136152bfb9bbSJeremy L Thompson p += 1; 136252bfb9bbSJeremy L Thompson else 136352bfb9bbSJeremy L Thompson break; 136452bfb9bbSJeremy L Thompson } 136552bfb9bbSJeremy L Thompson if (q == n-1) break; // Finished reducing 136652bfb9bbSJeremy L Thompson 136752bfb9bbSJeremy L Thompson // Reduce tridiagonal portion 136852bfb9bbSJeremy L Thompson CeedScalar tnn = matT[(n-1-q)+n*(n-1-q)], 136952bfb9bbSJeremy L Thompson tnnm1 = matT[(n-2-q)+n*(n-1-q)]; 137052bfb9bbSJeremy L Thompson CeedScalar d = (matT[(n-2-q)+n*(n-2-q)] - tnn)/2; 137152bfb9bbSJeremy L Thompson CeedScalar mu = tnn - tnnm1*tnnm1 / 137252bfb9bbSJeremy L Thompson (d + copysign(sqrt(d*d + tnnm1*tnnm1), d)); 137352bfb9bbSJeremy L Thompson CeedScalar x = matT[p+n*p] - mu; 137452bfb9bbSJeremy L Thompson CeedScalar z = matT[p+n*(p+1)]; 137552bfb9bbSJeremy L Thompson for (CeedInt k=p; k<n-1-q; k++) { 137652bfb9bbSJeremy L Thompson // Compute Givens rotation 137752bfb9bbSJeremy L Thompson CeedScalar c = 1, s = 0; 137852bfb9bbSJeremy L Thompson if (fabs(z) > tol) { 137952bfb9bbSJeremy L Thompson if (fabs(z) > fabs(x)) { 138052bfb9bbSJeremy L Thompson CeedScalar tau = -x/z; 138152bfb9bbSJeremy L Thompson s = 1/sqrt(1+tau*tau), c = s*tau; 138252bfb9bbSJeremy L Thompson } else { 138352bfb9bbSJeremy L Thompson CeedScalar tau = -z/x; 138452bfb9bbSJeremy L Thompson c = 1/sqrt(1+tau*tau), s = c*tau; 138552bfb9bbSJeremy L Thompson } 138652bfb9bbSJeremy L Thompson } 138752bfb9bbSJeremy L Thompson 138852bfb9bbSJeremy L Thompson // Apply Givens rotation to T 138952bfb9bbSJeremy L Thompson CeedGivensRotation(matT, c, s, CEED_NOTRANSPOSE, k, k+1, n, n); 139052bfb9bbSJeremy L Thompson CeedGivensRotation(matT, c, s, CEED_TRANSPOSE, k, k+1, n, n); 139152bfb9bbSJeremy L Thompson 139252bfb9bbSJeremy L Thompson // Apply Givens rotation to Q 139352bfb9bbSJeremy L Thompson CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k+1, n, n); 139452bfb9bbSJeremy L Thompson 139552bfb9bbSJeremy L Thompson // Update x, z 139652bfb9bbSJeremy L Thompson if (k < n-q-2) { 139752bfb9bbSJeremy L Thompson x = matT[k+n*(k+1)]; 139852bfb9bbSJeremy L Thompson z = matT[k+n*(k+2)]; 139952bfb9bbSJeremy L Thompson } 140052bfb9bbSJeremy L Thompson } 140152bfb9bbSJeremy L Thompson itr++; 140252bfb9bbSJeremy L Thompson } 140352bfb9bbSJeremy L Thompson // Save eigenvalues 140452bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) 140552bfb9bbSJeremy L Thompson lambda[i] = matT[i+n*i]; 140652bfb9bbSJeremy L Thompson 140752bfb9bbSJeremy L Thompson // Check convergence 140852bfb9bbSJeremy L Thompson if (itr == maxitr && q < n-1) 1409c042f62fSJeremy L Thompson // LCOV_EXCL_START 1410e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_MINOR, 1411e15f9bd0SJeremy L Thompson "Symmetric QR failed to converge"); 1412c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 1413e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 141452bfb9bbSJeremy L Thompson } 141552bfb9bbSJeremy L Thompson 141652bfb9bbSJeremy L Thompson /** 141752bfb9bbSJeremy L Thompson @brief Return Simultaneous Diagonalization of two matrices. This solves the 141852bfb9bbSJeremy L Thompson generalized eigenvalue problem A x = lambda B x, where A and B 141952bfb9bbSJeremy L Thompson are symmetric and B is positive definite. We generate the matrix X 142052bfb9bbSJeremy L Thompson and vector Lambda such that X^T A X = Lambda and X^T B X = I. This 142152bfb9bbSJeremy L Thompson is equivalent to the LAPACK routine 'sygv' with TYPE = 1. 142252bfb9bbSJeremy L Thompson 142377645d7bSjeremylt @param ceed A Ceed context for error handling 1424d1d35e2fSjeremylt @param[in] mat_A Row-major matrix to be factorized with eigenvalues 1425d1d35e2fSjeremylt @param[in] mat_B Row-major matrix to be factorized to identity 142652bfb9bbSJeremy L Thompson @param[out] x Row-major orthogonal matrix 1427460bf743SValeria Barra @param[out] lambda Vector of length n of generalized eigenvalues 142852bfb9bbSJeremy L Thompson @param n Number of rows/columns 142952bfb9bbSJeremy L Thompson 143052bfb9bbSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 143152bfb9bbSJeremy L Thompson 143252bfb9bbSJeremy L Thompson @ref Utility 143352bfb9bbSJeremy L Thompson **/ 1434d1d35e2fSjeremylt int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *mat_A, 1435d1d35e2fSjeremylt CeedScalar *mat_B, CeedScalar *x, 143652bfb9bbSJeremy L Thompson CeedScalar *lambda, CeedInt n) { 143752bfb9bbSJeremy L Thompson int ierr; 1438d1d35e2fSjeremylt CeedScalar mat_C[n*n], matG[n*n], vecD[n]; 143952bfb9bbSJeremy L Thompson 144052bfb9bbSJeremy L Thompson // Compute B = G D G^T 1441d1d35e2fSjeremylt memcpy(matG, mat_B, n*n*sizeof(mat_B[0])); 144252bfb9bbSJeremy L Thompson ierr = CeedSymmetricSchurDecomposition(ceed, matG, vecD, n); CeedChk(ierr); 1443fb551037Sjeremylt for (CeedInt i=0; i<n; i++) 1444fb551037Sjeremylt vecD[i] = sqrt(vecD[i]); 144552bfb9bbSJeremy L Thompson 1446fb551037Sjeremylt // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T 1447fb551037Sjeremylt // = D^-1/2 G^T A G D^-1/2 144852bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) 144952bfb9bbSJeremy L Thompson for (CeedInt j=0; j<n; j++) 1450d1d35e2fSjeremylt mat_C[j+i*n] = matG[i+j*n] / vecD[i]; 1451d1d35e2fSjeremylt ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)mat_C, 1452d1d35e2fSjeremylt (const CeedScalar *)mat_A, x, n, n, n); 14539289e5bfSjeremylt CeedChk(ierr); 145452bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) 145552bfb9bbSJeremy L Thompson for (CeedInt j=0; j<n; j++) 1456fb551037Sjeremylt matG[j+i*n] = matG[j+i*n] / vecD[j]; 14579289e5bfSjeremylt ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)x, 1458d1d35e2fSjeremylt (const CeedScalar *)matG, mat_C, n, n, n); 14599289e5bfSjeremylt CeedChk(ierr); 146052bfb9bbSJeremy L Thompson 146152bfb9bbSJeremy L Thompson // Compute Q^T C Q = lambda 1462d1d35e2fSjeremylt ierr = CeedSymmetricSchurDecomposition(ceed, mat_C, lambda, n); CeedChk(ierr); 146352bfb9bbSJeremy L Thompson 1464fb551037Sjeremylt // Set x = (G D^1/2)^-T Q 1465fb551037Sjeremylt // = G D^-1/2 Q 14669289e5bfSjeremylt ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)matG, 1467d1d35e2fSjeremylt (const CeedScalar *)mat_C, x, n, n, n); 14689289e5bfSjeremylt CeedChk(ierr); 1469e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 147052bfb9bbSJeremy L Thompson } 147152bfb9bbSJeremy L Thompson 1472d7b241e6Sjeremylt /// @} 1473