1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3 // 4 // SPDX-License-Identifier: BSD-2-Clause 5 // 6 // This file is part of CEED: http://github.com/ceed 7 8 #include <ceed-impl.h> 9 #include <ceed.h> 10 #include <ceed/backend.h> 11 #include <math.h> 12 #include <stdbool.h> 13 #include <stdio.h> 14 #include <string.h> 15 16 /// @file 17 /// Implementation of CeedBasis interfaces 18 19 /// @cond DOXYGEN_SKIP 20 static struct CeedBasis_private ceed_basis_collocated; 21 /// @endcond 22 23 /// @addtogroup CeedBasisUser 24 /// @{ 25 26 /// Indicate that the quadrature points are collocated with the nodes 27 const CeedBasis CEED_BASIS_COLLOCATED = &ceed_basis_collocated; 28 29 /// @} 30 31 /// ---------------------------------------------------------------------------- 32 /// CeedBasis Library Internal Functions 33 /// ---------------------------------------------------------------------------- 34 /// @addtogroup CeedBasisDeveloper 35 /// @{ 36 37 /** 38 @brief Compute Householder reflection 39 40 Computes A = (I - b v v^T) A, where A is an mxn matrix indexed as A[i*row + j*col] 41 42 @param[in,out] A Matrix to apply Householder reflection to, in place 43 @param[in] v Householder vector 44 @param[in] b Scaling factor 45 @param[in] m Number of rows in A 46 @param[in] n Number of columns in A 47 @param[in] row Row stride 48 @param[in] col Col stride 49 50 @return An error code: 0 - success, otherwise - failure 51 52 @ref Developer 53 **/ 54 static int CeedHouseholderReflect(CeedScalar *A, const CeedScalar *v, CeedScalar b, CeedInt m, CeedInt n, CeedInt row, CeedInt col) { 55 for (CeedInt j = 0; j < n; j++) { 56 CeedScalar w = A[0 * row + j * col]; 57 for (CeedInt i = 1; i < m; i++) w += v[i] * A[i * row + j * col]; 58 A[0 * row + j * col] -= b * w; 59 for (CeedInt i = 1; i < m; i++) A[i * row + j * col] -= b * w * v[i]; 60 } 61 return CEED_ERROR_SUCCESS; 62 } 63 64 /** 65 @brief Compute Givens rotation 66 67 Computes A = G A (or G^T A in transpose mode), where A is an mxn matrix indexed as A[i*n + j*m] 68 69 @param[in,out] A Row major matrix to apply Givens rotation to, in place 70 @param[in] c Cosine factor 71 @param[in] s Sine factor 72 @param[in] t_mode @ref CEED_NOTRANSPOSE to rotate the basis counter-clockwise, which has the effect of rotating columns of A clockwise; 73 @ref CEED_TRANSPOSE for the opposite rotation 74 @param[in] i First row/column to apply rotation 75 @param[in] k Second row/column to apply rotation 76 @param[in] m Number of rows in A 77 @param[in] n Number of columns in A 78 79 @return An error code: 0 - success, otherwise - failure 80 81 @ref Developer 82 **/ 83 static int CeedGivensRotation(CeedScalar *A, CeedScalar c, CeedScalar s, CeedTransposeMode t_mode, CeedInt i, CeedInt k, CeedInt m, CeedInt n) { 84 CeedInt stride_j = 1, stride_ik = m, num_its = n; 85 if (t_mode == CEED_NOTRANSPOSE) { 86 stride_j = n; 87 stride_ik = 1; 88 num_its = m; 89 } 90 91 // Apply rotation 92 for (CeedInt j = 0; j < num_its; j++) { 93 CeedScalar tau1 = A[i * stride_ik + j * stride_j], tau2 = A[k * stride_ik + j * stride_j]; 94 A[i * stride_ik + j * stride_j] = c * tau1 - s * tau2; 95 A[k * stride_ik + j * stride_j] = s * tau1 + c * tau2; 96 } 97 return CEED_ERROR_SUCCESS; 98 } 99 100 /** 101 @brief View an array stored in a CeedBasis 102 103 @param[in] name Name of array 104 @param[in] fp_fmt Printing format 105 @param[in] m Number of rows in array 106 @param[in] n Number of columns in array 107 @param[in] a Array to be viewed 108 @param[in] stream Stream to view to, e.g., stdout 109 110 @return An error code: 0 - success, otherwise - failure 111 112 @ref Developer 113 **/ 114 static int CeedScalarView(const char *name, const char *fp_fmt, CeedInt m, CeedInt n, const CeedScalar *a, FILE *stream) { 115 for (CeedInt i = 0; i < m; i++) { 116 if (m > 1) fprintf(stream, "%12s[%" CeedInt_FMT "]:", name, i); 117 else fprintf(stream, "%12s:", name); 118 for (CeedInt j = 0; j < n; j++) fprintf(stream, fp_fmt, fabs(a[i * n + j]) > 1E-14 ? a[i * n + j] : 0); 119 fputs("\n", stream); 120 } 121 return CEED_ERROR_SUCCESS; 122 } 123 124 /** 125 @brief Create the interpolation and gradient matrices for projection from the nodes of `basis_from` to the nodes of `basis_to`. 126 127 The interpolation is given by `interp_project = interp_to^+ * interp_from`, where the pesudoinverse `interp_to^+` is given by QR factorization. 128 The gradient is given by `grad_project = interp_to^+ * grad_from`. 129 Note: `basis_from` and `basis_to` must have compatible quadrature spaces. 130 131 @param[in] basis_from CeedBasis to project from 132 @param[in] basis_to CeedBasis to project to 133 @param[out] interp_project Address of the variable where the newly created interpolation matrix will be stored. 134 @param[out] grad_project Address of the variable where the newly created gradient matrix will be stored. 135 136 @return An error code: 0 - success, otherwise - failure 137 138 @ref Developer 139 **/ 140 static int CeedBasisCreateProjectionMatrices(CeedBasis basis_from, CeedBasis basis_to, CeedScalar **interp_project, CeedScalar **grad_project) { 141 Ceed ceed; 142 CeedCall(CeedBasisGetCeed(basis_to, &ceed)); 143 144 // Check for compatible quadrature spaces 145 CeedInt Q_to, Q_from; 146 CeedCall(CeedBasisGetNumQuadraturePoints(basis_to, &Q_to)); 147 CeedCall(CeedBasisGetNumQuadraturePoints(basis_from, &Q_from)); 148 if (Q_to != Q_from) { 149 // LCOV_EXCL_START 150 return CeedError(ceed, CEED_ERROR_DIMENSION, "Bases must have compatible quadrature spaces"); 151 // LCOV_EXCL_STOP 152 } 153 154 // Check for matching tensor or non-tensor 155 CeedInt P_to, P_from, Q = Q_to; 156 bool is_tensor_to, is_tensor_from; 157 CeedCall(CeedBasisIsTensor(basis_to, &is_tensor_to)); 158 CeedCall(CeedBasisIsTensor(basis_from, &is_tensor_from)); 159 if (is_tensor_to && is_tensor_from) { 160 CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_to)); 161 CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_from)); 162 CeedCall(CeedBasisGetNumQuadraturePoints1D(basis_from, &Q)); 163 } else if (!is_tensor_to && !is_tensor_from) { 164 CeedCall(CeedBasisGetNumNodes(basis_to, &P_to)); 165 CeedCall(CeedBasisGetNumNodes(basis_from, &P_from)); 166 } else { 167 // LCOV_EXCL_START 168 return CeedError(ceed, CEED_ERROR_MINOR, "Bases must both be tensor or non-tensor"); 169 // LCOV_EXCL_STOP 170 } 171 172 // Get source matrices 173 CeedInt dim; 174 CeedScalar *interp_to, *interp_from, *tau; 175 CeedCall(CeedBasisGetDimension(basis_to, &dim)); 176 CeedCall(CeedMalloc(Q * P_from, &interp_from)); 177 CeedCall(CeedMalloc(Q * P_to, &interp_to)); 178 CeedCall(CeedCalloc(P_to * P_from, interp_project)); 179 CeedCall(CeedCalloc(P_to * P_from * (is_tensor_to ? 1 : dim), grad_project)); 180 CeedCall(CeedMalloc(Q, &tau)); 181 const CeedScalar *interp_to_source = NULL, *interp_from_source = NULL, *grad_from_source; 182 if (is_tensor_to) { 183 CeedCall(CeedBasisGetInterp1D(basis_to, &interp_to_source)); 184 CeedCall(CeedBasisGetInterp1D(basis_from, &interp_from_source)); 185 CeedCall(CeedBasisGetGrad1D(basis_from, &grad_from_source)); 186 } else { 187 CeedCall(CeedBasisGetInterp(basis_to, &interp_to_source)); 188 CeedCall(CeedBasisGetInterp(basis_from, &interp_from_source)); 189 CeedCall(CeedBasisGetGrad(basis_from, &grad_from_source)); 190 } 191 192 // Build matrices 193 CeedInt num_matrices = 1 + (is_tensor_to ? 1 : dim); 194 CeedScalar *input_from[num_matrices], *output_project[num_matrices]; 195 input_from[0] = (CeedScalar *)interp_from_source; 196 output_project[0] = *interp_project; 197 for (CeedInt m = 1; m < num_matrices; m++) { 198 input_from[m] = (CeedScalar *)&grad_from_source[(m - 1) * Q * P_from]; 199 output_project[m] = &((*grad_project)[(m - 1) * P_to * P_from]); 200 } 201 for (CeedInt m = 0; m < num_matrices; m++) { 202 // -- QR Factorization, interp_to = Q R 203 memcpy(interp_to, interp_to_source, Q * P_to * sizeof(interp_to_source[0])); 204 CeedCall(CeedQRFactorization(ceed, interp_to, tau, Q, P_to)); 205 206 // -- Apply Qtranspose, interp_to = Qtranspose interp_from 207 memcpy(interp_from, input_from[m], Q * P_from * sizeof(input_from[m][0])); 208 CeedCall(CeedHouseholderApplyQ(interp_from, interp_to, tau, CEED_TRANSPOSE, Q, P_from, P_to, P_from, 1)); 209 210 // -- Apply Rinv, interp_project = Rinv interp_c 211 for (CeedInt j = 0; j < P_from; j++) { // Column j 212 output_project[m][j + P_from * (P_to - 1)] = interp_from[j + P_from * (P_to - 1)] / interp_to[P_to * P_to - 1]; 213 for (CeedInt i = P_to - 2; i >= 0; i--) { // Row i 214 output_project[m][j + P_from * i] = interp_from[j + P_from * i]; 215 for (CeedInt k = i + 1; k < P_to; k++) { 216 output_project[m][j + P_from * i] -= interp_to[k + P_to * i] * output_project[m][j + P_from * k]; 217 } 218 output_project[m][j + P_from * i] /= interp_to[i + P_to * i]; 219 } 220 } 221 } 222 223 // Cleanup 224 CeedCall(CeedFree(&tau)); 225 CeedCall(CeedFree(&interp_to)); 226 CeedCall(CeedFree(&interp_from)); 227 228 return CEED_ERROR_SUCCESS; 229 } 230 231 /// @} 232 233 /// ---------------------------------------------------------------------------- 234 /// Ceed Backend API 235 /// ---------------------------------------------------------------------------- 236 /// @addtogroup CeedBasisBackend 237 /// @{ 238 239 /** 240 @brief Return collocated grad matrix 241 242 @param[in] basis CeedBasis 243 @param[out] collo_grad_1d Row-major (Q_1d * Q_1d) matrix expressing derivatives of basis functions at quadrature points 244 245 @return An error code: 0 - success, otherwise - failure 246 247 @ref Backend 248 **/ 249 int CeedBasisGetCollocatedGrad(CeedBasis basis, CeedScalar *collo_grad_1d) { 250 int i, j, k; 251 Ceed ceed; 252 CeedInt P_1d = (basis)->P_1d, Q_1d = (basis)->Q_1d; 253 CeedScalar *interp_1d, *grad_1d, *tau; 254 255 CeedCall(CeedMalloc(Q_1d * P_1d, &interp_1d)); 256 CeedCall(CeedMalloc(Q_1d * P_1d, &grad_1d)); 257 CeedCall(CeedMalloc(Q_1d, &tau)); 258 memcpy(interp_1d, (basis)->interp_1d, Q_1d * P_1d * sizeof(basis)->interp_1d[0]); 259 memcpy(grad_1d, (basis)->grad_1d, Q_1d * P_1d * sizeof(basis)->interp_1d[0]); 260 261 // QR Factorization, interp_1d = Q R 262 CeedCall(CeedBasisGetCeed(basis, &ceed)); 263 CeedCall(CeedQRFactorization(ceed, interp_1d, tau, Q_1d, P_1d)); 264 // Note: This function is for backend use, so all errors are terminal and we do not need to clean up memory on failure. 265 266 // Apply Rinv, collo_grad_1d = grad_1d Rinv 267 for (i = 0; i < Q_1d; i++) { // Row i 268 collo_grad_1d[Q_1d * i] = grad_1d[P_1d * i] / interp_1d[0]; 269 for (j = 1; j < P_1d; j++) { // Column j 270 collo_grad_1d[j + Q_1d * i] = grad_1d[j + P_1d * i]; 271 for (k = 0; k < j; k++) collo_grad_1d[j + Q_1d * i] -= interp_1d[j + P_1d * k] * collo_grad_1d[k + Q_1d * i]; 272 collo_grad_1d[j + Q_1d * i] /= interp_1d[j + P_1d * j]; 273 } 274 for (j = P_1d; j < Q_1d; j++) collo_grad_1d[j + Q_1d * i] = 0; 275 } 276 277 // Apply Qtranspose, collo_grad = collo_grad Q_transpose 278 CeedCall(CeedHouseholderApplyQ(collo_grad_1d, interp_1d, tau, CEED_NOTRANSPOSE, Q_1d, Q_1d, P_1d, 1, Q_1d)); 279 280 CeedCall(CeedFree(&interp_1d)); 281 CeedCall(CeedFree(&grad_1d)); 282 CeedCall(CeedFree(&tau)); 283 return CEED_ERROR_SUCCESS; 284 } 285 286 /** 287 @brief Get tensor status for given CeedBasis 288 289 @param[in] basis CeedBasis 290 @param[out] is_tensor Variable to store tensor status 291 292 @return An error code: 0 - success, otherwise - failure 293 294 @ref Backend 295 **/ 296 int CeedBasisIsTensor(CeedBasis basis, bool *is_tensor) { 297 *is_tensor = basis->tensor_basis; 298 return CEED_ERROR_SUCCESS; 299 } 300 301 /** 302 @brief Get backend data of a CeedBasis 303 304 @param[in] basis CeedBasis 305 @param[out] data Variable to store data 306 307 @return An error code: 0 - success, otherwise - failure 308 309 @ref Backend 310 **/ 311 int CeedBasisGetData(CeedBasis basis, void *data) { 312 *(void **)data = basis->data; 313 return CEED_ERROR_SUCCESS; 314 } 315 316 /** 317 @brief Set backend data of a CeedBasis 318 319 @param[in,out] basis CeedBasis 320 @param[in] data Data to set 321 322 @return An error code: 0 - success, otherwise - failure 323 324 @ref Backend 325 **/ 326 int CeedBasisSetData(CeedBasis basis, void *data) { 327 basis->data = data; 328 return CEED_ERROR_SUCCESS; 329 } 330 331 /** 332 @brief Increment the reference counter for a CeedBasis 333 334 @param[in,out] basis Basis to increment the reference counter 335 336 @return An error code: 0 - success, otherwise - failure 337 338 @ref Backend 339 **/ 340 int CeedBasisReference(CeedBasis basis) { 341 basis->ref_count++; 342 return CEED_ERROR_SUCCESS; 343 } 344 345 /** 346 @brief Get number of Q-vector components for given CeedBasis 347 348 @param[in] basis CeedBasis 349 @param[in] eval_mode \ref CEED_EVAL_INTERP to use interpolated values, 350 \ref CEED_EVAL_GRAD to use gradients, 351 \ref CEED_EVAL_DIV to use divergence, 352 \ref CEED_EVAL_CURL to use curl. 353 @param[out] q_comp Variable to store number of Q-vector components of basis 354 355 @return An error code: 0 - success, otherwise - failure 356 357 @ref Backend 358 **/ 359 int CeedBasisGetNumQuadratureComponents(CeedBasis basis, CeedEvalMode eval_mode, CeedInt *q_comp) { 360 switch (eval_mode) { 361 case CEED_EVAL_INTERP: 362 *q_comp = (basis->fe_space == CEED_FE_SPACE_H1) ? 1 : basis->dim; 363 break; 364 case CEED_EVAL_GRAD: 365 *q_comp = basis->dim; 366 break; 367 case CEED_EVAL_DIV: 368 *q_comp = 1; 369 break; 370 case CEED_EVAL_CURL: 371 *q_comp = (basis->dim < 3) ? 1 : basis->dim; 372 break; 373 case CEED_EVAL_NONE: 374 case CEED_EVAL_WEIGHT: 375 *q_comp = 0; 376 break; 377 } 378 return CEED_ERROR_SUCCESS; 379 } 380 381 /** 382 @brief Estimate number of FLOPs required to apply CeedBasis in t_mode and eval_mode 383 384 @param[in] basis Basis to estimate FLOPs for 385 @param[in] t_mode Apply basis or transpose 386 @param[in] eval_mode Basis evaluation mode 387 @param[out] flops Address of variable to hold FLOPs estimate 388 389 @ref Backend 390 **/ 391 int CeedBasisGetFlopsEstimate(CeedBasis basis, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedSize *flops) { 392 bool is_tensor; 393 394 CeedCall(CeedBasisIsTensor(basis, &is_tensor)); 395 if (is_tensor) { 396 CeedInt dim, num_comp, P_1d, Q_1d; 397 CeedCall(CeedBasisGetDimension(basis, &dim)); 398 CeedCall(CeedBasisGetNumComponents(basis, &num_comp)); 399 CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d)); 400 CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 401 if (t_mode == CEED_TRANSPOSE) { 402 P_1d = Q_1d; 403 Q_1d = P_1d; 404 } 405 CeedInt tensor_flops = 0, pre = num_comp * CeedIntPow(P_1d, dim - 1), post = 1; 406 for (CeedInt d = 0; d < dim; d++) { 407 tensor_flops += 2 * pre * P_1d * post * Q_1d; 408 pre /= P_1d; 409 post *= Q_1d; 410 } 411 switch (eval_mode) { 412 case CEED_EVAL_NONE: 413 *flops = 0; 414 break; 415 case CEED_EVAL_INTERP: 416 *flops = tensor_flops; 417 break; 418 case CEED_EVAL_GRAD: 419 *flops = tensor_flops * 2; 420 break; 421 case CEED_EVAL_DIV: 422 // LCOV_EXCL_START 423 return CeedError(basis->ceed, CEED_ERROR_INCOMPATIBLE, "Tensor CEED_EVAL_DIV not supported"); 424 break; 425 case CEED_EVAL_CURL: 426 return CeedError(basis->ceed, CEED_ERROR_INCOMPATIBLE, "Tensor CEED_EVAL_CURL not supported"); 427 break; 428 // LCOV_EXCL_STOP 429 case CEED_EVAL_WEIGHT: 430 *flops = dim * CeedIntPow(Q_1d, dim); 431 break; 432 } 433 } else { 434 CeedInt dim, num_comp, q_comp, num_nodes, num_qpts; 435 CeedCall(CeedBasisGetDimension(basis, &dim)); 436 CeedCall(CeedBasisGetNumComponents(basis, &num_comp)); 437 CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp)); 438 CeedCall(CeedBasisGetNumNodes(basis, &num_nodes)); 439 CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts)); 440 switch (eval_mode) { 441 case CEED_EVAL_NONE: 442 *flops = 0; 443 break; 444 case CEED_EVAL_INTERP: 445 case CEED_EVAL_GRAD: 446 case CEED_EVAL_DIV: 447 case CEED_EVAL_CURL: 448 *flops = num_nodes * num_qpts * num_comp * q_comp; 449 break; 450 case CEED_EVAL_WEIGHT: 451 *flops = 0; 452 break; 453 } 454 } 455 456 return CEED_ERROR_SUCCESS; 457 } 458 459 /** 460 @brief Get CeedFESpace for a CeedBasis 461 462 @param[in] basis CeedBasis 463 @param[out] fe_space Variable to store CeedFESpace 464 465 @return An error code: 0 - success, otherwise - failure 466 467 @ref Backend 468 **/ 469 int CeedBasisGetFESpace(CeedBasis basis, CeedFESpace *fe_space) { 470 *fe_space = basis->fe_space; 471 return CEED_ERROR_SUCCESS; 472 } 473 474 /** 475 @brief Get dimension for given CeedElemTopology 476 477 @param[in] topo CeedElemTopology 478 @param[out] dim Variable to store dimension of topology 479 480 @return An error code: 0 - success, otherwise - failure 481 482 @ref Backend 483 **/ 484 int CeedBasisGetTopologyDimension(CeedElemTopology topo, CeedInt *dim) { 485 *dim = (CeedInt)topo >> 16; 486 return CEED_ERROR_SUCCESS; 487 } 488 489 /** 490 @brief Get CeedTensorContract of a CeedBasis 491 492 @param[in] basis CeedBasis 493 @param[out] contract Variable to store CeedTensorContract 494 495 @return An error code: 0 - success, otherwise - failure 496 497 @ref Backend 498 **/ 499 int CeedBasisGetTensorContract(CeedBasis basis, CeedTensorContract *contract) { 500 *contract = basis->contract; 501 return CEED_ERROR_SUCCESS; 502 } 503 504 /** 505 @brief Set CeedTensorContract of a CeedBasis 506 507 @param[in,out] basis CeedBasis 508 @param[in] contract CeedTensorContract to set 509 510 @return An error code: 0 - success, otherwise - failure 511 512 @ref Backend 513 **/ 514 int CeedBasisSetTensorContract(CeedBasis basis, CeedTensorContract contract) { 515 basis->contract = contract; 516 CeedCall(CeedTensorContractReference(contract)); 517 return CEED_ERROR_SUCCESS; 518 } 519 520 /** 521 @brief Return a reference implementation of matrix multiplication C = A B. 522 523 Note: This is a reference implementation for CPU CeedScalar pointers that is not intended for high performance. 524 525 @param[in] ceed Ceed context for error handling 526 @param[in] mat_A Row-major matrix A 527 @param[in] mat_B Row-major matrix B 528 @param[out] mat_C Row-major output matrix C 529 @param[in] m Number of rows of C 530 @param[in] n Number of columns of C 531 @param[in] kk Number of columns of A/rows of B 532 533 @return An error code: 0 - success, otherwise - failure 534 535 @ref Utility 536 **/ 537 int CeedMatrixMatrixMultiply(Ceed ceed, const CeedScalar *mat_A, const CeedScalar *mat_B, CeedScalar *mat_C, CeedInt m, CeedInt n, CeedInt kk) { 538 for (CeedInt i = 0; i < m; i++) { 539 for (CeedInt j = 0; j < n; j++) { 540 CeedScalar sum = 0; 541 for (CeedInt k = 0; k < kk; k++) sum += mat_A[k + i * kk] * mat_B[j + k * n]; 542 mat_C[j + i * n] = sum; 543 } 544 } 545 return CEED_ERROR_SUCCESS; 546 } 547 548 /** 549 @brief Return QR Factorization of a matrix 550 551 @param[in] ceed Ceed context for error handling 552 @param[in,out] mat Row-major matrix to be factorized in place 553 @param[in,out] tau Vector of length m of scaling factors 554 @param[in] m Number of rows 555 @param[in] n Number of columns 556 557 @return An error code: 0 - success, otherwise - failure 558 559 @ref Utility 560 **/ 561 int CeedQRFactorization(Ceed ceed, CeedScalar *mat, CeedScalar *tau, CeedInt m, CeedInt n) { 562 CeedScalar v[m]; 563 564 // Check matrix shape 565 if (n > m) { 566 // LCOV_EXCL_START 567 return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute QR factorization with n > m"); 568 // LCOV_EXCL_STOP 569 } 570 571 for (CeedInt i = 0; i < n; i++) { 572 if (i >= m - 1) { // last row of matrix, no reflection needed 573 tau[i] = 0.; 574 break; 575 } 576 // Calculate Householder vector, magnitude 577 CeedScalar sigma = 0.0; 578 v[i] = mat[i + n * i]; 579 for (CeedInt j = i + 1; j < m; j++) { 580 v[j] = mat[i + n * j]; 581 sigma += v[j] * v[j]; 582 } 583 CeedScalar norm = sqrt(v[i] * v[i] + sigma); // norm of v[i:m] 584 CeedScalar R_ii = -copysign(norm, v[i]); 585 v[i] -= R_ii; 586 // norm of v[i:m] after modification above and scaling below 587 // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 588 // tau = 2 / (norm*norm) 589 tau[i] = 2 * v[i] * v[i] / (v[i] * v[i] + sigma); 590 for (CeedInt j = i + 1; j < m; j++) v[j] /= v[i]; 591 592 // Apply Householder reflector to lower right panel 593 CeedHouseholderReflect(&mat[i * n + i + 1], &v[i], tau[i], m - i, n - i - 1, n, 1); 594 // Save v 595 mat[i + n * i] = R_ii; 596 for (CeedInt j = i + 1; j < m; j++) mat[i + n * j] = v[j]; 597 } 598 return CEED_ERROR_SUCCESS; 599 } 600 601 /** 602 @brief Apply Householder Q matrix 603 604 Compute mat_A = mat_Q mat_A, where mat_Q is mxm and mat_A is mxn. 605 606 @param[in,out] mat_A Matrix to apply Householder Q to, in place 607 @param[in] mat_Q Householder Q matrix 608 @param[in] tau Householder scaling factors 609 @param[in] t_mode Transpose mode for application 610 @param[in] m Number of rows in A 611 @param[in] n Number of columns in A 612 @param[in] k Number of elementary reflectors in Q, k<m 613 @param[in] row Row stride in A 614 @param[in] col Col stride in A 615 616 @return An error code: 0 - success, otherwise - failure 617 618 @ref Utility 619 **/ 620 int CeedHouseholderApplyQ(CeedScalar *mat_A, const CeedScalar *mat_Q, const CeedScalar *tau, CeedTransposeMode t_mode, CeedInt m, CeedInt n, 621 CeedInt k, CeedInt row, CeedInt col) { 622 CeedScalar *v; 623 CeedCall(CeedMalloc(m, &v)); 624 for (CeedInt ii = 0; ii < k; ii++) { 625 CeedInt i = t_mode == CEED_TRANSPOSE ? ii : k - 1 - ii; 626 for (CeedInt j = i + 1; j < m; j++) v[j] = mat_Q[j * k + i]; 627 // Apply Householder reflector (I - tau v v^T) collo_grad_1d^T 628 CeedCall(CeedHouseholderReflect(&mat_A[i * row], &v[i], tau[i], m - i, n, row, col)); 629 } 630 CeedCall(CeedFree(&v)); 631 return CEED_ERROR_SUCCESS; 632 } 633 634 /** 635 @brief Return symmetric Schur decomposition of the symmetric matrix mat via symmetric QR factorization 636 637 @param[in] ceed Ceed context for error handling 638 @param[in,out] mat Row-major matrix to be factorized in place 639 @param[out] lambda Vector of length n of eigenvalues 640 @param[in] n Number of rows/columns 641 642 @return An error code: 0 - success, otherwise - failure 643 644 @ref Utility 645 **/ 646 CeedPragmaOptimizeOff int CeedSymmetricSchurDecomposition(Ceed ceed, CeedScalar *mat, CeedScalar *lambda, CeedInt n) { 647 // Check bounds for clang-tidy 648 if (n < 2) { 649 // LCOV_EXCL_START 650 return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Cannot compute symmetric Schur decomposition of scalars"); 651 // LCOV_EXCL_STOP 652 } 653 654 CeedScalar v[n - 1], tau[n - 1], mat_T[n * n]; 655 656 // Copy mat to mat_T and set mat to I 657 memcpy(mat_T, mat, n * n * sizeof(mat[0])); 658 for (CeedInt i = 0; i < n; i++) { 659 for (CeedInt j = 0; j < n; j++) mat[j + n * i] = (i == j) ? 1 : 0; 660 } 661 662 // Reduce to tridiagonal 663 for (CeedInt i = 0; i < n - 1; i++) { 664 // Calculate Householder vector, magnitude 665 CeedScalar sigma = 0.0; 666 v[i] = mat_T[i + n * (i + 1)]; 667 for (CeedInt j = i + 1; j < n - 1; j++) { 668 v[j] = mat_T[i + n * (j + 1)]; 669 sigma += v[j] * v[j]; 670 } 671 CeedScalar norm = sqrt(v[i] * v[i] + sigma); // norm of v[i:n-1] 672 CeedScalar R_ii = -copysign(norm, v[i]); 673 v[i] -= R_ii; 674 // norm of v[i:m] after modification above and scaling below 675 // norm = sqrt(v[i]*v[i] + sigma) / v[i]; 676 // tau = 2 / (norm*norm) 677 tau[i] = i == n - 2 ? 2 : 2 * v[i] * v[i] / (v[i] * v[i] + sigma); 678 for (CeedInt j = i + 1; j < n - 1; j++) v[j] /= v[i]; 679 680 // Update sub and super diagonal 681 for (CeedInt j = i + 2; j < n; j++) { 682 mat_T[i + n * j] = 0; 683 mat_T[j + n * i] = 0; 684 } 685 // Apply symmetric Householder reflector to lower right panel 686 CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1); 687 CeedHouseholderReflect(&mat_T[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), 1, n); 688 689 // Save v 690 mat_T[i + n * (i + 1)] = R_ii; 691 mat_T[(i + 1) + n * i] = R_ii; 692 for (CeedInt j = i + 1; j < n - 1; j++) { 693 mat_T[i + n * (j + 1)] = v[j]; 694 } 695 } 696 // Backwards accumulation of Q 697 for (CeedInt i = n - 2; i >= 0; i--) { 698 if (tau[i] > 0.0) { 699 v[i] = 1; 700 for (CeedInt j = i + 1; j < n - 1; j++) { 701 v[j] = mat_T[i + n * (j + 1)]; 702 mat_T[i + n * (j + 1)] = 0; 703 } 704 CeedHouseholderReflect(&mat[(i + 1) + n * (i + 1)], &v[i], tau[i], n - (i + 1), n - (i + 1), n, 1); 705 } 706 } 707 708 // Reduce sub and super diagonal 709 CeedInt p = 0, q = 0, itr = 0, max_itr = n * n * n * n; 710 CeedScalar tol = CEED_EPSILON; 711 712 while (itr < max_itr) { 713 // Update p, q, size of reduced portions of diagonal 714 p = 0; 715 q = 0; 716 for (CeedInt i = n - 2; i >= 0; i--) { 717 if (fabs(mat_T[i + n * (i + 1)]) < tol) q += 1; 718 else break; 719 } 720 for (CeedInt i = 0; i < n - q - 1; i++) { 721 if (fabs(mat_T[i + n * (i + 1)]) < tol) p += 1; 722 else break; 723 } 724 if (q == n - 1) break; // Finished reducing 725 726 // Reduce tridiagonal portion 727 CeedScalar t_nn = mat_T[(n - 1 - q) + n * (n - 1 - q)], t_nnm1 = mat_T[(n - 2 - q) + n * (n - 1 - q)]; 728 CeedScalar d = (mat_T[(n - 2 - q) + n * (n - 2 - q)] - t_nn) / 2; 729 CeedScalar mu = t_nn - t_nnm1 * t_nnm1 / (d + copysign(sqrt(d * d + t_nnm1 * t_nnm1), d)); 730 CeedScalar x = mat_T[p + n * p] - mu; 731 CeedScalar z = mat_T[p + n * (p + 1)]; 732 for (CeedInt k = p; k < n - q - 1; k++) { 733 // Compute Givens rotation 734 CeedScalar c = 1, s = 0; 735 if (fabs(z) > tol) { 736 if (fabs(z) > fabs(x)) { 737 CeedScalar tau = -x / z; 738 s = 1 / sqrt(1 + tau * tau), c = s * tau; 739 } else { 740 CeedScalar tau = -z / x; 741 c = 1 / sqrt(1 + tau * tau), s = c * tau; 742 } 743 } 744 745 // Apply Givens rotation to T 746 CeedGivensRotation(mat_T, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n); 747 CeedGivensRotation(mat_T, c, s, CEED_TRANSPOSE, k, k + 1, n, n); 748 749 // Apply Givens rotation to Q 750 CeedGivensRotation(mat, c, s, CEED_NOTRANSPOSE, k, k + 1, n, n); 751 752 // Update x, z 753 if (k < n - q - 2) { 754 x = mat_T[k + n * (k + 1)]; 755 z = mat_T[k + n * (k + 2)]; 756 } 757 } 758 itr++; 759 } 760 761 // Save eigenvalues 762 for (CeedInt i = 0; i < n; i++) lambda[i] = mat_T[i + n * i]; 763 764 // Check convergence 765 if (itr == max_itr && q < n - 1) { 766 // LCOV_EXCL_START 767 return CeedError(ceed, CEED_ERROR_MINOR, "Symmetric QR failed to converge"); 768 // LCOV_EXCL_STOP 769 } 770 return CEED_ERROR_SUCCESS; 771 } 772 CeedPragmaOptimizeOn; 773 774 /** 775 @brief Return Simultaneous Diagonalization of two matrices. 776 777 This solves the generalized eigenvalue problem A x = lambda B x, where A and B are symmetric and B is positive definite. 778 We generate the matrix X and vector Lambda such that X^T A X = Lambda and X^T B X = I. 779 This is equivalent to the LAPACK routine 'sygv' with TYPE = 1. 780 781 @param[in] ceed Ceed context for error handling 782 @param[in] mat_A Row-major matrix to be factorized with eigenvalues 783 @param[in] mat_B Row-major matrix to be factorized to identity 784 @param[out] mat_X Row-major orthogonal matrix 785 @param[out] lambda Vector of length n of generalized eigenvalues 786 @param[in] n Number of rows/columns 787 788 @return An error code: 0 - success, otherwise - failure 789 790 @ref Utility 791 **/ 792 CeedPragmaOptimizeOff int CeedSimultaneousDiagonalization(Ceed ceed, CeedScalar *mat_A, CeedScalar *mat_B, CeedScalar *mat_X, CeedScalar *lambda, 793 CeedInt n) { 794 CeedScalar *mat_C, *mat_G, *vec_D; 795 CeedCall(CeedCalloc(n * n, &mat_C)); 796 CeedCall(CeedCalloc(n * n, &mat_G)); 797 CeedCall(CeedCalloc(n, &vec_D)); 798 799 // Compute B = G D G^T 800 memcpy(mat_G, mat_B, n * n * sizeof(mat_B[0])); 801 CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_G, vec_D, n)); 802 803 // Sort eigenvalues 804 for (CeedInt i = n - 1; i >= 0; i--) { 805 for (CeedInt j = 0; j < i; j++) { 806 if (fabs(vec_D[j]) > fabs(vec_D[j + 1])) { 807 CeedScalar temp; 808 temp = vec_D[j]; 809 vec_D[j] = vec_D[j + 1]; 810 vec_D[j + 1] = temp; 811 for (CeedInt k = 0; k < n; k++) { 812 temp = mat_G[k * n + j]; 813 mat_G[k * n + j] = mat_G[k * n + j + 1]; 814 mat_G[k * n + j + 1] = temp; 815 } 816 } 817 } 818 } 819 820 // Compute C = (G D^1/2)^-1 A (G D^1/2)^-T 821 // = D^-1/2 G^T A G D^-1/2 822 // -- D = D^-1/2 823 for (CeedInt i = 0; i < n; i++) vec_D[i] = 1. / sqrt(vec_D[i]); 824 // -- G = G D^-1/2 825 // -- C = D^-1/2 G^T 826 for (CeedInt i = 0; i < n; i++) { 827 for (CeedInt j = 0; j < n; j++) { 828 mat_G[i * n + j] *= vec_D[j]; 829 mat_C[j * n + i] = mat_G[i * n + j]; 830 } 831 } 832 // -- X = (D^-1/2 G^T) A 833 CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_C, (const CeedScalar *)mat_A, mat_X, n, n, n)); 834 // -- C = (D^-1/2 G^T A) (G D^-1/2) 835 CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_X, (const CeedScalar *)mat_G, mat_C, n, n, n)); 836 837 // Compute Q^T C Q = lambda 838 CeedCall(CeedSymmetricSchurDecomposition(ceed, mat_C, lambda, n)); 839 840 // Sort eigenvalues 841 for (CeedInt i = n - 1; i >= 0; i--) { 842 for (CeedInt j = 0; j < i; j++) { 843 if (fabs(lambda[j]) > fabs(lambda[j + 1])) { 844 CeedScalar temp; 845 temp = lambda[j]; 846 lambda[j] = lambda[j + 1]; 847 lambda[j + 1] = temp; 848 for (CeedInt k = 0; k < n; k++) { 849 temp = mat_C[k * n + j]; 850 mat_C[k * n + j] = mat_C[k * n + j + 1]; 851 mat_C[k * n + j + 1] = temp; 852 } 853 } 854 } 855 } 856 857 // Set X = (G D^1/2)^-T Q 858 // = G D^-1/2 Q 859 CeedCall(CeedMatrixMatrixMultiply(ceed, (const CeedScalar *)mat_G, (const CeedScalar *)mat_C, mat_X, n, n, n)); 860 861 // Cleanup 862 CeedCall(CeedFree(&mat_C)); 863 CeedCall(CeedFree(&mat_G)); 864 CeedCall(CeedFree(&vec_D)); 865 return CEED_ERROR_SUCCESS; 866 } 867 CeedPragmaOptimizeOn; 868 869 /// @} 870 871 /// ---------------------------------------------------------------------------- 872 /// CeedBasis Public API 873 /// ---------------------------------------------------------------------------- 874 /// @addtogroup CeedBasisUser 875 /// @{ 876 877 /** 878 @brief Create a tensor-product basis for H^1 discretizations 879 880 @param[in] ceed Ceed object where the CeedBasis will be created 881 @param[in] dim Topological dimension 882 @param[in] num_comp Number of field components (1 for scalar fields) 883 @param[in] P_1d Number of nodes in one dimension 884 @param[in] Q_1d Number of quadrature points in one dimension 885 @param[in] interp_1d Row-major (Q_1d * P_1d) matrix expressing the values of nodal basis functions at quadrature points 886 @param[in] grad_1d Row-major (Q_1d * P_1d) matrix expressing derivatives of nodal basis functions at quadrature points 887 @param[in] q_ref_1d Array of length Q_1d holding the locations of quadrature points on the 1D reference element [-1, 1] 888 @param[in] q_weight_1d Array of length Q_1d holding the quadrature weights on the reference element 889 @param[out] basis Address of the variable where the newly created CeedBasis will be stored. 890 891 @return An error code: 0 - success, otherwise - failure 892 893 @ref User 894 **/ 895 int CeedBasisCreateTensorH1(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, 896 const CeedScalar *grad_1d, const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis *basis) { 897 if (!ceed->BasisCreateTensorH1) { 898 Ceed delegate; 899 CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis")); 900 901 if (!delegate) { 902 // LCOV_EXCL_START 903 return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisCreateTensorH1"); 904 // LCOV_EXCL_STOP 905 } 906 907 CeedCall(CeedBasisCreateTensorH1(delegate, dim, num_comp, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis)); 908 return CEED_ERROR_SUCCESS; 909 } 910 911 if (dim < 1) { 912 // LCOV_EXCL_START 913 return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis dimension must be a positive value"); 914 // LCOV_EXCL_STOP 915 } 916 917 if (num_comp < 1) { 918 // LCOV_EXCL_START 919 return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component"); 920 // LCOV_EXCL_STOP 921 } 922 923 if (P_1d < 1) { 924 // LCOV_EXCL_START 925 return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node"); 926 // LCOV_EXCL_STOP 927 } 928 929 if (Q_1d < 1) { 930 // LCOV_EXCL_START 931 return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point"); 932 // LCOV_EXCL_STOP 933 } 934 935 CeedElemTopology topo = dim == 1 ? CEED_TOPOLOGY_LINE : dim == 2 ? CEED_TOPOLOGY_QUAD : CEED_TOPOLOGY_HEX; 936 937 CeedCall(CeedCalloc(1, basis)); 938 (*basis)->ceed = ceed; 939 CeedCall(CeedReference(ceed)); 940 (*basis)->ref_count = 1; 941 (*basis)->tensor_basis = 1; 942 (*basis)->dim = dim; 943 (*basis)->topo = topo; 944 (*basis)->num_comp = num_comp; 945 (*basis)->P_1d = P_1d; 946 (*basis)->Q_1d = Q_1d; 947 (*basis)->P = CeedIntPow(P_1d, dim); 948 (*basis)->Q = CeedIntPow(Q_1d, dim); 949 (*basis)->fe_space = CEED_FE_SPACE_H1; 950 CeedCall(CeedCalloc(Q_1d, &(*basis)->q_ref_1d)); 951 CeedCall(CeedCalloc(Q_1d, &(*basis)->q_weight_1d)); 952 if (q_ref_1d) memcpy((*basis)->q_ref_1d, q_ref_1d, Q_1d * sizeof(q_ref_1d[0])); 953 if (q_weight_1d) memcpy((*basis)->q_weight_1d, q_weight_1d, Q_1d * sizeof(q_weight_1d[0])); 954 CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->interp_1d)); 955 CeedCall(CeedCalloc(Q_1d * P_1d, &(*basis)->grad_1d)); 956 if (interp_1d) memcpy((*basis)->interp_1d, interp_1d, Q_1d * P_1d * sizeof(interp_1d[0])); 957 if (grad_1d) memcpy((*basis)->grad_1d, grad_1d, Q_1d * P_1d * sizeof(grad_1d[0])); 958 CeedCall(ceed->BasisCreateTensorH1(dim, P_1d, Q_1d, interp_1d, grad_1d, q_ref_1d, q_weight_1d, *basis)); 959 return CEED_ERROR_SUCCESS; 960 } 961 962 /** 963 @brief Create a tensor-product Lagrange basis 964 965 @param[in] ceed Ceed object where the CeedBasis will be created 966 @param[in] dim Topological dimension of element 967 @param[in] num_comp Number of field components (1 for scalar fields) 968 @param[in] P Number of Gauss-Lobatto nodes in one dimension. 969 The polynomial degree of the resulting Q_k element is k=P-1. 970 @param[in] Q Number of quadrature points in one dimension. 971 @param[in] quad_mode Distribution of the Q quadrature points (affects order of accuracy for the quadrature) 972 @param[out] basis Address of the variable where the newly created CeedBasis will be stored. 973 974 @return An error code: 0 - success, otherwise - failure 975 976 @ref User 977 **/ 978 int CeedBasisCreateTensorH1Lagrange(Ceed ceed, CeedInt dim, CeedInt num_comp, CeedInt P, CeedInt Q, CeedQuadMode quad_mode, CeedBasis *basis) { 979 // Allocate 980 int ierr = CEED_ERROR_SUCCESS, i, j, k; 981 CeedScalar c1, c2, c3, c4, dx, *nodes, *interp_1d, *grad_1d, *q_ref_1d, *q_weight_1d; 982 983 if (dim < 1) { 984 // LCOV_EXCL_START 985 return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis dimension must be a positive value"); 986 // LCOV_EXCL_STOP 987 } 988 989 if (num_comp < 1) { 990 // LCOV_EXCL_START 991 return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component"); 992 // LCOV_EXCL_STOP 993 } 994 995 if (P < 1) { 996 // LCOV_EXCL_START 997 return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node"); 998 // LCOV_EXCL_STOP 999 } 1000 1001 if (Q < 1) { 1002 // LCOV_EXCL_START 1003 return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point"); 1004 // LCOV_EXCL_STOP 1005 } 1006 1007 // Get Nodes and Weights 1008 CeedCall(CeedCalloc(P * Q, &interp_1d)); 1009 CeedCall(CeedCalloc(P * Q, &grad_1d)); 1010 CeedCall(CeedCalloc(P, &nodes)); 1011 CeedCall(CeedCalloc(Q, &q_ref_1d)); 1012 CeedCall(CeedCalloc(Q, &q_weight_1d)); 1013 if (CeedLobattoQuadrature(P, nodes, NULL) != CEED_ERROR_SUCCESS) goto cleanup; 1014 switch (quad_mode) { 1015 case CEED_GAUSS: 1016 ierr = CeedGaussQuadrature(Q, q_ref_1d, q_weight_1d); 1017 break; 1018 case CEED_GAUSS_LOBATTO: 1019 ierr = CeedLobattoQuadrature(Q, q_ref_1d, q_weight_1d); 1020 break; 1021 } 1022 if (ierr != CEED_ERROR_SUCCESS) goto cleanup; 1023 1024 // Build B, D matrix 1025 // Fornberg, 1998 1026 for (i = 0; i < Q; i++) { 1027 c1 = 1.0; 1028 c3 = nodes[0] - q_ref_1d[i]; 1029 interp_1d[i * P + 0] = 1.0; 1030 for (j = 1; j < P; j++) { 1031 c2 = 1.0; 1032 c4 = c3; 1033 c3 = nodes[j] - q_ref_1d[i]; 1034 for (k = 0; k < j; k++) { 1035 dx = nodes[j] - nodes[k]; 1036 c2 *= dx; 1037 if (k == j - 1) { 1038 grad_1d[i * P + j] = c1 * (interp_1d[i * P + k] - c4 * grad_1d[i * P + k]) / c2; 1039 interp_1d[i * P + j] = -c1 * c4 * interp_1d[i * P + k] / c2; 1040 } 1041 grad_1d[i * P + k] = (c3 * grad_1d[i * P + k] - interp_1d[i * P + k]) / dx; 1042 interp_1d[i * P + k] = c3 * interp_1d[i * P + k] / dx; 1043 } 1044 c1 = c2; 1045 } 1046 } 1047 // Pass to CeedBasisCreateTensorH1 1048 CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P, Q, interp_1d, grad_1d, q_ref_1d, q_weight_1d, basis)); 1049 cleanup: 1050 CeedCall(CeedFree(&interp_1d)); 1051 CeedCall(CeedFree(&grad_1d)); 1052 CeedCall(CeedFree(&nodes)); 1053 CeedCall(CeedFree(&q_ref_1d)); 1054 CeedCall(CeedFree(&q_weight_1d)); 1055 return CEED_ERROR_SUCCESS; 1056 } 1057 1058 /** 1059 @brief Create a non tensor-product basis for H^1 discretizations 1060 1061 @param[in] ceed Ceed object where the CeedBasis will be created 1062 @param[in] topo Topology of element, e.g. hypercube, simplex, ect 1063 @param[in] num_comp Number of field components (1 for scalar fields) 1064 @param[in] num_nodes Total number of nodes 1065 @param[in] num_qpts Total number of quadrature points 1066 @param[in] interp Row-major (num_qpts * num_nodes) matrix expressing the values of nodal basis functions at quadrature points 1067 @param[in] grad Row-major (dim * num_qpts * num_nodes) matrix expressing derivatives of nodal basis functions at quadrature points 1068 @param[in] q_ref Array of length num_qpts * dim holding the locations of quadrature points on the reference element 1069 @param[in] q_weight Array of length num_qpts holding the quadrature weights on the reference element 1070 @param[out] basis Address of the variable where the newly created CeedBasis will be stored. 1071 1072 @return An error code: 0 - success, otherwise - failure 1073 1074 @ref User 1075 **/ 1076 int CeedBasisCreateH1(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 1077 const CeedScalar *grad, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) { 1078 CeedInt P = num_nodes, Q = num_qpts, dim = 0; 1079 1080 if (!ceed->BasisCreateH1) { 1081 Ceed delegate; 1082 CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis")); 1083 1084 if (!delegate) { 1085 // LCOV_EXCL_START 1086 return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisCreateH1"); 1087 // LCOV_EXCL_STOP 1088 } 1089 1090 CeedCall(CeedBasisCreateH1(delegate, topo, num_comp, num_nodes, num_qpts, interp, grad, q_ref, q_weight, basis)); 1091 return CEED_ERROR_SUCCESS; 1092 } 1093 1094 if (num_comp < 1) { 1095 // LCOV_EXCL_START 1096 return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component"); 1097 // LCOV_EXCL_STOP 1098 } 1099 1100 if (num_nodes < 1) { 1101 // LCOV_EXCL_START 1102 return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node"); 1103 // LCOV_EXCL_STOP 1104 } 1105 1106 if (num_qpts < 1) { 1107 // LCOV_EXCL_START 1108 return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point"); 1109 // LCOV_EXCL_STOP 1110 } 1111 1112 CeedCall(CeedCalloc(1, basis)); 1113 1114 CeedCall(CeedBasisGetTopologyDimension(topo, &dim)); 1115 1116 (*basis)->ceed = ceed; 1117 CeedCall(CeedReference(ceed)); 1118 (*basis)->ref_count = 1; 1119 (*basis)->tensor_basis = 0; 1120 (*basis)->dim = dim; 1121 (*basis)->topo = topo; 1122 (*basis)->num_comp = num_comp; 1123 (*basis)->P = P; 1124 (*basis)->Q = Q; 1125 (*basis)->fe_space = CEED_FE_SPACE_H1; 1126 CeedCall(CeedCalloc(Q * dim, &(*basis)->q_ref_1d)); 1127 CeedCall(CeedCalloc(Q, &(*basis)->q_weight_1d)); 1128 if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0])); 1129 if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0])); 1130 CeedCall(CeedCalloc(Q * P, &(*basis)->interp)); 1131 CeedCall(CeedCalloc(dim * Q * P, &(*basis)->grad)); 1132 if (interp) memcpy((*basis)->interp, interp, Q * P * sizeof(interp[0])); 1133 if (grad) memcpy((*basis)->grad, grad, dim * Q * P * sizeof(grad[0])); 1134 CeedCall(ceed->BasisCreateH1(topo, dim, P, Q, interp, grad, q_ref, q_weight, *basis)); 1135 return CEED_ERROR_SUCCESS; 1136 } 1137 1138 /** 1139 @brief Create a non tensor-product basis for \f$H(\mathrm{div})\f$ discretizations 1140 1141 @param[in] ceed Ceed object where the CeedBasis will be created 1142 @param[in] topo Topology of element (`CEED_TOPOLOGY_QUAD`, `CEED_TOPOLOGY_PRISM`, etc.), dimension of which is used in some array sizes below 1143 @param[in] num_comp Number of components (usually 1 for vectors in H(div) bases) 1144 @param[in] num_nodes Total number of nodes (dofs per element) 1145 @param[in] num_qpts Total number of quadrature points 1146 @param[in] interp Row-major (dim * num_qpts * num_nodes) matrix expressing the values of basis functions at quadrature points 1147 @param[in] div Row-major (num_qpts * num_nodes) matrix expressing divergence of basis functions at quadrature points 1148 @param[in] q_ref Array of length num_qpts * dim holding the locations of quadrature points on the reference element 1149 @param[in] q_weight Array of length num_qpts holding the quadrature weights on the reference element 1150 @param[out] basis Address of the variable where the newly created CeedBasis will be stored. 1151 1152 @return An error code: 0 - success, otherwise - failure 1153 1154 @ref User 1155 **/ 1156 int CeedBasisCreateHdiv(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 1157 const CeedScalar *div, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) { 1158 CeedInt Q = num_qpts, P = num_nodes, dim = 0; 1159 1160 if (!ceed->BasisCreateHdiv) { 1161 Ceed delegate; 1162 CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis")); 1163 1164 if (!delegate) { 1165 // LCOV_EXCL_START 1166 return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateHdiv"); 1167 // LCOV_EXCL_STOP 1168 } 1169 1170 CeedCall(CeedBasisCreateHdiv(delegate, topo, num_comp, num_nodes, num_qpts, interp, div, q_ref, q_weight, basis)); 1171 return CEED_ERROR_SUCCESS; 1172 } 1173 1174 if (num_comp < 1) { 1175 // LCOV_EXCL_START 1176 return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component"); 1177 // LCOV_EXCL_STOP 1178 } 1179 1180 if (num_nodes < 1) { 1181 // LCOV_EXCL_START 1182 return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node"); 1183 // LCOV_EXCL_STOP 1184 } 1185 1186 if (num_qpts < 1) { 1187 // LCOV_EXCL_START 1188 return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point"); 1189 // LCOV_EXCL_STOP 1190 } 1191 1192 CeedCall(CeedCalloc(1, basis)); 1193 1194 CeedCall(CeedBasisGetTopologyDimension(topo, &dim)); 1195 1196 (*basis)->ceed = ceed; 1197 CeedCall(CeedReference(ceed)); 1198 (*basis)->ref_count = 1; 1199 (*basis)->tensor_basis = 0; 1200 (*basis)->dim = dim; 1201 (*basis)->topo = topo; 1202 (*basis)->num_comp = num_comp; 1203 (*basis)->P = P; 1204 (*basis)->Q = Q; 1205 (*basis)->fe_space = CEED_FE_SPACE_HDIV; 1206 CeedCall(CeedMalloc(Q * dim, &(*basis)->q_ref_1d)); 1207 CeedCall(CeedMalloc(Q, &(*basis)->q_weight_1d)); 1208 if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0])); 1209 if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0])); 1210 CeedCall(CeedMalloc(dim * Q * P, &(*basis)->interp)); 1211 CeedCall(CeedMalloc(Q * P, &(*basis)->div)); 1212 if (interp) memcpy((*basis)->interp, interp, dim * Q * P * sizeof(interp[0])); 1213 if (div) memcpy((*basis)->div, div, Q * P * sizeof(div[0])); 1214 CeedCall(ceed->BasisCreateHdiv(topo, dim, P, Q, interp, div, q_ref, q_weight, *basis)); 1215 return CEED_ERROR_SUCCESS; 1216 } 1217 1218 /** 1219 @brief Create a non tensor-product basis for H(curl) discretizations 1220 1221 @param[in] ceed Ceed object where the CeedBasis will be created 1222 @param[in] topo Topology of element (`CEED_TOPOLOGY_QUAD`, `CEED_TOPOLOGY_PRISM`, etc.), dimension of which is used in some array sizes below 1223 @param[in] num_comp Number of components (usually 1 for vectors in H(curl) bases) 1224 @param[in] num_nodes Total number of nodes (dofs per element) 1225 @param[in] num_qpts Total number of quadrature points 1226 @param[in] interp Row-major (dim * num_qpts * num_nodes) matrix expressing the values of basis functions at quadrature points 1227 @param[in] curl Row-major (curl_comp * num_qpts * num_nodes, curl_comp = 1 if dim < 3 else dim) matrix expressing curl of basis functions at 1228 quadrature points 1229 @param[in] q_ref Array of length num_qpts * dim holding the locations of quadrature points on the reference element 1230 @param[in] q_weight Array of length num_qpts holding the quadrature weights on the reference element 1231 @param[out] basis Address of the variable where the newly created CeedBasis will be stored. 1232 1233 @return An error code: 0 - success, otherwise - failure 1234 1235 @ref User 1236 **/ 1237 int CeedBasisCreateHcurl(Ceed ceed, CeedElemTopology topo, CeedInt num_comp, CeedInt num_nodes, CeedInt num_qpts, const CeedScalar *interp, 1238 const CeedScalar *curl, const CeedScalar *q_ref, const CeedScalar *q_weight, CeedBasis *basis) { 1239 CeedInt Q = num_qpts, P = num_nodes, dim = 0, curl_comp = 0; 1240 1241 if (!ceed->BasisCreateHdiv) { 1242 Ceed delegate; 1243 CeedCall(CeedGetObjectDelegate(ceed, &delegate, "Basis")); 1244 1245 if (!delegate) { 1246 // LCOV_EXCL_START 1247 return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement BasisCreateHcurl"); 1248 // LCOV_EXCL_STOP 1249 } 1250 1251 CeedCall(CeedBasisCreateHcurl(delegate, topo, num_comp, num_nodes, num_qpts, interp, curl, q_ref, q_weight, basis)); 1252 return CEED_ERROR_SUCCESS; 1253 } 1254 1255 if (num_comp < 1) { 1256 // LCOV_EXCL_START 1257 return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 component"); 1258 // LCOV_EXCL_STOP 1259 } 1260 1261 if (num_nodes < 1) { 1262 // LCOV_EXCL_START 1263 return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 node"); 1264 // LCOV_EXCL_STOP 1265 } 1266 1267 if (num_qpts < 1) { 1268 // LCOV_EXCL_START 1269 return CeedError(ceed, CEED_ERROR_DIMENSION, "Basis must have at least 1 quadrature point"); 1270 // LCOV_EXCL_STOP 1271 } 1272 1273 CeedCall(CeedCalloc(1, basis)); 1274 1275 CeedCall(CeedBasisGetTopologyDimension(topo, &dim)); 1276 curl_comp = (dim < 3) ? 1 : dim; 1277 1278 (*basis)->ceed = ceed; 1279 CeedCall(CeedReference(ceed)); 1280 (*basis)->ref_count = 1; 1281 (*basis)->tensor_basis = 0; 1282 (*basis)->dim = dim; 1283 (*basis)->topo = topo; 1284 (*basis)->num_comp = num_comp; 1285 (*basis)->P = P; 1286 (*basis)->Q = Q; 1287 (*basis)->fe_space = CEED_FE_SPACE_HCURL; 1288 CeedCall(CeedMalloc(Q * dim, &(*basis)->q_ref_1d)); 1289 CeedCall(CeedMalloc(Q, &(*basis)->q_weight_1d)); 1290 if (q_ref) memcpy((*basis)->q_ref_1d, q_ref, Q * dim * sizeof(q_ref[0])); 1291 if (q_weight) memcpy((*basis)->q_weight_1d, q_weight, Q * sizeof(q_weight[0])); 1292 CeedCall(CeedMalloc(dim * Q * P, &(*basis)->interp)); 1293 CeedCall(CeedMalloc(curl_comp * Q * P, &(*basis)->curl)); 1294 if (interp) memcpy((*basis)->interp, interp, dim * Q * P * sizeof(interp[0])); 1295 if (curl) memcpy((*basis)->curl, curl, curl_comp * Q * P * sizeof(curl[0])); 1296 CeedCall(ceed->BasisCreateHcurl(topo, dim, P, Q, interp, curl, q_ref, q_weight, *basis)); 1297 return CEED_ERROR_SUCCESS; 1298 } 1299 1300 /** 1301 @brief Create a CeedBasis for projection from the nodes of `basis_from` to the nodes of `basis_to`. 1302 1303 Only `CEED_EVAL_INTERP` and `CEED_EVAL_GRAD` will be valid for the new basis, `basis_project`. 1304 The interpolation is given by `interp_project = interp_to^+ * interp_from`, where the pesudoinverse `interp_to^+` is given by QR 1305 factorization. The gradient is given by `grad_project = interp_to^+ * grad_from`. Note: `basis_from` and `basis_to` must have compatible quadrature 1306 spaces. 1307 Note: `basis_project` will have the same number of components as `basis_from`, regardless of the number of components that `basis_to` has. If 1308 `basis_from` has 3 components and `basis_to` has 5 components, then `basis_project` will have 3 components. 1309 1310 @param[in] basis_from CeedBasis to prolong from 1311 @param[in] basis_to CeedBasis to prolong to 1312 @param[out] basis_project Address of the variable where the newly created CeedBasis will be stored. 1313 1314 @return An error code: 0 - success, otherwise - failure 1315 1316 @ref User 1317 **/ 1318 int CeedBasisCreateProjection(CeedBasis basis_from, CeedBasis basis_to, CeedBasis *basis_project) { 1319 Ceed ceed; 1320 CeedCall(CeedBasisGetCeed(basis_to, &ceed)); 1321 1322 // Create projection matrix 1323 CeedScalar *interp_project, *grad_project; 1324 CeedCall(CeedBasisCreateProjectionMatrices(basis_from, basis_to, &interp_project, &grad_project)); 1325 1326 // Build basis 1327 bool is_tensor; 1328 CeedInt dim, num_comp; 1329 CeedScalar *q_ref, *q_weight; 1330 CeedCall(CeedBasisIsTensor(basis_to, &is_tensor)); 1331 CeedCall(CeedBasisGetDimension(basis_to, &dim)); 1332 CeedCall(CeedBasisGetNumComponents(basis_from, &num_comp)); 1333 if (is_tensor) { 1334 CeedInt P_1d_to, P_1d_from; 1335 CeedCall(CeedBasisGetNumNodes1D(basis_from, &P_1d_from)); 1336 CeedCall(CeedBasisGetNumNodes1D(basis_to, &P_1d_to)); 1337 CeedCall(CeedCalloc(P_1d_to, &q_ref)); 1338 CeedCall(CeedCalloc(P_1d_to, &q_weight)); 1339 CeedCall(CeedBasisCreateTensorH1(ceed, dim, num_comp, P_1d_from, P_1d_to, interp_project, grad_project, q_ref, q_weight, basis_project)); 1340 } else { 1341 CeedElemTopology topo; 1342 CeedCall(CeedBasisGetTopology(basis_to, &topo)); 1343 CeedInt num_nodes_to, num_nodes_from; 1344 CeedCall(CeedBasisGetNumNodes(basis_from, &num_nodes_from)); 1345 CeedCall(CeedBasisGetNumNodes(basis_to, &num_nodes_to)); 1346 CeedCall(CeedCalloc(num_nodes_to * dim, &q_ref)); 1347 CeedCall(CeedCalloc(num_nodes_to, &q_weight)); 1348 CeedCall(CeedBasisCreateH1(ceed, topo, num_comp, num_nodes_from, num_nodes_to, interp_project, grad_project, q_ref, q_weight, basis_project)); 1349 } 1350 1351 // Cleanup 1352 CeedCall(CeedFree(&interp_project)); 1353 CeedCall(CeedFree(&grad_project)); 1354 CeedCall(CeedFree(&q_ref)); 1355 CeedCall(CeedFree(&q_weight)); 1356 1357 return CEED_ERROR_SUCCESS; 1358 } 1359 1360 /** 1361 @brief Copy the pointer to a CeedBasis. 1362 1363 Note: If the value of `basis_copy` passed into this function is non-NULL, then it is assumed that `basis_copy` is a pointer to a CeedBasis. 1364 This CeedBasis will be destroyed if `basis_copy` is the only reference to this CeedBasis. 1365 1366 @param[in] basis CeedBasis to copy reference to 1367 @param[in,out] basis_copy Variable to store copied reference 1368 1369 @return An error code: 0 - success, otherwise - failure 1370 1371 @ref User 1372 **/ 1373 int CeedBasisReferenceCopy(CeedBasis basis, CeedBasis *basis_copy) { 1374 if (basis != CEED_BASIS_COLLOCATED) CeedCall(CeedBasisReference(basis)); 1375 CeedCall(CeedBasisDestroy(basis_copy)); 1376 *basis_copy = basis; 1377 return CEED_ERROR_SUCCESS; 1378 } 1379 1380 /** 1381 @brief View a CeedBasis 1382 1383 @param[in] basis CeedBasis to view 1384 @param[in] stream Stream to view to, e.g., stdout 1385 1386 @return An error code: 0 - success, otherwise - failure 1387 1388 @ref User 1389 **/ 1390 int CeedBasisView(CeedBasis basis, FILE *stream) { 1391 CeedElemTopology topo = basis->topo; 1392 CeedFESpace fe_space = basis->fe_space; 1393 CeedInt q_comp = 0; 1394 1395 // Print FE space and element topology of the basis 1396 if (basis->tensor_basis) { 1397 fprintf(stream, "CeedBasis (%s on a %s element): dim=%" CeedInt_FMT " P=%" CeedInt_FMT " Q=%" CeedInt_FMT "\n", CeedFESpaces[fe_space], 1398 CeedElemTopologies[topo], basis->dim, basis->P_1d, basis->Q_1d); 1399 } else { 1400 fprintf(stream, "CeedBasis (%s on a %s element): dim=%" CeedInt_FMT " P=%" CeedInt_FMT " Q=%" CeedInt_FMT "\n", CeedFESpaces[fe_space], 1401 CeedElemTopologies[topo], basis->dim, basis->P, basis->Q); 1402 } 1403 // Print quadrature data, interpolation/gradient/divergence/curl of the basis 1404 if (basis->tensor_basis) { // tensor basis 1405 CeedCall(CeedScalarView("qref1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_ref_1d, stream)); 1406 CeedCall(CeedScalarView("qweight1d", "\t% 12.8f", 1, basis->Q_1d, basis->q_weight_1d, stream)); 1407 CeedCall(CeedScalarView("interp1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, basis->interp_1d, stream)); 1408 CeedCall(CeedScalarView("grad1d", "\t% 12.8f", basis->Q_1d, basis->P_1d, basis->grad_1d, stream)); 1409 } else { // non-tensor basis 1410 CeedCall(CeedScalarView("qref", "\t% 12.8f", 1, basis->Q * basis->dim, basis->q_ref_1d, stream)); 1411 CeedCall(CeedScalarView("qweight", "\t% 12.8f", 1, basis->Q, basis->q_weight_1d, stream)); 1412 CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp)); 1413 CeedCall(CeedScalarView("interp", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->interp, stream)); 1414 if (basis->grad) { 1415 CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp)); 1416 CeedCall(CeedScalarView("grad", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->grad, stream)); 1417 } 1418 if (basis->div) { 1419 CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_DIV, &q_comp)); 1420 CeedCall(CeedScalarView("div", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->div, stream)); 1421 } 1422 if (basis->curl) { 1423 CeedCall(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_CURL, &q_comp)); 1424 CeedCall(CeedScalarView("curl", "\t% 12.8f", q_comp * basis->Q, basis->P, basis->curl, stream)); 1425 } 1426 } 1427 return CEED_ERROR_SUCCESS; 1428 } 1429 1430 /** 1431 @brief Apply basis evaluation from nodes to quadrature points or vice versa 1432 1433 @param[in] basis CeedBasis to evaluate 1434 @param[in] num_elem The number of elements to apply the basis evaluation to; 1435 the backend will specify the ordering in CeedElemRestrictionCreateBlocked() 1436 @param[in] t_mode \ref CEED_NOTRANSPOSE to evaluate from nodes to quadrature points; 1437 \ref CEED_TRANSPOSE to apply the transpose, mapping from quadrature points to nodes 1438 @param[in] eval_mode \ref CEED_EVAL_NONE to use values directly, 1439 \ref CEED_EVAL_INTERP to use interpolated values, 1440 \ref CEED_EVAL_GRAD to use gradients, 1441 \ref CEED_EVAL_DIV to use divergence, 1442 \ref CEED_EVAL_CURL to use curl, 1443 \ref CEED_EVAL_WEIGHT to use quadrature weights. 1444 @param[in] u Input CeedVector 1445 @param[out] v Output CeedVector 1446 1447 @return An error code: 0 - success, otherwise - failure 1448 1449 @ref User 1450 **/ 1451 int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, CeedVector v) { 1452 CeedSize u_length = 0, v_length; 1453 CeedInt dim, num_comp, q_comp, num_nodes, num_qpts; 1454 CeedCall(CeedBasisGetDimension(basis, &dim)); 1455 CeedCall(CeedBasisGetNumComponents(basis, &num_comp)); 1456 CeedCall(CeedBasisGetNumQuadratureComponents(basis, eval_mode, &q_comp)); 1457 CeedCall(CeedBasisGetNumNodes(basis, &num_nodes)); 1458 CeedCall(CeedBasisGetNumQuadraturePoints(basis, &num_qpts)); 1459 CeedCall(CeedVectorGetLength(v, &v_length)); 1460 if (u) { 1461 CeedCall(CeedVectorGetLength(u, &u_length)); 1462 } 1463 1464 if (!basis->Apply) { 1465 // LCOV_EXCL_START 1466 return CeedError(basis->ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support BasisApply"); 1467 // LCOV_EXCL_STOP 1468 } 1469 1470 // Check compatibility of topological and geometrical dimensions 1471 if ((t_mode == CEED_TRANSPOSE && (v_length % num_nodes != 0 || u_length % num_qpts != 0)) || 1472 (t_mode == CEED_NOTRANSPOSE && (u_length % num_nodes != 0 || v_length % num_qpts != 0))) { 1473 // LCOV_EXCL_START 1474 return CeedError(basis->ceed, CEED_ERROR_DIMENSION, "Length of input/output vectors incompatible with basis dimensions"); 1475 // LCOV_EXCL_STOP 1476 } 1477 1478 // Check vector lengths to prevent out of bounds issues 1479 bool bad_dims = false; 1480 switch (eval_mode) { 1481 case CEED_EVAL_NONE: 1482 case CEED_EVAL_INTERP: 1483 case CEED_EVAL_GRAD: 1484 case CEED_EVAL_DIV: 1485 case CEED_EVAL_CURL: 1486 bad_dims = ((t_mode == CEED_TRANSPOSE && (u_length < num_elem * num_comp * num_qpts * q_comp || v_length < num_elem * num_comp * num_nodes)) || 1487 (t_mode == CEED_NOTRANSPOSE && (v_length < num_elem * num_qpts * num_comp * q_comp || u_length < num_elem * num_comp * num_nodes))); 1488 break; 1489 case CEED_EVAL_WEIGHT: 1490 bad_dims = v_length < num_elem * num_qpts; 1491 break; 1492 } 1493 if (bad_dims) { 1494 // LCOV_EXCL_START 1495 return CeedError(basis->ceed, CEED_ERROR_DIMENSION, "Input/output vectors too short for basis and evaluation mode"); 1496 // LCOV_EXCL_STOP 1497 } 1498 1499 CeedCall(basis->Apply(basis, num_elem, t_mode, eval_mode, u, v)); 1500 return CEED_ERROR_SUCCESS; 1501 } 1502 1503 /** 1504 @brief Get Ceed associated with a CeedBasis 1505 1506 @param[in] basis CeedBasis 1507 @param[out] ceed Variable to store Ceed 1508 1509 @return An error code: 0 - success, otherwise - failure 1510 1511 @ref Advanced 1512 **/ 1513 int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed) { 1514 *ceed = basis->ceed; 1515 return CEED_ERROR_SUCCESS; 1516 } 1517 1518 /** 1519 @brief Get dimension for given CeedBasis 1520 1521 @param[in] basis CeedBasis 1522 @param[out] dim Variable to store dimension of basis 1523 1524 @return An error code: 0 - success, otherwise - failure 1525 1526 @ref Advanced 1527 **/ 1528 int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim) { 1529 *dim = basis->dim; 1530 return CEED_ERROR_SUCCESS; 1531 } 1532 1533 /** 1534 @brief Get topology for given CeedBasis 1535 1536 @param[in] basis CeedBasis 1537 @param[out] topo Variable to store topology of basis 1538 1539 @return An error code: 0 - success, otherwise - failure 1540 1541 @ref Advanced 1542 **/ 1543 int CeedBasisGetTopology(CeedBasis basis, CeedElemTopology *topo) { 1544 *topo = basis->topo; 1545 return CEED_ERROR_SUCCESS; 1546 } 1547 1548 /** 1549 @brief Get number of components for given CeedBasis 1550 1551 @param[in] basis CeedBasis 1552 @param[out] num_comp Variable to store number of components of basis 1553 1554 @return An error code: 0 - success, otherwise - failure 1555 1556 @ref Advanced 1557 **/ 1558 int CeedBasisGetNumComponents(CeedBasis basis, CeedInt *num_comp) { 1559 *num_comp = basis->num_comp; 1560 return CEED_ERROR_SUCCESS; 1561 } 1562 1563 /** 1564 @brief Get total number of nodes (in dim dimensions) of a CeedBasis 1565 1566 @param[in] basis CeedBasis 1567 @param[out] P Variable to store number of nodes 1568 1569 @return An error code: 0 - success, otherwise - failure 1570 1571 @ref Utility 1572 **/ 1573 int CeedBasisGetNumNodes(CeedBasis basis, CeedInt *P) { 1574 *P = basis->P; 1575 return CEED_ERROR_SUCCESS; 1576 } 1577 1578 /** 1579 @brief Get total number of nodes (in 1 dimension) of a CeedBasis 1580 1581 @param[in] basis CeedBasis 1582 @param[out] P_1d Variable to store number of nodes 1583 1584 @return An error code: 0 - success, otherwise - failure 1585 1586 @ref Advanced 1587 **/ 1588 int CeedBasisGetNumNodes1D(CeedBasis basis, CeedInt *P_1d) { 1589 if (!basis->tensor_basis) { 1590 // LCOV_EXCL_START 1591 return CeedError(basis->ceed, CEED_ERROR_MINOR, "Cannot supply P_1d for non-tensor basis"); 1592 // LCOV_EXCL_STOP 1593 } 1594 1595 *P_1d = basis->P_1d; 1596 return CEED_ERROR_SUCCESS; 1597 } 1598 1599 /** 1600 @brief Get total number of quadrature points (in dim dimensions) of a CeedBasis 1601 1602 @param[in] basis CeedBasis 1603 @param[out] Q Variable to store number of quadrature points 1604 1605 @return An error code: 0 - success, otherwise - failure 1606 1607 @ref Utility 1608 **/ 1609 int CeedBasisGetNumQuadraturePoints(CeedBasis basis, CeedInt *Q) { 1610 *Q = basis->Q; 1611 return CEED_ERROR_SUCCESS; 1612 } 1613 1614 /** 1615 @brief Get total number of quadrature points (in 1 dimension) of a CeedBasis 1616 1617 @param[in] basis CeedBasis 1618 @param[out] Q_1d Variable to store number of quadrature points 1619 1620 @return An error code: 0 - success, otherwise - failure 1621 1622 @ref Advanced 1623 **/ 1624 int CeedBasisGetNumQuadraturePoints1D(CeedBasis basis, CeedInt *Q_1d) { 1625 if (!basis->tensor_basis) { 1626 // LCOV_EXCL_START 1627 return CeedError(basis->ceed, CEED_ERROR_MINOR, "Cannot supply Q_1d for non-tensor basis"); 1628 // LCOV_EXCL_STOP 1629 } 1630 1631 *Q_1d = basis->Q_1d; 1632 return CEED_ERROR_SUCCESS; 1633 } 1634 1635 /** 1636 @brief Get reference coordinates of quadrature points (in dim dimensions) of a CeedBasis 1637 1638 @param[in] basis CeedBasis 1639 @param[out] q_ref Variable to store reference coordinates of quadrature points 1640 1641 @return An error code: 0 - success, otherwise - failure 1642 1643 @ref Advanced 1644 **/ 1645 int CeedBasisGetQRef(CeedBasis basis, const CeedScalar **q_ref) { 1646 *q_ref = basis->q_ref_1d; 1647 return CEED_ERROR_SUCCESS; 1648 } 1649 1650 /** 1651 @brief Get quadrature weights of quadrature points (in dim dimensions) of a CeedBasis 1652 1653 @param[in] basis CeedBasis 1654 @param[out] q_weight Variable to store quadrature weights 1655 1656 @return An error code: 0 - success, otherwise - failure 1657 1658 @ref Advanced 1659 **/ 1660 int CeedBasisGetQWeights(CeedBasis basis, const CeedScalar **q_weight) { 1661 *q_weight = basis->q_weight_1d; 1662 return CEED_ERROR_SUCCESS; 1663 } 1664 1665 /** 1666 @brief Get interpolation matrix of a CeedBasis 1667 1668 @param[in] basis CeedBasis 1669 @param[out] interp Variable to store interpolation matrix 1670 1671 @return An error code: 0 - success, otherwise - failure 1672 1673 @ref Advanced 1674 **/ 1675 int CeedBasisGetInterp(CeedBasis basis, const CeedScalar **interp) { 1676 if (!basis->interp && basis->tensor_basis) { 1677 // Allocate 1678 CeedCall(CeedMalloc(basis->Q * basis->P, &basis->interp)); 1679 1680 // Initialize 1681 for (CeedInt i = 0; i < basis->Q * basis->P; i++) basis->interp[i] = 1.0; 1682 1683 // Calculate 1684 for (CeedInt d = 0; d < basis->dim; d++) { 1685 for (CeedInt qpt = 0; qpt < basis->Q; qpt++) { 1686 for (CeedInt node = 0; node < basis->P; node++) { 1687 CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d; 1688 CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d; 1689 basis->interp[qpt * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p]; 1690 } 1691 } 1692 } 1693 } 1694 *interp = basis->interp; 1695 return CEED_ERROR_SUCCESS; 1696 } 1697 1698 /** 1699 @brief Get 1D interpolation matrix of a tensor product CeedBasis 1700 1701 @param[in] basis CeedBasis 1702 @param[out] interp_1d Variable to store interpolation matrix 1703 1704 @return An error code: 0 - success, otherwise - failure 1705 1706 @ref Backend 1707 **/ 1708 int CeedBasisGetInterp1D(CeedBasis basis, const CeedScalar **interp_1d) { 1709 if (!basis->tensor_basis) { 1710 // LCOV_EXCL_START 1711 return CeedError(basis->ceed, CEED_ERROR_MINOR, "CeedBasis is not a tensor product basis."); 1712 // LCOV_EXCL_STOP 1713 } 1714 1715 *interp_1d = basis->interp_1d; 1716 return CEED_ERROR_SUCCESS; 1717 } 1718 1719 /** 1720 @brief Get gradient matrix of a CeedBasis 1721 1722 @param[in] basis CeedBasis 1723 @param[out] grad Variable to store gradient matrix 1724 1725 @return An error code: 0 - success, otherwise - failure 1726 1727 @ref Advanced 1728 **/ 1729 int CeedBasisGetGrad(CeedBasis basis, const CeedScalar **grad) { 1730 if (!basis->grad && basis->tensor_basis) { 1731 // Allocate 1732 CeedCall(CeedMalloc(basis->dim * basis->Q * basis->P, &basis->grad)); 1733 1734 // Initialize 1735 for (CeedInt i = 0; i < basis->dim * basis->Q * basis->P; i++) basis->grad[i] = 1.0; 1736 1737 // Calculate 1738 for (CeedInt d = 0; d < basis->dim; d++) { 1739 for (CeedInt i = 0; i < basis->dim; i++) { 1740 for (CeedInt qpt = 0; qpt < basis->Q; qpt++) { 1741 for (CeedInt node = 0; node < basis->P; node++) { 1742 CeedInt p = (node / CeedIntPow(basis->P_1d, d)) % basis->P_1d; 1743 CeedInt q = (qpt / CeedIntPow(basis->Q_1d, d)) % basis->Q_1d; 1744 if (i == d) basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->grad_1d[q * basis->P_1d + p]; 1745 else basis->grad[(i * basis->Q + qpt) * (basis->P) + node] *= basis->interp_1d[q * basis->P_1d + p]; 1746 } 1747 } 1748 } 1749 } 1750 } 1751 *grad = basis->grad; 1752 return CEED_ERROR_SUCCESS; 1753 } 1754 1755 /** 1756 @brief Get 1D gradient matrix of a tensor product CeedBasis 1757 1758 @param[in] basis CeedBasis 1759 @param[out] grad_1d Variable to store gradient matrix 1760 1761 @return An error code: 0 - success, otherwise - failure 1762 1763 @ref Advanced 1764 **/ 1765 int CeedBasisGetGrad1D(CeedBasis basis, const CeedScalar **grad_1d) { 1766 if (!basis->tensor_basis) { 1767 // LCOV_EXCL_START 1768 return CeedError(basis->ceed, CEED_ERROR_MINOR, "CeedBasis is not a tensor product basis."); 1769 // LCOV_EXCL_STOP 1770 } 1771 1772 *grad_1d = basis->grad_1d; 1773 return CEED_ERROR_SUCCESS; 1774 } 1775 1776 /** 1777 @brief Get divergence matrix of a CeedBasis 1778 1779 @param[in] basis CeedBasis 1780 @param[out] div Variable to store divergence matrix 1781 1782 @return An error code: 0 - success, otherwise - failure 1783 1784 @ref Advanced 1785 **/ 1786 int CeedBasisGetDiv(CeedBasis basis, const CeedScalar **div) { 1787 if (!basis->div) { 1788 // LCOV_EXCL_START 1789 return CeedError(basis->ceed, CEED_ERROR_MINOR, "CeedBasis does not have divergence matrix."); 1790 // LCOV_EXCL_STOP 1791 } 1792 1793 *div = basis->div; 1794 return CEED_ERROR_SUCCESS; 1795 } 1796 1797 /** 1798 @brief Get curl matrix of a CeedBasis 1799 1800 @param[in] basis CeedBasis 1801 @param[out] curl Variable to store curl matrix 1802 1803 @return An error code: 0 - success, otherwise - failure 1804 1805 @ref Advanced 1806 **/ 1807 int CeedBasisGetCurl(CeedBasis basis, const CeedScalar **curl) { 1808 if (!basis->curl) { 1809 // LCOV_EXCL_START 1810 return CeedError(basis->ceed, CEED_ERROR_MINOR, "CeedBasis does not have curl matrix."); 1811 // LCOV_EXCL_STOP 1812 } 1813 1814 *curl = basis->curl; 1815 return CEED_ERROR_SUCCESS; 1816 } 1817 1818 /** 1819 @brief Destroy a CeedBasis 1820 1821 @param[in,out] basis CeedBasis to destroy 1822 1823 @return An error code: 0 - success, otherwise - failure 1824 1825 @ref User 1826 **/ 1827 int CeedBasisDestroy(CeedBasis *basis) { 1828 if (!*basis || *basis == CEED_BASIS_COLLOCATED || --(*basis)->ref_count > 0) { 1829 *basis = NULL; 1830 return CEED_ERROR_SUCCESS; 1831 } 1832 if ((*basis)->Destroy) CeedCall((*basis)->Destroy(*basis)); 1833 if ((*basis)->contract) CeedCall(CeedTensorContractDestroy(&(*basis)->contract)); 1834 CeedCall(CeedFree(&(*basis)->q_ref_1d)); 1835 CeedCall(CeedFree(&(*basis)->q_weight_1d)); 1836 CeedCall(CeedFree(&(*basis)->interp)); 1837 CeedCall(CeedFree(&(*basis)->interp_1d)); 1838 CeedCall(CeedFree(&(*basis)->grad)); 1839 CeedCall(CeedFree(&(*basis)->grad_1d)); 1840 CeedCall(CeedFree(&(*basis)->div)); 1841 CeedCall(CeedFree(&(*basis)->curl)); 1842 CeedCall(CeedDestroy(&(*basis)->ceed)); 1843 CeedCall(CeedFree(basis)); 1844 return CEED_ERROR_SUCCESS; 1845 } 1846 1847 /** 1848 @brief Construct a Gauss-Legendre quadrature 1849 1850 @param[in] Q Number of quadrature points (integrates polynomials of degree 2*Q-1 exactly) 1851 @param[out] q_ref_1d Array of length Q to hold the abscissa on [-1, 1] 1852 @param[out] q_weight_1d Array of length Q to hold the weights 1853 1854 @return An error code: 0 - success, otherwise - failure 1855 1856 @ref Utility 1857 **/ 1858 int CeedGaussQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) { 1859 // Allocate 1860 CeedScalar P0, P1, P2, dP2, xi, wi, PI = 4.0 * atan(1.0); 1861 // Build q_ref_1d, q_weight_1d 1862 for (CeedInt i = 0; i <= Q / 2; i++) { 1863 // Guess 1864 xi = cos(PI * (CeedScalar)(2 * i + 1) / ((CeedScalar)(2 * Q))); 1865 // Pn(xi) 1866 P0 = 1.0; 1867 P1 = xi; 1868 P2 = 0.0; 1869 for (CeedInt j = 2; j <= Q; j++) { 1870 P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j)); 1871 P0 = P1; 1872 P1 = P2; 1873 } 1874 // First Newton Step 1875 dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0); 1876 xi = xi - P2 / dP2; 1877 // Newton to convergence 1878 for (CeedInt k = 0; k < 100 && fabs(P2) > 10 * CEED_EPSILON; k++) { 1879 P0 = 1.0; 1880 P1 = xi; 1881 for (CeedInt j = 2; j <= Q; j++) { 1882 P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j)); 1883 P0 = P1; 1884 P1 = P2; 1885 } 1886 dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0); 1887 xi = xi - P2 / dP2; 1888 } 1889 // Save xi, wi 1890 wi = 2.0 / ((1.0 - xi * xi) * dP2 * dP2); 1891 q_weight_1d[i] = wi; 1892 q_weight_1d[Q - 1 - i] = wi; 1893 q_ref_1d[i] = -xi; 1894 q_ref_1d[Q - 1 - i] = xi; 1895 } 1896 return CEED_ERROR_SUCCESS; 1897 } 1898 1899 /** 1900 @brief Construct a Gauss-Legendre-Lobatto quadrature 1901 1902 @param[in] Q Number of quadrature points (integrates polynomials of degree 2*Q-3 exactly) 1903 @param[out] q_ref_1d Array of length Q to hold the abscissa on [-1, 1] 1904 @param[out] q_weight_1d Array of length Q to hold the weights 1905 1906 @return An error code: 0 - success, otherwise - failure 1907 1908 @ref Utility 1909 **/ 1910 int CeedLobattoQuadrature(CeedInt Q, CeedScalar *q_ref_1d, CeedScalar *q_weight_1d) { 1911 // Allocate 1912 CeedScalar P0, P1, P2, dP2, d2P2, xi, wi, PI = 4.0 * atan(1.0); 1913 // Build q_ref_1d, q_weight_1d 1914 // Set endpoints 1915 if (Q < 2) { 1916 // LCOV_EXCL_START 1917 return CeedError(NULL, CEED_ERROR_DIMENSION, "Cannot create Lobatto quadrature with Q=%" CeedInt_FMT " < 2 points", Q); 1918 // LCOV_EXCL_STOP 1919 } 1920 wi = 2.0 / ((CeedScalar)(Q * (Q - 1))); 1921 if (q_weight_1d) { 1922 q_weight_1d[0] = wi; 1923 q_weight_1d[Q - 1] = wi; 1924 } 1925 q_ref_1d[0] = -1.0; 1926 q_ref_1d[Q - 1] = 1.0; 1927 // Interior 1928 for (CeedInt i = 1; i <= (Q - 1) / 2; i++) { 1929 // Guess 1930 xi = cos(PI * (CeedScalar)(i) / (CeedScalar)(Q - 1)); 1931 // Pn(xi) 1932 P0 = 1.0; 1933 P1 = xi; 1934 P2 = 0.0; 1935 for (CeedInt j = 2; j < Q; j++) { 1936 P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j)); 1937 P0 = P1; 1938 P1 = P2; 1939 } 1940 // First Newton step 1941 dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0); 1942 d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi); 1943 xi = xi - dP2 / d2P2; 1944 // Newton to convergence 1945 for (CeedInt k = 0; k < 100 && fabs(dP2) > 10 * CEED_EPSILON; k++) { 1946 P0 = 1.0; 1947 P1 = xi; 1948 for (CeedInt j = 2; j < Q; j++) { 1949 P2 = (((CeedScalar)(2 * j - 1)) * xi * P1 - ((CeedScalar)(j - 1)) * P0) / ((CeedScalar)(j)); 1950 P0 = P1; 1951 P1 = P2; 1952 } 1953 dP2 = (xi * P2 - P0) * (CeedScalar)Q / (xi * xi - 1.0); 1954 d2P2 = (2 * xi * dP2 - (CeedScalar)(Q * (Q - 1)) * P2) / (1.0 - xi * xi); 1955 xi = xi - dP2 / d2P2; 1956 } 1957 // Save xi, wi 1958 wi = 2.0 / (((CeedScalar)(Q * (Q - 1))) * P2 * P2); 1959 if (q_weight_1d) { 1960 q_weight_1d[i] = wi; 1961 q_weight_1d[Q - 1 - i] = wi; 1962 } 1963 q_ref_1d[i] = -xi; 1964 q_ref_1d[Q - 1 - i] = xi; 1965 } 1966 return CEED_ERROR_SUCCESS; 1967 } 1968 1969 /// @} 1970