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