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