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