1 #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/ 2 #include <petsc-private/dmimpl.h> /*I "petscdm.h" I*/ 3 4 #undef __FUNCT__ 5 #define __FUNCT__ "DMTSDestroy" 6 static PetscErrorCode DMTSDestroy(DMTS *kdm) 7 { 8 PetscErrorCode ierr; 9 10 PetscFunctionBegin; 11 if (!*kdm) PetscFunctionReturn(0); 12 PetscValidHeaderSpecific((*kdm),DMTS_CLASSID,1); 13 if (--((PetscObject)(*kdm))->refct > 0) {*kdm = 0; PetscFunctionReturn(0);} 14 if ((*kdm)->ops->destroy) {ierr = ((*kdm)->ops->destroy)(*kdm);CHKERRQ(ierr);} 15 ierr = PetscHeaderDestroy(kdm);CHKERRQ(ierr); 16 PetscFunctionReturn(0); 17 } 18 19 #undef __FUNCT__ 20 #define __FUNCT__ "DMTSLoad" 21 PetscErrorCode DMTSLoad(DMTS kdm,PetscViewer viewer) 22 { 23 PetscErrorCode ierr; 24 25 PetscFunctionBegin; 26 ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ifunction,1,PETSC_FUNCTION);CHKERRQ(ierr); 27 ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ifunctionview,1,PETSC_FUNCTION);CHKERRQ(ierr); 28 ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ifunctionload,1,PETSC_FUNCTION);CHKERRQ(ierr); 29 if (kdm->ops->ifunctionload) { 30 ierr = (*kdm->ops->ifunctionload)(&kdm->ifunctionctx,viewer);CHKERRQ(ierr); 31 } 32 ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ijacobian,1,PETSC_FUNCTION);CHKERRQ(ierr); 33 ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ijacobianview,1,PETSC_FUNCTION);CHKERRQ(ierr); 34 ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ijacobianload,1,PETSC_FUNCTION);CHKERRQ(ierr); 35 if (kdm->ops->ijacobianload) { 36 ierr = (*kdm->ops->ijacobianload)(&kdm->ijacobianctx,viewer);CHKERRQ(ierr); 37 } 38 PetscFunctionReturn(0); 39 } 40 41 #undef __FUNCT__ 42 #define __FUNCT__ "DMTSView" 43 PetscErrorCode DMTSView(DMTS kdm,PetscViewer viewer) 44 { 45 PetscErrorCode ierr; 46 PetscBool isascii,isbinary; 47 48 PetscFunctionBegin; 49 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 50 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 51 if (isascii) { 52 const char *fname; 53 54 ierr = PetscFPTFind(kdm->ops->ifunction,&fname);CHKERRQ(ierr); 55 if (fname) { 56 ierr = PetscViewerASCIIPrintf(viewer," IFunction used by TS: %s\n",fname);CHKERRQ(ierr); 57 } 58 ierr = PetscFPTFind(kdm->ops->ijacobian,&fname);CHKERRQ(ierr); 59 if (fname) { 60 ierr = PetscViewerASCIIPrintf(viewer," IJacobian function used by TS: %s\n",fname);CHKERRQ(ierr); 61 } 62 } else if (isbinary) { 63 ierr = PetscViewerBinaryWrite(viewer,kdm->ops->ifunction,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr); 64 ierr = PetscViewerBinaryWrite(viewer,kdm->ops->ifunctionview,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr); 65 ierr = PetscViewerBinaryWrite(viewer,kdm->ops->ifunctionload,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr); 66 if (kdm->ops->ifunctionview) { 67 ierr = (*kdm->ops->ifunctionview)(kdm->ifunctionctx,viewer);CHKERRQ(ierr); 68 } 69 ierr = PetscViewerBinaryWrite(viewer,kdm->ops->ijacobian,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr); 70 ierr = PetscViewerBinaryWrite(viewer,kdm->ops->ijacobianview,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr); 71 ierr = PetscViewerBinaryWrite(viewer,kdm->ops->ijacobianload,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr); 72 if (kdm->ops->ijacobianview) { 73 ierr = (*kdm->ops->ijacobianview)(kdm->ijacobianctx,viewer);CHKERRQ(ierr); 74 } 75 } 76 PetscFunctionReturn(0); 77 } 78 79 #undef __FUNCT__ 80 #define __FUNCT__ "DMTSCreate" 81 static PetscErrorCode DMTSCreate(MPI_Comm comm,DMTS *kdm) 82 { 83 PetscErrorCode ierr; 84 85 PetscFunctionBegin; 86 #ifndef PETSC_USE_DYNAMIC_LIBRARIES 87 ierr = TSInitializePackage(PETSC_NULL);CHKERRQ(ierr); 88 #endif 89 ierr = PetscHeaderCreate(*kdm, _p_DMTS, struct _DMTSOps, DMTS_CLASSID, -1, "DMTS", "DMTS", "DMTS", comm, DMTSDestroy, DMTSView);CHKERRQ(ierr); 90 ierr = PetscMemzero((*kdm)->ops, sizeof(struct _DMTSOps));CHKERRQ(ierr); 91 PetscFunctionReturn(0); 92 } 93 94 #undef __FUNCT__ 95 #define __FUNCT__ "DMCoarsenHook_DMTS" 96 /* Attaches the DMTS to the coarse level. 97 * Under what conditions should we copy versus duplicate? 98 */ 99 static PetscErrorCode DMCoarsenHook_DMTS(DM dm,DM dmc,void *ctx) 100 { 101 PetscErrorCode ierr; 102 103 PetscFunctionBegin; 104 ierr = DMCopyDMTS(dm,dmc);CHKERRQ(ierr); 105 PetscFunctionReturn(0); 106 } 107 108 #undef __FUNCT__ 109 #define __FUNCT__ "DMRestrictHook_DMTS" 110 /* This could restrict auxiliary information to the coarse level. 111 */ 112 static PetscErrorCode DMRestrictHook_DMTS(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx) 113 { 114 115 PetscFunctionBegin; 116 PetscFunctionReturn(0); 117 } 118 119 #undef __FUNCT__ 120 #define __FUNCT__ "DMTSCopy" 121 /*@C 122 DMTSCopy - copies the information in a DMTS to another DMTS 123 124 Not Collective 125 126 Input Argument: 127 + kdm - Original DMTS 128 - nkdm - DMTS to receive the data, should have been created with DMTSCreate() 129 130 Level: developer 131 132 .seealso: DMTSCreate(), DMTSDestroy() 133 @*/ 134 PetscErrorCode DMTSCopy(DMTS kdm,DMTS nkdm) 135 { 136 PetscErrorCode ierr; 137 138 PetscFunctionBegin; 139 PetscValidHeaderSpecific(kdm,DMTS_CLASSID,1); 140 PetscValidHeaderSpecific(nkdm,DMTS_CLASSID,2); 141 nkdm->ops->rhsfunction = kdm->ops->rhsfunction; 142 nkdm->ops->rhsjacobian = kdm->ops->rhsjacobian; 143 nkdm->ops->ifunction = kdm->ops->ifunction; 144 nkdm->ops->ijacobian = kdm->ops->ijacobian; 145 nkdm->ops->solution = kdm->ops->solution; 146 nkdm->ops->destroy = kdm->ops->destroy; 147 nkdm->ops->duplicate = kdm->ops->duplicate; 148 149 nkdm->rhsfunctionctx = kdm->rhsfunctionctx; 150 nkdm->rhsjacobianctx = kdm->rhsjacobianctx; 151 nkdm->ifunctionctx = kdm->ifunctionctx; 152 nkdm->ijacobianctx = kdm->ijacobianctx; 153 nkdm->solutionctx = kdm->solutionctx; 154 155 nkdm->data = kdm->data; 156 157 /* 158 nkdm->fortran_func_pointers[0] = kdm->fortran_func_pointers[0]; 159 nkdm->fortran_func_pointers[1] = kdm->fortran_func_pointers[1]; 160 nkdm->fortran_func_pointers[2] = kdm->fortran_func_pointers[2]; 161 */ 162 163 /* implementation specific copy hooks */ 164 if (kdm->ops->duplicate) {ierr = (*kdm->ops->duplicate)(kdm,nkdm);CHKERRQ(ierr);} 165 PetscFunctionReturn(0); 166 } 167 168 #undef __FUNCT__ 169 #define __FUNCT__ "DMGetDMTS" 170 /*@C 171 DMGetDMTS - get read-only private DMTS context from a DM 172 173 Not Collective 174 175 Input Argument: 176 . dm - DM to be used with TS 177 178 Output Argument: 179 . tsdm - private DMTS context 180 181 Level: developer 182 183 Notes: 184 Use DMGetDMTSWrite() if write access is needed. The DMTSSetXXX API should be used wherever possible. 185 186 .seealso: DMGetDMTSWrite() 187 @*/ 188 PetscErrorCode DMGetDMTS(DM dm,DMTS *tsdm) 189 { 190 PetscErrorCode ierr; 191 192 PetscFunctionBegin; 193 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 194 *tsdm = (DMTS) dm->dmts; 195 if (!*tsdm) { 196 ierr = PetscInfo(dm,"Creating new DMTS\n");CHKERRQ(ierr); 197 ierr = DMTSCreate(((PetscObject)dm)->comm,tsdm);CHKERRQ(ierr); 198 dm->dmts = (PetscObject) *tsdm; 199 ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,PETSC_NULL);CHKERRQ(ierr); 200 } 201 PetscFunctionReturn(0); 202 } 203 204 #undef __FUNCT__ 205 #define __FUNCT__ "DMGetDMTSWrite" 206 /*@C 207 DMGetDMTSWrite - get write access to private DMTS context from a DM 208 209 Not Collective 210 211 Input Argument: 212 . dm - DM to be used with TS 213 214 Output Argument: 215 . tsdm - private DMTS context 216 217 Level: developer 218 219 .seealso: DMGetDMTS() 220 @*/ 221 PetscErrorCode DMGetDMTSWrite(DM dm,DMTS *tsdm) 222 { 223 PetscErrorCode ierr; 224 DMTS sdm; 225 226 PetscFunctionBegin; 227 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 228 ierr = DMGetDMTS(dm,&sdm);CHKERRQ(ierr); 229 if (!sdm->originaldm) sdm->originaldm = dm; 230 if (sdm->originaldm != dm) { /* Copy on write */ 231 DMTS oldsdm = sdm; 232 ierr = PetscInfo(dm,"Copying DMTS due to write\n");CHKERRQ(ierr); 233 ierr = DMTSCreate(((PetscObject)dm)->comm,&sdm);CHKERRQ(ierr); 234 ierr = DMTSCopy(oldsdm,sdm);CHKERRQ(ierr); 235 ierr = DMTSDestroy((DMTS*)&dm->dmts);CHKERRQ(ierr); 236 dm->dmts = (PetscObject) sdm; 237 } 238 *tsdm = sdm; 239 PetscFunctionReturn(0); 240 } 241 242 #undef __FUNCT__ 243 #define __FUNCT__ "DMCopyDMTS" 244 /*@C 245 DMCopyDMTS - copies a DM context to a new DM 246 247 Logically Collective 248 249 Input Arguments: 250 + dmsrc - DM to obtain context from 251 - dmdest - DM to add context to 252 253 Level: developer 254 255 Note: 256 The context is copied by reference. This function does not ensure that a context exists. 257 258 .seealso: DMGetDMTS(), TSSetDM() 259 @*/ 260 PetscErrorCode DMCopyDMTS(DM dmsrc,DM dmdest) 261 { 262 PetscErrorCode ierr; 263 264 PetscFunctionBegin; 265 PetscValidHeaderSpecific(dmsrc,DM_CLASSID,1); 266 PetscValidHeaderSpecific(dmdest,DM_CLASSID,2); 267 ierr = DMTSDestroy((DMTS*)&dmdest->dmts);CHKERRQ(ierr); 268 dmdest->dmts = dmsrc->dmts; 269 ierr = PetscObjectReference(dmdest->dmts);CHKERRQ(ierr); 270 ierr = DMCoarsenHookAdd(dmdest,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,PETSC_NULL);CHKERRQ(ierr); 271 PetscFunctionReturn(0); 272 } 273 274 #undef __FUNCT__ 275 #define __FUNCT__ "DMTSSetIFunction" 276 /*@C 277 DMTSSetIFunction - set TS implicit function evaluation function 278 279 Not Collective 280 281 Input Arguments: 282 + dm - DM to be used with TS 283 . func - function evaluation function, see TSSetIFunction() for calling sequence 284 - ctx - context for residual evaluation 285 286 Level: advanced 287 288 Note: 289 TSSetFunction() is normally used, but it calls this function internally because the user context is actually 290 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 291 not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 292 293 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian() 294 @*/ 295 PetscErrorCode DMTSSetIFunction(DM dm,TSIFunction func,void *ctx) 296 { 297 PetscErrorCode ierr; 298 DMTS tsdm; 299 300 PetscFunctionBegin; 301 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 302 ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr); 303 if (func) tsdm->ops->ifunction = func; 304 if (ctx) tsdm->ifunctionctx = ctx; 305 PetscFunctionReturn(0); 306 } 307 308 #undef __FUNCT__ 309 #define __FUNCT__ "DMTSGetIFunction" 310 /*@C 311 DMTSGetIFunction - get TS implicit residual evaluation function 312 313 Not Collective 314 315 Input Argument: 316 . dm - DM to be used with TS 317 318 Output Arguments: 319 + func - function evaluation function, see TSSetIFunction() for calling sequence 320 - ctx - context for residual evaluation 321 322 Level: advanced 323 324 Note: 325 TSGetFunction() is normally used, but it calls this function internally because the user context is actually 326 associated with the DM. 327 328 .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction() 329 @*/ 330 PetscErrorCode DMTSGetIFunction(DM dm,TSIFunction *func,void **ctx) 331 { 332 PetscErrorCode ierr; 333 DMTS tsdm; 334 335 PetscFunctionBegin; 336 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 337 ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr); 338 if (func) *func = tsdm->ops->ifunction; 339 if (ctx) *ctx = tsdm->ifunctionctx; 340 PetscFunctionReturn(0); 341 } 342 343 344 #undef __FUNCT__ 345 #define __FUNCT__ "DMTSSetRHSFunction" 346 /*@C 347 DMTSSetRHSFunction - set TS explicit residual evaluation function 348 349 Not Collective 350 351 Input Arguments: 352 + dm - DM to be used with TS 353 . func - RHS function evaluation function, see TSSetRHSFunction() for calling sequence 354 - ctx - context for residual evaluation 355 356 Level: advanced 357 358 Note: 359 TSSetRSHFunction() is normally used, but it calls this function internally because the user context is actually 360 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 361 not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 362 363 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian() 364 @*/ 365 PetscErrorCode DMTSSetRHSFunction(DM dm,TSRHSFunction func,void *ctx) 366 { 367 PetscErrorCode ierr; 368 DMTS tsdm; 369 370 PetscFunctionBegin; 371 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 372 ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr); 373 if (func) tsdm->ops->rhsfunction = func; 374 if (ctx) tsdm->rhsfunctionctx = ctx; 375 PetscFunctionReturn(0); 376 } 377 378 #undef __FUNCT__ 379 #define __FUNCT__ "DMTSGetSolutionFunction" 380 /*@C 381 DMTSGetSolutionFunction - gets the TS solution evaluation function 382 383 Not Collective 384 385 Input Arguments: 386 . dm - DM to be used with TS 387 388 Output Parameters: 389 + func - solution function evaluation function, see TSSetSolution() for calling sequence 390 - ctx - context for solution evaluation 391 392 Level: advanced 393 394 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian() 395 @*/ 396 PetscErrorCode DMTSGetSolutionFunction(DM dm,TSSolutionFunction *func,void **ctx) 397 { 398 PetscErrorCode ierr; 399 DMTS tsdm; 400 401 PetscFunctionBegin; 402 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 403 ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr); 404 if (func) *func = tsdm->ops->solution; 405 if (ctx) *ctx = tsdm->solutionctx; 406 PetscFunctionReturn(0); 407 } 408 409 #undef __FUNCT__ 410 #define __FUNCT__ "DMTSSetSolutionFunction" 411 /*@C 412 DMTSSetSolutionFunction - set TS solution evaluation function 413 414 Not Collective 415 416 Input Arguments: 417 + dm - DM to be used with TS 418 . func - solution function evaluation function, see TSSetSolution() for calling sequence 419 - ctx - context for solution evaluation 420 421 Level: advanced 422 423 Note: 424 TSSetSolutionFunction() is normally used, but it calls this function internally because the user context is actually 425 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 426 not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 427 428 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian() 429 @*/ 430 PetscErrorCode DMTSSetSolutionFunction(DM dm,TSSolutionFunction func,void *ctx) 431 { 432 PetscErrorCode ierr; 433 DMTS tsdm; 434 435 PetscFunctionBegin; 436 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 437 ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr); 438 if (func) tsdm->ops->solution = func; 439 if (ctx) tsdm->solutionctx = ctx; 440 PetscFunctionReturn(0); 441 } 442 443 #undef __FUNCT__ 444 #define __FUNCT__ "DMTSGetRHSFunction" 445 /*@C 446 DMTSGetRHSFunction - get TS explicit residual evaluation function 447 448 Not Collective 449 450 Input Argument: 451 . dm - DM to be used with TS 452 453 Output Arguments: 454 + func - residual evaluation function, see TSSetRHSFunction() for calling sequence 455 - ctx - context for residual evaluation 456 457 Level: advanced 458 459 Note: 460 TSGetFunction() is normally used, but it calls this function internally because the user context is actually 461 associated with the DM. 462 463 .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction() 464 @*/ 465 PetscErrorCode DMTSGetRHSFunction(DM dm,TSRHSFunction *func,void **ctx) 466 { 467 PetscErrorCode ierr; 468 DMTS tsdm; 469 470 PetscFunctionBegin; 471 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 472 ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr); 473 if (func) *func = tsdm->ops->rhsfunction; 474 if (ctx) *ctx = tsdm->rhsfunctionctx; 475 PetscFunctionReturn(0); 476 } 477 478 #undef __FUNCT__ 479 #define __FUNCT__ "DMTSSetIJacobian" 480 /*@C 481 DMTSSetIJacobian - set TS Jacobian evaluation function 482 483 Not Collective 484 485 Input Argument: 486 + dm - DM to be used with TS 487 . func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence 488 - ctx - context for residual evaluation 489 490 Level: advanced 491 492 Note: 493 TSSetJacobian() is normally used, but it calls this function internally because the user context is actually 494 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 495 not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 496 497 .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian() 498 @*/ 499 PetscErrorCode DMTSSetIJacobian(DM dm,TSIJacobian func,void *ctx) 500 { 501 PetscErrorCode ierr; 502 DMTS sdm; 503 504 PetscFunctionBegin; 505 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 506 ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr); 507 if (func) sdm->ops->ijacobian = func; 508 if (ctx) sdm->ijacobianctx = ctx; 509 PetscFunctionReturn(0); 510 } 511 512 #undef __FUNCT__ 513 #define __FUNCT__ "DMTSGetIJacobian" 514 /*@C 515 DMTSGetIJacobian - get TS Jacobian evaluation function 516 517 Not Collective 518 519 Input Argument: 520 . dm - DM to be used with TS 521 522 Output Arguments: 523 + func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence 524 - ctx - context for residual evaluation 525 526 Level: advanced 527 528 Note: 529 TSGetJacobian() is normally used, but it calls this function internally because the user context is actually 530 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 531 not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 532 533 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian() 534 @*/ 535 PetscErrorCode DMTSGetIJacobian(DM dm,TSIJacobian *func,void **ctx) 536 { 537 PetscErrorCode ierr; 538 DMTS tsdm; 539 540 PetscFunctionBegin; 541 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 542 ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr); 543 if (func) *func = tsdm->ops->ijacobian; 544 if (ctx) *ctx = tsdm->ijacobianctx; 545 PetscFunctionReturn(0); 546 } 547 548 549 #undef __FUNCT__ 550 #define __FUNCT__ "DMTSSetRHSJacobian" 551 /*@C 552 DMTSSetRHSJacobian - set TS Jacobian evaluation function 553 554 Not Collective 555 556 Input Argument: 557 + dm - DM to be used with TS 558 . func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence 559 - ctx - context for residual evaluation 560 561 Level: advanced 562 563 Note: 564 TSSetJacobian() is normally used, but it calls this function internally because the user context is actually 565 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 566 not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 567 568 .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian() 569 @*/ 570 PetscErrorCode DMTSSetRHSJacobian(DM dm,TSRHSJacobian func,void *ctx) 571 { 572 PetscErrorCode ierr; 573 DMTS tsdm; 574 575 PetscFunctionBegin; 576 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 577 ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr); 578 if (func) tsdm->ops->rhsjacobian = func; 579 if (ctx) tsdm->rhsjacobianctx = ctx; 580 PetscFunctionReturn(0); 581 } 582 583 #undef __FUNCT__ 584 #define __FUNCT__ "DMTSGetRHSJacobian" 585 /*@C 586 DMTSGetRHSJacobian - get TS Jacobian evaluation function 587 588 Not Collective 589 590 Input Argument: 591 . dm - DM to be used with TS 592 593 Output Arguments: 594 + func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence 595 - ctx - context for residual evaluation 596 597 Level: advanced 598 599 Note: 600 TSGetJacobian() is normally used, but it calls this function internally because the user context is actually 601 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 602 not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 603 604 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian() 605 @*/ 606 PetscErrorCode DMTSGetRHSJacobian(DM dm,TSRHSJacobian *func,void **ctx) 607 { 608 PetscErrorCode ierr; 609 DMTS tsdm; 610 611 PetscFunctionBegin; 612 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 613 ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr); 614 if (func) *func = tsdm->ops->rhsjacobian; 615 if (ctx) *ctx = tsdm->rhsjacobianctx; 616 PetscFunctionReturn(0); 617 } 618 619 #undef __FUNCT__ 620 #define __FUNCT__ "DMTSSetIFunctionSerialize" 621 /*@C 622 DMTSSetIFunctionSerialize - sets functions used to view and load a IFunction context 623 624 Not Collective 625 626 Input Arguments: 627 + dm - DM to be used with TS 628 . view - viewer function 629 - load - loading function 630 631 Level: advanced 632 633 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian() 634 @*/ 635 PetscErrorCode DMTSSetIFunctionSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer)) 636 { 637 PetscErrorCode ierr; 638 DMTS tsdm; 639 640 PetscFunctionBegin; 641 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 642 ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr); 643 tsdm->ops->ifunctionview = view; 644 tsdm->ops->ifunctionload = load; 645 PetscFunctionReturn(0); 646 } 647 648 #undef __FUNCT__ 649 #define __FUNCT__ "DMTSSetIJacobianSerialize" 650 /*@C 651 DMTSSetIJacobianSerialize - sets functions used to view and load a IJacobian context 652 653 Not Collective 654 655 Input Arguments: 656 + dm - DM to be used with TS 657 . view - viewer function 658 - load - loading function 659 660 Level: advanced 661 662 .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian() 663 @*/ 664 PetscErrorCode DMTSSetIJacobianSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer)) 665 { 666 PetscErrorCode ierr; 667 DMTS tsdm; 668 669 PetscFunctionBegin; 670 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 671 ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr); 672 tsdm->ops->ijacobianview = view; 673 tsdm->ops->ijacobianload = load; 674 PetscFunctionReturn(0); 675 } 676