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