1 #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 2 #include <petsc/private/dmimpl.h> /*I "petscdm.h" I*/ 3 4 static PetscErrorCode DMSNESUnsetFunctionContext_DMSNES(DMSNES sdm) 5 { 6 PetscFunctionBegin; 7 PetscCall(PetscObjectCompose((PetscObject)sdm, "function ctx", NULL)); 8 sdm->functionctxcontainer = NULL; 9 PetscFunctionReturn(PETSC_SUCCESS); 10 } 11 12 static PetscErrorCode DMSNESUnsetJacobianContext_DMSNES(DMSNES sdm) 13 { 14 PetscFunctionBegin; 15 PetscCall(PetscObjectCompose((PetscObject)sdm, "jacobian ctx", NULL)); 16 sdm->jacobianctxcontainer = NULL; 17 PetscFunctionReturn(PETSC_SUCCESS); 18 } 19 20 static PetscErrorCode DMSNESDestroy(DMSNES *kdm) 21 { 22 PetscFunctionBegin; 23 if (!*kdm) PetscFunctionReturn(PETSC_SUCCESS); 24 PetscValidHeaderSpecific((*kdm), DMSNES_CLASSID, 1); 25 if (--((PetscObject)(*kdm))->refct > 0) { 26 *kdm = NULL; 27 PetscFunctionReturn(PETSC_SUCCESS); 28 } 29 PetscCall(DMSNESUnsetFunctionContext_DMSNES(*kdm)); 30 PetscCall(DMSNESUnsetJacobianContext_DMSNES(*kdm)); 31 if ((*kdm)->ops->destroy) PetscCall(((*kdm)->ops->destroy)(*kdm)); 32 PetscCall(PetscHeaderDestroy(kdm)); 33 PetscFunctionReturn(PETSC_SUCCESS); 34 } 35 36 PetscErrorCode DMSNESLoad(DMSNES kdm, PetscViewer viewer) 37 { 38 PetscFunctionBegin; 39 PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->computefunction, 1, NULL, PETSC_FUNCTION)); 40 PetscCall(PetscViewerBinaryRead(viewer, &kdm->ops->computejacobian, 1, NULL, PETSC_FUNCTION)); 41 PetscFunctionReturn(PETSC_SUCCESS); 42 } 43 44 PetscErrorCode DMSNESView(DMSNES kdm, PetscViewer viewer) 45 { 46 PetscBool isascii, isbinary; 47 48 PetscFunctionBegin; 49 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 50 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 51 if (isascii) { 52 #if defined(PETSC_SERIALIZE_FUNCTIONS) 53 const char *fname; 54 55 PetscCall(PetscFPTFind(kdm->ops->computefunction, &fname)); 56 if (fname) PetscCall(PetscViewerASCIIPrintf(viewer, "Function used by SNES: %s\n", fname)); 57 PetscCall(PetscFPTFind(kdm->ops->computejacobian, &fname)); 58 if (fname) PetscCall(PetscViewerASCIIPrintf(viewer, "Jacobian function used by SNES: %s\n", fname)); 59 #endif 60 } else if (isbinary) { 61 struct { 62 PetscErrorCode (*func)(SNES, Vec, Vec, void *); 63 } funcstruct; 64 struct { 65 PetscErrorCode (*jac)(SNES, Vec, Mat, Mat, void *); 66 } jacstruct; 67 funcstruct.func = kdm->ops->computefunction; 68 jacstruct.jac = kdm->ops->computejacobian; 69 PetscCall(PetscViewerBinaryWrite(viewer, &funcstruct, 1, PETSC_FUNCTION)); 70 PetscCall(PetscViewerBinaryWrite(viewer, &jacstruct, 1, PETSC_FUNCTION)); 71 } 72 PetscFunctionReturn(PETSC_SUCCESS); 73 } 74 75 static PetscErrorCode DMSNESCreate(MPI_Comm comm, DMSNES *kdm) 76 { 77 PetscFunctionBegin; 78 PetscCall(SNESInitializePackage()); 79 PetscCall(PetscHeaderCreate(*kdm, DMSNES_CLASSID, "DMSNES", "DMSNES", "DMSNES", comm, DMSNESDestroy, DMSNESView)); 80 PetscFunctionReturn(PETSC_SUCCESS); 81 } 82 83 /* Attaches the DMSNES to the coarse level. 84 * Under what conditions should we copy versus duplicate? 85 */ 86 static PetscErrorCode DMCoarsenHook_DMSNES(DM dm, DM dmc, void *ctx) 87 { 88 PetscFunctionBegin; 89 PetscCall(DMCopyDMSNES(dm, dmc)); 90 PetscFunctionReturn(PETSC_SUCCESS); 91 } 92 93 /* This could restrict auxiliary information to the coarse level. 94 */ 95 static PetscErrorCode DMRestrictHook_DMSNES(DM dm, Mat Restrict, Vec rscale, Mat Inject, DM dmc, void *ctx) 96 { 97 PetscFunctionBegin; 98 PetscFunctionReturn(PETSC_SUCCESS); 99 } 100 101 /* Attaches the DMSNES to the subdomain. */ 102 static PetscErrorCode DMSubDomainHook_DMSNES(DM dm, DM subdm, void *ctx) 103 { 104 PetscFunctionBegin; 105 PetscCall(DMCopyDMSNES(dm, subdm)); 106 PetscFunctionReturn(PETSC_SUCCESS); 107 } 108 109 /* This could restrict auxiliary information to the coarse level. 110 */ 111 static PetscErrorCode DMSubDomainRestrictHook_DMSNES(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx) 112 { 113 PetscFunctionBegin; 114 PetscFunctionReturn(PETSC_SUCCESS); 115 } 116 117 static PetscErrorCode DMRefineHook_DMSNES(DM dm, DM dmf, void *ctx) 118 { 119 PetscFunctionBegin; 120 PetscCall(DMCopyDMSNES(dm, dmf)); 121 PetscFunctionReturn(PETSC_SUCCESS); 122 } 123 124 /* This could restrict auxiliary information to the coarse level. 125 */ 126 static PetscErrorCode DMInterpolateHook_DMSNES(DM dm, Mat Interp, DM dmf, void *ctx) 127 { 128 PetscFunctionBegin; 129 PetscFunctionReturn(PETSC_SUCCESS); 130 } 131 132 /*@C 133 DMSNESCopy - copies the information in a `DMSNES` to another `DMSNES` 134 135 Not Collective 136 137 Input Parameters: 138 + kdm - Original `DMSNES` 139 - nkdm - `DMSNES` to receive the data, should have been created with `DMSNESCreate()` 140 141 Level: developer 142 143 .seealso: `DMSNES`, `DMSNESCreate()`, `DMSNESDestroy()` 144 @*/ 145 PetscErrorCode DMSNESCopy(DMSNES kdm, DMSNES nkdm) 146 { 147 PetscFunctionBegin; 148 PetscValidHeaderSpecific(kdm, DMSNES_CLASSID, 1); 149 PetscValidHeaderSpecific(nkdm, DMSNES_CLASSID, 2); 150 nkdm->ops->computefunction = kdm->ops->computefunction; 151 nkdm->ops->computejacobian = kdm->ops->computejacobian; 152 nkdm->ops->computegs = kdm->ops->computegs; 153 nkdm->ops->computeobjective = kdm->ops->computeobjective; 154 nkdm->ops->computepjacobian = kdm->ops->computepjacobian; 155 nkdm->ops->computepfunction = kdm->ops->computepfunction; 156 nkdm->ops->destroy = kdm->ops->destroy; 157 nkdm->ops->duplicate = kdm->ops->duplicate; 158 159 nkdm->gsctx = kdm->gsctx; 160 nkdm->pctx = kdm->pctx; 161 nkdm->objectivectx = kdm->objectivectx; 162 nkdm->originaldm = kdm->originaldm; 163 nkdm->functionctxcontainer = kdm->functionctxcontainer; 164 nkdm->jacobianctxcontainer = kdm->jacobianctxcontainer; 165 if (nkdm->functionctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "function ctx", (PetscObject)nkdm->functionctxcontainer)); 166 if (nkdm->jacobianctxcontainer) PetscCall(PetscObjectCompose((PetscObject)nkdm, "jacobian ctx", (PetscObject)nkdm->jacobianctxcontainer)); 167 168 /* 169 nkdm->fortran_func_pointers[0] = kdm->fortran_func_pointers[0]; 170 nkdm->fortran_func_pointers[1] = kdm->fortran_func_pointers[1]; 171 nkdm->fortran_func_pointers[2] = kdm->fortran_func_pointers[2]; 172 */ 173 174 /* implementation specific copy hooks */ 175 PetscTryTypeMethod(kdm, duplicate, nkdm); 176 PetscFunctionReturn(PETSC_SUCCESS); 177 } 178 179 /*@C 180 DMGetDMSNES - get read-only private `DMSNES` context from a `DM` 181 182 Not Collective 183 184 Input Parameter: 185 . dm - `DM` to be used with `SNES` 186 187 Output Parameter: 188 . snesdm - private `DMSNES` context 189 190 Level: developer 191 192 Note: 193 Use `DMGetDMSNESWrite()` if write access is needed. The DMSNESSetXXX API should be used wherever possible. 194 195 .seealso: `DMGetDMSNESWrite()` 196 @*/ 197 PetscErrorCode DMGetDMSNES(DM dm, DMSNES *snesdm) 198 { 199 PetscFunctionBegin; 200 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 201 *snesdm = (DMSNES)dm->dmsnes; 202 if (!*snesdm) { 203 PetscCall(PetscInfo(dm, "Creating new DMSNES\n")); 204 PetscCall(DMSNESCreate(PetscObjectComm((PetscObject)dm), snesdm)); 205 206 dm->dmsnes = (PetscObject)*snesdm; 207 (*snesdm)->originaldm = dm; 208 PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_DMSNES, DMRestrictHook_DMSNES, NULL)); 209 PetscCall(DMRefineHookAdd(dm, DMRefineHook_DMSNES, DMInterpolateHook_DMSNES, NULL)); 210 PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_DMSNES, DMSubDomainRestrictHook_DMSNES, NULL)); 211 } 212 PetscFunctionReturn(PETSC_SUCCESS); 213 } 214 215 /*@C 216 DMGetDMSNESWrite - get write access to private `DMSNES` context from a `DM` 217 218 Not Collective 219 220 Input Parameter: 221 . dm - `DM` to be used with `SNES` 222 223 Output Parameter: 224 . snesdm - private `DMSNES` context 225 226 Level: developer 227 228 .seealso: `DMGetDMSNES()` 229 @*/ 230 PetscErrorCode DMGetDMSNESWrite(DM dm, DMSNES *snesdm) 231 { 232 DMSNES sdm; 233 234 PetscFunctionBegin; 235 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 236 PetscCall(DMGetDMSNES(dm, &sdm)); 237 PetscCheck(sdm->originaldm, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "DMSNES has a NULL originaldm"); 238 if (sdm->originaldm != dm) { /* Copy on write */ 239 DMSNES oldsdm = sdm; 240 PetscCall(PetscInfo(dm, "Copying DMSNES due to write\n")); 241 PetscCall(DMSNESCreate(PetscObjectComm((PetscObject)dm), &sdm)); 242 PetscCall(DMSNESCopy(oldsdm, sdm)); 243 PetscCall(DMSNESDestroy((DMSNES *)&dm->dmsnes)); 244 dm->dmsnes = (PetscObject)sdm; 245 sdm->originaldm = dm; 246 } 247 *snesdm = sdm; 248 PetscFunctionReturn(PETSC_SUCCESS); 249 } 250 251 /*@C 252 DMCopyDMSNES - copies a `DMSNES` context to a new `DM` 253 254 Logically Collective 255 256 Input Parameters: 257 + dmsrc - `DM` to obtain context from 258 - dmdest - `DM` to add context to 259 260 Level: developer 261 262 Note: 263 The context is copied by reference. This function does not ensure that a context exists. 264 265 .seealso: `DMGetDMSNES()`, `SNESSetDM()` 266 @*/ 267 PetscErrorCode DMCopyDMSNES(DM dmsrc, DM dmdest) 268 { 269 PetscFunctionBegin; 270 PetscValidHeaderSpecific(dmsrc, DM_CLASSID, 1); 271 PetscValidHeaderSpecific(dmdest, DM_CLASSID, 2); 272 if (!dmdest->dmsnes) PetscCall(DMSNESCreate(PetscObjectComm((PetscObject)dmdest), (DMSNES *)&dmdest->dmsnes)); 273 PetscCall(DMSNESCopy((DMSNES)dmsrc->dmsnes, (DMSNES)dmdest->dmsnes)); 274 PetscCall(DMCoarsenHookAdd(dmdest, DMCoarsenHook_DMSNES, NULL, NULL)); 275 PetscCall(DMRefineHookAdd(dmdest, DMRefineHook_DMSNES, NULL, NULL)); 276 PetscCall(DMSubDomainHookAdd(dmdest, DMSubDomainHook_DMSNES, DMSubDomainRestrictHook_DMSNES, NULL)); 277 PetscFunctionReturn(PETSC_SUCCESS); 278 } 279 280 /*@C 281 DMSNESSetFunction - set `SNES` residual evaluation function 282 283 Not Collective 284 285 Input Parameters: 286 + dm - DM to be used with `SNES` 287 . f - residual evaluation function; see `SNESFunction` for details 288 - ctx - context for residual evaluation 289 290 Level: advanced 291 292 Note: 293 `SNESSetFunction()` is normally used, but it calls this function internally because the user context is actually 294 associated with the `DM`. This makes the interface consistent regardless of whether the user interacts with a `DM` or 295 not. 296 297 Developer Notes: 298 If `DM` took a more central role at some later date, this could become the primary method of setting the residual. 299 300 .seealso: `DMSNESSetContext()`, `SNESSetFunction()`, `DMSNESSetJacobian()`, `SNESFunction` 301 @*/ 302 PetscErrorCode DMSNESSetFunction(DM dm, PetscErrorCode (*f)(SNES, Vec, Vec, void *), void *ctx) 303 { 304 DMSNES sdm; 305 306 PetscFunctionBegin; 307 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 308 PetscCall(DMGetDMSNESWrite(dm, &sdm)); 309 if (f) sdm->ops->computefunction = f; 310 if (ctx) { 311 PetscContainer ctxcontainer; 312 PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)sdm), &ctxcontainer)); 313 PetscCall(PetscContainerSetPointer(ctxcontainer, ctx)); 314 PetscCall(PetscObjectCompose((PetscObject)sdm, "function ctx", (PetscObject)ctxcontainer)); 315 sdm->functionctxcontainer = ctxcontainer; 316 PetscCall(PetscContainerDestroy(&ctxcontainer)); 317 } 318 PetscFunctionReturn(PETSC_SUCCESS); 319 } 320 321 /*@C 322 DMSNESSetFunctionContextDestroy - set `SNES` residual evaluation context destroy function 323 324 Not Collective 325 326 Input Parameters: 327 + dm - `DM` to be used with `SNES` 328 - f - residual evaluation context destroy function 329 330 Level: advanced 331 332 .seealso: `DMSNESSetFunction()`, `SNESSetFunction()` 333 @*/ 334 PetscErrorCode DMSNESSetFunctionContextDestroy(DM dm, PetscErrorCode (*f)(void *)) 335 { 336 DMSNES sdm; 337 338 PetscFunctionBegin; 339 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 340 PetscCall(DMGetDMSNESWrite(dm, &sdm)); 341 if (sdm->functionctxcontainer) PetscCall(PetscContainerSetUserDestroy(sdm->functionctxcontainer, f)); 342 PetscFunctionReturn(PETSC_SUCCESS); 343 } 344 345 PetscErrorCode DMSNESUnsetFunctionContext_Internal(DM dm) 346 { 347 DMSNES sdm; 348 349 PetscFunctionBegin; 350 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 351 PetscCall(DMGetDMSNESWrite(dm, &sdm)); 352 PetscCall(DMSNESUnsetFunctionContext_DMSNES(sdm)); 353 PetscFunctionReturn(PETSC_SUCCESS); 354 } 355 356 /*@C 357 DMSNESSetMFFunction - set `SNES` residual evaluation function used in applying the matrix-free Jacobian with -snes_mf_operator 358 359 Logically Collective 360 361 Input Parameters: 362 + dm - `DM` to be used with `SNES` 363 . func - residual evaluation function; see `SNESFunction` for details 364 - ctx - optional function context 365 366 Calling sequence of `func`: 367 $ PetscErrorCode func(SNES snes, Vec x, Vec f, void *ctx); 368 + snes - the solver object 369 . x - the input vector 370 . f - the output vector 371 - ctx - the optional function context 372 373 Level: advanced 374 375 .seealso: `DMSNESSetContext()`, `SNESSetFunction()`, `DMSNESSetJacobian()`, `SNESFunction`, `DMSNESSetFunction()` 376 @*/ 377 PetscErrorCode DMSNESSetMFFunction(DM dm, PetscErrorCode (*func)(SNES snes, Vec x, Vec f, void *ctx), void *ctx) 378 { 379 DMSNES sdm; 380 381 PetscFunctionBegin; 382 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 383 if (func || ctx) PetscCall(DMGetDMSNESWrite(dm, &sdm)); 384 if (func) sdm->ops->computemffunction = func; 385 if (ctx) sdm->mffunctionctx = ctx; 386 PetscFunctionReturn(PETSC_SUCCESS); 387 } 388 389 /*@C 390 DMSNESGetFunction - get `SNES` residual evaluation function 391 392 Not Collective 393 394 Input Parameter: 395 . dm - `DM` to be used with `SNES` 396 397 Output Parameters: 398 + f - residual evaluation function; see `SNESFunction` for details 399 - ctx - context for residual evaluation 400 401 Level: advanced 402 403 Note: 404 `SNESGetFunction()` is normally used, but it calls this function internally because the user context is actually 405 associated with the `DM`. 406 407 .seealso: `DMSNESSetContext()`, `DMSNESSetFunction()`, `SNESSetFunction()`, `SNESFunction` 408 @*/ 409 PetscErrorCode DMSNESGetFunction(DM dm, PetscErrorCode (**f)(SNES, Vec, Vec, void *), void **ctx) 410 { 411 DMSNES sdm; 412 413 PetscFunctionBegin; 414 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 415 PetscCall(DMGetDMSNES(dm, &sdm)); 416 if (f) *f = sdm->ops->computefunction; 417 if (ctx) { 418 if (sdm->functionctxcontainer) PetscCall(PetscContainerGetPointer(sdm->functionctxcontainer, ctx)); 419 else *ctx = NULL; 420 } 421 PetscFunctionReturn(PETSC_SUCCESS); 422 } 423 424 /*@C 425 DMSNESSetObjective - set `SNES` objective evaluation function 426 427 Not Collective 428 429 Input Parameters: 430 + dm - `DM` to be used with `SNES` 431 . obj - objective evaluation function; see `SNESObjectiveFunction` for details 432 - ctx - context for residual evaluation 433 434 Level: advanced 435 436 .seealso: `DMSNESSetContext()`, `SNESGetObjective()`, `DMSNESSetFunction()` 437 @*/ 438 PetscErrorCode DMSNESSetObjective(DM dm, PetscErrorCode (*obj)(SNES, Vec, PetscReal *, void *), void *ctx) 439 { 440 DMSNES sdm; 441 442 PetscFunctionBegin; 443 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 444 if (obj || ctx) PetscCall(DMGetDMSNESWrite(dm, &sdm)); 445 if (obj) sdm->ops->computeobjective = obj; 446 if (ctx) sdm->objectivectx = ctx; 447 PetscFunctionReturn(PETSC_SUCCESS); 448 } 449 450 /*@C 451 DMSNESGetObjective - get `SNES` objective evaluation function 452 453 Not Collective 454 455 Input Parameter: 456 . dm - `DM` to be used with `SNES` 457 458 Output Parameters: 459 + obj - residual evaluation function; see `SNESObjectiveFunction` for details 460 - ctx - context for residual evaluation 461 462 Level: advanced 463 464 Note: 465 `SNESGetFunction()` is normally used, but it calls this function internally because the user context is actually 466 associated with the `DM`. 467 468 .seealso: `DMSNESSetContext()`, `DMSNESSetObjective()`, `SNESSetFunction()` 469 @*/ 470 PetscErrorCode DMSNESGetObjective(DM dm, PetscErrorCode (**obj)(SNES, Vec, PetscReal *, void *), void **ctx) 471 { 472 DMSNES sdm; 473 474 PetscFunctionBegin; 475 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 476 PetscCall(DMGetDMSNES(dm, &sdm)); 477 if (obj) *obj = sdm->ops->computeobjective; 478 if (ctx) *ctx = sdm->objectivectx; 479 PetscFunctionReturn(PETSC_SUCCESS); 480 } 481 482 /*@C 483 DMSNESSetNGS - set `SNES` Gauss-Seidel relaxation function 484 485 Not Collective 486 487 Input Parameters: 488 + dm - `DM` to be used with `SNES` 489 . f - relaxation function, see `SNESGSFunction` 490 - ctx - context for residual evaluation 491 492 Level: advanced 493 494 Note: 495 `SNESSetNGS()` is normally used, but it calls this function internally because the user context is actually 496 associated with the `DM`. This makes the interface consistent regardless of whether the user interacts with a `DM` or 497 not. 498 499 Developer Notes: 500 If `DM` took a more central role at some later date, this could become the primary method of supplying the smoother 501 502 .seealso: `DMSNESSetContext()`, `SNESSetFunction()`, `DMSNESSetJacobian()`, `DMSNESSetFunction()`, `SNESGSFunction` 503 @*/ 504 PetscErrorCode DMSNESSetNGS(DM dm, PetscErrorCode (*f)(SNES, Vec, Vec, void *), void *ctx) 505 { 506 DMSNES sdm; 507 508 PetscFunctionBegin; 509 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 510 if (f || ctx) PetscCall(DMGetDMSNESWrite(dm, &sdm)); 511 if (f) sdm->ops->computegs = f; 512 if (ctx) sdm->gsctx = ctx; 513 PetscFunctionReturn(PETSC_SUCCESS); 514 } 515 516 /*@C 517 DMSNESGetNGS - get `SNES` Gauss-Seidel relaxation function 518 519 Not Collective 520 521 Input Parameter: 522 . dm - `DM` to be used with `SNES` 523 524 Output Parameters: 525 + f - relaxation function which performs Gauss-Seidel sweeps, see `SNESSetNGS()` 526 - ctx - context for residual evaluation 527 528 Level: advanced 529 530 Note: 531 `SNESGetNGS()` is normally used, but it calls this function internally because the user context is actually 532 associated with the `DM`. 533 534 Developer Notes: 535 This makes the interface consistent regardless of whether the user interacts with a `DM` or 536 not. If `DM` took a more central role at some later date, this could become the primary method of setting the residual. 537 538 .seealso: `DMSNESSetContext()`, `SNESGetNGS()`, `DMSNESGetJacobian()`, `DMSNESGetFunction()` 539 @*/ 540 PetscErrorCode DMSNESGetNGS(DM dm, PetscErrorCode (**f)(SNES, Vec, Vec, void *), void **ctx) 541 { 542 DMSNES sdm; 543 544 PetscFunctionBegin; 545 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 546 PetscCall(DMGetDMSNES(dm, &sdm)); 547 if (f) *f = sdm->ops->computegs; 548 if (ctx) *ctx = sdm->gsctx; 549 PetscFunctionReturn(PETSC_SUCCESS); 550 } 551 552 /*@C 553 DMSNESSetJacobian - set `SNES` Jacobian evaluation function 554 555 Not Collective 556 557 Input Parameters: 558 + dm - `DM` to be used with `SNES` 559 . J - Jacobian evaluation function 560 - ctx - context for residual evaluation 561 562 Level: advanced 563 564 Note: 565 `SNESSetJacobian()` is normally used, but it calls this function internally because the user context is actually 566 associated with the `DM`. 567 568 Developer Notes: 569 This makes the interface consistent regardless of whether the user interacts with a `DM` or 570 not. If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian. 571 572 .seealso: `DMSNESSetContext()`, `SNESSetFunction()`, `DMSNESGetJacobian()`, `SNESSetJacobian()`, `SNESJacobianFunction` 573 @*/ 574 PetscErrorCode DMSNESSetJacobian(DM dm, PetscErrorCode (*J)(SNES, Vec, Mat, Mat, void *), void *ctx) 575 { 576 DMSNES sdm; 577 578 PetscFunctionBegin; 579 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 580 if (J || ctx) PetscCall(DMGetDMSNESWrite(dm, &sdm)); 581 if (J) sdm->ops->computejacobian = J; 582 if (ctx) { 583 PetscContainer ctxcontainer; 584 PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)sdm), &ctxcontainer)); 585 PetscCall(PetscContainerSetPointer(ctxcontainer, ctx)); 586 PetscCall(PetscObjectCompose((PetscObject)sdm, "jacobian ctx", (PetscObject)ctxcontainer)); 587 sdm->jacobianctxcontainer = ctxcontainer; 588 PetscCall(PetscContainerDestroy(&ctxcontainer)); 589 } 590 PetscFunctionReturn(PETSC_SUCCESS); 591 } 592 593 /*@C 594 DMSNESSetJacobianContextDestroy - set `SNES` Jacobian evaluation context destroy function 595 596 Not Collective 597 598 Input Parameters: 599 + dm - `DM` to be used with `SNES` 600 - f - Jacobian evaluation context destroy function 601 602 Level: advanced 603 604 .seealso: `DMSNESSetJacobian()` 605 @*/ 606 PetscErrorCode DMSNESSetJacobianContextDestroy(DM dm, PetscErrorCode (*f)(void *)) 607 { 608 DMSNES sdm; 609 610 PetscFunctionBegin; 611 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 612 PetscCall(DMGetDMSNESWrite(dm, &sdm)); 613 if (sdm->jacobianctxcontainer) PetscCall(PetscContainerSetUserDestroy(sdm->jacobianctxcontainer, f)); 614 PetscFunctionReturn(PETSC_SUCCESS); 615 } 616 617 PetscErrorCode DMSNESUnsetJacobianContext_Internal(DM dm) 618 { 619 DMSNES sdm; 620 621 PetscFunctionBegin; 622 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 623 PetscCall(DMGetDMSNESWrite(dm, &sdm)); 624 PetscCall(DMSNESUnsetJacobianContext_DMSNES(sdm)); 625 PetscFunctionReturn(PETSC_SUCCESS); 626 } 627 628 /*@C 629 DMSNESGetJacobian - get `SNES` Jacobian evaluation function 630 631 Not Collective 632 633 Input Parameter: 634 . dm - `DM` to be used with `SNES` 635 636 Output Parameters: 637 + J - Jacobian evaluation function; for all calling sequence see `SNESJacobianFunction` 638 - ctx - context for residual evaluation 639 640 Level: advanced 641 642 Note: 643 `SNESGetJacobian()` is normally used, but it calls this function internally because the user context is actually 644 associated with the `DM`. This makes the interface consistent regardless of whether the user interacts with a `DM` or 645 not. 646 647 Developer Notes: 648 If `DM` took a more central role at some later date, this could become the primary method of setting the Jacobian. 649 650 .seealso: `DMSNESSetContext()`, `SNESSetFunction()`, `DMSNESSetJacobian()`, `SNESJacobianFunction` 651 @*/ 652 PetscErrorCode DMSNESGetJacobian(DM dm, PetscErrorCode (**J)(SNES, Vec, Mat, Mat, void *), void **ctx) 653 { 654 DMSNES sdm; 655 656 PetscFunctionBegin; 657 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 658 PetscCall(DMGetDMSNES(dm, &sdm)); 659 if (J) *J = sdm->ops->computejacobian; 660 if (ctx) { 661 if (sdm->jacobianctxcontainer) PetscCall(PetscContainerGetPointer(sdm->jacobianctxcontainer, ctx)); 662 else *ctx = NULL; 663 } 664 PetscFunctionReturn(PETSC_SUCCESS); 665 } 666 667 /*@C 668 DMSNESSetPicard - set SNES Picard iteration matrix and RHS evaluation functions. 669 670 Not Collective 671 672 Input Parameters: 673 + dm - `DM` to be used with `SNES` 674 . b - RHS evaluation function 675 . J - Picard matrix evaluation function 676 - ctx - context for residual evaluation 677 678 Level: advanced 679 680 .seealso: `SNESSetPicard()`, `DMSNESSetFunction()`, `DMSNESSetJacobian()` 681 @*/ 682 PetscErrorCode DMSNESSetPicard(DM dm, PetscErrorCode (*b)(SNES, Vec, Vec, void *), PetscErrorCode (*J)(SNES, Vec, Mat, Mat, void *), void *ctx) 683 { 684 DMSNES sdm; 685 686 PetscFunctionBegin; 687 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 688 PetscCall(DMGetDMSNES(dm, &sdm)); 689 if (b) sdm->ops->computepfunction = b; 690 if (J) sdm->ops->computepjacobian = J; 691 if (ctx) sdm->pctx = ctx; 692 PetscFunctionReturn(PETSC_SUCCESS); 693 } 694 695 /*@C 696 DMSNESGetPicard - get `SNES` Picard iteration evaluation functions 697 698 Not Collective 699 700 Input Parameter: 701 . dm - `DM` to be used with `SNES` 702 703 Output Parameters: 704 + b - RHS evaluation function; see `SNESFunction` for details 705 . J - RHS evaluation function; see `SNESJacobianFunction` for detailsa 706 - ctx - context for residual evaluation 707 708 Level: advanced 709 710 .seealso: `DMSNESSetContext()`, `SNESSetFunction()`, `DMSNESSetJacobian()` 711 @*/ 712 PetscErrorCode DMSNESGetPicard(DM dm, PetscErrorCode (**b)(SNES, Vec, Vec, void *), PetscErrorCode (**J)(SNES, Vec, Mat, Mat, void *), void **ctx) 713 { 714 DMSNES sdm; 715 716 PetscFunctionBegin; 717 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 718 PetscCall(DMGetDMSNES(dm, &sdm)); 719 if (b) *b = sdm->ops->computepfunction; 720 if (J) *J = sdm->ops->computepjacobian; 721 if (ctx) *ctx = sdm->pctx; 722 PetscFunctionReturn(PETSC_SUCCESS); 723 } 724