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