1 #define PETSCSNES_DLL 2 3 #include "src/snes/snesimpl.h" /*I "petscsnes.h" I*/ 4 5 PetscTruth SNESRegisterAllCalled = PETSC_FALSE; 6 PetscFList SNESList = PETSC_NULL; 7 8 /* Logging support */ 9 PetscCookie PETSCSNES_DLLEXPORT SNES_COOKIE = 0; 10 PetscEvent SNES_Solve = 0, SNES_LineSearch = 0, SNES_FunctionEval = 0, SNES_JacobianEval = 0; 11 12 #undef __FUNCT__ 13 #define __FUNCT__ "SNESView" 14 /*@C 15 SNESView - Prints the SNES data structure. 16 17 Collective on SNES 18 19 Input Parameters: 20 + SNES - the SNES context 21 - viewer - visualization context 22 23 Options Database Key: 24 . -snes_view - Calls SNESView() at end of SNESSolve() 25 26 Notes: 27 The available visualization contexts include 28 + PETSC_VIEWER_STDOUT_SELF - standard output (default) 29 - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 30 output where only the first processor opens 31 the file. All other processors send their 32 data to the first processor to print. 33 34 The user can open an alternative visualization context with 35 PetscViewerASCIIOpen() - output to a specified file. 36 37 Level: beginner 38 39 .keywords: SNES, view 40 41 .seealso: PetscViewerASCIIOpen() 42 @*/ 43 PetscErrorCode PETSCSNES_DLLEXPORT SNESView(SNES snes,PetscViewer viewer) 44 { 45 SNES_KSP_EW_ConvCtx *kctx; 46 PetscErrorCode ierr; 47 KSP ksp; 48 SNESType type; 49 PetscTruth iascii,isstring; 50 51 PetscFunctionBegin; 52 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 53 if (!viewer) viewer = PETSC_VIEWER_STDOUT_(snes->comm); 54 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2); 55 PetscCheckSameComm(snes,1,viewer,2); 56 57 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 58 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 59 if (iascii) { 60 if (snes->prefix) { 61 ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:(%s)\n",snes->prefix);CHKERRQ(ierr); 62 } else { 63 ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:\n");CHKERRQ(ierr); 64 } 65 ierr = SNESGetType(snes,&type);CHKERRQ(ierr); 66 if (type) { 67 ierr = PetscViewerASCIIPrintf(viewer," type: %s\n",type);CHKERRQ(ierr); 68 } else { 69 ierr = PetscViewerASCIIPrintf(viewer," type: not set yet\n");CHKERRQ(ierr); 70 } 71 if (snes->view) { 72 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 73 ierr = (*snes->view)(snes,viewer);CHKERRQ(ierr); 74 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 75 } 76 ierr = PetscViewerASCIIPrintf(viewer," maximum iterations=%D, maximum function evaluations=%D\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr); 77 ierr = PetscViewerASCIIPrintf(viewer," tolerances: relative=%G, absolute=%G, solution=%G\n", 78 snes->rtol,snes->abstol,snes->xtol);CHKERRQ(ierr); 79 ierr = PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",snes->linear_its);CHKERRQ(ierr); 80 ierr = PetscViewerASCIIPrintf(viewer," total number of function evaluations=%D\n",snes->nfuncs);CHKERRQ(ierr); 81 if (snes->ksp_ewconv) { 82 kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 83 if (kctx) { 84 ierr = PetscViewerASCIIPrintf(viewer," Eisenstat-Walker computation of KSP relative tolerance (version %D)\n",kctx->version);CHKERRQ(ierr); 85 ierr = PetscViewerASCIIPrintf(viewer," rtol_0=%G, rtol_max=%G, threshold=%G\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr); 86 ierr = PetscViewerASCIIPrintf(viewer," gamma=%G, alpha=%G, alpha2=%G\n",kctx->gamma,kctx->alpha,kctx->alpha2);CHKERRQ(ierr); 87 } 88 } 89 } else if (isstring) { 90 ierr = SNESGetType(snes,&type);CHKERRQ(ierr); 91 ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr); 92 } 93 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 94 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 95 ierr = KSPView(ksp,viewer);CHKERRQ(ierr); 96 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 97 PetscFunctionReturn(0); 98 } 99 100 /* 101 We retain a list of functions that also take SNES command 102 line options. These are called at the end SNESSetFromOptions() 103 */ 104 #define MAXSETFROMOPTIONS 5 105 static PetscInt numberofsetfromoptions; 106 static PetscErrorCode (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES); 107 108 #undef __FUNCT__ 109 #define __FUNCT__ "SNESAddOptionsChecker" 110 /*@C 111 SNESAddOptionsChecker - Adds an additional function to check for SNES options. 112 113 Not Collective 114 115 Input Parameter: 116 . snescheck - function that checks for options 117 118 Level: developer 119 120 .seealso: SNESSetFromOptions() 121 @*/ 122 PetscErrorCode PETSCSNES_DLLEXPORT SNESAddOptionsChecker(PetscErrorCode (*snescheck)(SNES)) 123 { 124 PetscFunctionBegin; 125 if (numberofsetfromoptions >= MAXSETFROMOPTIONS) { 126 SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %D allowed", MAXSETFROMOPTIONS); 127 } 128 othersetfromoptions[numberofsetfromoptions++] = snescheck; 129 PetscFunctionReturn(0); 130 } 131 132 #undef __FUNCT__ 133 #define __FUNCT__ "SNESSetFromOptions" 134 /*@ 135 SNESSetFromOptions - Sets various SNES and KSP parameters from user options. 136 137 Collective on SNES 138 139 Input Parameter: 140 . snes - the SNES context 141 142 Options Database Keys: 143 + -snes_type <type> - ls, tr, umls, umtr, test 144 . -snes_stol - convergence tolerance in terms of the norm 145 of the change in the solution between steps 146 . -snes_atol <abstol> - absolute tolerance of residual norm 147 . -snes_rtol <rtol> - relative decrease in tolerance norm from initial 148 . -snes_max_it <max_it> - maximum number of iterations 149 . -snes_max_funcs <max_funcs> - maximum number of function evaluations 150 . -snes_max_fail <max_fail> - maximum number of failures 151 . -snes_trtol <trtol> - trust region tolerance 152 . -snes_no_convergence_test - skip convergence test in nonlinear 153 solver; hence iterations will continue until max_it 154 or some other criterion is reached. Saves expense 155 of convergence test 156 . -snes_monitor <optional filename> - prints residual norm at each iteration. if no 157 filename given prints to stdout 158 . -snes_vecmonitor - plots solution at each iteration 159 . -snes_vecmonitor_update - plots update to solution at each iteration 160 . -snes_xmonitor - plots residual norm at each iteration 161 . -snes_fd - use finite differences to compute Jacobian; very slow, only for testing 162 . -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration 163 - -snes_print_converged_reason - print the reason for convergence/divergence after each solve 164 165 Options Database for Eisenstat-Walker method: 166 + -snes_ksp_ew_conv - use Eisenstat-Walker method for determining linear system convergence 167 . -snes_ksp_ew_version ver - version of Eisenstat-Walker method 168 . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0 169 . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax 170 . -snes_ksp_ew_gamma <gamma> - Sets gamma 171 . -snes_ksp_ew_alpha <alpha> - Sets alpha 172 . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2 173 - -snes_ksp_ew_threshold <threshold> - Sets threshold 174 175 Notes: 176 To see all options, run your program with the -help option or consult 177 the users manual. 178 179 Level: beginner 180 181 .keywords: SNES, nonlinear, set, options, database 182 183 .seealso: SNESSetOptionsPrefix() 184 @*/ 185 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFromOptions(SNES snes) 186 { 187 KSP ksp; 188 SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 189 PetscTruth flg; 190 PetscErrorCode ierr; 191 PetscInt i; 192 const char *deft; 193 char type[256], monfilename[PETSC_MAX_PATH_LEN]; 194 PetscViewer monviewer; 195 196 PetscFunctionBegin; 197 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 198 199 ierr = PetscOptionsBegin(snes->comm,snes->prefix,"Nonlinear solver (SNES) options","SNES");CHKERRQ(ierr); 200 if (snes->type_name) { 201 deft = snes->type_name; 202 } else { 203 deft = SNESLS; 204 } 205 206 if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 207 ierr = PetscOptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr); 208 if (flg) { 209 ierr = SNESSetType(snes,type);CHKERRQ(ierr); 210 } else if (!snes->type_name) { 211 ierr = SNESSetType(snes,deft);CHKERRQ(ierr); 212 } 213 ierr = PetscOptionsName("-snes_view","Print detailed information on solver used","SNESView",0);CHKERRQ(ierr); 214 215 ierr = PetscOptionsReal("-snes_stol","Stop if step length less then","SNESSetTolerances",snes->xtol,&snes->xtol,0);CHKERRQ(ierr); 216 ierr = PetscOptionsReal("-snes_atol","Stop if function norm less then","SNESSetTolerances",snes->abstol,&snes->abstol,0);CHKERRQ(ierr); 217 218 ierr = PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less then","SNESSetTolerances",snes->rtol,&snes->rtol,0);CHKERRQ(ierr); 219 ierr = PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,PETSC_NULL);CHKERRQ(ierr); 220 ierr = PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr); 221 ierr = PetscOptionsInt("-snes_max_fail","Maximum failures","SNESSetTolerances",snes->maxFailures,&snes->maxFailures,PETSC_NULL);CHKERRQ(ierr); 222 ierr = PetscOptionsName("-snes_converged_reason","Print reason for converged or diverged","SNESSolve",&flg);CHKERRQ(ierr); 223 if (flg) { 224 snes->printreason = PETSC_TRUE; 225 } 226 227 ierr = PetscOptionsTruth("-snes_ksp_ew_conv","Use Eisentat-Walker linear system convergence test","SNES_KSP_SetParametersEW",snes->ksp_ewconv,&snes->ksp_ewconv,PETSC_NULL);CHKERRQ(ierr); 228 229 ierr = PetscOptionsInt("-snes_ksp_ew_version","Version 1 or 2","SNES_KSP_SetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr); 230 ierr = PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNES_KSP_SetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr); 231 ierr = PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNES_KSP_SetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr); 232 ierr = PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNES_KSP_SetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr); 233 ierr = PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNES_KSP_SetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr); 234 ierr = PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNES_KSP_SetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr); 235 ierr = PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNES_KSP_SetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr); 236 237 ierr = PetscOptionsName("-snes_no_convergence_test","Don't test for convergence","None",&flg);CHKERRQ(ierr); 238 if (flg) {snes->converged = 0;} 239 ierr = PetscOptionsName("-snes_cancelmonitors","Remove all monitors","SNESClearMonitor",&flg);CHKERRQ(ierr); 240 if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);} 241 242 ierr = PetscOptionsString("-snes_monitor","Monitor norm of function","SNESSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 243 if (flg) { 244 ierr = PetscViewerASCIIOpen(snes->comm,monfilename,&monviewer);CHKERRQ(ierr); 245 ierr = SNESSetMonitor(snes,SNESDefaultMonitor,monviewer,(PetscErrorCode (*)(void*))PetscViewerDestroy);CHKERRQ(ierr); 246 } 247 248 ierr = PetscOptionsString("-snes_ratiomonitor","Monitor ratios of norms of function","SNESSetRatioMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 249 if (flg) { 250 ierr = PetscViewerASCIIOpen(snes->comm,monfilename,&monviewer);CHKERRQ(ierr); 251 ierr = SNESSetRatioMonitor(snes,monviewer);CHKERRQ(ierr); 252 } 253 254 ierr = PetscOptionsString("-snes_smonitor","Monitor norm of function (fewer digits)","SNESSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 255 if (flg) { 256 ierr = PetscViewerASCIIOpen(snes->comm,monfilename,&monviewer);CHKERRQ(ierr); 257 ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,monviewer,(PetscErrorCode (*)(void*))PetscViewerDestroy);CHKERRQ(ierr); 258 } 259 260 ierr = PetscOptionsName("-snes_vecmonitor","Plot solution at each iteration","SNESVecViewMonitor",&flg);CHKERRQ(ierr); 261 if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewMonitor,0,0);CHKERRQ(ierr);} 262 ierr = PetscOptionsName("-snes_vecmonitor_update","Plot correction at each iteration","SNESVecViewUpdateMonitor",&flg);CHKERRQ(ierr); 263 if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewUpdateMonitor,0,0);CHKERRQ(ierr);} 264 ierr = PetscOptionsName("-snes_vecmonitor_residual","Plot residual at each iteration","SNESVecViewResidualMonitor",&flg);CHKERRQ(ierr); 265 if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewResidualMonitor,0,0);CHKERRQ(ierr);} 266 ierr = PetscOptionsName("-snes_xmonitor","Plot function norm at each iteration","SNESLGMonitor",&flg);CHKERRQ(ierr); 267 if (flg) {ierr = SNESSetMonitor(snes,SNESLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);} 268 269 ierr = PetscOptionsName("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",&flg);CHKERRQ(ierr); 270 if (flg) { 271 ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr); 272 ierr = PetscInfo(snes,"Setting default finite difference Jacobian matrix\n");CHKERRQ(ierr); 273 } 274 275 for(i = 0; i < numberofsetfromoptions; i++) { 276 ierr = (*othersetfromoptions[i])(snes);CHKERRQ(ierr); 277 } 278 279 if (snes->setfromoptions) { 280 ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr); 281 } 282 ierr = PetscOptionsEnd();CHKERRQ(ierr); 283 284 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 285 ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr); 286 287 PetscFunctionReturn(0); 288 } 289 290 291 #undef __FUNCT__ 292 #define __FUNCT__ "SNESSetApplicationContext" 293 /*@ 294 SNESSetApplicationContext - Sets the optional user-defined context for 295 the nonlinear solvers. 296 297 Collective on SNES 298 299 Input Parameters: 300 + snes - the SNES context 301 - usrP - optional user context 302 303 Level: intermediate 304 305 .keywords: SNES, nonlinear, set, application, context 306 307 .seealso: SNESGetApplicationContext() 308 @*/ 309 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetApplicationContext(SNES snes,void *usrP) 310 { 311 PetscFunctionBegin; 312 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 313 snes->user = usrP; 314 PetscFunctionReturn(0); 315 } 316 317 #undef __FUNCT__ 318 #define __FUNCT__ "SNESGetApplicationContext" 319 /*@C 320 SNESGetApplicationContext - Gets the user-defined context for the 321 nonlinear solvers. 322 323 Not Collective 324 325 Input Parameter: 326 . snes - SNES context 327 328 Output Parameter: 329 . usrP - user context 330 331 Level: intermediate 332 333 .keywords: SNES, nonlinear, get, application, context 334 335 .seealso: SNESSetApplicationContext() 336 @*/ 337 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetApplicationContext(SNES snes,void **usrP) 338 { 339 PetscFunctionBegin; 340 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 341 *usrP = snes->user; 342 PetscFunctionReturn(0); 343 } 344 345 #undef __FUNCT__ 346 #define __FUNCT__ "SNESGetIterationNumber" 347 /*@ 348 SNESGetIterationNumber - Gets the number of nonlinear iterations completed 349 at this time. 350 351 Not Collective 352 353 Input Parameter: 354 . snes - SNES context 355 356 Output Parameter: 357 . iter - iteration number 358 359 Notes: 360 For example, during the computation of iteration 2 this would return 1. 361 362 This is useful for using lagged Jacobians (where one does not recompute the 363 Jacobian at each SNES iteration). For example, the code 364 .vb 365 ierr = SNESGetIterationNumber(snes,&it); 366 if (!(it % 2)) { 367 [compute Jacobian here] 368 } 369 .ve 370 can be used in your ComputeJacobian() function to cause the Jacobian to be 371 recomputed every second SNES iteration. 372 373 Level: intermediate 374 375 .keywords: SNES, nonlinear, get, iteration, number, 376 377 .seealso: SNESGetFunctionNorm(), SNESGetNumberLinearIterations() 378 @*/ 379 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetIterationNumber(SNES snes,PetscInt* iter) 380 { 381 PetscFunctionBegin; 382 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 383 PetscValidIntPointer(iter,2); 384 *iter = snes->iter; 385 PetscFunctionReturn(0); 386 } 387 388 #undef __FUNCT__ 389 #define __FUNCT__ "SNESGetFunctionNorm" 390 /*@ 391 SNESGetFunctionNorm - Gets the norm of the current function that was set 392 with SNESSSetFunction(). 393 394 Collective on SNES 395 396 Input Parameter: 397 . snes - SNES context 398 399 Output Parameter: 400 . fnorm - 2-norm of function 401 402 Level: intermediate 403 404 .keywords: SNES, nonlinear, get, function, norm 405 406 .seealso: SNESGetFunction(), SNESGetIterationNumber(), SNESGetNumberLinearIterations() 407 @*/ 408 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunctionNorm(SNES snes,PetscScalar *fnorm) 409 { 410 PetscFunctionBegin; 411 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 412 PetscValidScalarPointer(fnorm,2); 413 *fnorm = snes->norm; 414 PetscFunctionReturn(0); 415 } 416 417 #undef __FUNCT__ 418 #define __FUNCT__ "SNESGetNumberUnsuccessfulSteps" 419 /*@ 420 SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps 421 attempted by the nonlinear solver. 422 423 Not Collective 424 425 Input Parameter: 426 . snes - SNES context 427 428 Output Parameter: 429 . nfails - number of unsuccessful steps attempted 430 431 Notes: 432 This counter is reset to zero for each successive call to SNESSolve(). 433 434 Level: intermediate 435 436 .keywords: SNES, nonlinear, get, number, unsuccessful, steps 437 @*/ 438 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetNumberUnsuccessfulSteps(SNES snes,PetscInt* nfails) 439 { 440 PetscFunctionBegin; 441 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 442 PetscValidIntPointer(nfails,2); 443 *nfails = snes->numFailures; 444 PetscFunctionReturn(0); 445 } 446 447 #undef __FUNCT__ 448 #define __FUNCT__ "SNESSetMaximumUnsuccessfulSteps" 449 /*@ 450 SNESSetMaximumUnsuccessfulSteps - Sets the maximum number of unsuccessful steps 451 attempted by the nonlinear solver before it gives up. 452 453 Not Collective 454 455 Input Parameters: 456 + snes - SNES context 457 - maxFails - maximum of unsuccessful steps 458 459 Level: intermediate 460 461 .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps 462 @*/ 463 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMaximumUnsuccessfulSteps(SNES snes, PetscInt maxFails) 464 { 465 PetscFunctionBegin; 466 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 467 snes->maxFailures = maxFails; 468 PetscFunctionReturn(0); 469 } 470 471 #undef __FUNCT__ 472 #define __FUNCT__ "SNESGetMaximumUnsuccessfulSteps" 473 /*@ 474 SNESGetMaximumUnsuccessfulSteps - Gets the maximum number of unsuccessful steps 475 attempted by the nonlinear solver before it gives up. 476 477 Not Collective 478 479 Input Parameter: 480 . snes - SNES context 481 482 Output Parameter: 483 . maxFails - maximum of unsuccessful steps 484 485 Level: intermediate 486 487 .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps 488 @*/ 489 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetMaximumUnsuccessfulSteps(SNES snes, PetscInt *maxFails) 490 { 491 PetscFunctionBegin; 492 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 493 PetscValidIntPointer(maxFails,2); 494 *maxFails = snes->maxFailures; 495 PetscFunctionReturn(0); 496 } 497 498 #undef __FUNCT__ 499 #define __FUNCT__ "SNESGetNumberLinearIterations" 500 /*@ 501 SNESGetNumberLinearIterations - Gets the total number of linear iterations 502 used by the nonlinear solver. 503 504 Not Collective 505 506 Input Parameter: 507 . snes - SNES context 508 509 Output Parameter: 510 . lits - number of linear iterations 511 512 Notes: 513 This counter is reset to zero for each successive call to SNESSolve(). 514 515 Level: intermediate 516 517 .keywords: SNES, nonlinear, get, number, linear, iterations 518 519 .seealso: SNESGetIterationNumber(), SNESGetFunctionNorm() 520 @*/ 521 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetNumberLinearIterations(SNES snes,PetscInt* lits) 522 { 523 PetscFunctionBegin; 524 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 525 PetscValidIntPointer(lits,2); 526 *lits = snes->linear_its; 527 PetscFunctionReturn(0); 528 } 529 530 #undef __FUNCT__ 531 #define __FUNCT__ "SNESGetKSP" 532 /*@ 533 SNESGetKSP - Returns the KSP context for a SNES solver. 534 535 Not Collective, but if SNES object is parallel, then KSP object is parallel 536 537 Input Parameter: 538 . snes - the SNES context 539 540 Output Parameter: 541 . ksp - the KSP context 542 543 Notes: 544 The user can then directly manipulate the KSP context to set various 545 options, etc. Likewise, the user can then extract and manipulate the 546 KSP and PC contexts as well. 547 548 Level: beginner 549 550 .keywords: SNES, nonlinear, get, KSP, context 551 552 .seealso: KSPGetPC() 553 @*/ 554 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetKSP(SNES snes,KSP *ksp) 555 { 556 PetscFunctionBegin; 557 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 558 PetscValidPointer(ksp,2); 559 *ksp = snes->ksp; 560 PetscFunctionReturn(0); 561 } 562 563 #undef __FUNCT__ 564 #define __FUNCT__ "SNESPublish_Petsc" 565 static PetscErrorCode SNESPublish_Petsc(PetscObject obj) 566 { 567 PetscFunctionBegin; 568 PetscFunctionReturn(0); 569 } 570 571 /* -----------------------------------------------------------*/ 572 #undef __FUNCT__ 573 #define __FUNCT__ "SNESCreate" 574 /*@ 575 SNESCreate - Creates a nonlinear solver context. 576 577 Collective on MPI_Comm 578 579 Input Parameters: 580 + comm - MPI communicator 581 582 Output Parameter: 583 . outsnes - the new SNES context 584 585 Options Database Keys: 586 + -snes_mf - Activates default matrix-free Jacobian-vector products, 587 and no preconditioning matrix 588 . -snes_mf_operator - Activates default matrix-free Jacobian-vector 589 products, and a user-provided preconditioning matrix 590 as set by SNESSetJacobian() 591 - -snes_fd - Uses (slow!) finite differences to compute Jacobian 592 593 Level: beginner 594 595 .keywords: SNES, nonlinear, create, context 596 597 .seealso: SNESSolve(), SNESDestroy(), SNES 598 @*/ 599 PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate(MPI_Comm comm,SNES *outsnes) 600 { 601 PetscErrorCode ierr; 602 SNES snes; 603 SNES_KSP_EW_ConvCtx *kctx; 604 605 PetscFunctionBegin; 606 PetscValidPointer(outsnes,2); 607 *outsnes = PETSC_NULL; 608 #ifndef PETSC_USE_DYNAMIC_LIBRARIES 609 ierr = SNESInitializePackage(PETSC_NULL);CHKERRQ(ierr); 610 #endif 611 612 ierr = PetscHeaderCreate(snes,_p_SNES,PetscInt,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView);CHKERRQ(ierr); 613 snes->bops->publish = SNESPublish_Petsc; 614 snes->max_its = 50; 615 snes->max_funcs = 10000; 616 snes->norm = 0.0; 617 snes->rtol = 1.e-8; 618 snes->ttol = 0.0; 619 snes->abstol = 1.e-50; 620 snes->xtol = 1.e-8; 621 snes->deltatol = 1.e-12; 622 snes->nfuncs = 0; 623 snes->numFailures = 0; 624 snes->maxFailures = 1; 625 snes->linear_its = 0; 626 snes->numbermonitors = 0; 627 snes->data = 0; 628 snes->view = 0; 629 snes->setupcalled = PETSC_FALSE; 630 snes->ksp_ewconv = PETSC_FALSE; 631 snes->vwork = 0; 632 snes->nwork = 0; 633 snes->conv_hist_len = 0; 634 snes->conv_hist_max = 0; 635 snes->conv_hist = PETSC_NULL; 636 snes->conv_hist_its = PETSC_NULL; 637 snes->conv_hist_reset = PETSC_TRUE; 638 snes->reason = SNES_CONVERGED_ITERATING; 639 640 /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 641 ierr = PetscNew(SNES_KSP_EW_ConvCtx,&kctx);CHKERRQ(ierr); 642 ierr = PetscLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx));CHKERRQ(ierr); 643 snes->kspconvctx = (void*)kctx; 644 kctx->version = 2; 645 kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 646 this was too large for some test cases */ 647 kctx->rtol_last = 0; 648 kctx->rtol_max = .9; 649 kctx->gamma = 1.0; 650 kctx->alpha2 = .5*(1.0 + sqrt(5.0)); 651 kctx->alpha = kctx->alpha2; 652 kctx->threshold = .1; 653 kctx->lresid_last = 0; 654 kctx->norm_last = 0; 655 656 ierr = KSPCreate(comm,&snes->ksp);CHKERRQ(ierr); 657 ierr = PetscLogObjectParent(snes,snes->ksp);CHKERRQ(ierr); 658 659 *outsnes = snes; 660 ierr = PetscPublishAll(snes);CHKERRQ(ierr); 661 PetscFunctionReturn(0); 662 } 663 664 #undef __FUNCT__ 665 #define __FUNCT__ "SNESSetFunction" 666 /*@C 667 SNESSetFunction - Sets the function evaluation routine and function 668 vector for use by the SNES routines in solving systems of nonlinear 669 equations. 670 671 Collective on SNES 672 673 Input Parameters: 674 + snes - the SNES context 675 . func - function evaluation routine 676 . r - vector to store function value 677 - ctx - [optional] user-defined context for private data for the 678 function evaluation routine (may be PETSC_NULL) 679 680 Calling sequence of func: 681 $ func (SNES snes,Vec x,Vec f,void *ctx); 682 683 . f - function vector 684 - ctx - optional user-defined function context 685 686 Notes: 687 The Newton-like methods typically solve linear systems of the form 688 $ f'(x) x = -f(x), 689 where f'(x) denotes the Jacobian matrix and f(x) is the function. 690 691 Level: beginner 692 693 .keywords: SNES, nonlinear, set, function 694 695 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 696 @*/ 697 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx) 698 { 699 PetscFunctionBegin; 700 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 701 PetscValidHeaderSpecific(r,VEC_COOKIE,2); 702 PetscCheckSameComm(snes,1,r,2); 703 704 snes->computefunction = func; 705 snes->vec_func = snes->vec_func_always = r; 706 snes->funP = ctx; 707 PetscFunctionReturn(0); 708 } 709 710 /* --------------------------------------------------------------- */ 711 #undef __FUNCT__ 712 #define __FUNCT__ "SNESSetRhs" 713 /*@C 714 SNESSetRhs - Sets the vector for solving F(x) = rhs. If rhs is not set 715 it assumes a zero right hand side. 716 717 Collective on SNES 718 719 Input Parameters: 720 + snes - the SNES context 721 - rhs - the right hand side vector or PETSC_NULL for a zero right hand side 722 723 Level: intermediate 724 725 .keywords: SNES, nonlinear, set, function, right hand side 726 727 .seealso: SNESGetRhs(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction() 728 @*/ 729 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetRhs(SNES snes,Vec rhs) 730 { 731 PetscErrorCode ierr; 732 733 PetscFunctionBegin; 734 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 735 if (rhs) { 736 PetscValidHeaderSpecific(rhs,VEC_COOKIE,2); 737 PetscCheckSameComm(snes,1,rhs,2); 738 ierr = PetscObjectReference((PetscObject)rhs);CHKERRQ(ierr); 739 } 740 if (snes->afine) { 741 ierr = VecDestroy(snes->afine);CHKERRQ(ierr); 742 } 743 snes->afine = rhs; 744 PetscFunctionReturn(0); 745 } 746 747 #undef __FUNCT__ 748 #define __FUNCT__ "SNESGetRhs" 749 /*@C 750 SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set 751 it assumes a zero right hand side. 752 753 Collective on SNES 754 755 Input Parameter: 756 . snes - the SNES context 757 758 Output Parameter: 759 . rhs - the right hand side vector or PETSC_NULL for a zero right hand side 760 761 Level: intermediate 762 763 .keywords: SNES, nonlinear, get, function, right hand side 764 765 .seealso: SNESSetRhs(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction() 766 @*/ 767 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetRhs(SNES snes,Vec *rhs) 768 { 769 PetscFunctionBegin; 770 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 771 PetscValidPointer(rhs,2); 772 *rhs = snes->afine; 773 PetscFunctionReturn(0); 774 } 775 776 #undef __FUNCT__ 777 #define __FUNCT__ "SNESComputeFunction" 778 /*@ 779 SNESComputeFunction - Calls the function that has been set with 780 SNESSetFunction(). 781 782 Collective on SNES 783 784 Input Parameters: 785 + snes - the SNES context 786 - x - input vector 787 788 Output Parameter: 789 . y - function vector, as set by SNESSetFunction() 790 791 Notes: 792 SNESComputeFunction() is typically used within nonlinear solvers 793 implementations, so most users would not generally call this routine 794 themselves. 795 796 Level: developer 797 798 .keywords: SNES, nonlinear, compute, function 799 800 .seealso: SNESSetFunction(), SNESGetFunction() 801 @*/ 802 PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeFunction(SNES snes,Vec x,Vec y) 803 { 804 PetscErrorCode ierr; 805 806 PetscFunctionBegin; 807 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 808 PetscValidHeaderSpecific(x,VEC_COOKIE,2); 809 PetscValidHeaderSpecific(y,VEC_COOKIE,3); 810 PetscCheckSameComm(snes,1,x,2); 811 PetscCheckSameComm(snes,1,y,3); 812 813 ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 814 if (snes->computefunction) { 815 PetscStackPush("SNES user function"); 816 CHKMEMQ; 817 ierr = (*snes->computefunction)(snes,x,y,snes->funP); 818 CHKMEMQ; 819 PetscStackPop; 820 if (PetscExceptionValue(ierr)) { 821 PetscErrorCode pierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(pierr); 822 } 823 CHKERRQ(ierr); 824 } else if (snes->afine) { 825 ierr = MatMult(snes->jacobian, x, y);CHKERRQ(ierr); 826 } else { 827 SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() before SNESComputeFunction(), likely called from SNESSolve()."); 828 } 829 if (snes->afine) { 830 ierr = VecAXPY(y,-1.0,snes->afine);CHKERRQ(ierr); 831 } 832 snes->nfuncs++; 833 ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 834 PetscFunctionReturn(0); 835 } 836 837 #undef __FUNCT__ 838 #define __FUNCT__ "SNESComputeJacobian" 839 /*@ 840 SNESComputeJacobian - Computes the Jacobian matrix that has been 841 set with SNESSetJacobian(). 842 843 Collective on SNES and Mat 844 845 Input Parameters: 846 + snes - the SNES context 847 - x - input vector 848 849 Output Parameters: 850 + A - Jacobian matrix 851 . B - optional preconditioning matrix 852 - flag - flag indicating matrix structure (one of, SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER) 853 854 Notes: 855 Most users should not need to explicitly call this routine, as it 856 is used internally within the nonlinear solvers. 857 858 See KSPSetOperators() for important information about setting the 859 flag parameter. 860 861 Level: developer 862 863 .keywords: SNES, compute, Jacobian, matrix 864 865 .seealso: SNESSetJacobian(), KSPSetOperators(), MatStructure 866 @*/ 867 PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 868 { 869 PetscErrorCode ierr; 870 871 PetscFunctionBegin; 872 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 873 PetscValidHeaderSpecific(X,VEC_COOKIE,2); 874 PetscValidPointer(flg,5); 875 PetscCheckSameComm(snes,1,X,2); 876 if (!snes->computejacobian) PetscFunctionReturn(0); 877 ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 878 *flg = DIFFERENT_NONZERO_PATTERN; 879 PetscStackPush("SNES user Jacobian function"); 880 CHKMEMQ; 881 ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr); 882 CHKMEMQ; 883 PetscStackPop; 884 ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 885 /* make sure user returned a correct Jacobian and preconditioner */ 886 PetscValidHeaderSpecific(*A,MAT_COOKIE,3); 887 PetscValidHeaderSpecific(*B,MAT_COOKIE,4); 888 PetscFunctionReturn(0); 889 } 890 891 #undef __FUNCT__ 892 #define __FUNCT__ "SNESSetJacobian" 893 /*@C 894 SNESSetJacobian - Sets the function to compute Jacobian as well as the 895 location to store the matrix. 896 897 Collective on SNES and Mat 898 899 Input Parameters: 900 + snes - the SNES context 901 . A - Jacobian matrix 902 . B - preconditioner matrix (usually same as the Jacobian) 903 . func - Jacobian evaluation routine 904 - ctx - [optional] user-defined context for private data for the 905 Jacobian evaluation routine (may be PETSC_NULL) 906 907 Calling sequence of func: 908 $ func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 909 910 + x - input vector 911 . A - Jacobian matrix 912 . B - preconditioner matrix, usually the same as A 913 . flag - flag indicating information about the preconditioner matrix 914 structure (same as flag in KSPSetOperators()), one of SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER 915 - ctx - [optional] user-defined Jacobian context 916 917 Notes: 918 See KSPSetOperators() for important information about setting the flag 919 output parameter in the routine func(). Be sure to read this information! 920 921 The routine func() takes Mat * as the matrix arguments rather than Mat. 922 This allows the Jacobian evaluation routine to replace A and/or B with a 923 completely new new matrix structure (not just different matrix elements) 924 when appropriate, for instance, if the nonzero structure is changing 925 throughout the global iterations. 926 927 Level: beginner 928 929 .keywords: SNES, nonlinear, set, Jacobian, matrix 930 931 .seealso: KSPSetOperators(), SNESSetFunction(), MatSNESMFComputeJacobian(), SNESDefaultComputeJacobianColor(), MatStructure 932 @*/ 933 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 934 { 935 PetscErrorCode ierr; 936 937 PetscFunctionBegin; 938 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 939 if (A) PetscValidHeaderSpecific(A,MAT_COOKIE,2); 940 if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3); 941 if (A) PetscCheckSameComm(snes,1,A,2); 942 if (B) PetscCheckSameComm(snes,1,B,2); 943 if (func) snes->computejacobian = func; 944 if (ctx) snes->jacP = ctx; 945 if (A) { 946 if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 947 snes->jacobian = A; 948 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 949 } 950 if (B) { 951 if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 952 snes->jacobian_pre = B; 953 ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 954 } 955 ierr = KSPSetOperators(snes->ksp,A,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 956 PetscFunctionReturn(0); 957 } 958 959 #undef __FUNCT__ 960 #define __FUNCT__ "SNESGetJacobian" 961 /*@C 962 SNESGetJacobian - Returns the Jacobian matrix and optionally the user 963 provided context for evaluating the Jacobian. 964 965 Not Collective, but Mat object will be parallel if SNES object is 966 967 Input Parameter: 968 . snes - the nonlinear solver context 969 970 Output Parameters: 971 + A - location to stash Jacobian matrix (or PETSC_NULL) 972 . B - location to stash preconditioner matrix (or PETSC_NULL) 973 . func - location to put Jacobian function (or PETSC_NULL) 974 - ctx - location to stash Jacobian ctx (or PETSC_NULL) 975 976 Level: advanced 977 978 .seealso: SNESSetJacobian(), SNESComputeJacobian() 979 @*/ 980 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetJacobian(SNES snes,Mat *A,Mat *B,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx) 981 { 982 PetscFunctionBegin; 983 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 984 if (A) *A = snes->jacobian; 985 if (B) *B = snes->jacobian_pre; 986 if (func) *func = snes->computejacobian; 987 if (ctx) *ctx = snes->jacP; 988 PetscFunctionReturn(0); 989 } 990 991 /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 992 EXTERN PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*); 993 994 #undef __FUNCT__ 995 #define __FUNCT__ "SNESSetUp" 996 /*@ 997 SNESSetUp - Sets up the internal data structures for the later use 998 of a nonlinear solver. 999 1000 Collective on SNES 1001 1002 Input Parameters: 1003 . snes - the SNES context 1004 1005 Notes: 1006 For basic use of the SNES solvers the user need not explicitly call 1007 SNESSetUp(), since these actions will automatically occur during 1008 the call to SNESSolve(). However, if one wishes to control this 1009 phase separately, SNESSetUp() should be called after SNESCreate() 1010 and optional routines of the form SNESSetXXX(), but before SNESSolve(). 1011 1012 Level: advanced 1013 1014 .keywords: SNES, nonlinear, setup 1015 1016 .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 1017 @*/ 1018 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUp(SNES snes) 1019 { 1020 PetscErrorCode ierr; 1021 PetscTruth flg, iseqtr; 1022 1023 PetscFunctionBegin; 1024 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1025 if (snes->setupcalled) PetscFunctionReturn(0); 1026 1027 ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr); 1028 /* 1029 This version replaces the user provided Jacobian matrix with a 1030 matrix-free version but still employs the user-provided preconditioner matrix 1031 */ 1032 if (flg) { 1033 Mat J; 1034 ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1035 ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 1036 ierr = PetscInfo(snes,"Setting default matrix-free operator routines\n");CHKERRQ(ierr); 1037 ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 1038 ierr = MatDestroy(J);CHKERRQ(ierr); 1039 } 1040 1041 #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) && !defined(PETSC_USE_MAT_SINGLE) && !defined(PETSC_USE_LONG_DOUBLE) && !defined(PETSC_USE_INT) 1042 ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr); 1043 if (flg) { 1044 Mat J; 1045 ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1046 ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 1047 ierr = MatDestroy(J);CHKERRQ(ierr); 1048 } 1049 #endif 1050 1051 ierr = PetscOptionsHasName(snes->prefix,"-snes_mf",&flg);CHKERRQ(ierr); 1052 /* 1053 This version replaces both the user-provided Jacobian and the user- 1054 provided preconditioner matrix with the default matrix free version. 1055 */ 1056 if (flg) { 1057 Mat J; 1058 KSP ksp; 1059 PC pc; 1060 1061 ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1062 ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 1063 ierr = PetscInfo(snes,"Setting default matrix-free operator and preconditioner routines;\nThat is no preconditioner is being used.\n");CHKERRQ(ierr); 1064 ierr = SNESSetJacobian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr); 1065 ierr = MatDestroy(J);CHKERRQ(ierr); 1066 1067 /* force no preconditioner */ 1068 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1069 ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 1070 ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr); 1071 if (!flg) { 1072 ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 1073 } 1074 } 1075 1076 if (!snes->vec_func && !snes->afine) { 1077 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 1078 } 1079 if (!snes->computefunction && !snes->afine) { 1080 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 1081 } 1082 if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option"); 1083 if (snes->vec_func == snes->vec_sol) { 1084 SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector"); 1085 } 1086 1087 /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 1088 ierr = PetscTypeCompare((PetscObject)snes,SNESTR,&iseqtr);CHKERRQ(ierr); 1089 if (snes->ksp_ewconv && !iseqtr) { 1090 KSP ksp; 1091 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1092 ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,snes);CHKERRQ(ierr); 1093 } 1094 1095 if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);} 1096 snes->setupcalled = PETSC_TRUE; 1097 PetscFunctionReturn(0); 1098 } 1099 1100 #undef __FUNCT__ 1101 #define __FUNCT__ "SNESDestroy" 1102 /*@ 1103 SNESDestroy - Destroys the nonlinear solver context that was created 1104 with SNESCreate(). 1105 1106 Collective on SNES 1107 1108 Input Parameter: 1109 . snes - the SNES context 1110 1111 Level: beginner 1112 1113 .keywords: SNES, nonlinear, destroy 1114 1115 .seealso: SNESCreate(), SNESSolve() 1116 @*/ 1117 PetscErrorCode PETSCSNES_DLLEXPORT SNESDestroy(SNES snes) 1118 { 1119 PetscErrorCode ierr; 1120 1121 PetscFunctionBegin; 1122 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1123 if (--snes->refct > 0) PetscFunctionReturn(0); 1124 1125 /* if memory was published with AMS then destroy it */ 1126 ierr = PetscObjectDepublish(snes);CHKERRQ(ierr); 1127 1128 if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);} 1129 ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr); 1130 if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 1131 if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 1132 if (snes->afine) {ierr = VecDestroy(snes->afine);CHKERRQ(ierr);} 1133 ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr); 1134 if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);} 1135 ierr = SNESClearMonitor(snes);CHKERRQ(ierr); 1136 ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr); 1137 PetscFunctionReturn(0); 1138 } 1139 1140 /* ----------- Routines to set solver parameters ---------- */ 1141 1142 #undef __FUNCT__ 1143 #define __FUNCT__ "SNESSetTolerances" 1144 /*@ 1145 SNESSetTolerances - Sets various parameters used in convergence tests. 1146 1147 Collective on SNES 1148 1149 Input Parameters: 1150 + snes - the SNES context 1151 . abstol - absolute convergence tolerance 1152 . rtol - relative convergence tolerance 1153 . stol - convergence tolerance in terms of the norm 1154 of the change in the solution between steps 1155 . maxit - maximum number of iterations 1156 - maxf - maximum number of function evaluations 1157 1158 Options Database Keys: 1159 + -snes_atol <abstol> - Sets abstol 1160 . -snes_rtol <rtol> - Sets rtol 1161 . -snes_stol <stol> - Sets stol 1162 . -snes_max_it <maxit> - Sets maxit 1163 - -snes_max_funcs <maxf> - Sets maxf 1164 1165 Notes: 1166 The default maximum number of iterations is 50. 1167 The default maximum number of function evaluations is 1000. 1168 1169 Level: intermediate 1170 1171 .keywords: SNES, nonlinear, set, convergence, tolerances 1172 1173 .seealso: SNESSetTrustRegionTolerance() 1174 @*/ 1175 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf) 1176 { 1177 PetscFunctionBegin; 1178 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1179 if (abstol != PETSC_DEFAULT) snes->abstol = abstol; 1180 if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1181 if (stol != PETSC_DEFAULT) snes->xtol = stol; 1182 if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1183 if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 1184 PetscFunctionReturn(0); 1185 } 1186 1187 #undef __FUNCT__ 1188 #define __FUNCT__ "SNESGetTolerances" 1189 /*@ 1190 SNESGetTolerances - Gets various parameters used in convergence tests. 1191 1192 Not Collective 1193 1194 Input Parameters: 1195 + snes - the SNES context 1196 . abstol - absolute convergence tolerance 1197 . rtol - relative convergence tolerance 1198 . stol - convergence tolerance in terms of the norm 1199 of the change in the solution between steps 1200 . maxit - maximum number of iterations 1201 - maxf - maximum number of function evaluations 1202 1203 Notes: 1204 The user can specify PETSC_NULL for any parameter that is not needed. 1205 1206 Level: intermediate 1207 1208 .keywords: SNES, nonlinear, get, convergence, tolerances 1209 1210 .seealso: SNESSetTolerances() 1211 @*/ 1212 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetTolerances(SNES snes,PetscReal *abstol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf) 1213 { 1214 PetscFunctionBegin; 1215 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1216 if (abstol) *abstol = snes->abstol; 1217 if (rtol) *rtol = snes->rtol; 1218 if (stol) *stol = snes->xtol; 1219 if (maxit) *maxit = snes->max_its; 1220 if (maxf) *maxf = snes->max_funcs; 1221 PetscFunctionReturn(0); 1222 } 1223 1224 #undef __FUNCT__ 1225 #define __FUNCT__ "SNESSetTrustRegionTolerance" 1226 /*@ 1227 SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 1228 1229 Collective on SNES 1230 1231 Input Parameters: 1232 + snes - the SNES context 1233 - tol - tolerance 1234 1235 Options Database Key: 1236 . -snes_trtol <tol> - Sets tol 1237 1238 Level: intermediate 1239 1240 .keywords: SNES, nonlinear, set, trust region, tolerance 1241 1242 .seealso: SNESSetTolerances() 1243 @*/ 1244 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTrustRegionTolerance(SNES snes,PetscReal tol) 1245 { 1246 PetscFunctionBegin; 1247 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1248 snes->deltatol = tol; 1249 PetscFunctionReturn(0); 1250 } 1251 1252 /* 1253 Duplicate the lg monitors for SNES from KSP; for some reason with 1254 dynamic libraries things don't work under Sun4 if we just use 1255 macros instead of functions 1256 */ 1257 #undef __FUNCT__ 1258 #define __FUNCT__ "SNESLGMonitor" 1259 PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitor(SNES snes,PetscInt it,PetscReal norm,void *ctx) 1260 { 1261 PetscErrorCode ierr; 1262 1263 PetscFunctionBegin; 1264 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1265 ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr); 1266 PetscFunctionReturn(0); 1267 } 1268 1269 #undef __FUNCT__ 1270 #define __FUNCT__ "SNESLGMonitorCreate" 1271 PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1272 { 1273 PetscErrorCode ierr; 1274 1275 PetscFunctionBegin; 1276 ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr); 1277 PetscFunctionReturn(0); 1278 } 1279 1280 #undef __FUNCT__ 1281 #define __FUNCT__ "SNESLGMonitorDestroy" 1282 PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitorDestroy(PetscDrawLG draw) 1283 { 1284 PetscErrorCode ierr; 1285 1286 PetscFunctionBegin; 1287 ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr); 1288 PetscFunctionReturn(0); 1289 } 1290 1291 /* ------------ Routines to set performance monitoring options ----------- */ 1292 1293 #undef __FUNCT__ 1294 #define __FUNCT__ "SNESSetMonitor" 1295 /*@C 1296 SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every 1297 iteration of the nonlinear solver to display the iteration's 1298 progress. 1299 1300 Collective on SNES 1301 1302 Input Parameters: 1303 + snes - the SNES context 1304 . func - monitoring routine 1305 . mctx - [optional] user-defined context for private data for the 1306 monitor routine (use PETSC_NULL if no context is desired) 1307 - monitordestroy - [optional] routine that frees monitor context 1308 (may be PETSC_NULL) 1309 1310 Calling sequence of func: 1311 $ int func(SNES snes,PetscInt its, PetscReal norm,void *mctx) 1312 1313 + snes - the SNES context 1314 . its - iteration number 1315 . norm - 2-norm function value (may be estimated) 1316 - mctx - [optional] monitoring context 1317 1318 Options Database Keys: 1319 + -snes_monitor - sets SNESDefaultMonitor() 1320 . -snes_xmonitor - sets line graph monitor, 1321 uses SNESLGMonitorCreate() 1322 _ -snes_cancelmonitors - cancels all monitors that have 1323 been hardwired into a code by 1324 calls to SNESSetMonitor(), but 1325 does not cancel those set via 1326 the options database. 1327 1328 Notes: 1329 Several different monitoring routines may be set by calling 1330 SNESSetMonitor() multiple times; all will be called in the 1331 order in which they were set. 1332 1333 Level: intermediate 1334 1335 .keywords: SNES, nonlinear, set, monitor 1336 1337 .seealso: SNESDefaultMonitor(), SNESClearMonitor() 1338 @*/ 1339 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMonitor(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void*)) 1340 { 1341 PetscFunctionBegin; 1342 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1343 if (snes->numbermonitors >= MAXSNESMONITORS) { 1344 SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1345 } 1346 snes->monitor[snes->numbermonitors] = func; 1347 snes->monitordestroy[snes->numbermonitors] = monitordestroy; 1348 snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 1349 PetscFunctionReturn(0); 1350 } 1351 1352 #undef __FUNCT__ 1353 #define __FUNCT__ "SNESClearMonitor" 1354 /*@C 1355 SNESClearMonitor - Clears all the monitor functions for a SNES object. 1356 1357 Collective on SNES 1358 1359 Input Parameters: 1360 . snes - the SNES context 1361 1362 Options Database Key: 1363 . -snes_cancelmonitors - cancels all monitors that have been hardwired 1364 into a code by calls to SNESSetMonitor(), but does not cancel those 1365 set via the options database 1366 1367 Notes: 1368 There is no way to clear one specific monitor from a SNES object. 1369 1370 Level: intermediate 1371 1372 .keywords: SNES, nonlinear, set, monitor 1373 1374 .seealso: SNESDefaultMonitor(), SNESSetMonitor() 1375 @*/ 1376 PetscErrorCode PETSCSNES_DLLEXPORT SNESClearMonitor(SNES snes) 1377 { 1378 PetscErrorCode ierr; 1379 PetscInt i; 1380 1381 PetscFunctionBegin; 1382 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1383 for (i=0; i<snes->numbermonitors; i++) { 1384 if (snes->monitordestroy[i]) { 1385 ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr); 1386 } 1387 } 1388 snes->numbermonitors = 0; 1389 PetscFunctionReturn(0); 1390 } 1391 1392 #undef __FUNCT__ 1393 #define __FUNCT__ "SNESSetConvergenceTest" 1394 /*@C 1395 SNESSetConvergenceTest - Sets the function that is to be used 1396 to test for convergence of the nonlinear iterative solution. 1397 1398 Collective on SNES 1399 1400 Input Parameters: 1401 + snes - the SNES context 1402 . func - routine to test for convergence 1403 - cctx - [optional] context for private data for the convergence routine 1404 (may be PETSC_NULL) 1405 1406 Calling sequence of func: 1407 $ PetscErrorCode func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx) 1408 1409 + snes - the SNES context 1410 . cctx - [optional] convergence context 1411 . reason - reason for convergence/divergence 1412 . xnorm - 2-norm of current iterate 1413 . gnorm - 2-norm of current step 1414 - f - 2-norm of function 1415 1416 Level: advanced 1417 1418 .keywords: SNES, nonlinear, set, convergence, test 1419 1420 .seealso: SNESConverged_LS(), SNESConverged_TR() 1421 @*/ 1422 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx) 1423 { 1424 PetscFunctionBegin; 1425 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1426 (snes)->converged = func; 1427 (snes)->cnvP = cctx; 1428 PetscFunctionReturn(0); 1429 } 1430 1431 #undef __FUNCT__ 1432 #define __FUNCT__ "SNESGetConvergedReason" 1433 /*@ 1434 SNESGetConvergedReason - Gets the reason the SNES iteration was stopped. 1435 1436 Not Collective 1437 1438 Input Parameter: 1439 . snes - the SNES context 1440 1441 Output Parameter: 1442 . reason - negative value indicates diverged, positive value converged, see petscsnes.h or the 1443 manual pages for the individual convergence tests for complete lists 1444 1445 Level: intermediate 1446 1447 Notes: Can only be called after the call the SNESSolve() is complete. 1448 1449 .keywords: SNES, nonlinear, set, convergence, test 1450 1451 .seealso: SNESSetConvergenceTest(), SNESConverged_LS(), SNESConverged_TR(), SNESConvergedReason 1452 @*/ 1453 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason) 1454 { 1455 PetscFunctionBegin; 1456 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1457 PetscValidPointer(reason,2); 1458 *reason = snes->reason; 1459 PetscFunctionReturn(0); 1460 } 1461 1462 #undef __FUNCT__ 1463 #define __FUNCT__ "SNESSetConvergenceHistory" 1464 /*@ 1465 SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1466 1467 Collective on SNES 1468 1469 Input Parameters: 1470 + snes - iterative context obtained from SNESCreate() 1471 . a - array to hold history 1472 . its - integer array holds the number of linear iterations for each solve. 1473 . na - size of a and its 1474 - reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero, 1475 else it continues storing new values for new nonlinear solves after the old ones 1476 1477 Notes: 1478 If set, this array will contain the function norms computed 1479 at each step. 1480 1481 This routine is useful, e.g., when running a code for purposes 1482 of accurate performance monitoring, when no I/O should be done 1483 during the section of code that is being timed. 1484 1485 Level: intermediate 1486 1487 .keywords: SNES, set, convergence, history 1488 1489 .seealso: SNESGetConvergenceHistory() 1490 1491 @*/ 1492 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt *its,PetscInt na,PetscTruth reset) 1493 { 1494 PetscFunctionBegin; 1495 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1496 if (na) PetscValidScalarPointer(a,2); 1497 snes->conv_hist = a; 1498 snes->conv_hist_its = its; 1499 snes->conv_hist_max = na; 1500 snes->conv_hist_reset = reset; 1501 PetscFunctionReturn(0); 1502 } 1503 1504 #undef __FUNCT__ 1505 #define __FUNCT__ "SNESGetConvergenceHistory" 1506 /*@C 1507 SNESGetConvergenceHistory - Gets the array used to hold the convergence history. 1508 1509 Collective on SNES 1510 1511 Input Parameter: 1512 . snes - iterative context obtained from SNESCreate() 1513 1514 Output Parameters: 1515 . a - array to hold history 1516 . its - integer array holds the number of linear iterations (or 1517 negative if not converged) for each solve. 1518 - na - size of a and its 1519 1520 Notes: 1521 The calling sequence for this routine in Fortran is 1522 $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr) 1523 1524 This routine is useful, e.g., when running a code for purposes 1525 of accurate performance monitoring, when no I/O should be done 1526 during the section of code that is being timed. 1527 1528 Level: intermediate 1529 1530 .keywords: SNES, get, convergence, history 1531 1532 .seealso: SNESSetConvergencHistory() 1533 1534 @*/ 1535 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na) 1536 { 1537 PetscFunctionBegin; 1538 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1539 if (a) *a = snes->conv_hist; 1540 if (its) *its = snes->conv_hist_its; 1541 if (na) *na = snes->conv_hist_len; 1542 PetscFunctionReturn(0); 1543 } 1544 1545 #undef __FUNCT__ 1546 #define __FUNCT__ "SNESSetUpdate" 1547 /*@C 1548 SNESSetUpdate - Sets the general-purpose update function called 1549 at the beginning o every iteration of the nonlinear solve. Specifically 1550 it is called just before the Jacobian is "evaluated". 1551 1552 Collective on SNES 1553 1554 Input Parameters: 1555 . snes - The nonlinear solver context 1556 . func - The function 1557 1558 Calling sequence of func: 1559 . func (SNES snes, PetscInt step); 1560 1561 . step - The current step of the iteration 1562 1563 Level: intermediate 1564 1565 .keywords: SNES, update 1566 1567 .seealso SNESDefaultUpdate(), SNESSetRhsBC(), SNESSetSolutionBC() 1568 @*/ 1569 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt)) 1570 { 1571 PetscFunctionBegin; 1572 PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 1573 snes->update = func; 1574 PetscFunctionReturn(0); 1575 } 1576 1577 #undef __FUNCT__ 1578 #define __FUNCT__ "SNESDefaultUpdate" 1579 /*@ 1580 SNESDefaultUpdate - The default update function which does nothing. 1581 1582 Not collective 1583 1584 Input Parameters: 1585 . snes - The nonlinear solver context 1586 . step - The current step of the iteration 1587 1588 Level: intermediate 1589 1590 .keywords: SNES, update 1591 .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultSolutionBC() 1592 @*/ 1593 PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultUpdate(SNES snes, PetscInt step) 1594 { 1595 PetscFunctionBegin; 1596 PetscFunctionReturn(0); 1597 } 1598 1599 #undef __FUNCT__ 1600 #define __FUNCT__ "SNESScaleStep_Private" 1601 /* 1602 SNESScaleStep_Private - Scales a step so that its length is less than the 1603 positive parameter delta. 1604 1605 Input Parameters: 1606 + snes - the SNES context 1607 . y - approximate solution of linear system 1608 . fnorm - 2-norm of current function 1609 - delta - trust region size 1610 1611 Output Parameters: 1612 + gpnorm - predicted function norm at the new point, assuming local 1613 linearization. The value is zero if the step lies within the trust 1614 region, and exceeds zero otherwise. 1615 - ynorm - 2-norm of the step 1616 1617 Note: 1618 For non-trust region methods such as SNESLS, the parameter delta 1619 is set to be the maximum allowable step size. 1620 1621 .keywords: SNES, nonlinear, scale, step 1622 */ 1623 PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm) 1624 { 1625 PetscReal nrm; 1626 PetscScalar cnorm; 1627 PetscErrorCode ierr; 1628 1629 PetscFunctionBegin; 1630 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1631 PetscValidHeaderSpecific(y,VEC_COOKIE,2); 1632 PetscCheckSameComm(snes,1,y,2); 1633 1634 ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr); 1635 if (nrm > *delta) { 1636 nrm = *delta/nrm; 1637 *gpnorm = (1.0 - nrm)*(*fnorm); 1638 cnorm = nrm; 1639 ierr = VecScale(y,cnorm);CHKERRQ(ierr); 1640 *ynorm = *delta; 1641 } else { 1642 *gpnorm = 0.0; 1643 *ynorm = nrm; 1644 } 1645 PetscFunctionReturn(0); 1646 } 1647 1648 #undef __FUNCT__ 1649 #define __FUNCT__ "SNESSolve" 1650 /*@ 1651 SNESSolve - Solves a nonlinear system F(x) = b. 1652 Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX(). 1653 1654 Collective on SNES 1655 1656 Input Parameters: 1657 + snes - the SNES context 1658 . b - the constant part of the equation, or PETSC_NULL to use zero. 1659 - x - the solution vector, or PETSC_NULL if it was set with SNESSetSolution() 1660 1661 Notes: 1662 The user should initialize the vector,x, with the initial guess 1663 for the nonlinear solve prior to calling SNESSolve. In particular, 1664 to employ an initial guess of zero, the user should explicitly set 1665 this vector to zero by calling VecSet(). 1666 1667 Level: beginner 1668 1669 .keywords: SNES, nonlinear, solve 1670 1671 .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetRhs(), SNESSetSolution() 1672 @*/ 1673 PetscErrorCode PETSCSNES_DLLEXPORT SNESSolve(SNES snes,Vec b,Vec x) 1674 { 1675 PetscErrorCode ierr; 1676 PetscTruth flg; 1677 char filename[PETSC_MAX_PATH_LEN]; 1678 PetscViewer viewer; 1679 1680 PetscFunctionBegin; 1681 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1682 if (!snes->solve) SETERRQ(PETSC_ERR_ORDER,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()"); 1683 1684 if (b) { 1685 ierr = SNESSetRhs(snes, b); CHKERRQ(ierr); 1686 if (!snes->vec_func) { 1687 Vec r; 1688 1689 ierr = VecDuplicate(b, &r); CHKERRQ(ierr); 1690 ierr = SNESSetFunction(snes, r, PETSC_NULL, PETSC_NULL); CHKERRQ(ierr); 1691 } 1692 } 1693 if (x) { 1694 PetscValidHeaderSpecific(x,VEC_COOKIE,3); 1695 PetscCheckSameComm(snes,1,x,3); 1696 } else { 1697 ierr = SNESGetSolution(snes, &x); CHKERRQ(ierr); 1698 if (!x) { 1699 ierr = VecDuplicate(snes->vec_func_always, &x); CHKERRQ(ierr); 1700 } 1701 } 1702 snes->vec_sol = snes->vec_sol_always = x; 1703 if (!snes->setupcalled) { 1704 ierr = SNESSetUp(snes);CHKERRQ(ierr); 1705 } 1706 if (snes->conv_hist_reset) snes->conv_hist_len = 0; 1707 ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 1708 snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0; 1709 1710 ierr = PetscExceptionTry1((*(snes)->solve)(snes),PETSC_ERR_ARG_DOMAIN); 1711 if (PetscExceptionValue(ierr)) { 1712 /* this means that a caller above me has also tryed this exception so I don't handle it here, pass it up */ 1713 PetscErrorCode pierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(pierr); 1714 } else if (PetscExceptionCaught(ierr,PETSC_ERR_ARG_DOMAIN)) { 1715 /* translate exception into SNES not converged reason */ 1716 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1717 ierr = 0; 1718 } 1719 CHKERRQ(ierr); 1720 1721 ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 1722 ierr = PetscOptionsGetString(snes->prefix,"-snes_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 1723 if (flg && !PetscPreLoadingOn) { 1724 ierr = PetscViewerASCIIOpen(snes->comm,filename,&viewer);CHKERRQ(ierr); 1725 ierr = SNESView(snes,viewer);CHKERRQ(ierr); 1726 ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr); 1727 } 1728 1729 ierr = PetscOptionsHasName(snes->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr); 1730 if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); } 1731 if (snes->printreason) { 1732 if (snes->reason > 0) { 1733 ierr = PetscPrintf(snes->comm,"Nonlinear solve converged due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr); 1734 } else { 1735 ierr = PetscPrintf(snes->comm,"Nonlinear solve did not converge due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr); 1736 } 1737 } 1738 1739 PetscFunctionReturn(0); 1740 } 1741 1742 /* --------- Internal routines for SNES Package --------- */ 1743 1744 #undef __FUNCT__ 1745 #define __FUNCT__ "SNESSetType" 1746 /*@C 1747 SNESSetType - Sets the method for the nonlinear solver. 1748 1749 Collective on SNES 1750 1751 Input Parameters: 1752 + snes - the SNES context 1753 - type - a known method 1754 1755 Options Database Key: 1756 . -snes_type <type> - Sets the method; use -help for a list 1757 of available methods (for instance, ls or tr) 1758 1759 Notes: 1760 See "petsc/include/petscsnes.h" for available methods (for instance) 1761 + SNESLS - Newton's method with line search 1762 (systems of nonlinear equations) 1763 . SNESTR - Newton's method with trust region 1764 (systems of nonlinear equations) 1765 1766 Normally, it is best to use the SNESSetFromOptions() command and then 1767 set the SNES solver type from the options database rather than by using 1768 this routine. Using the options database provides the user with 1769 maximum flexibility in evaluating the many nonlinear solvers. 1770 The SNESSetType() routine is provided for those situations where it 1771 is necessary to set the nonlinear solver independently of the command 1772 line or options database. This might be the case, for example, when 1773 the choice of solver changes during the execution of the program, 1774 and the user's application is taking responsibility for choosing the 1775 appropriate method. 1776 1777 Level: intermediate 1778 1779 .keywords: SNES, set, type 1780 1781 .seealso: SNESType, SNESCreate() 1782 1783 @*/ 1784 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetType(SNES snes,SNESType type) 1785 { 1786 PetscErrorCode ierr,(*r)(SNES); 1787 PetscTruth match; 1788 1789 PetscFunctionBegin; 1790 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1791 PetscValidCharPointer(type,2); 1792 1793 ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr); 1794 if (match) PetscFunctionReturn(0); 1795 1796 if (snes->setupcalled) { 1797 snes->setupcalled = PETSC_FALSE; 1798 ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr); 1799 snes->data = 0; 1800 } 1801 1802 /* Get the function pointers for the iterative method requested */ 1803 if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 1804 ierr = PetscFListFind(snes->comm,SNESList,type,(void (**)(void)) &r);CHKERRQ(ierr); 1805 if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type); 1806 ierr = PetscFree(snes->data);CHKERRQ(ierr); 1807 snes->data = 0; 1808 ierr = (*r)(snes);CHKERRQ(ierr); 1809 ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr); 1810 PetscFunctionReturn(0); 1811 } 1812 1813 1814 /* --------------------------------------------------------------------- */ 1815 #undef __FUNCT__ 1816 #define __FUNCT__ "SNESRegisterDestroy" 1817 /*@ 1818 SNESRegisterDestroy - Frees the list of nonlinear solvers that were 1819 registered by SNESRegisterDynamic(). 1820 1821 Not Collective 1822 1823 Level: advanced 1824 1825 .keywords: SNES, nonlinear, register, destroy 1826 1827 .seealso: SNESRegisterAll(), SNESRegisterAll() 1828 @*/ 1829 PetscErrorCode PETSCSNES_DLLEXPORT SNESRegisterDestroy(void) 1830 { 1831 PetscErrorCode ierr; 1832 1833 PetscFunctionBegin; 1834 if (SNESList) { 1835 ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr); 1836 SNESList = 0; 1837 } 1838 SNESRegisterAllCalled = PETSC_FALSE; 1839 PetscFunctionReturn(0); 1840 } 1841 1842 #undef __FUNCT__ 1843 #define __FUNCT__ "SNESGetType" 1844 /*@C 1845 SNESGetType - Gets the SNES method type and name (as a string). 1846 1847 Not Collective 1848 1849 Input Parameter: 1850 . snes - nonlinear solver context 1851 1852 Output Parameter: 1853 . type - SNES method (a character string) 1854 1855 Level: intermediate 1856 1857 .keywords: SNES, nonlinear, get, type, name 1858 @*/ 1859 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetType(SNES snes,SNESType *type) 1860 { 1861 PetscFunctionBegin; 1862 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1863 PetscValidPointer(type,2); 1864 *type = snes->type_name; 1865 PetscFunctionReturn(0); 1866 } 1867 1868 #undef __FUNCT__ 1869 #define __FUNCT__ "SNESGetSolution" 1870 /*@ 1871 SNESGetSolution - Returns the vector where the approximate solution is 1872 stored. 1873 1874 Not Collective, but Vec is parallel if SNES is parallel 1875 1876 Input Parameter: 1877 . snes - the SNES context 1878 1879 Output Parameter: 1880 . x - the solution 1881 1882 Level: intermediate 1883 1884 .keywords: SNES, nonlinear, get, solution 1885 1886 .seealso: SNESSetSolution(), SNESGetFunction(), SNESGetSolutionUpdate() 1887 @*/ 1888 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolution(SNES snes,Vec *x) 1889 { 1890 PetscFunctionBegin; 1891 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1892 PetscValidPointer(x,2); 1893 *x = snes->vec_sol_always; 1894 PetscFunctionReturn(0); 1895 } 1896 1897 #undef __FUNCT__ 1898 #define __FUNCT__ "SNESSetSolution" 1899 /*@ 1900 SNESSetSolution - Sets the vector where the approximate solution is stored. 1901 1902 Not Collective, but Vec is parallel if SNES is parallel 1903 1904 Input Parameters: 1905 + snes - the SNES context 1906 - x - the solution 1907 1908 Output Parameter: 1909 1910 Level: intermediate 1911 1912 Notes: this is not normally used, rather one simply calls SNESSolve() with 1913 the appropriate solution vector. 1914 1915 .keywords: SNES, nonlinear, set, solution 1916 1917 .seealso: SNESGetSolution(), SNESGetFunction(), SNESGetSolutionUpdate() 1918 @*/ 1919 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetSolution(SNES snes,Vec x) 1920 { 1921 PetscFunctionBegin; 1922 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1923 PetscValidHeaderSpecific(x,VEC_COOKIE,2); 1924 PetscCheckSameComm(snes,1,x,2); 1925 snes->vec_sol_always = x; 1926 PetscFunctionReturn(0); 1927 } 1928 1929 #undef __FUNCT__ 1930 #define __FUNCT__ "SNESGetSolutionUpdate" 1931 /*@ 1932 SNESGetSolutionUpdate - Returns the vector where the solution update is 1933 stored. 1934 1935 Not Collective, but Vec is parallel if SNES is parallel 1936 1937 Input Parameter: 1938 . snes - the SNES context 1939 1940 Output Parameter: 1941 . x - the solution update 1942 1943 Level: advanced 1944 1945 .keywords: SNES, nonlinear, get, solution, update 1946 1947 .seealso: SNESGetSolution(), SNESGetFunction 1948 @*/ 1949 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolutionUpdate(SNES snes,Vec *x) 1950 { 1951 PetscFunctionBegin; 1952 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1953 PetscValidPointer(x,2); 1954 *x = snes->vec_sol_update_always; 1955 PetscFunctionReturn(0); 1956 } 1957 1958 #undef __FUNCT__ 1959 #define __FUNCT__ "SNESGetFunction" 1960 /*@C 1961 SNESGetFunction - Returns the vector where the function is stored. 1962 1963 Not Collective, but Vec is parallel if SNES is parallel 1964 1965 Input Parameter: 1966 . snes - the SNES context 1967 1968 Output Parameter: 1969 + r - the function (or PETSC_NULL) 1970 . func - the function (or PETSC_NULL) 1971 - ctx - the function context (or PETSC_NULL) 1972 1973 Level: advanced 1974 1975 .keywords: SNES, nonlinear, get, function 1976 1977 .seealso: SNESSetFunction(), SNESGetSolution() 1978 @*/ 1979 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx) 1980 { 1981 PetscFunctionBegin; 1982 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1983 if (r) *r = snes->vec_func_always; 1984 if (func) *func = snes->computefunction; 1985 if (ctx) *ctx = snes->funP; 1986 PetscFunctionReturn(0); 1987 } 1988 1989 #undef __FUNCT__ 1990 #define __FUNCT__ "SNESSetOptionsPrefix" 1991 /*@C 1992 SNESSetOptionsPrefix - Sets the prefix used for searching for all 1993 SNES options in the database. 1994 1995 Collective on SNES 1996 1997 Input Parameter: 1998 + snes - the SNES context 1999 - prefix - the prefix to prepend to all option names 2000 2001 Notes: 2002 A hyphen (-) must NOT be given at the beginning of the prefix name. 2003 The first character of all runtime options is AUTOMATICALLY the hyphen. 2004 2005 Level: advanced 2006 2007 .keywords: SNES, set, options, prefix, database 2008 2009 .seealso: SNESSetFromOptions() 2010 @*/ 2011 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetOptionsPrefix(SNES snes,const char prefix[]) 2012 { 2013 PetscErrorCode ierr; 2014 2015 PetscFunctionBegin; 2016 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2017 ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2018 ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 2019 PetscFunctionReturn(0); 2020 } 2021 2022 #undef __FUNCT__ 2023 #define __FUNCT__ "SNESAppendOptionsPrefix" 2024 /*@C 2025 SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 2026 SNES options in the database. 2027 2028 Collective on SNES 2029 2030 Input Parameters: 2031 + snes - the SNES context 2032 - prefix - the prefix to prepend to all option names 2033 2034 Notes: 2035 A hyphen (-) must NOT be given at the beginning of the prefix name. 2036 The first character of all runtime options is AUTOMATICALLY the hyphen. 2037 2038 Level: advanced 2039 2040 .keywords: SNES, append, options, prefix, database 2041 2042 .seealso: SNESGetOptionsPrefix() 2043 @*/ 2044 PetscErrorCode PETSCSNES_DLLEXPORT SNESAppendOptionsPrefix(SNES snes,const char prefix[]) 2045 { 2046 PetscErrorCode ierr; 2047 2048 PetscFunctionBegin; 2049 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2050 ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2051 ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 2052 PetscFunctionReturn(0); 2053 } 2054 2055 #undef __FUNCT__ 2056 #define __FUNCT__ "SNESGetOptionsPrefix" 2057 /*@C 2058 SNESGetOptionsPrefix - Sets the prefix used for searching for all 2059 SNES options in the database. 2060 2061 Not Collective 2062 2063 Input Parameter: 2064 . snes - the SNES context 2065 2066 Output Parameter: 2067 . prefix - pointer to the prefix string used 2068 2069 Notes: On the fortran side, the user should pass in a string 'prifix' of 2070 sufficient length to hold the prefix. 2071 2072 Level: advanced 2073 2074 .keywords: SNES, get, options, prefix, database 2075 2076 .seealso: SNESAppendOptionsPrefix() 2077 @*/ 2078 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetOptionsPrefix(SNES snes,const char *prefix[]) 2079 { 2080 PetscErrorCode ierr; 2081 2082 PetscFunctionBegin; 2083 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2084 ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2085 PetscFunctionReturn(0); 2086 } 2087 2088 2089 #undef __FUNCT__ 2090 #define __FUNCT__ "SNESRegister" 2091 /*@C 2092 SNESRegister - See SNESRegisterDynamic() 2093 2094 Level: advanced 2095 @*/ 2096 PetscErrorCode PETSCSNES_DLLEXPORT SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES)) 2097 { 2098 char fullname[PETSC_MAX_PATH_LEN]; 2099 PetscErrorCode ierr; 2100 2101 PetscFunctionBegin; 2102 ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 2103 ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 2104 PetscFunctionReturn(0); 2105 } 2106 2107 #undef __FUNCT__ 2108 #define __FUNCT__ "SNESTestLocalMin" 2109 PetscErrorCode PETSCSNES_DLLEXPORT SNESTestLocalMin(SNES snes) 2110 { 2111 PetscErrorCode ierr; 2112 PetscInt N,i,j; 2113 Vec u,uh,fh; 2114 PetscScalar value; 2115 PetscReal norm; 2116 2117 PetscFunctionBegin; 2118 ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr); 2119 ierr = VecDuplicate(u,&uh);CHKERRQ(ierr); 2120 ierr = VecDuplicate(u,&fh);CHKERRQ(ierr); 2121 2122 /* currently only works for sequential */ 2123 ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n"); 2124 ierr = VecGetSize(u,&N);CHKERRQ(ierr); 2125 for (i=0; i<N; i++) { 2126 ierr = VecCopy(u,uh);CHKERRQ(ierr); 2127 ierr = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr); 2128 for (j=-10; j<11; j++) { 2129 value = PetscSign(j)*exp(PetscAbs(j)-10.0); 2130 ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2131 ierr = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr); 2132 ierr = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr); 2133 ierr = PetscPrintf(PETSC_COMM_WORLD," j norm %D %18.16e\n",j,norm);CHKERRQ(ierr); 2134 value = -value; 2135 ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2136 } 2137 } 2138 ierr = VecDestroy(uh);CHKERRQ(ierr); 2139 ierr = VecDestroy(fh);CHKERRQ(ierr); 2140 PetscFunctionReturn(0); 2141 } 2142