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 `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 char *field_name; 129 CeedInt size; 130 CeedEvalMode eval_mode; 131 132 CeedCall(CeedQFunctionFieldGetName(field, &field_name)); 133 CeedCall(CeedQFunctionFieldGetSize(field, &size)); 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 216 CeedCall(CeedQFunctionGetCeed(qf, &ceed)); 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 CeedCall(CeedCheckFilePath(ceed, qf->user_source, &is_absolute_path)); 253 if (is_absolute_path) { 254 absolute_path = (char *)qf->user_source; 255 } else { 256 CeedCall(CeedGetJitAbsolutePath(ceed, qf->user_source, &absolute_path)); 257 } 258 259 size_t source_len = strlen(absolute_path) - kernel_name_len - 1; 260 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 `CeedQFunction` 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 `CeedQFunction` source file. 276 277 Note: Caller is responsible for freeing the string buffer with @ref 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 return CEED_ERROR_SUCCESS; 295 } 296 297 /** 298 @brief Get the User Function for a `CeedQFunction` 299 300 @param[in] qf `CeedQFunction` 301 @param[out] f Variable to store user function 302 303 @return An error code: 0 - success, otherwise - failure 304 305 @ref Backend 306 **/ 307 int CeedQFunctionGetUserFunction(CeedQFunction qf, CeedQFunctionUser *f) { 308 *f = qf->function; 309 return CEED_ERROR_SUCCESS; 310 } 311 312 /** 313 @brief Get global context for a `CeedQFunction`. 314 315 Note: For `CeedQFunction` from the Fortran interface, this function will return the Fortran context `CeedQFunctionContext`. 316 317 @param[in] qf CeedQFunction 318 @param[out] ctx Variable to store CeedQFunctionContext 319 320 @return An error code: 0 - success, otherwise - failure 321 322 @ref Backend 323 **/ 324 int CeedQFunctionGetContext(CeedQFunction qf, CeedQFunctionContext *ctx) { 325 *ctx = qf->ctx; 326 return CEED_ERROR_SUCCESS; 327 } 328 329 /** 330 @brief Get context data of a `CeedQFunction` 331 332 @param[in] qf `CeedQFunction` 333 @param[in] mem_type Memory type on which to access the data. 334 If the backend uses a different memory type, this will perform a copy. 335 @param[out] data Data on memory type mem_type 336 337 @return An error code: 0 - success, otherwise - failure 338 339 @ref Backend 340 **/ 341 int CeedQFunctionGetContextData(CeedQFunction qf, CeedMemType mem_type, void *data) { 342 bool is_writable; 343 CeedQFunctionContext ctx; 344 345 CeedCall(CeedQFunctionGetContext(qf, &ctx)); 346 if (ctx) { 347 CeedCall(CeedQFunctionIsContextWritable(qf, &is_writable)); 348 if (is_writable) { 349 CeedCall(CeedQFunctionContextGetData(ctx, mem_type, data)); 350 } else { 351 CeedCall(CeedQFunctionContextGetDataRead(ctx, mem_type, data)); 352 } 353 } else { 354 *(void **)data = NULL; 355 } 356 return CEED_ERROR_SUCCESS; 357 } 358 359 /** 360 @brief Restore context data of a `CeedQFunction` 361 362 @param[in] qf `CeedQFunction` 363 @param[in,out] data Data to restore 364 365 @return An error code: 0 - success, otherwise - failure 366 367 @ref Backend 368 **/ 369 int CeedQFunctionRestoreContextData(CeedQFunction qf, void *data) { 370 bool is_writable; 371 CeedQFunctionContext ctx; 372 373 CeedCall(CeedQFunctionGetContext(qf, &ctx)); 374 if (ctx) { 375 CeedCall(CeedQFunctionIsContextWritable(qf, &is_writable)); 376 if (is_writable) { 377 CeedCall(CeedQFunctionContextRestoreData(ctx, data)); 378 } else { 379 CeedCall(CeedQFunctionContextRestoreDataRead(ctx, data)); 380 } 381 } 382 *(void **)data = NULL; 383 return CEED_ERROR_SUCCESS; 384 } 385 386 /** 387 @brief Get true user context for a `CeedQFunction` 388 389 Note: For all `CeedQFunction` this function will return the user `CeedQFunctionContext` and not interface context `CeedQFunctionContext`, if any such object exists. 390 391 @param[in] qf `CeedQFunction` 392 @param[out] ctx Variable to store `CeedQFunctionContext` 393 394 @return An error code: 0 - success, otherwise - failure 395 @ref Backend 396 **/ 397 int CeedQFunctionGetInnerContext(CeedQFunction qf, CeedQFunctionContext *ctx) { 398 if (qf->is_fortran) { 399 CeedFortranContext fortran_ctx = NULL; 400 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 *(void **)data = NULL; 464 return CEED_ERROR_SUCCESS; 465 } 466 467 /** 468 @brief Determine if `CeedQFunction` is identity 469 470 @param[in] qf `CeedQFunction` 471 @param[out] is_identity Variable to store identity status 472 473 @return An error code: 0 - success, otherwise - failure 474 475 @ref Backend 476 **/ 477 int CeedQFunctionIsIdentity(CeedQFunction qf, bool *is_identity) { 478 *is_identity = qf->is_identity; 479 return CEED_ERROR_SUCCESS; 480 } 481 482 /** 483 @brief Determine if `CeedQFunctionContext` is writable 484 485 @param[in] qf `CeedQFunction` 486 @param[out] is_writable Variable to store context writeable status 487 488 @return An error code: 0 - success, otherwise - failure 489 490 @ref Backend 491 **/ 492 int CeedQFunctionIsContextWritable(CeedQFunction qf, bool *is_writable) { 493 *is_writable = qf->is_context_writable; 494 return CEED_ERROR_SUCCESS; 495 } 496 497 /** 498 @brief Get backend data of a `CeedQFunction` 499 500 @param[in] qf `CeedQFunction` 501 @param[out] data Variable to store data 502 503 @return An error code: 0 - success, otherwise - failure 504 505 @ref Backend 506 **/ 507 int CeedQFunctionGetData(CeedQFunction qf, void *data) { 508 *(void **)data = qf->data; 509 return CEED_ERROR_SUCCESS; 510 } 511 512 /** 513 @brief Set backend data of a `CeedQFunction` 514 515 @param[in,out] qf `CeedQFunction` 516 @param[in] data Data to set 517 518 @return An error code: 0 - success, otherwise - failure 519 520 @ref Backend 521 **/ 522 int CeedQFunctionSetData(CeedQFunction qf, void *data) { 523 qf->data = data; 524 return CEED_ERROR_SUCCESS; 525 } 526 527 /** 528 @brief Increment the reference counter for a `CeedQFunction` 529 530 @param[in,out] qf `CeedQFunction` to increment the reference counter 531 532 @return An error code: 0 - success, otherwise - failure 533 534 @ref Backend 535 **/ 536 int CeedQFunctionReference(CeedQFunction qf) { 537 qf->ref_count++; 538 return CEED_ERROR_SUCCESS; 539 } 540 541 /** 542 @brief Estimate number of FLOPs per quadrature required to apply `CeedQFunction` 543 544 @param[in] qf `CeedQFunction` to estimate FLOPs for 545 @param[out] flops Address of variable to hold FLOPs estimate 546 547 @ref Backend 548 **/ 549 int CeedQFunctionGetFlopsEstimate(CeedQFunction qf, CeedSize *flops) { 550 CeedCheck(qf->user_flop_estimate > -1, qf->ceed, CEED_ERROR_INCOMPLETE, "Must set FLOPs estimate with CeedQFunctionSetUserFlopsEstimate"); 551 *flops = qf->user_flop_estimate; 552 return CEED_ERROR_SUCCESS; 553 } 554 555 /// @} 556 557 /// ---------------------------------------------------------------------------- 558 /// CeedQFunction Public API 559 /// ---------------------------------------------------------------------------- 560 /// @addtogroup CeedQFunctionUser 561 /// @{ 562 563 /** 564 @brief Create a `CeedQFunction` for evaluating interior (volumetric) terms 565 566 @param[in] ceed `Ceed` object used to create the `CeedQFunction` 567 @param[in] vec_length Vector length. 568 Caller must ensure that number of quadrature points is a multiple of `vec_length`. 569 @param[in] f Function pointer to evaluate action at quadrature points. 570 See `CeedQFunctionUser`. 571 @param[in] source Absolute path to source of `CeedQFunctionUser`, "\abs_path\file.h:function_name". 572 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.). 573 The entire contents of this file and all locally included files are used during JiT compilation for GPU backends. 574 All source files must be at the provided filepath at runtime for JiT to function. 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 `f` 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 589 CeedCall(CeedGetObjectDelegate(ceed, &delegate, "QFunction")); 590 CeedCheck(delegate, ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support CeedQFunctionCreateInterior"); 591 CeedCall(CeedQFunctionCreateInterior(delegate, vec_length, f, source, qf)); 592 return CEED_ERROR_SUCCESS; 593 } 594 595 CeedCheck(!strlen(source) || strrchr(source, ':'), ceed, CEED_ERROR_INCOMPLETE, 596 "Provided path to source does not include function name. Provided: \"%s\"\nRequired: \"\\abs_path\\file.h:function_name\"", source); 597 598 CeedCall(CeedCalloc(1, qf)); 599 CeedCall(CeedReferenceCopy(ceed, &(*qf)->ceed)); 600 (*qf)->ref_count = 1; 601 (*qf)->vec_length = vec_length; 602 (*qf)->is_identity = false; 603 (*qf)->is_context_writable = true; 604 (*qf)->function = f; 605 (*qf)->user_flop_estimate = -1; 606 if (strlen(source)) { 607 size_t user_source_len = strlen(source); 608 609 CeedCall(CeedCalloc(user_source_len + 1, &user_source_copy)); 610 memcpy(user_source_copy, source, user_source_len); 611 (*qf)->user_source = user_source_copy; 612 } 613 CeedCall(CeedCalloc(CEED_FIELD_MAX, &(*qf)->input_fields)); 614 CeedCall(CeedCalloc(CEED_FIELD_MAX, &(*qf)->output_fields)); 615 CeedCall(ceed->QFunctionCreate(*qf)); 616 return CEED_ERROR_SUCCESS; 617 } 618 619 /** 620 @brief Create a `CeedQFunction` for evaluating interior (volumetric) terms by name 621 622 @param[in] ceed `Ceed` object used to create the `CeedQFunction` 623 @param[in] name Name of `CeedQFunction` to use from gallery 624 @param[out] qf Address of the variable where the newly created `CeedQFunction` will be stored 625 626 @return An error code: 0 - success, otherwise - failure 627 628 @ref User 629 **/ 630 int CeedQFunctionCreateInteriorByName(Ceed ceed, const char *name, CeedQFunction *qf) { 631 size_t match_len = 0, match_index = UINT_MAX; 632 633 CeedCall(CeedQFunctionRegisterAll()); 634 // Find matching backend 635 CeedCheck(name, ceed, CEED_ERROR_INCOMPLETE, "No CeedQFunction name provided"); 636 for (size_t i = 0; i < num_qfunctions; i++) { 637 size_t n; 638 const char *curr_name = gallery_qfunctions[i].name; 639 for (n = 0; curr_name[n] && curr_name[n] == name[n]; n++) { 640 } 641 if (n > match_len) { 642 match_len = n; 643 match_index = i; 644 } 645 } 646 CeedCheck(match_len > 0, ceed, CEED_ERROR_UNSUPPORTED, "No suitable gallery CeedQFunction"); 647 648 // Create QFunction 649 CeedCall(CeedQFunctionCreateInterior(ceed, gallery_qfunctions[match_index].vec_length, gallery_qfunctions[match_index].f, 650 gallery_qfunctions[match_index].source, qf)); 651 652 // QFunction specific setup 653 CeedCall(gallery_qfunctions[match_index].init(ceed, name, *qf)); 654 655 // Copy name 656 CeedCall(CeedStringAllocCopy(name, (char **)&(*qf)->gallery_name)); 657 (*qf)->is_gallery = true; 658 return CEED_ERROR_SUCCESS; 659 } 660 661 /** 662 @brief Create an identity `CeedQFunction`. 663 664 Inputs are written into outputs in the order given. 665 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. 666 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. 667 668 @param[in] ceed `Ceed` object used to create the `CeedQFunction` 669 @param[in] size Size of the `CeedQFunction` fields 670 @param[in] in_mode @ref CeedEvalMode for input to `CeedQFunction` 671 @param[in] out_mode @ref CeedEvalMode for output to `CeedQFunction` 672 @param[out] qf Address of the variable where the newly created `CeedQFunction` will be stored 673 674 @return An error code: 0 - success, otherwise - failure 675 676 @ref User 677 **/ 678 int CeedQFunctionCreateIdentity(Ceed ceed, CeedInt size, CeedEvalMode in_mode, CeedEvalMode out_mode, CeedQFunction *qf) { 679 CeedQFunctionContext ctx; 680 CeedContextFieldLabel size_label; 681 682 CeedCall(CeedQFunctionCreateInteriorByName(ceed, "Identity", qf)); 683 CeedCall(CeedQFunctionAddInput(*qf, "input", size, in_mode)); 684 CeedCall(CeedQFunctionAddOutput(*qf, "output", size, out_mode)); 685 686 (*qf)->is_identity = true; 687 688 CeedCall(CeedQFunctionGetContext(*qf, &ctx)); 689 CeedCall(CeedQFunctionContextGetFieldLabel(ctx, "size", &size_label)); 690 CeedCall(CeedQFunctionContextSetInt32(ctx, size_label, &size)); 691 return CEED_ERROR_SUCCESS; 692 } 693 694 /** 695 @brief Copy the pointer to a `CeedQFunction`. 696 697 Both pointers should be destroyed with @ref CeedQFunctionDestroy(). 698 699 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`. 700 This `CeedQFunction` will be destroyed if `*qf_copy` is the only reference to this `CeedQFunction`. 701 702 @param[in] qf `CeedQFunction` to copy reference to 703 @param[out] qf_copy Variable to store copied reference 704 705 @return An error code: 0 - success, otherwise - failure 706 707 @ref User 708 **/ 709 int CeedQFunctionReferenceCopy(CeedQFunction qf, CeedQFunction *qf_copy) { 710 CeedCall(CeedQFunctionReference(qf)); 711 CeedCall(CeedQFunctionDestroy(qf_copy)); 712 *qf_copy = qf; 713 return CEED_ERROR_SUCCESS; 714 } 715 716 /** 717 @brief Add a `CeedQFunction` input 718 719 @param[in,out] qf `CeedQFunction` 720 @param[in] field_name Name of `CeedQFunction` field 721 @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. 722 @param[in] eval_mode @ref CEED_EVAL_NONE to use values directly, 723 @ref CEED_EVAL_INTERP to use interpolated values, 724 @ref CEED_EVAL_GRAD to use gradients, 725 @ref CEED_EVAL_DIV to use divergence, 726 @ref CEED_EVAL_CURL to use curl 727 728 @return An error code: 0 - success, otherwise - failure 729 730 @ref User 731 **/ 732 int CeedQFunctionAddInput(CeedQFunction qf, const char *field_name, CeedInt size, CeedEvalMode eval_mode) { 733 CeedCheck(!qf->is_immutable, qf->ceed, CEED_ERROR_MAJOR, "QFunction cannot be changed after set as immutable"); 734 CeedCheck(eval_mode != CEED_EVAL_WEIGHT || size == 1, qf->ceed, CEED_ERROR_DIMENSION, "CEED_EVAL_WEIGHT should have size 1"); 735 for (CeedInt i = 0; i < qf->num_input_fields; i++) { 736 CeedCheck(strcmp(field_name, qf->input_fields[i]->field_name), qf->ceed, CEED_ERROR_MINOR, "QFunction field names must be unique"); 737 } 738 for (CeedInt i = 0; i < qf->num_output_fields; i++) { 739 CeedCheck(strcmp(field_name, qf->output_fields[i]->field_name), qf->ceed, CEED_ERROR_MINOR, "QFunction field names must be unique"); 740 } 741 CeedCall(CeedQFunctionFieldSet(&qf->input_fields[qf->num_input_fields], field_name, size, eval_mode)); 742 qf->num_input_fields++; 743 return CEED_ERROR_SUCCESS; 744 } 745 746 /** 747 @brief Add a `CeedQFunction` output 748 749 @param[in,out] qf `CeedQFunction` 750 @param[in] field_name Name of `CeedQFunction` field 751 @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. 752 @param[in] eval_mode @ref CEED_EVAL_NONE to use values directly, 753 @ref CEED_EVAL_INTERP to use interpolated values, 754 @ref CEED_EVAL_GRAD to use gradients, 755 @ref CEED_EVAL_DIV to use divergence, 756 @ref CEED_EVAL_CURL to use curl. 757 758 @return An error code: 0 - success, otherwise - failure 759 760 @ref User 761 **/ 762 int CeedQFunctionAddOutput(CeedQFunction qf, const char *field_name, CeedInt size, CeedEvalMode eval_mode) { 763 CeedCheck(!qf->is_immutable, qf->ceed, CEED_ERROR_MAJOR, "CeedQFunction cannot be changed after set as immutable"); 764 CeedCheck(eval_mode != CEED_EVAL_WEIGHT, qf->ceed, CEED_ERROR_DIMENSION, "Cannot create CeedQFunction output with CEED_EVAL_WEIGHT"); 765 for (CeedInt i = 0; i < qf->num_input_fields; i++) { 766 CeedCheck(strcmp(field_name, qf->input_fields[i]->field_name), qf->ceed, CEED_ERROR_MINOR, "CeedQFunction field names must be unique"); 767 } 768 for (CeedInt i = 0; i < qf->num_output_fields; i++) { 769 CeedCheck(strcmp(field_name, qf->output_fields[i]->field_name), qf->ceed, CEED_ERROR_MINOR, "CeedQFunction field names must be unique"); 770 } 771 CeedCall(CeedQFunctionFieldSet(&qf->output_fields[qf->num_output_fields], field_name, size, eval_mode)); 772 qf->num_output_fields++; 773 return CEED_ERROR_SUCCESS; 774 } 775 776 /** 777 @brief Get the `CeedQFunctionField` of a `CeedQFunction` 778 779 Note: Calling this function asserts that setup is complete and sets the `CeedQFunction` as immutable. 780 781 @param[in] qf `CeedQFunction` 782 @param[out] num_input_fields Variable to store number of input fields 783 @param[out] input_fields Variable to store input fields 784 @param[out] num_output_fields Variable to store number of output fields 785 @param[out] output_fields Variable to store output fields 786 787 @return An error code: 0 - success, otherwise - failure 788 789 @ref Advanced 790 **/ 791 int CeedQFunctionGetFields(CeedQFunction qf, CeedInt *num_input_fields, CeedQFunctionField **input_fields, CeedInt *num_output_fields, 792 CeedQFunctionField **output_fields) { 793 qf->is_immutable = true; 794 if (num_input_fields) *num_input_fields = qf->num_input_fields; 795 if (input_fields) *input_fields = qf->input_fields; 796 if (num_output_fields) *num_output_fields = qf->num_output_fields; 797 if (output_fields) *output_fields = qf->output_fields; 798 return CEED_ERROR_SUCCESS; 799 } 800 801 /** 802 @brief Get the name of a `CeedQFunctionField` 803 804 @param[in] qf_field `CeedQFunctionField` 805 @param[out] field_name Variable to store the field name 806 807 @return An error code: 0 - success, otherwise - failure 808 809 @ref Advanced 810 **/ 811 int CeedQFunctionFieldGetName(CeedQFunctionField qf_field, char **field_name) { 812 *field_name = (char *)qf_field->field_name; 813 return CEED_ERROR_SUCCESS; 814 } 815 816 /** 817 @brief Get the number of components of a `CeedQFunctionField` 818 819 @param[in] qf_field `CeedQFunctionField` 820 @param[out] size Variable to store the size of the field 821 822 @return An error code: 0 - success, otherwise - failure 823 824 @ref Advanced 825 **/ 826 int CeedQFunctionFieldGetSize(CeedQFunctionField qf_field, CeedInt *size) { 827 *size = qf_field->size; 828 return CEED_ERROR_SUCCESS; 829 } 830 831 /** 832 @brief Get the @ref CeedEvalMode of a `CeedQFunctionField` 833 834 @param[in] qf_field `CeedQFunctionField` 835 @param[out] eval_mode Variable to store the field evaluation mode 836 837 @return An error code: 0 - success, otherwise - failure 838 839 @ref Advanced 840 **/ 841 int CeedQFunctionFieldGetEvalMode(CeedQFunctionField qf_field, CeedEvalMode *eval_mode) { 842 *eval_mode = qf_field->eval_mode; 843 return CEED_ERROR_SUCCESS; 844 } 845 846 /** 847 @brief Set global context for a `CeedQFunction` 848 849 @param[in,out] qf `CeedQFunction` 850 @param[in] ctx Context data to set 851 852 @return An error code: 0 - success, otherwise - failure 853 854 @ref User 855 **/ 856 int CeedQFunctionSetContext(CeedQFunction qf, CeedQFunctionContext ctx) { 857 CeedCall(CeedQFunctionContextDestroy(&qf->ctx)); 858 qf->ctx = ctx; 859 if (ctx) CeedCall(CeedQFunctionContextReference(ctx)); 860 return CEED_ERROR_SUCCESS; 861 } 862 863 /** 864 @brief Set writability of `CeedQFunctionContext` when calling the `CeedQFunctionUser`. 865 866 The default value is `is_writable == true`. 867 868 Setting `is_writable == true` indicates the `CeedQFunctionUser` writes into the `CeedQFunctionContext` and requires memory synchronization after calling @ref CeedQFunctionApply(). 869 870 Setting 'is_writable == false' asserts that `CeedQFunctionUser` does not modify the `CeedQFunctionContext`. 871 Violating this assertion may lead to inconsistent data. 872 873 Setting `is_writable == false` may offer a performance improvement on GPU backends. 874 875 @param[in,out] qf `CeedQFunction` 876 @param[in] is_writable Boolean flag for writability status 877 878 @return An error code: 0 - success, otherwise - failure 879 880 @ref User 881 **/ 882 int CeedQFunctionSetContextWritable(CeedQFunction qf, bool is_writable) { 883 qf->is_context_writable = is_writable; 884 return CEED_ERROR_SUCCESS; 885 } 886 887 /** 888 @brief Set estimated number of FLOPs per quadrature required to apply `CeedQFunction` 889 890 @param[in] qf `CeedQFunction` to estimate FLOPs for 891 @param[out] flops FLOPs per quadrature point estimate 892 893 @ref Backend 894 **/ 895 int CeedQFunctionSetUserFlopsEstimate(CeedQFunction qf, CeedSize flops) { 896 CeedCheck(flops >= 0, qf->ceed, CEED_ERROR_INCOMPATIBLE, "Must set non-negative FLOPs estimate"); 897 qf->user_flop_estimate = flops; 898 return CEED_ERROR_SUCCESS; 899 } 900 901 /** 902 @brief View a `CeedQFunction` 903 904 @param[in] qf `CeedQFunction` to view 905 @param[in] stream Stream to write; typically `stdout` or a file 906 907 @return Error code: 0 - success, otherwise - failure 908 909 @ref User 910 **/ 911 int CeedQFunctionView(CeedQFunction qf, FILE *stream) { 912 char *kernel_name; 913 914 CeedCall(CeedQFunctionGetKernelName(qf, &kernel_name)); 915 fprintf(stream, "%sCeedQFunction - %s\n", qf->is_gallery ? "Gallery " : "User ", qf->is_gallery ? qf->gallery_name : kernel_name); 916 917 fprintf(stream, " %" CeedInt_FMT " input field%s:\n", qf->num_input_fields, qf->num_input_fields > 1 ? "s" : ""); 918 for (CeedInt i = 0; i < qf->num_input_fields; i++) { 919 CeedCall(CeedQFunctionFieldView(qf->input_fields[i], i, 1, stream)); 920 } 921 922 fprintf(stream, " %" CeedInt_FMT " output field%s:\n", qf->num_output_fields, qf->num_output_fields > 1 ? "s" : ""); 923 for (CeedInt i = 0; i < qf->num_output_fields; i++) { 924 CeedCall(CeedQFunctionFieldView(qf->output_fields[i], i, 0, stream)); 925 } 926 return CEED_ERROR_SUCCESS; 927 } 928 929 /** 930 @brief Get the `Ceed` associated with a `CeedQFunction` 931 932 @param[in] qf `CeedQFunction` 933 @param[out] ceed Variable to store`Ceed` 934 935 @return An error code: 0 - success, otherwise - failure 936 937 @ref Advanced 938 **/ 939 int CeedQFunctionGetCeed(CeedQFunction qf, Ceed *ceed) { 940 *ceed = qf->ceed; 941 return CEED_ERROR_SUCCESS; 942 } 943 944 /** 945 @brief Apply the action of a `CeedQFunction` 946 947 Note: Calling this function asserts that setup is complete and sets the `CeedQFunction` as immutable. 948 949 @param[in] qf `CeedQFunction` 950 @param[in] Q Number of quadrature points 951 @param[in] u Array of input `CeedVector` 952 @param[out] v Array of output `CeedVector` 953 954 @return An error code: 0 - success, otherwise - failure 955 956 @ref User 957 **/ 958 int CeedQFunctionApply(CeedQFunction qf, CeedInt Q, CeedVector *u, CeedVector *v) { 959 CeedCheck(qf->Apply, qf->ceed, CEED_ERROR_UNSUPPORTED, "Backend does not support CeedQFunctionApply"); 960 CeedCheck(Q % qf->vec_length == 0, qf->ceed, CEED_ERROR_DIMENSION, 961 "Number of quadrature points %" CeedInt_FMT " must be a multiple of %" CeedInt_FMT, Q, qf->vec_length); 962 qf->is_immutable = true; 963 CeedCall(qf->Apply(qf, Q, u, v)); 964 return CEED_ERROR_SUCCESS; 965 } 966 967 /** 968 @brief Destroy a `CeedQFunction` 969 970 @param[in,out] qf `CeedQFunction` to destroy 971 972 @return An error code: 0 - success, otherwise - failure 973 974 @ref User 975 **/ 976 int CeedQFunctionDestroy(CeedQFunction *qf) { 977 if (!*qf || --(*qf)->ref_count > 0) { 978 *qf = NULL; 979 return CEED_ERROR_SUCCESS; 980 } 981 // Backend destroy 982 if ((*qf)->Destroy) { 983 CeedCall((*qf)->Destroy(*qf)); 984 } 985 // Free fields 986 for (CeedInt i = 0; i < (*qf)->num_input_fields; i++) { 987 CeedCall(CeedFree(&(*(*qf)->input_fields[i]).field_name)); 988 CeedCall(CeedFree(&(*qf)->input_fields[i])); 989 } 990 for (CeedInt i = 0; i < (*qf)->num_output_fields; i++) { 991 CeedCall(CeedFree(&(*(*qf)->output_fields[i]).field_name)); 992 CeedCall(CeedFree(&(*qf)->output_fields[i])); 993 } 994 CeedCall(CeedFree(&(*qf)->input_fields)); 995 CeedCall(CeedFree(&(*qf)->output_fields)); 996 997 // User context data object 998 CeedCall(CeedQFunctionContextDestroy(&(*qf)->ctx)); 999 1000 CeedCall(CeedFree(&(*qf)->user_source)); 1001 CeedCall(CeedFree(&(*qf)->source_path)); 1002 CeedCall(CeedFree(&(*qf)->gallery_name)); 1003 CeedCall(CeedFree(&(*qf)->kernel_name)); 1004 CeedCall(CeedDestroy(&(*qf)->ceed)); 1005 CeedCall(CeedFree(qf)); 1006 return CEED_ERROR_SUCCESS; 1007 } 1008 1009 /// @} 1010