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,NULL,PETSC_FUNCTION);CHKERRQ(ierr); 27 ierr = PetscViewerBinaryRead(viewer,&kdm->ops->computejacobian,1,NULL,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 } funcstruct; 58 struct { 59 PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*); 60 } jacstruct; 61 funcstruct.func = kdm->ops->computefunction; 62 jacstruct.jac = kdm->ops->computejacobian; 63 ierr = PetscViewerBinaryWrite(viewer,&funcstruct,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr); 64 ierr = PetscViewerBinaryWrite(viewer,&jacstruct,1,PETSC_FUNCTION,PETSC_FALSE);CHKERRQ(ierr); 65 } 66 PetscFunctionReturn(0); 67 } 68 69 #undef __FUNCT__ 70 #define __FUNCT__ "DMSNESCreate" 71 static PetscErrorCode DMSNESCreate(MPI_Comm comm,DMSNES *kdm) 72 { 73 PetscErrorCode ierr; 74 75 PetscFunctionBegin; 76 ierr = SNESInitializePackage();CHKERRQ(ierr); 77 ierr = PetscHeaderCreate(*kdm, _p_DMSNES, struct _DMSNESOps, DMSNES_CLASSID, "DMSNES", "DMSNES", "DMSNES", comm, DMSNESDestroy, DMSNESView);CHKERRQ(ierr); 78 ierr = PetscMemzero((*kdm)->ops, sizeof(struct _DMSNESOps));CHKERRQ(ierr); 79 PetscFunctionReturn(0); 80 } 81 82 #undef __FUNCT__ 83 #define __FUNCT__ "DMCoarsenHook_DMSNES" 84 /* Attaches the DMSNES to the coarse level. 85 * Under what conditions should we copy versus duplicate? 86 */ 87 static PetscErrorCode DMCoarsenHook_DMSNES(DM dm,DM dmc,void *ctx) 88 { 89 PetscErrorCode ierr; 90 91 PetscFunctionBegin; 92 ierr = DMCopyDMSNES(dm,dmc);CHKERRQ(ierr); 93 PetscFunctionReturn(0); 94 } 95 96 #undef __FUNCT__ 97 #define __FUNCT__ "DMRestrictHook_DMSNES" 98 /* This could restrict auxiliary information to the coarse level. 99 */ 100 static PetscErrorCode DMRestrictHook_DMSNES(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx) 101 { 102 103 PetscFunctionBegin; 104 PetscFunctionReturn(0); 105 } 106 107 #undef __FUNCT__ 108 #define __FUNCT__ "DMSubDomainHook_DMSNES" 109 /* Attaches the DMSNES to the subdomain. */ 110 static PetscErrorCode DMSubDomainHook_DMSNES(DM dm,DM subdm,void *ctx) 111 { 112 PetscErrorCode ierr; 113 114 PetscFunctionBegin; 115 ierr = DMCopyDMSNES(dm,subdm);CHKERRQ(ierr); 116 PetscFunctionReturn(0); 117 } 118 119 #undef __FUNCT__ 120 #define __FUNCT__ "DMSubDomainRestrictHook_DMSNES" 121 /* This could restrict auxiliary information to the coarse level. 122 */ 123 static PetscErrorCode DMSubDomainRestrictHook_DMSNES(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 124 { 125 126 PetscFunctionBegin; 127 PetscFunctionReturn(0); 128 } 129 130 #undef __FUNCT__ 131 #define __FUNCT__ "DMRefineHook_DMSNES" 132 static PetscErrorCode DMRefineHook_DMSNES(DM dm,DM dmf,void *ctx) 133 { 134 PetscErrorCode ierr; 135 136 PetscFunctionBegin; 137 ierr = DMCopyDMSNES(dm,dmf);CHKERRQ(ierr); 138 PetscFunctionReturn(0); 139 } 140 141 #undef __FUNCT__ 142 #define __FUNCT__ "DMInterpolateHook_DMSNES" 143 /* This could restrict auxiliary information to the coarse level. 144 */ 145 static PetscErrorCode DMInterpolateHook_DMSNES(DM dm,Mat Interp,DM dmf,void *ctx) 146 { 147 148 PetscFunctionBegin; 149 PetscFunctionReturn(0); 150 } 151 152 #undef __FUNCT__ 153 #define __FUNCT__ "DMSNESCopy" 154 /*@C 155 DMSNESCopy - copies the information in a DMSNES to another DMSNES 156 157 Not Collective 158 159 Input Argument: 160 + kdm - Original DMSNES 161 - nkdm - DMSNES to receive the data, should have been created with DMSNESCreate() 162 163 Level: developer 164 165 .seealso: DMSNESCreate(), DMSNESDestroy() 166 @*/ 167 PetscErrorCode DMSNESCopy(DMSNES kdm,DMSNES nkdm) 168 { 169 PetscErrorCode ierr; 170 171 PetscFunctionBegin; 172 PetscValidHeaderSpecific(kdm,DMSNES_CLASSID,1); 173 PetscValidHeaderSpecific(nkdm,DMSNES_CLASSID,2); 174 nkdm->ops->computefunction = kdm->ops->computefunction; 175 nkdm->ops->computejacobian = kdm->ops->computejacobian; 176 nkdm->ops->computegs = kdm->ops->computegs; 177 nkdm->ops->computeobjective = kdm->ops->computeobjective; 178 nkdm->ops->computepjacobian = kdm->ops->computepjacobian; 179 nkdm->ops->computepfunction = kdm->ops->computepfunction; 180 nkdm->ops->destroy = kdm->ops->destroy; 181 nkdm->ops->duplicate = kdm->ops->duplicate; 182 183 nkdm->functionctx = kdm->functionctx; 184 nkdm->gsctx = kdm->gsctx; 185 nkdm->pctx = kdm->pctx; 186 nkdm->jacobianctx = kdm->jacobianctx; 187 nkdm->objectivectx = kdm->objectivectx; 188 nkdm->data = kdm->data; 189 190 /* 191 nkdm->fortran_func_pointers[0] = kdm->fortran_func_pointers[0]; 192 nkdm->fortran_func_pointers[1] = kdm->fortran_func_pointers[1]; 193 nkdm->fortran_func_pointers[2] = kdm->fortran_func_pointers[2]; 194 */ 195 196 /* implementation specific copy hooks */ 197 if (kdm->ops->duplicate) {ierr = (*kdm->ops->duplicate)(kdm,nkdm);CHKERRQ(ierr);} 198 PetscFunctionReturn(0); 199 } 200 201 #undef __FUNCT__ 202 #define __FUNCT__ "DMGetDMSNES" 203 /*@C 204 DMGetDMSNES - get read-only private DMSNES context from a DM 205 206 Not Collective 207 208 Input Argument: 209 . dm - DM to be used with SNES 210 211 Output Argument: 212 . snesdm - private DMSNES context 213 214 Level: developer 215 216 Notes: 217 Use DMGetDMSNESWrite() if write access is needed. The DMSNESSetXXX API should be used wherever possible. 218 219 .seealso: DMGetDMSNESWrite() 220 @*/ 221 PetscErrorCode DMGetDMSNES(DM dm,DMSNES *snesdm) 222 { 223 PetscErrorCode ierr; 224 225 PetscFunctionBegin; 226 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 227 *snesdm = (DMSNES) dm->dmsnes; 228 if (!*snesdm) { 229 ierr = PetscInfo(dm,"Creating new DMSNES\n");CHKERRQ(ierr); 230 ierr = DMSNESCreate(PetscObjectComm((PetscObject)dm),snesdm);CHKERRQ(ierr); 231 232 dm->dmsnes = (PetscObject) *snesdm; 233 234 ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_DMSNES,DMRestrictHook_DMSNES,NULL);CHKERRQ(ierr); 235 ierr = DMRefineHookAdd(dm,DMRefineHook_DMSNES,DMInterpolateHook_DMSNES,NULL);CHKERRQ(ierr); 236 ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_DMSNES,DMSubDomainRestrictHook_DMSNES,NULL);CHKERRQ(ierr); 237 } 238 PetscFunctionReturn(0); 239 } 240 241 #undef __FUNCT__ 242 #define __FUNCT__ "DMGetDMSNESWrite" 243 /*@C 244 DMGetDMSNESWrite - get write access to private DMSNES context from a DM 245 246 Not Collective 247 248 Input Argument: 249 . dm - DM to be used with SNES 250 251 Output Argument: 252 . snesdm - private DMSNES context 253 254 Level: developer 255 256 .seealso: DMGetDMSNES() 257 @*/ 258 PetscErrorCode DMGetDMSNESWrite(DM dm,DMSNES *snesdm) 259 { 260 PetscErrorCode ierr; 261 DMSNES sdm; 262 263 PetscFunctionBegin; 264 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 265 ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr); 266 if (!sdm->originaldm) sdm->originaldm = dm; 267 if (sdm->originaldm != dm) { /* Copy on write */ 268 DMSNES oldsdm = sdm; 269 ierr = PetscInfo(dm,"Copying DMSNES due to write\n");CHKERRQ(ierr); 270 ierr = DMSNESCreate(PetscObjectComm((PetscObject)dm),&sdm);CHKERRQ(ierr); 271 ierr = DMSNESCopy(oldsdm,sdm);CHKERRQ(ierr); 272 ierr = DMSNESDestroy((DMSNES*)&dm->dmsnes);CHKERRQ(ierr); 273 dm->dmsnes = (PetscObject)sdm; 274 } 275 *snesdm = sdm; 276 PetscFunctionReturn(0); 277 } 278 279 #undef __FUNCT__ 280 #define __FUNCT__ "DMCopyDMSNES" 281 /*@C 282 DMCopyDMSNES - copies a DM context to a new DM 283 284 Logically Collective 285 286 Input Arguments: 287 + dmsrc - DM to obtain context from 288 - dmdest - DM to add context to 289 290 Level: developer 291 292 Note: 293 The context is copied by reference. This function does not ensure that a context exists. 294 295 .seealso: DMGetDMSNES(), SNESSetDM() 296 @*/ 297 PetscErrorCode DMCopyDMSNES(DM dmsrc,DM dmdest) 298 { 299 PetscErrorCode ierr; 300 301 PetscFunctionBegin; 302 PetscValidHeaderSpecific(dmsrc,DM_CLASSID,1); 303 PetscValidHeaderSpecific(dmdest,DM_CLASSID,2); 304 if (!dmdest->dmsnes) {ierr = DMSNESCreate(PetscObjectComm((PetscObject) dmdest), (DMSNES *) &dmdest->dmsnes);CHKERRQ(ierr);} 305 ierr = DMSNESCopy((DMSNES) dmsrc->dmsnes, (DMSNES) 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 . f - residual evaluation function; see SNESFunction for details 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 (*f)(SNES,Vec,Vec,void*),void *ctx) 334 { 335 PetscErrorCode ierr; 336 DMSNES sdm; 337 338 PetscFunctionBegin; 339 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 340 if (f || ctx) { 341 ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 342 } 343 if (f) sdm->ops->computefunction = f; 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 + f - residual evaluation function; see SNESFunction for details 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 (**f)(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 (f) *f = 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 . obj - objective evaluation function; see SNESObjectiveFunction for details 393 - ctx - context for residual evaluation 394 395 Level: advanced 396 397 .seealso: DMSNESSetContext(), SNESGetObjective(), DMSNESSetFunction() 398 @*/ 399 PetscErrorCode DMSNESSetObjective(DM dm,PetscErrorCode (*obj)(SNES,Vec,PetscReal*,void*),void *ctx) 400 { 401 PetscErrorCode ierr; 402 DMSNES sdm; 403 404 PetscFunctionBegin; 405 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 406 if (obj || ctx) { 407 ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 408 } 409 if (obj) sdm->ops->computeobjective = obj; 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 + obj- residual evaluation function; see SNESObjectiveFunction for details 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 (**obj)(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 (obj) *obj = sdm->ops->computeobjective; 445 if (ctx) *ctx = sdm->objectivectx; 446 PetscFunctionReturn(0); 447 } 448 449 #undef __FUNCT__ 450 #define __FUNCT__ "DMSNESSetNGS" 451 /*@C 452 DMSNESSetNGS - set SNES Gauss-Seidel relaxation function 453 454 Not Collective 455 456 Input Argument: 457 + dm - DM to be used with SNES 458 . f - relaxation function, see SNESGSFunction 459 - ctx - context for residual evaluation 460 461 Level: advanced 462 463 Note: 464 SNESSetNGS() 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 DMSNESSetNGS(DM dm,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx) 471 { 472 PetscErrorCode ierr; 473 DMSNES sdm; 474 475 PetscFunctionBegin; 476 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 477 if (f || ctx) { 478 ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 479 } 480 if (f) sdm->ops->computegs = f; 481 if (ctx) sdm->gsctx = ctx; 482 PetscFunctionReturn(0); 483 } 484 485 #undef __FUNCT__ 486 #define __FUNCT__ "DMSNESGetNGS" 487 /*@C 488 DMSNESGetNGS - 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 + f - relaxation function which performs Gauss-Seidel sweeps, see SNESGSFunction 497 - ctx - context for residual evaluation 498 499 Level: advanced 500 501 Note: 502 SNESGetNGS() 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(), SNESGetNGS(), DMSNESGetJacobian(), DMSNESGetFunction(), SNESNGSFunction 507 @*/ 508 PetscErrorCode DMSNESGetNGS(DM dm,PetscErrorCode (**f)(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 (f) *f = 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 . J - 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 (*J)(SNES,Vec,Mat,Mat,void*),void *ctx) 543 { 544 PetscErrorCode ierr; 545 DMSNES sdm; 546 547 PetscFunctionBegin; 548 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 549 if (J || ctx) { 550 ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 551 } 552 if (J) sdm->ops->computejacobian = J; 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 + J - Jacobian evaluation function; see SNESJacobianFunction for all calling sequence 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 (**J)(SNES,Vec,Mat,Mat,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 (J) *J = 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 . b - RHS evaluation function 603 . J - 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 (*b)(SNES,Vec,Vec,void*),PetscErrorCode (*J)(SNES,Vec,Mat,Mat,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 (b) sdm->ops->computepfunction = b; 619 if (J) sdm->ops->computepjacobian = J; 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 + b - RHS evaluation function; see SNESFunction for details 636 . J - RHS evaluation function; see SNESJacobianFunction for detailsa 637 - ctx - context for residual evaluation 638 639 Level: advanced 640 641 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian() 642 @*/ 643 PetscErrorCode DMSNESGetPicard(DM dm,PetscErrorCode (**b)(SNES,Vec,Vec,void*),PetscErrorCode (**J)(SNES,Vec,Mat,Mat,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 (b) *b = sdm->ops->computepfunction; 652 if (J) *J = sdm->ops->computepjacobian; 653 if (ctx) *ctx = sdm->pctx; 654 PetscFunctionReturn(0); 655 } 656