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