1 #define PETSCSNES_DLL 2 3 #include "include/private/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; 10 PetscLogEvent SNES_Solve, SNES_LineSearch, SNES_FunctionEval, SNES_JacobianEval; 11 12 #undef __FUNCT__ 13 #define __FUNCT__ "SNESSetFunctionDomainError" 14 /*@ 15 SNESSetFunctionDomainError - tells SNES that the input vector to your FormFunction is not 16 in the functions domain. For example, negative pressure. 17 18 Collective on SNES 19 20 Input Parameters: 21 . SNES - the SNES context 22 23 Level: advanced 24 25 .keywords: SNES, view 26 27 .seealso: SNESCreate(), SNESSetFunction() 28 @*/ 29 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFunctionDomainError(SNES snes) 30 { 31 PetscFunctionBegin; 32 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 33 snes->domainerror = PETSC_TRUE; 34 PetscFunctionReturn(0); 35 } 36 37 #undef __FUNCT__ 38 #define __FUNCT__ "SNESView" 39 /*@C 40 SNESView - Prints the SNES data structure. 41 42 Collective on SNES 43 44 Input Parameters: 45 + SNES - the SNES context 46 - viewer - visualization context 47 48 Options Database Key: 49 . -snes_view - Calls SNESView() at end of SNESSolve() 50 51 Notes: 52 The available visualization contexts include 53 + PETSC_VIEWER_STDOUT_SELF - standard output (default) 54 - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 55 output where only the first processor opens 56 the file. All other processors send their 57 data to the first processor to print. 58 59 The user can open an alternative visualization context with 60 PetscViewerASCIIOpen() - output to a specified file. 61 62 Level: beginner 63 64 .keywords: SNES, view 65 66 .seealso: PetscViewerASCIIOpen() 67 @*/ 68 PetscErrorCode PETSCSNES_DLLEXPORT SNESView(SNES snes,PetscViewer viewer) 69 { 70 SNESKSPEW *kctx; 71 PetscErrorCode ierr; 72 KSP ksp; 73 const SNESType type; 74 PetscTruth iascii,isstring; 75 76 PetscFunctionBegin; 77 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 78 if (!viewer) { 79 ierr = PetscViewerASCIIGetStdout(((PetscObject)snes)->comm,&viewer);CHKERRQ(ierr); 80 } 81 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2); 82 PetscCheckSameComm(snes,1,viewer,2); 83 84 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 85 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 86 if (iascii) { 87 if (((PetscObject)snes)->prefix) { 88 ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:(%s)\n",((PetscObject)snes)->prefix);CHKERRQ(ierr); 89 } else { 90 ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:\n");CHKERRQ(ierr); 91 } 92 ierr = SNESGetType(snes,&type);CHKERRQ(ierr); 93 if (type) { 94 ierr = PetscViewerASCIIPrintf(viewer," type: %s\n",type);CHKERRQ(ierr); 95 } else { 96 ierr = PetscViewerASCIIPrintf(viewer," type: not set yet\n");CHKERRQ(ierr); 97 } 98 if (snes->ops->view) { 99 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 100 ierr = (*snes->ops->view)(snes,viewer);CHKERRQ(ierr); 101 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 102 } 103 ierr = PetscViewerASCIIPrintf(viewer," maximum iterations=%D, maximum function evaluations=%D\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr); 104 ierr = PetscViewerASCIIPrintf(viewer," tolerances: relative=%G, absolute=%G, solution=%G\n", 105 snes->rtol,snes->abstol,snes->xtol);CHKERRQ(ierr); 106 ierr = PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",snes->linear_its);CHKERRQ(ierr); 107 ierr = PetscViewerASCIIPrintf(viewer," total number of function evaluations=%D\n",snes->nfuncs);CHKERRQ(ierr); 108 if (snes->ksp_ewconv) { 109 kctx = (SNESKSPEW *)snes->kspconvctx; 110 if (kctx) { 111 ierr = PetscViewerASCIIPrintf(viewer," Eisenstat-Walker computation of KSP relative tolerance (version %D)\n",kctx->version);CHKERRQ(ierr); 112 ierr = PetscViewerASCIIPrintf(viewer," rtol_0=%G, rtol_max=%G, threshold=%G\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr); 113 ierr = PetscViewerASCIIPrintf(viewer," gamma=%G, alpha=%G, alpha2=%G\n",kctx->gamma,kctx->alpha,kctx->alpha2);CHKERRQ(ierr); 114 } 115 } 116 } else if (isstring) { 117 ierr = SNESGetType(snes,&type);CHKERRQ(ierr); 118 ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr); 119 } 120 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 121 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 122 ierr = KSPView(ksp,viewer);CHKERRQ(ierr); 123 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 124 PetscFunctionReturn(0); 125 } 126 127 /* 128 We retain a list of functions that also take SNES command 129 line options. These are called at the end SNESSetFromOptions() 130 */ 131 #define MAXSETFROMOPTIONS 5 132 static PetscInt numberofsetfromoptions; 133 static PetscErrorCode (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES); 134 135 #undef __FUNCT__ 136 #define __FUNCT__ "SNESAddOptionsChecker" 137 /*@C 138 SNESAddOptionsChecker - Adds an additional function to check for SNES options. 139 140 Not Collective 141 142 Input Parameter: 143 . snescheck - function that checks for options 144 145 Level: developer 146 147 .seealso: SNESSetFromOptions() 148 @*/ 149 PetscErrorCode PETSCSNES_DLLEXPORT SNESAddOptionsChecker(PetscErrorCode (*snescheck)(SNES)) 150 { 151 PetscFunctionBegin; 152 if (numberofsetfromoptions >= MAXSETFROMOPTIONS) { 153 SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %D allowed", MAXSETFROMOPTIONS); 154 } 155 othersetfromoptions[numberofsetfromoptions++] = snescheck; 156 PetscFunctionReturn(0); 157 } 158 159 #undef __FUNCT__ 160 #define __FUNCT__ "SNESSetFromOptions" 161 /*@ 162 SNESSetFromOptions - Sets various SNES and KSP parameters from user options. 163 164 Collective on SNES 165 166 Input Parameter: 167 . snes - the SNES context 168 169 Options Database Keys: 170 + -snes_type <type> - ls, tr, umls, umtr, test 171 . -snes_stol - convergence tolerance in terms of the norm 172 of the change in the solution between steps 173 . -snes_atol <abstol> - absolute tolerance of residual norm 174 . -snes_rtol <rtol> - relative decrease in tolerance norm from initial 175 . -snes_max_it <max_it> - maximum number of iterations 176 . -snes_max_funcs <max_funcs> - maximum number of function evaluations 177 . -snes_max_fail <max_fail> - maximum number of failures 178 . -snes_max_linear_solve_fail - number of linear solver failures before SNESSolve() stops 179 . -snes_lag_preconditioner <lag> - how often preconditioner is rebuilt (use -1 to never rebuild) 180 . -snes_lag_jacobian <lag> - how often Jacobian is rebuilt (use -1 to never rebuild) 181 . -snes_trtol <trtol> - trust region tolerance 182 . -snes_no_convergence_test - skip convergence test in nonlinear 183 solver; hence iterations will continue until max_it 184 or some other criterion is reached. Saves expense 185 of convergence test 186 . -snes_monitor <optional filename> - prints residual norm at each iteration. if no 187 filename given prints to stdout 188 . -snes_monitor_solution - plots solution at each iteration 189 . -snes_monitor_residual - plots residual (not its norm) at each iteration 190 . -snes_monitor_solution_update - plots update to solution at each iteration 191 . -snes_monitor_draw - plots residual norm at each iteration 192 . -snes_fd - use finite differences to compute Jacobian; very slow, only for testing 193 . -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration 194 - -snes_converged_reason - print the reason for convergence/divergence after each solve 195 196 Options Database for Eisenstat-Walker method: 197 + -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence 198 . -snes_ksp_ew_version ver - version of Eisenstat-Walker method 199 . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0 200 . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax 201 . -snes_ksp_ew_gamma <gamma> - Sets gamma 202 . -snes_ksp_ew_alpha <alpha> - Sets alpha 203 . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2 204 - -snes_ksp_ew_threshold <threshold> - Sets threshold 205 206 Notes: 207 To see all options, run your program with the -help option or consult 208 the users manual. 209 210 Level: beginner 211 212 .keywords: SNES, nonlinear, set, options, database 213 214 .seealso: SNESSetOptionsPrefix() 215 @*/ 216 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFromOptions(SNES snes) 217 { 218 PetscTruth flg; 219 PetscInt i,indx,lag; 220 const char *deft = SNESLS; 221 const char *convtests[] = {"default","skip"}; 222 SNESKSPEW *kctx = NULL; 223 char type[256], monfilename[PETSC_MAX_PATH_LEN]; 224 PetscViewerASCIIMonitor monviewer; 225 PetscErrorCode ierr; 226 227 PetscFunctionBegin; 228 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 229 230 ierr = PetscOptionsBegin(((PetscObject)snes)->comm,((PetscObject)snes)->prefix,"Nonlinear solver (SNES) options","SNES");CHKERRQ(ierr); 231 if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 232 if (((PetscObject)snes)->type_name) { deft = ((PetscObject)snes)->type_name; } 233 ierr = PetscOptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr); 234 if (flg) { 235 ierr = SNESSetType(snes,type);CHKERRQ(ierr); 236 } else if (!((PetscObject)snes)->type_name) { 237 ierr = SNESSetType(snes,deft);CHKERRQ(ierr); 238 } 239 ierr = PetscOptionsName("-snes_view","Print detailed information on solver used","SNESView",0);CHKERRQ(ierr); 240 241 ierr = PetscOptionsReal("-snes_stol","Stop if step length less then","SNESSetTolerances",snes->xtol,&snes->xtol,0);CHKERRQ(ierr); 242 ierr = PetscOptionsReal("-snes_atol","Stop if function norm less then","SNESSetTolerances",snes->abstol,&snes->abstol,0);CHKERRQ(ierr); 243 244 ierr = PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less then","SNESSetTolerances",snes->rtol,&snes->rtol,0);CHKERRQ(ierr); 245 ierr = PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,PETSC_NULL);CHKERRQ(ierr); 246 ierr = PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr); 247 ierr = PetscOptionsInt("-snes_max_fail","Maximum failures","SNESSetTolerances",snes->maxFailures,&snes->maxFailures,PETSC_NULL);CHKERRQ(ierr); 248 ierr = PetscOptionsInt("-snes_max_linear_solve_fail","Maximum failures in linear solves allowed","SNESSetMaxLinearSolveFailures",snes->maxLinearSolveFailures,&snes->maxLinearSolveFailures,PETSC_NULL);CHKERRQ(ierr); 249 250 ierr = PetscOptionsInt("-snes_lag_preconditioner","How often to rebuild preconditioner","SNESSetLagPreconditioner",snes->lagpreconditioner,&lag,&flg);CHKERRQ(ierr); 251 if (flg) { 252 ierr = SNESSetLagPreconditioner(snes,lag);CHKERRQ(ierr); 253 } 254 ierr = PetscOptionsInt("-snes_lag_jacobian","How often to rebuild Jacobian","SNESSetLagJacobian",snes->lagjacobian,&lag,&flg);CHKERRQ(ierr); 255 if (flg) { 256 ierr = SNESSetLagJacobian(snes,lag);CHKERRQ(ierr); 257 } 258 259 ierr = PetscOptionsEList("-snes_convergence_test","Convergence test","SNESSetConvergenceTest",convtests,2,"default",&indx,&flg);CHKERRQ(ierr); 260 if (flg) { 261 switch (indx) { 262 case 0: ierr = SNESSetConvergenceTest(snes,SNESDefaultConverged,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); break; 263 case 1: ierr = SNESSetConvergenceTest(snes,SNESSkipConverged,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); break; 264 } 265 } 266 267 ierr = PetscOptionsName("-snes_converged_reason","Print reason for converged or diverged","SNESSolve",&flg);CHKERRQ(ierr); 268 if (flg) { 269 snes->printreason = PETSC_TRUE; 270 } 271 272 kctx = (SNESKSPEW *)snes->kspconvctx; 273 274 ierr = PetscOptionsTruth("-snes_ksp_ew","Use Eisentat-Walker linear system convergence test","SNESKSPSetUseEW",snes->ksp_ewconv,&snes->ksp_ewconv,PETSC_NULL);CHKERRQ(ierr); 275 276 ierr = PetscOptionsInt("-snes_ksp_ew_version","Version 1, 2 or 3","SNESKSPSetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr); 277 ierr = PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNESKSPSetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr); 278 ierr = PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNESKSPSetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr); 279 ierr = PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNESKSPSetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr); 280 ierr = PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNESKSPSetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr); 281 ierr = PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNESKSPSetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr); 282 ierr = PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNESKSPSetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr); 283 284 ierr = PetscOptionsName("-snes_monitor_cancel","Remove all monitors","SNESMonitorCancel",&flg);CHKERRQ(ierr); 285 if (flg) {ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);} 286 287 ierr = PetscOptionsString("-snes_monitor","Monitor norm of function","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 288 if (flg) { 289 ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,monfilename,((PetscObject)snes)->tablevel,&monviewer);CHKERRQ(ierr); 290 ierr = SNESMonitorSet(snes,SNESMonitorDefault,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 291 } 292 293 ierr = PetscOptionsString("-snes_ratiomonitor","Monitor ratios of norms of function","SNESMonitorSetRatio","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 294 if (flg) { 295 ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,monfilename,((PetscObject)snes)->tablevel,&monviewer);CHKERRQ(ierr); 296 ierr = SNESMonitorSetRatio(snes,monviewer);CHKERRQ(ierr); 297 } 298 299 ierr = PetscOptionsString("-snes_monitor_short","Monitor norm of function (fewer digits)","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 300 if (flg) { 301 ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,monfilename,((PetscObject)snes)->tablevel,&monviewer);CHKERRQ(ierr); 302 ierr = SNESMonitorSet(snes,SNESMonitorDefaultShort,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 303 } 304 305 ierr = PetscOptionsName("-snes_monitor_solution","Plot solution at each iteration","SNESMonitorSolution",&flg);CHKERRQ(ierr); 306 if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorSolution,0,0);CHKERRQ(ierr);} 307 ierr = PetscOptionsName("-snes_monitor_solution_update","Plot correction at each iteration","SNESMonitorSolutionUpdate",&flg);CHKERRQ(ierr); 308 if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorSolutionUpdate,0,0);CHKERRQ(ierr);} 309 ierr = PetscOptionsName("-snes_monitor_residual","Plot residual at each iteration","SNESMonitorResidual",&flg);CHKERRQ(ierr); 310 if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorResidual,0,0);CHKERRQ(ierr);} 311 ierr = PetscOptionsName("-snes_monitor_draw","Plot function norm at each iteration","SNESMonitorLG",&flg);CHKERRQ(ierr); 312 if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);} 313 314 ierr = PetscOptionsName("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",&flg);CHKERRQ(ierr); 315 if (flg) { 316 ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr); 317 ierr = PetscInfo(snes,"Setting default finite difference Jacobian matrix\n");CHKERRQ(ierr); 318 } 319 320 for(i = 0; i < numberofsetfromoptions; i++) { 321 ierr = (*othersetfromoptions[i])(snes);CHKERRQ(ierr); 322 } 323 324 if (snes->ops->setfromoptions) { 325 ierr = (*snes->ops->setfromoptions)(snes);CHKERRQ(ierr); 326 } 327 ierr = PetscOptionsEnd();CHKERRQ(ierr); 328 329 330 if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);} 331 ierr = KSPSetFromOptions(snes->ksp);CHKERRQ(ierr); 332 333 PetscFunctionReturn(0); 334 } 335 336 337 #undef __FUNCT__ 338 #define __FUNCT__ "SNESSetApplicationContext" 339 /*@ 340 SNESSetApplicationContext - Sets the optional user-defined context for 341 the nonlinear solvers. 342 343 Collective on SNES 344 345 Input Parameters: 346 + snes - the SNES context 347 - usrP - optional user context 348 349 Level: intermediate 350 351 .keywords: SNES, nonlinear, set, application, context 352 353 .seealso: SNESGetApplicationContext() 354 @*/ 355 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetApplicationContext(SNES snes,void *usrP) 356 { 357 PetscFunctionBegin; 358 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 359 snes->user = usrP; 360 PetscFunctionReturn(0); 361 } 362 363 #undef __FUNCT__ 364 #define __FUNCT__ "SNESGetApplicationContext" 365 /*@C 366 SNESGetApplicationContext - Gets the user-defined context for the 367 nonlinear solvers. 368 369 Not Collective 370 371 Input Parameter: 372 . snes - SNES context 373 374 Output Parameter: 375 . usrP - user context 376 377 Level: intermediate 378 379 .keywords: SNES, nonlinear, get, application, context 380 381 .seealso: SNESSetApplicationContext() 382 @*/ 383 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetApplicationContext(SNES snes,void **usrP) 384 { 385 PetscFunctionBegin; 386 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 387 *usrP = snes->user; 388 PetscFunctionReturn(0); 389 } 390 391 #undef __FUNCT__ 392 #define __FUNCT__ "SNESGetIterationNumber" 393 /*@ 394 SNESGetIterationNumber - Gets the number of nonlinear iterations completed 395 at this time. 396 397 Not Collective 398 399 Input Parameter: 400 . snes - SNES context 401 402 Output Parameter: 403 . iter - iteration number 404 405 Notes: 406 For example, during the computation of iteration 2 this would return 1. 407 408 This is useful for using lagged Jacobians (where one does not recompute the 409 Jacobian at each SNES iteration). For example, the code 410 .vb 411 ierr = SNESGetIterationNumber(snes,&it); 412 if (!(it % 2)) { 413 [compute Jacobian here] 414 } 415 .ve 416 can be used in your ComputeJacobian() function to cause the Jacobian to be 417 recomputed every second SNES iteration. 418 419 Level: intermediate 420 421 .keywords: SNES, nonlinear, get, iteration, number, 422 423 .seealso: SNESGetFunctionNorm(), SNESGetLinearSolveIterations() 424 @*/ 425 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetIterationNumber(SNES snes,PetscInt* iter) 426 { 427 PetscFunctionBegin; 428 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 429 PetscValidIntPointer(iter,2); 430 *iter = snes->iter; 431 PetscFunctionReturn(0); 432 } 433 434 #undef __FUNCT__ 435 #define __FUNCT__ "SNESGetFunctionNorm" 436 /*@ 437 SNESGetFunctionNorm - Gets the norm of the current function that was set 438 with SNESSSetFunction(). 439 440 Collective on SNES 441 442 Input Parameter: 443 . snes - SNES context 444 445 Output Parameter: 446 . fnorm - 2-norm of function 447 448 Level: intermediate 449 450 .keywords: SNES, nonlinear, get, function, norm 451 452 .seealso: SNESGetFunction(), SNESGetIterationNumber(), SNESGetLinearSolveIterations() 453 @*/ 454 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunctionNorm(SNES snes,PetscReal *fnorm) 455 { 456 PetscFunctionBegin; 457 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 458 PetscValidScalarPointer(fnorm,2); 459 *fnorm = snes->norm; 460 PetscFunctionReturn(0); 461 } 462 463 #undef __FUNCT__ 464 #define __FUNCT__ "SNESGetNonlinearStepFailures" 465 /*@ 466 SNESGetNonlinearStepFailures - Gets the number of unsuccessful steps 467 attempted by the nonlinear solver. 468 469 Not Collective 470 471 Input Parameter: 472 . snes - SNES context 473 474 Output Parameter: 475 . nfails - number of unsuccessful steps attempted 476 477 Notes: 478 This counter is reset to zero for each successive call to SNESSolve(). 479 480 Level: intermediate 481 482 .keywords: SNES, nonlinear, get, number, unsuccessful, steps 483 484 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), 485 SNESSetMaxNonlinearStepFailures(), SNESGetMaxNonlinearStepFailures() 486 @*/ 487 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetNonlinearStepFailures(SNES snes,PetscInt* nfails) 488 { 489 PetscFunctionBegin; 490 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 491 PetscValidIntPointer(nfails,2); 492 *nfails = snes->numFailures; 493 PetscFunctionReturn(0); 494 } 495 496 #undef __FUNCT__ 497 #define __FUNCT__ "SNESSetMaxNonlinearStepFailures" 498 /*@ 499 SNESSetMaxNonlinearStepFailures - Sets the maximum number of unsuccessful steps 500 attempted by the nonlinear solver before it gives up. 501 502 Not Collective 503 504 Input Parameters: 505 + snes - SNES context 506 - maxFails - maximum of unsuccessful steps 507 508 Level: intermediate 509 510 .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps 511 512 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), 513 SNESGetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures() 514 @*/ 515 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMaxNonlinearStepFailures(SNES snes, PetscInt maxFails) 516 { 517 PetscFunctionBegin; 518 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 519 snes->maxFailures = maxFails; 520 PetscFunctionReturn(0); 521 } 522 523 #undef __FUNCT__ 524 #define __FUNCT__ "SNESGetMaxNonlinearStepFailures" 525 /*@ 526 SNESGetMaxNonlinearStepFailures - Gets the maximum number of unsuccessful steps 527 attempted by the nonlinear solver before it gives up. 528 529 Not Collective 530 531 Input Parameter: 532 . snes - SNES context 533 534 Output Parameter: 535 . maxFails - maximum of unsuccessful steps 536 537 Level: intermediate 538 539 .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps 540 541 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), 542 SNESSetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures() 543 544 @*/ 545 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetMaxNonlinearStepFailures(SNES snes, PetscInt *maxFails) 546 { 547 PetscFunctionBegin; 548 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 549 PetscValidIntPointer(maxFails,2); 550 *maxFails = snes->maxFailures; 551 PetscFunctionReturn(0); 552 } 553 554 #undef __FUNCT__ 555 #define __FUNCT__ "SNESGetNumberFunctionEvals" 556 /*@ 557 SNESGetNumberFunctionEvals - Gets the number of user provided function evaluations 558 done by SNES. 559 560 Not Collective 561 562 Input Parameter: 563 . snes - SNES context 564 565 Output Parameter: 566 . nfuncs - number of evaluations 567 568 Level: intermediate 569 570 .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps 571 572 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures() 573 @*/ 574 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetNumberFunctionEvals(SNES snes, PetscInt *nfuncs) 575 { 576 PetscFunctionBegin; 577 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 578 PetscValidIntPointer(nfuncs,2); 579 *nfuncs = snes->nfuncs; 580 PetscFunctionReturn(0); 581 } 582 583 #undef __FUNCT__ 584 #define __FUNCT__ "SNESGetLinearSolveFailures" 585 /*@ 586 SNESGetLinearSolveFailures - Gets the number of failed (non-converged) 587 linear solvers. 588 589 Not Collective 590 591 Input Parameter: 592 . snes - SNES context 593 594 Output Parameter: 595 . nfails - number of failed solves 596 597 Notes: 598 This counter is reset to zero for each successive call to SNESSolve(). 599 600 Level: intermediate 601 602 .keywords: SNES, nonlinear, get, number, unsuccessful, steps 603 604 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures() 605 @*/ 606 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLinearSolveFailures(SNES snes,PetscInt* nfails) 607 { 608 PetscFunctionBegin; 609 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 610 PetscValidIntPointer(nfails,2); 611 *nfails = snes->numLinearSolveFailures; 612 PetscFunctionReturn(0); 613 } 614 615 #undef __FUNCT__ 616 #define __FUNCT__ "SNESSetMaxLinearSolveFailures" 617 /*@ 618 SNESSetMaxLinearSolveFailures - the number of failed linear solve attempts 619 allowed before SNES returns with a diverged reason of SNES_DIVERGED_LINEAR_SOLVE 620 621 Collective on SNES 622 623 Input Parameters: 624 + snes - SNES context 625 - maxFails - maximum allowed linear solve failures 626 627 Level: intermediate 628 629 Notes: By default this is 0; that is SNES returns on the first failed linear solve 630 631 .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps 632 633 .seealso: SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations() 634 @*/ 635 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMaxLinearSolveFailures(SNES snes, PetscInt maxFails) 636 { 637 PetscFunctionBegin; 638 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 639 snes->maxLinearSolveFailures = maxFails; 640 PetscFunctionReturn(0); 641 } 642 643 #undef __FUNCT__ 644 #define __FUNCT__ "SNESGetMaxLinearSolveFailures" 645 /*@ 646 SNESGetMaxLinearSolveFailures - gets the maximum number of linear solve failures that 647 are allowed before SNES terminates 648 649 Not Collective 650 651 Input Parameter: 652 . snes - SNES context 653 654 Output Parameter: 655 . maxFails - maximum of unsuccessful solves allowed 656 657 Level: intermediate 658 659 Notes: By default this is 1; that is SNES returns on the first failed linear solve 660 661 .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps 662 663 .seealso: SNESGetLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), 664 @*/ 665 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetMaxLinearSolveFailures(SNES snes, PetscInt *maxFails) 666 { 667 PetscFunctionBegin; 668 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 669 PetscValidIntPointer(maxFails,2); 670 *maxFails = snes->maxLinearSolveFailures; 671 PetscFunctionReturn(0); 672 } 673 674 #undef __FUNCT__ 675 #define __FUNCT__ "SNESGetLinearSolveIterations" 676 /*@ 677 SNESGetLinearSolveIterations - Gets the total number of linear iterations 678 used by the nonlinear solver. 679 680 Not Collective 681 682 Input Parameter: 683 . snes - SNES context 684 685 Output Parameter: 686 . lits - number of linear iterations 687 688 Notes: 689 This counter is reset to zero for each successive call to SNESSolve(). 690 691 Level: intermediate 692 693 .keywords: SNES, nonlinear, get, number, linear, iterations 694 695 .seealso: SNESGetIterationNumber(), SNESGetFunctionNorm()S, NESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures() 696 @*/ 697 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLinearSolveIterations(SNES snes,PetscInt* lits) 698 { 699 PetscFunctionBegin; 700 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 701 PetscValidIntPointer(lits,2); 702 *lits = snes->linear_its; 703 PetscFunctionReturn(0); 704 } 705 706 #undef __FUNCT__ 707 #define __FUNCT__ "SNESGetKSP" 708 /*@ 709 SNESGetKSP - Returns the KSP context for a SNES solver. 710 711 Not Collective, but if SNES object is parallel, then KSP object is parallel 712 713 Input Parameter: 714 . snes - the SNES context 715 716 Output Parameter: 717 . ksp - the KSP context 718 719 Notes: 720 The user can then directly manipulate the KSP context to set various 721 options, etc. Likewise, the user can then extract and manipulate the 722 PC contexts as well. 723 724 Level: beginner 725 726 .keywords: SNES, nonlinear, get, KSP, context 727 728 .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP() 729 @*/ 730 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetKSP(SNES snes,KSP *ksp) 731 { 732 PetscErrorCode ierr; 733 734 PetscFunctionBegin; 735 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 736 PetscValidPointer(ksp,2); 737 738 if (!snes->ksp) { 739 ierr = KSPCreate(((PetscObject)snes)->comm,&snes->ksp);CHKERRQ(ierr); 740 ierr = PetscObjectIncrementTabLevel((PetscObject)snes->ksp,(PetscObject)snes,1);CHKERRQ(ierr); 741 ierr = PetscLogObjectParent(snes,snes->ksp);CHKERRQ(ierr); 742 } 743 *ksp = snes->ksp; 744 PetscFunctionReturn(0); 745 } 746 747 #undef __FUNCT__ 748 #define __FUNCT__ "SNESSetKSP" 749 /*@ 750 SNESSetKSP - Sets a KSP context for the SNES object to use 751 752 Not Collective, but the SNES and KSP objects must live on the same MPI_Comm 753 754 Input Parameters: 755 + snes - the SNES context 756 - ksp - the KSP context 757 758 Notes: 759 The SNES object already has its KSP object, you can obtain with SNESGetKSP() 760 so this routine is rarely needed. 761 762 The KSP object that is already in the SNES object has its reference count 763 decreased by one. 764 765 Level: developer 766 767 .keywords: SNES, nonlinear, get, KSP, context 768 769 .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP() 770 @*/ 771 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetKSP(SNES snes,KSP ksp) 772 { 773 PetscErrorCode ierr; 774 775 PetscFunctionBegin; 776 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 777 PetscValidHeaderSpecific(ksp,KSP_COOKIE,2); 778 PetscCheckSameComm(snes,1,ksp,2); 779 ierr = PetscObjectReference((PetscObject)ksp);CHKERRQ(ierr); 780 if (snes->ksp) {ierr = PetscObjectDereference((PetscObject)snes->ksp);CHKERRQ(ierr);} 781 snes->ksp = ksp; 782 PetscFunctionReturn(0); 783 } 784 785 #if 0 786 #undef __FUNCT__ 787 #define __FUNCT__ "SNESPublish_Petsc" 788 static PetscErrorCode SNESPublish_Petsc(PetscObject obj) 789 { 790 PetscFunctionBegin; 791 PetscFunctionReturn(0); 792 } 793 #endif 794 795 /* -----------------------------------------------------------*/ 796 #undef __FUNCT__ 797 #define __FUNCT__ "SNESCreate" 798 /*@ 799 SNESCreate - Creates a nonlinear solver context. 800 801 Collective on MPI_Comm 802 803 Input Parameters: 804 . comm - MPI communicator 805 806 Output Parameter: 807 . outsnes - the new SNES context 808 809 Options Database Keys: 810 + -snes_mf - Activates default matrix-free Jacobian-vector products, 811 and no preconditioning matrix 812 . -snes_mf_operator - Activates default matrix-free Jacobian-vector 813 products, and a user-provided preconditioning matrix 814 as set by SNESSetJacobian() 815 - -snes_fd - Uses (slow!) finite differences to compute Jacobian 816 817 Level: beginner 818 819 .keywords: SNES, nonlinear, create, context 820 821 .seealso: SNESSolve(), SNESDestroy(), SNES, SNESSetLagPreconditioner() 822 823 @*/ 824 PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate(MPI_Comm comm,SNES *outsnes) 825 { 826 PetscErrorCode ierr; 827 SNES snes; 828 SNESKSPEW *kctx; 829 830 PetscFunctionBegin; 831 PetscValidPointer(outsnes,2); 832 *outsnes = PETSC_NULL; 833 #ifndef PETSC_USE_DYNAMIC_LIBRARIES 834 ierr = SNESInitializePackage(PETSC_NULL);CHKERRQ(ierr); 835 #endif 836 837 ierr = PetscHeaderCreate(snes,_p_SNES,struct _SNESOps,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView);CHKERRQ(ierr); 838 839 snes->ops->converged = SNESDefaultConverged; 840 snes->max_its = 50; 841 snes->max_funcs = 10000; 842 snes->norm = 0.0; 843 snes->rtol = 1.e-8; 844 snes->ttol = 0.0; 845 snes->abstol = 1.e-50; 846 snes->xtol = 1.e-8; 847 snes->deltatol = 1.e-12; 848 snes->nfuncs = 0; 849 snes->numFailures = 0; 850 snes->maxFailures = 1; 851 snes->linear_its = 0; 852 snes->lagjacobian = 1; 853 snes->lagpreconditioner = 1; 854 snes->numbermonitors = 0; 855 snes->data = 0; 856 snes->setupcalled = PETSC_FALSE; 857 snes->ksp_ewconv = PETSC_FALSE; 858 snes->vwork = 0; 859 snes->nwork = 0; 860 snes->conv_hist_len = 0; 861 snes->conv_hist_max = 0; 862 snes->conv_hist = PETSC_NULL; 863 snes->conv_hist_its = PETSC_NULL; 864 snes->conv_hist_reset = PETSC_TRUE; 865 snes->reason = SNES_CONVERGED_ITERATING; 866 867 snes->numLinearSolveFailures = 0; 868 snes->maxLinearSolveFailures = 1; 869 870 /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 871 ierr = PetscNewLog(snes,SNESKSPEW,&kctx);CHKERRQ(ierr); 872 snes->kspconvctx = (void*)kctx; 873 kctx->version = 2; 874 kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 875 this was too large for some test cases */ 876 kctx->rtol_last = 0; 877 kctx->rtol_max = .9; 878 kctx->gamma = 1.0; 879 kctx->alpha = .5*(1.0 + sqrt(5.0)); 880 kctx->alpha2 = kctx->alpha; 881 kctx->threshold = .1; 882 kctx->lresid_last = 0; 883 kctx->norm_last = 0; 884 885 *outsnes = snes; 886 ierr = PetscPublishAll(snes);CHKERRQ(ierr); 887 PetscFunctionReturn(0); 888 } 889 890 #undef __FUNCT__ 891 #define __FUNCT__ "SNESSetFunction" 892 /*@C 893 SNESSetFunction - Sets the function evaluation routine and function 894 vector for use by the SNES routines in solving systems of nonlinear 895 equations. 896 897 Collective on SNES 898 899 Input Parameters: 900 + snes - the SNES context 901 . r - vector to store function value 902 . func - function evaluation routine 903 - ctx - [optional] user-defined context for private data for the 904 function evaluation routine (may be PETSC_NULL) 905 906 Calling sequence of func: 907 $ func (SNES snes,Vec x,Vec f,void *ctx); 908 909 . f - function vector 910 - ctx - optional user-defined function context 911 912 Notes: 913 The Newton-like methods typically solve linear systems of the form 914 $ f'(x) x = -f(x), 915 where f'(x) denotes the Jacobian matrix and f(x) is the function. 916 917 Level: beginner 918 919 .keywords: SNES, nonlinear, set, function 920 921 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 922 @*/ 923 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx) 924 { 925 PetscErrorCode ierr; 926 PetscFunctionBegin; 927 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 928 PetscValidHeaderSpecific(r,VEC_COOKIE,2); 929 PetscCheckSameComm(snes,1,r,2); 930 ierr = PetscObjectReference((PetscObject)r);CHKERRQ(ierr); 931 if (snes->vec_func) { ierr = VecDestroy(snes->vec_func);CHKERRQ(ierr); } 932 snes->ops->computefunction = func; 933 snes->vec_func = r; 934 snes->funP = ctx; 935 PetscFunctionReturn(0); 936 } 937 938 /* --------------------------------------------------------------- */ 939 #undef __FUNCT__ 940 #define __FUNCT__ "SNESGetRhs" 941 /*@C 942 SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set 943 it assumes a zero right hand side. 944 945 Collective on SNES 946 947 Input Parameter: 948 . snes - the SNES context 949 950 Output Parameter: 951 . rhs - the right hand side vector or PETSC_NULL if the right hand side vector is null 952 953 Level: intermediate 954 955 .keywords: SNES, nonlinear, get, function, right hand side 956 957 .seealso: SNESGetSolution(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction() 958 @*/ 959 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetRhs(SNES snes,Vec *rhs) 960 { 961 PetscFunctionBegin; 962 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 963 PetscValidPointer(rhs,2); 964 *rhs = snes->vec_rhs; 965 PetscFunctionReturn(0); 966 } 967 968 #undef __FUNCT__ 969 #define __FUNCT__ "SNESComputeFunction" 970 /*@ 971 SNESComputeFunction - Calls the function that has been set with 972 SNESSetFunction(). 973 974 Collective on SNES 975 976 Input Parameters: 977 + snes - the SNES context 978 - x - input vector 979 980 Output Parameter: 981 . y - function vector, as set by SNESSetFunction() 982 983 Notes: 984 SNESComputeFunction() is typically used within nonlinear solvers 985 implementations, so most users would not generally call this routine 986 themselves. 987 988 Level: developer 989 990 .keywords: SNES, nonlinear, compute, function 991 992 .seealso: SNESSetFunction(), SNESGetFunction() 993 @*/ 994 PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeFunction(SNES snes,Vec x,Vec y) 995 { 996 PetscErrorCode ierr; 997 998 PetscFunctionBegin; 999 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1000 PetscValidHeaderSpecific(x,VEC_COOKIE,2); 1001 PetscValidHeaderSpecific(y,VEC_COOKIE,3); 1002 PetscCheckSameComm(snes,1,x,2); 1003 PetscCheckSameComm(snes,1,y,3); 1004 1005 ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 1006 if (snes->ops->computefunction) { 1007 PetscStackPush("SNES user function"); 1008 CHKMEMQ; 1009 ierr = (*snes->ops->computefunction)(snes,x,y,snes->funP); 1010 CHKMEMQ; 1011 PetscStackPop; 1012 if (PetscExceptionValue(ierr)) { 1013 PetscErrorCode pierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(pierr); 1014 } 1015 CHKERRQ(ierr); 1016 } else if (snes->vec_rhs) { 1017 ierr = MatMult(snes->jacobian, x, y);CHKERRQ(ierr); 1018 } else { 1019 SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() before SNESComputeFunction(), likely called from SNESSolve()."); 1020 } 1021 if (snes->vec_rhs) { 1022 ierr = VecAXPY(y,-1.0,snes->vec_rhs);CHKERRQ(ierr); 1023 } 1024 snes->nfuncs++; 1025 ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 1026 PetscFunctionReturn(0); 1027 } 1028 1029 #undef __FUNCT__ 1030 #define __FUNCT__ "SNESComputeJacobian" 1031 /*@ 1032 SNESComputeJacobian - Computes the Jacobian matrix that has been 1033 set with SNESSetJacobian(). 1034 1035 Collective on SNES and Mat 1036 1037 Input Parameters: 1038 + snes - the SNES context 1039 - x - input vector 1040 1041 Output Parameters: 1042 + A - Jacobian matrix 1043 . B - optional preconditioning matrix 1044 - flag - flag indicating matrix structure (one of, SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER) 1045 1046 Options Database Keys: 1047 + -snes_lag_preconditioner <lag> 1048 - -snes_lag_jacobian <lag> 1049 1050 Notes: 1051 Most users should not need to explicitly call this routine, as it 1052 is used internally within the nonlinear solvers. 1053 1054 See KSPSetOperators() for important information about setting the 1055 flag parameter. 1056 1057 Level: developer 1058 1059 .keywords: SNES, compute, Jacobian, matrix 1060 1061 .seealso: SNESSetJacobian(), KSPSetOperators(), MatStructure, SNESSetLagPreconditioner(), SNESSetLagJacobian() 1062 @*/ 1063 PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 1064 { 1065 PetscErrorCode ierr; 1066 PetscTruth flag; 1067 1068 PetscFunctionBegin; 1069 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1070 PetscValidHeaderSpecific(X,VEC_COOKIE,2); 1071 PetscValidPointer(flg,5); 1072 PetscCheckSameComm(snes,1,X,2); 1073 if (!snes->ops->computejacobian) PetscFunctionReturn(0); 1074 1075 /* make sure that MatAssemblyBegin/End() is called on A matrix if it is matrix free */ 1076 1077 if (snes->lagjacobian == -2) { 1078 snes->lagjacobian = -1; 1079 ierr = PetscInfo(snes,"Recomputing Jacobian/preconditioner because lag is -2 (means compute Jacobian, but then never again) \n");CHKERRQ(ierr); 1080 } else if (snes->lagjacobian == -1) { 1081 *flg = SAME_PRECONDITIONER; 1082 ierr = PetscInfo(snes,"Reusing Jacobian/preconditioner because lag is -1\n");CHKERRQ(ierr); 1083 ierr = PetscTypeCompare((PetscObject)*A,MATMFFD,&flag);CHKERRQ(ierr); 1084 if (flag) { 1085 ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1086 ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1087 } 1088 PetscFunctionReturn(0); 1089 } else if (snes->lagjacobian > 1 && snes->iter % snes->lagjacobian) { 1090 *flg = SAME_PRECONDITIONER; 1091 ierr = PetscInfo2(snes,"Reusing Jacobian/preconditioner because lag is %D and SNES iteration is %D\n",snes->lagjacobian,snes->iter);CHKERRQ(ierr); 1092 ierr = PetscTypeCompare((PetscObject)*A,MATMFFD,&flag);CHKERRQ(ierr); 1093 if (flag) { 1094 ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1095 ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1096 } 1097 PetscFunctionReturn(0); 1098 } 1099 1100 *flg = DIFFERENT_NONZERO_PATTERN; 1101 ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 1102 PetscStackPush("SNES user Jacobian function"); 1103 CHKMEMQ; 1104 ierr = (*snes->ops->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr); 1105 CHKMEMQ; 1106 PetscStackPop; 1107 ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 1108 1109 if (snes->lagpreconditioner == -1) { 1110 *flg = SAME_PRECONDITIONER; 1111 ierr = PetscInfo(snes,"Reusing preconditioner because lag is -1\n");CHKERRQ(ierr); 1112 } else if (snes->lagpreconditioner > 1 && snes->iter % snes->lagpreconditioner) { 1113 *flg = SAME_PRECONDITIONER; 1114 ierr = PetscInfo2(snes,"Reusing preconditioner because lag is %D and SNES iteration is %D\n",snes->lagpreconditioner,snes->iter);CHKERRQ(ierr); 1115 } 1116 1117 /* make sure user returned a correct Jacobian and preconditioner */ 1118 /* PetscValidHeaderSpecific(*A,MAT_COOKIE,3); 1119 PetscValidHeaderSpecific(*B,MAT_COOKIE,4); */ 1120 PetscFunctionReturn(0); 1121 } 1122 1123 #undef __FUNCT__ 1124 #define __FUNCT__ "SNESSetJacobian" 1125 /*@C 1126 SNESSetJacobian - Sets the function to compute Jacobian as well as the 1127 location to store the matrix. 1128 1129 Collective on SNES and Mat 1130 1131 Input Parameters: 1132 + snes - the SNES context 1133 . A - Jacobian matrix 1134 . B - preconditioner matrix (usually same as the Jacobian) 1135 . func - Jacobian evaluation routine 1136 - ctx - [optional] user-defined context for private data for the 1137 Jacobian evaluation routine (may be PETSC_NULL) 1138 1139 Calling sequence of func: 1140 $ func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 1141 1142 + x - input vector 1143 . A - Jacobian matrix 1144 . B - preconditioner matrix, usually the same as A 1145 . flag - flag indicating information about the preconditioner matrix 1146 structure (same as flag in KSPSetOperators()), one of SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER 1147 - ctx - [optional] user-defined Jacobian context 1148 1149 Notes: 1150 See KSPSetOperators() for important information about setting the flag 1151 output parameter in the routine func(). Be sure to read this information! 1152 1153 The routine func() takes Mat * as the matrix arguments rather than Mat. 1154 This allows the Jacobian evaluation routine to replace A and/or B with a 1155 completely new new matrix structure (not just different matrix elements) 1156 when appropriate, for instance, if the nonzero structure is changing 1157 throughout the global iterations. 1158 1159 If the A matrix and B matrix are different you must call MatAssemblyBegin/End() on 1160 each matrix. 1161 1162 Level: beginner 1163 1164 .keywords: SNES, nonlinear, set, Jacobian, matrix 1165 1166 .seealso: KSPSetOperators(), SNESSetFunction(), MatMFFDComputeJacobian(), SNESDefaultComputeJacobianColor(), MatStructure 1167 @*/ 1168 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 1169 { 1170 PetscErrorCode ierr; 1171 1172 PetscFunctionBegin; 1173 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1174 if (A) PetscValidHeaderSpecific(A,MAT_COOKIE,2); 1175 if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3); 1176 if (A) PetscCheckSameComm(snes,1,A,2); 1177 if (B) PetscCheckSameComm(snes,1,B,2); 1178 if (func) snes->ops->computejacobian = func; 1179 if (ctx) snes->jacP = ctx; 1180 if (A) { 1181 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 1182 if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 1183 snes->jacobian = A; 1184 } 1185 if (B) { 1186 ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 1187 if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 1188 snes->jacobian_pre = B; 1189 } 1190 PetscFunctionReturn(0); 1191 } 1192 1193 #undef __FUNCT__ 1194 #define __FUNCT__ "SNESGetJacobian" 1195 /*@C 1196 SNESGetJacobian - Returns the Jacobian matrix and optionally the user 1197 provided context for evaluating the Jacobian. 1198 1199 Not Collective, but Mat object will be parallel if SNES object is 1200 1201 Input Parameter: 1202 . snes - the nonlinear solver context 1203 1204 Output Parameters: 1205 + A - location to stash Jacobian matrix (or PETSC_NULL) 1206 . B - location to stash preconditioner matrix (or PETSC_NULL) 1207 . func - location to put Jacobian function (or PETSC_NULL) 1208 - ctx - location to stash Jacobian ctx (or PETSC_NULL) 1209 1210 Level: advanced 1211 1212 .seealso: SNESSetJacobian(), SNESComputeJacobian() 1213 @*/ 1214 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetJacobian(SNES snes,Mat *A,Mat *B,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx) 1215 { 1216 PetscFunctionBegin; 1217 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1218 if (A) *A = snes->jacobian; 1219 if (B) *B = snes->jacobian_pre; 1220 if (func) *func = snes->ops->computejacobian; 1221 if (ctx) *ctx = snes->jacP; 1222 PetscFunctionReturn(0); 1223 } 1224 1225 /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 1226 EXTERN PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*); 1227 1228 #undef __FUNCT__ 1229 #define __FUNCT__ "SNESSetUp" 1230 /*@ 1231 SNESSetUp - Sets up the internal data structures for the later use 1232 of a nonlinear solver. 1233 1234 Collective on SNES 1235 1236 Input Parameters: 1237 . snes - the SNES context 1238 1239 Notes: 1240 For basic use of the SNES solvers the user need not explicitly call 1241 SNESSetUp(), since these actions will automatically occur during 1242 the call to SNESSolve(). However, if one wishes to control this 1243 phase separately, SNESSetUp() should be called after SNESCreate() 1244 and optional routines of the form SNESSetXXX(), but before SNESSolve(). 1245 1246 Level: advanced 1247 1248 .keywords: SNES, nonlinear, setup 1249 1250 .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 1251 @*/ 1252 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUp(SNES snes) 1253 { 1254 PetscErrorCode ierr; 1255 PetscTruth flg; 1256 1257 PetscFunctionBegin; 1258 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1259 if (snes->setupcalled) PetscFunctionReturn(0); 1260 1261 if (!((PetscObject)snes)->type_name) { 1262 ierr = SNESSetType(snes,SNESLS);CHKERRQ(ierr); 1263 } 1264 1265 ierr = PetscOptionsHasName(((PetscObject)snes)->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr); 1266 /* 1267 This version replaces the user provided Jacobian matrix with a 1268 matrix-free version but still employs the user-provided preconditioner matrix 1269 */ 1270 if (flg) { 1271 Mat J; 1272 ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr); 1273 ierr = MatMFFDSetFromOptions(J);CHKERRQ(ierr); 1274 ierr = PetscInfo(snes,"Setting default matrix-free operator routines\n");CHKERRQ(ierr); 1275 ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 1276 ierr = MatDestroy(J);CHKERRQ(ierr); 1277 } 1278 1279 #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) && !defined(PETSC_USE_MAT_SINGLE) && !defined(PETSC_USE_LONG_DOUBLE) && !defined(PETSC_USE_INT) 1280 ierr = PetscOptionsHasName(((PetscObject)snes)->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr); 1281 if (flg) { 1282 Mat J; 1283 ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1284 ierr = PetscInfo(snes,"Setting default matrix-free operator routines (version 2)\n");CHKERRQ(ierr); 1285 ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 1286 ierr = MatDestroy(J);CHKERRQ(ierr); 1287 } 1288 #endif 1289 1290 ierr = PetscOptionsHasName(((PetscObject)snes)->prefix,"-snes_mf",&flg);CHKERRQ(ierr); 1291 /* 1292 This version replaces both the user-provided Jacobian and the user- 1293 provided preconditioner matrix with the default matrix free version. 1294 */ 1295 if (flg) { 1296 Mat J; 1297 KSP ksp; 1298 PC pc; 1299 /* create and set matrix-free operator */ 1300 ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr); 1301 ierr = MatMFFDSetFromOptions(J);CHKERRQ(ierr); 1302 ierr = PetscInfo(snes,"Setting default matrix-free operator routines\n");CHKERRQ(ierr); 1303 ierr = SNESSetJacobian(snes,J,J,MatMFFDComputeJacobian,snes->funP);CHKERRQ(ierr); 1304 ierr = MatDestroy(J);CHKERRQ(ierr); 1305 /* force no preconditioner */ 1306 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1307 ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 1308 ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr); 1309 if (!flg) { 1310 ierr = PetscInfo(snes,"Setting default matrix-free preconditioner routines;\nThat is no preconditioner is being used\n");CHKERRQ(ierr); 1311 ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 1312 } 1313 } 1314 1315 if (!snes->vec_func && !snes->vec_rhs) { 1316 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 1317 } 1318 if (!snes->ops->computefunction && !snes->vec_rhs) { 1319 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 1320 } 1321 if (!snes->jacobian) { 1322 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option"); 1323 } 1324 if (snes->vec_func == snes->vec_sol) { 1325 SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector"); 1326 } 1327 1328 if (snes->ops->setup) { 1329 ierr = (*snes->ops->setup)(snes);CHKERRQ(ierr); 1330 } 1331 snes->setupcalled = PETSC_TRUE; 1332 PetscFunctionReturn(0); 1333 } 1334 1335 #undef __FUNCT__ 1336 #define __FUNCT__ "SNESDestroy" 1337 /*@ 1338 SNESDestroy - Destroys the nonlinear solver context that was created 1339 with SNESCreate(). 1340 1341 Collective on SNES 1342 1343 Input Parameter: 1344 . snes - the SNES context 1345 1346 Level: beginner 1347 1348 .keywords: SNES, nonlinear, destroy 1349 1350 .seealso: SNESCreate(), SNESSolve() 1351 @*/ 1352 PetscErrorCode PETSCSNES_DLLEXPORT SNESDestroy(SNES snes) 1353 { 1354 PetscErrorCode ierr; 1355 1356 PetscFunctionBegin; 1357 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1358 if (--((PetscObject)snes)->refct > 0) PetscFunctionReturn(0); 1359 1360 /* if memory was published with AMS then destroy it */ 1361 ierr = PetscObjectDepublish(snes);CHKERRQ(ierr); 1362 if (snes->ops->destroy) {ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr);} 1363 1364 if (snes->vec_rhs) {ierr = VecDestroy(snes->vec_rhs);CHKERRQ(ierr);} 1365 if (snes->vec_sol) {ierr = VecDestroy(snes->vec_sol);CHKERRQ(ierr);} 1366 if (snes->vec_func) {ierr = VecDestroy(snes->vec_func);CHKERRQ(ierr);} 1367 if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 1368 if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 1369 if (snes->ksp) {ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr);} 1370 ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr); 1371 if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);} 1372 ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 1373 if (snes->ops->convergeddestroy) {ierr = (*snes->ops->convergeddestroy)(snes->cnvP);CHKERRQ(ierr);} 1374 ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr); 1375 PetscFunctionReturn(0); 1376 } 1377 1378 /* ----------- Routines to set solver parameters ---------- */ 1379 1380 #undef __FUNCT__ 1381 #define __FUNCT__ "SNESSetLagPreconditioner" 1382 /*@ 1383 SNESSetLagPreconditioner - Determines when the preconditioner is rebuilt in the nonlinear solve. 1384 1385 Collective on SNES 1386 1387 Input Parameters: 1388 + snes - the SNES context 1389 - lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time 1390 the Jacobian is built etc. 1391 1392 Options Database Keys: 1393 . -snes_lag_preconditioner <lag> 1394 1395 Notes: 1396 The default is 1 1397 The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 1398 If -1 is used before the very first nonlinear solve the preconditioner is still built because there is no previous preconditioner to use 1399 1400 Level: intermediate 1401 1402 .keywords: SNES, nonlinear, set, convergence, tolerances 1403 1404 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian() 1405 1406 @*/ 1407 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLagPreconditioner(SNES snes,PetscInt lag) 1408 { 1409 PetscFunctionBegin; 1410 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1411 if (lag < -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -1, 1 or greater"); 1412 if (!lag) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0"); 1413 snes->lagpreconditioner = lag; 1414 PetscFunctionReturn(0); 1415 } 1416 1417 #undef __FUNCT__ 1418 #define __FUNCT__ "SNESGetLagPreconditioner" 1419 /*@ 1420 SNESGetLagPreconditioner - Indicates how often the preconditioner is rebuilt 1421 1422 Collective on SNES 1423 1424 Input Parameter: 1425 . snes - the SNES context 1426 1427 Output Parameter: 1428 . lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time 1429 the Jacobian is built etc. 1430 1431 Options Database Keys: 1432 . -snes_lag_preconditioner <lag> 1433 1434 Notes: 1435 The default is 1 1436 The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 1437 1438 Level: intermediate 1439 1440 .keywords: SNES, nonlinear, set, convergence, tolerances 1441 1442 .seealso: SNESSetTrustRegionTolerance(), SNESSetLagPreconditioner() 1443 1444 @*/ 1445 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLagPreconditioner(SNES snes,PetscInt *lag) 1446 { 1447 PetscFunctionBegin; 1448 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1449 *lag = snes->lagpreconditioner; 1450 PetscFunctionReturn(0); 1451 } 1452 1453 #undef __FUNCT__ 1454 #define __FUNCT__ "SNESSetLagJacobian" 1455 /*@ 1456 SNESSetLagJacobian - Determines when the Jacobian is rebuilt in the nonlinear solve. See SNESSetLagPreconditioner() for determining how 1457 often the preconditioner is rebuilt. 1458 1459 Collective on SNES 1460 1461 Input Parameters: 1462 + snes - the SNES context 1463 - lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time 1464 the Jacobian is built etc. -2 means rebuild at next chance but then never again 1465 1466 Options Database Keys: 1467 . -snes_lag_jacobian <lag> 1468 1469 Notes: 1470 The default is 1 1471 The Jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 1472 If -1 is used before the very first nonlinear solve the CODE WILL FAIL! because no Jacobian is used, use -2 to indicate you want it recomputed 1473 at the next Newton step but never again (unless it is reset to another value) 1474 1475 Level: intermediate 1476 1477 .keywords: SNES, nonlinear, set, convergence, tolerances 1478 1479 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagPreconditioner(), SNESGetLagJacobian() 1480 1481 @*/ 1482 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLagJacobian(SNES snes,PetscInt lag) 1483 { 1484 PetscFunctionBegin; 1485 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1486 if (lag < -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater"); 1487 if (!lag) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0"); 1488 snes->lagjacobian = lag; 1489 PetscFunctionReturn(0); 1490 } 1491 1492 #undef __FUNCT__ 1493 #define __FUNCT__ "SNESGetLagJacobian" 1494 /*@ 1495 SNESGetLagJacobian - Indicates how often the Jacobian is rebuilt. See SNESGetLagPreconditioner() to determine when the preconditioner is rebuilt 1496 1497 Collective on SNES 1498 1499 Input Parameter: 1500 . snes - the SNES context 1501 1502 Output Parameter: 1503 . lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time 1504 the Jacobian is built etc. 1505 1506 Options Database Keys: 1507 . -snes_lag_jacobian <lag> 1508 1509 Notes: 1510 The default is 1 1511 The jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 1512 1513 Level: intermediate 1514 1515 .keywords: SNES, nonlinear, set, convergence, tolerances 1516 1517 .seealso: SNESSetTrustRegionTolerance(), SNESSetLagJacobian(), SNESSetLagPreconditioner(), SNESGetLagPreconditioner() 1518 1519 @*/ 1520 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLagJacobian(SNES snes,PetscInt *lag) 1521 { 1522 PetscFunctionBegin; 1523 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1524 *lag = snes->lagjacobian; 1525 PetscFunctionReturn(0); 1526 } 1527 1528 #undef __FUNCT__ 1529 #define __FUNCT__ "SNESSetTolerances" 1530 /*@ 1531 SNESSetTolerances - Sets various parameters used in convergence tests. 1532 1533 Collective on SNES 1534 1535 Input Parameters: 1536 + snes - the SNES context 1537 . abstol - absolute convergence tolerance 1538 . rtol - relative convergence tolerance 1539 . stol - convergence tolerance in terms of the norm 1540 of the change in the solution between steps 1541 . maxit - maximum number of iterations 1542 - maxf - maximum number of function evaluations 1543 1544 Options Database Keys: 1545 + -snes_atol <abstol> - Sets abstol 1546 . -snes_rtol <rtol> - Sets rtol 1547 . -snes_stol <stol> - Sets stol 1548 . -snes_max_it <maxit> - Sets maxit 1549 - -snes_max_funcs <maxf> - Sets maxf 1550 1551 Notes: 1552 The default maximum number of iterations is 50. 1553 The default maximum number of function evaluations is 1000. 1554 1555 Level: intermediate 1556 1557 .keywords: SNES, nonlinear, set, convergence, tolerances 1558 1559 .seealso: SNESSetTrustRegionTolerance() 1560 @*/ 1561 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf) 1562 { 1563 PetscFunctionBegin; 1564 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1565 if (abstol != PETSC_DEFAULT) snes->abstol = abstol; 1566 if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1567 if (stol != PETSC_DEFAULT) snes->xtol = stol; 1568 if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1569 if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 1570 PetscFunctionReturn(0); 1571 } 1572 1573 #undef __FUNCT__ 1574 #define __FUNCT__ "SNESGetTolerances" 1575 /*@ 1576 SNESGetTolerances - Gets various parameters used in convergence tests. 1577 1578 Not Collective 1579 1580 Input Parameters: 1581 + snes - the SNES context 1582 . atol - absolute convergence tolerance 1583 . rtol - relative convergence tolerance 1584 . stol - convergence tolerance in terms of the norm 1585 of the change in the solution between steps 1586 . maxit - maximum number of iterations 1587 - maxf - maximum number of function evaluations 1588 1589 Notes: 1590 The user can specify PETSC_NULL for any parameter that is not needed. 1591 1592 Level: intermediate 1593 1594 .keywords: SNES, nonlinear, get, convergence, tolerances 1595 1596 .seealso: SNESSetTolerances() 1597 @*/ 1598 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf) 1599 { 1600 PetscFunctionBegin; 1601 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1602 if (atol) *atol = snes->abstol; 1603 if (rtol) *rtol = snes->rtol; 1604 if (stol) *stol = snes->xtol; 1605 if (maxit) *maxit = snes->max_its; 1606 if (maxf) *maxf = snes->max_funcs; 1607 PetscFunctionReturn(0); 1608 } 1609 1610 #undef __FUNCT__ 1611 #define __FUNCT__ "SNESSetTrustRegionTolerance" 1612 /*@ 1613 SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 1614 1615 Collective on SNES 1616 1617 Input Parameters: 1618 + snes - the SNES context 1619 - tol - tolerance 1620 1621 Options Database Key: 1622 . -snes_trtol <tol> - Sets tol 1623 1624 Level: intermediate 1625 1626 .keywords: SNES, nonlinear, set, trust region, tolerance 1627 1628 .seealso: SNESSetTolerances() 1629 @*/ 1630 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTrustRegionTolerance(SNES snes,PetscReal tol) 1631 { 1632 PetscFunctionBegin; 1633 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1634 snes->deltatol = tol; 1635 PetscFunctionReturn(0); 1636 } 1637 1638 /* 1639 Duplicate the lg monitors for SNES from KSP; for some reason with 1640 dynamic libraries things don't work under Sun4 if we just use 1641 macros instead of functions 1642 */ 1643 #undef __FUNCT__ 1644 #define __FUNCT__ "SNESMonitorLG" 1645 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLG(SNES snes,PetscInt it,PetscReal norm,void *ctx) 1646 { 1647 PetscErrorCode ierr; 1648 1649 PetscFunctionBegin; 1650 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1651 ierr = KSPMonitorLG((KSP)snes,it,norm,ctx);CHKERRQ(ierr); 1652 PetscFunctionReturn(0); 1653 } 1654 1655 #undef __FUNCT__ 1656 #define __FUNCT__ "SNESMonitorLGCreate" 1657 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1658 { 1659 PetscErrorCode ierr; 1660 1661 PetscFunctionBegin; 1662 ierr = KSPMonitorLGCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr); 1663 PetscFunctionReturn(0); 1664 } 1665 1666 #undef __FUNCT__ 1667 #define __FUNCT__ "SNESMonitorLGDestroy" 1668 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGDestroy(PetscDrawLG draw) 1669 { 1670 PetscErrorCode ierr; 1671 1672 PetscFunctionBegin; 1673 ierr = KSPMonitorLGDestroy(draw);CHKERRQ(ierr); 1674 PetscFunctionReturn(0); 1675 } 1676 1677 /* ------------ Routines to set performance monitoring options ----------- */ 1678 1679 #undef __FUNCT__ 1680 #define __FUNCT__ "SNESMonitorSet" 1681 /*@C 1682 SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every 1683 iteration of the nonlinear solver to display the iteration's 1684 progress. 1685 1686 Collective on SNES 1687 1688 Input Parameters: 1689 + snes - the SNES context 1690 . func - monitoring routine 1691 . mctx - [optional] user-defined context for private data for the 1692 monitor routine (use PETSC_NULL if no context is desired) 1693 - monitordestroy - [optional] routine that frees monitor context 1694 (may be PETSC_NULL) 1695 1696 Calling sequence of func: 1697 $ int func(SNES snes,PetscInt its, PetscReal norm,void *mctx) 1698 1699 + snes - the SNES context 1700 . its - iteration number 1701 . norm - 2-norm function value (may be estimated) 1702 - mctx - [optional] monitoring context 1703 1704 Options Database Keys: 1705 + -snes_monitor - sets SNESMonitorDefault() 1706 . -snes_monitor_draw - sets line graph monitor, 1707 uses SNESMonitorLGCreate() 1708 _ -snes_monitor_cancel - cancels all monitors that have 1709 been hardwired into a code by 1710 calls to SNESMonitorSet(), but 1711 does not cancel those set via 1712 the options database. 1713 1714 Notes: 1715 Several different monitoring routines may be set by calling 1716 SNESMonitorSet() multiple times; all will be called in the 1717 order in which they were set. 1718 1719 Fortran notes: Only a single monitor function can be set for each SNES object 1720 1721 Level: intermediate 1722 1723 .keywords: SNES, nonlinear, set, monitor 1724 1725 .seealso: SNESMonitorDefault(), SNESMonitorCancel() 1726 @*/ 1727 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorSet(SNES snes,PetscErrorCode (*monitor)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void*)) 1728 { 1729 PetscInt i; 1730 1731 PetscFunctionBegin; 1732 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1733 if (snes->numbermonitors >= MAXSNESMONITORS) { 1734 SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1735 } 1736 for (i=0; i<snes->numbermonitors;i++) { 1737 if (monitor == snes->monitor[i] && monitordestroy == snes->monitordestroy[i] && mctx == snes->monitorcontext[i]) PetscFunctionReturn(0); 1738 1739 /* check if both default monitors that share common ASCII viewer */ 1740 if (monitor == snes->monitor[i] && monitor == SNESMonitorDefault) { 1741 if (mctx && snes->monitorcontext[i]) { 1742 PetscErrorCode ierr; 1743 PetscViewerASCIIMonitor viewer1 = (PetscViewerASCIIMonitor) mctx; 1744 PetscViewerASCIIMonitor viewer2 = (PetscViewerASCIIMonitor) snes->monitorcontext[i]; 1745 if (viewer1->viewer == viewer2->viewer) { 1746 ierr = (*monitordestroy)(mctx);CHKERRQ(ierr); 1747 PetscFunctionReturn(0); 1748 } 1749 } 1750 } 1751 } 1752 snes->monitor[snes->numbermonitors] = monitor; 1753 snes->monitordestroy[snes->numbermonitors] = monitordestroy; 1754 snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 1755 PetscFunctionReturn(0); 1756 } 1757 1758 #undef __FUNCT__ 1759 #define __FUNCT__ "SNESMonitorCancel" 1760 /*@C 1761 SNESMonitorCancel - Clears all the monitor functions for a SNES object. 1762 1763 Collective on SNES 1764 1765 Input Parameters: 1766 . snes - the SNES context 1767 1768 Options Database Key: 1769 . -snes_monitor_cancel - cancels all monitors that have been hardwired 1770 into a code by calls to SNESMonitorSet(), but does not cancel those 1771 set via the options database 1772 1773 Notes: 1774 There is no way to clear one specific monitor from a SNES object. 1775 1776 Level: intermediate 1777 1778 .keywords: SNES, nonlinear, set, monitor 1779 1780 .seealso: SNESMonitorDefault(), SNESMonitorSet() 1781 @*/ 1782 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorCancel(SNES snes) 1783 { 1784 PetscErrorCode ierr; 1785 PetscInt i; 1786 1787 PetscFunctionBegin; 1788 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1789 for (i=0; i<snes->numbermonitors; i++) { 1790 if (snes->monitordestroy[i]) { 1791 ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr); 1792 } 1793 } 1794 snes->numbermonitors = 0; 1795 PetscFunctionReturn(0); 1796 } 1797 1798 #undef __FUNCT__ 1799 #define __FUNCT__ "SNESSetConvergenceTest" 1800 /*@C 1801 SNESSetConvergenceTest - Sets the function that is to be used 1802 to test for convergence of the nonlinear iterative solution. 1803 1804 Collective on SNES 1805 1806 Input Parameters: 1807 + snes - the SNES context 1808 . func - routine to test for convergence 1809 . cctx - [optional] context for private data for the convergence routine (may be PETSC_NULL) 1810 - destroy - [optional] destructor for the context (may be PETSC_NULL; PETSC_NULL_FUNCTION in Fortran) 1811 1812 Calling sequence of func: 1813 $ PetscErrorCode func (SNES snes,PetscInt it,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx) 1814 1815 + snes - the SNES context 1816 . it - current iteration (0 is the first and is before any Newton step) 1817 . cctx - [optional] convergence context 1818 . reason - reason for convergence/divergence 1819 . xnorm - 2-norm of current iterate 1820 . gnorm - 2-norm of current step 1821 - f - 2-norm of function 1822 1823 Level: advanced 1824 1825 .keywords: SNES, nonlinear, set, convergence, test 1826 1827 .seealso: SNESDefaultConverged(), SNESSkipConverged() 1828 @*/ 1829 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*)) 1830 { 1831 PetscErrorCode ierr; 1832 1833 PetscFunctionBegin; 1834 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1835 if (!func) func = SNESSkipConverged; 1836 if (snes->ops->convergeddestroy) { 1837 ierr = (*snes->ops->convergeddestroy)(snes->cnvP);CHKERRQ(ierr); 1838 } 1839 snes->ops->converged = func; 1840 snes->ops->convergeddestroy = destroy; 1841 snes->cnvP = cctx; 1842 PetscFunctionReturn(0); 1843 } 1844 1845 #undef __FUNCT__ 1846 #define __FUNCT__ "SNESGetConvergedReason" 1847 /*@ 1848 SNESGetConvergedReason - Gets the reason the SNES iteration was stopped. 1849 1850 Not Collective 1851 1852 Input Parameter: 1853 . snes - the SNES context 1854 1855 Output Parameter: 1856 . reason - negative value indicates diverged, positive value converged, see SNESConvergedReason or the 1857 manual pages for the individual convergence tests for complete lists 1858 1859 Level: intermediate 1860 1861 Notes: Can only be called after the call the SNESSolve() is complete. 1862 1863 .keywords: SNES, nonlinear, set, convergence, test 1864 1865 .seealso: SNESSetConvergenceTest(), SNESConvergedReason 1866 @*/ 1867 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason) 1868 { 1869 PetscFunctionBegin; 1870 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1871 PetscValidPointer(reason,2); 1872 *reason = snes->reason; 1873 PetscFunctionReturn(0); 1874 } 1875 1876 #undef __FUNCT__ 1877 #define __FUNCT__ "SNESSetConvergenceHistory" 1878 /*@ 1879 SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1880 1881 Collective on SNES 1882 1883 Input Parameters: 1884 + snes - iterative context obtained from SNESCreate() 1885 . a - array to hold history 1886 . its - integer array holds the number of linear iterations for each solve. 1887 . na - size of a and its 1888 - reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero, 1889 else it continues storing new values for new nonlinear solves after the old ones 1890 1891 Notes: 1892 If set, this array will contain the function norms computed 1893 at each step. 1894 1895 This routine is useful, e.g., when running a code for purposes 1896 of accurate performance monitoring, when no I/O should be done 1897 during the section of code that is being timed. 1898 1899 Level: intermediate 1900 1901 .keywords: SNES, set, convergence, history 1902 1903 .seealso: SNESGetConvergenceHistory() 1904 1905 @*/ 1906 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscTruth reset) 1907 { 1908 PetscFunctionBegin; 1909 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1910 if (na) PetscValidScalarPointer(a,2); 1911 if (its) PetscValidIntPointer(its,3); 1912 snes->conv_hist = a; 1913 snes->conv_hist_its = its; 1914 snes->conv_hist_max = na; 1915 snes->conv_hist_len = 0; 1916 snes->conv_hist_reset = reset; 1917 PetscFunctionReturn(0); 1918 } 1919 1920 #undef __FUNCT__ 1921 #define __FUNCT__ "SNESGetConvergenceHistory" 1922 /*@C 1923 SNESGetConvergenceHistory - Gets the array used to hold the convergence history. 1924 1925 Collective on SNES 1926 1927 Input Parameter: 1928 . snes - iterative context obtained from SNESCreate() 1929 1930 Output Parameters: 1931 . a - array to hold history 1932 . its - integer array holds the number of linear iterations (or 1933 negative if not converged) for each solve. 1934 - na - size of a and its 1935 1936 Notes: 1937 The calling sequence for this routine in Fortran is 1938 $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr) 1939 1940 This routine is useful, e.g., when running a code for purposes 1941 of accurate performance monitoring, when no I/O should be done 1942 during the section of code that is being timed. 1943 1944 Level: intermediate 1945 1946 .keywords: SNES, get, convergence, history 1947 1948 .seealso: SNESSetConvergencHistory() 1949 1950 @*/ 1951 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na) 1952 { 1953 PetscFunctionBegin; 1954 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1955 if (a) *a = snes->conv_hist; 1956 if (its) *its = snes->conv_hist_its; 1957 if (na) *na = snes->conv_hist_len; 1958 PetscFunctionReturn(0); 1959 } 1960 1961 #undef __FUNCT__ 1962 #define __FUNCT__ "SNESSetUpdate" 1963 /*@C 1964 SNESSetUpdate - Sets the general-purpose update function called 1965 at the beginning o every iteration of the nonlinear solve. Specifically 1966 it is called just before the Jacobian is "evaluated". 1967 1968 Collective on SNES 1969 1970 Input Parameters: 1971 . snes - The nonlinear solver context 1972 . func - The function 1973 1974 Calling sequence of func: 1975 . func (SNES snes, PetscInt step); 1976 1977 . step - The current step of the iteration 1978 1979 Level: intermediate 1980 1981 .keywords: SNES, update 1982 1983 .seealso SNESDefaultUpdate(), SNESSetJacobian(), SNESSolve() 1984 @*/ 1985 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt)) 1986 { 1987 PetscFunctionBegin; 1988 PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 1989 snes->ops->update = func; 1990 PetscFunctionReturn(0); 1991 } 1992 1993 #undef __FUNCT__ 1994 #define __FUNCT__ "SNESDefaultUpdate" 1995 /*@ 1996 SNESDefaultUpdate - The default update function which does nothing. 1997 1998 Not collective 1999 2000 Input Parameters: 2001 . snes - The nonlinear solver context 2002 . step - The current step of the iteration 2003 2004 Level: intermediate 2005 2006 .keywords: SNES, update 2007 .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultShortolutionBC() 2008 @*/ 2009 PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultUpdate(SNES snes, PetscInt step) 2010 { 2011 PetscFunctionBegin; 2012 PetscFunctionReturn(0); 2013 } 2014 2015 #undef __FUNCT__ 2016 #define __FUNCT__ "SNESScaleStep_Private" 2017 /* 2018 SNESScaleStep_Private - Scales a step so that its length is less than the 2019 positive parameter delta. 2020 2021 Input Parameters: 2022 + snes - the SNES context 2023 . y - approximate solution of linear system 2024 . fnorm - 2-norm of current function 2025 - delta - trust region size 2026 2027 Output Parameters: 2028 + gpnorm - predicted function norm at the new point, assuming local 2029 linearization. The value is zero if the step lies within the trust 2030 region, and exceeds zero otherwise. 2031 - ynorm - 2-norm of the step 2032 2033 Note: 2034 For non-trust region methods such as SNESLS, the parameter delta 2035 is set to be the maximum allowable step size. 2036 2037 .keywords: SNES, nonlinear, scale, step 2038 */ 2039 PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm) 2040 { 2041 PetscReal nrm; 2042 PetscScalar cnorm; 2043 PetscErrorCode ierr; 2044 2045 PetscFunctionBegin; 2046 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2047 PetscValidHeaderSpecific(y,VEC_COOKIE,2); 2048 PetscCheckSameComm(snes,1,y,2); 2049 2050 ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr); 2051 if (nrm > *delta) { 2052 nrm = *delta/nrm; 2053 *gpnorm = (1.0 - nrm)*(*fnorm); 2054 cnorm = nrm; 2055 ierr = VecScale(y,cnorm);CHKERRQ(ierr); 2056 *ynorm = *delta; 2057 } else { 2058 *gpnorm = 0.0; 2059 *ynorm = nrm; 2060 } 2061 PetscFunctionReturn(0); 2062 } 2063 2064 #undef __FUNCT__ 2065 #define __FUNCT__ "SNESSolve" 2066 /*@C 2067 SNESSolve - Solves a nonlinear system F(x) = b. 2068 Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX(). 2069 2070 Collective on SNES 2071 2072 Input Parameters: 2073 + snes - the SNES context 2074 . b - the constant part of the equation, or PETSC_NULL to use zero. 2075 - x - the solution vector. 2076 2077 Notes: 2078 The user should initialize the vector,x, with the initial guess 2079 for the nonlinear solve prior to calling SNESSolve. In particular, 2080 to employ an initial guess of zero, the user should explicitly set 2081 this vector to zero by calling VecSet(). 2082 2083 Level: beginner 2084 2085 .keywords: SNES, nonlinear, solve 2086 2087 .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian() 2088 @*/ 2089 PetscErrorCode PETSCSNES_DLLEXPORT SNESSolve(SNES snes,Vec b,Vec x) 2090 { 2091 PetscErrorCode ierr; 2092 PetscTruth flg; 2093 char filename[PETSC_MAX_PATH_LEN]; 2094 PetscViewer viewer; 2095 2096 PetscFunctionBegin; 2097 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2098 PetscValidHeaderSpecific(x,VEC_COOKIE,3); 2099 PetscCheckSameComm(snes,1,x,3); 2100 if (b) PetscValidHeaderSpecific(b,VEC_COOKIE,2); 2101 if (b) PetscCheckSameComm(snes,1,b,2); 2102 2103 /* set solution vector */ 2104 ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 2105 if (snes->vec_sol) { ierr = VecDestroy(snes->vec_sol);CHKERRQ(ierr); } 2106 snes->vec_sol = x; 2107 /* set afine vector if provided */ 2108 if (b) { ierr = PetscObjectReference((PetscObject)b);CHKERRQ(ierr); } 2109 if (snes->vec_rhs) { ierr = VecDestroy(snes->vec_rhs);CHKERRQ(ierr); } 2110 snes->vec_rhs = b; 2111 2112 if (!snes->vec_func && snes->vec_rhs) { 2113 ierr = VecDuplicate(b, &snes->vec_func); CHKERRQ(ierr); 2114 } 2115 2116 ierr = SNESSetUp(snes);CHKERRQ(ierr); 2117 2118 if (snes->conv_hist_reset) snes->conv_hist_len = 0; 2119 snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0; 2120 2121 ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 2122 ierr = (*snes->ops->solve)(snes);CHKERRQ(ierr); 2123 ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 2124 if (snes->domainerror){ 2125 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 2126 snes->domainerror = PETSC_FALSE; 2127 } 2128 2129 if (!snes->reason) { 2130 SETERRQ(PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason"); 2131 } 2132 2133 ierr = PetscOptionsGetString(((PetscObject)snes)->prefix,"-snes_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 2134 if (flg && !PetscPreLoadingOn) { 2135 ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,filename,&viewer);CHKERRQ(ierr); 2136 ierr = SNESView(snes,viewer);CHKERRQ(ierr); 2137 ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr); 2138 } 2139 2140 ierr = PetscOptionsHasName(((PetscObject)snes)->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr); 2141 if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); } 2142 if (snes->printreason) { 2143 if (snes->reason > 0) { 2144 ierr = PetscPrintf(((PetscObject)snes)->comm,"Nonlinear solve converged due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr); 2145 } else { 2146 ierr = PetscPrintf(((PetscObject)snes)->comm,"Nonlinear solve did not converge due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr); 2147 } 2148 } 2149 2150 PetscFunctionReturn(0); 2151 } 2152 2153 /* --------- Internal routines for SNES Package --------- */ 2154 2155 #undef __FUNCT__ 2156 #define __FUNCT__ "SNESSetType" 2157 /*@C 2158 SNESSetType - Sets the method for the nonlinear solver. 2159 2160 Collective on SNES 2161 2162 Input Parameters: 2163 + snes - the SNES context 2164 - type - a known method 2165 2166 Options Database Key: 2167 . -snes_type <type> - Sets the method; use -help for a list 2168 of available methods (for instance, ls or tr) 2169 2170 Notes: 2171 See "petsc/include/petscsnes.h" for available methods (for instance) 2172 + SNESLS - Newton's method with line search 2173 (systems of nonlinear equations) 2174 . SNESTR - Newton's method with trust region 2175 (systems of nonlinear equations) 2176 2177 Normally, it is best to use the SNESSetFromOptions() command and then 2178 set the SNES solver type from the options database rather than by using 2179 this routine. Using the options database provides the user with 2180 maximum flexibility in evaluating the many nonlinear solvers. 2181 The SNESSetType() routine is provided for those situations where it 2182 is necessary to set the nonlinear solver independently of the command 2183 line or options database. This might be the case, for example, when 2184 the choice of solver changes during the execution of the program, 2185 and the user's application is taking responsibility for choosing the 2186 appropriate method. 2187 2188 Level: intermediate 2189 2190 .keywords: SNES, set, type 2191 2192 .seealso: SNESType, SNESCreate() 2193 2194 @*/ 2195 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetType(SNES snes,const SNESType type) 2196 { 2197 PetscErrorCode ierr,(*r)(SNES); 2198 PetscTruth match; 2199 2200 PetscFunctionBegin; 2201 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2202 PetscValidCharPointer(type,2); 2203 2204 ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr); 2205 if (match) PetscFunctionReturn(0); 2206 2207 ierr = PetscFListFind(SNESList,((PetscObject)snes)->comm,type,(void (**)(void)) &r);CHKERRQ(ierr); 2208 if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type); 2209 /* Destroy the previous private SNES context */ 2210 if (snes->ops->destroy) { ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr); } 2211 /* Reinitialize function pointers in SNESOps structure */ 2212 snes->ops->setup = 0; 2213 snes->ops->solve = 0; 2214 snes->ops->view = 0; 2215 snes->ops->setfromoptions = 0; 2216 snes->ops->destroy = 0; 2217 /* Call the SNESCreate_XXX routine for this particular Nonlinear solver */ 2218 snes->setupcalled = PETSC_FALSE; 2219 ierr = (*r)(snes);CHKERRQ(ierr); 2220 ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr); 2221 PetscFunctionReturn(0); 2222 } 2223 2224 2225 /* --------------------------------------------------------------------- */ 2226 #undef __FUNCT__ 2227 #define __FUNCT__ "SNESRegisterDestroy" 2228 /*@ 2229 SNESRegisterDestroy - Frees the list of nonlinear solvers that were 2230 registered by SNESRegisterDynamic(). 2231 2232 Not Collective 2233 2234 Level: advanced 2235 2236 .keywords: SNES, nonlinear, register, destroy 2237 2238 .seealso: SNESRegisterAll(), SNESRegisterAll() 2239 @*/ 2240 PetscErrorCode PETSCSNES_DLLEXPORT SNESRegisterDestroy(void) 2241 { 2242 PetscErrorCode ierr; 2243 2244 PetscFunctionBegin; 2245 ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr); 2246 SNESRegisterAllCalled = PETSC_FALSE; 2247 PetscFunctionReturn(0); 2248 } 2249 2250 #undef __FUNCT__ 2251 #define __FUNCT__ "SNESGetType" 2252 /*@C 2253 SNESGetType - Gets the SNES method type and name (as a string). 2254 2255 Not Collective 2256 2257 Input Parameter: 2258 . snes - nonlinear solver context 2259 2260 Output Parameter: 2261 . type - SNES method (a character string) 2262 2263 Level: intermediate 2264 2265 .keywords: SNES, nonlinear, get, type, name 2266 @*/ 2267 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetType(SNES snes,const SNESType *type) 2268 { 2269 PetscFunctionBegin; 2270 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2271 PetscValidPointer(type,2); 2272 *type = ((PetscObject)snes)->type_name; 2273 PetscFunctionReturn(0); 2274 } 2275 2276 #undef __FUNCT__ 2277 #define __FUNCT__ "SNESGetSolution" 2278 /*@ 2279 SNESGetSolution - Returns the vector where the approximate solution is 2280 stored. 2281 2282 Not Collective, but Vec is parallel if SNES is parallel 2283 2284 Input Parameter: 2285 . snes - the SNES context 2286 2287 Output Parameter: 2288 . x - the solution 2289 2290 Level: intermediate 2291 2292 .keywords: SNES, nonlinear, get, solution 2293 2294 .seealso: SNESGetSolutionUpdate(), SNESGetFunction() 2295 @*/ 2296 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolution(SNES snes,Vec *x) 2297 { 2298 PetscFunctionBegin; 2299 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2300 PetscValidPointer(x,2); 2301 *x = snes->vec_sol; 2302 PetscFunctionReturn(0); 2303 } 2304 2305 #undef __FUNCT__ 2306 #define __FUNCT__ "SNESGetSolutionUpdate" 2307 /*@ 2308 SNESGetSolutionUpdate - Returns the vector where the solution update is 2309 stored. 2310 2311 Not Collective, but Vec is parallel if SNES is parallel 2312 2313 Input Parameter: 2314 . snes - the SNES context 2315 2316 Output Parameter: 2317 . x - the solution update 2318 2319 Level: advanced 2320 2321 .keywords: SNES, nonlinear, get, solution, update 2322 2323 .seealso: SNESGetSolution(), SNESGetFunction() 2324 @*/ 2325 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolutionUpdate(SNES snes,Vec *x) 2326 { 2327 PetscFunctionBegin; 2328 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2329 PetscValidPointer(x,2); 2330 *x = snes->vec_sol_update; 2331 PetscFunctionReturn(0); 2332 } 2333 2334 #undef __FUNCT__ 2335 #define __FUNCT__ "SNESGetFunction" 2336 /*@C 2337 SNESGetFunction - Returns the vector where the function is stored. 2338 2339 Not Collective, but Vec is parallel if SNES is parallel 2340 2341 Input Parameter: 2342 . snes - the SNES context 2343 2344 Output Parameter: 2345 + r - the function (or PETSC_NULL) 2346 . func - the function (or PETSC_NULL) 2347 - ctx - the function context (or PETSC_NULL) 2348 2349 Level: advanced 2350 2351 .keywords: SNES, nonlinear, get, function 2352 2353 .seealso: SNESSetFunction(), SNESGetSolution() 2354 @*/ 2355 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx) 2356 { 2357 PetscFunctionBegin; 2358 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2359 if (r) *r = snes->vec_func; 2360 if (func) *func = snes->ops->computefunction; 2361 if (ctx) *ctx = snes->funP; 2362 PetscFunctionReturn(0); 2363 } 2364 2365 #undef __FUNCT__ 2366 #define __FUNCT__ "SNESSetOptionsPrefix" 2367 /*@C 2368 SNESSetOptionsPrefix - Sets the prefix used for searching for all 2369 SNES options in the database. 2370 2371 Collective on SNES 2372 2373 Input Parameter: 2374 + snes - the SNES context 2375 - prefix - the prefix to prepend to all option names 2376 2377 Notes: 2378 A hyphen (-) must NOT be given at the beginning of the prefix name. 2379 The first character of all runtime options is AUTOMATICALLY the hyphen. 2380 2381 Level: advanced 2382 2383 .keywords: SNES, set, options, prefix, database 2384 2385 .seealso: SNESSetFromOptions() 2386 @*/ 2387 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetOptionsPrefix(SNES snes,const char prefix[]) 2388 { 2389 PetscErrorCode ierr; 2390 2391 PetscFunctionBegin; 2392 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2393 ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2394 if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);} 2395 ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 2396 PetscFunctionReturn(0); 2397 } 2398 2399 #undef __FUNCT__ 2400 #define __FUNCT__ "SNESAppendOptionsPrefix" 2401 /*@C 2402 SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 2403 SNES options in the database. 2404 2405 Collective on SNES 2406 2407 Input Parameters: 2408 + snes - the SNES context 2409 - prefix - the prefix to prepend to all option names 2410 2411 Notes: 2412 A hyphen (-) must NOT be given at the beginning of the prefix name. 2413 The first character of all runtime options is AUTOMATICALLY the hyphen. 2414 2415 Level: advanced 2416 2417 .keywords: SNES, append, options, prefix, database 2418 2419 .seealso: SNESGetOptionsPrefix() 2420 @*/ 2421 PetscErrorCode PETSCSNES_DLLEXPORT SNESAppendOptionsPrefix(SNES snes,const char prefix[]) 2422 { 2423 PetscErrorCode ierr; 2424 2425 PetscFunctionBegin; 2426 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2427 ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2428 if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);} 2429 ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 2430 PetscFunctionReturn(0); 2431 } 2432 2433 #undef __FUNCT__ 2434 #define __FUNCT__ "SNESGetOptionsPrefix" 2435 /*@C 2436 SNESGetOptionsPrefix - Sets the prefix used for searching for all 2437 SNES options in the database. 2438 2439 Not Collective 2440 2441 Input Parameter: 2442 . snes - the SNES context 2443 2444 Output Parameter: 2445 . prefix - pointer to the prefix string used 2446 2447 Notes: On the fortran side, the user should pass in a string 'prefix' of 2448 sufficient length to hold the prefix. 2449 2450 Level: advanced 2451 2452 .keywords: SNES, get, options, prefix, database 2453 2454 .seealso: SNESAppendOptionsPrefix() 2455 @*/ 2456 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetOptionsPrefix(SNES snes,const char *prefix[]) 2457 { 2458 PetscErrorCode ierr; 2459 2460 PetscFunctionBegin; 2461 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2462 ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2463 PetscFunctionReturn(0); 2464 } 2465 2466 2467 #undef __FUNCT__ 2468 #define __FUNCT__ "SNESRegister" 2469 /*@C 2470 SNESRegister - See SNESRegisterDynamic() 2471 2472 Level: advanced 2473 @*/ 2474 PetscErrorCode PETSCSNES_DLLEXPORT SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES)) 2475 { 2476 char fullname[PETSC_MAX_PATH_LEN]; 2477 PetscErrorCode ierr; 2478 2479 PetscFunctionBegin; 2480 ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 2481 ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 2482 PetscFunctionReturn(0); 2483 } 2484 2485 #undef __FUNCT__ 2486 #define __FUNCT__ "SNESTestLocalMin" 2487 PetscErrorCode PETSCSNES_DLLEXPORT SNESTestLocalMin(SNES snes) 2488 { 2489 PetscErrorCode ierr; 2490 PetscInt N,i,j; 2491 Vec u,uh,fh; 2492 PetscScalar value; 2493 PetscReal norm; 2494 2495 PetscFunctionBegin; 2496 ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr); 2497 ierr = VecDuplicate(u,&uh);CHKERRQ(ierr); 2498 ierr = VecDuplicate(u,&fh);CHKERRQ(ierr); 2499 2500 /* currently only works for sequential */ 2501 ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n"); 2502 ierr = VecGetSize(u,&N);CHKERRQ(ierr); 2503 for (i=0; i<N; i++) { 2504 ierr = VecCopy(u,uh);CHKERRQ(ierr); 2505 ierr = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr); 2506 for (j=-10; j<11; j++) { 2507 value = PetscSign(j)*exp(PetscAbs(j)-10.0); 2508 ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2509 ierr = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr); 2510 ierr = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr); 2511 ierr = PetscPrintf(PETSC_COMM_WORLD," j norm %D %18.16e\n",j,norm);CHKERRQ(ierr); 2512 value = -value; 2513 ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2514 } 2515 } 2516 ierr = VecDestroy(uh);CHKERRQ(ierr); 2517 ierr = VecDestroy(fh);CHKERRQ(ierr); 2518 PetscFunctionReturn(0); 2519 } 2520 2521 #undef __FUNCT__ 2522 #define __FUNCT__ "SNESKSPSetUseEW" 2523 /*@ 2524 SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for 2525 computing relative tolerance for linear solvers within an inexact 2526 Newton method. 2527 2528 Collective on SNES 2529 2530 Input Parameters: 2531 + snes - SNES context 2532 - flag - PETSC_TRUE or PETSC_FALSE 2533 2534 Notes: 2535 Currently, the default is to use a constant relative tolerance for 2536 the inner linear solvers. Alternatively, one can use the 2537 Eisenstat-Walker method, where the relative convergence tolerance 2538 is reset at each Newton iteration according progress of the nonlinear 2539 solver. 2540 2541 Level: advanced 2542 2543 Reference: 2544 S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an 2545 inexact Newton method", SISC 17 (1), pp.16-32, 1996. 2546 2547 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton 2548 2549 .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW() 2550 @*/ 2551 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPSetUseEW(SNES snes,PetscTruth flag) 2552 { 2553 PetscFunctionBegin; 2554 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2555 snes->ksp_ewconv = flag; 2556 PetscFunctionReturn(0); 2557 } 2558 2559 #undef __FUNCT__ 2560 #define __FUNCT__ "SNESKSPGetUseEW" 2561 /*@ 2562 SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method 2563 for computing relative tolerance for linear solvers within an 2564 inexact Newton method. 2565 2566 Not Collective 2567 2568 Input Parameter: 2569 . snes - SNES context 2570 2571 Output Parameter: 2572 . flag - PETSC_TRUE or PETSC_FALSE 2573 2574 Notes: 2575 Currently, the default is to use a constant relative tolerance for 2576 the inner linear solvers. Alternatively, one can use the 2577 Eisenstat-Walker method, where the relative convergence tolerance 2578 is reset at each Newton iteration according progress of the nonlinear 2579 solver. 2580 2581 Level: advanced 2582 2583 Reference: 2584 S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an 2585 inexact Newton method", SISC 17 (1), pp.16-32, 1996. 2586 2587 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton 2588 2589 .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW() 2590 @*/ 2591 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPGetUseEW(SNES snes, PetscTruth *flag) 2592 { 2593 PetscFunctionBegin; 2594 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2595 PetscValidPointer(flag,2); 2596 *flag = snes->ksp_ewconv; 2597 PetscFunctionReturn(0); 2598 } 2599 2600 #undef __FUNCT__ 2601 #define __FUNCT__ "SNESKSPSetParametersEW" 2602 /*@ 2603 SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker 2604 convergence criteria for the linear solvers within an inexact 2605 Newton method. 2606 2607 Collective on SNES 2608 2609 Input Parameters: 2610 + snes - SNES context 2611 . version - version 1, 2 (default is 2) or 3 2612 . rtol_0 - initial relative tolerance (0 <= rtol_0 < 1) 2613 . rtol_max - maximum relative tolerance (0 <= rtol_max < 1) 2614 . gamma - multiplicative factor for version 2 rtol computation 2615 (0 <= gamma2 <= 1) 2616 . alpha - power for version 2 rtol computation (1 < alpha <= 2) 2617 . alpha2 - power for safeguard 2618 - threshold - threshold for imposing safeguard (0 < threshold < 1) 2619 2620 Note: 2621 Version 3 was contributed by Luis Chacon, June 2006. 2622 2623 Use PETSC_DEFAULT to retain the default for any of the parameters. 2624 2625 Level: advanced 2626 2627 Reference: 2628 S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an 2629 inexact Newton method", Utah State University Math. Stat. Dept. Res. 2630 Report 6/94/75, June, 1994, to appear in SIAM J. Sci. Comput. 2631 2632 .keywords: SNES, KSP, Eisenstat, Walker, set, parameters 2633 2634 .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW() 2635 @*/ 2636 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max, 2637 PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold) 2638 { 2639 SNESKSPEW *kctx; 2640 PetscFunctionBegin; 2641 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2642 kctx = (SNESKSPEW*)snes->kspconvctx; 2643 if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing"); 2644 2645 if (version != PETSC_DEFAULT) kctx->version = version; 2646 if (rtol_0 != PETSC_DEFAULT) kctx->rtol_0 = rtol_0; 2647 if (rtol_max != PETSC_DEFAULT) kctx->rtol_max = rtol_max; 2648 if (gamma != PETSC_DEFAULT) kctx->gamma = gamma; 2649 if (alpha != PETSC_DEFAULT) kctx->alpha = alpha; 2650 if (alpha2 != PETSC_DEFAULT) kctx->alpha2 = alpha2; 2651 if (threshold != PETSC_DEFAULT) kctx->threshold = threshold; 2652 2653 if (kctx->version < 1 || kctx->version > 3) { 2654 SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 and 3 are supported: %D",kctx->version); 2655 } 2656 if (kctx->rtol_0 < 0.0 || kctx->rtol_0 >= 1.0) { 2657 SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_0 < 1.0: %G",kctx->rtol_0); 2658 } 2659 if (kctx->rtol_max < 0.0 || kctx->rtol_max >= 1.0) { 2660 SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_max (%G) < 1.0\n",kctx->rtol_max); 2661 } 2662 if (kctx->gamma < 0.0 || kctx->gamma > 1.0) { 2663 SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= gamma (%G) <= 1.0\n",kctx->gamma); 2664 } 2665 if (kctx->alpha <= 1.0 || kctx->alpha > 2.0) { 2666 SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"1.0 < alpha (%G) <= 2.0\n",kctx->alpha); 2667 } 2668 if (kctx->threshold <= 0.0 || kctx->threshold >= 1.0) { 2669 SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 < threshold (%G) < 1.0\n",kctx->threshold); 2670 } 2671 PetscFunctionReturn(0); 2672 } 2673 2674 #undef __FUNCT__ 2675 #define __FUNCT__ "SNESKSPGetParametersEW" 2676 /*@ 2677 SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker 2678 convergence criteria for the linear solvers within an inexact 2679 Newton method. 2680 2681 Not Collective 2682 2683 Input Parameters: 2684 snes - SNES context 2685 2686 Output Parameters: 2687 + version - version 1, 2 (default is 2) or 3 2688 . rtol_0 - initial relative tolerance (0 <= rtol_0 < 1) 2689 . rtol_max - maximum relative tolerance (0 <= rtol_max < 1) 2690 . gamma - multiplicative factor for version 2 rtol computation 2691 (0 <= gamma2 <= 1) 2692 . alpha - power for version 2 rtol computation (1 < alpha <= 2) 2693 . alpha2 - power for safeguard 2694 - threshold - threshold for imposing safeguard (0 < threshold < 1) 2695 2696 Level: advanced 2697 2698 .keywords: SNES, KSP, Eisenstat, Walker, get, parameters 2699 2700 .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW() 2701 @*/ 2702 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max, 2703 PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold) 2704 { 2705 SNESKSPEW *kctx; 2706 PetscFunctionBegin; 2707 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2708 kctx = (SNESKSPEW*)snes->kspconvctx; 2709 if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing"); 2710 if(version) *version = kctx->version; 2711 if(rtol_0) *rtol_0 = kctx->rtol_0; 2712 if(rtol_max) *rtol_max = kctx->rtol_max; 2713 if(gamma) *gamma = kctx->gamma; 2714 if(alpha) *alpha = kctx->alpha; 2715 if(alpha2) *alpha2 = kctx->alpha2; 2716 if(threshold) *threshold = kctx->threshold; 2717 PetscFunctionReturn(0); 2718 } 2719 2720 #undef __FUNCT__ 2721 #define __FUNCT__ "SNESKSPEW_PreSolve" 2722 static PetscErrorCode SNESKSPEW_PreSolve(SNES snes, KSP ksp, Vec b, Vec x) 2723 { 2724 PetscErrorCode ierr; 2725 SNESKSPEW *kctx = (SNESKSPEW*)snes->kspconvctx; 2726 PetscReal rtol=PETSC_DEFAULT,stol; 2727 2728 PetscFunctionBegin; 2729 if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists"); 2730 if (!snes->iter) { /* first time in, so use the original user rtol */ 2731 rtol = kctx->rtol_0; 2732 } else { 2733 if (kctx->version == 1) { 2734 rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last; 2735 if (rtol < 0.0) rtol = -rtol; 2736 stol = pow(kctx->rtol_last,kctx->alpha2); 2737 if (stol > kctx->threshold) rtol = PetscMax(rtol,stol); 2738 } else if (kctx->version == 2) { 2739 rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha); 2740 stol = kctx->gamma * pow(kctx->rtol_last,kctx->alpha); 2741 if (stol > kctx->threshold) rtol = PetscMax(rtol,stol); 2742 } else if (kctx->version == 3) {/* contributed by Luis Chacon, June 2006. */ 2743 rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha); 2744 /* safeguard: avoid sharp decrease of rtol */ 2745 stol = kctx->gamma*pow(kctx->rtol_last,kctx->alpha); 2746 stol = PetscMax(rtol,stol); 2747 rtol = PetscMin(kctx->rtol_0,stol); 2748 /* safeguard: avoid oversolving */ 2749 stol = kctx->gamma*(snes->ttol)/snes->norm; 2750 stol = PetscMax(rtol,stol); 2751 rtol = PetscMin(kctx->rtol_0,stol); 2752 } else SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 or 3 are supported: %D",kctx->version); 2753 } 2754 /* safeguard: avoid rtol greater than one */ 2755 rtol = PetscMin(rtol,kctx->rtol_max); 2756 ierr = KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 2757 ierr = PetscInfo3(snes,"iter %D, Eisenstat-Walker (version %D) KSP rtol=%G\n",snes->iter,kctx->version,rtol);CHKERRQ(ierr); 2758 PetscFunctionReturn(0); 2759 } 2760 2761 #undef __FUNCT__ 2762 #define __FUNCT__ "SNESKSPEW_PostSolve" 2763 static PetscErrorCode SNESKSPEW_PostSolve(SNES snes, KSP ksp, Vec b, Vec x) 2764 { 2765 PetscErrorCode ierr; 2766 SNESKSPEW *kctx = (SNESKSPEW*)snes->kspconvctx; 2767 PCSide pcside; 2768 Vec lres; 2769 2770 PetscFunctionBegin; 2771 if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists"); 2772 ierr = KSPGetTolerances(ksp,&kctx->rtol_last,0,0,0);CHKERRQ(ierr); 2773 ierr = SNESGetFunctionNorm(snes,&kctx->norm_last);CHKERRQ(ierr); 2774 if (kctx->version == 1) { 2775 ierr = KSPGetPreconditionerSide(ksp,&pcside);CHKERRQ(ierr); 2776 if (pcside == PC_RIGHT) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */ 2777 /* KSP residual is true linear residual */ 2778 ierr = KSPGetResidualNorm(ksp,&kctx->lresid_last);CHKERRQ(ierr); 2779 } else { 2780 /* KSP residual is preconditioned residual */ 2781 /* compute true linear residual norm */ 2782 ierr = VecDuplicate(b,&lres);CHKERRQ(ierr); 2783 ierr = MatMult(snes->jacobian,x,lres);CHKERRQ(ierr); 2784 ierr = VecAYPX(lres,-1.0,b);CHKERRQ(ierr); 2785 ierr = VecNorm(lres,NORM_2,&kctx->lresid_last);CHKERRQ(ierr); 2786 ierr = VecDestroy(lres);CHKERRQ(ierr); 2787 } 2788 } 2789 PetscFunctionReturn(0); 2790 } 2791 2792 #undef __FUNCT__ 2793 #define __FUNCT__ "SNES_KSPSolve" 2794 PetscErrorCode SNES_KSPSolve(SNES snes, KSP ksp, Vec b, Vec x) 2795 { 2796 PetscErrorCode ierr; 2797 2798 PetscFunctionBegin; 2799 if (snes->ksp_ewconv) { ierr = SNESKSPEW_PreSolve(snes,ksp,b,x);CHKERRQ(ierr); } 2800 ierr = KSPSolve(ksp,b,x);CHKERRQ(ierr); 2801 if (snes->ksp_ewconv) { ierr = SNESKSPEW_PostSolve(snes,ksp,b,x);CHKERRQ(ierr); } 2802 PetscFunctionReturn(0); 2803 } 2804