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