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