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