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