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