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 86*d1d35e2fSjeremylt @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, 98*d1d35e2fSjeremylt 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++) { 104*d1d35e2fSjeremylt 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]; 107*d1d35e2fSjeremylt // 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 123*d1d35e2fSjeremylt @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, 136*d1d35e2fSjeremylt CeedTransposeMode t_mode, CeedInt i, CeedInt k, 1377a982d89SJeremy L. Thompson CeedInt m, CeedInt n) { 138*d1d35e2fSjeremylt CeedInt stride_j = 1, stride_ik = m, num_its = n; 139*d1d35e2fSjeremylt if (t_mode == CEED_NOTRANSPOSE) { 140*d1d35e2fSjeremylt stride_j = n; stride_ik = 1; num_its = m; 1417a982d89SJeremy L. Thompson } 1427a982d89SJeremy L. Thompson 1437a982d89SJeremy L. Thompson // Apply rotation 144*d1d35e2fSjeremylt for (CeedInt j=0; j<num_its; j++) { 145*d1d35e2fSjeremylt CeedScalar tau1 = A[i*stride_ik+j*stride_j], tau2 = A[k*stride_ik+j*stride_j]; 146*d1d35e2fSjeremylt A[i*stride_ik+j*stride_j] = c*tau1 - s*tau2; 147*d1d35e2fSjeremylt 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 156*d1d35e2fSjeremylt @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 **/ 166*d1d35e2fSjeremylt 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++) 174*d1d35e2fSjeremylt 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 192*d1d35e2fSjeremylt @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 **/ 199*d1d35e2fSjeremylt int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *collo_grad_1d) { 2007a982d89SJeremy L. Thompson int i, j, k; 2017a982d89SJeremy L. Thompson Ceed ceed; 202*d1d35e2fSjeremylt CeedInt ierr, P_1d=(basis)->P_1d, Q_1d=(basis)->Q_1d; 203*d1d35e2fSjeremylt CeedScalar *interp_1d, *grad_1d, tau[Q_1d]; 2047a982d89SJeremy L. Thompson 205*d1d35e2fSjeremylt ierr = CeedMalloc(Q_1d*P_1d, &interp_1d); CeedChk(ierr); 206*d1d35e2fSjeremylt ierr = CeedMalloc(Q_1d*P_1d, &grad_1d); CeedChk(ierr); 207*d1d35e2fSjeremylt memcpy(interp_1d, (basis)->interp_1d, Q_1d*P_1d*sizeof(basis)->interp_1d[0]); 208*d1d35e2fSjeremylt memcpy(grad_1d, (basis)->grad_1d, Q_1d*P_1d*sizeof(basis)->interp_1d[0]); 2097a982d89SJeremy L. Thompson 210*d1d35e2fSjeremylt // QR Factorization, interp_1d = Q R 2117a982d89SJeremy L. Thompson ierr = CeedBasisGetCeed(basis, &ceed); CeedChk(ierr); 212*d1d35e2fSjeremylt 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 216*d1d35e2fSjeremylt // Apply Rinv, collo_grad_1d = grad_1d Rinv 217*d1d35e2fSjeremylt for (i=0; i<Q_1d; i++) { // Row i 218*d1d35e2fSjeremylt collo_grad_1d[Q_1d*i] = grad_1d[P_1d*i]/interp_1d[0]; 219*d1d35e2fSjeremylt for (j=1; j<P_1d; j++) { // Column j 220*d1d35e2fSjeremylt collo_grad_1d[j+Q_1d*i] = grad_1d[j+P_1d*i]; 2217a982d89SJeremy L. Thompson for (k=0; k<j; k++) 222*d1d35e2fSjeremylt collo_grad_1d[j+Q_1d*i] -= interp_1d[j+P_1d*k]*collo_grad_1d[k+Q_1d*i]; 223*d1d35e2fSjeremylt collo_grad_1d[j+Q_1d*i] /= interp_1d[j+P_1d*j]; 2247a982d89SJeremy L. Thompson } 225*d1d35e2fSjeremylt for (j=P_1d; j<Q_1d; j++) 226*d1d35e2fSjeremylt collo_grad_1d[j+Q_1d*i] = 0; 2277a982d89SJeremy L. Thompson } 2287a982d89SJeremy L. Thompson 229*d1d35e2fSjeremylt // Apply Qtranspose, collograd = collo_grad Q_transpose 230*d1d35e2fSjeremylt ierr = CeedHouseholderApplyQ(collo_grad_1d, interp_1d, tau, CEED_NOTRANSPOSE, 231*d1d35e2fSjeremylt Q_1d, Q_1d, P_1d, 1, Q_1d); CeedChk(ierr); 2327a982d89SJeremy L. Thompson 233*d1d35e2fSjeremylt ierr = CeedFree(&interp_1d); CeedChk(ierr); 234*d1d35e2fSjeremylt 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 257*d1d35e2fSjeremylt @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 **/ 263*d1d35e2fSjeremylt int CeedBasisIsTensor(CeedBasis basis, bool *is_tensor) { 264*d1d35e2fSjeremylt *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 /** 2997a982d89SJeremy L. Thompson @brief Get dimension for given CeedElemTopology 3007a982d89SJeremy L. Thompson 3017a982d89SJeremy L. Thompson @param topo CeedElemTopology 3027a982d89SJeremy L. Thompson @param[out] dim Variable to store dimension of topology 3037a982d89SJeremy L. Thompson 3047a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3057a982d89SJeremy L. Thompson 3067a982d89SJeremy L. Thompson @ref Backend 3077a982d89SJeremy L. Thompson **/ 3087a982d89SJeremy L. Thompson int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) { 3097a982d89SJeremy L. Thompson *dim = (CeedInt) topo >> 16; 310e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3117a982d89SJeremy L. Thompson } 3127a982d89SJeremy L. Thompson 3137a982d89SJeremy L. Thompson /** 3147a982d89SJeremy L. Thompson @brief Get CeedTensorContract of a CeedBasis 3157a982d89SJeremy L. Thompson 3167a982d89SJeremy L. Thompson @param basis CeedBasis 3177a982d89SJeremy L. Thompson @param[out] contract Variable to store CeedTensorContract 3187a982d89SJeremy L. Thompson 3197a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3207a982d89SJeremy L. Thompson 3217a982d89SJeremy L. Thompson @ref Backend 3227a982d89SJeremy L. Thompson **/ 3237a982d89SJeremy L. Thompson int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) { 3247a982d89SJeremy L. Thompson *contract = basis->contract; 325e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3267a982d89SJeremy L. Thompson } 3277a982d89SJeremy L. Thompson 3287a982d89SJeremy L. Thompson /** 3297a982d89SJeremy L. Thompson @brief Set CeedTensorContract of a CeedBasis 3307a982d89SJeremy L. Thompson 3317a982d89SJeremy L. Thompson @param[out] basis CeedBasis 3327a982d89SJeremy L. Thompson @param contract CeedTensorContract to set 3337a982d89SJeremy L. Thompson 3347a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3357a982d89SJeremy L. Thompson 3367a982d89SJeremy L. Thompson @ref Backend 3377a982d89SJeremy L. Thompson **/ 3387a982d89SJeremy L. Thompson int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract *contract) { 3397a982d89SJeremy L. Thompson basis->contract = *contract; 340e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3417a982d89SJeremy L. Thompson } 3427a982d89SJeremy L. Thompson 3437a982d89SJeremy L. Thompson /** 3447a982d89SJeremy L. Thompson @brief Return a reference implementation of matrix multiplication C = A B. 3457a982d89SJeremy L. Thompson Note, this is a reference implementation for CPU CeedScalar pointers 3467a982d89SJeremy L. Thompson that is not intended for high performance. 3477a982d89SJeremy L. Thompson 3487a982d89SJeremy L. Thompson @param ceed A Ceed context for error handling 349*d1d35e2fSjeremylt @param[in] mat_A Row-major matrix A 350*d1d35e2fSjeremylt @param[in] mat_B Row-major matrix B 351*d1d35e2fSjeremylt @param[out] mat_C Row-major output matrix C 3527a982d89SJeremy L. Thompson @param m Number of rows of C 3537a982d89SJeremy L. Thompson @param n Number of columns of C 3547a982d89SJeremy L. Thompson @param kk Number of columns of A/rows of B 3557a982d89SJeremy L. Thompson 3567a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3577a982d89SJeremy L. Thompson 3587a982d89SJeremy L. Thompson @ref Utility 3597a982d89SJeremy L. Thompson **/ 360*d1d35e2fSjeremylt int CeedMatrixMultiply(Ceed ceed, const CeedScalar *mat_A, 361*d1d35e2fSjeremylt const CeedScalar *mat_B, CeedScalar *mat_C, CeedInt m, 3627a982d89SJeremy L. Thompson CeedInt n, CeedInt kk) { 3637a982d89SJeremy L. Thompson for (CeedInt i=0; i<m; i++) 3647a982d89SJeremy L. Thompson for (CeedInt j=0; j<n; j++) { 3657a982d89SJeremy L. Thompson CeedScalar sum = 0; 3667a982d89SJeremy L. Thompson for (CeedInt k=0; k<kk; k++) 367*d1d35e2fSjeremylt sum += mat_A[k+i*kk]*mat_B[j+k*n]; 368*d1d35e2fSjeremylt mat_C[j+i*n] = sum; 3697a982d89SJeremy L. Thompson } 370e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 3717a982d89SJeremy L. Thompson } 3727a982d89SJeremy L. Thompson 3737a982d89SJeremy L. Thompson /// @} 3747a982d89SJeremy L. Thompson 3757a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 3767a982d89SJeremy L. Thompson /// CeedBasis Public API 3777a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 3787a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser 379d7b241e6Sjeremylt /// @{ 380d7b241e6Sjeremylt 381b11c1e72Sjeremylt /** 38295bb1877Svaleriabarra @brief Create a tensor-product basis for H^1 discretizations 383b11c1e72Sjeremylt 384b11c1e72Sjeremylt @param ceed A Ceed object where the CeedBasis will be created 385b11c1e72Sjeremylt @param dim Topological dimension 386*d1d35e2fSjeremylt @param num_comp Number of field components (1 for scalar fields) 387*d1d35e2fSjeremylt @param P_1d Number of nodes in one dimension 388*d1d35e2fSjeremylt @param Q_1d Number of quadrature points in one dimension 389*d1d35e2fSjeremylt @param interp_1d Row-major (Q_1d * P_1d) matrix expressing the values of nodal 390b11c1e72Sjeremylt basis functions at quadrature points 391*d1d35e2fSjeremylt @param grad_1d Row-major (Q_1d * P_1d) matrix expressing derivatives of nodal 392b11c1e72Sjeremylt basis functions at quadrature points 393*d1d35e2fSjeremylt @param q_ref_1d Array of length Q_1d holding the locations of quadrature points 394b11c1e72Sjeremylt on the 1D reference element [-1, 1] 395*d1d35e2fSjeremylt @param q_weight_1d Array of length Q_1d holding the quadrature weights on the 396b11c1e72Sjeremylt reference element 397b11c1e72Sjeremylt @param[out] basis Address of the variable where the newly created 398b11c1e72Sjeremylt CeedBasis will be stored. 399b11c1e72Sjeremylt 400b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 401dfdf5a53Sjeremylt 4027a982d89SJeremy L. Thompson @ref User 403b11c1e72Sjeremylt **/ 404*d1d35e2fSjeremylt int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt num_comp, 405*d1d35e2fSjeremylt CeedInt P_1d, CeedInt Q_1d, 406*d1d35e2fSjeremylt const CeedScalar *interp_1d, 407*d1d35e2fSjeremylt const CeedScalar *grad_1d, const CeedScalar *q_ref_1d, 408*d1d35e2fSjeremylt const CeedScalar *q_weight_1d, CeedBasis *basis) { 409d7b241e6Sjeremylt int ierr; 410d7b241e6Sjeremylt 4115fe0d4faSjeremylt if (!ceed->BasisCreateTensorH1) { 4125fe0d4faSjeremylt Ceed delegate; 413aefd8378Sjeremylt ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr); 4145fe0d4faSjeremylt 4155fe0d4faSjeremylt if (!delegate) 416c042f62fSJeremy L Thompson // LCOV_EXCL_START 417e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_UNSUPPORTED, 418e15f9bd0SJeremy L Thompson "Backend does not support BasisCreateTensorH1"); 419c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 4205fe0d4faSjeremylt 421*d1d35e2fSjeremylt ierr = CeedBasisCreateTensorH1(delegate, dim, num_comp, P_1d, 422*d1d35e2fSjeremylt Q_1d, interp_1d, grad_1d, q_ref_1d, 423*d1d35e2fSjeremylt q_weight_1d, basis); CeedChk(ierr); 424e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 4255fe0d4faSjeremylt } 426e15f9bd0SJeremy L Thompson 427e15f9bd0SJeremy L Thompson if (dim<1) 428e15f9bd0SJeremy L Thompson // LCOV_EXCL_START 429e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, 430e15f9bd0SJeremy L Thompson "Basis dimension must be a positive value"); 431e15f9bd0SJeremy L Thompson // LCOV_EXCL_STOP 432*d1d35e2fSjeremylt CeedElemTopology topo = dim == 1 ? CEED_LINE 433*d1d35e2fSjeremylt : dim == 2 ? CEED_QUAD 434*d1d35e2fSjeremylt : CEED_HEX; 435e15f9bd0SJeremy L Thompson 436d7b241e6Sjeremylt ierr = CeedCalloc(1, basis); CeedChk(ierr); 437d7b241e6Sjeremylt (*basis)->ceed = ceed; 438*d1d35e2fSjeremylt ceed->ref_count++; 439*d1d35e2fSjeremylt (*basis)->ref_count = 1; 440*d1d35e2fSjeremylt (*basis)->tensor_basis = 1; 441d7b241e6Sjeremylt (*basis)->dim = dim; 442d99fa3c5SJeremy L Thompson (*basis)->topo = topo; 443*d1d35e2fSjeremylt (*basis)->num_comp = num_comp; 444*d1d35e2fSjeremylt (*basis)->P_1d = P_1d; 445*d1d35e2fSjeremylt (*basis)->Q_1d = Q_1d; 446*d1d35e2fSjeremylt (*basis)->P = CeedIntPow(P_1d, dim); 447*d1d35e2fSjeremylt (*basis)->Q = CeedIntPow(Q_1d, dim); 448*d1d35e2fSjeremylt ierr = CeedMalloc(Q_1d,&(*basis)->q_ref_1d); CeedChk(ierr); 449*d1d35e2fSjeremylt ierr = CeedMalloc(Q_1d,&(*basis)->q_weight_1d); CeedChk(ierr); 450*d1d35e2fSjeremylt memcpy((*basis)->q_ref_1d, q_ref_1d, Q_1d*sizeof(q_ref_1d[0])); 451*d1d35e2fSjeremylt memcpy((*basis)->q_weight_1d, q_weight_1d, Q_1d*sizeof(q_weight_1d[0])); 452*d1d35e2fSjeremylt ierr = CeedMalloc(Q_1d*P_1d,&(*basis)->interp_1d); CeedChk(ierr); 453*d1d35e2fSjeremylt ierr = CeedMalloc(Q_1d*P_1d,&(*basis)->grad_1d); CeedChk(ierr); 454*d1d35e2fSjeremylt memcpy((*basis)->interp_1d, interp_1d, Q_1d*P_1d*sizeof(interp_1d[0])); 455*d1d35e2fSjeremylt memcpy((*basis)->grad_1d, grad_1d, Q_1d*P_1d*sizeof(grad_1d[0])); 456*d1d35e2fSjeremylt ierr = ceed->BasisCreateTensorH1(dim, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, 457*d1d35e2fSjeremylt q_weight_1d, *basis); CeedChk(ierr); 458e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 459d7b241e6Sjeremylt } 460d7b241e6Sjeremylt 461b11c1e72Sjeremylt /** 46295bb1877Svaleriabarra @brief Create a tensor-product Lagrange basis 463b11c1e72Sjeremylt 464b11c1e72Sjeremylt @param ceed A Ceed object where the CeedBasis will be created 465b11c1e72Sjeremylt @param dim Topological dimension of element 466*d1d35e2fSjeremylt @param num_comp Number of field components (1 for scalar fields) 467b11c1e72Sjeremylt @param P Number of Gauss-Lobatto nodes in one dimension. The 468b11c1e72Sjeremylt polynomial degree of the resulting Q_k element is k=P-1. 469b11c1e72Sjeremylt @param Q Number of quadrature points in one dimension. 470*d1d35e2fSjeremylt @param quad_mode Distribution of the Q quadrature points (affects order of 471b11c1e72Sjeremylt accuracy for the quadrature) 472b11c1e72Sjeremylt @param[out] basis Address of the variable where the newly created 473b11c1e72Sjeremylt CeedBasis will be stored. 474b11c1e72Sjeremylt 475b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 476dfdf5a53Sjeremylt 4777a982d89SJeremy L. Thompson @ref User 478b11c1e72Sjeremylt **/ 479*d1d35e2fSjeremylt int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt num_comp, 480*d1d35e2fSjeremylt CeedInt P, CeedInt Q, CeedQuadMode quad_mode, 481692c2638Sjeremylt CeedBasis *basis) { 482d7b241e6Sjeremylt // Allocate 483e15f9bd0SJeremy L Thompson int ierr, ierr2, i, j, k; 484*d1d35e2fSjeremylt CeedScalar c1, c2, c3, c4, dx, *nodes, *interp_1d, *grad_1d, *q_ref_1d, 485*d1d35e2fSjeremylt *q_weight_1d; 4864d537eeaSYohann 4874d537eeaSYohann if (dim<1) 488c042f62fSJeremy L Thompson // LCOV_EXCL_START 489e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_DIMENSION, 490e15f9bd0SJeremy L Thompson "Basis dimension must be a positive value"); 491c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 4924d537eeaSYohann 493e15f9bd0SJeremy L Thompson // Get Nodes and Weights 494*d1d35e2fSjeremylt ierr = CeedCalloc(P*Q, &interp_1d); CeedChk(ierr); 495*d1d35e2fSjeremylt ierr = CeedCalloc(P*Q, &grad_1d); CeedChk(ierr); 496d7b241e6Sjeremylt ierr = CeedCalloc(P, &nodes); CeedChk(ierr); 497*d1d35e2fSjeremylt ierr = CeedCalloc(Q, &q_ref_1d); CeedChk(ierr); 498*d1d35e2fSjeremylt ierr = CeedCalloc(Q, &q_weight_1d); CeedChk(ierr); 499e15f9bd0SJeremy L Thompson ierr = CeedLobattoQuadrature(P, nodes, NULL); 500e15f9bd0SJeremy L Thompson if (ierr) { goto cleanup; } CeedChk(ierr); 501*d1d35e2fSjeremylt switch (quad_mode) { 502d7b241e6Sjeremylt case CEED_GAUSS: 503*d1d35e2fSjeremylt ierr = CeedGaussQuadrature(Q, q_ref_1d, q_weight_1d); 504d7b241e6Sjeremylt break; 505d7b241e6Sjeremylt case CEED_GAUSS_LOBATTO: 506*d1d35e2fSjeremylt ierr = CeedLobattoQuadrature(Q, q_ref_1d, q_weight_1d); 507d7b241e6Sjeremylt break; 508d7b241e6Sjeremylt } 509e15f9bd0SJeremy L Thompson if (ierr) { goto cleanup; } CeedChk(ierr); 510e15f9bd0SJeremy L Thompson 511d7b241e6Sjeremylt // Build B, D matrix 512d7b241e6Sjeremylt // Fornberg, 1998 513d7b241e6Sjeremylt for (i = 0; i < Q; i++) { 514d7b241e6Sjeremylt c1 = 1.0; 515*d1d35e2fSjeremylt c3 = nodes[0] - q_ref_1d[i]; 516*d1d35e2fSjeremylt interp_1d[i*P+0] = 1.0; 517d7b241e6Sjeremylt for (j = 1; j < P; j++) { 518d7b241e6Sjeremylt c2 = 1.0; 519d7b241e6Sjeremylt c4 = c3; 520*d1d35e2fSjeremylt c3 = nodes[j] - q_ref_1d[i]; 521d7b241e6Sjeremylt for (k = 0; k < j; k++) { 522d7b241e6Sjeremylt dx = nodes[j] - nodes[k]; 523d7b241e6Sjeremylt c2 *= dx; 524d7b241e6Sjeremylt if (k == j - 1) { 525*d1d35e2fSjeremylt grad_1d[i*P + j] = c1*(interp_1d[i*P + k] - c4*grad_1d[i*P + k]) / c2; 526*d1d35e2fSjeremylt interp_1d[i*P + j] = - c1*c4*interp_1d[i*P + k] / c2; 527d7b241e6Sjeremylt } 528*d1d35e2fSjeremylt grad_1d[i*P + k] = (c3*grad_1d[i*P + k] - interp_1d[i*P + k]) / dx; 529*d1d35e2fSjeremylt interp_1d[i*P + k] = c3*interp_1d[i*P + k] / dx; 530d7b241e6Sjeremylt } 531d7b241e6Sjeremylt c1 = c2; 532d7b241e6Sjeremylt } 533d7b241e6Sjeremylt } 534d7b241e6Sjeremylt // // Pass to CeedBasisCreateTensorH1 535*d1d35e2fSjeremylt ierr = CeedBasisCreateTensorH1(ceed, dim, num_comp, P, Q, interp_1d, grad_1d, 536*d1d35e2fSjeremylt q_ref_1d, 537*d1d35e2fSjeremylt q_weight_1d, basis); CeedChk(ierr); 538e15f9bd0SJeremy L Thompson cleanup: 539*d1d35e2fSjeremylt ierr2 = CeedFree(&interp_1d); CeedChk(ierr2); 540*d1d35e2fSjeremylt ierr2 = CeedFree(&grad_1d); CeedChk(ierr2); 541e15f9bd0SJeremy L Thompson ierr2 = CeedFree(&nodes); CeedChk(ierr2); 542*d1d35e2fSjeremylt ierr2 = CeedFree(&q_ref_1d); CeedChk(ierr2); 543*d1d35e2fSjeremylt ierr2 = CeedFree(&q_weight_1d); CeedChk(ierr2); 544e15f9bd0SJeremy L Thompson CeedChk(ierr); 545e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 546d7b241e6Sjeremylt } 547d7b241e6Sjeremylt 548b11c1e72Sjeremylt /** 54995bb1877Svaleriabarra @brief Create a non tensor-product basis for H^1 discretizations 550a8de75f0Sjeremylt 551a8de75f0Sjeremylt @param ceed A Ceed object where the CeedBasis will be created 552a8de75f0Sjeremylt @param topo Topology of element, e.g. hypercube, simplex, ect 553*d1d35e2fSjeremylt @param num_comp Number of field components (1 for scalar fields) 554*d1d35e2fSjeremylt @param num_nodes Total number of nodes 555*d1d35e2fSjeremylt @param num_qpts Total number of quadrature points 556*d1d35e2fSjeremylt @param interp Row-major (num_qpts * num_nodes) matrix expressing the values of 5578795c945Sjeremylt nodal basis functions at quadrature points 558*d1d35e2fSjeremylt @param grad Row-major (num_qpts * dim * num_nodes) matrix expressing 5598795c945Sjeremylt derivatives of nodal basis functions at quadrature points 560*d1d35e2fSjeremylt @param q_ref Array of length num_qpts holding the locations of quadrature 5618795c945Sjeremylt points on the reference element [-1, 1] 562*d1d35e2fSjeremylt @param q_weight Array of length num_qpts holding the quadrature weights on the 563a8de75f0Sjeremylt reference element 564a8de75f0Sjeremylt @param[out] basis Address of the variable where the newly created 565a8de75f0Sjeremylt CeedBasis will be stored. 566a8de75f0Sjeremylt 567a8de75f0Sjeremylt @return An error code: 0 - success, otherwise - failure 568a8de75f0Sjeremylt 5697a982d89SJeremy L. Thompson @ref User 570a8de75f0Sjeremylt **/ 571*d1d35e2fSjeremylt int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, 572*d1d35e2fSjeremylt CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 573*d1d35e2fSjeremylt const CeedScalar *grad, const CeedScalar *q_ref, 574*d1d35e2fSjeremylt const CeedScalar *q_weight, CeedBasis *basis) { 575a8de75f0Sjeremylt int ierr; 576*d1d35e2fSjeremylt CeedInt P = num_nodes, Q = num_qpts, dim = 0; 577a8de75f0Sjeremylt 5785fe0d4faSjeremylt if (!ceed->BasisCreateH1) { 5795fe0d4faSjeremylt Ceed delegate; 580aefd8378Sjeremylt ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr); 5815fe0d4faSjeremylt 5825fe0d4faSjeremylt if (!delegate) 583c042f62fSJeremy L Thompson // LCOV_EXCL_START 584e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_UNSUPPORTED, 585e15f9bd0SJeremy L Thompson "Backend does not support BasisCreateH1"); 586c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 5875fe0d4faSjeremylt 588*d1d35e2fSjeremylt ierr = CeedBasisCreateH1(delegate, topo, num_comp, num_nodes, 589*d1d35e2fSjeremylt num_qpts, interp, grad, q_ref, 590*d1d35e2fSjeremylt q_weight, basis); CeedChk(ierr); 591e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 5925fe0d4faSjeremylt } 5935fe0d4faSjeremylt 594a8de75f0Sjeremylt ierr = CeedCalloc(1,basis); CeedChk(ierr); 595a8de75f0Sjeremylt 596a8de75f0Sjeremylt ierr = CeedBasisGetTopologyDimension(topo, &dim); CeedChk(ierr); 597a8de75f0Sjeremylt 598a8de75f0Sjeremylt (*basis)->ceed = ceed; 599*d1d35e2fSjeremylt ceed->ref_count++; 600*d1d35e2fSjeremylt (*basis)->ref_count = 1; 601*d1d35e2fSjeremylt (*basis)->tensor_basis = 0; 602a8de75f0Sjeremylt (*basis)->dim = dim; 603d99fa3c5SJeremy L Thompson (*basis)->topo = topo; 604*d1d35e2fSjeremylt (*basis)->num_comp = num_comp; 605a8de75f0Sjeremylt (*basis)->P = P; 606a8de75f0Sjeremylt (*basis)->Q = Q; 607*d1d35e2fSjeremylt ierr = CeedMalloc(Q*dim,&(*basis)->q_ref_1d); CeedChk(ierr); 608*d1d35e2fSjeremylt ierr = CeedMalloc(Q,&(*basis)->q_weight_1d); CeedChk(ierr); 609*d1d35e2fSjeremylt memcpy((*basis)->q_ref_1d, q_ref, Q*dim*sizeof(q_ref[0])); 610*d1d35e2fSjeremylt memcpy((*basis)->q_weight_1d, q_weight, Q*sizeof(q_weight[0])); 61100f91b2bSjeremylt ierr = CeedMalloc(Q*P, &(*basis)->interp); CeedChk(ierr); 61200f91b2bSjeremylt ierr = CeedMalloc(dim*Q*P, &(*basis)->grad); CeedChk(ierr); 61300f91b2bSjeremylt memcpy((*basis)->interp, interp, Q*P*sizeof(interp[0])); 61400f91b2bSjeremylt memcpy((*basis)->grad, grad, dim*Q*P*sizeof(grad[0])); 615*d1d35e2fSjeremylt ierr = ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, q_ref, 616*d1d35e2fSjeremylt q_weight, *basis); CeedChk(ierr); 617e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 618a8de75f0Sjeremylt } 619a8de75f0Sjeremylt 620a8de75f0Sjeremylt /** 6217a982d89SJeremy L. Thompson @brief View a CeedBasis 6227a982d89SJeremy L. Thompson 6237a982d89SJeremy L. Thompson @param basis CeedBasis to view 6247a982d89SJeremy L. Thompson @param stream Stream to view to, e.g., stdout 6257a982d89SJeremy L. Thompson 6267a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 6277a982d89SJeremy L. Thompson 6287a982d89SJeremy L. Thompson @ref User 6297a982d89SJeremy L. Thompson **/ 6307a982d89SJeremy L. Thompson int CeedBasisView(CeedBasis basis, FILE *stream) { 6317a982d89SJeremy L. Thompson int ierr; 6327a982d89SJeremy L. Thompson 633*d1d35e2fSjeremylt if (basis->tensor_basis) { 634*d1d35e2fSjeremylt fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P_1d, 635*d1d35e2fSjeremylt basis->Q_1d); 636*d1d35e2fSjeremylt ierr = CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_ref_1d, 6377a982d89SJeremy L. Thompson stream); CeedChk(ierr); 638*d1d35e2fSjeremylt ierr = CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q_1d, 639*d1d35e2fSjeremylt basis->q_weight_1d, stream); CeedChk(ierr); 640*d1d35e2fSjeremylt ierr = CeedScalarView("interp1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, 641*d1d35e2fSjeremylt basis->interp_1d, stream); CeedChk(ierr); 642*d1d35e2fSjeremylt ierr = CeedScalarView("grad1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, 643*d1d35e2fSjeremylt basis->grad_1d, stream); CeedChk(ierr); 6447a982d89SJeremy L. Thompson } else { 6457a982d89SJeremy L. Thompson fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P, 6467a982d89SJeremy L. Thompson basis->Q); 6477a982d89SJeremy L. Thompson ierr = CeedScalarView("qref", "\t% 12.8f", 1, basis->Q*basis->dim, 648*d1d35e2fSjeremylt basis->q_ref_1d, 6497a982d89SJeremy L. Thompson stream); CeedChk(ierr); 650*d1d35e2fSjeremylt ierr = CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->q_weight_1d, 6517a982d89SJeremy L. Thompson stream); CeedChk(ierr); 6527a982d89SJeremy L. Thompson ierr = CeedScalarView("interp", "\t% 12.8f", basis->Q, basis->P, 6537a982d89SJeremy L. Thompson basis->interp, stream); CeedChk(ierr); 6547a982d89SJeremy L. Thompson ierr = CeedScalarView("grad", "\t% 12.8f", basis->dim*basis->Q, basis->P, 6557a982d89SJeremy L. Thompson basis->grad, stream); CeedChk(ierr); 6567a982d89SJeremy L. Thompson } 657e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 6587a982d89SJeremy L. Thompson } 6597a982d89SJeremy L. Thompson 6607a982d89SJeremy L. Thompson /** 6617a982d89SJeremy L. Thompson @brief Apply basis evaluation from nodes to quadrature points or vice versa 6627a982d89SJeremy L. Thompson 6637a982d89SJeremy L. Thompson @param basis CeedBasis to evaluate 664*d1d35e2fSjeremylt @param num_elem The number of elements to apply the basis evaluation to; 6657a982d89SJeremy L. Thompson the backend will specify the ordering in 6664cc79fe7SJed Brown CeedElemRestrictionCreateBlocked() 667*d1d35e2fSjeremylt @param t_mode \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature 6687a982d89SJeremy L. Thompson points, \ref CEED_TRANSPOSE to apply the transpose, mapping 6697a982d89SJeremy L. Thompson from quadrature points to nodes 670*d1d35e2fSjeremylt @param eval_mode \ref CEED_EVAL_NONE to use values directly, 6717a982d89SJeremy L. Thompson \ref CEED_EVAL_INTERP to use interpolated values, 6727a982d89SJeremy L. Thompson \ref CEED_EVAL_GRAD to use gradients, 6737a982d89SJeremy L. Thompson \ref CEED_EVAL_WEIGHT to use quadrature weights. 6747a982d89SJeremy L. Thompson @param[in] u Input CeedVector 6757a982d89SJeremy L. Thompson @param[out] v Output CeedVector 6767a982d89SJeremy L. Thompson 6777a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 6787a982d89SJeremy L. Thompson 6797a982d89SJeremy L. Thompson @ref User 6807a982d89SJeremy L. Thompson **/ 681*d1d35e2fSjeremylt int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode, 682*d1d35e2fSjeremylt CeedEvalMode eval_mode, CeedVector u, CeedVector v) { 6837a982d89SJeremy L. Thompson int ierr; 684*d1d35e2fSjeremylt CeedInt u_length = 0, v_length, dim, num_comp, num_nodes, num_qpts; 685e15f9bd0SJeremy L Thompson ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr); 686*d1d35e2fSjeremylt ierr = CeedBasisGetNumComponents(basis, &num_comp); CeedChk(ierr); 687*d1d35e2fSjeremylt ierr = CeedBasisGetNumNodes(basis, &num_nodes); CeedChk(ierr); 688*d1d35e2fSjeremylt ierr = CeedBasisGetNumQuadraturePoints(basis, &num_qpts); CeedChk(ierr); 689*d1d35e2fSjeremylt ierr = CeedVectorGetLength(v, &v_length); CeedChk(ierr); 6907a982d89SJeremy L. Thompson if (u) { 691*d1d35e2fSjeremylt ierr = CeedVectorGetLength(u, &u_length); CeedChk(ierr); 6927a982d89SJeremy L. Thompson } 6937a982d89SJeremy L. Thompson 694e15f9bd0SJeremy L Thompson if (!basis->Apply) 695e15f9bd0SJeremy L Thompson // LCOV_EXCL_START 696e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_UNSUPPORTED, 697e15f9bd0SJeremy L Thompson "Backend does not support BasisApply"); 698e15f9bd0SJeremy L Thompson // LCOV_EXCL_STOP 699e15f9bd0SJeremy L Thompson 700e15f9bd0SJeremy L Thompson // Check compatibility of topological and geometrical dimensions 701*d1d35e2fSjeremylt if ((t_mode == CEED_TRANSPOSE && (v_length%num_nodes != 0 || 702*d1d35e2fSjeremylt u_length%num_qpts != 0)) || 703*d1d35e2fSjeremylt (t_mode == CEED_NOTRANSPOSE && (u_length%num_nodes != 0 || 704*d1d35e2fSjeremylt v_length%num_qpts != 0))) 705e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_DIMENSION, 706e15f9bd0SJeremy L Thompson "Length of input/output vectors " 7077a982d89SJeremy L. Thompson "incompatible with basis dimensions"); 7087a982d89SJeremy L. Thompson 709e15f9bd0SJeremy L Thompson // Check vector lengths to prevent out of bounds issues 710*d1d35e2fSjeremylt bool bad_dims = false; 711*d1d35e2fSjeremylt switch (eval_mode) { 712e15f9bd0SJeremy L Thompson case CEED_EVAL_NONE: 713*d1d35e2fSjeremylt case CEED_EVAL_INTERP: bad_dims = 714*d1d35e2fSjeremylt ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts || 715*d1d35e2fSjeremylt v_length < num_elem*num_comp*num_nodes)) || 716*d1d35e2fSjeremylt (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp || 717*d1d35e2fSjeremylt u_length < num_elem*num_comp*num_nodes))); 718e15f9bd0SJeremy L Thompson break; 719*d1d35e2fSjeremylt case CEED_EVAL_GRAD: bad_dims = 720*d1d35e2fSjeremylt ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts*dim || 721*d1d35e2fSjeremylt v_length < num_elem*num_comp*num_nodes)) || 722*d1d35e2fSjeremylt (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp*dim || 723*d1d35e2fSjeremylt u_length < num_elem*num_comp*num_nodes))); 724e15f9bd0SJeremy L Thompson break; 725e15f9bd0SJeremy L Thompson case CEED_EVAL_WEIGHT: 726*d1d35e2fSjeremylt bad_dims = v_length < num_elem*num_qpts; 727e15f9bd0SJeremy L Thompson break; 728e15f9bd0SJeremy L Thompson // LCOV_EXCL_START 729*d1d35e2fSjeremylt case CEED_EVAL_DIV: bad_dims = 730*d1d35e2fSjeremylt ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts || 731*d1d35e2fSjeremylt v_length < num_elem*num_comp*num_nodes)) || 732*d1d35e2fSjeremylt (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp || 733*d1d35e2fSjeremylt u_length < num_elem*num_comp*num_nodes))); 734e15f9bd0SJeremy L Thompson break; 735*d1d35e2fSjeremylt case CEED_EVAL_CURL: bad_dims = 736*d1d35e2fSjeremylt ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts || 737*d1d35e2fSjeremylt v_length < num_elem*num_comp*num_nodes)) || 738*d1d35e2fSjeremylt (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp || 739*d1d35e2fSjeremylt u_length < num_elem*num_comp*num_nodes))); 740e15f9bd0SJeremy L Thompson break; 741e15f9bd0SJeremy L Thompson // LCOV_EXCL_STOP 742e15f9bd0SJeremy L Thompson } 743*d1d35e2fSjeremylt if (bad_dims) 744e15f9bd0SJeremy L Thompson // LCOV_EXCL_START 745e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_DIMENSION, 746*d1d35e2fSjeremylt "Input/output vectors too short for basis and evaluation mode"); 747e15f9bd0SJeremy L Thompson // LCOV_EXCL_STOP 748e15f9bd0SJeremy L Thompson 749*d1d35e2fSjeremylt ierr = basis->Apply(basis, num_elem, t_mode, eval_mode, u, v); CeedChk(ierr); 750e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 7517a982d89SJeremy L. Thompson } 7527a982d89SJeremy L. Thompson 7537a982d89SJeremy L. Thompson /** 7549d007619Sjeremylt @brief Get dimension for given CeedBasis 7559d007619Sjeremylt 7569d007619Sjeremylt @param basis CeedBasis 7579d007619Sjeremylt @param[out] dim Variable to store dimension of basis 7589d007619Sjeremylt 7599d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 7609d007619Sjeremylt 7619d007619Sjeremylt @ref Backend 7629d007619Sjeremylt **/ 7639d007619Sjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) { 7649d007619Sjeremylt *dim = basis->dim; 765e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 7669d007619Sjeremylt } 7679d007619Sjeremylt 7689d007619Sjeremylt /** 769d99fa3c5SJeremy L Thompson @brief Get topology for given CeedBasis 770d99fa3c5SJeremy L Thompson 771d99fa3c5SJeremy L Thompson @param basis CeedBasis 772d99fa3c5SJeremy L Thompson @param[out] topo Variable to store topology of basis 773d99fa3c5SJeremy L Thompson 774d99fa3c5SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 775d99fa3c5SJeremy L Thompson 776d99fa3c5SJeremy L Thompson @ref Backend 777d99fa3c5SJeremy L Thompson **/ 778d99fa3c5SJeremy L Thompson int CeedBasisGetTopology(CeedBasis basis, CeedElemTopology *topo) { 779d99fa3c5SJeremy L Thompson *topo = basis->topo; 780e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 781d99fa3c5SJeremy L Thompson } 782d99fa3c5SJeremy L Thompson 783d99fa3c5SJeremy L Thompson /** 7849d007619Sjeremylt @brief Get number of components for given CeedBasis 7859d007619Sjeremylt 7869d007619Sjeremylt @param basis CeedBasis 787*d1d35e2fSjeremylt @param[out] num_comp Variable to store number of components of basis 7889d007619Sjeremylt 7899d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 7909d007619Sjeremylt 7919d007619Sjeremylt @ref Backend 7929d007619Sjeremylt **/ 793*d1d35e2fSjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *num_comp) { 794*d1d35e2fSjeremylt *num_comp = basis->num_comp; 795e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 7969d007619Sjeremylt } 7979d007619Sjeremylt 7989d007619Sjeremylt /** 7999d007619Sjeremylt @brief Get total number of nodes (in dim dimensions) of a CeedBasis 8009d007619Sjeremylt 8019d007619Sjeremylt @param basis CeedBasis 8029d007619Sjeremylt @param[out] P Variable to store number of nodes 8039d007619Sjeremylt 8049d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 8059d007619Sjeremylt 8069d007619Sjeremylt @ref Utility 8079d007619Sjeremylt **/ 8089d007619Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) { 8099d007619Sjeremylt *P = basis->P; 810e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 8119d007619Sjeremylt } 8129d007619Sjeremylt 8139d007619Sjeremylt /** 8149d007619Sjeremylt @brief Get total number of nodes (in 1 dimension) of a CeedBasis 8159d007619Sjeremylt 8169d007619Sjeremylt @param basis CeedBasis 817*d1d35e2fSjeremylt @param[out] P_1d Variable to store number of nodes 8189d007619Sjeremylt 8199d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 8209d007619Sjeremylt 8219d007619Sjeremylt @ref Backend 8229d007619Sjeremylt **/ 823*d1d35e2fSjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P_1d) { 824*d1d35e2fSjeremylt if (!basis->tensor_basis) 8259d007619Sjeremylt // LCOV_EXCL_START 826e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_MINOR, 827*d1d35e2fSjeremylt "Cannot supply P_1d for non-tensor basis"); 8289d007619Sjeremylt // LCOV_EXCL_STOP 8299d007619Sjeremylt 830*d1d35e2fSjeremylt *P_1d = basis->P_1d; 831e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 8329d007619Sjeremylt } 8339d007619Sjeremylt 8349d007619Sjeremylt /** 8359d007619Sjeremylt @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis 8369d007619Sjeremylt 8379d007619Sjeremylt @param basis CeedBasis 8389d007619Sjeremylt @param[out] Q Variable to store number of quadrature points 8399d007619Sjeremylt 8409d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 8419d007619Sjeremylt 8429d007619Sjeremylt @ref Utility 8439d007619Sjeremylt **/ 8449d007619Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) { 8459d007619Sjeremylt *Q = basis->Q; 846e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 8479d007619Sjeremylt } 8489d007619Sjeremylt 8499d007619Sjeremylt /** 8509d007619Sjeremylt @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis 8519d007619Sjeremylt 8529d007619Sjeremylt @param basis CeedBasis 853*d1d35e2fSjeremylt @param[out] Q_1d Variable to store number of quadrature points 8549d007619Sjeremylt 8559d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 8569d007619Sjeremylt 8579d007619Sjeremylt @ref Backend 8589d007619Sjeremylt **/ 859*d1d35e2fSjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q_1d) { 860*d1d35e2fSjeremylt if (!basis->tensor_basis) 8619d007619Sjeremylt // LCOV_EXCL_START 862e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_MINOR, 863*d1d35e2fSjeremylt "Cannot supply Q_1d for non-tensor basis"); 8649d007619Sjeremylt // LCOV_EXCL_STOP 8659d007619Sjeremylt 866*d1d35e2fSjeremylt *Q_1d = basis->Q_1d; 867e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 8689d007619Sjeremylt } 8699d007619Sjeremylt 8709d007619Sjeremylt /** 8719d007619Sjeremylt @brief Get reference coordinates of quadrature points (in dim dimensions) 8729d007619Sjeremylt of a CeedBasis 8739d007619Sjeremylt 8749d007619Sjeremylt @param basis CeedBasis 875*d1d35e2fSjeremylt @param[out] q_ref Variable to store reference coordinates of quadrature points 8769d007619Sjeremylt 8779d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 8789d007619Sjeremylt 8799d007619Sjeremylt @ref Backend 8809d007619Sjeremylt **/ 881*d1d35e2fSjeremylt int CeedBasisGetQRef(CeedBasis basis, const CeedScalar **q_ref) { 882*d1d35e2fSjeremylt *q_ref = basis->q_ref_1d; 883e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 8849d007619Sjeremylt } 8859d007619Sjeremylt 8869d007619Sjeremylt /** 8879d007619Sjeremylt @brief Get quadrature weights of quadrature points (in dim dimensions) 8889d007619Sjeremylt of a CeedBasis 8899d007619Sjeremylt 8909d007619Sjeremylt @param basis CeedBasis 891*d1d35e2fSjeremylt @param[out] q_weight Variable to store quadrature weights 8929d007619Sjeremylt 8939d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 8949d007619Sjeremylt 8959d007619Sjeremylt @ref Backend 8969d007619Sjeremylt **/ 897*d1d35e2fSjeremylt int CeedBasisGetQWeights(CeedBasis basis, const CeedScalar **q_weight) { 898*d1d35e2fSjeremylt *q_weight = basis->q_weight_1d; 899e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 9009d007619Sjeremylt } 9019d007619Sjeremylt 9029d007619Sjeremylt /** 9039d007619Sjeremylt @brief Get interpolation matrix of a CeedBasis 9049d007619Sjeremylt 9059d007619Sjeremylt @param basis CeedBasis 9069d007619Sjeremylt @param[out] interp Variable to store interpolation matrix 9079d007619Sjeremylt 9089d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 9099d007619Sjeremylt 9109d007619Sjeremylt @ref Backend 9119d007619Sjeremylt **/ 9126c58de82SJeremy L Thompson int CeedBasisGetInterp(CeedBasis basis, const CeedScalar **interp) { 913*d1d35e2fSjeremylt if (!basis->interp && basis->tensor_basis) { 9149d007619Sjeremylt // Allocate 9159d007619Sjeremylt int ierr; 9169d007619Sjeremylt ierr = CeedMalloc(basis->Q*basis->P, &basis->interp); CeedChk(ierr); 9179d007619Sjeremylt 9189d007619Sjeremylt // Initialize 9199d007619Sjeremylt for (CeedInt i=0; i<basis->Q*basis->P; i++) 9209d007619Sjeremylt basis->interp[i] = 1.0; 9219d007619Sjeremylt 9229d007619Sjeremylt // Calculate 9239d007619Sjeremylt for (CeedInt d=0; d<basis->dim; d++) 9249d007619Sjeremylt for (CeedInt qpt=0; qpt<basis->Q; qpt++) 9259d007619Sjeremylt for (CeedInt node=0; node<basis->P; node++) { 926*d1d35e2fSjeremylt CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d; 927*d1d35e2fSjeremylt CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d; 928*d1d35e2fSjeremylt basis->interp[qpt*(basis->P)+node] *= basis->interp_1d[q*basis->P_1d+p]; 9299d007619Sjeremylt } 9309d007619Sjeremylt } 9319d007619Sjeremylt *interp = basis->interp; 932e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 9339d007619Sjeremylt } 9349d007619Sjeremylt 9359d007619Sjeremylt /** 9369d007619Sjeremylt @brief Get 1D interpolation matrix of a tensor product CeedBasis 9379d007619Sjeremylt 9389d007619Sjeremylt @param basis CeedBasis 939*d1d35e2fSjeremylt @param[out] interp_1d Variable to store interpolation matrix 9409d007619Sjeremylt 9419d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 9429d007619Sjeremylt 9439d007619Sjeremylt @ref Backend 9449d007619Sjeremylt **/ 945*d1d35e2fSjeremylt int CeedBasisGetInterp1D(CeedBasis basis, const CeedScalar **interp_1d) { 946*d1d35e2fSjeremylt if (!basis->tensor_basis) 9479d007619Sjeremylt // LCOV_EXCL_START 948e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_MINOR, 949e15f9bd0SJeremy L Thompson "CeedBasis is not a tensor product basis."); 9509d007619Sjeremylt // LCOV_EXCL_STOP 9519d007619Sjeremylt 952*d1d35e2fSjeremylt *interp_1d = basis->interp_1d; 953e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 9549d007619Sjeremylt } 9559d007619Sjeremylt 9569d007619Sjeremylt /** 9579d007619Sjeremylt @brief Get gradient matrix of a CeedBasis 9589d007619Sjeremylt 9599d007619Sjeremylt @param basis CeedBasis 9609d007619Sjeremylt @param[out] grad Variable to store gradient matrix 9619d007619Sjeremylt 9629d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 9639d007619Sjeremylt 9649d007619Sjeremylt @ref Backend 9659d007619Sjeremylt **/ 9666c58de82SJeremy L Thompson int CeedBasisGetGrad(CeedBasis basis, const CeedScalar **grad) { 967*d1d35e2fSjeremylt if (!basis->grad && basis->tensor_basis) { 9689d007619Sjeremylt // Allocate 9699d007619Sjeremylt int ierr; 9709d007619Sjeremylt ierr = CeedMalloc(basis->dim*basis->Q*basis->P, &basis->grad); 9719d007619Sjeremylt CeedChk(ierr); 9729d007619Sjeremylt 9739d007619Sjeremylt // Initialize 9749d007619Sjeremylt for (CeedInt i=0; i<basis->dim*basis->Q*basis->P; i++) 9759d007619Sjeremylt basis->grad[i] = 1.0; 9769d007619Sjeremylt 9779d007619Sjeremylt // Calculate 9789d007619Sjeremylt for (CeedInt d=0; d<basis->dim; d++) 9799d007619Sjeremylt for (CeedInt i=0; i<basis->dim; i++) 9809d007619Sjeremylt for (CeedInt qpt=0; qpt<basis->Q; qpt++) 9819d007619Sjeremylt for (CeedInt node=0; node<basis->P; node++) { 982*d1d35e2fSjeremylt CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d; 983*d1d35e2fSjeremylt CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d; 9849d007619Sjeremylt if (i == d) 9859d007619Sjeremylt basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *= 986*d1d35e2fSjeremylt basis->grad_1d[q*basis->P_1d+p]; 9879d007619Sjeremylt else 9889d007619Sjeremylt basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *= 989*d1d35e2fSjeremylt basis->interp_1d[q*basis->P_1d+p]; 9909d007619Sjeremylt } 9919d007619Sjeremylt } 9929d007619Sjeremylt *grad = basis->grad; 993e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 9949d007619Sjeremylt } 9959d007619Sjeremylt 9969d007619Sjeremylt /** 9979d007619Sjeremylt @brief Get 1D gradient matrix of a tensor product CeedBasis 9989d007619Sjeremylt 9999d007619Sjeremylt @param basis CeedBasis 1000*d1d35e2fSjeremylt @param[out] grad_1d Variable to store gradient matrix 10019d007619Sjeremylt 10029d007619Sjeremylt @return An error code: 0 - success, otherwise - failure 10039d007619Sjeremylt 10049d007619Sjeremylt @ref Backend 10059d007619Sjeremylt **/ 1006*d1d35e2fSjeremylt int CeedBasisGetGrad1D(CeedBasis basis, const CeedScalar **grad_1d) { 1007*d1d35e2fSjeremylt if (!basis->tensor_basis) 10089d007619Sjeremylt // LCOV_EXCL_START 1009e15f9bd0SJeremy L Thompson return CeedError(basis->ceed, CEED_ERROR_MINOR, 1010e15f9bd0SJeremy L Thompson "CeedBasis is not a tensor product basis."); 10119d007619Sjeremylt // LCOV_EXCL_STOP 10129d007619Sjeremylt 1013*d1d35e2fSjeremylt *grad_1d = basis->grad_1d; 1014e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 10159d007619Sjeremylt } 10169d007619Sjeremylt 10179d007619Sjeremylt /** 10187a982d89SJeremy L. Thompson @brief Destroy a CeedBasis 10197a982d89SJeremy L. Thompson 10207a982d89SJeremy L. Thompson @param basis CeedBasis to destroy 10217a982d89SJeremy L. Thompson 10227a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 10237a982d89SJeremy L. Thompson 10247a982d89SJeremy L. Thompson @ref User 10257a982d89SJeremy L. Thompson **/ 10267a982d89SJeremy L. Thompson int CeedBasisDestroy(CeedBasis *basis) { 10277a982d89SJeremy L. Thompson int ierr; 10287a982d89SJeremy L. Thompson 1029*d1d35e2fSjeremylt if (!*basis || --(*basis)->ref_count > 0) return CEED_ERROR_SUCCESS; 10307a982d89SJeremy L. Thompson if ((*basis)->Destroy) { 10317a982d89SJeremy L. Thompson ierr = (*basis)->Destroy(*basis); CeedChk(ierr); 10327a982d89SJeremy L. Thompson } 10337a982d89SJeremy L. Thompson ierr = CeedFree(&(*basis)->interp); CeedChk(ierr); 1034*d1d35e2fSjeremylt ierr = CeedFree(&(*basis)->interp_1d); CeedChk(ierr); 10357a982d89SJeremy L. Thompson ierr = CeedFree(&(*basis)->grad); CeedChk(ierr); 1036*d1d35e2fSjeremylt ierr = CeedFree(&(*basis)->grad_1d); CeedChk(ierr); 1037*d1d35e2fSjeremylt ierr = CeedFree(&(*basis)->q_ref_1d); CeedChk(ierr); 1038*d1d35e2fSjeremylt ierr = CeedFree(&(*basis)->q_weight_1d); CeedChk(ierr); 10397a982d89SJeremy L. Thompson ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr); 10407a982d89SJeremy L. Thompson ierr = CeedFree(basis); CeedChk(ierr); 1041e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 10427a982d89SJeremy L. Thompson } 10437a982d89SJeremy L. Thompson 10447a982d89SJeremy L. Thompson /** 1045b11c1e72Sjeremylt @brief Construct a Gauss-Legendre quadrature 1046b11c1e72Sjeremylt 1047b11c1e72Sjeremylt @param Q Number of quadrature points (integrates polynomials of 1048b11c1e72Sjeremylt degree 2*Q-1 exactly) 1049*d1d35e2fSjeremylt @param[out] q_ref_1d Array of length Q to hold the abscissa on [-1, 1] 1050*d1d35e2fSjeremylt @param[out] q_weight_1d Array of length Q to hold the weights 1051b11c1e72Sjeremylt 1052b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 1053dfdf5a53Sjeremylt 1054dfdf5a53Sjeremylt @ref Utility 1055b11c1e72Sjeremylt **/ 1056*d1d35e2fSjeremylt int CeedGaussQuadrature(CeedInt Q, CeedScalar *q_ref_1d, 1057*d1d35e2fSjeremylt CeedScalar *q_weight_1d) { 1058d7b241e6Sjeremylt // Allocate 1059d7b241e6Sjeremylt CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0*atan(1.0); 1060*d1d35e2fSjeremylt // Build q_ref_1d, q_weight_1d 1061d7b241e6Sjeremylt for (int i = 0; i <= Q/2; i++) { 1062d7b241e6Sjeremylt // Guess 1063d7b241e6Sjeremylt xi = cos(PI*(CeedScalar)(2*i+1)/((CeedScalar)(2*Q))); 1064d7b241e6Sjeremylt // Pn(xi) 1065d7b241e6Sjeremylt P0 = 1.0; 1066d7b241e6Sjeremylt P1 = xi; 1067d7b241e6Sjeremylt P2 = 0.0; 1068d7b241e6Sjeremylt for (int j = 2; j <= Q; j++) { 1069d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 1070d7b241e6Sjeremylt P0 = P1; 1071d7b241e6Sjeremylt P1 = P2; 1072d7b241e6Sjeremylt } 1073d7b241e6Sjeremylt // First Newton Step 1074d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 1075d7b241e6Sjeremylt xi = xi-P2/dP2; 1076d7b241e6Sjeremylt // Newton to convergence 10770e4d4210Sjeremylt for (int k=0; k<100 && fabs(P2)>10*CEED_EPSILON; k++) { 1078d7b241e6Sjeremylt P0 = 1.0; 1079d7b241e6Sjeremylt P1 = xi; 1080d7b241e6Sjeremylt for (int j = 2; j <= Q; j++) { 1081d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 1082d7b241e6Sjeremylt P0 = P1; 1083d7b241e6Sjeremylt P1 = P2; 1084d7b241e6Sjeremylt } 1085d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 1086d7b241e6Sjeremylt xi = xi-P2/dP2; 1087d7b241e6Sjeremylt } 1088d7b241e6Sjeremylt // Save xi, wi 1089d7b241e6Sjeremylt wi = 2.0/((1.0-xi*xi)*dP2*dP2); 1090*d1d35e2fSjeremylt q_weight_1d[i] = wi; 1091*d1d35e2fSjeremylt q_weight_1d[Q-1-i] = wi; 1092*d1d35e2fSjeremylt q_ref_1d[i] = -xi; 1093*d1d35e2fSjeremylt q_ref_1d[Q-1-i]= xi; 1094d7b241e6Sjeremylt } 1095e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1096d7b241e6Sjeremylt } 1097d7b241e6Sjeremylt 1098b11c1e72Sjeremylt /** 1099b11c1e72Sjeremylt @brief Construct a Gauss-Legendre-Lobatto quadrature 1100b11c1e72Sjeremylt 1101b11c1e72Sjeremylt @param Q Number of quadrature points (integrates polynomials of 1102b11c1e72Sjeremylt degree 2*Q-3 exactly) 1103*d1d35e2fSjeremylt @param[out] q_ref_1d Array of length Q to hold the abscissa on [-1, 1] 1104*d1d35e2fSjeremylt @param[out] q_weight_1d Array of length Q to hold the weights 1105b11c1e72Sjeremylt 1106b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 1107dfdf5a53Sjeremylt 1108dfdf5a53Sjeremylt @ref Utility 1109b11c1e72Sjeremylt **/ 1110*d1d35e2fSjeremylt int CeedLobattoQuadrature(CeedInt Q, CeedScalar *q_ref_1d, 1111*d1d35e2fSjeremylt CeedScalar *q_weight_1d) { 1112d7b241e6Sjeremylt // Allocate 1113d7b241e6Sjeremylt CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0*atan(1.0); 1114*d1d35e2fSjeremylt // Build q_ref_1d, q_weight_1d 1115d7b241e6Sjeremylt // Set endpoints 111630a100c3SJed Brown if (Q < 2) 1117b0d62198Sjeremylt // LCOV_EXCL_START 1118e15f9bd0SJeremy L Thompson return CeedError(NULL, CEED_ERROR_DIMENSION, 11197ed52d01Sjeremylt "Cannot create Lobatto quadrature with Q=%d < 2 points", Q); 1120b0d62198Sjeremylt // LCOV_EXCL_STOP 1121d7b241e6Sjeremylt wi = 2.0/((CeedScalar)(Q*(Q-1))); 1122*d1d35e2fSjeremylt if (q_weight_1d) { 1123*d1d35e2fSjeremylt q_weight_1d[0] = wi; 1124*d1d35e2fSjeremylt q_weight_1d[Q-1] = wi; 1125d7b241e6Sjeremylt } 1126*d1d35e2fSjeremylt q_ref_1d[0] = -1.0; 1127*d1d35e2fSjeremylt q_ref_1d[Q-1] = 1.0; 1128d7b241e6Sjeremylt // Interior 1129d7b241e6Sjeremylt for (int i = 1; i <= (Q-1)/2; i++) { 1130d7b241e6Sjeremylt // Guess 1131d7b241e6Sjeremylt xi = cos(PI*(CeedScalar)(i)/(CeedScalar)(Q-1)); 1132d7b241e6Sjeremylt // Pn(xi) 1133d7b241e6Sjeremylt P0 = 1.0; 1134d7b241e6Sjeremylt P1 = xi; 1135d7b241e6Sjeremylt P2 = 0.0; 1136d7b241e6Sjeremylt for (int j = 2; j < Q; j++) { 1137d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 1138d7b241e6Sjeremylt P0 = P1; 1139d7b241e6Sjeremylt P1 = P2; 1140d7b241e6Sjeremylt } 1141d7b241e6Sjeremylt // First Newton step 1142d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 1143d7b241e6Sjeremylt d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi); 1144d7b241e6Sjeremylt xi = xi-dP2/d2P2; 1145d7b241e6Sjeremylt // Newton to convergence 11460e4d4210Sjeremylt for (int k=0; k<100 && fabs(dP2)>10*CEED_EPSILON; k++) { 1147d7b241e6Sjeremylt P0 = 1.0; 1148d7b241e6Sjeremylt P1 = xi; 1149d7b241e6Sjeremylt for (int j = 2; j < Q; j++) { 1150d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 1151d7b241e6Sjeremylt P0 = P1; 1152d7b241e6Sjeremylt P1 = P2; 1153d7b241e6Sjeremylt } 1154d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 1155d7b241e6Sjeremylt d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi); 1156d7b241e6Sjeremylt xi = xi-dP2/d2P2; 1157d7b241e6Sjeremylt } 1158d7b241e6Sjeremylt // Save xi, wi 1159d7b241e6Sjeremylt wi = 2.0/(((CeedScalar)(Q*(Q-1)))*P2*P2); 1160*d1d35e2fSjeremylt if (q_weight_1d) { 1161*d1d35e2fSjeremylt q_weight_1d[i] = wi; 1162*d1d35e2fSjeremylt q_weight_1d[Q-1-i] = wi; 1163d7b241e6Sjeremylt } 1164*d1d35e2fSjeremylt q_ref_1d[i] = -xi; 1165*d1d35e2fSjeremylt q_ref_1d[Q-1-i]= xi; 1166d7b241e6Sjeremylt } 1167e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1168d7b241e6Sjeremylt } 1169d7b241e6Sjeremylt 1170dfdf5a53Sjeremylt /** 117195bb1877Svaleriabarra @brief Return QR Factorization of a matrix 1172b11c1e72Sjeremylt 117377645d7bSjeremylt @param ceed A Ceed context for error handling 117452bfb9bbSJeremy L Thompson @param[in,out] mat Row-major matrix to be factorized in place 117552bfb9bbSJeremy L Thompson @param[in,out] tau Vector of length m of scaling factors 1176b11c1e72Sjeremylt @param m Number of rows 1177b11c1e72Sjeremylt @param n Number of columns 1178b11c1e72Sjeremylt 1179b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 1180dfdf5a53Sjeremylt 1181dfdf5a53Sjeremylt @ref Utility 1182b11c1e72Sjeremylt **/ 1183a7bd39daSjeremylt int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau, 1184d7b241e6Sjeremylt CeedInt m, CeedInt n) { 1185d7b241e6Sjeremylt CeedScalar v[m]; 1186d7b241e6Sjeremylt 1187a7bd39daSjeremylt // Check m >= n 1188a7bd39daSjeremylt if (n > m) 1189c042f62fSJeremy L Thompson // LCOV_EXCL_START 1190e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_UNSUPPORTED, 1191e15f9bd0SJeremy L Thompson "Cannot compute QR factorization with n > m"); 1192c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 1193a7bd39daSjeremylt 119452bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) { 1195d7b241e6Sjeremylt // Calculate Householder vector, magnitude 1196d7b241e6Sjeremylt CeedScalar sigma = 0.0; 1197d7b241e6Sjeremylt v[i] = mat[i+n*i]; 119852bfb9bbSJeremy L Thompson for (CeedInt j=i+1; j<m; j++) { 1199d7b241e6Sjeremylt v[j] = mat[i+n*j]; 1200d7b241e6Sjeremylt sigma += v[j] * v[j]; 1201d7b241e6Sjeremylt } 1202d7b241e6Sjeremylt CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m] 1203d7b241e6Sjeremylt CeedScalar Rii = -copysign(norm, v[i]); 1204d7b241e6Sjeremylt v[i] -= Rii; 1205d7b241e6Sjeremylt // norm of v[i:m] after modification above and scaling below 1206d7b241e6Sjeremylt // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 1207d7b241e6Sjeremylt // tau = 2 / (norm*norm) 1208d7b241e6Sjeremylt tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma); 1209fb551037Sjeremylt 12101d102b48SJeremy L Thompson for (CeedInt j=i+1; j<m; j++) 12111d102b48SJeremy L Thompson v[j] /= v[i]; 1212d7b241e6Sjeremylt 1213d7b241e6Sjeremylt // Apply Householder reflector to lower right panel 1214d7b241e6Sjeremylt CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1); 1215d7b241e6Sjeremylt // Save v 1216d7b241e6Sjeremylt mat[i+n*i] = Rii; 12171d102b48SJeremy L Thompson for (CeedInt j=i+1; j<m; j++) 1218d7b241e6Sjeremylt mat[i+n*j] = v[j]; 1219d7b241e6Sjeremylt } 1220e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1221d7b241e6Sjeremylt } 1222d7b241e6Sjeremylt 1223b11c1e72Sjeremylt /** 122452bfb9bbSJeremy L Thompson @brief Return symmetric Schur decomposition of the symmetric matrix mat via 122552bfb9bbSJeremy L Thompson symmetric QR factorization 122652bfb9bbSJeremy L Thompson 122777645d7bSjeremylt @param ceed A Ceed context for error handling 122852bfb9bbSJeremy L Thompson @param[in,out] mat Row-major matrix to be factorized in place 1229460bf743SValeria Barra @param[out] lambda Vector of length n of eigenvalues 123052bfb9bbSJeremy L Thompson @param n Number of rows/columns 123152bfb9bbSJeremy L Thompson 123252bfb9bbSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 123352bfb9bbSJeremy L Thompson 123452bfb9bbSJeremy L Thompson @ref Utility 123552bfb9bbSJeremy L Thompson **/ 123652bfb9bbSJeremy L Thompson int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat, 123752bfb9bbSJeremy L Thompson CeedScalar *lambda, CeedInt n) { 123852bfb9bbSJeremy L Thompson // Check bounds for clang-tidy 123952bfb9bbSJeremy L Thompson if (n<2) 1240c042f62fSJeremy L Thompson // LCOV_EXCL_START 1241e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_UNSUPPORTED, 1242c042f62fSJeremy L Thompson "Cannot compute symmetric Schur decomposition of scalars"); 1243c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 124452bfb9bbSJeremy L Thompson 124552bfb9bbSJeremy L Thompson CeedScalar v[n-1], tau[n-1], matT[n*n]; 124652bfb9bbSJeremy L Thompson 124752bfb9bbSJeremy L Thompson // Copy mat to matT and set mat to I 124852bfb9bbSJeremy L Thompson memcpy(matT, mat, n*n*sizeof(mat[0])); 124952bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) 125052bfb9bbSJeremy L Thompson for (CeedInt j=0; j<n; j++) 125152bfb9bbSJeremy L Thompson mat[j+n*i] = (i==j) ? 1 : 0; 125252bfb9bbSJeremy L Thompson 125352bfb9bbSJeremy L Thompson // Reduce to tridiagonal 125452bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n-1; i++) { 125552bfb9bbSJeremy L Thompson // Calculate Householder vector, magnitude 125652bfb9bbSJeremy L Thompson CeedScalar sigma = 0.0; 125752bfb9bbSJeremy L Thompson v[i] = matT[i+n*(i+1)]; 125852bfb9bbSJeremy L Thompson for (CeedInt j=i+1; j<n-1; j++) { 125952bfb9bbSJeremy L Thompson v[j] = matT[i+n*(j+1)]; 126052bfb9bbSJeremy L Thompson sigma += v[j] * v[j]; 126152bfb9bbSJeremy L Thompson } 126252bfb9bbSJeremy L Thompson CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:n-1] 126352bfb9bbSJeremy L Thompson CeedScalar Rii = -copysign(norm, v[i]); 126452bfb9bbSJeremy L Thompson v[i] -= Rii; 126552bfb9bbSJeremy L Thompson // norm of v[i:m] after modification above and scaling below 126652bfb9bbSJeremy L Thompson // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 126752bfb9bbSJeremy L Thompson // tau = 2 / (norm*norm) 12680e4d4210Sjeremylt if (sigma > 10*CEED_EPSILON) 126952bfb9bbSJeremy L Thompson tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma); 1270fb551037Sjeremylt else 1271fb551037Sjeremylt tau[i] = 0; 1272fb551037Sjeremylt 1273fb551037Sjeremylt for (CeedInt j=i+1; j<n-1; j++) 1274fb551037Sjeremylt v[j] /= v[i]; 127552bfb9bbSJeremy L Thompson 127652bfb9bbSJeremy L Thompson // Update sub and super diagonal 127752bfb9bbSJeremy L Thompson matT[i+n*(i+1)] = Rii; 127852bfb9bbSJeremy L Thompson matT[(i+1)+n*i] = Rii; 127952bfb9bbSJeremy L Thompson for (CeedInt j=i+2; j<n; j++) { 128052bfb9bbSJeremy L Thompson matT[i+n*j] = 0; matT[j+n*i] = 0; 128152bfb9bbSJeremy L Thompson } 128252bfb9bbSJeremy L Thompson // Apply symmetric Householder reflector to lower right panel 128352bfb9bbSJeremy L Thompson CeedHouseholderReflect(&matT[(i+1)+n*(i+1)], &v[i], tau[i], 128452bfb9bbSJeremy L Thompson n-(i+1), n-(i+1), n, 1); 128552bfb9bbSJeremy L Thompson CeedHouseholderReflect(&matT[(i+1)+n*(i+1)], &v[i], tau[i], 128652bfb9bbSJeremy L Thompson n-(i+1), n-(i+1), 1, n); 128752bfb9bbSJeremy L Thompson // Save v 128852bfb9bbSJeremy L Thompson for (CeedInt j=i+1; j<n-1; j++) { 128952bfb9bbSJeremy L Thompson matT[i+n*(j+1)] = v[j]; 129052bfb9bbSJeremy L Thompson } 129152bfb9bbSJeremy L Thompson } 129252bfb9bbSJeremy L Thompson // Backwards accumulation of Q 129352bfb9bbSJeremy L Thompson for (CeedInt i=n-2; i>=0; i--) { 129452bfb9bbSJeremy L Thompson v[i] = 1; 129552bfb9bbSJeremy L Thompson for (CeedInt j=i+1; j<n-1; j++) { 129652bfb9bbSJeremy L Thompson v[j] = matT[i+n*(j+1)]; 129752bfb9bbSJeremy L Thompson matT[i+n*(j+1)] = 0; 129852bfb9bbSJeremy L Thompson } 129952bfb9bbSJeremy L Thompson CeedHouseholderReflect(&mat[(i+1)+n*(i+1)], &v[i], tau[i], 130052bfb9bbSJeremy L Thompson n-(i+1), n-(i+1), n, 1); 130152bfb9bbSJeremy L Thompson } 130252bfb9bbSJeremy L Thompson 130352bfb9bbSJeremy L Thompson // Reduce sub and super diagonal 130452bfb9bbSJeremy L Thompson CeedInt p = 0, q = 0, itr = 0, maxitr = n*n*n; 13050e4d4210Sjeremylt CeedScalar tol = 10*CEED_EPSILON; 130652bfb9bbSJeremy L Thompson 130752bfb9bbSJeremy L Thompson while (q < n && itr < maxitr) { 130852bfb9bbSJeremy L Thompson // Update p, q, size of reduced portions of diagonal 130952bfb9bbSJeremy L Thompson p = 0; q = 0; 131052bfb9bbSJeremy L Thompson for (CeedInt i=n-2; i>=0; i--) { 131152bfb9bbSJeremy L Thompson if (fabs(matT[i+n*(i+1)]) < tol) 131252bfb9bbSJeremy L Thompson q += 1; 131352bfb9bbSJeremy L Thompson else 131452bfb9bbSJeremy L Thompson break; 131552bfb9bbSJeremy L Thompson } 131652bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n-1-q; i++) { 131752bfb9bbSJeremy L Thompson if (fabs(matT[i+n*(i+1)]) < tol) 131852bfb9bbSJeremy L Thompson p += 1; 131952bfb9bbSJeremy L Thompson else 132052bfb9bbSJeremy L Thompson break; 132152bfb9bbSJeremy L Thompson } 132252bfb9bbSJeremy L Thompson if (q == n-1) break; // Finished reducing 132352bfb9bbSJeremy L Thompson 132452bfb9bbSJeremy L Thompson // Reduce tridiagonal portion 132552bfb9bbSJeremy L Thompson CeedScalar tnn = matT[(n-1-q)+n*(n-1-q)], 132652bfb9bbSJeremy L Thompson tnnm1 = matT[(n-2-q)+n*(n-1-q)]; 132752bfb9bbSJeremy L Thompson CeedScalar d = (matT[(n-2-q)+n*(n-2-q)] - tnn)/2; 132852bfb9bbSJeremy L Thompson CeedScalar mu = tnn - tnnm1*tnnm1 / 132952bfb9bbSJeremy L Thompson (d + copysign(sqrt(d*d + tnnm1*tnnm1), d)); 133052bfb9bbSJeremy L Thompson CeedScalar x = matT[p+n*p] - mu; 133152bfb9bbSJeremy L Thompson CeedScalar z = matT[p+n*(p+1)]; 133252bfb9bbSJeremy L Thompson for (CeedInt k=p; k<n-1-q; k++) { 133352bfb9bbSJeremy L Thompson // Compute Givens rotation 133452bfb9bbSJeremy L Thompson CeedScalar c = 1, s = 0; 133552bfb9bbSJeremy L Thompson if (fabs(z) > tol) { 133652bfb9bbSJeremy L Thompson if (fabs(z) > fabs(x)) { 133752bfb9bbSJeremy L Thompson CeedScalar tau = -x/z; 133852bfb9bbSJeremy L Thompson s = 1/sqrt(1+tau*tau), c = s*tau; 133952bfb9bbSJeremy L Thompson } else { 134052bfb9bbSJeremy L Thompson CeedScalar tau = -z/x; 134152bfb9bbSJeremy L Thompson c = 1/sqrt(1+tau*tau), s = c*tau; 134252bfb9bbSJeremy L Thompson } 134352bfb9bbSJeremy L Thompson } 134452bfb9bbSJeremy L Thompson 134552bfb9bbSJeremy L Thompson // Apply Givens rotation to T 134652bfb9bbSJeremy L Thompson CeedGivensRotation(matT, c, s, CEED_NOTRANSPOSE, k, k+1, n, n); 134752bfb9bbSJeremy L Thompson CeedGivensRotation(matT, c, s, CEED_TRANSPOSE, k, k+1, n, n); 134852bfb9bbSJeremy L Thompson 134952bfb9bbSJeremy L Thompson // Apply Givens rotation to Q 135052bfb9bbSJeremy L Thompson CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k+1, n, n); 135152bfb9bbSJeremy L Thompson 135252bfb9bbSJeremy L Thompson // Update x, z 135352bfb9bbSJeremy L Thompson if (k < n-q-2) { 135452bfb9bbSJeremy L Thompson x = matT[k+n*(k+1)]; 135552bfb9bbSJeremy L Thompson z = matT[k+n*(k+2)]; 135652bfb9bbSJeremy L Thompson } 135752bfb9bbSJeremy L Thompson } 135852bfb9bbSJeremy L Thompson itr++; 135952bfb9bbSJeremy L Thompson } 136052bfb9bbSJeremy L Thompson // Save eigenvalues 136152bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) 136252bfb9bbSJeremy L Thompson lambda[i] = matT[i+n*i]; 136352bfb9bbSJeremy L Thompson 136452bfb9bbSJeremy L Thompson // Check convergence 136552bfb9bbSJeremy L Thompson if (itr == maxitr && q < n-1) 1366c042f62fSJeremy L Thompson // LCOV_EXCL_START 1367e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_MINOR, 1368e15f9bd0SJeremy L Thompson "Symmetric QR failed to converge"); 1369c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 1370e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 137152bfb9bbSJeremy L Thompson } 137252bfb9bbSJeremy L Thompson 137352bfb9bbSJeremy L Thompson /** 137452bfb9bbSJeremy L Thompson @brief Return Simultaneous Diagonalization of two matrices. This solves the 137552bfb9bbSJeremy L Thompson generalized eigenvalue problem A x = lambda B x, where A and B 137652bfb9bbSJeremy L Thompson are symmetric and B is positive definite. We generate the matrix X 137752bfb9bbSJeremy L Thompson and vector Lambda such that X^T A X = Lambda and X^T B X = I. This 137852bfb9bbSJeremy L Thompson is equivalent to the LAPACK routine 'sygv' with TYPE = 1. 137952bfb9bbSJeremy L Thompson 138077645d7bSjeremylt @param ceed A Ceed context for error handling 1381*d1d35e2fSjeremylt @param[in] mat_A Row-major matrix to be factorized with eigenvalues 1382*d1d35e2fSjeremylt @param[in] mat_B Row-major matrix to be factorized to identity 138352bfb9bbSJeremy L Thompson @param[out] x Row-major orthogonal matrix 1384460bf743SValeria Barra @param[out] lambda Vector of length n of generalized eigenvalues 138552bfb9bbSJeremy L Thompson @param n Number of rows/columns 138652bfb9bbSJeremy L Thompson 138752bfb9bbSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 138852bfb9bbSJeremy L Thompson 138952bfb9bbSJeremy L Thompson @ref Utility 139052bfb9bbSJeremy L Thompson **/ 1391*d1d35e2fSjeremylt int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *mat_A, 1392*d1d35e2fSjeremylt CeedScalar *mat_B, CeedScalar *x, 139352bfb9bbSJeremy L Thompson CeedScalar *lambda, CeedInt n) { 139452bfb9bbSJeremy L Thompson int ierr; 1395*d1d35e2fSjeremylt CeedScalar mat_C[n*n], matG[n*n], vecD[n]; 139652bfb9bbSJeremy L Thompson 139752bfb9bbSJeremy L Thompson // Compute B = G D G^T 1398*d1d35e2fSjeremylt memcpy(matG, mat_B, n*n*sizeof(mat_B[0])); 139952bfb9bbSJeremy L Thompson ierr = CeedSymmetricSchurDecomposition(ceed, matG, vecD, n); CeedChk(ierr); 1400fb551037Sjeremylt for (CeedInt i=0; i<n; i++) 1401fb551037Sjeremylt vecD[i] = sqrt(vecD[i]); 140252bfb9bbSJeremy L Thompson 1403fb551037Sjeremylt // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T 1404fb551037Sjeremylt // = D^-1/2 G^T A G D^-1/2 140552bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) 140652bfb9bbSJeremy L Thompson for (CeedInt j=0; j<n; j++) 1407*d1d35e2fSjeremylt mat_C[j+i*n] = matG[i+j*n] / vecD[i]; 1408*d1d35e2fSjeremylt ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)mat_C, 1409*d1d35e2fSjeremylt (const CeedScalar *)mat_A, x, n, n, n); 14109289e5bfSjeremylt CeedChk(ierr); 141152bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) 141252bfb9bbSJeremy L Thompson for (CeedInt j=0; j<n; j++) 1413fb551037Sjeremylt matG[j+i*n] = matG[j+i*n] / vecD[j]; 14149289e5bfSjeremylt ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)x, 1415*d1d35e2fSjeremylt (const CeedScalar *)matG, mat_C, n, n, n); 14169289e5bfSjeremylt CeedChk(ierr); 141752bfb9bbSJeremy L Thompson 141852bfb9bbSJeremy L Thompson // Compute Q^T C Q = lambda 1419*d1d35e2fSjeremylt ierr = CeedSymmetricSchurDecomposition(ceed, mat_C, lambda, n); CeedChk(ierr); 142052bfb9bbSJeremy L Thompson 1421fb551037Sjeremylt // Set x = (G D^1/2)^-T Q 1422fb551037Sjeremylt // = G D^-1/2 Q 14239289e5bfSjeremylt ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)matG, 1424*d1d35e2fSjeremylt (const CeedScalar *)mat_C, x, n, n, n); 14259289e5bfSjeremylt CeedChk(ierr); 1426e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 142752bfb9bbSJeremy L Thompson } 142852bfb9bbSJeremy L Thompson 1429d7b241e6Sjeremylt /// @} 1430