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