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