1 // Copyright (c) 2017-2025, 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 <ceed/jit-tools.h> 12 #include <limits.h> 13 #include <stdbool.h> 14 #include <stdio.h> 15 #include <string.h> 16 17 /// @file 18 /// Implementation of public CeedQFunction interfaces 19 20 /// @cond DOXYGEN_SKIP 21 static struct CeedQFunction_private ceed_qfunction_none; 22 /// @endcond 23 24 /// @addtogroup CeedQFunctionUser 25 /// @{ 26 27 // Indicate that no `CeedQFunction` is provided by the user 28 const CeedQFunction CEED_QFUNCTION_NONE = &ceed_qfunction_none; 29 30 /// @} 31 32 /// @cond DOXYGEN_SKIP 33 static struct { 34 char name[CEED_MAX_RESOURCE_LEN]; 35 char source[CEED_MAX_RESOURCE_LEN]; 36 CeedInt vec_length; 37 CeedQFunctionUser f; 38 int (*init)(Ceed ceed, const char *name, CeedQFunction qf); 39 } gallery_qfunctions[1024]; 40 static size_t num_qfunctions; 41 /// @endcond 42 43 /// ---------------------------------------------------------------------------- 44 /// CeedQFunction Library Internal Functions 45 /// ---------------------------------------------------------------------------- 46 /// @addtogroup CeedQFunctionDeveloper 47 /// @{ 48 49 /** 50 @brief Register a gallery `CeedQFunction` 51 52 @param[in] name Name for this backend to respond to 53 @param[in] source Absolute path to source of `CeedQFunction`, "\path\CEED_DIR\gallery\folder\file.h:function_name" 54 @param[in] vec_length Vector length. 55 Caller must ensure that number of quadrature points is a multiple of `vec_length`. 56 @param[in] f Function pointer to evaluate action at quadrature points. 57 See `CeedQFunctionUser`. 58 @param[in] init Initialization function called by @ref CeedQFunctionCreateInteriorByName() when the `CeedQFunction` is selected. 59 60 @return An error code: 0 - success, otherwise - failure 61 62 @ref Developer 63 **/ 64 int CeedQFunctionRegister(const char *name, const char *source, CeedInt vec_length, CeedQFunctionUser f, 65 int (*init)(Ceed, const char *, CeedQFunction)) { 66 const char *relative_file_path; 67 int ierr = 0; 68 69 CeedDebugEnv("Gallery Register: %s", name); 70 CeedCall(CeedGetJitRelativePath(source, &relative_file_path)); 71 CeedPragmaCritical(CeedQFunctionRegister) { 72 if (num_qfunctions < sizeof(gallery_qfunctions) / sizeof(gallery_qfunctions[0])) { 73 strncpy(gallery_qfunctions[num_qfunctions].name, name, CEED_MAX_RESOURCE_LEN); 74 gallery_qfunctions[num_qfunctions].name[CEED_MAX_RESOURCE_LEN - 1] = 0; 75 strncpy(gallery_qfunctions[num_qfunctions].source, relative_file_path, CEED_MAX_RESOURCE_LEN); 76 gallery_qfunctions[num_qfunctions].source[CEED_MAX_RESOURCE_LEN - 1] = 0; 77 gallery_qfunctions[num_qfunctions].vec_length = vec_length; 78 gallery_qfunctions[num_qfunctions].f = f; 79 gallery_qfunctions[num_qfunctions].init = init; 80 num_qfunctions++; 81 } else { 82 ierr = 1; 83 } 84 } 85 CeedCheck(ierr == 0, NULL, CEED_ERROR_MAJOR, "Too many gallery CeedQFunctions"); 86 return CEED_ERROR_SUCCESS; 87 } 88 89 /** 90 @brief Set a `CeedQFunction` field, used by @ref CeedQFunctionAddInput() and @ref CeedQFunctionAddOutput() 91 92 @param[out] f `CeedQFunctionField` 93 @param[in] field_name Name of `CeedQFunction` field 94 @param[in] size Size of `CeedQFunction` field, (`num_comp * 1`) for @ref CEED_EVAL_NONE and @ref CEED_EVAL_WEIGHT, (`num_comp * 1`) for @ref CEED_EVAL_INTERP for an \f$H^1\f$ space or (`num_comp * dim`) for an \f$H(\mathrm{div})\f$ or \f$H(\mathrm{curl})\f$ space, (`num_comp * dim`) for @ref CEED_EVAL_GRAD, or (num_comp * 1) for @ref CEED_EVAL_DIV, and (`num_comp * curl_dim`) with `curl_dim = 1` if `dim < 3` and `curl_dim = dim` for @ref CEED_EVAL_CURL. 95 @param[in] eval_mode @ref CEED_EVAL_NONE to use values directly, 96 @ref CEED_EVAL_WEIGHT to use quadrature weights, 97 @ref CEED_EVAL_INTERP to use interpolated values, 98 @ref CEED_EVAL_GRAD to use gradients, 99 @ref CEED_EVAL_DIV to use divergence, 100 @ref CEED_EVAL_CURL to use curl 101 102 @return An error code: 0 - success, otherwise - failure 103 104 @ref Developer 105 **/ 106 static int CeedQFunctionFieldSet(CeedQFunctionField *f, const char *field_name, CeedInt size, CeedEvalMode eval_mode) { 107 CeedCall(CeedCalloc(1, f)); 108 CeedCall(CeedStringAllocCopy(field_name, (char **)&(*f)->field_name)); 109 (*f)->size = size; 110 (*f)->eval_mode = eval_mode; 111 return CEED_ERROR_SUCCESS; 112 } 113 114 /** 115 @brief View a field of a `CeedQFunction` 116 117 @param[in] field `CeedQFunction` field to view 118 @param[in] field_number Number of field being viewed 119 @param[in] in true for input field, false for output 120 @param[in] tabs Tabs to append before each new line 121 @param[in] stream Stream to view to, e.g., `stdout` 122 123 @return An error code: 0 - success, otherwise - failure 124 125 @ref Utility 126 **/ 127 static int CeedQFunctionFieldView(CeedQFunctionField field, CeedInt field_number, bool in, const char *tabs, FILE *stream) { 128 const char *inout = in ? "Input" : "Output"; 129 const char *field_name; 130 CeedInt size; 131 CeedEvalMode eval_mode; 132 133 CeedCall(CeedQFunctionFieldGetData(field, &field_name, &size, &eval_mode)); 134 fprintf(stream, 135 "%s %s field %" CeedInt_FMT 136 ":\n%s" 137 " Name: \"%s\"\n%s" 138 " Size: %" CeedInt_FMT 139 "\n%s" 140 " EvalMode: \"%s\"\n", 141 tabs, inout, field_number, tabs, field_name, tabs, size, tabs, CeedEvalModes[eval_mode]); 142 return CEED_ERROR_SUCCESS; 143 } 144 145 /** 146 @brief Set flag to determine if Fortran interface is used 147 148 @param[in,out] qf CeedQFunction 149 @param[in] status Boolean value to set as Fortran status 150 151 @return An error code: 0 - success, otherwise - failure 152 153 @ref Backend 154 **/ 155 int CeedQFunctionSetFortranStatus(CeedQFunction qf, bool status) { 156 qf->is_fortran = status; 157 return CEED_ERROR_SUCCESS; 158 } 159 160 /// @} 161 162 /// ---------------------------------------------------------------------------- 163 /// CeedQFunction Backend API 164 /// ---------------------------------------------------------------------------- 165 /// @addtogroup CeedQFunctionBackend 166 /// @{ 167 168 /** 169 @brief Get the vector length of a `CeedQFunction` 170 171 @param[in] qf `CeedQFunction` 172 @param[out] vec_length Variable to store vector length 173 174 @return An error code: 0 - success, otherwise - failure 175 176 @ref Backend 177 **/ 178 int CeedQFunctionGetVectorLength(CeedQFunction qf, CeedInt *vec_length) { 179 *vec_length = qf->vec_length; 180 return CEED_ERROR_SUCCESS; 181 } 182 183 /** 184 @brief Get the number of inputs and outputs to a `CeedQFunction` 185 186 @param[in] qf `CeedQFunction` 187 @param[out] num_input Variable to store number of input fields 188 @param[out] num_output Variable to store number of output fields 189 190 @return An error code: 0 - success, otherwise - failure 191 192 @ref Backend 193 **/ 194 int CeedQFunctionGetNumArgs(CeedQFunction qf, CeedInt *num_input, CeedInt *num_output) { 195 if (num_input) *num_input = qf->num_input_fields; 196 if (num_output) *num_output = qf->num_output_fields; 197 return CEED_ERROR_SUCCESS; 198 } 199 200 /** 201 @brief Get the name of the `CeedQFunction`. 202 Use the `name` if created via @ref `CeedQFunctionCreateInteriorByName`, otherwise return the kernel name via @ref `CeedQFunctionGetKernelName`. 203 204 @param[in] qf `CeedQFunction` 205 @param[out] kernel_name Variable to store `CeedQFunction` name 206 207 @return An error code: 0 - success, otherwise - failure 208 209 @ref Backend 210 **/ 211 int CeedQFunctionGetName(CeedQFunction qf, const char **name) { 212 if (qf->is_gallery) { 213 *name = qf->gallery_name; 214 } else { 215 CeedCall(CeedQFunctionGetKernelName(qf, name)); 216 } 217 return CEED_ERROR_SUCCESS; 218 } 219 220 /** 221 @brief Get the name of the user function for a `CeedQFunction` 222 223 @param[in] qf `CeedQFunction` 224 @param[out] kernel_name Variable to store string holding kernel name 225 226 @return An error code: 0 - success, otherwise - failure 227 228 @ref Backend 229 **/ 230 int CeedQFunctionGetKernelName(CeedQFunction qf, const char **kernel_name) { 231 if (!qf->kernel_name) { 232 char *kernel_name_copy; 233 234 if (qf->user_source) { 235 const char *kernel_name = strrchr(qf->user_source, ':') + 1; 236 size_t kernel_name_len = strlen(kernel_name); 237 238 CeedCall(CeedCalloc(kernel_name_len + 1, &kernel_name_copy)); 239 memcpy(kernel_name_copy, kernel_name, kernel_name_len); 240 } else { 241 CeedCall(CeedCalloc(1, &kernel_name_copy)); 242 } 243 qf->kernel_name = kernel_name_copy; 244 } 245 246 *kernel_name = qf->kernel_name; 247 return CEED_ERROR_SUCCESS; 248 } 249 250 /** 251 @brief Get the source path string for a `CeedQFunction` 252 253 @param[in] qf `CeedQFunction` 254 @param[out] source_path Variable to store source path string 255 256 @return An error code: 0 - success, otherwise - failure 257 258 @ref Backend 259 **/ 260 int CeedQFunctionGetSourcePath(CeedQFunction qf, const char **source_path) { 261 if (!qf->source_path && qf->user_source) { 262 Ceed ceed; 263 bool is_absolute_path; 264 char *source_path_copy; 265 const char *absolute_path; 266 const char *kernel_name = strrchr(qf->user_source, ':') + 1; 267 size_t kernel_name_len = strlen(kernel_name); 268 269 CeedCall(CeedQFunctionGetCeed(qf, &ceed)); 270 CeedCall(CeedCheckFilePath(ceed, qf->user_source, &is_absolute_path)); 271 if (is_absolute_path) { 272 absolute_path = (char *)qf->user_source; 273 } else { 274 CeedCall(CeedGetJitAbsolutePath(ceed, qf->user_source, &absolute_path)); 275 } 276 CeedCall(CeedDestroy(&ceed)); 277 278 size_t source_len = strlen(absolute_path) - kernel_name_len - 1; 279 280 CeedCall(CeedCalloc(source_len + 1, &source_path_copy)); 281 memcpy(source_path_copy, absolute_path, source_len); 282 qf->source_path = source_path_copy; 283 284 if (!is_absolute_path) CeedCall(CeedFree(&absolute_path)); 285 } 286 287 *source_path = (char *)qf->source_path; 288 return CEED_ERROR_SUCCESS; 289 } 290 291 /** 292 @brief Initialize and load `CeedQFunction` source file into string buffer, including full text of local files in place of `#include "local.h"`. 293 294 The `buffer` is set to `NULL` if there is no `CeedQFunction` source file. 295 296 Note: This function may as well return a mutable buffer, but all current uses 297 do not modify it. (This is just a downside of `const` semantics with output 298 arguments instead of returns.) 299 300 Note: Caller is responsible for freeing the string buffer with @ref CeedFree(). 301 302 @param[in] qf `CeedQFunction` 303 @param[out] source_buffer String buffer for source file contents 304 305 @return An error code: 0 - success, otherwise - failure 306 307 @ref Backend 308 **/ 309 int CeedQFunctionLoadSourceToBuffer(CeedQFunction qf, const char **source_buffer) { 310 const char *source_path; 311 312 CeedCall(CeedQFunctionGetSourcePath(qf, &source_path)); 313 *source_buffer = NULL; 314 if (source_path) { 315 Ceed ceed; 316 char *buffer = NULL; 317 318 CeedCall(CeedQFunctionGetCeed(qf, &ceed)); 319 CeedCall(CeedLoadSourceToBuffer(ceed, source_path, &buffer)); 320 CeedCall(CeedDestroy(&ceed)); 321 *source_buffer = buffer; 322 } 323 return CEED_ERROR_SUCCESS; 324 } 325 326 /** 327 @brief Get the User Function for a `CeedQFunction` 328 329 @param[in] qf `CeedQFunction` 330 @param[out] f Variable to store user function 331 332 @return An error code: 0 - success, otherwise - failure 333 334 @ref Backend 335 **/ 336 int CeedQFunctionGetUserFunction(CeedQFunction qf, CeedQFunctionUser *f) { 337 *f = qf->function; 338 return CEED_ERROR_SUCCESS; 339 } 340 341 /** 342 @brief Get the number of tabs to indent for @ref CeedQFunctionView() output 343 344 @param[in] qf `CeedQFunction` to get the number of view tabs 345 @param[out] num_tabs Number of view tabs 346 347 @return Error code: 0 - success, otherwise - failure 348 349 @ref Backend 350 **/ 351 int CeedQFunctionGetNumViewTabs(CeedQFunction qf, CeedInt *num_tabs) { 352 *num_tabs = qf->num_tabs; 353 return CEED_ERROR_SUCCESS; 354 } 355 356 /** 357 @brief Get global context for a `CeedQFunction`. 358 359 Note: For `CeedQFunction` from the Fortran interface, this function will return the Fortran context `CeedQFunctionContext`. 360 361 @param[in] qf CeedQFunction 362 @param[out] ctx Variable to store CeedQFunctionContext 363 364 @return An error code: 0 - success, otherwise - failure 365 366 @ref Backend 367 **/ 368 int CeedQFunctionGetContext(CeedQFunction qf, CeedQFunctionContext *ctx) { 369 *ctx = NULL; 370 if (qf->ctx) CeedCall(CeedQFunctionContextReferenceCopy(qf->ctx, ctx)); 371 return CEED_ERROR_SUCCESS; 372 } 373 374 /** 375 @brief Get context data of a `CeedQFunction` 376 377 @param[in] qf `CeedQFunction` 378 @param[in] mem_type Memory type on which to access the data. 379 If the backend uses a different memory type, this will perform a copy. 380 @param[out] data Data on memory type mem_type 381 382 @return An error code: 0 - success, otherwise - failure 383 384 @ref Backend 385 **/ 386 int CeedQFunctionGetContextData(CeedQFunction qf, CeedMemType mem_type, void *data) { 387 bool is_writable; 388 CeedQFunctionContext ctx; 389 390 CeedCall(CeedQFunctionGetContext(qf, &ctx)); 391 if (ctx) { 392 CeedCall(CeedQFunctionIsContextWritable(qf, &is_writable)); 393 if (is_writable) { 394 CeedCall(CeedQFunctionContextGetData(ctx, mem_type, data)); 395 } else { 396 CeedCall(CeedQFunctionContextGetDataRead(ctx, mem_type, data)); 397 } 398 } else { 399 *(void **)data = NULL; 400 } 401 CeedCall(CeedQFunctionContextDestroy(&ctx)); 402 return CEED_ERROR_SUCCESS; 403 } 404 405 /** 406 @brief Restore context data of a `CeedQFunction` 407 408 @param[in] qf `CeedQFunction` 409 @param[in,out] data Data to restore 410 411 @return An error code: 0 - success, otherwise - failure 412 413 @ref Backend 414 **/ 415 int CeedQFunctionRestoreContextData(CeedQFunction qf, void *data) { 416 bool is_writable; 417 CeedQFunctionContext ctx; 418 419 CeedCall(CeedQFunctionGetContext(qf, &ctx)); 420 if (ctx) { 421 CeedCall(CeedQFunctionIsContextWritable(qf, &is_writable)); 422 if (is_writable) { 423 CeedCall(CeedQFunctionContextRestoreData(ctx, data)); 424 } else { 425 CeedCall(CeedQFunctionContextRestoreDataRead(ctx, data)); 426 } 427 } 428 CeedCall(CeedQFunctionContextDestroy(&ctx)); 429 return CEED_ERROR_SUCCESS; 430 } 431 432 /** 433 @brief Get true user context for a `CeedQFunction` 434 435 Note: For all `CeedQFunction` this function will return the user `CeedQFunctionContext` and not interface context `CeedQFunctionContext`, if any such object exists. 436 437 @param[in] qf `CeedQFunction` 438 @param[out] ctx Variable to store `CeedQFunctionContext` 439 440 @return An error code: 0 - success, otherwise - failure 441 @ref Backend 442 **/ 443 int CeedQFunctionGetInnerContext(CeedQFunction qf, CeedQFunctionContext *ctx) { 444 CeedQFunctionContext qf_ctx; 445 446 CeedCall(CeedQFunctionGetContext(qf, &qf_ctx)); 447 if (qf->is_fortran) { 448 CeedFortranContext fortran_ctx = NULL; 449 450 CeedCall(CeedQFunctionContextGetData(qf_ctx, CEED_MEM_HOST, &fortran_ctx)); 451 *ctx = fortran_ctx->inner_ctx; 452 CeedCall(CeedQFunctionContextRestoreData(qf_ctx, &fortran_ctx)); 453 } else { 454 *ctx = qf_ctx; 455 } 456 CeedCall(CeedQFunctionContextDestroy(&qf_ctx)); 457 return CEED_ERROR_SUCCESS; 458 } 459 460 /** 461 @brief Get inner context data of a `CeedQFunction` 462 463 @param[in] qf `CeedQFunction` 464 @param[in] mem_type Memory type on which to access the data. 465 If the backend uses a different memory type, this will perform a copy. 466 @param[out] data Data on memory type mem_type 467 468 @return An error code: 0 - success, otherwise - failure 469 470 @ref Backend 471 **/ 472 int CeedQFunctionGetInnerContextData(CeedQFunction qf, CeedMemType mem_type, void *data) { 473 bool is_writable; 474 CeedQFunctionContext ctx; 475 476 CeedCall(CeedQFunctionGetInnerContext(qf, &ctx)); 477 if (ctx) { 478 CeedCall(CeedQFunctionIsContextWritable(qf, &is_writable)); 479 if (is_writable) { 480 CeedCall(CeedQFunctionContextGetData(ctx, mem_type, data)); 481 } else { 482 CeedCall(CeedQFunctionContextGetDataRead(ctx, mem_type, data)); 483 } 484 } else { 485 *(void **)data = NULL; 486 } 487 return CEED_ERROR_SUCCESS; 488 } 489 490 /** 491 @brief Restore inner context data of a `CeedQFunction` 492 493 @param[in] qf `CeedQFunction` 494 @param[in,out] data Data to restore 495 496 @return An error code: 0 - success, otherwise - failure 497 498 @ref Backend 499 **/ 500 int CeedQFunctionRestoreInnerContextData(CeedQFunction qf, void *data) { 501 bool is_writable; 502 CeedQFunctionContext ctx; 503 504 CeedCall(CeedQFunctionGetInnerContext(qf, &ctx)); 505 if (ctx) { 506 CeedCall(CeedQFunctionIsContextWritable(qf, &is_writable)); 507 if (is_writable) { 508 CeedCall(CeedQFunctionContextRestoreData(ctx, data)); 509 } else { 510 CeedCall(CeedQFunctionContextRestoreDataRead(ctx, data)); 511 } 512 } 513 *(void **)data = NULL; 514 return CEED_ERROR_SUCCESS; 515 } 516 517 /** 518 @brief Determine if `CeedQFunction` is identity 519 520 @param[in] qf `CeedQFunction` 521 @param[out] is_identity Variable to store identity status 522 523 @return An error code: 0 - success, otherwise - failure 524 525 @ref Backend 526 **/ 527 int CeedQFunctionIsIdentity(CeedQFunction qf, bool *is_identity) { 528 *is_identity = qf->is_identity; 529 return CEED_ERROR_SUCCESS; 530 } 531 532 /** 533 @brief Determine if `CeedQFunctionContext` is writable 534 535 @param[in] qf `CeedQFunction` 536 @param[out] is_writable Variable to store context writeable status 537 538 @return An error code: 0 - success, otherwise - failure 539 540 @ref Backend 541 **/ 542 int CeedQFunctionIsContextWritable(CeedQFunction qf, bool *is_writable) { 543 *is_writable = qf->is_context_writable; 544 return CEED_ERROR_SUCCESS; 545 } 546 547 /** 548 @brief Get backend data of a `CeedQFunction` 549 550 @param[in] qf `CeedQFunction` 551 @param[out] data Variable to store data 552 553 @return An error code: 0 - success, otherwise - failure 554 555 @ref Backend 556 **/ 557 int CeedQFunctionGetData(CeedQFunction qf, void *data) { 558 *(void **)data = qf->data; 559 return CEED_ERROR_SUCCESS; 560 } 561 562 /** 563 @brief Set backend data of a `CeedQFunction` 564 565 @param[in,out] qf `CeedQFunction` 566 @param[in] data Data to set 567 568 @return An error code: 0 - success, otherwise - failure 569 570 @ref Backend 571 **/ 572 int CeedQFunctionSetData(CeedQFunction qf, void *data) { 573 qf->data = data; 574 return CEED_ERROR_SUCCESS; 575 } 576 577 /** 578 @brief Get a boolean value indicating if the `CeedQFunction` is immutable 579 580 @param[in] qf `CeedOperator` 581 @param[out] is_immutable Variable to store immutability status 582 583 @return An error code: 0 - success, otherwise - failure 584 585 @ref Backend 586 **/ 587 int CeedQFunctionIsImmutable(CeedQFunction qf, bool *is_immutable) { 588 *is_immutable = qf->is_immutable; 589 return CEED_ERROR_SUCCESS; 590 } 591 592 /** 593 @brief Set the immutable flag of a `CeedQFunction` to `true` 594 595 @param[in,out] qf `CeedQFunction` 596 597 @return An error code: 0 - success, otherwise - failure 598 599 @ref Backend 600 **/ 601 int CeedQFunctionSetImmutable(CeedQFunction qf) { 602 qf->is_immutable = true; 603 return CEED_ERROR_SUCCESS; 604 } 605 606 /** 607 @brief Increment the reference counter for a `CeedQFunction` 608 609 @param[in,out] qf `CeedQFunction` to increment the reference counter 610 611 @return An error code: 0 - success, otherwise - failure 612 613 @ref Backend 614 **/ 615 int CeedQFunctionReference(CeedQFunction qf) { 616 qf->ref_count++; 617 return CEED_ERROR_SUCCESS; 618 } 619 620 /** 621 @brief Estimate number of FLOPs per quadrature required to apply `CeedQFunction` 622 623 @param[in] qf `CeedQFunction` to estimate FLOPs for 624 @param[out] flops Address of variable to hold FLOPs estimate 625 626 @ref Backend 627 **/ 628 int CeedQFunctionGetFlopsEstimate(CeedQFunction qf, CeedSize *flops) { 629 *flops = qf->user_flop_estimate; 630 return CEED_ERROR_SUCCESS; 631 } 632 633 /// @} 634 635 /// ---------------------------------------------------------------------------- 636 /// CeedQFunction Public API 637 /// ---------------------------------------------------------------------------- 638 /// @addtogroup CeedQFunctionUser 639 /// @{ 640 641 /** 642 @brief Create a `CeedQFunction` for evaluating interior (volumetric) terms 643 644 @param[in] ceed `Ceed` object used to create the `CeedQFunction` 645 @param[in] vec_length Vector length. 646 Caller must ensure that number of quadrature points is a multiple of `vec_length`. 647 @param[in] f Function pointer to evaluate action at quadrature points. 648 See `CeedQFunctionUser`. 649 @param[in] source Absolute path to source of `CeedQFunctionUser`, "\abs_path\file.h:function_name". 650 The entire source file must only contain constructs supported by all targeted backends (i.e. CUDA for `/gpu/cuda`, OpenCL/SYCL for `/gpu/sycl`, etc.). 651 The entire contents of this file and all locally included files are used during JiT compilation for GPU backends. 652 The header `ceed/types.h` is preferred over `ceed.h` or `ceed/ceed.h` for `CeedQFunction` source files. 653 The macro `CEED_RUNNING_JIT_PASS` is set during JiT and can be used to guard include statements that JiT compilers cannot use, such as `math.h` or `std*.h`. 654 All source files must be at the provided filepath at runtime for JiT to function. 655 @param[out] qf Address of the variable where the newly created `CeedQFunction` will be stored 656 657 @return An error code: 0 - success, otherwise - failure 658 659 See \ref CeedQFunctionUser for details on the call-back function `f` arguments. 660 661 @ref User 662 **/ 663 int CeedQFunctionCreateInterior(Ceed ceed, CeedInt vec_length, CeedQFunctionUser f, const char *source, CeedQFunction *qf) { 664 char *user_source_copy; 665 666 if (!ceed->QFunctionCreate) { 667 Ceed delegate; 668 669 CeedCall(CeedGetObjectDelegate(ceed, &delegate, "QFunction")); 670 CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement CeedQFunctionCreateInterior"); 671 CeedCall(CeedQFunctionCreateInterior(delegate, vec_length, f, source, qf)); 672 CeedCall(CeedDestroy(&delegate)); 673 return CEED_ERROR_SUCCESS; 674 } 675 676 CeedCheck(!strlen(source) || strrchr(source, ':'), ceed, CEED_ERROR_INCOMPLETE, 677 "Provided path to source does not include function name. Provided: \"%s\"\nRequired: \"\\abs_path\\file.h:function_name\"", source); 678 679 CeedCall(CeedCalloc(1, qf)); 680 CeedCall(CeedReferenceCopy(ceed, &(*qf)->ceed)); 681 (*qf)->ref_count = 1; 682 (*qf)->vec_length = vec_length; 683 (*qf)->is_identity = false; 684 (*qf)->is_context_writable = true; 685 (*qf)->function = f; 686 (*qf)->user_flop_estimate = -1; 687 if (strlen(source)) { 688 size_t user_source_len = strlen(source); 689 690 CeedCall(CeedCalloc(user_source_len + 1, &user_source_copy)); 691 memcpy(user_source_copy, source, user_source_len); 692 (*qf)->user_source = user_source_copy; 693 } 694 CeedCall(CeedCalloc(CEED_FIELD_MAX, &(*qf)->input_fields)); 695 CeedCall(CeedCalloc(CEED_FIELD_MAX, &(*qf)->output_fields)); 696 CeedCall(ceed->QFunctionCreate(*qf)); 697 return CEED_ERROR_SUCCESS; 698 } 699 700 /** 701 @brief Create a `CeedQFunction` for evaluating interior (volumetric) terms by name 702 703 @param[in] ceed `Ceed` object used to create the `CeedQFunction` 704 @param[in] name Name of `CeedQFunction` to use from gallery 705 @param[out] qf Address of the variable where the newly created `CeedQFunction` will be stored 706 707 @return An error code: 0 - success, otherwise - failure 708 709 @ref User 710 **/ 711 int CeedQFunctionCreateInteriorByName(Ceed ceed, const char *name, CeedQFunction *qf) { 712 size_t match_len = 0, match_index = UINT_MAX; 713 714 CeedCall(CeedQFunctionRegisterAll()); 715 // Find matching backend 716 CeedCheck(name, ceed, CEED_ERROR_INCOMPLETE, "No CeedQFunction name provided"); 717 for (size_t i = 0; i < num_qfunctions; i++) { 718 size_t n; 719 const char *curr_name = gallery_qfunctions[i].name; 720 for (n = 0; curr_name[n] && curr_name[n] == name[n]; n++) { 721 } 722 if (n > match_len) { 723 match_len = n; 724 match_index = i; 725 } 726 } 727 CeedCheck(match_len > 0, ceed, CEED_ERROR_UNSUPPORTED, "No suitable gallery CeedQFunction"); 728 729 // Create QFunction 730 CeedCall(CeedQFunctionCreateInterior(ceed, gallery_qfunctions[match_index].vec_length, gallery_qfunctions[match_index].f, 731 gallery_qfunctions[match_index].source, qf)); 732 733 // QFunction specific setup 734 CeedCall(gallery_qfunctions[match_index].init(ceed, name, *qf)); 735 736 // Copy name 737 CeedCall(CeedStringAllocCopy(name, (char **)&(*qf)->gallery_name)); 738 (*qf)->is_gallery = true; 739 return CEED_ERROR_SUCCESS; 740 } 741 742 /** 743 @brief Create an identity `CeedQFunction`. 744 745 Inputs are written into outputs in the order given. 746 This is useful for `CeedOperator that can be represented with only the action of a `CeedElemRestriction` and `CeedBasis`, such as restriction and prolongation operators for p-multigrid. 747 Backends may optimize `CeedOperator` with this `CeedQFunction` to avoid the copy of input data to output fields by using the same memory location for both. 748 749 @param[in] ceed `Ceed` object used to create the `CeedQFunction` 750 @param[in] size Size of the `CeedQFunction` fields 751 @param[in] in_mode @ref CeedEvalMode for input to `CeedQFunction` 752 @param[in] out_mode @ref CeedEvalMode for output to `CeedQFunction` 753 @param[out] qf Address of the variable where the newly created `CeedQFunction` will be stored 754 755 @return An error code: 0 - success, otherwise - failure 756 757 @ref User 758 **/ 759 int CeedQFunctionCreateIdentity(Ceed ceed, CeedInt size, CeedEvalMode in_mode, CeedEvalMode out_mode, CeedQFunction *qf) { 760 CeedQFunctionContext ctx; 761 CeedContextFieldLabel size_label; 762 763 CeedCall(CeedQFunctionCreateInteriorByName(ceed, "Identity", qf)); 764 CeedCall(CeedQFunctionAddInput(*qf, "input", size, in_mode)); 765 CeedCall(CeedQFunctionAddOutput(*qf, "output", size, out_mode)); 766 767 (*qf)->is_identity = true; 768 769 CeedCall(CeedQFunctionGetContext(*qf, &ctx)); 770 CeedCall(CeedQFunctionContextGetFieldLabel(ctx, "size", &size_label)); 771 CeedCall(CeedQFunctionContextSetInt32(ctx, size_label, &size)); 772 CeedCall(CeedQFunctionContextDestroy(&ctx)); 773 return CEED_ERROR_SUCCESS; 774 } 775 776 /** 777 @brief Copy the pointer to a `CeedQFunction`. 778 779 Both pointers should be destroyed with @ref CeedQFunctionDestroy(). 780 781 Note: If the value of `*qf_copy` passed to this function is non-NULL, then it is assumed that `*qf_copy` is a pointer to a `CeedQFunction`. 782 This `CeedQFunction` will be destroyed if `*qf_copy` is the only reference to this `CeedQFunction`. 783 784 @param[in] qf `CeedQFunction` to copy reference to 785 @param[out] qf_copy Variable to store copied reference 786 787 @return An error code: 0 - success, otherwise - failure 788 789 @ref User 790 **/ 791 int CeedQFunctionReferenceCopy(CeedQFunction qf, CeedQFunction *qf_copy) { 792 CeedCall(CeedQFunctionReference(qf)); 793 CeedCall(CeedQFunctionDestroy(qf_copy)); 794 *qf_copy = qf; 795 return CEED_ERROR_SUCCESS; 796 } 797 798 /** 799 @brief Add a `CeedQFunction` input 800 801 @param[in,out] qf `CeedQFunction` 802 @param[in] field_name Name of `CeedQFunction` field 803 @param[in] size Size of `CeedQFunction` field, 804 (`num_comp * 1`) for @ref CEED_EVAL_NONE, 805 (`num_comp * 1`) for @ref CEED_EVAL_INTERP for an \f$H^1\f$ space or (`num_comp * dim`) for an \f$H(\mathrm{div})\f$ or \f$H(\mathrm{curl})\f$ space, 806 (`num_comp * dim`) for @ref CEED_EVAL_GRAD, 807 (`num_comp * 1`) for @ref CEED_EVAL_DIV, and 808 (`num_comp * curl_dim`) with `curl_dim = 1` if `dim < 3` otherwise `curl_dim = dim` for @ref CEED_EVAL_CURL. 809 @param[in] eval_mode @ref CEED_EVAL_NONE to use values directly, 810 @ref CEED_EVAL_INTERP to use interpolated values, 811 @ref CEED_EVAL_GRAD to use gradients, 812 @ref CEED_EVAL_DIV to use divergence, 813 @ref CEED_EVAL_CURL to use curl 814 815 Note: In the user `CeedQFunctionUser`, the `in` argument list the fields in the order given by the calls to `CeedQFunctionAddInput`. 816 817 @return An error code: 0 - success, otherwise - failure 818 819 @ref User 820 **/ 821 int CeedQFunctionAddInput(CeedQFunction qf, const char *field_name, CeedInt size, CeedEvalMode eval_mode) { 822 bool is_immutable; 823 824 CeedCall(CeedQFunctionIsImmutable(qf, &is_immutable)); 825 CeedCheck(!is_immutable, CeedQFunctionReturnCeed(qf), CEED_ERROR_MAJOR, "QFunction cannot be changed after set as immutable"); 826 CeedCheck(eval_mode != CEED_EVAL_WEIGHT || size == 1, CeedQFunctionReturnCeed(qf), CEED_ERROR_DIMENSION, "CEED_EVAL_WEIGHT should have size 1"); 827 for (CeedInt i = 0; i < qf->num_input_fields; i++) { 828 CeedCheck(strcmp(field_name, qf->input_fields[i]->field_name), CeedQFunctionReturnCeed(qf), CEED_ERROR_MINOR, 829 "CeedQFunction field names must be unique. Duplicate name: %s", field_name); 830 } 831 for (CeedInt i = 0; i < qf->num_output_fields; i++) { 832 CeedCheck(strcmp(field_name, qf->output_fields[i]->field_name), CeedQFunctionReturnCeed(qf), CEED_ERROR_MINOR, 833 "CeedQFunction field names must be unique. Duplicate name: %s", field_name); 834 } 835 CeedCall(CeedQFunctionFieldSet(&qf->input_fields[qf->num_input_fields], field_name, size, eval_mode)); 836 qf->num_input_fields++; 837 return CEED_ERROR_SUCCESS; 838 } 839 840 /** 841 @brief Add a `CeedQFunction` output 842 843 @param[in,out] qf `CeedQFunction` 844 @param[in] field_name Name of `CeedQFunction` field 845 @param[in] size Size of `CeedQFunction` field, 846 (`num_comp * 1`) for @ref CEED_EVAL_NONE, 847 (`num_comp * 1`) for @ref CEED_EVAL_INTERP for an \f$H^1\f$ space or (`num_comp * dim`) for an \f$H(\mathrm{div})\f$ or \f$H(\mathrm{curl})\f$ space, 848 (`num_comp * dim`) for @ref CEED_EVAL_GRAD, 849 (`num_comp * 1`) for @ref CEED_EVAL_DIV, and 850 (`num_comp * curl_dim`) with `curl_dim = 1` if `dim < 3` otherwise `curl_dim = dim` for @ref CEED_EVAL_CURL. 851 @param[in] eval_mode @ref CEED_EVAL_NONE to use values directly, 852 @ref CEED_EVAL_INTERP to use interpolated values, 853 @ref CEED_EVAL_GRAD to use gradients, 854 @ref CEED_EVAL_DIV to use divergence, 855 @ref CEED_EVAL_CURL to use curl. 856 857 Note: In the user `CeedQFunctionUser`, the `out` argument list the fields in the order given by the calls to `CeedQFunctionAddOutput`. 858 859 @return An error code: 0 - success, otherwise - failure 860 861 @ref User 862 **/ 863 int CeedQFunctionAddOutput(CeedQFunction qf, const char *field_name, CeedInt size, CeedEvalMode eval_mode) { 864 bool is_immutable; 865 866 CeedCall(CeedQFunctionIsImmutable(qf, &is_immutable)); 867 CeedCheck(!is_immutable, CeedQFunctionReturnCeed(qf), CEED_ERROR_MAJOR, "CeedQFunction cannot be changed after set as immutable"); 868 CeedCheck(eval_mode != CEED_EVAL_WEIGHT, CeedQFunctionReturnCeed(qf), CEED_ERROR_DIMENSION, 869 "Cannot create CeedQFunction output with CEED_EVAL_WEIGHT"); 870 for (CeedInt i = 0; i < qf->num_input_fields; i++) { 871 CeedCheck(strcmp(field_name, qf->input_fields[i]->field_name), CeedQFunctionReturnCeed(qf), CEED_ERROR_MINOR, 872 "CeedQFunction field names must be unique"); 873 } 874 for (CeedInt i = 0; i < qf->num_output_fields; i++) { 875 CeedCheck(strcmp(field_name, qf->output_fields[i]->field_name), CeedQFunctionReturnCeed(qf), CEED_ERROR_MINOR, 876 "CeedQFunction field names must be unique"); 877 } 878 CeedCall(CeedQFunctionFieldSet(&qf->output_fields[qf->num_output_fields], field_name, size, eval_mode)); 879 qf->num_output_fields++; 880 return CEED_ERROR_SUCCESS; 881 } 882 883 /** 884 @brief Get the `CeedQFunctionField` of a `CeedQFunction` 885 886 Note: Calling this function asserts that setup is complete and sets the `CeedQFunction` as immutable. 887 888 @param[in] qf `CeedQFunction` 889 @param[out] num_input_fields Variable to store number of input fields 890 @param[out] input_fields Variable to store input fields 891 @param[out] num_output_fields Variable to store number of output fields 892 @param[out] output_fields Variable to store output fields 893 894 @return An error code: 0 - success, otherwise - failure 895 896 @ref Advanced 897 **/ 898 int CeedQFunctionGetFields(CeedQFunction qf, CeedInt *num_input_fields, CeedQFunctionField **input_fields, CeedInt *num_output_fields, 899 CeedQFunctionField **output_fields) { 900 CeedCall(CeedQFunctionSetImmutable(qf)); 901 if (num_input_fields) *num_input_fields = qf->num_input_fields; 902 if (input_fields) *input_fields = qf->input_fields; 903 if (num_output_fields) *num_output_fields = qf->num_output_fields; 904 if (output_fields) *output_fields = qf->output_fields; 905 return CEED_ERROR_SUCCESS; 906 } 907 908 /** 909 @brief Get the name of a `CeedQFunctionField` 910 911 @param[in] qf_field `CeedQFunctionField` 912 @param[out] field_name Variable to store the field name 913 914 @return An error code: 0 - success, otherwise - failure 915 916 @ref Advanced 917 **/ 918 int CeedQFunctionFieldGetName(CeedQFunctionField qf_field, const char **field_name) { 919 *field_name = qf_field->field_name; 920 return CEED_ERROR_SUCCESS; 921 } 922 923 /** 924 @brief Get the number of components of a `CeedQFunctionField` 925 926 @param[in] qf_field `CeedQFunctionField` 927 @param[out] size Variable to store the size of the field 928 929 @return An error code: 0 - success, otherwise - failure 930 931 @ref Advanced 932 **/ 933 int CeedQFunctionFieldGetSize(CeedQFunctionField qf_field, CeedInt *size) { 934 *size = qf_field->size; 935 return CEED_ERROR_SUCCESS; 936 } 937 938 /** 939 @brief Get the @ref CeedEvalMode of a `CeedQFunctionField` 940 941 @param[in] qf_field `CeedQFunctionField` 942 @param[out] eval_mode Variable to store the field evaluation mode 943 944 @return An error code: 0 - success, otherwise - failure 945 946 @ref Advanced 947 **/ 948 int CeedQFunctionFieldGetEvalMode(CeedQFunctionField qf_field, CeedEvalMode *eval_mode) { 949 *eval_mode = qf_field->eval_mode; 950 return CEED_ERROR_SUCCESS; 951 } 952 953 /** 954 @brief Get the data of a `CeedQFunctionField`. 955 956 Any arguments set as `NULL` are ignored. 957 958 @param[in] qf_field `CeedQFunctionField` 959 @param[out] field_name Variable to store the field name 960 @param[out] size Variable to store the size of the field 961 @param[out] eval_mode Variable to store the field evaluation mode 962 963 @return An error code: 0 - success, otherwise - failure 964 965 @ref Advanced 966 **/ 967 int CeedQFunctionFieldGetData(CeedQFunctionField qf_field, const char **field_name, CeedInt *size, CeedEvalMode *eval_mode) { 968 if (field_name) CeedCall(CeedQFunctionFieldGetName(qf_field, field_name)); 969 if (size) CeedCall(CeedQFunctionFieldGetSize(qf_field, size)); 970 if (eval_mode) CeedCall(CeedQFunctionFieldGetEvalMode(qf_field, eval_mode)); 971 return CEED_ERROR_SUCCESS; 972 } 973 974 /** 975 @brief Set global context for a `CeedQFunction` 976 977 @param[in,out] qf `CeedQFunction` 978 @param[in] ctx Context data to set 979 980 @return An error code: 0 - success, otherwise - failure 981 982 @ref User 983 **/ 984 int CeedQFunctionSetContext(CeedQFunction qf, CeedQFunctionContext ctx) { 985 CeedCall(CeedQFunctionContextDestroy(&qf->ctx)); 986 qf->ctx = ctx; 987 if (ctx) CeedCall(CeedQFunctionContextReference(ctx)); 988 return CEED_ERROR_SUCCESS; 989 } 990 991 /** 992 @brief Set writability of `CeedQFunctionContext` when calling the `CeedQFunctionUser`. 993 994 The default value is `is_writable == true`. 995 996 Setting `is_writable == true` indicates the `CeedQFunctionUser` writes into the `CeedQFunctionContext` and requires memory synchronization after calling @ref CeedQFunctionApply(). 997 998 Setting 'is_writable == false' asserts that `CeedQFunctionUser` does not modify the `CeedQFunctionContext`. 999 Violating this assertion may lead to inconsistent data. 1000 1001 Setting `is_writable == false` may offer a performance improvement on GPU backends. 1002 1003 @param[in,out] qf `CeedQFunction` 1004 @param[in] is_writable Boolean flag for writability status 1005 1006 @return An error code: 0 - success, otherwise - failure 1007 1008 @ref User 1009 **/ 1010 int CeedQFunctionSetContextWritable(CeedQFunction qf, bool is_writable) { 1011 qf->is_context_writable = is_writable; 1012 return CEED_ERROR_SUCCESS; 1013 } 1014 1015 /** 1016 @brief Set estimated number of FLOPs per quadrature required to apply `CeedQFunction` 1017 1018 @param[in] qf `CeedQFunction` to estimate FLOPs for 1019 @param[out] flops FLOPs per quadrature point estimate 1020 1021 @ref Backend 1022 **/ 1023 int CeedQFunctionSetUserFlopsEstimate(CeedQFunction qf, CeedSize flops) { 1024 CeedCheck(flops >= 0, CeedQFunctionReturnCeed(qf), CEED_ERROR_INCOMPATIBLE, "Must set non-negative FLOPs estimate"); 1025 qf->user_flop_estimate = flops; 1026 return CEED_ERROR_SUCCESS; 1027 } 1028 1029 /** 1030 @brief Set the number of tabs to indent for @ref CeedQFunctionView() output 1031 1032 @param[in] qf `CeedQFunction` to set the number of view tabs 1033 @param[in] num_tabs Number of view tabs to set 1034 1035 @return Error code: 0 - success, otherwise - failure 1036 1037 @ref User 1038 **/ 1039 int CeedQFunctionSetNumViewTabs(CeedQFunction qf, CeedInt num_tabs) { 1040 CeedCheck(num_tabs >= 0, CeedQFunctionReturnCeed(qf), CEED_ERROR_MINOR, "Number of view tabs must be non-negative"); 1041 qf->num_tabs = num_tabs; 1042 return CEED_ERROR_SUCCESS; 1043 } 1044 1045 /** 1046 @brief View a `CeedQFunction` 1047 1048 @param[in] qf `CeedQFunction` to view 1049 @param[in] stream Stream to write; typically `stdout` or a file 1050 1051 @return Error code: 0 - success, otherwise - failure 1052 1053 @ref User 1054 **/ 1055 int CeedQFunctionView(CeedQFunction qf, FILE *stream) { 1056 char *tabs = NULL; 1057 const char *name; 1058 1059 { 1060 CeedInt num_tabs = 0; 1061 1062 CeedCall(CeedQFunctionGetNumViewTabs(qf, &num_tabs)); 1063 CeedCall(CeedCalloc(CEED_TAB_WIDTH * num_tabs + 1, &tabs)); 1064 for (CeedInt i = 0; i < CEED_TAB_WIDTH * num_tabs; i++) tabs[i] = ' '; 1065 } 1066 1067 CeedCall(CeedQFunctionGetName(qf, &name)); 1068 fprintf(stream, "%s%sCeedQFunction - %s\n", tabs, qf->is_gallery ? "Gallery " : "User ", name); 1069 1070 fprintf(stream, "%s %" CeedInt_FMT " input field%s:\n", tabs, qf->num_input_fields, qf->num_input_fields > 1 ? "s" : ""); 1071 for (CeedInt i = 0; i < qf->num_input_fields; i++) { 1072 CeedCall(CeedQFunctionFieldView(qf->input_fields[i], i, 1, tabs, stream)); 1073 } 1074 1075 fprintf(stream, "%s %" CeedInt_FMT " output field%s:\n", tabs, qf->num_output_fields, qf->num_output_fields > 1 ? "s" : ""); 1076 for (CeedInt i = 0; i < qf->num_output_fields; i++) { 1077 CeedCall(CeedQFunctionFieldView(qf->output_fields[i], i, 0, tabs, stream)); 1078 } 1079 CeedCall(CeedFree(&tabs)); 1080 return CEED_ERROR_SUCCESS; 1081 } 1082 1083 /** 1084 @brief Get the `Ceed` associated with a `CeedQFunction` 1085 1086 @param[in] qf `CeedQFunction` 1087 @param[out] ceed Variable to store`Ceed` 1088 1089 @return An error code: 0 - success, otherwise - failure 1090 1091 @ref Advanced 1092 **/ 1093 int CeedQFunctionGetCeed(CeedQFunction qf, Ceed *ceed) { 1094 *ceed = NULL; 1095 CeedCall(CeedReferenceCopy(CeedQFunctionReturnCeed(qf), ceed)); 1096 return CEED_ERROR_SUCCESS; 1097 } 1098 1099 /** 1100 @brief Return the `Ceed` associated with a `CeedQFunction` 1101 1102 @param[in] qf `CeedQFunction` 1103 1104 @return `Ceed` associated with the `qf` 1105 1106 @ref Advanced 1107 **/ 1108 Ceed CeedQFunctionReturnCeed(CeedQFunction qf) { return qf->ceed; } 1109 1110 /** 1111 @brief Apply the action of a `CeedQFunction` 1112 1113 Note: Calling this function asserts that setup is complete and sets the `CeedQFunction` as immutable. 1114 1115 @param[in] qf `CeedQFunction` 1116 @param[in] Q Number of quadrature points 1117 @param[in] u Array of input `CeedVector` 1118 @param[out] v Array of output `CeedVector` 1119 1120 @return An error code: 0 - success, otherwise - failure 1121 1122 @ref User 1123 **/ 1124 int CeedQFunctionApply(CeedQFunction qf, CeedInt Q, CeedVector *u, CeedVector *v) { 1125 CeedInt vec_length; 1126 1127 CeedCheck(qf->Apply, CeedQFunctionReturnCeed(qf), CEED_ERROR_UNSUPPORTED, "Backend does not support CeedQFunctionApply"); 1128 CeedCall(CeedQFunctionGetVectorLength(qf, &vec_length)); 1129 CeedCheck(Q % vec_length == 0, CeedQFunctionReturnCeed(qf), CEED_ERROR_DIMENSION, 1130 "Number of quadrature points %" CeedInt_FMT " must be a multiple of %" CeedInt_FMT, Q, qf->vec_length); 1131 CeedCall(CeedQFunctionSetImmutable(qf)); 1132 CeedCall(qf->Apply(qf, Q, u, v)); 1133 return CEED_ERROR_SUCCESS; 1134 } 1135 1136 /** 1137 @brief Destroy a `CeedQFunction` 1138 1139 @param[in,out] qf `CeedQFunction` to destroy 1140 1141 @return An error code: 0 - success, otherwise - failure 1142 1143 @ref User 1144 **/ 1145 int CeedQFunctionDestroy(CeedQFunction *qf) { 1146 if (!*qf || --(*qf)->ref_count > 0) { 1147 *qf = NULL; 1148 return CEED_ERROR_SUCCESS; 1149 } 1150 // Backend destroy 1151 if ((*qf)->Destroy) { 1152 CeedCall((*qf)->Destroy(*qf)); 1153 } 1154 // Free fields 1155 for (CeedInt i = 0; i < (*qf)->num_input_fields; i++) { 1156 CeedCall(CeedFree(&(*(*qf)->input_fields[i]).field_name)); 1157 CeedCall(CeedFree(&(*qf)->input_fields[i])); 1158 } 1159 for (CeedInt i = 0; i < (*qf)->num_output_fields; i++) { 1160 CeedCall(CeedFree(&(*(*qf)->output_fields[i]).field_name)); 1161 CeedCall(CeedFree(&(*qf)->output_fields[i])); 1162 } 1163 CeedCall(CeedFree(&(*qf)->input_fields)); 1164 CeedCall(CeedFree(&(*qf)->output_fields)); 1165 1166 // User context data object 1167 CeedCall(CeedQFunctionContextDestroy(&(*qf)->ctx)); 1168 1169 CeedCall(CeedFree(&(*qf)->user_source)); 1170 CeedCall(CeedFree(&(*qf)->source_path)); 1171 CeedCall(CeedFree(&(*qf)->gallery_name)); 1172 CeedCall(CeedFree(&(*qf)->kernel_name)); 1173 CeedCall(CeedDestroy(&(*qf)->ceed)); 1174 CeedCall(CeedFree(qf)); 1175 return CEED_ERROR_SUCCESS; 1176 } 1177 1178 /// @} 1179