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