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 24d7b241e6Sjeremylt /// @cond DOXYGEN_SKIP 25783c99b3SValeria Barra static struct CeedBasis_private ceed_basis_collocated; 26d7b241e6Sjeremylt /// @endcond 27d7b241e6Sjeremylt 28d7b241e6Sjeremylt /// @file 29d7b241e6Sjeremylt /// Implementation of public CeedBasis interfaces 30d7b241e6Sjeremylt /// 31dfdf5a53Sjeremylt /// @addtogroup CeedBasis 32d7b241e6Sjeremylt /// @{ 33d7b241e6Sjeremylt 34b11c1e72Sjeremylt /** 3595bb1877Svaleriabarra @brief Create a tensor-product basis for H^1 discretizations 36b11c1e72Sjeremylt 37b11c1e72Sjeremylt @param ceed A Ceed object where the CeedBasis will be created 38b11c1e72Sjeremylt @param dim Topological dimension 39b11c1e72Sjeremylt @param ncomp Number of field components (1 for scalar fields) 40b11c1e72Sjeremylt @param P1d Number of nodes in one dimension 41b11c1e72Sjeremylt @param Q1d Number of quadrature points in one dimension 4295bb1877Svaleriabarra @param interp1d Row-major (Q1d * P1d) matrix expressing the values of nodal 43b11c1e72Sjeremylt basis functions at quadrature points 4495bb1877Svaleriabarra @param grad1d Row-major (Q1d * P1d) matrix expressing derivatives of nodal 45b11c1e72Sjeremylt basis functions at quadrature points 46b11c1e72Sjeremylt @param qref1d Array of length Q1d holding the locations of quadrature points 47b11c1e72Sjeremylt on the 1D reference element [-1, 1] 48b11c1e72Sjeremylt @param qweight1d Array of length Q1d holding the quadrature weights on the 49b11c1e72Sjeremylt reference element 50b11c1e72Sjeremylt @param[out] basis Address of the variable where the newly created 51b11c1e72Sjeremylt CeedBasis will be stored. 52b11c1e72Sjeremylt 53b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 54dfdf5a53Sjeremylt 55dfdf5a53Sjeremylt @ref Basic 56b11c1e72Sjeremylt **/ 57d7b241e6Sjeremylt int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt ncomp, CeedInt P1d, 58d7b241e6Sjeremylt CeedInt Q1d, const CeedScalar *interp1d, 59d7b241e6Sjeremylt const CeedScalar *grad1d, const CeedScalar *qref1d, 60d7b241e6Sjeremylt const CeedScalar *qweight1d, CeedBasis *basis) { 61d7b241e6Sjeremylt int ierr; 62d7b241e6Sjeremylt 634d537eeaSYohann if (dim<1) 64c042f62fSJeremy L Thompson // LCOV_EXCL_START 654d537eeaSYohann return CeedError(ceed, 1, "Basis dimension must be a positive value"); 66c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 674d537eeaSYohann 685fe0d4faSjeremylt if (!ceed->BasisCreateTensorH1) { 695fe0d4faSjeremylt Ceed delegate; 70aefd8378Sjeremylt ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr); 715fe0d4faSjeremylt 725fe0d4faSjeremylt if (!delegate) 73c042f62fSJeremy L Thompson // LCOV_EXCL_START 74d7b241e6Sjeremylt return CeedError(ceed, 1, "Backend does not support BasisCreateTensorH1"); 75c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 765fe0d4faSjeremylt 775fe0d4faSjeremylt ierr = CeedBasisCreateTensorH1(delegate, dim, ncomp, P1d, 785fe0d4faSjeremylt Q1d, interp1d, grad1d, qref1d, 795fe0d4faSjeremylt qweight1d, basis); CeedChk(ierr); 805fe0d4faSjeremylt return 0; 815fe0d4faSjeremylt } 82d7b241e6Sjeremylt ierr = CeedCalloc(1,basis); CeedChk(ierr); 83d7b241e6Sjeremylt (*basis)->ceed = ceed; 84d7b241e6Sjeremylt ceed->refcount++; 85d7b241e6Sjeremylt (*basis)->refcount = 1; 86a8de75f0Sjeremylt (*basis)->tensorbasis = 1; 87d7b241e6Sjeremylt (*basis)->dim = dim; 88d7b241e6Sjeremylt (*basis)->ncomp = ncomp; 89d7b241e6Sjeremylt (*basis)->P1d = P1d; 90d7b241e6Sjeremylt (*basis)->Q1d = Q1d; 91a8de75f0Sjeremylt (*basis)->P = CeedIntPow(P1d, dim); 92a8de75f0Sjeremylt (*basis)->Q = CeedIntPow(Q1d, dim); 93d7b241e6Sjeremylt ierr = CeedMalloc(Q1d,&(*basis)->qref1d); CeedChk(ierr); 94d7b241e6Sjeremylt ierr = CeedMalloc(Q1d,&(*basis)->qweight1d); CeedChk(ierr); 95d7b241e6Sjeremylt memcpy((*basis)->qref1d, qref1d, Q1d*sizeof(qref1d[0])); 96d7b241e6Sjeremylt memcpy((*basis)->qweight1d, qweight1d, Q1d*sizeof(qweight1d[0])); 97d7b241e6Sjeremylt ierr = CeedMalloc(Q1d*P1d,&(*basis)->interp1d); CeedChk(ierr); 98d7b241e6Sjeremylt ierr = CeedMalloc(Q1d*P1d,&(*basis)->grad1d); CeedChk(ierr); 99d7b241e6Sjeremylt memcpy((*basis)->interp1d, interp1d, Q1d*P1d*sizeof(interp1d[0])); 10009486605Sjeremylt memcpy((*basis)->grad1d, grad1d, Q1d*P1d*sizeof(grad1d[0])); 101667bc5fcSjeremylt ierr = ceed->BasisCreateTensorH1(dim, P1d, Q1d, interp1d, grad1d, qref1d, 102d7b241e6Sjeremylt qweight1d, *basis); CeedChk(ierr); 103d7b241e6Sjeremylt return 0; 104d7b241e6Sjeremylt } 105d7b241e6Sjeremylt 106b11c1e72Sjeremylt /** 10795bb1877Svaleriabarra @brief Create a tensor-product Lagrange basis 108b11c1e72Sjeremylt 109b11c1e72Sjeremylt @param ceed A Ceed object where the CeedBasis will be created 110b11c1e72Sjeremylt @param dim Topological dimension of element 11195bb1877Svaleriabarra @param ncomp Number of field components (1 for scalar fields) 112b11c1e72Sjeremylt @param P Number of Gauss-Lobatto nodes in one dimension. The 113b11c1e72Sjeremylt polynomial degree of the resulting Q_k element is k=P-1. 114b11c1e72Sjeremylt @param Q Number of quadrature points in one dimension. 115b11c1e72Sjeremylt @param qmode Distribution of the Q quadrature points (affects order of 116b11c1e72Sjeremylt accuracy for the quadrature) 117b11c1e72Sjeremylt @param[out] basis Address of the variable where the newly created 118b11c1e72Sjeremylt CeedBasis will be stored. 119b11c1e72Sjeremylt 120b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 121dfdf5a53Sjeremylt 122dfdf5a53Sjeremylt @ref Basic 123b11c1e72Sjeremylt **/ 124d7b241e6Sjeremylt int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt ncomp, 125692c2638Sjeremylt CeedInt P, CeedInt Q, CeedQuadMode qmode, 126692c2638Sjeremylt CeedBasis *basis) { 127d7b241e6Sjeremylt // Allocate 128d7b241e6Sjeremylt int ierr, i, j, k; 129d7b241e6Sjeremylt CeedScalar c1, c2, c3, c4, dx, *nodes, *interp1d, *grad1d, *qref1d, *qweight1d; 1304d537eeaSYohann 1314d537eeaSYohann if (dim<1) 132c042f62fSJeremy L Thompson // LCOV_EXCL_START 1334d537eeaSYohann return CeedError(ceed, 1, "Basis dimension must be a positive value"); 134c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 1354d537eeaSYohann 136d7b241e6Sjeremylt ierr = CeedCalloc(P*Q, &interp1d); CeedChk(ierr); 137d7b241e6Sjeremylt ierr = CeedCalloc(P*Q, &grad1d); CeedChk(ierr); 138d7b241e6Sjeremylt ierr = CeedCalloc(P, &nodes); CeedChk(ierr); 139d7b241e6Sjeremylt ierr = CeedCalloc(Q, &qref1d); CeedChk(ierr); 140d7b241e6Sjeremylt ierr = CeedCalloc(Q, &qweight1d); CeedChk(ierr); 141d7b241e6Sjeremylt // Get Nodes and Weights 142d7b241e6Sjeremylt ierr = CeedLobattoQuadrature(P, nodes, NULL); CeedChk(ierr); 143d7b241e6Sjeremylt switch (qmode) { 144d7b241e6Sjeremylt case CEED_GAUSS: 145d7b241e6Sjeremylt ierr = CeedGaussQuadrature(Q, qref1d, qweight1d); CeedChk(ierr); 146d7b241e6Sjeremylt break; 147d7b241e6Sjeremylt case CEED_GAUSS_LOBATTO: 148d7b241e6Sjeremylt ierr = CeedLobattoQuadrature(Q, qref1d, qweight1d); CeedChk(ierr); 149d7b241e6Sjeremylt break; 150d7b241e6Sjeremylt } 151d7b241e6Sjeremylt // Build B, D matrix 152d7b241e6Sjeremylt // Fornberg, 1998 153d7b241e6Sjeremylt for (i = 0; i < Q; i++) { 154d7b241e6Sjeremylt c1 = 1.0; 155d7b241e6Sjeremylt c3 = nodes[0] - qref1d[i]; 156d7b241e6Sjeremylt interp1d[i*P+0] = 1.0; 157d7b241e6Sjeremylt for (j = 1; j < P; j++) { 158d7b241e6Sjeremylt c2 = 1.0; 159d7b241e6Sjeremylt c4 = c3; 160d7b241e6Sjeremylt c3 = nodes[j] - qref1d[i]; 161d7b241e6Sjeremylt for (k = 0; k < j; k++) { 162d7b241e6Sjeremylt dx = nodes[j] - nodes[k]; 163d7b241e6Sjeremylt c2 *= dx; 164d7b241e6Sjeremylt if (k == j - 1) { 165d7b241e6Sjeremylt grad1d[i*P + j] = c1*(interp1d[i*P + k] - c4*grad1d[i*P + k]) / c2; 166d7b241e6Sjeremylt interp1d[i*P + j] = - c1*c4*interp1d[i*P + k] / c2; 167d7b241e6Sjeremylt } 168d7b241e6Sjeremylt grad1d[i*P + k] = (c3*grad1d[i*P + k] - interp1d[i*P + k]) / dx; 169d7b241e6Sjeremylt interp1d[i*P + k] = c3*interp1d[i*P + k] / dx; 170d7b241e6Sjeremylt } 171d7b241e6Sjeremylt c1 = c2; 172d7b241e6Sjeremylt } 173d7b241e6Sjeremylt } 174d7b241e6Sjeremylt // // Pass to CeedBasisCreateTensorH1 175d7b241e6Sjeremylt ierr = CeedBasisCreateTensorH1(ceed, dim, ncomp, P, Q, interp1d, grad1d, qref1d, 176d7b241e6Sjeremylt qweight1d, basis); CeedChk(ierr); 177d7b241e6Sjeremylt ierr = CeedFree(&interp1d); CeedChk(ierr); 178d7b241e6Sjeremylt ierr = CeedFree(&grad1d); CeedChk(ierr); 179d7b241e6Sjeremylt ierr = CeedFree(&nodes); CeedChk(ierr); 180d7b241e6Sjeremylt ierr = CeedFree(&qref1d); CeedChk(ierr); 181d7b241e6Sjeremylt ierr = CeedFree(&qweight1d); CeedChk(ierr); 182d7b241e6Sjeremylt return 0; 183d7b241e6Sjeremylt } 184d7b241e6Sjeremylt 185b11c1e72Sjeremylt /** 18695bb1877Svaleriabarra @brief Create a non tensor-product basis for H^1 discretizations 187a8de75f0Sjeremylt 188a8de75f0Sjeremylt @param ceed A Ceed object where the CeedBasis will be created 189a8de75f0Sjeremylt @param topo Topology of element, e.g. hypercube, simplex, ect 190a8de75f0Sjeremylt @param ncomp Number of field components (1 for scalar fields) 1918795c945Sjeremylt @param nnodes Total number of nodes 192a8de75f0Sjeremylt @param nqpts Total number of quadrature points 19395bb1877Svaleriabarra @param interp Row-major (nqpts * nnodes) matrix expressing the values of 1948795c945Sjeremylt nodal basis functions at quadrature points 19595bb1877Svaleriabarra @param grad Row-major (nqpts * dim * nnodes) matrix expressing 1968795c945Sjeremylt derivatives of nodal basis functions at quadrature points 1978795c945Sjeremylt @param qref Array of length nqpts holding the locations of quadrature 1988795c945Sjeremylt points on the reference element [-1, 1] 199a8de75f0Sjeremylt @param qweight Array of length nqpts holding the quadrature weights on the 200a8de75f0Sjeremylt reference element 201a8de75f0Sjeremylt @param[out] basis Address of the variable where the newly created 202a8de75f0Sjeremylt CeedBasis will be stored. 203a8de75f0Sjeremylt 204a8de75f0Sjeremylt @return An error code: 0 - success, otherwise - failure 205a8de75f0Sjeremylt 206a8de75f0Sjeremylt @ref Basic 207a8de75f0Sjeremylt **/ 208a8de75f0Sjeremylt int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt ncomp, 209692c2638Sjeremylt CeedInt nnodes, CeedInt nqpts, const CeedScalar *interp, 210a8de75f0Sjeremylt const CeedScalar *grad, const CeedScalar *qref, 211a8de75f0Sjeremylt const CeedScalar *qweight, CeedBasis *basis) { 212a8de75f0Sjeremylt int ierr; 2138795c945Sjeremylt CeedInt P = nnodes, Q = nqpts, dim = 0; 214a8de75f0Sjeremylt 2155fe0d4faSjeremylt if (!ceed->BasisCreateH1) { 2165fe0d4faSjeremylt Ceed delegate; 217aefd8378Sjeremylt ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr); 2185fe0d4faSjeremylt 2195fe0d4faSjeremylt if (!delegate) 220c042f62fSJeremy L Thompson // LCOV_EXCL_START 221a8de75f0Sjeremylt return CeedError(ceed, 1, "Backend does not support BasisCreateH1"); 222c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 2235fe0d4faSjeremylt 2248795c945Sjeremylt ierr = CeedBasisCreateH1(delegate, topo, ncomp, nnodes, 2255fe0d4faSjeremylt nqpts, interp, grad, qref, 2265fe0d4faSjeremylt qweight, basis); CeedChk(ierr); 2275fe0d4faSjeremylt return 0; 2285fe0d4faSjeremylt } 2295fe0d4faSjeremylt 230a8de75f0Sjeremylt ierr = CeedCalloc(1,basis); CeedChk(ierr); 231a8de75f0Sjeremylt 232a8de75f0Sjeremylt ierr = CeedBasisGetTopologyDimension(topo, &dim); CeedChk(ierr); 233a8de75f0Sjeremylt 234a8de75f0Sjeremylt (*basis)->ceed = ceed; 235a8de75f0Sjeremylt ceed->refcount++; 236a8de75f0Sjeremylt (*basis)->refcount = 1; 237a8de75f0Sjeremylt (*basis)->tensorbasis = 0; 238a8de75f0Sjeremylt (*basis)->dim = dim; 239a8de75f0Sjeremylt (*basis)->ncomp = ncomp; 240a8de75f0Sjeremylt (*basis)->P = P; 241a8de75f0Sjeremylt (*basis)->Q = Q; 242a8de75f0Sjeremylt ierr = CeedMalloc(Q*dim,&(*basis)->qref1d); CeedChk(ierr); 243a8de75f0Sjeremylt ierr = CeedMalloc(Q,&(*basis)->qweight1d); CeedChk(ierr); 244a8de75f0Sjeremylt memcpy((*basis)->qref1d, qref, Q*dim*sizeof(qref[0])); 245a8de75f0Sjeremylt memcpy((*basis)->qweight1d, qweight, Q*sizeof(qweight[0])); 24600f91b2bSjeremylt ierr = CeedMalloc(Q*P, &(*basis)->interp); CeedChk(ierr); 24700f91b2bSjeremylt ierr = CeedMalloc(dim*Q*P, &(*basis)->grad); CeedChk(ierr); 24800f91b2bSjeremylt memcpy((*basis)->interp, interp, Q*P*sizeof(interp[0])); 24900f91b2bSjeremylt memcpy((*basis)->grad, grad, dim*Q*P*sizeof(grad[0])); 250667bc5fcSjeremylt ierr = ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, qref, 251a8de75f0Sjeremylt qweight, *basis); CeedChk(ierr); 252a8de75f0Sjeremylt return 0; 253a8de75f0Sjeremylt } 254a8de75f0Sjeremylt 255a8de75f0Sjeremylt /** 256b11c1e72Sjeremylt @brief Construct a Gauss-Legendre quadrature 257b11c1e72Sjeremylt 258b11c1e72Sjeremylt @param Q Number of quadrature points (integrates polynomials of 259b11c1e72Sjeremylt degree 2*Q-1 exactly) 260b11c1e72Sjeremylt @param[out] qref1d Array of length Q to hold the abscissa on [-1, 1] 261b11c1e72Sjeremylt @param[out] qweight1d Array of length Q to hold the weights 262b11c1e72Sjeremylt 263b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 264dfdf5a53Sjeremylt 265dfdf5a53Sjeremylt @ref Utility 266b11c1e72Sjeremylt **/ 267d7b241e6Sjeremylt int CeedGaussQuadrature(CeedInt Q, CeedScalar *qref1d, CeedScalar *qweight1d) { 268d7b241e6Sjeremylt // Allocate 269d7b241e6Sjeremylt CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0*atan(1.0); 270d7b241e6Sjeremylt // Build qref1d, qweight1d 271d7b241e6Sjeremylt for (int i = 0; i <= Q/2; i++) { 272d7b241e6Sjeremylt // Guess 273d7b241e6Sjeremylt xi = cos(PI*(CeedScalar)(2*i+1)/((CeedScalar)(2*Q))); 274d7b241e6Sjeremylt // Pn(xi) 275d7b241e6Sjeremylt P0 = 1.0; 276d7b241e6Sjeremylt P1 = xi; 277d7b241e6Sjeremylt P2 = 0.0; 278d7b241e6Sjeremylt for (int j = 2; j <= Q; j++) { 279d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 280d7b241e6Sjeremylt P0 = P1; 281d7b241e6Sjeremylt P1 = P2; 282d7b241e6Sjeremylt } 283d7b241e6Sjeremylt // First Newton Step 284d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 285d7b241e6Sjeremylt xi = xi-P2/dP2; 286d7b241e6Sjeremylt // Newton to convergence 2870e4d4210Sjeremylt for (int k=0; k<100 && fabs(P2)>10*CEED_EPSILON; k++) { 288d7b241e6Sjeremylt P0 = 1.0; 289d7b241e6Sjeremylt P1 = xi; 290d7b241e6Sjeremylt for (int j = 2; j <= Q; j++) { 291d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 292d7b241e6Sjeremylt P0 = P1; 293d7b241e6Sjeremylt P1 = P2; 294d7b241e6Sjeremylt } 295d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 296d7b241e6Sjeremylt xi = xi-P2/dP2; 297d7b241e6Sjeremylt } 298d7b241e6Sjeremylt // Save xi, wi 299d7b241e6Sjeremylt wi = 2.0/((1.0-xi*xi)*dP2*dP2); 300d7b241e6Sjeremylt qweight1d[i] = wi; 301d7b241e6Sjeremylt qweight1d[Q-1-i] = wi; 302d7b241e6Sjeremylt qref1d[i] = -xi; 303d7b241e6Sjeremylt qref1d[Q-1-i]= xi; 304d7b241e6Sjeremylt } 305d7b241e6Sjeremylt return 0; 306d7b241e6Sjeremylt } 307d7b241e6Sjeremylt 308b11c1e72Sjeremylt /** 309b11c1e72Sjeremylt @brief Construct a Gauss-Legendre-Lobatto quadrature 310b11c1e72Sjeremylt 311b11c1e72Sjeremylt @param Q Number of quadrature points (integrates polynomials of 312b11c1e72Sjeremylt degree 2*Q-3 exactly) 313b11c1e72Sjeremylt @param[out] qref1d Array of length Q to hold the abscissa on [-1, 1] 314b11c1e72Sjeremylt @param[out] qweight1d Array of length Q to hold the weights 315b11c1e72Sjeremylt 316b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 317dfdf5a53Sjeremylt 318dfdf5a53Sjeremylt @ref Utility 319b11c1e72Sjeremylt **/ 320d7b241e6Sjeremylt int CeedLobattoQuadrature(CeedInt Q, CeedScalar *qref1d, 321d7b241e6Sjeremylt CeedScalar *qweight1d) { 322d7b241e6Sjeremylt // Allocate 323d7b241e6Sjeremylt CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0*atan(1.0); 324d7b241e6Sjeremylt // Build qref1d, qweight1d 325d7b241e6Sjeremylt // Set endpoints 326d7b241e6Sjeremylt wi = 2.0/((CeedScalar)(Q*(Q-1))); 327d7b241e6Sjeremylt if (qweight1d) { 328d7b241e6Sjeremylt qweight1d[0] = wi; 329d7b241e6Sjeremylt qweight1d[Q-1] = wi; 330d7b241e6Sjeremylt } 331d7b241e6Sjeremylt qref1d[0] = -1.0; 332d7b241e6Sjeremylt qref1d[Q-1] = 1.0; 333d7b241e6Sjeremylt // Interior 334d7b241e6Sjeremylt for (int i = 1; i <= (Q-1)/2; i++) { 335d7b241e6Sjeremylt // Guess 336d7b241e6Sjeremylt xi = cos(PI*(CeedScalar)(i)/(CeedScalar)(Q-1)); 337d7b241e6Sjeremylt // Pn(xi) 338d7b241e6Sjeremylt P0 = 1.0; 339d7b241e6Sjeremylt P1 = xi; 340d7b241e6Sjeremylt P2 = 0.0; 341d7b241e6Sjeremylt for (int j = 2; j < Q; j++) { 342d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 343d7b241e6Sjeremylt P0 = P1; 344d7b241e6Sjeremylt P1 = P2; 345d7b241e6Sjeremylt } 346d7b241e6Sjeremylt // First Newton step 347d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 348d7b241e6Sjeremylt d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi); 349d7b241e6Sjeremylt xi = xi-dP2/d2P2; 350d7b241e6Sjeremylt // Newton to convergence 3510e4d4210Sjeremylt for (int k=0; k<100 && fabs(dP2)>10*CEED_EPSILON; k++) { 352d7b241e6Sjeremylt P0 = 1.0; 353d7b241e6Sjeremylt P1 = xi; 354d7b241e6Sjeremylt for (int j = 2; j < Q; j++) { 355d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 356d7b241e6Sjeremylt P0 = P1; 357d7b241e6Sjeremylt P1 = P2; 358d7b241e6Sjeremylt } 359d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 360d7b241e6Sjeremylt d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi); 361d7b241e6Sjeremylt xi = xi-dP2/d2P2; 362d7b241e6Sjeremylt } 363d7b241e6Sjeremylt // Save xi, wi 364d7b241e6Sjeremylt wi = 2.0/(((CeedScalar)(Q*(Q-1)))*P2*P2); 365d7b241e6Sjeremylt if (qweight1d) { 366d7b241e6Sjeremylt qweight1d[i] = wi; 367d7b241e6Sjeremylt qweight1d[Q-1-i] = wi; 368d7b241e6Sjeremylt } 369d7b241e6Sjeremylt qref1d[i] = -xi; 370d7b241e6Sjeremylt qref1d[Q-1-i]= xi; 371d7b241e6Sjeremylt } 372d7b241e6Sjeremylt return 0; 373d7b241e6Sjeremylt } 374d7b241e6Sjeremylt 375dfdf5a53Sjeremylt /** 376dfdf5a53Sjeremylt @brief View an array stored in a CeedBasis 377dfdf5a53Sjeremylt 378dfdf5a53Sjeremylt @param name Name of array 379dfdf5a53Sjeremylt @param fpformat Printing format 380dfdf5a53Sjeremylt @param m Number of rows in array 381dfdf5a53Sjeremylt @param n Number of columns in array 382dfdf5a53Sjeremylt @param a Array to be viewed 383dfdf5a53Sjeremylt @param stream Stream to view to, e.g., stdout 384dfdf5a53Sjeremylt 385dfdf5a53Sjeremylt @return An error code: 0 - success, otherwise - failure 386dfdf5a53Sjeremylt 387dfdf5a53Sjeremylt @ref Utility 388dfdf5a53Sjeremylt **/ 389d7b241e6Sjeremylt static int CeedScalarView(const char *name, const char *fpformat, CeedInt m, 390d7b241e6Sjeremylt CeedInt n, const CeedScalar *a, FILE *stream) { 391d7b241e6Sjeremylt for (int i=0; i<m; i++) { 3921d102b48SJeremy L Thompson if (m > 1) 3931d102b48SJeremy L Thompson fprintf(stream, "%12s[%d]:", name, i); 3941d102b48SJeremy L Thompson else 3951d102b48SJeremy L Thompson fprintf(stream, "%12s:", name); 3961d102b48SJeremy L Thompson for (int j=0; j<n; j++) 397d7b241e6Sjeremylt fprintf(stream, fpformat, fabs(a[i*n+j]) > 1E-14 ? a[i*n+j] : 0); 398d7b241e6Sjeremylt fputs("\n", stream); 399d7b241e6Sjeremylt } 400d7b241e6Sjeremylt return 0; 401d7b241e6Sjeremylt } 402d7b241e6Sjeremylt 403b11c1e72Sjeremylt /** 404b11c1e72Sjeremylt @brief View a CeedBasis 405b11c1e72Sjeremylt 406b11c1e72Sjeremylt @param basis CeedBasis to view 407b11c1e72Sjeremylt @param stream Stream to view to, e.g., stdout 408b11c1e72Sjeremylt 409b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 410dfdf5a53Sjeremylt 411dfdf5a53Sjeremylt @ref Utility 412b11c1e72Sjeremylt **/ 413d7b241e6Sjeremylt int CeedBasisView(CeedBasis basis, FILE *stream) { 414d7b241e6Sjeremylt int ierr; 415d7b241e6Sjeremylt 416a8de75f0Sjeremylt if (basis->tensorbasis) { 417d7b241e6Sjeremylt fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P1d, 418d7b241e6Sjeremylt basis->Q1d); 419d7b241e6Sjeremylt ierr = CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q1d, basis->qref1d, 420d7b241e6Sjeremylt stream); CeedChk(ierr); 4218795c945Sjeremylt ierr = CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q1d, 4228795c945Sjeremylt basis->qweight1d, stream); CeedChk(ierr); 423d7b241e6Sjeremylt ierr = CeedScalarView("interp1d", "\t% 12.8f", basis->Q1d, basis->P1d, 424d7b241e6Sjeremylt basis->interp1d, stream); CeedChk(ierr); 425d7b241e6Sjeremylt ierr = CeedScalarView("grad1d", "\t% 12.8f", basis->Q1d, basis->P1d, 426d7b241e6Sjeremylt basis->grad1d, stream); CeedChk(ierr); 427a8de75f0Sjeremylt } else { 428a8de75f0Sjeremylt fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P, 429a8de75f0Sjeremylt basis->Q); 430a8de75f0Sjeremylt ierr = CeedScalarView("qref", "\t% 12.8f", 1, basis->Q*basis->dim, 431a8de75f0Sjeremylt basis->qref1d, 432a8de75f0Sjeremylt stream); CeedChk(ierr); 433a8de75f0Sjeremylt ierr = CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->qweight1d, 434a8de75f0Sjeremylt stream); CeedChk(ierr); 435a8de75f0Sjeremylt ierr = CeedScalarView("interp", "\t% 12.8f", basis->Q, basis->P, 43600f91b2bSjeremylt basis->interp, stream); CeedChk(ierr); 437a8de75f0Sjeremylt ierr = CeedScalarView("grad", "\t% 12.8f", basis->dim*basis->Q, basis->P, 43800f91b2bSjeremylt basis->grad, stream); CeedChk(ierr); 439a8de75f0Sjeremylt } 440d7b241e6Sjeremylt return 0; 441d7b241e6Sjeremylt } 442d7b241e6Sjeremylt 443dfdf5a53Sjeremylt /** 44452bfb9bbSJeremy L Thompson @brief Compute Householder reflection 445dfdf5a53Sjeremylt 446dfdf5a53Sjeremylt Computes A = (I - b v v^T) A 447dfdf5a53Sjeremylt where A is an mxn matrix indexed as A[i*row + j*col] 448dfdf5a53Sjeremylt 44952bfb9bbSJeremy L Thompson @param[in,out] A Matrix to apply Householder reflection to, in place 450dfdf5a53Sjeremylt @param v Householder vector 451dfdf5a53Sjeremylt @param b Scaling factor 452dfdf5a53Sjeremylt @param m Number of rows in A 453dfdf5a53Sjeremylt @param n Number of columns in A 45452bfb9bbSJeremy L Thompson @param row Row stride 45552bfb9bbSJeremy L Thompson @param col Col stride 456dfdf5a53Sjeremylt 457dfdf5a53Sjeremylt @return An error code: 0 - success, otherwise - failure 458dfdf5a53Sjeremylt 459dfdf5a53Sjeremylt @ref Developer 460dfdf5a53Sjeremylt **/ 461d7b241e6Sjeremylt static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v, 462d7b241e6Sjeremylt CeedScalar b, CeedInt m, CeedInt n, 463d7b241e6Sjeremylt CeedInt row, CeedInt col) { 464d7b241e6Sjeremylt for (CeedInt j=0; j<n; j++) { 465d7b241e6Sjeremylt CeedScalar w = A[0*row + j*col]; 4661d102b48SJeremy L Thompson for (CeedInt i=1; i<m; i++) 4671d102b48SJeremy L Thompson w += v[i] * A[i*row + j*col]; 468d7b241e6Sjeremylt A[0*row + j*col] -= b * w; 4691d102b48SJeremy L Thompson for (CeedInt i=1; i<m; i++) 4701d102b48SJeremy L Thompson A[i*row + j*col] -= b * w * v[i]; 471d7b241e6Sjeremylt } 472d7b241e6Sjeremylt return 0; 473d7b241e6Sjeremylt } 474d7b241e6Sjeremylt 475dfdf5a53Sjeremylt /** 476dfdf5a53Sjeremylt @brief Apply Householder Q matrix 477dfdf5a53Sjeremylt 47852bfb9bbSJeremy L Thompson Compute A = Q A where Q is mxm and A is mxn. 479dfdf5a53Sjeremylt 48052bfb9bbSJeremy L Thompson @param[in,out] A Matrix to apply Householder Q to, in place 481dfdf5a53Sjeremylt @param Q Householder Q matrix 482dfdf5a53Sjeremylt @param tau Householder scaling factors 483dfdf5a53Sjeremylt @param tmode Transpose mode for application 484dfdf5a53Sjeremylt @param m Number of rows in A 485dfdf5a53Sjeremylt @param n Number of columns in A 48652bfb9bbSJeremy L Thompson @param k Number of elementary reflectors in Q, k<m 48752bfb9bbSJeremy L Thompson @param row Row stride in A 48852bfb9bbSJeremy L Thompson @param col Col stride in A 489dfdf5a53Sjeremylt 490dfdf5a53Sjeremylt @return An error code: 0 - success, otherwise - failure 491dfdf5a53Sjeremylt 492dfdf5a53Sjeremylt @ref Developer 493dfdf5a53Sjeremylt **/ 494d7b241e6Sjeremylt static int CeedHouseholderApplyQ(CeedScalar *A, const CeedScalar *Q, 495d7b241e6Sjeremylt const CeedScalar *tau, CeedTransposeMode tmode, 496d7b241e6Sjeremylt CeedInt m, CeedInt n, CeedInt k, 497d7b241e6Sjeremylt CeedInt row, CeedInt col) { 498d7b241e6Sjeremylt CeedScalar v[m]; 499d7b241e6Sjeremylt for (CeedInt ii=0; ii<k; ii++) { 500d7b241e6Sjeremylt CeedInt i = tmode == CEED_TRANSPOSE ? ii : k-1-ii; 50152bfb9bbSJeremy L Thompson for (CeedInt j=i+1; j<m; j++) 502d7b241e6Sjeremylt v[j] = Q[j*k+i]; 503dc7d240cSValeria Barra // Apply Householder reflector (I - tau v v^T) collograd1d^T 504d7b241e6Sjeremylt CeedHouseholderReflect(&A[i*row], &v[i], tau[i], m-i, n, row, col); 505d7b241e6Sjeremylt } 506d7b241e6Sjeremylt return 0; 507d7b241e6Sjeremylt } 508d7b241e6Sjeremylt 509b11c1e72Sjeremylt /** 51052bfb9bbSJeremy L Thompson @brief Compute Givens rotation 51152bfb9bbSJeremy L Thompson 51252bfb9bbSJeremy L Thompson Computes A = G A (or G^T A in transpose mode) 51352bfb9bbSJeremy L Thompson where A is an mxn matrix indexed as A[i*n + j*m] 51452bfb9bbSJeremy L Thompson 51552bfb9bbSJeremy L Thompson @param[in,out] A Row major matrix to apply Givens rotation to, in place 51652bfb9bbSJeremy L Thompson @param c Cosine factor 51752bfb9bbSJeremy L Thompson @param s Sine factor 51852bfb9bbSJeremy L Thompson @param i First row/column to apply rotation 51952bfb9bbSJeremy L Thompson @param k Second row/column to apply rotation 52052bfb9bbSJeremy L Thompson @param m Number of rows in A 52152bfb9bbSJeremy L Thompson @param n Number of columns in A 52252bfb9bbSJeremy L Thompson 52352bfb9bbSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 52452bfb9bbSJeremy L Thompson 52552bfb9bbSJeremy L Thompson @ref Developer 52652bfb9bbSJeremy L Thompson **/ 52752bfb9bbSJeremy L Thompson static int CeedGivensRotation(CeedScalar *A, CeedScalar c, CeedScalar s, 52852bfb9bbSJeremy L Thompson CeedTransposeMode tmode, CeedInt i, CeedInt k, 52952bfb9bbSJeremy L Thompson CeedInt m, CeedInt n) { 53052bfb9bbSJeremy L Thompson CeedInt stridej = 1, strideik = m, numits = n; 53152bfb9bbSJeremy L Thompson if (tmode == CEED_NOTRANSPOSE) { 53252bfb9bbSJeremy L Thompson stridej = n; strideik = 1; numits = m; 53352bfb9bbSJeremy L Thompson } 53452bfb9bbSJeremy L Thompson 53552bfb9bbSJeremy L Thompson // Apply rotation 53652bfb9bbSJeremy L Thompson for (CeedInt j=0; j<numits; j++) { 53752bfb9bbSJeremy L Thompson CeedScalar tau1 = A[i*strideik+j*stridej], tau2 = A[k*strideik+j*stridej]; 53852bfb9bbSJeremy L Thompson A[i*strideik+j*stridej] = c*tau1 - s*tau2; 53952bfb9bbSJeremy L Thompson A[k*strideik+j*stridej] = s*tau1 + c*tau2; 54052bfb9bbSJeremy L Thompson } 54152bfb9bbSJeremy L Thompson 54252bfb9bbSJeremy L Thompson return 0; 54352bfb9bbSJeremy L Thompson } 54452bfb9bbSJeremy L Thompson 54552bfb9bbSJeremy L Thompson /** 54695bb1877Svaleriabarra @brief Return QR Factorization of a matrix 547b11c1e72Sjeremylt 548*77645d7bSjeremylt @param ceed A Ceed context for error handling 54952bfb9bbSJeremy L Thompson @param[in,out] mat Row-major matrix to be factorized in place 55052bfb9bbSJeremy L Thompson @param[in,out] tau Vector of length m of scaling factors 551b11c1e72Sjeremylt @param m Number of rows 552b11c1e72Sjeremylt @param n Number of columns 553b11c1e72Sjeremylt 554b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 555dfdf5a53Sjeremylt 556dfdf5a53Sjeremylt @ref Utility 557b11c1e72Sjeremylt **/ 558a7bd39daSjeremylt int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau, 559d7b241e6Sjeremylt CeedInt m, CeedInt n) { 560d7b241e6Sjeremylt CeedScalar v[m]; 561d7b241e6Sjeremylt 562a7bd39daSjeremylt // Check m >= n 563a7bd39daSjeremylt if (n > m) 564c042f62fSJeremy L Thompson // LCOV_EXCL_START 565a7bd39daSjeremylt return CeedError(ceed, 1, "Cannot compute QR factorization with n > m"); 566c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 567a7bd39daSjeremylt 56852bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) { 569d7b241e6Sjeremylt // Calculate Householder vector, magnitude 570d7b241e6Sjeremylt CeedScalar sigma = 0.0; 571d7b241e6Sjeremylt v[i] = mat[i+n*i]; 57252bfb9bbSJeremy L Thompson for (CeedInt j=i+1; j<m; j++) { 573d7b241e6Sjeremylt v[j] = mat[i+n*j]; 574d7b241e6Sjeremylt sigma += v[j] * v[j]; 575d7b241e6Sjeremylt } 576d7b241e6Sjeremylt CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m] 577d7b241e6Sjeremylt CeedScalar Rii = -copysign(norm, v[i]); 578d7b241e6Sjeremylt v[i] -= Rii; 579d7b241e6Sjeremylt // norm of v[i:m] after modification above and scaling below 580d7b241e6Sjeremylt // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 581d7b241e6Sjeremylt // tau = 2 / (norm*norm) 582d7b241e6Sjeremylt tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma); 583fb551037Sjeremylt 5841d102b48SJeremy L Thompson for (CeedInt j=i+1; j<m; j++) 5851d102b48SJeremy L Thompson v[j] /= v[i]; 586d7b241e6Sjeremylt 587d7b241e6Sjeremylt // Apply Householder reflector to lower right panel 588d7b241e6Sjeremylt CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1); 589d7b241e6Sjeremylt // Save v 590d7b241e6Sjeremylt mat[i+n*i] = Rii; 5911d102b48SJeremy L Thompson for (CeedInt j=i+1; j<m; j++) 592d7b241e6Sjeremylt mat[i+n*j] = v[j]; 593d7b241e6Sjeremylt } 594d7b241e6Sjeremylt 595d7b241e6Sjeremylt return 0; 596d7b241e6Sjeremylt } 597d7b241e6Sjeremylt 598b11c1e72Sjeremylt /** 59952bfb9bbSJeremy L Thompson @brief Return symmetric Schur decomposition of the symmetric matrix mat via 60052bfb9bbSJeremy L Thompson symmetric QR factorization 60152bfb9bbSJeremy L Thompson 602*77645d7bSjeremylt @param ceed A Ceed context for error handling 60352bfb9bbSJeremy L Thompson @param[in,out] mat Row-major matrix to be factorized in place 604460bf743SValeria Barra @param[out] lambda Vector of length n of eigenvalues 60552bfb9bbSJeremy L Thompson @param n Number of rows/columns 60652bfb9bbSJeremy L Thompson 60752bfb9bbSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 60852bfb9bbSJeremy L Thompson 60952bfb9bbSJeremy L Thompson @ref Utility 61052bfb9bbSJeremy L Thompson **/ 61152bfb9bbSJeremy L Thompson int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat, 61252bfb9bbSJeremy L Thompson CeedScalar *lambda, CeedInt n) { 61352bfb9bbSJeremy L Thompson // Check bounds for clang-tidy 61452bfb9bbSJeremy L Thompson if (n<2) 615c042f62fSJeremy L Thompson // LCOV_EXCL_START 616c042f62fSJeremy L Thompson return CeedError(ceed, 1, 617c042f62fSJeremy L Thompson "Cannot compute symmetric Schur decomposition of scalars"); 618c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 61952bfb9bbSJeremy L Thompson 62052bfb9bbSJeremy L Thompson CeedScalar v[n-1], tau[n-1], matT[n*n]; 62152bfb9bbSJeremy L Thompson 62252bfb9bbSJeremy L Thompson // Copy mat to matT and set mat to I 62352bfb9bbSJeremy L Thompson memcpy(matT, mat, n*n*sizeof(mat[0])); 62452bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) 62552bfb9bbSJeremy L Thompson for (CeedInt j=0; j<n; j++) 62652bfb9bbSJeremy L Thompson mat[j+n*i] = (i==j) ? 1 : 0; 62752bfb9bbSJeremy L Thompson 62852bfb9bbSJeremy L Thompson // Reduce to tridiagonal 62952bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n-1; i++) { 63052bfb9bbSJeremy L Thompson // Calculate Householder vector, magnitude 63152bfb9bbSJeremy L Thompson CeedScalar sigma = 0.0; 63252bfb9bbSJeremy L Thompson v[i] = matT[i+n*(i+1)]; 63352bfb9bbSJeremy L Thompson for (CeedInt j=i+1; j<n-1; j++) { 63452bfb9bbSJeremy L Thompson v[j] = matT[i+n*(j+1)]; 63552bfb9bbSJeremy L Thompson sigma += v[j] * v[j]; 63652bfb9bbSJeremy L Thompson } 63752bfb9bbSJeremy L Thompson CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:n-1] 63852bfb9bbSJeremy L Thompson CeedScalar Rii = -copysign(norm, v[i]); 63952bfb9bbSJeremy L Thompson v[i] -= Rii; 64052bfb9bbSJeremy L Thompson // norm of v[i:m] after modification above and scaling below 64152bfb9bbSJeremy L Thompson // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 64252bfb9bbSJeremy L Thompson // tau = 2 / (norm*norm) 6430e4d4210Sjeremylt if (sigma > 10*CEED_EPSILON) 64452bfb9bbSJeremy L Thompson tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma); 645fb551037Sjeremylt else 646fb551037Sjeremylt tau[i] = 0; 647fb551037Sjeremylt 648fb551037Sjeremylt for (CeedInt j=i+1; j<n-1; j++) 649fb551037Sjeremylt v[j] /= v[i]; 65052bfb9bbSJeremy L Thompson 65152bfb9bbSJeremy L Thompson // Update sub and super diagonal 65252bfb9bbSJeremy L Thompson matT[i+n*(i+1)] = Rii; 65352bfb9bbSJeremy L Thompson matT[(i+1)+n*i] = Rii; 65452bfb9bbSJeremy L Thompson for (CeedInt j=i+2; j<n; j++) { 65552bfb9bbSJeremy L Thompson matT[i+n*j] = 0; matT[j+n*i] = 0; 65652bfb9bbSJeremy L Thompson } 65752bfb9bbSJeremy L Thompson // Apply symmetric Householder reflector to lower right panel 65852bfb9bbSJeremy L Thompson CeedHouseholderReflect(&matT[(i+1)+n*(i+1)], &v[i], tau[i], 65952bfb9bbSJeremy L Thompson n-(i+1), n-(i+1), n, 1); 66052bfb9bbSJeremy L Thompson CeedHouseholderReflect(&matT[(i+1)+n*(i+1)], &v[i], tau[i], 66152bfb9bbSJeremy L Thompson n-(i+1), n-(i+1), 1, n); 66252bfb9bbSJeremy L Thompson // Save v 66352bfb9bbSJeremy L Thompson for (CeedInt j=i+1; j<n-1; j++) { 66452bfb9bbSJeremy L Thompson matT[i+n*(j+1)] = v[j]; 66552bfb9bbSJeremy L Thompson } 66652bfb9bbSJeremy L Thompson } 66752bfb9bbSJeremy L Thompson // Backwards accumulation of Q 66852bfb9bbSJeremy L Thompson for (CeedInt i=n-2; i>=0; i--) { 66952bfb9bbSJeremy L Thompson v[i] = 1; 67052bfb9bbSJeremy L Thompson for (CeedInt j=i+1; j<n-1; j++) { 67152bfb9bbSJeremy L Thompson v[j] = matT[i+n*(j+1)]; 67252bfb9bbSJeremy L Thompson matT[i+n*(j+1)] = 0; 67352bfb9bbSJeremy L Thompson } 67452bfb9bbSJeremy L Thompson CeedHouseholderReflect(&mat[(i+1)+n*(i+1)], &v[i], tau[i], 67552bfb9bbSJeremy L Thompson n-(i+1), n-(i+1), n, 1); 67652bfb9bbSJeremy L Thompson } 67752bfb9bbSJeremy L Thompson 67852bfb9bbSJeremy L Thompson // Reduce sub and super diagonal 67952bfb9bbSJeremy L Thompson CeedInt p = 0, q = 0, itr = 0, maxitr = n*n*n; 6800e4d4210Sjeremylt CeedScalar tol = 10*CEED_EPSILON; 68152bfb9bbSJeremy L Thompson 68252bfb9bbSJeremy L Thompson while (q < n && itr < maxitr) { 68352bfb9bbSJeremy L Thompson // Update p, q, size of reduced portions of diagonal 68452bfb9bbSJeremy L Thompson p = 0; q = 0; 68552bfb9bbSJeremy L Thompson for (CeedInt i=n-2; i>=0; i--) { 68652bfb9bbSJeremy L Thompson if (fabs(matT[i+n*(i+1)]) < tol) 68752bfb9bbSJeremy L Thompson q += 1; 68852bfb9bbSJeremy L Thompson else 68952bfb9bbSJeremy L Thompson break; 69052bfb9bbSJeremy L Thompson } 69152bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n-1-q; i++) { 69252bfb9bbSJeremy L Thompson if (fabs(matT[i+n*(i+1)]) < tol) 69352bfb9bbSJeremy L Thompson p += 1; 69452bfb9bbSJeremy L Thompson else 69552bfb9bbSJeremy L Thompson break; 69652bfb9bbSJeremy L Thompson } 69752bfb9bbSJeremy L Thompson if (q == n-1) break; // Finished reducing 69852bfb9bbSJeremy L Thompson 69952bfb9bbSJeremy L Thompson // Reduce tridiagonal portion 70052bfb9bbSJeremy L Thompson CeedScalar tnn = matT[(n-1-q)+n*(n-1-q)], 70152bfb9bbSJeremy L Thompson tnnm1 = matT[(n-2-q)+n*(n-1-q)]; 70252bfb9bbSJeremy L Thompson CeedScalar d = (matT[(n-2-q)+n*(n-2-q)] - tnn)/2; 70352bfb9bbSJeremy L Thompson CeedScalar mu = tnn - tnnm1*tnnm1 / 70452bfb9bbSJeremy L Thompson (d + copysign(sqrt(d*d + tnnm1*tnnm1), d)); 70552bfb9bbSJeremy L Thompson CeedScalar x = matT[p+n*p] - mu; 70652bfb9bbSJeremy L Thompson CeedScalar z = matT[p+n*(p+1)]; 70752bfb9bbSJeremy L Thompson for (CeedInt k=p; k<n-1-q; k++) { 70852bfb9bbSJeremy L Thompson // Compute Givens rotation 70952bfb9bbSJeremy L Thompson CeedScalar c = 1, s = 0; 71052bfb9bbSJeremy L Thompson if (fabs(z) > tol) { 71152bfb9bbSJeremy L Thompson if (fabs(z) > fabs(x)) { 71252bfb9bbSJeremy L Thompson CeedScalar tau = -x/z; 71352bfb9bbSJeremy L Thompson s = 1/sqrt(1+tau*tau), c = s*tau; 71452bfb9bbSJeremy L Thompson } else { 71552bfb9bbSJeremy L Thompson CeedScalar tau = -z/x; 71652bfb9bbSJeremy L Thompson c = 1/sqrt(1+tau*tau), s = c*tau; 71752bfb9bbSJeremy L Thompson } 71852bfb9bbSJeremy L Thompson } 71952bfb9bbSJeremy L Thompson 72052bfb9bbSJeremy L Thompson // Apply Givens rotation to T 72152bfb9bbSJeremy L Thompson CeedGivensRotation(matT, c, s, CEED_NOTRANSPOSE, k, k+1, n, n); 72252bfb9bbSJeremy L Thompson CeedGivensRotation(matT, c, s, CEED_TRANSPOSE, k, k+1, n, n); 72352bfb9bbSJeremy L Thompson 72452bfb9bbSJeremy L Thompson // Apply Givens rotation to Q 72552bfb9bbSJeremy L Thompson CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k+1, n, n); 72652bfb9bbSJeremy L Thompson 72752bfb9bbSJeremy L Thompson // Update x, z 72852bfb9bbSJeremy L Thompson if (k < n-q-2) { 72952bfb9bbSJeremy L Thompson x = matT[k+n*(k+1)]; 73052bfb9bbSJeremy L Thompson z = matT[k+n*(k+2)]; 73152bfb9bbSJeremy L Thompson } 73252bfb9bbSJeremy L Thompson } 73352bfb9bbSJeremy L Thompson itr++; 73452bfb9bbSJeremy L Thompson } 73552bfb9bbSJeremy L Thompson // Save eigenvalues 73652bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) 73752bfb9bbSJeremy L Thompson lambda[i] = matT[i+n*i]; 73852bfb9bbSJeremy L Thompson 73952bfb9bbSJeremy L Thompson // Check convergence 74052bfb9bbSJeremy L Thompson if (itr == maxitr && q < n-1) 741c042f62fSJeremy L Thompson // LCOV_EXCL_START 74252bfb9bbSJeremy L Thompson return CeedError(ceed, 1, "Symmetric QR failed to converge"); 743c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 74452bfb9bbSJeremy L Thompson 74552bfb9bbSJeremy L Thompson return 0; 74652bfb9bbSJeremy L Thompson } 74752bfb9bbSJeremy L Thompson 74852bfb9bbSJeremy L Thompson /** 7499289e5bfSjeremylt @brief Return a reference implementation of matrix multiplication C = A B. 7509289e5bfSjeremylt Note, this is a reference implementation for CPU CeedScalar pointers 7519289e5bfSjeremylt that is not intended for high performance. 75252bfb9bbSJeremy L Thompson 753*77645d7bSjeremylt @param ceed A Ceed context for error handling 75452bfb9bbSJeremy L Thompson @param[in] matA Row-major matrix A 75552bfb9bbSJeremy L Thompson @param[in] matB Row-major matrix B 75652bfb9bbSJeremy L Thompson @param[out] matC Row-major output matrix C 75752bfb9bbSJeremy L Thompson @param m Number of rows of C 75852bfb9bbSJeremy L Thompson @param n Number of columns of C 75952bfb9bbSJeremy L Thompson @param kk Number of columns of A/rows of B 76052bfb9bbSJeremy L Thompson 76152bfb9bbSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 76252bfb9bbSJeremy L Thompson 76352bfb9bbSJeremy L Thompson @ref Utility 76452bfb9bbSJeremy L Thompson **/ 7659289e5bfSjeremylt int CeedMatrixMultiply(Ceed ceed, const CeedScalar *matA, 7669289e5bfSjeremylt const CeedScalar *matB, CeedScalar *matC, CeedInt m, 7679289e5bfSjeremylt CeedInt n, CeedInt kk) { 76852bfb9bbSJeremy L Thompson for (CeedInt i=0; i<m; i++) 76952bfb9bbSJeremy L Thompson for (CeedInt j=0; j<n; j++) { 77052bfb9bbSJeremy L Thompson CeedScalar sum = 0; 77152bfb9bbSJeremy L Thompson for (CeedInt k=0; k<kk; k++) 77252bfb9bbSJeremy L Thompson sum += matA[k+i*kk]*matB[j+k*n]; 77352bfb9bbSJeremy L Thompson matC[j+i*n] = sum; 77452bfb9bbSJeremy L Thompson } 77552bfb9bbSJeremy L Thompson return 0; 77652bfb9bbSJeremy L Thompson } 77752bfb9bbSJeremy L Thompson 77852bfb9bbSJeremy L Thompson /** 77952bfb9bbSJeremy L Thompson @brief Return Simultaneous Diagonalization of two matrices. This solves the 78052bfb9bbSJeremy L Thompson generalized eigenvalue problem A x = lambda B x, where A and B 78152bfb9bbSJeremy L Thompson are symmetric and B is positive definite. We generate the matrix X 78252bfb9bbSJeremy L Thompson and vector Lambda such that X^T A X = Lambda and X^T B X = I. This 78352bfb9bbSJeremy L Thompson is equivalent to the LAPACK routine 'sygv' with TYPE = 1. 78452bfb9bbSJeremy L Thompson 785*77645d7bSjeremylt @param ceed A Ceed context for error handling 78652bfb9bbSJeremy L Thompson @param[in] matA Row-major matrix to be factorized with eigenvalues 78752bfb9bbSJeremy L Thompson @param[in] matB Row-major matrix to be factorized to identity 78852bfb9bbSJeremy L Thompson @param[out] x Row-major orthogonal matrix 789460bf743SValeria Barra @param[out] lambda Vector of length n of generalized eigenvalues 79052bfb9bbSJeremy L Thompson @param n Number of rows/columns 79152bfb9bbSJeremy L Thompson 79252bfb9bbSJeremy L Thompson @return An error code: 0 - success, otherwise - failure 79352bfb9bbSJeremy L Thompson 79452bfb9bbSJeremy L Thompson @ref Utility 79552bfb9bbSJeremy L Thompson **/ 79652bfb9bbSJeremy L Thompson int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *matA, 79752bfb9bbSJeremy L Thompson CeedScalar *matB, CeedScalar *x, 79852bfb9bbSJeremy L Thompson CeedScalar *lambda, CeedInt n) { 79952bfb9bbSJeremy L Thompson int ierr; 80052bfb9bbSJeremy L Thompson CeedScalar matC[n*n], matG[n*n], vecD[n]; 80152bfb9bbSJeremy L Thompson 80252bfb9bbSJeremy L Thompson // Compute B = G D G^T 80352bfb9bbSJeremy L Thompson memcpy(matG, matB, n*n*sizeof(matB[0])); 80452bfb9bbSJeremy L Thompson ierr = CeedSymmetricSchurDecomposition(ceed, matG, vecD, n); CeedChk(ierr); 805fb551037Sjeremylt for (CeedInt i=0; i<n; i++) 806fb551037Sjeremylt vecD[i] = sqrt(vecD[i]); 80752bfb9bbSJeremy L Thompson 808fb551037Sjeremylt // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T 809fb551037Sjeremylt // = D^-1/2 G^T A G D^-1/2 81052bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) 81152bfb9bbSJeremy L Thompson for (CeedInt j=0; j<n; j++) 812fb551037Sjeremylt matC[j+i*n] = matG[i+j*n] / vecD[i]; 8139289e5bfSjeremylt ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)matC, 8149289e5bfSjeremylt (const CeedScalar *)matA, x, n, n, n); 8159289e5bfSjeremylt CeedChk(ierr); 81652bfb9bbSJeremy L Thompson for (CeedInt i=0; i<n; i++) 81752bfb9bbSJeremy L Thompson for (CeedInt j=0; j<n; j++) 818fb551037Sjeremylt matG[j+i*n] = matG[j+i*n] / vecD[j]; 8199289e5bfSjeremylt ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)x, 8209289e5bfSjeremylt (const CeedScalar *)matG, matC, n, n, n); 8219289e5bfSjeremylt CeedChk(ierr); 82252bfb9bbSJeremy L Thompson 82352bfb9bbSJeremy L Thompson // Compute Q^T C Q = lambda 82452bfb9bbSJeremy L Thompson ierr = CeedSymmetricSchurDecomposition(ceed, matC, lambda, n); CeedChk(ierr); 82552bfb9bbSJeremy L Thompson 826fb551037Sjeremylt // Set x = (G D^1/2)^-T Q 827fb551037Sjeremylt // = G D^-1/2 Q 8289289e5bfSjeremylt ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)matG, 8299289e5bfSjeremylt (const CeedScalar *)matC, x, n, n, n); 8309289e5bfSjeremylt CeedChk(ierr); 83152bfb9bbSJeremy L Thompson 83252bfb9bbSJeremy L Thompson return 0; 83352bfb9bbSJeremy L Thompson } 83452bfb9bbSJeremy L Thompson 83552bfb9bbSJeremy L Thompson /** 836783c99b3SValeria Barra @brief Return collocated grad matrix 837b11c1e72Sjeremylt 838b11c1e72Sjeremylt @param basis CeedBasis 83995bb1877Svaleriabarra @param[out] collograd1d Row-major (Q1d * Q1d) matrix expressing derivatives of 840b11c1e72Sjeremylt basis functions at quadrature points 841b11c1e72Sjeremylt 842b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 843dfdf5a53Sjeremylt 844dfdf5a53Sjeremylt @ref Advanced 845b11c1e72Sjeremylt **/ 846dc7d240cSValeria Barra int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *collograd1d) { 847d7b241e6Sjeremylt int i, j, k; 848a7bd39daSjeremylt Ceed ceed; 849d7b241e6Sjeremylt CeedInt ierr, P1d=(basis)->P1d, Q1d=(basis)->Q1d; 850d7b241e6Sjeremylt CeedScalar *interp1d, *grad1d, tau[Q1d]; 851d7b241e6Sjeremylt 852d7b241e6Sjeremylt ierr = CeedMalloc(Q1d*P1d, &interp1d); CeedChk(ierr); 853d7b241e6Sjeremylt ierr = CeedMalloc(Q1d*P1d, &grad1d); CeedChk(ierr); 854d7b241e6Sjeremylt memcpy(interp1d, (basis)->interp1d, Q1d*P1d*sizeof(basis)->interp1d[0]); 855d7b241e6Sjeremylt memcpy(grad1d, (basis)->grad1d, Q1d*P1d*sizeof(basis)->interp1d[0]); 856d7b241e6Sjeremylt 857d7b241e6Sjeremylt // QR Factorization, interp1d = Q R 858a7bd39daSjeremylt ierr = CeedBasisGetCeed(basis, &ceed); CeedChk(ierr); 859a7bd39daSjeremylt ierr = CeedQRFactorization(ceed, interp1d, tau, Q1d, P1d); CeedChk(ierr); 860d7b241e6Sjeremylt 861dc7d240cSValeria Barra // Apply Rinv, collograd1d = grad1d Rinv 862d7b241e6Sjeremylt for (i=0; i<Q1d; i++) { // Row i 863dc7d240cSValeria Barra collograd1d[Q1d*i] = grad1d[P1d*i]/interp1d[0]; 864d7b241e6Sjeremylt for (j=1; j<P1d; j++) { // Column j 865dc7d240cSValeria Barra collograd1d[j+Q1d*i] = grad1d[j+P1d*i]; 8661d102b48SJeremy L Thompson for (k=0; k<j; k++) 867dc7d240cSValeria Barra collograd1d[j+Q1d*i] -= interp1d[j+P1d*k]*collograd1d[k+Q1d*i]; 868dc7d240cSValeria Barra collograd1d[j+Q1d*i] /= interp1d[j+P1d*j]; 869d7b241e6Sjeremylt } 8701d102b48SJeremy L Thompson for (j=P1d; j<Q1d; j++) 871dc7d240cSValeria Barra collograd1d[j+Q1d*i] = 0; 872d7b241e6Sjeremylt } 873d7b241e6Sjeremylt 874dc7d240cSValeria Barra // Apply Qtranspose, collograd = collograd Qtranspose 875dc7d240cSValeria Barra CeedHouseholderApplyQ(collograd1d, interp1d, tau, CEED_NOTRANSPOSE, 876d7b241e6Sjeremylt Q1d, Q1d, P1d, 1, Q1d); 877d7b241e6Sjeremylt 878d7b241e6Sjeremylt ierr = CeedFree(&interp1d); CeedChk(ierr); 879d7b241e6Sjeremylt ierr = CeedFree(&grad1d); CeedChk(ierr); 880d7b241e6Sjeremylt 881d7b241e6Sjeremylt return 0; 882d7b241e6Sjeremylt } 883d7b241e6Sjeremylt 884b11c1e72Sjeremylt /** 88595bb1877Svaleriabarra @brief Apply basis evaluation from nodes to quadrature points or vice versa 886b11c1e72Sjeremylt 887b11c1e72Sjeremylt @param basis CeedBasis to evaluate 888b11c1e72Sjeremylt @param nelem The number of elements to apply the basis evaluation to; 889b11c1e72Sjeremylt the backend will specify the ordering in 890b11c1e72Sjeremylt ElemRestrictionCreateBlocked 891b11c1e72Sjeremylt @param tmode \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature 892b11c1e72Sjeremylt points, \ref CEED_TRANSPOSE to apply the transpose, mapping 893b11c1e72Sjeremylt from quadrature points to nodes 894bb64980dSjeremylt @param emode \ref CEED_EVAL_NONE to use values directly, 895bb64980dSjeremylt \ref CEED_EVAL_INTERP to use interpolated values, 896bb64980dSjeremylt \ref CEED_EVAL_GRAD to use gradients, 897bb64980dSjeremylt \ref CEED_EVAL_WEIGHT to use quadrature weights. 89834138859Sjeremylt @param[in] u Input CeedVector 89934138859Sjeremylt @param[out] v Output CeedVector 900b11c1e72Sjeremylt 901b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 902dfdf5a53Sjeremylt 903dfdf5a53Sjeremylt @ref Advanced 904b11c1e72Sjeremylt **/ 905d7b241e6Sjeremylt int CeedBasisApply(CeedBasis basis, CeedInt nelem, CeedTransposeMode tmode, 906aedaa0e5Sjeremylt CeedEvalMode emode, CeedVector u, CeedVector v) { 907d7b241e6Sjeremylt int ierr; 9088795c945Sjeremylt CeedInt ulength = 0, vlength, nnodes, nqpt; 909c042f62fSJeremy L Thompson if (!basis->Apply) 910c042f62fSJeremy L Thompson // LCOV_EXCL_START 911c042f62fSJeremy L Thompson return CeedError(basis->ceed, 1, "Backend does not support BasisApply"); 912c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 913c042f62fSJeremy L Thompson 914c042f62fSJeremy L Thompson // Check compatibility of topological and geometrical dimensions 9158795c945Sjeremylt ierr = CeedBasisGetNumNodes(basis, &nnodes); CeedChk(ierr); 916b502e64cSValeria Barra ierr = CeedBasisGetNumQuadraturePoints(basis, &nqpt); CeedChk(ierr); 917b502e64cSValeria Barra ierr = CeedVectorGetLength(v, &vlength); CeedChk(ierr); 918b502e64cSValeria Barra 919b502e64cSValeria Barra if (u) { 920b502e64cSValeria Barra ierr = CeedVectorGetLength(u, &ulength); CeedChk(ierr); 921b502e64cSValeria Barra } 922b502e64cSValeria Barra 9231d102b48SJeremy L Thompson if ((tmode == CEED_TRANSPOSE && (vlength%nnodes != 0 || ulength%nqpt != 0)) || 9248795c945Sjeremylt (tmode == CEED_NOTRANSPOSE && (ulength%nnodes != 0 || vlength%nqpt != 0))) 9251d102b48SJeremy L Thompson return CeedError(basis->ceed, 1, "Length of input/output vectors " 9261d102b48SJeremy L Thompson "incompatible with basis dimensions"); 927b502e64cSValeria Barra 928d7b241e6Sjeremylt ierr = basis->Apply(basis, nelem, tmode, emode, u, v); CeedChk(ierr); 929d7b241e6Sjeremylt return 0; 930d7b241e6Sjeremylt } 931d7b241e6Sjeremylt 932b11c1e72Sjeremylt /** 9334ce2993fSjeremylt @brief Get Ceed associated with a CeedBasis 934b11c1e72Sjeremylt 935b11c1e72Sjeremylt @param basis CeedBasis 9364ce2993fSjeremylt @param[out] ceed Variable to store Ceed 9374ce2993fSjeremylt 9384ce2993fSjeremylt @return An error code: 0 - success, otherwise - failure 9394ce2993fSjeremylt 94023617272Sjeremylt @ref Advanced 9414ce2993fSjeremylt **/ 9424ce2993fSjeremylt int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) { 9434ce2993fSjeremylt *ceed = basis->ceed; 9444ce2993fSjeremylt return 0; 9454ce2993fSjeremylt }; 9464ce2993fSjeremylt 9474ce2993fSjeremylt /** 9484ce2993fSjeremylt @brief Get dimension for given CeedBasis 9494ce2993fSjeremylt 9504ce2993fSjeremylt @param basis CeedBasis 9514ce2993fSjeremylt @param[out] dim Variable to store dimension of basis 9524ce2993fSjeremylt 9534ce2993fSjeremylt @return An error code: 0 - success, otherwise - failure 9544ce2993fSjeremylt 95523617272Sjeremylt @ref Advanced 9564ce2993fSjeremylt **/ 9574ce2993fSjeremylt int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) { 9584ce2993fSjeremylt *dim = basis->dim; 9594ce2993fSjeremylt return 0; 9604ce2993fSjeremylt }; 9614ce2993fSjeremylt 9624ce2993fSjeremylt /** 9634ce2993fSjeremylt @brief Get tensor status for given CeedBasis 9644ce2993fSjeremylt 9654ce2993fSjeremylt @param basis CeedBasis 9664ce2993fSjeremylt @param[out] tensor Variable to store tensor status 9674ce2993fSjeremylt 9684ce2993fSjeremylt @return An error code: 0 - success, otherwise - failure 9694ce2993fSjeremylt 97023617272Sjeremylt @ref Advanced 9714ce2993fSjeremylt **/ 9724ce2993fSjeremylt int CeedBasisGetTensorStatus(CeedBasis basis, bool *tensor) { 9734ce2993fSjeremylt *tensor = basis->tensorbasis; 9744ce2993fSjeremylt return 0; 9754ce2993fSjeremylt }; 9764ce2993fSjeremylt 9774ce2993fSjeremylt /** 9784ce2993fSjeremylt @brief Get number of components for given CeedBasis 9794ce2993fSjeremylt 9804ce2993fSjeremylt @param basis CeedBasis 981288c0443SJeremy L Thompson @param[out] numcomp Variable to store number of components of basis 9824ce2993fSjeremylt 9834ce2993fSjeremylt @return An error code: 0 - success, otherwise - failure 9844ce2993fSjeremylt 98523617272Sjeremylt @ref Advanced 9864ce2993fSjeremylt **/ 9874ce2993fSjeremylt int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *numcomp) { 9884ce2993fSjeremylt *numcomp = basis->ncomp; 9894ce2993fSjeremylt return 0; 9904ce2993fSjeremylt }; 9914ce2993fSjeremylt 9924ce2993fSjeremylt /** 9934ce2993fSjeremylt @brief Get total number of nodes (in 1 dimension) of a CeedBasis 9944ce2993fSjeremylt 9954ce2993fSjeremylt @param basis CeedBasis 9964ce2993fSjeremylt @param[out] P1d Variable to store number of nodes 9974ce2993fSjeremylt 9984ce2993fSjeremylt @return An error code: 0 - success, otherwise - failure 9994ce2993fSjeremylt 100023617272Sjeremylt @ref Advanced 10014ce2993fSjeremylt **/ 10024ce2993fSjeremylt int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P1d) { 1003c042f62fSJeremy L Thompson if (!basis->tensorbasis) 1004c042f62fSJeremy L Thompson // LCOV_EXCL_START 1005c042f62fSJeremy L Thompson return CeedError(basis->ceed, 1, "Cannot supply P1d for non-tensor basis"); 1006c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 1007c042f62fSJeremy L Thompson 10084ce2993fSjeremylt *P1d = basis->P1d; 10094ce2993fSjeremylt return 0; 10104ce2993fSjeremylt } 10114ce2993fSjeremylt 10124ce2993fSjeremylt /** 10134ce2993fSjeremylt @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis 10144ce2993fSjeremylt 10154ce2993fSjeremylt @param basis CeedBasis 10164ce2993fSjeremylt @param[out] Q1d Variable to store number of quadrature points 10174ce2993fSjeremylt 10184ce2993fSjeremylt @return An error code: 0 - success, otherwise - failure 10194ce2993fSjeremylt 102023617272Sjeremylt @ref Advanced 10214ce2993fSjeremylt **/ 10224ce2993fSjeremylt int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q1d) { 1023c042f62fSJeremy L Thompson if (!basis->tensorbasis) 1024c042f62fSJeremy L Thompson // LCOV_EXCL_START 1025c042f62fSJeremy L Thompson return CeedError(basis->ceed, 1, "Cannot supply Q1d for non-tensor basis"); 1026c042f62fSJeremy L Thompson // LCOV_EXCL_STOP 1027c042f62fSJeremy L Thompson 10284ce2993fSjeremylt *Q1d = basis->Q1d; 10294ce2993fSjeremylt return 0; 10304ce2993fSjeremylt } 10314ce2993fSjeremylt 10324ce2993fSjeremylt /** 10334ce2993fSjeremylt @brief Get total number of nodes (in dim dimensions) of a CeedBasis 10344ce2993fSjeremylt 10354ce2993fSjeremylt @param basis CeedBasis 10364ce2993fSjeremylt @param[out] P Variable to store number of nodes 1037b11c1e72Sjeremylt 1038b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 1039dfdf5a53Sjeremylt 1040dfdf5a53Sjeremylt @ref Utility 1041b11c1e72Sjeremylt **/ 1042d7b241e6Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) { 1043a8de75f0Sjeremylt *P = basis->P; 1044d7b241e6Sjeremylt return 0; 1045d7b241e6Sjeremylt } 1046d7b241e6Sjeremylt 1047b11c1e72Sjeremylt /** 10484ce2993fSjeremylt @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis 1049b11c1e72Sjeremylt 1050b11c1e72Sjeremylt @param basis CeedBasis 10514ce2993fSjeremylt @param[out] Q Variable to store number of quadrature points 1052b11c1e72Sjeremylt 1053b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 1054dfdf5a53Sjeremylt 1055dfdf5a53Sjeremylt @ref Utility 1056b11c1e72Sjeremylt **/ 1057d7b241e6Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) { 1058a8de75f0Sjeremylt *Q = basis->Q; 1059d7b241e6Sjeremylt return 0; 1060d7b241e6Sjeremylt } 1061d7b241e6Sjeremylt 1062b11c1e72Sjeremylt /** 10638c91a0c9SJeremy L Thompson @brief Get reference coordinates of quadrature points (in dim dimensions) 10644ce2993fSjeremylt of a CeedBasis 10654ce2993fSjeremylt 10664ce2993fSjeremylt @param basis CeedBasis 10678c91a0c9SJeremy L Thompson @param[out] qref Variable to store reference coordinates of quadrature points 10684ce2993fSjeremylt 10694ce2993fSjeremylt @return An error code: 0 - success, otherwise - failure 10704ce2993fSjeremylt 107123617272Sjeremylt @ref Advanced 10724ce2993fSjeremylt **/ 10734ce2993fSjeremylt int CeedBasisGetQRef(CeedBasis basis, CeedScalar **qref) { 10744ce2993fSjeremylt *qref = basis->qref1d; 10754ce2993fSjeremylt return 0; 10764ce2993fSjeremylt } 10774ce2993fSjeremylt 10784ce2993fSjeremylt /** 10794ce2993fSjeremylt @brief Get quadrature weights of quadrature points (in dim dimensions) 10804ce2993fSjeremylt of a CeedBasis 10814ce2993fSjeremylt 10824ce2993fSjeremylt @param basis CeedBasis 10834ce2993fSjeremylt @param[out] qweight Variable to store quadrature weights 10844ce2993fSjeremylt 10854ce2993fSjeremylt @return An error code: 0 - success, otherwise - failure 10864ce2993fSjeremylt 108723617272Sjeremylt @ref Advanced 10884ce2993fSjeremylt **/ 10894ce2993fSjeremylt int CeedBasisGetQWeights(CeedBasis basis, CeedScalar **qweight) { 10904ce2993fSjeremylt *qweight = basis->qweight1d; 10914ce2993fSjeremylt return 0; 10924ce2993fSjeremylt } 10934ce2993fSjeremylt 10944ce2993fSjeremylt /** 10954ce2993fSjeremylt @brief Get interpolation matrix of a CeedBasis 10964ce2993fSjeremylt 10974ce2993fSjeremylt @param basis CeedBasis 1098288c0443SJeremy L Thompson @param[out] interp Variable to store interpolation matrix 10994ce2993fSjeremylt 11004ce2993fSjeremylt @return An error code: 0 - success, otherwise - failure 11014ce2993fSjeremylt 110223617272Sjeremylt @ref Advanced 11034ce2993fSjeremylt **/ 11044ce2993fSjeremylt int CeedBasisGetInterp(CeedBasis basis, CeedScalar **interp) { 110500f91b2bSjeremylt if (!basis->interp && basis->tensorbasis) { 110600f91b2bSjeremylt // Allocate 110700f91b2bSjeremylt int ierr; 110800f91b2bSjeremylt ierr = CeedMalloc(basis->Q*basis->P, &basis->interp); CeedChk(ierr); 110900f91b2bSjeremylt 111000f91b2bSjeremylt // Initialize 111100f91b2bSjeremylt for (CeedInt i=0; i<basis->Q*basis->P; i++) 111200f91b2bSjeremylt basis->interp[i] = 1.0; 111300f91b2bSjeremylt 111400f91b2bSjeremylt // Calculate 111500f91b2bSjeremylt for (CeedInt d=0; d<basis->dim; d++) 111600f91b2bSjeremylt for (CeedInt qpt=0; qpt<basis->Q; qpt++) 111700f91b2bSjeremylt for (CeedInt node=0; node<basis->P; node++) { 111800f91b2bSjeremylt CeedInt p = (node / CeedIntPow(basis->P1d, d)) % basis->P1d; 111900f91b2bSjeremylt CeedInt q = (qpt / CeedIntPow(basis->Q1d, d)) % basis->Q1d; 112000f91b2bSjeremylt basis->interp[qpt*(basis->P)+node] *= basis->interp1d[q*basis->P1d+p]; 112100f91b2bSjeremylt } 112200f91b2bSjeremylt } 112300f91b2bSjeremylt 112400f91b2bSjeremylt *interp = basis->interp; 112500f91b2bSjeremylt 112600f91b2bSjeremylt return 0; 112700f91b2bSjeremylt } 112800f91b2bSjeremylt 112900f91b2bSjeremylt /** 113000f91b2bSjeremylt @brief Get 1D interpolation matrix of a tensor product CeedBasis 113100f91b2bSjeremylt 113200f91b2bSjeremylt @param basis CeedBasis 113300f91b2bSjeremylt @param[out] interp1d Variable to store interpolation matrix 113400f91b2bSjeremylt 113500f91b2bSjeremylt @return An error code: 0 - success, otherwise - failure 113600f91b2bSjeremylt 113700f91b2bSjeremylt @ref Advanced 113800f91b2bSjeremylt **/ 113900f91b2bSjeremylt int CeedBasisGetInterp1D(CeedBasis basis, CeedScalar **interp1d) { 114000f91b2bSjeremylt if (!basis->tensorbasis) 114100f91b2bSjeremylt // LCOV_EXCL_START 114200f91b2bSjeremylt return CeedError(basis->ceed, 1, "CeedBasis is not a tensor product basis."); 114300f91b2bSjeremylt // LCOV_EXCL_STOP 114400f91b2bSjeremylt 114500f91b2bSjeremylt *interp1d = basis->interp1d; 114600f91b2bSjeremylt 11474ce2993fSjeremylt return 0; 11484ce2993fSjeremylt } 11494ce2993fSjeremylt 11504ce2993fSjeremylt /** 11514ce2993fSjeremylt @brief Get gradient matrix of a CeedBasis 11524ce2993fSjeremylt 11534ce2993fSjeremylt @param basis CeedBasis 1154288c0443SJeremy L Thompson @param[out] grad Variable to store gradient matrix 11554ce2993fSjeremylt 11564ce2993fSjeremylt @return An error code: 0 - success, otherwise - failure 11574ce2993fSjeremylt 115823617272Sjeremylt @ref Advanced 11594ce2993fSjeremylt **/ 11604ce2993fSjeremylt int CeedBasisGetGrad(CeedBasis basis, CeedScalar **grad) { 116100f91b2bSjeremylt if (!basis->grad && basis->tensorbasis) { 116200f91b2bSjeremylt // Allocate 116300f91b2bSjeremylt int ierr; 116400f91b2bSjeremylt ierr = CeedMalloc(basis->dim*basis->Q*basis->P, &basis->grad); 116500f91b2bSjeremylt CeedChk(ierr); 116600f91b2bSjeremylt 116700f91b2bSjeremylt // Initialize 116800f91b2bSjeremylt for (CeedInt i=0; i<basis->dim*basis->Q*basis->P; i++) 116900f91b2bSjeremylt basis->grad[i] = 1.0; 117000f91b2bSjeremylt 117100f91b2bSjeremylt // Calculate 117200f91b2bSjeremylt for (CeedInt d=0; d<basis->dim; d++) 117300f91b2bSjeremylt for (CeedInt i=0; i<basis->dim; i++) 117400f91b2bSjeremylt for (CeedInt qpt=0; qpt<basis->Q; qpt++) 117500f91b2bSjeremylt for (CeedInt node=0; node<basis->P; node++) { 117600f91b2bSjeremylt CeedInt p = (node / CeedIntPow(basis->P1d, d)) % basis->P1d; 117700f91b2bSjeremylt CeedInt q = (qpt / CeedIntPow(basis->Q1d, d)) % basis->Q1d; 117800f91b2bSjeremylt if (i == d) 117900f91b2bSjeremylt basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *= 118000f91b2bSjeremylt basis->grad1d[q*basis->P1d+p]; 118100f91b2bSjeremylt else 118200f91b2bSjeremylt basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *= 118300f91b2bSjeremylt basis->interp1d[q*basis->P1d+p]; 118400f91b2bSjeremylt } 118500f91b2bSjeremylt } 118600f91b2bSjeremylt 118700f91b2bSjeremylt *grad = basis->grad; 118800f91b2bSjeremylt 11894ce2993fSjeremylt return 0; 11904ce2993fSjeremylt } 11914ce2993fSjeremylt 11924ce2993fSjeremylt /** 119300f91b2bSjeremylt @brief Get 1D gradient matrix of a tensor product CeedBasis 1194b7ec98d8SJeremy L Thompson 1195fb551037Sjeremylt @param basis CeedBasis 119600f91b2bSjeremylt @param[out] grad1d Variable to store gradient matrix 1197b7ec98d8SJeremy L Thompson 1198b7ec98d8SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 1199b7ec98d8SJeremy L Thompson 1200b7ec98d8SJeremy L Thompson @ref Advanced 1201b7ec98d8SJeremy L Thompson **/ 120200f91b2bSjeremylt int CeedBasisGetGrad1D(CeedBasis basis, CeedScalar **grad1d) { 120300f91b2bSjeremylt if (!basis->tensorbasis) 1204b7ec98d8SJeremy L Thompson // LCOV_EXCL_START 120500f91b2bSjeremylt return CeedError(basis->ceed, 1, "CeedBasis is not a tensor product basis."); 1206b7ec98d8SJeremy L Thompson // LCOV_EXCL_STOP 120700f91b2bSjeremylt 120800f91b2bSjeremylt *grad1d = basis->grad1d; 120900f91b2bSjeremylt 1210b7ec98d8SJeremy L Thompson return 0; 1211b7ec98d8SJeremy L Thompson } 1212b7ec98d8SJeremy L Thompson 1213b7ec98d8SJeremy L Thompson /** 12144ce2993fSjeremylt @brief Get backend data of a CeedBasis 12154ce2993fSjeremylt 12164ce2993fSjeremylt @param basis CeedBasis 12174ce2993fSjeremylt @param[out] data Variable to store data 12184ce2993fSjeremylt 12194ce2993fSjeremylt @return An error code: 0 - success, otherwise - failure 12204ce2993fSjeremylt 122123617272Sjeremylt @ref Advanced 12224ce2993fSjeremylt **/ 12234ce2993fSjeremylt int CeedBasisGetData(CeedBasis basis, void **data) { 12244ce2993fSjeremylt *data = basis->data; 12254ce2993fSjeremylt return 0; 12264ce2993fSjeremylt } 12274ce2993fSjeremylt 12284ce2993fSjeremylt /** 1229fe2413ffSjeremylt @brief Set backend data of a CeedBasis 1230fe2413ffSjeremylt 1231fe2413ffSjeremylt @param[out] basis CeedBasis 1232fe2413ffSjeremylt @param data Data to set 1233fe2413ffSjeremylt 1234fe2413ffSjeremylt @return An error code: 0 - success, otherwise - failure 1235fe2413ffSjeremylt 1236fe2413ffSjeremylt @ref Advanced 1237fe2413ffSjeremylt **/ 1238fe2413ffSjeremylt int CeedBasisSetData(CeedBasis basis, void **data) { 1239fe2413ffSjeremylt basis->data = *data; 1240fe2413ffSjeremylt return 0; 1241fe2413ffSjeremylt } 1242fe2413ffSjeremylt 1243fe2413ffSjeremylt /** 12442f86a920SJeremy L Thompson @brief Get CeedTensorContract of a CeedBasis 12452f86a920SJeremy L Thompson 12462f86a920SJeremy L Thompson @param basis CeedBasis 12472f86a920SJeremy L Thompson @param[out] contract Variable to store CeedTensorContract 12482f86a920SJeremy L Thompson 12492f86a920SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 12502f86a920SJeremy L Thompson 12512f86a920SJeremy L Thompson @ref Advanced 12522f86a920SJeremy L Thompson **/ 12531d102b48SJeremy L Thompson int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) { 12542f86a920SJeremy L Thompson *contract = basis->contract; 12552f86a920SJeremy L Thompson return 0; 12562f86a920SJeremy L Thompson } 12572f86a920SJeremy L Thompson 12582f86a920SJeremy L Thompson /** 12592f86a920SJeremy L Thompson @brief Set CeedTensorContract of a CeedBasis 12602f86a920SJeremy L Thompson 12612f86a920SJeremy L Thompson @param[out] basis CeedBasis 12622f86a920SJeremy L Thompson @param contract CeedTensorContract to set 12632f86a920SJeremy L Thompson 12642f86a920SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 12652f86a920SJeremy L Thompson 12662f86a920SJeremy L Thompson @ref Advanced 12672f86a920SJeremy L Thompson **/ 12681d102b48SJeremy L Thompson int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract *contract) { 12692f86a920SJeremy L Thompson basis->contract = *contract; 12702f86a920SJeremy L Thompson return 0; 12712f86a920SJeremy L Thompson } 12722f86a920SJeremy L Thompson 12732f86a920SJeremy L Thompson /** 1274a8de75f0Sjeremylt @brief Get dimension for given CeedElemTopology 1275a8de75f0Sjeremylt 1276a8de75f0Sjeremylt @param topo CeedElemTopology 12774ce2993fSjeremylt @param[out] dim Variable to store dimension of topology 1278a8de75f0Sjeremylt 1279a8de75f0Sjeremylt @return An error code: 0 - success, otherwise - failure 1280a8de75f0Sjeremylt 128123617272Sjeremylt @ref Advanced 1282a8de75f0Sjeremylt **/ 1283a8de75f0Sjeremylt int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) { 1284a8de75f0Sjeremylt *dim = (CeedInt) topo >> 16; 1285a8de75f0Sjeremylt return 0; 1286a8de75f0Sjeremylt }; 1287a8de75f0Sjeremylt 1288a8de75f0Sjeremylt /** 1289b11c1e72Sjeremylt @brief Destroy a CeedBasis 1290b11c1e72Sjeremylt 1291b11c1e72Sjeremylt @param basis CeedBasis to destroy 1292b11c1e72Sjeremylt 1293b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 1294dfdf5a53Sjeremylt 1295dfdf5a53Sjeremylt @ref Basic 1296b11c1e72Sjeremylt **/ 1297d7b241e6Sjeremylt int CeedBasisDestroy(CeedBasis *basis) { 1298d7b241e6Sjeremylt int ierr; 1299d7b241e6Sjeremylt 13001d102b48SJeremy L Thompson if (!*basis || --(*basis)->refcount > 0) 13011d102b48SJeremy L Thompson return 0; 1302d7b241e6Sjeremylt if ((*basis)->Destroy) { 1303d7b241e6Sjeremylt ierr = (*basis)->Destroy(*basis); CeedChk(ierr); 1304d7b241e6Sjeremylt } 130500f91b2bSjeremylt ierr = CeedFree(&(*basis)->interp); CeedChk(ierr); 1306d7b241e6Sjeremylt ierr = CeedFree(&(*basis)->interp1d); CeedChk(ierr); 130700f91b2bSjeremylt ierr = CeedFree(&(*basis)->grad); CeedChk(ierr); 1308d7b241e6Sjeremylt ierr = CeedFree(&(*basis)->grad1d); CeedChk(ierr); 1309d7b241e6Sjeremylt ierr = CeedFree(&(*basis)->qref1d); CeedChk(ierr); 1310d7b241e6Sjeremylt ierr = CeedFree(&(*basis)->qweight1d); CeedChk(ierr); 1311d7b241e6Sjeremylt ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr); 1312d7b241e6Sjeremylt ierr = CeedFree(basis); CeedChk(ierr); 1313d7b241e6Sjeremylt return 0; 1314d7b241e6Sjeremylt } 1315d7b241e6Sjeremylt 131633e6becaSjeremylt /// @cond DOXYGEN_SKIP 13178795c945Sjeremylt // Indicate that the quadrature points are collocated with the nodes 1318783c99b3SValeria Barra CeedBasis CEED_BASIS_COLLOCATED = &ceed_basis_collocated; 131933e6becaSjeremylt /// @endcond 1320d7b241e6Sjeremylt /// @} 1321