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 @param[out] qf Address of the variable where the newly created 444 CeedQFunction will be stored 445 446 @return An error code: 0 - success, otherwise - failure 447 448 See \ref CeedQFunctionUser for details on the call-back function @a f's 449 arguments. 450 451 @ref User 452 **/ 453 int CeedQFunctionCreateInterior(Ceed ceed, CeedInt vlength, CeedQFunctionUser f, 454 const char *source, CeedQFunction *qf) { 455 int ierr; 456 char *source_copy; 457 458 if (!ceed->QFunctionCreate) { 459 Ceed delegate; 460 ierr = CeedGetObjectDelegate(ceed, &delegate, "QFunction"); CeedChk(ierr); 461 462 if (!delegate) 463 // LCOV_EXCL_START 464 return CeedError(ceed, 1, "Backend does not support QFunctionCreate"); 465 // LCOV_EXCL_STOP 466 467 ierr = CeedQFunctionCreateInterior(delegate, vlength, f, source, qf); 468 CeedChk(ierr); 469 return 0; 470 } 471 472 ierr = CeedCalloc(1, qf); CeedChk(ierr); 473 (*qf)->ceed = ceed; 474 ceed->refcount++; 475 (*qf)->refcount = 1; 476 (*qf)->vlength = vlength; 477 (*qf)->identity = 0; 478 (*qf)->function = f; 479 size_t slen = strlen(source) + 1; 480 ierr = CeedMalloc(slen, &source_copy); CeedChk(ierr); 481 memcpy(source_copy, source, slen); 482 (*qf)->sourcepath = source_copy; 483 ierr = CeedCalloc(16, &(*qf)->inputfields); CeedChk(ierr); 484 ierr = CeedCalloc(16, &(*qf)->outputfields); CeedChk(ierr); 485 ierr = ceed->QFunctionCreate(*qf); CeedChk(ierr); 486 return 0; 487 } 488 489 /** 490 @brief Create a CeedQFunction for evaluating interior (volumetric) terms by name. 491 492 @param ceed A Ceed object where the CeedQFunction will be created 493 @param name Name of QFunction to use from gallery 494 @param[out] qf Address of the variable where the newly created 495 CeedQFunction will be stored 496 497 @return An error code: 0 - success, otherwise - failure 498 499 @ref User 500 **/ 501 int CeedQFunctionCreateInteriorByName(Ceed ceed, const char *name, 502 CeedQFunction *qf) { 503 int ierr; 504 size_t matchlen = 0, matchidx = UINT_MAX; 505 char *name_copy; 506 507 // Find matching backend 508 if (!name) return CeedError(NULL, 1, "No QFunction name provided"); 509 for (size_t i=0; i<num_qfunctions; i++) { 510 size_t n; 511 const char *currname = qfunctions[i].name; 512 for (n = 0; currname[n] && currname[n] == name[n]; n++) {} 513 if (n > matchlen) { 514 matchlen = n; 515 matchidx = i; 516 } 517 } 518 if (!matchlen) 519 // LCOV_EXCL_START 520 return CeedError(NULL, 1, "No suitable gallery QFunction"); 521 // LCOV_EXCL_STOP 522 523 // Create QFunction 524 ierr = CeedQFunctionCreateInterior(ceed, qfunctions[matchidx].vlength, 525 qfunctions[matchidx].f, 526 qfunctions[matchidx].source, qf); 527 CeedChk(ierr); 528 529 // QFunction specific setup 530 ierr = qfunctions[matchidx].init(ceed, name, *qf); CeedChk(ierr); 531 532 // Copy name 533 size_t slen = strlen(name) + 1; 534 ierr = CeedMalloc(slen, &name_copy); CeedChk(ierr); 535 memcpy(name_copy, name, slen); 536 (*qf)->qfname = name_copy; 537 538 return 0; 539 } 540 541 /** 542 @brief Create an identity CeedQFunction. Inputs are written into outputs in 543 the order given. This is useful for CeedOperators that can be 544 represented with only the action of a CeedRestriction and CeedBasis, 545 such as restriction and prolongation operators for p-multigrid. 546 Backends may optimize CeedOperators with this CeedQFunction to avoid 547 the copy of input data to output fields by using the same memory 548 location for both. 549 550 @param ceed A Ceed object where the CeedQFunction will be created 551 @param[in] size Size of the qfunction fields 552 @param[in] inmode CeedEvalMode for input to CeedQFunction 553 @param[in] outmode CeedEvalMode for output to CeedQFunction 554 @param[out] qf Address of the variable where the newly created 555 CeedQFunction will be stored 556 557 @return An error code: 0 - success, otherwise - failure 558 559 @ref User 560 **/ 561 int CeedQFunctionCreateIdentity(Ceed ceed, CeedInt size, CeedEvalMode inmode, 562 CeedEvalMode outmode, CeedQFunction *qf) { 563 int ierr; 564 565 if (inmode == CEED_EVAL_NONE && outmode == CEED_EVAL_NONE) 566 // LCOV_EXCL_START 567 return CeedError(ceed, 1, "CEED_EVAL_NONE for a both the input and " 568 "output does not make sense with an identity QFunction"); 569 // LCOV_EXCL_STOP 570 571 ierr = CeedQFunctionCreateInteriorByName(ceed, "Identity", qf); CeedChk(ierr); 572 ierr = CeedQFunctionAddInput(*qf, "input", 1, inmode); CeedChk(ierr); 573 ierr = CeedQFunctionAddOutput(*qf, "output", 1, outmode); CeedChk(ierr); 574 575 (*qf)->identity = 1; 576 if (size > 1) { 577 CeedInt *ctx; 578 ierr = CeedCalloc(1, &ctx); CeedChk(ierr); 579 ctx[0] = size; 580 ierr = CeedQFunctionSetContext(*qf, ctx, sizeof(ctx)); CeedChk(ierr); 581 (*qf)->inputfields[0]->size = size; 582 (*qf)->outputfields[0]->size = size; 583 } 584 585 return 0; 586 } 587 588 /** 589 @brief Add a CeedQFunction input 590 591 @param qf CeedQFunction 592 @param fieldname Name of QFunction field 593 @param size Size of QFunction field, (ncomp * dim) for @ref CEED_EVAL_GRAD or 594 (ncomp * 1) for @ref CEED_EVAL_NONE and @ref CEED_EVAL_INTERP 595 @param emode \ref CEED_EVAL_NONE to use values directly, 596 \ref CEED_EVAL_INTERP to use interpolated values, 597 \ref CEED_EVAL_GRAD to use gradients. 598 599 @return An error code: 0 - success, otherwise - failure 600 601 @ref User 602 **/ 603 int CeedQFunctionAddInput(CeedQFunction qf, const char *fieldname, CeedInt size, 604 CeedEvalMode emode) { 605 int ierr = CeedQFunctionFieldSet(&qf->inputfields[qf->numinputfields], 606 fieldname, size, emode); 607 CeedChk(ierr); 608 qf->numinputfields++; 609 return 0; 610 } 611 612 /** 613 @brief Add a CeedQFunction output 614 615 @param qf CeedQFunction 616 @param fieldname Name of QFunction field 617 @param size Size of QFunction field, (ncomp * dim) for @ref CEED_EVAL_GRAD or 618 (ncomp * 1) for @ref CEED_EVAL_NONE and @ref CEED_EVAL_INTERP 619 @param emode \ref CEED_EVAL_NONE to use values directly, 620 \ref CEED_EVAL_INTERP to use interpolated values, 621 \ref CEED_EVAL_GRAD to use gradients. 622 623 @return An error code: 0 - success, otherwise - failure 624 625 @ref User 626 **/ 627 int CeedQFunctionAddOutput(CeedQFunction qf, const char *fieldname, 628 CeedInt size, CeedEvalMode emode) { 629 if (emode == CEED_EVAL_WEIGHT) 630 // LCOV_EXCL_START 631 return CeedError(qf->ceed, 1, "Cannot create QFunction output with " 632 "CEED_EVAL_WEIGHT"); 633 // LCOV_EXCL_STOP 634 int ierr = CeedQFunctionFieldSet(&qf->outputfields[qf->numoutputfields], 635 fieldname, size, emode); 636 CeedChk(ierr); 637 qf->numoutputfields++; 638 return 0; 639 } 640 641 /** 642 @brief Set global context for a CeedQFunction 643 644 @param qf CeedQFunction 645 @param ctx Context data to set 646 @param ctxsize Size of context data values 647 648 @return An error code: 0 - success, otherwise - failure 649 650 @ref User 651 **/ 652 int CeedQFunctionSetContext(CeedQFunction qf, void *ctx, size_t ctxsize) { 653 qf->ctx = ctx; 654 qf->ctxsize = ctxsize; 655 return 0; 656 } 657 658 /** 659 @brief View a CeedQFunction 660 661 @param[in] qf CeedQFunction to view 662 @param[in] stream Stream to write; typically stdout/stderr or a file 663 664 @return Error code: 0 - success, otherwise - failure 665 666 @ref User 667 **/ 668 int CeedQFunctionView(CeedQFunction qf, FILE *stream) { 669 int ierr; 670 671 fprintf(stream, "%sCeedQFunction %s\n", 672 qf->qfname ? "Gallery " : "User ", qf->qfname ? qf->qfname : ""); 673 674 fprintf(stream, " %d Input Field%s:\n", qf->numinputfields, 675 qf->numinputfields>1 ? "s" : ""); 676 for (CeedInt i=0; i<qf->numinputfields; i++) { 677 ierr = CeedQFunctionFieldView(qf->inputfields[i], i, 1, stream); 678 CeedChk(ierr); 679 } 680 681 fprintf(stream, " %d Output Field%s:\n", qf->numoutputfields, 682 qf->numoutputfields>1 ? "s" : ""); 683 for (CeedInt i=0; i<qf->numoutputfields; i++) { 684 ierr = CeedQFunctionFieldView(qf->outputfields[i], i, 0, stream); 685 CeedChk(ierr); 686 } 687 return 0; 688 } 689 690 /** 691 @brief Apply the action of a CeedQFunction 692 693 @param qf CeedQFunction 694 @param Q Number of quadrature points 695 @param[in] u Array of input CeedVectors 696 @param[out] v Array of output CeedVectors 697 698 @return An error code: 0 - success, otherwise - failure 699 700 @ref User 701 **/ 702 int CeedQFunctionApply(CeedQFunction qf, CeedInt Q, 703 CeedVector *u, CeedVector *v) { 704 int ierr; 705 if (!qf->Apply) 706 // LCOV_EXCL_START 707 return CeedError(qf->ceed, 1, "Backend does not support QFunctionApply"); 708 // LCOV_EXCL_STOP 709 if (Q % qf->vlength) 710 // LCOV_EXCL_START 711 return CeedError(qf->ceed, 2, "Number of quadrature points %d must be a " 712 "multiple of %d", Q, qf->vlength); 713 // LCOV_EXCL_STOP 714 ierr = qf->Apply(qf, Q, u, v); CeedChk(ierr); 715 return 0; 716 } 717 718 /** 719 @brief Destroy a CeedQFunction 720 721 @param qf CeedQFunction to destroy 722 723 @return An error code: 0 - success, otherwise - failure 724 725 @ref User 726 **/ 727 int CeedQFunctionDestroy(CeedQFunction *qf) { 728 int ierr; 729 730 if (!*qf || --(*qf)->refcount > 0) 731 return 0; 732 // Backend destroy 733 if ((*qf)->Destroy) { 734 ierr = (*qf)->Destroy(*qf); CeedChk(ierr); 735 } 736 // Free fields 737 for (int i=0; i<(*qf)->numinputfields; i++) { 738 ierr = CeedFree(&(*(*qf)->inputfields[i]).fieldname); CeedChk(ierr); 739 ierr = CeedFree(&(*qf)->inputfields[i]); CeedChk(ierr); 740 } 741 for (int i=0; i<(*qf)->numoutputfields; i++) { 742 ierr = CeedFree(&(*(*qf)->outputfields[i]).fieldname); CeedChk(ierr); 743 ierr = CeedFree(&(*qf)->outputfields[i]); CeedChk(ierr); 744 } 745 ierr = CeedFree(&(*qf)->inputfields); CeedChk(ierr); 746 ierr = CeedFree(&(*qf)->outputfields); CeedChk(ierr); 747 // Free ctx if identity 748 if ((*qf)->identity) { 749 ierr = CeedFree(&(*qf)->ctx); CeedChk(ierr); 750 } 751 752 ierr = CeedFree(&(*qf)->sourcepath); CeedChk(ierr); 753 ierr = CeedFree(&(*qf)->qfname); CeedChk(ierr); 754 ierr = CeedDestroy(&(*qf)->ceed); CeedChk(ierr); 755 ierr = CeedFree(qf); CeedChk(ierr); 756 return 0; 757 } 758 759 /// @} 760