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) fprintf(stream, "%12s[%d]:", name, i); 394 else fprintf(stream, "%12s:", name); 395 for (int j=0; j<n; j++) { 396 fprintf(stream, fpformat, fabs(a[i*n+j]) > 1E-14 ? a[i*n+j] : 0); 397 } 398 fputs("\n", stream); 399 } 400 return 0; 401 } 402 403 /** 404 @brief View a CeedBasis 405 406 @param basis CeedBasis to view 407 @param stream Stream to view to, e.g., stdout 408 409 @return An error code: 0 - success, otherwise - failure 410 411 @ref Utility 412 **/ 413 int CeedBasisView(CeedBasis basis, FILE *stream) { 414 int ierr; 415 416 if (basis->tensorbasis) { 417 fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P1d, 418 basis->Q1d); 419 ierr = CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q1d, basis->qref1d, 420 stream); CeedChk(ierr); 421 ierr = CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q1d, 422 basis->qweight1d, stream); CeedChk(ierr); 423 ierr = CeedScalarView("interp1d", "\t% 12.8f", basis->Q1d, basis->P1d, 424 basis->interp1d, stream); CeedChk(ierr); 425 ierr = CeedScalarView("grad1d", "\t% 12.8f", basis->Q1d, basis->P1d, 426 basis->grad1d, stream); CeedChk(ierr); 427 } else { 428 fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P, 429 basis->Q); 430 ierr = CeedScalarView("qref", "\t% 12.8f", 1, basis->Q*basis->dim, 431 basis->qref1d, 432 stream); CeedChk(ierr); 433 ierr = CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->qweight1d, 434 stream); CeedChk(ierr); 435 ierr = CeedScalarView("interp", "\t% 12.8f", basis->Q, basis->P, 436 basis->interp1d, stream); CeedChk(ierr); 437 ierr = CeedScalarView("grad", "\t% 12.8f", basis->dim*basis->Q, basis->P, 438 basis->grad1d, stream); CeedChk(ierr); 439 } 440 return 0; 441 } 442 443 /** 444 @brief Compute Householder reflection 445 446 Computes A = (I - b v v^T) A 447 where A is an mxn matrix indexed as A[i*row + j*col] 448 449 @param[in,out] A Matrix to apply Householder reflection to, in place 450 @param v Householder vector 451 @param b Scaling factor 452 @param m Number of rows in A 453 @param n Number of columns in A 454 @param row Row stride 455 @param col Col stride 456 457 @return An error code: 0 - success, otherwise - failure 458 459 @ref Developer 460 **/ 461 static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v, 462 CeedScalar b, CeedInt m, CeedInt n, 463 CeedInt row, CeedInt col) { 464 for (CeedInt j=0; j<n; j++) { 465 CeedScalar w = A[0*row + j*col]; 466 for (CeedInt i=1; i<m; i++) w += v[i] * A[i*row + j*col]; 467 A[0*row + j*col] -= b * w; 468 for (CeedInt i=1; i<m; i++) A[i*row + j*col] -= b * w * v[i]; 469 } 470 return 0; 471 } 472 473 /** 474 @brief Apply Householder Q matrix 475 476 Compute A = Q A where Q is mxm and A is mxn. 477 478 @param[in,out] A Matrix to apply Householder Q to, in place 479 @param Q Householder Q matrix 480 @param tau Householder scaling factors 481 @param tmode Transpose mode for application 482 @param m Number of rows in A 483 @param n Number of columns in A 484 @param k Number of elementary reflectors in Q, k<m 485 @param row Row stride in A 486 @param col Col stride in A 487 488 @return An error code: 0 - success, otherwise - failure 489 490 @ref Developer 491 **/ 492 static int CeedHouseholderApplyQ(CeedScalar *A, const CeedScalar *Q, 493 const CeedScalar *tau, CeedTransposeMode tmode, 494 CeedInt m, CeedInt n, CeedInt k, 495 CeedInt row, CeedInt col) { 496 CeedScalar v[m]; 497 for (CeedInt ii=0; ii<k; ii++) { 498 CeedInt i = tmode == CEED_TRANSPOSE ? ii : k-1-ii; 499 for (CeedInt j=i+1; j<m; j++) 500 v[j] = Q[j*k+i]; 501 // Apply Householder reflector (I - tau v v^T) colograd1d^T 502 CeedHouseholderReflect(&A[i*row], &v[i], tau[i], m-i, n, row, col); 503 } 504 return 0; 505 } 506 507 /** 508 @brief Compute Givens rotation 509 510 Computes A = G A (or G^T A in transpose mode) 511 where A is an mxn matrix indexed as A[i*n + j*m] 512 513 @param[in,out] A Row major matrix to apply Givens rotation to, in place 514 @param c Cosine factor 515 @param s Sine factor 516 @param i First row/column to apply rotation 517 @param k Second row/column to apply rotation 518 @param m Number of rows in A 519 @param n Number of columns in A 520 521 @return An error code: 0 - success, otherwise - failure 522 523 @ref Developer 524 **/ 525 static int CeedGivensRotation(CeedScalar *A, CeedScalar c, CeedScalar s, 526 CeedTransposeMode tmode, CeedInt i, CeedInt k, 527 CeedInt m, CeedInt n) { 528 CeedInt stridej = 1, strideik = m, numits = n; 529 if (tmode == CEED_NOTRANSPOSE) { 530 stridej = n; strideik = 1; numits = m; 531 } 532 533 // Apply rotation 534 for (CeedInt j=0; j<numits; j++) { 535 CeedScalar tau1 = A[i*strideik+j*stridej], tau2 = A[k*strideik+j*stridej]; 536 A[i*strideik+j*stridej] = c*tau1 - s*tau2; 537 A[k*strideik+j*stridej] = s*tau1 + c*tau2; 538 } 539 540 return 0; 541 } 542 543 /** 544 @brief Return QR Factorization of matrix 545 546 @param ceed A Ceed object currently in use 547 @param[in,out] mat Row-major matrix to be factorized in place 548 @param[in,out] tau Vector of length m of scaling factors 549 @param m Number of rows 550 @param n Number of columns 551 552 @return An error code: 0 - success, otherwise - failure 553 554 @ref Utility 555 **/ 556 int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau, 557 CeedInt m, CeedInt n) { 558 CeedScalar v[m]; 559 560 // Check m >= n 561 if (n > m) 562 // LCOV_EXCL_START 563 return CeedError(ceed, 1, "Cannot compute QR factorization with n > m"); 564 // LCOV_EXCL_STOP 565 566 for (CeedInt i=0; i<n; i++) { 567 // Calculate Householder vector, magnitude 568 CeedScalar sigma = 0.0; 569 v[i] = mat[i+n*i]; 570 for (CeedInt j=i+1; j<m; j++) { 571 v[j] = mat[i+n*j]; 572 sigma += v[j] * v[j]; 573 } 574 CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m] 575 CeedScalar Rii = -copysign(norm, v[i]); 576 v[i] -= Rii; 577 // norm of v[i:m] after modification above and scaling below 578 // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 579 // tau = 2 / (norm*norm) 580 tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma); 581 for (CeedInt j=i+1; j<m; j++) v[j] /= v[i]; 582 583 // Apply Householder reflector to lower right panel 584 CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1); 585 // Save v 586 mat[i+n*i] = Rii; 587 for (CeedInt j=i+1; j<m; j++) { 588 mat[i+n*j] = v[j]; 589 } 590 } 591 592 return 0; 593 } 594 595 /** 596 @brief Return symmetric Schur decomposition of the symmetric matrix mat via 597 symmetric QR factorization 598 599 @param[in,out] mat Row-major matrix to be factorized in place 600 @param[out] lambda Vector of length m of eigenvalues 601 @param n Number of rows/columns 602 603 @return An error code: 0 - success, otherwise - failure 604 605 @ref Utility 606 **/ 607 int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat, 608 CeedScalar *lambda, CeedInt n) { 609 // Check bounds for clang-tidy 610 if (n<2) 611 // LCOV_EXCL_START 612 return CeedError(ceed, 1, 613 "Cannot compute symmetric Schur decomposition of scalars"); 614 // LCOV_EXCL_STOP 615 616 CeedScalar v[n-1], tau[n-1], matT[n*n]; 617 618 // Copy mat to matT and set mat to I 619 memcpy(matT, mat, n*n*sizeof(mat[0])); 620 for (CeedInt i=0; i<n; i++) 621 for (CeedInt j=0; j<n; j++) 622 mat[j+n*i] = (i==j) ? 1 : 0; 623 624 // Reduce to tridiagonal 625 for (CeedInt i=0; i<n-1; i++) { 626 // Calculate Householder vector, magnitude 627 CeedScalar sigma = 0.0; 628 v[i] = matT[i+n*(i+1)]; 629 for (CeedInt j=i+1; j<n-1; j++) { 630 v[j] = matT[i+n*(j+1)]; 631 sigma += v[j] * v[j]; 632 } 633 CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:n-1] 634 CeedScalar Rii = -copysign(norm, v[i]); 635 v[i] -= Rii; 636 // norm of v[i:m] after modification above and scaling below 637 // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 638 // tau = 2 / (norm*norm) 639 tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma); 640 for (CeedInt j=i+1; j<n-1; j++) v[j] /= v[i]; 641 642 // Update sub and super diagonal 643 matT[i+n*(i+1)] = Rii; 644 matT[(i+1)+n*i] = Rii; 645 for (CeedInt j=i+2; j<n; j++) { 646 matT[i+n*j] = 0; matT[j+n*i] = 0; 647 } 648 // Apply symmetric Householder reflector to lower right panel 649 CeedHouseholderReflect(&matT[(i+1)+n*(i+1)], &v[i], tau[i], 650 n-(i+1), n-(i+1), n, 1); 651 CeedHouseholderReflect(&matT[(i+1)+n*(i+1)], &v[i], tau[i], 652 n-(i+1), n-(i+1), 1, n); 653 // Save v 654 for (CeedInt j=i+1; j<n-1; j++) { 655 matT[i+n*(j+1)] = v[j]; 656 } 657 } 658 // Backwards accumulation of Q 659 for (CeedInt i=n-2; i>=0; i--) { 660 v[i] = 1; 661 for (CeedInt j=i+1; j<n-1; j++) { 662 v[j] = matT[i+n*(j+1)]; 663 matT[i+n*(j+1)] = 0; 664 } 665 CeedHouseholderReflect(&mat[(i+1)+n*(i+1)], &v[i], tau[i], 666 n-(i+1), n-(i+1), n, 1); 667 } 668 669 // Reduce sub and super diagonal 670 CeedInt p = 0, q = 0, itr = 0, maxitr = n*n*n; 671 CeedScalar tol = 1e-15; 672 673 while (q < n && itr < maxitr) { 674 // Update p, q, size of reduced portions of diagonal 675 p = 0; q = 0; 676 for (CeedInt i=n-2; i>=0; i--) { 677 if (fabs(matT[i+n*(i+1)]) < tol) 678 q += 1; 679 else 680 break; 681 } 682 for (CeedInt i=0; i<n-1-q; i++) { 683 if (fabs(matT[i+n*(i+1)]) < tol) 684 p += 1; 685 else 686 break; 687 } 688 if (q == n-1) break; // Finished reducing 689 690 // Reduce tridiagonal portion 691 CeedScalar tnn = matT[(n-1-q)+n*(n-1-q)], 692 tnnm1 = matT[(n-2-q)+n*(n-1-q)]; 693 CeedScalar d = (matT[(n-2-q)+n*(n-2-q)] - tnn)/2; 694 CeedScalar mu = tnn - tnnm1*tnnm1 / 695 (d + copysign(sqrt(d*d + tnnm1*tnnm1), d)); 696 CeedScalar x = matT[p+n*p] - mu; 697 CeedScalar z = matT[p+n*(p+1)]; 698 for (CeedInt k=p; k<n-1-q; k++) { 699 // Compute Givens rotation 700 CeedScalar c = 1, s = 0; 701 if (fabs(z) > tol) { 702 if (fabs(z) > fabs(x)) { 703 CeedScalar tau = -x/z; 704 s = 1/sqrt(1+tau*tau), c = s*tau; 705 } else { 706 CeedScalar tau = -z/x; 707 c = 1/sqrt(1+tau*tau), s = c*tau; 708 } 709 } 710 711 // Apply Givens rotation to T 712 CeedGivensRotation(matT, c, s, CEED_NOTRANSPOSE, k, k+1, n, n); 713 CeedGivensRotation(matT, c, s, CEED_TRANSPOSE, k, k+1, n, n); 714 715 // Apply Givens rotation to Q 716 CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k+1, n, n); 717 718 // Update x, z 719 if (k < n-q-2) { 720 x = matT[k+n*(k+1)]; 721 z = matT[k+n*(k+2)]; 722 } 723 } 724 itr++; 725 } 726 // Save eigenvalues 727 for (CeedInt i=0; i<n; i++) 728 lambda[i] = matT[i+n*i]; 729 730 // Check convergence 731 if (itr == maxitr && q < n-1) 732 // LCOV_EXCL_START 733 return CeedError(ceed, 1, "Symmetric QR failed to converge"); 734 // LCOV_EXCL_STOP 735 736 return 0; 737 } 738 739 /** 740 @brief Return C = A B 741 742 @param[in] matA Row-major matrix A 743 @param[in] matB Row-major matrix B 744 @param[out] matC Row-major output matrix C 745 @param m Number of rows of C 746 @param n Number of columns of C 747 @param kk Number of columns of A/rows of B 748 749 @return An error code: 0 - success, otherwise - failure 750 751 @ref Utility 752 **/ 753 static int CeedMatrixMultiply(Ceed ceed, CeedScalar *matA, CeedScalar *matB, 754 CeedScalar *matC, CeedInt m, CeedInt n, 755 CeedInt kk) { 756 for (CeedInt i=0; i<m; i++) 757 for (CeedInt j=0; j<n; j++) { 758 CeedScalar sum = 0; 759 for (CeedInt k=0; k<kk; k++) 760 sum += matA[k+i*kk]*matB[j+k*n]; 761 matC[j+i*n] = sum; 762 } 763 return 0; 764 } 765 766 /** 767 @brief Return Simultaneous Diagonalization of two matrices. This solves the 768 generalized eigenvalue problem A x = lambda B x, where A and B 769 are symmetric and B is positive definite. We generate the matrix X 770 and vector Lambda such that X^T A X = Lambda and X^T B X = I. This 771 is equivalent to the LAPACK routine 'sygv' with TYPE = 1. 772 773 @param[in] matA Row-major matrix to be factorized with eigenvalues 774 @param[in] matB Row-major matrix to be factorized to identity 775 @param[out] x Row-major orthogonal matrix 776 @param[out] lambda Vector of length m of generalized eigenvalues 777 @param n Number of rows/columns 778 779 @return An error code: 0 - success, otherwise - failure 780 781 @ref Utility 782 **/ 783 int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *matA, 784 CeedScalar *matB, CeedScalar *x, 785 CeedScalar *lambda, CeedInt n) { 786 int ierr; 787 CeedScalar matC[n*n], matG[n*n], vecD[n]; 788 789 // Compute B = G D G^T 790 memcpy(matG, matB, n*n*sizeof(matB[0])); 791 ierr = CeedSymmetricSchurDecomposition(ceed, matG, vecD, n); CeedChk(ierr); 792 for (CeedInt i=0; i<n; i++) vecD[i] = sqrt(vecD[i]); 793 794 // Compute C = (G D^-1/2)^-1 A (G D^-1/2)^-T 795 // = D^1/2 G^T A D^1/2 G 796 for (CeedInt i=0; i<n; i++) 797 for (CeedInt j=0; j<n; j++) 798 matC[j+i*n] = vecD[i] * matG[i+j*n]; 799 CeedMatrixMultiply(ceed, matC, matA, x, n, n, n); 800 for (CeedInt i=0; i<n; i++) 801 for (CeedInt j=0; j<n; j++) 802 matG[j+i*n] = vecD[i] * matG[j+i*n]; 803 CeedMatrixMultiply(ceed, x, matG, matC, n, n, n); 804 805 // Compute Q^T C Q = lambda 806 ierr = CeedSymmetricSchurDecomposition(ceed, matC, lambda, n); CeedChk(ierr); 807 808 // Set x = (G D^-1/2)^-T Q 809 // = D^1/2 G Q 810 CeedMatrixMultiply(ceed, matG, matC, x, n, n, n); 811 812 return 0; 813 } 814 815 /** 816 @brief Return collocated grad matrix 817 818 @param basis CeedBasis 819 @param[out] colograd1d Row-major Q1d × Q1d matrix expressing derivatives of 820 basis functions at quadrature points 821 822 @return An error code: 0 - success, otherwise - failure 823 824 @ref Advanced 825 **/ 826 int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *colograd1d) { 827 int i, j, k; 828 Ceed ceed; 829 CeedInt ierr, P1d=(basis)->P1d, Q1d=(basis)->Q1d; 830 CeedScalar *interp1d, *grad1d, tau[Q1d]; 831 832 ierr = CeedMalloc(Q1d*P1d, &interp1d); CeedChk(ierr); 833 ierr = CeedMalloc(Q1d*P1d, &grad1d); CeedChk(ierr); 834 memcpy(interp1d, (basis)->interp1d, Q1d*P1d*sizeof(basis)->interp1d[0]); 835 memcpy(grad1d, (basis)->grad1d, Q1d*P1d*sizeof(basis)->interp1d[0]); 836 837 // QR Factorization, interp1d = Q R 838 ierr = CeedBasisGetCeed(basis, &ceed); CeedChk(ierr); 839 ierr = CeedQRFactorization(ceed, interp1d, tau, Q1d, P1d); CeedChk(ierr); 840 841 // Apply Rinv, colograd1d = grad1d Rinv 842 for (i=0; i<Q1d; i++) { // Row i 843 colograd1d[Q1d*i] = grad1d[P1d*i]/interp1d[0]; 844 for (j=1; j<P1d; j++) { // Column j 845 colograd1d[j+Q1d*i] = grad1d[j+P1d*i]; 846 for (k=0; k<j; k++) { 847 colograd1d[j+Q1d*i] -= interp1d[j+P1d*k]*colograd1d[k+Q1d*i]; 848 } 849 colograd1d[j+Q1d*i] /= interp1d[j+P1d*j]; 850 } 851 for (j=P1d; j<Q1d; j++) { 852 colograd1d[j+Q1d*i] = 0; 853 } 854 } 855 856 // Apply Qtranspose, colograd = colograd Qtranspose 857 CeedHouseholderApplyQ(colograd1d, interp1d, tau, CEED_NOTRANSPOSE, 858 Q1d, Q1d, P1d, 1, Q1d); 859 860 ierr = CeedFree(&interp1d); CeedChk(ierr); 861 ierr = CeedFree(&grad1d); CeedChk(ierr); 862 863 return 0; 864 } 865 866 /** 867 @brief Apply basis evaluation from nodes to quadrature points or vice-versa 868 869 @param basis CeedBasis to evaluate 870 @param nelem The number of elements to apply the basis evaluation to; 871 the backend will specify the ordering in 872 ElemRestrictionCreateBlocked 873 @param tmode \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature 874 points, \ref CEED_TRANSPOSE to apply the transpose, mapping 875 from quadrature points to nodes 876 @param emode \ref CEED_EVAL_INTERP to obtain interpolated values, 877 \ref CEED_EVAL_GRAD to obtain gradients. 878 @param[in] u Input array 879 @param[out] v Output array 880 881 @return An error code: 0 - success, otherwise - failure 882 883 @ref Advanced 884 **/ 885 int CeedBasisApply(CeedBasis basis, CeedInt nelem, CeedTransposeMode tmode, 886 CeedEvalMode emode, CeedVector u, CeedVector v) { 887 int ierr; 888 CeedInt ulength = 0, vlength, nnodes, nqpt; 889 if (!basis->Apply) 890 // LCOV_EXCL_START 891 return CeedError(basis->ceed, 1, "Backend does not support BasisApply"); 892 // LCOV_EXCL_STOP 893 894 // Check compatibility of topological and geometrical dimensions 895 ierr = CeedBasisGetNumNodes(basis, &nnodes); CeedChk(ierr); 896 ierr = CeedBasisGetNumQuadraturePoints(basis, &nqpt); CeedChk(ierr); 897 ierr = CeedVectorGetLength(v, &vlength); CeedChk(ierr); 898 899 if (u) { 900 ierr = CeedVectorGetLength(u, &ulength); CeedChk(ierr); 901 } 902 903 if ((tmode == CEED_TRANSPOSE && (vlength % nnodes != 0 904 || ulength % nqpt != 0)) 905 || 906 (tmode == CEED_NOTRANSPOSE && (ulength % nnodes != 0 || vlength % nqpt != 0))) 907 return CeedError(basis->ceed, 1, 908 "Length of input/output vectors incompatible with basis dimensions"); 909 910 ierr = basis->Apply(basis, nelem, tmode, emode, u, v); CeedChk(ierr); 911 return 0; 912 } 913 914 /** 915 @brief Get Ceed associated with a CeedBasis 916 917 @param basis CeedBasis 918 @param[out] ceed Variable to store Ceed 919 920 @return An error code: 0 - success, otherwise - failure 921 922 @ref Advanced 923 **/ 924 int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) { 925 *ceed = basis->ceed; 926 927 return 0; 928 }; 929 930 /** 931 @brief Get dimension for given CeedBasis 932 933 @param basis CeedBasis 934 @param[out] dim Variable to store dimension of basis 935 936 @return An error code: 0 - success, otherwise - failure 937 938 @ref Advanced 939 **/ 940 int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) { 941 *dim = basis->dim; 942 943 return 0; 944 }; 945 946 /** 947 @brief Get tensor status for given CeedBasis 948 949 @param basis CeedBasis 950 @param[out] tensor Variable to store tensor status 951 952 @return An error code: 0 - success, otherwise - failure 953 954 @ref Advanced 955 **/ 956 int CeedBasisGetTensorStatus(CeedBasis basis, bool *tensor) { 957 *tensor = basis->tensorbasis; 958 959 return 0; 960 }; 961 962 /** 963 @brief Get number of components for given CeedBasis 964 965 @param basis CeedBasis 966 @param[out] numcomp Variable to store number of components of basis 967 968 @return An error code: 0 - success, otherwise - failure 969 970 @ref Advanced 971 **/ 972 int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *numcomp) { 973 *numcomp = basis->ncomp; 974 975 return 0; 976 }; 977 978 /** 979 @brief Get total number of nodes (in 1 dimension) of a CeedBasis 980 981 @param basis CeedBasis 982 @param[out] P1d Variable to store number of nodes 983 984 @return An error code: 0 - success, otherwise - failure 985 986 @ref Advanced 987 **/ 988 int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P1d) { 989 if (!basis->tensorbasis) 990 // LCOV_EXCL_START 991 return CeedError(basis->ceed, 1, "Cannot supply P1d for non-tensor basis"); 992 // LCOV_EXCL_STOP 993 994 *P1d = basis->P1d; 995 return 0; 996 } 997 998 /** 999 @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis 1000 1001 @param basis CeedBasis 1002 @param[out] Q1d Variable to store number of quadrature points 1003 1004 @return An error code: 0 - success, otherwise - failure 1005 1006 @ref Advanced 1007 **/ 1008 int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q1d) { 1009 if (!basis->tensorbasis) 1010 // LCOV_EXCL_START 1011 return CeedError(basis->ceed, 1, "Cannot supply Q1d for non-tensor basis"); 1012 // LCOV_EXCL_STOP 1013 1014 *Q1d = basis->Q1d; 1015 return 0; 1016 } 1017 1018 /** 1019 @brief Get total number of nodes (in dim dimensions) of a CeedBasis 1020 1021 @param basis CeedBasis 1022 @param[out] P Variable to store number of nodes 1023 1024 @return An error code: 0 - success, otherwise - failure 1025 1026 @ref Utility 1027 **/ 1028 int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) { 1029 *P = basis->P; 1030 return 0; 1031 } 1032 1033 /** 1034 @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis 1035 1036 @param basis CeedBasis 1037 @param[out] Q Variable to store number of quadrature points 1038 1039 @return An error code: 0 - success, otherwise - failure 1040 1041 @ref Utility 1042 **/ 1043 int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) { 1044 *Q = basis->Q; 1045 return 0; 1046 } 1047 1048 /** 1049 @brief Get reference coordinates of quadrature points (in dim dimensions) 1050 of a CeedBasis 1051 1052 @param basis CeedBasis 1053 @param[out] qref Variable to store reference coordinates of quadrature points 1054 1055 @return An error code: 0 - success, otherwise - failure 1056 1057 @ref Advanced 1058 **/ 1059 int CeedBasisGetQRef(CeedBasis basis, CeedScalar* *qref) { 1060 *qref = basis->qref1d; 1061 return 0; 1062 } 1063 1064 /** 1065 @brief Get quadrature weights of quadrature points (in dim dimensions) 1066 of a CeedBasis 1067 1068 @param basis CeedBasis 1069 @param[out] qweight Variable to store quadrature weights 1070 1071 @return An error code: 0 - success, otherwise - failure 1072 1073 @ref Advanced 1074 **/ 1075 int CeedBasisGetQWeights(CeedBasis basis, CeedScalar* *qweight) { 1076 *qweight = basis->qweight1d; 1077 return 0; 1078 } 1079 1080 /** 1081 @brief Get interpolation matrix of a CeedBasis 1082 1083 @param basis CeedBasis 1084 @param[out] interp Variable to store interpolation matrix 1085 1086 @return An error code: 0 - success, otherwise - failure 1087 1088 @ref Advanced 1089 **/ 1090 int CeedBasisGetInterp(CeedBasis basis, CeedScalar* *interp) { 1091 *interp = basis->interp1d; 1092 return 0; 1093 } 1094 1095 /** 1096 @brief Get gradient matrix of a CeedBasis 1097 1098 @param basis CeedBasis 1099 @param[out] grad Variable to store gradient matrix 1100 1101 @return An error code: 0 - success, otherwise - failure 1102 1103 @ref Advanced 1104 **/ 1105 int CeedBasisGetGrad(CeedBasis basis, CeedScalar* *grad) { 1106 *grad = basis->grad1d; 1107 return 0; 1108 } 1109 1110 /** 1111 @brief Get backend data of a CeedBasis 1112 1113 @param basis CeedBasis 1114 @param[out] data Variable to store data 1115 1116 @return An error code: 0 - success, otherwise - failure 1117 1118 @ref Advanced 1119 **/ 1120 int CeedBasisGetData(CeedBasis basis, void* *data) { 1121 *data = basis->data; 1122 return 0; 1123 } 1124 1125 /** 1126 @brief Set backend data of a CeedBasis 1127 1128 @param[out] basis CeedBasis 1129 @param data Data to set 1130 1131 @return An error code: 0 - success, otherwise - failure 1132 1133 @ref Advanced 1134 **/ 1135 int CeedBasisSetData(CeedBasis basis, void* *data) { 1136 basis->data = *data; 1137 return 0; 1138 } 1139 1140 /** 1141 @brief Get CeedTensorContract of a CeedBasis 1142 1143 @param basis CeedBasis 1144 @param[out] contract Variable to store CeedTensorContract 1145 1146 @return An error code: 0 - success, otherwise - failure 1147 1148 @ref Advanced 1149 **/ 1150 int CeedBasisGetTensorContract(CeedBasis basis, 1151 CeedTensorContract *contract) { 1152 *contract = basis->contract; 1153 return 0; 1154 } 1155 1156 /** 1157 @brief Set CeedTensorContract of a CeedBasis 1158 1159 @param[out] basis CeedBasis 1160 @param contract CeedTensorContract to set 1161 1162 @return An error code: 0 - success, otherwise - failure 1163 1164 @ref Advanced 1165 **/ 1166 int CeedBasisSetTensorContract(CeedBasis basis, 1167 CeedTensorContract *contract) { 1168 basis->contract = *contract; 1169 return 0; 1170 } 1171 1172 /** 1173 @brief Get dimension for given CeedElemTopology 1174 1175 @param topo CeedElemTopology 1176 @param[out] dim Variable to store dimension of topology 1177 1178 @return An error code: 0 - success, otherwise - failure 1179 1180 @ref Advanced 1181 **/ 1182 int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) { 1183 *dim = (CeedInt) topo >> 16; 1184 1185 return 0; 1186 }; 1187 1188 /** 1189 @brief Destroy a CeedBasis 1190 1191 @param basis CeedBasis to destroy 1192 1193 @return An error code: 0 - success, otherwise - failure 1194 1195 @ref Basic 1196 **/ 1197 int CeedBasisDestroy(CeedBasis *basis) { 1198 int ierr; 1199 1200 if (!*basis || --(*basis)->refcount > 0) return 0; 1201 if ((*basis)->Destroy) { 1202 ierr = (*basis)->Destroy(*basis); CeedChk(ierr); 1203 } 1204 ierr = CeedFree(&(*basis)->interp1d); CeedChk(ierr); 1205 ierr = CeedFree(&(*basis)->grad1d); CeedChk(ierr); 1206 ierr = CeedFree(&(*basis)->qref1d); CeedChk(ierr); 1207 ierr = CeedFree(&(*basis)->qweight1d); CeedChk(ierr); 1208 ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr); 1209 ierr = CeedFree(basis); CeedChk(ierr); 1210 return 0; 1211 } 1212 1213 /// @cond DOXYGEN_SKIP 1214 // Indicate that the quadrature points are collocated with the nodes 1215 CeedBasis CEED_BASIS_COLLOCATED = &ceed_basis_collocated; 1216 /// @endcond 1217 /// @} 1218