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