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(PETSC_SUCCESS); 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(PETSC_SUCCESS); 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(PETSC_SUCCESS); 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(PETSC_SUCCESS); 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(PETSC_SUCCESS); 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(PETSC_SUCCESS); 50 } 51 52 static PetscErrorCode DMTSDestroy(DMTS *kdm) 53 { 54 PetscFunctionBegin; 55 if (!*kdm) PetscFunctionReturn(PETSC_SUCCESS); 56 PetscValidHeaderSpecific((*kdm), DMTS_CLASSID, 1); 57 if (--((PetscObject)(*kdm))->refct > 0) { 58 *kdm = NULL; 59 PetscFunctionReturn(PETSC_SUCCESS); 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(PETSC_SUCCESS); 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(PETSC_SUCCESS); 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(PETSC_SUCCESS); 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(PETSC_SUCCESS); 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(PETSC_SUCCESS); 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(PETSC_SUCCESS); 184 } 185 186 static PetscErrorCode DMSubDomainHook_DMTS(DM dm, DM subdm, void *ctx) 187 { 188 PetscFunctionBegin; 189 PetscCall(DMCopyDMTS(dm, subdm)); 190 PetscFunctionReturn(PETSC_SUCCESS); 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(PETSC_SUCCESS); 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: [](ch_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(PETSC_SUCCESS); 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: [](ch_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(PETSC_SUCCESS); 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: [](ch_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(PETSC_SUCCESS); 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: [](ch_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(PETSC_SUCCESS); 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 Level: advanced 366 367 Note: 368 `TSSetIFunction()` is normally used, but it calls this function internally because the user context is actually 369 associated with the `DM`. This makes the interface consistent regardless of whether the user interacts with a `DM` or 370 not. If `DM` took a more central role at some later date, this could become the primary method of setting the residual. 371 372 .seealso: [](ch_ts), `TS`, `DM`, `TSIFunction` 373 @*/ 374 PetscErrorCode DMTSSetIFunction(DM dm, TSIFunction func, void *ctx) 375 { 376 DMTS tsdm; 377 378 PetscFunctionBegin; 379 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 380 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 381 if (func) tsdm->ops->ifunction = func; 382 if (ctx) { 383 PetscContainer ctxcontainer; 384 PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer)); 385 PetscCall(PetscContainerSetPointer(ctxcontainer, ctx)); 386 PetscCall(PetscObjectCompose((PetscObject)tsdm, "ifunction ctx", (PetscObject)ctxcontainer)); 387 tsdm->ifunctionctxcontainer = ctxcontainer; 388 PetscCall(PetscContainerDestroy(&ctxcontainer)); 389 } 390 PetscFunctionReturn(PETSC_SUCCESS); 391 } 392 393 /*@C 394 DMTSSetIFunctionContextDestroy - set `TS` implicit evaluation context destroy function 395 396 Not Collective 397 398 Input Parameters: 399 + dm - `DM` to be used with `TS` 400 - f - implicit evaluation context destroy function 401 402 Level: advanced 403 404 .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetIFunction()`, `TSSetIFunction()` 405 @*/ 406 PetscErrorCode DMTSSetIFunctionContextDestroy(DM dm, PetscErrorCode (*f)(void *)) 407 { 408 DMTS tsdm; 409 410 PetscFunctionBegin; 411 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 412 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 413 if (tsdm->ifunctionctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->ifunctionctxcontainer, f)); 414 PetscFunctionReturn(PETSC_SUCCESS); 415 } 416 417 PetscErrorCode DMTSUnsetIFunctionContext_Internal(DM dm) 418 { 419 DMTS tsdm; 420 421 PetscFunctionBegin; 422 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 423 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 424 PetscCall(DMTSUnsetIFunctionContext_DMTS(tsdm)); 425 PetscFunctionReturn(PETSC_SUCCESS); 426 } 427 428 /*@C 429 DMTSGetIFunction - get `TS` implicit residual evaluation function 430 431 Not Collective 432 433 Input Parameter: 434 . dm - `DM` to be used with `TS` 435 436 Output Parameters: 437 + func - function evaluation function, for calling sequence see `TSSetIFunction()` 438 - ctx - context for residual evaluation 439 440 Level: advanced 441 442 Note: 443 `TSGetIFunction()` is normally used, but it calls this function internally because the user context is actually 444 associated with the `DM`. 445 446 .seealso: [](ch_ts), `TS`, `DM`, `DMTSSetIFunction()` 447 @*/ 448 PetscErrorCode DMTSGetIFunction(DM dm, TSIFunction *func, void **ctx) 449 { 450 DMTS tsdm; 451 452 PetscFunctionBegin; 453 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 454 PetscCall(DMGetDMTS(dm, &tsdm)); 455 if (func) *func = tsdm->ops->ifunction; 456 if (ctx) { 457 if (tsdm->ifunctionctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->ifunctionctxcontainer, ctx)); 458 else *ctx = NULL; 459 } 460 PetscFunctionReturn(PETSC_SUCCESS); 461 } 462 463 /*@C 464 DMTSSetI2Function - set `TS` implicit function evaluation function for 2nd order systems 465 466 Not Collective 467 468 Input Parameters: 469 + dm - `DM` to be used with `TS` 470 . fun - function evaluation routine 471 - ctx - context for residual evaluation 472 473 Level: advanced 474 475 Note: 476 `TSSetI2Function()` is normally used, but it calls this function internally because the user context is actually 477 associated with the `DM`. 478 479 .seealso: [](ch_ts), `DM`, `TS`, `TSSetI2Function()` 480 @*/ 481 PetscErrorCode DMTSSetI2Function(DM dm, TSI2Function fun, void *ctx) 482 { 483 DMTS tsdm; 484 485 PetscFunctionBegin; 486 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 487 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 488 if (fun) tsdm->ops->i2function = fun; 489 if (ctx) { 490 PetscContainer ctxcontainer; 491 PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer)); 492 PetscCall(PetscContainerSetPointer(ctxcontainer, ctx)); 493 PetscCall(PetscObjectCompose((PetscObject)tsdm, "i2function ctx", (PetscObject)ctxcontainer)); 494 tsdm->i2functionctxcontainer = ctxcontainer; 495 PetscCall(PetscContainerDestroy(&ctxcontainer)); 496 } 497 PetscFunctionReturn(PETSC_SUCCESS); 498 } 499 500 /*@C 501 DMTSSetI2FunctionContextDestroy - set `TS` implicit evaluation for 2nd order systems context destroy 502 503 Not Collective 504 505 Input Parameters: 506 + dm - `DM` to be used with `TS` 507 - f - implicit evaluation context destroy function 508 509 Level: advanced 510 511 Note: 512 `TSSetI2FunctionContextDestroy()` is normally used, but it calls this function internally because the user context is actually 513 associated with the `DM`. 514 515 .seealso: [](ch_ts), `TSSetI2FunctionContextDestroy()`, `DMTSSetI2Function()`, `TSSetI2Function()` 516 @*/ 517 PetscErrorCode DMTSSetI2FunctionContextDestroy(DM dm, PetscErrorCode (*f)(void *)) 518 { 519 DMTS tsdm; 520 521 PetscFunctionBegin; 522 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 523 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 524 if (tsdm->i2functionctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->i2functionctxcontainer, f)); 525 PetscFunctionReturn(PETSC_SUCCESS); 526 } 527 528 PetscErrorCode DMTSUnsetI2FunctionContext_Internal(DM dm) 529 { 530 DMTS tsdm; 531 532 PetscFunctionBegin; 533 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 534 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 535 PetscCall(DMTSUnsetI2FunctionContext_DMTS(tsdm)); 536 PetscFunctionReturn(PETSC_SUCCESS); 537 } 538 539 /*@C 540 DMTSGetI2Function - get `TS` implicit residual evaluation function for 2nd order systems 541 542 Not Collective 543 544 Input Parameter: 545 . dm - `DM` to be used with `TS` 546 547 Output Parameters: 548 + fun - function evaluation function, for calling sequence see `TSSetI2Function()` 549 - ctx - context for residual evaluation 550 551 Level: advanced 552 553 Note: 554 `TSGetI2Function()` is normally used, but it calls this function internally because the user context is actually 555 associated with the `DM`. 556 557 .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetI2Function()`, `TSGetI2Function()` 558 @*/ 559 PetscErrorCode DMTSGetI2Function(DM dm, TSI2Function *fun, void **ctx) 560 { 561 DMTS tsdm; 562 563 PetscFunctionBegin; 564 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 565 PetscCall(DMGetDMTS(dm, &tsdm)); 566 if (fun) *fun = tsdm->ops->i2function; 567 if (ctx) { 568 if (tsdm->i2functionctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->i2functionctxcontainer, ctx)); 569 else *ctx = NULL; 570 } 571 PetscFunctionReturn(PETSC_SUCCESS); 572 } 573 574 /*@C 575 DMTSSetI2Jacobian - set `TS` implicit Jacobian evaluation function for 2nd order systems 576 577 Not Collective 578 579 Input Parameters: 580 + dm - `DM` to be used with `TS` 581 . jac - Jacobian evaluation routine 582 - ctx - context for Jacobian evaluation 583 584 Level: advanced 585 586 Note: 587 `TSSetI2Jacobian()` is normally used, but it calls this function internally because the user context is actually 588 associated with the `DM`. 589 590 .seealso: [](ch_ts), `DM`, `TS`, `TSI2Jacobian`, `TSSetI2Jacobian()` 591 @*/ 592 PetscErrorCode DMTSSetI2Jacobian(DM dm, TSI2Jacobian jac, void *ctx) 593 { 594 DMTS tsdm; 595 596 PetscFunctionBegin; 597 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 598 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 599 if (jac) tsdm->ops->i2jacobian = jac; 600 if (ctx) { 601 PetscContainer ctxcontainer; 602 PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer)); 603 PetscCall(PetscContainerSetPointer(ctxcontainer, ctx)); 604 PetscCall(PetscObjectCompose((PetscObject)tsdm, "i2jacobian ctx", (PetscObject)ctxcontainer)); 605 tsdm->i2jacobianctxcontainer = ctxcontainer; 606 PetscCall(PetscContainerDestroy(&ctxcontainer)); 607 } 608 PetscFunctionReturn(PETSC_SUCCESS); 609 } 610 611 /*@C 612 DMTSSetI2JacobianContextDestroy - set `TS` implicit Jacobian evaluation for 2nd order systems context destroy function 613 614 Not Collective 615 616 Input Parameters: 617 + dm - `DM` to be used with `TS` 618 - f - implicit Jacobian evaluation context destroy function 619 620 Level: advanced 621 622 .seealso: [](ch_ts), `DM`, `TS`, `TSSetI2JacobianContextDestroy()`, `DMTSSetI2Jacobian()`, `TSSetI2Jacobian()` 623 @*/ 624 PetscErrorCode DMTSSetI2JacobianContextDestroy(DM dm, PetscErrorCode (*f)(void *)) 625 { 626 DMTS tsdm; 627 628 PetscFunctionBegin; 629 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 630 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 631 if (tsdm->i2jacobianctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->i2jacobianctxcontainer, f)); 632 PetscFunctionReturn(PETSC_SUCCESS); 633 } 634 635 PetscErrorCode DMTSUnsetI2JacobianContext_Internal(DM dm) 636 { 637 DMTS tsdm; 638 639 PetscFunctionBegin; 640 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 641 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 642 PetscCall(DMTSUnsetI2JacobianContext_DMTS(tsdm)); 643 PetscFunctionReturn(PETSC_SUCCESS); 644 } 645 646 /*@C 647 DMTSGetI2Jacobian - get `TS` implicit Jacobian evaluation function for 2nd order systems 648 649 Not Collective 650 651 Input Parameter: 652 . dm - `DM` to be used with `TS` 653 654 Output Parameters: 655 + jac - Jacobian evaluation function, for calling sequence see `TSSetI2Jacobian()` 656 - ctx - context for Jacobian evaluation 657 658 Level: advanced 659 660 Note: 661 `TSGetI2Jacobian()` is normally used, but it calls this function internally because the user context is actually 662 associated with the `DM`. 663 664 .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetI2Jacobian()`, `TSGetI2Jacobian()` 665 @*/ 666 PetscErrorCode DMTSGetI2Jacobian(DM dm, TSI2Jacobian *jac, void **ctx) 667 { 668 DMTS tsdm; 669 670 PetscFunctionBegin; 671 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 672 PetscCall(DMGetDMTS(dm, &tsdm)); 673 if (jac) *jac = tsdm->ops->i2jacobian; 674 if (ctx) { 675 if (tsdm->i2jacobianctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->i2jacobianctxcontainer, ctx)); 676 else *ctx = NULL; 677 } 678 PetscFunctionReturn(PETSC_SUCCESS); 679 } 680 681 /*@C 682 DMTSSetRHSFunction - set `TS` explicit residual evaluation function 683 684 Not Collective 685 686 Input Parameters: 687 + dm - `DM` to be used with `TS` 688 . func - RHS function evaluation routine 689 - ctx - context for residual evaluation 690 691 Level: advanced 692 693 Note: 694 `TSSetRHSFunction()` is normally used, but it calls this function internally because the user context is actually 695 associated with the `DM`. This makes the interface consistent regardless of whether the user interacts with a `DM` or 696 not. If `DM` took a more central role at some later date, this could become the primary method of setting the residual. 697 698 .seealso: [](ch_ts), `DM`, `TS`, `TSRHSFunction` 699 @*/ 700 PetscErrorCode DMTSSetRHSFunction(DM dm, TSRHSFunction func, void *ctx) 701 { 702 DMTS tsdm; 703 704 PetscFunctionBegin; 705 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 706 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 707 if (func) tsdm->ops->rhsfunction = func; 708 if (ctx) { 709 PetscContainer ctxcontainer; 710 PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer)); 711 PetscCall(PetscContainerSetPointer(ctxcontainer, ctx)); 712 PetscCall(PetscObjectCompose((PetscObject)tsdm, "rhs function ctx", (PetscObject)ctxcontainer)); 713 tsdm->rhsfunctionctxcontainer = ctxcontainer; 714 PetscCall(PetscContainerDestroy(&ctxcontainer)); 715 } 716 PetscFunctionReturn(PETSC_SUCCESS); 717 } 718 719 /*@C 720 DMTSSetRHSFunctionContextDestroy - set `TS` explicit residual evaluation context destroy function 721 722 Not Collective 723 724 Input Parameters: 725 + dm - `DM` to be used with `TS` 726 - f - explicit evaluation context destroy function 727 728 Level: advanced 729 730 Note: 731 `TSSetRHSFunctionContextDestroy()` is normally used, but it calls this function internally because the user context is actually 732 associated with the `DM`. This makes the interface consistent regardless of whether the user interacts with a `DM` or 733 not. 734 735 Developer Notes: 736 If `DM` took a more central role at some later date, this could become the primary method of setting the residual. 737 738 .seealso: [](ch_ts), `TSSetRHSFunctionContextDestroy()`, `DMTSSetRHSFunction()`, `TSSetRHSFunction()` 739 @*/ 740 PetscErrorCode DMTSSetRHSFunctionContextDestroy(DM dm, PetscErrorCode (*f)(void *)) 741 { 742 DMTS tsdm; 743 744 PetscFunctionBegin; 745 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 746 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 747 if (tsdm->rhsfunctionctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->rhsfunctionctxcontainer, f)); 748 PetscFunctionReturn(PETSC_SUCCESS); 749 } 750 751 PetscErrorCode DMTSUnsetRHSFunctionContext_Internal(DM dm) 752 { 753 DMTS tsdm; 754 755 PetscFunctionBegin; 756 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 757 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 758 PetscCall(DMTSUnsetRHSFunctionContext_DMTS(tsdm)); 759 tsdm->rhsfunctionctxcontainer = NULL; 760 PetscFunctionReturn(PETSC_SUCCESS); 761 } 762 763 /*@C 764 DMTSSetTransientVariable - sets function to transform from state to transient variables 765 766 Logically Collective 767 768 Input Parameters: 769 + dm - `DM` to be used with `TS` 770 . tvar - a function that transforms to transient variables 771 - ctx - a context for tvar 772 773 Level: advanced 774 775 Notes: 776 This is typically used to transform from primitive to conservative variables so that a time integrator (e.g., `TSBDF`) 777 can be conservative. In this context, primitive variables P are used to model the state (e.g., because they lead to 778 well-conditioned formulations even in limiting cases such as low-Mach or zero porosity). The transient variable is 779 C(P), specified by calling this function. An IFunction thus receives arguments (P, Cdot) and the IJacobian must be 780 evaluated via the chain rule, as in 781 782 dF/dP + shift * dF/dCdot dC/dP. 783 784 .seealso: [](ch_ts), `TS`, `TSBDF`, `TSSetTransientVariable()`, `DMTSGetTransientVariable()`, `DMTSSetIFunction()`, `DMTSSetIJacobian()` 785 @*/ 786 PetscErrorCode DMTSSetTransientVariable(DM dm, TSTransientVariable tvar, void *ctx) 787 { 788 DMTS dmts; 789 790 PetscFunctionBegin; 791 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 792 PetscCall(DMGetDMTSWrite(dm, &dmts)); 793 dmts->ops->transientvar = tvar; 794 dmts->transientvarctx = ctx; 795 PetscFunctionReturn(PETSC_SUCCESS); 796 } 797 798 /*@C 799 DMTSGetTransientVariable - gets function to transform from state to transient variables set with `DMTSSetTransientVariable()` 800 801 Logically Collective 802 803 Input Parameter: 804 . dm - `DM` to be used with `TS` 805 806 Output Parameters: 807 + tvar - a function that transforms to transient variables 808 - ctx - a context for tvar 809 810 Level: advanced 811 812 .seealso: [](ch_ts), `DM`, `DMTSSetTransientVariable()`, `DMTSGetIFunction()`, `DMTSGetIJacobian()` 813 @*/ 814 PetscErrorCode DMTSGetTransientVariable(DM dm, TSTransientVariable *tvar, void *ctx) 815 { 816 DMTS dmts; 817 818 PetscFunctionBegin; 819 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 820 PetscCall(DMGetDMTS(dm, &dmts)); 821 if (tvar) *tvar = dmts->ops->transientvar; 822 if (ctx) *(void **)ctx = dmts->transientvarctx; 823 PetscFunctionReturn(PETSC_SUCCESS); 824 } 825 826 /*@C 827 DMTSGetSolutionFunction - gets the `TS` solution evaluation function 828 829 Not Collective 830 831 Input Parameter: 832 . dm - `DM` to be used with `TS` 833 834 Output Parameters: 835 + func - solution function evaluation function, for calling sequence see `TSSetSolution()` 836 - ctx - context for solution evaluation 837 838 Level: advanced 839 840 .seealso: [](ch_ts), `TS`, `DM`, `DMTSSetSolutionFunction()` 841 @*/ 842 PetscErrorCode DMTSGetSolutionFunction(DM dm, TSSolutionFunction *func, void **ctx) 843 { 844 DMTS tsdm; 845 846 PetscFunctionBegin; 847 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 848 PetscCall(DMGetDMTS(dm, &tsdm)); 849 if (func) *func = tsdm->ops->solution; 850 if (ctx) *ctx = tsdm->solutionctx; 851 PetscFunctionReturn(PETSC_SUCCESS); 852 } 853 854 /*@C 855 DMTSSetSolutionFunction - set `TS` solution evaluation function 856 857 Not Collective 858 859 Input Parameters: 860 + dm - `DM` to be used with `TS` 861 . func - solution function evaluation routine 862 - ctx - context for solution evaluation 863 864 Level: advanced 865 866 Note: 867 `TSSetSolutionFunction()` is normally used, but it calls this function internally because the user context is actually 868 associated with the `DM`. This makes the interface consistent regardless of whether the user interacts with a `DM` or 869 not. If `DM` took a more central role at some later date, this could become the primary method of setting the residual. 870 871 .seealso: [](ch_ts), `DM`, `TS`, `DMTSGetSolutionFunction()` 872 @*/ 873 PetscErrorCode DMTSSetSolutionFunction(DM dm, TSSolutionFunction func, void *ctx) 874 { 875 DMTS tsdm; 876 877 PetscFunctionBegin; 878 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 879 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 880 if (func) tsdm->ops->solution = func; 881 if (ctx) tsdm->solutionctx = ctx; 882 PetscFunctionReturn(PETSC_SUCCESS); 883 } 884 885 /*@C 886 DMTSSetForcingFunction - set `TS` forcing function evaluation function 887 888 Not Collective 889 890 Input Parameters: 891 + dm - `DM` to be used with `TS` 892 . func - forcing function evaluation routine 893 - ctx - context for solution evaluation 894 895 Level: advanced 896 897 Note: 898 `TSSetForcingFunction()` is normally used, but it calls this function internally because the user context is actually 899 associated with the `DM`. This makes the interface consistent regardless of whether the user interacts with a `DM` or 900 not. If `DM` took a more central role at some later date, this could become the primary method of setting the residual. 901 902 .seealso: [](ch_ts), `DM`, `TS`, `TSForcingFunction`, `TSSetForcingFunction()`, `DMTSGetForcingFunction()` 903 @*/ 904 PetscErrorCode DMTSSetForcingFunction(DM dm, TSForcingFunction func, void *ctx) 905 { 906 DMTS tsdm; 907 908 PetscFunctionBegin; 909 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 910 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 911 if (func) tsdm->ops->forcing = func; 912 if (ctx) tsdm->forcingctx = ctx; 913 PetscFunctionReturn(PETSC_SUCCESS); 914 } 915 916 /*@C 917 DMTSGetForcingFunction - get `TS` forcing function evaluation function 918 919 Not Collective 920 921 Input Parameter: 922 . dm - `DM` to be used with `TS` 923 924 Output Parameters: 925 + f - forcing function evaluation function; see `TSForcingFunction` for details 926 - ctx - context for solution evaluation 927 928 Level: advanced 929 930 Note: 931 `TSSetForcingFunction()` is normally used, but it calls this function internally because the user context is actually 932 associated with the `DM`. This makes the interface consistent regardless of whether the user interacts with a `DM` or 933 not. If `DM` took a more central role at some later date, this could become the primary method of setting the residual. 934 935 .seealso: [](ch_ts), `TS`, `DM`, `TSSetForcingFunction()` 936 @*/ 937 PetscErrorCode DMTSGetForcingFunction(DM dm, TSForcingFunction *f, void **ctx) 938 { 939 DMTS tsdm; 940 941 PetscFunctionBegin; 942 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 943 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 944 if (f) *f = tsdm->ops->forcing; 945 if (ctx) *ctx = tsdm->forcingctx; 946 PetscFunctionReturn(PETSC_SUCCESS); 947 } 948 949 /*@C 950 DMTSGetRHSFunction - get `TS` explicit residual evaluation function 951 952 Not Collective 953 954 Input Parameter: 955 . dm - `DM` to be used with `TS` 956 957 Output Parameters: 958 + func - residual evaluation function, for calling sequence see `TSSetRHSFunction()` 959 - ctx - context for residual evaluation 960 961 Level: advanced 962 963 Note: 964 `TSGetRHSFunction()` is normally used, but it calls this function internally because the user context is actually 965 associated with the DM. 966 967 .seealso: [](ch_ts), `DM`, `TS` 968 @*/ 969 PetscErrorCode DMTSGetRHSFunction(DM dm, TSRHSFunction *func, void **ctx) 970 { 971 DMTS tsdm; 972 973 PetscFunctionBegin; 974 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 975 PetscCall(DMGetDMTS(dm, &tsdm)); 976 if (func) *func = tsdm->ops->rhsfunction; 977 if (ctx) { 978 if (tsdm->rhsfunctionctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->rhsfunctionctxcontainer, ctx)); 979 else *ctx = NULL; 980 } 981 PetscFunctionReturn(PETSC_SUCCESS); 982 } 983 984 /*@C 985 DMTSSetIJacobian - set `TS` Jacobian evaluation function 986 987 Not Collective 988 989 Input Parameters: 990 + dm - `DM` to be used with `TS` 991 . func - Jacobian evaluation routine 992 - ctx - context for residual evaluation 993 994 Level: advanced 995 996 Note: 997 `TSSetIJacobian()` is normally used, but it calls this function internally because the user context is actually 998 associated with the `DM`. This makes the interface consistent regardless of whether the user interacts with a `DM` or 999 not. If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian. 1000 1001 .seealso: [](ch_ts), `TS`, `DM`, `TSIJacobian`, `DMTSGetIJacobian()` 1002 @*/ 1003 PetscErrorCode DMTSSetIJacobian(DM dm, TSIJacobian func, void *ctx) 1004 { 1005 DMTS tsdm; 1006 1007 PetscFunctionBegin; 1008 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1009 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 1010 if (func) tsdm->ops->ijacobian = func; 1011 if (ctx) { 1012 PetscContainer ctxcontainer; 1013 PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer)); 1014 PetscCall(PetscContainerSetPointer(ctxcontainer, ctx)); 1015 PetscCall(PetscObjectCompose((PetscObject)tsdm, "ijacobian ctx", (PetscObject)ctxcontainer)); 1016 tsdm->ijacobianctxcontainer = ctxcontainer; 1017 PetscCall(PetscContainerDestroy(&ctxcontainer)); 1018 } 1019 PetscFunctionReturn(PETSC_SUCCESS); 1020 } 1021 1022 /*@C 1023 DMTSSetIJacobianContextDestroy - set `TS` Jacobian evaluation context destroy function 1024 1025 Not Collective 1026 1027 Input Parameters: 1028 + dm - `DM` to be used with `TS` 1029 - f - Jacobian evaluation context destroy function 1030 1031 Level: advanced 1032 1033 Note: 1034 `TSSetIJacobianContextDestroy()` is normally used, but it calls this function internally because the user context is actually 1035 associated with the `DM`. This makes the interface consistent regardless of whether the user interacts with a `DM` or 1036 not. 1037 1038 Developer Notes: 1039 If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian. 1040 1041 .seealso: [](ch_ts), `TSSetIJacobianContextDestroy()`, `TSSetI2JacobianContextDestroy()`, `DMTSSetIJacobian()`, `TSSetIJacobian()` 1042 @*/ 1043 PetscErrorCode DMTSSetIJacobianContextDestroy(DM dm, PetscErrorCode (*f)(void *)) 1044 { 1045 DMTS tsdm; 1046 1047 PetscFunctionBegin; 1048 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1049 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 1050 if (tsdm->ijacobianctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->ijacobianctxcontainer, f)); 1051 PetscFunctionReturn(PETSC_SUCCESS); 1052 } 1053 1054 PetscErrorCode DMTSUnsetIJacobianContext_Internal(DM dm) 1055 { 1056 DMTS tsdm; 1057 1058 PetscFunctionBegin; 1059 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1060 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 1061 PetscCall(DMTSUnsetIJacobianContext_DMTS(tsdm)); 1062 PetscFunctionReturn(PETSC_SUCCESS); 1063 } 1064 1065 /*@C 1066 DMTSGetIJacobian - get `TS` Jacobian evaluation function 1067 1068 Not Collective 1069 1070 Input Parameter: 1071 . dm - `DM` to be used with `TS` 1072 1073 Output Parameters: 1074 + func - Jacobian evaluation function, for calling sequence see `TSSetIJacobian()` 1075 - ctx - context for residual evaluation 1076 1077 Level: advanced 1078 1079 Note: 1080 `TSGetIJacobian()` is normally used, but it calls this function internally because the user context is actually 1081 associated with the `DM`. This makes the interface consistent regardless of whether the user interacts with a `DM` or 1082 not. If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian. 1083 1084 .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetIJacobian()` 1085 @*/ 1086 PetscErrorCode DMTSGetIJacobian(DM dm, TSIJacobian *func, void **ctx) 1087 { 1088 DMTS tsdm; 1089 1090 PetscFunctionBegin; 1091 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1092 PetscCall(DMGetDMTS(dm, &tsdm)); 1093 if (func) *func = tsdm->ops->ijacobian; 1094 if (ctx) { 1095 if (tsdm->ijacobianctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->ijacobianctxcontainer, ctx)); 1096 else *ctx = NULL; 1097 } 1098 PetscFunctionReturn(PETSC_SUCCESS); 1099 } 1100 1101 /*@C 1102 DMTSSetRHSJacobian - set `TS` Jacobian evaluation function 1103 1104 Not Collective 1105 1106 Input Parameters: 1107 + dm - `DM` to be used with `TS` 1108 . func - Jacobian evaluation routine 1109 - ctx - context for residual evaluation 1110 1111 Level: advanced 1112 1113 Note: 1114 `TSSetRHSJacobian()` is normally used, but it calls this function internally because the user context is actually 1115 associated with the `DM`. This makes the interface consistent regardless of whether the user interacts with a `DM` or 1116 not. 1117 1118 Developer Notes: 1119 If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian. 1120 1121 .seealso: [](ch_ts), `TSRHSJacobian`, `DMTSGetRHSJacobian()` 1122 @*/ 1123 PetscErrorCode DMTSSetRHSJacobian(DM dm, TSRHSJacobian func, void *ctx) 1124 { 1125 DMTS tsdm; 1126 1127 PetscFunctionBegin; 1128 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1129 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 1130 if (func) tsdm->ops->rhsjacobian = func; 1131 if (ctx) { 1132 PetscContainer ctxcontainer; 1133 PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)tsdm), &ctxcontainer)); 1134 PetscCall(PetscContainerSetPointer(ctxcontainer, ctx)); 1135 PetscCall(PetscObjectCompose((PetscObject)tsdm, "rhs jacobian ctx", (PetscObject)ctxcontainer)); 1136 tsdm->rhsjacobianctxcontainer = ctxcontainer; 1137 PetscCall(PetscContainerDestroy(&ctxcontainer)); 1138 } 1139 PetscFunctionReturn(PETSC_SUCCESS); 1140 } 1141 1142 /*@C 1143 DMTSSetRHSJacobianContextDestroy - set `TS` Jacobian evaluation context destroy function 1144 1145 Not Collective 1146 1147 Input Parameters: 1148 + dm - `DM` to be used with `TS` 1149 - f - Jacobian evaluation context destroy function 1150 1151 Level: advanced 1152 1153 Note: 1154 The user usually calls `TSSetRHSJacobianContextDestroy()` which calls this routine 1155 1156 .seealso: [](ch_ts), `TS`, `TSSetRHSJacobianContextDestroy()`, `DMTSSetRHSJacobian()`, `TSSetRHSJacobian()` 1157 @*/ 1158 PetscErrorCode DMTSSetRHSJacobianContextDestroy(DM dm, PetscErrorCode (*f)(void *)) 1159 { 1160 DMTS tsdm; 1161 1162 PetscFunctionBegin; 1163 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1164 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 1165 if (tsdm->rhsjacobianctxcontainer) PetscCall(PetscContainerSetUserDestroy(tsdm->rhsjacobianctxcontainer, f)); 1166 PetscFunctionReturn(PETSC_SUCCESS); 1167 } 1168 1169 PetscErrorCode DMTSUnsetRHSJacobianContext_Internal(DM dm) 1170 { 1171 DMTS tsdm; 1172 1173 PetscFunctionBegin; 1174 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1175 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 1176 PetscCall(DMTSUnsetRHSJacobianContext_DMTS(tsdm)); 1177 PetscFunctionReturn(PETSC_SUCCESS); 1178 } 1179 1180 /*@C 1181 DMTSGetRHSJacobian - get `TS` Jacobian evaluation function 1182 1183 Not Collective 1184 1185 Input Parameter: 1186 . dm - `DM` to be used with `TS` 1187 1188 Output Parameters: 1189 + func - Jacobian evaluation function, for calling sequence see `TSSetRHSJacobian()` 1190 - ctx - context for residual evaluation 1191 1192 Level: advanced 1193 1194 Note: 1195 `TSGetRHSJacobian()` is normally used, but it calls this function internally because the user context is actually 1196 associated with the `DM`. This makes the interface consistent regardless of whether the user interacts with a `DM` or 1197 not. If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian. 1198 1199 .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetRHSJacobian()` 1200 @*/ 1201 PetscErrorCode DMTSGetRHSJacobian(DM dm, TSRHSJacobian *func, void **ctx) 1202 { 1203 DMTS tsdm; 1204 1205 PetscFunctionBegin; 1206 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1207 PetscCall(DMGetDMTS(dm, &tsdm)); 1208 if (func) *func = tsdm->ops->rhsjacobian; 1209 if (ctx) { 1210 if (tsdm->rhsjacobianctxcontainer) PetscCall(PetscContainerGetPointer(tsdm->rhsjacobianctxcontainer, ctx)); 1211 else *ctx = NULL; 1212 } 1213 PetscFunctionReturn(PETSC_SUCCESS); 1214 } 1215 1216 /*@C 1217 DMTSSetIFunctionSerialize - sets functions used to view and load a IFunction context 1218 1219 Not Collective 1220 1221 Input Parameters: 1222 + dm - `DM` to be used with `TS` 1223 . view - viewer function 1224 - load - loading function 1225 1226 Level: advanced 1227 1228 .seealso: [](ch_ts), `DM`, `TS` 1229 @*/ 1230 PetscErrorCode DMTSSetIFunctionSerialize(DM dm, PetscErrorCode (*view)(void *, PetscViewer), PetscErrorCode (*load)(void **, PetscViewer)) 1231 { 1232 DMTS tsdm; 1233 1234 PetscFunctionBegin; 1235 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1236 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 1237 tsdm->ops->ifunctionview = view; 1238 tsdm->ops->ifunctionload = load; 1239 PetscFunctionReturn(PETSC_SUCCESS); 1240 } 1241 1242 /*@C 1243 DMTSSetIJacobianSerialize - sets functions used to view and load a IJacobian context 1244 1245 Not Collective 1246 1247 Input Parameters: 1248 + dm - `DM` to be used with `TS` 1249 . view - viewer function 1250 - load - loading function 1251 1252 Level: advanced 1253 1254 .seealso: [](ch_ts), `DM`, `TS` 1255 @*/ 1256 PetscErrorCode DMTSSetIJacobianSerialize(DM dm, PetscErrorCode (*view)(void *, PetscViewer), PetscErrorCode (*load)(void **, PetscViewer)) 1257 { 1258 DMTS tsdm; 1259 1260 PetscFunctionBegin; 1261 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1262 PetscCall(DMGetDMTSWrite(dm, &tsdm)); 1263 tsdm->ops->ijacobianview = view; 1264 tsdm->ops->ijacobianload = load; 1265 PetscFunctionReturn(PETSC_SUCCESS); 1266 } 1267