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