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