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. If DM took a more central role at some later date, this could become the primary method of setting the residual. 334 335 .seealso: `DMSNESSetFunction()`, `SNESSetFunction()` 336 @*/ 337 PetscErrorCode DMSNESSetFunctionContextDestroy(DM dm,PetscErrorCode (*f)(void*)) 338 { 339 DMSNES sdm; 340 341 PetscFunctionBegin; 342 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 343 PetscCall(DMGetDMSNESWrite(dm,&sdm)); 344 if (sdm->functionctxcontainer) PetscCall(PetscContainerSetUserDestroy(sdm->functionctxcontainer,f)); 345 PetscFunctionReturn(0); 346 } 347 348 PetscErrorCode DMSNESUnsetFunctionContext_Internal(DM dm) 349 { 350 DMSNES sdm; 351 352 PetscFunctionBegin; 353 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 354 PetscCall(DMGetDMSNESWrite(dm,&sdm)); 355 PetscCall(DMSNESUnsetFunctionContext_DMSNES(sdm)); 356 PetscFunctionReturn(0); 357 } 358 359 /*@C 360 DMSNESSetMFFunction - set SNES residual evaluation function used in applying the matrix-free Jacobian with -snes_mf_operator 361 362 Logically Collective on dm 363 364 Input Parameters: 365 + dm - DM to be used with SNES 366 - f - residual evaluation function; see SNESFunction for details 367 368 Level: advanced 369 370 .seealso: `DMSNESSetContext()`, `SNESSetFunction()`, `DMSNESSetJacobian()`, `SNESFunction`, `DMSNESSetFunction()` 371 @*/ 372 PetscErrorCode DMSNESSetMFFunction(DM dm,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx) 373 { 374 DMSNES sdm; 375 376 PetscFunctionBegin; 377 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 378 if (f || ctx) { 379 PetscCall(DMGetDMSNESWrite(dm,&sdm)); 380 } 381 if (f) sdm->ops->computemffunction = f; 382 if (ctx) sdm->mffunctionctx = ctx; 383 PetscFunctionReturn(0); 384 } 385 386 /*@C 387 DMSNESGetFunction - get SNES residual evaluation function 388 389 Not Collective 390 391 Input Parameter: 392 . dm - DM to be used with SNES 393 394 Output Parameters: 395 + f - residual evaluation function; see SNESFunction for details 396 - ctx - context for residual evaluation 397 398 Level: advanced 399 400 Note: 401 SNESGetFunction() is normally used, but it calls this function internally because the user context is actually 402 associated with the DM. 403 404 .seealso: `DMSNESSetContext()`, `DMSNESSetFunction()`, `SNESSetFunction()`, `SNESFunction` 405 @*/ 406 PetscErrorCode DMSNESGetFunction(DM dm,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx) 407 { 408 DMSNES sdm; 409 410 PetscFunctionBegin; 411 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 412 PetscCall(DMGetDMSNES(dm,&sdm)); 413 if (f) *f = sdm->ops->computefunction; 414 if (ctx) { 415 if (sdm->functionctxcontainer) PetscCall(PetscContainerGetPointer(sdm->functionctxcontainer,ctx)); 416 else *ctx = NULL; 417 } 418 PetscFunctionReturn(0); 419 } 420 421 /*@C 422 DMSNESSetObjective - set SNES objective evaluation function 423 424 Not Collective 425 426 Input Parameters: 427 + dm - DM to be used with SNES 428 . obj - objective evaluation function; see SNESObjectiveFunction for details 429 - ctx - context for residual evaluation 430 431 Level: advanced 432 433 .seealso: `DMSNESSetContext()`, `SNESGetObjective()`, `DMSNESSetFunction()` 434 @*/ 435 PetscErrorCode DMSNESSetObjective(DM dm,PetscErrorCode (*obj)(SNES,Vec,PetscReal*,void*),void *ctx) 436 { 437 DMSNES sdm; 438 439 PetscFunctionBegin; 440 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 441 if (obj || ctx) { 442 PetscCall(DMGetDMSNESWrite(dm,&sdm)); 443 } 444 if (obj) sdm->ops->computeobjective = obj; 445 if (ctx) sdm->objectivectx = ctx; 446 PetscFunctionReturn(0); 447 } 448 449 /*@C 450 DMSNESGetObjective - get SNES objective evaluation function 451 452 Not Collective 453 454 Input Parameter: 455 . dm - DM to be used with SNES 456 457 Output Parameters: 458 + obj- residual evaluation function; see SNESObjectiveFunction for details 459 - ctx - context for residual evaluation 460 461 Level: advanced 462 463 Note: 464 SNESGetFunction() is normally used, but it calls this function internally because the user context is actually 465 associated with the DM. 466 467 .seealso: `DMSNESSetContext()`, `DMSNESSetObjective()`, `SNESSetFunction()` 468 @*/ 469 PetscErrorCode DMSNESGetObjective(DM dm,PetscErrorCode (**obj)(SNES,Vec,PetscReal*,void*),void **ctx) 470 { 471 DMSNES sdm; 472 473 PetscFunctionBegin; 474 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 475 PetscCall(DMGetDMSNES(dm,&sdm)); 476 if (obj) *obj = sdm->ops->computeobjective; 477 if (ctx) *ctx = sdm->objectivectx; 478 PetscFunctionReturn(0); 479 } 480 481 /*@C 482 DMSNESSetNGS - set SNES Gauss-Seidel relaxation function 483 484 Not Collective 485 486 Input Parameters: 487 + dm - DM to be used with SNES 488 . f - relaxation function, see SNESGSFunction 489 - ctx - context for residual evaluation 490 491 Level: advanced 492 493 Note: 494 SNESSetNGS() 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 residual. 497 498 .seealso: `DMSNESSetContext()`, `SNESSetFunction()`, `DMSNESSetJacobian()`, `DMSNESSetFunction()`, `SNESGSFunction` 499 @*/ 500 PetscErrorCode DMSNESSetNGS(DM dm,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx) 501 { 502 DMSNES sdm; 503 504 PetscFunctionBegin; 505 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 506 if (f || ctx) { 507 PetscCall(DMGetDMSNESWrite(dm,&sdm)); 508 } 509 if (f) sdm->ops->computegs = f; 510 if (ctx) sdm->gsctx = ctx; 511 PetscFunctionReturn(0); 512 } 513 514 /*@C 515 DMSNESGetNGS - get SNES Gauss-Seidel relaxation function 516 517 Not Collective 518 519 Input Parameter: 520 . dm - DM to be used with SNES 521 522 Output Parameters: 523 + f - relaxation function which performs Gauss-Seidel sweeps, see SNESGSFunction 524 - ctx - context for residual evaluation 525 526 Level: advanced 527 528 Note: 529 SNESGetNGS() is normally used, but it calls this function internally because the user context is actually 530 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 531 not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 532 533 .seealso: `DMSNESSetContext()`, `SNESGetNGS()`, `DMSNESGetJacobian()`, `DMSNESGetFunction()`, `SNESNGSFunction` 534 @*/ 535 PetscErrorCode DMSNESGetNGS(DM dm,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx) 536 { 537 DMSNES sdm; 538 539 PetscFunctionBegin; 540 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 541 PetscCall(DMGetDMSNES(dm,&sdm)); 542 if (f) *f = sdm->ops->computegs; 543 if (ctx) *ctx = sdm->gsctx; 544 PetscFunctionReturn(0); 545 } 546 547 /*@C 548 DMSNESSetJacobian - set SNES Jacobian evaluation function 549 550 Not Collective 551 552 Input Parameters: 553 + dm - DM to be used with SNES 554 . J - Jacobian evaluation function 555 - ctx - context for residual evaluation 556 557 Level: advanced 558 559 Note: 560 SNESSetJacobian() is normally used, but it calls this function internally because the user context is actually 561 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 562 not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 563 564 .seealso: `DMSNESSetContext()`, `SNESSetFunction()`, `DMSNESGetJacobian()`, `SNESSetJacobian()`, `SNESJacobianFunction` 565 @*/ 566 PetscErrorCode DMSNESSetJacobian(DM dm,PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx) 567 { 568 DMSNES sdm; 569 570 PetscFunctionBegin; 571 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 572 if (J || ctx) PetscCall(DMGetDMSNESWrite(dm,&sdm)); 573 if (J) sdm->ops->computejacobian = J; 574 if (ctx) { 575 PetscContainer ctxcontainer; 576 PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)sdm),&ctxcontainer)); 577 PetscCall(PetscContainerSetPointer(ctxcontainer,ctx)); 578 PetscCall(PetscObjectCompose((PetscObject)sdm,"jacobian ctx",(PetscObject)ctxcontainer)); 579 sdm->jacobianctxcontainer = ctxcontainer; 580 PetscCall(PetscContainerDestroy(&ctxcontainer)); 581 } 582 PetscFunctionReturn(0); 583 } 584 585 /*@C 586 DMSNESSetJacobianContextDestroy - set SNES Jacobian evaluation context destroy function 587 588 Not Collective 589 590 Input Parameters: 591 + dm - DM to be used with SNES 592 . f - Jacobian evaluation contex destroy function 593 594 Level: advanced 595 596 Note: 597 SNESSetJacobianContextDestroy() is normally used, but it calls this function internally because the user context is actually 598 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 599 not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 600 601 .seealso: `DMSNESSetJacobian()`, `SNESSetJacobianContextDestroyFunction()` 602 @*/ 603 PetscErrorCode DMSNESSetJacobianContextDestroy(DM dm,PetscErrorCode (*f)(void*)) 604 { 605 DMSNES sdm; 606 607 PetscFunctionBegin; 608 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 609 PetscCall(DMGetDMSNESWrite(dm,&sdm)); 610 if (sdm->jacobianctxcontainer) PetscCall(PetscContainerSetUserDestroy(sdm->jacobianctxcontainer,f)); 611 PetscFunctionReturn(0); 612 } 613 614 PetscErrorCode DMSNESUnsetJacobianContext_Internal(DM dm) 615 { 616 DMSNES sdm; 617 618 PetscFunctionBegin; 619 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 620 PetscCall(DMGetDMSNESWrite(dm,&sdm)); 621 PetscCall(DMSNESUnsetJacobianContext_DMSNES(sdm)); 622 PetscFunctionReturn(0); 623 } 624 625 /*@C 626 DMSNESGetJacobian - get SNES Jacobian evaluation function 627 628 Not Collective 629 630 Input Parameter: 631 . dm - DM to be used with SNES 632 633 Output Parameters: 634 + J - Jacobian evaluation function; see SNESJacobianFunction for all calling sequence 635 - ctx - context for residual evaluation 636 637 Level: advanced 638 639 Note: 640 SNESGetJacobian() is normally used, but it calls this function internally because the user context is actually 641 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 642 not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 643 644 .seealso: `DMSNESSetContext()`, `SNESSetFunction()`, `DMSNESSetJacobian()`, `SNESJacobianFunction` 645 @*/ 646 PetscErrorCode DMSNESGetJacobian(DM dm,PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx) 647 { 648 DMSNES sdm; 649 650 PetscFunctionBegin; 651 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 652 PetscCall(DMGetDMSNES(dm,&sdm)); 653 if (J) *J = sdm->ops->computejacobian; 654 if (ctx) { 655 if (sdm->jacobianctxcontainer) PetscCall(PetscContainerGetPointer(sdm->jacobianctxcontainer,ctx)); 656 else *ctx = NULL; 657 } 658 PetscFunctionReturn(0); 659 } 660 661 /*@C 662 DMSNESSetPicard - set SNES Picard iteration matrix and RHS evaluation functions. 663 664 Not Collective 665 666 Input Parameters: 667 + dm - DM to be used with SNES 668 . b - RHS evaluation function 669 . J - Picard matrix evaluation function 670 - ctx - context for residual evaluation 671 672 Level: advanced 673 674 .seealso: `SNESSetPicard()`, `DMSNESSetFunction()`, `DMSNESSetJacobian()` 675 @*/ 676 PetscErrorCode DMSNESSetPicard(DM dm,PetscErrorCode (*b)(SNES,Vec,Vec,void*),PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx) 677 { 678 DMSNES sdm; 679 680 PetscFunctionBegin; 681 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 682 PetscCall(DMGetDMSNES(dm,&sdm)); 683 if (b) sdm->ops->computepfunction = b; 684 if (J) sdm->ops->computepjacobian = J; 685 if (ctx) sdm->pctx = ctx; 686 PetscFunctionReturn(0); 687 } 688 689 /*@C 690 DMSNESGetPicard - get SNES Picard iteration evaluation functions 691 692 Not Collective 693 694 Input Parameter: 695 . dm - DM to be used with SNES 696 697 Output Parameters: 698 + b - RHS evaluation function; see SNESFunction for details 699 . J - RHS evaluation function; see SNESJacobianFunction for detailsa 700 - ctx - context for residual evaluation 701 702 Level: advanced 703 704 .seealso: `DMSNESSetContext()`, `SNESSetFunction()`, `DMSNESSetJacobian()` 705 @*/ 706 PetscErrorCode DMSNESGetPicard(DM dm,PetscErrorCode (**b)(SNES,Vec,Vec,void*),PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx) 707 { 708 DMSNES sdm; 709 710 PetscFunctionBegin; 711 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 712 PetscCall(DMGetDMSNES(dm,&sdm)); 713 if (b) *b = sdm->ops->computepfunction; 714 if (J) *J = sdm->ops->computepjacobian; 715 if (ctx) *ctx = sdm->pctx; 716 PetscFunctionReturn(0); 717 } 718