1 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 2 #include <petsc/private/dmimpl.h> 3 4 static PetscErrorCode DMTSDestroy(DMTS *kdm) 5 { 6 PetscErrorCode ierr; 7 8 PetscFunctionBegin; 9 if (!*kdm) PetscFunctionReturn(0); 10 PetscValidHeaderSpecific((*kdm),DMTS_CLASSID,1); 11 if (--((PetscObject)(*kdm))->refct > 0) {*kdm = 0; PetscFunctionReturn(0);} 12 if ((*kdm)->ops->destroy) {ierr = ((*kdm)->ops->destroy)(*kdm);CHKERRQ(ierr);} 13 ierr = PetscHeaderDestroy(kdm);CHKERRQ(ierr); 14 PetscFunctionReturn(0); 15 } 16 17 PetscErrorCode DMTSLoad(DMTS kdm,PetscViewer viewer) 18 { 19 PetscErrorCode ierr; 20 21 PetscFunctionBegin; 22 ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ifunction,1,NULL,PETSC_FUNCTION);CHKERRQ(ierr); 23 ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ifunctionview,1,NULL,PETSC_FUNCTION);CHKERRQ(ierr); 24 ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ifunctionload,1,NULL,PETSC_FUNCTION);CHKERRQ(ierr); 25 if (kdm->ops->ifunctionload) { 26 ierr = (*kdm->ops->ifunctionload)(&kdm->ifunctionctx,viewer);CHKERRQ(ierr); 27 } 28 ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ijacobian,1,NULL,PETSC_FUNCTION);CHKERRQ(ierr); 29 ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ijacobianview,1,NULL,PETSC_FUNCTION);CHKERRQ(ierr); 30 ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ijacobianload,1,NULL,PETSC_FUNCTION);CHKERRQ(ierr); 31 if (kdm->ops->ijacobianload) { 32 ierr = (*kdm->ops->ijacobianload)(&kdm->ijacobianctx,viewer);CHKERRQ(ierr); 33 } 34 PetscFunctionReturn(0); 35 } 36 37 PetscErrorCode DMTSView(DMTS kdm,PetscViewer viewer) 38 { 39 PetscErrorCode ierr; 40 PetscBool isascii,isbinary; 41 42 PetscFunctionBegin; 43 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 44 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 45 if (isascii) { 46 #if defined(PETSC_SERIALIZE_FUNCTIONS) 47 const char *fname; 48 49 ierr = PetscFPTFind(kdm->ops->ifunction,&fname);CHKERRQ(ierr); 50 if (fname) { 51 ierr = PetscViewerASCIIPrintf(viewer," IFunction used by TS: %s\n",fname);CHKERRQ(ierr); 52 } 53 ierr = PetscFPTFind(kdm->ops->ijacobian,&fname);CHKERRQ(ierr); 54 if (fname) { 55 ierr = PetscViewerASCIIPrintf(viewer," IJacobian function used by TS: %s\n",fname);CHKERRQ(ierr); 56 } 57 #endif 58 } else if (isbinary) { 59 struct { 60 TSIFunction ifunction; 61 } funcstruct; 62 struct { 63 PetscErrorCode (*ifunctionview)(void*,PetscViewer); 64 } funcviewstruct; 65 struct { 66 PetscErrorCode (*ifunctionload)(void**,PetscViewer); 67 } funcloadstruct; 68 struct { 69 TSIJacobian ijacobian; 70 } jacstruct; 71 struct { 72 PetscErrorCode (*ijacobianview)(void*,PetscViewer); 73 } jacviewstruct; 74 struct { 75 PetscErrorCode (*ijacobianload)(void**,PetscViewer); 76 } jacloadstruct; 77 78 funcstruct.ifunction = kdm->ops->ifunction; 79 funcviewstruct.ifunctionview = kdm->ops->ifunctionview; 80 funcloadstruct.ifunctionload = kdm->ops->ifunctionload; 81 ierr = PetscViewerBinaryWrite(viewer,&funcstruct,1,PETSC_FUNCTION);CHKERRQ(ierr); 82 ierr = PetscViewerBinaryWrite(viewer,&funcviewstruct,1,PETSC_FUNCTION);CHKERRQ(ierr); 83 ierr = PetscViewerBinaryWrite(viewer,&funcloadstruct,1,PETSC_FUNCTION);CHKERRQ(ierr); 84 if (kdm->ops->ifunctionview) { 85 ierr = (*kdm->ops->ifunctionview)(kdm->ifunctionctx,viewer);CHKERRQ(ierr); 86 } 87 jacstruct.ijacobian = kdm->ops->ijacobian; 88 jacviewstruct.ijacobianview = kdm->ops->ijacobianview; 89 jacloadstruct.ijacobianload = kdm->ops->ijacobianload; 90 ierr = PetscViewerBinaryWrite(viewer,&jacstruct,1,PETSC_FUNCTION);CHKERRQ(ierr); 91 ierr = PetscViewerBinaryWrite(viewer,&jacviewstruct,1,PETSC_FUNCTION);CHKERRQ(ierr); 92 ierr = PetscViewerBinaryWrite(viewer,&jacloadstruct,1,PETSC_FUNCTION);CHKERRQ(ierr); 93 if (kdm->ops->ijacobianview) { 94 ierr = (*kdm->ops->ijacobianview)(kdm->ijacobianctx,viewer);CHKERRQ(ierr); 95 } 96 } 97 PetscFunctionReturn(0); 98 } 99 100 static PetscErrorCode DMTSCreate(MPI_Comm comm,DMTS *kdm) 101 { 102 PetscErrorCode ierr; 103 104 PetscFunctionBegin; 105 ierr = TSInitializePackage();CHKERRQ(ierr); 106 ierr = PetscHeaderCreate(*kdm, DMTS_CLASSID, "DMTS", "DMTS", "DMTS", comm, DMTSDestroy, DMTSView);CHKERRQ(ierr); 107 PetscFunctionReturn(0); 108 } 109 110 /* Attaches the DMTS to the coarse level. 111 * Under what conditions should we copy versus duplicate? 112 */ 113 static PetscErrorCode DMCoarsenHook_DMTS(DM dm,DM dmc,void *ctx) 114 { 115 PetscErrorCode ierr; 116 117 PetscFunctionBegin; 118 ierr = DMCopyDMTS(dm,dmc);CHKERRQ(ierr); 119 PetscFunctionReturn(0); 120 } 121 122 /* This could restrict auxiliary information to the coarse level. 123 */ 124 static PetscErrorCode DMRestrictHook_DMTS(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx) 125 { 126 127 PetscFunctionBegin; 128 PetscFunctionReturn(0); 129 } 130 131 static PetscErrorCode DMSubDomainHook_DMTS(DM dm,DM subdm,void *ctx) 132 { 133 PetscErrorCode ierr; 134 135 PetscFunctionBegin; 136 ierr = DMCopyDMTS(dm,subdm);CHKERRQ(ierr); 137 PetscFunctionReturn(0); 138 } 139 140 /* This could restrict auxiliary information to the coarse level. 141 */ 142 static PetscErrorCode DMSubDomainRestrictHook_DMTS(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 143 { 144 PetscFunctionBegin; 145 PetscFunctionReturn(0); 146 } 147 148 /*@C 149 DMTSCopy - copies the information in a DMTS to another DMTS 150 151 Not Collective 152 153 Input Argument: 154 + kdm - Original DMTS 155 - nkdm - DMTS to receive the data, should have been created with DMTSCreate() 156 157 Level: developer 158 159 .seealso: DMTSCreate(), DMTSDestroy() 160 @*/ 161 PetscErrorCode DMTSCopy(DMTS kdm,DMTS nkdm) 162 { 163 PetscErrorCode ierr; 164 165 PetscFunctionBegin; 166 PetscValidHeaderSpecific(kdm,DMTS_CLASSID,1); 167 PetscValidHeaderSpecific(nkdm,DMTS_CLASSID,2); 168 nkdm->ops->rhsfunction = kdm->ops->rhsfunction; 169 nkdm->ops->rhsjacobian = kdm->ops->rhsjacobian; 170 nkdm->ops->ifunction = kdm->ops->ifunction; 171 nkdm->ops->ijacobian = kdm->ops->ijacobian; 172 nkdm->ops->i2function = kdm->ops->i2function; 173 nkdm->ops->i2jacobian = kdm->ops->i2jacobian; 174 nkdm->ops->solution = kdm->ops->solution; 175 nkdm->ops->destroy = kdm->ops->destroy; 176 nkdm->ops->duplicate = kdm->ops->duplicate; 177 178 nkdm->rhsfunctionctx = kdm->rhsfunctionctx; 179 nkdm->rhsjacobianctx = kdm->rhsjacobianctx; 180 nkdm->ifunctionctx = kdm->ifunctionctx; 181 nkdm->ijacobianctx = kdm->ijacobianctx; 182 nkdm->i2functionctx = kdm->i2functionctx; 183 nkdm->i2jacobianctx = kdm->i2jacobianctx; 184 nkdm->solutionctx = kdm->solutionctx; 185 186 nkdm->data = kdm->data; 187 188 /* 189 nkdm->fortran_func_pointers[0] = kdm->fortran_func_pointers[0]; 190 nkdm->fortran_func_pointers[1] = kdm->fortran_func_pointers[1]; 191 nkdm->fortran_func_pointers[2] = kdm->fortran_func_pointers[2]; 192 */ 193 194 /* implementation specific copy hooks */ 195 if (kdm->ops->duplicate) {ierr = (*kdm->ops->duplicate)(kdm,nkdm);CHKERRQ(ierr);} 196 PetscFunctionReturn(0); 197 } 198 199 /*@C 200 DMGetDMTS - get read-only private DMTS context from a DM 201 202 Not Collective 203 204 Input Argument: 205 . dm - DM to be used with TS 206 207 Output Argument: 208 . tsdm - private DMTS context 209 210 Level: developer 211 212 Notes: 213 Use DMGetDMTSWrite() if write access is needed. The DMTSSetXXX API should be used wherever possible. 214 215 .seealso: DMGetDMTSWrite() 216 @*/ 217 PetscErrorCode DMGetDMTS(DM dm,DMTS *tsdm) 218 { 219 PetscErrorCode ierr; 220 221 PetscFunctionBegin; 222 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 223 *tsdm = (DMTS) dm->dmts; 224 if (!*tsdm) { 225 ierr = PetscInfo(dm,"Creating new DMTS\n");CHKERRQ(ierr); 226 ierr = DMTSCreate(PetscObjectComm((PetscObject)dm),tsdm);CHKERRQ(ierr); 227 dm->dmts = (PetscObject) *tsdm; 228 (*tsdm)->originaldm = dm; 229 ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,NULL);CHKERRQ(ierr); 230 ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_DMTS,DMSubDomainRestrictHook_DMTS,NULL);CHKERRQ(ierr); 231 } 232 PetscFunctionReturn(0); 233 } 234 235 /*@C 236 DMGetDMTSWrite - get write access to private DMTS context from a DM 237 238 Not Collective 239 240 Input Argument: 241 . dm - DM to be used with TS 242 243 Output Argument: 244 . tsdm - private DMTS context 245 246 Level: developer 247 248 .seealso: DMGetDMTS() 249 @*/ 250 PetscErrorCode DMGetDMTSWrite(DM dm,DMTS *tsdm) 251 { 252 PetscErrorCode ierr; 253 DMTS sdm; 254 255 PetscFunctionBegin; 256 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 257 ierr = DMGetDMTS(dm,&sdm);CHKERRQ(ierr); 258 if (!sdm->originaldm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"DMTS has a NULL originaldm"); 259 if (sdm->originaldm != dm) { /* Copy on write */ 260 DMTS oldsdm = sdm; 261 ierr = PetscInfo(dm,"Copying DMTS due to write\n");CHKERRQ(ierr); 262 ierr = DMTSCreate(PetscObjectComm((PetscObject)dm),&sdm);CHKERRQ(ierr); 263 ierr = DMTSCopy(oldsdm,sdm);CHKERRQ(ierr); 264 ierr = DMTSDestroy((DMTS*)&dm->dmts);CHKERRQ(ierr); 265 dm->dmts = (PetscObject) sdm; 266 sdm->originaldm = dm; 267 } 268 *tsdm = sdm; 269 PetscFunctionReturn(0); 270 } 271 272 /*@C 273 DMCopyDMTS - copies a DM context to a new DM 274 275 Logically Collective 276 277 Input Arguments: 278 + dmsrc - DM to obtain context from 279 - dmdest - DM to add context to 280 281 Level: developer 282 283 Note: 284 The context is copied by reference. This function does not ensure that a context exists. 285 286 .seealso: DMGetDMTS(), TSSetDM() 287 @*/ 288 PetscErrorCode DMCopyDMTS(DM dmsrc,DM dmdest) 289 { 290 PetscErrorCode ierr; 291 292 PetscFunctionBegin; 293 PetscValidHeaderSpecific(dmsrc,DM_CLASSID,1); 294 PetscValidHeaderSpecific(dmdest,DM_CLASSID,2); 295 ierr = DMTSDestroy((DMTS*)&dmdest->dmts);CHKERRQ(ierr); 296 dmdest->dmts = dmsrc->dmts; 297 ierr = PetscObjectReference(dmdest->dmts);CHKERRQ(ierr); 298 ierr = DMCoarsenHookAdd(dmdest,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,NULL);CHKERRQ(ierr); 299 ierr = DMSubDomainHookAdd(dmdest,DMSubDomainHook_DMTS,DMSubDomainRestrictHook_DMTS,NULL);CHKERRQ(ierr); 300 PetscFunctionReturn(0); 301 } 302 303 /*@C 304 DMTSSetIFunction - set TS implicit function evaluation function 305 306 Not Collective 307 308 Input Arguments: 309 + dm - DM to be used with TS 310 . func - function evaluation function, see TSSetIFunction() for calling sequence 311 - ctx - context for residual evaluation 312 313 Level: advanced 314 315 Note: 316 TSSetFunction() is normally used, but it calls this function internally because the user context is actually 317 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 318 not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 319 320 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian() 321 @*/ 322 PetscErrorCode DMTSSetIFunction(DM dm,TSIFunction func,void *ctx) 323 { 324 PetscErrorCode ierr; 325 DMTS tsdm; 326 327 PetscFunctionBegin; 328 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 329 ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr); 330 if (func) tsdm->ops->ifunction = func; 331 if (ctx) tsdm->ifunctionctx = ctx; 332 PetscFunctionReturn(0); 333 } 334 335 /*@C 336 DMTSGetIFunction - get TS implicit residual evaluation function 337 338 Not Collective 339 340 Input Argument: 341 . dm - DM to be used with TS 342 343 Output Arguments: 344 + func - function evaluation function, see TSSetIFunction() for calling sequence 345 - ctx - context for residual evaluation 346 347 Level: advanced 348 349 Note: 350 TSGetFunction() is normally used, but it calls this function internally because the user context is actually 351 associated with the DM. 352 353 .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction() 354 @*/ 355 PetscErrorCode DMTSGetIFunction(DM dm,TSIFunction *func,void **ctx) 356 { 357 PetscErrorCode ierr; 358 DMTS tsdm; 359 360 PetscFunctionBegin; 361 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 362 ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr); 363 if (func) *func = tsdm->ops->ifunction; 364 if (ctx) *ctx = tsdm->ifunctionctx; 365 PetscFunctionReturn(0); 366 } 367 368 /*@C 369 DMTSSetI2Function - set TS implicit function evaluation function for 2nd order systems 370 371 Not Collective 372 373 Input Arguments: 374 + dm - DM to be used with TS 375 . fun - function evaluation function, see TSSetI2Function() for calling sequence 376 - ctx - context for residual evaluation 377 378 Level: advanced 379 380 Note: 381 TSSetI2Function() is normally used, but it calls this function internally because the user context is actually 382 associated with the DM. 383 384 .seealso: TSSetI2Function() 385 @*/ 386 PetscErrorCode DMTSSetI2Function(DM dm,TSI2Function fun,void *ctx) 387 { 388 DMTS tsdm; 389 PetscErrorCode ierr; 390 391 PetscFunctionBegin; 392 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 393 ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr); 394 if (fun) tsdm->ops->i2function = fun; 395 if (ctx) tsdm->i2functionctx = ctx; 396 PetscFunctionReturn(0); 397 } 398 399 /*@C 400 DMTSGetI2Function - get TS implicit residual evaluation function for 2nd order systems 401 402 Not Collective 403 404 Input Argument: 405 . dm - DM to be used with TS 406 407 Output Arguments: 408 + fun - function evaluation function, see TSSetI2Function() for calling sequence 409 - ctx - context for residual evaluation 410 411 Level: advanced 412 413 Note: 414 TSGetI2Function() is normally used, but it calls this function internally because the user context is actually 415 associated with the DM. 416 417 .seealso: DMTSSetI2Function(),TSGetI2Function() 418 @*/ 419 PetscErrorCode DMTSGetI2Function(DM dm,TSI2Function *fun,void **ctx) 420 { 421 DMTS tsdm; 422 PetscErrorCode ierr; 423 424 PetscFunctionBegin; 425 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 426 ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr); 427 if (fun) *fun = tsdm->ops->i2function; 428 if (ctx) *ctx = tsdm->i2functionctx; 429 PetscFunctionReturn(0); 430 } 431 432 /*@C 433 DMTSSetI2Jacobian - set TS implicit Jacobian evaluation function for 2nd order systems 434 435 Not Collective 436 437 Input Arguments: 438 + dm - DM to be used with TS 439 . fun - Jacobian evaluation function, see TSSetI2Jacobian() for calling sequence 440 - ctx - context for Jacobian evaluation 441 442 Level: advanced 443 444 Note: 445 TSSetI2Jacobian() is normally used, but it calls this function internally because the user context is actually 446 associated with the DM. 447 448 .seealso: TSSetI2Jacobian() 449 @*/ 450 PetscErrorCode DMTSSetI2Jacobian(DM dm,TSI2Jacobian jac,void *ctx) 451 { 452 DMTS tsdm; 453 PetscErrorCode ierr; 454 455 PetscFunctionBegin; 456 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 457 ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr); 458 if (jac) tsdm->ops->i2jacobian = jac; 459 if (ctx) tsdm->i2jacobianctx = ctx; 460 PetscFunctionReturn(0); 461 } 462 463 /*@C 464 DMTSGetI2Jacobian - get TS implicit Jacobian evaluation function for 2nd order systems 465 466 Not Collective 467 468 Input Argument: 469 . dm - DM to be used with TS 470 471 Output Arguments: 472 + jac - Jacobian evaluation function, see TSSetI2Jacobian() for calling sequence 473 - ctx - context for Jacobian evaluation 474 475 Level: advanced 476 477 Note: 478 TSGetI2Jacobian() is normally used, but it calls this function internally because the user context is actually 479 associated with the DM. 480 481 .seealso: DMTSSetI2Jacobian(),TSGetI2Jacobian() 482 @*/ 483 PetscErrorCode DMTSGetI2Jacobian(DM dm,TSI2Jacobian *jac,void **ctx) 484 { 485 DMTS tsdm; 486 PetscErrorCode ierr; 487 488 PetscFunctionBegin; 489 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 490 ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr); 491 if (jac) *jac = tsdm->ops->i2jacobian; 492 if (ctx) *ctx = tsdm->i2jacobianctx; 493 PetscFunctionReturn(0); 494 } 495 496 /*@C 497 DMTSSetRHSFunction - set TS explicit residual evaluation function 498 499 Not Collective 500 501 Input Arguments: 502 + dm - DM to be used with TS 503 . func - RHS function evaluation function, see TSSetRHSFunction() for calling sequence 504 - ctx - context for residual evaluation 505 506 Level: advanced 507 508 Note: 509 TSSetRSHFunction() is normally used, but it calls this function internally because the user context is actually 510 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 511 not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 512 513 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian() 514 @*/ 515 PetscErrorCode DMTSSetRHSFunction(DM dm,TSRHSFunction func,void *ctx) 516 { 517 PetscErrorCode ierr; 518 DMTS tsdm; 519 520 PetscFunctionBegin; 521 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 522 ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr); 523 if (func) tsdm->ops->rhsfunction = func; 524 if (ctx) tsdm->rhsfunctionctx = ctx; 525 PetscFunctionReturn(0); 526 } 527 528 /*@C 529 DMTSSetTransientVariable - sets function to transform from state to transient variables 530 531 Logically Collective 532 533 Input Arguments: 534 + dm - DM to be used with TS 535 . tvar - a function that transforms in-place to transient variables 536 - ctx - a context for tvar 537 538 Level: advanced 539 540 Notes: 541 This is typically used to transform from primitive to conservative variables so that a time integrator (e.g., TSBDF) 542 can be conservative. In this context, primitive variables P are used to model the state (e.g., because they lead to 543 well-conditioned formulations even in limiting cases such as low-Mach or zero porosity). The transient variable is 544 C(P), specified by calling this function. An IFunction thus receives arguments (P, Cdot) and the IJacobian must be 545 evaluated via the chain rule, as in 546 547 dF/dP + shift * dF/dCdot dC/dP. 548 549 .seealso: TSSetTransientVariable(), DMTSGetTransientVariable(), DMTSSetIFunction(), DMTSSetIJacobian() 550 @*/ 551 PetscErrorCode DMTSSetTransientVariable(DM dm,TSTransientVariable tvar,void *ctx) 552 { 553 PetscErrorCode ierr; 554 DMTS dmts; 555 556 PetscFunctionBegin; 557 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 558 ierr = DMGetDMTSWrite(dm,&dmts);CHKERRQ(ierr); 559 dmts->ops->transientvar = tvar; 560 dmts->transientvarctx = ctx; 561 PetscFunctionReturn(0); 562 } 563 564 /*@C 565 DMTSGetTransientVariable - gets function to transform from state to transient variables 566 567 Logically Collective 568 569 Input Arguments: 570 . dm - DM to be used with TS 571 572 Output Arguments: 573 + tvar - a function that transforms in-place to transient variables 574 - ctx - a context for tvar 575 576 Level: advanced 577 578 .seealso: DMTSSetTransientVariable(), DMTSGetIFunction(), DMTSGetIJacobian() 579 @*/ 580 PetscErrorCode DMTSGetTransientVariable(DM dm,TSTransientVariable *tvar,void *ctx) 581 { 582 PetscErrorCode ierr; 583 DMTS dmts; 584 585 PetscFunctionBegin; 586 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 587 ierr = DMGetDMTS(dm,&dmts);CHKERRQ(ierr); 588 if (tvar) *tvar = dmts->ops->transientvar; 589 if (ctx) *(void**)ctx = dmts->transientvarctx; 590 PetscFunctionReturn(0); 591 } 592 593 /*@C 594 DMTSGetSolutionFunction - gets the TS solution evaluation function 595 596 Not Collective 597 598 Input Arguments: 599 . dm - DM to be used with TS 600 601 Output Parameters: 602 + func - solution function evaluation function, see TSSetSolution() for calling sequence 603 - ctx - context for solution evaluation 604 605 Level: advanced 606 607 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), DMTSSetSolutionFunction() 608 @*/ 609 PetscErrorCode DMTSGetSolutionFunction(DM dm,TSSolutionFunction *func,void **ctx) 610 { 611 PetscErrorCode ierr; 612 DMTS tsdm; 613 614 PetscFunctionBegin; 615 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 616 ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr); 617 if (func) *func = tsdm->ops->solution; 618 if (ctx) *ctx = tsdm->solutionctx; 619 PetscFunctionReturn(0); 620 } 621 622 /*@C 623 DMTSSetSolutionFunction - set TS solution evaluation function 624 625 Not Collective 626 627 Input Arguments: 628 + dm - DM to be used with TS 629 . func - solution function evaluation function, see TSSetSolution() for calling sequence 630 - ctx - context for solution evaluation 631 632 Level: advanced 633 634 Note: 635 TSSetSolutionFunction() is normally used, but it calls this function internally because the user context is actually 636 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 637 not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 638 639 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), DMTSGetSolutionFunction() 640 @*/ 641 PetscErrorCode DMTSSetSolutionFunction(DM dm,TSSolutionFunction func,void *ctx) 642 { 643 PetscErrorCode ierr; 644 DMTS tsdm; 645 646 PetscFunctionBegin; 647 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 648 ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr); 649 if (func) tsdm->ops->solution = func; 650 if (ctx) tsdm->solutionctx = ctx; 651 PetscFunctionReturn(0); 652 } 653 654 /*@C 655 DMTSSetForcingFunction - set TS forcing function evaluation function 656 657 Not Collective 658 659 Input Arguments: 660 + dm - DM to be used with TS 661 . f - forcing function evaluation function; see TSForcingFunction 662 - ctx - context for solution evaluation 663 664 Level: advanced 665 666 Note: 667 TSSetForcingFunction() is normally used, but it calls this function internally because the user context is actually 668 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 669 not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 670 671 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), TSSetForcingFunction(), DMTSGetForcingFunction() 672 @*/ 673 PetscErrorCode DMTSSetForcingFunction(DM dm,TSForcingFunction f,void *ctx) 674 { 675 PetscErrorCode ierr; 676 DMTS tsdm; 677 678 PetscFunctionBegin; 679 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 680 ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr); 681 if (f) tsdm->ops->forcing = f; 682 if (ctx) tsdm->forcingctx = ctx; 683 PetscFunctionReturn(0); 684 } 685 686 687 /*@C 688 DMTSGetForcingFunction - get TS forcing function evaluation function 689 690 Not Collective 691 692 Input Argument: 693 . dm - DM to be used with TS 694 695 Output Arguments: 696 + f - forcing function evaluation function; see TSForcingFunction for details 697 - ctx - context for solution evaluation 698 699 Level: advanced 700 701 Note: 702 TSSetForcingFunction() is normally used, but it calls this function internally because the user context is actually 703 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 704 not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 705 706 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), TSSetForcingFunction(), DMTSGetForcingFunction() 707 @*/ 708 PetscErrorCode DMTSGetForcingFunction(DM dm,TSForcingFunction *f,void **ctx) 709 { 710 PetscErrorCode ierr; 711 DMTS tsdm; 712 713 PetscFunctionBegin; 714 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 715 ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr); 716 if (f) *f = tsdm->ops->forcing; 717 if (ctx) *ctx = tsdm->forcingctx; 718 PetscFunctionReturn(0); 719 } 720 721 /*@C 722 DMTSGetRHSFunction - get TS explicit residual evaluation function 723 724 Not Collective 725 726 Input Argument: 727 . dm - DM to be used with TS 728 729 Output Arguments: 730 + func - residual evaluation function, see TSSetRHSFunction() for calling sequence 731 - ctx - context for residual evaluation 732 733 Level: advanced 734 735 Note: 736 TSGetFunction() is normally used, but it calls this function internally because the user context is actually 737 associated with the DM. 738 739 .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction() 740 @*/ 741 PetscErrorCode DMTSGetRHSFunction(DM dm,TSRHSFunction *func,void **ctx) 742 { 743 PetscErrorCode ierr; 744 DMTS tsdm; 745 746 PetscFunctionBegin; 747 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 748 ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr); 749 if (func) *func = tsdm->ops->rhsfunction; 750 if (ctx) *ctx = tsdm->rhsfunctionctx; 751 PetscFunctionReturn(0); 752 } 753 754 /*@C 755 DMTSSetIJacobian - set TS Jacobian evaluation function 756 757 Not Collective 758 759 Input Argument: 760 + dm - DM to be used with TS 761 . func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence 762 - ctx - context for residual evaluation 763 764 Level: advanced 765 766 Note: 767 TSSetJacobian() is normally used, but it calls this function internally because the user context is actually 768 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 769 not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 770 771 .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian() 772 @*/ 773 PetscErrorCode DMTSSetIJacobian(DM dm,TSIJacobian func,void *ctx) 774 { 775 PetscErrorCode ierr; 776 DMTS sdm; 777 778 PetscFunctionBegin; 779 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 780 ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr); 781 if (func) sdm->ops->ijacobian = func; 782 if (ctx) sdm->ijacobianctx = ctx; 783 PetscFunctionReturn(0); 784 } 785 786 /*@C 787 DMTSGetIJacobian - get TS Jacobian evaluation function 788 789 Not Collective 790 791 Input Argument: 792 . dm - DM to be used with TS 793 794 Output Arguments: 795 + func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence 796 - ctx - context for residual evaluation 797 798 Level: advanced 799 800 Note: 801 TSGetJacobian() is normally used, but it calls this function internally because the user context is actually 802 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 803 not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 804 805 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian() 806 @*/ 807 PetscErrorCode DMTSGetIJacobian(DM dm,TSIJacobian *func,void **ctx) 808 { 809 PetscErrorCode ierr; 810 DMTS tsdm; 811 812 PetscFunctionBegin; 813 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 814 ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr); 815 if (func) *func = tsdm->ops->ijacobian; 816 if (ctx) *ctx = tsdm->ijacobianctx; 817 PetscFunctionReturn(0); 818 } 819 820 821 /*@C 822 DMTSSetRHSJacobian - set TS Jacobian evaluation function 823 824 Not Collective 825 826 Input Argument: 827 + dm - DM to be used with TS 828 . func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence 829 - ctx - context for residual evaluation 830 831 Level: advanced 832 833 Note: 834 TSSetJacobian() is normally used, but it calls this function internally because the user context is actually 835 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 836 not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 837 838 .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian() 839 @*/ 840 PetscErrorCode DMTSSetRHSJacobian(DM dm,TSRHSJacobian func,void *ctx) 841 { 842 PetscErrorCode ierr; 843 DMTS tsdm; 844 845 PetscFunctionBegin; 846 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 847 ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr); 848 if (func) tsdm->ops->rhsjacobian = func; 849 if (ctx) tsdm->rhsjacobianctx = ctx; 850 PetscFunctionReturn(0); 851 } 852 853 /*@C 854 DMTSGetRHSJacobian - get TS Jacobian evaluation function 855 856 Not Collective 857 858 Input Argument: 859 . dm - DM to be used with TS 860 861 Output Arguments: 862 + func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence 863 - ctx - context for residual evaluation 864 865 Level: advanced 866 867 Note: 868 TSGetJacobian() is normally used, but it calls this function internally because the user context is actually 869 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 870 not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 871 872 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian() 873 @*/ 874 PetscErrorCode DMTSGetRHSJacobian(DM dm,TSRHSJacobian *func,void **ctx) 875 { 876 PetscErrorCode ierr; 877 DMTS tsdm; 878 879 PetscFunctionBegin; 880 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 881 ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr); 882 if (func) *func = tsdm->ops->rhsjacobian; 883 if (ctx) *ctx = tsdm->rhsjacobianctx; 884 PetscFunctionReturn(0); 885 } 886 887 /*@C 888 DMTSSetIFunctionSerialize - sets functions used to view and load a IFunction context 889 890 Not Collective 891 892 Input Arguments: 893 + dm - DM to be used with TS 894 . view - viewer function 895 - load - loading function 896 897 Level: advanced 898 899 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian() 900 @*/ 901 PetscErrorCode DMTSSetIFunctionSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer)) 902 { 903 PetscErrorCode ierr; 904 DMTS tsdm; 905 906 PetscFunctionBegin; 907 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 908 ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr); 909 tsdm->ops->ifunctionview = view; 910 tsdm->ops->ifunctionload = load; 911 PetscFunctionReturn(0); 912 } 913 914 /*@C 915 DMTSSetIJacobianSerialize - sets functions used to view and load a IJacobian context 916 917 Not Collective 918 919 Input Arguments: 920 + dm - DM to be used with TS 921 . view - viewer function 922 - load - loading function 923 924 Level: advanced 925 926 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian() 927 @*/ 928 PetscErrorCode DMTSSetIJacobianSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer)) 929 { 930 PetscErrorCode ierr; 931 DMTS tsdm; 932 933 PetscFunctionBegin; 934 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 935 ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr); 936 tsdm->ops->ijacobianview = view; 937 tsdm->ops->ijacobianload = load; 938 PetscFunctionReturn(0); 939 } 940