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 17d7b241e6Sjeremylt #include <ceed-impl.h> 18d863ab9bSjeremylt #include <ceed-backend.h> 19d7b241e6Sjeremylt #include <math.h> 20d7b241e6Sjeremylt #include <stdio.h> 21d7b241e6Sjeremylt #include <stdlib.h> 22d7b241e6Sjeremylt #include <string.h> 23d7b241e6Sjeremylt 247a982d89SJeremy L. Thompson /// @file 257a982d89SJeremy L. Thompson /// Implementation of CeedBasis interfaces 267a982d89SJeremy L. Thompson 27d7b241e6Sjeremylt /// @cond DOXYGEN_SKIP 28783c99b3SValeria Barra static struct CeedBasis_private ceed_basis_collocated; 29d7b241e6Sjeremylt /// @endcond 30d7b241e6Sjeremylt 317a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser 327a982d89SJeremy L. Thompson /// @{ 337a982d89SJeremy L. Thompson 347a982d89SJeremy L. Thompson /// Indicate that the quadrature points are collocated with the nodes 357a982d89SJeremy L. Thompson const CeedBasis CEED_BASIS_COLLOCATED = &ceed_basis_collocated; 367a982d89SJeremy L. Thompson 377a982d89SJeremy L. Thompson /// @} 387a982d89SJeremy L. Thompson 397a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 407a982d89SJeremy L. Thompson /// CeedBasis Library Internal Functions 417a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 427a982d89SJeremy L. Thompson /// @addtogroup CeedBasisDeveloper 437a982d89SJeremy L. Thompson /// @{ 447a982d89SJeremy L. Thompson 457a982d89SJeremy L. Thompson /** 467a982d89SJeremy L. Thompson @brief Compute Householder reflection 477a982d89SJeremy L. Thompson 487a982d89SJeremy L. Thompson Computes A = (I - b v v^T) A 497a982d89SJeremy L. Thompson where A is an mxn matrix indexed as A[i*row + j*col] 507a982d89SJeremy L. Thompson 517a982d89SJeremy L. Thompson @param[in,out] A Matrix to apply Householder reflection to, in place 527a982d89SJeremy L. Thompson @param v Householder vector 537a982d89SJeremy L. Thompson @param b Scaling factor 547a982d89SJeremy L. Thompson @param m Number of rows in A 557a982d89SJeremy L. Thompson @param n Number of columns in A 567a982d89SJeremy L. Thompson @param row Row stride 577a982d89SJeremy L. Thompson @param col Col stride 587a982d89SJeremy L. Thompson 597a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 607a982d89SJeremy L. Thompson 617a982d89SJeremy L. Thompson @ref Developer 627a982d89SJeremy L. Thompson **/ 637a982d89SJeremy L. Thompson static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v, 647a982d89SJeremy L. Thompson CeedScalar b, CeedInt m, CeedInt n, 657a982d89SJeremy L. Thompson CeedInt row, CeedInt col) { 667a982d89SJeremy L. Thompson for (CeedInt j=0; j<n; j++) { 677a982d89SJeremy L. Thompson CeedScalar w = A[0*row + j*col]; 687a982d89SJeremy L. Thompson for (CeedInt i=1; i<m; i++) 697a982d89SJeremy L. Thompson w += v[i] * A[i*row + j*col]; 707a982d89SJeremy L. Thompson A[0*row + j*col] -= b * w; 717a982d89SJeremy L. Thompson for (CeedInt i=1; i<m; i++) 727a982d89SJeremy L. Thompson A[i*row + j*col] -= b * w * v[i]; 737a982d89SJeremy L. Thompson } 747a982d89SJeremy L. Thompson return 0; 757a982d89SJeremy L. Thompson } 767a982d89SJeremy L. Thompson 777a982d89SJeremy L. Thompson /** 787a982d89SJeremy L. Thompson @brief Apply Householder Q matrix 797a982d89SJeremy L. Thompson 807a982d89SJeremy L. Thompson Compute A = Q A where Q is mxm and A is mxn. 817a982d89SJeremy L. Thompson 827a982d89SJeremy L. Thompson @param[in,out] A Matrix to apply Householder Q to, in place 837a982d89SJeremy L. Thompson @param Q Householder Q matrix 847a982d89SJeremy L. Thompson @param tau Householder scaling factors 857a982d89SJeremy L. Thompson @param tmode Transpose mode for application 867a982d89SJeremy L. Thompson @param m Number of rows in A 877a982d89SJeremy L. Thompson @param n Number of columns in A 887a982d89SJeremy L. Thompson @param k Number of elementary reflectors in Q, k<m 897a982d89SJeremy L. Thompson @param row Row stride in A 907a982d89SJeremy L. Thompson @param col Col stride in A 917a982d89SJeremy L. Thompson 927a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 937a982d89SJeremy L. Thompson 947a982d89SJeremy L. Thompson @ref Developer 957a982d89SJeremy L. Thompson **/ 967a982d89SJeremy L. Thompson static int CeedHouseholderApplyQ(CeedScalar *A, const CeedScalar *Q, 977a982d89SJeremy L. Thompson const CeedScalar *tau, CeedTransposeMode tmode, 987a982d89SJeremy L. Thompson CeedInt m, CeedInt n, CeedInt k, 997a982d89SJeremy L. Thompson CeedInt row, CeedInt col) { 1007a982d89SJeremy L. Thompson CeedScalar v[m]; 1017a982d89SJeremy L. Thompson for (CeedInt ii=0; ii<k; ii++) { 1027a982d89SJeremy L. Thompson CeedInt i = tmode == CEED_TRANSPOSE ? ii : k-1-ii; 1037a982d89SJeremy L. Thompson for (CeedInt j=i+1; j<m; j++) 1047a982d89SJeremy L. Thompson v[j] = Q[j*k+i]; 1057a982d89SJeremy L. Thompson // Apply Householder reflector (I - tau v v^T) collograd1d^T 1067a982d89SJeremy L. Thompson CeedHouseholderReflect(&A[i*row], &v[i], tau[i], m-i, n, row, col); 1077a982d89SJeremy L. Thompson } 1087a982d89SJeremy L. Thompson return 0; 1097a982d89SJeremy L. Thompson } 1107a982d89SJeremy L. Thompson 1117a982d89SJeremy L. Thompson /** 1127a982d89SJeremy L. Thompson @brief Compute Givens rotation 1137a982d89SJeremy L. Thompson 1147a982d89SJeremy L. Thompson Computes A = G A (or G^T A in transpose mode) 1157a982d89SJeremy L. Thompson where A is an mxn matrix indexed as A[i*n + j*m] 1167a982d89SJeremy L. Thompson 1177a982d89SJeremy L. Thompson @param[in,out] A Row major matrix to apply Givens rotation to, in place 1187a982d89SJeremy L. Thompson @param c Cosine factor 1197a982d89SJeremy L. Thompson @param s Sine factor 1204c4400c7SValeria Barra @param tmode CEED_NOTRANSPOSE to rotate the basis counter-clockwise, 1214c4400c7SValeria Barra which has the effect of rotating columns of A clockwise; 1224c4400c7SValeria Barra CEED_TRANSPOSE for the opposite rotation 1237a982d89SJeremy L. Thompson @param i First row/column to apply rotation 1247a982d89SJeremy L. Thompson @param k Second row/column to apply rotation 1257a982d89SJeremy L. Thompson @param m Number of rows in A 1267a982d89SJeremy L. Thompson @param n Number of columns in A 1277a982d89SJeremy L. Thompson 1287a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 1297a982d89SJeremy L. Thompson 1307a982d89SJeremy L. Thompson @ref Developer 1317a982d89SJeremy L. Thompson **/ 1327a982d89SJeremy L. Thompson static int CeedGivensRotation(CeedScalar *A, CeedScalar c, CeedScalar s, 1337a982d89SJeremy L. Thompson CeedTransposeMode tmode, CeedInt i, CeedInt k, 1347a982d89SJeremy L. Thompson CeedInt m, CeedInt n) { 1357a982d89SJeremy L. Thompson CeedInt stridej = 1, strideik = m, numits = n; 1367a982d89SJeremy L. Thompson if (tmode == CEED_NOTRANSPOSE) { 1377a982d89SJeremy L. Thompson stridej = n; strideik = 1; numits = m; 1387a982d89SJeremy L. Thompson } 1397a982d89SJeremy L. Thompson 1407a982d89SJeremy L. Thompson // Apply rotation 1417a982d89SJeremy L. Thompson for (CeedInt j=0; j<numits; j++) { 1427a982d89SJeremy L. Thompson CeedScalar tau1 = A[i*strideik+j*stridej], tau2 = A[k*strideik+j*stridej]; 1437a982d89SJeremy L. Thompson A[i*strideik+j*stridej] = c*tau1 - s*tau2; 1447a982d89SJeremy L. Thompson A[k*strideik+j*stridej] = s*tau1 + c*tau2; 1457a982d89SJeremy L. Thompson } 1467a982d89SJeremy L. Thompson 1477a982d89SJeremy L. Thompson return 0; 1487a982d89SJeremy L. Thompson } 1497a982d89SJeremy L. Thompson 1507a982d89SJeremy L. Thompson /** 1517a982d89SJeremy L. Thompson @brief View an array stored in a CeedBasis 1527a982d89SJeremy L. Thompson 153*0a0da059Sjeremylt @param[in] name Name of array 154*0a0da059Sjeremylt @param[in] fpformat Printing format 155*0a0da059Sjeremylt @param[in] m Number of rows in array 156*0a0da059Sjeremylt @param[in] n Number of columns in array 157*0a0da059Sjeremylt @param[in] a Array to be viewed 158*0a0da059Sjeremylt @param[in] stream Stream to view to, e.g., stdout 1597a982d89SJeremy L. Thompson 1607a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 1617a982d89SJeremy L. Thompson 1627a982d89SJeremy L. Thompson @ref Developer 1637a982d89SJeremy L. Thompson **/ 1647a982d89SJeremy L. Thompson static int CeedScalarView(const char *name, const char *fpformat, CeedInt m, 1657a982d89SJeremy L. Thompson CeedInt n, const CeedScalar *a, FILE *stream) { 1667a982d89SJeremy L. Thompson for (int i=0; i<m; i++) { 1677a982d89SJeremy L. Thompson if (m > 1) 1687a982d89SJeremy L. Thompson fprintf(stream, "%12s[%d]:", name, i); 1697a982d89SJeremy L. Thompson else 1707a982d89SJeremy L. Thompson fprintf(stream, "%12s:", name); 1717a982d89SJeremy L. Thompson for (int j=0; j<n; j++) 1727a982d89SJeremy L. Thompson fprintf(stream, fpformat, fabs(a[i*n+j]) > 1E-14 ? a[i*n+j] : 0); 1737a982d89SJeremy L. Thompson fputs("\n", stream); 1747a982d89SJeremy L. Thompson } 1757a982d89SJeremy L. Thompson return 0; 1767a982d89SJeremy L. Thompson } 1777a982d89SJeremy L. Thompson 1787a982d89SJeremy L. Thompson /// @} 1797a982d89SJeremy L. Thompson 1807a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 1817a982d89SJeremy L. Thompson /// Ceed Backend API 1827a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 1837a982d89SJeremy L. Thompson /// @addtogroup CeedBasisBackend 1847a982d89SJeremy L. Thompson /// @{ 1857a982d89SJeremy L. Thompson 1867a982d89SJeremy L. Thompson /** 1877a982d89SJeremy L. Thompson @brief Return collocated grad matrix 1887a982d89SJeremy L. Thompson 1897a982d89SJeremy L. Thompson @param basis CeedBasis 1907a982d89SJeremy L. Thompson @param[out] collograd1d Row-major (Q1d * Q1d) matrix expressing derivatives of 1917a982d89SJeremy L. Thompson basis functions at quadrature points 1927a982d89SJeremy L. Thompson 1937a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 1947a982d89SJeremy L. Thompson 1957a982d89SJeremy L. Thompson @ref Backend 1967a982d89SJeremy L. Thompson **/ 1977a982d89SJeremy L. Thompson int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *collograd1d) { 1987a982d89SJeremy L. Thompson int i, j, k; 1997a982d89SJeremy L. Thompson Ceed ceed; 2007a982d89SJeremy L. Thompson CeedInt ierr, P1d=(basis)->P1d, Q1d=(basis)->Q1d; 2017a982d89SJeremy L. Thompson CeedScalar *interp1d, *grad1d, tau[Q1d]; 2027a982d89SJeremy L. Thompson 2037a982d89SJeremy L. Thompson ierr = CeedMalloc(Q1d*P1d, &interp1d); CeedChk(ierr); 2047a982d89SJeremy L. Thompson ierr = CeedMalloc(Q1d*P1d, &grad1d); CeedChk(ierr); 2057a982d89SJeremy L. Thompson memcpy(interp1d, (basis)->interp1d, Q1d*P1d*sizeof(basis)->interp1d[0]); 2067a982d89SJeremy L. Thompson memcpy(grad1d, (basis)->grad1d, Q1d*P1d*sizeof(basis)->interp1d[0]); 2077a982d89SJeremy L. Thompson 2087a982d89SJeremy L. Thompson // QR Factorization, interp1d = Q R 2097a982d89SJeremy L. Thompson ierr = CeedBasisGetCeed(basis, &ceed); CeedChk(ierr); 2107a982d89SJeremy L. Thompson ierr = CeedQRFactorization(ceed, interp1d, tau, Q1d, P1d); CeedChk(ierr); 2117a982d89SJeremy L. Thompson 2127a982d89SJeremy L. Thompson // Apply Rinv, collograd1d = grad1d Rinv 2137a982d89SJeremy L. Thompson for (i=0; i<Q1d; i++) { // Row i 2147a982d89SJeremy L. Thompson collograd1d[Q1d*i] = grad1d[P1d*i]/interp1d[0]; 2157a982d89SJeremy L. Thompson for (j=1; j<P1d; j++) { // Column j 2167a982d89SJeremy L. Thompson collograd1d[j+Q1d*i] = grad1d[j+P1d*i]; 2177a982d89SJeremy L. Thompson for (k=0; k<j; k++) 2187a982d89SJeremy L. Thompson collograd1d[j+Q1d*i] -= interp1d[j+P1d*k]*collograd1d[k+Q1d*i]; 2197a982d89SJeremy L. Thompson collograd1d[j+Q1d*i] /= interp1d[j+P1d*j]; 2207a982d89SJeremy L. Thompson } 2217a982d89SJeremy L. Thompson for (j=P1d; j<Q1d; j++) 2227a982d89SJeremy L. Thompson collograd1d[j+Q1d*i] = 0; 2237a982d89SJeremy L. Thompson } 2247a982d89SJeremy L. Thompson 2257a982d89SJeremy L. Thompson // Apply Qtranspose, collograd = collograd Qtranspose 2267a982d89SJeremy L. Thompson CeedHouseholderApplyQ(collograd1d, interp1d, tau, CEED_NOTRANSPOSE, 2277a982d89SJeremy L. Thompson Q1d, Q1d, P1d, 1, Q1d); 2287a982d89SJeremy L. Thompson 2297a982d89SJeremy L. Thompson ierr = CeedFree(&interp1d); CeedChk(ierr); 2307a982d89SJeremy L. Thompson ierr = CeedFree(&grad1d); CeedChk(ierr); 2317a982d89SJeremy L. Thompson 2327a982d89SJeremy L. Thompson return 0; 2337a982d89SJeremy L. Thompson } 2347a982d89SJeremy L. Thompson 2357a982d89SJeremy L. Thompson /** 2367a982d89SJeremy L. Thompson @brief Get Ceed associated with a CeedBasis 2377a982d89SJeremy L. Thompson 2387a982d89SJeremy L. Thompson @param basis CeedBasis 2397a982d89SJeremy L. Thompson @param[out] ceed Variable to store Ceed 2407a982d89SJeremy L. Thompson 2417a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 2427a982d89SJeremy L. Thompson 2437a982d89SJeremy L. Thompson @ref Backend 2447a982d89SJeremy L. Thompson **/ 2457a982d89SJeremy L. Thompson int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) { 2467a982d89SJeremy L. Thompson *ceed = basis->ceed; 2477a982d89SJeremy L. Thompson return 0; 2487a982d89SJeremy L. Thompson } 2497a982d89SJeremy L. Thompson 2507a982d89SJeremy L. Thompson /** 2517a982d89SJeremy L. Thompson @brief Get tensor status for given CeedBasis 2527a982d89SJeremy L. Thompson 2537a982d89SJeremy L. Thompson @param basis CeedBasis 2547a982d89SJeremy L. Thompson @param[out] tensor Variable to store tensor status 2557a982d89SJeremy L. Thompson 2567a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 2577a982d89SJeremy L. Thompson 2587a982d89SJeremy L. Thompson @ref Backend 2597a982d89SJeremy L. Thompson **/ 2607a982d89SJeremy L. Thompson int CeedBasisGetTensorStatus(CeedBasis basis, bool *tensor) { 2617a982d89SJeremy L. Thompson *tensor = basis->tensorbasis; 2627a982d89SJeremy L. Thompson return 0; 2637a982d89SJeremy L. Thompson } 2647a982d89SJeremy L. Thompson 2657a982d89SJeremy L. Thompson /** 2667a982d89SJeremy L. Thompson @brief Get dimension for given CeedBasis 2677a982d89SJeremy L. Thompson 2687a982d89SJeremy L. Thompson @param basis CeedBasis 2697a982d89SJeremy L. Thompson @param[out] dim Variable to store dimension of basis 2707a982d89SJeremy L. Thompson 2717a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 2727a982d89SJeremy L. Thompson 2737a982d89SJeremy L. Thompson @ref Backend 2747a982d89SJeremy L. Thompson **/ 2757a982d89SJeremy L. Thompson int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) { 2767a982d89SJeremy L. Thompson *dim = basis->dim; 2777a982d89SJeremy L. Thompson return 0; 2787a982d89SJeremy L. Thompson } 2797a982d89SJeremy L. Thompson 2807a982d89SJeremy L. Thompson /** 2817a982d89SJeremy L. Thompson @brief Get number of components for given CeedBasis 2827a982d89SJeremy L. Thompson 2837a982d89SJeremy L. Thompson @param basis CeedBasis 2847a982d89SJeremy L. Thompson @param[out] numcomp Variable to store number of components of basis 2857a982d89SJeremy L. Thompson 2867a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 2877a982d89SJeremy L. Thompson 2887a982d89SJeremy L. Thompson @ref Backend 2897a982d89SJeremy L. Thompson **/ 2907a982d89SJeremy L. Thompson int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *numcomp) { 2917a982d89SJeremy L. Thompson *numcomp = basis->ncomp; 2927a982d89SJeremy L. Thompson return 0; 2937a982d89SJeremy L. Thompson } 2947a982d89SJeremy L. Thompson 2957a982d89SJeremy L. Thompson /** 2967a982d89SJeremy L. Thompson @brief Get total number of nodes (in 1 dimension) of a CeedBasis 2977a982d89SJeremy L. Thompson 2987a982d89SJeremy L. Thompson @param basis CeedBasis 2997a982d89SJeremy L. Thompson @param[out] P1d Variable to store number of nodes 3007a982d89SJeremy L. Thompson 3017a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3027a982d89SJeremy L. Thompson 3037a982d89SJeremy L. Thompson @ref Backend 3047a982d89SJeremy L. Thompson **/ 3057a982d89SJeremy L. Thompson int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P1d) { 3067a982d89SJeremy L. Thompson if (!basis->tensorbasis) 3077a982d89SJeremy L. Thompson // LCOV_EXCL_START 3087a982d89SJeremy L. Thompson return CeedError(basis->ceed, 1, "Cannot supply P1d for non-tensor basis"); 3097a982d89SJeremy L. Thompson // LCOV_EXCL_STOP 3107a982d89SJeremy L. Thompson 3117a982d89SJeremy L. Thompson *P1d = basis->P1d; 3127a982d89SJeremy L. Thompson return 0; 3137a982d89SJeremy L. Thompson } 3147a982d89SJeremy L. Thompson 3157a982d89SJeremy L. Thompson /** 3167a982d89SJeremy L. Thompson @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis 3177a982d89SJeremy L. Thompson 3187a982d89SJeremy L. Thompson @param basis CeedBasis 3197a982d89SJeremy L. Thompson @param[out] Q1d Variable to store number of quadrature points 3207a982d89SJeremy L. Thompson 3217a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3227a982d89SJeremy L. Thompson 3237a982d89SJeremy L. Thompson @ref Backend 3247a982d89SJeremy L. Thompson **/ 3257a982d89SJeremy L. Thompson int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q1d) { 3267a982d89SJeremy L. Thompson if (!basis->tensorbasis) 3277a982d89SJeremy L. Thompson // LCOV_EXCL_START 3287a982d89SJeremy L. Thompson return CeedError(basis->ceed, 1, "Cannot supply Q1d for non-tensor basis"); 3297a982d89SJeremy L. Thompson // LCOV_EXCL_STOP 3307a982d89SJeremy L. Thompson 3317a982d89SJeremy L. Thompson *Q1d = basis->Q1d; 3327a982d89SJeremy L. Thompson return 0; 3337a982d89SJeremy L. Thompson } 3347a982d89SJeremy L. Thompson 3357a982d89SJeremy L. Thompson /** 3367a982d89SJeremy L. Thompson @brief Get reference coordinates of quadrature points (in dim dimensions) 3377a982d89SJeremy L. Thompson of a CeedBasis 3387a982d89SJeremy L. Thompson 3397a982d89SJeremy L. Thompson @param basis CeedBasis 3407a982d89SJeremy L. Thompson @param[out] qref Variable to store reference coordinates of quadrature points 3417a982d89SJeremy L. Thompson 3427a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3437a982d89SJeremy L. Thompson 3447a982d89SJeremy L. Thompson @ref Backend 3457a982d89SJeremy L. Thompson **/ 3467a982d89SJeremy L. Thompson int CeedBasisGetQRef(CeedBasis basis, CeedScalar **qref) { 3477a982d89SJeremy L. Thompson *qref = basis->qref1d; 3487a982d89SJeremy L. Thompson return 0; 3497a982d89SJeremy L. Thompson } 3507a982d89SJeremy L. Thompson 3517a982d89SJeremy L. Thompson /** 3527a982d89SJeremy L. Thompson @brief Get quadrature weights of quadrature points (in dim dimensions) 3537a982d89SJeremy L. Thompson of a CeedBasis 3547a982d89SJeremy L. Thompson 3557a982d89SJeremy L. Thompson @param basis CeedBasis 3567a982d89SJeremy L. Thompson @param[out] qweight Variable to store quadrature weights 3577a982d89SJeremy L. Thompson 3587a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3597a982d89SJeremy L. Thompson 3607a982d89SJeremy L. Thompson @ref Backend 3617a982d89SJeremy L. Thompson **/ 3627a982d89SJeremy L. Thompson int CeedBasisGetQWeights(CeedBasis basis, CeedScalar **qweight) { 3637a982d89SJeremy L. Thompson *qweight = basis->qweight1d; 3647a982d89SJeremy L. Thompson return 0; 3657a982d89SJeremy L. Thompson } 3667a982d89SJeremy L. Thompson 3677a982d89SJeremy L. Thompson /** 3687a982d89SJeremy L. Thompson @brief Get interpolation matrix of a CeedBasis 3697a982d89SJeremy L. Thompson 3707a982d89SJeremy L. Thompson @param basis CeedBasis 3717a982d89SJeremy L. Thompson @param[out] interp Variable to store interpolation matrix 3727a982d89SJeremy L. Thompson 3737a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 3747a982d89SJeremy L. Thompson 3757a982d89SJeremy L. Thompson @ref Backend 3767a982d89SJeremy L. Thompson **/ 3777a982d89SJeremy L. Thompson int CeedBasisGetInterp(CeedBasis basis, CeedScalar **interp) { 3787a982d89SJeremy L. Thompson if (!basis->interp && basis->tensorbasis) { 3797a982d89SJeremy L. Thompson // Allocate 3807a982d89SJeremy L. Thompson int ierr; 3817a982d89SJeremy L. Thompson ierr = CeedMalloc(basis->Q*basis->P, &basis->interp); CeedChk(ierr); 3827a982d89SJeremy L. Thompson 3837a982d89SJeremy L. Thompson // Initialize 3847a982d89SJeremy L. Thompson for (CeedInt i=0; i<basis->Q*basis->P; i++) 3857a982d89SJeremy L. Thompson basis->interp[i] = 1.0; 3867a982d89SJeremy L. Thompson 3877a982d89SJeremy L. Thompson // Calculate 3887a982d89SJeremy L. Thompson for (CeedInt d=0; d<basis->dim; d++) 3897a982d89SJeremy L. Thompson for (CeedInt qpt=0; qpt<basis->Q; qpt++) 3907a982d89SJeremy L. Thompson for (CeedInt node=0; node<basis->P; node++) { 3917a982d89SJeremy L. Thompson CeedInt p = (node / CeedIntPow(basis->P1d, d)) % basis->P1d; 3927a982d89SJeremy L. Thompson CeedInt q = (qpt / CeedIntPow(basis->Q1d, d)) % basis->Q1d; 3937a982d89SJeremy L. Thompson basis->interp[qpt*(basis->P)+node] *= basis->interp1d[q*basis->P1d+p]; 3947a982d89SJeremy L. Thompson } 3957a982d89SJeremy L. Thompson } 3967a982d89SJeremy L. Thompson 3977a982d89SJeremy L. Thompson *interp = basis->interp; 3987a982d89SJeremy L. Thompson 3997a982d89SJeremy L. Thompson return 0; 4007a982d89SJeremy L. Thompson } 4017a982d89SJeremy L. Thompson 4027a982d89SJeremy L. Thompson /** 4037a982d89SJeremy L. Thompson @brief Get 1D interpolation matrix of a tensor product CeedBasis 4047a982d89SJeremy L. Thompson 4057a982d89SJeremy L. Thompson @param basis CeedBasis 4067a982d89SJeremy L. Thompson @param[out] interp1d Variable to store interpolation matrix 4077a982d89SJeremy L. Thompson 4087a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 4097a982d89SJeremy L. Thompson 4107a982d89SJeremy L. Thompson @ref Backend 4117a982d89SJeremy L. Thompson **/ 4127a982d89SJeremy L. Thompson int CeedBasisGetInterp1D(CeedBasis basis, CeedScalar **interp1d) { 4137a982d89SJeremy L. Thompson if (!basis->tensorbasis) 4147a982d89SJeremy L. Thompson // LCOV_EXCL_START 4157a982d89SJeremy L. Thompson return CeedError(basis->ceed, 1, "CeedBasis is not a tensor product basis."); 4167a982d89SJeremy L. Thompson // LCOV_EXCL_STOP 4177a982d89SJeremy L. Thompson 4187a982d89SJeremy L. Thompson *interp1d = basis->interp1d; 4197a982d89SJeremy L. Thompson 4207a982d89SJeremy L. Thompson return 0; 4217a982d89SJeremy L. Thompson } 4227a982d89SJeremy L. Thompson 4237a982d89SJeremy L. Thompson /** 4247a982d89SJeremy L. Thompson @brief Get gradient matrix of a CeedBasis 4257a982d89SJeremy L. Thompson 4267a982d89SJeremy L. Thompson @param basis CeedBasis 4277a982d89SJeremy L. Thompson @param[out] grad Variable to store gradient matrix 4287a982d89SJeremy L. Thompson 4297a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 4307a982d89SJeremy L. Thompson 4317a982d89SJeremy L. Thompson @ref Backend 4327a982d89SJeremy L. Thompson **/ 4337a982d89SJeremy L. Thompson int CeedBasisGetGrad(CeedBasis basis, CeedScalar **grad) { 4347a982d89SJeremy L. Thompson if (!basis->grad && basis->tensorbasis) { 4357a982d89SJeremy L. Thompson // Allocate 4367a982d89SJeremy L. Thompson int ierr; 4377a982d89SJeremy L. Thompson ierr = CeedMalloc(basis->dim*basis->Q*basis->P, &basis->grad); 4387a982d89SJeremy L. Thompson CeedChk(ierr); 4397a982d89SJeremy L. Thompson 4407a982d89SJeremy L. Thompson // Initialize 4417a982d89SJeremy L. Thompson for (CeedInt i=0; i<basis->dim*basis->Q*basis->P; i++) 4427a982d89SJeremy L. Thompson basis->grad[i] = 1.0; 4437a982d89SJeremy L. Thompson 4447a982d89SJeremy L. Thompson // Calculate 4457a982d89SJeremy L. Thompson for (CeedInt d=0; d<basis->dim; d++) 4467a982d89SJeremy L. Thompson for (CeedInt i=0; i<basis->dim; i++) 4477a982d89SJeremy L. Thompson for (CeedInt qpt=0; qpt<basis->Q; qpt++) 4487a982d89SJeremy L. Thompson for (CeedInt node=0; node<basis->P; node++) { 4497a982d89SJeremy L. Thompson CeedInt p = (node / CeedIntPow(basis->P1d, d)) % basis->P1d; 4507a982d89SJeremy L. Thompson CeedInt q = (qpt / CeedIntPow(basis->Q1d, d)) % basis->Q1d; 4517a982d89SJeremy L. Thompson if (i == d) 4527a982d89SJeremy L. Thompson basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *= 4537a982d89SJeremy L. Thompson basis->grad1d[q*basis->P1d+p]; 4547a982d89SJeremy L. Thompson else 4557a982d89SJeremy L. Thompson basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *= 4567a982d89SJeremy L. Thompson basis->interp1d[q*basis->P1d+p]; 4577a982d89SJeremy L. Thompson } 4587a982d89SJeremy L. Thompson } 4597a982d89SJeremy L. Thompson 4607a982d89SJeremy L. Thompson *grad = basis->grad; 4617a982d89SJeremy L. Thompson 4627a982d89SJeremy L. Thompson return 0; 4637a982d89SJeremy L. Thompson } 4647a982d89SJeremy L. Thompson 4657a982d89SJeremy L. Thompson /** 4667a982d89SJeremy L. Thompson @brief Get 1D gradient matrix of a tensor product CeedBasis 4677a982d89SJeremy L. Thompson 4687a982d89SJeremy L. Thompson @param basis CeedBasis 4697a982d89SJeremy L. Thompson @param[out] grad1d Variable to store gradient matrix 4707a982d89SJeremy L. Thompson 4717a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 4727a982d89SJeremy L. Thompson 4737a982d89SJeremy L. Thompson @ref Backend 4747a982d89SJeremy L. Thompson **/ 4757a982d89SJeremy L. Thompson int CeedBasisGetGrad1D(CeedBasis basis, CeedScalar **grad1d) { 4767a982d89SJeremy L. Thompson if (!basis->tensorbasis) 4777a982d89SJeremy L. Thompson // LCOV_EXCL_START 4787a982d89SJeremy L. Thompson return CeedError(basis->ceed, 1, "CeedBasis is not a tensor product basis."); 4797a982d89SJeremy L. Thompson // LCOV_EXCL_STOP 4807a982d89SJeremy L. Thompson 4817a982d89SJeremy L. Thompson *grad1d = basis->grad1d; 4827a982d89SJeremy L. Thompson 4837a982d89SJeremy L. Thompson return 0; 4847a982d89SJeremy L. Thompson } 4857a982d89SJeremy L. Thompson 4867a982d89SJeremy L. Thompson /** 4877a982d89SJeremy L. Thompson @brief Get backend data of a CeedBasis 4887a982d89SJeremy L. Thompson 4897a982d89SJeremy L. Thompson @param basis CeedBasis 4907a982d89SJeremy L. Thompson @param[out] data Variable to store data 4917a982d89SJeremy L. Thompson 4927a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 4937a982d89SJeremy L. Thompson 4947a982d89SJeremy L. Thompson @ref Backend 4957a982d89SJeremy L. Thompson **/ 4967a982d89SJeremy L. Thompson int CeedBasisGetData(CeedBasis basis, void **data) { 4977a982d89SJeremy L. Thompson *data = basis->data; 4987a982d89SJeremy L. Thompson return 0; 4997a982d89SJeremy L. Thompson } 5007a982d89SJeremy L. Thompson 5017a982d89SJeremy L. Thompson /** 5027a982d89SJeremy L. Thompson @brief Set backend data of a CeedBasis 5037a982d89SJeremy L. Thompson 5047a982d89SJeremy L. Thompson @param[out] basis CeedBasis 5057a982d89SJeremy L. Thompson @param data Data to set 5067a982d89SJeremy L. Thompson 5077a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 5087a982d89SJeremy L. Thompson 5097a982d89SJeremy L. Thompson @ref Backend 5107a982d89SJeremy L. Thompson **/ 5117a982d89SJeremy L. Thompson int CeedBasisSetData(CeedBasis basis, void **data) { 5127a982d89SJeremy L. Thompson basis->data = *data; 5137a982d89SJeremy L. Thompson return 0; 5147a982d89SJeremy L. Thompson } 5157a982d89SJeremy L. Thompson 5167a982d89SJeremy L. Thompson /** 5177a982d89SJeremy L. Thompson @brief Get dimension for given CeedElemTopology 5187a982d89SJeremy L. Thompson 5197a982d89SJeremy L. Thompson @param topo CeedElemTopology 5207a982d89SJeremy L. Thompson @param[out] dim Variable to store dimension of topology 5217a982d89SJeremy L. Thompson 5227a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 5237a982d89SJeremy L. Thompson 5247a982d89SJeremy L. Thompson @ref Backend 5257a982d89SJeremy L. Thompson **/ 5267a982d89SJeremy L. Thompson int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) { 5277a982d89SJeremy L. Thompson *dim = (CeedInt) topo >> 16; 5287a982d89SJeremy L. Thompson return 0; 5297a982d89SJeremy L. Thompson } 5307a982d89SJeremy L. Thompson 5317a982d89SJeremy L. Thompson /** 5327a982d89SJeremy L. Thompson @brief Get CeedTensorContract of a CeedBasis 5337a982d89SJeremy L. Thompson 5347a982d89SJeremy L. Thompson @param basis CeedBasis 5357a982d89SJeremy L. Thompson @param[out] contract Variable to store CeedTensorContract 5367a982d89SJeremy L. Thompson 5377a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 5387a982d89SJeremy L. Thompson 5397a982d89SJeremy L. Thompson @ref Backend 5407a982d89SJeremy L. Thompson **/ 5417a982d89SJeremy L. Thompson int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) { 5427a982d89SJeremy L. Thompson *contract = basis->contract; 5437a982d89SJeremy L. Thompson return 0; 5447a982d89SJeremy L. Thompson } 5457a982d89SJeremy L. Thompson 5467a982d89SJeremy L. Thompson /** 5477a982d89SJeremy L. Thompson @brief Set CeedTensorContract of a CeedBasis 5487a982d89SJeremy L. Thompson 5497a982d89SJeremy L. Thompson @param[out] basis CeedBasis 5507a982d89SJeremy L. Thompson @param contract CeedTensorContract to set 5517a982d89SJeremy L. Thompson 5527a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 5537a982d89SJeremy L. Thompson 5547a982d89SJeremy L. Thompson @ref Backend 5557a982d89SJeremy L. Thompson **/ 5567a982d89SJeremy L. Thompson int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract *contract) { 5577a982d89SJeremy L. Thompson basis->contract = *contract; 5587a982d89SJeremy L. Thompson return 0; 5597a982d89SJeremy L. Thompson } 5607a982d89SJeremy L. Thompson 5617a982d89SJeremy L. Thompson /** 5627a982d89SJeremy L. Thompson @brief Return a reference implementation of matrix multiplication C = A B. 5637a982d89SJeremy L. Thompson Note, this is a reference implementation for CPU CeedScalar pointers 5647a982d89SJeremy L. Thompson that is not intended for high performance. 5657a982d89SJeremy L. Thompson 5667a982d89SJeremy L. Thompson @param ceed A Ceed context for error handling 5677a982d89SJeremy L. Thompson @param[in] matA Row-major matrix A 5687a982d89SJeremy L. Thompson @param[in] matB Row-major matrix B 5697a982d89SJeremy L. Thompson @param[out] matC Row-major output matrix C 5707a982d89SJeremy L. Thompson @param m Number of rows of C 5717a982d89SJeremy L. Thompson @param n Number of columns of C 5727a982d89SJeremy L. Thompson @param kk Number of columns of A/rows of B 5737a982d89SJeremy L. Thompson 5747a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 5757a982d89SJeremy L. Thompson 5767a982d89SJeremy L. Thompson @ref Utility 5777a982d89SJeremy L. Thompson **/ 5787a982d89SJeremy L. Thompson int CeedMatrixMultiply(Ceed ceed, const CeedScalar *matA, 5797a982d89SJeremy L. Thompson const CeedScalar *matB, CeedScalar *matC, CeedInt m, 5807a982d89SJeremy L. Thompson CeedInt n, CeedInt kk) { 5817a982d89SJeremy L. Thompson for (CeedInt i=0; i<m; i++) 5827a982d89SJeremy L. Thompson for (CeedInt j=0; j<n; j++) { 5837a982d89SJeremy L. Thompson CeedScalar sum = 0; 5847a982d89SJeremy L. Thompson for (CeedInt k=0; k<kk; k++) 5857a982d89SJeremy L. Thompson sum += matA[k+i*kk]*matB[j+k*n]; 5867a982d89SJeremy L. Thompson matC[j+i*n] = sum; 5877a982d89SJeremy L. Thompson } 5887a982d89SJeremy L. Thompson return 0; 5897a982d89SJeremy L. Thompson } 5907a982d89SJeremy L. Thompson 5917a982d89SJeremy L. Thompson /// @} 5927a982d89SJeremy L. Thompson 5937a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 5947a982d89SJeremy L. Thompson /// CeedBasis Public API 5957a982d89SJeremy L. Thompson /// ---------------------------------------------------------------------------- 5967a982d89SJeremy L. Thompson /// @addtogroup CeedBasisUser 597d7b241e6Sjeremylt /// @{ 598d7b241e6Sjeremylt 599b11c1e72Sjeremylt /** 60095bb1877Svaleriabarra @brief Create a tensor-product basis for H^1 discretizations 601b11c1e72Sjeremylt 602b11c1e72Sjeremylt @param ceed A Ceed object where the CeedBasis will be created 603b11c1e72Sjeremylt @param dim Topological dimension 604b11c1e72Sjeremylt @param ncomp Number of field components (1 for scalar fields) 605b11c1e72Sjeremylt @param P1d Number of nodes in one dimension 606b11c1e72Sjeremylt @param Q1d Number of quadrature points in one dimension 60795bb1877Svaleriabarra @param interp1d Row-major (Q1d * P1d) matrix expressing the values of nodal 608b11c1e72Sjeremylt basis functions at quadrature points 60995bb1877Svaleriabarra @param grad1d Row-major (Q1d * P1d) matrix expressing derivatives of nodal 610b11c1e72Sjeremylt basis functions at quadrature points 611b11c1e72Sjeremylt @param qref1d Array of length Q1d holding the locations of quadrature points 612b11c1e72Sjeremylt on the 1D reference element [-1, 1] 613b11c1e72Sjeremylt @param qweight1d Array of length Q1d holding the quadrature weights on the 614b11c1e72Sjeremylt reference element 615b11c1e72Sjeremylt @param[out] basis Address of the variable where the newly created 616b11c1e72Sjeremylt CeedBasis will be stored. 617b11c1e72Sjeremylt 618b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 619dfdf5a53Sjeremylt 6207a982d89SJeremy L. Thompson @ref User 621b11c1e72Sjeremylt **/ 622d7b241e6Sjeremylt int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt ncomp, CeedInt P1d, 623d7b241e6Sjeremylt CeedInt Q1d, const CeedScalar *interp1d, 624d7b241e6Sjeremylt const CeedScalar *grad1d, const CeedScalar *qref1d, 625d7b241e6Sjeremylt const CeedScalar *qweight1d, CeedBasis *basis) { 626d7b241e6Sjeremylt int ierr; 627d7b241e6Sjeremylt 6284d537eeaSYohann if (dim<1) 629c042f62fSJeremy L Thompson // LCOV_EXCL_START 6304d537eeaSYohann return CeedError(ceed, 1, "Basis dimension must be a positive value"); 631c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 6324d537eeaSYohann 6335fe0d4faSjeremylt if (!ceed->BasisCreateTensorH1) { 6345fe0d4faSjeremylt Ceed delegate; 635aefd8378Sjeremylt ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr); 6365fe0d4faSjeremylt 6375fe0d4faSjeremylt if (!delegate) 638c042f62fSJeremy L Thompson // LCOV_EXCL_START 639d7b241e6Sjeremylt return CeedError(ceed, 1, "Backend does not support BasisCreateTensorH1"); 640c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 6415fe0d4faSjeremylt 6425fe0d4faSjeremylt ierr = CeedBasisCreateTensorH1(delegate, dim, ncomp, P1d, 6435fe0d4faSjeremylt Q1d, interp1d, grad1d, qref1d, 6445fe0d4faSjeremylt qweight1d, basis); CeedChk(ierr); 6455fe0d4faSjeremylt return 0; 6465fe0d4faSjeremylt } 647d7b241e6Sjeremylt ierr = CeedCalloc(1,basis); CeedChk(ierr); 648d7b241e6Sjeremylt (*basis)->ceed = ceed; 649d7b241e6Sjeremylt ceed->refcount++; 650d7b241e6Sjeremylt (*basis)->refcount = 1; 651a8de75f0Sjeremylt (*basis)->tensorbasis = 1; 652d7b241e6Sjeremylt (*basis)->dim = dim; 653d7b241e6Sjeremylt (*basis)->ncomp = ncomp; 654d7b241e6Sjeremylt (*basis)->P1d = P1d; 655d7b241e6Sjeremylt (*basis)->Q1d = Q1d; 656a8de75f0Sjeremylt (*basis)->P = CeedIntPow(P1d, dim); 657a8de75f0Sjeremylt (*basis)->Q = CeedIntPow(Q1d, dim); 658d7b241e6Sjeremylt ierr = CeedMalloc(Q1d,&(*basis)->qref1d); CeedChk(ierr); 659d7b241e6Sjeremylt ierr = CeedMalloc(Q1d,&(*basis)->qweight1d); CeedChk(ierr); 660d7b241e6Sjeremylt memcpy((*basis)->qref1d, qref1d, Q1d*sizeof(qref1d[0])); 661d7b241e6Sjeremylt memcpy((*basis)->qweight1d, qweight1d, Q1d*sizeof(qweight1d[0])); 662d7b241e6Sjeremylt ierr = CeedMalloc(Q1d*P1d,&(*basis)->interp1d); CeedChk(ierr); 663d7b241e6Sjeremylt ierr = CeedMalloc(Q1d*P1d,&(*basis)->grad1d); CeedChk(ierr); 664d7b241e6Sjeremylt memcpy((*basis)->interp1d, interp1d, Q1d*P1d*sizeof(interp1d[0])); 66509486605Sjeremylt memcpy((*basis)->grad1d, grad1d, Q1d*P1d*sizeof(grad1d[0])); 666667bc5fcSjeremylt ierr = ceed->BasisCreateTensorH1(dim, P1d, Q1d, interp1d, grad1d, qref1d, 667d7b241e6Sjeremylt qweight1d, *basis); CeedChk(ierr); 668d7b241e6Sjeremylt return 0; 669d7b241e6Sjeremylt } 670d7b241e6Sjeremylt 671b11c1e72Sjeremylt /** 67295bb1877Svaleriabarra @brief Create a tensor-product Lagrange basis 673b11c1e72Sjeremylt 674b11c1e72Sjeremylt @param ceed A Ceed object where the CeedBasis will be created 675b11c1e72Sjeremylt @param dim Topological dimension of element 67695bb1877Svaleriabarra @param ncomp Number of field components (1 for scalar fields) 677b11c1e72Sjeremylt @param P Number of Gauss-Lobatto nodes in one dimension. The 678b11c1e72Sjeremylt polynomial degree of the resulting Q_k element is k=P-1. 679b11c1e72Sjeremylt @param Q Number of quadrature points in one dimension. 680b11c1e72Sjeremylt @param qmode Distribution of the Q quadrature points (affects order of 681b11c1e72Sjeremylt accuracy for the quadrature) 682b11c1e72Sjeremylt @param[out] basis Address of the variable where the newly created 683b11c1e72Sjeremylt CeedBasis will be stored. 684b11c1e72Sjeremylt 685b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 686dfdf5a53Sjeremylt 6877a982d89SJeremy L. Thompson @ref User 688b11c1e72Sjeremylt **/ 689d7b241e6Sjeremylt int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt ncomp, 690692c2638Sjeremylt CeedInt P, CeedInt Q, CeedQuadMode qmode, 691692c2638Sjeremylt CeedBasis *basis) { 692d7b241e6Sjeremylt // Allocate 693d7b241e6Sjeremylt int ierr, i, j, k; 694d7b241e6Sjeremylt CeedScalar c1, c2, c3, c4, dx, *nodes, *interp1d, *grad1d, *qref1d, *qweight1d; 6954d537eeaSYohann 6964d537eeaSYohann if (dim<1) 697c042f62fSJeremy L Thompson // LCOV_EXCL_START 6984d537eeaSYohann return CeedError(ceed, 1, "Basis dimension must be a positive value"); 699c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 7004d537eeaSYohann 701d7b241e6Sjeremylt ierr = CeedCalloc(P*Q, &interp1d); CeedChk(ierr); 702d7b241e6Sjeremylt ierr = CeedCalloc(P*Q, &grad1d); CeedChk(ierr); 703d7b241e6Sjeremylt ierr = CeedCalloc(P, &nodes); CeedChk(ierr); 704d7b241e6Sjeremylt ierr = CeedCalloc(Q, &qref1d); CeedChk(ierr); 705d7b241e6Sjeremylt ierr = CeedCalloc(Q, &qweight1d); CeedChk(ierr); 706d7b241e6Sjeremylt // Get Nodes and Weights 707d7b241e6Sjeremylt ierr = CeedLobattoQuadrature(P, nodes, NULL); CeedChk(ierr); 708d7b241e6Sjeremylt switch (qmode) { 709d7b241e6Sjeremylt case CEED_GAUSS: 710d7b241e6Sjeremylt ierr = CeedGaussQuadrature(Q, qref1d, qweight1d); CeedChk(ierr); 711d7b241e6Sjeremylt break; 712d7b241e6Sjeremylt case CEED_GAUSS_LOBATTO: 713d7b241e6Sjeremylt ierr = CeedLobattoQuadrature(Q, qref1d, qweight1d); CeedChk(ierr); 714d7b241e6Sjeremylt break; 715d7b241e6Sjeremylt } 716d7b241e6Sjeremylt // Build B, D matrix 717d7b241e6Sjeremylt // Fornberg, 1998 718d7b241e6Sjeremylt for (i = 0; i < Q; i++) { 719d7b241e6Sjeremylt c1 = 1.0; 720d7b241e6Sjeremylt c3 = nodes[0] - qref1d[i]; 721d7b241e6Sjeremylt interp1d[i*P+0] = 1.0; 722d7b241e6Sjeremylt for (j = 1; j < P; j++) { 723d7b241e6Sjeremylt c2 = 1.0; 724d7b241e6Sjeremylt c4 = c3; 725d7b241e6Sjeremylt c3 = nodes[j] - qref1d[i]; 726d7b241e6Sjeremylt for (k = 0; k < j; k++) { 727d7b241e6Sjeremylt dx = nodes[j] - nodes[k]; 728d7b241e6Sjeremylt c2 *= dx; 729d7b241e6Sjeremylt if (k == j - 1) { 730d7b241e6Sjeremylt grad1d[i*P + j] = c1*(interp1d[i*P + k] - c4*grad1d[i*P + k]) / c2; 731d7b241e6Sjeremylt interp1d[i*P + j] = - c1*c4*interp1d[i*P + k] / c2; 732d7b241e6Sjeremylt } 733d7b241e6Sjeremylt grad1d[i*P + k] = (c3*grad1d[i*P + k] - interp1d[i*P + k]) / dx; 734d7b241e6Sjeremylt interp1d[i*P + k] = c3*interp1d[i*P + k] / dx; 735d7b241e6Sjeremylt } 736d7b241e6Sjeremylt c1 = c2; 737d7b241e6Sjeremylt } 738d7b241e6Sjeremylt } 739d7b241e6Sjeremylt // // Pass to CeedBasisCreateTensorH1 740d7b241e6Sjeremylt ierr = CeedBasisCreateTensorH1(ceed, dim, ncomp, P, Q, interp1d, grad1d, qref1d, 741d7b241e6Sjeremylt qweight1d, basis); CeedChk(ierr); 742d7b241e6Sjeremylt ierr = CeedFree(&interp1d); CeedChk(ierr); 743d7b241e6Sjeremylt ierr = CeedFree(&grad1d); CeedChk(ierr); 744d7b241e6Sjeremylt ierr = CeedFree(&nodes); CeedChk(ierr); 745d7b241e6Sjeremylt ierr = CeedFree(&qref1d); CeedChk(ierr); 746d7b241e6Sjeremylt ierr = CeedFree(&qweight1d); CeedChk(ierr); 747d7b241e6Sjeremylt return 0; 748d7b241e6Sjeremylt } 749d7b241e6Sjeremylt 750b11c1e72Sjeremylt /** 75195bb1877Svaleriabarra @brief Create a non tensor-product basis for H^1 discretizations 752a8de75f0Sjeremylt 753a8de75f0Sjeremylt @param ceed A Ceed object where the CeedBasis will be created 754a8de75f0Sjeremylt @param topo Topology of element, e.g. hypercube, simplex, ect 755a8de75f0Sjeremylt @param ncomp Number of field components (1 for scalar fields) 7568795c945Sjeremylt @param nnodes Total number of nodes 757a8de75f0Sjeremylt @param nqpts Total number of quadrature points 75895bb1877Svaleriabarra @param interp Row-major (nqpts * nnodes) matrix expressing the values of 7598795c945Sjeremylt nodal basis functions at quadrature points 76095bb1877Svaleriabarra @param grad Row-major (nqpts * dim * nnodes) matrix expressing 7618795c945Sjeremylt derivatives of nodal basis functions at quadrature points 7628795c945Sjeremylt @param qref Array of length nqpts holding the locations of quadrature 7638795c945Sjeremylt points on the reference element [-1, 1] 764a8de75f0Sjeremylt @param qweight Array of length nqpts holding the quadrature weights on the 765a8de75f0Sjeremylt reference element 766a8de75f0Sjeremylt @param[out] basis Address of the variable where the newly created 767a8de75f0Sjeremylt CeedBasis will be stored. 768a8de75f0Sjeremylt 769a8de75f0Sjeremylt @return An error code: 0 - success, otherwise - failure 770a8de75f0Sjeremylt 7717a982d89SJeremy L. Thompson @ref User 772a8de75f0Sjeremylt **/ 773a8de75f0Sjeremylt int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt ncomp, 774692c2638Sjeremylt CeedInt nnodes, CeedInt nqpts, const CeedScalar *interp, 775a8de75f0Sjeremylt const CeedScalar *grad, const CeedScalar *qref, 776a8de75f0Sjeremylt const CeedScalar *qweight, CeedBasis *basis) { 777a8de75f0Sjeremylt int ierr; 7788795c945Sjeremylt CeedInt P = nnodes, Q = nqpts, dim = 0; 779a8de75f0Sjeremylt 7805fe0d4faSjeremylt if (!ceed->BasisCreateH1) { 7815fe0d4faSjeremylt Ceed delegate; 782aefd8378Sjeremylt ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr); 7835fe0d4faSjeremylt 7845fe0d4faSjeremylt if (!delegate) 785c042f62fSJeremy L Thompson // LCOV_EXCL_START 786a8de75f0Sjeremylt return CeedError(ceed, 1, "Backend does not support BasisCreateH1"); 787c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 7885fe0d4faSjeremylt 7898795c945Sjeremylt ierr = CeedBasisCreateH1(delegate, topo, ncomp, nnodes, 7905fe0d4faSjeremylt nqpts, interp, grad, qref, 7915fe0d4faSjeremylt qweight, basis); CeedChk(ierr); 7925fe0d4faSjeremylt return 0; 7935fe0d4faSjeremylt } 7945fe0d4faSjeremylt 795a8de75f0Sjeremylt ierr = CeedCalloc(1,basis); CeedChk(ierr); 796a8de75f0Sjeremylt 797a8de75f0Sjeremylt ierr = CeedBasisGetTopologyDimension(topo, &dim); CeedChk(ierr); 798a8de75f0Sjeremylt 799a8de75f0Sjeremylt (*basis)->ceed = ceed; 800a8de75f0Sjeremylt ceed->refcount++; 801a8de75f0Sjeremylt (*basis)->refcount = 1; 802a8de75f0Sjeremylt (*basis)->tensorbasis = 0; 803a8de75f0Sjeremylt (*basis)->dim = dim; 804a8de75f0Sjeremylt (*basis)->ncomp = ncomp; 805a8de75f0Sjeremylt (*basis)->P = P; 806a8de75f0Sjeremylt (*basis)->Q = Q; 807a8de75f0Sjeremylt ierr = CeedMalloc(Q*dim,&(*basis)->qref1d); CeedChk(ierr); 808a8de75f0Sjeremylt ierr = CeedMalloc(Q,&(*basis)->qweight1d); CeedChk(ierr); 809a8de75f0Sjeremylt memcpy((*basis)->qref1d, qref, Q*dim*sizeof(qref[0])); 810a8de75f0Sjeremylt memcpy((*basis)->qweight1d, qweight, Q*sizeof(qweight[0])); 81100f91b2bSjeremylt ierr = CeedMalloc(Q*P, &(*basis)->interp); CeedChk(ierr); 81200f91b2bSjeremylt ierr = CeedMalloc(dim*Q*P, &(*basis)->grad); CeedChk(ierr); 81300f91b2bSjeremylt memcpy((*basis)->interp, interp, Q*P*sizeof(interp[0])); 81400f91b2bSjeremylt memcpy((*basis)->grad, grad, dim*Q*P*sizeof(grad[0])); 815667bc5fcSjeremylt ierr = ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, qref, 816a8de75f0Sjeremylt qweight, *basis); CeedChk(ierr); 817a8de75f0Sjeremylt return 0; 818a8de75f0Sjeremylt } 819a8de75f0Sjeremylt 820a8de75f0Sjeremylt /** 8217a982d89SJeremy L. Thompson @brief View a CeedBasis 8227a982d89SJeremy L. Thompson 8237a982d89SJeremy L. Thompson @param basis CeedBasis to view 8247a982d89SJeremy L. Thompson @param stream Stream to view to, e.g., stdout 8257a982d89SJeremy L. Thompson 8267a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 8277a982d89SJeremy L. Thompson 8287a982d89SJeremy L. Thompson @ref User 8297a982d89SJeremy L. Thompson **/ 8307a982d89SJeremy L. Thompson int CeedBasisView(CeedBasis basis, FILE *stream) { 8317a982d89SJeremy L. Thompson int ierr; 8327a982d89SJeremy L. Thompson 8337a982d89SJeremy L. Thompson if (basis->tensorbasis) { 8347a982d89SJeremy L. Thompson fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P1d, 8357a982d89SJeremy L. Thompson basis->Q1d); 8367a982d89SJeremy L. Thompson ierr = CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q1d, basis->qref1d, 8377a982d89SJeremy L. Thompson stream); CeedChk(ierr); 8387a982d89SJeremy L. Thompson ierr = CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q1d, 8397a982d89SJeremy L. Thompson basis->qweight1d, stream); CeedChk(ierr); 8407a982d89SJeremy L. Thompson ierr = CeedScalarView("interp1d", "\t% 12.8f", basis->Q1d, basis->P1d, 8417a982d89SJeremy L. Thompson basis->interp1d, stream); CeedChk(ierr); 8427a982d89SJeremy L. Thompson ierr = CeedScalarView("grad1d", "\t% 12.8f", basis->Q1d, basis->P1d, 8437a982d89SJeremy L. Thompson basis->grad1d, stream); CeedChk(ierr); 8447a982d89SJeremy L. Thompson } else { 8457a982d89SJeremy L. Thompson fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P, 8467a982d89SJeremy L. Thompson basis->Q); 8477a982d89SJeremy L. Thompson ierr = CeedScalarView("qref", "\t% 12.8f", 1, basis->Q*basis->dim, 8487a982d89SJeremy L. Thompson basis->qref1d, 8497a982d89SJeremy L. Thompson stream); CeedChk(ierr); 8507a982d89SJeremy L. Thompson ierr = CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->qweight1d, 8517a982d89SJeremy L. Thompson stream); CeedChk(ierr); 8527a982d89SJeremy L. Thompson ierr = CeedScalarView("interp", "\t% 12.8f", basis->Q, basis->P, 8537a982d89SJeremy L. Thompson basis->interp, stream); CeedChk(ierr); 8547a982d89SJeremy L. Thompson ierr = CeedScalarView("grad", "\t% 12.8f", basis->dim*basis->Q, basis->P, 8557a982d89SJeremy L. Thompson basis->grad, stream); CeedChk(ierr); 8567a982d89SJeremy L. Thompson } 8577a982d89SJeremy L. Thompson return 0; 8587a982d89SJeremy L. Thompson } 8597a982d89SJeremy L. Thompson 8607a982d89SJeremy L. Thompson /** 8617a982d89SJeremy L. Thompson @brief Get total number of nodes (in dim dimensions) of a CeedBasis 8627a982d89SJeremy L. Thompson 8637a982d89SJeremy L. Thompson @param basis CeedBasis 8647a982d89SJeremy L. Thompson @param[out] P Variable to store number of nodes 8657a982d89SJeremy L. Thompson 8667a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 8677a982d89SJeremy L. Thompson 8687a982d89SJeremy L. Thompson @ref Utility 8697a982d89SJeremy L. Thompson **/ 8707a982d89SJeremy L. Thompson int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) { 8717a982d89SJeremy L. Thompson *P = basis->P; 8727a982d89SJeremy L. Thompson return 0; 8737a982d89SJeremy L. Thompson } 8747a982d89SJeremy L. Thompson 8757a982d89SJeremy L. Thompson /** 8767a982d89SJeremy L. Thompson @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis 8777a982d89SJeremy L. Thompson 8787a982d89SJeremy L. Thompson @param basis CeedBasis 8797a982d89SJeremy L. Thompson @param[out] Q Variable to store number of quadrature points 8807a982d89SJeremy L. Thompson 8817a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 8827a982d89SJeremy L. Thompson 8837a982d89SJeremy L. Thompson @ref Utility 8847a982d89SJeremy L. Thompson **/ 8857a982d89SJeremy L. Thompson int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) { 8867a982d89SJeremy L. Thompson *Q = basis->Q; 8877a982d89SJeremy L. Thompson return 0; 8887a982d89SJeremy L. Thompson } 8897a982d89SJeremy L. Thompson 8907a982d89SJeremy L. Thompson /** 8917a982d89SJeremy L. Thompson @brief Apply basis evaluation from nodes to quadrature points or vice versa 8927a982d89SJeremy L. Thompson 8937a982d89SJeremy L. Thompson @param basis CeedBasis to evaluate 8947a982d89SJeremy L. Thompson @param nelem The number of elements to apply the basis evaluation to; 8957a982d89SJeremy L. Thompson the backend will specify the ordering in 8967a982d89SJeremy L. Thompson ElemRestrictionCreateBlocked 8977a982d89SJeremy L. Thompson @param tmode \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature 8987a982d89SJeremy L. Thompson points, \ref CEED_TRANSPOSE to apply the transpose, mapping 8997a982d89SJeremy L. Thompson from quadrature points to nodes 9007a982d89SJeremy L. Thompson @param emode \ref CEED_EVAL_NONE to use values directly, 9017a982d89SJeremy L. Thompson \ref CEED_EVAL_INTERP to use interpolated values, 9027a982d89SJeremy L. Thompson \ref CEED_EVAL_GRAD to use gradients, 9037a982d89SJeremy L. Thompson \ref CEED_EVAL_WEIGHT to use quadrature weights. 9047a982d89SJeremy L. Thompson @param[in] u Input CeedVector 9057a982d89SJeremy L. Thompson @param[out] v Output CeedVector 9067a982d89SJeremy L. Thompson 9077a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 9087a982d89SJeremy L. Thompson 9097a982d89SJeremy L. Thompson @ref User 9107a982d89SJeremy L. Thompson **/ 9117a982d89SJeremy L. Thompson int CeedBasisApply(CeedBasis basis, CeedInt nelem, CeedTransposeMode tmode, 9127a982d89SJeremy L. Thompson CeedEvalMode emode, CeedVector u, CeedVector v) { 9137a982d89SJeremy L. Thompson int ierr; 9147a982d89SJeremy L. Thompson CeedInt ulength = 0, vlength, nnodes, nqpt; 9157a982d89SJeremy L. Thompson if (!basis->Apply) 9167a982d89SJeremy L. Thompson // LCOV_EXCL_START 9177a982d89SJeremy L. Thompson return CeedError(basis->ceed, 1, "Backend does not support BasisApply"); 9187a982d89SJeremy L. Thompson // LCOV_EXCL_STOP 9197a982d89SJeremy L. Thompson 9207a982d89SJeremy L. Thompson // Check compatibility of topological and geometrical dimensions 9217a982d89SJeremy L. Thompson ierr = CeedBasisGetNumNodes(basis, &nnodes); CeedChk(ierr); 9227a982d89SJeremy L. Thompson ierr = CeedBasisGetNumQuadraturePoints(basis, &nqpt); CeedChk(ierr); 9237a982d89SJeremy L. Thompson ierr = CeedVectorGetLength(v, &vlength); CeedChk(ierr); 9247a982d89SJeremy L. Thompson 9257a982d89SJeremy L. Thompson if (u) { 9267a982d89SJeremy L. Thompson ierr = CeedVectorGetLength(u, &ulength); CeedChk(ierr); 9277a982d89SJeremy L. Thompson } 9287a982d89SJeremy L. Thompson 9297a982d89SJeremy L. Thompson if ((tmode == CEED_TRANSPOSE && (vlength%nnodes != 0 || ulength%nqpt != 0)) || 9307a982d89SJeremy L. Thompson (tmode == CEED_NOTRANSPOSE && (ulength%nnodes != 0 || vlength%nqpt != 0))) 9317a982d89SJeremy L. Thompson return CeedError(basis->ceed, 1, "Length of input/output vectors " 9327a982d89SJeremy L. Thompson "incompatible with basis dimensions"); 9337a982d89SJeremy L. Thompson 9347a982d89SJeremy L. Thompson ierr = basis->Apply(basis, nelem, tmode, emode, u, v); CeedChk(ierr); 9357a982d89SJeremy L. Thompson return 0; 9367a982d89SJeremy L. Thompson } 9377a982d89SJeremy L. Thompson 9387a982d89SJeremy L. Thompson /** 9397a982d89SJeremy L. Thompson @brief Destroy a CeedBasis 9407a982d89SJeremy L. Thompson 9417a982d89SJeremy L. Thompson @param basis CeedBasis to destroy 9427a982d89SJeremy L. Thompson 9437a982d89SJeremy L. Thompson @return An error code: 0 - success, otherwise - failure 9447a982d89SJeremy L. Thompson 9457a982d89SJeremy L. Thompson @ref User 9467a982d89SJeremy L. Thompson **/ 9477a982d89SJeremy L. Thompson int CeedBasisDestroy(CeedBasis *basis) { 9487a982d89SJeremy L. Thompson int ierr; 9497a982d89SJeremy L. Thompson 9507a982d89SJeremy L. Thompson if (!*basis || --(*basis)->refcount > 0) 9517a982d89SJeremy L. Thompson return 0; 9527a982d89SJeremy L. Thompson if ((*basis)->Destroy) { 9537a982d89SJeremy L. Thompson ierr = (*basis)->Destroy(*basis); CeedChk(ierr); 9547a982d89SJeremy L. Thompson } 9557a982d89SJeremy L. Thompson ierr = CeedFree(&(*basis)->interp); CeedChk(ierr); 9567a982d89SJeremy L. Thompson ierr = CeedFree(&(*basis)->interp1d); CeedChk(ierr); 9577a982d89SJeremy L. Thompson ierr = CeedFree(&(*basis)->grad); CeedChk(ierr); 9587a982d89SJeremy L. Thompson ierr = CeedFree(&(*basis)->grad1d); CeedChk(ierr); 9597a982d89SJeremy L. Thompson ierr = CeedFree(&(*basis)->qref1d); CeedChk(ierr); 9607a982d89SJeremy L. Thompson ierr = CeedFree(&(*basis)->qweight1d); CeedChk(ierr); 9617a982d89SJeremy L. Thompson ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr); 9627a982d89SJeremy L. Thompson ierr = CeedFree(basis); CeedChk(ierr); 9637a982d89SJeremy L. Thompson return 0; 9647a982d89SJeremy L. Thompson } 9657a982d89SJeremy L. Thompson 9667a982d89SJeremy L. Thompson /** 967b11c1e72Sjeremylt @brief Construct a Gauss-Legendre quadrature 968b11c1e72Sjeremylt 969b11c1e72Sjeremylt @param Q Number of quadrature points (integrates polynomials of 970b11c1e72Sjeremylt degree 2*Q-1 exactly) 971b11c1e72Sjeremylt @param[out] qref1d Array of length Q to hold the abscissa on [-1, 1] 972b11c1e72Sjeremylt @param[out] qweight1d Array of length Q to hold the weights 973b11c1e72Sjeremylt 974b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 975dfdf5a53Sjeremylt 976dfdf5a53Sjeremylt @ref Utility 977b11c1e72Sjeremylt **/ 978d7b241e6Sjeremylt int CeedGaussQuadrature(CeedInt Q, CeedScalar *qref1d, CeedScalar *qweight1d) { 979d7b241e6Sjeremylt // Allocate 980d7b241e6Sjeremylt CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0*atan(1.0); 981d7b241e6Sjeremylt // Build qref1d, qweight1d 982d7b241e6Sjeremylt for (int i = 0; i <= Q/2; i++) { 983d7b241e6Sjeremylt // Guess 984d7b241e6Sjeremylt xi = cos(PI*(CeedScalar)(2*i+1)/((CeedScalar)(2*Q))); 985d7b241e6Sjeremylt // Pn(xi) 986d7b241e6Sjeremylt P0 = 1.0; 987d7b241e6Sjeremylt P1 = xi; 988d7b241e6Sjeremylt P2 = 0.0; 989d7b241e6Sjeremylt for (int j = 2; j <= Q; j++) { 990d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 991d7b241e6Sjeremylt P0 = P1; 992d7b241e6Sjeremylt P1 = P2; 993d7b241e6Sjeremylt } 994d7b241e6Sjeremylt // First Newton Step 995d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 996d7b241e6Sjeremylt xi = xi-P2/dP2; 997d7b241e6Sjeremylt // Newton to convergence 9980e4d4210Sjeremylt for (int k=0; k<100 && fabs(P2)>10*CEED_EPSILON; k++) { 999d7b241e6Sjeremylt P0 = 1.0; 1000d7b241e6Sjeremylt P1 = xi; 1001d7b241e6Sjeremylt for (int j = 2; j <= Q; j++) { 1002d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 1003d7b241e6Sjeremylt P0 = P1; 1004d7b241e6Sjeremylt P1 = P2; 1005d7b241e6Sjeremylt } 1006d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 1007d7b241e6Sjeremylt xi = xi-P2/dP2; 1008d7b241e6Sjeremylt } 1009d7b241e6Sjeremylt // Save xi, wi 1010d7b241e6Sjeremylt wi = 2.0/((1.0-xi*xi)*dP2*dP2); 1011d7b241e6Sjeremylt qweight1d[i] = wi; 1012d7b241e6Sjeremylt qweight1d[Q-1-i] = wi; 1013d7b241e6Sjeremylt qref1d[i] = -xi; 1014d7b241e6Sjeremylt qref1d[Q-1-i]= xi; 1015d7b241e6Sjeremylt } 1016d7b241e6Sjeremylt return 0; 1017d7b241e6Sjeremylt } 1018d7b241e6Sjeremylt 1019b11c1e72Sjeremylt /** 1020b11c1e72Sjeremylt @brief Construct a Gauss-Legendre-Lobatto quadrature 1021b11c1e72Sjeremylt 1022b11c1e72Sjeremylt @param Q Number of quadrature points (integrates polynomials of 1023b11c1e72Sjeremylt degree 2*Q-3 exactly) 1024b11c1e72Sjeremylt @param[out] qref1d Array of length Q to hold the abscissa on [-1, 1] 1025b11c1e72Sjeremylt @param[out] qweight1d Array of length Q to hold the weights 1026b11c1e72Sjeremylt 1027b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 1028dfdf5a53Sjeremylt 1029dfdf5a53Sjeremylt @ref Utility 1030b11c1e72Sjeremylt **/ 1031d7b241e6Sjeremylt int CeedLobattoQuadrature(CeedInt Q, CeedScalar *qref1d, 1032d7b241e6Sjeremylt CeedScalar *qweight1d) { 1033d7b241e6Sjeremylt // Allocate 1034d7b241e6Sjeremylt CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0*atan(1.0); 1035d7b241e6Sjeremylt // Build qref1d, qweight1d 1036d7b241e6Sjeremylt // Set endpoints 103730a100c3SJed Brown if (Q < 2) 10387ed52d01Sjeremylt return CeedError(NULL, 1, 10397ed52d01Sjeremylt "Cannot create Lobatto quadrature with Q=%d < 2 points", Q); 1040d7b241e6Sjeremylt wi = 2.0/((CeedScalar)(Q*(Q-1))); 1041d7b241e6Sjeremylt if (qweight1d) { 1042d7b241e6Sjeremylt qweight1d[0] = wi; 1043d7b241e6Sjeremylt qweight1d[Q-1] = wi; 1044d7b241e6Sjeremylt } 1045d7b241e6Sjeremylt qref1d[0] = -1.0; 1046d7b241e6Sjeremylt qref1d[Q-1] = 1.0; 1047d7b241e6Sjeremylt // Interior 1048d7b241e6Sjeremylt for (int i = 1; i <= (Q-1)/2; i++) { 1049d7b241e6Sjeremylt // Guess 1050d7b241e6Sjeremylt xi = cos(PI*(CeedScalar)(i)/(CeedScalar)(Q-1)); 1051d7b241e6Sjeremylt // Pn(xi) 1052d7b241e6Sjeremylt P0 = 1.0; 1053d7b241e6Sjeremylt P1 = xi; 1054d7b241e6Sjeremylt P2 = 0.0; 1055d7b241e6Sjeremylt for (int j = 2; j < Q; j++) { 1056d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 1057d7b241e6Sjeremylt P0 = P1; 1058d7b241e6Sjeremylt P1 = P2; 1059d7b241e6Sjeremylt } 1060d7b241e6Sjeremylt // First Newton step 1061d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 1062d7b241e6Sjeremylt d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi); 1063d7b241e6Sjeremylt xi = xi-dP2/d2P2; 1064d7b241e6Sjeremylt // Newton to convergence 10650e4d4210Sjeremylt for (int k=0; k<100 && fabs(dP2)>10*CEED_EPSILON; k++) { 1066d7b241e6Sjeremylt P0 = 1.0; 1067d7b241e6Sjeremylt P1 = xi; 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 dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 1074d7b241e6Sjeremylt d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi); 1075d7b241e6Sjeremylt xi = xi-dP2/d2P2; 1076d7b241e6Sjeremylt } 1077d7b241e6Sjeremylt // Save xi, wi 1078d7b241e6Sjeremylt wi = 2.0/(((CeedScalar)(Q*(Q-1)))*P2*P2); 1079d7b241e6Sjeremylt if (qweight1d) { 1080d7b241e6Sjeremylt qweight1d[i] = wi; 1081d7b241e6Sjeremylt qweight1d[Q-1-i] = wi; 1082d7b241e6Sjeremylt } 1083d7b241e6Sjeremylt qref1d[i] = -xi; 1084d7b241e6Sjeremylt qref1d[Q-1-i]= xi; 1085d7b241e6Sjeremylt } 1086d7b241e6Sjeremylt return 0; 1087d7b241e6Sjeremylt } 1088d7b241e6Sjeremylt 1089dfdf5a53Sjeremylt /** 109095bb1877Svaleriabarra @brief Return QR Factorization of a matrix 1091b11c1e72Sjeremylt 109277645d7bSjeremylt @param ceed A Ceed context for error handling 109352bfb9bbSJeremy L Thompson @param[in,out] mat Row-major matrix to be factorized in place 109452bfb9bbSJeremy L Thompson @param[in,out] tau Vector of length m of scaling factors 1095b11c1e72Sjeremylt @param m Number of rows 1096b11c1e72Sjeremylt @param n Number of columns 1097b11c1e72Sjeremylt 1098b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 1099dfdf5a53Sjeremylt 1100dfdf5a53Sjeremylt @ref Utility 1101b11c1e72Sjeremylt **/ 1102a7bd39daSjeremylt int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau, 1103d7b241e6Sjeremylt CeedInt m, CeedInt n) { 1104d7b241e6Sjeremylt CeedScalar v[m]; 1105d7b241e6Sjeremylt 1106a7bd39daSjeremylt // Check m >= n 1107a7bd39daSjeremylt if (n > m) 1108c042f62fSJeremy L Thompson // LCOV_EXCL_START 1109a7bd39daSjeremylt return CeedError(ceed, 1, "Cannot compute QR factorization with n > m"); 1110c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 1111a7bd39daSjeremylt 111252bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) { 1113d7b241e6Sjeremylt // Calculate Householder vector, magnitude 1114d7b241e6Sjeremylt CeedScalar sigma = 0.0; 1115d7b241e6Sjeremylt v[i] = mat[i+n*i]; 111652bfb9bbSJeremy L Thompson for (CeedInt j=i+1; j<m; j++) { 1117d7b241e6Sjeremylt v[j] = mat[i+n*j]; 1118d7b241e6Sjeremylt sigma += v[j] * v[j]; 1119d7b241e6Sjeremylt } 1120d7b241e6Sjeremylt CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m] 1121d7b241e6Sjeremylt CeedScalar Rii = -copysign(norm, v[i]); 1122d7b241e6Sjeremylt v[i] -= Rii; 1123d7b241e6Sjeremylt // norm of v[i:m] after modification above and scaling below 1124d7b241e6Sjeremylt // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 1125d7b241e6Sjeremylt // tau = 2 / (norm*norm) 1126d7b241e6Sjeremylt tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma); 1127fb551037Sjeremylt 11281d102b48SJeremy L Thompson for (CeedInt j=i+1; j<m; j++) 11291d102b48SJeremy L Thompson v[j] /= v[i]; 1130d7b241e6Sjeremylt 1131d7b241e6Sjeremylt // Apply Householder reflector to lower right panel 1132d7b241e6Sjeremylt CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1); 1133d7b241e6Sjeremylt // Save v 1134d7b241e6Sjeremylt mat[i+n*i] = Rii; 11351d102b48SJeremy L Thompson for (CeedInt j=i+1; j<m; j++) 1136d7b241e6Sjeremylt mat[i+n*j] = v[j]; 1137d7b241e6Sjeremylt } 1138d7b241e6Sjeremylt 1139d7b241e6Sjeremylt return 0; 1140d7b241e6Sjeremylt } 1141d7b241e6Sjeremylt 1142b11c1e72Sjeremylt /** 114352bfb9bbSJeremy L Thompson @brief Return symmetric Schur decomposition of the symmetric matrix mat via 114452bfb9bbSJeremy L Thompson symmetric QR factorization 114552bfb9bbSJeremy L Thompson 114677645d7bSjeremylt @param ceed A Ceed context for error handling 114752bfb9bbSJeremy L Thompson @param[in,out] mat Row-major matrix to be factorized in place 1148460bf743SValeria Barra @param[out] lambda Vector of length n of eigenvalues 114952bfb9bbSJeremy L Thompson @param n Number of rows/columns 115052bfb9bbSJeremy L Thompson 115152bfb9bbSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 115252bfb9bbSJeremy L Thompson 115352bfb9bbSJeremy L Thompson @ref Utility 115452bfb9bbSJeremy L Thompson **/ 115552bfb9bbSJeremy L Thompson int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat, 115652bfb9bbSJeremy L Thompson CeedScalar *lambda, CeedInt n) { 115752bfb9bbSJeremy L Thompson // Check bounds for clang-tidy 115852bfb9bbSJeremy L Thompson if (n<2) 1159c042f62fSJeremy L Thompson // LCOV_EXCL_START 1160c042f62fSJeremy L Thompson return CeedError(ceed, 1, 1161c042f62fSJeremy L Thompson "Cannot compute symmetric Schur decomposition of scalars"); 1162c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 116352bfb9bbSJeremy L Thompson 116452bfb9bbSJeremy L Thompson CeedScalar v[n-1], tau[n-1], matT[n*n]; 116552bfb9bbSJeremy L Thompson 116652bfb9bbSJeremy L Thompson // Copy mat to matT and set mat to I 116752bfb9bbSJeremy L Thompson memcpy(matT, mat, n*n*sizeof(mat[0])); 116852bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) 116952bfb9bbSJeremy L Thompson for (CeedInt j=0; j<n; j++) 117052bfb9bbSJeremy L Thompson mat[j+n*i] = (i==j) ? 1 : 0; 117152bfb9bbSJeremy L Thompson 117252bfb9bbSJeremy L Thompson // Reduce to tridiagonal 117352bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n-1; i++) { 117452bfb9bbSJeremy L Thompson // Calculate Householder vector, magnitude 117552bfb9bbSJeremy L Thompson CeedScalar sigma = 0.0; 117652bfb9bbSJeremy L Thompson v[i] = matT[i+n*(i+1)]; 117752bfb9bbSJeremy L Thompson for (CeedInt j=i+1; j<n-1; j++) { 117852bfb9bbSJeremy L Thompson v[j] = matT[i+n*(j+1)]; 117952bfb9bbSJeremy L Thompson sigma += v[j] * v[j]; 118052bfb9bbSJeremy L Thompson } 118152bfb9bbSJeremy L Thompson CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:n-1] 118252bfb9bbSJeremy L Thompson CeedScalar Rii = -copysign(norm, v[i]); 118352bfb9bbSJeremy L Thompson v[i] -= Rii; 118452bfb9bbSJeremy L Thompson // norm of v[i:m] after modification above and scaling below 118552bfb9bbSJeremy L Thompson // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 118652bfb9bbSJeremy L Thompson // tau = 2 / (norm*norm) 11870e4d4210Sjeremylt if (sigma > 10*CEED_EPSILON) 118852bfb9bbSJeremy L Thompson tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma); 1189fb551037Sjeremylt else 1190fb551037Sjeremylt tau[i] = 0; 1191fb551037Sjeremylt 1192fb551037Sjeremylt for (CeedInt j=i+1; j<n-1; j++) 1193fb551037Sjeremylt v[j] /= v[i]; 119452bfb9bbSJeremy L Thompson 119552bfb9bbSJeremy L Thompson // Update sub and super diagonal 119652bfb9bbSJeremy L Thompson matT[i+n*(i+1)] = Rii; 119752bfb9bbSJeremy L Thompson matT[(i+1)+n*i] = Rii; 119852bfb9bbSJeremy L Thompson for (CeedInt j=i+2; j<n; j++) { 119952bfb9bbSJeremy L Thompson matT[i+n*j] = 0; matT[j+n*i] = 0; 120052bfb9bbSJeremy L Thompson } 120152bfb9bbSJeremy L Thompson // Apply symmetric Householder reflector to lower right panel 120252bfb9bbSJeremy L Thompson CeedHouseholderReflect(&matT[(i+1)+n*(i+1)], &v[i], tau[i], 120352bfb9bbSJeremy L Thompson n-(i+1), n-(i+1), n, 1); 120452bfb9bbSJeremy L Thompson CeedHouseholderReflect(&matT[(i+1)+n*(i+1)], &v[i], tau[i], 120552bfb9bbSJeremy L Thompson n-(i+1), n-(i+1), 1, n); 120652bfb9bbSJeremy L Thompson // Save v 120752bfb9bbSJeremy L Thompson for (CeedInt j=i+1; j<n-1; j++) { 120852bfb9bbSJeremy L Thompson matT[i+n*(j+1)] = v[j]; 120952bfb9bbSJeremy L Thompson } 121052bfb9bbSJeremy L Thompson } 121152bfb9bbSJeremy L Thompson // Backwards accumulation of Q 121252bfb9bbSJeremy L Thompson for (CeedInt i=n-2; i>=0; i--) { 121352bfb9bbSJeremy L Thompson v[i] = 1; 121452bfb9bbSJeremy L Thompson for (CeedInt j=i+1; j<n-1; j++) { 121552bfb9bbSJeremy L Thompson v[j] = matT[i+n*(j+1)]; 121652bfb9bbSJeremy L Thompson matT[i+n*(j+1)] = 0; 121752bfb9bbSJeremy L Thompson } 121852bfb9bbSJeremy L Thompson CeedHouseholderReflect(&mat[(i+1)+n*(i+1)], &v[i], tau[i], 121952bfb9bbSJeremy L Thompson n-(i+1), n-(i+1), n, 1); 122052bfb9bbSJeremy L Thompson } 122152bfb9bbSJeremy L Thompson 122252bfb9bbSJeremy L Thompson // Reduce sub and super diagonal 122352bfb9bbSJeremy L Thompson CeedInt p = 0, q = 0, itr = 0, maxitr = n*n*n; 12240e4d4210Sjeremylt CeedScalar tol = 10*CEED_EPSILON; 122552bfb9bbSJeremy L Thompson 122652bfb9bbSJeremy L Thompson while (q < n && itr < maxitr) { 122752bfb9bbSJeremy L Thompson // Update p, q, size of reduced portions of diagonal 122852bfb9bbSJeremy L Thompson p = 0; q = 0; 122952bfb9bbSJeremy L Thompson for (CeedInt i=n-2; i>=0; i--) { 123052bfb9bbSJeremy L Thompson if (fabs(matT[i+n*(i+1)]) < tol) 123152bfb9bbSJeremy L Thompson q += 1; 123252bfb9bbSJeremy L Thompson else 123352bfb9bbSJeremy L Thompson break; 123452bfb9bbSJeremy L Thompson } 123552bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n-1-q; i++) { 123652bfb9bbSJeremy L Thompson if (fabs(matT[i+n*(i+1)]) < tol) 123752bfb9bbSJeremy L Thompson p += 1; 123852bfb9bbSJeremy L Thompson else 123952bfb9bbSJeremy L Thompson break; 124052bfb9bbSJeremy L Thompson } 124152bfb9bbSJeremy L Thompson if (q == n-1) break; // Finished reducing 124252bfb9bbSJeremy L Thompson 124352bfb9bbSJeremy L Thompson // Reduce tridiagonal portion 124452bfb9bbSJeremy L Thompson CeedScalar tnn = matT[(n-1-q)+n*(n-1-q)], 124552bfb9bbSJeremy L Thompson tnnm1 = matT[(n-2-q)+n*(n-1-q)]; 124652bfb9bbSJeremy L Thompson CeedScalar d = (matT[(n-2-q)+n*(n-2-q)] - tnn)/2; 124752bfb9bbSJeremy L Thompson CeedScalar mu = tnn - tnnm1*tnnm1 / 124852bfb9bbSJeremy L Thompson (d + copysign(sqrt(d*d + tnnm1*tnnm1), d)); 124952bfb9bbSJeremy L Thompson CeedScalar x = matT[p+n*p] - mu; 125052bfb9bbSJeremy L Thompson CeedScalar z = matT[p+n*(p+1)]; 125152bfb9bbSJeremy L Thompson for (CeedInt k=p; k<n-1-q; k++) { 125252bfb9bbSJeremy L Thompson // Compute Givens rotation 125352bfb9bbSJeremy L Thompson CeedScalar c = 1, s = 0; 125452bfb9bbSJeremy L Thompson if (fabs(z) > tol) { 125552bfb9bbSJeremy L Thompson if (fabs(z) > fabs(x)) { 125652bfb9bbSJeremy L Thompson CeedScalar tau = -x/z; 125752bfb9bbSJeremy L Thompson s = 1/sqrt(1+tau*tau), c = s*tau; 125852bfb9bbSJeremy L Thompson } else { 125952bfb9bbSJeremy L Thompson CeedScalar tau = -z/x; 126052bfb9bbSJeremy L Thompson c = 1/sqrt(1+tau*tau), s = c*tau; 126152bfb9bbSJeremy L Thompson } 126252bfb9bbSJeremy L Thompson } 126352bfb9bbSJeremy L Thompson 126452bfb9bbSJeremy L Thompson // Apply Givens rotation to T 126552bfb9bbSJeremy L Thompson CeedGivensRotation(matT, c, s, CEED_NOTRANSPOSE, k, k+1, n, n); 126652bfb9bbSJeremy L Thompson CeedGivensRotation(matT, c, s, CEED_TRANSPOSE, k, k+1, n, n); 126752bfb9bbSJeremy L Thompson 126852bfb9bbSJeremy L Thompson // Apply Givens rotation to Q 126952bfb9bbSJeremy L Thompson CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k+1, n, n); 127052bfb9bbSJeremy L Thompson 127152bfb9bbSJeremy L Thompson // Update x, z 127252bfb9bbSJeremy L Thompson if (k < n-q-2) { 127352bfb9bbSJeremy L Thompson x = matT[k+n*(k+1)]; 127452bfb9bbSJeremy L Thompson z = matT[k+n*(k+2)]; 127552bfb9bbSJeremy L Thompson } 127652bfb9bbSJeremy L Thompson } 127752bfb9bbSJeremy L Thompson itr++; 127852bfb9bbSJeremy L Thompson } 127952bfb9bbSJeremy L Thompson // Save eigenvalues 128052bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) 128152bfb9bbSJeremy L Thompson lambda[i] = matT[i+n*i]; 128252bfb9bbSJeremy L Thompson 128352bfb9bbSJeremy L Thompson // Check convergence 128452bfb9bbSJeremy L Thompson if (itr == maxitr && q < n-1) 1285c042f62fSJeremy L Thompson // LCOV_EXCL_START 128652bfb9bbSJeremy L Thompson return CeedError(ceed, 1, "Symmetric QR failed to converge"); 1287c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 128852bfb9bbSJeremy L Thompson 128952bfb9bbSJeremy L Thompson return 0; 129052bfb9bbSJeremy L Thompson } 129152bfb9bbSJeremy L Thompson 129252bfb9bbSJeremy L Thompson /** 129352bfb9bbSJeremy L Thompson @brief Return Simultaneous Diagonalization of two matrices. This solves the 129452bfb9bbSJeremy L Thompson generalized eigenvalue problem A x = lambda B x, where A and B 129552bfb9bbSJeremy L Thompson are symmetric and B is positive definite. We generate the matrix X 129652bfb9bbSJeremy L Thompson and vector Lambda such that X^T A X = Lambda and X^T B X = I. This 129752bfb9bbSJeremy L Thompson is equivalent to the LAPACK routine 'sygv' with TYPE = 1. 129852bfb9bbSJeremy L Thompson 129977645d7bSjeremylt @param ceed A Ceed context for error handling 130052bfb9bbSJeremy L Thompson @param[in] matA Row-major matrix to be factorized with eigenvalues 130152bfb9bbSJeremy L Thompson @param[in] matB Row-major matrix to be factorized to identity 130252bfb9bbSJeremy L Thompson @param[out] x Row-major orthogonal matrix 1303460bf743SValeria Barra @param[out] lambda Vector of length n of generalized eigenvalues 130452bfb9bbSJeremy L Thompson @param n Number of rows/columns 130552bfb9bbSJeremy L Thompson 130652bfb9bbSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 130752bfb9bbSJeremy L Thompson 130852bfb9bbSJeremy L Thompson @ref Utility 130952bfb9bbSJeremy L Thompson **/ 131052bfb9bbSJeremy L Thompson int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *matA, 131152bfb9bbSJeremy L Thompson CeedScalar *matB, CeedScalar *x, 131252bfb9bbSJeremy L Thompson CeedScalar *lambda, CeedInt n) { 131352bfb9bbSJeremy L Thompson int ierr; 131452bfb9bbSJeremy L Thompson CeedScalar matC[n*n], matG[n*n], vecD[n]; 131552bfb9bbSJeremy L Thompson 131652bfb9bbSJeremy L Thompson // Compute B = G D G^T 131752bfb9bbSJeremy L Thompson memcpy(matG, matB, n*n*sizeof(matB[0])); 131852bfb9bbSJeremy L Thompson ierr = CeedSymmetricSchurDecomposition(ceed, matG, vecD, n); CeedChk(ierr); 1319fb551037Sjeremylt for (CeedInt i=0; i<n; i++) 1320fb551037Sjeremylt vecD[i] = sqrt(vecD[i]); 132152bfb9bbSJeremy L Thompson 1322fb551037Sjeremylt // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T 1323fb551037Sjeremylt // = D^-1/2 G^T A G D^-1/2 132452bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) 132552bfb9bbSJeremy L Thompson for (CeedInt j=0; j<n; j++) 1326fb551037Sjeremylt matC[j+i*n] = matG[i+j*n] / vecD[i]; 13279289e5bfSjeremylt ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)matC, 13289289e5bfSjeremylt (const CeedScalar *)matA, x, n, n, n); 13299289e5bfSjeremylt CeedChk(ierr); 133052bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) 133152bfb9bbSJeremy L Thompson for (CeedInt j=0; j<n; j++) 1332fb551037Sjeremylt matG[j+i*n] = matG[j+i*n] / vecD[j]; 13339289e5bfSjeremylt ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)x, 13349289e5bfSjeremylt (const CeedScalar *)matG, matC, n, n, n); 13359289e5bfSjeremylt CeedChk(ierr); 133652bfb9bbSJeremy L Thompson 133752bfb9bbSJeremy L Thompson // Compute Q^T C Q = lambda 133852bfb9bbSJeremy L Thompson ierr = CeedSymmetricSchurDecomposition(ceed, matC, lambda, n); CeedChk(ierr); 133952bfb9bbSJeremy L Thompson 1340fb551037Sjeremylt // Set x = (G D^1/2)^-T Q 1341fb551037Sjeremylt // = G D^-1/2 Q 13429289e5bfSjeremylt ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)matG, 13439289e5bfSjeremylt (const CeedScalar *)matC, x, n, n, n); 13449289e5bfSjeremylt CeedChk(ierr); 134552bfb9bbSJeremy L Thompson 134652bfb9bbSJeremy L Thompson return 0; 134752bfb9bbSJeremy L Thompson } 134852bfb9bbSJeremy L Thompson 1349d7b241e6Sjeremylt /// @} 1350