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