1 #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 2 #include <petscdmshell.h> 3 #include <petscdraw.h> 4 #include <petscds.h> 5 #include <petscdmadaptor.h> 6 #include <petscconvest.h> 7 8 PetscBool SNESRegisterAllCalled = PETSC_FALSE; 9 PetscFunctionList SNESList = NULL; 10 11 /* Logging support */ 12 PetscClassId SNES_CLASSID, DMSNES_CLASSID; 13 PetscLogEvent SNES_Solve, SNES_Setup, SNES_FunctionEval, SNES_JacobianEval, SNES_NGSEval, SNES_NGSFuncEval, SNES_NPCSolve, SNES_ObjectiveEval; 14 15 /*@ 16 SNESSetErrorIfNotConverged - Causes SNESSolve() to generate an error if the solver has not converged. 17 18 Logically Collective on SNES 19 20 Input Parameters: 21 + snes - iterative context obtained from SNESCreate() 22 - flg - PETSC_TRUE indicates you want the error generated 23 24 Options database keys: 25 . -snes_error_if_not_converged : this takes an optional truth value (0/1/no/yes/true/false) 26 27 Level: intermediate 28 29 Notes: 30 Normally PETSc continues if a linear solver fails to converge, you can call SNESGetConvergedReason() after a SNESSolve() 31 to determine if it has converged. 32 33 .seealso: SNESGetErrorIfNotConverged(), KSPGetErrorIfNotConverged(), KSPSetErrorIFNotConverged() 34 @*/ 35 PetscErrorCode SNESSetErrorIfNotConverged(SNES snes,PetscBool flg) 36 { 37 PetscFunctionBegin; 38 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 39 PetscValidLogicalCollectiveBool(snes,flg,2); 40 snes->errorifnotconverged = flg; 41 PetscFunctionReturn(0); 42 } 43 44 /*@ 45 SNESGetErrorIfNotConverged - Will SNESSolve() generate an error if the solver does not converge? 46 47 Not Collective 48 49 Input Parameter: 50 . snes - iterative context obtained from SNESCreate() 51 52 Output Parameter: 53 . flag - PETSC_TRUE if it will generate an error, else PETSC_FALSE 54 55 Level: intermediate 56 57 .seealso: SNESSetErrorIfNotConverged(), KSPGetErrorIfNotConverged(), KSPSetErrorIFNotConverged() 58 @*/ 59 PetscErrorCode SNESGetErrorIfNotConverged(SNES snes,PetscBool *flag) 60 { 61 PetscFunctionBegin; 62 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 63 PetscValidPointer(flag,2); 64 *flag = snes->errorifnotconverged; 65 PetscFunctionReturn(0); 66 } 67 68 /*@ 69 SNESSetAlwaysComputesFinalResidual - does the SNES always compute the residual at the final solution? 70 71 Logically Collective on SNES 72 73 Input Parameters: 74 + snes - the shell SNES 75 - flg - is the residual computed? 76 77 Level: advanced 78 79 .seealso: SNESGetAlwaysComputesFinalResidual() 80 @*/ 81 PetscErrorCode SNESSetAlwaysComputesFinalResidual(SNES snes, PetscBool flg) 82 { 83 PetscFunctionBegin; 84 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 85 snes->alwayscomputesfinalresidual = flg; 86 PetscFunctionReturn(0); 87 } 88 89 /*@ 90 SNESGetAlwaysComputesFinalResidual - does the SNES always compute the residual at the final solution? 91 92 Logically Collective on SNES 93 94 Input Parameter: 95 . snes - the shell SNES 96 97 Output Parameter: 98 . flg - is the residual computed? 99 100 Level: advanced 101 102 .seealso: SNESSetAlwaysComputesFinalResidual() 103 @*/ 104 PetscErrorCode SNESGetAlwaysComputesFinalResidual(SNES snes, PetscBool *flg) 105 { 106 PetscFunctionBegin; 107 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 108 *flg = snes->alwayscomputesfinalresidual; 109 PetscFunctionReturn(0); 110 } 111 112 /*@ 113 SNESSetFunctionDomainError - tells SNES that the input vector to your SNESFunction is not 114 in the functions domain. For example, negative pressure. 115 116 Logically Collective on SNES 117 118 Input Parameters: 119 . snes - the SNES context 120 121 Level: advanced 122 123 .seealso: SNESCreate(), SNESSetFunction(), SNESFunction 124 @*/ 125 PetscErrorCode SNESSetFunctionDomainError(SNES snes) 126 { 127 PetscFunctionBegin; 128 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 129 if (snes->errorifnotconverged) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"User code indicates input vector is not in the function domain"); 130 snes->domainerror = PETSC_TRUE; 131 PetscFunctionReturn(0); 132 } 133 134 /*@ 135 SNESSetJacobianDomainError - tells SNES that computeJacobian does not make sense any more. For example there is a negative element transformation. 136 137 Logically Collective on SNES 138 139 Input Parameters: 140 . snes - the SNES context 141 142 Level: advanced 143 144 .seealso: SNESCreate(), SNESSetFunction(), SNESFunction(), SNESSetFunctionDomainError() 145 @*/ 146 PetscErrorCode SNESSetJacobianDomainError(SNES snes) 147 { 148 PetscFunctionBegin; 149 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 150 if (snes->errorifnotconverged) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"User code indicates computeJacobian does not make sense"); 151 snes->jacobiandomainerror = PETSC_TRUE; 152 PetscFunctionReturn(0); 153 } 154 155 /*@ 156 SNESSetCheckJacobianDomainError - if or not to check jacobian domain error after each Jacobian evaluation. By default, we check Jacobian domain error 157 in the debug mode, and do not check it in the optimized mode. 158 159 Logically Collective on SNES 160 161 Input Parameters: 162 . snes - the SNES context 163 . flg - indicates if or not to check jacobian domain error after each Jacobian evaluation 164 165 Level: advanced 166 167 .seealso: SNESCreate(), SNESSetFunction(), SNESFunction(), SNESSetFunctionDomainError(), SNESGetCheckJacobianDomainError() 168 @*/ 169 PetscErrorCode SNESSetCheckJacobianDomainError(SNES snes, PetscBool flg) 170 { 171 PetscFunctionBegin; 172 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 173 snes->checkjacdomainerror = flg; 174 PetscFunctionReturn(0); 175 } 176 177 /*@ 178 SNESGetCheckJacobianDomainError - Get an indicator whether or not we are checking Jacobian domain errors after each Jacobian evaluation. 179 180 Logically Collective on SNES 181 182 Input Parameters: 183 . snes - the SNES context 184 185 Output Parameters: 186 . flg - PETSC_FALSE indicates that we don't check jacobian domain errors after each Jacobian evaluation 187 188 Level: advanced 189 190 .seealso: SNESCreate(), SNESSetFunction(), SNESFunction(), SNESSetFunctionDomainError(), SNESSetCheckJacobianDomainError() 191 @*/ 192 PetscErrorCode SNESGetCheckJacobianDomainError(SNES snes, PetscBool *flg) 193 { 194 PetscFunctionBegin; 195 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 196 PetscValidPointer(flg, 2); 197 *flg = snes->checkjacdomainerror; 198 PetscFunctionReturn(0); 199 } 200 201 /*@ 202 SNESGetFunctionDomainError - Gets the status of the domain error after a call to SNESComputeFunction; 203 204 Logically Collective on SNES 205 206 Input Parameters: 207 . snes - the SNES context 208 209 Output Parameters: 210 . domainerror - Set to PETSC_TRUE if there's a domain error; PETSC_FALSE otherwise. 211 212 Level: advanced 213 214 .seealso: SNESSetFunctionDomainError(), SNESComputeFunction() 215 @*/ 216 PetscErrorCode SNESGetFunctionDomainError(SNES snes, PetscBool *domainerror) 217 { 218 PetscFunctionBegin; 219 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 220 PetscValidPointer(domainerror, 2); 221 *domainerror = snes->domainerror; 222 PetscFunctionReturn(0); 223 } 224 225 /*@ 226 SNESGetJacobianDomainError - Gets the status of the Jacobian domain error after a call to SNESComputeJacobian; 227 228 Logically Collective on SNES 229 230 Input Parameters: 231 . snes - the SNES context 232 233 Output Parameters: 234 . domainerror - Set to PETSC_TRUE if there's a jacobian domain error; PETSC_FALSE otherwise. 235 236 Level: advanced 237 238 .seealso: SNESSetFunctionDomainError(), SNESComputeFunction(),SNESGetFunctionDomainError() 239 @*/ 240 PetscErrorCode SNESGetJacobianDomainError(SNES snes, PetscBool *domainerror) 241 { 242 PetscFunctionBegin; 243 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 244 PetscValidPointer(domainerror, 2); 245 *domainerror = snes->jacobiandomainerror; 246 PetscFunctionReturn(0); 247 } 248 249 /*@C 250 SNESLoad - Loads a SNES that has been stored in binary with SNESView(). 251 252 Collective on PetscViewer 253 254 Input Parameters: 255 + newdm - the newly loaded SNES, this needs to have been created with SNESCreate() or 256 some related function before a call to SNESLoad(). 257 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() 258 259 Level: intermediate 260 261 Notes: 262 The type is determined by the data in the file, any type set into the SNES before this call is ignored. 263 264 Notes for advanced users: 265 Most users should not need to know the details of the binary storage 266 format, since SNESLoad() and TSView() completely hide these details. 267 But for anyone who's interested, the standard binary matrix storage 268 format is 269 .vb 270 has not yet been determined 271 .ve 272 273 .seealso: PetscViewerBinaryOpen(), SNESView(), MatLoad(), VecLoad() 274 @*/ 275 PetscErrorCode SNESLoad(SNES snes, PetscViewer viewer) 276 { 277 PetscErrorCode ierr; 278 PetscBool isbinary; 279 PetscInt classid; 280 char type[256]; 281 KSP ksp; 282 DM dm; 283 DMSNES dmsnes; 284 285 PetscFunctionBegin; 286 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 287 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 288 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 289 if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()"); 290 291 ierr = PetscViewerBinaryRead(viewer,&classid,1,NULL,PETSC_INT);CHKERRQ(ierr); 292 if (classid != SNES_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONG,"Not SNES next in file"); 293 ierr = PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);CHKERRQ(ierr); 294 ierr = SNESSetType(snes, type);CHKERRQ(ierr); 295 if (snes->ops->load) { 296 ierr = (*snes->ops->load)(snes,viewer);CHKERRQ(ierr); 297 } 298 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 299 ierr = DMGetDMSNES(dm,&dmsnes);CHKERRQ(ierr); 300 ierr = DMSNESLoad(dmsnes,viewer);CHKERRQ(ierr); 301 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 302 ierr = KSPLoad(ksp,viewer);CHKERRQ(ierr); 303 PetscFunctionReturn(0); 304 } 305 306 #include <petscdraw.h> 307 #if defined(PETSC_HAVE_SAWS) 308 #include <petscviewersaws.h> 309 #endif 310 311 /*@C 312 SNESView - Prints the SNES data structure. 313 314 Collective on SNES 315 316 Input Parameters: 317 + SNES - the SNES context 318 - viewer - visualization context 319 320 Options Database Key: 321 . -snes_view - Calls SNESView() at end of SNESSolve() 322 323 Notes: 324 The available visualization contexts include 325 + PETSC_VIEWER_STDOUT_SELF - standard output (default) 326 - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 327 output where only the first processor opens 328 the file. All other processors send their 329 data to the first processor to print. 330 331 The user can open an alternative visualization context with 332 PetscViewerASCIIOpen() - output to a specified file. 333 334 Level: beginner 335 336 .seealso: PetscViewerASCIIOpen() 337 @*/ 338 PetscErrorCode SNESView(SNES snes,PetscViewer viewer) 339 { 340 SNESKSPEW *kctx; 341 PetscErrorCode ierr; 342 KSP ksp; 343 SNESLineSearch linesearch; 344 PetscBool iascii,isstring,isbinary,isdraw; 345 DMSNES dmsnes; 346 #if defined(PETSC_HAVE_SAWS) 347 PetscBool issaws; 348 #endif 349 350 PetscFunctionBegin; 351 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 352 if (!viewer) { 353 ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&viewer);CHKERRQ(ierr); 354 } 355 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 356 PetscCheckSameComm(snes,1,viewer,2); 357 358 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 359 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 360 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 361 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 362 #if defined(PETSC_HAVE_SAWS) 363 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSAWS,&issaws);CHKERRQ(ierr); 364 #endif 365 if (iascii) { 366 SNESNormSchedule normschedule; 367 DM dm; 368 PetscErrorCode (*cJ)(SNES,Vec,Mat,Mat,void*); 369 void *ctx; 370 371 ierr = PetscObjectPrintClassNamePrefixType((PetscObject)snes,viewer);CHKERRQ(ierr); 372 if (!snes->setupcalled) { 373 ierr = PetscViewerASCIIPrintf(viewer," SNES has not been set up so information may be incomplete\n");CHKERRQ(ierr); 374 } 375 if (snes->ops->view) { 376 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 377 ierr = (*snes->ops->view)(snes,viewer);CHKERRQ(ierr); 378 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 379 } 380 ierr = PetscViewerASCIIPrintf(viewer," maximum iterations=%D, maximum function evaluations=%D\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr); 381 ierr = PetscViewerASCIIPrintf(viewer," tolerances: relative=%g, absolute=%g, solution=%g\n",(double)snes->rtol,(double)snes->abstol,(double)snes->stol);CHKERRQ(ierr); 382 if (snes->usesksp) { 383 ierr = PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",snes->linear_its);CHKERRQ(ierr); 384 } 385 ierr = PetscViewerASCIIPrintf(viewer," total number of function evaluations=%D\n",snes->nfuncs);CHKERRQ(ierr); 386 ierr = SNESGetNormSchedule(snes, &normschedule);CHKERRQ(ierr); 387 if (normschedule > 0) {ierr = PetscViewerASCIIPrintf(viewer," norm schedule %s\n",SNESNormSchedules[normschedule]);CHKERRQ(ierr);} 388 if (snes->gridsequence) { 389 ierr = PetscViewerASCIIPrintf(viewer," total number of grid sequence refinements=%D\n",snes->gridsequence);CHKERRQ(ierr); 390 } 391 if (snes->ksp_ewconv) { 392 kctx = (SNESKSPEW*)snes->kspconvctx; 393 if (kctx) { 394 ierr = PetscViewerASCIIPrintf(viewer," Eisenstat-Walker computation of KSP relative tolerance (version %D)\n",kctx->version);CHKERRQ(ierr); 395 ierr = PetscViewerASCIIPrintf(viewer," rtol_0=%g, rtol_max=%g, threshold=%g\n",(double)kctx->rtol_0,(double)kctx->rtol_max,(double)kctx->threshold);CHKERRQ(ierr); 396 ierr = PetscViewerASCIIPrintf(viewer," gamma=%g, alpha=%g, alpha2=%g\n",(double)kctx->gamma,(double)kctx->alpha,(double)kctx->alpha2);CHKERRQ(ierr); 397 } 398 } 399 if (snes->lagpreconditioner == -1) { 400 ierr = PetscViewerASCIIPrintf(viewer," Preconditioned is never rebuilt\n");CHKERRQ(ierr); 401 } else if (snes->lagpreconditioner > 1) { 402 ierr = PetscViewerASCIIPrintf(viewer," Preconditioned is rebuilt every %D new Jacobians\n",snes->lagpreconditioner);CHKERRQ(ierr); 403 } 404 if (snes->lagjacobian == -1) { 405 ierr = PetscViewerASCIIPrintf(viewer," Jacobian is never rebuilt\n");CHKERRQ(ierr); 406 } else if (snes->lagjacobian > 1) { 407 ierr = PetscViewerASCIIPrintf(viewer," Jacobian is rebuilt every %D SNES iterations\n",snes->lagjacobian);CHKERRQ(ierr); 408 } 409 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 410 ierr = DMSNESGetJacobian(dm,&cJ,&ctx);CHKERRQ(ierr); 411 if (cJ == SNESComputeJacobianDefault) { 412 ierr = PetscViewerASCIIPrintf(viewer," Jacobian is built using finite differences one column at a time\n");CHKERRQ(ierr); 413 } else if (cJ == SNESComputeJacobianDefaultColor) { 414 ierr = PetscViewerASCIIPrintf(viewer," Jacobian is built using finite differences with coloring\n");CHKERRQ(ierr); 415 } 416 } else if (isstring) { 417 const char *type; 418 ierr = SNESGetType(snes,&type);CHKERRQ(ierr); 419 ierr = PetscViewerStringSPrintf(viewer," SNESType: %-7.7s",type);CHKERRQ(ierr); 420 if (snes->ops->view) {ierr = (*snes->ops->view)(snes,viewer);CHKERRQ(ierr);} 421 } else if (isbinary) { 422 PetscInt classid = SNES_FILE_CLASSID; 423 MPI_Comm comm; 424 PetscMPIInt rank; 425 char type[256]; 426 427 ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 428 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 429 if (!rank) { 430 ierr = PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr); 431 ierr = PetscStrncpy(type,((PetscObject)snes)->type_name,sizeof(type));CHKERRQ(ierr); 432 ierr = PetscViewerBinaryWrite(viewer,type,sizeof(type),PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr); 433 } 434 if (snes->ops->view) { 435 ierr = (*snes->ops->view)(snes,viewer);CHKERRQ(ierr); 436 } 437 } else if (isdraw) { 438 PetscDraw draw; 439 char str[36]; 440 PetscReal x,y,bottom,h; 441 442 ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 443 ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 444 ierr = PetscStrncpy(str,"SNES: ",sizeof(str));CHKERRQ(ierr); 445 ierr = PetscStrlcat(str,((PetscObject)snes)->type_name,sizeof(str));CHKERRQ(ierr); 446 ierr = PetscDrawStringBoxed(draw,x,y,PETSC_DRAW_BLUE,PETSC_DRAW_BLACK,str,NULL,&h);CHKERRQ(ierr); 447 bottom = y - h; 448 ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr); 449 if (snes->ops->view) { 450 ierr = (*snes->ops->view)(snes,viewer);CHKERRQ(ierr); 451 } 452 #if defined(PETSC_HAVE_SAWS) 453 } else if (issaws) { 454 PetscMPIInt rank; 455 const char *name; 456 457 ierr = PetscObjectGetName((PetscObject)snes,&name);CHKERRQ(ierr); 458 ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 459 if (!((PetscObject)snes)->amsmem && !rank) { 460 char dir[1024]; 461 462 ierr = PetscObjectViewSAWs((PetscObject)snes,viewer);CHKERRQ(ierr); 463 ierr = PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/its",name);CHKERRQ(ierr); 464 PetscStackCallSAWs(SAWs_Register,(dir,&snes->iter,1,SAWs_READ,SAWs_INT)); 465 if (!snes->conv_hist) { 466 ierr = SNESSetConvergenceHistory(snes,NULL,NULL,PETSC_DECIDE,PETSC_TRUE);CHKERRQ(ierr); 467 } 468 ierr = PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/conv_hist",name);CHKERRQ(ierr); 469 PetscStackCallSAWs(SAWs_Register,(dir,snes->conv_hist,10,SAWs_READ,SAWs_DOUBLE)); 470 } 471 #endif 472 } 473 if (snes->linesearch) { 474 ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 475 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 476 ierr = SNESLineSearchView(linesearch, viewer);CHKERRQ(ierr); 477 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 478 } 479 if (snes->npc && snes->usesnpc) { 480 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 481 ierr = SNESView(snes->npc, viewer);CHKERRQ(ierr); 482 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 483 } 484 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 485 ierr = DMGetDMSNES(snes->dm,&dmsnes);CHKERRQ(ierr); 486 ierr = DMSNESView(dmsnes, viewer);CHKERRQ(ierr); 487 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 488 if (snes->usesksp) { 489 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 490 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 491 ierr = KSPView(ksp,viewer);CHKERRQ(ierr); 492 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 493 } 494 if (isdraw) { 495 PetscDraw draw; 496 ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 497 ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 498 } 499 PetscFunctionReturn(0); 500 } 501 502 /* 503 We retain a list of functions that also take SNES command 504 line options. These are called at the end SNESSetFromOptions() 505 */ 506 #define MAXSETFROMOPTIONS 5 507 static PetscInt numberofsetfromoptions; 508 static PetscErrorCode (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES); 509 510 /*@C 511 SNESAddOptionsChecker - Adds an additional function to check for SNES options. 512 513 Not Collective 514 515 Input Parameter: 516 . snescheck - function that checks for options 517 518 Level: developer 519 520 .seealso: SNESSetFromOptions() 521 @*/ 522 PetscErrorCode SNESAddOptionsChecker(PetscErrorCode (*snescheck)(SNES)) 523 { 524 PetscFunctionBegin; 525 if (numberofsetfromoptions >= MAXSETFROMOPTIONS) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %D allowed", MAXSETFROMOPTIONS); 526 othersetfromoptions[numberofsetfromoptions++] = snescheck; 527 PetscFunctionReturn(0); 528 } 529 530 PETSC_INTERN PetscErrorCode SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*); 531 532 static PetscErrorCode SNESSetUpMatrixFree_Private(SNES snes, PetscBool hasOperator, PetscInt version) 533 { 534 Mat J; 535 KSP ksp; 536 PC pc; 537 PetscBool match; 538 PetscErrorCode ierr; 539 MatNullSpace nullsp; 540 541 PetscFunctionBegin; 542 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 543 544 if (!snes->vec_func && (snes->jacobian || snes->jacobian_pre)) { 545 Mat A = snes->jacobian, B = snes->jacobian_pre; 546 ierr = MatCreateVecs(A ? A : B, NULL,&snes->vec_func);CHKERRQ(ierr); 547 } 548 549 if (version == 1) { 550 ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr); 551 ierr = MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);CHKERRQ(ierr); 552 ierr = MatSetFromOptions(J);CHKERRQ(ierr); 553 } else if (version == 2) { 554 if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"SNESSetFunction() must be called first"); 555 #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128) && !defined(PETSC_USE_REAL___FP16) 556 ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_func,&J);CHKERRQ(ierr); 557 #else 558 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "matrix-free operator rutines (version 2)"); 559 #endif 560 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "matrix-free operator rutines, only version 1 and 2"); 561 562 /* attach any user provided null space that was on Amat to the newly created matrix free matrix */ 563 if (snes->jacobian) { 564 ierr = MatGetNullSpace(snes->jacobian,&nullsp);CHKERRQ(ierr); 565 if (nullsp) { 566 ierr = MatSetNullSpace(J,nullsp);CHKERRQ(ierr); 567 } 568 } 569 570 ierr = PetscInfo1(snes,"Setting default matrix-free operator routines (version %D)\n", version);CHKERRQ(ierr); 571 if (hasOperator) { 572 573 /* This version replaces the user provided Jacobian matrix with a 574 matrix-free version but still employs the user-provided preconditioner matrix. */ 575 ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 576 } else { 577 /* This version replaces both the user-provided Jacobian and the user- 578 provided preconditioner Jacobian with the default matrix free version. */ 579 if ((snes->npcside== PC_LEFT) && snes->npc) { 580 if (!snes->jacobian){ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);} 581 } else { 582 ierr = SNESSetJacobian(snes,J,J,MatMFFDComputeJacobian,0);CHKERRQ(ierr); 583 } 584 /* Force no preconditioner */ 585 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 586 ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 587 ierr = PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&match);CHKERRQ(ierr); 588 if (!match) { 589 ierr = PetscInfo(snes,"Setting default matrix-free preconditioner routines\nThat is no preconditioner is being used\n");CHKERRQ(ierr); 590 ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 591 } 592 } 593 ierr = MatDestroy(&J);CHKERRQ(ierr); 594 PetscFunctionReturn(0); 595 } 596 597 static PetscErrorCode DMRestrictHook_SNESVecSol(DM dmfine,Mat Restrict,Vec Rscale,Mat Inject,DM dmcoarse,void *ctx) 598 { 599 SNES snes = (SNES)ctx; 600 PetscErrorCode ierr; 601 Vec Xfine,Xfine_named = NULL,Xcoarse; 602 603 PetscFunctionBegin; 604 if (PetscLogPrintInfo) { 605 PetscInt finelevel,coarselevel,fineclevel,coarseclevel; 606 ierr = DMGetRefineLevel(dmfine,&finelevel);CHKERRQ(ierr); 607 ierr = DMGetCoarsenLevel(dmfine,&fineclevel);CHKERRQ(ierr); 608 ierr = DMGetRefineLevel(dmcoarse,&coarselevel);CHKERRQ(ierr); 609 ierr = DMGetCoarsenLevel(dmcoarse,&coarseclevel);CHKERRQ(ierr); 610 ierr = PetscInfo4(dmfine,"Restricting SNES solution vector from level %D-%D to level %D-%D\n",finelevel,fineclevel,coarselevel,coarseclevel);CHKERRQ(ierr); 611 } 612 if (dmfine == snes->dm) Xfine = snes->vec_sol; 613 else { 614 ierr = DMGetNamedGlobalVector(dmfine,"SNESVecSol",&Xfine_named);CHKERRQ(ierr); 615 Xfine = Xfine_named; 616 } 617 ierr = DMGetNamedGlobalVector(dmcoarse,"SNESVecSol",&Xcoarse);CHKERRQ(ierr); 618 if (Inject) { 619 ierr = MatRestrict(Inject,Xfine,Xcoarse);CHKERRQ(ierr); 620 } else { 621 ierr = MatRestrict(Restrict,Xfine,Xcoarse);CHKERRQ(ierr); 622 ierr = VecPointwiseMult(Xcoarse,Xcoarse,Rscale);CHKERRQ(ierr); 623 } 624 ierr = DMRestoreNamedGlobalVector(dmcoarse,"SNESVecSol",&Xcoarse);CHKERRQ(ierr); 625 if (Xfine_named) {ierr = DMRestoreNamedGlobalVector(dmfine,"SNESVecSol",&Xfine_named);CHKERRQ(ierr);} 626 PetscFunctionReturn(0); 627 } 628 629 static PetscErrorCode DMCoarsenHook_SNESVecSol(DM dm,DM dmc,void *ctx) 630 { 631 PetscErrorCode ierr; 632 633 PetscFunctionBegin; 634 ierr = DMCoarsenHookAdd(dmc,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,ctx);CHKERRQ(ierr); 635 PetscFunctionReturn(0); 636 } 637 638 /* This may be called to rediscretize the operator on levels of linear multigrid. The DM shuffle is so the user can 639 * safely call SNESGetDM() in their residual evaluation routine. */ 640 static PetscErrorCode KSPComputeOperators_SNES(KSP ksp,Mat A,Mat B,void *ctx) 641 { 642 SNES snes = (SNES)ctx; 643 PetscErrorCode ierr; 644 Vec X,Xnamed = NULL; 645 DM dmsave; 646 void *ctxsave; 647 PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*) = NULL; 648 649 PetscFunctionBegin; 650 dmsave = snes->dm; 651 ierr = KSPGetDM(ksp,&snes->dm);CHKERRQ(ierr); 652 if (dmsave == snes->dm) X = snes->vec_sol; /* We are on the finest level */ 653 else { /* We are on a coarser level, this vec was initialized using a DM restrict hook */ 654 ierr = DMGetNamedGlobalVector(snes->dm,"SNESVecSol",&Xnamed);CHKERRQ(ierr); 655 X = Xnamed; 656 ierr = SNESGetJacobian(snes,NULL,NULL,&jac,&ctxsave);CHKERRQ(ierr); 657 /* If the DM's don't match up, the MatFDColoring context needed for the jacobian won't match up either -- fixit. */ 658 if (jac == SNESComputeJacobianDefaultColor) { 659 ierr = SNESSetJacobian(snes,NULL,NULL,SNESComputeJacobianDefaultColor,0);CHKERRQ(ierr); 660 } 661 } 662 /* Make sure KSP DM has the Jacobian computation routine */ 663 { 664 DMSNES sdm; 665 666 ierr = DMGetDMSNES(snes->dm, &sdm);CHKERRQ(ierr); 667 if (!sdm->ops->computejacobian) { 668 ierr = DMCopyDMSNES(dmsave, snes->dm);CHKERRQ(ierr); 669 } 670 } 671 /* Compute the operators */ 672 ierr = SNESComputeJacobian(snes,X,A,B);CHKERRQ(ierr); 673 /* Put the previous context back */ 674 if (snes->dm != dmsave && jac == SNESComputeJacobianDefaultColor) { 675 ierr = SNESSetJacobian(snes,NULL,NULL,jac,ctxsave);CHKERRQ(ierr); 676 } 677 678 if (Xnamed) {ierr = DMRestoreNamedGlobalVector(snes->dm,"SNESVecSol",&Xnamed);CHKERRQ(ierr);} 679 snes->dm = dmsave; 680 PetscFunctionReturn(0); 681 } 682 683 /*@ 684 SNESSetUpMatrices - ensures that matrices are available for SNES, to be called by SNESSetUp_XXX() 685 686 Collective 687 688 Input Arguments: 689 . snes - snes to configure 690 691 Level: developer 692 693 .seealso: SNESSetUp() 694 @*/ 695 PetscErrorCode SNESSetUpMatrices(SNES snes) 696 { 697 PetscErrorCode ierr; 698 DM dm; 699 DMSNES sdm; 700 701 PetscFunctionBegin; 702 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 703 ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr); 704 if (!sdm->ops->computejacobian) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"DMSNES not properly configured"); 705 else if (!snes->jacobian && snes->mf) { 706 Mat J; 707 void *functx; 708 ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr); 709 ierr = MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);CHKERRQ(ierr); 710 ierr = MatSetFromOptions(J);CHKERRQ(ierr); 711 ierr = SNESGetFunction(snes,NULL,NULL,&functx);CHKERRQ(ierr); 712 ierr = SNESSetJacobian(snes,J,J,0,0);CHKERRQ(ierr); 713 ierr = MatDestroy(&J);CHKERRQ(ierr); 714 } else if (snes->mf_operator && !snes->jacobian_pre && !snes->jacobian) { 715 Mat J,B; 716 ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr); 717 ierr = MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);CHKERRQ(ierr); 718 ierr = MatSetFromOptions(J);CHKERRQ(ierr); 719 ierr = DMCreateMatrix(snes->dm,&B);CHKERRQ(ierr); 720 /* sdm->computejacobian was already set to reach here */ 721 ierr = SNESSetJacobian(snes,J,B,NULL,NULL);CHKERRQ(ierr); 722 ierr = MatDestroy(&J);CHKERRQ(ierr); 723 ierr = MatDestroy(&B);CHKERRQ(ierr); 724 } else if (!snes->jacobian_pre) { 725 PetscErrorCode (*nspconstr)(DM, PetscInt, MatNullSpace *); 726 PetscDS prob; 727 Mat J, B; 728 MatNullSpace nullspace = NULL; 729 PetscBool hasPrec = PETSC_FALSE; 730 PetscInt Nf; 731 732 J = snes->jacobian; 733 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 734 if (prob) {ierr = PetscDSHasJacobianPreconditioner(prob, &hasPrec);CHKERRQ(ierr);} 735 if (J) {ierr = PetscObjectReference((PetscObject) J);CHKERRQ(ierr);} 736 else if (hasPrec) {ierr = DMCreateMatrix(snes->dm, &J);CHKERRQ(ierr);} 737 ierr = DMCreateMatrix(snes->dm, &B);CHKERRQ(ierr); 738 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 739 ierr = DMGetNullSpaceConstructor(snes->dm, Nf, &nspconstr);CHKERRQ(ierr); 740 if (nspconstr) (*nspconstr)(snes->dm, -1, &nullspace); 741 ierr = MatSetNullSpace(B, nullspace);CHKERRQ(ierr); 742 ierr = MatNullSpaceDestroy(&nullspace);CHKERRQ(ierr); 743 ierr = SNESSetJacobian(snes, J ? J : B, B, NULL, NULL);CHKERRQ(ierr); 744 ierr = MatDestroy(&J);CHKERRQ(ierr); 745 ierr = MatDestroy(&B);CHKERRQ(ierr); 746 } 747 { 748 KSP ksp; 749 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 750 ierr = KSPSetComputeOperators(ksp,KSPComputeOperators_SNES,snes);CHKERRQ(ierr); 751 ierr = DMCoarsenHookAdd(snes->dm,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,snes);CHKERRQ(ierr); 752 } 753 PetscFunctionReturn(0); 754 } 755 756 /*@C 757 SNESMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user 758 759 Collective on SNES 760 761 Input Parameters: 762 + snes - SNES object you wish to monitor 763 . name - the monitor type one is seeking 764 . help - message indicating what monitoring is done 765 . manual - manual page for the monitor 766 . monitor - the monitor function 767 - monitorsetup - a function that is called once ONLY if the user selected this monitor that may set additional features of the SNES or PetscViewer objects 768 769 Level: developer 770 771 .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), 772 PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 773 PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 774 PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 775 PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 776 PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 777 PetscOptionsFList(), PetscOptionsEList() 778 @*/ 779 PetscErrorCode SNESMonitorSetFromOptions(SNES snes,const char name[],const char help[], const char manual[],PetscErrorCode (*monitor)(SNES,PetscInt,PetscReal,PetscViewerAndFormat*),PetscErrorCode (*monitorsetup)(SNES,PetscViewerAndFormat*)) 780 { 781 PetscErrorCode ierr; 782 PetscViewer viewer; 783 PetscViewerFormat format; 784 PetscBool flg; 785 786 PetscFunctionBegin; 787 ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject) snes)->options,((PetscObject)snes)->prefix,name,&viewer,&format,&flg);CHKERRQ(ierr); 788 if (flg) { 789 PetscViewerAndFormat *vf; 790 ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr); 791 ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr); 792 if (monitorsetup) { 793 ierr = (*monitorsetup)(snes,vf);CHKERRQ(ierr); 794 } 795 ierr = SNESMonitorSet(snes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);CHKERRQ(ierr); 796 } 797 PetscFunctionReturn(0); 798 } 799 800 /*@ 801 SNESSetFromOptions - Sets various SNES and KSP parameters from user options. 802 803 Collective on SNES 804 805 Input Parameter: 806 . snes - the SNES context 807 808 Options Database Keys: 809 + -snes_type <type> - newtonls, newtontr, ngmres, ncg, nrichardson, qn, vi, fas, SNESType for complete list 810 . -snes_stol - convergence tolerance in terms of the norm 811 of the change in the solution between steps 812 . -snes_atol <abstol> - absolute tolerance of residual norm 813 . -snes_rtol <rtol> - relative decrease in tolerance norm from initial 814 . -snes_divergence_tolerance <divtol> - if the residual goes above divtol*rnorm0, exit with divergence 815 . -snes_force_iteration <force> - force SNESSolve() to take at least one iteration 816 . -snes_max_it <max_it> - maximum number of iterations 817 . -snes_max_funcs <max_funcs> - maximum number of function evaluations 818 . -snes_max_fail <max_fail> - maximum number of line search failures allowed before stopping, default is none 819 . -snes_max_linear_solve_fail - number of linear solver failures before SNESSolve() stops 820 . -snes_lag_preconditioner <lag> - how often preconditioner is rebuilt (use -1 to never rebuild) 821 . -snes_lag_jacobian <lag> - how often Jacobian is rebuilt (use -1 to never rebuild) 822 . -snes_trtol <trtol> - trust region tolerance 823 . -snes_no_convergence_test - skip convergence test in nonlinear 824 solver; hence iterations will continue until max_it 825 or some other criterion is reached. Saves expense 826 of convergence test 827 . -snes_monitor [ascii][:filename][:viewer format] - prints residual norm at each iteration. if no filename given prints to stdout 828 . -snes_monitor_solution [ascii binary draw][:filename][:viewer format] - plots solution at each iteration 829 . -snes_monitor_residual [ascii binary draw][:filename][:viewer format] - plots residual (not its norm) at each iteration 830 . -snes_monitor_solution_update [ascii binary draw][:filename][:viewer format] - plots update to solution at each iteration 831 . -snes_monitor_lg_residualnorm - plots residual norm at each iteration 832 . -snes_monitor_lg_range - plots residual norm at each iteration 833 . -snes_fd - use finite differences to compute Jacobian; very slow, only for testing 834 . -snes_fd_color - use finite differences with coloring to compute Jacobian 835 . -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration 836 - -snes_converged_reason - print the reason for convergence/divergence after each solve 837 838 Options Database for Eisenstat-Walker method: 839 + -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence 840 . -snes_ksp_ew_version ver - version of Eisenstat-Walker method 841 . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0 842 . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax 843 . -snes_ksp_ew_gamma <gamma> - Sets gamma 844 . -snes_ksp_ew_alpha <alpha> - Sets alpha 845 . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2 846 - -snes_ksp_ew_threshold <threshold> - Sets threshold 847 848 Notes: 849 To see all options, run your program with the -help option or consult 850 Users-Manual: ch_snes 851 852 Level: beginner 853 854 .seealso: SNESSetOptionsPrefix(), SNESResetFromOptions() 855 @*/ 856 PetscErrorCode SNESSetFromOptions(SNES snes) 857 { 858 PetscBool flg,pcset,persist,set; 859 PetscInt i,indx,lag,grids; 860 const char *deft = SNESNEWTONLS; 861 const char *convtests[] = {"default","skip"}; 862 SNESKSPEW *kctx = NULL; 863 char type[256], monfilename[PETSC_MAX_PATH_LEN]; 864 PetscErrorCode ierr; 865 PCSide pcside; 866 const char *optionsprefix; 867 868 PetscFunctionBegin; 869 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 870 ierr = SNESRegisterAll();CHKERRQ(ierr); 871 ierr = PetscObjectOptionsBegin((PetscObject)snes);CHKERRQ(ierr); 872 if (((PetscObject)snes)->type_name) deft = ((PetscObject)snes)->type_name; 873 ierr = PetscOptionsFList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr); 874 if (flg) { 875 ierr = SNESSetType(snes,type);CHKERRQ(ierr); 876 } else if (!((PetscObject)snes)->type_name) { 877 ierr = SNESSetType(snes,deft);CHKERRQ(ierr); 878 } 879 ierr = PetscOptionsReal("-snes_stol","Stop if step length less than","SNESSetTolerances",snes->stol,&snes->stol,NULL);CHKERRQ(ierr); 880 ierr = PetscOptionsReal("-snes_atol","Stop if function norm less than","SNESSetTolerances",snes->abstol,&snes->abstol,NULL);CHKERRQ(ierr); 881 882 ierr = PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less than","SNESSetTolerances",snes->rtol,&snes->rtol,NULL);CHKERRQ(ierr); 883 ierr = PetscOptionsReal("-snes_divergence_tolerance","Stop if residual norm increases by this factor","SNESSetDivergenceTolerance",snes->divtol,&snes->divtol,NULL);CHKERRQ(ierr); 884 ierr = PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,NULL);CHKERRQ(ierr); 885 ierr = PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,NULL);CHKERRQ(ierr); 886 ierr = PetscOptionsInt("-snes_max_fail","Maximum nonlinear step failures","SNESSetMaxNonlinearStepFailures",snes->maxFailures,&snes->maxFailures,NULL);CHKERRQ(ierr); 887 ierr = PetscOptionsInt("-snes_max_linear_solve_fail","Maximum failures in linear solves allowed","SNESSetMaxLinearSolveFailures",snes->maxLinearSolveFailures,&snes->maxLinearSolveFailures,NULL);CHKERRQ(ierr); 888 ierr = PetscOptionsBool("-snes_error_if_not_converged","Generate error if solver does not converge","SNESSetErrorIfNotConverged",snes->errorifnotconverged,&snes->errorifnotconverged,NULL);CHKERRQ(ierr); 889 ierr = PetscOptionsBool("-snes_force_iteration","Force SNESSolve() to take at least one iteration","SNESSetForceIteration",snes->forceiteration,&snes->forceiteration,NULL);CHKERRQ(ierr); 890 ierr = PetscOptionsBool("-snes_check_jacobian_domain_error","Check Jacobian domain error after Jacobian evaluation","SNESCheckJacobianDomainError",snes->checkjacdomainerror,&snes->checkjacdomainerror,NULL);CHKERRQ(ierr); 891 892 ierr = PetscOptionsInt("-snes_lag_preconditioner","How often to rebuild preconditioner","SNESSetLagPreconditioner",snes->lagpreconditioner,&lag,&flg);CHKERRQ(ierr); 893 if (flg) { 894 ierr = SNESSetLagPreconditioner(snes,lag);CHKERRQ(ierr); 895 } 896 ierr = PetscOptionsBool("-snes_lag_preconditioner_persists","Preconditioner lagging through multiple solves","SNESSetLagPreconditionerPersists",snes->lagjac_persist,&persist,&flg);CHKERRQ(ierr); 897 if (flg) { 898 ierr = SNESSetLagPreconditionerPersists(snes,persist);CHKERRQ(ierr); 899 } 900 ierr = PetscOptionsInt("-snes_lag_jacobian","How often to rebuild Jacobian","SNESSetLagJacobian",snes->lagjacobian,&lag,&flg);CHKERRQ(ierr); 901 if (flg) { 902 ierr = SNESSetLagJacobian(snes,lag);CHKERRQ(ierr); 903 } 904 ierr = PetscOptionsBool("-snes_lag_jacobian_persists","Jacobian lagging through multiple solves","SNESSetLagJacobianPersists",snes->lagjac_persist,&persist,&flg);CHKERRQ(ierr); 905 if (flg) { 906 ierr = SNESSetLagJacobianPersists(snes,persist);CHKERRQ(ierr); 907 } 908 909 ierr = PetscOptionsInt("-snes_grid_sequence","Use grid sequencing to generate initial guess","SNESSetGridSequence",snes->gridsequence,&grids,&flg);CHKERRQ(ierr); 910 if (flg) { 911 ierr = SNESSetGridSequence(snes,grids);CHKERRQ(ierr); 912 } 913 914 ierr = PetscOptionsEList("-snes_convergence_test","Convergence test","SNESSetConvergenceTest",convtests,2,"default",&indx,&flg);CHKERRQ(ierr); 915 if (flg) { 916 switch (indx) { 917 case 0: ierr = SNESSetConvergenceTest(snes,SNESConvergedDefault,NULL,NULL);CHKERRQ(ierr); break; 918 case 1: ierr = SNESSetConvergenceTest(snes,SNESConvergedSkip,NULL,NULL);CHKERRQ(ierr); break; 919 } 920 } 921 922 ierr = PetscOptionsEList("-snes_norm_schedule","SNES Norm schedule","SNESSetNormSchedule",SNESNormSchedules,5,"function",&indx,&flg);CHKERRQ(ierr); 923 if (flg) { ierr = SNESSetNormSchedule(snes,(SNESNormSchedule)indx);CHKERRQ(ierr); } 924 925 ierr = PetscOptionsEList("-snes_function_type","SNES Norm schedule","SNESSetFunctionType",SNESFunctionTypes,2,"unpreconditioned",&indx,&flg);CHKERRQ(ierr); 926 if (flg) { ierr = SNESSetFunctionType(snes,(SNESFunctionType)indx);CHKERRQ(ierr); } 927 928 kctx = (SNESKSPEW*)snes->kspconvctx; 929 930 ierr = PetscOptionsBool("-snes_ksp_ew","Use Eisentat-Walker linear system convergence test","SNESKSPSetUseEW",snes->ksp_ewconv,&snes->ksp_ewconv,NULL);CHKERRQ(ierr); 931 932 ierr = PetscOptionsInt("-snes_ksp_ew_version","Version 1, 2 or 3","SNESKSPSetParametersEW",kctx->version,&kctx->version,NULL);CHKERRQ(ierr); 933 ierr = PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNESKSPSetParametersEW",kctx->rtol_0,&kctx->rtol_0,NULL);CHKERRQ(ierr); 934 ierr = PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNESKSPSetParametersEW",kctx->rtol_max,&kctx->rtol_max,NULL);CHKERRQ(ierr); 935 ierr = PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNESKSPSetParametersEW",kctx->gamma,&kctx->gamma,NULL);CHKERRQ(ierr); 936 ierr = PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNESKSPSetParametersEW",kctx->alpha,&kctx->alpha,NULL);CHKERRQ(ierr); 937 ierr = PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNESKSPSetParametersEW",kctx->alpha2,&kctx->alpha2,NULL);CHKERRQ(ierr); 938 ierr = PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNESKSPSetParametersEW",kctx->threshold,&kctx->threshold,NULL);CHKERRQ(ierr); 939 940 flg = PETSC_FALSE; 941 ierr = PetscOptionsBool("-snes_monitor_cancel","Remove all monitors","SNESMonitorCancel",flg,&flg,&set);CHKERRQ(ierr); 942 if (set && flg) {ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);} 943 944 ierr = SNESMonitorSetFromOptions(snes,"-snes_monitor","Monitor norm of function","SNESMonitorDefault",SNESMonitorDefault,NULL);CHKERRQ(ierr); 945 ierr = SNESMonitorSetFromOptions(snes,"-snes_monitor_short","Monitor norm of function with fewer digits","SNESMonitorDefaultShort",SNESMonitorDefaultShort,NULL);CHKERRQ(ierr); 946 ierr = SNESMonitorSetFromOptions(snes,"-snes_monitor_range","Monitor range of elements of function","SNESMonitorRange",SNESMonitorRange,NULL);CHKERRQ(ierr); 947 948 ierr = SNESMonitorSetFromOptions(snes,"-snes_monitor_ratio","Monitor ratios of the norm of function for consecutive steps","SNESMonitorRatio",SNESMonitorRatio,SNESMonitorRatioSetUp);CHKERRQ(ierr); 949 ierr = SNESMonitorSetFromOptions(snes,"-snes_monitor_field","Monitor norm of function (split into fields)","SNESMonitorDefaultField",SNESMonitorDefaultField,NULL);CHKERRQ(ierr); 950 ierr = SNESMonitorSetFromOptions(snes,"-snes_monitor_solution","View solution at each iteration","SNESMonitorSolution",SNESMonitorSolution,NULL);CHKERRQ(ierr); 951 ierr = SNESMonitorSetFromOptions(snes,"-snes_monitor_solution_update","View correction at each iteration","SNESMonitorSolutionUpdate",SNESMonitorSolutionUpdate,NULL);CHKERRQ(ierr); 952 ierr = SNESMonitorSetFromOptions(snes,"-snes_monitor_residual","View residual at each iteration","SNESMonitorResidual",SNESMonitorResidual,NULL);CHKERRQ(ierr); 953 ierr = SNESMonitorSetFromOptions(snes,"-snes_monitor_jacupdate_spectrum","Print the change in the spectrum of the Jacobian","SNESMonitorJacUpdateSpectrum",SNESMonitorJacUpdateSpectrum,NULL);CHKERRQ(ierr); 954 ierr = SNESMonitorSetFromOptions(snes,"-snes_monitor_fields","Monitor norm of function per field","SNESMonitorSet",SNESMonitorFields,NULL);CHKERRQ(ierr); 955 956 ierr = PetscOptionsString("-snes_monitor_python","Use Python function","SNESMonitorSet",0,monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 957 if (flg) {ierr = PetscPythonMonitorSet((PetscObject)snes,monfilename);CHKERRQ(ierr);} 958 959 960 flg = PETSC_FALSE; 961 ierr = PetscOptionsBool("-snes_monitor_lg_residualnorm","Plot function norm at each iteration","SNESMonitorLGResidualNorm",flg,&flg,NULL);CHKERRQ(ierr); 962 if (flg) { 963 PetscDrawLG ctx; 964 965 ierr = SNESMonitorLGCreate(PetscObjectComm((PetscObject)snes),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&ctx);CHKERRQ(ierr); 966 ierr = SNESMonitorSet(snes,SNESMonitorLGResidualNorm,ctx,(PetscErrorCode (*)(void**))PetscDrawLGDestroy);CHKERRQ(ierr); 967 } 968 flg = PETSC_FALSE; 969 ierr = PetscOptionsBool("-snes_monitor_lg_range","Plot function range at each iteration","SNESMonitorLGRange",flg,&flg,NULL);CHKERRQ(ierr); 970 if (flg) { 971 PetscViewer ctx; 972 973 ierr = PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&ctx);CHKERRQ(ierr); 974 ierr = SNESMonitorSet(snes,SNESMonitorLGRange,ctx,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 975 } 976 977 978 979 flg = PETSC_FALSE; 980 ierr = PetscOptionsBool("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESComputeJacobianDefault",flg,&flg,NULL);CHKERRQ(ierr); 981 if (flg) { 982 void *functx; 983 DM dm; 984 DMSNES sdm; 985 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 986 ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr); 987 sdm->jacobianctx = NULL; 988 ierr = SNESGetFunction(snes,NULL,NULL,&functx);CHKERRQ(ierr); 989 ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESComputeJacobianDefault,functx);CHKERRQ(ierr); 990 ierr = PetscInfo(snes,"Setting default finite difference Jacobian matrix\n");CHKERRQ(ierr); 991 } 992 993 flg = PETSC_FALSE; 994 ierr = PetscOptionsBool("-snes_fd_function","Use finite differences (slow) to compute function from user objective","SNESObjectiveComputeFunctionDefaultFD",flg,&flg,NULL);CHKERRQ(ierr); 995 if (flg) { 996 ierr = SNESSetFunction(snes,NULL,SNESObjectiveComputeFunctionDefaultFD,NULL);CHKERRQ(ierr); 997 } 998 999 flg = PETSC_FALSE; 1000 ierr = PetscOptionsBool("-snes_fd_color","Use finite differences with coloring to compute Jacobian","SNESComputeJacobianDefaultColor",flg,&flg,NULL);CHKERRQ(ierr); 1001 if (flg) { 1002 DM dm; 1003 DMSNES sdm; 1004 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1005 ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr); 1006 sdm->jacobianctx = NULL; 1007 ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESComputeJacobianDefaultColor,0);CHKERRQ(ierr); 1008 ierr = PetscInfo(snes,"Setting default finite difference coloring Jacobian matrix\n");CHKERRQ(ierr); 1009 } 1010 1011 flg = PETSC_FALSE; 1012 ierr = PetscOptionsBool("-snes_mf_operator","Use a Matrix-Free Jacobian with user-provided preconditioner matrix","SNESSetUseMatrixFree",PETSC_FALSE,&snes->mf_operator,&flg);CHKERRQ(ierr); 1013 if (flg && snes->mf_operator) { 1014 snes->mf_operator = PETSC_TRUE; 1015 snes->mf = PETSC_TRUE; 1016 } 1017 flg = PETSC_FALSE; 1018 ierr = PetscOptionsBool("-snes_mf","Use a Matrix-Free Jacobian with no preconditioner matrix","SNESSetUseMatrixFree",PETSC_FALSE,&snes->mf,&flg);CHKERRQ(ierr); 1019 if (!flg && snes->mf_operator) snes->mf = PETSC_TRUE; 1020 ierr = PetscOptionsInt("-snes_mf_version","Matrix-Free routines version 1 or 2","None",snes->mf_version,&snes->mf_version,0);CHKERRQ(ierr); 1021 1022 flg = PETSC_FALSE; 1023 ierr = SNESGetNPCSide(snes,&pcside);CHKERRQ(ierr); 1024 ierr = PetscOptionsEnum("-snes_npc_side","SNES nonlinear preconditioner side","SNESSetNPCSide",PCSides,(PetscEnum)pcside,(PetscEnum*)&pcside,&flg);CHKERRQ(ierr); 1025 if (flg) {ierr = SNESSetNPCSide(snes,pcside);CHKERRQ(ierr);} 1026 1027 #if defined(PETSC_HAVE_SAWS) 1028 /* 1029 Publish convergence information using SAWs 1030 */ 1031 flg = PETSC_FALSE; 1032 ierr = PetscOptionsBool("-snes_monitor_saws","Publish SNES progress using SAWs","SNESMonitorSet",flg,&flg,NULL);CHKERRQ(ierr); 1033 if (flg) { 1034 void *ctx; 1035 ierr = SNESMonitorSAWsCreate(snes,&ctx);CHKERRQ(ierr); 1036 ierr = SNESMonitorSet(snes,SNESMonitorSAWs,ctx,SNESMonitorSAWsDestroy);CHKERRQ(ierr); 1037 } 1038 #endif 1039 #if defined(PETSC_HAVE_SAWS) 1040 { 1041 PetscBool set; 1042 flg = PETSC_FALSE; 1043 ierr = PetscOptionsBool("-snes_saws_block","Block for SAWs at end of SNESSolve","PetscObjectSAWsBlock",((PetscObject)snes)->amspublishblock,&flg,&set);CHKERRQ(ierr); 1044 if (set) { 1045 ierr = PetscObjectSAWsSetBlock((PetscObject)snes,flg);CHKERRQ(ierr); 1046 } 1047 } 1048 #endif 1049 1050 for (i = 0; i < numberofsetfromoptions; i++) { 1051 ierr = (*othersetfromoptions[i])(snes);CHKERRQ(ierr); 1052 } 1053 1054 if (snes->ops->setfromoptions) { 1055 ierr = (*snes->ops->setfromoptions)(PetscOptionsObject,snes);CHKERRQ(ierr); 1056 } 1057 1058 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 1059 ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)snes);CHKERRQ(ierr); 1060 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1061 1062 if (snes->linesearch) { 1063 ierr = SNESGetLineSearch(snes, &snes->linesearch);CHKERRQ(ierr); 1064 ierr = SNESLineSearchSetFromOptions(snes->linesearch);CHKERRQ(ierr); 1065 } 1066 1067 if (snes->usesksp) { 1068 if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);} 1069 ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 1070 ierr = KSPSetFromOptions(snes->ksp);CHKERRQ(ierr); 1071 } 1072 1073 /* if someone has set the SNES NPC type, create it. */ 1074 ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 1075 ierr = PetscOptionsHasName(((PetscObject)snes)->options,optionsprefix, "-npc_snes_type", &pcset);CHKERRQ(ierr); 1076 if (pcset && (!snes->npc)) { 1077 ierr = SNESGetNPC(snes, &snes->npc);CHKERRQ(ierr); 1078 } 1079 snes->setfromoptionscalled++; 1080 PetscFunctionReturn(0); 1081 } 1082 1083 /*@ 1084 SNESResetFromOptions - Sets various SNES and KSP parameters from user options ONLY if the SNES was previously set from options 1085 1086 Collective on SNES 1087 1088 Input Parameter: 1089 . snes - the SNES context 1090 1091 Level: beginner 1092 1093 .seealso: SNESSetFromOptions(), SNESSetOptionsPrefix() 1094 @*/ 1095 PetscErrorCode SNESResetFromOptions(SNES snes) 1096 { 1097 PetscErrorCode ierr; 1098 1099 PetscFunctionBegin; 1100 if (snes->setfromoptionscalled) {ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);} 1101 PetscFunctionReturn(0); 1102 } 1103 1104 /*@C 1105 SNESSetComputeApplicationContext - Sets an optional function to compute a user-defined context for 1106 the nonlinear solvers. 1107 1108 Logically Collective on SNES 1109 1110 Input Parameters: 1111 + snes - the SNES context 1112 . compute - function to compute the context 1113 - destroy - function to destroy the context 1114 1115 Level: intermediate 1116 1117 Notes: 1118 This function is currently not available from Fortran. 1119 1120 .seealso: SNESGetApplicationContext(), SNESSetComputeApplicationContext(), SNESGetApplicationContext() 1121 @*/ 1122 PetscErrorCode SNESSetComputeApplicationContext(SNES snes,PetscErrorCode (*compute)(SNES,void**),PetscErrorCode (*destroy)(void**)) 1123 { 1124 PetscFunctionBegin; 1125 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1126 snes->ops->usercompute = compute; 1127 snes->ops->userdestroy = destroy; 1128 PetscFunctionReturn(0); 1129 } 1130 1131 /*@ 1132 SNESSetApplicationContext - Sets the optional user-defined context for 1133 the nonlinear solvers. 1134 1135 Logically Collective on SNES 1136 1137 Input Parameters: 1138 + snes - the SNES context 1139 - usrP - optional user context 1140 1141 Level: intermediate 1142 1143 Fortran Notes: 1144 To use this from Fortran you must write a Fortran interface definition for this 1145 function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument. 1146 1147 .seealso: SNESGetApplicationContext() 1148 @*/ 1149 PetscErrorCode SNESSetApplicationContext(SNES snes,void *usrP) 1150 { 1151 PetscErrorCode ierr; 1152 KSP ksp; 1153 1154 PetscFunctionBegin; 1155 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1156 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1157 ierr = KSPSetApplicationContext(ksp,usrP);CHKERRQ(ierr); 1158 snes->user = usrP; 1159 PetscFunctionReturn(0); 1160 } 1161 1162 /*@ 1163 SNESGetApplicationContext - Gets the user-defined context for the 1164 nonlinear solvers. 1165 1166 Not Collective 1167 1168 Input Parameter: 1169 . snes - SNES context 1170 1171 Output Parameter: 1172 . usrP - user context 1173 1174 Fortran Notes: 1175 To use this from Fortran you must write a Fortran interface definition for this 1176 function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument. 1177 1178 Level: intermediate 1179 1180 .seealso: SNESSetApplicationContext() 1181 @*/ 1182 PetscErrorCode SNESGetApplicationContext(SNES snes,void *usrP) 1183 { 1184 PetscFunctionBegin; 1185 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1186 *(void**)usrP = snes->user; 1187 PetscFunctionReturn(0); 1188 } 1189 1190 /*@ 1191 SNESSetUseMatrixFree - indicates that SNES should use matrix free finite difference matrix vector products internally to apply 1192 the Jacobian. 1193 1194 Collective on SNES 1195 1196 Input Parameters: 1197 + snes - SNES context 1198 . mf - use matrix-free for both the Amat and Pmat used by SNESSetJacobian(), both the Amat and Pmat set in SNESSetJacobian() will be ignored 1199 - mf_operator - use matrix-free only for the Amat used by SNESSetJacobian(), this means the user provided Pmat will continue to be used 1200 1201 Options Database: 1202 + -snes_mf - use matrix free for both the mat and pmat operator 1203 - -snes_mf_operator - use matrix free only for the mat operator 1204 1205 Level: intermediate 1206 1207 .seealso: SNESGetUseMatrixFree(), MatCreateSNESMF() 1208 @*/ 1209 PetscErrorCode SNESSetUseMatrixFree(SNES snes,PetscBool mf_operator,PetscBool mf) 1210 { 1211 PetscFunctionBegin; 1212 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1213 PetscValidLogicalCollectiveBool(snes,mf_operator,2); 1214 PetscValidLogicalCollectiveBool(snes,mf,3); 1215 if (mf && !mf_operator) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"If using mf must also use mf_operator"); 1216 snes->mf = mf; 1217 snes->mf_operator = mf_operator; 1218 PetscFunctionReturn(0); 1219 } 1220 1221 /*@ 1222 SNESGetUseMatrixFree - indicates if the SNES uses matrix free finite difference matrix vector products to apply 1223 the Jacobian. 1224 1225 Collective on SNES 1226 1227 Input Parameter: 1228 . snes - SNES context 1229 1230 Output Parameters: 1231 + mf - use matrix-free for both the Amat and Pmat used by SNESSetJacobian(), both the Amat and Pmat set in SNESSetJacobian() will be ignored 1232 - mf_operator - use matrix-free only for the Amat used by SNESSetJacobian(), this means the user provided Pmat will continue to be used 1233 1234 Options Database: 1235 + -snes_mf - use matrix free for both the mat and pmat operator 1236 - -snes_mf_operator - use matrix free only for the mat operator 1237 1238 Level: intermediate 1239 1240 .seealso: SNESSetUseMatrixFree(), MatCreateSNESMF() 1241 @*/ 1242 PetscErrorCode SNESGetUseMatrixFree(SNES snes,PetscBool *mf_operator,PetscBool *mf) 1243 { 1244 PetscFunctionBegin; 1245 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1246 if (mf) *mf = snes->mf; 1247 if (mf_operator) *mf_operator = snes->mf_operator; 1248 PetscFunctionReturn(0); 1249 } 1250 1251 /*@ 1252 SNESGetIterationNumber - Gets the number of nonlinear iterations completed 1253 at this time. 1254 1255 Not Collective 1256 1257 Input Parameter: 1258 . snes - SNES context 1259 1260 Output Parameter: 1261 . iter - iteration number 1262 1263 Notes: 1264 For example, during the computation of iteration 2 this would return 1. 1265 1266 This is useful for using lagged Jacobians (where one does not recompute the 1267 Jacobian at each SNES iteration). For example, the code 1268 .vb 1269 ierr = SNESGetIterationNumber(snes,&it); 1270 if (!(it % 2)) { 1271 [compute Jacobian here] 1272 } 1273 .ve 1274 can be used in your ComputeJacobian() function to cause the Jacobian to be 1275 recomputed every second SNES iteration. 1276 1277 After the SNES solve is complete this will return the number of nonlinear iterations used. 1278 1279 Level: intermediate 1280 1281 .seealso: SNESGetLinearSolveIterations() 1282 @*/ 1283 PetscErrorCode SNESGetIterationNumber(SNES snes,PetscInt *iter) 1284 { 1285 PetscFunctionBegin; 1286 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1287 PetscValidIntPointer(iter,2); 1288 *iter = snes->iter; 1289 PetscFunctionReturn(0); 1290 } 1291 1292 /*@ 1293 SNESSetIterationNumber - Sets the current iteration number. 1294 1295 Not Collective 1296 1297 Input Parameter: 1298 . snes - SNES context 1299 . iter - iteration number 1300 1301 Level: developer 1302 1303 .seealso: SNESGetLinearSolveIterations() 1304 @*/ 1305 PetscErrorCode SNESSetIterationNumber(SNES snes,PetscInt iter) 1306 { 1307 PetscErrorCode ierr; 1308 1309 PetscFunctionBegin; 1310 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1311 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 1312 snes->iter = iter; 1313 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 1314 PetscFunctionReturn(0); 1315 } 1316 1317 /*@ 1318 SNESGetNonlinearStepFailures - Gets the number of unsuccessful steps 1319 attempted by the nonlinear solver. 1320 1321 Not Collective 1322 1323 Input Parameter: 1324 . snes - SNES context 1325 1326 Output Parameter: 1327 . nfails - number of unsuccessful steps attempted 1328 1329 Notes: 1330 This counter is reset to zero for each successive call to SNESSolve(). 1331 1332 Level: intermediate 1333 1334 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), 1335 SNESSetMaxNonlinearStepFailures(), SNESGetMaxNonlinearStepFailures() 1336 @*/ 1337 PetscErrorCode SNESGetNonlinearStepFailures(SNES snes,PetscInt *nfails) 1338 { 1339 PetscFunctionBegin; 1340 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1341 PetscValidIntPointer(nfails,2); 1342 *nfails = snes->numFailures; 1343 PetscFunctionReturn(0); 1344 } 1345 1346 /*@ 1347 SNESSetMaxNonlinearStepFailures - Sets the maximum number of unsuccessful steps 1348 attempted by the nonlinear solver before it gives up. 1349 1350 Not Collective 1351 1352 Input Parameters: 1353 + snes - SNES context 1354 - maxFails - maximum of unsuccessful steps 1355 1356 Level: intermediate 1357 1358 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), 1359 SNESGetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures() 1360 @*/ 1361 PetscErrorCode SNESSetMaxNonlinearStepFailures(SNES snes, PetscInt maxFails) 1362 { 1363 PetscFunctionBegin; 1364 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1365 snes->maxFailures = maxFails; 1366 PetscFunctionReturn(0); 1367 } 1368 1369 /*@ 1370 SNESGetMaxNonlinearStepFailures - Gets the maximum number of unsuccessful steps 1371 attempted by the nonlinear solver before it gives up. 1372 1373 Not Collective 1374 1375 Input Parameter: 1376 . snes - SNES context 1377 1378 Output Parameter: 1379 . maxFails - maximum of unsuccessful steps 1380 1381 Level: intermediate 1382 1383 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), 1384 SNESSetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures() 1385 1386 @*/ 1387 PetscErrorCode SNESGetMaxNonlinearStepFailures(SNES snes, PetscInt *maxFails) 1388 { 1389 PetscFunctionBegin; 1390 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1391 PetscValidIntPointer(maxFails,2); 1392 *maxFails = snes->maxFailures; 1393 PetscFunctionReturn(0); 1394 } 1395 1396 /*@ 1397 SNESGetNumberFunctionEvals - Gets the number of user provided function evaluations 1398 done by SNES. 1399 1400 Not Collective 1401 1402 Input Parameter: 1403 . snes - SNES context 1404 1405 Output Parameter: 1406 . nfuncs - number of evaluations 1407 1408 Level: intermediate 1409 1410 Notes: 1411 Reset every time SNESSolve is called unless SNESSetCountersReset() is used. 1412 1413 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), SNESSetCountersReset() 1414 @*/ 1415 PetscErrorCode SNESGetNumberFunctionEvals(SNES snes, PetscInt *nfuncs) 1416 { 1417 PetscFunctionBegin; 1418 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1419 PetscValidIntPointer(nfuncs,2); 1420 *nfuncs = snes->nfuncs; 1421 PetscFunctionReturn(0); 1422 } 1423 1424 /*@ 1425 SNESGetLinearSolveFailures - Gets the number of failed (non-converged) 1426 linear solvers. 1427 1428 Not Collective 1429 1430 Input Parameter: 1431 . snes - SNES context 1432 1433 Output Parameter: 1434 . nfails - number of failed solves 1435 1436 Level: intermediate 1437 1438 Options Database Keys: 1439 . -snes_max_linear_solve_fail <num> - The number of failures before the solve is terminated 1440 1441 Notes: 1442 This counter is reset to zero for each successive call to SNESSolve(). 1443 1444 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures() 1445 @*/ 1446 PetscErrorCode SNESGetLinearSolveFailures(SNES snes,PetscInt *nfails) 1447 { 1448 PetscFunctionBegin; 1449 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1450 PetscValidIntPointer(nfails,2); 1451 *nfails = snes->numLinearSolveFailures; 1452 PetscFunctionReturn(0); 1453 } 1454 1455 /*@ 1456 SNESSetMaxLinearSolveFailures - the number of failed linear solve attempts 1457 allowed before SNES returns with a diverged reason of SNES_DIVERGED_LINEAR_SOLVE 1458 1459 Logically Collective on SNES 1460 1461 Input Parameters: 1462 + snes - SNES context 1463 - maxFails - maximum allowed linear solve failures 1464 1465 Level: intermediate 1466 1467 Options Database Keys: 1468 . -snes_max_linear_solve_fail <num> - The number of failures before the solve is terminated 1469 1470 Notes: 1471 By default this is 0; that is SNES returns on the first failed linear solve 1472 1473 .seealso: SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations() 1474 @*/ 1475 PetscErrorCode SNESSetMaxLinearSolveFailures(SNES snes, PetscInt maxFails) 1476 { 1477 PetscFunctionBegin; 1478 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1479 PetscValidLogicalCollectiveInt(snes,maxFails,2); 1480 snes->maxLinearSolveFailures = maxFails; 1481 PetscFunctionReturn(0); 1482 } 1483 1484 /*@ 1485 SNESGetMaxLinearSolveFailures - gets the maximum number of linear solve failures that 1486 are allowed before SNES terminates 1487 1488 Not Collective 1489 1490 Input Parameter: 1491 . snes - SNES context 1492 1493 Output Parameter: 1494 . maxFails - maximum of unsuccessful solves allowed 1495 1496 Level: intermediate 1497 1498 Notes: 1499 By default this is 1; that is SNES returns on the first failed linear solve 1500 1501 .seealso: SNESGetLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), 1502 @*/ 1503 PetscErrorCode SNESGetMaxLinearSolveFailures(SNES snes, PetscInt *maxFails) 1504 { 1505 PetscFunctionBegin; 1506 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1507 PetscValidIntPointer(maxFails,2); 1508 *maxFails = snes->maxLinearSolveFailures; 1509 PetscFunctionReturn(0); 1510 } 1511 1512 /*@ 1513 SNESGetLinearSolveIterations - Gets the total number of linear iterations 1514 used by the nonlinear solver. 1515 1516 Not Collective 1517 1518 Input Parameter: 1519 . snes - SNES context 1520 1521 Output Parameter: 1522 . lits - number of linear iterations 1523 1524 Notes: 1525 This counter is reset to zero for each successive call to SNESSolve() unless SNESSetCountersReset() is used. 1526 1527 If the linear solver fails inside the SNESSolve() the iterations for that call to the linear solver are not included. If you wish to count them 1528 then call KSPGetIterationNumber() after the failed solve. 1529 1530 Level: intermediate 1531 1532 .seealso: SNESGetIterationNumber(), SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESSetCountersReset() 1533 @*/ 1534 PetscErrorCode SNESGetLinearSolveIterations(SNES snes,PetscInt *lits) 1535 { 1536 PetscFunctionBegin; 1537 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1538 PetscValidIntPointer(lits,2); 1539 *lits = snes->linear_its; 1540 PetscFunctionReturn(0); 1541 } 1542 1543 /*@ 1544 SNESSetCountersReset - Sets whether or not the counters for linear iterations and function evaluations 1545 are reset every time SNESSolve() is called. 1546 1547 Logically Collective on SNES 1548 1549 Input Parameter: 1550 + snes - SNES context 1551 - reset - whether to reset the counters or not 1552 1553 Notes: 1554 This defaults to PETSC_TRUE 1555 1556 Level: developer 1557 1558 .seealso: SNESGetNumberFunctionEvals(), SNESGetLinearSolveIterations(), SNESGetNPC() 1559 @*/ 1560 PetscErrorCode SNESSetCountersReset(SNES snes,PetscBool reset) 1561 { 1562 PetscFunctionBegin; 1563 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1564 PetscValidLogicalCollectiveBool(snes,reset,2); 1565 snes->counters_reset = reset; 1566 PetscFunctionReturn(0); 1567 } 1568 1569 1570 /*@ 1571 SNESSetKSP - Sets a KSP context for the SNES object to use 1572 1573 Not Collective, but the SNES and KSP objects must live on the same MPI_Comm 1574 1575 Input Parameters: 1576 + snes - the SNES context 1577 - ksp - the KSP context 1578 1579 Notes: 1580 The SNES object already has its KSP object, you can obtain with SNESGetKSP() 1581 so this routine is rarely needed. 1582 1583 The KSP object that is already in the SNES object has its reference count 1584 decreased by one. 1585 1586 Level: developer 1587 1588 .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP() 1589 @*/ 1590 PetscErrorCode SNESSetKSP(SNES snes,KSP ksp) 1591 { 1592 PetscErrorCode ierr; 1593 1594 PetscFunctionBegin; 1595 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1596 PetscValidHeaderSpecific(ksp,KSP_CLASSID,2); 1597 PetscCheckSameComm(snes,1,ksp,2); 1598 ierr = PetscObjectReference((PetscObject)ksp);CHKERRQ(ierr); 1599 if (snes->ksp) {ierr = PetscObjectDereference((PetscObject)snes->ksp);CHKERRQ(ierr);} 1600 snes->ksp = ksp; 1601 PetscFunctionReturn(0); 1602 } 1603 1604 /* -----------------------------------------------------------*/ 1605 /*@ 1606 SNESCreate - Creates a nonlinear solver context. 1607 1608 Collective on MPI_Comm 1609 1610 Input Parameters: 1611 . comm - MPI communicator 1612 1613 Output Parameter: 1614 . outsnes - the new SNES context 1615 1616 Options Database Keys: 1617 + -snes_mf - Activates default matrix-free Jacobian-vector products, 1618 and no preconditioning matrix 1619 . -snes_mf_operator - Activates default matrix-free Jacobian-vector 1620 products, and a user-provided preconditioning matrix 1621 as set by SNESSetJacobian() 1622 - -snes_fd - Uses (slow!) finite differences to compute Jacobian 1623 1624 Level: beginner 1625 1626 Developer Notes: 1627 SNES always creates a KSP object even though many SNES methods do not use it. This is 1628 unfortunate and should be fixed at some point. The flag snes->usesksp indicates if the 1629 particular method does use KSP and regulates if the information about the KSP is printed 1630 in SNESView(). TSSetFromOptions() does call SNESSetFromOptions() which can lead to users being confused 1631 by help messages about meaningless SNES options. 1632 1633 SNES always creates the snes->kspconvctx even though it is used by only one type. This should 1634 be fixed. 1635 1636 .seealso: SNESSolve(), SNESDestroy(), SNES, SNESSetLagPreconditioner() 1637 1638 @*/ 1639 PetscErrorCode SNESCreate(MPI_Comm comm,SNES *outsnes) 1640 { 1641 PetscErrorCode ierr; 1642 SNES snes; 1643 SNESKSPEW *kctx; 1644 1645 PetscFunctionBegin; 1646 PetscValidPointer(outsnes,2); 1647 *outsnes = NULL; 1648 ierr = SNESInitializePackage();CHKERRQ(ierr); 1649 1650 ierr = PetscHeaderCreate(snes,SNES_CLASSID,"SNES","Nonlinear solver","SNES",comm,SNESDestroy,SNESView);CHKERRQ(ierr); 1651 1652 snes->ops->converged = SNESConvergedDefault; 1653 snes->usesksp = PETSC_TRUE; 1654 snes->tolerancesset = PETSC_FALSE; 1655 snes->max_its = 50; 1656 snes->max_funcs = 10000; 1657 snes->norm = 0.0; 1658 snes->xnorm = 0.0; 1659 snes->ynorm = 0.0; 1660 snes->normschedule = SNES_NORM_ALWAYS; 1661 snes->functype = SNES_FUNCTION_DEFAULT; 1662 #if defined(PETSC_USE_REAL_SINGLE) 1663 snes->rtol = 1.e-5; 1664 #else 1665 snes->rtol = 1.e-8; 1666 #endif 1667 snes->ttol = 0.0; 1668 #if defined(PETSC_USE_REAL_SINGLE) 1669 snes->abstol = 1.e-25; 1670 #else 1671 snes->abstol = 1.e-50; 1672 #endif 1673 #if defined(PETSC_USE_REAL_SINGLE) 1674 snes->stol = 1.e-5; 1675 #else 1676 snes->stol = 1.e-8; 1677 #endif 1678 #if defined(PETSC_USE_REAL_SINGLE) 1679 snes->deltatol = 1.e-6; 1680 #else 1681 snes->deltatol = 1.e-12; 1682 #endif 1683 snes->divtol = 1.e4; 1684 snes->rnorm0 = 0; 1685 snes->nfuncs = 0; 1686 snes->numFailures = 0; 1687 snes->maxFailures = 1; 1688 snes->linear_its = 0; 1689 snes->lagjacobian = 1; 1690 snes->jac_iter = 0; 1691 snes->lagjac_persist = PETSC_FALSE; 1692 snes->lagpreconditioner = 1; 1693 snes->pre_iter = 0; 1694 snes->lagpre_persist = PETSC_FALSE; 1695 snes->numbermonitors = 0; 1696 snes->data = 0; 1697 snes->setupcalled = PETSC_FALSE; 1698 snes->ksp_ewconv = PETSC_FALSE; 1699 snes->nwork = 0; 1700 snes->work = 0; 1701 snes->nvwork = 0; 1702 snes->vwork = 0; 1703 snes->conv_hist_len = 0; 1704 snes->conv_hist_max = 0; 1705 snes->conv_hist = NULL; 1706 snes->conv_hist_its = NULL; 1707 snes->conv_hist_reset = PETSC_TRUE; 1708 snes->counters_reset = PETSC_TRUE; 1709 snes->vec_func_init_set = PETSC_FALSE; 1710 snes->reason = SNES_CONVERGED_ITERATING; 1711 snes->npcside = PC_RIGHT; 1712 snes->setfromoptionscalled = 0; 1713 1714 snes->mf = PETSC_FALSE; 1715 snes->mf_operator = PETSC_FALSE; 1716 snes->mf_version = 1; 1717 1718 snes->numLinearSolveFailures = 0; 1719 snes->maxLinearSolveFailures = 1; 1720 1721 snes->vizerotolerance = 1.e-8; 1722 #if defined(PETSC_USE_DEBUG) 1723 snes->checkjacdomainerror = PETSC_TRUE; 1724 #else 1725 snes->checkjacdomainerror = PETSC_FALSE; 1726 #endif 1727 1728 /* Set this to true if the implementation of SNESSolve_XXX does compute the residual at the final solution. */ 1729 snes->alwayscomputesfinalresidual = PETSC_FALSE; 1730 1731 /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 1732 ierr = PetscNewLog(snes,&kctx);CHKERRQ(ierr); 1733 1734 snes->kspconvctx = (void*)kctx; 1735 kctx->version = 2; 1736 kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 1737 this was too large for some test cases */ 1738 kctx->rtol_last = 0.0; 1739 kctx->rtol_max = .9; 1740 kctx->gamma = 1.0; 1741 kctx->alpha = .5*(1.0 + PetscSqrtReal(5.0)); 1742 kctx->alpha2 = kctx->alpha; 1743 kctx->threshold = .1; 1744 kctx->lresid_last = 0.0; 1745 kctx->norm_last = 0.0; 1746 1747 *outsnes = snes; 1748 PetscFunctionReturn(0); 1749 } 1750 1751 /*MC 1752 SNESFunction - Functional form used to convey the nonlinear function to be solved by SNES 1753 1754 Synopsis: 1755 #include "petscsnes.h" 1756 PetscErrorCode SNESFunction(SNES snes,Vec x,Vec f,void *ctx); 1757 1758 Input Parameters: 1759 + snes - the SNES context 1760 . x - state at which to evaluate residual 1761 - ctx - optional user-defined function context, passed in with SNESSetFunction() 1762 1763 Output Parameter: 1764 . f - vector to put residual (function value) 1765 1766 Level: intermediate 1767 1768 .seealso: SNESSetFunction(), SNESGetFunction() 1769 M*/ 1770 1771 /*@C 1772 SNESSetFunction - Sets the function evaluation routine and function 1773 vector for use by the SNES routines in solving systems of nonlinear 1774 equations. 1775 1776 Logically Collective on SNES 1777 1778 Input Parameters: 1779 + snes - the SNES context 1780 . r - vector to store function value 1781 . f - function evaluation routine; see SNESFunction for calling sequence details 1782 - ctx - [optional] user-defined context for private data for the 1783 function evaluation routine (may be NULL) 1784 1785 Notes: 1786 The Newton-like methods typically solve linear systems of the form 1787 $ f'(x) x = -f(x), 1788 where f'(x) denotes the Jacobian matrix and f(x) is the function. 1789 1790 Level: beginner 1791 1792 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetPicard(), SNESFunction 1793 @*/ 1794 PetscErrorCode SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx) 1795 { 1796 PetscErrorCode ierr; 1797 DM dm; 1798 1799 PetscFunctionBegin; 1800 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1801 if (r) { 1802 PetscValidHeaderSpecific(r,VEC_CLASSID,2); 1803 PetscCheckSameComm(snes,1,r,2); 1804 ierr = PetscObjectReference((PetscObject)r);CHKERRQ(ierr); 1805 ierr = VecDestroy(&snes->vec_func);CHKERRQ(ierr); 1806 1807 snes->vec_func = r; 1808 } 1809 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1810 ierr = DMSNESSetFunction(dm,f,ctx);CHKERRQ(ierr); 1811 PetscFunctionReturn(0); 1812 } 1813 1814 1815 /*@C 1816 SNESSetInitialFunction - Sets the function vector to be used as the 1817 function norm at the initialization of the method. In some 1818 instances, the user has precomputed the function before calling 1819 SNESSolve. This function allows one to avoid a redundant call 1820 to SNESComputeFunction in that case. 1821 1822 Logically Collective on SNES 1823 1824 Input Parameters: 1825 + snes - the SNES context 1826 - f - vector to store function value 1827 1828 Notes: 1829 This should not be modified during the solution procedure. 1830 1831 This is used extensively in the SNESFAS hierarchy and in nonlinear preconditioning. 1832 1833 Level: developer 1834 1835 .seealso: SNESSetFunction(), SNESComputeFunction(), SNESSetInitialFunctionNorm() 1836 @*/ 1837 PetscErrorCode SNESSetInitialFunction(SNES snes, Vec f) 1838 { 1839 PetscErrorCode ierr; 1840 Vec vec_func; 1841 1842 PetscFunctionBegin; 1843 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1844 PetscValidHeaderSpecific(f,VEC_CLASSID,2); 1845 PetscCheckSameComm(snes,1,f,2); 1846 if (snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 1847 snes->vec_func_init_set = PETSC_FALSE; 1848 PetscFunctionReturn(0); 1849 } 1850 ierr = SNESGetFunction(snes,&vec_func,NULL,NULL);CHKERRQ(ierr); 1851 ierr = VecCopy(f, vec_func);CHKERRQ(ierr); 1852 1853 snes->vec_func_init_set = PETSC_TRUE; 1854 PetscFunctionReturn(0); 1855 } 1856 1857 /*@ 1858 SNESSetNormSchedule - Sets the SNESNormSchedule used in covergence and monitoring 1859 of the SNES method. 1860 1861 Logically Collective on SNES 1862 1863 Input Parameters: 1864 + snes - the SNES context 1865 - normschedule - the frequency of norm computation 1866 1867 Options Database Key: 1868 . -snes_norm_schedule <none, always, initialonly, finalonly, initalfinalonly> 1869 1870 Notes: 1871 Only certain SNES methods support certain SNESNormSchedules. Most require evaluation 1872 of the nonlinear function and the taking of its norm at every iteration to 1873 even ensure convergence at all. However, methods such as custom Gauss-Seidel methods 1874 (SNESNGS) and the like do not require the norm of the function to be computed, and therfore 1875 may either be monitored for convergence or not. As these are often used as nonlinear 1876 preconditioners, monitoring the norm of their error is not a useful enterprise within 1877 their solution. 1878 1879 Level: developer 1880 1881 .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule 1882 @*/ 1883 PetscErrorCode SNESSetNormSchedule(SNES snes, SNESNormSchedule normschedule) 1884 { 1885 PetscFunctionBegin; 1886 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1887 snes->normschedule = normschedule; 1888 PetscFunctionReturn(0); 1889 } 1890 1891 1892 /*@ 1893 SNESGetNormSchedule - Gets the SNESNormSchedule used in covergence and monitoring 1894 of the SNES method. 1895 1896 Logically Collective on SNES 1897 1898 Input Parameters: 1899 + snes - the SNES context 1900 - normschedule - the type of the norm used 1901 1902 Level: advanced 1903 1904 .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule 1905 @*/ 1906 PetscErrorCode SNESGetNormSchedule(SNES snes, SNESNormSchedule *normschedule) 1907 { 1908 PetscFunctionBegin; 1909 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1910 *normschedule = snes->normschedule; 1911 PetscFunctionReturn(0); 1912 } 1913 1914 1915 /*@ 1916 SNESSetFunctionNorm - Sets the last computed residual norm. 1917 1918 Logically Collective on SNES 1919 1920 Input Parameters: 1921 + snes - the SNES context 1922 1923 - normschedule - the frequency of norm computation 1924 1925 Level: developer 1926 1927 .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule 1928 @*/ 1929 PetscErrorCode SNESSetFunctionNorm(SNES snes, PetscReal norm) 1930 { 1931 PetscFunctionBegin; 1932 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1933 snes->norm = norm; 1934 PetscFunctionReturn(0); 1935 } 1936 1937 /*@ 1938 SNESGetFunctionNorm - Gets the last computed norm of the residual 1939 1940 Not Collective 1941 1942 Input Parameter: 1943 . snes - the SNES context 1944 1945 Output Parameter: 1946 . norm - the last computed residual norm 1947 1948 Level: developer 1949 1950 .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule 1951 @*/ 1952 PetscErrorCode SNESGetFunctionNorm(SNES snes, PetscReal *norm) 1953 { 1954 PetscFunctionBegin; 1955 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1956 PetscValidPointer(norm, 2); 1957 *norm = snes->norm; 1958 PetscFunctionReturn(0); 1959 } 1960 1961 /*@ 1962 SNESGetUpdateNorm - Gets the last computed norm of the Newton update 1963 1964 Not Collective 1965 1966 Input Parameter: 1967 . snes - the SNES context 1968 1969 Output Parameter: 1970 . ynorm - the last computed update norm 1971 1972 Level: developer 1973 1974 .seealso: SNESSetNormSchedule(), SNESComputeFunction(), SNESGetFunctionNorm() 1975 @*/ 1976 PetscErrorCode SNESGetUpdateNorm(SNES snes, PetscReal *ynorm) 1977 { 1978 PetscFunctionBegin; 1979 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1980 PetscValidPointer(ynorm, 2); 1981 *ynorm = snes->ynorm; 1982 PetscFunctionReturn(0); 1983 } 1984 1985 /*@ 1986 SNESGetSolutionNorm - Gets the last computed norm of the solution 1987 1988 Not Collective 1989 1990 Input Parameter: 1991 . snes - the SNES context 1992 1993 Output Parameter: 1994 . xnorm - the last computed solution norm 1995 1996 Level: developer 1997 1998 .seealso: SNESSetNormSchedule(), SNESComputeFunction(), SNESGetFunctionNorm(), SNESGetUpdateNorm() 1999 @*/ 2000 PetscErrorCode SNESGetSolutionNorm(SNES snes, PetscReal *xnorm) 2001 { 2002 PetscFunctionBegin; 2003 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2004 PetscValidPointer(xnorm, 2); 2005 *xnorm = snes->xnorm; 2006 PetscFunctionReturn(0); 2007 } 2008 2009 /*@C 2010 SNESSetFunctionType - Sets the SNESNormSchedule used in covergence and monitoring 2011 of the SNES method. 2012 2013 Logically Collective on SNES 2014 2015 Input Parameters: 2016 + snes - the SNES context 2017 - normschedule - the frequency of norm computation 2018 2019 Notes: 2020 Only certain SNES methods support certain SNESNormSchedules. Most require evaluation 2021 of the nonlinear function and the taking of its norm at every iteration to 2022 even ensure convergence at all. However, methods such as custom Gauss-Seidel methods 2023 (SNESNGS) and the like do not require the norm of the function to be computed, and therfore 2024 may either be monitored for convergence or not. As these are often used as nonlinear 2025 preconditioners, monitoring the norm of their error is not a useful enterprise within 2026 their solution. 2027 2028 Level: developer 2029 2030 .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule 2031 @*/ 2032 PetscErrorCode SNESSetFunctionType(SNES snes, SNESFunctionType type) 2033 { 2034 PetscFunctionBegin; 2035 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2036 snes->functype = type; 2037 PetscFunctionReturn(0); 2038 } 2039 2040 2041 /*@C 2042 SNESGetFunctionType - Gets the SNESNormSchedule used in covergence and monitoring 2043 of the SNES method. 2044 2045 Logically Collective on SNES 2046 2047 Input Parameters: 2048 + snes - the SNES context 2049 - normschedule - the type of the norm used 2050 2051 Level: advanced 2052 2053 .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule 2054 @*/ 2055 PetscErrorCode SNESGetFunctionType(SNES snes, SNESFunctionType *type) 2056 { 2057 PetscFunctionBegin; 2058 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2059 *type = snes->functype; 2060 PetscFunctionReturn(0); 2061 } 2062 2063 /*MC 2064 SNESNGSFunction - function used to convey a Gauss-Seidel sweep on the nonlinear function 2065 2066 Synopsis: 2067 #include <petscsnes.h> 2068 $ SNESNGSFunction(SNES snes,Vec x,Vec b,void *ctx); 2069 2070 + X - solution vector 2071 . B - RHS vector 2072 - ctx - optional user-defined Gauss-Seidel context 2073 2074 Level: intermediate 2075 2076 .seealso: SNESSetNGS(), SNESGetNGS() 2077 M*/ 2078 2079 /*@C 2080 SNESSetNGS - Sets the user nonlinear Gauss-Seidel routine for 2081 use with composed nonlinear solvers. 2082 2083 Input Parameters: 2084 + snes - the SNES context 2085 . f - function evaluation routine to apply Gauss-Seidel see SNESNGSFunction 2086 - ctx - [optional] user-defined context for private data for the 2087 smoother evaluation routine (may be NULL) 2088 2089 Notes: 2090 The NGS routines are used by the composed nonlinear solver to generate 2091 a problem appropriate update to the solution, particularly FAS. 2092 2093 Level: intermediate 2094 2095 .seealso: SNESGetFunction(), SNESComputeNGS() 2096 @*/ 2097 PetscErrorCode SNESSetNGS(SNES snes,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx) 2098 { 2099 PetscErrorCode ierr; 2100 DM dm; 2101 2102 PetscFunctionBegin; 2103 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2104 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 2105 ierr = DMSNESSetNGS(dm,f,ctx);CHKERRQ(ierr); 2106 PetscFunctionReturn(0); 2107 } 2108 2109 PetscErrorCode SNESPicardComputeFunction(SNES snes,Vec x,Vec f,void *ctx) 2110 { 2111 PetscErrorCode ierr; 2112 DM dm; 2113 DMSNES sdm; 2114 2115 PetscFunctionBegin; 2116 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 2117 ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr); 2118 if (!sdm->ops->computepfunction) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetPicard() to provide Picard function."); 2119 if (!sdm->ops->computepjacobian) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetPicard() to provide Picard Jacobian."); 2120 /* A(x)*x - b(x) */ 2121 PetscStackPush("SNES Picard user function"); 2122 ierr = (*sdm->ops->computepfunction)(snes,x,f,sdm->pctx);CHKERRQ(ierr); 2123 PetscStackPop; 2124 PetscStackPush("SNES Picard user Jacobian"); 2125 ierr = (*sdm->ops->computepjacobian)(snes,x,snes->jacobian,snes->jacobian_pre,sdm->pctx);CHKERRQ(ierr); 2126 PetscStackPop; 2127 ierr = VecScale(f,-1.0);CHKERRQ(ierr); 2128 ierr = MatMultAdd(snes->jacobian,x,f,f);CHKERRQ(ierr); 2129 PetscFunctionReturn(0); 2130 } 2131 2132 PetscErrorCode SNESPicardComputeJacobian(SNES snes,Vec x1,Mat J,Mat B,void *ctx) 2133 { 2134 PetscFunctionBegin; 2135 /* the jacobian matrix should be pre-filled in SNESPicardComputeFunction */ 2136 PetscFunctionReturn(0); 2137 } 2138 2139 /*@C 2140 SNESSetPicard - Use SNES to solve the semilinear-system A(x) x = b(x) via a Picard type iteration (Picard linearization) 2141 2142 Logically Collective on SNES 2143 2144 Input Parameters: 2145 + snes - the SNES context 2146 . r - vector to store function value 2147 . b - function evaluation routine 2148 . Amat - matrix with which A(x) x - b(x) is to be computed 2149 . Pmat - matrix from which preconditioner is computed (usually the same as Amat) 2150 . J - function to compute matrix value, see SNESJacobianFunction for details on its calling sequence 2151 - ctx - [optional] user-defined context for private data for the 2152 function evaluation routine (may be NULL) 2153 2154 Notes: 2155 We do not recomemend using this routine. It is far better to provide the nonlinear function F() and some approximation to the Jacobian and use 2156 an approximate Newton solver. This interface is provided to allow porting/testing a previous Picard based code in PETSc before converting it to approximate Newton. 2157 2158 One can call SNESSetPicard() or SNESSetFunction() (and possibly SNESSetJacobian()) but cannot call both 2159 2160 $ Solves the equation A(x) x = b(x) via the defect correction algorithm A(x^{n}) (x^{n+1} - x^{n}) = b(x^{n}) - A(x^{n})x^{n} 2161 $ Note that when an exact solver is used this corresponds to the "classic" Picard A(x^{n}) x^{n+1} = b(x^{n}) iteration. 2162 2163 Run with -snes_mf_operator to solve the system with Newton's method using A(x^{n}) to construct the preconditioner. 2164 2165 We implement the defect correction form of the Picard iteration because it converges much more generally when inexact linear solvers are used then 2166 the direct Picard iteration A(x^n) x^{n+1} = b(x^n) 2167 2168 There is some controversity over the definition of a Picard iteration for nonlinear systems but almost everyone agrees that it involves a linear solve and some 2169 believe it is the iteration A(x^{n}) x^{n+1} = b(x^{n}) hence we use the name Picard. If anyone has an authoritative reference that defines the Picard iteration 2170 different please contact us at petsc-dev@mcs.anl.gov and we'll have an entirely new argument :-). 2171 2172 Level: intermediate 2173 2174 .seealso: SNESGetFunction(), SNESSetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESGetPicard(), SNESLineSearchPreCheckPicard(), SNESJacobianFunction 2175 @*/ 2176 PetscErrorCode SNESSetPicard(SNES snes,Vec r,PetscErrorCode (*b)(SNES,Vec,Vec,void*),Mat Amat, Mat Pmat, PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx) 2177 { 2178 PetscErrorCode ierr; 2179 DM dm; 2180 2181 PetscFunctionBegin; 2182 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2183 ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr); 2184 ierr = DMSNESSetPicard(dm,b,J,ctx);CHKERRQ(ierr); 2185 ierr = SNESSetFunction(snes,r,SNESPicardComputeFunction,ctx);CHKERRQ(ierr); 2186 ierr = SNESSetJacobian(snes,Amat,Pmat,SNESPicardComputeJacobian,ctx);CHKERRQ(ierr); 2187 PetscFunctionReturn(0); 2188 } 2189 2190 /*@C 2191 SNESGetPicard - Returns the context for the Picard iteration 2192 2193 Not Collective, but Vec is parallel if SNES is parallel. Collective if Vec is requested, but has not been created yet. 2194 2195 Input Parameter: 2196 . snes - the SNES context 2197 2198 Output Parameter: 2199 + r - the function (or NULL) 2200 . f - the function (or NULL); see SNESFunction for calling sequence details 2201 . Amat - the matrix used to defined the operation A(x) x - b(x) (or NULL) 2202 . Pmat - the matrix from which the preconditioner will be constructed (or NULL) 2203 . J - the function for matrix evaluation (or NULL); see SNESJacobianFunction for calling sequence details 2204 - ctx - the function context (or NULL) 2205 2206 Level: advanced 2207 2208 .seealso: SNESSetPicard(), SNESGetFunction(), SNESGetJacobian(), SNESGetDM(), SNESFunction, SNESJacobianFunction 2209 @*/ 2210 PetscErrorCode SNESGetPicard(SNES snes,Vec *r,PetscErrorCode (**f)(SNES,Vec,Vec,void*),Mat *Amat, Mat *Pmat, PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx) 2211 { 2212 PetscErrorCode ierr; 2213 DM dm; 2214 2215 PetscFunctionBegin; 2216 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2217 ierr = SNESGetFunction(snes,r,NULL,NULL);CHKERRQ(ierr); 2218 ierr = SNESGetJacobian(snes,Amat,Pmat,NULL,NULL);CHKERRQ(ierr); 2219 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 2220 ierr = DMSNESGetPicard(dm,f,J,ctx);CHKERRQ(ierr); 2221 PetscFunctionReturn(0); 2222 } 2223 2224 /*@C 2225 SNESSetComputeInitialGuess - Sets a routine used to compute an initial guess for the problem 2226 2227 Logically Collective on SNES 2228 2229 Input Parameters: 2230 + snes - the SNES context 2231 . func - function evaluation routine 2232 - ctx - [optional] user-defined context for private data for the 2233 function evaluation routine (may be NULL) 2234 2235 Calling sequence of func: 2236 $ func (SNES snes,Vec x,void *ctx); 2237 2238 . f - function vector 2239 - ctx - optional user-defined function context 2240 2241 Level: intermediate 2242 2243 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 2244 @*/ 2245 PetscErrorCode SNESSetComputeInitialGuess(SNES snes,PetscErrorCode (*func)(SNES,Vec,void*),void *ctx) 2246 { 2247 PetscFunctionBegin; 2248 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2249 if (func) snes->ops->computeinitialguess = func; 2250 if (ctx) snes->initialguessP = ctx; 2251 PetscFunctionReturn(0); 2252 } 2253 2254 /* --------------------------------------------------------------- */ 2255 /*@C 2256 SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set 2257 it assumes a zero right hand side. 2258 2259 Logically Collective on SNES 2260 2261 Input Parameter: 2262 . snes - the SNES context 2263 2264 Output Parameter: 2265 . rhs - the right hand side vector or NULL if the right hand side vector is null 2266 2267 Level: intermediate 2268 2269 .seealso: SNESGetSolution(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction() 2270 @*/ 2271 PetscErrorCode SNESGetRhs(SNES snes,Vec *rhs) 2272 { 2273 PetscFunctionBegin; 2274 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2275 PetscValidPointer(rhs,2); 2276 *rhs = snes->vec_rhs; 2277 PetscFunctionReturn(0); 2278 } 2279 2280 /*@ 2281 SNESComputeFunction - Calls the function that has been set with SNESSetFunction(). 2282 2283 Collective on SNES 2284 2285 Input Parameters: 2286 + snes - the SNES context 2287 - x - input vector 2288 2289 Output Parameter: 2290 . y - function vector, as set by SNESSetFunction() 2291 2292 Notes: 2293 SNESComputeFunction() is typically used within nonlinear solvers 2294 implementations, so most users would not generally call this routine 2295 themselves. 2296 2297 Level: developer 2298 2299 .seealso: SNESSetFunction(), SNESGetFunction() 2300 @*/ 2301 PetscErrorCode SNESComputeFunction(SNES snes,Vec x,Vec y) 2302 { 2303 PetscErrorCode ierr; 2304 DM dm; 2305 DMSNES sdm; 2306 2307 PetscFunctionBegin; 2308 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2309 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 2310 PetscValidHeaderSpecific(y,VEC_CLASSID,3); 2311 PetscCheckSameComm(snes,1,x,2); 2312 PetscCheckSameComm(snes,1,y,3); 2313 ierr = VecValidValues(x,2,PETSC_TRUE);CHKERRQ(ierr); 2314 2315 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 2316 ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr); 2317 if (sdm->ops->computefunction) { 2318 if (sdm->ops->computefunction != SNESObjectiveComputeFunctionDefaultFD) { 2319 ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 2320 } 2321 ierr = VecLockReadPush(x);CHKERRQ(ierr); 2322 PetscStackPush("SNES user function"); 2323 ierr = (*sdm->ops->computefunction)(snes,x,y,sdm->functionctx);CHKERRQ(ierr); 2324 PetscStackPop; 2325 ierr = VecLockReadPop(x);CHKERRQ(ierr); 2326 if (sdm->ops->computefunction != SNESObjectiveComputeFunctionDefaultFD) { 2327 ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 2328 } 2329 } else if (snes->vec_rhs) { 2330 ierr = MatMult(snes->jacobian, x, y);CHKERRQ(ierr); 2331 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() or SNESSetDM() before SNESComputeFunction(), likely called from SNESSolve()."); 2332 if (snes->vec_rhs) { 2333 ierr = VecAXPY(y,-1.0,snes->vec_rhs);CHKERRQ(ierr); 2334 } 2335 snes->nfuncs++; 2336 /* 2337 domainerror might not be set on all processes; so we tag vector locally with Inf and the next inner product or norm will 2338 propagate the value to all processes 2339 */ 2340 if (snes->domainerror) { 2341 ierr = VecSetInf(y);CHKERRQ(ierr); 2342 } 2343 PetscFunctionReturn(0); 2344 } 2345 2346 /*@ 2347 SNESComputeNGS - Calls the Gauss-Seidel function that has been set with SNESSetNGS(). 2348 2349 Collective on SNES 2350 2351 Input Parameters: 2352 + snes - the SNES context 2353 . x - input vector 2354 - b - rhs vector 2355 2356 Output Parameter: 2357 . x - new solution vector 2358 2359 Notes: 2360 SNESComputeNGS() is typically used within composed nonlinear solver 2361 implementations, so most users would not generally call this routine 2362 themselves. 2363 2364 Level: developer 2365 2366 .seealso: SNESSetNGS(), SNESComputeFunction() 2367 @*/ 2368 PetscErrorCode SNESComputeNGS(SNES snes,Vec b,Vec x) 2369 { 2370 PetscErrorCode ierr; 2371 DM dm; 2372 DMSNES sdm; 2373 2374 PetscFunctionBegin; 2375 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2376 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 2377 if (b) PetscValidHeaderSpecific(b,VEC_CLASSID,3); 2378 PetscCheckSameComm(snes,1,x,2); 2379 if (b) PetscCheckSameComm(snes,1,b,3); 2380 if (b) {ierr = VecValidValues(b,2,PETSC_TRUE);CHKERRQ(ierr);} 2381 ierr = PetscLogEventBegin(SNES_NGSEval,snes,x,b,0);CHKERRQ(ierr); 2382 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 2383 ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr); 2384 if (sdm->ops->computegs) { 2385 if (b) {ierr = VecLockReadPush(b);CHKERRQ(ierr);} 2386 PetscStackPush("SNES user NGS"); 2387 ierr = (*sdm->ops->computegs)(snes,x,b,sdm->gsctx);CHKERRQ(ierr); 2388 PetscStackPop; 2389 if (b) {ierr = VecLockReadPop(b);CHKERRQ(ierr);} 2390 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetNGS() before SNESComputeNGS(), likely called from SNESSolve()."); 2391 ierr = PetscLogEventEnd(SNES_NGSEval,snes,x,b,0);CHKERRQ(ierr); 2392 PetscFunctionReturn(0); 2393 } 2394 2395 PetscErrorCode SNESTestJacobian(SNES snes) 2396 { 2397 Mat A,B,C,D,jacobian; 2398 Vec x = snes->vec_sol,f = snes->vec_func; 2399 PetscErrorCode ierr; 2400 PetscReal nrm,gnorm; 2401 PetscReal threshold = 1.e-5; 2402 PetscInt m,n,M,N; 2403 void *functx; 2404 PetscBool complete_print = PETSC_FALSE,threshold_print = PETSC_FALSE,test = PETSC_FALSE,flg; 2405 PetscViewer viewer,mviewer; 2406 MPI_Comm comm; 2407 PetscInt tabs; 2408 static PetscBool directionsprinted = PETSC_FALSE; 2409 PetscViewerFormat format; 2410 2411 PetscFunctionBegin; 2412 ierr = PetscObjectOptionsBegin((PetscObject)snes);CHKERRQ(ierr); 2413 ierr = PetscOptionsName("-snes_test_jacobian","Compare hand-coded and finite difference Jacobians","None",&test);CHKERRQ(ierr); 2414 ierr = PetscOptionsReal("-snes_test_jacobian", "Threshold for element difference between hand-coded and finite difference being meaningful", "None", threshold, &threshold,NULL);CHKERRQ(ierr); 2415 ierr = PetscOptionsViewer("-snes_test_jacobian_view","View difference between hand-coded and finite difference Jacobians element entries","None",&mviewer,&format,&complete_print);CHKERRQ(ierr); 2416 if (!complete_print) { 2417 ierr = PetscOptionsViewer("-snes_test_jacobian_display","Display difference between hand-coded and finite difference Jacobians","None",&mviewer,&format,&complete_print);CHKERRQ(ierr); 2418 } 2419 /* for compatibility with PETSc 3.9 and older. */ 2420 ierr = PetscOptionsReal("-snes_test_jacobian_display_threshold", "Display difference between hand-coded and finite difference Jacobians which exceed input threshold", "None", threshold, &threshold, &threshold_print);CHKERRQ(ierr); 2421 ierr = PetscOptionsEnd();CHKERRQ(ierr); 2422 if (!test) PetscFunctionReturn(0); 2423 2424 ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 2425 ierr = PetscViewerASCIIGetStdout(comm,&viewer);CHKERRQ(ierr); 2426 ierr = PetscViewerASCIIGetTab(viewer, &tabs);CHKERRQ(ierr); 2427 ierr = PetscViewerASCIISetTab(viewer, ((PetscObject)snes)->tablevel);CHKERRQ(ierr); 2428 ierr = PetscViewerASCIIPrintf(viewer," ---------- Testing Jacobian -------------\n");CHKERRQ(ierr); 2429 if (!complete_print && !directionsprinted) { 2430 ierr = PetscViewerASCIIPrintf(viewer," Run with -snes_test_jacobian_view and optionally -snes_test_jacobian <threshold> to show difference\n");CHKERRQ(ierr); 2431 ierr = PetscViewerASCIIPrintf(viewer," of hand-coded and finite difference Jacobian entries greater than <threshold>.\n");CHKERRQ(ierr); 2432 } 2433 if (!directionsprinted) { 2434 ierr = PetscViewerASCIIPrintf(viewer," Testing hand-coded Jacobian, if (for double precision runs) ||J - Jfd||_F/||J||_F is\n");CHKERRQ(ierr); 2435 ierr = PetscViewerASCIIPrintf(viewer," O(1.e-8), the hand-coded Jacobian is probably correct.\n");CHKERRQ(ierr); 2436 directionsprinted = PETSC_TRUE; 2437 } 2438 if (complete_print) { 2439 ierr = PetscViewerPushFormat(mviewer,format);CHKERRQ(ierr); 2440 } 2441 2442 /* evaluate the function at this point because SNESComputeJacobianDefault() assumes that the function has been evaluated and put into snes->vec_func */ 2443 ierr = SNESComputeFunction(snes,x,f);CHKERRQ(ierr); 2444 2445 ierr = PetscObjectTypeCompare((PetscObject)snes->jacobian,MATMFFD,&flg);CHKERRQ(ierr); 2446 if (!flg) jacobian = snes->jacobian; 2447 else jacobian = snes->jacobian_pre; 2448 2449 while (jacobian) { 2450 ierr = PetscObjectTypeCompareAny((PetscObject)jacobian,&flg,MATSEQAIJ,MATMPIAIJ,MATSEQDENSE,MATMPIDENSE,MATSEQBAIJ,MATMPIBAIJ,MATSEQSBAIJ,MATMPIBAIJ,"");CHKERRQ(ierr); 2451 if (flg) { 2452 A = jacobian; 2453 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 2454 } else { 2455 ierr = MatComputeOperator(jacobian,MATAIJ,&A);CHKERRQ(ierr); 2456 } 2457 2458 ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 2459 ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 2460 ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 2461 ierr = MatSetSizes(B,m,n,M,N);CHKERRQ(ierr); 2462 ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 2463 ierr = MatSetUp(B);CHKERRQ(ierr); 2464 ierr = MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 2465 2466 ierr = SNESGetFunction(snes,NULL,NULL,&functx);CHKERRQ(ierr); 2467 ierr = SNESComputeJacobianDefault(snes,x,B,B,functx);CHKERRQ(ierr); 2468 2469 ierr = MatDuplicate(B,MAT_COPY_VALUES,&D);CHKERRQ(ierr); 2470 ierr = MatAYPX(D,-1.0,A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 2471 ierr = MatNorm(D,NORM_FROBENIUS,&nrm);CHKERRQ(ierr); 2472 ierr = MatNorm(A,NORM_FROBENIUS,&gnorm);CHKERRQ(ierr); 2473 ierr = MatDestroy(&D);CHKERRQ(ierr); 2474 if (!gnorm) gnorm = 1; /* just in case */ 2475 ierr = PetscViewerASCIIPrintf(viewer," ||J - Jfd||_F/||J||_F = %g, ||J - Jfd||_F = %g\n",(double)(nrm/gnorm),(double)nrm);CHKERRQ(ierr); 2476 2477 if (complete_print) { 2478 ierr = PetscViewerASCIIPrintf(viewer," Hand-coded Jacobian ----------\n");CHKERRQ(ierr); 2479 ierr = MatView(jacobian,mviewer);CHKERRQ(ierr); 2480 ierr = PetscViewerASCIIPrintf(viewer," Finite difference Jacobian ----------\n");CHKERRQ(ierr); 2481 ierr = MatView(B,mviewer);CHKERRQ(ierr); 2482 } 2483 2484 if (threshold_print || complete_print) { 2485 PetscInt Istart, Iend, *ccols, bncols, cncols, j, row; 2486 PetscScalar *cvals; 2487 const PetscInt *bcols; 2488 const PetscScalar *bvals; 2489 2490 ierr = MatAYPX(B,-1.0,A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 2491 ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr); 2492 ierr = MatSetSizes(C,m,n,M,N);CHKERRQ(ierr); 2493 ierr = MatSetType(C,((PetscObject)A)->type_name);CHKERRQ(ierr); 2494 ierr = MatSetUp(C);CHKERRQ(ierr); 2495 ierr = MatSetOption(C,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 2496 ierr = MatGetOwnershipRange(B,&Istart,&Iend);CHKERRQ(ierr); 2497 2498 for (row = Istart; row < Iend; row++) { 2499 ierr = MatGetRow(B,row,&bncols,&bcols,&bvals);CHKERRQ(ierr); 2500 ierr = PetscMalloc2(bncols,&ccols,bncols,&cvals);CHKERRQ(ierr); 2501 for (j = 0, cncols = 0; j < bncols; j++) { 2502 if (PetscAbsScalar(bvals[j]) > threshold) { 2503 ccols[cncols] = bcols[j]; 2504 cvals[cncols] = bvals[j]; 2505 cncols += 1; 2506 } 2507 } 2508 if (cncols) { 2509 ierr = MatSetValues(C,1,&row,cncols,ccols,cvals,INSERT_VALUES);CHKERRQ(ierr); 2510 } 2511 ierr = MatRestoreRow(B,row,&bncols,&bcols,&bvals);CHKERRQ(ierr); 2512 ierr = PetscFree2(ccols,cvals);CHKERRQ(ierr); 2513 } 2514 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2515 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2516 ierr = PetscViewerASCIIPrintf(viewer," Hand-coded minus finite-difference Jacobian with tolerance %g ----------\n",(double)threshold);CHKERRQ(ierr); 2517 ierr = MatView(C,complete_print ? mviewer : viewer);CHKERRQ(ierr); 2518 ierr = MatDestroy(&C);CHKERRQ(ierr); 2519 } 2520 ierr = MatDestroy(&A);CHKERRQ(ierr); 2521 ierr = MatDestroy(&B);CHKERRQ(ierr); 2522 2523 if (jacobian != snes->jacobian_pre) { 2524 jacobian = snes->jacobian_pre; 2525 ierr = PetscViewerASCIIPrintf(viewer," ---------- Testing Jacobian for preconditioner -------------\n");CHKERRQ(ierr); 2526 } 2527 else jacobian = NULL; 2528 } 2529 if (complete_print) { 2530 ierr = PetscViewerPopFormat(mviewer);CHKERRQ(ierr); 2531 } 2532 if (mviewer) { ierr = PetscViewerDestroy(&mviewer);CHKERRQ(ierr); } 2533 ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr); 2534 PetscFunctionReturn(0); 2535 } 2536 2537 /*@ 2538 SNESComputeJacobian - Computes the Jacobian matrix that has been set with SNESSetJacobian(). 2539 2540 Collective on SNES and Mat 2541 2542 Input Parameters: 2543 + snes - the SNES context 2544 - x - input vector 2545 2546 Output Parameters: 2547 + A - Jacobian matrix 2548 - B - optional preconditioning matrix 2549 2550 Options Database Keys: 2551 + -snes_lag_preconditioner <lag> 2552 . -snes_lag_jacobian <lag> 2553 . -snes_test_jacobian - compare the user provided Jacobian with one compute via finite differences to check for errors 2554 . -snes_test_jacobian_display - display the user provided Jacobian, the finite difference Jacobian and the difference between them to help users detect the location of errors in the user provided Jacobian 2555 . -snes_test_jacobian_display_threshold <numerical value> - display entries in the difference between the user provided Jacobian and finite difference Jacobian that are greater than a certain value to help users detect errors 2556 . -snes_compare_explicit - Compare the computed Jacobian to the finite difference Jacobian and output the differences 2557 . -snes_compare_explicit_draw - Compare the computed Jacobian to the finite difference Jacobian and draw the result 2558 . -snes_compare_explicit_contour - Compare the computed Jacobian to the finite difference Jacobian and draw a contour plot with the result 2559 . -snes_compare_operator - Make the comparison options above use the operator instead of the preconditioning matrix 2560 . -snes_compare_coloring - Compute the finite difference Jacobian using coloring and display norms of difference 2561 . -snes_compare_coloring_display - Compute the finite differece Jacobian using coloring and display verbose differences 2562 . -snes_compare_coloring_threshold - Display only those matrix entries that differ by more than a given threshold 2563 . -snes_compare_coloring_threshold_atol - Absolute tolerance for difference in matrix entries to be displayed by -snes_compare_coloring_threshold 2564 . -snes_compare_coloring_threshold_rtol - Relative tolerance for difference in matrix entries to be displayed by -snes_compare_coloring_threshold 2565 . -snes_compare_coloring_draw - Compute the finite differece Jacobian using coloring and draw differences 2566 - -snes_compare_coloring_draw_contour - Compute the finite differece Jacobian using coloring and show contours of matrices and differences 2567 2568 2569 Notes: 2570 Most users should not need to explicitly call this routine, as it 2571 is used internally within the nonlinear solvers. 2572 2573 Developer Notes: 2574 This has duplicative ways of checking the accuracy of the user provided Jacobian (see the options above). This is for historical reasons, the routine SNESTestJacobian() use to used 2575 for with the SNESType of test that has been removed. 2576 2577 Level: developer 2578 2579 .seealso: SNESSetJacobian(), KSPSetOperators(), MatStructure, SNESSetLagPreconditioner(), SNESSetLagJacobian() 2580 @*/ 2581 PetscErrorCode SNESComputeJacobian(SNES snes,Vec X,Mat A,Mat B) 2582 { 2583 PetscErrorCode ierr; 2584 PetscBool flag; 2585 DM dm; 2586 DMSNES sdm; 2587 KSP ksp; 2588 2589 PetscFunctionBegin; 2590 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2591 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 2592 PetscCheckSameComm(snes,1,X,2); 2593 ierr = VecValidValues(X,2,PETSC_TRUE);CHKERRQ(ierr); 2594 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 2595 ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr); 2596 2597 if (!sdm->ops->computejacobian) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_USER,"Must call SNESSetJacobian(), DMSNESSetJacobian(), DMDASNESSetJacobianLocal(), etc"); 2598 2599 /* make sure that MatAssemblyBegin/End() is called on A matrix if it is matrix free */ 2600 2601 if (snes->lagjacobian == -2) { 2602 snes->lagjacobian = -1; 2603 2604 ierr = PetscInfo(snes,"Recomputing Jacobian/preconditioner because lag is -2 (means compute Jacobian, but then never again) \n");CHKERRQ(ierr); 2605 } else if (snes->lagjacobian == -1) { 2606 ierr = PetscInfo(snes,"Reusing Jacobian/preconditioner because lag is -1\n");CHKERRQ(ierr); 2607 ierr = PetscObjectTypeCompare((PetscObject)A,MATMFFD,&flag);CHKERRQ(ierr); 2608 if (flag) { 2609 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2610 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2611 } 2612 PetscFunctionReturn(0); 2613 } else if (snes->lagjacobian > 1 && (snes->iter + snes->jac_iter) % snes->lagjacobian) { 2614 ierr = PetscInfo2(snes,"Reusing Jacobian/preconditioner because lag is %D and SNES iteration is %D\n",snes->lagjacobian,snes->iter);CHKERRQ(ierr); 2615 ierr = PetscObjectTypeCompare((PetscObject)A,MATMFFD,&flag);CHKERRQ(ierr); 2616 if (flag) { 2617 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2618 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2619 } 2620 PetscFunctionReturn(0); 2621 } 2622 if (snes->npc && snes->npcside== PC_LEFT) { 2623 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2624 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2625 PetscFunctionReturn(0); 2626 } 2627 2628 ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,A,B);CHKERRQ(ierr); 2629 ierr = VecLockReadPush(X);CHKERRQ(ierr); 2630 PetscStackPush("SNES user Jacobian function"); 2631 ierr = (*sdm->ops->computejacobian)(snes,X,A,B,sdm->jacobianctx);CHKERRQ(ierr); 2632 PetscStackPop; 2633 ierr = VecLockReadPop(X);CHKERRQ(ierr); 2634 ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,A,B);CHKERRQ(ierr); 2635 2636 /* the next line ensures that snes->ksp exists */ 2637 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 2638 if (snes->lagpreconditioner == -2) { 2639 ierr = PetscInfo(snes,"Rebuilding preconditioner exactly once since lag is -2\n");CHKERRQ(ierr); 2640 ierr = KSPSetReusePreconditioner(snes->ksp,PETSC_FALSE);CHKERRQ(ierr); 2641 snes->lagpreconditioner = -1; 2642 } else if (snes->lagpreconditioner == -1) { 2643 ierr = PetscInfo(snes,"Reusing preconditioner because lag is -1\n");CHKERRQ(ierr); 2644 ierr = KSPSetReusePreconditioner(snes->ksp,PETSC_TRUE);CHKERRQ(ierr); 2645 } else if (snes->lagpreconditioner > 1 && (snes->iter + snes->pre_iter) % snes->lagpreconditioner) { 2646 ierr = PetscInfo2(snes,"Reusing preconditioner because lag is %D and SNES iteration is %D\n",snes->lagpreconditioner,snes->iter);CHKERRQ(ierr); 2647 ierr = KSPSetReusePreconditioner(snes->ksp,PETSC_TRUE);CHKERRQ(ierr); 2648 } else { 2649 ierr = PetscInfo(snes,"Rebuilding preconditioner\n");CHKERRQ(ierr); 2650 ierr = KSPSetReusePreconditioner(snes->ksp,PETSC_FALSE);CHKERRQ(ierr); 2651 } 2652 2653 ierr = SNESTestJacobian(snes);CHKERRQ(ierr); 2654 /* make sure user returned a correct Jacobian and preconditioner */ 2655 /* PetscValidHeaderSpecific(A,MAT_CLASSID,3); 2656 PetscValidHeaderSpecific(B,MAT_CLASSID,4); */ 2657 { 2658 PetscBool flag = PETSC_FALSE,flag_draw = PETSC_FALSE,flag_contour = PETSC_FALSE,flag_operator = PETSC_FALSE; 2659 ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject) snes)->options,((PetscObject)snes)->prefix,"-snes_compare_explicit",NULL,NULL,&flag);CHKERRQ(ierr); 2660 ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject) snes)->options,((PetscObject)snes)->prefix,"-snes_compare_explicit_draw",NULL,NULL,&flag_draw);CHKERRQ(ierr); 2661 ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject) snes)->options,((PetscObject)snes)->prefix,"-snes_compare_explicit_draw_contour",NULL,NULL,&flag_contour);CHKERRQ(ierr); 2662 ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject) snes)->options,((PetscObject)snes)->prefix,"-snes_compare_operator",NULL,NULL,&flag_operator);CHKERRQ(ierr); 2663 if (flag || flag_draw || flag_contour) { 2664 Mat Bexp_mine = NULL,Bexp,FDexp; 2665 PetscViewer vdraw,vstdout; 2666 PetscBool flg; 2667 if (flag_operator) { 2668 ierr = MatComputeOperator(A,MATAIJ,&Bexp_mine);CHKERRQ(ierr); 2669 Bexp = Bexp_mine; 2670 } else { 2671 /* See if the preconditioning matrix can be viewed and added directly */ 2672 ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQAIJ,MATMPIAIJ,MATSEQDENSE,MATMPIDENSE,MATSEQBAIJ,MATMPIBAIJ,MATSEQSBAIJ,MATMPIBAIJ,"");CHKERRQ(ierr); 2673 if (flg) Bexp = B; 2674 else { 2675 /* If the "preconditioning" matrix is itself MATSHELL or some other type without direct support */ 2676 ierr = MatComputeOperator(B,MATAIJ,&Bexp_mine);CHKERRQ(ierr); 2677 Bexp = Bexp_mine; 2678 } 2679 } 2680 ierr = MatConvert(Bexp,MATSAME,MAT_INITIAL_MATRIX,&FDexp);CHKERRQ(ierr); 2681 ierr = SNESComputeJacobianDefault(snes,X,FDexp,FDexp,NULL);CHKERRQ(ierr); 2682 ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&vstdout);CHKERRQ(ierr); 2683 if (flag_draw || flag_contour) { 2684 ierr = PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),0,"Explicit Jacobians",PETSC_DECIDE,PETSC_DECIDE,300,300,&vdraw);CHKERRQ(ierr); 2685 if (flag_contour) {ierr = PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);CHKERRQ(ierr);} 2686 } else vdraw = NULL; 2687 ierr = PetscViewerASCIIPrintf(vstdout,"Explicit %s\n",flag_operator ? "Jacobian" : "preconditioning Jacobian");CHKERRQ(ierr); 2688 if (flag) {ierr = MatView(Bexp,vstdout);CHKERRQ(ierr);} 2689 if (vdraw) {ierr = MatView(Bexp,vdraw);CHKERRQ(ierr);} 2690 ierr = PetscViewerASCIIPrintf(vstdout,"Finite difference Jacobian\n");CHKERRQ(ierr); 2691 if (flag) {ierr = MatView(FDexp,vstdout);CHKERRQ(ierr);} 2692 if (vdraw) {ierr = MatView(FDexp,vdraw);CHKERRQ(ierr);} 2693 ierr = MatAYPX(FDexp,-1.0,Bexp,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 2694 ierr = PetscViewerASCIIPrintf(vstdout,"User-provided matrix minus finite difference Jacobian\n");CHKERRQ(ierr); 2695 if (flag) {ierr = MatView(FDexp,vstdout);CHKERRQ(ierr);} 2696 if (vdraw) { /* Always use contour for the difference */ 2697 ierr = PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);CHKERRQ(ierr); 2698 ierr = MatView(FDexp,vdraw);CHKERRQ(ierr); 2699 ierr = PetscViewerPopFormat(vdraw);CHKERRQ(ierr); 2700 } 2701 if (flag_contour) {ierr = PetscViewerPopFormat(vdraw);CHKERRQ(ierr);} 2702 ierr = PetscViewerDestroy(&vdraw);CHKERRQ(ierr); 2703 ierr = MatDestroy(&Bexp_mine);CHKERRQ(ierr); 2704 ierr = MatDestroy(&FDexp);CHKERRQ(ierr); 2705 } 2706 } 2707 { 2708 PetscBool flag = PETSC_FALSE,flag_display = PETSC_FALSE,flag_draw = PETSC_FALSE,flag_contour = PETSC_FALSE,flag_threshold = PETSC_FALSE; 2709 PetscReal threshold_atol = PETSC_SQRT_MACHINE_EPSILON,threshold_rtol = 10*PETSC_SQRT_MACHINE_EPSILON; 2710 ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring",NULL,NULL,&flag);CHKERRQ(ierr); 2711 ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring_display",NULL,NULL,&flag_display);CHKERRQ(ierr); 2712 ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring_draw",NULL,NULL,&flag_draw);CHKERRQ(ierr); 2713 ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring_draw_contour",NULL,NULL,&flag_contour);CHKERRQ(ierr); 2714 ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold",NULL,NULL,&flag_threshold);CHKERRQ(ierr); 2715 if (flag_threshold) { 2716 ierr = PetscOptionsGetReal(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold_rtol",&threshold_rtol,NULL);CHKERRQ(ierr); 2717 ierr = PetscOptionsGetReal(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold_atol",&threshold_atol,NULL);CHKERRQ(ierr); 2718 } 2719 if (flag || flag_display || flag_draw || flag_contour || flag_threshold) { 2720 Mat Bfd; 2721 PetscViewer vdraw,vstdout; 2722 MatColoring coloring; 2723 ISColoring iscoloring; 2724 MatFDColoring matfdcoloring; 2725 PetscErrorCode (*func)(SNES,Vec,Vec,void*); 2726 void *funcctx; 2727 PetscReal norm1,norm2,normmax; 2728 2729 ierr = MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&Bfd);CHKERRQ(ierr); 2730 ierr = MatColoringCreate(Bfd,&coloring);CHKERRQ(ierr); 2731 ierr = MatColoringSetType(coloring,MATCOLORINGSL);CHKERRQ(ierr); 2732 ierr = MatColoringSetFromOptions(coloring);CHKERRQ(ierr); 2733 ierr = MatColoringApply(coloring,&iscoloring);CHKERRQ(ierr); 2734 ierr = MatColoringDestroy(&coloring);CHKERRQ(ierr); 2735 ierr = MatFDColoringCreate(Bfd,iscoloring,&matfdcoloring);CHKERRQ(ierr); 2736 ierr = MatFDColoringSetFromOptions(matfdcoloring);CHKERRQ(ierr); 2737 ierr = MatFDColoringSetUp(Bfd,iscoloring,matfdcoloring);CHKERRQ(ierr); 2738 ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 2739 2740 /* This method of getting the function is currently unreliable since it doesn't work for DM local functions. */ 2741 ierr = SNESGetFunction(snes,NULL,&func,&funcctx);CHKERRQ(ierr); 2742 ierr = MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))func,funcctx);CHKERRQ(ierr); 2743 ierr = PetscObjectSetOptionsPrefix((PetscObject)matfdcoloring,((PetscObject)snes)->prefix);CHKERRQ(ierr); 2744 ierr = PetscObjectAppendOptionsPrefix((PetscObject)matfdcoloring,"coloring_");CHKERRQ(ierr); 2745 ierr = MatFDColoringSetFromOptions(matfdcoloring);CHKERRQ(ierr); 2746 ierr = MatFDColoringApply(Bfd,matfdcoloring,X,snes);CHKERRQ(ierr); 2747 ierr = MatFDColoringDestroy(&matfdcoloring);CHKERRQ(ierr); 2748 2749 ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&vstdout);CHKERRQ(ierr); 2750 if (flag_draw || flag_contour) { 2751 ierr = PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),0,"Colored Jacobians",PETSC_DECIDE,PETSC_DECIDE,300,300,&vdraw);CHKERRQ(ierr); 2752 if (flag_contour) {ierr = PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);CHKERRQ(ierr);} 2753 } else vdraw = NULL; 2754 ierr = PetscViewerASCIIPrintf(vstdout,"Explicit preconditioning Jacobian\n");CHKERRQ(ierr); 2755 if (flag_display) {ierr = MatView(B,vstdout);CHKERRQ(ierr);} 2756 if (vdraw) {ierr = MatView(B,vdraw);CHKERRQ(ierr);} 2757 ierr = PetscViewerASCIIPrintf(vstdout,"Colored Finite difference Jacobian\n");CHKERRQ(ierr); 2758 if (flag_display) {ierr = MatView(Bfd,vstdout);CHKERRQ(ierr);} 2759 if (vdraw) {ierr = MatView(Bfd,vdraw);CHKERRQ(ierr);} 2760 ierr = MatAYPX(Bfd,-1.0,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 2761 ierr = MatNorm(Bfd,NORM_1,&norm1);CHKERRQ(ierr); 2762 ierr = MatNorm(Bfd,NORM_FROBENIUS,&norm2);CHKERRQ(ierr); 2763 ierr = MatNorm(Bfd,NORM_MAX,&normmax);CHKERRQ(ierr); 2764 ierr = PetscViewerASCIIPrintf(vstdout,"User-provided matrix minus finite difference Jacobian, norm1=%g normFrob=%g normmax=%g\n",(double)norm1,(double)norm2,(double)normmax);CHKERRQ(ierr); 2765 if (flag_display) {ierr = MatView(Bfd,vstdout);CHKERRQ(ierr);} 2766 if (vdraw) { /* Always use contour for the difference */ 2767 ierr = PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);CHKERRQ(ierr); 2768 ierr = MatView(Bfd,vdraw);CHKERRQ(ierr); 2769 ierr = PetscViewerPopFormat(vdraw);CHKERRQ(ierr); 2770 } 2771 if (flag_contour) {ierr = PetscViewerPopFormat(vdraw);CHKERRQ(ierr);} 2772 2773 if (flag_threshold) { 2774 PetscInt bs,rstart,rend,i; 2775 ierr = MatGetBlockSize(B,&bs);CHKERRQ(ierr); 2776 ierr = MatGetOwnershipRange(B,&rstart,&rend);CHKERRQ(ierr); 2777 for (i=rstart; i<rend; i++) { 2778 const PetscScalar *ba,*ca; 2779 const PetscInt *bj,*cj; 2780 PetscInt bn,cn,j,maxentrycol = -1,maxdiffcol = -1,maxrdiffcol = -1; 2781 PetscReal maxentry = 0,maxdiff = 0,maxrdiff = 0; 2782 ierr = MatGetRow(B,i,&bn,&bj,&ba);CHKERRQ(ierr); 2783 ierr = MatGetRow(Bfd,i,&cn,&cj,&ca);CHKERRQ(ierr); 2784 if (bn != cn) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_PLIB,"Unexpected different nonzero pattern in -snes_compare_coloring_threshold"); 2785 for (j=0; j<bn; j++) { 2786 PetscReal rdiff = PetscAbsScalar(ca[j]) / (threshold_atol + threshold_rtol*PetscAbsScalar(ba[j])); 2787 if (PetscAbsScalar(ba[j]) > PetscAbs(maxentry)) { 2788 maxentrycol = bj[j]; 2789 maxentry = PetscRealPart(ba[j]); 2790 } 2791 if (PetscAbsScalar(ca[j]) > PetscAbs(maxdiff)) { 2792 maxdiffcol = bj[j]; 2793 maxdiff = PetscRealPart(ca[j]); 2794 } 2795 if (rdiff > maxrdiff) { 2796 maxrdiffcol = bj[j]; 2797 maxrdiff = rdiff; 2798 } 2799 } 2800 if (maxrdiff > 1) { 2801 ierr = PetscViewerASCIIPrintf(vstdout,"row %D (maxentry=%g at %D, maxdiff=%g at %D, maxrdiff=%g at %D):",i,(double)maxentry,maxentrycol,(double)maxdiff,maxdiffcol,(double)maxrdiff,maxrdiffcol);CHKERRQ(ierr); 2802 for (j=0; j<bn; j++) { 2803 PetscReal rdiff; 2804 rdiff = PetscAbsScalar(ca[j]) / (threshold_atol + threshold_rtol*PetscAbsScalar(ba[j])); 2805 if (rdiff > 1) { 2806 ierr = PetscViewerASCIIPrintf(vstdout," (%D,%g:%g)",bj[j],(double)PetscRealPart(ba[j]),(double)PetscRealPart(ca[j]));CHKERRQ(ierr); 2807 } 2808 } 2809 ierr = PetscViewerASCIIPrintf(vstdout,"\n",i,maxentry,maxdiff,maxrdiff);CHKERRQ(ierr); 2810 } 2811 ierr = MatRestoreRow(B,i,&bn,&bj,&ba);CHKERRQ(ierr); 2812 ierr = MatRestoreRow(Bfd,i,&cn,&cj,&ca);CHKERRQ(ierr); 2813 } 2814 } 2815 ierr = PetscViewerDestroy(&vdraw);CHKERRQ(ierr); 2816 ierr = MatDestroy(&Bfd);CHKERRQ(ierr); 2817 } 2818 } 2819 PetscFunctionReturn(0); 2820 } 2821 2822 /*MC 2823 SNESJacobianFunction - Function used to convey the nonlinear Jacobian of the function to be solved by SNES 2824 2825 Synopsis: 2826 #include "petscsnes.h" 2827 PetscErrorCode SNESJacobianFunction(SNES snes,Vec x,Mat Amat,Mat Pmat,void *ctx); 2828 2829 + x - input vector 2830 . Amat - the matrix that defines the (approximate) Jacobian 2831 . Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat. 2832 - ctx - [optional] user-defined Jacobian context 2833 2834 Level: intermediate 2835 2836 .seealso: SNESSetFunction(), SNESGetFunction(), SNESSetJacobian(), SNESGetJacobian() 2837 M*/ 2838 2839 /*@C 2840 SNESSetJacobian - Sets the function to compute Jacobian as well as the 2841 location to store the matrix. 2842 2843 Logically Collective on SNES and Mat 2844 2845 Input Parameters: 2846 + snes - the SNES context 2847 . Amat - the matrix that defines the (approximate) Jacobian 2848 . Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat. 2849 . J - Jacobian evaluation routine (if NULL then SNES retains any previously set value), see SNESJacobianFunction for details 2850 - ctx - [optional] user-defined context for private data for the 2851 Jacobian evaluation routine (may be NULL) (if NULL then SNES retains any previously set value) 2852 2853 Notes: 2854 If the Amat matrix and Pmat matrix are different you must call MatAssemblyBegin/End() on 2855 each matrix. 2856 2857 If you know the operator Amat has a null space you can use MatSetNullSpace() and MatSetTransposeNullSpace() to supply the null 2858 space to Amat and the KSP solvers will automatically use that null space as needed during the solution process. 2859 2860 If using SNESComputeJacobianDefaultColor() to assemble a Jacobian, the ctx argument 2861 must be a MatFDColoring. 2862 2863 Other defect-correction schemes can be used by computing a different matrix in place of the Jacobian. One common 2864 example is to use the "Picard linearization" which only differentiates through the highest order parts of each term. 2865 2866 Level: beginner 2867 2868 .seealso: KSPSetOperators(), SNESSetFunction(), MatMFFDComputeJacobian(), SNESComputeJacobianDefaultColor(), MatStructure, J, 2869 SNESSetPicard(), SNESJacobianFunction 2870 @*/ 2871 PetscErrorCode SNESSetJacobian(SNES snes,Mat Amat,Mat Pmat,PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx) 2872 { 2873 PetscErrorCode ierr; 2874 DM dm; 2875 2876 PetscFunctionBegin; 2877 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2878 if (Amat) PetscValidHeaderSpecific(Amat,MAT_CLASSID,2); 2879 if (Pmat) PetscValidHeaderSpecific(Pmat,MAT_CLASSID,3); 2880 if (Amat) PetscCheckSameComm(snes,1,Amat,2); 2881 if (Pmat) PetscCheckSameComm(snes,1,Pmat,3); 2882 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 2883 ierr = DMSNESSetJacobian(dm,J,ctx);CHKERRQ(ierr); 2884 if (Amat) { 2885 ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr); 2886 ierr = MatDestroy(&snes->jacobian);CHKERRQ(ierr); 2887 2888 snes->jacobian = Amat; 2889 } 2890 if (Pmat) { 2891 ierr = PetscObjectReference((PetscObject)Pmat);CHKERRQ(ierr); 2892 ierr = MatDestroy(&snes->jacobian_pre);CHKERRQ(ierr); 2893 2894 snes->jacobian_pre = Pmat; 2895 } 2896 PetscFunctionReturn(0); 2897 } 2898 2899 /*@C 2900 SNESGetJacobian - Returns the Jacobian matrix and optionally the user 2901 provided context for evaluating the Jacobian. 2902 2903 Not Collective, but Mat object will be parallel if SNES object is 2904 2905 Input Parameter: 2906 . snes - the nonlinear solver context 2907 2908 Output Parameters: 2909 + Amat - location to stash (approximate) Jacobian matrix (or NULL) 2910 . Pmat - location to stash matrix used to compute the preconditioner (or NULL) 2911 . J - location to put Jacobian function (or NULL), see SNESJacobianFunction for details on its calling sequence 2912 - ctx - location to stash Jacobian ctx (or NULL) 2913 2914 Level: advanced 2915 2916 .seealso: SNESSetJacobian(), SNESComputeJacobian(), SNESJacobianFunction, SNESGetFunction() 2917 @*/ 2918 PetscErrorCode SNESGetJacobian(SNES snes,Mat *Amat,Mat *Pmat,PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx) 2919 { 2920 PetscErrorCode ierr; 2921 DM dm; 2922 DMSNES sdm; 2923 2924 PetscFunctionBegin; 2925 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2926 if (Amat) *Amat = snes->jacobian; 2927 if (Pmat) *Pmat = snes->jacobian_pre; 2928 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 2929 ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr); 2930 if (J) *J = sdm->ops->computejacobian; 2931 if (ctx) *ctx = sdm->jacobianctx; 2932 PetscFunctionReturn(0); 2933 } 2934 2935 /*@ 2936 SNESSetUp - Sets up the internal data structures for the later use 2937 of a nonlinear solver. 2938 2939 Collective on SNES 2940 2941 Input Parameters: 2942 . snes - the SNES context 2943 2944 Notes: 2945 For basic use of the SNES solvers the user need not explicitly call 2946 SNESSetUp(), since these actions will automatically occur during 2947 the call to SNESSolve(). However, if one wishes to control this 2948 phase separately, SNESSetUp() should be called after SNESCreate() 2949 and optional routines of the form SNESSetXXX(), but before SNESSolve(). 2950 2951 Level: advanced 2952 2953 .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 2954 @*/ 2955 PetscErrorCode SNESSetUp(SNES snes) 2956 { 2957 PetscErrorCode ierr; 2958 DM dm; 2959 DMSNES sdm; 2960 SNESLineSearch linesearch, pclinesearch; 2961 void *lsprectx,*lspostctx; 2962 PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*); 2963 PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*); 2964 PetscErrorCode (*func)(SNES,Vec,Vec,void*); 2965 Vec f,fpc; 2966 void *funcctx; 2967 PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*); 2968 void *jacctx,*appctx; 2969 Mat j,jpre; 2970 2971 PetscFunctionBegin; 2972 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2973 if (snes->setupcalled) PetscFunctionReturn(0); 2974 ierr = PetscLogEventBegin(SNES_Setup,snes,0,0,0);CHKERRQ(ierr); 2975 2976 if (!((PetscObject)snes)->type_name) { 2977 ierr = SNESSetType(snes,SNESNEWTONLS);CHKERRQ(ierr); 2978 } 2979 2980 ierr = SNESGetFunction(snes,&snes->vec_func,NULL,NULL);CHKERRQ(ierr); 2981 2982 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 2983 ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr); 2984 if (!sdm->ops->computefunction) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Function never provided to SNES object"); 2985 if (!sdm->ops->computejacobian) { 2986 ierr = DMSNESSetJacobian(dm,SNESComputeJacobianDefaultColor,NULL);CHKERRQ(ierr); 2987 } 2988 if (!snes->vec_func) { 2989 ierr = DMCreateGlobalVector(dm,&snes->vec_func);CHKERRQ(ierr); 2990 } 2991 2992 if (!snes->ksp) { 2993 ierr = SNESGetKSP(snes, &snes->ksp);CHKERRQ(ierr); 2994 } 2995 2996 if (snes->linesearch) { 2997 ierr = SNESGetLineSearch(snes, &snes->linesearch);CHKERRQ(ierr); 2998 ierr = SNESLineSearchSetFunction(snes->linesearch,SNESComputeFunction);CHKERRQ(ierr); 2999 } 3000 3001 if (snes->npc && (snes->npcside== PC_LEFT)) { 3002 snes->mf = PETSC_TRUE; 3003 snes->mf_operator = PETSC_FALSE; 3004 } 3005 3006 if (snes->npc) { 3007 /* copy the DM over */ 3008 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 3009 ierr = SNESSetDM(snes->npc,dm);CHKERRQ(ierr); 3010 3011 ierr = SNESGetFunction(snes,&f,&func,&funcctx);CHKERRQ(ierr); 3012 ierr = VecDuplicate(f,&fpc);CHKERRQ(ierr); 3013 ierr = SNESSetFunction(snes->npc,fpc,func,funcctx);CHKERRQ(ierr); 3014 ierr = SNESGetJacobian(snes,&j,&jpre,&jac,&jacctx);CHKERRQ(ierr); 3015 ierr = SNESSetJacobian(snes->npc,j,jpre,jac,jacctx);CHKERRQ(ierr); 3016 ierr = SNESGetApplicationContext(snes,&appctx);CHKERRQ(ierr); 3017 ierr = SNESSetApplicationContext(snes->npc,appctx);CHKERRQ(ierr); 3018 ierr = VecDestroy(&fpc);CHKERRQ(ierr); 3019 3020 /* copy the function pointers over */ 3021 ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes,(PetscObject)snes->npc);CHKERRQ(ierr); 3022 3023 /* default to 1 iteration */ 3024 ierr = SNESSetTolerances(snes->npc,0.0,0.0,0.0,1,snes->npc->max_funcs);CHKERRQ(ierr); 3025 if (snes->npcside==PC_RIGHT) { 3026 ierr = SNESSetNormSchedule(snes->npc,SNES_NORM_FINAL_ONLY);CHKERRQ(ierr); 3027 } else { 3028 ierr = SNESSetNormSchedule(snes->npc,SNES_NORM_NONE);CHKERRQ(ierr); 3029 } 3030 ierr = SNESSetFromOptions(snes->npc);CHKERRQ(ierr); 3031 3032 /* copy the line search context over */ 3033 if (snes->linesearch && snes->npc->linesearch) { 3034 ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); 3035 ierr = SNESGetLineSearch(snes->npc,&pclinesearch);CHKERRQ(ierr); 3036 ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr); 3037 ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr); 3038 ierr = SNESLineSearchSetPreCheck(pclinesearch,precheck,lsprectx);CHKERRQ(ierr); 3039 ierr = SNESLineSearchSetPostCheck(pclinesearch,postcheck,lspostctx);CHKERRQ(ierr); 3040 ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)pclinesearch);CHKERRQ(ierr); 3041 } 3042 } 3043 if (snes->mf) { 3044 ierr = SNESSetUpMatrixFree_Private(snes, snes->mf_operator, snes->mf_version);CHKERRQ(ierr); 3045 } 3046 if (snes->ops->usercompute && !snes->user) { 3047 ierr = (*snes->ops->usercompute)(snes,(void**)&snes->user);CHKERRQ(ierr); 3048 } 3049 3050 snes->jac_iter = 0; 3051 snes->pre_iter = 0; 3052 3053 if (snes->ops->setup) { 3054 ierr = (*snes->ops->setup)(snes);CHKERRQ(ierr); 3055 } 3056 3057 if (snes->npc && (snes->npcside== PC_LEFT)) { 3058 if (snes->functype == SNES_FUNCTION_PRECONDITIONED) { 3059 if (snes->linesearch){ 3060 ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); 3061 ierr = SNESLineSearchSetFunction(linesearch,SNESComputeFunctionDefaultNPC);CHKERRQ(ierr); 3062 } 3063 } 3064 } 3065 ierr = PetscLogEventEnd(SNES_Setup,snes,0,0,0);CHKERRQ(ierr); 3066 snes->setupcalled = PETSC_TRUE; 3067 PetscFunctionReturn(0); 3068 } 3069 3070 /*@ 3071 SNESReset - Resets a SNES context to the snessetupcalled = 0 state and removes any allocated Vecs and Mats 3072 3073 Collective on SNES 3074 3075 Input Parameter: 3076 . snes - iterative context obtained from SNESCreate() 3077 3078 Level: intermediate 3079 3080 Notes: 3081 Also calls the application context destroy routine set with SNESSetComputeApplicationContext() 3082 3083 .seealso: SNESCreate(), SNESSetUp(), SNESSolve() 3084 @*/ 3085 PetscErrorCode SNESReset(SNES snes) 3086 { 3087 PetscErrorCode ierr; 3088 3089 PetscFunctionBegin; 3090 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3091 if (snes->ops->userdestroy && snes->user) { 3092 ierr = (*snes->ops->userdestroy)((void**)&snes->user);CHKERRQ(ierr); 3093 snes->user = NULL; 3094 } 3095 if (snes->npc) { 3096 ierr = SNESReset(snes->npc);CHKERRQ(ierr); 3097 } 3098 3099 if (snes->ops->reset) { 3100 ierr = (*snes->ops->reset)(snes);CHKERRQ(ierr); 3101 } 3102 if (snes->ksp) { 3103 ierr = KSPReset(snes->ksp);CHKERRQ(ierr); 3104 } 3105 3106 if (snes->linesearch) { 3107 ierr = SNESLineSearchReset(snes->linesearch);CHKERRQ(ierr); 3108 } 3109 3110 ierr = VecDestroy(&snes->vec_rhs);CHKERRQ(ierr); 3111 ierr = VecDestroy(&snes->vec_sol);CHKERRQ(ierr); 3112 ierr = VecDestroy(&snes->vec_sol_update);CHKERRQ(ierr); 3113 ierr = VecDestroy(&snes->vec_func);CHKERRQ(ierr); 3114 ierr = MatDestroy(&snes->jacobian);CHKERRQ(ierr); 3115 ierr = MatDestroy(&snes->jacobian_pre);CHKERRQ(ierr); 3116 ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr); 3117 ierr = VecDestroyVecs(snes->nvwork,&snes->vwork);CHKERRQ(ierr); 3118 3119 snes->alwayscomputesfinalresidual = PETSC_FALSE; 3120 3121 snes->nwork = snes->nvwork = 0; 3122 snes->setupcalled = PETSC_FALSE; 3123 PetscFunctionReturn(0); 3124 } 3125 3126 /*@ 3127 SNESDestroy - Destroys the nonlinear solver context that was created 3128 with SNESCreate(). 3129 3130 Collective on SNES 3131 3132 Input Parameter: 3133 . snes - the SNES context 3134 3135 Level: beginner 3136 3137 .seealso: SNESCreate(), SNESSolve() 3138 @*/ 3139 PetscErrorCode SNESDestroy(SNES *snes) 3140 { 3141 PetscErrorCode ierr; 3142 3143 PetscFunctionBegin; 3144 if (!*snes) PetscFunctionReturn(0); 3145 PetscValidHeaderSpecific((*snes),SNES_CLASSID,1); 3146 if (--((PetscObject)(*snes))->refct > 0) {*snes = 0; PetscFunctionReturn(0);} 3147 3148 ierr = SNESReset((*snes));CHKERRQ(ierr); 3149 ierr = SNESDestroy(&(*snes)->npc);CHKERRQ(ierr); 3150 3151 /* if memory was published with SAWs then destroy it */ 3152 ierr = PetscObjectSAWsViewOff((PetscObject)*snes);CHKERRQ(ierr); 3153 if ((*snes)->ops->destroy) {ierr = (*((*snes))->ops->destroy)((*snes));CHKERRQ(ierr);} 3154 3155 if ((*snes)->dm) {ierr = DMCoarsenHookRemove((*snes)->dm,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,*snes);CHKERRQ(ierr);} 3156 ierr = DMDestroy(&(*snes)->dm);CHKERRQ(ierr); 3157 ierr = KSPDestroy(&(*snes)->ksp);CHKERRQ(ierr); 3158 ierr = SNESLineSearchDestroy(&(*snes)->linesearch);CHKERRQ(ierr); 3159 3160 ierr = PetscFree((*snes)->kspconvctx);CHKERRQ(ierr); 3161 if ((*snes)->ops->convergeddestroy) { 3162 ierr = (*(*snes)->ops->convergeddestroy)((*snes)->cnvP);CHKERRQ(ierr); 3163 } 3164 if ((*snes)->conv_malloc) { 3165 ierr = PetscFree((*snes)->conv_hist);CHKERRQ(ierr); 3166 ierr = PetscFree((*snes)->conv_hist_its);CHKERRQ(ierr); 3167 } 3168 ierr = SNESMonitorCancel((*snes));CHKERRQ(ierr); 3169 ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr); 3170 PetscFunctionReturn(0); 3171 } 3172 3173 /* ----------- Routines to set solver parameters ---------- */ 3174 3175 /*@ 3176 SNESSetLagPreconditioner - Determines when the preconditioner is rebuilt in the nonlinear solve. 3177 3178 Logically Collective on SNES 3179 3180 Input Parameters: 3181 + snes - the SNES context 3182 - lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time 3183 the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that 3184 3185 Options Database Keys: 3186 . -snes_lag_preconditioner <lag> 3187 3188 Notes: 3189 The default is 1 3190 The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 3191 If -1 is used before the very first nonlinear solve the preconditioner is still built because there is no previous preconditioner to use 3192 3193 Level: intermediate 3194 3195 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian() 3196 3197 @*/ 3198 PetscErrorCode SNESSetLagPreconditioner(SNES snes,PetscInt lag) 3199 { 3200 PetscFunctionBegin; 3201 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3202 if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater"); 3203 if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0"); 3204 PetscValidLogicalCollectiveInt(snes,lag,2); 3205 snes->lagpreconditioner = lag; 3206 PetscFunctionReturn(0); 3207 } 3208 3209 /*@ 3210 SNESSetGridSequence - sets the number of steps of grid sequencing that SNES does 3211 3212 Logically Collective on SNES 3213 3214 Input Parameters: 3215 + snes - the SNES context 3216 - steps - the number of refinements to do, defaults to 0 3217 3218 Options Database Keys: 3219 . -snes_grid_sequence <steps> 3220 3221 Level: intermediate 3222 3223 Notes: 3224 Use SNESGetSolution() to extract the fine grid solution after grid sequencing. 3225 3226 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESGetGridSequence() 3227 3228 @*/ 3229 PetscErrorCode SNESSetGridSequence(SNES snes,PetscInt steps) 3230 { 3231 PetscFunctionBegin; 3232 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3233 PetscValidLogicalCollectiveInt(snes,steps,2); 3234 snes->gridsequence = steps; 3235 PetscFunctionReturn(0); 3236 } 3237 3238 /*@ 3239 SNESGetGridSequence - gets the number of steps of grid sequencing that SNES does 3240 3241 Logically Collective on SNES 3242 3243 Input Parameter: 3244 . snes - the SNES context 3245 3246 Output Parameter: 3247 . steps - the number of refinements to do, defaults to 0 3248 3249 Options Database Keys: 3250 . -snes_grid_sequence <steps> 3251 3252 Level: intermediate 3253 3254 Notes: 3255 Use SNESGetSolution() to extract the fine grid solution after grid sequencing. 3256 3257 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESSetGridSequence() 3258 3259 @*/ 3260 PetscErrorCode SNESGetGridSequence(SNES snes,PetscInt *steps) 3261 { 3262 PetscFunctionBegin; 3263 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3264 *steps = snes->gridsequence; 3265 PetscFunctionReturn(0); 3266 } 3267 3268 /*@ 3269 SNESGetLagPreconditioner - Indicates how often the preconditioner is rebuilt 3270 3271 Not Collective 3272 3273 Input Parameter: 3274 . snes - the SNES context 3275 3276 Output Parameter: 3277 . lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time 3278 the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that 3279 3280 Options Database Keys: 3281 . -snes_lag_preconditioner <lag> 3282 3283 Notes: 3284 The default is 1 3285 The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 3286 3287 Level: intermediate 3288 3289 .seealso: SNESSetTrustRegionTolerance(), SNESSetLagPreconditioner() 3290 3291 @*/ 3292 PetscErrorCode SNESGetLagPreconditioner(SNES snes,PetscInt *lag) 3293 { 3294 PetscFunctionBegin; 3295 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3296 *lag = snes->lagpreconditioner; 3297 PetscFunctionReturn(0); 3298 } 3299 3300 /*@ 3301 SNESSetLagJacobian - Determines when the Jacobian is rebuilt in the nonlinear solve. See SNESSetLagPreconditioner() for determining how 3302 often the preconditioner is rebuilt. 3303 3304 Logically Collective on SNES 3305 3306 Input Parameters: 3307 + snes - the SNES context 3308 - lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time 3309 the Jacobian is built etc. -2 means rebuild at next chance but then never again 3310 3311 Options Database Keys: 3312 . -snes_lag_jacobian <lag> 3313 3314 Notes: 3315 The default is 1 3316 The Jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 3317 If -1 is used before the very first nonlinear solve the CODE WILL FAIL! because no Jacobian is used, use -2 to indicate you want it recomputed 3318 at the next Newton step but never again (unless it is reset to another value) 3319 3320 Level: intermediate 3321 3322 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagPreconditioner(), SNESGetLagJacobian() 3323 3324 @*/ 3325 PetscErrorCode SNESSetLagJacobian(SNES snes,PetscInt lag) 3326 { 3327 PetscFunctionBegin; 3328 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3329 if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater"); 3330 if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0"); 3331 PetscValidLogicalCollectiveInt(snes,lag,2); 3332 snes->lagjacobian = lag; 3333 PetscFunctionReturn(0); 3334 } 3335 3336 /*@ 3337 SNESGetLagJacobian - Indicates how often the Jacobian is rebuilt. See SNESGetLagPreconditioner() to determine when the preconditioner is rebuilt 3338 3339 Not Collective 3340 3341 Input Parameter: 3342 . snes - the SNES context 3343 3344 Output Parameter: 3345 . lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time 3346 the Jacobian is built etc. 3347 3348 Options Database Keys: 3349 . -snes_lag_jacobian <lag> 3350 3351 Notes: 3352 The default is 1 3353 The jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 3354 3355 Level: intermediate 3356 3357 .seealso: SNESSetTrustRegionTolerance(), SNESSetLagJacobian(), SNESSetLagPreconditioner(), SNESGetLagPreconditioner() 3358 3359 @*/ 3360 PetscErrorCode SNESGetLagJacobian(SNES snes,PetscInt *lag) 3361 { 3362 PetscFunctionBegin; 3363 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3364 *lag = snes->lagjacobian; 3365 PetscFunctionReturn(0); 3366 } 3367 3368 /*@ 3369 SNESSetLagJacobianPersists - Set whether or not the Jacobian lagging persists through multiple solves 3370 3371 Logically collective on SNES 3372 3373 Input Parameter: 3374 + snes - the SNES context 3375 - flg - jacobian lagging persists if true 3376 3377 Options Database Keys: 3378 . -snes_lag_jacobian_persists <flg> 3379 3380 Notes: 3381 This is useful both for nonlinear preconditioning, where it's appropriate to have the Jacobian be stale by 3382 several solves, and for implicit time-stepping, where Jacobian lagging in the inner nonlinear solve over several 3383 timesteps may present huge efficiency gains. 3384 3385 Level: developer 3386 3387 .seealso: SNESSetLagPreconditionerPersists(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESGetNPC() 3388 3389 @*/ 3390 PetscErrorCode SNESSetLagJacobianPersists(SNES snes,PetscBool flg) 3391 { 3392 PetscFunctionBegin; 3393 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3394 PetscValidLogicalCollectiveBool(snes,flg,2); 3395 snes->lagjac_persist = flg; 3396 PetscFunctionReturn(0); 3397 } 3398 3399 /*@ 3400 SNESSetLagPreconditionerPersists - Set whether or not the preconditioner lagging persists through multiple solves 3401 3402 Logically Collective on SNES 3403 3404 Input Parameter: 3405 + snes - the SNES context 3406 - flg - preconditioner lagging persists if true 3407 3408 Options Database Keys: 3409 . -snes_lag_jacobian_persists <flg> 3410 3411 Notes: 3412 This is useful both for nonlinear preconditioning, where it's appropriate to have the preconditioner be stale 3413 by several solves, and for implicit time-stepping, where preconditioner lagging in the inner nonlinear solve over 3414 several timesteps may present huge efficiency gains. 3415 3416 Level: developer 3417 3418 .seealso: SNESSetLagJacobianPersists(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESGetNPC() 3419 3420 @*/ 3421 PetscErrorCode SNESSetLagPreconditionerPersists(SNES snes,PetscBool flg) 3422 { 3423 PetscFunctionBegin; 3424 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3425 PetscValidLogicalCollectiveBool(snes,flg,2); 3426 snes->lagpre_persist = flg; 3427 PetscFunctionReturn(0); 3428 } 3429 3430 /*@ 3431 SNESSetForceIteration - force SNESSolve() to take at least one iteration regardless of the initial residual norm 3432 3433 Logically Collective on SNES 3434 3435 Input Parameters: 3436 + snes - the SNES context 3437 - force - PETSC_TRUE require at least one iteration 3438 3439 Options Database Keys: 3440 . -snes_force_iteration <force> - Sets forcing an iteration 3441 3442 Notes: 3443 This is used sometimes with TS to prevent TS from detecting a false steady state solution 3444 3445 Level: intermediate 3446 3447 .seealso: SNESSetTrustRegionTolerance(), SNESSetDivergenceTolerance() 3448 @*/ 3449 PetscErrorCode SNESSetForceIteration(SNES snes,PetscBool force) 3450 { 3451 PetscFunctionBegin; 3452 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3453 snes->forceiteration = force; 3454 PetscFunctionReturn(0); 3455 } 3456 3457 /*@ 3458 SNESGetForceIteration - Whether or not to force SNESSolve() take at least one iteration regardless of the initial residual norm 3459 3460 Logically Collective on SNES 3461 3462 Input Parameters: 3463 . snes - the SNES context 3464 3465 Output Parameter: 3466 . force - PETSC_TRUE requires at least one iteration. 3467 3468 Level: intermediate 3469 3470 .seealso: SNESSetForceIteration(), SNESSetTrustRegionTolerance(), SNESSetDivergenceTolerance() 3471 @*/ 3472 PetscErrorCode SNESGetForceIteration(SNES snes,PetscBool *force) 3473 { 3474 PetscFunctionBegin; 3475 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3476 *force = snes->forceiteration; 3477 PetscFunctionReturn(0); 3478 } 3479 3480 /*@ 3481 SNESSetTolerances - Sets various parameters used in convergence tests. 3482 3483 Logically Collective on SNES 3484 3485 Input Parameters: 3486 + snes - the SNES context 3487 . abstol - absolute convergence tolerance 3488 . rtol - relative convergence tolerance 3489 . stol - convergence tolerance in terms of the norm of the change in the solution between steps, || delta x || < stol*|| x || 3490 . maxit - maximum number of iterations 3491 - maxf - maximum number of function evaluations (-1 indicates no limit) 3492 3493 Options Database Keys: 3494 + -snes_atol <abstol> - Sets abstol 3495 . -snes_rtol <rtol> - Sets rtol 3496 . -snes_stol <stol> - Sets stol 3497 . -snes_max_it <maxit> - Sets maxit 3498 - -snes_max_funcs <maxf> - Sets maxf 3499 3500 Notes: 3501 The default maximum number of iterations is 50. 3502 The default maximum number of function evaluations is 1000. 3503 3504 Level: intermediate 3505 3506 .seealso: SNESSetTrustRegionTolerance(), SNESSetDivergenceTolerance(), SNESSetForceIteration() 3507 @*/ 3508 PetscErrorCode SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf) 3509 { 3510 PetscFunctionBegin; 3511 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3512 PetscValidLogicalCollectiveReal(snes,abstol,2); 3513 PetscValidLogicalCollectiveReal(snes,rtol,3); 3514 PetscValidLogicalCollectiveReal(snes,stol,4); 3515 PetscValidLogicalCollectiveInt(snes,maxit,5); 3516 PetscValidLogicalCollectiveInt(snes,maxf,6); 3517 3518 if (abstol != PETSC_DEFAULT) { 3519 if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol); 3520 snes->abstol = abstol; 3521 } 3522 if (rtol != PETSC_DEFAULT) { 3523 if (rtol < 0.0 || 1.0 <= rtol) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Relative tolerance %g must be non-negative and less than 1.0",(double)rtol); 3524 snes->rtol = rtol; 3525 } 3526 if (stol != PETSC_DEFAULT) { 3527 if (stol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Step tolerance %g must be non-negative",(double)stol); 3528 snes->stol = stol; 3529 } 3530 if (maxit != PETSC_DEFAULT) { 3531 if (maxit < 0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxit); 3532 snes->max_its = maxit; 3533 } 3534 if (maxf != PETSC_DEFAULT) { 3535 if (maxf < -1) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of function evaluations %D must be -1 or nonnegative",maxf); 3536 snes->max_funcs = maxf; 3537 } 3538 snes->tolerancesset = PETSC_TRUE; 3539 PetscFunctionReturn(0); 3540 } 3541 3542 /*@ 3543 SNESSetDivergenceTolerance - Sets the divergence tolerance used for the SNES divergence test. 3544 3545 Logically Collective on SNES 3546 3547 Input Parameters: 3548 + snes - the SNES context 3549 - divtol - the divergence tolerance. Use -1 to deactivate the test. 3550 3551 Options Database Keys: 3552 + -snes_divergence_tolerance <divtol> - Sets divtol 3553 3554 Notes: 3555 The default divergence tolerance is 1e4. 3556 3557 Level: intermediate 3558 3559 .seealso: SNESSetTolerances(), SNESGetDivergenceTolerance 3560 @*/ 3561 PetscErrorCode SNESSetDivergenceTolerance(SNES snes,PetscReal divtol) 3562 { 3563 PetscFunctionBegin; 3564 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3565 PetscValidLogicalCollectiveReal(snes,divtol,2); 3566 3567 if (divtol != PETSC_DEFAULT) { 3568 snes->divtol = divtol; 3569 } 3570 else { 3571 snes->divtol = 1.0e4; 3572 } 3573 PetscFunctionReturn(0); 3574 } 3575 3576 /*@ 3577 SNESGetTolerances - Gets various parameters used in convergence tests. 3578 3579 Not Collective 3580 3581 Input Parameters: 3582 + snes - the SNES context 3583 . atol - absolute convergence tolerance 3584 . rtol - relative convergence tolerance 3585 . stol - convergence tolerance in terms of the norm 3586 of the change in the solution between steps 3587 . maxit - maximum number of iterations 3588 - maxf - maximum number of function evaluations 3589 3590 Notes: 3591 The user can specify NULL for any parameter that is not needed. 3592 3593 Level: intermediate 3594 3595 .seealso: SNESSetTolerances() 3596 @*/ 3597 PetscErrorCode SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf) 3598 { 3599 PetscFunctionBegin; 3600 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3601 if (atol) *atol = snes->abstol; 3602 if (rtol) *rtol = snes->rtol; 3603 if (stol) *stol = snes->stol; 3604 if (maxit) *maxit = snes->max_its; 3605 if (maxf) *maxf = snes->max_funcs; 3606 PetscFunctionReturn(0); 3607 } 3608 3609 /*@ 3610 SNESGetDivergenceTolerance - Gets divergence tolerance used in divergence test. 3611 3612 Not Collective 3613 3614 Input Parameters: 3615 + snes - the SNES context 3616 - divtol - divergence tolerance 3617 3618 Level: intermediate 3619 3620 .seealso: SNESSetDivergenceTolerance() 3621 @*/ 3622 PetscErrorCode SNESGetDivergenceTolerance(SNES snes,PetscReal *divtol) 3623 { 3624 PetscFunctionBegin; 3625 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3626 if (divtol) *divtol = snes->divtol; 3627 PetscFunctionReturn(0); 3628 } 3629 3630 /*@ 3631 SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 3632 3633 Logically Collective on SNES 3634 3635 Input Parameters: 3636 + snes - the SNES context 3637 - tol - tolerance 3638 3639 Options Database Key: 3640 . -snes_trtol <tol> - Sets tol 3641 3642 Level: intermediate 3643 3644 .seealso: SNESSetTolerances() 3645 @*/ 3646 PetscErrorCode SNESSetTrustRegionTolerance(SNES snes,PetscReal tol) 3647 { 3648 PetscFunctionBegin; 3649 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3650 PetscValidLogicalCollectiveReal(snes,tol,2); 3651 snes->deltatol = tol; 3652 PetscFunctionReturn(0); 3653 } 3654 3655 /* 3656 Duplicate the lg monitors for SNES from KSP; for some reason with 3657 dynamic libraries things don't work under Sun4 if we just use 3658 macros instead of functions 3659 */ 3660 PetscErrorCode SNESMonitorLGResidualNorm(SNES snes,PetscInt it,PetscReal norm,void *ctx) 3661 { 3662 PetscErrorCode ierr; 3663 3664 PetscFunctionBegin; 3665 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3666 ierr = KSPMonitorLGResidualNorm((KSP)snes,it,norm,ctx);CHKERRQ(ierr); 3667 PetscFunctionReturn(0); 3668 } 3669 3670 PetscErrorCode SNESMonitorLGCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *lgctx) 3671 { 3672 PetscErrorCode ierr; 3673 3674 PetscFunctionBegin; 3675 ierr = KSPMonitorLGResidualNormCreate(comm,host,label,x,y,m,n,lgctx);CHKERRQ(ierr); 3676 PetscFunctionReturn(0); 3677 } 3678 3679 PETSC_INTERN PetscErrorCode SNESMonitorRange_Private(SNES,PetscInt,PetscReal*); 3680 3681 PetscErrorCode SNESMonitorLGRange(SNES snes,PetscInt n,PetscReal rnorm,void *monctx) 3682 { 3683 PetscDrawLG lg; 3684 PetscErrorCode ierr; 3685 PetscReal x,y,per; 3686 PetscViewer v = (PetscViewer)monctx; 3687 static PetscReal prev; /* should be in the context */ 3688 PetscDraw draw; 3689 3690 PetscFunctionBegin; 3691 PetscValidHeaderSpecific(v,PETSC_VIEWER_CLASSID,4); 3692 ierr = PetscViewerDrawGetDrawLG(v,0,&lg);CHKERRQ(ierr); 3693 if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 3694 ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr); 3695 ierr = PetscDrawSetTitle(draw,"Residual norm");CHKERRQ(ierr); 3696 x = (PetscReal)n; 3697 if (rnorm > 0.0) y = PetscLog10Real(rnorm); 3698 else y = -15.0; 3699 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 3700 if (n < 20 || !(n % 5) || snes->reason) { 3701 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 3702 ierr = PetscDrawLGSave(lg);CHKERRQ(ierr); 3703 } 3704 3705 ierr = PetscViewerDrawGetDrawLG(v,1,&lg);CHKERRQ(ierr); 3706 if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 3707 ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr); 3708 ierr = PetscDrawSetTitle(draw,"% elemts > .2*max elemt");CHKERRQ(ierr); 3709 ierr = SNESMonitorRange_Private(snes,n,&per);CHKERRQ(ierr); 3710 x = (PetscReal)n; 3711 y = 100.0*per; 3712 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 3713 if (n < 20 || !(n % 5) || snes->reason) { 3714 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 3715 ierr = PetscDrawLGSave(lg);CHKERRQ(ierr); 3716 } 3717 3718 ierr = PetscViewerDrawGetDrawLG(v,2,&lg);CHKERRQ(ierr); 3719 if (!n) {prev = rnorm;ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 3720 ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr); 3721 ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm");CHKERRQ(ierr); 3722 x = (PetscReal)n; 3723 y = (prev - rnorm)/prev; 3724 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 3725 if (n < 20 || !(n % 5) || snes->reason) { 3726 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 3727 ierr = PetscDrawLGSave(lg);CHKERRQ(ierr); 3728 } 3729 3730 ierr = PetscViewerDrawGetDrawLG(v,3,&lg);CHKERRQ(ierr); 3731 if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 3732 ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr); 3733 ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");CHKERRQ(ierr); 3734 x = (PetscReal)n; 3735 y = (prev - rnorm)/(prev*per); 3736 if (n > 2) { /*skip initial crazy value */ 3737 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 3738 } 3739 if (n < 20 || !(n % 5) || snes->reason) { 3740 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 3741 ierr = PetscDrawLGSave(lg);CHKERRQ(ierr); 3742 } 3743 prev = rnorm; 3744 PetscFunctionReturn(0); 3745 } 3746 3747 /*@ 3748 SNESMonitor - runs the user provided monitor routines, if they exist 3749 3750 Collective on SNES 3751 3752 Input Parameters: 3753 + snes - nonlinear solver context obtained from SNESCreate() 3754 . iter - iteration number 3755 - rnorm - relative norm of the residual 3756 3757 Notes: 3758 This routine is called by the SNES implementations. 3759 It does not typically need to be called by the user. 3760 3761 Level: developer 3762 3763 .seealso: SNESMonitorSet() 3764 @*/ 3765 PetscErrorCode SNESMonitor(SNES snes,PetscInt iter,PetscReal rnorm) 3766 { 3767 PetscErrorCode ierr; 3768 PetscInt i,n = snes->numbermonitors; 3769 3770 PetscFunctionBegin; 3771 ierr = VecLockReadPush(snes->vec_sol);CHKERRQ(ierr); 3772 for (i=0; i<n; i++) { 3773 ierr = (*snes->monitor[i])(snes,iter,rnorm,snes->monitorcontext[i]);CHKERRQ(ierr); 3774 } 3775 ierr = VecLockReadPop(snes->vec_sol);CHKERRQ(ierr); 3776 PetscFunctionReturn(0); 3777 } 3778 3779 /* ------------ Routines to set performance monitoring options ----------- */ 3780 3781 /*MC 3782 SNESMonitorFunction - functional form passed to SNESMonitorSet() to monitor convergence of nonlinear solver 3783 3784 Synopsis: 3785 #include <petscsnes.h> 3786 $ PetscErrorCode SNESMonitorFunction(SNES snes,PetscInt its, PetscReal norm,void *mctx) 3787 3788 + snes - the SNES context 3789 . its - iteration number 3790 . norm - 2-norm function value (may be estimated) 3791 - mctx - [optional] monitoring context 3792 3793 Level: advanced 3794 3795 .seealso: SNESMonitorSet(), SNESMonitorGet() 3796 M*/ 3797 3798 /*@C 3799 SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every 3800 iteration of the nonlinear solver to display the iteration's 3801 progress. 3802 3803 Logically Collective on SNES 3804 3805 Input Parameters: 3806 + snes - the SNES context 3807 . f - the monitor function, see SNESMonitorFunction for the calling sequence 3808 . mctx - [optional] user-defined context for private data for the 3809 monitor routine (use NULL if no context is desired) 3810 - monitordestroy - [optional] routine that frees monitor context 3811 (may be NULL) 3812 3813 Options Database Keys: 3814 + -snes_monitor - sets SNESMonitorDefault() 3815 . -snes_monitor_lg_residualnorm - sets line graph monitor, 3816 uses SNESMonitorLGCreate() 3817 - -snes_monitor_cancel - cancels all monitors that have 3818 been hardwired into a code by 3819 calls to SNESMonitorSet(), but 3820 does not cancel those set via 3821 the options database. 3822 3823 Notes: 3824 Several different monitoring routines may be set by calling 3825 SNESMonitorSet() multiple times; all will be called in the 3826 order in which they were set. 3827 3828 Fortran Notes: 3829 Only a single monitor function can be set for each SNES object 3830 3831 Level: intermediate 3832 3833 .seealso: SNESMonitorDefault(), SNESMonitorCancel(), SNESMonitorFunction 3834 @*/ 3835 PetscErrorCode SNESMonitorSet(SNES snes,PetscErrorCode (*f)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**)) 3836 { 3837 PetscInt i; 3838 PetscErrorCode ierr; 3839 PetscBool identical; 3840 3841 PetscFunctionBegin; 3842 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3843 for (i=0; i<snes->numbermonitors;i++) { 3844 ierr = PetscMonitorCompare((PetscErrorCode (*)(void))f,mctx,monitordestroy,(PetscErrorCode (*)(void))snes->monitor[i],snes->monitorcontext[i],snes->monitordestroy[i],&identical);CHKERRQ(ierr); 3845 if (identical) PetscFunctionReturn(0); 3846 } 3847 if (snes->numbermonitors >= MAXSNESMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 3848 snes->monitor[snes->numbermonitors] = f; 3849 snes->monitordestroy[snes->numbermonitors] = monitordestroy; 3850 snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 3851 PetscFunctionReturn(0); 3852 } 3853 3854 /*@ 3855 SNESMonitorCancel - Clears all the monitor functions for a SNES object. 3856 3857 Logically Collective on SNES 3858 3859 Input Parameters: 3860 . snes - the SNES context 3861 3862 Options Database Key: 3863 . -snes_monitor_cancel - cancels all monitors that have been hardwired 3864 into a code by calls to SNESMonitorSet(), but does not cancel those 3865 set via the options database 3866 3867 Notes: 3868 There is no way to clear one specific monitor from a SNES object. 3869 3870 Level: intermediate 3871 3872 .seealso: SNESMonitorDefault(), SNESMonitorSet() 3873 @*/ 3874 PetscErrorCode SNESMonitorCancel(SNES snes) 3875 { 3876 PetscErrorCode ierr; 3877 PetscInt i; 3878 3879 PetscFunctionBegin; 3880 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3881 for (i=0; i<snes->numbermonitors; i++) { 3882 if (snes->monitordestroy[i]) { 3883 ierr = (*snes->monitordestroy[i])(&snes->monitorcontext[i]);CHKERRQ(ierr); 3884 } 3885 } 3886 snes->numbermonitors = 0; 3887 PetscFunctionReturn(0); 3888 } 3889 3890 /*MC 3891 SNESConvergenceTestFunction - functional form used for testing of convergence of nonlinear solver 3892 3893 Synopsis: 3894 #include <petscsnes.h> 3895 $ PetscErrorCode SNESConvergenceTest(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx) 3896 3897 + snes - the SNES context 3898 . it - current iteration (0 is the first and is before any Newton step) 3899 . cctx - [optional] convergence context 3900 . reason - reason for convergence/divergence 3901 . xnorm - 2-norm of current iterate 3902 . gnorm - 2-norm of current step 3903 - f - 2-norm of function 3904 3905 Level: intermediate 3906 3907 .seealso: SNESSetConvergenceTest(), SNESGetConvergenceTest() 3908 M*/ 3909 3910 /*@C 3911 SNESSetConvergenceTest - Sets the function that is to be used 3912 to test for convergence of the nonlinear iterative solution. 3913 3914 Logically Collective on SNES 3915 3916 Input Parameters: 3917 + snes - the SNES context 3918 . SNESConvergenceTestFunction - routine to test for convergence 3919 . cctx - [optional] context for private data for the convergence routine (may be NULL) 3920 - destroy - [optional] destructor for the context (may be NULL; PETSC_NULL_FUNCTION in Fortran) 3921 3922 Level: advanced 3923 3924 .seealso: SNESConvergedDefault(), SNESConvergedSkip(), SNESConvergenceTestFunction 3925 @*/ 3926 PetscErrorCode SNESSetConvergenceTest(SNES snes,PetscErrorCode (*SNESConvergenceTestFunction)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*)) 3927 { 3928 PetscErrorCode ierr; 3929 3930 PetscFunctionBegin; 3931 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3932 if (!SNESConvergenceTestFunction) SNESConvergenceTestFunction = SNESConvergedSkip; 3933 if (snes->ops->convergeddestroy) { 3934 ierr = (*snes->ops->convergeddestroy)(snes->cnvP);CHKERRQ(ierr); 3935 } 3936 snes->ops->converged = SNESConvergenceTestFunction; 3937 snes->ops->convergeddestroy = destroy; 3938 snes->cnvP = cctx; 3939 PetscFunctionReturn(0); 3940 } 3941 3942 /*@ 3943 SNESGetConvergedReason - Gets the reason the SNES iteration was stopped. 3944 3945 Not Collective 3946 3947 Input Parameter: 3948 . snes - the SNES context 3949 3950 Output Parameter: 3951 . reason - negative value indicates diverged, positive value converged, see SNESConvergedReason or the 3952 manual pages for the individual convergence tests for complete lists 3953 3954 Options Database: 3955 . -snes_converged_reason - prints the reason to standard out 3956 3957 Level: intermediate 3958 3959 Notes: 3960 Should only be called after the call the SNESSolve() is complete, if it is called earlier it returns the value SNES__CONVERGED_ITERATING. 3961 3962 .seealso: SNESSetConvergenceTest(), SNESSetConvergedReason(), SNESConvergedReason 3963 @*/ 3964 PetscErrorCode SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason) 3965 { 3966 PetscFunctionBegin; 3967 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3968 PetscValidPointer(reason,2); 3969 *reason = snes->reason; 3970 PetscFunctionReturn(0); 3971 } 3972 3973 /*@ 3974 SNESSetConvergedReason - Sets the reason the SNES iteration was stopped. 3975 3976 Not Collective 3977 3978 Input Parameters: 3979 + snes - the SNES context 3980 - reason - negative value indicates diverged, positive value converged, see SNESConvergedReason or the 3981 manual pages for the individual convergence tests for complete lists 3982 3983 Level: intermediate 3984 3985 .seealso: SNESGetConvergedReason(), SNESSetConvergenceTest(), SNESConvergedReason 3986 @*/ 3987 PetscErrorCode SNESSetConvergedReason(SNES snes,SNESConvergedReason reason) 3988 { 3989 PetscFunctionBegin; 3990 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3991 snes->reason = reason; 3992 PetscFunctionReturn(0); 3993 } 3994 3995 /*@ 3996 SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 3997 3998 Logically Collective on SNES 3999 4000 Input Parameters: 4001 + snes - iterative context obtained from SNESCreate() 4002 . a - array to hold history, this array will contain the function norms computed at each step 4003 . its - integer array holds the number of linear iterations for each solve. 4004 . na - size of a and its 4005 - reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero, 4006 else it continues storing new values for new nonlinear solves after the old ones 4007 4008 Notes: 4009 If 'a' and 'its' are NULL then space is allocated for the history. If 'na' PETSC_DECIDE or PETSC_DEFAULT then a 4010 default array of length 10000 is allocated. 4011 4012 This routine is useful, e.g., when running a code for purposes 4013 of accurate performance monitoring, when no I/O should be done 4014 during the section of code that is being timed. 4015 4016 Level: intermediate 4017 4018 .seealso: SNESGetConvergenceHistory() 4019 4020 @*/ 4021 PetscErrorCode SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscBool reset) 4022 { 4023 PetscErrorCode ierr; 4024 4025 PetscFunctionBegin; 4026 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 4027 if (a) PetscValidScalarPointer(a,2); 4028 if (its) PetscValidIntPointer(its,3); 4029 if (!a) { 4030 if (na == PETSC_DECIDE || na == PETSC_DEFAULT) na = 1000; 4031 ierr = PetscCalloc1(na,&a);CHKERRQ(ierr); 4032 ierr = PetscCalloc1(na,&its);CHKERRQ(ierr); 4033 4034 snes->conv_malloc = PETSC_TRUE; 4035 } 4036 snes->conv_hist = a; 4037 snes->conv_hist_its = its; 4038 snes->conv_hist_max = na; 4039 snes->conv_hist_len = 0; 4040 snes->conv_hist_reset = reset; 4041 PetscFunctionReturn(0); 4042 } 4043 4044 #if defined(PETSC_HAVE_MATLAB_ENGINE) 4045 #include <engine.h> /* MATLAB include file */ 4046 #include <mex.h> /* MATLAB include file */ 4047 4048 PETSC_EXTERN mxArray *SNESGetConvergenceHistoryMatlab(SNES snes) 4049 { 4050 mxArray *mat; 4051 PetscInt i; 4052 PetscReal *ar; 4053 4054 PetscFunctionBegin; 4055 mat = mxCreateDoubleMatrix(snes->conv_hist_len,1,mxREAL); 4056 ar = (PetscReal*) mxGetData(mat); 4057 for (i=0; i<snes->conv_hist_len; i++) ar[i] = snes->conv_hist[i]; 4058 PetscFunctionReturn(mat); 4059 } 4060 #endif 4061 4062 /*@C 4063 SNESGetConvergenceHistory - Gets the array used to hold the convergence history. 4064 4065 Not Collective 4066 4067 Input Parameter: 4068 . snes - iterative context obtained from SNESCreate() 4069 4070 Output Parameters: 4071 . a - array to hold history 4072 . its - integer array holds the number of linear iterations (or 4073 negative if not converged) for each solve. 4074 - na - size of a and its 4075 4076 Notes: 4077 The calling sequence for this routine in Fortran is 4078 $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr) 4079 4080 This routine is useful, e.g., when running a code for purposes 4081 of accurate performance monitoring, when no I/O should be done 4082 during the section of code that is being timed. 4083 4084 Level: intermediate 4085 4086 .seealso: SNESSetConvergencHistory() 4087 4088 @*/ 4089 PetscErrorCode SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na) 4090 { 4091 PetscFunctionBegin; 4092 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 4093 if (a) *a = snes->conv_hist; 4094 if (its) *its = snes->conv_hist_its; 4095 if (na) *na = snes->conv_hist_len; 4096 PetscFunctionReturn(0); 4097 } 4098 4099 /*@C 4100 SNESSetUpdate - Sets the general-purpose update function called 4101 at the beginning of every iteration of the nonlinear solve. Specifically 4102 it is called just before the Jacobian is "evaluated". 4103 4104 Logically Collective on SNES 4105 4106 Input Parameters: 4107 . snes - The nonlinear solver context 4108 . func - The function 4109 4110 Calling sequence of func: 4111 . func (SNES snes, PetscInt step); 4112 4113 . step - The current step of the iteration 4114 4115 Level: advanced 4116 4117 Note: This is NOT what one uses to update the ghost points before a function evaluation, that should be done at the beginning of your FormFunction() 4118 This is not used by most users. 4119 4120 .seealso SNESSetJacobian(), SNESSolve() 4121 @*/ 4122 PetscErrorCode SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt)) 4123 { 4124 PetscFunctionBegin; 4125 PetscValidHeaderSpecific(snes, SNES_CLASSID,1); 4126 snes->ops->update = func; 4127 PetscFunctionReturn(0); 4128 } 4129 4130 /* 4131 SNESScaleStep_Private - Scales a step so that its length is less than the 4132 positive parameter delta. 4133 4134 Input Parameters: 4135 + snes - the SNES context 4136 . y - approximate solution of linear system 4137 . fnorm - 2-norm of current function 4138 - delta - trust region size 4139 4140 Output Parameters: 4141 + gpnorm - predicted function norm at the new point, assuming local 4142 linearization. The value is zero if the step lies within the trust 4143 region, and exceeds zero otherwise. 4144 - ynorm - 2-norm of the step 4145 4146 Note: 4147 For non-trust region methods such as SNESNEWTONLS, the parameter delta 4148 is set to be the maximum allowable step size. 4149 4150 */ 4151 PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm) 4152 { 4153 PetscReal nrm; 4154 PetscScalar cnorm; 4155 PetscErrorCode ierr; 4156 4157 PetscFunctionBegin; 4158 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 4159 PetscValidHeaderSpecific(y,VEC_CLASSID,2); 4160 PetscCheckSameComm(snes,1,y,2); 4161 4162 ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr); 4163 if (nrm > *delta) { 4164 nrm = *delta/nrm; 4165 *gpnorm = (1.0 - nrm)*(*fnorm); 4166 cnorm = nrm; 4167 ierr = VecScale(y,cnorm);CHKERRQ(ierr); 4168 *ynorm = *delta; 4169 } else { 4170 *gpnorm = 0.0; 4171 *ynorm = nrm; 4172 } 4173 PetscFunctionReturn(0); 4174 } 4175 4176 /*@ 4177 SNESReasonView - Displays the reason a SNES solve converged or diverged to a viewer 4178 4179 Collective on SNES 4180 4181 Parameter: 4182 + snes - iterative context obtained from SNESCreate() 4183 - viewer - the viewer to display the reason 4184 4185 4186 Options Database Keys: 4187 . -snes_converged_reason - print reason for converged or diverged, also prints number of iterations 4188 4189 Level: beginner 4190 4191 .seealso: SNESCreate(), SNESSetUp(), SNESDestroy(), SNESSetTolerances(), SNESConvergedDefault() 4192 4193 @*/ 4194 PetscErrorCode SNESReasonView(SNES snes,PetscViewer viewer) 4195 { 4196 PetscViewerFormat format; 4197 PetscBool isAscii; 4198 PetscErrorCode ierr; 4199 4200 PetscFunctionBegin; 4201 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii);CHKERRQ(ierr); 4202 if (isAscii) { 4203 ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 4204 ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 4205 if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 4206 DM dm; 4207 Vec u; 4208 PetscDS prob; 4209 PetscInt Nf, f; 4210 PetscErrorCode (**exactSol)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *); 4211 void **exactCtx; 4212 PetscReal error; 4213 4214 ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr); 4215 ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr); 4216 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 4217 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 4218 ierr = PetscMalloc2(Nf, &exactSol, Nf, &exactCtx);CHKERRQ(ierr); 4219 for (f = 0; f < Nf; ++f) {ierr = PetscDSGetExactSolution(prob, f, &exactSol[f], &exactCtx[f]);CHKERRQ(ierr);} 4220 ierr = DMComputeL2Diff(dm, 0.0, exactSol, exactCtx, u, &error);CHKERRQ(ierr); 4221 ierr = PetscFree2(exactSol, exactCtx);CHKERRQ(ierr); 4222 if (error < 1.0e-11) {ierr = PetscViewerASCIIPrintf(viewer, "L_2 Error: < 1.0e-11\n");CHKERRQ(ierr);} 4223 else {ierr = PetscViewerASCIIPrintf(viewer, "L_2 Error: %g\n", error);CHKERRQ(ierr);} 4224 } 4225 if (snes->reason > 0) { 4226 if (((PetscObject) snes)->prefix) { 4227 ierr = PetscViewerASCIIPrintf(viewer,"Nonlinear %s solve converged due to %s iterations %D\n",((PetscObject) snes)->prefix,SNESConvergedReasons[snes->reason],snes->iter);CHKERRQ(ierr); 4228 } else { 4229 ierr = PetscViewerASCIIPrintf(viewer,"Nonlinear solve converged due to %s iterations %D\n",SNESConvergedReasons[snes->reason],snes->iter);CHKERRQ(ierr); 4230 } 4231 } else { 4232 if (((PetscObject) snes)->prefix) { 4233 ierr = PetscViewerASCIIPrintf(viewer,"Nonlinear %s solve did not converge due to %s iterations %D\n",((PetscObject) snes)->prefix,SNESConvergedReasons[snes->reason],snes->iter);CHKERRQ(ierr); 4234 } else { 4235 ierr = PetscViewerASCIIPrintf(viewer,"Nonlinear solve did not converge due to %s iterations %D\n",SNESConvergedReasons[snes->reason],snes->iter);CHKERRQ(ierr); 4236 } 4237 } 4238 ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 4239 } 4240 PetscFunctionReturn(0); 4241 } 4242 4243 /*@C 4244 SNESReasonViewFromOptions - Processes command line options to determine if/how a SNESReason is to be viewed. 4245 4246 Collective on SNES 4247 4248 Input Parameters: 4249 . snes - the SNES object 4250 4251 Level: intermediate 4252 4253 @*/ 4254 PetscErrorCode SNESReasonViewFromOptions(SNES snes) 4255 { 4256 PetscErrorCode ierr; 4257 PetscViewer viewer; 4258 PetscBool flg; 4259 static PetscBool incall = PETSC_FALSE; 4260 PetscViewerFormat format; 4261 4262 PetscFunctionBegin; 4263 if (incall) PetscFunctionReturn(0); 4264 incall = PETSC_TRUE; 4265 ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_converged_reason",&viewer,&format,&flg);CHKERRQ(ierr); 4266 if (flg) { 4267 ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr); 4268 ierr = SNESReasonView(snes,viewer);CHKERRQ(ierr); 4269 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 4270 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 4271 } 4272 incall = PETSC_FALSE; 4273 PetscFunctionReturn(0); 4274 } 4275 4276 /*@ 4277 SNESSolve - Solves a nonlinear system F(x) = b. 4278 Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX(). 4279 4280 Collective on SNES 4281 4282 Input Parameters: 4283 + snes - the SNES context 4284 . b - the constant part of the equation F(x) = b, or NULL to use zero. 4285 - x - the solution vector. 4286 4287 Notes: 4288 The user should initialize the vector,x, with the initial guess 4289 for the nonlinear solve prior to calling SNESSolve. In particular, 4290 to employ an initial guess of zero, the user should explicitly set 4291 this vector to zero by calling VecSet(). 4292 4293 Level: beginner 4294 4295 .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetGridSequence(), SNESGetSolution() 4296 @*/ 4297 PetscErrorCode SNESSolve(SNES snes,Vec b,Vec x) 4298 { 4299 PetscErrorCode ierr; 4300 PetscBool flg; 4301 PetscInt grid; 4302 Vec xcreated = NULL; 4303 DM dm; 4304 4305 PetscFunctionBegin; 4306 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 4307 if (x) PetscValidHeaderSpecific(x,VEC_CLASSID,3); 4308 if (x) PetscCheckSameComm(snes,1,x,3); 4309 if (b) PetscValidHeaderSpecific(b,VEC_CLASSID,2); 4310 if (b) PetscCheckSameComm(snes,1,b,2); 4311 4312 /* High level operations using the nonlinear solver */ 4313 { 4314 PetscViewer viewer; 4315 PetscViewerFormat format; 4316 PetscInt num; 4317 PetscBool flg; 4318 static PetscBool incall = PETSC_FALSE; 4319 4320 if (!incall) { 4321 /* Estimate the convergence rate of the discretization */ 4322 ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) snes),((PetscObject)snes)->options, ((PetscObject) snes)->prefix, "-snes_convergence_estimate", &viewer, &format, &flg);CHKERRQ(ierr); 4323 if (flg) { 4324 PetscConvEst conv; 4325 DM dm; 4326 PetscReal *alpha; /* Convergence rate of the solution error for each field in the L_2 norm */ 4327 PetscInt Nf; 4328 4329 incall = PETSC_TRUE; 4330 ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr); 4331 ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 4332 ierr = PetscMalloc1(Nf, &alpha);CHKERRQ(ierr); 4333 ierr = PetscConvEstCreate(PetscObjectComm((PetscObject) snes), &conv);CHKERRQ(ierr); 4334 ierr = PetscConvEstSetSolver(conv, snes);CHKERRQ(ierr); 4335 ierr = PetscConvEstSetFromOptions(conv);CHKERRQ(ierr); 4336 ierr = PetscConvEstSetUp(conv);CHKERRQ(ierr); 4337 ierr = PetscConvEstGetConvRate(conv, alpha);CHKERRQ(ierr); 4338 ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr); 4339 ierr = PetscConvEstRateView(conv, alpha, viewer);CHKERRQ(ierr); 4340 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 4341 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 4342 ierr = PetscConvEstDestroy(&conv);CHKERRQ(ierr); 4343 ierr = PetscFree(alpha);CHKERRQ(ierr); 4344 incall = PETSC_FALSE; 4345 } 4346 /* Adaptively refine the initial grid */ 4347 num = 1; 4348 ierr = PetscOptionsGetInt(NULL, ((PetscObject) snes)->prefix, "-snes_adapt_initial", &num, &flg);CHKERRQ(ierr); 4349 if (flg) { 4350 DMAdaptor adaptor; 4351 4352 incall = PETSC_TRUE; 4353 ierr = DMAdaptorCreate(PETSC_COMM_WORLD, &adaptor);CHKERRQ(ierr); 4354 ierr = DMAdaptorSetSolver(adaptor, snes);CHKERRQ(ierr); 4355 ierr = DMAdaptorSetSequenceLength(adaptor, num);CHKERRQ(ierr); 4356 ierr = DMAdaptorSetFromOptions(adaptor);CHKERRQ(ierr); 4357 ierr = DMAdaptorSetUp(adaptor);CHKERRQ(ierr); 4358 ierr = DMAdaptorAdapt(adaptor, x, DM_ADAPTATION_INITIAL, &dm, &x);CHKERRQ(ierr); 4359 ierr = DMAdaptorDestroy(&adaptor);CHKERRQ(ierr); 4360 incall = PETSC_FALSE; 4361 } 4362 /* Use grid sequencing to adapt */ 4363 num = 0; 4364 ierr = PetscOptionsGetInt(NULL, ((PetscObject) snes)->prefix, "-snes_adapt_sequence", &num, NULL);CHKERRQ(ierr); 4365 if (num) { 4366 DMAdaptor adaptor; 4367 4368 incall = PETSC_TRUE; 4369 ierr = DMAdaptorCreate(PETSC_COMM_WORLD, &adaptor);CHKERRQ(ierr); 4370 ierr = DMAdaptorSetSolver(adaptor, snes);CHKERRQ(ierr); 4371 ierr = DMAdaptorSetSequenceLength(adaptor, num);CHKERRQ(ierr); 4372 ierr = DMAdaptorSetFromOptions(adaptor);CHKERRQ(ierr); 4373 ierr = DMAdaptorSetUp(adaptor);CHKERRQ(ierr); 4374 ierr = DMAdaptorAdapt(adaptor, x, DM_ADAPTATION_SEQUENTIAL, &dm, &x);CHKERRQ(ierr); 4375 ierr = DMAdaptorDestroy(&adaptor);CHKERRQ(ierr); 4376 incall = PETSC_FALSE; 4377 } 4378 } 4379 } 4380 if (!x) { 4381 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 4382 ierr = DMCreateGlobalVector(dm,&xcreated);CHKERRQ(ierr); 4383 x = xcreated; 4384 } 4385 ierr = SNESViewFromOptions(snes,NULL,"-snes_view_pre");CHKERRQ(ierr); 4386 4387 for (grid=0; grid<snes->gridsequence; grid++) {ierr = PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)));CHKERRQ(ierr);} 4388 for (grid=0; grid<snes->gridsequence+1; grid++) { 4389 4390 /* set solution vector */ 4391 if (!grid) {ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);} 4392 ierr = VecDestroy(&snes->vec_sol);CHKERRQ(ierr); 4393 snes->vec_sol = x; 4394 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 4395 4396 /* set affine vector if provided */ 4397 if (b) { ierr = PetscObjectReference((PetscObject)b);CHKERRQ(ierr); } 4398 ierr = VecDestroy(&snes->vec_rhs);CHKERRQ(ierr); 4399 snes->vec_rhs = b; 4400 4401 if (snes->vec_rhs && (snes->vec_func == snes->vec_rhs)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Right hand side vector cannot be function vector"); 4402 if (snes->vec_func == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector"); 4403 if (snes->vec_rhs == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be right hand side vector"); 4404 if (!snes->vec_sol_update /* && snes->vec_sol */) { 4405 ierr = VecDuplicate(snes->vec_sol,&snes->vec_sol_update);CHKERRQ(ierr); 4406 ierr = PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->vec_sol_update);CHKERRQ(ierr); 4407 } 4408 ierr = DMShellSetGlobalVector(dm,snes->vec_sol);CHKERRQ(ierr); 4409 ierr = SNESSetUp(snes);CHKERRQ(ierr); 4410 4411 if (!grid) { 4412 if (snes->ops->computeinitialguess) { 4413 ierr = (*snes->ops->computeinitialguess)(snes,snes->vec_sol,snes->initialguessP);CHKERRQ(ierr); 4414 } 4415 } 4416 4417 if (snes->conv_hist_reset) snes->conv_hist_len = 0; 4418 if (snes->counters_reset) {snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;} 4419 4420 ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 4421 ierr = (*snes->ops->solve)(snes);CHKERRQ(ierr); 4422 ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 4423 if (!snes->reason) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason"); 4424 snes->domainerror = PETSC_FALSE; /* clear the flag if it has been set */ 4425 4426 if (snes->lagjac_persist) snes->jac_iter += snes->iter; 4427 if (snes->lagpre_persist) snes->pre_iter += snes->iter; 4428 4429 ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_test_local_min",NULL,NULL,&flg);CHKERRQ(ierr); 4430 if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); } 4431 ierr = SNESReasonViewFromOptions(snes);CHKERRQ(ierr); 4432 4433 if (snes->errorifnotconverged && snes->reason < 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_NOT_CONVERGED,"SNESSolve has not converged"); 4434 if (snes->reason < 0) break; 4435 if (grid < snes->gridsequence) { 4436 DM fine; 4437 Vec xnew; 4438 Mat interp; 4439 4440 ierr = DMRefine(snes->dm,PetscObjectComm((PetscObject)snes),&fine);CHKERRQ(ierr); 4441 if (!fine) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"DMRefine() did not perform any refinement, cannot continue grid sequencing"); 4442 ierr = DMCreateInterpolation(snes->dm,fine,&interp,NULL);CHKERRQ(ierr); 4443 ierr = DMCreateGlobalVector(fine,&xnew);CHKERRQ(ierr); 4444 ierr = MatInterpolate(interp,x,xnew);CHKERRQ(ierr); 4445 ierr = DMInterpolate(snes->dm,interp,fine);CHKERRQ(ierr); 4446 ierr = MatDestroy(&interp);CHKERRQ(ierr); 4447 x = xnew; 4448 4449 ierr = SNESReset(snes);CHKERRQ(ierr); 4450 ierr = SNESSetDM(snes,fine);CHKERRQ(ierr); 4451 ierr = SNESResetFromOptions(snes);CHKERRQ(ierr); 4452 ierr = DMDestroy(&fine);CHKERRQ(ierr); 4453 ierr = PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)));CHKERRQ(ierr); 4454 } 4455 } 4456 ierr = SNESViewFromOptions(snes,NULL,"-snes_view");CHKERRQ(ierr); 4457 ierr = VecViewFromOptions(snes->vec_sol,(PetscObject)snes,"-snes_view_solution");CHKERRQ(ierr); 4458 4459 ierr = VecDestroy(&xcreated);CHKERRQ(ierr); 4460 ierr = PetscObjectSAWsBlock((PetscObject)snes);CHKERRQ(ierr); 4461 PetscFunctionReturn(0); 4462 } 4463 4464 /* --------- Internal routines for SNES Package --------- */ 4465 4466 /*@C 4467 SNESSetType - Sets the method for the nonlinear solver. 4468 4469 Collective on SNES 4470 4471 Input Parameters: 4472 + snes - the SNES context 4473 - type - a known method 4474 4475 Options Database Key: 4476 . -snes_type <type> - Sets the method; use -help for a list 4477 of available methods (for instance, newtonls or newtontr) 4478 4479 Notes: 4480 See "petsc/include/petscsnes.h" for available methods (for instance) 4481 + SNESNEWTONLS - Newton's method with line search 4482 (systems of nonlinear equations) 4483 . SNESNEWTONTR - Newton's method with trust region 4484 (systems of nonlinear equations) 4485 4486 Normally, it is best to use the SNESSetFromOptions() command and then 4487 set the SNES solver type from the options database rather than by using 4488 this routine. Using the options database provides the user with 4489 maximum flexibility in evaluating the many nonlinear solvers. 4490 The SNESSetType() routine is provided for those situations where it 4491 is necessary to set the nonlinear solver independently of the command 4492 line or options database. This might be the case, for example, when 4493 the choice of solver changes during the execution of the program, 4494 and the user's application is taking responsibility for choosing the 4495 appropriate method. 4496 4497 Developer Notes: 4498 SNESRegister() adds a constructor for a new SNESType to SNESList, SNESSetType() locates 4499 the constructor in that list and calls it to create the spexific object. 4500 4501 Level: intermediate 4502 4503 .seealso: SNESType, SNESCreate(), SNESDestroy(), SNESGetType(), SNESSetFromOptions() 4504 4505 @*/ 4506 PetscErrorCode SNESSetType(SNES snes,SNESType type) 4507 { 4508 PetscErrorCode ierr,(*r)(SNES); 4509 PetscBool match; 4510 4511 PetscFunctionBegin; 4512 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 4513 PetscValidCharPointer(type,2); 4514 4515 ierr = PetscObjectTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr); 4516 if (match) PetscFunctionReturn(0); 4517 4518 ierr = PetscFunctionListFind(SNESList,type,&r);CHKERRQ(ierr); 4519 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type); 4520 /* Destroy the previous private SNES context */ 4521 if (snes->ops->destroy) { 4522 ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr); 4523 snes->ops->destroy = NULL; 4524 } 4525 /* Reinitialize function pointers in SNESOps structure */ 4526 snes->ops->setup = 0; 4527 snes->ops->solve = 0; 4528 snes->ops->view = 0; 4529 snes->ops->setfromoptions = 0; 4530 snes->ops->destroy = 0; 4531 ierr = SNESLineSearchDestroy(&snes->linesearch);CHKERRQ(ierr); 4532 /* Call the SNESCreate_XXX routine for this particular Nonlinear solver */ 4533 snes->setupcalled = PETSC_FALSE; 4534 4535 ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr); 4536 ierr = (*r)(snes);CHKERRQ(ierr); 4537 PetscFunctionReturn(0); 4538 } 4539 4540 /*@C 4541 SNESGetType - Gets the SNES method type and name (as a string). 4542 4543 Not Collective 4544 4545 Input Parameter: 4546 . snes - nonlinear solver context 4547 4548 Output Parameter: 4549 . type - SNES method (a character string) 4550 4551 Level: intermediate 4552 4553 @*/ 4554 PetscErrorCode SNESGetType(SNES snes,SNESType *type) 4555 { 4556 PetscFunctionBegin; 4557 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 4558 PetscValidPointer(type,2); 4559 *type = ((PetscObject)snes)->type_name; 4560 PetscFunctionReturn(0); 4561 } 4562 4563 /*@ 4564 SNESSetSolution - Sets the solution vector for use by the SNES routines. 4565 4566 Logically Collective on SNES and Vec 4567 4568 Input Parameters: 4569 + snes - the SNES context obtained from SNESCreate() 4570 - u - the solution vector 4571 4572 Level: beginner 4573 4574 @*/ 4575 PetscErrorCode SNESSetSolution(SNES snes, Vec u) 4576 { 4577 DM dm; 4578 PetscErrorCode ierr; 4579 4580 PetscFunctionBegin; 4581 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 4582 PetscValidHeaderSpecific(u, VEC_CLASSID, 2); 4583 ierr = PetscObjectReference((PetscObject) u);CHKERRQ(ierr); 4584 ierr = VecDestroy(&snes->vec_sol);CHKERRQ(ierr); 4585 4586 snes->vec_sol = u; 4587 4588 ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr); 4589 ierr = DMShellSetGlobalVector(dm, u);CHKERRQ(ierr); 4590 PetscFunctionReturn(0); 4591 } 4592 4593 /*@ 4594 SNESGetSolution - Returns the vector where the approximate solution is 4595 stored. This is the fine grid solution when using SNESSetGridSequence(). 4596 4597 Not Collective, but Vec is parallel if SNES is parallel 4598 4599 Input Parameter: 4600 . snes - the SNES context 4601 4602 Output Parameter: 4603 . x - the solution 4604 4605 Level: intermediate 4606 4607 .seealso: SNESGetSolutionUpdate(), SNESGetFunction() 4608 @*/ 4609 PetscErrorCode SNESGetSolution(SNES snes,Vec *x) 4610 { 4611 PetscFunctionBegin; 4612 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 4613 PetscValidPointer(x,2); 4614 *x = snes->vec_sol; 4615 PetscFunctionReturn(0); 4616 } 4617 4618 /*@ 4619 SNESGetSolutionUpdate - Returns the vector where the solution update is 4620 stored. 4621 4622 Not Collective, but Vec is parallel if SNES is parallel 4623 4624 Input Parameter: 4625 . snes - the SNES context 4626 4627 Output Parameter: 4628 . x - the solution update 4629 4630 Level: advanced 4631 4632 .seealso: SNESGetSolution(), SNESGetFunction() 4633 @*/ 4634 PetscErrorCode SNESGetSolutionUpdate(SNES snes,Vec *x) 4635 { 4636 PetscFunctionBegin; 4637 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 4638 PetscValidPointer(x,2); 4639 *x = snes->vec_sol_update; 4640 PetscFunctionReturn(0); 4641 } 4642 4643 /*@C 4644 SNESGetFunction - Returns the vector where the function is stored. 4645 4646 Not Collective, but Vec is parallel if SNES is parallel. Collective if Vec is requested, but has not been created yet. 4647 4648 Input Parameter: 4649 . snes - the SNES context 4650 4651 Output Parameter: 4652 + r - the vector that is used to store residuals (or NULL if you don't want it) 4653 . f - the function (or NULL if you don't want it); see SNESFunction for calling sequence details 4654 - ctx - the function context (or NULL if you don't want it) 4655 4656 Level: advanced 4657 4658 Notes: The vector r DOES NOT, in general contain the current value of the SNES nonlinear function 4659 4660 .seealso: SNESSetFunction(), SNESGetSolution(), SNESFunction 4661 @*/ 4662 PetscErrorCode SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx) 4663 { 4664 PetscErrorCode ierr; 4665 DM dm; 4666 4667 PetscFunctionBegin; 4668 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 4669 if (r) { 4670 if (!snes->vec_func) { 4671 if (snes->vec_rhs) { 4672 ierr = VecDuplicate(snes->vec_rhs,&snes->vec_func);CHKERRQ(ierr); 4673 } else if (snes->vec_sol) { 4674 ierr = VecDuplicate(snes->vec_sol,&snes->vec_func);CHKERRQ(ierr); 4675 } else if (snes->dm) { 4676 ierr = DMCreateGlobalVector(snes->dm,&snes->vec_func);CHKERRQ(ierr); 4677 } 4678 } 4679 *r = snes->vec_func; 4680 } 4681 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 4682 ierr = DMSNESGetFunction(dm,f,ctx);CHKERRQ(ierr); 4683 PetscFunctionReturn(0); 4684 } 4685 4686 /*@C 4687 SNESGetNGS - Returns the NGS function and context. 4688 4689 Input Parameter: 4690 . snes - the SNES context 4691 4692 Output Parameter: 4693 + f - the function (or NULL) see SNESNGSFunction for details 4694 - ctx - the function context (or NULL) 4695 4696 Level: advanced 4697 4698 .seealso: SNESSetNGS(), SNESGetFunction() 4699 @*/ 4700 4701 PetscErrorCode SNESGetNGS (SNES snes, PetscErrorCode (**f)(SNES, Vec, Vec, void*), void ** ctx) 4702 { 4703 PetscErrorCode ierr; 4704 DM dm; 4705 4706 PetscFunctionBegin; 4707 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 4708 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 4709 ierr = DMSNESGetNGS(dm,f,ctx);CHKERRQ(ierr); 4710 PetscFunctionReturn(0); 4711 } 4712 4713 /*@C 4714 SNESSetOptionsPrefix - Sets the prefix used for searching for all 4715 SNES options in the database. 4716 4717 Logically Collective on SNES 4718 4719 Input Parameter: 4720 + snes - the SNES context 4721 - prefix - the prefix to prepend to all option names 4722 4723 Notes: 4724 A hyphen (-) must NOT be given at the beginning of the prefix name. 4725 The first character of all runtime options is AUTOMATICALLY the hyphen. 4726 4727 Level: advanced 4728 4729 .seealso: SNESSetFromOptions() 4730 @*/ 4731 PetscErrorCode SNESSetOptionsPrefix(SNES snes,const char prefix[]) 4732 { 4733 PetscErrorCode ierr; 4734 4735 PetscFunctionBegin; 4736 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 4737 ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 4738 if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);} 4739 if (snes->linesearch) { 4740 ierr = SNESGetLineSearch(snes,&snes->linesearch);CHKERRQ(ierr); 4741 ierr = PetscObjectSetOptionsPrefix((PetscObject)snes->linesearch,prefix);CHKERRQ(ierr); 4742 } 4743 ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 4744 PetscFunctionReturn(0); 4745 } 4746 4747 /*@C 4748 SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 4749 SNES options in the database. 4750 4751 Logically Collective on SNES 4752 4753 Input Parameters: 4754 + snes - the SNES context 4755 - prefix - the prefix to prepend to all option names 4756 4757 Notes: 4758 A hyphen (-) must NOT be given at the beginning of the prefix name. 4759 The first character of all runtime options is AUTOMATICALLY the hyphen. 4760 4761 Level: advanced 4762 4763 .seealso: SNESGetOptionsPrefix() 4764 @*/ 4765 PetscErrorCode SNESAppendOptionsPrefix(SNES snes,const char prefix[]) 4766 { 4767 PetscErrorCode ierr; 4768 4769 PetscFunctionBegin; 4770 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 4771 ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 4772 if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);} 4773 if (snes->linesearch) { 4774 ierr = SNESGetLineSearch(snes,&snes->linesearch);CHKERRQ(ierr); 4775 ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes->linesearch,prefix);CHKERRQ(ierr); 4776 } 4777 ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 4778 PetscFunctionReturn(0); 4779 } 4780 4781 /*@C 4782 SNESGetOptionsPrefix - Sets the prefix used for searching for all 4783 SNES options in the database. 4784 4785 Not Collective 4786 4787 Input Parameter: 4788 . snes - the SNES context 4789 4790 Output Parameter: 4791 . prefix - pointer to the prefix string used 4792 4793 Notes: 4794 On the fortran side, the user should pass in a string 'prefix' of 4795 sufficient length to hold the prefix. 4796 4797 Level: advanced 4798 4799 .seealso: SNESAppendOptionsPrefix() 4800 @*/ 4801 PetscErrorCode SNESGetOptionsPrefix(SNES snes,const char *prefix[]) 4802 { 4803 PetscErrorCode ierr; 4804 4805 PetscFunctionBegin; 4806 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 4807 ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 4808 PetscFunctionReturn(0); 4809 } 4810 4811 4812 /*@C 4813 SNESRegister - Adds a method to the nonlinear solver package. 4814 4815 Not collective 4816 4817 Input Parameters: 4818 + name_solver - name of a new user-defined solver 4819 - routine_create - routine to create method context 4820 4821 Notes: 4822 SNESRegister() may be called multiple times to add several user-defined solvers. 4823 4824 Sample usage: 4825 .vb 4826 SNESRegister("my_solver",MySolverCreate); 4827 .ve 4828 4829 Then, your solver can be chosen with the procedural interface via 4830 $ SNESSetType(snes,"my_solver") 4831 or at runtime via the option 4832 $ -snes_type my_solver 4833 4834 Level: advanced 4835 4836 Note: If your function is not being put into a shared library then use SNESRegister() instead 4837 4838 .seealso: SNESRegisterAll(), SNESRegisterDestroy() 4839 4840 Level: advanced 4841 @*/ 4842 PetscErrorCode SNESRegister(const char sname[],PetscErrorCode (*function)(SNES)) 4843 { 4844 PetscErrorCode ierr; 4845 4846 PetscFunctionBegin; 4847 ierr = SNESInitializePackage();CHKERRQ(ierr); 4848 ierr = PetscFunctionListAdd(&SNESList,sname,function);CHKERRQ(ierr); 4849 PetscFunctionReturn(0); 4850 } 4851 4852 PetscErrorCode SNESTestLocalMin(SNES snes) 4853 { 4854 PetscErrorCode ierr; 4855 PetscInt N,i,j; 4856 Vec u,uh,fh; 4857 PetscScalar value; 4858 PetscReal norm; 4859 4860 PetscFunctionBegin; 4861 ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr); 4862 ierr = VecDuplicate(u,&uh);CHKERRQ(ierr); 4863 ierr = VecDuplicate(u,&fh);CHKERRQ(ierr); 4864 4865 /* currently only works for sequential */ 4866 ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");CHKERRQ(ierr); 4867 ierr = VecGetSize(u,&N);CHKERRQ(ierr); 4868 for (i=0; i<N; i++) { 4869 ierr = VecCopy(u,uh);CHKERRQ(ierr); 4870 ierr = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr); 4871 for (j=-10; j<11; j++) { 4872 value = PetscSign(j)*PetscExpReal(PetscAbs(j)-10.0); 4873 ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 4874 ierr = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr); 4875 ierr = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr); 4876 ierr = PetscPrintf(PETSC_COMM_WORLD," j norm %D %18.16e\n",j,norm);CHKERRQ(ierr); 4877 value = -value; 4878 ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 4879 } 4880 } 4881 ierr = VecDestroy(&uh);CHKERRQ(ierr); 4882 ierr = VecDestroy(&fh);CHKERRQ(ierr); 4883 PetscFunctionReturn(0); 4884 } 4885 4886 /*@ 4887 SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for 4888 computing relative tolerance for linear solvers within an inexact 4889 Newton method. 4890 4891 Logically Collective on SNES 4892 4893 Input Parameters: 4894 + snes - SNES context 4895 - flag - PETSC_TRUE or PETSC_FALSE 4896 4897 Options Database: 4898 + -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence 4899 . -snes_ksp_ew_version ver - version of Eisenstat-Walker method 4900 . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0 4901 . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax 4902 . -snes_ksp_ew_gamma <gamma> - Sets gamma 4903 . -snes_ksp_ew_alpha <alpha> - Sets alpha 4904 . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2 4905 - -snes_ksp_ew_threshold <threshold> - Sets threshold 4906 4907 Notes: 4908 Currently, the default is to use a constant relative tolerance for 4909 the inner linear solvers. Alternatively, one can use the 4910 Eisenstat-Walker method, where the relative convergence tolerance 4911 is reset at each Newton iteration according progress of the nonlinear 4912 solver. 4913 4914 Level: advanced 4915 4916 Reference: 4917 S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an 4918 inexact Newton method", SISC 17 (1), pp.16-32, 1996. 4919 4920 .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW() 4921 @*/ 4922 PetscErrorCode SNESKSPSetUseEW(SNES snes,PetscBool flag) 4923 { 4924 PetscFunctionBegin; 4925 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 4926 PetscValidLogicalCollectiveBool(snes,flag,2); 4927 snes->ksp_ewconv = flag; 4928 PetscFunctionReturn(0); 4929 } 4930 4931 /*@ 4932 SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method 4933 for computing relative tolerance for linear solvers within an 4934 inexact Newton method. 4935 4936 Not Collective 4937 4938 Input Parameter: 4939 . snes - SNES context 4940 4941 Output Parameter: 4942 . flag - PETSC_TRUE or PETSC_FALSE 4943 4944 Notes: 4945 Currently, the default is to use a constant relative tolerance for 4946 the inner linear solvers. Alternatively, one can use the 4947 Eisenstat-Walker method, where the relative convergence tolerance 4948 is reset at each Newton iteration according progress of the nonlinear 4949 solver. 4950 4951 Level: advanced 4952 4953 Reference: 4954 S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an 4955 inexact Newton method", SISC 17 (1), pp.16-32, 1996. 4956 4957 .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW() 4958 @*/ 4959 PetscErrorCode SNESKSPGetUseEW(SNES snes, PetscBool *flag) 4960 { 4961 PetscFunctionBegin; 4962 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 4963 PetscValidPointer(flag,2); 4964 *flag = snes->ksp_ewconv; 4965 PetscFunctionReturn(0); 4966 } 4967 4968 /*@ 4969 SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker 4970 convergence criteria for the linear solvers within an inexact 4971 Newton method. 4972 4973 Logically Collective on SNES 4974 4975 Input Parameters: 4976 + snes - SNES context 4977 . version - version 1, 2 (default is 2) or 3 4978 . rtol_0 - initial relative tolerance (0 <= rtol_0 < 1) 4979 . rtol_max - maximum relative tolerance (0 <= rtol_max < 1) 4980 . gamma - multiplicative factor for version 2 rtol computation 4981 (0 <= gamma2 <= 1) 4982 . alpha - power for version 2 rtol computation (1 < alpha <= 2) 4983 . alpha2 - power for safeguard 4984 - threshold - threshold for imposing safeguard (0 < threshold < 1) 4985 4986 Note: 4987 Version 3 was contributed by Luis Chacon, June 2006. 4988 4989 Use PETSC_DEFAULT to retain the default for any of the parameters. 4990 4991 Level: advanced 4992 4993 Reference: 4994 S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an 4995 inexact Newton method", Utah State University Math. Stat. Dept. Res. 4996 Report 6/94/75, June, 1994, to appear in SIAM J. Sci. Comput. 4997 4998 .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW() 4999 @*/ 5000 PetscErrorCode SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max,PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold) 5001 { 5002 SNESKSPEW *kctx; 5003 5004 PetscFunctionBegin; 5005 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 5006 kctx = (SNESKSPEW*)snes->kspconvctx; 5007 if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing"); 5008 PetscValidLogicalCollectiveInt(snes,version,2); 5009 PetscValidLogicalCollectiveReal(snes,rtol_0,3); 5010 PetscValidLogicalCollectiveReal(snes,rtol_max,4); 5011 PetscValidLogicalCollectiveReal(snes,gamma,5); 5012 PetscValidLogicalCollectiveReal(snes,alpha,6); 5013 PetscValidLogicalCollectiveReal(snes,alpha2,7); 5014 PetscValidLogicalCollectiveReal(snes,threshold,8); 5015 5016 if (version != PETSC_DEFAULT) kctx->version = version; 5017 if (rtol_0 != PETSC_DEFAULT) kctx->rtol_0 = rtol_0; 5018 if (rtol_max != PETSC_DEFAULT) kctx->rtol_max = rtol_max; 5019 if (gamma != PETSC_DEFAULT) kctx->gamma = gamma; 5020 if (alpha != PETSC_DEFAULT) kctx->alpha = alpha; 5021 if (alpha2 != PETSC_DEFAULT) kctx->alpha2 = alpha2; 5022 if (threshold != PETSC_DEFAULT) kctx->threshold = threshold; 5023 5024 if (kctx->version < 1 || kctx->version > 3) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 and 3 are supported: %D",kctx->version); 5025 if (kctx->rtol_0 < 0.0 || kctx->rtol_0 >= 1.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_0 < 1.0: %g",(double)kctx->rtol_0); 5026 if (kctx->rtol_max < 0.0 || kctx->rtol_max >= 1.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_max (%g) < 1.0\n",(double)kctx->rtol_max); 5027 if (kctx->gamma < 0.0 || kctx->gamma > 1.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= gamma (%g) <= 1.0\n",(double)kctx->gamma); 5028 if (kctx->alpha <= 1.0 || kctx->alpha > 2.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"1.0 < alpha (%g) <= 2.0\n",(double)kctx->alpha); 5029 if (kctx->threshold <= 0.0 || kctx->threshold >= 1.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 < threshold (%g) < 1.0\n",(double)kctx->threshold); 5030 PetscFunctionReturn(0); 5031 } 5032 5033 /*@ 5034 SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker 5035 convergence criteria for the linear solvers within an inexact 5036 Newton method. 5037 5038 Not Collective 5039 5040 Input Parameters: 5041 snes - SNES context 5042 5043 Output Parameters: 5044 + version - version 1, 2 (default is 2) or 3 5045 . rtol_0 - initial relative tolerance (0 <= rtol_0 < 1) 5046 . rtol_max - maximum relative tolerance (0 <= rtol_max < 1) 5047 . gamma - multiplicative factor for version 2 rtol computation (0 <= gamma2 <= 1) 5048 . alpha - power for version 2 rtol computation (1 < alpha <= 2) 5049 . alpha2 - power for safeguard 5050 - threshold - threshold for imposing safeguard (0 < threshold < 1) 5051 5052 Level: advanced 5053 5054 .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW() 5055 @*/ 5056 PetscErrorCode SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max,PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold) 5057 { 5058 SNESKSPEW *kctx; 5059 5060 PetscFunctionBegin; 5061 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 5062 kctx = (SNESKSPEW*)snes->kspconvctx; 5063 if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing"); 5064 if (version) *version = kctx->version; 5065 if (rtol_0) *rtol_0 = kctx->rtol_0; 5066 if (rtol_max) *rtol_max = kctx->rtol_max; 5067 if (gamma) *gamma = kctx->gamma; 5068 if (alpha) *alpha = kctx->alpha; 5069 if (alpha2) *alpha2 = kctx->alpha2; 5070 if (threshold) *threshold = kctx->threshold; 5071 PetscFunctionReturn(0); 5072 } 5073 5074 PetscErrorCode KSPPreSolve_SNESEW(KSP ksp, Vec b, Vec x, SNES snes) 5075 { 5076 PetscErrorCode ierr; 5077 SNESKSPEW *kctx = (SNESKSPEW*)snes->kspconvctx; 5078 PetscReal rtol = PETSC_DEFAULT,stol; 5079 5080 PetscFunctionBegin; 5081 if (!snes->ksp_ewconv) PetscFunctionReturn(0); 5082 if (!snes->iter) { 5083 rtol = kctx->rtol_0; /* first time in, so use the original user rtol */ 5084 ierr = VecNorm(snes->vec_func,NORM_2,&kctx->norm_first);CHKERRQ(ierr); 5085 } 5086 else { 5087 if (kctx->version == 1) { 5088 rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last; 5089 if (rtol < 0.0) rtol = -rtol; 5090 stol = PetscPowReal(kctx->rtol_last,kctx->alpha2); 5091 if (stol > kctx->threshold) rtol = PetscMax(rtol,stol); 5092 } else if (kctx->version == 2) { 5093 rtol = kctx->gamma * PetscPowReal(snes->norm/kctx->norm_last,kctx->alpha); 5094 stol = kctx->gamma * PetscPowReal(kctx->rtol_last,kctx->alpha); 5095 if (stol > kctx->threshold) rtol = PetscMax(rtol,stol); 5096 } else if (kctx->version == 3) { /* contributed by Luis Chacon, June 2006. */ 5097 rtol = kctx->gamma * PetscPowReal(snes->norm/kctx->norm_last,kctx->alpha); 5098 /* safeguard: avoid sharp decrease of rtol */ 5099 stol = kctx->gamma*PetscPowReal(kctx->rtol_last,kctx->alpha); 5100 stol = PetscMax(rtol,stol); 5101 rtol = PetscMin(kctx->rtol_0,stol); 5102 /* safeguard: avoid oversolving */ 5103 stol = kctx->gamma*(kctx->norm_first*snes->rtol)/snes->norm; 5104 stol = PetscMax(rtol,stol); 5105 rtol = PetscMin(kctx->rtol_0,stol); 5106 } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 or 3 are supported: %D",kctx->version); 5107 } 5108 /* safeguard: avoid rtol greater than one */ 5109 rtol = PetscMin(rtol,kctx->rtol_max); 5110 ierr = KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 5111 ierr = PetscInfo3(snes,"iter %D, Eisenstat-Walker (version %D) KSP rtol=%g\n",snes->iter,kctx->version,(double)rtol);CHKERRQ(ierr); 5112 PetscFunctionReturn(0); 5113 } 5114 5115 PetscErrorCode KSPPostSolve_SNESEW(KSP ksp, Vec b, Vec x, SNES snes) 5116 { 5117 PetscErrorCode ierr; 5118 SNESKSPEW *kctx = (SNESKSPEW*)snes->kspconvctx; 5119 PCSide pcside; 5120 Vec lres; 5121 5122 PetscFunctionBegin; 5123 if (!snes->ksp_ewconv) PetscFunctionReturn(0); 5124 ierr = KSPGetTolerances(ksp,&kctx->rtol_last,0,0,0);CHKERRQ(ierr); 5125 kctx->norm_last = snes->norm; 5126 if (kctx->version == 1) { 5127 PC pc; 5128 PetscBool isNone; 5129 5130 ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 5131 ierr = PetscObjectTypeCompare((PetscObject) pc, PCNONE, &isNone);CHKERRQ(ierr); 5132 ierr = KSPGetPCSide(ksp,&pcside);CHKERRQ(ierr); 5133 if (pcside == PC_RIGHT || isNone) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */ 5134 /* KSP residual is true linear residual */ 5135 ierr = KSPGetResidualNorm(ksp,&kctx->lresid_last);CHKERRQ(ierr); 5136 } else { 5137 /* KSP residual is preconditioned residual */ 5138 /* compute true linear residual norm */ 5139 ierr = VecDuplicate(b,&lres);CHKERRQ(ierr); 5140 ierr = MatMult(snes->jacobian,x,lres);CHKERRQ(ierr); 5141 ierr = VecAYPX(lres,-1.0,b);CHKERRQ(ierr); 5142 ierr = VecNorm(lres,NORM_2,&kctx->lresid_last);CHKERRQ(ierr); 5143 ierr = VecDestroy(&lres);CHKERRQ(ierr); 5144 } 5145 } 5146 PetscFunctionReturn(0); 5147 } 5148 5149 /*@ 5150 SNESGetKSP - Returns the KSP context for a SNES solver. 5151 5152 Not Collective, but if SNES object is parallel, then KSP object is parallel 5153 5154 Input Parameter: 5155 . snes - the SNES context 5156 5157 Output Parameter: 5158 . ksp - the KSP context 5159 5160 Notes: 5161 The user can then directly manipulate the KSP context to set various 5162 options, etc. Likewise, the user can then extract and manipulate the 5163 PC contexts as well. 5164 5165 Level: beginner 5166 5167 .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP() 5168 @*/ 5169 PetscErrorCode SNESGetKSP(SNES snes,KSP *ksp) 5170 { 5171 PetscErrorCode ierr; 5172 5173 PetscFunctionBegin; 5174 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 5175 PetscValidPointer(ksp,2); 5176 5177 if (!snes->ksp) { 5178 PetscBool monitor = PETSC_FALSE; 5179 5180 ierr = KSPCreate(PetscObjectComm((PetscObject)snes),&snes->ksp);CHKERRQ(ierr); 5181 ierr = PetscObjectIncrementTabLevel((PetscObject)snes->ksp,(PetscObject)snes,1);CHKERRQ(ierr); 5182 ierr = PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->ksp);CHKERRQ(ierr); 5183 5184 ierr = KSPSetPreSolve(snes->ksp,(PetscErrorCode (*)(KSP,Vec,Vec,void*))KSPPreSolve_SNESEW,snes);CHKERRQ(ierr); 5185 ierr = KSPSetPostSolve(snes->ksp,(PetscErrorCode (*)(KSP,Vec,Vec,void*))KSPPostSolve_SNESEW,snes);CHKERRQ(ierr); 5186 5187 ierr = PetscOptionsGetBool(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-ksp_monitor_snes",&monitor,NULL);CHKERRQ(ierr); 5188 if (monitor) { 5189 ierr = KSPMonitorSet(snes->ksp,KSPMonitorSNES,snes,NULL);CHKERRQ(ierr); 5190 } 5191 monitor = PETSC_FALSE; 5192 ierr = PetscOptionsGetBool(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-ksp_monitor_snes_lg",&monitor,NULL);CHKERRQ(ierr); 5193 if (monitor) { 5194 PetscObject *objs; 5195 ierr = KSPMonitorSNESLGResidualNormCreate(PetscObjectComm((PetscObject)snes),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,600,600,&objs);CHKERRQ(ierr); 5196 objs[0] = (PetscObject) snes; 5197 ierr = KSPMonitorSet(snes->ksp,(PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*))KSPMonitorSNESLGResidualNorm,objs,(PetscErrorCode (*)(void**))KSPMonitorSNESLGResidualNormDestroy);CHKERRQ(ierr); 5198 } 5199 ierr = PetscObjectSetOptions((PetscObject)snes->ksp,((PetscObject)snes)->options);CHKERRQ(ierr); 5200 } 5201 *ksp = snes->ksp; 5202 PetscFunctionReturn(0); 5203 } 5204 5205 5206 #include <petsc/private/dmimpl.h> 5207 /*@ 5208 SNESSetDM - Sets the DM that may be used by some nonlinear solvers or their underlying preconditioners 5209 5210 Logically Collective on SNES 5211 5212 Input Parameters: 5213 + snes - the nonlinear solver context 5214 - dm - the dm, cannot be NULL 5215 5216 Notes: 5217 A DM can only be used for solving one problem at a time because information about the problem is stored on the DM, 5218 even when not using interfaces like DMSNESSetFunction(). Use DMClone() to get a distinct DM when solving different 5219 problems using the same function space. 5220 5221 Level: intermediate 5222 5223 .seealso: SNESGetDM(), KSPSetDM(), KSPGetDM() 5224 @*/ 5225 PetscErrorCode SNESSetDM(SNES snes,DM dm) 5226 { 5227 PetscErrorCode ierr; 5228 KSP ksp; 5229 DMSNES sdm; 5230 5231 PetscFunctionBegin; 5232 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 5233 PetscValidHeaderSpecific(dm,DM_CLASSID,2); 5234 ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr); 5235 if (snes->dm) { /* Move the DMSNES context over to the new DM unless the new DM already has one */ 5236 if (snes->dm->dmsnes && !dm->dmsnes) { 5237 ierr = DMCopyDMSNES(snes->dm,dm);CHKERRQ(ierr); 5238 ierr = DMGetDMSNES(snes->dm,&sdm);CHKERRQ(ierr); 5239 if (sdm->originaldm == snes->dm) sdm->originaldm = dm; /* Grant write privileges to the replacement DM */ 5240 } 5241 ierr = DMCoarsenHookRemove(snes->dm,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,snes);CHKERRQ(ierr); 5242 ierr = DMDestroy(&snes->dm);CHKERRQ(ierr); 5243 } 5244 snes->dm = dm; 5245 snes->dmAuto = PETSC_FALSE; 5246 5247 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 5248 ierr = KSPSetDM(ksp,dm);CHKERRQ(ierr); 5249 ierr = KSPSetDMActive(ksp,PETSC_FALSE);CHKERRQ(ierr); 5250 if (snes->npc) { 5251 ierr = SNESSetDM(snes->npc, snes->dm);CHKERRQ(ierr); 5252 ierr = SNESSetNPCSide(snes,snes->npcside);CHKERRQ(ierr); 5253 } 5254 PetscFunctionReturn(0); 5255 } 5256 5257 /*@ 5258 SNESGetDM - Gets the DM that may be used by some preconditioners 5259 5260 Not Collective but DM obtained is parallel on SNES 5261 5262 Input Parameter: 5263 . snes - the preconditioner context 5264 5265 Output Parameter: 5266 . dm - the dm 5267 5268 Level: intermediate 5269 5270 .seealso: SNESSetDM(), KSPSetDM(), KSPGetDM() 5271 @*/ 5272 PetscErrorCode SNESGetDM(SNES snes,DM *dm) 5273 { 5274 PetscErrorCode ierr; 5275 5276 PetscFunctionBegin; 5277 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 5278 if (!snes->dm) { 5279 ierr = DMShellCreate(PetscObjectComm((PetscObject)snes),&snes->dm);CHKERRQ(ierr); 5280 snes->dmAuto = PETSC_TRUE; 5281 } 5282 *dm = snes->dm; 5283 PetscFunctionReturn(0); 5284 } 5285 5286 /*@ 5287 SNESSetNPC - Sets the nonlinear preconditioner to be used. 5288 5289 Collective on SNES 5290 5291 Input Parameters: 5292 + snes - iterative context obtained from SNESCreate() 5293 - pc - the preconditioner object 5294 5295 Notes: 5296 Use SNESGetNPC() to retrieve the preconditioner context (for example, 5297 to configure it using the API). 5298 5299 Level: developer 5300 5301 .seealso: SNESGetNPC(), SNESHasNPC() 5302 @*/ 5303 PetscErrorCode SNESSetNPC(SNES snes, SNES pc) 5304 { 5305 PetscErrorCode ierr; 5306 5307 PetscFunctionBegin; 5308 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 5309 PetscValidHeaderSpecific(pc, SNES_CLASSID, 2); 5310 PetscCheckSameComm(snes, 1, pc, 2); 5311 ierr = PetscObjectReference((PetscObject) pc);CHKERRQ(ierr); 5312 ierr = SNESDestroy(&snes->npc);CHKERRQ(ierr); 5313 snes->npc = pc; 5314 ierr = PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->npc);CHKERRQ(ierr); 5315 PetscFunctionReturn(0); 5316 } 5317 5318 /*@ 5319 SNESGetNPC - Creates a nonlinear preconditioning solver (SNES) to be used to precondition the nonlinear solver. 5320 5321 Not Collective; but any changes to the obtained SNES object must be applied collectively 5322 5323 Input Parameter: 5324 . snes - iterative context obtained from SNESCreate() 5325 5326 Output Parameter: 5327 . pc - preconditioner context 5328 5329 Notes: 5330 If a SNES was previously set with SNESSetNPC() then that SNES is returned. 5331 5332 The (preconditioner) SNES returned automatically inherits the same nonlinear function and Jacobian supplied to the original 5333 SNES during SNESSetUp() 5334 5335 Level: developer 5336 5337 .seealso: SNESSetNPC(), SNESHasNPC(), SNES, SNESCreate() 5338 @*/ 5339 PetscErrorCode SNESGetNPC(SNES snes, SNES *pc) 5340 { 5341 PetscErrorCode ierr; 5342 const char *optionsprefix; 5343 5344 PetscFunctionBegin; 5345 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 5346 PetscValidPointer(pc, 2); 5347 if (!snes->npc) { 5348 ierr = SNESCreate(PetscObjectComm((PetscObject)snes),&snes->npc);CHKERRQ(ierr); 5349 ierr = PetscObjectIncrementTabLevel((PetscObject)snes->npc,(PetscObject)snes,1);CHKERRQ(ierr); 5350 ierr = PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->npc);CHKERRQ(ierr); 5351 ierr = SNESGetOptionsPrefix(snes,&optionsprefix);CHKERRQ(ierr); 5352 ierr = SNESSetOptionsPrefix(snes->npc,optionsprefix);CHKERRQ(ierr); 5353 ierr = SNESAppendOptionsPrefix(snes->npc,"npc_");CHKERRQ(ierr); 5354 ierr = SNESSetCountersReset(snes->npc,PETSC_FALSE);CHKERRQ(ierr); 5355 } 5356 *pc = snes->npc; 5357 PetscFunctionReturn(0); 5358 } 5359 5360 /*@ 5361 SNESHasNPC - Returns whether a nonlinear preconditioner exists 5362 5363 Not Collective 5364 5365 Input Parameter: 5366 . snes - iterative context obtained from SNESCreate() 5367 5368 Output Parameter: 5369 . has_npc - whether the SNES has an NPC or not 5370 5371 Level: developer 5372 5373 .seealso: SNESSetNPC(), SNESGetNPC() 5374 @*/ 5375 PetscErrorCode SNESHasNPC(SNES snes, PetscBool *has_npc) 5376 { 5377 PetscFunctionBegin; 5378 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 5379 *has_npc = (PetscBool) (snes->npc ? PETSC_TRUE : PETSC_FALSE); 5380 PetscFunctionReturn(0); 5381 } 5382 5383 /*@ 5384 SNESSetNPCSide - Sets the preconditioning side. 5385 5386 Logically Collective on SNES 5387 5388 Input Parameter: 5389 . snes - iterative context obtained from SNESCreate() 5390 5391 Output Parameter: 5392 . side - the preconditioning side, where side is one of 5393 .vb 5394 PC_LEFT - left preconditioning 5395 PC_RIGHT - right preconditioning (default for most nonlinear solvers) 5396 .ve 5397 5398 Options Database Keys: 5399 . -snes_pc_side <right,left> 5400 5401 Notes: 5402 SNESNRICHARDSON and SNESNCG only support left preconditioning. 5403 5404 Level: intermediate 5405 5406 .seealso: SNESGetNPCSide(), KSPSetPCSide() 5407 @*/ 5408 PetscErrorCode SNESSetNPCSide(SNES snes,PCSide side) 5409 { 5410 PetscFunctionBegin; 5411 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 5412 PetscValidLogicalCollectiveEnum(snes,side,2); 5413 snes->npcside= side; 5414 PetscFunctionReturn(0); 5415 } 5416 5417 /*@ 5418 SNESGetNPCSide - Gets the preconditioning side. 5419 5420 Not Collective 5421 5422 Input Parameter: 5423 . snes - iterative context obtained from SNESCreate() 5424 5425 Output Parameter: 5426 . side - the preconditioning side, where side is one of 5427 .vb 5428 PC_LEFT - left preconditioning 5429 PC_RIGHT - right preconditioning (default for most nonlinear solvers) 5430 .ve 5431 5432 Level: intermediate 5433 5434 .seealso: SNESSetNPCSide(), KSPGetPCSide() 5435 @*/ 5436 PetscErrorCode SNESGetNPCSide(SNES snes,PCSide *side) 5437 { 5438 PetscFunctionBegin; 5439 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 5440 PetscValidPointer(side,2); 5441 *side = snes->npcside; 5442 PetscFunctionReturn(0); 5443 } 5444 5445 /*@ 5446 SNESSetLineSearch - Sets the linesearch on the SNES instance. 5447 5448 Collective on SNES 5449 5450 Input Parameters: 5451 + snes - iterative context obtained from SNESCreate() 5452 - linesearch - the linesearch object 5453 5454 Notes: 5455 Use SNESGetLineSearch() to retrieve the preconditioner context (for example, 5456 to configure it using the API). 5457 5458 Level: developer 5459 5460 .seealso: SNESGetLineSearch() 5461 @*/ 5462 PetscErrorCode SNESSetLineSearch(SNES snes, SNESLineSearch linesearch) 5463 { 5464 PetscErrorCode ierr; 5465 5466 PetscFunctionBegin; 5467 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 5468 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 2); 5469 PetscCheckSameComm(snes, 1, linesearch, 2); 5470 ierr = PetscObjectReference((PetscObject) linesearch);CHKERRQ(ierr); 5471 ierr = SNESLineSearchDestroy(&snes->linesearch);CHKERRQ(ierr); 5472 5473 snes->linesearch = linesearch; 5474 5475 ierr = PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->linesearch);CHKERRQ(ierr); 5476 PetscFunctionReturn(0); 5477 } 5478 5479 /*@ 5480 SNESGetLineSearch - Returns a pointer to the line search context set with SNESSetLineSearch() 5481 or creates a default line search instance associated with the SNES and returns it. 5482 5483 Not Collective 5484 5485 Input Parameter: 5486 . snes - iterative context obtained from SNESCreate() 5487 5488 Output Parameter: 5489 . linesearch - linesearch context 5490 5491 Level: beginner 5492 5493 .seealso: SNESSetLineSearch(), SNESLineSearchCreate() 5494 @*/ 5495 PetscErrorCode SNESGetLineSearch(SNES snes, SNESLineSearch *linesearch) 5496 { 5497 PetscErrorCode ierr; 5498 const char *optionsprefix; 5499 5500 PetscFunctionBegin; 5501 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 5502 PetscValidPointer(linesearch, 2); 5503 if (!snes->linesearch) { 5504 ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 5505 ierr = SNESLineSearchCreate(PetscObjectComm((PetscObject)snes), &snes->linesearch);CHKERRQ(ierr); 5506 ierr = SNESLineSearchSetSNES(snes->linesearch, snes);CHKERRQ(ierr); 5507 ierr = SNESLineSearchAppendOptionsPrefix(snes->linesearch, optionsprefix);CHKERRQ(ierr); 5508 ierr = PetscObjectIncrementTabLevel((PetscObject) snes->linesearch, (PetscObject) snes, 1);CHKERRQ(ierr); 5509 ierr = PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->linesearch);CHKERRQ(ierr); 5510 } 5511 *linesearch = snes->linesearch; 5512 PetscFunctionReturn(0); 5513 } 5514 5515 #if defined(PETSC_HAVE_MATLAB_ENGINE) 5516 #include <mex.h> 5517 5518 typedef struct {char *funcname; mxArray *ctx;} SNESMatlabContext; 5519 5520 /* 5521 SNESComputeFunction_Matlab - Calls the function that has been set with SNESSetFunctionMatlab(). 5522 5523 Collective on SNES 5524 5525 Input Parameters: 5526 + snes - the SNES context 5527 - x - input vector 5528 5529 Output Parameter: 5530 . y - function vector, as set by SNESSetFunction() 5531 5532 Notes: 5533 SNESComputeFunction() is typically used within nonlinear solvers 5534 implementations, so most users would not generally call this routine 5535 themselves. 5536 5537 Level: developer 5538 5539 .seealso: SNESSetFunction(), SNESGetFunction() 5540 */ 5541 PetscErrorCode SNESComputeFunction_Matlab(SNES snes,Vec x,Vec y, void *ctx) 5542 { 5543 PetscErrorCode ierr; 5544 SNESMatlabContext *sctx = (SNESMatlabContext*)ctx; 5545 int nlhs = 1,nrhs = 5; 5546 mxArray *plhs[1],*prhs[5]; 5547 long long int lx = 0,ly = 0,ls = 0; 5548 5549 PetscFunctionBegin; 5550 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 5551 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 5552 PetscValidHeaderSpecific(y,VEC_CLASSID,3); 5553 PetscCheckSameComm(snes,1,x,2); 5554 PetscCheckSameComm(snes,1,y,3); 5555 5556 /* call Matlab function in ctx with arguments x and y */ 5557 5558 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 5559 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 5560 ierr = PetscMemcpy(&ly,&y,sizeof(x));CHKERRQ(ierr); 5561 prhs[0] = mxCreateDoubleScalar((double)ls); 5562 prhs[1] = mxCreateDoubleScalar((double)lx); 5563 prhs[2] = mxCreateDoubleScalar((double)ly); 5564 prhs[3] = mxCreateString(sctx->funcname); 5565 prhs[4] = sctx->ctx; 5566 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeFunctionInternal");CHKERRQ(ierr); 5567 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 5568 mxDestroyArray(prhs[0]); 5569 mxDestroyArray(prhs[1]); 5570 mxDestroyArray(prhs[2]); 5571 mxDestroyArray(prhs[3]); 5572 mxDestroyArray(plhs[0]); 5573 PetscFunctionReturn(0); 5574 } 5575 5576 /* 5577 SNESSetFunctionMatlab - Sets the function evaluation routine and function 5578 vector for use by the SNES routines in solving systems of nonlinear 5579 equations from MATLAB. Here the function is a string containing the name of a MATLAB function 5580 5581 Logically Collective on SNES 5582 5583 Input Parameters: 5584 + snes - the SNES context 5585 . r - vector to store function value 5586 - f - function evaluation routine 5587 5588 Notes: 5589 The Newton-like methods typically solve linear systems of the form 5590 $ f'(x) x = -f(x), 5591 where f'(x) denotes the Jacobian matrix and f(x) is the function. 5592 5593 Level: beginner 5594 5595 Developer Note: This bleeds the allocated memory SNESMatlabContext *sctx; 5596 5597 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction() 5598 */ 5599 PetscErrorCode SNESSetFunctionMatlab(SNES snes,Vec r,const char *f,mxArray *ctx) 5600 { 5601 PetscErrorCode ierr; 5602 SNESMatlabContext *sctx; 5603 5604 PetscFunctionBegin; 5605 /* currently sctx is memory bleed */ 5606 ierr = PetscNew(&sctx);CHKERRQ(ierr); 5607 ierr = PetscStrallocpy(f,&sctx->funcname);CHKERRQ(ierr); 5608 /* 5609 This should work, but it doesn't 5610 sctx->ctx = ctx; 5611 mexMakeArrayPersistent(sctx->ctx); 5612 */ 5613 sctx->ctx = mxDuplicateArray(ctx); 5614 ierr = SNESSetFunction(snes,r,SNESComputeFunction_Matlab,sctx);CHKERRQ(ierr); 5615 PetscFunctionReturn(0); 5616 } 5617 5618 /* 5619 SNESComputeJacobian_Matlab - Calls the function that has been set with SNESSetJacobianMatlab(). 5620 5621 Collective on SNES 5622 5623 Input Parameters: 5624 + snes - the SNES context 5625 . x - input vector 5626 . A, B - the matrices 5627 - ctx - user context 5628 5629 Level: developer 5630 5631 .seealso: SNESSetFunction(), SNESGetFunction() 5632 @*/ 5633 PetscErrorCode SNESComputeJacobian_Matlab(SNES snes,Vec x,Mat A,Mat B,void *ctx) 5634 { 5635 PetscErrorCode ierr; 5636 SNESMatlabContext *sctx = (SNESMatlabContext*)ctx; 5637 int nlhs = 2,nrhs = 6; 5638 mxArray *plhs[2],*prhs[6]; 5639 long long int lx = 0,lA = 0,ls = 0, lB = 0; 5640 5641 PetscFunctionBegin; 5642 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 5643 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 5644 5645 /* call Matlab function in ctx with arguments x and y */ 5646 5647 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 5648 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 5649 ierr = PetscMemcpy(&lA,A,sizeof(x));CHKERRQ(ierr); 5650 ierr = PetscMemcpy(&lB,B,sizeof(x));CHKERRQ(ierr); 5651 prhs[0] = mxCreateDoubleScalar((double)ls); 5652 prhs[1] = mxCreateDoubleScalar((double)lx); 5653 prhs[2] = mxCreateDoubleScalar((double)lA); 5654 prhs[3] = mxCreateDoubleScalar((double)lB); 5655 prhs[4] = mxCreateString(sctx->funcname); 5656 prhs[5] = sctx->ctx; 5657 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeJacobianInternal");CHKERRQ(ierr); 5658 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 5659 mxDestroyArray(prhs[0]); 5660 mxDestroyArray(prhs[1]); 5661 mxDestroyArray(prhs[2]); 5662 mxDestroyArray(prhs[3]); 5663 mxDestroyArray(prhs[4]); 5664 mxDestroyArray(plhs[0]); 5665 mxDestroyArray(plhs[1]); 5666 PetscFunctionReturn(0); 5667 } 5668 5669 /* 5670 SNESSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices 5671 vector for use by the SNES routines in solving systems of nonlinear 5672 equations from MATLAB. Here the function is a string containing the name of a MATLAB function 5673 5674 Logically Collective on SNES 5675 5676 Input Parameters: 5677 + snes - the SNES context 5678 . A,B - Jacobian matrices 5679 . J - function evaluation routine 5680 - ctx - user context 5681 5682 Level: developer 5683 5684 Developer Note: This bleeds the allocated memory SNESMatlabContext *sctx; 5685 5686 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction(), J 5687 */ 5688 PetscErrorCode SNESSetJacobianMatlab(SNES snes,Mat A,Mat B,const char *J,mxArray *ctx) 5689 { 5690 PetscErrorCode ierr; 5691 SNESMatlabContext *sctx; 5692 5693 PetscFunctionBegin; 5694 /* currently sctx is memory bleed */ 5695 ierr = PetscNew(&sctx);CHKERRQ(ierr); 5696 ierr = PetscStrallocpy(J,&sctx->funcname);CHKERRQ(ierr); 5697 /* 5698 This should work, but it doesn't 5699 sctx->ctx = ctx; 5700 mexMakeArrayPersistent(sctx->ctx); 5701 */ 5702 sctx->ctx = mxDuplicateArray(ctx); 5703 ierr = SNESSetJacobian(snes,A,B,SNESComputeJacobian_Matlab,sctx);CHKERRQ(ierr); 5704 PetscFunctionReturn(0); 5705 } 5706 5707 /* 5708 SNESMonitor_Matlab - Calls the function that has been set with SNESMonitorSetMatlab(). 5709 5710 Collective on SNES 5711 5712 .seealso: SNESSetFunction(), SNESGetFunction() 5713 @*/ 5714 PetscErrorCode SNESMonitor_Matlab(SNES snes,PetscInt it, PetscReal fnorm, void *ctx) 5715 { 5716 PetscErrorCode ierr; 5717 SNESMatlabContext *sctx = (SNESMatlabContext*)ctx; 5718 int nlhs = 1,nrhs = 6; 5719 mxArray *plhs[1],*prhs[6]; 5720 long long int lx = 0,ls = 0; 5721 Vec x = snes->vec_sol; 5722 5723 PetscFunctionBegin; 5724 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 5725 5726 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 5727 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 5728 prhs[0] = mxCreateDoubleScalar((double)ls); 5729 prhs[1] = mxCreateDoubleScalar((double)it); 5730 prhs[2] = mxCreateDoubleScalar((double)fnorm); 5731 prhs[3] = mxCreateDoubleScalar((double)lx); 5732 prhs[4] = mxCreateString(sctx->funcname); 5733 prhs[5] = sctx->ctx; 5734 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESMonitorInternal");CHKERRQ(ierr); 5735 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 5736 mxDestroyArray(prhs[0]); 5737 mxDestroyArray(prhs[1]); 5738 mxDestroyArray(prhs[2]); 5739 mxDestroyArray(prhs[3]); 5740 mxDestroyArray(prhs[4]); 5741 mxDestroyArray(plhs[0]); 5742 PetscFunctionReturn(0); 5743 } 5744 5745 /* 5746 SNESMonitorSetMatlab - Sets the monitor function from MATLAB 5747 5748 Level: developer 5749 5750 Developer Note: This bleeds the allocated memory SNESMatlabContext *sctx; 5751 5752 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction() 5753 */ 5754 PetscErrorCode SNESMonitorSetMatlab(SNES snes,const char *f,mxArray *ctx) 5755 { 5756 PetscErrorCode ierr; 5757 SNESMatlabContext *sctx; 5758 5759 PetscFunctionBegin; 5760 /* currently sctx is memory bleed */ 5761 ierr = PetscNew(&sctx);CHKERRQ(ierr); 5762 ierr = PetscStrallocpy(f,&sctx->funcname);CHKERRQ(ierr); 5763 /* 5764 This should work, but it doesn't 5765 sctx->ctx = ctx; 5766 mexMakeArrayPersistent(sctx->ctx); 5767 */ 5768 sctx->ctx = mxDuplicateArray(ctx); 5769 ierr = SNESMonitorSet(snes,SNESMonitor_Matlab,sctx,NULL);CHKERRQ(ierr); 5770 PetscFunctionReturn(0); 5771 } 5772 5773 #endif 5774