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