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