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 /** 34 @brief Create a tensor product basis for H^1 discretizations 35 36 @param ceed A Ceed object where the CeedBasis will be created 37 @param dim Topological dimension 38 @param ncomp Number of field components (1 for scalar fields) 39 @param P1d Number of nodes in one dimension 40 @param Q1d Number of quadrature points in one dimension 41 @param interp1d Row-major Q1d × P1d matrix expressing the values of nodal 42 basis functions at quadrature points 43 @param grad1d Row-major Q1d × P1d matrix expressing derivatives of nodal 44 basis functions at quadrature points 45 @param qref1d Array of length Q1d holding the locations of quadrature points 46 on the 1D reference element [-1, 1] 47 @param qweight1d Array of length Q1d holding the quadrature weights on the 48 reference element 49 @param[out] basis Address of the variable where the newly created 50 CeedBasis will be stored. 51 52 @return An error code: 0 - success, otherwise - failure 53 **/ 54 int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt ncomp, CeedInt P1d, 55 CeedInt Q1d, const CeedScalar *interp1d, 56 const CeedScalar *grad1d, const CeedScalar *qref1d, 57 const CeedScalar *qweight1d, CeedBasis *basis) { 58 int ierr; 59 60 if (!ceed->BasisCreateTensorH1) 61 return CeedError(ceed, 1, "Backend does not support BasisCreateTensorH1"); 62 ierr = CeedCalloc(1,basis); CeedChk(ierr); 63 (*basis)->ceed = ceed; 64 ceed->refcount++; 65 (*basis)->refcount = 1; 66 (*basis)->dim = dim; 67 (*basis)->ncomp = ncomp; 68 (*basis)->P1d = P1d; 69 (*basis)->Q1d = Q1d; 70 ierr = CeedMalloc(Q1d,&(*basis)->qref1d); CeedChk(ierr); 71 ierr = CeedMalloc(Q1d,&(*basis)->qweight1d); CeedChk(ierr); 72 memcpy((*basis)->qref1d, qref1d, Q1d*sizeof(qref1d[0])); 73 memcpy((*basis)->qweight1d, qweight1d, Q1d*sizeof(qweight1d[0])); 74 ierr = CeedMalloc(Q1d*P1d,&(*basis)->interp1d); CeedChk(ierr); 75 ierr = CeedMalloc(Q1d*P1d,&(*basis)->grad1d); CeedChk(ierr); 76 memcpy((*basis)->interp1d, interp1d, Q1d*P1d*sizeof(interp1d[0])); 77 memcpy((*basis)->grad1d, grad1d, Q1d*P1d*sizeof(interp1d[0])); 78 ierr = ceed->BasisCreateTensorH1(ceed, dim, P1d, Q1d, interp1d, grad1d, qref1d, 79 qweight1d, *basis); CeedChk(ierr); 80 return 0; 81 } 82 83 /** 84 @brief Create a tensor product Lagrange basis 85 86 @param ceed A Ceed object where the CeedBasis will be created 87 @param dim Topological dimension of element 88 @param ncomp Number of field components 89 @param P Number of Gauss-Lobatto nodes in one dimension. The 90 polynomial degree of the resulting Q_k element is k=P-1. 91 @param Q Number of quadrature points in one dimension. 92 @param qmode Distribution of the Q quadrature points (affects order of 93 accuracy for the quadrature) 94 @param[out] basis Address of the variable where the newly created 95 CeedBasis will be stored. 96 97 @return An error code: 0 - success, otherwise - failure 98 **/ 99 int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt ncomp, 100 CeedInt P, CeedInt Q, 101 CeedQuadMode qmode, CeedBasis *basis) { 102 // Allocate 103 int ierr, i, j, k; 104 CeedScalar c1, c2, c3, c4, dx, *nodes, *interp1d, *grad1d, *qref1d, *qweight1d; 105 ierr = CeedCalloc(P*Q, &interp1d); CeedChk(ierr); 106 ierr = CeedCalloc(P*Q, &grad1d); CeedChk(ierr); 107 ierr = CeedCalloc(P, &nodes); CeedChk(ierr); 108 ierr = CeedCalloc(Q, &qref1d); CeedChk(ierr); 109 ierr = CeedCalloc(Q, &qweight1d); CeedChk(ierr); 110 // Get Nodes and Weights 111 ierr = CeedLobattoQuadrature(P, nodes, NULL); CeedChk(ierr); 112 switch (qmode) { 113 case CEED_GAUSS: 114 ierr = CeedGaussQuadrature(Q, qref1d, qweight1d); CeedChk(ierr); 115 break; 116 case CEED_GAUSS_LOBATTO: 117 ierr = CeedLobattoQuadrature(Q, qref1d, qweight1d); CeedChk(ierr); 118 break; 119 } 120 // Build B, D matrix 121 // Fornberg, 1998 122 for (i = 0; i < Q; i++) { 123 c1 = 1.0; 124 c3 = nodes[0] - qref1d[i]; 125 interp1d[i*P+0] = 1.0; 126 for (j = 1; j < P; j++) { 127 c2 = 1.0; 128 c4 = c3; 129 c3 = nodes[j] - qref1d[i]; 130 for (k = 0; k < j; k++) { 131 dx = nodes[j] - nodes[k]; 132 c2 *= dx; 133 if (k == j - 1) { 134 grad1d[i*P + j] = c1*(interp1d[i*P + k] - c4*grad1d[i*P + k]) / c2; 135 interp1d[i*P + j] = - c1*c4*interp1d[i*P + k] / c2; 136 } 137 grad1d[i*P + k] = (c3*grad1d[i*P + k] - interp1d[i*P + k]) / dx; 138 interp1d[i*P + k] = c3*interp1d[i*P + k] / dx; 139 } 140 c1 = c2; 141 } 142 } 143 // // Pass to CeedBasisCreateTensorH1 144 ierr = CeedBasisCreateTensorH1(ceed, dim, ncomp, P, Q, interp1d, grad1d, qref1d, 145 qweight1d, basis); CeedChk(ierr); 146 ierr = CeedFree(&interp1d); CeedChk(ierr); 147 ierr = CeedFree(&grad1d); CeedChk(ierr); 148 ierr = CeedFree(&nodes); CeedChk(ierr); 149 ierr = CeedFree(&qref1d); CeedChk(ierr); 150 ierr = CeedFree(&qweight1d); CeedChk(ierr); 151 return 0; 152 } 153 154 /** 155 @brief Construct a Gauss-Legendre quadrature 156 157 @param Q Number of quadrature points (integrates polynomials of 158 degree 2*Q-1 exactly) 159 @param[out] qref1d Array of length Q to hold the abscissa on [-1, 1] 160 @param[out] qweight1d Array of length Q to hold the weights 161 162 @return An error code: 0 - success, otherwise - failure 163 **/ 164 int CeedGaussQuadrature(CeedInt Q, CeedScalar *qref1d, CeedScalar *qweight1d) { 165 // Allocate 166 CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0*atan(1.0); 167 // Build qref1d, qweight1d 168 for (int i = 0; i <= Q/2; i++) { 169 // Guess 170 xi = cos(PI*(CeedScalar)(2*i+1)/((CeedScalar)(2*Q))); 171 // Pn(xi) 172 P0 = 1.0; 173 P1 = xi; 174 P2 = 0.0; 175 for (int j = 2; j <= Q; j++) { 176 P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 177 P0 = P1; 178 P1 = P2; 179 } 180 // First Newton Step 181 dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 182 xi = xi-P2/dP2; 183 // Newton to convergence 184 for (int k=0; k<100 && fabs(P2)>1e-15; k++) { 185 P0 = 1.0; 186 P1 = xi; 187 for (int j = 2; j <= Q; j++) { 188 P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 189 P0 = P1; 190 P1 = P2; 191 } 192 dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 193 xi = xi-P2/dP2; 194 } 195 // Save xi, wi 196 wi = 2.0/((1.0-xi*xi)*dP2*dP2); 197 qweight1d[i] = wi; 198 qweight1d[Q-1-i] = wi; 199 qref1d[i] = -xi; 200 qref1d[Q-1-i]= xi; 201 } 202 return 0; 203 } 204 205 /** 206 @brief Construct a Gauss-Legendre-Lobatto quadrature 207 208 @param Q Number of quadrature points (integrates polynomials of 209 degree 2*Q-3 exactly) 210 @param[out] qref1d Array of length Q to hold the abscissa on [-1, 1] 211 @param[out] qweight1d Array of length Q to hold the weights 212 213 @return An error code: 0 - success, otherwise - failure 214 **/ 215 int CeedLobattoQuadrature(CeedInt Q, CeedScalar *qref1d, 216 CeedScalar *qweight1d) { 217 // Allocate 218 CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0*atan(1.0); 219 // Build qref1d, qweight1d 220 // Set endpoints 221 wi = 2.0/((CeedScalar)(Q*(Q-1))); 222 if (qweight1d) { 223 qweight1d[0] = wi; 224 qweight1d[Q-1] = wi; 225 } 226 qref1d[0] = -1.0; 227 qref1d[Q-1] = 1.0; 228 // Interior 229 for (int i = 1; i <= (Q-1)/2; i++) { 230 // Guess 231 xi = cos(PI*(CeedScalar)(i)/(CeedScalar)(Q-1)); 232 // Pn(xi) 233 P0 = 1.0; 234 P1 = xi; 235 P2 = 0.0; 236 for (int j = 2; j < Q; j++) { 237 P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 238 P0 = P1; 239 P1 = P2; 240 } 241 // First Newton step 242 dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 243 d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi); 244 xi = xi-dP2/d2P2; 245 // Newton to convergence 246 for (int k=0; k<100 && fabs(dP2)>1e-15; k++) { 247 P0 = 1.0; 248 P1 = xi; 249 for (int j = 2; j < Q; j++) { 250 P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 251 P0 = P1; 252 P1 = P2; 253 } 254 dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 255 d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi); 256 xi = xi-dP2/d2P2; 257 } 258 // Save xi, wi 259 wi = 2.0/(((CeedScalar)(Q*(Q-1)))*P2*P2); 260 if (qweight1d) { 261 qweight1d[i] = wi; 262 qweight1d[Q-1-i] = wi; 263 } 264 qref1d[i] = -xi; 265 qref1d[Q-1-i]= xi; 266 } 267 return 0; 268 } 269 270 static int CeedScalarView(const char *name, const char *fpformat, CeedInt m, 271 CeedInt n, const CeedScalar *a, FILE *stream) { 272 for (int i=0; i<m; i++) { 273 if (m > 1) fprintf(stream, "%12s[%d]:", name, i); 274 else fprintf(stream, "%12s:", name); 275 for (int j=0; j<n; j++) { 276 fprintf(stream, fpformat, fabs(a[i*n+j]) > 1E-14 ? a[i*n+j] : 0); 277 } 278 fputs("\n", stream); 279 } 280 return 0; 281 } 282 283 /** 284 @brief View a CeedBasis 285 286 @param basis CeedBasis to view 287 @param stream Stream to view to, e.g., stdout 288 289 @return An error code: 0 - success, otherwise - failure 290 **/ 291 int CeedBasisView(CeedBasis basis, FILE *stream) { 292 int ierr; 293 294 fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P1d, 295 basis->Q1d); 296 ierr = CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q1d, basis->qref1d, 297 stream); CeedChk(ierr); 298 ierr = CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q1d, basis->qweight1d, 299 stream); CeedChk(ierr); 300 ierr = CeedScalarView("interp1d", "\t% 12.8f", basis->Q1d, basis->P1d, 301 basis->interp1d, stream); CeedChk(ierr); 302 ierr = CeedScalarView("grad1d", "\t% 12.8f", basis->Q1d, basis->P1d, 303 basis->grad1d, stream); CeedChk(ierr); 304 return 0; 305 } 306 307 // Computes A = (I - b v v^T) A 308 // where A is an mxn matrix indexed as A[i*row + j*col] 309 static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v, 310 CeedScalar b, CeedInt m, CeedInt n, 311 CeedInt row, CeedInt col) { 312 for (CeedInt j=0; j<n; j++) { 313 CeedScalar w = A[0*row + j*col]; 314 for (CeedInt i=1; i<m; i++) w += v[i] * A[i*row + j*col]; 315 A[0*row + j*col] -= b * w; 316 for (CeedInt i=1; i<m; i++) A[i*row + j*col] -= b * w * v[i]; 317 } 318 return 0; 319 } 320 321 // Compute A = Q A where Q is mxk and A is mxn. k<m 322 static int CeedHouseholderApplyQ(CeedScalar *A, const CeedScalar *Q, 323 const CeedScalar *tau, CeedTransposeMode tmode, 324 CeedInt m, CeedInt n, CeedInt k, 325 CeedInt row, CeedInt col) { 326 CeedScalar v[m]; 327 for (CeedInt ii=0; ii<k; ii++) { 328 CeedInt i = tmode == CEED_TRANSPOSE ? ii : k-1-ii; 329 for (CeedInt j=i+1; j<m; j++) { 330 v[j] = Q[j*k+i]; 331 } 332 // Apply Householder reflector (I - tau v v^T) colograd1d^T 333 CeedHouseholderReflect(&A[i*row], &v[i], tau[i], m-i, n, row, col); 334 } 335 return 0; 336 } 337 338 /** 339 @brief Return QR Factorization of matrix 340 341 @param[out] mat Row-major matrix to be factorized in place 342 @param[out] tau Vector of length m of scaling fators 343 @param m Number of rows 344 @param n Number of columns 345 346 @return An error code: 0 - success, otherwise - failure 347 **/ 348 int CeedQRFactorization(CeedScalar *mat, CeedScalar *tau, 349 CeedInt m, CeedInt n) { 350 CeedInt i, j; 351 CeedScalar v[m]; 352 353 for (i=0; i<n; i++) { 354 // Calculate Householder vector, magnitude 355 CeedScalar sigma = 0.0; 356 v[i] = mat[i+n*i]; 357 for (j=i+1; j<m; j++) { 358 v[j] = mat[i+n*j]; 359 sigma += v[j] * v[j]; 360 } 361 CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m] 362 CeedScalar Rii = -copysign(norm, v[i]); 363 v[i] -= Rii; 364 // norm of v[i:m] after modification above and scaling below 365 // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 366 // tau = 2 / (norm*norm) 367 tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma); 368 for (j=i+1; j<m; j++) v[j] /= v[i]; 369 370 // Apply Householder reflector to lower right panel 371 CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1); 372 // Save v 373 mat[i+n*i] = Rii; 374 for (j=i+1; j<m; j++) { 375 mat[i+n*j] = v[j]; 376 } 377 } 378 379 return 0; 380 } 381 382 /** 383 @brief Return colocated grad matrix 384 385 @param basis CeedBasis 386 @param[out] colograd1d Row-major Q1d × Q1d matrix expressing derivatives of 387 basis functions at quadrature points 388 389 @return An error code: 0 - success, otherwise - failure 390 **/ 391 int CeedBasisGetColocatedGrad(CeedBasis basis, CeedScalar *colograd1d) { 392 int i, j, k; 393 CeedInt ierr, P1d=(basis)->P1d, Q1d=(basis)->Q1d; 394 CeedScalar *interp1d, *grad1d, tau[Q1d]; 395 396 ierr = CeedMalloc(Q1d*P1d, &interp1d); CeedChk(ierr); 397 ierr = CeedMalloc(Q1d*P1d, &grad1d); CeedChk(ierr); 398 memcpy(interp1d, (basis)->interp1d, Q1d*P1d*sizeof(basis)->interp1d[0]); 399 memcpy(grad1d, (basis)->grad1d, Q1d*P1d*sizeof(basis)->interp1d[0]); 400 401 // QR Factorization, interp1d = Q R 402 ierr = CeedQRFactorization(interp1d, tau, Q1d, P1d); CeedChk(ierr); 403 404 // Apply Rinv, colograd1d = grad1d Rinv 405 for (i=0; i<Q1d; i++) { // Row i 406 colograd1d[Q1d*i] = grad1d[P1d*i]/interp1d[0]; 407 for (j=1; j<P1d; j++) { // Column j 408 colograd1d[j+Q1d*i] = grad1d[j+P1d*i]; 409 for (k=0; k<j; k++) { 410 colograd1d[j+Q1d*i] -= interp1d[j+P1d*k]*colograd1d[k+Q1d*i]; 411 } 412 colograd1d[j+Q1d*i] /= interp1d[j+P1d*j]; 413 } 414 for (j=P1d; j<Q1d; j++) { 415 colograd1d[j+Q1d*i] = 0; 416 } 417 } 418 419 // Apply Qtranspose, colograd = colograd Qtranspose 420 CeedHouseholderApplyQ(colograd1d, interp1d, tau, CEED_NOTRANSPOSE, 421 Q1d, Q1d, P1d, 1, Q1d); 422 423 ierr = CeedFree(&interp1d); CeedChk(ierr); 424 ierr = CeedFree(&grad1d); CeedChk(ierr); 425 426 return 0; 427 } 428 429 /** 430 @brief Apply basis evaluation from nodes to quadrature points or vice-versa 431 432 @param basis CeedBasis to evaluate 433 @param nelem The number of elements to apply the basis evaluation to; 434 the backend will specify the ordering in 435 ElemRestrictionCreateBlocked 436 @param tmode \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature 437 points, \ref CEED_TRANSPOSE to apply the transpose, mapping 438 from quadrature points to nodes 439 @param emode \ref CEED_EVAL_INTERP to obtain interpolated values, 440 \ref CEED_EVAL_GRAD to obtain gradients. 441 @param[in] u Input array 442 @param[out] v Output array 443 444 @return An error code: 0 - success, otherwise - failure 445 **/ 446 int CeedBasisApply(CeedBasis basis, CeedInt nelem, CeedTransposeMode tmode, 447 CeedEvalMode emode, const CeedScalar *u, CeedScalar *v) { 448 int ierr; 449 if (!basis->Apply) return CeedError(basis->ceed, 1, 450 "Backend does not support BasisApply"); 451 ierr = basis->Apply(basis, nelem, tmode, emode, u, v); CeedChk(ierr); 452 return 0; 453 } 454 455 /** 456 @brief Get total number of nodes (in dim dimensions) 457 458 @param basis CeedBasis 459 @param[out] P Number of nodes 460 461 @return An error code: 0 - success, otherwise - failure 462 **/ 463 int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) { 464 *P = CeedPowInt(basis->P1d, basis->dim); 465 return 0; 466 } 467 468 /** 469 @brief Get total number of quadrature points (in dim dimensions) 470 471 @param basis CeedBasis 472 @param[out] Q Number of quadrature points 473 474 @return An error code: 0 - success, otherwise - failure 475 **/ 476 int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) { 477 *Q = CeedPowInt(basis->Q1d, basis->dim); 478 return 0; 479 } 480 481 /** 482 @brief Destroy a CeedBasis 483 484 @param basis CeedBasis to destroy 485 486 @return An error code: 0 - success, otherwise - failure 487 **/ 488 int CeedBasisDestroy(CeedBasis *basis) { 489 int ierr; 490 491 if (!*basis || --(*basis)->refcount > 0) return 0; 492 if ((*basis)->Destroy) { 493 ierr = (*basis)->Destroy(*basis); CeedChk(ierr); 494 } 495 ierr = CeedFree(&(*basis)->interp1d); CeedChk(ierr); 496 ierr = CeedFree(&(*basis)->grad1d); CeedChk(ierr); 497 ierr = CeedFree(&(*basis)->qref1d); CeedChk(ierr); 498 ierr = CeedFree(&(*basis)->qweight1d); CeedChk(ierr); 499 ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr); 500 ierr = CeedFree(basis); CeedChk(ierr); 501 return 0; 502 } 503 504 /// @cond DOXYGEN_SKIP 505 // Indicate that the quadrature points are colocated with the dofs 506 CeedBasis CEED_BASIS_COLOCATED = &ceed_basis_colocated; 507 /// @endcond 508 /// @} 509