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