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 // LCOV_EXCL_START 65 return CeedError(ceed, 1, "Basis dimension must be a positive value"); 66 // LCOV_EXCL_STOP 67 68 if (!ceed->BasisCreateTensorH1) { 69 Ceed delegate; 70 ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr); 71 72 if (!delegate) 73 // LCOV_EXCL_START 74 return CeedError(ceed, 1, "Backend does not support BasisCreateTensorH1"); 75 // LCOV_EXCL_STOP 76 77 ierr = CeedBasisCreateTensorH1(delegate, dim, ncomp, P1d, 78 Q1d, interp1d, grad1d, qref1d, 79 qweight1d, basis); CeedChk(ierr); 80 return 0; 81 } 82 ierr = CeedCalloc(1,basis); CeedChk(ierr); 83 (*basis)->ceed = ceed; 84 ceed->refcount++; 85 (*basis)->refcount = 1; 86 (*basis)->tensorbasis = 1; 87 (*basis)->dim = dim; 88 (*basis)->ncomp = ncomp; 89 (*basis)->P1d = P1d; 90 (*basis)->Q1d = Q1d; 91 (*basis)->P = CeedIntPow(P1d, dim); 92 (*basis)->Q = CeedIntPow(Q1d, dim); 93 ierr = CeedMalloc(Q1d,&(*basis)->qref1d); CeedChk(ierr); 94 ierr = CeedMalloc(Q1d,&(*basis)->qweight1d); CeedChk(ierr); 95 memcpy((*basis)->qref1d, qref1d, Q1d*sizeof(qref1d[0])); 96 memcpy((*basis)->qweight1d, qweight1d, Q1d*sizeof(qweight1d[0])); 97 ierr = CeedMalloc(Q1d*P1d,&(*basis)->interp1d); CeedChk(ierr); 98 ierr = CeedMalloc(Q1d*P1d,&(*basis)->grad1d); CeedChk(ierr); 99 memcpy((*basis)->interp1d, interp1d, Q1d*P1d*sizeof(interp1d[0])); 100 memcpy((*basis)->grad1d, grad1d, Q1d*P1d*sizeof(grad1d[0])); 101 ierr = ceed->BasisCreateTensorH1(dim, P1d, Q1d, interp1d, grad1d, qref1d, 102 qweight1d, *basis); CeedChk(ierr); 103 return 0; 104 } 105 106 /** 107 @brief Create a tensor product Lagrange basis 108 109 @param ceed A Ceed object where the CeedBasis will be created 110 @param dim Topological dimension of element 111 @param ncomp Number of field components 112 @param P Number of Gauss-Lobatto nodes in one dimension. The 113 polynomial degree of the resulting Q_k element is k=P-1. 114 @param Q Number of quadrature points in one dimension. 115 @param qmode Distribution of the Q quadrature points (affects order of 116 accuracy for the quadrature) 117 @param[out] basis Address of the variable where the newly created 118 CeedBasis will be stored. 119 120 @return An error code: 0 - success, otherwise - failure 121 122 @ref Basic 123 **/ 124 int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt ncomp, 125 CeedInt P, CeedInt Q, 126 CeedQuadMode qmode, CeedBasis *basis) { 127 // Allocate 128 int ierr, i, j, k; 129 CeedScalar c1, c2, c3, c4, dx, *nodes, *interp1d, *grad1d, *qref1d, *qweight1d; 130 131 if (dim<1) 132 // LCOV_EXCL_START 133 return CeedError(ceed, 1, "Basis dimension must be a positive value"); 134 // LCOV_EXCL_STOP 135 136 ierr = CeedCalloc(P*Q, &interp1d); CeedChk(ierr); 137 ierr = CeedCalloc(P*Q, &grad1d); CeedChk(ierr); 138 ierr = CeedCalloc(P, &nodes); CeedChk(ierr); 139 ierr = CeedCalloc(Q, &qref1d); CeedChk(ierr); 140 ierr = CeedCalloc(Q, &qweight1d); CeedChk(ierr); 141 // Get Nodes and Weights 142 ierr = CeedLobattoQuadrature(P, nodes, NULL); CeedChk(ierr); 143 switch (qmode) { 144 case CEED_GAUSS: 145 ierr = CeedGaussQuadrature(Q, qref1d, qweight1d); CeedChk(ierr); 146 break; 147 case CEED_GAUSS_LOBATTO: 148 ierr = CeedLobattoQuadrature(Q, qref1d, qweight1d); CeedChk(ierr); 149 break; 150 } 151 // Build B, D matrix 152 // Fornberg, 1998 153 for (i = 0; i < Q; i++) { 154 c1 = 1.0; 155 c3 = nodes[0] - qref1d[i]; 156 interp1d[i*P+0] = 1.0; 157 for (j = 1; j < P; j++) { 158 c2 = 1.0; 159 c4 = c3; 160 c3 = nodes[j] - qref1d[i]; 161 for (k = 0; k < j; k++) { 162 dx = nodes[j] - nodes[k]; 163 c2 *= dx; 164 if (k == j - 1) { 165 grad1d[i*P + j] = c1*(interp1d[i*P + k] - c4*grad1d[i*P + k]) / c2; 166 interp1d[i*P + j] = - c1*c4*interp1d[i*P + k] / c2; 167 } 168 grad1d[i*P + k] = (c3*grad1d[i*P + k] - interp1d[i*P + k]) / dx; 169 interp1d[i*P + k] = c3*interp1d[i*P + k] / dx; 170 } 171 c1 = c2; 172 } 173 } 174 // // Pass to CeedBasisCreateTensorH1 175 ierr = CeedBasisCreateTensorH1(ceed, dim, ncomp, P, Q, interp1d, grad1d, qref1d, 176 qweight1d, basis); CeedChk(ierr); 177 ierr = CeedFree(&interp1d); CeedChk(ierr); 178 ierr = CeedFree(&grad1d); CeedChk(ierr); 179 ierr = CeedFree(&nodes); CeedChk(ierr); 180 ierr = CeedFree(&qref1d); CeedChk(ierr); 181 ierr = CeedFree(&qweight1d); CeedChk(ierr); 182 return 0; 183 } 184 185 /** 186 @brief Create a non tensor product basis for H^1 discretizations 187 188 @param ceed A Ceed object where the CeedBasis will be created 189 @param topo Topology of element, e.g. hypercube, simplex, ect 190 @param ncomp Number of field components (1 for scalar fields) 191 @param nnodes Total number of nodes 192 @param nqpts Total number of quadrature points 193 @param interp Row-major nqpts × nnodes matrix expressing the values of 194 nodal basis functions at quadrature points 195 @param grad Row-major (nqpts x dim) × nnodes matrix expressing 196 derivatives of nodal basis functions at quadrature points 197 @param qref Array of length nqpts holding the locations of quadrature 198 points on the reference element [-1, 1] 199 @param qweight Array of length nqpts holding the quadrature weights on the 200 reference element 201 @param[out] basis Address of the variable where the newly created 202 CeedBasis will be stored. 203 204 @return An error code: 0 - success, otherwise - failure 205 206 @ref Basic 207 **/ 208 int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt ncomp, 209 CeedInt nnodes, CeedInt nqpts, 210 const CeedScalar *interp, 211 const CeedScalar *grad, const CeedScalar *qref, 212 const CeedScalar *qweight, CeedBasis *basis) { 213 int ierr; 214 CeedInt P = nnodes, Q = nqpts, dim = 0; 215 216 if (!ceed->BasisCreateH1) { 217 Ceed delegate; 218 ierr = CeedGetObjectDelegate(ceed, &delegate, "Basis"); CeedChk(ierr); 219 220 if (!delegate) 221 // LCOV_EXCL_START 222 return CeedError(ceed, 1, "Backend does not support BasisCreateH1"); 223 // LCOV_EXCL_STOP 224 225 ierr = CeedBasisCreateH1(delegate, topo, ncomp, nnodes, 226 nqpts, interp, grad, qref, 227 qweight, basis); CeedChk(ierr); 228 return 0; 229 } 230 231 ierr = CeedCalloc(1,basis); CeedChk(ierr); 232 233 ierr = CeedBasisGetTopologyDimension(topo, &dim); CeedChk(ierr); 234 235 (*basis)->ceed = ceed; 236 ceed->refcount++; 237 (*basis)->refcount = 1; 238 (*basis)->tensorbasis = 0; 239 (*basis)->dim = dim; 240 (*basis)->ncomp = ncomp; 241 (*basis)->P = P; 242 (*basis)->Q = Q; 243 ierr = CeedMalloc(Q*dim,&(*basis)->qref1d); CeedChk(ierr); 244 ierr = CeedMalloc(Q,&(*basis)->qweight1d); CeedChk(ierr); 245 memcpy((*basis)->qref1d, qref, Q*dim*sizeof(qref[0])); 246 memcpy((*basis)->qweight1d, qweight, Q*sizeof(qweight[0])); 247 ierr = CeedMalloc(Q*P,&(*basis)->interp1d); CeedChk(ierr); 248 ierr = CeedMalloc(dim*Q*P,&(*basis)->grad1d); CeedChk(ierr); 249 memcpy((*basis)->interp1d, interp, Q*P*sizeof(interp[0])); 250 memcpy((*basis)->grad1d, grad, dim*Q*P*sizeof(grad[0])); 251 ierr = ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, qref, 252 qweight, *basis); CeedChk(ierr); 253 return 0; 254 } 255 256 /** 257 @brief Construct a Gauss-Legendre quadrature 258 259 @param Q Number of quadrature points (integrates polynomials of 260 degree 2*Q-1 exactly) 261 @param[out] qref1d Array of length Q to hold the abscissa on [-1, 1] 262 @param[out] qweight1d Array of length Q to hold the weights 263 264 @return An error code: 0 - success, otherwise - failure 265 266 @ref Utility 267 **/ 268 int CeedGaussQuadrature(CeedInt Q, CeedScalar *qref1d, CeedScalar *qweight1d) { 269 // Allocate 270 CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0*atan(1.0); 271 // Build qref1d, qweight1d 272 for (int i = 0; i <= Q/2; i++) { 273 // Guess 274 xi = cos(PI*(CeedScalar)(2*i+1)/((CeedScalar)(2*Q))); 275 // Pn(xi) 276 P0 = 1.0; 277 P1 = xi; 278 P2 = 0.0; 279 for (int j = 2; j <= Q; j++) { 280 P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 281 P0 = P1; 282 P1 = P2; 283 } 284 // First Newton Step 285 dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 286 xi = xi-P2/dP2; 287 // Newton to convergence 288 for (int k=0; k<100 && fabs(P2)>1e-15; k++) { 289 P0 = 1.0; 290 P1 = xi; 291 for (int j = 2; j <= Q; j++) { 292 P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 293 P0 = P1; 294 P1 = P2; 295 } 296 dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 297 xi = xi-P2/dP2; 298 } 299 // Save xi, wi 300 wi = 2.0/((1.0-xi*xi)*dP2*dP2); 301 qweight1d[i] = wi; 302 qweight1d[Q-1-i] = wi; 303 qref1d[i] = -xi; 304 qref1d[Q-1-i]= xi; 305 } 306 return 0; 307 } 308 309 /** 310 @brief Construct a Gauss-Legendre-Lobatto quadrature 311 312 @param Q Number of quadrature points (integrates polynomials of 313 degree 2*Q-3 exactly) 314 @param[out] qref1d Array of length Q to hold the abscissa on [-1, 1] 315 @param[out] qweight1d Array of length Q to hold the weights 316 317 @return An error code: 0 - success, otherwise - failure 318 319 @ref Utility 320 **/ 321 int CeedLobattoQuadrature(CeedInt Q, CeedScalar *qref1d, 322 CeedScalar *qweight1d) { 323 // Allocate 324 CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0*atan(1.0); 325 // Build qref1d, qweight1d 326 // Set endpoints 327 wi = 2.0/((CeedScalar)(Q*(Q-1))); 328 if (qweight1d) { 329 qweight1d[0] = wi; 330 qweight1d[Q-1] = wi; 331 } 332 qref1d[0] = -1.0; 333 qref1d[Q-1] = 1.0; 334 // Interior 335 for (int i = 1; i <= (Q-1)/2; i++) { 336 // Guess 337 xi = cos(PI*(CeedScalar)(i)/(CeedScalar)(Q-1)); 338 // Pn(xi) 339 P0 = 1.0; 340 P1 = xi; 341 P2 = 0.0; 342 for (int j = 2; j < Q; j++) { 343 P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 344 P0 = P1; 345 P1 = P2; 346 } 347 // First Newton step 348 dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 349 d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi); 350 xi = xi-dP2/d2P2; 351 // Newton to convergence 352 for (int k=0; k<100 && fabs(dP2)>1e-15; k++) { 353 P0 = 1.0; 354 P1 = xi; 355 for (int j = 2; j < Q; j++) { 356 P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 357 P0 = P1; 358 P1 = P2; 359 } 360 dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 361 d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi); 362 xi = xi-dP2/d2P2; 363 } 364 // Save xi, wi 365 wi = 2.0/(((CeedScalar)(Q*(Q-1)))*P2*P2); 366 if (qweight1d) { 367 qweight1d[i] = wi; 368 qweight1d[Q-1-i] = wi; 369 } 370 qref1d[i] = -xi; 371 qref1d[Q-1-i]= xi; 372 } 373 return 0; 374 } 375 376 /** 377 @brief View an array stored in a CeedBasis 378 379 @param name Name of array 380 @param fpformat Printing format 381 @param m Number of rows in array 382 @param n Number of columns in array 383 @param a Array to be viewed 384 @param stream Stream to view to, e.g., stdout 385 386 @return An error code: 0 - success, otherwise - failure 387 388 @ref Utility 389 **/ 390 static int CeedScalarView(const char *name, const char *fpformat, CeedInt m, 391 CeedInt n, const CeedScalar *a, FILE *stream) { 392 for (int i=0; i<m; i++) { 393 if (m > 1) 394 fprintf(stream, "%12s[%d]:", name, i); 395 else 396 fprintf(stream, "%12s:", name); 397 for (int j=0; j<n; j++) 398 fprintf(stream, fpformat, fabs(a[i*n+j]) > 1E-14 ? a[i*n+j] : 0); 399 fputs("\n", stream); 400 } 401 return 0; 402 } 403 404 /** 405 @brief View a CeedBasis 406 407 @param basis CeedBasis to view 408 @param stream Stream to view to, e.g., stdout 409 410 @return An error code: 0 - success, otherwise - failure 411 412 @ref Utility 413 **/ 414 int CeedBasisView(CeedBasis basis, FILE *stream) { 415 int ierr; 416 417 if (basis->tensorbasis) { 418 fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P1d, 419 basis->Q1d); 420 ierr = CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q1d, basis->qref1d, 421 stream); CeedChk(ierr); 422 ierr = CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q1d, 423 basis->qweight1d, stream); CeedChk(ierr); 424 ierr = CeedScalarView("interp1d", "\t% 12.8f", basis->Q1d, basis->P1d, 425 basis->interp1d, stream); CeedChk(ierr); 426 ierr = CeedScalarView("grad1d", "\t% 12.8f", basis->Q1d, basis->P1d, 427 basis->grad1d, stream); CeedChk(ierr); 428 } else { 429 fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P, 430 basis->Q); 431 ierr = CeedScalarView("qref", "\t% 12.8f", 1, basis->Q*basis->dim, 432 basis->qref1d, 433 stream); CeedChk(ierr); 434 ierr = CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->qweight1d, 435 stream); CeedChk(ierr); 436 ierr = CeedScalarView("interp", "\t% 12.8f", basis->Q, basis->P, 437 basis->interp1d, stream); CeedChk(ierr); 438 ierr = CeedScalarView("grad", "\t% 12.8f", basis->dim*basis->Q, basis->P, 439 basis->grad1d, stream); CeedChk(ierr); 440 } 441 return 0; 442 } 443 444 /** 445 @brief Compute Householder reflection 446 447 Computes A = (I - b v v^T) A 448 where A is an mxn matrix indexed as A[i*row + j*col] 449 450 @param[in,out] A Matrix to apply Householder reflection to, in place 451 @param v Householder vector 452 @param b Scaling factor 453 @param m Number of rows in A 454 @param n Number of columns in A 455 @param row Row stride 456 @param col Col stride 457 458 @return An error code: 0 - success, otherwise - failure 459 460 @ref Developer 461 **/ 462 static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v, 463 CeedScalar b, CeedInt m, CeedInt n, 464 CeedInt row, CeedInt col) { 465 for (CeedInt j=0; j<n; j++) { 466 CeedScalar w = A[0*row + j*col]; 467 for (CeedInt i=1; i<m; i++) 468 w += v[i] * A[i*row + j*col]; 469 A[0*row + j*col] -= b * w; 470 for (CeedInt i=1; i<m; i++) 471 A[i*row + j*col] -= b * w * v[i]; 472 } 473 return 0; 474 } 475 476 /** 477 @brief Apply Householder Q matrix 478 479 Compute A = Q A where Q is mxm and A is mxn. 480 481 @param[in,out] A Matrix to apply Householder Q to, in place 482 @param Q Householder Q matrix 483 @param tau Householder scaling factors 484 @param tmode Transpose mode for application 485 @param m Number of rows in A 486 @param n Number of columns in A 487 @param k Number of elementary reflectors in Q, k<m 488 @param row Row stride in A 489 @param col Col stride in A 490 491 @return An error code: 0 - success, otherwise - failure 492 493 @ref Developer 494 **/ 495 static int CeedHouseholderApplyQ(CeedScalar *A, const CeedScalar *Q, 496 const CeedScalar *tau, CeedTransposeMode tmode, 497 CeedInt m, CeedInt n, CeedInt k, 498 CeedInt row, CeedInt col) { 499 CeedScalar v[m]; 500 for (CeedInt ii=0; ii<k; ii++) { 501 CeedInt i = tmode == CEED_TRANSPOSE ? ii : k-1-ii; 502 for (CeedInt j=i+1; j<m; j++) 503 v[j] = Q[j*k+i]; 504 // Apply Householder reflector (I - tau v v^T) colograd1d^T 505 CeedHouseholderReflect(&A[i*row], &v[i], tau[i], m-i, n, row, col); 506 } 507 return 0; 508 } 509 510 /** 511 @brief Compute Givens rotation 512 513 Computes A = G A (or G^T A in transpose mode) 514 where A is an mxn matrix indexed as A[i*n + j*m] 515 516 @param[in,out] A Row major matrix to apply Givens rotation to, in place 517 @param c Cosine factor 518 @param s Sine factor 519 @param i First row/column to apply rotation 520 @param k Second row/column to apply rotation 521 @param m Number of rows in A 522 @param n Number of columns in A 523 524 @return An error code: 0 - success, otherwise - failure 525 526 @ref Developer 527 **/ 528 static int CeedGivensRotation(CeedScalar *A, CeedScalar c, CeedScalar s, 529 CeedTransposeMode tmode, CeedInt i, CeedInt k, 530 CeedInt m, CeedInt n) { 531 CeedInt stridej = 1, strideik = m, numits = n; 532 if (tmode == CEED_NOTRANSPOSE) { 533 stridej = n; strideik = 1; numits = m; 534 } 535 536 // Apply rotation 537 for (CeedInt j=0; j<numits; j++) { 538 CeedScalar tau1 = A[i*strideik+j*stridej], tau2 = A[k*strideik+j*stridej]; 539 A[i*strideik+j*stridej] = c*tau1 - s*tau2; 540 A[k*strideik+j*stridej] = s*tau1 + c*tau2; 541 } 542 543 return 0; 544 } 545 546 /** 547 @brief Return QR Factorization of matrix 548 549 @param ceed A Ceed object currently in use 550 @param[in,out] mat Row-major matrix to be factorized in place 551 @param[in,out] tau Vector of length m of scaling factors 552 @param m Number of rows 553 @param n Number of columns 554 555 @return An error code: 0 - success, otherwise - failure 556 557 @ref Utility 558 **/ 559 int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau, 560 CeedInt m, CeedInt n) { 561 CeedScalar v[m]; 562 563 // Check m >= n 564 if (n > m) 565 // LCOV_EXCL_START 566 return CeedError(ceed, 1, "Cannot compute QR factorization with n > m"); 567 // LCOV_EXCL_STOP 568 569 for (CeedInt i=0; i<n; i++) { 570 // Calculate Householder vector, magnitude 571 CeedScalar sigma = 0.0; 572 v[i] = mat[i+n*i]; 573 for (CeedInt j=i+1; j<m; j++) { 574 v[j] = mat[i+n*j]; 575 sigma += v[j] * v[j]; 576 } 577 CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m] 578 CeedScalar Rii = -copysign(norm, v[i]); 579 v[i] -= Rii; 580 // norm of v[i:m] after modification above and scaling below 581 // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 582 // tau = 2 / (norm*norm) 583 tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma); 584 for (CeedInt j=i+1; j<m; j++) 585 v[j] /= v[i]; 586 587 // Apply Householder reflector to lower right panel 588 CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1); 589 // Save v 590 mat[i+n*i] = Rii; 591 for (CeedInt j=i+1; j<m; j++) 592 mat[i+n*j] = v[j]; 593 } 594 595 return 0; 596 } 597 598 /** 599 @brief Return symmetric Schur decomposition of the symmetric matrix mat via 600 symmetric QR factorization 601 602 @param[in,out] mat Row-major matrix to be factorized in place 603 @param[out] lambda Vector of length m of eigenvalues 604 @param n Number of rows/columns 605 606 @return An error code: 0 - success, otherwise - failure 607 608 @ref Utility 609 **/ 610 int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat, 611 CeedScalar *lambda, CeedInt n) { 612 // Check bounds for clang-tidy 613 if (n<2) 614 // LCOV_EXCL_START 615 return CeedError(ceed, 1, 616 "Cannot compute symmetric Schur decomposition of scalars"); 617 // LCOV_EXCL_STOP 618 619 CeedScalar v[n-1], tau[n-1], matT[n*n]; 620 621 // Copy mat to matT and set mat to I 622 memcpy(matT, mat, n*n*sizeof(mat[0])); 623 for (CeedInt i=0; i<n; i++) 624 for (CeedInt j=0; j<n; j++) 625 mat[j+n*i] = (i==j) ? 1 : 0; 626 627 // Reduce to tridiagonal 628 for (CeedInt i=0; i<n-1; i++) { 629 // Calculate Householder vector, magnitude 630 CeedScalar sigma = 0.0; 631 v[i] = matT[i+n*(i+1)]; 632 for (CeedInt j=i+1; j<n-1; j++) { 633 v[j] = matT[i+n*(j+1)]; 634 sigma += v[j] * v[j]; 635 } 636 CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:n-1] 637 CeedScalar Rii = -copysign(norm, v[i]); 638 v[i] -= Rii; 639 // norm of v[i:m] after modification above and scaling below 640 // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 641 // tau = 2 / (norm*norm) 642 tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma); 643 for (CeedInt j=i+1; j<n-1; j++) v[j] /= v[i]; 644 645 // Update sub and super diagonal 646 matT[i+n*(i+1)] = Rii; 647 matT[(i+1)+n*i] = Rii; 648 for (CeedInt j=i+2; j<n; j++) { 649 matT[i+n*j] = 0; matT[j+n*i] = 0; 650 } 651 // Apply symmetric Householder reflector to lower right panel 652 CeedHouseholderReflect(&matT[(i+1)+n*(i+1)], &v[i], tau[i], 653 n-(i+1), n-(i+1), n, 1); 654 CeedHouseholderReflect(&matT[(i+1)+n*(i+1)], &v[i], tau[i], 655 n-(i+1), n-(i+1), 1, n); 656 // Save v 657 for (CeedInt j=i+1; j<n-1; j++) { 658 matT[i+n*(j+1)] = v[j]; 659 } 660 } 661 // Backwards accumulation of Q 662 for (CeedInt i=n-2; i>=0; i--) { 663 v[i] = 1; 664 for (CeedInt j=i+1; j<n-1; j++) { 665 v[j] = matT[i+n*(j+1)]; 666 matT[i+n*(j+1)] = 0; 667 } 668 CeedHouseholderReflect(&mat[(i+1)+n*(i+1)], &v[i], tau[i], 669 n-(i+1), n-(i+1), n, 1); 670 } 671 672 // Reduce sub and super diagonal 673 CeedInt p = 0, q = 0, itr = 0, maxitr = n*n*n; 674 CeedScalar tol = 1e-15; 675 676 while (q < n && itr < maxitr) { 677 // Update p, q, size of reduced portions of diagonal 678 p = 0; q = 0; 679 for (CeedInt i=n-2; i>=0; i--) { 680 if (fabs(matT[i+n*(i+1)]) < tol) 681 q += 1; 682 else 683 break; 684 } 685 for (CeedInt i=0; i<n-1-q; i++) { 686 if (fabs(matT[i+n*(i+1)]) < tol) 687 p += 1; 688 else 689 break; 690 } 691 if (q == n-1) break; // Finished reducing 692 693 // Reduce tridiagonal portion 694 CeedScalar tnn = matT[(n-1-q)+n*(n-1-q)], 695 tnnm1 = matT[(n-2-q)+n*(n-1-q)]; 696 CeedScalar d = (matT[(n-2-q)+n*(n-2-q)] - tnn)/2; 697 CeedScalar mu = tnn - tnnm1*tnnm1 / 698 (d + copysign(sqrt(d*d + tnnm1*tnnm1), d)); 699 CeedScalar x = matT[p+n*p] - mu; 700 CeedScalar z = matT[p+n*(p+1)]; 701 for (CeedInt k=p; k<n-1-q; k++) { 702 // Compute Givens rotation 703 CeedScalar c = 1, s = 0; 704 if (fabs(z) > tol) { 705 if (fabs(z) > fabs(x)) { 706 CeedScalar tau = -x/z; 707 s = 1/sqrt(1+tau*tau), c = s*tau; 708 } else { 709 CeedScalar tau = -z/x; 710 c = 1/sqrt(1+tau*tau), s = c*tau; 711 } 712 } 713 714 // Apply Givens rotation to T 715 CeedGivensRotation(matT, c, s, CEED_NOTRANSPOSE, k, k+1, n, n); 716 CeedGivensRotation(matT, c, s, CEED_TRANSPOSE, k, k+1, n, n); 717 718 // Apply Givens rotation to Q 719 CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k+1, n, n); 720 721 // Update x, z 722 if (k < n-q-2) { 723 x = matT[k+n*(k+1)]; 724 z = matT[k+n*(k+2)]; 725 } 726 } 727 itr++; 728 } 729 // Save eigenvalues 730 for (CeedInt i=0; i<n; i++) 731 lambda[i] = matT[i+n*i]; 732 733 // Check convergence 734 if (itr == maxitr && q < n-1) 735 // LCOV_EXCL_START 736 return CeedError(ceed, 1, "Symmetric QR failed to converge"); 737 // LCOV_EXCL_STOP 738 739 return 0; 740 } 741 742 /** 743 @brief Return C = A B 744 745 @param[in] matA Row-major matrix A 746 @param[in] matB Row-major matrix B 747 @param[out] matC Row-major output matrix C 748 @param m Number of rows of C 749 @param n Number of columns of C 750 @param kk Number of columns of A/rows of B 751 752 @return An error code: 0 - success, otherwise - failure 753 754 @ref Utility 755 **/ 756 static int CeedMatrixMultiply(Ceed ceed, CeedScalar *matA, CeedScalar *matB, 757 CeedScalar *matC, CeedInt m, CeedInt n, 758 CeedInt kk) { 759 for (CeedInt i=0; i<m; i++) 760 for (CeedInt j=0; j<n; j++) { 761 CeedScalar sum = 0; 762 for (CeedInt k=0; k<kk; k++) 763 sum += matA[k+i*kk]*matB[j+k*n]; 764 matC[j+i*n] = sum; 765 } 766 return 0; 767 } 768 769 /** 770 @brief Return Simultaneous Diagonalization of two matrices. This solves the 771 generalized eigenvalue problem A x = lambda B x, where A and B 772 are symmetric and B is positive definite. We generate the matrix X 773 and vector Lambda such that X^T A X = Lambda and X^T B X = I. This 774 is equivalent to the LAPACK routine 'sygv' with TYPE = 1. 775 776 @param[in] matA Row-major matrix to be factorized with eigenvalues 777 @param[in] matB Row-major matrix to be factorized to identity 778 @param[out] x Row-major orthogonal matrix 779 @param[out] lambda Vector of length m of generalized eigenvalues 780 @param n Number of rows/columns 781 782 @return An error code: 0 - success, otherwise - failure 783 784 @ref Utility 785 **/ 786 int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *matA, 787 CeedScalar *matB, CeedScalar *x, 788 CeedScalar *lambda, CeedInt n) { 789 int ierr; 790 CeedScalar matC[n*n], matG[n*n], vecD[n]; 791 792 // Compute B = G D G^T 793 memcpy(matG, matB, n*n*sizeof(matB[0])); 794 ierr = CeedSymmetricSchurDecomposition(ceed, matG, vecD, n); CeedChk(ierr); 795 for (CeedInt i=0; i<n; i++) vecD[i] = sqrt(vecD[i]); 796 797 // Compute C = (G D^-1/2)^-1 A (G D^-1/2)^-T 798 // = D^1/2 G^T A D^1/2 G 799 for (CeedInt i=0; i<n; i++) 800 for (CeedInt j=0; j<n; j++) 801 matC[j+i*n] = vecD[i] * matG[i+j*n]; 802 CeedMatrixMultiply(ceed, matC, matA, x, n, n, n); 803 for (CeedInt i=0; i<n; i++) 804 for (CeedInt j=0; j<n; j++) 805 matG[j+i*n] = vecD[i] * matG[j+i*n]; 806 CeedMatrixMultiply(ceed, x, matG, matC, n, n, n); 807 808 // Compute Q^T C Q = lambda 809 ierr = CeedSymmetricSchurDecomposition(ceed, matC, lambda, n); CeedChk(ierr); 810 811 // Set x = (G D^-1/2)^-T Q 812 // = D^1/2 G Q 813 CeedMatrixMultiply(ceed, matG, matC, x, n, n, n); 814 815 return 0; 816 } 817 818 /** 819 @brief Return collocated grad matrix 820 821 @param basis CeedBasis 822 @param[out] colograd1d Row-major Q1d × Q1d matrix expressing derivatives of 823 basis functions at quadrature points 824 825 @return An error code: 0 - success, otherwise - failure 826 827 @ref Advanced 828 **/ 829 int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *colograd1d) { 830 int i, j, k; 831 Ceed ceed; 832 CeedInt ierr, P1d=(basis)->P1d, Q1d=(basis)->Q1d; 833 CeedScalar *interp1d, *grad1d, tau[Q1d]; 834 835 ierr = CeedMalloc(Q1d*P1d, &interp1d); CeedChk(ierr); 836 ierr = CeedMalloc(Q1d*P1d, &grad1d); CeedChk(ierr); 837 memcpy(interp1d, (basis)->interp1d, Q1d*P1d*sizeof(basis)->interp1d[0]); 838 memcpy(grad1d, (basis)->grad1d, Q1d*P1d*sizeof(basis)->interp1d[0]); 839 840 // QR Factorization, interp1d = Q R 841 ierr = CeedBasisGetCeed(basis, &ceed); CeedChk(ierr); 842 ierr = CeedQRFactorization(ceed, interp1d, tau, Q1d, P1d); CeedChk(ierr); 843 844 // Apply Rinv, colograd1d = grad1d Rinv 845 for (i=0; i<Q1d; i++) { // Row i 846 colograd1d[Q1d*i] = grad1d[P1d*i]/interp1d[0]; 847 for (j=1; j<P1d; j++) { // Column j 848 colograd1d[j+Q1d*i] = grad1d[j+P1d*i]; 849 for (k=0; k<j; k++) 850 colograd1d[j+Q1d*i] -= interp1d[j+P1d*k]*colograd1d[k+Q1d*i]; 851 colograd1d[j+Q1d*i] /= interp1d[j+P1d*j]; 852 } 853 for (j=P1d; j<Q1d; j++) 854 colograd1d[j+Q1d*i] = 0; 855 } 856 857 // Apply Qtranspose, colograd = colograd Qtranspose 858 CeedHouseholderApplyQ(colograd1d, interp1d, tau, CEED_NOTRANSPOSE, 859 Q1d, Q1d, P1d, 1, Q1d); 860 861 ierr = CeedFree(&interp1d); CeedChk(ierr); 862 ierr = CeedFree(&grad1d); CeedChk(ierr); 863 864 return 0; 865 } 866 867 /** 868 @brief Apply basis evaluation from nodes to quadrature points or vice-versa 869 870 @param basis CeedBasis to evaluate 871 @param nelem The number of elements to apply the basis evaluation to; 872 the backend will specify the ordering in 873 ElemRestrictionCreateBlocked 874 @param tmode \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature 875 points, \ref CEED_TRANSPOSE to apply the transpose, mapping 876 from quadrature points to nodes 877 @param emode \ref CEED_EVAL_INTERP to obtain interpolated values, 878 \ref CEED_EVAL_GRAD to obtain gradients. 879 @param[in] u Input array 880 @param[out] v Output array 881 882 @return An error code: 0 - success, otherwise - failure 883 884 @ref Advanced 885 **/ 886 int CeedBasisApply(CeedBasis basis, CeedInt nelem, CeedTransposeMode tmode, 887 CeedEvalMode emode, CeedVector u, CeedVector v) { 888 int ierr; 889 CeedInt ulength = 0, vlength, nnodes, nqpt; 890 if (!basis->Apply) 891 // LCOV_EXCL_START 892 return CeedError(basis->ceed, 1, "Backend does not support BasisApply"); 893 // LCOV_EXCL_STOP 894 895 // Check compatibility of topological and geometrical dimensions 896 ierr = CeedBasisGetNumNodes(basis, &nnodes); CeedChk(ierr); 897 ierr = CeedBasisGetNumQuadraturePoints(basis, &nqpt); CeedChk(ierr); 898 ierr = CeedVectorGetLength(v, &vlength); CeedChk(ierr); 899 900 if (u) { 901 ierr = CeedVectorGetLength(u, &ulength); CeedChk(ierr); 902 } 903 904 if ((tmode == CEED_TRANSPOSE && (vlength%nnodes != 0 || ulength%nqpt != 0)) || 905 (tmode == CEED_NOTRANSPOSE && (ulength%nnodes != 0 || vlength%nqpt != 0))) 906 return CeedError(basis->ceed, 1, "Length of input/output vectors " 907 "incompatible with basis dimensions"); 908 909 ierr = basis->Apply(basis, nelem, tmode, emode, u, v); CeedChk(ierr); 910 return 0; 911 } 912 913 /** 914 @brief Get Ceed associated with a CeedBasis 915 916 @param basis CeedBasis 917 @param[out] ceed Variable to store Ceed 918 919 @return An error code: 0 - success, otherwise - failure 920 921 @ref Advanced 922 **/ 923 int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) { 924 *ceed = basis->ceed; 925 return 0; 926 }; 927 928 /** 929 @brief Get dimension for given CeedBasis 930 931 @param basis CeedBasis 932 @param[out] dim Variable to store dimension of basis 933 934 @return An error code: 0 - success, otherwise - failure 935 936 @ref Advanced 937 **/ 938 int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) { 939 *dim = basis->dim; 940 return 0; 941 }; 942 943 /** 944 @brief Get tensor status for given CeedBasis 945 946 @param basis CeedBasis 947 @param[out] tensor Variable to store tensor status 948 949 @return An error code: 0 - success, otherwise - failure 950 951 @ref Advanced 952 **/ 953 int CeedBasisGetTensorStatus(CeedBasis basis, bool *tensor) { 954 *tensor = basis->tensorbasis; 955 return 0; 956 }; 957 958 /** 959 @brief Get number of components for given CeedBasis 960 961 @param basis CeedBasis 962 @param[out] numcomp Variable to store number of components of basis 963 964 @return An error code: 0 - success, otherwise - failure 965 966 @ref Advanced 967 **/ 968 int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *numcomp) { 969 *numcomp = basis->ncomp; 970 return 0; 971 }; 972 973 /** 974 @brief Get total number of nodes (in 1 dimension) of a CeedBasis 975 976 @param basis CeedBasis 977 @param[out] P1d Variable to store number of nodes 978 979 @return An error code: 0 - success, otherwise - failure 980 981 @ref Advanced 982 **/ 983 int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P1d) { 984 if (!basis->tensorbasis) 985 // LCOV_EXCL_START 986 return CeedError(basis->ceed, 1, "Cannot supply P1d for non-tensor basis"); 987 // LCOV_EXCL_STOP 988 989 *P1d = basis->P1d; 990 return 0; 991 } 992 993 /** 994 @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis 995 996 @param basis CeedBasis 997 @param[out] Q1d Variable to store number of quadrature points 998 999 @return An error code: 0 - success, otherwise - failure 1000 1001 @ref Advanced 1002 **/ 1003 int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q1d) { 1004 if (!basis->tensorbasis) 1005 // LCOV_EXCL_START 1006 return CeedError(basis->ceed, 1, "Cannot supply Q1d for non-tensor basis"); 1007 // LCOV_EXCL_STOP 1008 1009 *Q1d = basis->Q1d; 1010 return 0; 1011 } 1012 1013 /** 1014 @brief Get total number of nodes (in dim dimensions) of a CeedBasis 1015 1016 @param basis CeedBasis 1017 @param[out] P Variable to store number of nodes 1018 1019 @return An error code: 0 - success, otherwise - failure 1020 1021 @ref Utility 1022 **/ 1023 int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) { 1024 *P = basis->P; 1025 return 0; 1026 } 1027 1028 /** 1029 @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis 1030 1031 @param basis CeedBasis 1032 @param[out] Q Variable to store number of quadrature points 1033 1034 @return An error code: 0 - success, otherwise - failure 1035 1036 @ref Utility 1037 **/ 1038 int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) { 1039 *Q = basis->Q; 1040 return 0; 1041 } 1042 1043 /** 1044 @brief Get reference coordinates of quadrature points (in dim dimensions) 1045 of a CeedBasis 1046 1047 @param basis CeedBasis 1048 @param[out] qref Variable to store reference coordinates of quadrature points 1049 1050 @return An error code: 0 - success, otherwise - failure 1051 1052 @ref Advanced 1053 **/ 1054 int CeedBasisGetQRef(CeedBasis basis, CeedScalar* *qref) { 1055 *qref = basis->qref1d; 1056 return 0; 1057 } 1058 1059 /** 1060 @brief Get quadrature weights of quadrature points (in dim dimensions) 1061 of a CeedBasis 1062 1063 @param basis CeedBasis 1064 @param[out] qweight Variable to store quadrature weights 1065 1066 @return An error code: 0 - success, otherwise - failure 1067 1068 @ref Advanced 1069 **/ 1070 int CeedBasisGetQWeights(CeedBasis basis, CeedScalar* *qweight) { 1071 *qweight = basis->qweight1d; 1072 return 0; 1073 } 1074 1075 /** 1076 @brief Get interpolation matrix of a CeedBasis 1077 1078 @param basis CeedBasis 1079 @param[out] interp Variable to store interpolation matrix 1080 1081 @return An error code: 0 - success, otherwise - failure 1082 1083 @ref Advanced 1084 **/ 1085 int CeedBasisGetInterp(CeedBasis basis, CeedScalar* *interp) { 1086 *interp = basis->interp1d; 1087 return 0; 1088 } 1089 1090 /** 1091 @brief Get gradient matrix of a CeedBasis 1092 1093 @param basis CeedBasis 1094 @param[out] grad Variable to store gradient matrix 1095 1096 @return An error code: 0 - success, otherwise - failure 1097 1098 @ref Advanced 1099 **/ 1100 int CeedBasisGetGrad(CeedBasis basis, CeedScalar* *grad) { 1101 *grad = basis->grad1d; 1102 return 0; 1103 } 1104 1105 /** 1106 @brief Get value in CeedEvalMode matrix of a CeedBasis 1107 1108 @param basis CeedBasis 1109 @param[in] emode CeedEvalMode to retrieve value 1110 @param[in] node Node (column) to retrieve value 1111 @param[in] qpt Quadrature point (row) to retrieve value 1112 @param[in] dim Dimension to retrieve value for, for CEED_EVAL_GRAD 1113 @param[out] value Variable to store value 1114 1115 @return An error code: 0 - success, otherwise - failure 1116 1117 @ref Advanced 1118 **/ 1119 int CeedBasisGetValue(CeedBasis basis, CeedEvalMode emode, CeedInt qpt, 1120 CeedInt node, CeedInt dim, CeedScalar *value) { 1121 bool tensor = basis->tensorbasis; 1122 1123 switch (emode) { 1124 case CEED_EVAL_NONE: 1125 if (node == qpt) 1126 *value = 0.0; 1127 else 1128 *value = 1.0; 1129 break; 1130 case CEED_EVAL_INTERP: { 1131 CeedScalar *interp = basis->interp1d; 1132 1133 if (tensor) { 1134 CeedInt n, q; 1135 1136 *value = 1.0; 1137 for (CeedInt d=0; d<basis->dim; d++) { 1138 n = (node / CeedIntPow(basis->P1d, d)) % basis->P1d; 1139 q = (qpt / CeedIntPow(basis->Q1d, d)) % basis->Q1d; 1140 *value *= interp[q*(basis->P1d)+n]; 1141 } 1142 } else { 1143 *value = interp[qpt*(basis->P)+node]; 1144 } 1145 } break; 1146 case CEED_EVAL_GRAD: { 1147 CeedScalar *grad = basis->grad1d; 1148 1149 if (tensor) { 1150 CeedInt n, q; 1151 CeedScalar *interp = basis->interp1d; 1152 1153 *value = 1.0; 1154 for (CeedInt d=0; d<basis->dim; d++) { 1155 n = (node / CeedIntPow(basis->P1d, d)) % basis->P1d; 1156 q = (qpt / CeedIntPow(basis->Q1d, d)) % basis->Q1d; 1157 if (d == dim) 1158 *value *= grad[q*(basis->P1d)+n]; 1159 else 1160 *value *= interp[q*(basis->P1d)+n]; 1161 } 1162 } else { 1163 *value = grad[(dim*(basis->Q)+qpt)*(basis->P)+node]; 1164 } 1165 } break; 1166 case CEED_EVAL_WEIGHT: 1167 // LCOV_EXCL_START 1168 return CeedError(basis->ceed, 1, "CEED_EVAL_WEIGHT does not make sense in " 1169 "this context"); 1170 // LCOV_EXCL_STOP 1171 case CEED_EVAL_DIV: 1172 // LCOV_EXCL_START 1173 return CeedError(basis->ceed, 1, "CEED_EVAL_DIV not supported"); 1174 // LCOV_EXCL_STOP 1175 case CEED_EVAL_CURL: 1176 // LCOV_EXCL_START 1177 return CeedError(basis->ceed, 1, "CEED_EVAL_CURL not supported"); 1178 // LCOV_EXCL_STOP 1179 } 1180 return 0; 1181 } 1182 1183 /** 1184 @brief Get backend data of a CeedBasis 1185 1186 @param basis CeedBasis 1187 @param[out] data Variable to store data 1188 1189 @return An error code: 0 - success, otherwise - failure 1190 1191 @ref Advanced 1192 **/ 1193 int CeedBasisGetData(CeedBasis basis, void* *data) { 1194 *data = basis->data; 1195 return 0; 1196 } 1197 1198 /** 1199 @brief Set backend data of a CeedBasis 1200 1201 @param[out] basis CeedBasis 1202 @param data Data to set 1203 1204 @return An error code: 0 - success, otherwise - failure 1205 1206 @ref Advanced 1207 **/ 1208 int CeedBasisSetData(CeedBasis basis, void* *data) { 1209 basis->data = *data; 1210 return 0; 1211 } 1212 1213 /** 1214 @brief Get CeedTensorContract of a CeedBasis 1215 1216 @param basis CeedBasis 1217 @param[out] contract Variable to store CeedTensorContract 1218 1219 @return An error code: 0 - success, otherwise - failure 1220 1221 @ref Advanced 1222 **/ 1223 int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) { 1224 *contract = basis->contract; 1225 return 0; 1226 } 1227 1228 /** 1229 @brief Set CeedTensorContract of a CeedBasis 1230 1231 @param[out] basis CeedBasis 1232 @param contract CeedTensorContract to set 1233 1234 @return An error code: 0 - success, otherwise - failure 1235 1236 @ref Advanced 1237 **/ 1238 int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract *contract) { 1239 basis->contract = *contract; 1240 return 0; 1241 } 1242 1243 /** 1244 @brief Get dimension for given CeedElemTopology 1245 1246 @param topo CeedElemTopology 1247 @param[out] dim Variable to store dimension of topology 1248 1249 @return An error code: 0 - success, otherwise - failure 1250 1251 @ref Advanced 1252 **/ 1253 int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) { 1254 *dim = (CeedInt) topo >> 16; 1255 return 0; 1256 }; 1257 1258 /** 1259 @brief Destroy a CeedBasis 1260 1261 @param basis CeedBasis to destroy 1262 1263 @return An error code: 0 - success, otherwise - failure 1264 1265 @ref Basic 1266 **/ 1267 int CeedBasisDestroy(CeedBasis *basis) { 1268 int ierr; 1269 1270 if (!*basis || --(*basis)->refcount > 0) 1271 return 0; 1272 if ((*basis)->Destroy) { 1273 ierr = (*basis)->Destroy(*basis); CeedChk(ierr); 1274 } 1275 ierr = CeedFree(&(*basis)->interp1d); CeedChk(ierr); 1276 ierr = CeedFree(&(*basis)->grad1d); CeedChk(ierr); 1277 ierr = CeedFree(&(*basis)->qref1d); CeedChk(ierr); 1278 ierr = CeedFree(&(*basis)->qweight1d); CeedChk(ierr); 1279 ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr); 1280 ierr = CeedFree(basis); CeedChk(ierr); 1281 return 0; 1282 } 1283 1284 /// @cond DOXYGEN_SKIP 1285 // Indicate that the quadrature points are collocated with the nodes 1286 CeedBasis CEED_BASIS_COLLOCATED = &ceed_basis_collocated; 1287 /// @endcond 1288 /// @} 1289