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 <ceed-backend.h> 19 #include <math.h> 20 #include <stdio.h> 21 #include <stdlib.h> 22 #include <string.h> 23 24 /// @cond DOXYGEN_SKIP 25 static struct CeedBasis_private ceed_basis_collocated; 26 /// @endcond 27 28 /// @file 29 /// Implementation of public CeedBasis interfaces 30 /// 31 /// @addtogroup CeedBasis 32 /// @{ 33 34 /** 35 @brief Create a tensor product basis for H^1 discretizations 36 37 @param ceed A Ceed object where the CeedBasis will be created 38 @param dim Topological dimension 39 @param ncomp Number of field components (1 for scalar fields) 40 @param P1d Number of nodes in one dimension 41 @param Q1d Number of quadrature points in one dimension 42 @param interp1d Row-major Q1d × P1d matrix expressing the values of nodal 43 basis functions at quadrature points 44 @param grad1d Row-major Q1d × P1d matrix expressing derivatives of nodal 45 basis functions at quadrature points 46 @param qref1d Array of length Q1d holding the locations of quadrature points 47 on the 1D reference element [-1, 1] 48 @param qweight1d Array of length Q1d holding the quadrature weights on the 49 reference element 50 @param[out] basis Address of the variable where the newly created 51 CeedBasis will be stored. 52 53 @return An error code: 0 - success, otherwise - failure 54 55 @ref Basic 56 **/ 57 int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt ncomp, CeedInt P1d, 58 CeedInt Q1d, const CeedScalar *interp1d, 59 const CeedScalar *grad1d, const CeedScalar *qref1d, 60 const CeedScalar *qweight1d, CeedBasis *basis) { 61 int ierr; 62 63 if (dim<1) 64 return CeedError(ceed, 1, "Basis dimension must be a positive value"); 65 66 if (!ceed->BasisCreateTensorH1) { 67 Ceed delegate; 68 ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr); 69 70 if (!delegate) 71 return CeedError(ceed, 1, "Backend does not support BasisCreateTensorH1"); 72 73 ierr = CeedBasisCreateTensorH1(delegate, dim, ncomp, P1d, 74 Q1d, interp1d, grad1d, qref1d, 75 qweight1d, basis); CeedChk(ierr); 76 return 0; 77 } 78 ierr = CeedCalloc(1,basis); CeedChk(ierr); 79 (*basis)->ceed = ceed; 80 ceed->refcount++; 81 (*basis)->refcount = 1; 82 (*basis)->tensorbasis = 1; 83 (*basis)->dim = dim; 84 (*basis)->ncomp = ncomp; 85 (*basis)->P1d = P1d; 86 (*basis)->Q1d = Q1d; 87 (*basis)->P = CeedIntPow(P1d, dim); 88 (*basis)->Q = CeedIntPow(Q1d, dim); 89 ierr = CeedMalloc(Q1d,&(*basis)->qref1d); CeedChk(ierr); 90 ierr = CeedMalloc(Q1d,&(*basis)->qweight1d); CeedChk(ierr); 91 memcpy((*basis)->qref1d, qref1d, Q1d*sizeof(qref1d[0])); 92 memcpy((*basis)->qweight1d, qweight1d, Q1d*sizeof(qweight1d[0])); 93 ierr = CeedMalloc(Q1d*P1d,&(*basis)->interp1d); CeedChk(ierr); 94 ierr = CeedMalloc(Q1d*P1d,&(*basis)->grad1d); CeedChk(ierr); 95 memcpy((*basis)->interp1d, interp1d, Q1d*P1d*sizeof(interp1d[0])); 96 memcpy((*basis)->grad1d, grad1d, Q1d*P1d*sizeof(grad1d[0])); 97 ierr = ceed->BasisCreateTensorH1(dim, P1d, Q1d, interp1d, grad1d, qref1d, 98 qweight1d, *basis); CeedChk(ierr); 99 return 0; 100 } 101 102 /** 103 @brief Create a tensor product Lagrange basis 104 105 @param ceed A Ceed object where the CeedBasis will be created 106 @param dim Topological dimension of element 107 @param ncomp Number of field components 108 @param P Number of Gauss-Lobatto nodes in one dimension. The 109 polynomial degree of the resulting Q_k element is k=P-1. 110 @param Q Number of quadrature points in one dimension. 111 @param qmode Distribution of the Q quadrature points (affects order of 112 accuracy for the quadrature) 113 @param[out] basis Address of the variable where the newly created 114 CeedBasis will be stored. 115 116 @return An error code: 0 - success, otherwise - failure 117 118 @ref Basic 119 **/ 120 int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt ncomp, 121 CeedInt P, CeedInt Q, 122 CeedQuadMode qmode, CeedBasis *basis) { 123 // Allocate 124 int ierr, i, j, k; 125 CeedScalar c1, c2, c3, c4, dx, *nodes, *interp1d, *grad1d, *qref1d, *qweight1d; 126 127 if (dim<1) 128 return CeedError(ceed, 1, "Basis dimension must be a positive value"); 129 130 ierr = CeedCalloc(P*Q, &interp1d); CeedChk(ierr); 131 ierr = CeedCalloc(P*Q, &grad1d); CeedChk(ierr); 132 ierr = CeedCalloc(P, &nodes); CeedChk(ierr); 133 ierr = CeedCalloc(Q, &qref1d); CeedChk(ierr); 134 ierr = CeedCalloc(Q, &qweight1d); CeedChk(ierr); 135 // Get Nodes and Weights 136 ierr = CeedLobattoQuadrature(P, nodes, NULL); CeedChk(ierr); 137 switch (qmode) { 138 case CEED_GAUSS: 139 ierr = CeedGaussQuadrature(Q, qref1d, qweight1d); CeedChk(ierr); 140 break; 141 case CEED_GAUSS_LOBATTO: 142 ierr = CeedLobattoQuadrature(Q, qref1d, qweight1d); CeedChk(ierr); 143 break; 144 } 145 // Build B, D matrix 146 // Fornberg, 1998 147 for (i = 0; i < Q; i++) { 148 c1 = 1.0; 149 c3 = nodes[0] - qref1d[i]; 150 interp1d[i*P+0] = 1.0; 151 for (j = 1; j < P; j++) { 152 c2 = 1.0; 153 c4 = c3; 154 c3 = nodes[j] - qref1d[i]; 155 for (k = 0; k < j; k++) { 156 dx = nodes[j] - nodes[k]; 157 c2 *= dx; 158 if (k == j - 1) { 159 grad1d[i*P + j] = c1*(interp1d[i*P + k] - c4*grad1d[i*P + k]) / c2; 160 interp1d[i*P + j] = - c1*c4*interp1d[i*P + k] / c2; 161 } 162 grad1d[i*P + k] = (c3*grad1d[i*P + k] - interp1d[i*P + k]) / dx; 163 interp1d[i*P + k] = c3*interp1d[i*P + k] / dx; 164 } 165 c1 = c2; 166 } 167 } 168 // // Pass to CeedBasisCreateTensorH1 169 ierr = CeedBasisCreateTensorH1(ceed, dim, ncomp, P, Q, interp1d, grad1d, qref1d, 170 qweight1d, basis); CeedChk(ierr); 171 ierr = CeedFree(&interp1d); CeedChk(ierr); 172 ierr = CeedFree(&grad1d); CeedChk(ierr); 173 ierr = CeedFree(&nodes); CeedChk(ierr); 174 ierr = CeedFree(&qref1d); CeedChk(ierr); 175 ierr = CeedFree(&qweight1d); CeedChk(ierr); 176 return 0; 177 } 178 179 /** 180 @brief Create a non tensor product basis for H^1 discretizations 181 182 @param ceed A Ceed object where the CeedBasis will be created 183 @param topo Topology of element, e.g. hypercube, simplex, ect 184 @param ncomp Number of field components (1 for scalar fields) 185 @param nnodes Total number of nodes 186 @param nqpts Total number of quadrature points 187 @param interp Row-major nqpts × nnodes matrix expressing the values of 188 nodal basis functions at quadrature points 189 @param grad Row-major (nqpts x dim) × nnodes matrix expressing 190 derivatives of nodal basis functions at quadrature points 191 @param qref Array of length nqpts holding the locations of quadrature 192 points on the reference element [-1, 1] 193 @param qweight Array of length nqpts holding the quadrature weights on the 194 reference element 195 @param[out] basis Address of the variable where the newly created 196 CeedBasis will be stored. 197 198 @return An error code: 0 - success, otherwise - failure 199 200 @ref Basic 201 **/ 202 int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt ncomp, 203 CeedInt nnodes, CeedInt nqpts, 204 const CeedScalar *interp, 205 const CeedScalar *grad, const CeedScalar *qref, 206 const CeedScalar *qweight, CeedBasis *basis) { 207 int ierr; 208 CeedInt P = nnodes, Q = nqpts, dim = 0; 209 210 if (!ceed->BasisCreateH1) { 211 Ceed delegate; 212 ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr); 213 214 if (!delegate) 215 return CeedError(ceed, 1, "Backend does not support BasisCreateH1"); 216 217 ierr = CeedBasisCreateH1(delegate, topo, ncomp, nnodes, 218 nqpts, interp, grad, qref, 219 qweight, basis); CeedChk(ierr); 220 return 0; 221 } 222 223 ierr = CeedCalloc(1,basis); CeedChk(ierr); 224 225 ierr = CeedBasisGetTopologyDimension(topo, &dim); CeedChk(ierr); 226 227 (*basis)->ceed = ceed; 228 ceed->refcount++; 229 (*basis)->refcount = 1; 230 (*basis)->tensorbasis = 0; 231 (*basis)->dim = dim; 232 (*basis)->ncomp = ncomp; 233 (*basis)->P = P; 234 (*basis)->Q = Q; 235 ierr = CeedMalloc(Q*dim,&(*basis)->qref1d); CeedChk(ierr); 236 ierr = CeedMalloc(Q,&(*basis)->qweight1d); CeedChk(ierr); 237 memcpy((*basis)->qref1d, qref, Q*dim*sizeof(qref[0])); 238 memcpy((*basis)->qweight1d, qweight, Q*sizeof(qweight[0])); 239 ierr = CeedMalloc(Q*P,&(*basis)->interp1d); CeedChk(ierr); 240 ierr = CeedMalloc(dim*Q*P,&(*basis)->grad1d); CeedChk(ierr); 241 memcpy((*basis)->interp1d, interp, Q*P*sizeof(interp[0])); 242 memcpy((*basis)->grad1d, grad, dim*Q*P*sizeof(grad[0])); 243 ierr = ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, qref, 244 qweight, *basis); CeedChk(ierr); 245 return 0; 246 } 247 248 /** 249 @brief Construct a Gauss-Legendre quadrature 250 251 @param Q Number of quadrature points (integrates polynomials of 252 degree 2*Q-1 exactly) 253 @param[out] qref1d Array of length Q to hold the abscissa on [-1, 1] 254 @param[out] qweight1d Array of length Q to hold the weights 255 256 @return An error code: 0 - success, otherwise - failure 257 258 @ref Utility 259 **/ 260 int CeedGaussQuadrature(CeedInt Q, CeedScalar *qref1d, CeedScalar *qweight1d) { 261 // Allocate 262 CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0*atan(1.0); 263 // Build qref1d, qweight1d 264 for (int i = 0; i <= Q/2; i++) { 265 // Guess 266 xi = cos(PI*(CeedScalar)(2*i+1)/((CeedScalar)(2*Q))); 267 // Pn(xi) 268 P0 = 1.0; 269 P1 = xi; 270 P2 = 0.0; 271 for (int j = 2; j <= Q; j++) { 272 P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 273 P0 = P1; 274 P1 = P2; 275 } 276 // First Newton Step 277 dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 278 xi = xi-P2/dP2; 279 // Newton to convergence 280 for (int k=0; k<100 && fabs(P2)>1e-15; k++) { 281 P0 = 1.0; 282 P1 = xi; 283 for (int j = 2; j <= Q; j++) { 284 P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 285 P0 = P1; 286 P1 = P2; 287 } 288 dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 289 xi = xi-P2/dP2; 290 } 291 // Save xi, wi 292 wi = 2.0/((1.0-xi*xi)*dP2*dP2); 293 qweight1d[i] = wi; 294 qweight1d[Q-1-i] = wi; 295 qref1d[i] = -xi; 296 qref1d[Q-1-i]= xi; 297 } 298 return 0; 299 } 300 301 /** 302 @brief Construct a Gauss-Legendre-Lobatto quadrature 303 304 @param Q Number of quadrature points (integrates polynomials of 305 degree 2*Q-3 exactly) 306 @param[out] qref1d Array of length Q to hold the abscissa on [-1, 1] 307 @param[out] qweight1d Array of length Q to hold the weights 308 309 @return An error code: 0 - success, otherwise - failure 310 311 @ref Utility 312 **/ 313 int CeedLobattoQuadrature(CeedInt Q, CeedScalar *qref1d, 314 CeedScalar *qweight1d) { 315 // Allocate 316 CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0*atan(1.0); 317 // Build qref1d, qweight1d 318 // Set endpoints 319 wi = 2.0/((CeedScalar)(Q*(Q-1))); 320 if (qweight1d) { 321 qweight1d[0] = wi; 322 qweight1d[Q-1] = wi; 323 } 324 qref1d[0] = -1.0; 325 qref1d[Q-1] = 1.0; 326 // Interior 327 for (int i = 1; i <= (Q-1)/2; i++) { 328 // Guess 329 xi = cos(PI*(CeedScalar)(i)/(CeedScalar)(Q-1)); 330 // Pn(xi) 331 P0 = 1.0; 332 P1 = xi; 333 P2 = 0.0; 334 for (int j = 2; j < Q; j++) { 335 P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 336 P0 = P1; 337 P1 = P2; 338 } 339 // First Newton step 340 dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 341 d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi); 342 xi = xi-dP2/d2P2; 343 // Newton to convergence 344 for (int k=0; k<100 && fabs(dP2)>1e-15; k++) { 345 P0 = 1.0; 346 P1 = xi; 347 for (int j = 2; j < Q; j++) { 348 P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 349 P0 = P1; 350 P1 = P2; 351 } 352 dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 353 d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi); 354 xi = xi-dP2/d2P2; 355 } 356 // Save xi, wi 357 wi = 2.0/(((CeedScalar)(Q*(Q-1)))*P2*P2); 358 if (qweight1d) { 359 qweight1d[i] = wi; 360 qweight1d[Q-1-i] = wi; 361 } 362 qref1d[i] = -xi; 363 qref1d[Q-1-i]= xi; 364 } 365 return 0; 366 } 367 368 /** 369 @brief View an array stored in a CeedBasis 370 371 @param name Name of array 372 @param fpformat Printing format 373 @param m Number of rows in array 374 @param n Number of columns in array 375 @param a Array to be viewed 376 @param stream Stream to view to, e.g., stdout 377 378 @return An error code: 0 - success, otherwise - failure 379 380 @ref Utility 381 **/ 382 static int CeedScalarView(const char *name, const char *fpformat, CeedInt m, 383 CeedInt n, const CeedScalar *a, FILE *stream) { 384 for (int i=0; i<m; i++) { 385 if (m > 1) fprintf(stream, "%12s[%d]:", name, i); 386 else fprintf(stream, "%12s:", name); 387 for (int j=0; j<n; j++) { 388 fprintf(stream, fpformat, fabs(a[i*n+j]) > 1E-14 ? a[i*n+j] : 0); 389 } 390 fputs("\n", stream); 391 } 392 return 0; 393 } 394 395 /** 396 @brief View a CeedBasis 397 398 @param basis CeedBasis to view 399 @param stream Stream to view to, e.g., stdout 400 401 @return An error code: 0 - success, otherwise - failure 402 403 @ref Utility 404 **/ 405 int CeedBasisView(CeedBasis basis, FILE *stream) { 406 int ierr; 407 408 if (basis->tensorbasis) { 409 fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P1d, 410 basis->Q1d); 411 ierr = CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q1d, basis->qref1d, 412 stream); CeedChk(ierr); 413 ierr = CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q1d, 414 basis->qweight1d, stream); CeedChk(ierr); 415 ierr = CeedScalarView("interp1d", "\t% 12.8f", basis->Q1d, basis->P1d, 416 basis->interp1d, stream); CeedChk(ierr); 417 ierr = CeedScalarView("grad1d", "\t% 12.8f", basis->Q1d, basis->P1d, 418 basis->grad1d, stream); CeedChk(ierr); 419 } else { 420 fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P, 421 basis->Q); 422 ierr = CeedScalarView("qref", "\t% 12.8f", 1, basis->Q*basis->dim, 423 basis->qref1d, 424 stream); CeedChk(ierr); 425 ierr = CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->qweight1d, 426 stream); CeedChk(ierr); 427 ierr = CeedScalarView("interp", "\t% 12.8f", basis->Q, basis->P, 428 basis->interp1d, stream); CeedChk(ierr); 429 ierr = CeedScalarView("grad", "\t% 12.8f", basis->dim*basis->Q, basis->P, 430 basis->grad1d, stream); CeedChk(ierr); 431 } 432 return 0; 433 } 434 435 /** 436 @brief Compute Householder Reflection 437 438 Computes A = (I - b v v^T) A 439 where A is an mxn matrix indexed as A[i*row + j*col] 440 441 @param[out] A Matrix to apply Householder reflection to, in place 442 @param v Householder vector 443 @param b Scaling factor 444 @param m Number of rows in A 445 @param n Number of columns in A 446 @param row Col stride 447 @param col Row stride 448 449 @return An error code: 0 - success, otherwise - failure 450 451 @ref Developer 452 **/ 453 static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v, 454 CeedScalar b, CeedInt m, CeedInt n, 455 CeedInt row, CeedInt col) { 456 for (CeedInt j=0; j<n; j++) { 457 CeedScalar w = A[0*row + j*col]; 458 for (CeedInt i=1; i<m; i++) w += v[i] * A[i*row + j*col]; 459 A[0*row + j*col] -= b * w; 460 for (CeedInt i=1; i<m; i++) A[i*row + j*col] -= b * w * v[i]; 461 } 462 return 0; 463 } 464 465 /** 466 @brief Apply Householder Q matrix 467 468 Compute A = Q A where Q is mxk and A is mxn. k<m 469 470 @param[out] A Matrix to apply Householder Q to, in place 471 @param Q Householder Q matrix 472 @param tau Householder scaling factors 473 @param tmode Transpose mode for application 474 @param m Number of rows in A 475 @param n Number of columns in A 476 @param k Index of row targeted 477 @param row Col stride 478 @param col Row stride 479 480 @return An error code: 0 - success, otherwise - failure 481 482 @ref Developer 483 **/ 484 static int CeedHouseholderApplyQ(CeedScalar *A, const CeedScalar *Q, 485 const CeedScalar *tau, CeedTransposeMode tmode, 486 CeedInt m, CeedInt n, CeedInt k, 487 CeedInt row, CeedInt col) { 488 CeedScalar v[m]; 489 for (CeedInt ii=0; ii<k; ii++) { 490 CeedInt i = tmode == CEED_TRANSPOSE ? ii : k-1-ii; 491 for (CeedInt j=i+1; j<m; j++) { 492 v[j] = Q[j*k+i]; 493 } 494 // Apply Householder reflector (I - tau v v^T) colograd1d^T 495 CeedHouseholderReflect(&A[i*row], &v[i], tau[i], m-i, n, row, col); 496 } 497 return 0; 498 } 499 500 /** 501 @brief Return QR Factorization of matrix 502 503 @param[out] mat Row-major matrix to be factorized in place 504 @param[out] tau Vector of length m of scaling factors 505 @param m Number of rows 506 @param n Number of columns 507 508 @return An error code: 0 - success, otherwise - failure 509 510 @ref Utility 511 **/ 512 int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau, 513 CeedInt m, CeedInt n) { 514 CeedInt i, j; 515 CeedScalar v[m]; 516 517 // Check m >= n 518 if (n > m) 519 return CeedError(ceed, 1, "Cannot compute QR factorization with n > m"); 520 521 for (i=0; i<n; i++) { 522 // Calculate Householder vector, magnitude 523 CeedScalar sigma = 0.0; 524 v[i] = mat[i+n*i]; 525 for (j=i+1; j<m; j++) { 526 v[j] = mat[i+n*j]; 527 sigma += v[j] * v[j]; 528 } 529 CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m] 530 CeedScalar Rii = -copysign(norm, v[i]); 531 v[i] -= Rii; 532 // norm of v[i:m] after modification above and scaling below 533 // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 534 // tau = 2 / (norm*norm) 535 tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma); 536 for (j=i+1; j<m; j++) v[j] /= v[i]; 537 538 // Apply Householder reflector to lower right panel 539 CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1); 540 // Save v 541 mat[i+n*i] = Rii; 542 for (j=i+1; j<m; j++) { 543 mat[i+n*j] = v[j]; 544 } 545 } 546 547 return 0; 548 } 549 550 /** 551 @brief Return collocated grad matrix 552 553 @param basis CeedBasis 554 @param[out] colograd1d Row-major Q1d × Q1d matrix expressing derivatives of 555 basis functions at quadrature points 556 557 @return An error code: 0 - success, otherwise - failure 558 559 @ref Advanced 560 **/ 561 int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *colograd1d) { 562 int i, j, k; 563 Ceed ceed; 564 CeedInt ierr, P1d=(basis)->P1d, Q1d=(basis)->Q1d; 565 CeedScalar *interp1d, *grad1d, tau[Q1d]; 566 567 ierr = CeedMalloc(Q1d*P1d, &interp1d); CeedChk(ierr); 568 ierr = CeedMalloc(Q1d*P1d, &grad1d); CeedChk(ierr); 569 memcpy(interp1d, (basis)->interp1d, Q1d*P1d*sizeof(basis)->interp1d[0]); 570 memcpy(grad1d, (basis)->grad1d, Q1d*P1d*sizeof(basis)->interp1d[0]); 571 572 // QR Factorization, interp1d = Q R 573 ierr = CeedBasisGetCeed(basis, &ceed); CeedChk(ierr); 574 ierr = CeedQRFactorization(ceed, interp1d, tau, Q1d, P1d); CeedChk(ierr); 575 576 // Apply Rinv, colograd1d = grad1d Rinv 577 for (i=0; i<Q1d; i++) { // Row i 578 colograd1d[Q1d*i] = grad1d[P1d*i]/interp1d[0]; 579 for (j=1; j<P1d; j++) { // Column j 580 colograd1d[j+Q1d*i] = grad1d[j+P1d*i]; 581 for (k=0; k<j; k++) { 582 colograd1d[j+Q1d*i] -= interp1d[j+P1d*k]*colograd1d[k+Q1d*i]; 583 } 584 colograd1d[j+Q1d*i] /= interp1d[j+P1d*j]; 585 } 586 for (j=P1d; j<Q1d; j++) { 587 colograd1d[j+Q1d*i] = 0; 588 } 589 } 590 591 // Apply Qtranspose, colograd = colograd Qtranspose 592 CeedHouseholderApplyQ(colograd1d, interp1d, tau, CEED_NOTRANSPOSE, 593 Q1d, Q1d, P1d, 1, Q1d); 594 595 ierr = CeedFree(&interp1d); CeedChk(ierr); 596 ierr = CeedFree(&grad1d); CeedChk(ierr); 597 598 return 0; 599 } 600 601 /** 602 @brief Apply basis evaluation from nodes to quadrature points or vice-versa 603 604 @param basis CeedBasis to evaluate 605 @param nelem The number of elements to apply the basis evaluation to; 606 the backend will specify the ordering in 607 ElemRestrictionCreateBlocked 608 @param tmode \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature 609 points, \ref CEED_TRANSPOSE to apply the transpose, mapping 610 from quadrature points to nodes 611 @param emode \ref CEED_EVAL_INTERP to obtain interpolated values, 612 \ref CEED_EVAL_GRAD to obtain gradients. 613 @param[in] u Input array 614 @param[out] v Output array 615 616 @return An error code: 0 - success, otherwise - failure 617 618 @ref Advanced 619 **/ 620 int CeedBasisApply(CeedBasis basis, CeedInt nelem, CeedTransposeMode tmode, 621 CeedEvalMode emode, CeedVector u, CeedVector v) { 622 int ierr; 623 CeedInt ulength = 0, vlength, nnodes, nqpt; 624 if (!basis->Apply) return CeedError(basis->ceed, 1, 625 "Backend does not support BasisApply"); 626 // check compatibility of topological and geometrical dimensions 627 ierr = CeedBasisGetNumNodes(basis, &nnodes); CeedChk(ierr); 628 ierr = CeedBasisGetNumQuadraturePoints(basis, &nqpt); CeedChk(ierr); 629 ierr = CeedVectorGetLength(v, &vlength); CeedChk(ierr); 630 631 if (u) { 632 ierr = CeedVectorGetLength(u, &ulength); CeedChk(ierr); 633 } 634 635 if ((tmode == CEED_TRANSPOSE && (vlength % nnodes != 0 636 || ulength % nqpt != 0)) 637 || 638 (tmode == CEED_NOTRANSPOSE && (ulength % nnodes != 0 || vlength % nqpt != 0))) 639 return CeedError(basis->ceed, 1, 640 "Length of input/output vectors incompatible with basis dimensions"); 641 642 ierr = basis->Apply(basis, nelem, tmode, emode, u, v); CeedChk(ierr); 643 return 0; 644 } 645 646 /** 647 @brief Get Ceed associated with a CeedBasis 648 649 @param basis CeedBasis 650 @param[out] ceed Variable to store Ceed 651 652 @return An error code: 0 - success, otherwise - failure 653 654 @ref Advanced 655 **/ 656 int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) { 657 *ceed = basis->ceed; 658 659 return 0; 660 }; 661 662 /** 663 @brief Get dimension for given CeedBasis 664 665 @param basis CeedBasis 666 @param[out] dim Variable to store dimension of basis 667 668 @return An error code: 0 - success, otherwise - failure 669 670 @ref Advanced 671 **/ 672 int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) { 673 *dim = basis->dim; 674 675 return 0; 676 }; 677 678 /** 679 @brief Get tensor status for given CeedBasis 680 681 @param basis CeedBasis 682 @param[out] tensor Variable to store tensor status 683 684 @return An error code: 0 - success, otherwise - failure 685 686 @ref Advanced 687 **/ 688 int CeedBasisGetTensorStatus(CeedBasis basis, bool *tensor) { 689 *tensor = basis->tensorbasis; 690 691 return 0; 692 }; 693 694 /** 695 @brief Get number of components for given CeedBasis 696 697 @param basis CeedBasis 698 @param[out] dim Variable to store number of components of basis 699 700 @return An error code: 0 - success, otherwise - failure 701 702 @ref Advanced 703 **/ 704 int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *numcomp) { 705 *numcomp = basis->ncomp; 706 707 return 0; 708 }; 709 710 /** 711 @brief Get total number of nodes (in 1 dimension) of a CeedBasis 712 713 @param basis CeedBasis 714 @param[out] P1d Variable to store number of nodes 715 716 @return An error code: 0 - success, otherwise - failure 717 718 @ref Advanced 719 **/ 720 int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P1d) { 721 if (!basis->tensorbasis) return CeedError(basis->ceed, 1, 722 "Cannot supply P1d for non-tensor basis"); 723 *P1d = basis->P1d; 724 return 0; 725 } 726 727 /** 728 @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis 729 730 @param basis CeedBasis 731 @param[out] Q1d Variable to store number of quadrature points 732 733 @return An error code: 0 - success, otherwise - failure 734 735 @ref Advanced 736 **/ 737 int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q1d) { 738 if (!basis->tensorbasis) return CeedError(basis->ceed, 1, 739 "Cannot supply Q1d for non-tensor basis"); 740 *Q1d = basis->Q1d; 741 return 0; 742 } 743 744 /** 745 @brief Get total number of nodes (in dim dimensions) of a CeedBasis 746 747 @param basis CeedBasis 748 @param[out] P Variable to store number of nodes 749 750 @return An error code: 0 - success, otherwise - failure 751 752 @ref Utility 753 **/ 754 int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) { 755 *P = basis->P; 756 return 0; 757 } 758 759 /** 760 @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis 761 762 @param basis CeedBasis 763 @param[out] Q Variable to store number of quadrature points 764 765 @return An error code: 0 - success, otherwise - failure 766 767 @ref Utility 768 **/ 769 int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) { 770 *Q = basis->Q; 771 return 0; 772 } 773 774 /** 775 @brief Get reference coordinates of quadrature points (in dim dimensions) 776 of a CeedBasis 777 778 @param basis CeedBasis 779 @param[out] qref Variable to store reference coordinates of quadrature points 780 781 @return An error code: 0 - success, otherwise - failure 782 783 @ref Advanced 784 **/ 785 int CeedBasisGetQRef(CeedBasis basis, CeedScalar* *qref) { 786 *qref = basis->qref1d; 787 return 0; 788 } 789 790 /** 791 @brief Get quadrature weights of quadrature points (in dim dimensions) 792 of a CeedBasis 793 794 @param basis CeedBasis 795 @param[out] qweight Variable to store quadrature weights 796 797 @return An error code: 0 - success, otherwise - failure 798 799 @ref Advanced 800 **/ 801 int CeedBasisGetQWeights(CeedBasis basis, CeedScalar* *qweight) { 802 *qweight = basis->qweight1d; 803 return 0; 804 } 805 806 /** 807 @brief Get interpolation matrix of a CeedBasis 808 809 @param basis CeedBasis 810 @param[out] qref Variable to store interpolation matrix 811 812 @return An error code: 0 - success, otherwise - failure 813 814 @ref Advanced 815 **/ 816 int CeedBasisGetInterp(CeedBasis basis, CeedScalar* *interp) { 817 *interp = basis->interp1d; 818 return 0; 819 } 820 821 /** 822 @brief Get gradient matrix of a CeedBasis 823 824 @param basis CeedBasis 825 @param[out] qref Variable to store gradient matrix 826 827 @return An error code: 0 - success, otherwise - failure 828 829 @ref Advanced 830 **/ 831 int CeedBasisGetGrad(CeedBasis basis, CeedScalar* *grad) { 832 *grad = basis->grad1d; 833 return 0; 834 } 835 836 /** 837 @brief Get backend data of a CeedBasis 838 839 @param basis CeedBasis 840 @param[out] data Variable to store data 841 842 @return An error code: 0 - success, otherwise - failure 843 844 @ref Advanced 845 **/ 846 int CeedBasisGetData(CeedBasis basis, void* *data) { 847 *data = basis->data; 848 return 0; 849 } 850 851 /** 852 @brief Set backend data of a CeedBasis 853 854 @param[out] basis CeedBasis 855 @param data Data to set 856 857 @return An error code: 0 - success, otherwise - failure 858 859 @ref Advanced 860 **/ 861 int CeedBasisSetData(CeedBasis basis, void* *data) { 862 basis->data = *data; 863 return 0; 864 } 865 866 /** 867 @brief Get CeedTensorContract of a CeedBasis 868 869 @param basis CeedBasis 870 @param[out] contract Variable to store CeedTensorContract 871 872 @return An error code: 0 - success, otherwise - failure 873 874 @ref Advanced 875 **/ 876 int CeedBasisGetTensorContract(CeedBasis basis, 877 CeedTensorContract *contract) { 878 *contract = basis->contract; 879 return 0; 880 } 881 882 /** 883 @brief Set CeedTensorContract of a CeedBasis 884 885 @param[out] basis CeedBasis 886 @param contract CeedTensorContract to set 887 888 @return An error code: 0 - success, otherwise - failure 889 890 @ref Advanced 891 **/ 892 int CeedBasisSetTensorContract(CeedBasis basis, 893 CeedTensorContract *contract) { 894 basis->contract = *contract; 895 return 0; 896 } 897 898 /** 899 @brief Get dimension for given CeedElemTopology 900 901 @param topo CeedElemTopology 902 @param[out] dim Variable to store dimension of topology 903 904 @return An error code: 0 - success, otherwise - failure 905 906 @ref Advanced 907 **/ 908 int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) { 909 *dim = (CeedInt) topo >> 16; 910 911 return 0; 912 }; 913 914 /** 915 @brief Destroy a CeedBasis 916 917 @param basis CeedBasis to destroy 918 919 @return An error code: 0 - success, otherwise - failure 920 921 @ref Basic 922 **/ 923 int CeedBasisDestroy(CeedBasis *basis) { 924 int ierr; 925 926 if (!*basis || --(*basis)->refcount > 0) return 0; 927 if ((*basis)->Destroy) { 928 ierr = (*basis)->Destroy(*basis); CeedChk(ierr); 929 } 930 ierr = CeedFree(&(*basis)->interp1d); CeedChk(ierr); 931 ierr = CeedFree(&(*basis)->grad1d); CeedChk(ierr); 932 ierr = CeedFree(&(*basis)->qref1d); CeedChk(ierr); 933 ierr = CeedFree(&(*basis)->qweight1d); CeedChk(ierr); 934 ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr); 935 ierr = CeedFree(basis); CeedChk(ierr); 936 return 0; 937 } 938 939 /// @cond DOXYGEN_SKIP 940 // Indicate that the quadrature points are collocated with the nodes 941 CeedBasis CEED_BASIS_COLLOCATED = &ceed_basis_collocated; 942 /// @endcond 943 /// @} 944