1 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 2 #include <petsc/private/dmimpl.h> 3 4 static PetscErrorCode DMTSUnsetRHSFunctionContext_DMTS(DMTS tsdm) { 5 PetscFunctionBegin; 6 PetscCall(PetscObjectCompose((PetscObject)tsdm, "rhs function ctx", NULL)); 7 tsdm->rhsfunctionctxcontainer = NULL; 8 PetscFunctionReturn(0); 9 } 10 11 static PetscErrorCode DMTSUnsetRHSJacobianContext_DMTS(DMTS tsdm) { 12 PetscFunctionBegin; 13 PetscCall(PetscObjectCompose((PetscObject)tsdm, "rhs jacobian ctx", NULL)); 14 tsdm->rhsjacobianctxcontainer = NULL; 15 PetscFunctionReturn(0); 16 } 17 18 static PetscErrorCode DMTSUnsetIFunctionContext_DMTS(DMTS tsdm) { 19 PetscFunctionBegin; 20 PetscCall(PetscObjectCompose((PetscObject)tsdm, "ifunction ctx", NULL)); 21 tsdm->ifunctionctxcontainer = NULL; 22 PetscFunctionReturn(0); 23 } 24 25 static PetscErrorCode DMTSUnsetIJacobianContext_DMTS(DMTS tsdm) { 26 PetscFunctionBegin; 27 PetscCall(PetscObjectCompose((PetscObject)tsdm, "ijacobian ctx", NULL)); 28 tsdm->ijacobianctxcontainer = NULL; 29 PetscFunctionReturn(0); 30 } 31 32 static PetscErrorCode DMTSUnsetI2FunctionContext_DMTS(DMTS tsdm) { 33 PetscFunctionBegin; 34 PetscCall(PetscObjectCompose((PetscObject)tsdm, "i2function ctx", NULL)); 35 tsdm->i2functionctxcontainer = NULL; 36 PetscFunctionReturn(0); 37 } 38 39 static PetscErrorCode DMTSUnsetI2JacobianContext_DMTS(DMTS tsdm) { 40 PetscFunctionBegin; 41 PetscCall(PetscObjectCompose((PetscObject)tsdm, "i2jacobian ctx", NULL)); 42 tsdm->i2jacobianctxcontainer = NULL; 43 PetscFunctionReturn(0); 44 } 45 46 static PetscErrorCode DMTSDestroy(DMTS *kdm) { 47 PetscFunctionBegin; 48 if (!*kdm) PetscFunctionReturn(0); 49 PetscValidHeaderSpecific((*kdm), DMTS_CLASSID, 1); 50 if (--((PetscObject)(*kdm))->refct > 0) { 51 *kdm = NULL; 52 PetscFunctionReturn(0); 53 } 54 PetscCall(DMTSUnsetRHSFunctionContext_DMTS(*kdm)); 55 PetscCall(DMTSUnsetRHSJacobianContext_DMTS(*kdm)); 56 PetscCall(DMTSUnsetIFunctionContext_DMTS(*kdm)); 57 PetscCall(DMTSUnsetIJacobianContext_DMTS(*kdm)); 58 PetscCall(DMTSUnsetI2FunctionContext_DMTS(*kdm)); 59 PetscCall(DMTSUnsetI2JacobianContext_DMTS(*kdm)); 60 PetscTryTypeMethod(*kdm, destroy); 61 PetscCall(PetscHeaderDestroy(kdm)); 62 PetscFunctionReturn(0); 63 } 64 65 PetscErrorCode DMTSLoad(DMTS kdm, PetscViewer viewer) { 66 PetscFunctionBegin; 67 PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->ifunction, 1, NULL, PETSC_FUNCTION)); 68 PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->ifunctionview, 1, NULL, PETSC_FUNCTION)); 69 PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->ifunctionload, 1, NULL, PETSC_FUNCTION)); 70 if (kdm->ops->ifunctionload) { 71 void *ctx; 72 73 PetscCall(PetscContainerGetPointer(kdm->ifunctionctxcontainer, &ctx)); 74 PetscCall((*kdm->ops->ifunctionload)(&ctx, viewer)); 75 } 76 PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->ijacobian, 1, NULL, PETSC_FUNCTION)); 77 PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->ijacobianview, 1, NULL, PETSC_FUNCTION)); 78 PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->ijacobianload, 1, NULL, PETSC_FUNCTION)); 79 if (kdm->ops->ijacobianload) { 80 void *ctx; 81 82 PetscCall(PetscContainerGetPointer(kdm->ijacobianctxcontainer, &ctx)); 83 PetscCall((*kdm->ops->ijacobianload)(&ctx, viewer)); 84 } 85 PetscFunctionReturn(0); 86 } 87 88 PetscErrorCode DMTSView(DMTS kdm, PetscViewer viewer) { 89 PetscBool isascii, isbinary; 90 91 PetscFunctionBegin; 92 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 93 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 94 if (isascii) { 95 #if defined(PETSC_SERIALIZE_FUNCTIONS) 96 const char *fname; 97 98 PetscCall(PetscFPTFind(kdm->ops->ifunction, &fname)); 99 if (fname) PetscCall(PetscViewerASCIIPrintf(viewer, " IFunction used by TS: %s\n", fname)); 100 PetscCall(PetscFPTFind(kdm->ops->ijacobian, &fname)); 101 if (fname) PetscCall(PetscViewerASCIIPrintf(viewer, " IJacobian function used by TS: %s\n", fname)); 102 #endif 103 } else if (isbinary) { 104 struct { 105 TSIFunction ifunction; 106 } funcstruct; 107 struct { 108 PetscErrorCode (*ifunctionview)(void *, PetscViewer); 109 } funcviewstruct; 110 struct { 111 PetscErrorCode (*ifunctionload)(void **, PetscViewer); 112 } funcloadstruct; 113 struct { 114 TSIJacobian ijacobian; 115 } jacstruct; 116 struct { 117 PetscErrorCode (*ijacobianview)(void *, PetscViewer); 118 } jacviewstruct; 119 struct { 120 PetscErrorCode (*ijacobianload)(void **, PetscViewer); 121 } jacloadstruct; 122 123 funcstruct.ifunction = kdm->ops->ifunction; 124 funcviewstruct.ifunctionview = kdm->ops->ifunctionview; 125 funcloadstruct.ifunctionload = kdm->ops->ifunctionload; 126 PetscCall(PetscViewerBinaryWrite(viewer, &funcstruct, 1, PETSC_FUNCTION)); 127 PetscCall(PetscViewerBinaryWrite(viewer, &funcviewstruct, 1, PETSC_FUNCTION)); 128 PetscCall(PetscViewerBinaryWrite(viewer, &funcloadstruct, 1, PETSC_FUNCTION)); 129 if (kdm->ops->ifunctionview) { 130 void *ctx; 131 132 PetscCall(PetscContainerGetPointer(kdm->ifunctionctxcontainer, &ctx)); 133 PetscCall((*kdm->ops->ifunctionview)(ctx, viewer)); 134 } 135 jacstruct.ijacobian = kdm->ops->ijacobian; 136 jacviewstruct.ijacobianview = kdm->ops->ijacobianview; 137 jacloadstruct.ijacobianload = kdm->ops->ijacobianload; 138 PetscCall(PetscViewerBinaryWrite(viewer, &jacstruct, 1, PETSC_FUNCTION)); 139 PetscCall(PetscViewerBinaryWrite(viewer, &jacviewstruct, 1, PETSC_FUNCTION)); 140 PetscCall(PetscViewerBinaryWrite(viewer, &jacloadstruct, 1, PETSC_FUNCTION)); 141 if (kdm->ops->ijacobianview) { 142 void *ctx; 143 144 PetscCall(PetscContainerGetPointer(kdm->ijacobianctxcontainer, &ctx)); 145 PetscCall((*kdm->ops->ijacobianview)(ctx, viewer)); 146 } 147 } 148 PetscFunctionReturn(0); 149 } 150 151 static PetscErrorCode DMTSCreate(MPI_Comm comm, DMTS *kdm) { 152 PetscFunctionBegin; 153 PetscCall(TSInitializePackage()); 154 PetscCall(PetscHeaderCreate(*kdm, DMTS_CLASSID, "DMTS", "DMTS", "DMTS", comm, DMTSDestroy, DMTSView)); 155 PetscFunctionReturn(0); 156 } 157 158 /* Attaches the DMTS to the coarse level. 159 * Under what conditions should we copy versus duplicate? 160 */ 161 static PetscErrorCode DMCoarsenHook_DMTS(DM dm, DM dmc, void *ctx) { 162 PetscFunctionBegin; 163 PetscCall(DMCopyDMTS(dm, dmc)); 164 PetscFunctionReturn(0); 165 } 166 167 /* This could restrict auxiliary information to the coarse level. 168 */ 169 static PetscErrorCode DMRestrictHook_DMTS(DM dm, Mat Restrict, Vec rscale, Mat Inject, DM dmc, void *ctx) { 170 PetscFunctionBegin; 171 PetscFunctionReturn(0); 172 } 173 174 static PetscErrorCode DMSubDomainHook_DMTS(DM dm, DM subdm, void *ctx) { 175 PetscFunctionBegin; 176 PetscCall(DMCopyDMTS(dm, subdm)); 177 PetscFunctionReturn(0); 178 } 179 180 /* This could restrict auxiliary information to the coarse level. 181 */ 182 static PetscErrorCode DMSubDomainRestrictHook_DMTS(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx) { 183 PetscFunctionBegin; 184 PetscFunctionReturn(0); 185 } 186 187 /*@C 188 DMTSCopy - copies the information in a DMTS to another DMTS 189 190 Not Collective 191 192 Input Parameters: 193 + kdm - Original DMTS 194 - nkdm - DMTS to receive the data, should have been created with DMTSCreate() 195 196 Level: developer 197 198 .seealso: `DMTSCreate()`, `DMTSDestroy()` 199 @*/ 200 PetscErrorCode DMTSCopy(DMTS kdm, DMTS nkdm) { 201 PetscFunctionBegin; 202 PetscValidHeaderSpecific(kdm, DMTS_CLASSID, 1); 203 PetscValidHeaderSpecific(nkdm, DMTS_CLASSID, 2); 204 nkdm->ops->rhsfunction = kdm->ops->rhsfunction; 205 nkdm->ops->rhsjacobian = kdm->ops->rhsjacobian; 206 nkdm->ops->ifunction = kdm->ops->ifunction; 207 nkdm->ops->ijacobian = kdm->ops->ijacobian; 208 nkdm->ops->i2function = kdm->ops->i2function; 209 nkdm->ops->i2jacobian = kdm->ops->i2jacobian; 210 nkdm->ops->solution = kdm->ops->solution; 211 nkdm->ops->destroy = kdm->ops->destroy; 212 nkdm->ops->duplicate = kdm->ops->duplicate; 213 214 nkdm->solutionctx = kdm->solutionctx; 215 nkdm->rhsfunctionctxcontainer = kdm->rhsfunctionctxcontainer; 216 nkdm->rhsjacobianctxcontainer = kdm->rhsjacobianctxcontainer; 217 nkdm->ifunctionctxcontainer = kdm->ifunctionctxcontainer; 218 nkdm->ijacobianctxcontainer = kdm->ijacobianctxcontainer; 219 nkdm->i2functionctxcontainer = kdm->i2functionctxcontainer; 220 nkdm->i2jacobianctxcontainer = kdm->i2jacobianctxcontainer; 221 if (nkdm->rhsfunctionctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "rhs function ctx", (PetscObject)nkdm->rhsfunctionctxcontainer)); 222 if (nkdm->rhsjacobianctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "rhs jacobian ctx", (PetscObject)nkdm->rhsjacobianctxcontainer)); 223 if (nkdm->ifunctionctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "ifunction ctx", (PetscObject)nkdm->ifunctionctxcontainer)); 224 if (nkdm->ijacobianctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "ijacobian ctx", (PetscObject)nkdm->ijacobianctxcontainer)); 225 if (nkdm->i2functionctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "i2function ctx", (PetscObject)nkdm->i2functionctxcontainer)); 226 if (nkdm->i2jacobianctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "i2jacobian ctx", (PetscObject)nkdm->i2jacobianctxcontainer)); 227 228 nkdm->data = kdm->data; 229 230 /* 231 nkdm->fortran_func_pointers[0] = kdm->fortran_func_pointers[0]; 232 nkdm->fortran_func_pointers[1] = kdm->fortran_func_pointers[1]; 233 nkdm->fortran_func_pointers[2] = kdm->fortran_func_pointers[2]; 234 */ 235 236 /* implementation specific copy hooks */ 237 PetscTryTypeMethod(kdm, duplicate, nkdm); 238 PetscFunctionReturn(0); 239 } 240 241 /*@C 242 DMGetDMTS - get read-only private DMTS context from a DM 243 244 Not Collective 245 246 Input Parameter: 247 . dm - DM to be used with TS 248 249 Output Parameter: 250 . tsdm - private DMTS context 251 252 Level: developer 253 254 Notes: 255 Use DMGetDMTSWrite() if write access is needed. The DMTSSetXXX API should be used wherever possible. 256 257 .seealso: `DMGetDMTSWrite()` 258 @*/ 259 PetscErrorCode DMGetDMTS(DM dm, DMTS *tsdm) { 260 PetscFunctionBegin; 261 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 262 *tsdm = (DMTS)dm->dmts; 263 if (!*tsdm) { 264 PetscCall(PetscInfo(dm, "Creating new DMTS\n")); 265 PetscCall(DMTSCreate(PetscObjectComm((PetscObject)dm), tsdm)); 266 dm->dmts = (PetscObject)*tsdm; 267 (*tsdm)->originaldm = dm; 268 PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_DMTS, DMRestrictHook_DMTS, NULL)); 269 PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_DMTS, DMSubDomainRestrictHook_DMTS, NULL)); 270 } 271 PetscFunctionReturn(0); 272 } 273 274 /*@C 275 DMGetDMTSWrite - get write access to private DMTS context from a DM 276 277 Not Collective 278 279 Input Parameter: 280 . dm - DM to be used with TS 281 282 Output Parameter: 283 . tsdm - private DMTS context 284 285 Level: developer 286 287 .seealso: `DMGetDMTS()` 288 @*/ 289 PetscErrorCode DMGetDMTSWrite(DM dm, DMTS *tsdm) { 290 DMTS sdm; 291 292 PetscFunctionBegin; 293 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 294 PetscCall(DMGetDMTS(dm, &sdm)); 295 PetscCheck(sdm->originaldm, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "DMTS has a NULL originaldm"); 296 if (sdm->originaldm != dm) { /* Copy on write */ 297 DMTS oldsdm = sdm; 298 PetscCall(PetscInfo(dm, "Copying DMTS due to write\n")); 299 PetscCall(DMTSCreate(PetscObjectComm((PetscObject)dm), &sdm)); 300 PetscCall(DMTSCopy(oldsdm, sdm)); 301 PetscCall(DMTSDestroy((DMTS *)&dm->dmts)); 302 dm->dmts = (PetscObject)sdm; 303 sdm->originaldm = dm; 304 } 305 *tsdm = sdm; 306 PetscFunctionReturn(0); 307 } 308 309 /*@C 310 DMCopyDMTS - copies a DM context to a new DM 311 312 Logically Collective 313 314 Input Parameters: 315 + dmsrc - DM to obtain context from 316 - dmdest - DM to add context to 317 318 Level: developer 319 320 Note: 321 The context is copied by reference. This function does not ensure that a context exists. 322 323 .seealso: `DMGetDMTS()`, `TSSetDM()` 324 @*/ 325 PetscErrorCode DMCopyDMTS(DM dmsrc, DM dmdest) { 326 PetscFunctionBegin; 327 PetscValidHeaderSpecific(dmsrc, DM_CLASSID, 1); 328 PetscValidHeaderSpecific(dmdest, DM_CLASSID, 2); 329 PetscCall(DMTSDestroy((DMTS *)&dmdest->dmts)); 330 dmdest->dmts = dmsrc->dmts; 331 PetscCall(PetscObjectReference(dmdest->dmts)); 332 PetscCall(DMCoarsenHookAdd(dmdest, DMCoarsenHook_DMTS, DMRestrictHook_DMTS, NULL)); 333 PetscCall(DMSubDomainHookAdd(dmdest, DMSubDomainHook_DMTS, DMSubDomainRestrictHook_DMTS, NULL)); 334 PetscFunctionReturn(0); 335 } 336 337 /*@C 338 DMTSSetIFunction - set TS implicit function evaluation function 339 340 Not Collective 341 342 Input Parameters: 343 + dm - DM to be used with TS 344 . func - function evaluating f(t,u,u_t) 345 - ctx - context for residual evaluation 346 347 Calling sequence of func: 348 $ PetscErrorCode func(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx); 349 350 + t - time at step/stage being solved 351 . u - state vector 352 . u_t - time derivative of state vector 353 . F - function vector 354 - ctx - [optional] user-defined context for matrix evaluation routine 355 356 Level: advanced 357 358 Note: 359 TSSetFunction() is normally used, but it calls this function internally because the user context is actually 360 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 361 not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 362 363 .seealso: `DMTSSetContext()`, `TSSetIFunction()`, `DMTSSetJacobian()` 364 @*/ 365 PetscErrorCode DMTSSetIFunction(DM dm, TSIFunction func, void *ctx) { 366 DMTS tsdm; 367 368 PetscFunctionBegin; 369 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 370 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 371 if (func) tsdm->ops->ifunction = func; 372 if (ctx) { 373 PetscContainer ctxcontainer; 374 PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer)); 375 PetscCall(PetscContainerSetPointer(ctxcontainer, ctx)); 376 PetscCall(PetscObjectCompose((PetscObject)tsdm, "ifunction ctx", (PetscObject)ctxcontainer)); 377 tsdm->ifunctionctxcontainer = ctxcontainer; 378 PetscCall(PetscContainerDestroy(&ctxcontainer)); 379 } 380 PetscFunctionReturn(0); 381 } 382 383 /*@C 384 DMTSSetIFunctionContextDestroy - set `TS` implicit evaluation context destroy function 385 386 Not Collective 387 388 Input Parameters: 389 + dm - `DM` to be used with `TS` 390 - f - implicit evaluation context destroy function 391 392 Level: advanced 393 394 .seealso: `DMTSSetIFunction()`, `TSSetIFunction()` 395 @*/ 396 PetscErrorCode DMTSSetIFunctionContextDestroy(DM dm, PetscErrorCode (*f)(void *)) { 397 DMTS tsdm; 398 399 PetscFunctionBegin; 400 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 401 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 402 if (tsdm->ifunctionctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->ifunctionctxcontainer, f)); 403 PetscFunctionReturn(0); 404 } 405 406 PetscErrorCode DMTSUnsetIFunctionContext_Internal(DM dm) { 407 DMTS tsdm; 408 409 PetscFunctionBegin; 410 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 411 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 412 PetscCall(DMTSUnsetIFunctionContext_DMTS(tsdm)); 413 PetscFunctionReturn(0); 414 } 415 416 /*@C 417 DMTSGetIFunction - get TS implicit residual evaluation function 418 419 Not Collective 420 421 Input Parameter: 422 . dm - DM to be used with TS 423 424 Output Parameters: 425 + func - function evaluation function, see TSSetIFunction() for calling sequence 426 - ctx - context for residual evaluation 427 428 Level: advanced 429 430 Note: 431 TSGetFunction() is normally used, but it calls this function internally because the user context is actually 432 associated with the DM. 433 434 .seealso: `DMTSSetContext()`, `DMTSSetFunction()`, `TSSetFunction()` 435 @*/ 436 PetscErrorCode DMTSGetIFunction(DM dm, TSIFunction *func, void **ctx) { 437 DMTS tsdm; 438 439 PetscFunctionBegin; 440 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 441 PetscCall(DMGetDMTS(dm, &tsdm)); 442 if (func) *func = tsdm->ops->ifunction; 443 if (ctx) { 444 if (tsdm->ifunctionctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->ifunctionctxcontainer, ctx)); 445 else *ctx = NULL; 446 } 447 PetscFunctionReturn(0); 448 } 449 450 /*@C 451 DMTSSetI2Function - set TS implicit function evaluation function for 2nd order systems 452 453 Not Collective 454 455 Input Parameters: 456 + dm - DM to be used with TS 457 . fun - function evaluation routine 458 - ctx - context for residual evaluation 459 460 Calling sequence of fun: 461 $ PetscErrorCode fun(TS ts,PetscReal t,Vec U,Vec U_t,Vec U_tt,Vec F,ctx); 462 463 + t - time at step/stage being solved 464 . U - state vector 465 . U_t - time derivative of state vector 466 . U_tt - second time derivative of state vector 467 . F - function vector 468 - ctx - [optional] user-defined context for matrix evaluation routine (may be NULL) 469 470 Level: advanced 471 472 Note: 473 TSSetI2Function() is normally used, but it calls this function internally because the user context is actually 474 associated with the DM. 475 476 .seealso: `TSSetI2Function()` 477 @*/ 478 PetscErrorCode DMTSSetI2Function(DM dm, TSI2Function fun, void *ctx) { 479 DMTS tsdm; 480 481 PetscFunctionBegin; 482 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 483 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 484 if (fun) tsdm->ops->i2function = fun; 485 if (ctx) { 486 PetscContainer ctxcontainer; 487 PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer)); 488 PetscCall(PetscContainerSetPointer(ctxcontainer, ctx)); 489 PetscCall(PetscObjectCompose((PetscObject)tsdm, "i2function ctx", (PetscObject)ctxcontainer)); 490 tsdm->i2functionctxcontainer = ctxcontainer; 491 PetscCall(PetscContainerDestroy(&ctxcontainer)); 492 } 493 PetscFunctionReturn(0); 494 } 495 496 /*@C 497 DMTSSetI2FunctionContextDestroy - set `TS` implicit evaluation for 2nd order systems context destroy 498 499 Not Collective 500 501 Input Parameters: 502 + dm - `DM` to be used with `TS` 503 - f - implicit evaluation context destroy function 504 505 Level: advanced 506 507 Note: 508 `TSSetI2FunctionContextDestroy()` is normally used, but it calls this function internally because the user context is actually 509 associated with the `DM`. 510 511 .seealso: `TSSetI2FunctionContextDestroy()`, `DMTSSetI2Function()`, `TSSetI2Function()` 512 @*/ 513 PetscErrorCode DMTSSetI2FunctionContextDestroy(DM dm, PetscErrorCode (*f)(void *)) { 514 DMTS tsdm; 515 516 PetscFunctionBegin; 517 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 518 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 519 if (tsdm->i2functionctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->i2functionctxcontainer, f)); 520 PetscFunctionReturn(0); 521 } 522 523 PetscErrorCode DMTSUnsetI2FunctionContext_Internal(DM dm) { 524 DMTS tsdm; 525 526 PetscFunctionBegin; 527 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 528 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 529 PetscCall(DMTSUnsetI2FunctionContext_DMTS(tsdm)); 530 PetscFunctionReturn(0); 531 } 532 533 /*@C 534 DMTSGetI2Function - get TS implicit residual evaluation function for 2nd order systems 535 536 Not Collective 537 538 Input Parameter: 539 . dm - DM to be used with TS 540 541 Output Parameters: 542 + fun - function evaluation function, see TSSetI2Function() for calling sequence 543 - ctx - context for residual evaluation 544 545 Level: advanced 546 547 Note: 548 TSGetI2Function() is normally used, but it calls this function internally because the user context is actually 549 associated with the DM. 550 551 .seealso: `DMTSSetI2Function()`, `TSGetI2Function()` 552 @*/ 553 PetscErrorCode DMTSGetI2Function(DM dm, TSI2Function *fun, void **ctx) { 554 DMTS tsdm; 555 556 PetscFunctionBegin; 557 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 558 PetscCall(DMGetDMTS(dm, &tsdm)); 559 if (fun) *fun = tsdm->ops->i2function; 560 if (ctx) { 561 if (tsdm->i2functionctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->i2functionctxcontainer, ctx)); 562 else *ctx = NULL; 563 } 564 PetscFunctionReturn(0); 565 } 566 567 /*@C 568 DMTSSetI2Jacobian - set TS implicit Jacobian evaluation function for 2nd order systems 569 570 Not Collective 571 572 Input Parameters: 573 + dm - DM to be used with TS 574 . fun - Jacobian evaluation routine 575 - ctx - context for Jacobian evaluation 576 577 Calling sequence of jac: 578 $ PetscErrorCode jac(TS ts,PetscReal t,Vec U,Vec U_t,Vec U_tt,PetscReal v,PetscReal a,Mat J,Mat P,void *ctx); 579 580 + t - time at step/stage being solved 581 . U - state vector 582 . U_t - time derivative of state vector 583 . U_tt - second time derivative of state vector 584 . v - shift for U_t 585 . a - shift for U_tt 586 . J - Jacobian of G(U) = F(t,U,W+v*U,W'+a*U), equivalent to dF/dU + v*dF/dU_t + a*dF/dU_tt 587 . P - preconditioning matrix for J, may be same as J 588 - ctx - [optional] user-defined context for matrix evaluation routine 589 590 Level: advanced 591 592 Note: 593 TSSetI2Jacobian() is normally used, but it calls this function internally because the user context is actually 594 associated with the DM. 595 596 .seealso: `TSSetI2Jacobian()` 597 @*/ 598 PetscErrorCode DMTSSetI2Jacobian(DM dm, TSI2Jacobian jac, void *ctx) { 599 DMTS tsdm; 600 601 PetscFunctionBegin; 602 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 603 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 604 if (jac) tsdm->ops->i2jacobian = jac; 605 if (ctx) { 606 PetscContainer ctxcontainer; 607 PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer)); 608 PetscCall(PetscContainerSetPointer(ctxcontainer, ctx)); 609 PetscCall(PetscObjectCompose((PetscObject)tsdm, "i2jacobian ctx", (PetscObject)ctxcontainer)); 610 tsdm->i2jacobianctxcontainer = ctxcontainer; 611 PetscCall(PetscContainerDestroy(&ctxcontainer)); 612 } 613 PetscFunctionReturn(0); 614 } 615 616 /*@C 617 DMTSSetI2JacobianContextDestroy - set `TS` implicit Jacobian evaluation for 2nd order systems context destroy function 618 619 Not Collective 620 621 Input Parameters: 622 + dm - `DM` to be used with `TS` 623 - f - implicit Jacobian evaluation context destroy function 624 625 Level: advanced 626 627 .seealso: `TSSetI2JacobianContextDestroy()`, `DMTSSetI2Jacobian()`, `TSSetI2Jacobian()` 628 @*/ 629 PetscErrorCode DMTSSetI2JacobianContextDestroy(DM dm, PetscErrorCode (*f)(void *)) { 630 DMTS tsdm; 631 632 PetscFunctionBegin; 633 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 634 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 635 if (tsdm->i2jacobianctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->i2jacobianctxcontainer, f)); 636 PetscFunctionReturn(0); 637 } 638 639 PetscErrorCode DMTSUnsetI2JacobianContext_Internal(DM dm) { 640 DMTS tsdm; 641 642 PetscFunctionBegin; 643 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 644 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 645 PetscCall(DMTSUnsetI2JacobianContext_DMTS(tsdm)); 646 PetscFunctionReturn(0); 647 } 648 649 /*@C 650 DMTSGetI2Jacobian - get TS implicit Jacobian evaluation function for 2nd order systems 651 652 Not Collective 653 654 Input Parameter: 655 . dm - DM to be used with TS 656 657 Output Parameters: 658 + jac - Jacobian evaluation function, see TSSetI2Jacobian() for calling sequence 659 - ctx - context for Jacobian evaluation 660 661 Level: advanced 662 663 Note: 664 TSGetI2Jacobian() is normally used, but it calls this function internally because the user context is actually 665 associated with the DM. 666 667 .seealso: `DMTSSetI2Jacobian()`, `TSGetI2Jacobian()` 668 @*/ 669 PetscErrorCode DMTSGetI2Jacobian(DM dm, TSI2Jacobian *jac, void **ctx) { 670 DMTS tsdm; 671 672 PetscFunctionBegin; 673 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 674 PetscCall(DMGetDMTS(dm, &tsdm)); 675 if (jac) *jac = tsdm->ops->i2jacobian; 676 if (ctx) { 677 if (tsdm->i2jacobianctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->i2jacobianctxcontainer, ctx)); 678 else *ctx = NULL; 679 } 680 PetscFunctionReturn(0); 681 } 682 683 /*@C 684 DMTSSetRHSFunction - set TS explicit residual evaluation function 685 686 Not Collective 687 688 Input Parameters: 689 + dm - DM to be used with TS 690 . func - RHS function evaluation routine 691 - ctx - context for residual evaluation 692 693 Calling sequence of func: 694 $ PetscErrorCode func(TS ts,PetscReal t,Vec u,Vec F,void *ctx); 695 696 + ts - timestep context 697 . t - current timestep 698 . u - input vector 699 . F - function vector 700 - ctx - [optional] user-defined function context 701 702 Level: advanced 703 704 Note: 705 TSSetRHSFunction() is normally used, but it calls this function internally because the user context is actually 706 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 707 not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 708 709 .seealso: `DMTSSetContext()`, `TSSetRHSFunction()`, `DMTSSetJacobian()` 710 @*/ 711 PetscErrorCode DMTSSetRHSFunction(DM dm, TSRHSFunction func, void *ctx) { 712 DMTS tsdm; 713 714 PetscFunctionBegin; 715 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 716 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 717 if (func) tsdm->ops->rhsfunction = func; 718 if (ctx) { 719 PetscContainer ctxcontainer; 720 PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer)); 721 PetscCall(PetscContainerSetPointer(ctxcontainer, ctx)); 722 PetscCall(PetscObjectCompose((PetscObject)tsdm, "rhs function ctx", (PetscObject)ctxcontainer)); 723 tsdm->rhsfunctionctxcontainer = ctxcontainer; 724 PetscCall(PetscContainerDestroy(&ctxcontainer)); 725 } 726 PetscFunctionReturn(0); 727 } 728 729 /*@C 730 DMTSSetRHSFunctionContextDestroy - set `TS` explicit residual evaluation context destroy function 731 732 Not Collective 733 734 Input Parameters: 735 + dm - `DM` to be used with `TS` 736 - f - explicit evaluation context destroy function 737 738 Level: advanced 739 740 Note: 741 `TSSetRHSFunctionContextDestroy()` is normally used, but it calls this function internally because the user context is actually 742 associated with the `DM`. This makes the interface consistent regardless of whether the user interacts with a `DM` or 743 not. 744 745 Developer Note: 746 If `DM` took a more central role at some later date, this could become the primary method of setting the residual. 747 748 .seealso: `TSSetRHSFunctionContextDestroy()`, `DMTSSetRHSFunction()`, `TSSetRHSFunction()` 749 @*/ 750 PetscErrorCode DMTSSetRHSFunctionContextDestroy(DM dm, PetscErrorCode (*f)(void *)) { 751 DMTS tsdm; 752 753 PetscFunctionBegin; 754 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 755 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 756 if (tsdm->rhsfunctionctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->rhsfunctionctxcontainer, f)); 757 PetscFunctionReturn(0); 758 } 759 760 PetscErrorCode DMTSUnsetRHSFunctionContext_Internal(DM dm) { 761 DMTS tsdm; 762 763 PetscFunctionBegin; 764 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 765 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 766 PetscCall(DMTSUnsetRHSFunctionContext_DMTS(tsdm)); 767 tsdm->rhsfunctionctxcontainer = NULL; 768 PetscFunctionReturn(0); 769 } 770 771 /*@C 772 DMTSSetTransientVariable - sets function to transform from state to transient variables 773 774 Logically Collective 775 776 Input Parameters: 777 + dm - DM to be used with TS 778 . tvar - a function that transforms to transient variables 779 - ctx - a context for tvar 780 781 Calling sequence of tvar: 782 $ PetscErrorCode tvar(TS ts,Vec p,Vec c,void *ctx); 783 784 + ts - timestep context 785 . p - input vector (primitive form) 786 . c - output vector, transient variables (conservative form) 787 - ctx - [optional] user-defined function context 788 789 Level: advanced 790 791 Notes: 792 This is typically used to transform from primitive to conservative variables so that a time integrator (e.g., TSBDF) 793 can be conservative. In this context, primitive variables P are used to model the state (e.g., because they lead to 794 well-conditioned formulations even in limiting cases such as low-Mach or zero porosity). The transient variable is 795 C(P), specified by calling this function. An IFunction thus receives arguments (P, Cdot) and the IJacobian must be 796 evaluated via the chain rule, as in 797 798 dF/dP + shift * dF/dCdot dC/dP. 799 800 .seealso: `TSSetTransientVariable()`, `DMTSGetTransientVariable()`, `DMTSSetIFunction()`, `DMTSSetIJacobian()` 801 @*/ 802 PetscErrorCode DMTSSetTransientVariable(DM dm, TSTransientVariable tvar, void *ctx) { 803 DMTS dmts; 804 805 PetscFunctionBegin; 806 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 807 PetscCall(DMGetDMTSWrite(dm, &dmts)); 808 dmts->ops->transientvar = tvar; 809 dmts->transientvarctx = ctx; 810 PetscFunctionReturn(0); 811 } 812 813 /*@C 814 DMTSGetTransientVariable - gets function to transform from state to transient variables 815 816 Logically Collective 817 818 Input Parameter: 819 . dm - DM to be used with TS 820 821 Output Parameters: 822 + tvar - a function that transforms to transient variables 823 - ctx - a context for tvar 824 825 Level: advanced 826 827 .seealso: `DMTSSetTransientVariable()`, `DMTSGetIFunction()`, `DMTSGetIJacobian()` 828 @*/ 829 PetscErrorCode DMTSGetTransientVariable(DM dm, TSTransientVariable *tvar, void *ctx) { 830 DMTS dmts; 831 832 PetscFunctionBegin; 833 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 834 PetscCall(DMGetDMTS(dm, &dmts)); 835 if (tvar) *tvar = dmts->ops->transientvar; 836 if (ctx) *(void **)ctx = dmts->transientvarctx; 837 PetscFunctionReturn(0); 838 } 839 840 /*@C 841 DMTSGetSolutionFunction - gets the TS solution evaluation function 842 843 Not Collective 844 845 Input Parameter: 846 . dm - DM to be used with TS 847 848 Output Parameters: 849 + func - solution function evaluation function, see TSSetSolution() for calling sequence 850 - ctx - context for solution evaluation 851 852 Level: advanced 853 854 .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `DMTSSetSolutionFunction()` 855 @*/ 856 PetscErrorCode DMTSGetSolutionFunction(DM dm, TSSolutionFunction *func, void **ctx) { 857 DMTS tsdm; 858 859 PetscFunctionBegin; 860 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 861 PetscCall(DMGetDMTS(dm, &tsdm)); 862 if (func) *func = tsdm->ops->solution; 863 if (ctx) *ctx = tsdm->solutionctx; 864 PetscFunctionReturn(0); 865 } 866 867 /*@C 868 DMTSSetSolutionFunction - set TS solution evaluation function 869 870 Not Collective 871 872 Input Parameters: 873 + dm - DM to be used with TS 874 . func - solution function evaluation routine 875 - ctx - context for solution evaluation 876 877 Calling sequence of f: 878 $ PetscErrorCode f(TS ts,PetscReal t,Vec u,void *ctx); 879 880 + ts - timestep context 881 . t - current timestep 882 . u - output vector 883 - ctx - [optional] user-defined function context 884 885 Level: advanced 886 887 Note: 888 TSSetSolutionFunction() is normally used, but it calls this function internally because the user context is actually 889 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 890 not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 891 892 .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `DMTSGetSolutionFunction()` 893 @*/ 894 PetscErrorCode DMTSSetSolutionFunction(DM dm, TSSolutionFunction func, void *ctx) { 895 DMTS tsdm; 896 897 PetscFunctionBegin; 898 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 899 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 900 if (func) tsdm->ops->solution = func; 901 if (ctx) tsdm->solutionctx = ctx; 902 PetscFunctionReturn(0); 903 } 904 905 /*@C 906 DMTSSetForcingFunction - set TS forcing function evaluation function 907 908 Not Collective 909 910 Input Parameters: 911 + dm - DM to be used with TS 912 . f - forcing function evaluation routine 913 - ctx - context for solution evaluation 914 915 Calling sequence of func: 916 $ PetscErrorCode func (TS ts,PetscReal t,Vec f,void *ctx); 917 918 + ts - timestep context 919 . t - current timestep 920 . f - output vector 921 - ctx - [optional] user-defined function context 922 923 Level: advanced 924 925 Note: 926 TSSetForcingFunction() is normally used, but it calls this function internally because the user context is actually 927 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 928 not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 929 930 .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `TSSetForcingFunction()`, `DMTSGetForcingFunction()` 931 @*/ 932 PetscErrorCode DMTSSetForcingFunction(DM dm, TSForcingFunction f, void *ctx) { 933 DMTS tsdm; 934 935 PetscFunctionBegin; 936 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 937 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 938 if (f) tsdm->ops->forcing = f; 939 if (ctx) tsdm->forcingctx = ctx; 940 PetscFunctionReturn(0); 941 } 942 943 /*@C 944 DMTSGetForcingFunction - get TS forcing function evaluation function 945 946 Not Collective 947 948 Input Parameter: 949 . dm - DM to be used with TS 950 951 Output Parameters: 952 + f - forcing function evaluation function; see TSForcingFunction for details 953 - ctx - context for solution evaluation 954 955 Level: advanced 956 957 Note: 958 TSSetForcingFunction() is normally used, but it calls this function internally because the user context is actually 959 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 960 not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 961 962 .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()`, `TSSetForcingFunction()`, `DMTSGetForcingFunction()` 963 @*/ 964 PetscErrorCode DMTSGetForcingFunction(DM dm, TSForcingFunction *f, void **ctx) { 965 DMTS tsdm; 966 967 PetscFunctionBegin; 968 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 969 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 970 if (f) *f = tsdm->ops->forcing; 971 if (ctx) *ctx = tsdm->forcingctx; 972 PetscFunctionReturn(0); 973 } 974 975 /*@C 976 DMTSGetRHSFunction - get TS explicit residual evaluation function 977 978 Not Collective 979 980 Input Parameter: 981 . dm - DM to be used with TS 982 983 Output Parameters: 984 + func - residual evaluation function, see TSSetRHSFunction() for calling sequence 985 - ctx - context for residual evaluation 986 987 Level: advanced 988 989 Note: 990 TSGetFunction() is normally used, but it calls this function internally because the user context is actually 991 associated with the DM. 992 993 .seealso: `DMTSSetContext()`, `DMTSSetRHSFunction()`, `TSSetRHSFunction()` 994 @*/ 995 PetscErrorCode DMTSGetRHSFunction(DM dm, TSRHSFunction *func, void **ctx) { 996 DMTS tsdm; 997 998 PetscFunctionBegin; 999 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1000 PetscCall(DMGetDMTS(dm, &tsdm)); 1001 if (func) *func = tsdm->ops->rhsfunction; 1002 if (ctx) { 1003 if (tsdm->rhsfunctionctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->rhsfunctionctxcontainer, ctx)); 1004 else *ctx = NULL; 1005 } 1006 PetscFunctionReturn(0); 1007 } 1008 1009 /*@C 1010 DMTSSetIJacobian - set TS Jacobian evaluation function 1011 1012 Not Collective 1013 1014 Input Parameters: 1015 + dm - DM to be used with TS 1016 . func - Jacobian evaluation routine 1017 - ctx - context for residual evaluation 1018 1019 Calling sequence of f: 1020 $ PetscErrorCode f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat Amat,Mat Pmat,void *ctx); 1021 1022 + t - time at step/stage being solved 1023 . U - state vector 1024 . U_t - time derivative of state vector 1025 . a - shift 1026 . Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t 1027 . Pmat - matrix used for constructing preconditioner, usually the same as Amat 1028 - ctx - [optional] user-defined context for matrix evaluation routine 1029 1030 Level: advanced 1031 1032 Note: 1033 TSSetJacobian() is normally used, but it calls this function internally because the user context is actually 1034 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 1035 not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 1036 1037 .seealso: `DMTSSetContext()`, `TSSetRHSFunction()`, `DMTSGetJacobian()`, `TSSetIJacobian()`, `TSSetIFunction()` 1038 @*/ 1039 PetscErrorCode DMTSSetIJacobian(DM dm, TSIJacobian func, void *ctx) { 1040 DMTS tsdm; 1041 1042 PetscFunctionBegin; 1043 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1044 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 1045 if (func) tsdm->ops->ijacobian = func; 1046 if (ctx) { 1047 PetscContainer ctxcontainer; 1048 PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer)); 1049 PetscCall(PetscContainerSetPointer(ctxcontainer, ctx)); 1050 PetscCall(PetscObjectCompose((PetscObject)tsdm, "ijacobian ctx", (PetscObject)ctxcontainer)); 1051 tsdm->ijacobianctxcontainer = ctxcontainer; 1052 PetscCall(PetscContainerDestroy(&ctxcontainer)); 1053 } 1054 PetscFunctionReturn(0); 1055 } 1056 1057 /*@C 1058 DMTSSetIJacobianContextDestroy - set `TS` Jacobian evaluation context destroy function 1059 1060 Not Collective 1061 1062 Input Parameters: 1063 + dm - `DM` to be used with `TS` 1064 - f - Jacobian evaluation context destroy function 1065 1066 Level: advanced 1067 1068 Note: 1069 `TSSetIJacobianContextDestroy()` is normally used, but it calls this function internally because the user context is actually 1070 associated with the `DM`. This makes the interface consistent regardless of whether the user interacts with a `DM` or 1071 not. 1072 1073 Developer Note: 1074 If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian. 1075 1076 .seealso: `TSSetIJacobianContextDestroy()`, `TSSetI2JacobianContextDestroy()`, `DMTSSetIJacobian()`, `TSSetIJacobian()` 1077 @*/ 1078 PetscErrorCode DMTSSetIJacobianContextDestroy(DM dm, PetscErrorCode (*f)(void *)) { 1079 DMTS tsdm; 1080 1081 PetscFunctionBegin; 1082 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1083 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 1084 if (tsdm->ijacobianctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->ijacobianctxcontainer, f)); 1085 PetscFunctionReturn(0); 1086 } 1087 1088 PetscErrorCode DMTSUnsetIJacobianContext_Internal(DM dm) { 1089 DMTS tsdm; 1090 1091 PetscFunctionBegin; 1092 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1093 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 1094 PetscCall(DMTSUnsetIJacobianContext_DMTS(tsdm)); 1095 PetscFunctionReturn(0); 1096 } 1097 1098 /*@C 1099 DMTSGetIJacobian - get TS Jacobian evaluation function 1100 1101 Not Collective 1102 1103 Input Parameter: 1104 . dm - DM to be used with TS 1105 1106 Output Parameters: 1107 + func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence 1108 - ctx - context for residual evaluation 1109 1110 Level: advanced 1111 1112 Note: 1113 TSGetJacobian() is normally used, but it calls this function internally because the user context is actually 1114 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 1115 not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 1116 1117 .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()` 1118 @*/ 1119 PetscErrorCode DMTSGetIJacobian(DM dm, TSIJacobian *func, void **ctx) { 1120 DMTS tsdm; 1121 1122 PetscFunctionBegin; 1123 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1124 PetscCall(DMGetDMTS(dm, &tsdm)); 1125 if (func) *func = tsdm->ops->ijacobian; 1126 if (ctx) { 1127 if (tsdm->ijacobianctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->ijacobianctxcontainer, ctx)); 1128 else *ctx = NULL; 1129 } 1130 PetscFunctionReturn(0); 1131 } 1132 1133 /*@C 1134 DMTSSetRHSJacobian - set TS Jacobian evaluation function 1135 1136 Not Collective 1137 1138 Input Parameters: 1139 + dm - DM to be used with TS 1140 . func - Jacobian evaluation routine 1141 - ctx - context for residual evaluation 1142 1143 Calling sequence of func: 1144 $ PetscErrorCode func(TS ts,PetscReal t,Vec u,Mat A,Mat B,void *ctx); 1145 1146 + t - current timestep 1147 . u - input vector 1148 . Amat - (approximate) Jacobian matrix 1149 . Pmat - matrix from which preconditioner is to be constructed (usually the same as Amat) 1150 - ctx - [optional] user-defined context for matrix evaluation routine 1151 1152 Level: advanced 1153 1154 Note: 1155 TSSetJacobian() is normally used, but it calls this function internally because the user context is actually 1156 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 1157 not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 1158 1159 .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSGetJacobian()`, `TSSetRHSJacobian()` 1160 @*/ 1161 PetscErrorCode DMTSSetRHSJacobian(DM dm, TSRHSJacobian func, void *ctx) { 1162 DMTS tsdm; 1163 1164 PetscFunctionBegin; 1165 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1166 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 1167 if (func) tsdm->ops->rhsjacobian = func; 1168 if (ctx) { 1169 PetscContainer ctxcontainer; 1170 PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer)); 1171 PetscCall(PetscContainerSetPointer(ctxcontainer, ctx)); 1172 PetscCall(PetscObjectCompose((PetscObject)tsdm, "rhs jacobian ctx", (PetscObject)ctxcontainer)); 1173 tsdm->rhsjacobianctxcontainer = ctxcontainer; 1174 PetscCall(PetscContainerDestroy(&ctxcontainer)); 1175 } 1176 PetscFunctionReturn(0); 1177 } 1178 1179 /*@C 1180 DMTSSetRHSJacobianContextDestroy - set `TS` Jacobian evaluation context destroy function 1181 1182 Not Collective 1183 1184 Input Parameters: 1185 + dm - `DM` to be used with `TS` 1186 - f - Jacobian evaluation context destroy function 1187 1188 Level: advanced 1189 1190 Note: 1191 The user usually calls `TSSetRHSJacobianContextDestroy()` which calls this routine 1192 1193 Level: advanced 1194 1195 .seealso: `TSSetRHSJacobianContextDestroy()`, `DMTSSetRHSJacobian()`, `TSSetRHSJacobian()` 1196 @*/ 1197 PetscErrorCode DMTSSetRHSJacobianContextDestroy(DM dm, PetscErrorCode (*f)(void *)) { 1198 DMTS tsdm; 1199 1200 PetscFunctionBegin; 1201 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1202 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 1203 if (tsdm->rhsjacobianctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->rhsjacobianctxcontainer, f)); 1204 PetscFunctionReturn(0); 1205 } 1206 1207 PetscErrorCode DMTSUnsetRHSJacobianContext_Internal(DM dm) { 1208 DMTS tsdm; 1209 1210 PetscFunctionBegin; 1211 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1212 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 1213 PetscCall(DMTSUnsetRHSJacobianContext_DMTS(tsdm)); 1214 PetscFunctionReturn(0); 1215 } 1216 1217 /*@C 1218 DMTSGetRHSJacobian - get TS Jacobian evaluation function 1219 1220 Not Collective 1221 1222 Input Parameter: 1223 . dm - DM to be used with TS 1224 1225 Output Parameters: 1226 + func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence 1227 - ctx - context for residual evaluation 1228 1229 Level: advanced 1230 1231 Note: 1232 TSGetJacobian() is normally used, but it calls this function internally because the user context is actually 1233 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 1234 not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 1235 1236 .seealso: `DMTSSetContext()`, `TSSetRHSFunction()`, `DMTSSetRHSJacobian()`, `TSSetRHSJacobian()` 1237 @*/ 1238 PetscErrorCode DMTSGetRHSJacobian(DM dm, TSRHSJacobian *func, void **ctx) { 1239 DMTS tsdm; 1240 1241 PetscFunctionBegin; 1242 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1243 PetscCall(DMGetDMTS(dm, &tsdm)); 1244 if (func) *func = tsdm->ops->rhsjacobian; 1245 if (ctx) { 1246 if (tsdm->rhsjacobianctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->rhsjacobianctxcontainer, ctx)); 1247 else *ctx = NULL; 1248 } 1249 PetscFunctionReturn(0); 1250 } 1251 1252 /*@C 1253 DMTSSetIFunctionSerialize - sets functions used to view and load a IFunction context 1254 1255 Not Collective 1256 1257 Input Parameters: 1258 + dm - DM to be used with TS 1259 . view - viewer function 1260 - load - loading function 1261 1262 Level: advanced 1263 1264 .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()` 1265 @*/ 1266 PetscErrorCode DMTSSetIFunctionSerialize(DM dm, PetscErrorCode (*view)(void *, PetscViewer), PetscErrorCode (*load)(void **, PetscViewer)) { 1267 DMTS tsdm; 1268 1269 PetscFunctionBegin; 1270 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1271 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 1272 tsdm->ops->ifunctionview = view; 1273 tsdm->ops->ifunctionload = load; 1274 PetscFunctionReturn(0); 1275 } 1276 1277 /*@C 1278 DMTSSetIJacobianSerialize - sets functions used to view and load a IJacobian context 1279 1280 Not Collective 1281 1282 Input Parameters: 1283 + dm - DM to be used with TS 1284 . view - viewer function 1285 - load - loading function 1286 1287 Level: advanced 1288 1289 .seealso: `DMTSSetContext()`, `TSSetFunction()`, `DMTSSetJacobian()` 1290 @*/ 1291 PetscErrorCode DMTSSetIJacobianSerialize(DM dm, PetscErrorCode (*view)(void *, PetscViewer), PetscErrorCode (*load)(void **, PetscViewer)) { 1292 DMTS tsdm; 1293 1294 PetscFunctionBegin; 1295 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1296 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 1297 tsdm->ops->ijacobianview = view; 1298 tsdm->ops->ijacobianload = load; 1299 PetscFunctionReturn(0); 1300 } 1301