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