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