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