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