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> 18d7b241e6Sjeremylt #include <math.h> 19d7b241e6Sjeremylt #include <stdio.h> 20d7b241e6Sjeremylt #include <stdlib.h> 21d7b241e6Sjeremylt #include <string.h> 22d7b241e6Sjeremylt 23d7b241e6Sjeremylt /// @cond DOXYGEN_SKIP 24783c99b3SValeria Barra static struct CeedBasis_private ceed_basis_collocated; 25d7b241e6Sjeremylt /// @endcond 26d7b241e6Sjeremylt 27d7b241e6Sjeremylt /// @file 28d7b241e6Sjeremylt /// Implementation of public CeedBasis interfaces 29d7b241e6Sjeremylt /// 30dfdf5a53Sjeremylt /// @addtogroup CeedBasis 31d7b241e6Sjeremylt /// @{ 32d7b241e6Sjeremylt 33b11c1e72Sjeremylt /** 34b11c1e72Sjeremylt @brief Create a tensor product basis for H^1 discretizations 35b11c1e72Sjeremylt 36b11c1e72Sjeremylt @param ceed A Ceed object where the CeedBasis will be created 37b11c1e72Sjeremylt @param dim Topological dimension 38b11c1e72Sjeremylt @param ncomp Number of field components (1 for scalar fields) 39b11c1e72Sjeremylt @param P1d Number of nodes in one dimension 40b11c1e72Sjeremylt @param Q1d Number of quadrature points in one dimension 41b11c1e72Sjeremylt @param interp1d Row-major Q1d × P1d matrix expressing the values of nodal 42b11c1e72Sjeremylt basis functions at quadrature points 43b11c1e72Sjeremylt @param grad1d Row-major Q1d × P1d matrix expressing derivatives of nodal 44b11c1e72Sjeremylt basis functions at quadrature points 45b11c1e72Sjeremylt @param qref1d Array of length Q1d holding the locations of quadrature points 46b11c1e72Sjeremylt on the 1D reference element [-1, 1] 47b11c1e72Sjeremylt @param qweight1d Array of length Q1d holding the quadrature weights on the 48b11c1e72Sjeremylt reference element 49b11c1e72Sjeremylt @param[out] basis Address of the variable where the newly created 50b11c1e72Sjeremylt CeedBasis will be stored. 51b11c1e72Sjeremylt 52b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 53dfdf5a53Sjeremylt 54dfdf5a53Sjeremylt @ref Basic 55b11c1e72Sjeremylt **/ 56d7b241e6Sjeremylt int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt ncomp, CeedInt P1d, 57d7b241e6Sjeremylt CeedInt Q1d, const CeedScalar *interp1d, 58d7b241e6Sjeremylt const CeedScalar *grad1d, const CeedScalar *qref1d, 59d7b241e6Sjeremylt const CeedScalar *qweight1d, CeedBasis *basis) { 60d7b241e6Sjeremylt int ierr; 61d7b241e6Sjeremylt 62d7b241e6Sjeremylt if (!ceed->BasisCreateTensorH1) 63d7b241e6Sjeremylt return CeedError(ceed, 1, "Backend does not support BasisCreateTensorH1"); 64d7b241e6Sjeremylt ierr = CeedCalloc(1,basis); CeedChk(ierr); 65d7b241e6Sjeremylt (*basis)->ceed = ceed; 66d7b241e6Sjeremylt ceed->refcount++; 67d7b241e6Sjeremylt (*basis)->refcount = 1; 68*a8de75f0Sjeremylt (*basis)->tensorbasis = 1; 69d7b241e6Sjeremylt (*basis)->dim = dim; 70d7b241e6Sjeremylt (*basis)->ncomp = ncomp; 71d7b241e6Sjeremylt (*basis)->P1d = P1d; 72d7b241e6Sjeremylt (*basis)->Q1d = Q1d; 73*a8de75f0Sjeremylt (*basis)->P = CeedIntPow(P1d, dim); 74*a8de75f0Sjeremylt (*basis)->Q = CeedIntPow(Q1d, dim); 75d7b241e6Sjeremylt ierr = CeedMalloc(Q1d,&(*basis)->qref1d); CeedChk(ierr); 76d7b241e6Sjeremylt ierr = CeedMalloc(Q1d,&(*basis)->qweight1d); CeedChk(ierr); 77d7b241e6Sjeremylt memcpy((*basis)->qref1d, qref1d, Q1d*sizeof(qref1d[0])); 78d7b241e6Sjeremylt memcpy((*basis)->qweight1d, qweight1d, Q1d*sizeof(qweight1d[0])); 79d7b241e6Sjeremylt ierr = CeedMalloc(Q1d*P1d,&(*basis)->interp1d); CeedChk(ierr); 80d7b241e6Sjeremylt ierr = CeedMalloc(Q1d*P1d,&(*basis)->grad1d); CeedChk(ierr); 81d7b241e6Sjeremylt memcpy((*basis)->interp1d, interp1d, Q1d*P1d*sizeof(interp1d[0])); 8209486605Sjeremylt memcpy((*basis)->grad1d, grad1d, Q1d*P1d*sizeof(grad1d[0])); 83d7b241e6Sjeremylt ierr = ceed->BasisCreateTensorH1(ceed, dim, P1d, Q1d, interp1d, grad1d, qref1d, 84d7b241e6Sjeremylt qweight1d, *basis); CeedChk(ierr); 85d7b241e6Sjeremylt return 0; 86d7b241e6Sjeremylt } 87d7b241e6Sjeremylt 88b11c1e72Sjeremylt /** 89b11c1e72Sjeremylt @brief Create a tensor product Lagrange basis 90b11c1e72Sjeremylt 91b11c1e72Sjeremylt @param ceed A Ceed object where the CeedBasis will be created 92b11c1e72Sjeremylt @param dim Topological dimension of element 93b11c1e72Sjeremylt @param ncomp Number of field components 94b11c1e72Sjeremylt @param P Number of Gauss-Lobatto nodes in one dimension. The 95b11c1e72Sjeremylt polynomial degree of the resulting Q_k element is k=P-1. 96b11c1e72Sjeremylt @param Q Number of quadrature points in one dimension. 97b11c1e72Sjeremylt @param qmode Distribution of the Q quadrature points (affects order of 98b11c1e72Sjeremylt accuracy for the quadrature) 99b11c1e72Sjeremylt @param[out] basis Address of the variable where the newly created 100b11c1e72Sjeremylt CeedBasis will be stored. 101b11c1e72Sjeremylt 102b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 103dfdf5a53Sjeremylt 104dfdf5a53Sjeremylt @ref Basic 105b11c1e72Sjeremylt **/ 106d7b241e6Sjeremylt int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt ncomp, 107d7b241e6Sjeremylt CeedInt P, CeedInt Q, 108d7b241e6Sjeremylt CeedQuadMode qmode, CeedBasis *basis) { 109d7b241e6Sjeremylt // Allocate 110d7b241e6Sjeremylt int ierr, i, j, k; 111d7b241e6Sjeremylt CeedScalar c1, c2, c3, c4, dx, *nodes, *interp1d, *grad1d, *qref1d, *qweight1d; 112d7b241e6Sjeremylt ierr = CeedCalloc(P*Q, &interp1d); CeedChk(ierr); 113d7b241e6Sjeremylt ierr = CeedCalloc(P*Q, &grad1d); CeedChk(ierr); 114d7b241e6Sjeremylt ierr = CeedCalloc(P, &nodes); CeedChk(ierr); 115d7b241e6Sjeremylt ierr = CeedCalloc(Q, &qref1d); CeedChk(ierr); 116d7b241e6Sjeremylt ierr = CeedCalloc(Q, &qweight1d); CeedChk(ierr); 117d7b241e6Sjeremylt // Get Nodes and Weights 118d7b241e6Sjeremylt ierr = CeedLobattoQuadrature(P, nodes, NULL); CeedChk(ierr); 119d7b241e6Sjeremylt switch (qmode) { 120d7b241e6Sjeremylt case CEED_GAUSS: 121d7b241e6Sjeremylt ierr = CeedGaussQuadrature(Q, qref1d, qweight1d); CeedChk(ierr); 122d7b241e6Sjeremylt break; 123d7b241e6Sjeremylt case CEED_GAUSS_LOBATTO: 124d7b241e6Sjeremylt ierr = CeedLobattoQuadrature(Q, qref1d, qweight1d); CeedChk(ierr); 125d7b241e6Sjeremylt break; 126d7b241e6Sjeremylt } 127d7b241e6Sjeremylt // Build B, D matrix 128d7b241e6Sjeremylt // Fornberg, 1998 129d7b241e6Sjeremylt for (i = 0; i < Q; i++) { 130d7b241e6Sjeremylt c1 = 1.0; 131d7b241e6Sjeremylt c3 = nodes[0] - qref1d[i]; 132d7b241e6Sjeremylt interp1d[i*P+0] = 1.0; 133d7b241e6Sjeremylt for (j = 1; j < P; j++) { 134d7b241e6Sjeremylt c2 = 1.0; 135d7b241e6Sjeremylt c4 = c3; 136d7b241e6Sjeremylt c3 = nodes[j] - qref1d[i]; 137d7b241e6Sjeremylt for (k = 0; k < j; k++) { 138d7b241e6Sjeremylt dx = nodes[j] - nodes[k]; 139d7b241e6Sjeremylt c2 *= dx; 140d7b241e6Sjeremylt if (k == j - 1) { 141d7b241e6Sjeremylt grad1d[i*P + j] = c1*(interp1d[i*P + k] - c4*grad1d[i*P + k]) / c2; 142d7b241e6Sjeremylt interp1d[i*P + j] = - c1*c4*interp1d[i*P + k] / c2; 143d7b241e6Sjeremylt } 144d7b241e6Sjeremylt grad1d[i*P + k] = (c3*grad1d[i*P + k] - interp1d[i*P + k]) / dx; 145d7b241e6Sjeremylt interp1d[i*P + k] = c3*interp1d[i*P + k] / dx; 146d7b241e6Sjeremylt } 147d7b241e6Sjeremylt c1 = c2; 148d7b241e6Sjeremylt } 149d7b241e6Sjeremylt } 150d7b241e6Sjeremylt // // Pass to CeedBasisCreateTensorH1 151d7b241e6Sjeremylt ierr = CeedBasisCreateTensorH1(ceed, dim, ncomp, P, Q, interp1d, grad1d, qref1d, 152d7b241e6Sjeremylt qweight1d, basis); CeedChk(ierr); 153d7b241e6Sjeremylt ierr = CeedFree(&interp1d); CeedChk(ierr); 154d7b241e6Sjeremylt ierr = CeedFree(&grad1d); CeedChk(ierr); 155d7b241e6Sjeremylt ierr = CeedFree(&nodes); CeedChk(ierr); 156d7b241e6Sjeremylt ierr = CeedFree(&qref1d); CeedChk(ierr); 157d7b241e6Sjeremylt ierr = CeedFree(&qweight1d); CeedChk(ierr); 158d7b241e6Sjeremylt return 0; 159d7b241e6Sjeremylt } 160d7b241e6Sjeremylt 161b11c1e72Sjeremylt /** 162*a8de75f0Sjeremylt @brief Create a non tensor product basis for H^1 discretizations 163*a8de75f0Sjeremylt 164*a8de75f0Sjeremylt @param ceed A Ceed object where the CeedBasis will be created 165*a8de75f0Sjeremylt @param topo Topology of element, e.g. hypercube, simplex, ect 166*a8de75f0Sjeremylt @param ncomp Number of field components (1 for scalar fields) 167*a8de75f0Sjeremylt @param ndof Total number of nodes 168*a8de75f0Sjeremylt @param nqpts Total number of quadrature points 169*a8de75f0Sjeremylt @param interp Row-major nqpts × ndof matrix expressing the values of nodal 170*a8de75f0Sjeremylt basis functions at quadrature points 171*a8de75f0Sjeremylt @param grad Row-major (nqpts x dim) × ndof matrix expressing derivatives 172*a8de75f0Sjeremylt of nodal basis functions at quadrature points 173*a8de75f0Sjeremylt @param qref Array of length nqpts holding the locations of quadrature points 174*a8de75f0Sjeremylt on the reference element [-1, 1] 175*a8de75f0Sjeremylt @param qweight Array of length nqpts holding the quadrature weights on the 176*a8de75f0Sjeremylt reference element 177*a8de75f0Sjeremylt @param[out] basis Address of the variable where the newly created 178*a8de75f0Sjeremylt CeedBasis will be stored. 179*a8de75f0Sjeremylt 180*a8de75f0Sjeremylt @return An error code: 0 - success, otherwise - failure 181*a8de75f0Sjeremylt 182*a8de75f0Sjeremylt @ref Basic 183*a8de75f0Sjeremylt **/ 184*a8de75f0Sjeremylt int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt ncomp, 185*a8de75f0Sjeremylt CeedInt ndof, CeedInt nqpts, 186*a8de75f0Sjeremylt const CeedScalar *interp, 187*a8de75f0Sjeremylt const CeedScalar *grad, const CeedScalar *qref, 188*a8de75f0Sjeremylt const CeedScalar *qweight, CeedBasis *basis) { 189*a8de75f0Sjeremylt int ierr; 190*a8de75f0Sjeremylt CeedInt P = ndof, Q = nqpts, dim = 0; 191*a8de75f0Sjeremylt 192*a8de75f0Sjeremylt if (!ceed->BasisCreateH1) 193*a8de75f0Sjeremylt return CeedError(ceed, 1, "Backend does not support BasisCreateH1"); 194*a8de75f0Sjeremylt ierr = CeedCalloc(1,basis); CeedChk(ierr); 195*a8de75f0Sjeremylt 196*a8de75f0Sjeremylt ierr = CeedBasisGetTopologyDimension(topo, &dim); CeedChk(ierr); 197*a8de75f0Sjeremylt 198*a8de75f0Sjeremylt (*basis)->ceed = ceed; 199*a8de75f0Sjeremylt ceed->refcount++; 200*a8de75f0Sjeremylt (*basis)->refcount = 1; 201*a8de75f0Sjeremylt (*basis)->tensorbasis = 0; 202*a8de75f0Sjeremylt (*basis)->dim = dim; 203*a8de75f0Sjeremylt (*basis)->ncomp = ncomp; 204*a8de75f0Sjeremylt (*basis)->P = P; 205*a8de75f0Sjeremylt (*basis)->Q = Q; 206*a8de75f0Sjeremylt ierr = CeedMalloc(Q*dim,&(*basis)->qref1d); CeedChk(ierr); 207*a8de75f0Sjeremylt ierr = CeedMalloc(Q,&(*basis)->qweight1d); CeedChk(ierr); 208*a8de75f0Sjeremylt memcpy((*basis)->qref1d, qref, Q*dim*sizeof(qref[0])); 209*a8de75f0Sjeremylt memcpy((*basis)->qweight1d, qweight, Q*sizeof(qweight[0])); 210*a8de75f0Sjeremylt ierr = CeedMalloc(Q*P,&(*basis)->interp1d); CeedChk(ierr); 211*a8de75f0Sjeremylt ierr = CeedMalloc(dim*Q*P,&(*basis)->grad1d); CeedChk(ierr); 212*a8de75f0Sjeremylt memcpy((*basis)->interp1d, interp, Q*P*sizeof(interp[0])); 213*a8de75f0Sjeremylt memcpy((*basis)->grad1d, grad, dim*Q*P*sizeof(grad[0])); 214*a8de75f0Sjeremylt ierr = ceed->BasisCreateH1(ceed, topo, dim, P, Q, interp, grad, qref, 215*a8de75f0Sjeremylt qweight, *basis); CeedChk(ierr); 216*a8de75f0Sjeremylt return 0; 217*a8de75f0Sjeremylt } 218*a8de75f0Sjeremylt 219*a8de75f0Sjeremylt /** 220b11c1e72Sjeremylt @brief Construct a Gauss-Legendre quadrature 221b11c1e72Sjeremylt 222b11c1e72Sjeremylt @param Q Number of quadrature points (integrates polynomials of 223b11c1e72Sjeremylt degree 2*Q-1 exactly) 224b11c1e72Sjeremylt @param[out] qref1d Array of length Q to hold the abscissa on [-1, 1] 225b11c1e72Sjeremylt @param[out] qweight1d Array of length Q to hold the weights 226b11c1e72Sjeremylt 227b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 228dfdf5a53Sjeremylt 229dfdf5a53Sjeremylt @ref Utility 230b11c1e72Sjeremylt **/ 231d7b241e6Sjeremylt int CeedGaussQuadrature(CeedInt Q, CeedScalar *qref1d, CeedScalar *qweight1d) { 232d7b241e6Sjeremylt // Allocate 233d7b241e6Sjeremylt CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0*atan(1.0); 234d7b241e6Sjeremylt // Build qref1d, qweight1d 235d7b241e6Sjeremylt for (int i = 0; i <= Q/2; i++) { 236d7b241e6Sjeremylt // Guess 237d7b241e6Sjeremylt xi = cos(PI*(CeedScalar)(2*i+1)/((CeedScalar)(2*Q))); 238d7b241e6Sjeremylt // Pn(xi) 239d7b241e6Sjeremylt P0 = 1.0; 240d7b241e6Sjeremylt P1 = xi; 241d7b241e6Sjeremylt P2 = 0.0; 242d7b241e6Sjeremylt for (int j = 2; j <= Q; j++) { 243d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 244d7b241e6Sjeremylt P0 = P1; 245d7b241e6Sjeremylt P1 = P2; 246d7b241e6Sjeremylt } 247d7b241e6Sjeremylt // First Newton Step 248d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 249d7b241e6Sjeremylt xi = xi-P2/dP2; 250d7b241e6Sjeremylt // Newton to convergence 251d7b241e6Sjeremylt for (int k=0; k<100 && fabs(P2)>1e-15; k++) { 252d7b241e6Sjeremylt P0 = 1.0; 253d7b241e6Sjeremylt P1 = xi; 254d7b241e6Sjeremylt for (int j = 2; j <= Q; j++) { 255d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 256d7b241e6Sjeremylt P0 = P1; 257d7b241e6Sjeremylt P1 = P2; 258d7b241e6Sjeremylt } 259d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 260d7b241e6Sjeremylt xi = xi-P2/dP2; 261d7b241e6Sjeremylt } 262d7b241e6Sjeremylt // Save xi, wi 263d7b241e6Sjeremylt wi = 2.0/((1.0-xi*xi)*dP2*dP2); 264d7b241e6Sjeremylt qweight1d[i] = wi; 265d7b241e6Sjeremylt qweight1d[Q-1-i] = wi; 266d7b241e6Sjeremylt qref1d[i] = -xi; 267d7b241e6Sjeremylt qref1d[Q-1-i]= xi; 268d7b241e6Sjeremylt } 269d7b241e6Sjeremylt return 0; 270d7b241e6Sjeremylt } 271d7b241e6Sjeremylt 272b11c1e72Sjeremylt /** 273b11c1e72Sjeremylt @brief Construct a Gauss-Legendre-Lobatto quadrature 274b11c1e72Sjeremylt 275b11c1e72Sjeremylt @param Q Number of quadrature points (integrates polynomials of 276b11c1e72Sjeremylt degree 2*Q-3 exactly) 277b11c1e72Sjeremylt @param[out] qref1d Array of length Q to hold the abscissa on [-1, 1] 278b11c1e72Sjeremylt @param[out] qweight1d Array of length Q to hold the weights 279b11c1e72Sjeremylt 280b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 281dfdf5a53Sjeremylt 282dfdf5a53Sjeremylt @ref Utility 283b11c1e72Sjeremylt **/ 284d7b241e6Sjeremylt int CeedLobattoQuadrature(CeedInt Q, CeedScalar *qref1d, 285d7b241e6Sjeremylt CeedScalar *qweight1d) { 286d7b241e6Sjeremylt // Allocate 287d7b241e6Sjeremylt CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0*atan(1.0); 288d7b241e6Sjeremylt // Build qref1d, qweight1d 289d7b241e6Sjeremylt // Set endpoints 290d7b241e6Sjeremylt wi = 2.0/((CeedScalar)(Q*(Q-1))); 291d7b241e6Sjeremylt if (qweight1d) { 292d7b241e6Sjeremylt qweight1d[0] = wi; 293d7b241e6Sjeremylt qweight1d[Q-1] = wi; 294d7b241e6Sjeremylt } 295d7b241e6Sjeremylt qref1d[0] = -1.0; 296d7b241e6Sjeremylt qref1d[Q-1] = 1.0; 297d7b241e6Sjeremylt // Interior 298d7b241e6Sjeremylt for (int i = 1; i <= (Q-1)/2; i++) { 299d7b241e6Sjeremylt // Guess 300d7b241e6Sjeremylt xi = cos(PI*(CeedScalar)(i)/(CeedScalar)(Q-1)); 301d7b241e6Sjeremylt // Pn(xi) 302d7b241e6Sjeremylt P0 = 1.0; 303d7b241e6Sjeremylt P1 = xi; 304d7b241e6Sjeremylt P2 = 0.0; 305d7b241e6Sjeremylt for (int j = 2; j < Q; j++) { 306d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 307d7b241e6Sjeremylt P0 = P1; 308d7b241e6Sjeremylt P1 = P2; 309d7b241e6Sjeremylt } 310d7b241e6Sjeremylt // First Newton step 311d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 312d7b241e6Sjeremylt d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi); 313d7b241e6Sjeremylt xi = xi-dP2/d2P2; 314d7b241e6Sjeremylt // Newton to convergence 315d7b241e6Sjeremylt for (int k=0; k<100 && fabs(dP2)>1e-15; k++) { 316d7b241e6Sjeremylt P0 = 1.0; 317d7b241e6Sjeremylt P1 = xi; 318d7b241e6Sjeremylt for (int j = 2; j < Q; j++) { 319d7b241e6Sjeremylt P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 320d7b241e6Sjeremylt P0 = P1; 321d7b241e6Sjeremylt P1 = P2; 322d7b241e6Sjeremylt } 323d7b241e6Sjeremylt dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 324d7b241e6Sjeremylt d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi); 325d7b241e6Sjeremylt xi = xi-dP2/d2P2; 326d7b241e6Sjeremylt } 327d7b241e6Sjeremylt // Save xi, wi 328d7b241e6Sjeremylt wi = 2.0/(((CeedScalar)(Q*(Q-1)))*P2*P2); 329d7b241e6Sjeremylt if (qweight1d) { 330d7b241e6Sjeremylt qweight1d[i] = wi; 331d7b241e6Sjeremylt qweight1d[Q-1-i] = wi; 332d7b241e6Sjeremylt } 333d7b241e6Sjeremylt qref1d[i] = -xi; 334d7b241e6Sjeremylt qref1d[Q-1-i]= xi; 335d7b241e6Sjeremylt } 336d7b241e6Sjeremylt return 0; 337d7b241e6Sjeremylt } 338d7b241e6Sjeremylt 339dfdf5a53Sjeremylt /** 340dfdf5a53Sjeremylt @brief View an array stored in a CeedBasis 341dfdf5a53Sjeremylt 342dfdf5a53Sjeremylt @param name Name of array 343dfdf5a53Sjeremylt @param fpformat Printing format 344dfdf5a53Sjeremylt @param m Number of rows in array 345dfdf5a53Sjeremylt @param n Number of columns in array 346dfdf5a53Sjeremylt @param a Array to be viewed 347dfdf5a53Sjeremylt @param stream Stream to view to, e.g., stdout 348dfdf5a53Sjeremylt 349dfdf5a53Sjeremylt @return An error code: 0 - success, otherwise - failure 350dfdf5a53Sjeremylt 351dfdf5a53Sjeremylt @ref Utility 352dfdf5a53Sjeremylt **/ 353d7b241e6Sjeremylt static int CeedScalarView(const char *name, const char *fpformat, CeedInt m, 354d7b241e6Sjeremylt CeedInt n, const CeedScalar *a, FILE *stream) { 355d7b241e6Sjeremylt for (int i=0; i<m; i++) { 356d7b241e6Sjeremylt if (m > 1) fprintf(stream, "%12s[%d]:", name, i); 357d7b241e6Sjeremylt else fprintf(stream, "%12s:", name); 358d7b241e6Sjeremylt for (int j=0; j<n; j++) { 359d7b241e6Sjeremylt fprintf(stream, fpformat, fabs(a[i*n+j]) > 1E-14 ? a[i*n+j] : 0); 360d7b241e6Sjeremylt } 361d7b241e6Sjeremylt fputs("\n", stream); 362d7b241e6Sjeremylt } 363d7b241e6Sjeremylt return 0; 364d7b241e6Sjeremylt } 365d7b241e6Sjeremylt 366b11c1e72Sjeremylt /** 367b11c1e72Sjeremylt @brief View a CeedBasis 368b11c1e72Sjeremylt 369b11c1e72Sjeremylt @param basis CeedBasis to view 370b11c1e72Sjeremylt @param stream Stream to view to, e.g., stdout 371b11c1e72Sjeremylt 372b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 373dfdf5a53Sjeremylt 374dfdf5a53Sjeremylt @ref Utility 375b11c1e72Sjeremylt **/ 376d7b241e6Sjeremylt int CeedBasisView(CeedBasis basis, FILE *stream) { 377d7b241e6Sjeremylt int ierr; 378d7b241e6Sjeremylt 379*a8de75f0Sjeremylt if (basis->tensorbasis) { 380d7b241e6Sjeremylt fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P1d, 381d7b241e6Sjeremylt basis->Q1d); 382d7b241e6Sjeremylt ierr = CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q1d, basis->qref1d, 383d7b241e6Sjeremylt stream); CeedChk(ierr); 384d7b241e6Sjeremylt ierr = CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q1d, basis->qweight1d, 385d7b241e6Sjeremylt stream); CeedChk(ierr); 386d7b241e6Sjeremylt ierr = CeedScalarView("interp1d", "\t% 12.8f", basis->Q1d, basis->P1d, 387d7b241e6Sjeremylt basis->interp1d, stream); CeedChk(ierr); 388d7b241e6Sjeremylt ierr = CeedScalarView("grad1d", "\t% 12.8f", basis->Q1d, basis->P1d, 389d7b241e6Sjeremylt basis->grad1d, stream); CeedChk(ierr); 390*a8de75f0Sjeremylt } else { 391*a8de75f0Sjeremylt fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P, 392*a8de75f0Sjeremylt basis->Q); 393*a8de75f0Sjeremylt ierr = CeedScalarView("qref", "\t% 12.8f", 1, basis->Q*basis->dim, 394*a8de75f0Sjeremylt basis->qref1d, 395*a8de75f0Sjeremylt stream); CeedChk(ierr); 396*a8de75f0Sjeremylt ierr = CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->qweight1d, 397*a8de75f0Sjeremylt stream); CeedChk(ierr); 398*a8de75f0Sjeremylt ierr = CeedScalarView("interp", "\t% 12.8f", basis->Q, basis->P, 399*a8de75f0Sjeremylt basis->interp1d, stream); CeedChk(ierr); 400*a8de75f0Sjeremylt ierr = CeedScalarView("grad", "\t% 12.8f", basis->dim*basis->Q, basis->P, 401*a8de75f0Sjeremylt basis->grad1d, stream); CeedChk(ierr); 402*a8de75f0Sjeremylt } 403d7b241e6Sjeremylt return 0; 404d7b241e6Sjeremylt } 405d7b241e6Sjeremylt 406dfdf5a53Sjeremylt /** 407dfdf5a53Sjeremylt @brief Compute Householder Reflection 408dfdf5a53Sjeremylt 409dfdf5a53Sjeremylt Computes A = (I - b v v^T) A 410dfdf5a53Sjeremylt where A is an mxn matrix indexed as A[i*row + j*col] 411dfdf5a53Sjeremylt 412dfdf5a53Sjeremylt @param[out] A Matrix to apply Householder reflection to, in place 413dfdf5a53Sjeremylt @param v Householder vector 414dfdf5a53Sjeremylt @param b Scaling factor 415dfdf5a53Sjeremylt @param m Number of rows in A 416dfdf5a53Sjeremylt @param n Number of columns in A 417dfdf5a53Sjeremylt @param row Col stride 418dfdf5a53Sjeremylt @param col Row stride 419dfdf5a53Sjeremylt 420dfdf5a53Sjeremylt @return An error code: 0 - success, otherwise - failure 421dfdf5a53Sjeremylt 422dfdf5a53Sjeremylt @ref Developer 423dfdf5a53Sjeremylt **/ 424d7b241e6Sjeremylt static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v, 425d7b241e6Sjeremylt CeedScalar b, CeedInt m, CeedInt n, 426d7b241e6Sjeremylt CeedInt row, CeedInt col) { 427d7b241e6Sjeremylt for (CeedInt j=0; j<n; j++) { 428d7b241e6Sjeremylt CeedScalar w = A[0*row + j*col]; 429d7b241e6Sjeremylt for (CeedInt i=1; i<m; i++) w += v[i] * A[i*row + j*col]; 430d7b241e6Sjeremylt A[0*row + j*col] -= b * w; 431d7b241e6Sjeremylt for (CeedInt i=1; i<m; i++) A[i*row + j*col] -= b * w * v[i]; 432d7b241e6Sjeremylt } 433d7b241e6Sjeremylt return 0; 434d7b241e6Sjeremylt } 435d7b241e6Sjeremylt 436dfdf5a53Sjeremylt /** 437dfdf5a53Sjeremylt @brief Apply Householder Q matrix 438dfdf5a53Sjeremylt 439dfdf5a53Sjeremylt Compute A = Q A where Q is mxk and A is mxn. k<m 440dfdf5a53Sjeremylt 441dfdf5a53Sjeremylt @param[out] A Matrix to apply Householder Q to, in place 442dfdf5a53Sjeremylt @param Q Householder Q matrix 443dfdf5a53Sjeremylt @param tau Householder scaling factors 444dfdf5a53Sjeremylt @param tmode Transpose mode for application 445dfdf5a53Sjeremylt @param m Number of rows in A 446dfdf5a53Sjeremylt @param n Number of columns in A 447dfdf5a53Sjeremylt @param k Index of row targeted 448dfdf5a53Sjeremylt @param row Col stride 449dfdf5a53Sjeremylt @param col Row stride 450dfdf5a53Sjeremylt 451dfdf5a53Sjeremylt @return An error code: 0 - success, otherwise - failure 452dfdf5a53Sjeremylt 453dfdf5a53Sjeremylt @ref Developer 454dfdf5a53Sjeremylt **/ 455d7b241e6Sjeremylt static int CeedHouseholderApplyQ(CeedScalar *A, const CeedScalar *Q, 456d7b241e6Sjeremylt const CeedScalar *tau, CeedTransposeMode tmode, 457d7b241e6Sjeremylt CeedInt m, CeedInt n, CeedInt k, 458d7b241e6Sjeremylt CeedInt row, CeedInt col) { 459d7b241e6Sjeremylt CeedScalar v[m]; 460d7b241e6Sjeremylt for (CeedInt ii=0; ii<k; ii++) { 461d7b241e6Sjeremylt CeedInt i = tmode == CEED_TRANSPOSE ? ii : k-1-ii; 462d7b241e6Sjeremylt for (CeedInt j=i+1; j<m; j++) { 463d7b241e6Sjeremylt v[j] = Q[j*k+i]; 464d7b241e6Sjeremylt } 465d7b241e6Sjeremylt // Apply Householder reflector (I - tau v v^T) colograd1d^T 466d7b241e6Sjeremylt CeedHouseholderReflect(&A[i*row], &v[i], tau[i], m-i, n, row, col); 467d7b241e6Sjeremylt } 468d7b241e6Sjeremylt return 0; 469d7b241e6Sjeremylt } 470d7b241e6Sjeremylt 471b11c1e72Sjeremylt /** 472b11c1e72Sjeremylt @brief Return QR Factorization of matrix 473b11c1e72Sjeremylt 474b11c1e72Sjeremylt @param[out] mat Row-major matrix to be factorized in place 475b11c1e72Sjeremylt @param[out] tau Vector of length m of scaling fators 476b11c1e72Sjeremylt @param m Number of rows 477b11c1e72Sjeremylt @param n Number of columns 478b11c1e72Sjeremylt 479b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 480dfdf5a53Sjeremylt 481dfdf5a53Sjeremylt @ref Utility 482b11c1e72Sjeremylt **/ 483d7b241e6Sjeremylt int CeedQRFactorization(CeedScalar *mat, CeedScalar *tau, 484d7b241e6Sjeremylt CeedInt m, CeedInt n) { 485d7b241e6Sjeremylt CeedInt i, j; 486d7b241e6Sjeremylt CeedScalar v[m]; 487d7b241e6Sjeremylt 488d7b241e6Sjeremylt for (i=0; i<n; i++) { 489d7b241e6Sjeremylt // Calculate Householder vector, magnitude 490d7b241e6Sjeremylt CeedScalar sigma = 0.0; 491d7b241e6Sjeremylt v[i] = mat[i+n*i]; 492d7b241e6Sjeremylt for (j=i+1; j<m; j++) { 493d7b241e6Sjeremylt v[j] = mat[i+n*j]; 494d7b241e6Sjeremylt sigma += v[j] * v[j]; 495d7b241e6Sjeremylt } 496d7b241e6Sjeremylt CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m] 497d7b241e6Sjeremylt CeedScalar Rii = -copysign(norm, v[i]); 498d7b241e6Sjeremylt v[i] -= Rii; 499d7b241e6Sjeremylt // norm of v[i:m] after modification above and scaling below 500d7b241e6Sjeremylt // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 501d7b241e6Sjeremylt // tau = 2 / (norm*norm) 502d7b241e6Sjeremylt tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma); 503d7b241e6Sjeremylt for (j=i+1; j<m; j++) v[j] /= v[i]; 504d7b241e6Sjeremylt 505d7b241e6Sjeremylt // Apply Householder reflector to lower right panel 506d7b241e6Sjeremylt CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1); 507d7b241e6Sjeremylt // Save v 508d7b241e6Sjeremylt mat[i+n*i] = Rii; 509d7b241e6Sjeremylt for (j=i+1; j<m; j++) { 510d7b241e6Sjeremylt mat[i+n*j] = v[j]; 511d7b241e6Sjeremylt } 512d7b241e6Sjeremylt } 513d7b241e6Sjeremylt 514d7b241e6Sjeremylt return 0; 515d7b241e6Sjeremylt } 516d7b241e6Sjeremylt 517b11c1e72Sjeremylt /** 518783c99b3SValeria Barra @brief Return collocated grad matrix 519b11c1e72Sjeremylt 520b11c1e72Sjeremylt @param basis CeedBasis 521b11c1e72Sjeremylt @param[out] colograd1d Row-major Q1d × Q1d matrix expressing derivatives of 522b11c1e72Sjeremylt basis functions at quadrature points 523b11c1e72Sjeremylt 524b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 525dfdf5a53Sjeremylt 526dfdf5a53Sjeremylt @ref Advanced 527b11c1e72Sjeremylt **/ 528783c99b3SValeria Barra int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *colograd1d) { 529d7b241e6Sjeremylt int i, j, k; 530d7b241e6Sjeremylt CeedInt ierr, P1d=(basis)->P1d, Q1d=(basis)->Q1d; 531d7b241e6Sjeremylt CeedScalar *interp1d, *grad1d, tau[Q1d]; 532d7b241e6Sjeremylt 533d7b241e6Sjeremylt ierr = CeedMalloc(Q1d*P1d, &interp1d); CeedChk(ierr); 534d7b241e6Sjeremylt ierr = CeedMalloc(Q1d*P1d, &grad1d); CeedChk(ierr); 535d7b241e6Sjeremylt memcpy(interp1d, (basis)->interp1d, Q1d*P1d*sizeof(basis)->interp1d[0]); 536d7b241e6Sjeremylt memcpy(grad1d, (basis)->grad1d, Q1d*P1d*sizeof(basis)->interp1d[0]); 537d7b241e6Sjeremylt 538d7b241e6Sjeremylt // QR Factorization, interp1d = Q R 539d7b241e6Sjeremylt ierr = CeedQRFactorization(interp1d, tau, Q1d, P1d); CeedChk(ierr); 540d7b241e6Sjeremylt 541d7b241e6Sjeremylt // Apply Rinv, colograd1d = grad1d Rinv 542d7b241e6Sjeremylt for (i=0; i<Q1d; i++) { // Row i 543d7b241e6Sjeremylt colograd1d[Q1d*i] = grad1d[P1d*i]/interp1d[0]; 544d7b241e6Sjeremylt for (j=1; j<P1d; j++) { // Column j 545d7b241e6Sjeremylt colograd1d[j+Q1d*i] = grad1d[j+P1d*i]; 546d7b241e6Sjeremylt for (k=0; k<j; k++) { 547d7b241e6Sjeremylt colograd1d[j+Q1d*i] -= interp1d[j+P1d*k]*colograd1d[k+Q1d*i]; 548d7b241e6Sjeremylt } 549d7b241e6Sjeremylt colograd1d[j+Q1d*i] /= interp1d[j+P1d*j]; 550d7b241e6Sjeremylt } 551d7b241e6Sjeremylt for (j=P1d; j<Q1d; j++) { 552d7b241e6Sjeremylt colograd1d[j+Q1d*i] = 0; 553d7b241e6Sjeremylt } 554d7b241e6Sjeremylt } 555d7b241e6Sjeremylt 556d7b241e6Sjeremylt // Apply Qtranspose, colograd = colograd Qtranspose 557d7b241e6Sjeremylt CeedHouseholderApplyQ(colograd1d, interp1d, tau, CEED_NOTRANSPOSE, 558d7b241e6Sjeremylt Q1d, Q1d, P1d, 1, Q1d); 559d7b241e6Sjeremylt 560d7b241e6Sjeremylt ierr = CeedFree(&interp1d); CeedChk(ierr); 561d7b241e6Sjeremylt ierr = CeedFree(&grad1d); CeedChk(ierr); 562d7b241e6Sjeremylt 563d7b241e6Sjeremylt return 0; 564d7b241e6Sjeremylt } 565d7b241e6Sjeremylt 566b11c1e72Sjeremylt /** 567b11c1e72Sjeremylt @brief Apply basis evaluation from nodes to quadrature points or vice-versa 568b11c1e72Sjeremylt 569b11c1e72Sjeremylt @param basis CeedBasis to evaluate 570b11c1e72Sjeremylt @param nelem The number of elements to apply the basis evaluation to; 571b11c1e72Sjeremylt the backend will specify the ordering in 572b11c1e72Sjeremylt ElemRestrictionCreateBlocked 573b11c1e72Sjeremylt @param tmode \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature 574b11c1e72Sjeremylt points, \ref CEED_TRANSPOSE to apply the transpose, mapping 575b11c1e72Sjeremylt from quadrature points to nodes 576b11c1e72Sjeremylt @param emode \ref CEED_EVAL_INTERP to obtain interpolated values, 577b11c1e72Sjeremylt \ref CEED_EVAL_GRAD to obtain gradients. 578b11c1e72Sjeremylt @param[in] u Input array 579b11c1e72Sjeremylt @param[out] v Output array 580b11c1e72Sjeremylt 581b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 582dfdf5a53Sjeremylt 583dfdf5a53Sjeremylt @ref Advanced 584b11c1e72Sjeremylt **/ 585d7b241e6Sjeremylt int CeedBasisApply(CeedBasis basis, CeedInt nelem, CeedTransposeMode tmode, 586d7b241e6Sjeremylt CeedEvalMode emode, const CeedScalar *u, CeedScalar *v) { 587d7b241e6Sjeremylt int ierr; 588d7b241e6Sjeremylt if (!basis->Apply) return CeedError(basis->ceed, 1, 589d7b241e6Sjeremylt "Backend does not support BasisApply"); 590d7b241e6Sjeremylt ierr = basis->Apply(basis, nelem, tmode, emode, u, v); CeedChk(ierr); 591d7b241e6Sjeremylt return 0; 592d7b241e6Sjeremylt } 593d7b241e6Sjeremylt 594b11c1e72Sjeremylt /** 595b11c1e72Sjeremylt @brief Get total number of nodes (in dim dimensions) 596b11c1e72Sjeremylt 597b11c1e72Sjeremylt @param basis CeedBasis 598b11c1e72Sjeremylt @param[out] P Number of nodes 599b11c1e72Sjeremylt 600b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 601dfdf5a53Sjeremylt 602dfdf5a53Sjeremylt @ref Utility 603b11c1e72Sjeremylt **/ 604d7b241e6Sjeremylt int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) { 605*a8de75f0Sjeremylt *P = basis->P; 606d7b241e6Sjeremylt return 0; 607d7b241e6Sjeremylt } 608d7b241e6Sjeremylt 609b11c1e72Sjeremylt /** 610b11c1e72Sjeremylt @brief Get total number of quadrature points (in dim dimensions) 611b11c1e72Sjeremylt 612b11c1e72Sjeremylt @param basis CeedBasis 613b11c1e72Sjeremylt @param[out] Q Number of quadrature points 614b11c1e72Sjeremylt 615b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 616dfdf5a53Sjeremylt 617dfdf5a53Sjeremylt @ref Utility 618b11c1e72Sjeremylt **/ 619d7b241e6Sjeremylt int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) { 620*a8de75f0Sjeremylt *Q = basis->Q; 621d7b241e6Sjeremylt return 0; 622d7b241e6Sjeremylt } 623d7b241e6Sjeremylt 624b11c1e72Sjeremylt /** 625*a8de75f0Sjeremylt @brief Get dimension for given CeedElemTopology 626*a8de75f0Sjeremylt 627*a8de75f0Sjeremylt @param topo CeedElemTopology 628*a8de75f0Sjeremylt @param[out] dim Dimension of topology 629*a8de75f0Sjeremylt 630*a8de75f0Sjeremylt @return An error code: 0 - success, otherwise - failure 631*a8de75f0Sjeremylt 632*a8de75f0Sjeremylt @ref Utility 633*a8de75f0Sjeremylt **/ 634*a8de75f0Sjeremylt int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) { 635*a8de75f0Sjeremylt *dim = (CeedInt) topo >> 16; 636*a8de75f0Sjeremylt 637*a8de75f0Sjeremylt return 0; 638*a8de75f0Sjeremylt }; 639*a8de75f0Sjeremylt 640*a8de75f0Sjeremylt /** 641b11c1e72Sjeremylt @brief Destroy a CeedBasis 642b11c1e72Sjeremylt 643b11c1e72Sjeremylt @param basis CeedBasis to destroy 644b11c1e72Sjeremylt 645b11c1e72Sjeremylt @return An error code: 0 - success, otherwise - failure 646dfdf5a53Sjeremylt 647dfdf5a53Sjeremylt @ref Basic 648b11c1e72Sjeremylt **/ 649d7b241e6Sjeremylt int CeedBasisDestroy(CeedBasis *basis) { 650d7b241e6Sjeremylt int ierr; 651d7b241e6Sjeremylt 652d7b241e6Sjeremylt if (!*basis || --(*basis)->refcount > 0) return 0; 653d7b241e6Sjeremylt if ((*basis)->Destroy) { 654d7b241e6Sjeremylt ierr = (*basis)->Destroy(*basis); CeedChk(ierr); 655d7b241e6Sjeremylt } 656d7b241e6Sjeremylt ierr = CeedFree(&(*basis)->interp1d); CeedChk(ierr); 657d7b241e6Sjeremylt ierr = CeedFree(&(*basis)->grad1d); CeedChk(ierr); 658d7b241e6Sjeremylt ierr = CeedFree(&(*basis)->qref1d); CeedChk(ierr); 659d7b241e6Sjeremylt ierr = CeedFree(&(*basis)->qweight1d); CeedChk(ierr); 660d7b241e6Sjeremylt ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr); 661d7b241e6Sjeremylt ierr = CeedFree(basis); CeedChk(ierr); 662d7b241e6Sjeremylt return 0; 663d7b241e6Sjeremylt } 664d7b241e6Sjeremylt 66533e6becaSjeremylt /// @cond DOXYGEN_SKIP 666783c99b3SValeria Barra // Indicate that the quadrature points are collocated with the dofs 667783c99b3SValeria Barra CeedBasis CEED_BASIS_COLLOCATED = &ceed_basis_collocated; 66833e6becaSjeremylt /// @endcond 669d7b241e6Sjeremylt /// @} 670