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 CeedBasisIncrementRefCounter(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 = CeedTensorContractIncrementRefCounter(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 = CeedIncrementRefCounter(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 = CeedIncrementRefCounter(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 View a CeedBasis 638 639 @param basis CeedBasis to view 640 @param stream Stream to view to, e.g., stdout 641 642 @return An error code: 0 - success, otherwise - failure 643 644 @ref User 645 **/ 646 int CeedBasisView(CeedBasis basis, FILE *stream) { 647 int ierr; 648 649 if (basis->tensor_basis) { 650 fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P_1d, 651 basis->Q_1d); 652 ierr = CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_ref_1d, 653 stream); CeedChk(ierr); 654 ierr = CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q_1d, 655 basis->q_weight_1d, stream); CeedChk(ierr); 656 ierr = CeedScalarView("interp1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, 657 basis->interp_1d, stream); CeedChk(ierr); 658 ierr = CeedScalarView("grad1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, 659 basis->grad_1d, stream); CeedChk(ierr); 660 } else { 661 fprintf(stream, "CeedBasis: dim=%d P=%d Q=%d\n", basis->dim, basis->P, 662 basis->Q); 663 ierr = CeedScalarView("qref", "\t% 12.8f", 1, basis->Q*basis->dim, 664 basis->q_ref_1d, 665 stream); CeedChk(ierr); 666 ierr = CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->q_weight_1d, 667 stream); CeedChk(ierr); 668 ierr = CeedScalarView("interp", "\t% 12.8f", basis->Q, basis->P, 669 basis->interp, stream); CeedChk(ierr); 670 ierr = CeedScalarView("grad", "\t% 12.8f", basis->dim*basis->Q, basis->P, 671 basis->grad, stream); CeedChk(ierr); 672 } 673 return CEED_ERROR_SUCCESS; 674 } 675 676 /** 677 @brief Apply basis evaluation from nodes to quadrature points or vice versa 678 679 @param basis CeedBasis to evaluate 680 @param num_elem The number of elements to apply the basis evaluation to; 681 the backend will specify the ordering in 682 CeedElemRestrictionCreateBlocked() 683 @param t_mode \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature 684 points, \ref CEED_TRANSPOSE to apply the transpose, mapping 685 from quadrature points to nodes 686 @param eval_mode \ref CEED_EVAL_NONE to use values directly, 687 \ref CEED_EVAL_INTERP to use interpolated values, 688 \ref CEED_EVAL_GRAD to use gradients, 689 \ref CEED_EVAL_WEIGHT to use quadrature weights. 690 @param[in] u Input CeedVector 691 @param[out] v Output CeedVector 692 693 @return An error code: 0 - success, otherwise - failure 694 695 @ref User 696 **/ 697 int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode, 698 CeedEvalMode eval_mode, CeedVector u, CeedVector v) { 699 int ierr; 700 CeedInt u_length = 0, v_length, dim, num_comp, num_nodes, num_qpts; 701 ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr); 702 ierr = CeedBasisGetNumComponents(basis, &num_comp); CeedChk(ierr); 703 ierr = CeedBasisGetNumNodes(basis, &num_nodes); CeedChk(ierr); 704 ierr = CeedBasisGetNumQuadraturePoints(basis, &num_qpts); CeedChk(ierr); 705 ierr = CeedVectorGetLength(v, &v_length); CeedChk(ierr); 706 if (u) { 707 ierr = CeedVectorGetLength(u, &u_length); CeedChk(ierr); 708 } 709 710 if (!basis->Apply) 711 // LCOV_EXCL_START 712 return CeedError(basis->ceed, CEED_ERROR_UNSUPPORTED, 713 "Backend does not support BasisApply"); 714 // LCOV_EXCL_STOP 715 716 // Check compatibility of topological and geometrical dimensions 717 if ((t_mode == CEED_TRANSPOSE && (v_length%num_nodes != 0 || 718 u_length%num_qpts != 0)) || 719 (t_mode == CEED_NOTRANSPOSE && (u_length%num_nodes != 0 || 720 v_length%num_qpts != 0))) 721 return CeedError(basis->ceed, CEED_ERROR_DIMENSION, 722 "Length of input/output vectors " 723 "incompatible with basis dimensions"); 724 725 // Check vector lengths to prevent out of bounds issues 726 bool bad_dims = false; 727 switch (eval_mode) { 728 case CEED_EVAL_NONE: 729 case CEED_EVAL_INTERP: bad_dims = 730 ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts || 731 v_length < num_elem*num_comp*num_nodes)) || 732 (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp || 733 u_length < num_elem*num_comp*num_nodes))); 734 break; 735 case CEED_EVAL_GRAD: bad_dims = 736 ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts*dim || 737 v_length < num_elem*num_comp*num_nodes)) || 738 (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp*dim || 739 u_length < num_elem*num_comp*num_nodes))); 740 break; 741 case CEED_EVAL_WEIGHT: 742 bad_dims = v_length < num_elem*num_qpts; 743 break; 744 // LCOV_EXCL_START 745 case CEED_EVAL_DIV: bad_dims = 746 ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts || 747 v_length < num_elem*num_comp*num_nodes)) || 748 (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp || 749 u_length < num_elem*num_comp*num_nodes))); 750 break; 751 case CEED_EVAL_CURL: bad_dims = 752 ((t_mode == CEED_TRANSPOSE && (u_length < num_elem*num_comp*num_qpts || 753 v_length < num_elem*num_comp*num_nodes)) || 754 (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem*num_qpts*num_comp || 755 u_length < num_elem*num_comp*num_nodes))); 756 break; 757 // LCOV_EXCL_STOP 758 } 759 if (bad_dims) 760 // LCOV_EXCL_START 761 return CeedError(basis->ceed, CEED_ERROR_DIMENSION, 762 "Input/output vectors too short for basis and evaluation mode"); 763 // LCOV_EXCL_STOP 764 765 ierr = basis->Apply(basis, num_elem, t_mode, eval_mode, u, v); CeedChk(ierr); 766 return CEED_ERROR_SUCCESS; 767 } 768 769 /** 770 @brief Get dimension for given CeedBasis 771 772 @param basis CeedBasis 773 @param[out] dim Variable to store dimension of basis 774 775 @return An error code: 0 - success, otherwise - failure 776 777 @ref Backend 778 **/ 779 int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) { 780 *dim = basis->dim; 781 return CEED_ERROR_SUCCESS; 782 } 783 784 /** 785 @brief Get topology for given CeedBasis 786 787 @param basis CeedBasis 788 @param[out] topo Variable to store topology of basis 789 790 @return An error code: 0 - success, otherwise - failure 791 792 @ref Backend 793 **/ 794 int CeedBasisGetTopology(CeedBasis basis, CeedElemTopology *topo) { 795 *topo = basis->topo; 796 return CEED_ERROR_SUCCESS; 797 } 798 799 /** 800 @brief Get number of components for given CeedBasis 801 802 @param basis CeedBasis 803 @param[out] num_comp Variable to store number of components of basis 804 805 @return An error code: 0 - success, otherwise - failure 806 807 @ref Backend 808 **/ 809 int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *num_comp) { 810 *num_comp = basis->num_comp; 811 return CEED_ERROR_SUCCESS; 812 } 813 814 /** 815 @brief Get total number of nodes (in dim dimensions) of a CeedBasis 816 817 @param basis CeedBasis 818 @param[out] P Variable to store number of nodes 819 820 @return An error code: 0 - success, otherwise - failure 821 822 @ref Utility 823 **/ 824 int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) { 825 *P = basis->P; 826 return CEED_ERROR_SUCCESS; 827 } 828 829 /** 830 @brief Get total number of nodes (in 1 dimension) of a CeedBasis 831 832 @param basis CeedBasis 833 @param[out] P_1d Variable to store number of nodes 834 835 @return An error code: 0 - success, otherwise - failure 836 837 @ref Backend 838 **/ 839 int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P_1d) { 840 if (!basis->tensor_basis) 841 // LCOV_EXCL_START 842 return CeedError(basis->ceed, CEED_ERROR_MINOR, 843 "Cannot supply P_1d for non-tensor basis"); 844 // LCOV_EXCL_STOP 845 846 *P_1d = basis->P_1d; 847 return CEED_ERROR_SUCCESS; 848 } 849 850 /** 851 @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis 852 853 @param basis CeedBasis 854 @param[out] Q Variable to store number of quadrature points 855 856 @return An error code: 0 - success, otherwise - failure 857 858 @ref Utility 859 **/ 860 int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) { 861 *Q = basis->Q; 862 return CEED_ERROR_SUCCESS; 863 } 864 865 /** 866 @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis 867 868 @param basis CeedBasis 869 @param[out] Q_1d Variable to store number of quadrature points 870 871 @return An error code: 0 - success, otherwise - failure 872 873 @ref Backend 874 **/ 875 int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q_1d) { 876 if (!basis->tensor_basis) 877 // LCOV_EXCL_START 878 return CeedError(basis->ceed, CEED_ERROR_MINOR, 879 "Cannot supply Q_1d for non-tensor basis"); 880 // LCOV_EXCL_STOP 881 882 *Q_1d = basis->Q_1d; 883 return CEED_ERROR_SUCCESS; 884 } 885 886 /** 887 @brief Get reference coordinates of quadrature points (in dim dimensions) 888 of a CeedBasis 889 890 @param basis CeedBasis 891 @param[out] q_ref Variable to store reference coordinates of quadrature points 892 893 @return An error code: 0 - success, otherwise - failure 894 895 @ref Backend 896 **/ 897 int CeedBasisGetQRef(CeedBasis basis, const CeedScalar **q_ref) { 898 *q_ref = basis->q_ref_1d; 899 return CEED_ERROR_SUCCESS; 900 } 901 902 /** 903 @brief Get quadrature weights of quadrature points (in dim dimensions) 904 of a CeedBasis 905 906 @param basis CeedBasis 907 @param[out] q_weight Variable to store quadrature weights 908 909 @return An error code: 0 - success, otherwise - failure 910 911 @ref Backend 912 **/ 913 int CeedBasisGetQWeights(CeedBasis basis, const CeedScalar **q_weight) { 914 *q_weight = basis->q_weight_1d; 915 return CEED_ERROR_SUCCESS; 916 } 917 918 /** 919 @brief Get interpolation matrix of a CeedBasis 920 921 @param basis CeedBasis 922 @param[out] interp Variable to store interpolation matrix 923 924 @return An error code: 0 - success, otherwise - failure 925 926 @ref Backend 927 **/ 928 int CeedBasisGetInterp(CeedBasis basis, const CeedScalar **interp) { 929 if (!basis->interp && basis->tensor_basis) { 930 // Allocate 931 int ierr; 932 ierr = CeedMalloc(basis->Q*basis->P, &basis->interp); CeedChk(ierr); 933 934 // Initialize 935 for (CeedInt i=0; i<basis->Q*basis->P; i++) 936 basis->interp[i] = 1.0; 937 938 // Calculate 939 for (CeedInt d=0; d<basis->dim; d++) 940 for (CeedInt qpt=0; qpt<basis->Q; qpt++) 941 for (CeedInt node=0; node<basis->P; node++) { 942 CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d; 943 CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d; 944 basis->interp[qpt*(basis->P)+node] *= basis->interp_1d[q*basis->P_1d+p]; 945 } 946 } 947 *interp = basis->interp; 948 return CEED_ERROR_SUCCESS; 949 } 950 951 /** 952 @brief Get 1D interpolation matrix of a tensor product CeedBasis 953 954 @param basis CeedBasis 955 @param[out] interp_1d Variable to store interpolation matrix 956 957 @return An error code: 0 - success, otherwise - failure 958 959 @ref Backend 960 **/ 961 int CeedBasisGetInterp1D(CeedBasis basis, const CeedScalar **interp_1d) { 962 if (!basis->tensor_basis) 963 // LCOV_EXCL_START 964 return CeedError(basis->ceed, CEED_ERROR_MINOR, 965 "CeedBasis is not a tensor product basis."); 966 // LCOV_EXCL_STOP 967 968 *interp_1d = basis->interp_1d; 969 return CEED_ERROR_SUCCESS; 970 } 971 972 /** 973 @brief Get gradient matrix of a CeedBasis 974 975 @param basis CeedBasis 976 @param[out] grad Variable to store gradient matrix 977 978 @return An error code: 0 - success, otherwise - failure 979 980 @ref Backend 981 **/ 982 int CeedBasisGetGrad(CeedBasis basis, const CeedScalar **grad) { 983 if (!basis->grad && basis->tensor_basis) { 984 // Allocate 985 int ierr; 986 ierr = CeedMalloc(basis->dim*basis->Q*basis->P, &basis->grad); 987 CeedChk(ierr); 988 989 // Initialize 990 for (CeedInt i=0; i<basis->dim*basis->Q*basis->P; i++) 991 basis->grad[i] = 1.0; 992 993 // Calculate 994 for (CeedInt d=0; d<basis->dim; d++) 995 for (CeedInt i=0; i<basis->dim; i++) 996 for (CeedInt qpt=0; qpt<basis->Q; qpt++) 997 for (CeedInt node=0; node<basis->P; node++) { 998 CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d; 999 CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d; 1000 if (i == d) 1001 basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *= 1002 basis->grad_1d[q*basis->P_1d+p]; 1003 else 1004 basis->grad[(i*basis->Q+qpt)*(basis->P)+node] *= 1005 basis->interp_1d[q*basis->P_1d+p]; 1006 } 1007 } 1008 *grad = basis->grad; 1009 return CEED_ERROR_SUCCESS; 1010 } 1011 1012 /** 1013 @brief Get 1D gradient matrix of a tensor product CeedBasis 1014 1015 @param basis CeedBasis 1016 @param[out] grad_1d Variable to store gradient matrix 1017 1018 @return An error code: 0 - success, otherwise - failure 1019 1020 @ref Backend 1021 **/ 1022 int CeedBasisGetGrad1D(CeedBasis basis, const CeedScalar **grad_1d) { 1023 if (!basis->tensor_basis) 1024 // LCOV_EXCL_START 1025 return CeedError(basis->ceed, CEED_ERROR_MINOR, 1026 "CeedBasis is not a tensor product basis."); 1027 // LCOV_EXCL_STOP 1028 1029 *grad_1d = basis->grad_1d; 1030 return CEED_ERROR_SUCCESS; 1031 } 1032 1033 /** 1034 @brief Destroy a CeedBasis 1035 1036 @param basis CeedBasis to destroy 1037 1038 @return An error code: 0 - success, otherwise - failure 1039 1040 @ref User 1041 **/ 1042 int CeedBasisDestroy(CeedBasis *basis) { 1043 int ierr; 1044 1045 if (!*basis || --(*basis)->ref_count > 0) return CEED_ERROR_SUCCESS; 1046 if ((*basis)->Destroy) { 1047 ierr = (*basis)->Destroy(*basis); CeedChk(ierr); 1048 } 1049 if ((*basis)->contract) { 1050 ierr = CeedTensorContractDestroy(&(*basis)->contract); CeedChk(ierr); 1051 } 1052 ierr = CeedFree(&(*basis)->interp); CeedChk(ierr); 1053 ierr = CeedFree(&(*basis)->interp_1d); CeedChk(ierr); 1054 ierr = CeedFree(&(*basis)->grad); CeedChk(ierr); 1055 ierr = CeedFree(&(*basis)->grad_1d); CeedChk(ierr); 1056 ierr = CeedFree(&(*basis)->q_ref_1d); CeedChk(ierr); 1057 ierr = CeedFree(&(*basis)->q_weight_1d); CeedChk(ierr); 1058 ierr = CeedDestroy(&(*basis)->ceed); CeedChk(ierr); 1059 ierr = CeedFree(basis); CeedChk(ierr); 1060 return CEED_ERROR_SUCCESS; 1061 } 1062 1063 /** 1064 @brief Construct a Gauss-Legendre quadrature 1065 1066 @param Q Number of quadrature points (integrates polynomials of 1067 degree 2*Q-1 exactly) 1068 @param[out] q_ref_1d Array of length Q to hold the abscissa on [-1, 1] 1069 @param[out] q_weight_1d Array of length Q to hold the weights 1070 1071 @return An error code: 0 - success, otherwise - failure 1072 1073 @ref Utility 1074 **/ 1075 int CeedGaussQuadrature(CeedInt Q, CeedScalar *q_ref_1d, 1076 CeedScalar *q_weight_1d) { 1077 // Allocate 1078 CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0*atan(1.0); 1079 // Build q_ref_1d, q_weight_1d 1080 for (int i = 0; i <= Q/2; i++) { 1081 // Guess 1082 xi = cos(PI*(CeedScalar)(2*i+1)/((CeedScalar)(2*Q))); 1083 // Pn(xi) 1084 P0 = 1.0; 1085 P1 = xi; 1086 P2 = 0.0; 1087 for (int j = 2; j <= Q; j++) { 1088 P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 1089 P0 = P1; 1090 P1 = P2; 1091 } 1092 // First Newton Step 1093 dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 1094 xi = xi-P2/dP2; 1095 // Newton to convergence 1096 for (int k=0; k<100 && fabs(P2)>10*CEED_EPSILON; k++) { 1097 P0 = 1.0; 1098 P1 = xi; 1099 for (int j = 2; j <= Q; j++) { 1100 P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 1101 P0 = P1; 1102 P1 = P2; 1103 } 1104 dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 1105 xi = xi-P2/dP2; 1106 } 1107 // Save xi, wi 1108 wi = 2.0/((1.0-xi*xi)*dP2*dP2); 1109 q_weight_1d[i] = wi; 1110 q_weight_1d[Q-1-i] = wi; 1111 q_ref_1d[i] = -xi; 1112 q_ref_1d[Q-1-i]= xi; 1113 } 1114 return CEED_ERROR_SUCCESS; 1115 } 1116 1117 /** 1118 @brief Construct a Gauss-Legendre-Lobatto quadrature 1119 1120 @param Q Number of quadrature points (integrates polynomials of 1121 degree 2*Q-3 exactly) 1122 @param[out] q_ref_1d Array of length Q to hold the abscissa on [-1, 1] 1123 @param[out] q_weight_1d Array of length Q to hold the weights 1124 1125 @return An error code: 0 - success, otherwise - failure 1126 1127 @ref Utility 1128 **/ 1129 int CeedLobattoQuadrature(CeedInt Q, CeedScalar *q_ref_1d, 1130 CeedScalar *q_weight_1d) { 1131 // Allocate 1132 CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0*atan(1.0); 1133 // Build q_ref_1d, q_weight_1d 1134 // Set endpoints 1135 if (Q < 2) 1136 // LCOV_EXCL_START 1137 return CeedError(NULL, CEED_ERROR_DIMENSION, 1138 "Cannot create Lobatto quadrature with Q=%d < 2 points", Q); 1139 // LCOV_EXCL_STOP 1140 wi = 2.0/((CeedScalar)(Q*(Q-1))); 1141 if (q_weight_1d) { 1142 q_weight_1d[0] = wi; 1143 q_weight_1d[Q-1] = wi; 1144 } 1145 q_ref_1d[0] = -1.0; 1146 q_ref_1d[Q-1] = 1.0; 1147 // Interior 1148 for (int i = 1; i <= (Q-1)/2; i++) { 1149 // Guess 1150 xi = cos(PI*(CeedScalar)(i)/(CeedScalar)(Q-1)); 1151 // Pn(xi) 1152 P0 = 1.0; 1153 P1 = xi; 1154 P2 = 0.0; 1155 for (int j = 2; j < Q; j++) { 1156 P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 1157 P0 = P1; 1158 P1 = P2; 1159 } 1160 // First Newton step 1161 dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 1162 d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi); 1163 xi = xi-dP2/d2P2; 1164 // Newton to convergence 1165 for (int k=0; k<100 && fabs(dP2)>10*CEED_EPSILON; k++) { 1166 P0 = 1.0; 1167 P1 = xi; 1168 for (int j = 2; j < Q; j++) { 1169 P2 = (((CeedScalar)(2*j-1))*xi*P1-((CeedScalar)(j-1))*P0)/((CeedScalar)(j)); 1170 P0 = P1; 1171 P1 = P2; 1172 } 1173 dP2 = (xi*P2 - P0)*(CeedScalar)Q/(xi*xi-1.0); 1174 d2P2 = (2*xi*dP2 - (CeedScalar)(Q*(Q-1))*P2)/(1.0-xi*xi); 1175 xi = xi-dP2/d2P2; 1176 } 1177 // Save xi, wi 1178 wi = 2.0/(((CeedScalar)(Q*(Q-1)))*P2*P2); 1179 if (q_weight_1d) { 1180 q_weight_1d[i] = wi; 1181 q_weight_1d[Q-1-i] = wi; 1182 } 1183 q_ref_1d[i] = -xi; 1184 q_ref_1d[Q-1-i]= xi; 1185 } 1186 return CEED_ERROR_SUCCESS; 1187 } 1188 1189 /** 1190 @brief Return QR Factorization of a matrix 1191 1192 @param ceed A Ceed context for error handling 1193 @param[in,out] mat Row-major matrix to be factorized in place 1194 @param[in,out] tau Vector of length m of scaling factors 1195 @param m Number of rows 1196 @param n Number of columns 1197 1198 @return An error code: 0 - success, otherwise - failure 1199 1200 @ref Utility 1201 **/ 1202 int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau, 1203 CeedInt m, CeedInt n) { 1204 CeedScalar v[m]; 1205 1206 // Check m >= n 1207 if (n > m) 1208 // LCOV_EXCL_START 1209 return CeedError(ceed, CEED_ERROR_UNSUPPORTED, 1210 "Cannot compute QR factorization with n > m"); 1211 // LCOV_EXCL_STOP 1212 1213 for (CeedInt i=0; i<n; i++) { 1214 // Calculate Householder vector, magnitude 1215 CeedScalar sigma = 0.0; 1216 v[i] = mat[i+n*i]; 1217 for (CeedInt j=i+1; j<m; j++) { 1218 v[j] = mat[i+n*j]; 1219 sigma += v[j] * v[j]; 1220 } 1221 CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:m] 1222 CeedScalar Rii = -copysign(norm, v[i]); 1223 v[i] -= Rii; 1224 // norm of v[i:m] after modification above and scaling below 1225 // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 1226 // tau = 2 / (norm*norm) 1227 tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma); 1228 1229 for (CeedInt j=i+1; j<m; j++) 1230 v[j] /= v[i]; 1231 1232 // Apply Householder reflector to lower right panel 1233 CeedHouseholderReflect(&mat[i*n+i+1], &v[i], tau[i], m-i, n-i-1, n, 1); 1234 // Save v 1235 mat[i+n*i] = Rii; 1236 for (CeedInt j=i+1; j<m; j++) 1237 mat[i+n*j] = v[j]; 1238 } 1239 return CEED_ERROR_SUCCESS; 1240 } 1241 1242 /** 1243 @brief Return symmetric Schur decomposition of the symmetric matrix mat via 1244 symmetric QR factorization 1245 1246 @param ceed A Ceed context for error handling 1247 @param[in,out] mat Row-major matrix to be factorized in place 1248 @param[out] lambda Vector of length n of eigenvalues 1249 @param n Number of rows/columns 1250 1251 @return An error code: 0 - success, otherwise - failure 1252 1253 @ref Utility 1254 **/ 1255 int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat, 1256 CeedScalar *lambda, CeedInt n) { 1257 // Check bounds for clang-tidy 1258 if (n<2) 1259 // LCOV_EXCL_START 1260 return CeedError(ceed, CEED_ERROR_UNSUPPORTED, 1261 "Cannot compute symmetric Schur decomposition of scalars"); 1262 // LCOV_EXCL_STOP 1263 1264 CeedScalar v[n-1], tau[n-1], matT[n*n]; 1265 1266 // Copy mat to matT and set mat to I 1267 memcpy(matT, mat, n*n*sizeof(mat[0])); 1268 for (CeedInt i=0; i<n; i++) 1269 for (CeedInt j=0; j<n; j++) 1270 mat[j+n*i] = (i==j) ? 1 : 0; 1271 1272 // Reduce to tridiagonal 1273 for (CeedInt i=0; i<n-1; i++) { 1274 // Calculate Householder vector, magnitude 1275 CeedScalar sigma = 0.0; 1276 v[i] = matT[i+n*(i+1)]; 1277 for (CeedInt j=i+1; j<n-1; j++) { 1278 v[j] = matT[i+n*(j+1)]; 1279 sigma += v[j] * v[j]; 1280 } 1281 CeedScalar norm = sqrt(v[i]*v[i] + sigma); // norm of v[i:n-1] 1282 CeedScalar Rii = -copysign(norm, v[i]); 1283 v[i] -= Rii; 1284 // norm of v[i:m] after modification above and scaling below 1285 // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 1286 // tau = 2 / (norm*norm) 1287 if (sigma > 10*CEED_EPSILON) 1288 tau[i] = 2 * v[i]*v[i] / (v[i]*v[i] + sigma); 1289 else 1290 tau[i] = 0; 1291 1292 for (CeedInt j=i+1; j<n-1; j++) 1293 v[j] /= v[i]; 1294 1295 // Update sub and super diagonal 1296 matT[i+n*(i+1)] = Rii; 1297 matT[(i+1)+n*i] = Rii; 1298 for (CeedInt j=i+2; j<n; j++) { 1299 matT[i+n*j] = 0; matT[j+n*i] = 0; 1300 } 1301 // Apply symmetric Householder reflector to lower right panel 1302 CeedHouseholderReflect(&matT[(i+1)+n*(i+1)], &v[i], tau[i], 1303 n-(i+1), n-(i+1), n, 1); 1304 CeedHouseholderReflect(&matT[(i+1)+n*(i+1)], &v[i], tau[i], 1305 n-(i+1), n-(i+1), 1, n); 1306 // Save v 1307 for (CeedInt j=i+1; j<n-1; j++) { 1308 matT[i+n*(j+1)] = v[j]; 1309 } 1310 } 1311 // Backwards accumulation of Q 1312 for (CeedInt i=n-2; i>=0; i--) { 1313 v[i] = 1; 1314 for (CeedInt j=i+1; j<n-1; j++) { 1315 v[j] = matT[i+n*(j+1)]; 1316 matT[i+n*(j+1)] = 0; 1317 } 1318 CeedHouseholderReflect(&mat[(i+1)+n*(i+1)], &v[i], tau[i], 1319 n-(i+1), n-(i+1), n, 1); 1320 } 1321 1322 // Reduce sub and super diagonal 1323 CeedInt p = 0, q = 0, itr = 0, maxitr = n*n*n; 1324 CeedScalar tol = 10*CEED_EPSILON; 1325 1326 while (q < n && itr < maxitr) { 1327 // Update p, q, size of reduced portions of diagonal 1328 p = 0; q = 0; 1329 for (CeedInt i=n-2; i>=0; i--) { 1330 if (fabs(matT[i+n*(i+1)]) < tol) 1331 q += 1; 1332 else 1333 break; 1334 } 1335 for (CeedInt i=0; i<n-1-q; i++) { 1336 if (fabs(matT[i+n*(i+1)]) < tol) 1337 p += 1; 1338 else 1339 break; 1340 } 1341 if (q == n-1) break; // Finished reducing 1342 1343 // Reduce tridiagonal portion 1344 CeedScalar tnn = matT[(n-1-q)+n*(n-1-q)], 1345 tnnm1 = matT[(n-2-q)+n*(n-1-q)]; 1346 CeedScalar d = (matT[(n-2-q)+n*(n-2-q)] - tnn)/2; 1347 CeedScalar mu = tnn - tnnm1*tnnm1 / 1348 (d + copysign(sqrt(d*d + tnnm1*tnnm1), d)); 1349 CeedScalar x = matT[p+n*p] - mu; 1350 CeedScalar z = matT[p+n*(p+1)]; 1351 for (CeedInt k=p; k<n-1-q; k++) { 1352 // Compute Givens rotation 1353 CeedScalar c = 1, s = 0; 1354 if (fabs(z) > tol) { 1355 if (fabs(z) > fabs(x)) { 1356 CeedScalar tau = -x/z; 1357 s = 1/sqrt(1+tau*tau), c = s*tau; 1358 } else { 1359 CeedScalar tau = -z/x; 1360 c = 1/sqrt(1+tau*tau), s = c*tau; 1361 } 1362 } 1363 1364 // Apply Givens rotation to T 1365 CeedGivensRotation(matT, c, s, CEED_NOTRANSPOSE, k, k+1, n, n); 1366 CeedGivensRotation(matT, c, s, CEED_TRANSPOSE, k, k+1, n, n); 1367 1368 // Apply Givens rotation to Q 1369 CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k+1, n, n); 1370 1371 // Update x, z 1372 if (k < n-q-2) { 1373 x = matT[k+n*(k+1)]; 1374 z = matT[k+n*(k+2)]; 1375 } 1376 } 1377 itr++; 1378 } 1379 // Save eigenvalues 1380 for (CeedInt i=0; i<n; i++) 1381 lambda[i] = matT[i+n*i]; 1382 1383 // Check convergence 1384 if (itr == maxitr && q < n-1) 1385 // LCOV_EXCL_START 1386 return CeedError(ceed, CEED_ERROR_MINOR, 1387 "Symmetric QR failed to converge"); 1388 // LCOV_EXCL_STOP 1389 return CEED_ERROR_SUCCESS; 1390 } 1391 1392 /** 1393 @brief Return Simultaneous Diagonalization of two matrices. This solves the 1394 generalized eigenvalue problem A x = lambda B x, where A and B 1395 are symmetric and B is positive definite. We generate the matrix X 1396 and vector Lambda such that X^T A X = Lambda and X^T B X = I. This 1397 is equivalent to the LAPACK routine 'sygv' with TYPE = 1. 1398 1399 @param ceed A Ceed context for error handling 1400 @param[in] mat_A Row-major matrix to be factorized with eigenvalues 1401 @param[in] mat_B Row-major matrix to be factorized to identity 1402 @param[out] x Row-major orthogonal matrix 1403 @param[out] lambda Vector of length n of generalized eigenvalues 1404 @param n Number of rows/columns 1405 1406 @return An error code: 0 - success, otherwise - failure 1407 1408 @ref Utility 1409 **/ 1410 int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *mat_A, 1411 CeedScalar *mat_B, CeedScalar *x, 1412 CeedScalar *lambda, CeedInt n) { 1413 int ierr; 1414 CeedScalar mat_C[n*n], matG[n*n], vecD[n]; 1415 1416 // Compute B = G D G^T 1417 memcpy(matG, mat_B, n*n*sizeof(mat_B[0])); 1418 ierr = CeedSymmetricSchurDecomposition(ceed, matG, vecD, n); CeedChk(ierr); 1419 for (CeedInt i=0; i<n; i++) 1420 vecD[i] = sqrt(vecD[i]); 1421 1422 // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T 1423 // = D^-1/2 G^T A G D^-1/2 1424 for (CeedInt i=0; i<n; i++) 1425 for (CeedInt j=0; j<n; j++) 1426 mat_C[j+i*n] = matG[i+j*n] / vecD[i]; 1427 ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)mat_C, 1428 (const CeedScalar *)mat_A, x, n, n, n); 1429 CeedChk(ierr); 1430 for (CeedInt i=0; i<n; i++) 1431 for (CeedInt j=0; j<n; j++) 1432 matG[j+i*n] = matG[j+i*n] / vecD[j]; 1433 ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)x, 1434 (const CeedScalar *)matG, mat_C, n, n, n); 1435 CeedChk(ierr); 1436 1437 // Compute Q^T C Q = lambda 1438 ierr = CeedSymmetricSchurDecomposition(ceed, mat_C, lambda, n); CeedChk(ierr); 1439 1440 // Set x = (G D^1/2)^-T Q 1441 // = G D^-1/2 Q 1442 ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)matG, 1443 (const CeedScalar *)mat_C, x, n, n, n); 1444 CeedChk(ierr); 1445 return CEED_ERROR_SUCCESS; 1446 } 1447 1448 /// @} 1449