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