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