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