1 // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 2 // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 3 // reserved. See files LICENSE and NOTICE for details. 4 // 5 // This file is part of CEED, a collection of benchmarks, miniapps, software 6 // libraries and APIs for efficient high-order finite element and spectral 7 // element discretizations for exascale applications. For more information and 8 // source code availability see http://github.com/ceed. 9 // 10 // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 11 // a collaborative effort of two U.S. Department of Energy organizations (Office 12 // of Science and the National Nuclear Security Administration) responsible for 13 // the planning and preparation of a capable exascale ecosystem, including 14 // software, applications, hardware, advanced system engineering and early 15 // testbed platforms, in support of the nation's exascale computing imperative. 16 17 #include <ceed/ceed.h> 18 #include <ceed/backend.h> 19 #include <ceed-impl.h> 20 #include <limits.h> 21 #include <stdbool.h> 22 #include <stdio.h> 23 #include <string.h> 24 25 /// @file 26 /// Implementation of public CeedQFunction interfaces 27 28 /// @cond DOXYGEN_SKIP 29 static struct CeedQFunction_private ceed_qfunction_none; 30 /// @endcond 31 32 /// @addtogroup CeedQFunctionUser 33 /// @{ 34 35 // Indicate that no QFunction is provided by the user 36 const CeedQFunction CEED_QFUNCTION_NONE = &ceed_qfunction_none; 37 38 /// @} 39 40 /// @cond DOXYGEN_SKIP 41 static struct { 42 char name[CEED_MAX_RESOURCE_LEN]; 43 char source[CEED_MAX_RESOURCE_LEN]; 44 CeedInt vec_length; 45 CeedQFunctionUser f; 46 int (*init)(Ceed ceed, const char *name, CeedQFunction qf); 47 } gallery_qfunctions[1024]; 48 static size_t num_qfunctions; 49 /// @endcond 50 51 /// ---------------------------------------------------------------------------- 52 /// CeedQFunction Library Internal Functions 53 /// ---------------------------------------------------------------------------- 54 /// @addtogroup CeedQFunctionDeveloper 55 /// @{ 56 57 /** 58 @brief Register a gallery QFunction 59 60 @param name Name for this backend to respond to 61 @param source Absolute path to source of QFunction, 62 "\path\CEED_DIR\gallery\folder\file.h:function_name" 63 @param vec_length Vector length. Caller must ensure that number of quadrature 64 points is a multiple of vec_length. 65 @param f Function pointer to evaluate action at quadrature points. 66 See \ref CeedQFunctionUser. 67 @param init Initialization function called by CeedQFunctionInit() when the 68 QFunction is selected. 69 70 @return An error code: 0 - success, otherwise - failure 71 72 @ref Developer 73 **/ 74 int CeedQFunctionRegister(const char *name, const char *source, 75 CeedInt vec_length, CeedQFunctionUser f, 76 int (*init)(Ceed, const char *, CeedQFunction)) { 77 if (num_qfunctions >= sizeof(gallery_qfunctions) / sizeof( 78 gallery_qfunctions[0])) 79 // LCOV_EXCL_START 80 return CeedError(NULL, CEED_ERROR_MAJOR, "Too many gallery QFunctions"); 81 // LCOV_EXCL_STOP 82 83 strncpy(gallery_qfunctions[num_qfunctions].name, name, CEED_MAX_RESOURCE_LEN); 84 gallery_qfunctions[num_qfunctions].name[CEED_MAX_RESOURCE_LEN-1] = 0; 85 strncpy(gallery_qfunctions[num_qfunctions].source, source, 86 CEED_MAX_RESOURCE_LEN); 87 gallery_qfunctions[num_qfunctions].source[CEED_MAX_RESOURCE_LEN-1] = 0; 88 gallery_qfunctions[num_qfunctions].vec_length = vec_length; 89 gallery_qfunctions[num_qfunctions].f = f; 90 gallery_qfunctions[num_qfunctions].init = init; 91 num_qfunctions++; 92 return CEED_ERROR_SUCCESS; 93 } 94 95 /** 96 @brief Set a CeedQFunction field, used by CeedQFunctionAddInput/Output 97 98 @param f CeedQFunctionField 99 @param field_name Name of QFunction field 100 @param size Size of QFunction field, (num_comp * dim) for @ref CEED_EVAL_GRAD or 101 (num_comp * 1) for @ref CEED_EVAL_NONE, @ref CEED_EVAL_INTERP, and @ref CEED_EVAL_WEIGHT 102 @param eval_mode \ref CEED_EVAL_NONE to use values directly, 103 \ref CEED_EVAL_INTERP to use interpolated values, 104 \ref CEED_EVAL_GRAD to use gradients, 105 \ref CEED_EVAL_WEIGHT to use quadrature weights. 106 107 @return An error code: 0 - success, otherwise - failure 108 109 @ref Developer 110 **/ 111 static int CeedQFunctionFieldSet(CeedQFunctionField *f,const char *field_name, 112 CeedInt size, CeedEvalMode eval_mode) { 113 size_t len = strlen(field_name); 114 char *tmp; 115 int ierr; 116 117 ierr = CeedCalloc(1, f); CeedChk(ierr); 118 ierr = CeedCalloc(len+1, &tmp); CeedChk(ierr); 119 memcpy(tmp, field_name, len+1); 120 (*f)->field_name = tmp; 121 (*f)->size = size; 122 (*f)->eval_mode = eval_mode; 123 return CEED_ERROR_SUCCESS; 124 } 125 126 /** 127 @brief View a field of a CeedQFunction 128 129 @param[in] field QFunction field to view 130 @param[in] field_number Number of field being viewed 131 @param[in] in true for input field, false for output 132 @param[in] stream Stream to view to, e.g., stdout 133 134 @return An error code: 0 - success, otherwise - failure 135 136 @ref Utility 137 **/ 138 static int CeedQFunctionFieldView(CeedQFunctionField field, 139 CeedInt field_number, 140 bool in, FILE *stream) { 141 int ierr; 142 const char *inout = in ? "Input" : "Output"; 143 char *field_name; 144 ierr = CeedQFunctionFieldGetName(field, &field_name); CeedChk(ierr); 145 CeedInt size; 146 ierr = CeedQFunctionFieldGetSize(field, &size); CeedChk(ierr); 147 CeedEvalMode eval_mode; 148 ierr = CeedQFunctionFieldGetEvalMode(field, &eval_mode); CeedChk(ierr); 149 fprintf(stream, " %s Field [%d]:\n" 150 " Name: \"%s\"\n" 151 " Size: %d\n" 152 " EvalMode: \"%s\"\n", 153 inout, field_number, field_name, size, CeedEvalModes[eval_mode]); 154 return CEED_ERROR_SUCCESS; 155 } 156 157 /** 158 @brief Set flag to determine if Fortran interface is used 159 160 @param qf CeedQFunction 161 @param status Boolean value to set as Fortran status 162 163 @return An error code: 0 - success, otherwise - failure 164 165 @ref Backend 166 **/ 167 int CeedQFunctionSetFortranStatus(CeedQFunction qf, bool status) { 168 qf->fortran_status = status; 169 return CEED_ERROR_SUCCESS; 170 } 171 172 /// @} 173 174 /// ---------------------------------------------------------------------------- 175 /// CeedQFunction Backend API 176 /// ---------------------------------------------------------------------------- 177 /// @addtogroup CeedQFunctionBackend 178 /// @{ 179 180 /** 181 @brief Get the Ceed associated with a CeedQFunction 182 183 @param qf CeedQFunction 184 @param[out] ceed Variable to store Ceed 185 186 @return An error code: 0 - success, otherwise - failure 187 188 @ref Backend 189 **/ 190 int CeedQFunctionGetCeed(CeedQFunction qf, Ceed *ceed) { 191 *ceed = qf->ceed; 192 return CEED_ERROR_SUCCESS; 193 } 194 195 /** 196 @brief Get the vector length of a CeedQFunction 197 198 @param qf CeedQFunction 199 @param[out] vec_length Variable to store vector length 200 201 @return An error code: 0 - success, otherwise - failure 202 203 @ref Backend 204 **/ 205 int CeedQFunctionGetVectorLength(CeedQFunction qf, CeedInt *vec_length) { 206 *vec_length = qf->vec_length; 207 return CEED_ERROR_SUCCESS; 208 } 209 210 /** 211 @brief Get the number of inputs and outputs to a CeedQFunction 212 213 @param qf CeedQFunction 214 @param[out] num_input Variable to store number of input fields 215 @param[out] num_output Variable to store number of output fields 216 217 @return An error code: 0 - success, otherwise - failure 218 219 @ref Backend 220 **/ 221 int CeedQFunctionGetNumArgs(CeedQFunction qf, CeedInt *num_input, 222 CeedInt *num_output) { 223 if (num_input) *num_input = qf->num_input_fields; 224 if (num_output) *num_output = qf->num_output_fields; 225 return CEED_ERROR_SUCCESS; 226 } 227 228 /** 229 @brief Get the source path string for a CeedQFunction 230 231 @param qf CeedQFunction 232 @param[out] source Variable to store source path string 233 234 @return An error code: 0 - success, otherwise - failure 235 236 @ref Backend 237 **/ 238 int CeedQFunctionGetSourcePath(CeedQFunction qf, char **source) { 239 *source = (char *) qf->source_path; 240 return CEED_ERROR_SUCCESS; 241 } 242 243 /** 244 @brief Get the User Function for a CeedQFunction 245 246 @param qf CeedQFunction 247 @param[out] f Variable to store user function 248 249 @return An error code: 0 - success, otherwise - failure 250 251 @ref Backend 252 **/ 253 int CeedQFunctionGetUserFunction(CeedQFunction qf, CeedQFunctionUser *f) { 254 *f = qf->function; 255 return CEED_ERROR_SUCCESS; 256 } 257 258 /** 259 @brief Get global context for a CeedQFunction. 260 Note: For QFunctions from the Fortran interface, this 261 function will return the Fortran context 262 CeedQFunctionContext. 263 264 @param qf CeedQFunction 265 @param[out] ctx Variable to store CeedQFunctionContext 266 267 @return An error code: 0 - success, otherwise - failure 268 269 @ref Backend 270 **/ 271 int CeedQFunctionGetContext(CeedQFunction qf, CeedQFunctionContext *ctx) { 272 *ctx = qf->ctx; 273 return CEED_ERROR_SUCCESS; 274 } 275 276 /** 277 @brief Get true user context for a CeedQFunction 278 Note: For all QFunctions this function will return the user 279 CeedQFunctionContext and not interface context 280 CeedQFunctionContext, if any such object exists. 281 282 @param qf CeedQFunction 283 @param[out] ctx Variable to store CeedQFunctionContext 284 285 @return An error code: 0 - success, otherwise - failure 286 @ref Backend 287 **/ 288 int CeedQFunctionGetInnerContext(CeedQFunction qf, CeedQFunctionContext *ctx) { 289 int ierr; 290 if (qf->fortran_status) { 291 CeedFortranContext fortran_ctx = NULL; 292 ierr = CeedQFunctionContextGetData(qf->ctx, CEED_MEM_HOST, &fortran_ctx); 293 CeedChk(ierr); 294 *ctx = fortran_ctx->innerctx; 295 ierr = CeedQFunctionContextRestoreData(qf->ctx, (void *)&fortran_ctx); 296 CeedChk(ierr); 297 } else { 298 *ctx = qf->ctx; 299 } 300 return CEED_ERROR_SUCCESS; 301 } 302 303 /** 304 @brief Determine if QFunction is identity 305 306 @param qf CeedQFunction 307 @param[out] is_identity Variable to store identity status 308 309 @return An error code: 0 - success, otherwise - failure 310 311 @ref Backend 312 **/ 313 int CeedQFunctionIsIdentity(CeedQFunction qf, bool *is_identity) { 314 *is_identity = qf->identity; 315 return CEED_ERROR_SUCCESS; 316 } 317 318 /** 319 @brief Get backend data of a CeedQFunction 320 321 @param qf CeedQFunction 322 @param[out] data Variable to store data 323 324 @return An error code: 0 - success, otherwise - failure 325 326 @ref Backend 327 **/ 328 int CeedQFunctionGetData(CeedQFunction qf, void *data) { 329 *(void **)data = qf->data; 330 return CEED_ERROR_SUCCESS; 331 } 332 333 /** 334 @brief Set backend data of a CeedQFunction 335 336 @param[out] qf CeedQFunction 337 @param data Data to set 338 339 @return An error code: 0 - success, otherwise - failure 340 341 @ref Backend 342 **/ 343 int CeedQFunctionSetData(CeedQFunction qf, void *data) { 344 qf->data = data; 345 return CEED_ERROR_SUCCESS; 346 } 347 348 /** 349 @brief Increment the reference counter for a CeedQFunction 350 351 @param qf CeedQFunction to increment the reference counter 352 353 @return An error code: 0 - success, otherwise - failure 354 355 @ref Backend 356 **/ 357 int CeedQFunctionReference(CeedQFunction qf) { 358 qf->ref_count++; 359 return CEED_ERROR_SUCCESS; 360 } 361 362 /// @} 363 364 /// ---------------------------------------------------------------------------- 365 /// CeedQFunction Public API 366 /// ---------------------------------------------------------------------------- 367 /// @addtogroup CeedQFunctionUser 368 /// @{ 369 370 /** 371 @brief Create a CeedQFunction for evaluating interior (volumetric) terms. 372 373 @param ceed A Ceed object where the CeedQFunction will be created 374 @param vec_length Vector length. Caller must ensure that number of quadrature 375 points is a multiple of vec_length. 376 @param f Function pointer to evaluate action at quadrature points. 377 See \ref CeedQFunctionUser. 378 @param source Absolute path to source of QFunction, 379 "\abs_path\file.h:function_name". 380 For support across all backends, this source must only 381 contain constructs supported by C99, C++11, and CUDA. 382 @param[out] qf Address of the variable where the newly created 383 CeedQFunction will be stored 384 385 @return An error code: 0 - success, otherwise - failure 386 387 See \ref CeedQFunctionUser for details on the call-back function @a f's 388 arguments. 389 390 @ref User 391 **/ 392 int CeedQFunctionCreateInterior(Ceed ceed, CeedInt vec_length, 393 CeedQFunctionUser f, 394 const char *source, CeedQFunction *qf) { 395 int ierr; 396 char *source_copy; 397 398 if (!ceed->QFunctionCreate) { 399 Ceed delegate; 400 ierr = CeedGetObjectDelegate(ceed, &delegate, "QFunction"); CeedChk(ierr); 401 402 if (!delegate) 403 // LCOV_EXCL_START 404 return CeedError(ceed, CEED_ERROR_UNSUPPORTED, 405 "Backend does not support QFunctionCreate"); 406 // LCOV_EXCL_STOP 407 408 ierr = CeedQFunctionCreateInterior(delegate, vec_length, f, source, qf); 409 CeedChk(ierr); 410 return CEED_ERROR_SUCCESS; 411 } 412 413 ierr = CeedCalloc(1, qf); CeedChk(ierr); 414 (*qf)->ceed = ceed; 415 ierr = CeedReference(ceed); CeedChk(ierr); 416 (*qf)->ref_count = 1; 417 (*qf)->vec_length = vec_length; 418 (*qf)->identity = 0; 419 (*qf)->function = f; 420 size_t slen = strlen(source) + 1; 421 ierr = CeedMalloc(slen, &source_copy); CeedChk(ierr); 422 memcpy(source_copy, source, slen); 423 (*qf)->source_path = source_copy; 424 ierr = CeedCalloc(16, &(*qf)->input_fields); CeedChk(ierr); 425 ierr = CeedCalloc(16, &(*qf)->output_fields); CeedChk(ierr); 426 ierr = ceed->QFunctionCreate(*qf); CeedChk(ierr); 427 return CEED_ERROR_SUCCESS; 428 } 429 430 /** 431 @brief Create a CeedQFunction for evaluating interior (volumetric) terms by name. 432 433 @param ceed A Ceed object where the CeedQFunction will be created 434 @param name Name of QFunction to use from gallery 435 @param[out] qf Address of the variable where the newly created 436 CeedQFunction will be stored 437 438 @return An error code: 0 - success, otherwise - failure 439 440 @ref User 441 **/ 442 int CeedQFunctionCreateInteriorByName(Ceed ceed, const char *name, 443 CeedQFunction *qf) { 444 int ierr; 445 size_t match_len = 0, match_idx = UINT_MAX; 446 char *name_copy; 447 448 ierr = CeedQFunctionRegisterAll(); CeedChk(ierr); 449 // Find matching backend 450 if (!name) return CeedError(ceed, CEED_ERROR_INCOMPLETE, 451 "No QFunction name provided"); 452 for (size_t i=0; i<num_qfunctions; i++) { 453 size_t n; 454 const char *curr_name = gallery_qfunctions[i].name; 455 for (n = 0; curr_name[n] && curr_name[n] == name[n]; n++) {} 456 if (n > match_len) { 457 match_len = n; 458 match_idx = i; 459 } 460 } 461 if (!match_len) 462 // LCOV_EXCL_START 463 return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "No suitable gallery QFunction"); 464 // LCOV_EXCL_STOP 465 466 // Create QFunction 467 ierr = CeedQFunctionCreateInterior(ceed, 468 gallery_qfunctions[match_idx].vec_length, 469 gallery_qfunctions[match_idx].f, 470 gallery_qfunctions[match_idx].source, qf); 471 CeedChk(ierr); 472 473 // QFunction specific setup 474 ierr = gallery_qfunctions[match_idx].init(ceed, name, *qf); CeedChk(ierr); 475 476 // Copy name 477 size_t slen = strlen(name) + 1; 478 ierr = CeedMalloc(slen, &name_copy); CeedChk(ierr); 479 memcpy(name_copy, name, slen); 480 (*qf)->qf_name = name_copy; 481 return CEED_ERROR_SUCCESS; 482 } 483 484 /** 485 @brief Create an identity CeedQFunction. Inputs are written into outputs in 486 the order given. This is useful for CeedOperators that can be 487 represented with only the action of a CeedRestriction and CeedBasis, 488 such as restriction and prolongation operators for p-multigrid. 489 Backends may optimize CeedOperators with this CeedQFunction to avoid 490 the copy of input data to output fields by using the same memory 491 location for both. 492 493 @param ceed A Ceed object where the CeedQFunction will be created 494 @param[in] size Size of the QFunction fields 495 @param[in] in_mode CeedEvalMode for input to CeedQFunction 496 @param[in] out_mode CeedEvalMode for output to CeedQFunction 497 @param[out] qf Address of the variable where the newly created 498 CeedQFunction will be stored 499 500 @return An error code: 0 - success, otherwise - failure 501 502 @ref User 503 **/ 504 int CeedQFunctionCreateIdentity(Ceed ceed, CeedInt size, CeedEvalMode in_mode, 505 CeedEvalMode out_mode, CeedQFunction *qf) { 506 int ierr; 507 508 ierr = CeedQFunctionCreateInteriorByName(ceed, "Identity", qf); CeedChk(ierr); 509 ierr = CeedQFunctionAddInput(*qf, "input", size, in_mode); CeedChk(ierr); 510 ierr = CeedQFunctionAddOutput(*qf, "output", size, out_mode); CeedChk(ierr); 511 512 (*qf)->identity = 1; 513 CeedInt *size_data; 514 ierr = CeedCalloc(1, &size_data); CeedChk(ierr); 515 size_data[0] = size; 516 CeedQFunctionContext ctx; 517 ierr = CeedQFunctionContextCreate(ceed, &ctx); CeedChk(ierr); 518 ierr = CeedQFunctionContextSetData(ctx, CEED_MEM_HOST, CEED_OWN_POINTER, 519 sizeof(*size_data), (void *)size_data); 520 CeedChk(ierr); 521 ierr = CeedQFunctionSetContext(*qf, ctx); CeedChk(ierr); 522 ierr = CeedQFunctionContextDestroy(&ctx); CeedChk(ierr); 523 return CEED_ERROR_SUCCESS; 524 } 525 526 /** 527 @brief Copy the pointer to a CeedQFunction. Both pointers should 528 be destroyed with `CeedQFunctionDestroy()`; 529 Note: If `*qf_copy` is non-NULL, then it is assumed that 530 `*qf_copy` is a pointer to a CeedQFunction. This 531 CeedQFunction will be destroyed if `*qf_copy` is the only 532 reference to this CeedQFunction. 533 534 @param qf CeedQFunction to copy reference to 535 @param[out] qf_copy Variable to store copied reference 536 537 @return An error code: 0 - success, otherwise - failure 538 539 @ref User 540 **/ 541 int CeedQFunctionReferenceCopy(CeedQFunction qf, CeedQFunction *qf_copy) { 542 int ierr; 543 544 ierr = CeedQFunctionReference(qf); CeedChk(ierr); 545 ierr = CeedQFunctionDestroy(qf_copy); CeedChk(ierr); 546 *qf_copy = qf; 547 return CEED_ERROR_SUCCESS; 548 } 549 550 /** 551 @brief Add a CeedQFunction input 552 553 @param qf CeedQFunction 554 @param field_name Name of QFunction field 555 @param size Size of QFunction field, (num_comp * dim) for @ref CEED_EVAL_GRAD or 556 (num_comp * 1) for @ref CEED_EVAL_NONE and @ref CEED_EVAL_INTERP 557 @param eval_mode \ref CEED_EVAL_NONE to use values directly, 558 \ref CEED_EVAL_INTERP to use interpolated values, 559 \ref CEED_EVAL_GRAD to use gradients. 560 561 @return An error code: 0 - success, otherwise - failure 562 563 @ref User 564 **/ 565 int CeedQFunctionAddInput(CeedQFunction qf, const char *field_name, 566 CeedInt size, 567 CeedEvalMode eval_mode) { 568 if (qf->operators_set) 569 // LCOV_EXCL_START 570 return CeedError(qf->ceed, CEED_ERROR_MAJOR, 571 "QFunction cannot be changed when in use by an operator"); 572 // LCOV_EXCL_STOP 573 if ((eval_mode == CEED_EVAL_WEIGHT) && (size != 1)) 574 // LCOV_EXCL_START 575 return CeedError(qf->ceed, CEED_ERROR_DIMENSION, 576 "CEED_EVAL_WEIGHT should have size 1"); 577 // LCOV_EXCL_STOP 578 int ierr = CeedQFunctionFieldSet(&qf->input_fields[qf->num_input_fields], 579 field_name, size, eval_mode); 580 CeedChk(ierr); 581 qf->num_input_fields++; 582 return CEED_ERROR_SUCCESS; 583 } 584 585 /** 586 @brief Add a CeedQFunction output 587 588 @param qf CeedQFunction 589 @param field_name Name of QFunction field 590 @param size Size of QFunction field, (num_comp * dim) for @ref CEED_EVAL_GRAD or 591 (num_comp * 1) for @ref CEED_EVAL_NONE and @ref CEED_EVAL_INTERP 592 @param eval_mode \ref CEED_EVAL_NONE to use values directly, 593 \ref CEED_EVAL_INTERP to use interpolated values, 594 \ref CEED_EVAL_GRAD to use gradients. 595 596 @return An error code: 0 - success, otherwise - failure 597 598 @ref User 599 **/ 600 int CeedQFunctionAddOutput(CeedQFunction qf, const char *field_name, 601 CeedInt size, CeedEvalMode eval_mode) { 602 if (qf->operators_set) 603 // LCOV_EXCL_START 604 return CeedError(qf->ceed, CEED_ERROR_MAJOR, 605 "QFunction cannot be changed when in use by an operator"); 606 // LCOV_EXCL_STOP 607 if (eval_mode == CEED_EVAL_WEIGHT) 608 // LCOV_EXCL_START 609 return CeedError(qf->ceed, CEED_ERROR_DIMENSION, 610 "Cannot create QFunction output with " 611 "CEED_EVAL_WEIGHT"); 612 // LCOV_EXCL_STOP 613 int ierr = CeedQFunctionFieldSet(&qf->output_fields[qf->num_output_fields], 614 field_name, size, eval_mode); 615 CeedChk(ierr); 616 qf->num_output_fields++; 617 return CEED_ERROR_SUCCESS; 618 } 619 620 /** 621 @brief Get the CeedQFunctionFields of a CeedQFunction 622 623 @param qf CeedQFunction 624 @param[out] input_fields Variable to store input_fields 625 @param[out] output_fields Variable to store output_fields 626 627 @return An error code: 0 - success, otherwise - failure 628 629 @ref Backend 630 **/ 631 int CeedQFunctionGetFields(CeedQFunction qf, CeedInt *num_input_fields, 632 CeedQFunctionField **input_fields, 633 CeedInt *num_output_fields, 634 CeedQFunctionField **output_fields) { 635 if (num_input_fields) *num_input_fields = qf->num_input_fields; 636 if (input_fields) *input_fields = qf->input_fields; 637 if (num_output_fields) *num_output_fields = qf->num_output_fields; 638 if (output_fields) *output_fields = qf->output_fields; 639 return CEED_ERROR_SUCCESS; 640 } 641 642 /** 643 @brief Get the name of a CeedQFunctionField 644 645 @param qf_field CeedQFunctionField 646 @param[out] field_name Variable to store the field name 647 648 @return An error code: 0 - success, otherwise - failure 649 650 @ref Backend 651 **/ 652 int CeedQFunctionFieldGetName(CeedQFunctionField qf_field, char **field_name) { 653 *field_name = (char *)qf_field->field_name; 654 return CEED_ERROR_SUCCESS; 655 } 656 657 /** 658 @brief Get the number of components of a CeedQFunctionField 659 660 @param qf_field CeedQFunctionField 661 @param[out] size Variable to store the size of the field 662 663 @return An error code: 0 - success, otherwise - failure 664 665 @ref Backend 666 **/ 667 int CeedQFunctionFieldGetSize(CeedQFunctionField qf_field, CeedInt *size) { 668 *size = qf_field->size; 669 return CEED_ERROR_SUCCESS; 670 } 671 672 /** 673 @brief Get the CeedEvalMode of a CeedQFunctionField 674 675 @param qf_field CeedQFunctionField 676 @param[out] eval_mode Variable to store the field evaluation mode 677 678 @return An error code: 0 - success, otherwise - failure 679 680 @ref Backend 681 **/ 682 int CeedQFunctionFieldGetEvalMode(CeedQFunctionField qf_field, 683 CeedEvalMode *eval_mode) { 684 *eval_mode = qf_field->eval_mode; 685 return CEED_ERROR_SUCCESS; 686 } 687 688 /** 689 @brief Set global context for a CeedQFunction 690 691 @param qf CeedQFunction 692 @param ctx Context data to set 693 694 @return An error code: 0 - success, otherwise - failure 695 696 @ref User 697 **/ 698 int CeedQFunctionSetContext(CeedQFunction qf, CeedQFunctionContext ctx) { 699 int ierr; 700 qf->ctx = ctx; 701 ierr = CeedQFunctionContextReference(ctx); CeedChk(ierr); 702 return CEED_ERROR_SUCCESS; 703 } 704 705 /** 706 @brief View a CeedQFunction 707 708 @param[in] qf CeedQFunction to view 709 @param[in] stream Stream to write; typically stdout/stderr or a file 710 711 @return Error code: 0 - success, otherwise - failure 712 713 @ref User 714 **/ 715 int CeedQFunctionView(CeedQFunction qf, FILE *stream) { 716 int ierr; 717 718 fprintf(stream, "%sCeedQFunction %s\n", 719 qf->qf_name ? "Gallery " : "User ", qf->qf_name ? qf->qf_name : ""); 720 721 fprintf(stream, " %d Input Field%s:\n", qf->num_input_fields, 722 qf->num_input_fields>1 ? "s" : ""); 723 for (CeedInt i=0; i<qf->num_input_fields; i++) { 724 ierr = CeedQFunctionFieldView(qf->input_fields[i], i, 1, stream); 725 CeedChk(ierr); 726 } 727 728 fprintf(stream, " %d Output Field%s:\n", qf->num_output_fields, 729 qf->num_output_fields>1 ? "s" : ""); 730 for (CeedInt i=0; i<qf->num_output_fields; i++) { 731 ierr = CeedQFunctionFieldView(qf->output_fields[i], i, 0, stream); 732 CeedChk(ierr); 733 } 734 return CEED_ERROR_SUCCESS; 735 } 736 737 /** 738 @brief Apply the action of a CeedQFunction 739 740 @param qf CeedQFunction 741 @param Q Number of quadrature points 742 @param[in] u Array of input CeedVectors 743 @param[out] v Array of output CeedVectors 744 745 @return An error code: 0 - success, otherwise - failure 746 747 @ref User 748 **/ 749 int CeedQFunctionApply(CeedQFunction qf, CeedInt Q, 750 CeedVector *u, CeedVector *v) { 751 int ierr; 752 if (!qf->Apply) 753 // LCOV_EXCL_START 754 return CeedError(qf->ceed, CEED_ERROR_UNSUPPORTED, 755 "Backend does not support QFunctionApply"); 756 // LCOV_EXCL_STOP 757 if (Q % qf->vec_length) 758 // LCOV_EXCL_START 759 return CeedError(qf->ceed, CEED_ERROR_DIMENSION, 760 "Number of quadrature points %d must be a " 761 "multiple of %d", Q, qf->vec_length); 762 // LCOV_EXCL_STOP 763 ierr = qf->Apply(qf, Q, u, v); CeedChk(ierr); 764 return CEED_ERROR_SUCCESS; 765 } 766 767 /** 768 @brief Destroy a CeedQFunction 769 770 @param qf CeedQFunction to destroy 771 772 @return An error code: 0 - success, otherwise - failure 773 774 @ref User 775 **/ 776 int CeedQFunctionDestroy(CeedQFunction *qf) { 777 int ierr; 778 779 if (!*qf || --(*qf)->ref_count > 0) return CEED_ERROR_SUCCESS; 780 // Backend destroy 781 if ((*qf)->Destroy) { 782 ierr = (*qf)->Destroy(*qf); CeedChk(ierr); 783 } 784 // Free fields 785 for (int i=0; i<(*qf)->num_input_fields; i++) { 786 ierr = CeedFree(&(*(*qf)->input_fields[i]).field_name); CeedChk(ierr); 787 ierr = CeedFree(&(*qf)->input_fields[i]); CeedChk(ierr); 788 } 789 for (int i=0; i<(*qf)->num_output_fields; i++) { 790 ierr = CeedFree(&(*(*qf)->output_fields[i]).field_name); CeedChk(ierr); 791 ierr = CeedFree(&(*qf)->output_fields[i]); CeedChk(ierr); 792 } 793 ierr = CeedFree(&(*qf)->input_fields); CeedChk(ierr); 794 ierr = CeedFree(&(*qf)->output_fields); CeedChk(ierr); 795 796 // User context data object 797 ierr = CeedQFunctionContextDestroy(&(*qf)->ctx); CeedChk(ierr); 798 799 ierr = CeedFree(&(*qf)->source_path); CeedChk(ierr); 800 ierr = CeedFree(&(*qf)->qf_name); CeedChk(ierr); 801 ierr = CeedDestroy(&(*qf)->ceed); CeedChk(ierr); 802 ierr = CeedFree(qf); CeedChk(ierr); 803 return CEED_ERROR_SUCCESS; 804 } 805 806 /// @} 807