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