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