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