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