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