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__ "DMSNESSetGS" 367 /*@C 368 DMSNESSetGS - set SNES Gauss-Seidel relaxation function 369 370 Not Collective 371 372 Input Argument: 373 + dm - DM to be used with SNES 374 . func - relaxation function, see SNESSetGS() for calling sequence 375 - ctx - context for residual evaluation 376 377 Level: advanced 378 379 Note: 380 SNESSetGS() is normally used, but it calls this function internally because the user context is actually 381 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 382 not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 383 384 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian(), DMSNESSetFunction() 385 @*/ 386 PetscErrorCode DMSNESSetGS(DM dm,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx) 387 { 388 PetscErrorCode ierr; 389 SNESDM sdm; 390 391 PetscFunctionBegin; 392 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 393 ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr); 394 if (func) sdm->computegs = func; 395 if (ctx) sdm->gsctx = ctx; 396 PetscFunctionReturn(0); 397 } 398 399 #undef __FUNCT__ 400 #define __FUNCT__ "DMSNESGetGS" 401 /*@C 402 DMSNESGetGS - get SNES Gauss-Seidel relaxation function 403 404 Not Collective 405 406 Input Argument: 407 . dm - DM to be used with SNES 408 409 Output Arguments: 410 + func - relaxation function, see SNESSetGS() for calling sequence 411 - ctx - context for residual evaluation 412 413 Level: advanced 414 415 Note: 416 SNESGetGS() is normally used, but it calls this function internally because the user context is actually 417 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 418 not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 419 420 .seealso: DMSNESSetContext(), SNESGetGS(), DMSNESGetJacobian(), DMSNESGetFunction() 421 @*/ 422 PetscErrorCode DMSNESGetGS(DM dm,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx) 423 { 424 PetscErrorCode ierr; 425 SNESDM sdm; 426 427 PetscFunctionBegin; 428 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 429 ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 430 if (func) *func = sdm->computegs; 431 if (ctx) *ctx = sdm->gsctx; 432 PetscFunctionReturn(0); 433 } 434 435 #undef __FUNCT__ 436 #define __FUNCT__ "DMSNESSetJacobian" 437 /*@C 438 DMSNESSetJacobian - set SNES Jacobian evaluation function 439 440 Not Collective 441 442 Input Argument: 443 + dm - DM to be used with SNES 444 . func - Jacobian evaluation function, see SNESSetJacobian() for calling sequence 445 - ctx - context for residual evaluation 446 447 Level: advanced 448 449 Note: 450 SNESSetJacobian() is normally used, but it calls this function internally because the user context is actually 451 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 452 not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 453 454 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESGetJacobian(), SNESSetJacobian() 455 @*/ 456 PetscErrorCode DMSNESSetJacobian(DM dm,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 457 { 458 PetscErrorCode ierr; 459 SNESDM sdm; 460 461 PetscFunctionBegin; 462 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 463 ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr); 464 if (func) sdm->computejacobian = func; 465 if (ctx) sdm->jacobianctx = ctx; 466 PetscFunctionReturn(0); 467 } 468 469 #undef __FUNCT__ 470 #define __FUNCT__ "DMSNESGetJacobian" 471 /*@C 472 DMSNESGetJacobian - get SNES Jacobian evaluation function 473 474 Not Collective 475 476 Input Argument: 477 . dm - DM to be used with SNES 478 479 Output Arguments: 480 + func - Jacobian evaluation function, see SNESSetJacobian() for calling sequence 481 - ctx - context for residual evaluation 482 483 Level: advanced 484 485 Note: 486 SNESGetJacobian() is normally used, but it calls this function internally because the user context is actually 487 associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 488 not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 489 490 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian() 491 @*/ 492 PetscErrorCode DMSNESGetJacobian(DM dm,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx) 493 { 494 PetscErrorCode ierr; 495 SNESDM sdm; 496 497 PetscFunctionBegin; 498 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 499 ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 500 if (func) *func = sdm->computejacobian; 501 if (ctx) *ctx = sdm->jacobianctx; 502 PetscFunctionReturn(0); 503 } 504 505 #undef __FUNCT__ 506 #define __FUNCT__ "DMSNESSetPicard" 507 /*@C 508 DMSNESSetPicard - set SNES Picard iteration matrix and RHS evaluation functions. 509 510 Not Collective 511 512 Input Argument: 513 + dm - DM to be used with SNES 514 . func - RHS evaluation function, see SNESSetFunction() for calling sequence 515 . pjac - Picard matrix evaluation function, see SNESSetJacobian() for calling sequence 516 - ctx - context for residual evaluation 517 518 Level: advanced 519 520 .seealso: SNESSetPicard(), DMSNESSetFunction(), DMSNESSetJacobian() 521 @*/ 522 PetscErrorCode DMSNESSetPicard(DM dm,PetscErrorCode (*pfunc)(SNES,Vec,Vec,void*),PetscErrorCode (*pjac)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 523 { 524 PetscErrorCode ierr; 525 SNESDM sdm; 526 527 PetscFunctionBegin; 528 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 529 ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 530 if (pfunc) sdm->computepfunction = pfunc; 531 if (pjac) sdm->computepjacobian = pjac; 532 if (ctx) sdm->pctx = ctx; 533 PetscFunctionReturn(0); 534 } 535 536 537 #undef __FUNCT__ 538 #define __FUNCT__ "DMSNESGetPicard" 539 /*@C 540 DMSNESGetPicard - get SNES Picard iteration evaluation functions 541 542 Not Collective 543 544 Input Argument: 545 . dm - DM to be used with SNES 546 547 Output Arguments: 548 + pfunc - Jacobian evaluation function, see SNESSetJacobian() for calling sequence 549 . pjac - RHS evaluation function, see SNESSetFunction() for calling sequence 550 - ctx - context for residual evaluation 551 552 Level: advanced 553 554 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian() 555 @*/ 556 PetscErrorCode DMSNESGetPicard(DM dm,PetscErrorCode (**pfunc)(SNES,Vec,Vec,void*),PetscErrorCode (**pjac)(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 (pfunc) *pfunc = sdm->computepfunction; 565 if (pjac) *pjac = sdm->computepjacobian; 566 if (ctx) *ctx = sdm->pctx; 567 PetscFunctionReturn(0); 568 } 569 570 571 #undef __FUNCT__ 572 #define __FUNCT__ "SNESDefaultComputeFunction_DMLegacy" 573 static PetscErrorCode SNESDefaultComputeFunction_DMLegacy(SNES snes,Vec X,Vec F,void *ctx) 574 { 575 PetscErrorCode ierr; 576 DM dm; 577 578 PetscFunctionBegin; 579 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 580 ierr = DMComputeFunction(dm,X,F);CHKERRQ(ierr); 581 PetscFunctionReturn(0); 582 } 583 584 #undef __FUNCT__ 585 #define __FUNCT__ "SNESDefaultComputeJacobian_DMLegacy" 586 static PetscErrorCode SNESDefaultComputeJacobian_DMLegacy(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx) 587 { 588 PetscErrorCode ierr; 589 DM dm; 590 591 PetscFunctionBegin; 592 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 593 ierr = DMComputeJacobian(dm,X,*A,*B,mstr);CHKERRQ(ierr); 594 PetscFunctionReturn(0); 595 } 596 597 #undef __FUNCT__ 598 #define __FUNCT__ "DMSNESSetUpLegacy" 599 /* Sets up calling of legacy DM routines */ 600 PetscErrorCode DMSNESSetUpLegacy(DM dm) 601 { 602 PetscErrorCode ierr; 603 SNESDM sdm; 604 605 PetscFunctionBegin; 606 ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 607 if (!sdm->computefunction) {ierr = DMSNESSetFunction(dm,SNESDefaultComputeFunction_DMLegacy,PETSC_NULL);CHKERRQ(ierr);} 608 ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 609 if (!sdm->computejacobian) { 610 if (dm->ops->functionj) { 611 ierr = DMSNESSetJacobian(dm,SNESDefaultComputeJacobian_DMLegacy,PETSC_NULL);CHKERRQ(ierr); 612 } else { 613 ierr = DMSNESSetJacobian(dm,SNESDefaultComputeJacobianColor,PETSC_NULL);CHKERRQ(ierr); 614 } 615 } 616 PetscFunctionReturn(0); 617 } 618 619 /* block functions */ 620 621 #undef __FUNCT__ 622 #define __FUNCT__ "DMSNESSetBlockFunction" 623 /*@C 624 DMSNESSetBlockFunction - set SNES residual evaluation function 625 626 Not Collective 627 628 Input Arguments: 629 + dm - DM to be used with SNES 630 . func - residual evaluation function, see SNESSetFunction() for calling sequence 631 - ctx - context for residual evaluation 632 633 Level: developer 634 635 Note: 636 Mostly for use in DM implementations and transferred to a block function rather than being called from here. 637 638 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian() 639 @*/ 640 PetscErrorCode DMSNESSetBlockFunction(DM dm,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx) 641 { 642 PetscErrorCode ierr; 643 SNESDM sdm; 644 645 PetscFunctionBegin; 646 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 647 ierr = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr); 648 if (func) sdm->computeblockfunction = func; 649 if (ctx) sdm->blockfunctionctx = ctx; 650 PetscFunctionReturn(0); 651 } 652 653 #undef __FUNCT__ 654 #define __FUNCT__ "DMSNESGetBlockFunction" 655 /*@C 656 DMSNESGetBlockFunction - get SNES residual evaluation function 657 658 Not Collective 659 660 Input Argument: 661 . dm - DM to be used with SNES 662 663 Output Arguments: 664 + func - residual evaluation function, see SNESSetFunction() for calling sequence 665 - ctx - context for residual evaluation 666 667 Level: developer 668 669 .seealso: DMSNESSetContext(), DMSNESSetFunction(), SNESSetFunction() 670 @*/ 671 PetscErrorCode DMSNESGetBlockFunction(DM dm,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx) 672 { 673 PetscErrorCode ierr; 674 SNESDM sdm; 675 676 PetscFunctionBegin; 677 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 678 ierr = DMSNESGetContext(dm,&sdm);CHKERRQ(ierr); 679 if (func) *func = sdm->computeblockfunction; 680 if (ctx) *ctx = sdm->blockfunctionctx; 681 PetscFunctionReturn(0); 682 } 683 684 685 #undef __FUNCT__ 686 #define __FUNCT__ "DMSNESSetBlockJacobian" 687 /*@C 688 DMSNESSetJacobian - set SNES Jacobian evaluation function 689 690 Not Collective 691 692 Input Argument: 693 + dm - DM to be used with SNES 694 . func - Jacobian evaluation function, see SNESSetJacobian() for calling sequence 695 - ctx - context for residual evaluation 696 697 Level: advanced 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(), DMSNESGetJacobian(), SNESSetJacobian() 703 @*/ 704 PetscErrorCode DMSNESSetBlockJacobian(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 = DMSNESGetContextWrite(dm,&sdm);CHKERRQ(ierr); 712 if (func) sdm->computeblockjacobian = func; 713 if (ctx) sdm->blockjacobianctx = ctx; 714 PetscFunctionReturn(0); 715 } 716 717 #undef __FUNCT__ 718 #define __FUNCT__ "DMSNESGetBlockJacobian" 719 /*@C 720 DMSNESGetBlockJacobian - get SNES Jacobian evaluation function 721 722 Not Collective 723 724 Input Argument: 725 . dm - DM to be used with SNES 726 727 Output Arguments: 728 + func - Jacobian evaluation function, see SNESSetJacobian() for calling sequence 729 - ctx - context for residual evaluation 730 731 Level: advanced 732 733 .seealso: DMSNESSetContext(), SNESSetFunction(), DMSNESSetJacobian() 734 @*/ 735 PetscErrorCode DMSNESGetBlockJacobian(DM dm,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,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->computeblockjacobian; 744 if (ctx) *ctx = sdm->blockjacobianctx; 745 PetscFunctionReturn(0); 746 } 747