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 } 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 ierr = DMSNESDestroy((DMSNES*)&dmdest->dmsnes);CHKERRQ(ierr); 305 306 dmdest->dmsnes = dmsrc->dmsnes; 307 308 ierr = PetscObjectReference(dmdest->dmsnes);CHKERRQ(ierr); 309 ierr = DMCoarsenHookAdd(dmdest,DMCoarsenHook_DMSNES,NULL,NULL);CHKERRQ(ierr); 310 ierr = DMRefineHookAdd(dmdest,DMRefineHook_DMSNES,NULL,NULL);CHKERRQ(ierr); 311 ierr = DMSubDomainHookAdd(dmdest,DMSubDomainHook_DMSNES,DMSubDomainRestrictHook_DMSNES,NULL);CHKERRQ(ierr); 312 PetscFunctionReturn(0); 313 } 314 315 #undef __FUNCT__ 316 #define __FUNCT__ "DMSNESSetFunction" 317 /*@C 318 DMSNESSetFunction - set SNES residual evaluation function 319 320 Not Collective 321 322 Input Arguments: 323 + dm - DM to be used with SNES 324 . f - residual evaluation function; see SNESFunction for details 325 - ctx - context for residual evaluation 326 327 Level: advanced 328 329 Note: 330 SNESSetFunction() is normally used, but it calls this function internally because the user context is actually 331 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 332 not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 333 334 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), SNESFunction 335 @*/ 336 PetscErrorCode DMSNESSetFunction(DM dm,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx) 337 { 338 PetscErrorCode ierr; 339 DMSNES sdm; 340 341 PetscFunctionBegin; 342 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 343 if (f || ctx) { 344 ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 345 } 346 if (f) sdm->ops->computefunction = f; 347 if (ctx) sdm->functionctx = ctx; 348 PetscFunctionReturn(0); 349 } 350 351 #undef __FUNCT__ 352 #define __FUNCT__ "DMSNESGetFunction" 353 /*@C 354 DMSNESGetFunction - get SNES residual evaluation function 355 356 Not Collective 357 358 Input Argument: 359 . dm - DM to be used with SNES 360 361 Output Arguments: 362 + f - residual evaluation function; see SNESFunction for details 363 - ctx - context for residual evaluation 364 365 Level: advanced 366 367 Note: 368 SNESGetFunction() is normally used, but it calls this function internally because the user context is actually 369 associated with the DM. 370 371 .seealso: DMSNESSetContext(), DMSNESSetFunction(), SNESSetFunction(), SNESFunction 372 @*/ 373 PetscErrorCode DMSNESGetFunction(DM dm,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx) 374 { 375 PetscErrorCode ierr; 376 DMSNES sdm; 377 378 PetscFunctionBegin; 379 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 380 ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr); 381 if (f) *f = sdm->ops->computefunction; 382 if (ctx) *ctx = sdm->functionctx; 383 PetscFunctionReturn(0); 384 } 385 386 #undef __FUNCT__ 387 #define __FUNCT__ "DMSNESSetObjective" 388 /*@C 389 DMSNESSetObjective - set SNES objective evaluation function 390 391 Not Collective 392 393 Input Arguments: 394 + dm - DM to be used with SNES 395 . obj - objective evaluation function; see SNESObjectiveFunction for details 396 - ctx - context for residual evaluation 397 398 Level: advanced 399 400 .seealso: DMSNESSetContext(), SNESGetObjective(), DMSNESSetFunction() 401 @*/ 402 PetscErrorCode DMSNESSetObjective(DM dm,PetscErrorCode (*obj)(SNES,Vec,PetscReal*,void*),void *ctx) 403 { 404 PetscErrorCode ierr; 405 DMSNES sdm; 406 407 PetscFunctionBegin; 408 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 409 if (obj || ctx) { 410 ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 411 } 412 if (obj) sdm->ops->computeobjective = obj; 413 if (ctx) sdm->objectivectx = ctx; 414 PetscFunctionReturn(0); 415 } 416 417 #undef __FUNCT__ 418 #define __FUNCT__ "DMSNESGetObjective" 419 /*@C 420 DMSNESGetObjective - get SNES objective evaluation function 421 422 Not Collective 423 424 Input Argument: 425 . dm - DM to be used with SNES 426 427 Output Arguments: 428 + obj- residual evaluation function; see SNESObjectiveFunction for details 429 - ctx - context for residual evaluation 430 431 Level: advanced 432 433 Note: 434 SNESGetFunction() is normally used, but it calls this function internally because the user context is actually 435 associated with the DM. 436 437 .seealso: DMSNESSetContext(), DMSNESSetObjective(), SNESSetFunction() 438 @*/ 439 PetscErrorCode DMSNESGetObjective(DM dm,PetscErrorCode (**obj)(SNES,Vec,PetscReal*,void*),void **ctx) 440 { 441 PetscErrorCode ierr; 442 DMSNES sdm; 443 444 PetscFunctionBegin; 445 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 446 ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr); 447 if (obj) *obj = sdm->ops->computeobjective; 448 if (ctx) *ctx = sdm->objectivectx; 449 PetscFunctionReturn(0); 450 } 451 452 #undef __FUNCT__ 453 #define __FUNCT__ "DMSNESSetNGS" 454 /*@C 455 DMSNESSetNGS - set SNES Gauss-Seidel relaxation function 456 457 Not Collective 458 459 Input Argument: 460 + dm - DM to be used with SNES 461 . f - relaxation function, see SNESGSFunction 462 - ctx - context for residual evaluation 463 464 Level: advanced 465 466 Note: 467 SNESSetNGS() is normally used, but it calls this function internally because the user context is actually 468 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 469 not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 470 471 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), DMSNESSetFunction(), SNESGSFunction 472 @*/ 473 PetscErrorCode DMSNESSetNGS(DM dm,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx) 474 { 475 PetscErrorCode ierr; 476 DMSNES sdm; 477 478 PetscFunctionBegin; 479 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 480 if (f || ctx) { 481 ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 482 } 483 if (f) sdm->ops->computegs = f; 484 if (ctx) sdm->gsctx = ctx; 485 PetscFunctionReturn(0); 486 } 487 488 #undef __FUNCT__ 489 #define __FUNCT__ "DMSNESGetNGS" 490 /*@C 491 DMSNESGetNGS - get SNES Gauss-Seidel relaxation function 492 493 Not Collective 494 495 Input Argument: 496 . dm - DM to be used with SNES 497 498 Output Arguments: 499 + f - relaxation function which performs Gauss-Seidel sweeps, see SNESGSFunction 500 - ctx - context for residual evaluation 501 502 Level: advanced 503 504 Note: 505 SNESGetNGS() is normally used, but it calls this function internally because the user context is actually 506 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 507 not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 508 509 .seealso: DMSNESSetContext(), SNESGetNGS(), DMSNESGetJacobian(), DMSNESGetFunction(), SNESNGSFunction 510 @*/ 511 PetscErrorCode DMSNESGetNGS(DM dm,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx) 512 { 513 PetscErrorCode ierr; 514 DMSNES sdm; 515 516 PetscFunctionBegin; 517 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 518 ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr); 519 if (f) *f = sdm->ops->computegs; 520 if (ctx) *ctx = sdm->gsctx; 521 PetscFunctionReturn(0); 522 } 523 524 #undef __FUNCT__ 525 #define __FUNCT__ "DMSNESSetJacobian" 526 /*@C 527 DMSNESSetJacobian - set SNES Jacobian evaluation function 528 529 Not Collective 530 531 Input Argument: 532 + dm - DM to be used with SNES 533 . J - Jacobian evaluation function 534 - ctx - context for residual evaluation 535 536 Level: advanced 537 538 Note: 539 SNESSetJacobian() is normally used, but it calls this function internally because the user context is actually 540 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 541 not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 542 543 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESGetJacobian(), SNESSetJacobian(), SNESJacobianFunction 544 @*/ 545 PetscErrorCode DMSNESSetJacobian(DM dm,PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx) 546 { 547 PetscErrorCode ierr; 548 DMSNES sdm; 549 550 PetscFunctionBegin; 551 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 552 if (J || ctx) { 553 ierr = DMGetDMSNESWrite(dm,&sdm);CHKERRQ(ierr); 554 } 555 if (J) sdm->ops->computejacobian = J; 556 if (ctx) sdm->jacobianctx = ctx; 557 PetscFunctionReturn(0); 558 } 559 560 #undef __FUNCT__ 561 #define __FUNCT__ "DMSNESGetJacobian" 562 /*@C 563 DMSNESGetJacobian - get SNES Jacobian evaluation function 564 565 Not Collective 566 567 Input Argument: 568 . dm - DM to be used with SNES 569 570 Output Arguments: 571 + J - Jacobian evaluation function; see SNESJacobianFunction for all calling sequence 572 - ctx - context for residual evaluation 573 574 Level: advanced 575 576 Note: 577 SNESGetJacobian() is normally used, but it calls this function internally because the user context is actually 578 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 579 not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 580 581 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), SNESJacobianFunction 582 @*/ 583 PetscErrorCode DMSNESGetJacobian(DM dm,PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx) 584 { 585 PetscErrorCode ierr; 586 DMSNES sdm; 587 588 PetscFunctionBegin; 589 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 590 ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr); 591 if (J) *J = sdm->ops->computejacobian; 592 if (ctx) *ctx = sdm->jacobianctx; 593 PetscFunctionReturn(0); 594 } 595 596 #undef __FUNCT__ 597 #define __FUNCT__ "DMSNESSetPicard" 598 /*@C 599 DMSNESSetPicard - set SNES Picard iteration matrix and RHS evaluation functions. 600 601 Not Collective 602 603 Input Argument: 604 + dm - DM to be used with SNES 605 . b - RHS evaluation function 606 . J - Picard matrix evaluation function 607 - ctx - context for residual evaluation 608 609 Level: advanced 610 611 .seealso: SNESSetPicard(), DMSNESSetFunction(), DMSNESSetJacobian() 612 @*/ 613 PetscErrorCode DMSNESSetPicard(DM dm,PetscErrorCode (*b)(SNES,Vec,Vec,void*),PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx) 614 { 615 PetscErrorCode ierr; 616 DMSNES sdm; 617 618 PetscFunctionBegin; 619 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 620 ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr); 621 if (b) sdm->ops->computepfunction = b; 622 if (J) sdm->ops->computepjacobian = J; 623 if (ctx) sdm->pctx = ctx; 624 PetscFunctionReturn(0); 625 } 626 627 #undef __FUNCT__ 628 #define __FUNCT__ "DMSNESGetPicard" 629 /*@C 630 DMSNESGetPicard - get SNES Picard iteration evaluation functions 631 632 Not Collective 633 634 Input Argument: 635 . dm - DM to be used with SNES 636 637 Output Arguments: 638 + b - RHS evaluation function; see SNESFunction for details 639 . J - RHS evaluation function; see SNESJacobianFunction for detailsa 640 - ctx - context for residual evaluation 641 642 Level: advanced 643 644 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian() 645 @*/ 646 PetscErrorCode DMSNESGetPicard(DM dm,PetscErrorCode (**b)(SNES,Vec,Vec,void*),PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx) 647 { 648 PetscErrorCode ierr; 649 DMSNES sdm; 650 651 PetscFunctionBegin; 652 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 653 ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr); 654 if (b) *b = sdm->ops->computepfunction; 655 if (J) *J = sdm->ops->computepjacobian; 656 if (ctx) *ctx = sdm->pctx; 657 PetscFunctionReturn(0); 658 } 659