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