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