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 return CeedError(NULL, 1, "No suitable gallery QFunction"); 162 163 // Create QFunction 164 ierr = CeedQFunctionCreateInterior(ceed, qfunctions[matchidx].vlength, 165 qfunctions[matchidx].f, 166 qfunctions[matchidx].source, qf); 167 CeedChk(ierr); 168 169 // QFunction specific setup 170 ierr = qfunctions[matchidx].init(ceed, name, *qf); CeedChk(ierr); 171 172 return 0; 173 } 174 175 /** 176 @brief Create an identity CeedQFunction. Inputs are written into outputs in 177 the order given. This is useful for CeedOperators that can be 178 represented with only the action of a CeedRestriction and CeedBasis, 179 such as restriction and prolongation operators for p-multigrid. 180 Backends may optimize CeedOperators with this CeedQFunction to avoid 181 the copy of input data to output fields by using the same memory 182 location for both. 183 184 @param ceed A Ceed object where the CeedQFunction will be created 185 @param size Size of the qfunction fields 186 @param[out] qf Address of the variable where the newly created 187 CeedQFunction will be stored 188 189 @return An error code: 0 - success, otherwise - failure 190 191 @ref Basic 192 **/ 193 int CeedQFunctionCreateIdentity(Ceed ceed, CeedInt size, CeedQFunction *qf) { 194 int ierr; 195 196 ierr = CeedQFunctionCreateInteriorByName(ceed, "Identity", qf); CeedChk(ierr); 197 198 (*qf)->identity = 1; 199 if (size > 1) { 200 CeedInt *ctx; 201 ierr = CeedCalloc(1, &ctx); CeedChk(ierr); 202 ctx[0] = size; 203 ierr = CeedQFunctionSetContext(*qf, ctx, sizeof(ctx)); CeedChk(ierr); 204 (*qf)->inputfields[0]->size = size; 205 (*qf)->outputfields[0]->size = size; 206 } 207 208 return 0; 209 } 210 211 /** 212 @brief Set a CeedQFunction field, used by CeedQFunctionAddInput/Output 213 214 @param f CeedQFunctionField 215 @param fieldname Name of QFunction field 216 @param size Size of QFunction field, ncomp * (dim for CEED_EVAL_GRAD or 217 1 for CEED_EVAL_NONE and CEED_EVAL_INTERP) 218 @param emode \ref CEED_EVAL_NONE to use values directly, 219 \ref CEED_EVAL_INTERP to use interpolated values, 220 \ref CEED_EVAL_GRAD to use gradients. 221 222 @return An error code: 0 - success, otherwise - failure 223 224 @ref Developer 225 **/ 226 static int CeedQFunctionFieldSet(CeedQFunctionField *f,const char *fieldname, 227 CeedInt size, CeedEvalMode emode) { 228 size_t len = strlen(fieldname); 229 char *tmp; 230 int ierr; 231 ierr = CeedCalloc(1,f); CeedChk(ierr); 232 233 ierr = CeedCalloc(len+1, &tmp); CeedChk(ierr); 234 memcpy(tmp, fieldname, len+1); 235 (*f)->fieldname = tmp; 236 (*f)->size = size; 237 (*f)->emode = emode; 238 return 0; 239 } 240 241 /** 242 @brief Add a CeedQFunction input 243 244 @param qf CeedQFunction 245 @param fieldname Name of QFunction field 246 @param size Size of QFunction field, ncomp * (dim for CEED_EVAL_GRAD or 247 1 for CEED_EVAL_NONE and CEED_EVAL_INTERP) 248 @param emode \ref CEED_EVAL_NONE to use values directly, 249 \ref CEED_EVAL_INTERP to use interpolated values, 250 \ref CEED_EVAL_GRAD to use gradients. 251 252 @return An error code: 0 - success, otherwise - failure 253 254 @ref Basic 255 **/ 256 int CeedQFunctionAddInput(CeedQFunction qf, const char *fieldname, CeedInt size, 257 CeedEvalMode emode) { 258 int ierr = CeedQFunctionFieldSet(&qf->inputfields[qf->numinputfields], 259 fieldname, size, emode); 260 CeedChk(ierr); 261 qf->numinputfields++; 262 return 0; 263 } 264 265 /** 266 @brief Add a CeedQFunction output 267 268 @param qf CeedQFunction 269 @param fieldname Name of QFunction field 270 @param size Size of QFunction field, ncomp * (dim for CEED_EVAL_GRAD or 271 1 for CEED_EVAL_NONE and CEED_EVAL_INTERP) 272 @param emode \ref CEED_EVAL_NONE to use values directly, 273 \ref CEED_EVAL_INTERP to use interpolated values, 274 \ref CEED_EVAL_GRAD to use gradients. 275 276 @return An error code: 0 - success, otherwise - failure 277 278 @ref Basic 279 **/ 280 int CeedQFunctionAddOutput(CeedQFunction qf, const char *fieldname, 281 CeedInt size, CeedEvalMode emode) { 282 if (emode == CEED_EVAL_WEIGHT) 283 // LCOV_EXCL_START 284 return CeedError(qf->ceed, 1, "Cannot create QFunction output with " 285 "CEED_EVAL_WEIGHT"); 286 // LCOV_EXCL_STOP 287 int ierr = CeedQFunctionFieldSet(&qf->outputfields[qf->numoutputfields], 288 fieldname, size, emode); 289 CeedChk(ierr); 290 qf->numoutputfields++; 291 return 0; 292 } 293 294 /** 295 @brief Get the Ceed associated with a CeedQFunction 296 297 @param qf CeedQFunction 298 @param[out] ceed Variable to store Ceed 299 300 @return An error code: 0 - success, otherwise - failure 301 302 @ref Advanced 303 **/ 304 305 int CeedQFunctionGetCeed(CeedQFunction qf, Ceed *ceed) { 306 *ceed = qf->ceed; 307 return 0; 308 } 309 310 /** 311 @brief Get the vector length of a CeedQFunction 312 313 @param qf CeedQFunction 314 @param[out] vlength Variable to store vector length 315 316 @return An error code: 0 - success, otherwise - failure 317 318 @ref Advanced 319 **/ 320 321 int CeedQFunctionGetVectorLength(CeedQFunction qf, CeedInt *vlength) { 322 *vlength = qf->vlength; 323 return 0; 324 } 325 326 /** 327 @brief Get the number of inputs and outputs to a CeedQFunction 328 329 @param qf CeedQFunction 330 @param[out] numinput Variable to store number of input fields 331 @param[out] numoutput Variable to store number of output fields 332 333 @return An error code: 0 - success, otherwise - failure 334 335 @ref Advanced 336 **/ 337 338 int CeedQFunctionGetNumArgs(CeedQFunction qf, CeedInt *numinput, 339 CeedInt *numoutput) { 340 if (numinput) *numinput = qf->numinputfields; 341 if (numoutput) *numoutput = qf->numoutputfields; 342 return 0; 343 } 344 345 /** 346 @brief Get the source path string for a CeedQFunction 347 348 @param qf CeedQFunction 349 @param[out] source Variable to store source path string 350 351 @return An error code: 0 - success, otherwise - failure 352 353 @ref Advanced 354 **/ 355 356 int CeedQFunctionGetSourcePath(CeedQFunction qf, char* *source) { 357 *source = (char *) qf->sourcepath; 358 return 0; 359 } 360 361 /** 362 @brief Get the User Function for a CeedQFunction 363 364 @param qf CeedQFunction 365 @param[out] f Variable to store user function 366 367 @return An error code: 0 - success, otherwise - failure 368 369 @ref Advanced 370 **/ 371 372 int CeedQFunctionGetUserFunction(CeedQFunction qf, int (**f)()) { 373 *f = (int (*)())qf->function; 374 return 0; 375 } 376 377 /** 378 @brief Get global context size for a CeedQFunction 379 380 @param qf CeedQFunction 381 @param[out] ctxsize Variable to store size of context data values 382 383 @return An error code: 0 - success, otherwise - failure 384 385 @ref Advanced 386 **/ 387 388 int CeedQFunctionGetContextSize(CeedQFunction qf, size_t *ctxsize) { 389 if (qf->fortranstatus) { 390 fContext *fctx = qf->ctx; 391 *ctxsize = fctx->innerctxsize; 392 } else { 393 *ctxsize = qf->ctxsize; 394 } 395 return 0; 396 } 397 398 /** 399 @brief Get global context for a CeedQFunction 400 401 @param qf CeedQFunction 402 @param[out] ctx Variable to store context data values 403 404 @return An error code: 0 - success, otherwise - failure 405 406 @ref Advanced 407 **/ 408 409 int CeedQFunctionGetContext(CeedQFunction qf, void* *ctx) { 410 *ctx = qf->ctx; 411 return 0; 412 } 413 414 /** 415 @brief Determine if Fortran interface was used 416 417 @param qf CeedQFunction 418 @param[out] fortranstatus Variable to store Fortran status 419 420 @return An error code: 0 - success, otherwise - failure 421 422 @ref Advanced 423 **/ 424 425 int CeedQFunctionGetFortranStatus(CeedQFunction qf, bool *fortranstatus) { 426 *fortranstatus = qf->fortranstatus; 427 return 0; 428 } 429 430 /** 431 @brief Determine if QFunction is identity 432 433 @param qf CeedQFunction 434 @param[out] identity Variable to store identity status 435 436 @return An error code: 0 - success, otherwise - failure 437 438 @ref Advanced 439 **/ 440 441 int CeedQFunctionGetIdentityStatus(CeedQFunction qf, bool *identity) { 442 *identity = qf->identity; 443 return 0; 444 } 445 446 /** 447 @brief Get true user context for a CeedQFunction 448 449 @param qf CeedQFunction 450 @param[out] ctx Variable to store context data values 451 452 @return An error code: 0 - success, otherwise - failure 453 454 @ref Advanced 455 **/ 456 457 int CeedQFunctionGetInnerContext(CeedQFunction qf, void* *ctx) { 458 if (qf->fortranstatus) { 459 fContext *fctx = qf->ctx; 460 *ctx = fctx->innerctx; 461 } else { 462 *ctx = qf->ctx; 463 } 464 465 466 return 0; 467 } 468 469 /** 470 @brief Get backend data of a CeedQFunction 471 472 @param qf CeedQFunction 473 @param[out] data Variable to store data 474 475 @return An error code: 0 - success, otherwise - failure 476 477 @ref Advanced 478 **/ 479 480 int CeedQFunctionGetData(CeedQFunction qf, void* *data) { 481 *data = qf->data; 482 return 0; 483 } 484 485 /** 486 @brief Set backend data of a CeedQFunction 487 488 @param[out] qf CeedQFunction 489 @param data Data to set 490 491 @return An error code: 0 - success, otherwise - failure 492 493 @ref Advanced 494 **/ 495 496 int CeedQFunctionSetData(CeedQFunction qf, void* *data) { 497 qf->data = *data; 498 return 0; 499 } 500 501 /** 502 @brief Set global context for a CeedQFunction 503 504 @param qf CeedQFunction 505 @param ctx Context data to set 506 @param ctxsize Size of context data values 507 508 @return An error code: 0 - success, otherwise - failure 509 510 @ref Basic 511 **/ 512 int CeedQFunctionSetContext(CeedQFunction qf, void *ctx, size_t ctxsize) { 513 qf->ctx = ctx; 514 qf->ctxsize = ctxsize; 515 return 0; 516 } 517 518 /** 519 @brief Apply the action of a CeedQFunction 520 521 @param qf CeedQFunction 522 @param Q Number of quadrature points 523 @param[in] u Array of input data arrays 524 @param[out] v Array of output data arrays 525 526 @return An error code: 0 - success, otherwise - failure 527 528 @ref Advanced 529 **/ 530 int CeedQFunctionApply(CeedQFunction qf, CeedInt Q, 531 CeedVector *u, CeedVector *v) { 532 int ierr; 533 if (!qf->Apply) 534 // LCOV_EXCL_START 535 return CeedError(qf->ceed, 1, "Backend does not support QFunctionApply"); 536 // LCOV_EXCL_STOP 537 if (Q % qf->vlength) 538 // LCOV_EXCL_START 539 return CeedError(qf->ceed, 2, "Number of quadrature points %d must be a " 540 "multiple of %d", Q, qf->vlength); 541 // LCOV_EXCL_STOP 542 ierr = qf->Apply(qf, Q, u, v); CeedChk(ierr); 543 return 0; 544 } 545 546 /** 547 @brief Get the CeedQFunctionFields of a CeedQFunction 548 549 @param qf CeedQFunction 550 @param[out] inputfields Variable to store inputfields 551 @param[out] outputfields Variable to store outputfields 552 553 @return An error code: 0 - success, otherwise - failure 554 555 @ref Advanced 556 **/ 557 558 int CeedQFunctionGetFields(CeedQFunction qf, 559 CeedQFunctionField* *inputfields, 560 CeedQFunctionField* *outputfields) { 561 if (inputfields) 562 *inputfields = qf->inputfields; 563 if (outputfields) 564 *outputfields = qf->outputfields; 565 return 0; 566 } 567 568 /** 569 @brief Get the name of a CeedQFunctionField 570 571 @param qffield CeedQFunctionField 572 @param[out] fieldname Variable to store the field name 573 574 @return An error code: 0 - success, otherwise - failure 575 576 @ref Advanced 577 **/ 578 579 int CeedQFunctionFieldGetName(CeedQFunctionField qffield, 580 char* *fieldname) { 581 *fieldname = (char *)qffield->fieldname; 582 return 0; 583 } 584 585 /** 586 @brief Get the number of components of a CeedQFunctionField 587 588 @param qffield CeedQFunctionField 589 @param[out] size Variable to store the size of the field 590 591 @return An error code: 0 - success, otherwise - failure 592 593 @ref Advanced 594 **/ 595 596 int CeedQFunctionFieldGetSize(CeedQFunctionField qffield, CeedInt *size) { 597 *size = qffield->size; 598 return 0; 599 } 600 601 /** 602 @brief Get the CeedEvalMode of a CeedQFunctionField 603 604 @param qffield CeedQFunctionField 605 @param[out] emode Variable to store the field evaluation mode 606 607 @return An error code: 0 - success, otherwise - failure 608 609 @ref Advanced 610 **/ 611 612 int CeedQFunctionFieldGetEvalMode(CeedQFunctionField qffield, 613 CeedEvalMode *emode) { 614 *emode = qffield->emode; 615 return 0; 616 } 617 618 /** 619 @brief Destroy a CeedQFunction 620 621 @param qf CeedQFunction to destroy 622 623 @return An error code: 0 - success, otherwise - failure 624 625 @ref Basic 626 **/ 627 int CeedQFunctionDestroy(CeedQFunction *qf) { 628 int ierr; 629 630 if (!*qf || --(*qf)->refcount > 0) 631 return 0; 632 // Backend destroy 633 if ((*qf)->Destroy) { 634 ierr = (*qf)->Destroy(*qf); CeedChk(ierr); 635 } 636 // Free fields 637 for (int i=0; i<(*qf)->numinputfields; i++) { 638 ierr = CeedFree(&(*(*qf)->inputfields[i]).fieldname); CeedChk(ierr); 639 ierr = CeedFree(&(*qf)->inputfields[i]); CeedChk(ierr); 640 } 641 for (int i=0; i<(*qf)->numoutputfields; i++) { 642 ierr = CeedFree(&(*(*qf)->outputfields[i]).fieldname); CeedChk(ierr); 643 ierr = CeedFree(&(*qf)->outputfields[i]); CeedChk(ierr); 644 } 645 ierr = CeedFree(&(*qf)->inputfields); CeedChk(ierr); 646 ierr = CeedFree(&(*qf)->outputfields); CeedChk(ierr); 647 // Free ctx if identity 648 if ((*qf)->identity) { 649 ierr = CeedFree(&(*qf)->ctx); CeedChk(ierr); 650 } 651 652 ierr = CeedFree(&(*qf)->sourcepath); CeedChk(ierr); 653 ierr = CeedDestroy(&(*qf)->ceed); CeedChk(ierr); 654 ierr = CeedFree(qf); CeedChk(ierr); 655 return 0; 656 } 657 658 /// @} 659