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