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