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