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