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