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