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