1 // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 2 // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 3 // reserved. See files LICENSE and NOTICE for details. 4 // 5 // This file is part of CEED, a collection of benchmarks, miniapps, software 6 // libraries and APIs for efficient high-order finite element and spectral 7 // element discretizations for exascale applications. For more information and 8 // source code availability see http://github.com/ceed. 9 // 10 // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 11 // a collaborative effort of two U.S. Department of Energy organizations (Office 12 // of Science and the National Nuclear Security Administration) responsible for 13 // the planning and preparation of a capable exascale ecosystem, including 14 // software, applications, hardware, advanced system engineering and early 15 // testbed platforms, in support of the nation's exascale computing imperative. 16 17 #include <ceed-impl.h> 18 #include <math.h> 19 #include <stdio.h> 20 #include <stdlib.h> 21 #include <string.h> 22 23 /// @cond DOXYGEN_SKIP 24 static struct CeedBasis_private ceed_basis_colocated; 25 /// @endcond 26 27 /// @file 28 /// Implementation of public CeedBasis interfaces 29 /// 30 /// @defgroup CeedBasis CeedBasis: fully discrete finite element-like objects 31 /// @{ 32 33 /// Create a tensor product basis for H^1 discretizations 34 /// 35 /// @param ceed Ceed 36 /// @param dim Topological dimension 37 /// @param ncomp Number of field components (1 for scalar fields) 38 /// @param P1d Number of nodes in one dimension 39 /// @param Q1d Number of quadrature points in one dimension 40 /// @param interp1d Row-major Q1d × P1d matrix expressing the values of nodal 41 /// basis functions at quadrature points 42 /// @param grad1d Row-major Q1d × P1d matrix expressing derivatives of nodal 43 /// basis functions at quadrature points 44 /// @param qref1d Array of length Q1d holding the locations of quadrature points 45 /// on the 1D reference element [-1, 1] 46 /// @param qweight1d Array of length Q1d holding the quadrature weights on the 47 /// reference element 48 /// @param[out] basis New basis 49 /// 50 /// @sa CeedBasisCreateTensorH1Lagrange() 51 int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt ncomp, CeedInt P1d, 52 CeedInt Q1d, const CeedScalar *interp1d, 53 const CeedScalar *grad1d, const CeedScalar *qref1d, 54 const CeedScalar *qweight1d, CeedBasis *basis) { 55 int ierr; 56 57 if (!ceed->BasisCreateTensorH1) 58 return CeedError(ceed, 1, "Backend does not support BasisCreateTensorH1"); 59 ierr = CeedCalloc(1,basis); CeedChk(ierr); 60 (*basis)->ceed = ceed; 61 ceed->refcount++; 62 (*basis)->refcount = 1; 63 (*basis)->dim = dim; 64 (*basis)->ncomp = ncomp; 65 (*basis)->P1d = P1d; 66 (*basis)->Q1d = Q1d; 67 ierr = CeedMalloc(Q1d,&(*basis)->qref1d); CeedChk(ierr); 68 ierr = CeedMalloc(Q1d,&(*basis)->qweight1d); CeedChk(ierr); 69 memcpy((*basis)->qref1d, qref1d, Q1d*sizeof(qref1d[0])); 70 memcpy((*basis)->qweight1d, qweight1d, Q1d*sizeof(qweight1d[0])); 71 ierr = CeedMalloc(Q1d*P1d,&(*basis)->interp1d); CeedChk(ierr); 72 ierr = CeedMalloc(Q1d*P1d,&(*basis)->grad1d); CeedChk(ierr); 73 memcpy((*basis)->interp1d, interp1d, Q1d*P1d*sizeof(interp1d[0])); 74 memcpy((*basis)->grad1d, grad1d, Q1d*P1d*sizeof(interp1d[0])); 75 ierr = ceed->BasisCreateTensorH1(ceed, dim, P1d, Q1d, interp1d, grad1d, qref1d, 76 qweight1d, *basis); CeedChk(ierr); 77 return 0; 78 } 79 80 /// Create a tensor product Lagrange basis 81 /// 82 /// @param ceed Ceed 83 /// @param dim Topological dimension of element 84 /// @param ncomp Number of field components 85 /// @param P Number of Gauss-Lobatto nodes in one dimension. The polynomial degree 86 /// of the resulting Q_k element is k=P-1. 87 /// @param Q Number of quadrature points in one dimension. 88 /// @param qmode Distribution of the Q quadrature points (affects order of 89 /// accuracy for the quadrature) 90 /// @param[out] basis New basis 91 /// 92 /// @sa CeedBasisCreateTensorH1() 93 int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt ncomp, 94 CeedInt P, CeedInt Q, 95 CeedQuadMode qmode, CeedBasis *basis) { 96 // Allocate 97 int ierr, i, j, k; 98 CeedScalar c1, c2, c3, c4, dx, *nodes, *interp1d, *grad1d, *qref1d, *qweight1d; 99 ierr = CeedCalloc(P*Q, &interp1d); CeedChk(ierr); 100 ierr = CeedCalloc(P*Q, &grad1d); CeedChk(ierr); 101 ierr = CeedCalloc(P, &nodes); CeedChk(ierr); 102 ierr = CeedCalloc(Q, &qref1d); CeedChk(ierr); 103 ierr = CeedCalloc(Q, &qweight1d); CeedChk(ierr); 104 // Get Nodes and Weights 105 ierr = CeedLobattoQuadrature(P, nodes, NULL); CeedChk(ierr); 106 switch (qmode) { 107 case CEED_GAUSS: 108 ierr = CeedGaussQuadrature(Q, qref1d, qweight1d); CeedChk(ierr); 109 break; 110 case CEED_GAUSS_LOBATTO: 111 ierr = CeedLobattoQuadrature(Q, qref1d, qweight1d); CeedChk(ierr); 112 break; 113 } 114 // Build B, D matrix 115 // Fornberg, 1998 116 for (i = 0; i < Q; i++) { 117 c1 = 1.0; 118 c3 = nodes[0] - qref1d[i]; 119 interp1d[i*P+0] = 1.0; 120 for (j = 1; j < P; j++) { 121 c2 = 1.0; 122 c4 = c3; 123 c3 = nodes[j] - qref1d[i]; 124 for (k = 0; k < j; k++) { 125 dx = nodes[j] - nodes[k]; 126 c2 *= dx; 127 if (k == j - 1) { 128 grad1d[i*P + j] = c1*(interp1d[i*P + k] - c4*grad1d[i*P + k]) / c2; 129 interp1d[i*P + j] = - c1*c4*interp1d[i*P + k] / c2; 130 } 131 grad1d[i*P + k] = (c3*grad1d[i*P + k] - interp1d[i*P + k]) / dx; 132 interp1d[i*P + k] = c3*interp1d[i*P + k] / dx; 133 } 134 c1 = c2; 135 } 136 } 137 // // Pass to CeedBasisCreateTensorH1 138 ierr = CeedBasisCreateTensorH1(ceed, dim, ncomp, P, Q, interp1d, grad1d, qref1d, 139 qweight1d, basis); CeedChk(ierr); 140 ierr = CeedFree(&interp1d); CeedChk(ierr); 141 ierr = CeedFree(&grad1d); CeedChk(ierr); 142 ierr = CeedFree(&nodes); CeedChk(ierr); 143 ierr = CeedFree(&qref1d); CeedChk(ierr); 144 ierr = CeedFree(&qweight1d); CeedChk(ierr); 145 return 0; 146 } 147 148 /// Construct a Gauss-Legendre quadrature 149 /// 150 /// @param Q Number of quadrature points (integrates polynomials of degree 2*Q-1 exactly) 151 /// @param qref1d Array of length Q to hold the abscissa on [-1, 1] 152 /// @param qweight1d Array of length Q to hold the weights 153 /// @sa CeedLobattoQuadrature() 154 int CeedGaussQuadrature(CeedInt Q, CeedScalar *qref1d, CeedScalar *qweight1d) { 155 // Allocate 156 CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0*atan(1.0); 157 // Build qref1d, qweight1d 158 for (int i = 0; i <= Q/2; i++) { 159 // Guess 160 xi = cos(PI*(CeedScalar)(2*i+1)/((CeedScalar)(2*Q))); 161 // Pn(xi) 162 P0 = 1.0; 163 P1 = xi; 164 P2 = 0.0; 165 for (int j = 2; j <= Q; j++) { 166 P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 167 P0 = P1; 168 P1 = P2; 169 } 170 // First Newton Step 171 dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 172 xi = xi-P2/dP2; 173 // Newton to convergence 174 for (int k=0; k<100 && fabs(P2)>1e-15; k++) { 175 P0 = 1.0; 176 P1 = xi; 177 for (int j = 2; j <= Q; j++) { 178 P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 179 P0 = P1; 180 P1 = P2; 181 } 182 dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 183 xi = xi-P2/dP2; 184 } 185 // Save xi, wi 186 wi = 2.0/((1.0-xi*xi)*dP2*dP2); 187 qweight1d[i] = wi; 188 qweight1d[Q-1-i] = wi; 189 qref1d[i] = -xi; 190 qref1d[Q-1-i]= xi; 191 } 192 return 0; 193 } 194 195 /// Construct a Gauss-Legendre-Lobatto quadrature 196 /// 197 /// @param Q Number of quadrature points (integrates polynomials of degree 2*Q-3 exactly) 198 /// @param qref1d Array of length Q to hold the abscissa on [-1, 1] 199 /// @param qweight1d Array of length Q to hold the weights 200 /// @sa CeedGaussQuadrature() 201 int CeedLobattoQuadrature(CeedInt Q, CeedScalar *qref1d, 202 CeedScalar *qweight1d) { 203 // Allocate 204 CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0*atan(1.0); 205 // Build qref1d, qweight1d 206 // Set endpoints 207 wi = 2.0/((CeedScalar)(Q*(Q-1))); 208 if (qweight1d) { 209 qweight1d[0] = wi; 210 qweight1d[Q-1] = wi; 211 } 212 qref1d[0] = -1.0; 213 qref1d[Q-1] = 1.0; 214 // Interior 215 for (int i = 1; i <= (Q-1)/2; i++) { 216 // Guess 217 xi = cos(PI*(CeedScalar)(i)/(CeedScalar)(Q-1)); 218 // Pn(xi) 219 P0 = 1.0; 220 P1 = xi; 221 P2 = 0.0; 222 for (int j = 2; j < Q; j++) { 223 P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 224 P0 = P1; 225 P1 = P2; 226 } 227 // First Newton step 228 dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 229 d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi); 230 xi = xi-dP2/d2P2; 231 // Newton to convergence 232 for (int k=0; k<100 && fabs(dP2)>1e-15; k++) { 233 P0 = 1.0; 234 P1 = xi; 235 for (int j = 2; j < Q; j++) { 236 P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 237 P0 = P1; 238 P1 = P2; 239 } 240 dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 241 d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi); 242 xi = xi-dP2/d2P2; 243 } 244 // Save xi, wi 245 wi = 2.0/(((CeedScalar)(Q*(Q-1)))*P2*P2); 246 if (qweight1d) { 247 qweight1d[i] = wi; 248 qweight1d[Q-1-i] = wi; 249 } 250 qref1d[i] = -xi; 251 qref1d[Q-1-i]= xi; 252 } 253 return 0; 254 } 255 256 static int CeedScalarView(const char *name, const char *fpformat, CeedInt m, 257 CeedInt n, const CeedScalar *a, FILE *stream) { 258 for (int i=0; i<m; i++) { 259 if (m > 1) fprintf(stream, "%12s[%d]:", name, i); 260 else fprintf(stream, "%12s:", name); 261 for (int j=0; j<n; j++) { 262 fprintf(stream, fpformat, fabs(a[i*n+j]) > 1E-14 ? a[i*n+j] : 0); 263 } 264 fputs("\n", stream); 265 } 266 return 0; 267 } 268 269 /// View a basis 270 /// 271 /// @param basis Basis to view 272 /// @param stream Stream to view to, e.g., stdout 273 int CeedBasisView(CeedBasis basis, FILE *stream) { 274 int ierr; 275 276 fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P1d, 277 basis->Q1d); 278 ierr = CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q1d, basis->qref1d, 279 stream); CeedChk(ierr); 280 ierr = CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q1d, basis->qweight1d, 281 stream); CeedChk(ierr); 282 ierr = CeedScalarView("interp1d", "\t% 12.8f", basis->Q1d, basis->P1d, 283 basis->interp1d, stream); CeedChk(ierr); 284 ierr = CeedScalarView("grad1d", "\t% 12.8f", basis->Q1d, basis->P1d, 285 basis->grad1d, stream); CeedChk(ierr); 286 return 0; 287 } 288 289 // Computes A = (I - b v v^T) A 290 // where A is an mxn matrix indexed as A[i*row + j*col] 291 static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v, 292 CeedScalar b, CeedInt m, CeedInt n, 293 CeedInt row, CeedInt col) { 294 for (CeedInt j=0; j<n; j++) { 295 CeedScalar w = A[0*row + j*col]; 296 for (CeedInt i=1; i<m; i++) w += v[i] * A[i*row + j*col]; 297 A[0*row + j*col] -= b * w; 298 for (CeedInt i=1; i<m; i++) A[i*row + j*col] -= b * w * v[i]; 299 } 300 return 0; 301 } 302 303 // Compute A = Q A where Q is mxk and A is mxn. k<m 304 static int CeedHouseholderApplyQ(CeedScalar *A, const CeedScalar *Q, 305 const CeedScalar *tau, CeedTransposeMode tmode, 306 CeedInt m, CeedInt n, CeedInt k, 307 CeedInt row, CeedInt col) { 308 CeedScalar v[m]; 309 for (CeedInt ii=0; ii<k; ii++) { 310 CeedInt i = tmode == CEED_TRANSPOSE ? ii : k-1-ii; 311 for (CeedInt j=i+1; j<m; j++) { 312 v[j] = Q[j*k+i]; 313 } 314 // Apply Householder reflector (I - tau v v^T) colograd1d^T 315 CeedHouseholderReflect(&A[i*row], &v[i], tau[i], m-i, n, row, col); 316 } 317 return 0; 318 } 319 320 /// Return QR Factorization of matrix 321 /// @param mat Row-major matrix to be factorized in place 322 /// @param tau Vector of length m of scaling fators 323 /// @param m Number of rows 324 /// @param n Number of columns 325 int CeedQRFactorization(CeedScalar *mat, CeedScalar *tau, 326 CeedInt m, CeedInt n) { 327 CeedInt i, j; 328 CeedScalar v[m]; 329 330 for (i=0; i<n; i++) { 331 // Calculate Householder vector, magnitude 332 CeedScalar sigma = 0.0; 333 v[i] = mat[i+n*i]; 334 for (j=i+1; j<m; j++) { 335 v[j] = mat[i+n*j]; 336 sigma += v[j] * v[j]; 337 } 338 CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m] 339 CeedScalar Rii = -copysign(norm, v[i]); 340 v[i] -= Rii; 341 // norm of v[i:m] after modification above and scaling below 342 // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 343 // tau = 2 / (norm*norm) 344 tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma); 345 for (j=i+1; j<m; j++) v[j] /= v[i]; 346 347 // Apply Householder reflector to lower right panel 348 CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1); 349 // Save v 350 mat[i+n*i] = Rii; 351 for (j=i+1; j<m; j++) { 352 mat[i+n*j] = v[j]; 353 } 354 } 355 356 return 0; 357 } 358 359 /// Return colocated grad matrix 360 /// @param basis Basis 361 /// @param colograd1d Row-major Q1d × Q1d matrix expressing derivatives of 362 /// basis functions at quadrature points 363 int CeedBasisGetColocatedGrad(CeedBasis basis, CeedScalar *colograd1d) { 364 int i, j, k; 365 CeedInt ierr, P1d=(basis)->P1d, Q1d=(basis)->Q1d; 366 CeedScalar *interp1d, *grad1d, tau[Q1d]; 367 368 ierr = CeedMalloc(Q1d*P1d, &interp1d); CeedChk(ierr); 369 ierr = CeedMalloc(Q1d*P1d, &grad1d); CeedChk(ierr); 370 memcpy(interp1d, (basis)->interp1d, Q1d*P1d*sizeof(basis)->interp1d[0]); 371 memcpy(grad1d, (basis)->grad1d, Q1d*P1d*sizeof(basis)->interp1d[0]); 372 373 // QR Factorization, interp1d = Q R 374 ierr = CeedQRFactorization(interp1d, tau, Q1d, P1d); CeedChk(ierr); 375 376 // Apply Rinv, colograd1d = grad1d Rinv 377 for (i=0; i<Q1d; i++) { // Row i 378 colograd1d[Q1d*i] = grad1d[P1d*i]/interp1d[0]; 379 for (j=1; j<P1d; j++) { // Column j 380 colograd1d[j+Q1d*i] = grad1d[j+P1d*i]; 381 for (k=0; k<j; k++) { 382 colograd1d[j+Q1d*i] -= interp1d[j+P1d*k]*colograd1d[k+Q1d*i]; 383 } 384 colograd1d[j+Q1d*i] /= interp1d[j+P1d*j]; 385 } 386 for (j=P1d; j<Q1d; j++) { 387 colograd1d[j+Q1d*i] = 0; 388 } 389 } 390 391 // Apply Qtranspose, colograd = colograd Qtranspose 392 CeedHouseholderApplyQ(colograd1d, interp1d, tau, CEED_NOTRANSPOSE, 393 Q1d, Q1d, P1d, 1, Q1d); 394 395 ierr = CeedFree(&interp1d); CeedChk(ierr); 396 ierr = CeedFree(&grad1d); CeedChk(ierr); 397 398 return 0; 399 } 400 401 /// Apply basis evaluation from nodes to quadrature points or vice-versa 402 /// 403 /// @param basis Basis to evaluate 404 /// @param nelem the number of elements to apply the basis evaluation to; 405 /// the backend will specify the ordering in ElemRestrictionCreateBlocked 406 /// @param tmode \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature 407 /// points, \ref CEED_TRANSPOSE to apply the transpose, mapping from 408 /// quadrature points to nodes 409 /// @param emode \ref CEED_EVAL_INTERP to obtain interpolated values, 410 /// \ref CEED_EVAL_GRAD to obtain gradients. 411 /// @param u input vector 412 /// @param v output vector 413 int CeedBasisApply(CeedBasis basis, CeedInt nelem, CeedTransposeMode tmode, 414 CeedEvalMode emode, const CeedScalar *u, CeedScalar *v) { 415 int ierr; 416 if (!basis->Apply) return CeedError(basis->ceed, 1, 417 "Backend does not support BasisApply"); 418 ierr = basis->Apply(basis, nelem, tmode, emode, u, v); CeedChk(ierr); 419 return 0; 420 } 421 422 /// Get total number of nodes (in dim dimensions) 423 int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) { 424 *P = CeedPowInt(basis->P1d, basis->dim); 425 return 0; 426 } 427 428 /// Get total number of quadrature points (in dim dimensions) 429 int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) { 430 *Q = CeedPowInt(basis->Q1d, basis->dim); 431 return 0; 432 } 433 434 /// Destroy a CeedBasis 435 int CeedBasisDestroy(CeedBasis *basis) { 436 int ierr; 437 438 if (!*basis || --(*basis)->refcount > 0) return 0; 439 if ((*basis)->Destroy) { 440 ierr = (*basis)->Destroy(*basis); CeedChk(ierr); 441 } 442 ierr = CeedFree(&(*basis)->interp1d); CeedChk(ierr); 443 ierr = CeedFree(&(*basis)->grad1d); CeedChk(ierr); 444 ierr = CeedFree(&(*basis)->qref1d); CeedChk(ierr); 445 ierr = CeedFree(&(*basis)->qweight1d); CeedChk(ierr); 446 ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr); 447 ierr = CeedFree(basis); CeedChk(ierr); 448 return 0; 449 } 450 451 /// Indicate that the quadrature points are colocated with the dofs 452 CeedBasis CEED_BASIS_COLOCATED = &ceed_basis_colocated; 453 454 /// @} 455