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