1 // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 2 // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 3 // reserved. See files LICENSE and NOTICE for details. 4 // 5 // This file is part of CEED, a collection of benchmarks, miniapps, software 6 // libraries and APIs for efficient high-order finite element and spectral 7 // element discretizations for exascale applications. For more information and 8 // source code availability see http://github.com/ceed. 9 // 10 // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 11 // a collaborative effort of two U.S. Department of Energy organizations (Office 12 // of Science and the National Nuclear Security Administration) responsible for 13 // the planning and preparation of a capable exascale ecosystem, including 14 // software, applications, hardware, advanced system engineering and early 15 // testbed platforms, in support of the nation's exascale computing imperative. 16 17 #include <ceed-impl.h> 18 #include <ceed-backend.h> 19 #include <string.h> 20 #include <limits.h> 21 22 /// @cond DOXYGEN_SKIP 23 static struct { 24 char name[CEED_MAX_RESOURCE_LEN]; 25 char source[CEED_MAX_RESOURCE_LEN]; 26 CeedInt vlength; 27 CeedQFunctionUser f; 28 int (*init)(Ceed ceed, const char *name, CeedQFunction qf); 29 } qfunctions[1024]; 30 static size_t num_qfunctions; 31 /// @endcond 32 33 /// @file 34 /// Implementation of public CeedQFunction interfaces 35 /// 36 /// @addtogroup CeedQFunction 37 /// @{ 38 39 /** 40 @brief Create a CeedQFunction for evaluating interior (volumetric) terms. 41 42 @param ceed A Ceed object where the CeedQFunction will be created 43 @param vlength Vector length. Caller must ensure that number of quadrature 44 points is a multiple of vlength. 45 @param f Function pointer to evaluate action at quadrature points. 46 See \ref CeedQFunctionUser. 47 @param source Absolute path to source of QFunction, 48 "\abs_path\file.h:function_name" 49 @param[out] qf Address of the variable where the newly created 50 CeedQFunction will be stored 51 52 @return An error code: 0 - success, otherwise - failure 53 54 See \ref CeedQFunctionUser for details on the call-back function @a f's 55 arguments. 56 57 @ref Basic 58 **/ 59 int CeedQFunctionCreateInterior(Ceed ceed, CeedInt vlength, 60 CeedQFunctionUser f, 61 const char *source, CeedQFunction *qf) { 62 int ierr; 63 char *source_copy; 64 65 if (!ceed->QFunctionCreate) { 66 Ceed delegate; 67 ierr = CeedGetObjectDelegate(ceed, &delegate, "QFunction"); CeedChk(ierr); 68 69 if (!delegate) 70 // LCOV_EXCL_START 71 return CeedError(ceed, 1, "Backend does not support QFunctionCreate"); 72 // LCOV_EXCL_STOP 73 74 ierr = CeedQFunctionCreateInterior(delegate, vlength, f, source, qf); 75 CeedChk(ierr); 76 return 0; 77 } 78 79 ierr = CeedCalloc(1,qf); CeedChk(ierr); 80 (*qf)->ceed = ceed; 81 ceed->refcount++; 82 (*qf)->refcount = 1; 83 (*qf)->vlength = vlength; 84 (*qf)->identity = 0; 85 (*qf)->function = f; 86 size_t slen = strlen(source) + 1; 87 ierr = CeedMalloc(slen, &source_copy); CeedChk(ierr); 88 memcpy(source_copy, source, slen); 89 (*qf)->sourcepath = source_copy; 90 ierr = CeedCalloc(16, &(*qf)->inputfields); CeedChk(ierr); 91 ierr = CeedCalloc(16, &(*qf)->outputfields); CeedChk(ierr); 92 ierr = ceed->QFunctionCreate(*qf); CeedChk(ierr); 93 return 0; 94 } 95 96 /** 97 @brief Register a gallery QFunction 98 99 @param name Name for this backend to respond to 100 @param source Absolute path to source of QFunction, 101 "\path\CEED_DIR\gallery\folder\file.h:function_name" 102 @param vlength Vector length. Caller must ensure that number of quadrature 103 points is a multiple of vlength. 104 @param f Function pointer to evaluate action at quadrature points. 105 See \ref CeedQFunctionUser. 106 @param init Initialization function called by CeedQFunctionInit() when the 107 QFunction is selected. 108 109 @return An error code: 0 - success, otherwise - failure 110 111 @ref Advanced 112 **/ 113 int CeedQFunctionRegister(const char *name, const char *source, 114 CeedInt vlength, CeedQFunctionUser f, 115 int (*init)(Ceed, const char *, CeedQFunction)) { 116 if (num_qfunctions >= sizeof(qfunctions) / sizeof(qfunctions[0])) 117 // LCOV_EXCL_START 118 return CeedError(NULL, 1, "Too many gallery QFunctions"); 119 // LCOV_EXCL_STOP 120 121 strncpy(qfunctions[num_qfunctions].name, name, CEED_MAX_RESOURCE_LEN); 122 qfunctions[num_qfunctions].name[CEED_MAX_RESOURCE_LEN-1] = 0; 123 strncpy(qfunctions[num_qfunctions].source, source, CEED_MAX_RESOURCE_LEN); 124 qfunctions[num_qfunctions].source[CEED_MAX_RESOURCE_LEN-1] = 0; 125 qfunctions[num_qfunctions].vlength = vlength; 126 qfunctions[num_qfunctions].f = f; 127 qfunctions[num_qfunctions].init = init; 128 num_qfunctions++; 129 return 0; 130 } 131 132 /** 133 @brief Create a CeedQFunction for evaluating interior (volumetric) terms by name. 134 135 @param ceed A Ceed object where the CeedQFunction will be created 136 @param name Name of QFunction to use from gallery 137 @param[out] qf Address of the variable where the newly created 138 CeedQFunction will be stored 139 140 @return An error code: 0 - success, otherwise - failure 141 142 @ref Basic 143 **/ 144 int CeedQFunctionCreateInteriorByName(Ceed ceed, const char *name, 145 CeedQFunction *qf) { 146 int ierr; 147 size_t matchlen = 0, matchidx = UINT_MAX; 148 149 // Find matching backend 150 if (!name) return CeedError(NULL, 1, "No QFunction name provided"); 151 for (size_t i=0; i<num_qfunctions; i++) { 152 size_t n; 153 const char *currname = qfunctions[i].name; 154 for (n = 0; currname[n] && currname[n] == name[n]; n++) {} 155 if (n > matchlen) { 156 matchlen = n; 157 matchidx = i; 158 } 159 } 160 if (!matchlen) 161 // LCOV_EXCL_START 162 return CeedError(NULL, 1, "No suitable gallery QFunction"); 163 // LCOV_EXCL_STOP 164 165 // Create QFunction 166 ierr = CeedQFunctionCreateInterior(ceed, qfunctions[matchidx].vlength, 167 qfunctions[matchidx].f, 168 qfunctions[matchidx].source, qf); 169 CeedChk(ierr); 170 171 // QFunction specific setup 172 ierr = qfunctions[matchidx].init(ceed, name, *qf); CeedChk(ierr); 173 174 return 0; 175 } 176 177 /** 178 @brief Create an identity CeedQFunction. Inputs are written into outputs in 179 the order given. This is useful for CeedOperators that can be 180 represented with only the action of a CeedRestriction and CeedBasis, 181 such as restriction and prolongation operators for p-multigrid. 182 Backends may optimize CeedOperators with this CeedQFunction to avoid 183 the copy of input data to output fields by using the same memory 184 location for both. 185 186 @param ceed A Ceed object where the CeedQFunction will be created 187 @param size Size of the qfunction fields 188 @param[out] qf Address of the variable where the newly created 189 CeedQFunction will be stored 190 191 @return An error code: 0 - success, otherwise - failure 192 193 @ref Basic 194 **/ 195 int CeedQFunctionCreateIdentity(Ceed ceed, CeedInt size, CeedQFunction *qf) { 196 int ierr; 197 198 ierr = CeedQFunctionCreateInteriorByName(ceed, "Identity", qf); CeedChk(ierr); 199 200 (*qf)->identity = 1; 201 if (size > 1) { 202 CeedInt *ctx; 203 ierr = CeedCalloc(1, &ctx); CeedChk(ierr); 204 ctx[0] = size; 205 ierr = CeedQFunctionSetContext(*qf, ctx, sizeof(ctx)); CeedChk(ierr); 206 (*qf)->inputfields[0]->size = size; 207 (*qf)->outputfields[0]->size = size; 208 } 209 210 return 0; 211 } 212 213 /** 214 @brief Set a CeedQFunction field, used by CeedQFunctionAddInput/Output 215 216 @param f CeedQFunctionField 217 @param fieldname Name of QFunction field 218 @param size Size of QFunction field, ncomp * (dim for CEED_EVAL_GRAD or 219 1 for CEED_EVAL_NONE, CEED_EVAL_INTERP, and CEED_EVAL_WEIGHT) 220 @param emode \ref CEED_EVAL_NONE to use values directly, 221 \ref CEED_EVAL_INTERP to use interpolated values, 222 \ref CEED_EVAL_GRAD to use gradients, 223 \ref CEED_EVAL_WEIGHT to use quadrature weights. 224 225 @return An error code: 0 - success, otherwise - failure 226 227 @ref Developer 228 **/ 229 static int CeedQFunctionFieldSet(CeedQFunctionField *f,const char *fieldname, 230 CeedInt size, CeedEvalMode emode) { 231 size_t len = strlen(fieldname); 232 char *tmp; 233 int ierr; 234 ierr = CeedCalloc(1,f); CeedChk(ierr); 235 236 ierr = CeedCalloc(len+1, &tmp); CeedChk(ierr); 237 memcpy(tmp, fieldname, len+1); 238 (*f)->fieldname = tmp; 239 (*f)->size = size; 240 (*f)->emode = emode; 241 return 0; 242 } 243 244 /** 245 @brief Add a CeedQFunction input 246 247 @param qf CeedQFunction 248 @param fieldname Name of QFunction field 249 @param size Size of QFunction field, ncomp * (dim for CEED_EVAL_GRAD or 250 1 for CEED_EVAL_NONE and CEED_EVAL_INTERP) 251 @param emode \ref CEED_EVAL_NONE to use values directly, 252 \ref CEED_EVAL_INTERP to use interpolated values, 253 \ref CEED_EVAL_GRAD to use gradients. 254 255 @return An error code: 0 - success, otherwise - failure 256 257 @ref Basic 258 **/ 259 int CeedQFunctionAddInput(CeedQFunction qf, const char *fieldname, CeedInt size, 260 CeedEvalMode emode) { 261 int ierr = CeedQFunctionFieldSet(&qf->inputfields[qf->numinputfields], 262 fieldname, size, emode); 263 CeedChk(ierr); 264 qf->numinputfields++; 265 return 0; 266 } 267 268 /** 269 @brief Add a CeedQFunction output 270 271 @param qf CeedQFunction 272 @param fieldname Name of QFunction field 273 @param size Size of QFunction field, ncomp * (dim for CEED_EVAL_GRAD or 274 1 for CEED_EVAL_NONE and CEED_EVAL_INTERP) 275 @param emode \ref CEED_EVAL_NONE to use values directly, 276 \ref CEED_EVAL_INTERP to use interpolated values, 277 \ref CEED_EVAL_GRAD to use gradients. 278 279 @return An error code: 0 - success, otherwise - failure 280 281 @ref Basic 282 **/ 283 int CeedQFunctionAddOutput(CeedQFunction qf, const char *fieldname, 284 CeedInt size, CeedEvalMode emode) { 285 if (emode == CEED_EVAL_WEIGHT) 286 // LCOV_EXCL_START 287 return CeedError(qf->ceed, 1, "Cannot create QFunction output with " 288 "CEED_EVAL_WEIGHT"); 289 // LCOV_EXCL_STOP 290 int ierr = CeedQFunctionFieldSet(&qf->outputfields[qf->numoutputfields], 291 fieldname, size, emode); 292 CeedChk(ierr); 293 qf->numoutputfields++; 294 return 0; 295 } 296 297 /** 298 @brief Get the Ceed associated with a CeedQFunction 299 300 @param qf CeedQFunction 301 @param[out] ceed Variable to store Ceed 302 303 @return An error code: 0 - success, otherwise - failure 304 305 @ref Advanced 306 **/ 307 308 int CeedQFunctionGetCeed(CeedQFunction qf, Ceed *ceed) { 309 *ceed = qf->ceed; 310 return 0; 311 } 312 313 /** 314 @brief Get the vector length of a CeedQFunction 315 316 @param qf CeedQFunction 317 @param[out] vlength Variable to store vector length 318 319 @return An error code: 0 - success, otherwise - failure 320 321 @ref Advanced 322 **/ 323 324 int CeedQFunctionGetVectorLength(CeedQFunction qf, CeedInt *vlength) { 325 *vlength = qf->vlength; 326 return 0; 327 } 328 329 /** 330 @brief Get the number of inputs and outputs to a CeedQFunction 331 332 @param qf CeedQFunction 333 @param[out] numinput Variable to store number of input fields 334 @param[out] numoutput Variable to store number of output fields 335 336 @return An error code: 0 - success, otherwise - failure 337 338 @ref Advanced 339 **/ 340 341 int CeedQFunctionGetNumArgs(CeedQFunction qf, CeedInt *numinput, 342 CeedInt *numoutput) { 343 if (numinput) *numinput = qf->numinputfields; 344 if (numoutput) *numoutput = qf->numoutputfields; 345 return 0; 346 } 347 348 /** 349 @brief Get the source path string for a CeedQFunction 350 351 @param qf CeedQFunction 352 @param[out] source Variable to store source path string 353 354 @return An error code: 0 - success, otherwise - failure 355 356 @ref Advanced 357 **/ 358 359 int CeedQFunctionGetSourcePath(CeedQFunction qf, char* *source) { 360 *source = (char *) qf->sourcepath; 361 return 0; 362 } 363 364 /** 365 @brief Get the User Function for a CeedQFunction 366 367 @param qf CeedQFunction 368 @param[out] f Variable to store user function 369 370 @return An error code: 0 - success, otherwise - failure 371 372 @ref Advanced 373 **/ 374 375 int CeedQFunctionGetUserFunction(CeedQFunction qf, int (**f)()) { 376 *f = (int (*)())qf->function; 377 return 0; 378 } 379 380 /** 381 @brief Get global context size for a CeedQFunction 382 383 @param qf CeedQFunction 384 @param[out] ctxsize Variable to store size of context data values 385 386 @return An error code: 0 - success, otherwise - failure 387 388 @ref Advanced 389 **/ 390 391 int CeedQFunctionGetContextSize(CeedQFunction qf, size_t *ctxsize) { 392 if (qf->fortranstatus) { 393 fContext *fctx = qf->ctx; 394 *ctxsize = fctx->innerctxsize; 395 } else { 396 *ctxsize = qf->ctxsize; 397 } 398 return 0; 399 } 400 401 /** 402 @brief Get global context for a CeedQFunction 403 404 @param qf CeedQFunction 405 @param[out] ctx Variable to store context data values 406 407 @return An error code: 0 - success, otherwise - failure 408 409 @ref Advanced 410 **/ 411 412 int CeedQFunctionGetContext(CeedQFunction qf, void* *ctx) { 413 *ctx = qf->ctx; 414 return 0; 415 } 416 417 /** 418 @brief Determine if Fortran interface was used 419 420 @param qf CeedQFunction 421 @param[out] fortranstatus Variable to store Fortran status 422 423 @return An error code: 0 - success, otherwise - failure 424 425 @ref Advanced 426 **/ 427 428 int CeedQFunctionGetFortranStatus(CeedQFunction qf, bool *fortranstatus) { 429 *fortranstatus = qf->fortranstatus; 430 return 0; 431 } 432 433 /** 434 @brief Determine if QFunction is identity 435 436 @param qf CeedQFunction 437 @param[out] identity Variable to store identity status 438 439 @return An error code: 0 - success, otherwise - failure 440 441 @ref Advanced 442 **/ 443 444 int CeedQFunctionGetIdentityStatus(CeedQFunction qf, bool *identity) { 445 *identity = qf->identity; 446 return 0; 447 } 448 449 /** 450 @brief Get true user context for a CeedQFunction 451 452 @param qf CeedQFunction 453 @param[out] ctx Variable to store context data values 454 455 @return An error code: 0 - success, otherwise - failure 456 457 @ref Advanced 458 **/ 459 460 int CeedQFunctionGetInnerContext(CeedQFunction qf, void* *ctx) { 461 if (qf->fortranstatus) { 462 fContext *fctx = qf->ctx; 463 *ctx = fctx->innerctx; 464 } else { 465 *ctx = qf->ctx; 466 } 467 468 469 return 0; 470 } 471 472 /** 473 @brief Get backend data of a CeedQFunction 474 475 @param qf CeedQFunction 476 @param[out] data Variable to store data 477 478 @return An error code: 0 - success, otherwise - failure 479 480 @ref Advanced 481 **/ 482 483 int CeedQFunctionGetData(CeedQFunction qf, void* *data) { 484 *data = qf->data; 485 return 0; 486 } 487 488 /** 489 @brief Set backend data of a CeedQFunction 490 491 @param[out] qf CeedQFunction 492 @param data Data to set 493 494 @return An error code: 0 - success, otherwise - failure 495 496 @ref Advanced 497 **/ 498 499 int CeedQFunctionSetData(CeedQFunction qf, void* *data) { 500 qf->data = *data; 501 return 0; 502 } 503 504 /** 505 @brief Set global context for a CeedQFunction 506 507 @param qf CeedQFunction 508 @param ctx Context data to set 509 @param ctxsize Size of context data values 510 511 @return An error code: 0 - success, otherwise - failure 512 513 @ref Basic 514 **/ 515 int CeedQFunctionSetContext(CeedQFunction qf, void *ctx, size_t ctxsize) { 516 qf->ctx = ctx; 517 qf->ctxsize = ctxsize; 518 return 0; 519 } 520 521 /** 522 @brief Apply the action of a CeedQFunction 523 524 @param qf CeedQFunction 525 @param Q Number of quadrature points 526 @param[in] u Array of input data arrays 527 @param[out] v Array of output data arrays 528 529 @return An error code: 0 - success, otherwise - failure 530 531 @ref Advanced 532 **/ 533 int CeedQFunctionApply(CeedQFunction qf, CeedInt Q, 534 CeedVector *u, CeedVector *v) { 535 int ierr; 536 if (!qf->Apply) 537 // LCOV_EXCL_START 538 return CeedError(qf->ceed, 1, "Backend does not support QFunctionApply"); 539 // LCOV_EXCL_STOP 540 if (Q % qf->vlength) 541 // LCOV_EXCL_START 542 return CeedError(qf->ceed, 2, "Number of quadrature points %d must be a " 543 "multiple of %d", Q, qf->vlength); 544 // LCOV_EXCL_STOP 545 ierr = qf->Apply(qf, Q, u, v); CeedChk(ierr); 546 return 0; 547 } 548 549 /** 550 @brief Get the CeedQFunctionFields of a CeedQFunction 551 552 @param qf CeedQFunction 553 @param[out] inputfields Variable to store inputfields 554 @param[out] outputfields Variable to store outputfields 555 556 @return An error code: 0 - success, otherwise - failure 557 558 @ref Advanced 559 **/ 560 561 int CeedQFunctionGetFields(CeedQFunction qf, 562 CeedQFunctionField* *inputfields, 563 CeedQFunctionField* *outputfields) { 564 if (inputfields) 565 *inputfields = qf->inputfields; 566 if (outputfields) 567 *outputfields = qf->outputfields; 568 return 0; 569 } 570 571 /** 572 @brief Get the name of a CeedQFunctionField 573 574 @param qffield CeedQFunctionField 575 @param[out] fieldname Variable to store the field name 576 577 @return An error code: 0 - success, otherwise - failure 578 579 @ref Advanced 580 **/ 581 582 int CeedQFunctionFieldGetName(CeedQFunctionField qffield, 583 char* *fieldname) { 584 *fieldname = (char *)qffield->fieldname; 585 return 0; 586 } 587 588 /** 589 @brief Get the number of components of a CeedQFunctionField 590 591 @param qffield CeedQFunctionField 592 @param[out] size Variable to store the size of the field 593 594 @return An error code: 0 - success, otherwise - failure 595 596 @ref Advanced 597 **/ 598 599 int CeedQFunctionFieldGetSize(CeedQFunctionField qffield, CeedInt *size) { 600 *size = qffield->size; 601 return 0; 602 } 603 604 /** 605 @brief Get the CeedEvalMode of a CeedQFunctionField 606 607 @param qffield CeedQFunctionField 608 @param[out] emode Variable to store the field evaluation mode 609 610 @return An error code: 0 - success, otherwise - failure 611 612 @ref Advanced 613 **/ 614 615 int CeedQFunctionFieldGetEvalMode(CeedQFunctionField qffield, 616 CeedEvalMode *emode) { 617 *emode = qffield->emode; 618 return 0; 619 } 620 621 /** 622 @brief Destroy a CeedQFunction 623 624 @param qf CeedQFunction to destroy 625 626 @return An error code: 0 - success, otherwise - failure 627 628 @ref Basic 629 **/ 630 int CeedQFunctionDestroy(CeedQFunction *qf) { 631 int ierr; 632 633 if (!*qf || --(*qf)->refcount > 0) 634 return 0; 635 // Backend destroy 636 if ((*qf)->Destroy) { 637 ierr = (*qf)->Destroy(*qf); CeedChk(ierr); 638 } 639 // Free fields 640 for (int i=0; i<(*qf)->numinputfields; i++) { 641 ierr = CeedFree(&(*(*qf)->inputfields[i]).fieldname); CeedChk(ierr); 642 ierr = CeedFree(&(*qf)->inputfields[i]); CeedChk(ierr); 643 } 644 for (int i=0; i<(*qf)->numoutputfields; i++) { 645 ierr = CeedFree(&(*(*qf)->outputfields[i]).fieldname); CeedChk(ierr); 646 ierr = CeedFree(&(*qf)->outputfields[i]); CeedChk(ierr); 647 } 648 ierr = CeedFree(&(*qf)->inputfields); CeedChk(ierr); 649 ierr = CeedFree(&(*qf)->outputfields); CeedChk(ierr); 650 // Free ctx if identity 651 if ((*qf)->identity) { 652 ierr = CeedFree(&(*qf)->ctx); CeedChk(ierr); 653 } 654 655 ierr = CeedFree(&(*qf)->sourcepath); CeedChk(ierr); 656 ierr = CeedDestroy(&(*qf)->ceed); CeedChk(ierr); 657 ierr = CeedFree(qf); CeedChk(ierr); 658 return 0; 659 } 660 661 /// @} 662