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