1 2 #include <private/snesimpl.h> /*I "petscsnes.h" I*/ 3 4 PetscBool SNESRegisterAllCalled = PETSC_FALSE; 5 PetscFList SNESList = PETSC_NULL; 6 7 /* Logging support */ 8 PetscClassId SNES_CLASSID; 9 PetscLogEvent SNES_Solve, SNES_LineSearch, SNES_FunctionEval, SNES_JacobianEval; 10 11 #undef __FUNCT__ 12 #define __FUNCT__ "SNESSetErrorIfNotConverged" 13 /*@ 14 SNESSetErrorIfNotConverged - Causes SNESSolve() to generate an error if the solver has not converged. 15 16 Logically Collective on SNES 17 18 Input Parameters: 19 + snes - iterative context obtained from SNESCreate() 20 - flg - PETSC_TRUE indicates you want the error generated 21 22 Options database keys: 23 . -snes_error_if_not_converged : this takes an optional truth value (0/1/no/yes/true/false) 24 25 Level: intermediate 26 27 Notes: 28 Normally PETSc continues if a linear solver fails to converge, you can call SNESGetConvergedReason() after a SNESSolve() 29 to determine if it has converged. 30 31 .keywords: SNES, set, initial guess, nonzero 32 33 .seealso: SNESGetErrorIfNotConverged(), KSPGetErrorIfNotConverged(), KSPSetErrorIFNotConverged() 34 @*/ 35 PetscErrorCode SNESSetErrorIfNotConverged(SNES snes,PetscBool flg) 36 { 37 PetscFunctionBegin; 38 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 39 PetscValidLogicalCollectiveBool(snes,flg,2); 40 snes->errorifnotconverged = flg; 41 PetscFunctionReturn(0); 42 } 43 44 #undef __FUNCT__ 45 #define __FUNCT__ "SNESGetErrorIfNotConverged" 46 /*@ 47 SNESGetErrorIfNotConverged - Will SNESSolve() generate an error if the solver does not converge? 48 49 Not Collective 50 51 Input Parameter: 52 . snes - iterative context obtained from SNESCreate() 53 54 Output Parameter: 55 . flag - PETSC_TRUE if it will generate an error, else PETSC_FALSE 56 57 Level: intermediate 58 59 .keywords: SNES, set, initial guess, nonzero 60 61 .seealso: SNESSetErrorIfNotConverged(), KSPGetErrorIfNotConverged(), KSPSetErrorIFNotConverged() 62 @*/ 63 PetscErrorCode SNESGetErrorIfNotConverged(SNES snes,PetscBool *flag) 64 { 65 PetscFunctionBegin; 66 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 67 PetscValidPointer(flag,2); 68 *flag = snes->errorifnotconverged; 69 PetscFunctionReturn(0); 70 } 71 72 #undef __FUNCT__ 73 #define __FUNCT__ "SNESSetFunctionDomainError" 74 /*@ 75 SNESSetFunctionDomainError - tells SNES that the input vector to your FormFunction is not 76 in the functions domain. For example, negative pressure. 77 78 Logically Collective on SNES 79 80 Input Parameters: 81 . SNES - the SNES context 82 83 Level: advanced 84 85 .keywords: SNES, view 86 87 .seealso: SNESCreate(), SNESSetFunction() 88 @*/ 89 PetscErrorCode SNESSetFunctionDomainError(SNES snes) 90 { 91 PetscFunctionBegin; 92 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 93 snes->domainerror = PETSC_TRUE; 94 PetscFunctionReturn(0); 95 } 96 97 #undef __FUNCT__ 98 #define __FUNCT__ "SNESView" 99 /*@C 100 SNESView - Prints the SNES data structure. 101 102 Collective on SNES 103 104 Input Parameters: 105 + SNES - the SNES context 106 - viewer - visualization context 107 108 Options Database Key: 109 . -snes_view - Calls SNESView() at end of SNESSolve() 110 111 Notes: 112 The available visualization contexts include 113 + PETSC_VIEWER_STDOUT_SELF - standard output (default) 114 - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 115 output where only the first processor opens 116 the file. All other processors send their 117 data to the first processor to print. 118 119 The user can open an alternative visualization context with 120 PetscViewerASCIIOpen() - output to a specified file. 121 122 Level: beginner 123 124 .keywords: SNES, view 125 126 .seealso: PetscViewerASCIIOpen() 127 @*/ 128 PetscErrorCode SNESView(SNES snes,PetscViewer viewer) 129 { 130 SNESKSPEW *kctx; 131 PetscErrorCode ierr; 132 KSP ksp; 133 PetscBool iascii,isstring; 134 135 PetscFunctionBegin; 136 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 137 if (!viewer) { 138 ierr = PetscViewerASCIIGetStdout(((PetscObject)snes)->comm,&viewer);CHKERRQ(ierr); 139 } 140 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 141 PetscCheckSameComm(snes,1,viewer,2); 142 143 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 144 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 145 if (iascii) { 146 ierr = PetscObjectPrintClassNamePrefixType((PetscObject)snes,viewer,"SNES Object");CHKERRQ(ierr); 147 if (snes->ops->view) { 148 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 149 ierr = (*snes->ops->view)(snes,viewer);CHKERRQ(ierr); 150 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 151 } 152 ierr = PetscViewerASCIIPrintf(viewer," maximum iterations=%D, maximum function evaluations=%D\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr); 153 ierr = PetscViewerASCIIPrintf(viewer," tolerances: relative=%G, absolute=%G, solution=%G\n", 154 snes->rtol,snes->abstol,snes->xtol);CHKERRQ(ierr); 155 ierr = PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",snes->linear_its);CHKERRQ(ierr); 156 ierr = PetscViewerASCIIPrintf(viewer," total number of function evaluations=%D\n",snes->nfuncs);CHKERRQ(ierr); 157 if (snes->ksp_ewconv) { 158 kctx = (SNESKSPEW *)snes->kspconvctx; 159 if (kctx) { 160 ierr = PetscViewerASCIIPrintf(viewer," Eisenstat-Walker computation of KSP relative tolerance (version %D)\n",kctx->version);CHKERRQ(ierr); 161 ierr = PetscViewerASCIIPrintf(viewer," rtol_0=%G, rtol_max=%G, threshold=%G\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr); 162 ierr = PetscViewerASCIIPrintf(viewer," gamma=%G, alpha=%G, alpha2=%G\n",kctx->gamma,kctx->alpha,kctx->alpha2);CHKERRQ(ierr); 163 } 164 } 165 if (snes->lagpreconditioner == -1) { 166 ierr = PetscViewerASCIIPrintf(viewer," Preconditioned is never rebuilt\n");CHKERRQ(ierr); 167 } else if (snes->lagpreconditioner > 1) { 168 ierr = PetscViewerASCIIPrintf(viewer," Preconditioned is rebuilt every %D new Jacobians\n",snes->lagpreconditioner);CHKERRQ(ierr); 169 } 170 if (snes->lagjacobian == -1) { 171 ierr = PetscViewerASCIIPrintf(viewer," Jacobian is never rebuilt\n");CHKERRQ(ierr); 172 } else if (snes->lagjacobian > 1) { 173 ierr = PetscViewerASCIIPrintf(viewer," Jacobian is rebuilt every %D Newton iterations\n",snes->lagjacobian);CHKERRQ(ierr); 174 } 175 } else if (isstring) { 176 const char *type; 177 ierr = SNESGetType(snes,&type);CHKERRQ(ierr); 178 ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr); 179 } 180 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 181 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 182 ierr = KSPView(ksp,viewer);CHKERRQ(ierr); 183 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 184 PetscFunctionReturn(0); 185 } 186 187 /* 188 We retain a list of functions that also take SNES command 189 line options. These are called at the end SNESSetFromOptions() 190 */ 191 #define MAXSETFROMOPTIONS 5 192 static PetscInt numberofsetfromoptions; 193 static PetscErrorCode (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES); 194 195 #undef __FUNCT__ 196 #define __FUNCT__ "SNESAddOptionsChecker" 197 /*@C 198 SNESAddOptionsChecker - Adds an additional function to check for SNES options. 199 200 Not Collective 201 202 Input Parameter: 203 . snescheck - function that checks for options 204 205 Level: developer 206 207 .seealso: SNESSetFromOptions() 208 @*/ 209 PetscErrorCode SNESAddOptionsChecker(PetscErrorCode (*snescheck)(SNES)) 210 { 211 PetscFunctionBegin; 212 if (numberofsetfromoptions >= MAXSETFROMOPTIONS) { 213 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %D allowed", MAXSETFROMOPTIONS); 214 } 215 othersetfromoptions[numberofsetfromoptions++] = snescheck; 216 PetscFunctionReturn(0); 217 } 218 219 extern PetscErrorCode SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*); 220 221 #undef __FUNCT__ 222 #define __FUNCT__ "SNESSetUpMatrixFree_Private" 223 static PetscErrorCode SNESSetUpMatrixFree_Private(SNES snes, PetscBool hasOperator, PetscInt version) 224 { 225 Mat J; 226 KSP ksp; 227 PC pc; 228 PetscBool match; 229 PetscErrorCode ierr; 230 231 PetscFunctionBegin; 232 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 233 234 if(!snes->vec_func && (snes->jacobian || snes->jacobian_pre)) { 235 Mat A = snes->jacobian, B = snes->jacobian_pre; 236 ierr = MatGetVecs(A ? A : B, PETSC_NULL,&snes->vec_func);CHKERRQ(ierr); 237 } 238 239 if (version == 1) { 240 ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr); 241 ierr = MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);CHKERRQ(ierr); 242 ierr = MatSetFromOptions(J);CHKERRQ(ierr); 243 } else if (version == 2) { 244 if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"SNESSetFunction() must be called first"); 245 #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128) && !defined(PETSC_USE_REAL_LONG_DOUBLE) 246 ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_func,&J);CHKERRQ(ierr); 247 #else 248 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "matrix-free operator rutines (version 2)"); 249 #endif 250 } else { 251 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "matrix-free operator rutines, only version 1 and 2"); 252 } 253 254 ierr = PetscInfo1(snes,"Setting default matrix-free operator routines (version %D)\n", version);CHKERRQ(ierr); 255 if (hasOperator) { 256 /* This version replaces the user provided Jacobian matrix with a 257 matrix-free version but still employs the user-provided preconditioner matrix. */ 258 ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 259 } else { 260 /* This version replaces both the user-provided Jacobian and the user- 261 provided preconditioner matrix with the default matrix free version. */ 262 ierr = SNESSetJacobian(snes,J,J,MatMFFDComputeJacobian,snes->funP);CHKERRQ(ierr); 263 /* Force no preconditioner */ 264 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 265 ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 266 ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&match);CHKERRQ(ierr); 267 if (!match) { 268 ierr = PetscInfo(snes,"Setting default matrix-free preconditioner routines\nThat is no preconditioner is being used\n");CHKERRQ(ierr); 269 ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 270 } 271 } 272 ierr = MatDestroy(&J);CHKERRQ(ierr); 273 274 PetscFunctionReturn(0); 275 } 276 277 #undef __FUNCT__ 278 #define __FUNCT__ "SNESSetFromOptions" 279 /*@ 280 SNESSetFromOptions - Sets various SNES and KSP parameters from user options. 281 282 Collective on SNES 283 284 Input Parameter: 285 . snes - the SNES context 286 287 Options Database Keys: 288 + -snes_type <type> - ls, tr, umls, umtr, test 289 . -snes_stol - convergence tolerance in terms of the norm 290 of the change in the solution between steps 291 . -snes_atol <abstol> - absolute tolerance of residual norm 292 . -snes_rtol <rtol> - relative decrease in tolerance norm from initial 293 . -snes_max_it <max_it> - maximum number of iterations 294 . -snes_max_funcs <max_funcs> - maximum number of function evaluations 295 . -snes_max_fail <max_fail> - maximum number of failures 296 . -snes_max_linear_solve_fail - number of linear solver failures before SNESSolve() stops 297 . -snes_lag_preconditioner <lag> - how often preconditioner is rebuilt (use -1 to never rebuild) 298 . -snes_lag_jacobian <lag> - how often Jacobian is rebuilt (use -1 to never rebuild) 299 . -snes_trtol <trtol> - trust region tolerance 300 . -snes_no_convergence_test - skip convergence test in nonlinear 301 solver; hence iterations will continue until max_it 302 or some other criterion is reached. Saves expense 303 of convergence test 304 . -snes_monitor <optional filename> - prints residual norm at each iteration. if no 305 filename given prints to stdout 306 . -snes_monitor_solution - plots solution at each iteration 307 . -snes_monitor_residual - plots residual (not its norm) at each iteration 308 . -snes_monitor_solution_update - plots update to solution at each iteration 309 . -snes_monitor_draw - plots residual norm at each iteration 310 . -snes_fd - use finite differences to compute Jacobian; very slow, only for testing 311 . -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration 312 - -snes_converged_reason - print the reason for convergence/divergence after each solve 313 314 Options Database for Eisenstat-Walker method: 315 + -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence 316 . -snes_ksp_ew_version ver - version of Eisenstat-Walker method 317 . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0 318 . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax 319 . -snes_ksp_ew_gamma <gamma> - Sets gamma 320 . -snes_ksp_ew_alpha <alpha> - Sets alpha 321 . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2 322 - -snes_ksp_ew_threshold <threshold> - Sets threshold 323 324 Notes: 325 To see all options, run your program with the -help option or consult 326 the <A href="../../docs/manual.pdf#nameddest=ch_snes">SNES chapter of the users manual</A>. 327 328 Level: beginner 329 330 .keywords: SNES, nonlinear, set, options, database 331 332 .seealso: SNESSetOptionsPrefix() 333 @*/ 334 PetscErrorCode SNESSetFromOptions(SNES snes) 335 { 336 PetscBool flg,mf,mf_operator; 337 PetscInt i,indx,lag,mf_version,grids; 338 MatStructure matflag; 339 const char *deft = SNESLS; 340 const char *convtests[] = {"default","skip"}; 341 SNESKSPEW *kctx = NULL; 342 char type[256], monfilename[PETSC_MAX_PATH_LEN]; 343 PetscViewerASCIIMonitor monviewer; 344 PetscErrorCode ierr; 345 346 PetscFunctionBegin; 347 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 348 349 if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 350 ierr = PetscOptionsBegin(((PetscObject)snes)->comm,((PetscObject)snes)->prefix,"Nonlinear solver (SNES) options","SNES");CHKERRQ(ierr); 351 if (((PetscObject)snes)->type_name) { deft = ((PetscObject)snes)->type_name; } 352 ierr = PetscOptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr); 353 if (flg) { 354 ierr = SNESSetType(snes,type);CHKERRQ(ierr); 355 } else if (!((PetscObject)snes)->type_name) { 356 ierr = SNESSetType(snes,deft);CHKERRQ(ierr); 357 } 358 /* not used here, but called so will go into help messaage */ 359 ierr = PetscOptionsName("-snes_view","Print detailed information on solver used","SNESView",0);CHKERRQ(ierr); 360 361 ierr = PetscOptionsReal("-snes_stol","Stop if step length less than","SNESSetTolerances",snes->xtol,&snes->xtol,0);CHKERRQ(ierr); 362 ierr = PetscOptionsReal("-snes_atol","Stop if function norm less than","SNESSetTolerances",snes->abstol,&snes->abstol,0);CHKERRQ(ierr); 363 364 ierr = PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less than","SNESSetTolerances",snes->rtol,&snes->rtol,0);CHKERRQ(ierr); 365 ierr = PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,PETSC_NULL);CHKERRQ(ierr); 366 ierr = PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr); 367 ierr = PetscOptionsInt("-snes_max_fail","Maximum failures","SNESSetTolerances",snes->maxFailures,&snes->maxFailures,PETSC_NULL);CHKERRQ(ierr); 368 ierr = PetscOptionsInt("-snes_max_linear_solve_fail","Maximum failures in linear solves allowed","SNESSetMaxLinearSolveFailures",snes->maxLinearSolveFailures,&snes->maxLinearSolveFailures,PETSC_NULL);CHKERRQ(ierr); 369 ierr = PetscOptionsBool("-snes_error_if_not_converged","Generate error if solver does not converge","SNESSetErrorIfNotConverged",snes->errorifnotconverged,&snes->errorifnotconverged,PETSC_NULL);CHKERRQ(ierr); 370 371 ierr = PetscOptionsInt("-snes_lag_preconditioner","How often to rebuild preconditioner","SNESSetLagPreconditioner",snes->lagpreconditioner,&lag,&flg);CHKERRQ(ierr); 372 if (flg) { 373 ierr = SNESSetLagPreconditioner(snes,lag);CHKERRQ(ierr); 374 } 375 ierr = PetscOptionsInt("-snes_lag_jacobian","How often to rebuild Jacobian","SNESSetLagJacobian",snes->lagjacobian,&lag,&flg);CHKERRQ(ierr); 376 if (flg) { 377 ierr = SNESSetLagJacobian(snes,lag);CHKERRQ(ierr); 378 } 379 ierr = PetscOptionsInt("-snes_grid_sequence","Use grid sequencing to generate initial guess","SNESSetGridSequence",snes->gridsequence,&grids,&flg);CHKERRQ(ierr); 380 if (flg) { 381 ierr = SNESSetGridSequence(snes,grids);CHKERRQ(ierr); 382 } 383 384 ierr = PetscOptionsEList("-snes_convergence_test","Convergence test","SNESSetConvergenceTest",convtests,2,"default",&indx,&flg);CHKERRQ(ierr); 385 if (flg) { 386 switch (indx) { 387 case 0: ierr = SNESSetConvergenceTest(snes,SNESDefaultConverged,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); break; 388 case 1: ierr = SNESSetConvergenceTest(snes,SNESSkipConverged,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); break; 389 } 390 } 391 392 ierr = PetscOptionsBool("-snes_converged_reason","Print reason for converged or diverged","SNESSolve",snes->printreason,&snes->printreason,PETSC_NULL);CHKERRQ(ierr); 393 394 kctx = (SNESKSPEW *)snes->kspconvctx; 395 396 ierr = PetscOptionsBool("-snes_ksp_ew","Use Eisentat-Walker linear system convergence test","SNESKSPSetUseEW",snes->ksp_ewconv,&snes->ksp_ewconv,PETSC_NULL);CHKERRQ(ierr); 397 398 ierr = PetscOptionsInt("-snes_ksp_ew_version","Version 1, 2 or 3","SNESKSPSetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr); 399 ierr = PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNESKSPSetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr); 400 ierr = PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNESKSPSetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr); 401 ierr = PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNESKSPSetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr); 402 ierr = PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNESKSPSetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr); 403 ierr = PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNESKSPSetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr); 404 ierr = PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNESKSPSetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr); 405 406 flg = PETSC_FALSE; 407 ierr = PetscOptionsBool("-snes_monitor_cancel","Remove all monitors","SNESMonitorCancel",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 408 if (flg) {ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);} 409 410 ierr = PetscOptionsString("-snes_monitor","Monitor norm of function","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 411 if (flg) { 412 ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,monfilename,((PetscObject)snes)->tablevel,&monviewer);CHKERRQ(ierr); 413 ierr = SNESMonitorSet(snes,SNESMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 414 } 415 416 ierr = PetscOptionsString("-snes_monitor_range","Monitor range of elements of function","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 417 if (flg) { 418 ierr = SNESMonitorSet(snes,SNESMonitorRange,0,0);CHKERRQ(ierr); 419 } 420 421 ierr = PetscOptionsString("-snes_ratiomonitor","Monitor ratios of norms of function","SNESMonitorSetRatio","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 422 if (flg) { 423 ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,monfilename,((PetscObject)snes)->tablevel,&monviewer);CHKERRQ(ierr); 424 ierr = SNESMonitorSetRatio(snes,monviewer);CHKERRQ(ierr); 425 } 426 427 ierr = PetscOptionsString("-snes_monitor_short","Monitor norm of function (fewer digits)","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 428 if (flg) { 429 ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,monfilename,((PetscObject)snes)->tablevel,&monviewer);CHKERRQ(ierr); 430 ierr = SNESMonitorSet(snes,SNESMonitorDefaultShort,monviewer,(PetscErrorCode (*)(void**))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 431 } 432 433 flg = PETSC_FALSE; 434 ierr = PetscOptionsBool("-snes_monitor_solution","Plot solution at each iteration","SNESMonitorSolution",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 435 if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorSolution,0,0);CHKERRQ(ierr);} 436 flg = PETSC_FALSE; 437 ierr = PetscOptionsBool("-snes_monitor_solution_update","Plot correction at each iteration","SNESMonitorSolutionUpdate",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 438 if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorSolutionUpdate,0,0);CHKERRQ(ierr);} 439 flg = PETSC_FALSE; 440 ierr = PetscOptionsBool("-snes_monitor_residual","Plot residual at each iteration","SNESMonitorResidual",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 441 if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorResidual,0,0);CHKERRQ(ierr);} 442 flg = PETSC_FALSE; 443 ierr = PetscOptionsBool("-snes_monitor_draw","Plot function norm at each iteration","SNESMonitorLG",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 444 if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);} 445 flg = PETSC_FALSE; 446 ierr = PetscOptionsBool("-snes_monitor_range_draw","Plot function range at each iteration","SNESMonitorLG",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 447 if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorLGRange,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);} 448 449 flg = PETSC_FALSE; 450 ierr = PetscOptionsBool("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 451 if (flg) { 452 ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr); 453 ierr = PetscInfo(snes,"Setting default finite difference Jacobian matrix\n");CHKERRQ(ierr); 454 } 455 456 mf = mf_operator = PETSC_FALSE; 457 flg = PETSC_FALSE; 458 ierr = PetscOptionsBool("-snes_mf_operator","Use a Matrix-Free Jacobian with user-provided preconditioner matrix","MatCreateSNESMF",PETSC_FALSE,&mf_operator,&flg);CHKERRQ(ierr); 459 if (flg && mf_operator) mf = PETSC_TRUE; 460 flg = PETSC_FALSE; 461 ierr = PetscOptionsBool("-snes_mf","Use a Matrix-Free Jacobian with no preconditioner matrix","MatCreateSNESMF",PETSC_FALSE,&mf,&flg);CHKERRQ(ierr); 462 if (!flg && mf_operator) mf = PETSC_TRUE; 463 mf_version = 1; 464 ierr = PetscOptionsInt("-snes_mf_version","Matrix-Free routines version 1 or 2","None",mf_version,&mf_version,0);CHKERRQ(ierr); 465 466 for(i = 0; i < numberofsetfromoptions; i++) { 467 ierr = (*othersetfromoptions[i])(snes);CHKERRQ(ierr); 468 } 469 470 if (snes->ops->setfromoptions) { 471 ierr = (*snes->ops->setfromoptions)(snes);CHKERRQ(ierr); 472 } 473 474 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 475 ierr = PetscObjectProcessOptionsHandlers((PetscObject)snes);CHKERRQ(ierr); 476 ierr = PetscOptionsEnd();CHKERRQ(ierr); 477 478 if (mf) { ierr = SNESSetUpMatrixFree_Private(snes, mf_operator, mf_version);CHKERRQ(ierr); } 479 480 if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);} 481 ierr = KSPGetOperators(snes->ksp,PETSC_NULL,PETSC_NULL,&matflag); 482 ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,matflag);CHKERRQ(ierr); 483 ierr = KSPSetFromOptions(snes->ksp);CHKERRQ(ierr); 484 485 PetscFunctionReturn(0); 486 } 487 488 #undef __FUNCT__ 489 #define __FUNCT__ "SNESSetComputeApplicationContext" 490 /*@ 491 SNESSetComputeApplicationContext - Sets an optional function to compute a user-defined context for 492 the nonlinear solvers. 493 494 Logically Collective on SNES 495 496 Input Parameters: 497 + snes - the SNES context 498 . compute - function to compute the context 499 - destroy - function to destroy the context 500 501 Level: intermediate 502 503 .keywords: SNES, nonlinear, set, application, context 504 505 .seealso: SNESGetApplicationContext(), SNESSetComputeApplicationContext(), SNESGetApplicationContext() 506 @*/ 507 PetscErrorCode SNESSetComputeApplicationContext(SNES snes,PetscErrorCode (*compute)(SNES,void**),PetscErrorCode (*destroy)(void**)) 508 { 509 PetscFunctionBegin; 510 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 511 snes->ops->usercompute = compute; 512 snes->ops->userdestroy = destroy; 513 PetscFunctionReturn(0); 514 } 515 516 #undef __FUNCT__ 517 #define __FUNCT__ "SNESSetApplicationContext" 518 /*@ 519 SNESSetApplicationContext - Sets the optional user-defined context for 520 the nonlinear solvers. 521 522 Logically Collective on SNES 523 524 Input Parameters: 525 + snes - the SNES context 526 - usrP - optional user context 527 528 Level: intermediate 529 530 .keywords: SNES, nonlinear, set, application, context 531 532 .seealso: SNESGetApplicationContext(), SNESSetApplicationContext() 533 @*/ 534 PetscErrorCode SNESSetApplicationContext(SNES snes,void *usrP) 535 { 536 PetscErrorCode ierr; 537 538 PetscFunctionBegin; 539 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 540 ierr = KSPSetApplicationContext(snes->ksp,usrP);CHKERRQ(ierr); 541 snes->user = usrP; 542 PetscFunctionReturn(0); 543 } 544 545 #undef __FUNCT__ 546 #define __FUNCT__ "SNESGetApplicationContext" 547 /*@C 548 SNESGetApplicationContext - Gets the user-defined context for the 549 nonlinear solvers. 550 551 Not Collective 552 553 Input Parameter: 554 . snes - SNES context 555 556 Output Parameter: 557 . usrP - user context 558 559 Level: intermediate 560 561 .keywords: SNES, nonlinear, get, application, context 562 563 .seealso: SNESSetApplicationContext() 564 @*/ 565 PetscErrorCode SNESGetApplicationContext(SNES snes,void **usrP) 566 { 567 PetscFunctionBegin; 568 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 569 *usrP = snes->user; 570 PetscFunctionReturn(0); 571 } 572 573 #undef __FUNCT__ 574 #define __FUNCT__ "SNESGetIterationNumber" 575 /*@ 576 SNESGetIterationNumber - Gets the number of nonlinear iterations completed 577 at this time. 578 579 Not Collective 580 581 Input Parameter: 582 . snes - SNES context 583 584 Output Parameter: 585 . iter - iteration number 586 587 Notes: 588 For example, during the computation of iteration 2 this would return 1. 589 590 This is useful for using lagged Jacobians (where one does not recompute the 591 Jacobian at each SNES iteration). For example, the code 592 .vb 593 ierr = SNESGetIterationNumber(snes,&it); 594 if (!(it % 2)) { 595 [compute Jacobian here] 596 } 597 .ve 598 can be used in your ComputeJacobian() function to cause the Jacobian to be 599 recomputed every second SNES iteration. 600 601 Level: intermediate 602 603 .keywords: SNES, nonlinear, get, iteration, number, 604 605 .seealso: SNESGetFunctionNorm(), SNESGetLinearSolveIterations() 606 @*/ 607 PetscErrorCode SNESGetIterationNumber(SNES snes,PetscInt* iter) 608 { 609 PetscFunctionBegin; 610 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 611 PetscValidIntPointer(iter,2); 612 *iter = snes->iter; 613 PetscFunctionReturn(0); 614 } 615 616 #undef __FUNCT__ 617 #define __FUNCT__ "SNESGetFunctionNorm" 618 /*@ 619 SNESGetFunctionNorm - Gets the norm of the current function that was set 620 with SNESSSetFunction(). 621 622 Collective on SNES 623 624 Input Parameter: 625 . snes - SNES context 626 627 Output Parameter: 628 . fnorm - 2-norm of function 629 630 Level: intermediate 631 632 .keywords: SNES, nonlinear, get, function, norm 633 634 .seealso: SNESGetFunction(), SNESGetIterationNumber(), SNESGetLinearSolveIterations() 635 @*/ 636 PetscErrorCode SNESGetFunctionNorm(SNES snes,PetscReal *fnorm) 637 { 638 PetscFunctionBegin; 639 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 640 PetscValidScalarPointer(fnorm,2); 641 *fnorm = snes->norm; 642 PetscFunctionReturn(0); 643 } 644 645 #undef __FUNCT__ 646 #define __FUNCT__ "SNESGetNonlinearStepFailures" 647 /*@ 648 SNESGetNonlinearStepFailures - Gets the number of unsuccessful steps 649 attempted by the nonlinear solver. 650 651 Not Collective 652 653 Input Parameter: 654 . snes - SNES context 655 656 Output Parameter: 657 . nfails - number of unsuccessful steps attempted 658 659 Notes: 660 This counter is reset to zero for each successive call to SNESSolve(). 661 662 Level: intermediate 663 664 .keywords: SNES, nonlinear, get, number, unsuccessful, steps 665 666 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), 667 SNESSetMaxNonlinearStepFailures(), SNESGetMaxNonlinearStepFailures() 668 @*/ 669 PetscErrorCode SNESGetNonlinearStepFailures(SNES snes,PetscInt* nfails) 670 { 671 PetscFunctionBegin; 672 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 673 PetscValidIntPointer(nfails,2); 674 *nfails = snes->numFailures; 675 PetscFunctionReturn(0); 676 } 677 678 #undef __FUNCT__ 679 #define __FUNCT__ "SNESSetMaxNonlinearStepFailures" 680 /*@ 681 SNESSetMaxNonlinearStepFailures - Sets the maximum number of unsuccessful steps 682 attempted by the nonlinear solver before it gives up. 683 684 Not Collective 685 686 Input Parameters: 687 + snes - SNES context 688 - maxFails - maximum of unsuccessful steps 689 690 Level: intermediate 691 692 .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps 693 694 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), 695 SNESGetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures() 696 @*/ 697 PetscErrorCode SNESSetMaxNonlinearStepFailures(SNES snes, PetscInt maxFails) 698 { 699 PetscFunctionBegin; 700 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 701 snes->maxFailures = maxFails; 702 PetscFunctionReturn(0); 703 } 704 705 #undef __FUNCT__ 706 #define __FUNCT__ "SNESGetMaxNonlinearStepFailures" 707 /*@ 708 SNESGetMaxNonlinearStepFailures - Gets the maximum number of unsuccessful steps 709 attempted by the nonlinear solver before it gives up. 710 711 Not Collective 712 713 Input Parameter: 714 . snes - SNES context 715 716 Output Parameter: 717 . maxFails - maximum of unsuccessful steps 718 719 Level: intermediate 720 721 .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps 722 723 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), 724 SNESSetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures() 725 726 @*/ 727 PetscErrorCode SNESGetMaxNonlinearStepFailures(SNES snes, PetscInt *maxFails) 728 { 729 PetscFunctionBegin; 730 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 731 PetscValidIntPointer(maxFails,2); 732 *maxFails = snes->maxFailures; 733 PetscFunctionReturn(0); 734 } 735 736 #undef __FUNCT__ 737 #define __FUNCT__ "SNESGetNumberFunctionEvals" 738 /*@ 739 SNESGetNumberFunctionEvals - Gets the number of user provided function evaluations 740 done by SNES. 741 742 Not Collective 743 744 Input Parameter: 745 . snes - SNES context 746 747 Output Parameter: 748 . nfuncs - number of evaluations 749 750 Level: intermediate 751 752 .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps 753 754 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures() 755 @*/ 756 PetscErrorCode SNESGetNumberFunctionEvals(SNES snes, PetscInt *nfuncs) 757 { 758 PetscFunctionBegin; 759 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 760 PetscValidIntPointer(nfuncs,2); 761 *nfuncs = snes->nfuncs; 762 PetscFunctionReturn(0); 763 } 764 765 #undef __FUNCT__ 766 #define __FUNCT__ "SNESGetLinearSolveFailures" 767 /*@ 768 SNESGetLinearSolveFailures - Gets the number of failed (non-converged) 769 linear solvers. 770 771 Not Collective 772 773 Input Parameter: 774 . snes - SNES context 775 776 Output Parameter: 777 . nfails - number of failed solves 778 779 Notes: 780 This counter is reset to zero for each successive call to SNESSolve(). 781 782 Level: intermediate 783 784 .keywords: SNES, nonlinear, get, number, unsuccessful, steps 785 786 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures() 787 @*/ 788 PetscErrorCode SNESGetLinearSolveFailures(SNES snes,PetscInt* nfails) 789 { 790 PetscFunctionBegin; 791 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 792 PetscValidIntPointer(nfails,2); 793 *nfails = snes->numLinearSolveFailures; 794 PetscFunctionReturn(0); 795 } 796 797 #undef __FUNCT__ 798 #define __FUNCT__ "SNESSetMaxLinearSolveFailures" 799 /*@ 800 SNESSetMaxLinearSolveFailures - the number of failed linear solve attempts 801 allowed before SNES returns with a diverged reason of SNES_DIVERGED_LINEAR_SOLVE 802 803 Logically Collective on SNES 804 805 Input Parameters: 806 + snes - SNES context 807 - maxFails - maximum allowed linear solve failures 808 809 Level: intermediate 810 811 Notes: By default this is 0; that is SNES returns on the first failed linear solve 812 813 .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps 814 815 .seealso: SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations() 816 @*/ 817 PetscErrorCode SNESSetMaxLinearSolveFailures(SNES snes, PetscInt maxFails) 818 { 819 PetscFunctionBegin; 820 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 821 PetscValidLogicalCollectiveInt(snes,maxFails,2); 822 snes->maxLinearSolveFailures = maxFails; 823 PetscFunctionReturn(0); 824 } 825 826 #undef __FUNCT__ 827 #define __FUNCT__ "SNESGetMaxLinearSolveFailures" 828 /*@ 829 SNESGetMaxLinearSolveFailures - gets the maximum number of linear solve failures that 830 are allowed before SNES terminates 831 832 Not Collective 833 834 Input Parameter: 835 . snes - SNES context 836 837 Output Parameter: 838 . maxFails - maximum of unsuccessful solves allowed 839 840 Level: intermediate 841 842 Notes: By default this is 1; that is SNES returns on the first failed linear solve 843 844 .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps 845 846 .seealso: SNESGetLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), 847 @*/ 848 PetscErrorCode SNESGetMaxLinearSolveFailures(SNES snes, PetscInt *maxFails) 849 { 850 PetscFunctionBegin; 851 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 852 PetscValidIntPointer(maxFails,2); 853 *maxFails = snes->maxLinearSolveFailures; 854 PetscFunctionReturn(0); 855 } 856 857 #undef __FUNCT__ 858 #define __FUNCT__ "SNESGetLinearSolveIterations" 859 /*@ 860 SNESGetLinearSolveIterations - Gets the total number of linear iterations 861 used by the nonlinear solver. 862 863 Not Collective 864 865 Input Parameter: 866 . snes - SNES context 867 868 Output Parameter: 869 . lits - number of linear iterations 870 871 Notes: 872 This counter is reset to zero for each successive call to SNESSolve(). 873 874 Level: intermediate 875 876 .keywords: SNES, nonlinear, get, number, linear, iterations 877 878 .seealso: SNESGetIterationNumber(), SNESGetFunctionNorm(), SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures() 879 @*/ 880 PetscErrorCode SNESGetLinearSolveIterations(SNES snes,PetscInt* lits) 881 { 882 PetscFunctionBegin; 883 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 884 PetscValidIntPointer(lits,2); 885 *lits = snes->linear_its; 886 PetscFunctionReturn(0); 887 } 888 889 #undef __FUNCT__ 890 #define __FUNCT__ "SNESGetKSP" 891 /*@ 892 SNESGetKSP - Returns the KSP context for a SNES solver. 893 894 Not Collective, but if SNES object is parallel, then KSP object is parallel 895 896 Input Parameter: 897 . snes - the SNES context 898 899 Output Parameter: 900 . ksp - the KSP context 901 902 Notes: 903 The user can then directly manipulate the KSP context to set various 904 options, etc. Likewise, the user can then extract and manipulate the 905 PC contexts as well. 906 907 Level: beginner 908 909 .keywords: SNES, nonlinear, get, KSP, context 910 911 .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP() 912 @*/ 913 PetscErrorCode SNESGetKSP(SNES snes,KSP *ksp) 914 { 915 PetscErrorCode ierr; 916 917 PetscFunctionBegin; 918 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 919 PetscValidPointer(ksp,2); 920 921 if (!snes->ksp) { 922 ierr = KSPCreate(((PetscObject)snes)->comm,&snes->ksp);CHKERRQ(ierr); 923 ierr = PetscObjectIncrementTabLevel((PetscObject)snes->ksp,(PetscObject)snes,1);CHKERRQ(ierr); 924 ierr = PetscLogObjectParent(snes,snes->ksp);CHKERRQ(ierr); 925 } 926 *ksp = snes->ksp; 927 PetscFunctionReturn(0); 928 } 929 930 #undef __FUNCT__ 931 #define __FUNCT__ "SNESSetKSP" 932 /*@ 933 SNESSetKSP - Sets a KSP context for the SNES object to use 934 935 Not Collective, but the SNES and KSP objects must live on the same MPI_Comm 936 937 Input Parameters: 938 + snes - the SNES context 939 - ksp - the KSP context 940 941 Notes: 942 The SNES object already has its KSP object, you can obtain with SNESGetKSP() 943 so this routine is rarely needed. 944 945 The KSP object that is already in the SNES object has its reference count 946 decreased by one. 947 948 Level: developer 949 950 .keywords: SNES, nonlinear, get, KSP, context 951 952 .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP() 953 @*/ 954 PetscErrorCode SNESSetKSP(SNES snes,KSP ksp) 955 { 956 PetscErrorCode ierr; 957 958 PetscFunctionBegin; 959 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 960 PetscValidHeaderSpecific(ksp,KSP_CLASSID,2); 961 PetscCheckSameComm(snes,1,ksp,2); 962 ierr = PetscObjectReference((PetscObject)ksp);CHKERRQ(ierr); 963 if (snes->ksp) {ierr = PetscObjectDereference((PetscObject)snes->ksp);CHKERRQ(ierr);} 964 snes->ksp = ksp; 965 PetscFunctionReturn(0); 966 } 967 968 #if 0 969 #undef __FUNCT__ 970 #define __FUNCT__ "SNESPublish_Petsc" 971 static PetscErrorCode SNESPublish_Petsc(PetscObject obj) 972 { 973 PetscFunctionBegin; 974 PetscFunctionReturn(0); 975 } 976 #endif 977 978 /* -----------------------------------------------------------*/ 979 #undef __FUNCT__ 980 #define __FUNCT__ "SNESCreate" 981 /*@ 982 SNESCreate - Creates a nonlinear solver context. 983 984 Collective on MPI_Comm 985 986 Input Parameters: 987 . comm - MPI communicator 988 989 Output Parameter: 990 . outsnes - the new SNES context 991 992 Options Database Keys: 993 + -snes_mf - Activates default matrix-free Jacobian-vector products, 994 and no preconditioning matrix 995 . -snes_mf_operator - Activates default matrix-free Jacobian-vector 996 products, and a user-provided preconditioning matrix 997 as set by SNESSetJacobian() 998 - -snes_fd - Uses (slow!) finite differences to compute Jacobian 999 1000 Level: beginner 1001 1002 .keywords: SNES, nonlinear, create, context 1003 1004 .seealso: SNESSolve(), SNESDestroy(), SNES, SNESSetLagPreconditioner() 1005 1006 @*/ 1007 PetscErrorCode SNESCreate(MPI_Comm comm,SNES *outsnes) 1008 { 1009 PetscErrorCode ierr; 1010 SNES snes; 1011 SNESKSPEW *kctx; 1012 1013 PetscFunctionBegin; 1014 PetscValidPointer(outsnes,2); 1015 *outsnes = PETSC_NULL; 1016 #ifndef PETSC_USE_DYNAMIC_LIBRARIES 1017 ierr = SNESInitializePackage(PETSC_NULL);CHKERRQ(ierr); 1018 #endif 1019 1020 ierr = PetscHeaderCreate(snes,_p_SNES,struct _SNESOps,SNES_CLASSID,0,"SNES",comm,SNESDestroy,SNESView);CHKERRQ(ierr); 1021 1022 snes->ops->converged = SNESDefaultConverged; 1023 snes->max_its = 50; 1024 snes->max_funcs = 10000; 1025 snes->norm = 0.0; 1026 snes->rtol = 1.e-8; 1027 snes->ttol = 0.0; 1028 snes->abstol = 1.e-50; 1029 snes->xtol = 1.e-8; 1030 snes->deltatol = 1.e-12; 1031 snes->nfuncs = 0; 1032 snes->numFailures = 0; 1033 snes->maxFailures = 1; 1034 snes->linear_its = 0; 1035 snes->lagjacobian = 1; 1036 snes->lagpreconditioner = 1; 1037 snes->numbermonitors = 0; 1038 snes->data = 0; 1039 snes->setupcalled = PETSC_FALSE; 1040 snes->ksp_ewconv = PETSC_FALSE; 1041 snes->nwork = 0; 1042 snes->work = 0; 1043 snes->nvwork = 0; 1044 snes->vwork = 0; 1045 snes->conv_hist_len = 0; 1046 snes->conv_hist_max = 0; 1047 snes->conv_hist = PETSC_NULL; 1048 snes->conv_hist_its = PETSC_NULL; 1049 snes->conv_hist_reset = PETSC_TRUE; 1050 snes->reason = SNES_CONVERGED_ITERATING; 1051 1052 snes->numLinearSolveFailures = 0; 1053 snes->maxLinearSolveFailures = 1; 1054 1055 /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 1056 ierr = PetscNewLog(snes,SNESKSPEW,&kctx);CHKERRQ(ierr); 1057 snes->kspconvctx = (void*)kctx; 1058 kctx->version = 2; 1059 kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 1060 this was too large for some test cases */ 1061 kctx->rtol_last = 0.0; 1062 kctx->rtol_max = .9; 1063 kctx->gamma = 1.0; 1064 kctx->alpha = .5*(1.0 + sqrt(5.0)); 1065 kctx->alpha2 = kctx->alpha; 1066 kctx->threshold = .1; 1067 kctx->lresid_last = 0.0; 1068 kctx->norm_last = 0.0; 1069 1070 *outsnes = snes; 1071 PetscFunctionReturn(0); 1072 } 1073 1074 #undef __FUNCT__ 1075 #define __FUNCT__ "SNESSetFunction" 1076 /*@C 1077 SNESSetFunction - Sets the function evaluation routine and function 1078 vector for use by the SNES routines in solving systems of nonlinear 1079 equations. 1080 1081 Logically Collective on SNES 1082 1083 Input Parameters: 1084 + snes - the SNES context 1085 . r - vector to store function value 1086 . func - function evaluation routine 1087 - ctx - [optional] user-defined context for private data for the 1088 function evaluation routine (may be PETSC_NULL) 1089 1090 Calling sequence of func: 1091 $ func (SNES snes,Vec x,Vec f,void *ctx); 1092 1093 . f - function vector 1094 - ctx - optional user-defined function context 1095 1096 Notes: 1097 The Newton-like methods typically solve linear systems of the form 1098 $ f'(x) x = -f(x), 1099 where f'(x) denotes the Jacobian matrix and f(x) is the function. 1100 1101 Level: beginner 1102 1103 .keywords: SNES, nonlinear, set, function 1104 1105 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 1106 @*/ 1107 PetscErrorCode SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx) 1108 { 1109 PetscErrorCode ierr; 1110 PetscFunctionBegin; 1111 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1112 if (r) { 1113 PetscValidHeaderSpecific(r,VEC_CLASSID,2); 1114 PetscCheckSameComm(snes,1,r,2); 1115 ierr = PetscObjectReference((PetscObject)r);CHKERRQ(ierr); 1116 ierr = VecDestroy(&snes->vec_func);CHKERRQ(ierr); 1117 snes->vec_func = r; 1118 } else if (!snes->vec_func && snes->dm) { 1119 ierr = DMCreateGlobalVector(snes->dm,&snes->vec_func);CHKERRQ(ierr); 1120 } 1121 if (func) snes->ops->computefunction = func; 1122 if (ctx) snes->funP = ctx; 1123 PetscFunctionReturn(0); 1124 } 1125 1126 #undef __FUNCT__ 1127 #define __FUNCT__ "SNESSetComputeInitialGuess" 1128 /*@C 1129 SNESSetComputeInitialGuess - Sets a routine used to compute an initial guess for the problem 1130 1131 Logically Collective on SNES 1132 1133 Input Parameters: 1134 + snes - the SNES context 1135 . func - function evaluation routine 1136 - ctx - [optional] user-defined context for private data for the 1137 function evaluation routine (may be PETSC_NULL) 1138 1139 Calling sequence of func: 1140 $ func (SNES snes,Vec x,void *ctx); 1141 1142 . f - function vector 1143 - ctx - optional user-defined function context 1144 1145 Level: intermediate 1146 1147 .keywords: SNES, nonlinear, set, function 1148 1149 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 1150 @*/ 1151 PetscErrorCode SNESSetComputeInitialGuess(SNES snes,PetscErrorCode (*func)(SNES,Vec,void*),void *ctx) 1152 { 1153 PetscFunctionBegin; 1154 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1155 if (func) snes->ops->computeinitialguess = func; 1156 if (ctx) snes->initialguessP = ctx; 1157 PetscFunctionReturn(0); 1158 } 1159 1160 /* --------------------------------------------------------------- */ 1161 #undef __FUNCT__ 1162 #define __FUNCT__ "SNESGetRhs" 1163 /*@C 1164 SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set 1165 it assumes a zero right hand side. 1166 1167 Logically Collective on SNES 1168 1169 Input Parameter: 1170 . snes - the SNES context 1171 1172 Output Parameter: 1173 . rhs - the right hand side vector or PETSC_NULL if the right hand side vector is null 1174 1175 Level: intermediate 1176 1177 .keywords: SNES, nonlinear, get, function, right hand side 1178 1179 .seealso: SNESGetSolution(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction() 1180 @*/ 1181 PetscErrorCode SNESGetRhs(SNES snes,Vec *rhs) 1182 { 1183 PetscFunctionBegin; 1184 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1185 PetscValidPointer(rhs,2); 1186 *rhs = snes->vec_rhs; 1187 PetscFunctionReturn(0); 1188 } 1189 1190 #undef __FUNCT__ 1191 #define __FUNCT__ "SNESComputeFunction" 1192 /*@ 1193 SNESComputeFunction - Calls the function that has been set with 1194 SNESSetFunction(). 1195 1196 Collective on SNES 1197 1198 Input Parameters: 1199 + snes - the SNES context 1200 - x - input vector 1201 1202 Output Parameter: 1203 . y - function vector, as set by SNESSetFunction() 1204 1205 Notes: 1206 SNESComputeFunction() is typically used within nonlinear solvers 1207 implementations, so most users would not generally call this routine 1208 themselves. 1209 1210 Level: developer 1211 1212 .keywords: SNES, nonlinear, compute, function 1213 1214 .seealso: SNESSetFunction(), SNESGetFunction() 1215 @*/ 1216 PetscErrorCode SNESComputeFunction(SNES snes,Vec x,Vec y) 1217 { 1218 PetscErrorCode ierr; 1219 1220 PetscFunctionBegin; 1221 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1222 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 1223 PetscValidHeaderSpecific(y,VEC_CLASSID,3); 1224 PetscCheckSameComm(snes,1,x,2); 1225 PetscCheckSameComm(snes,1,y,3); 1226 1227 ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 1228 if (snes->ops->computefunction) { 1229 PetscStackPush("SNES user function"); 1230 ierr = (*snes->ops->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr); 1231 PetscStackPop; 1232 } else if (snes->vec_rhs) { 1233 ierr = MatMult(snes->jacobian, x, y);CHKERRQ(ierr); 1234 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() before SNESComputeFunction(), likely called from SNESSolve()."); 1235 if (snes->vec_rhs) { 1236 ierr = VecAXPY(y,-1.0,snes->vec_rhs);CHKERRQ(ierr); 1237 } 1238 snes->nfuncs++; 1239 ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 1240 PetscFunctionReturn(0); 1241 } 1242 1243 #undef __FUNCT__ 1244 #define __FUNCT__ "SNESComputeJacobian" 1245 /*@ 1246 SNESComputeJacobian - Computes the Jacobian matrix that has been 1247 set with SNESSetJacobian(). 1248 1249 Collective on SNES and Mat 1250 1251 Input Parameters: 1252 + snes - the SNES context 1253 - x - input vector 1254 1255 Output Parameters: 1256 + A - Jacobian matrix 1257 . B - optional preconditioning matrix 1258 - flag - flag indicating matrix structure (one of, SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER) 1259 1260 Options Database Keys: 1261 + -snes_lag_preconditioner <lag> 1262 - -snes_lag_jacobian <lag> 1263 1264 Notes: 1265 Most users should not need to explicitly call this routine, as it 1266 is used internally within the nonlinear solvers. 1267 1268 See KSPSetOperators() for important information about setting the 1269 flag parameter. 1270 1271 Level: developer 1272 1273 .keywords: SNES, compute, Jacobian, matrix 1274 1275 .seealso: SNESSetJacobian(), KSPSetOperators(), MatStructure, SNESSetLagPreconditioner(), SNESSetLagJacobian() 1276 @*/ 1277 PetscErrorCode SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 1278 { 1279 PetscErrorCode ierr; 1280 PetscBool flag; 1281 1282 PetscFunctionBegin; 1283 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1284 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 1285 PetscValidPointer(flg,5); 1286 PetscCheckSameComm(snes,1,X,2); 1287 if (!snes->ops->computejacobian) PetscFunctionReturn(0); 1288 1289 /* make sure that MatAssemblyBegin/End() is called on A matrix if it is matrix free */ 1290 1291 if (snes->lagjacobian == -2) { 1292 snes->lagjacobian = -1; 1293 ierr = PetscInfo(snes,"Recomputing Jacobian/preconditioner because lag is -2 (means compute Jacobian, but then never again) \n");CHKERRQ(ierr); 1294 } else if (snes->lagjacobian == -1) { 1295 *flg = SAME_PRECONDITIONER; 1296 ierr = PetscInfo(snes,"Reusing Jacobian/preconditioner because lag is -1\n");CHKERRQ(ierr); 1297 ierr = PetscTypeCompare((PetscObject)*A,MATMFFD,&flag);CHKERRQ(ierr); 1298 if (flag) { 1299 ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1300 ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1301 } 1302 PetscFunctionReturn(0); 1303 } else if (snes->lagjacobian > 1 && snes->iter % snes->lagjacobian) { 1304 *flg = SAME_PRECONDITIONER; 1305 ierr = PetscInfo2(snes,"Reusing Jacobian/preconditioner because lag is %D and SNES iteration is %D\n",snes->lagjacobian,snes->iter);CHKERRQ(ierr); 1306 ierr = PetscTypeCompare((PetscObject)*A,MATMFFD,&flag);CHKERRQ(ierr); 1307 if (flag) { 1308 ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1309 ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1310 } 1311 PetscFunctionReturn(0); 1312 } 1313 1314 *flg = DIFFERENT_NONZERO_PATTERN; 1315 ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 1316 PetscStackPush("SNES user Jacobian function"); 1317 ierr = (*snes->ops->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr); 1318 PetscStackPop; 1319 ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 1320 1321 if (snes->lagpreconditioner == -2) { 1322 ierr = PetscInfo(snes,"Rebuilding preconditioner exactly once since lag is -2\n");CHKERRQ(ierr); 1323 snes->lagpreconditioner = -1; 1324 } else if (snes->lagpreconditioner == -1) { 1325 *flg = SAME_PRECONDITIONER; 1326 ierr = PetscInfo(snes,"Reusing preconditioner because lag is -1\n");CHKERRQ(ierr); 1327 } else if (snes->lagpreconditioner > 1 && snes->iter % snes->lagpreconditioner) { 1328 *flg = SAME_PRECONDITIONER; 1329 ierr = PetscInfo2(snes,"Reusing preconditioner because lag is %D and SNES iteration is %D\n",snes->lagpreconditioner,snes->iter);CHKERRQ(ierr); 1330 } 1331 1332 /* make sure user returned a correct Jacobian and preconditioner */ 1333 /* PetscValidHeaderSpecific(*A,MAT_CLASSID,3); 1334 PetscValidHeaderSpecific(*B,MAT_CLASSID,4); */ 1335 PetscFunctionReturn(0); 1336 } 1337 1338 #undef __FUNCT__ 1339 #define __FUNCT__ "SNESSetJacobian" 1340 /*@C 1341 SNESSetJacobian - Sets the function to compute Jacobian as well as the 1342 location to store the matrix. 1343 1344 Logically Collective on SNES and Mat 1345 1346 Input Parameters: 1347 + snes - the SNES context 1348 . A - Jacobian matrix 1349 . B - preconditioner matrix (usually same as the Jacobian) 1350 . func - Jacobian evaluation routine (if PETSC_NULL then SNES retains any previously set value) 1351 - ctx - [optional] user-defined context for private data for the 1352 Jacobian evaluation routine (may be PETSC_NULL) (if PETSC_NULL then SNES retains any previously set value) 1353 1354 Calling sequence of func: 1355 $ func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 1356 1357 + x - input vector 1358 . A - Jacobian matrix 1359 . B - preconditioner matrix, usually the same as A 1360 . flag - flag indicating information about the preconditioner matrix 1361 structure (same as flag in KSPSetOperators()), one of SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER 1362 - ctx - [optional] user-defined Jacobian context 1363 1364 Notes: 1365 See KSPSetOperators() for important information about setting the flag 1366 output parameter in the routine func(). Be sure to read this information! 1367 1368 The routine func() takes Mat * as the matrix arguments rather than Mat. 1369 This allows the Jacobian evaluation routine to replace A and/or B with a 1370 completely new new matrix structure (not just different matrix elements) 1371 when appropriate, for instance, if the nonzero structure is changing 1372 throughout the global iterations. 1373 1374 If the A matrix and B matrix are different you must call MatAssemblyBegin/End() on 1375 each matrix. 1376 1377 If using SNESDefaultComputeJacobianColor() to assemble a Jacobian, the ctx argument 1378 must be a MatFDColoring. 1379 1380 Level: beginner 1381 1382 .keywords: SNES, nonlinear, set, Jacobian, matrix 1383 1384 .seealso: KSPSetOperators(), SNESSetFunction(), MatMFFDComputeJacobian(), SNESDefaultComputeJacobianColor(), MatStructure 1385 @*/ 1386 PetscErrorCode SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 1387 { 1388 PetscErrorCode ierr; 1389 1390 PetscFunctionBegin; 1391 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1392 if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2); 1393 if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3); 1394 if (A) PetscCheckSameComm(snes,1,A,2); 1395 if (B) PetscCheckSameComm(snes,1,B,3); 1396 if (func) snes->ops->computejacobian = func; 1397 if (ctx) snes->jacP = ctx; 1398 if (A) { 1399 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 1400 ierr = MatDestroy(&snes->jacobian);CHKERRQ(ierr); 1401 snes->jacobian = A; 1402 } 1403 if (B) { 1404 ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 1405 ierr = MatDestroy(&snes->jacobian_pre);CHKERRQ(ierr); 1406 snes->jacobian_pre = B; 1407 } 1408 PetscFunctionReturn(0); 1409 } 1410 1411 #undef __FUNCT__ 1412 #define __FUNCT__ "SNESGetJacobian" 1413 /*@C 1414 SNESGetJacobian - Returns the Jacobian matrix and optionally the user 1415 provided context for evaluating the Jacobian. 1416 1417 Not Collective, but Mat object will be parallel if SNES object is 1418 1419 Input Parameter: 1420 . snes - the nonlinear solver context 1421 1422 Output Parameters: 1423 + A - location to stash Jacobian matrix (or PETSC_NULL) 1424 . B - location to stash preconditioner matrix (or PETSC_NULL) 1425 . func - location to put Jacobian function (or PETSC_NULL) 1426 - ctx - location to stash Jacobian ctx (or PETSC_NULL) 1427 1428 Level: advanced 1429 1430 .seealso: SNESSetJacobian(), SNESComputeJacobian() 1431 @*/ 1432 PetscErrorCode SNESGetJacobian(SNES snes,Mat *A,Mat *B,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx) 1433 { 1434 PetscFunctionBegin; 1435 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1436 if (A) *A = snes->jacobian; 1437 if (B) *B = snes->jacobian_pre; 1438 if (func) *func = snes->ops->computejacobian; 1439 if (ctx) *ctx = snes->jacP; 1440 PetscFunctionReturn(0); 1441 } 1442 1443 /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 1444 1445 #undef __FUNCT__ 1446 #define __FUNCT__ "SNESSetUp" 1447 /*@ 1448 SNESSetUp - Sets up the internal data structures for the later use 1449 of a nonlinear solver. 1450 1451 Collective on SNES 1452 1453 Input Parameters: 1454 . snes - the SNES context 1455 1456 Notes: 1457 For basic use of the SNES solvers the user need not explicitly call 1458 SNESSetUp(), since these actions will automatically occur during 1459 the call to SNESSolve(). However, if one wishes to control this 1460 phase separately, SNESSetUp() should be called after SNESCreate() 1461 and optional routines of the form SNESSetXXX(), but before SNESSolve(). 1462 1463 Level: advanced 1464 1465 .keywords: SNES, nonlinear, setup 1466 1467 .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 1468 @*/ 1469 PetscErrorCode SNESSetUp(SNES snes) 1470 { 1471 PetscErrorCode ierr; 1472 1473 PetscFunctionBegin; 1474 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1475 if (snes->setupcalled) PetscFunctionReturn(0); 1476 1477 if (!((PetscObject)snes)->type_name) { 1478 ierr = SNESSetType(snes,SNESLS);CHKERRQ(ierr); 1479 } 1480 1481 if (!snes->vec_func && snes->dm && !snes->vec_rhs) { 1482 ierr = DMCreateGlobalVector(snes->dm,&snes->vec_func);CHKERRQ(ierr); 1483 } 1484 if (!snes->vec_func && !snes->vec_rhs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 1485 if (!snes->ops->computefunction && !snes->vec_rhs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 1486 if (snes->vec_func == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector"); 1487 if (snes->vec_rhs == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be right hand side vector"); 1488 1489 if (!snes->vec_func && snes->vec_rhs) { 1490 ierr = VecDuplicate(snes->vec_rhs, &snes->vec_func);CHKERRQ(ierr); 1491 } 1492 if (!snes->vec_sol_update /* && snes->vec_sol */) { 1493 ierr = VecDuplicate(snes->vec_sol,&snes->vec_sol_update);CHKERRQ(ierr); 1494 ierr = PetscLogObjectParent(snes,snes->vec_sol_update);CHKERRQ(ierr); 1495 } 1496 1497 if (!snes->ops->computejacobian && snes->dm) { 1498 Mat J; 1499 ISColoring coloring; 1500 MatFDColoring fd; 1501 1502 ierr = DMGetMatrix(snes->dm,MATAIJ,&J);CHKERRQ(ierr); 1503 ierr = DMGetColoring(snes->dm,IS_COLORING_GHOSTED,MATAIJ,&coloring);CHKERRQ(ierr); 1504 ierr = MatFDColoringCreate(J,coloring,&fd);CHKERRQ(ierr); 1505 ierr = MatFDColoringSetFunction(fd,(PetscErrorCode (*)(void))snes->ops->computefunction,snes->funP);CHKERRQ(ierr); 1506 ierr = SNESSetJacobian(snes,J,J,SNESDefaultComputeJacobianColor,fd);CHKERRQ(ierr); 1507 ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr); 1508 } else if (snes->dm && !snes->jacobian_pre){ 1509 Mat J; 1510 ierr = DMGetMatrix(snes->dm,MATAIJ,&J);CHKERRQ(ierr); 1511 ierr = SNESSetJacobian(snes,J,J,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 1512 ierr = MatDestroy(&J);CHKERRQ(ierr); 1513 } 1514 1515 if (!snes->ksp) {ierr = SNESGetKSP(snes, &snes->ksp);CHKERRQ(ierr);} 1516 1517 if (snes->ops->usercompute && !snes->user) { 1518 ierr = (*snes->ops->usercompute)(snes,(void**)&snes->user);CHKERRQ(ierr); 1519 } 1520 1521 if (snes->ops->setup) { 1522 ierr = (*snes->ops->setup)(snes);CHKERRQ(ierr); 1523 } 1524 1525 snes->setupcalled = PETSC_TRUE; 1526 PetscFunctionReturn(0); 1527 } 1528 1529 #undef __FUNCT__ 1530 #define __FUNCT__ "SNESReset" 1531 /*@ 1532 SNESReset - Resets a SNES context to the snessetupcalled = 0 state and removes any allocated Vecs and Mats 1533 1534 Collective on SNES 1535 1536 Input Parameter: 1537 . snes - iterative context obtained from SNESCreate() 1538 1539 Level: intermediate 1540 1541 Notes: Also calls the application context destroy routine set with SNESSetComputeApplicationContext() 1542 1543 .keywords: SNES, destroy 1544 1545 .seealso: SNESCreate(), SNESSetUp(), SNESSolve() 1546 @*/ 1547 PetscErrorCode SNESReset(SNES snes) 1548 { 1549 PetscErrorCode ierr; 1550 1551 PetscFunctionBegin; 1552 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1553 if (snes->ops->userdestroy && snes->user) { 1554 ierr = (*snes->ops->userdestroy)((void**)&snes->user);CHKERRQ(ierr); 1555 snes->user = PETSC_NULL; 1556 } 1557 1558 if (snes->ops->reset) { 1559 ierr = (*snes->ops->reset)(snes);CHKERRQ(ierr); 1560 } 1561 if (snes->ksp) {ierr = KSPReset(snes->ksp);CHKERRQ(ierr);} 1562 ierr = VecDestroy(&snes->vec_rhs);CHKERRQ(ierr); 1563 ierr = VecDestroy(&snes->vec_sol);CHKERRQ(ierr); 1564 ierr = VecDestroy(&snes->vec_sol_update);CHKERRQ(ierr); 1565 ierr = VecDestroy(&snes->vec_func);CHKERRQ(ierr); 1566 ierr = MatDestroy(&snes->jacobian);CHKERRQ(ierr); 1567 ierr = MatDestroy(&snes->jacobian_pre);CHKERRQ(ierr); 1568 if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);} 1569 if (snes->vwork) {ierr = VecDestroyVecs(snes->nvwork,&snes->vwork);CHKERRQ(ierr);} 1570 snes->nwork = snes->nvwork = 0; 1571 snes->setupcalled = PETSC_FALSE; 1572 PetscFunctionReturn(0); 1573 } 1574 1575 #undef __FUNCT__ 1576 #define __FUNCT__ "SNESDestroy" 1577 /*@ 1578 SNESDestroy - Destroys the nonlinear solver context that was created 1579 with SNESCreate(). 1580 1581 Collective on SNES 1582 1583 Input Parameter: 1584 . snes - the SNES context 1585 1586 Level: beginner 1587 1588 .keywords: SNES, nonlinear, destroy 1589 1590 .seealso: SNESCreate(), SNESSolve() 1591 @*/ 1592 PetscErrorCode SNESDestroy(SNES *snes) 1593 { 1594 PetscErrorCode ierr; 1595 1596 PetscFunctionBegin; 1597 if (!*snes) PetscFunctionReturn(0); 1598 PetscValidHeaderSpecific((*snes),SNES_CLASSID,1); 1599 if (--((PetscObject)(*snes))->refct > 0) {*snes = 0; PetscFunctionReturn(0);} 1600 1601 ierr = SNESReset((*snes));CHKERRQ(ierr); 1602 1603 /* if memory was published with AMS then destroy it */ 1604 ierr = PetscObjectDepublish((*snes));CHKERRQ(ierr); 1605 if ((*snes)->ops->destroy) {ierr = (*((*snes))->ops->destroy)((*snes));CHKERRQ(ierr);} 1606 1607 ierr = DMDestroy(&(*snes)->dm);CHKERRQ(ierr); 1608 ierr = KSPDestroy(&(*snes)->ksp);CHKERRQ(ierr); 1609 1610 ierr = PetscFree((*snes)->kspconvctx);CHKERRQ(ierr); 1611 if ((*snes)->ops->convergeddestroy) { 1612 ierr = (*(*snes)->ops->convergeddestroy)((*snes)->cnvP);CHKERRQ(ierr); 1613 } 1614 if ((*snes)->conv_malloc) { 1615 ierr = PetscFree((*snes)->conv_hist);CHKERRQ(ierr); 1616 ierr = PetscFree((*snes)->conv_hist_its);CHKERRQ(ierr); 1617 } 1618 ierr = SNESMonitorCancel((*snes));CHKERRQ(ierr); 1619 ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr); 1620 PetscFunctionReturn(0); 1621 } 1622 1623 /* ----------- Routines to set solver parameters ---------- */ 1624 1625 #undef __FUNCT__ 1626 #define __FUNCT__ "SNESSetLagPreconditioner" 1627 /*@ 1628 SNESSetLagPreconditioner - Determines when the preconditioner is rebuilt in the nonlinear solve. 1629 1630 Logically Collective on SNES 1631 1632 Input Parameters: 1633 + snes - the SNES context 1634 - lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time 1635 the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that 1636 1637 Options Database Keys: 1638 . -snes_lag_preconditioner <lag> 1639 1640 Notes: 1641 The default is 1 1642 The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 1643 If -1 is used before the very first nonlinear solve the preconditioner is still built because there is no previous preconditioner to use 1644 1645 Level: intermediate 1646 1647 .keywords: SNES, nonlinear, set, convergence, tolerances 1648 1649 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian() 1650 1651 @*/ 1652 PetscErrorCode SNESSetLagPreconditioner(SNES snes,PetscInt lag) 1653 { 1654 PetscFunctionBegin; 1655 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1656 if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater"); 1657 if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0"); 1658 PetscValidLogicalCollectiveInt(snes,lag,2); 1659 snes->lagpreconditioner = lag; 1660 PetscFunctionReturn(0); 1661 } 1662 1663 #undef __FUNCT__ 1664 #define __FUNCT__ "SNESSetGridSequence" 1665 /*@ 1666 SNESSetGridSequence - sets the number of steps of grid sequencing that SNES does 1667 1668 Logically Collective on SNES 1669 1670 Input Parameters: 1671 + snes - the SNES context 1672 - steps - the number of refinements to do, defaults to 0 1673 1674 Options Database Keys: 1675 . -snes_grid_sequence <steps> 1676 1677 Level: intermediate 1678 1679 .keywords: SNES, nonlinear, set, convergence, tolerances 1680 1681 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian() 1682 1683 @*/ 1684 PetscErrorCode SNESSetGridSequence(SNES snes,PetscInt steps) 1685 { 1686 PetscFunctionBegin; 1687 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1688 PetscValidLogicalCollectiveInt(snes,steps,2); 1689 snes->gridsequence = steps; 1690 PetscFunctionReturn(0); 1691 } 1692 1693 #undef __FUNCT__ 1694 #define __FUNCT__ "SNESGetLagPreconditioner" 1695 /*@ 1696 SNESGetLagPreconditioner - Indicates how often the preconditioner is rebuilt 1697 1698 Not Collective 1699 1700 Input Parameter: 1701 . snes - the SNES context 1702 1703 Output Parameter: 1704 . lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time 1705 the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that 1706 1707 Options Database Keys: 1708 . -snes_lag_preconditioner <lag> 1709 1710 Notes: 1711 The default is 1 1712 The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 1713 1714 Level: intermediate 1715 1716 .keywords: SNES, nonlinear, set, convergence, tolerances 1717 1718 .seealso: SNESSetTrustRegionTolerance(), SNESSetLagPreconditioner() 1719 1720 @*/ 1721 PetscErrorCode SNESGetLagPreconditioner(SNES snes,PetscInt *lag) 1722 { 1723 PetscFunctionBegin; 1724 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1725 *lag = snes->lagpreconditioner; 1726 PetscFunctionReturn(0); 1727 } 1728 1729 #undef __FUNCT__ 1730 #define __FUNCT__ "SNESSetLagJacobian" 1731 /*@ 1732 SNESSetLagJacobian - Determines when the Jacobian is rebuilt in the nonlinear solve. See SNESSetLagPreconditioner() for determining how 1733 often the preconditioner is rebuilt. 1734 1735 Logically Collective on SNES 1736 1737 Input Parameters: 1738 + snes - the SNES context 1739 - lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time 1740 the Jacobian is built etc. -2 means rebuild at next chance but then never again 1741 1742 Options Database Keys: 1743 . -snes_lag_jacobian <lag> 1744 1745 Notes: 1746 The default is 1 1747 The Jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 1748 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 1749 at the next Newton step but never again (unless it is reset to another value) 1750 1751 Level: intermediate 1752 1753 .keywords: SNES, nonlinear, set, convergence, tolerances 1754 1755 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagPreconditioner(), SNESGetLagJacobian() 1756 1757 @*/ 1758 PetscErrorCode SNESSetLagJacobian(SNES snes,PetscInt lag) 1759 { 1760 PetscFunctionBegin; 1761 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1762 if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater"); 1763 if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0"); 1764 PetscValidLogicalCollectiveInt(snes,lag,2); 1765 snes->lagjacobian = lag; 1766 PetscFunctionReturn(0); 1767 } 1768 1769 #undef __FUNCT__ 1770 #define __FUNCT__ "SNESGetLagJacobian" 1771 /*@ 1772 SNESGetLagJacobian - Indicates how often the Jacobian is rebuilt. See SNESGetLagPreconditioner() to determine when the preconditioner is rebuilt 1773 1774 Not Collective 1775 1776 Input Parameter: 1777 . snes - the SNES context 1778 1779 Output Parameter: 1780 . lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time 1781 the Jacobian is built etc. 1782 1783 Options Database Keys: 1784 . -snes_lag_jacobian <lag> 1785 1786 Notes: 1787 The default is 1 1788 The jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 1789 1790 Level: intermediate 1791 1792 .keywords: SNES, nonlinear, set, convergence, tolerances 1793 1794 .seealso: SNESSetTrustRegionTolerance(), SNESSetLagJacobian(), SNESSetLagPreconditioner(), SNESGetLagPreconditioner() 1795 1796 @*/ 1797 PetscErrorCode SNESGetLagJacobian(SNES snes,PetscInt *lag) 1798 { 1799 PetscFunctionBegin; 1800 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1801 *lag = snes->lagjacobian; 1802 PetscFunctionReturn(0); 1803 } 1804 1805 #undef __FUNCT__ 1806 #define __FUNCT__ "SNESSetTolerances" 1807 /*@ 1808 SNESSetTolerances - Sets various parameters used in convergence tests. 1809 1810 Logically Collective on SNES 1811 1812 Input Parameters: 1813 + snes - the SNES context 1814 . abstol - absolute convergence tolerance 1815 . rtol - relative convergence tolerance 1816 . stol - convergence tolerance in terms of the norm 1817 of the change in the solution between steps 1818 . maxit - maximum number of iterations 1819 - maxf - maximum number of function evaluations 1820 1821 Options Database Keys: 1822 + -snes_atol <abstol> - Sets abstol 1823 . -snes_rtol <rtol> - Sets rtol 1824 . -snes_stol <stol> - Sets stol 1825 . -snes_max_it <maxit> - Sets maxit 1826 - -snes_max_funcs <maxf> - Sets maxf 1827 1828 Notes: 1829 The default maximum number of iterations is 50. 1830 The default maximum number of function evaluations is 1000. 1831 1832 Level: intermediate 1833 1834 .keywords: SNES, nonlinear, set, convergence, tolerances 1835 1836 .seealso: SNESSetTrustRegionTolerance() 1837 @*/ 1838 PetscErrorCode SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf) 1839 { 1840 PetscFunctionBegin; 1841 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1842 PetscValidLogicalCollectiveReal(snes,abstol,2); 1843 PetscValidLogicalCollectiveReal(snes,rtol,3); 1844 PetscValidLogicalCollectiveReal(snes,stol,4); 1845 PetscValidLogicalCollectiveInt(snes,maxit,5); 1846 PetscValidLogicalCollectiveInt(snes,maxf,6); 1847 1848 if (abstol != PETSC_DEFAULT) snes->abstol = abstol; 1849 if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1850 if (stol != PETSC_DEFAULT) snes->xtol = stol; 1851 if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1852 if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 1853 PetscFunctionReturn(0); 1854 } 1855 1856 #undef __FUNCT__ 1857 #define __FUNCT__ "SNESGetTolerances" 1858 /*@ 1859 SNESGetTolerances - Gets various parameters used in convergence tests. 1860 1861 Not Collective 1862 1863 Input Parameters: 1864 + snes - the SNES context 1865 . atol - absolute convergence tolerance 1866 . rtol - relative convergence tolerance 1867 . stol - convergence tolerance in terms of the norm 1868 of the change in the solution between steps 1869 . maxit - maximum number of iterations 1870 - maxf - maximum number of function evaluations 1871 1872 Notes: 1873 The user can specify PETSC_NULL for any parameter that is not needed. 1874 1875 Level: intermediate 1876 1877 .keywords: SNES, nonlinear, get, convergence, tolerances 1878 1879 .seealso: SNESSetTolerances() 1880 @*/ 1881 PetscErrorCode SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf) 1882 { 1883 PetscFunctionBegin; 1884 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1885 if (atol) *atol = snes->abstol; 1886 if (rtol) *rtol = snes->rtol; 1887 if (stol) *stol = snes->xtol; 1888 if (maxit) *maxit = snes->max_its; 1889 if (maxf) *maxf = snes->max_funcs; 1890 PetscFunctionReturn(0); 1891 } 1892 1893 #undef __FUNCT__ 1894 #define __FUNCT__ "SNESSetTrustRegionTolerance" 1895 /*@ 1896 SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 1897 1898 Logically Collective on SNES 1899 1900 Input Parameters: 1901 + snes - the SNES context 1902 - tol - tolerance 1903 1904 Options Database Key: 1905 . -snes_trtol <tol> - Sets tol 1906 1907 Level: intermediate 1908 1909 .keywords: SNES, nonlinear, set, trust region, tolerance 1910 1911 .seealso: SNESSetTolerances() 1912 @*/ 1913 PetscErrorCode SNESSetTrustRegionTolerance(SNES snes,PetscReal tol) 1914 { 1915 PetscFunctionBegin; 1916 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1917 PetscValidLogicalCollectiveReal(snes,tol,2); 1918 snes->deltatol = tol; 1919 PetscFunctionReturn(0); 1920 } 1921 1922 /* 1923 Duplicate the lg monitors for SNES from KSP; for some reason with 1924 dynamic libraries things don't work under Sun4 if we just use 1925 macros instead of functions 1926 */ 1927 #undef __FUNCT__ 1928 #define __FUNCT__ "SNESMonitorLG" 1929 PetscErrorCode SNESMonitorLG(SNES snes,PetscInt it,PetscReal norm,void *ctx) 1930 { 1931 PetscErrorCode ierr; 1932 1933 PetscFunctionBegin; 1934 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1935 ierr = KSPMonitorLG((KSP)snes,it,norm,ctx);CHKERRQ(ierr); 1936 PetscFunctionReturn(0); 1937 } 1938 1939 #undef __FUNCT__ 1940 #define __FUNCT__ "SNESMonitorLGCreate" 1941 PetscErrorCode SNESMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1942 { 1943 PetscErrorCode ierr; 1944 1945 PetscFunctionBegin; 1946 ierr = KSPMonitorLGCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr); 1947 PetscFunctionReturn(0); 1948 } 1949 1950 #undef __FUNCT__ 1951 #define __FUNCT__ "SNESMonitorLGDestroy" 1952 PetscErrorCode SNESMonitorLGDestroy(PetscDrawLG *draw) 1953 { 1954 PetscErrorCode ierr; 1955 1956 PetscFunctionBegin; 1957 ierr = KSPMonitorLGDestroy(draw);CHKERRQ(ierr); 1958 PetscFunctionReturn(0); 1959 } 1960 1961 extern PetscErrorCode SNESMonitorRange_Private(SNES,PetscInt,PetscReal*); 1962 #undef __FUNCT__ 1963 #define __FUNCT__ "SNESMonitorLGRange" 1964 PetscErrorCode SNESMonitorLGRange(SNES snes,PetscInt n,PetscReal rnorm,void *monctx) 1965 { 1966 PetscDrawLG lg; 1967 PetscErrorCode ierr; 1968 PetscReal x,y,per; 1969 PetscViewer v = (PetscViewer)monctx; 1970 static PetscReal prev; /* should be in the context */ 1971 PetscDraw draw; 1972 PetscFunctionBegin; 1973 if (!monctx) { 1974 MPI_Comm comm; 1975 1976 ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 1977 v = PETSC_VIEWER_DRAW_(comm); 1978 } 1979 ierr = PetscViewerDrawGetDrawLG(v,0,&lg);CHKERRQ(ierr); 1980 if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 1981 ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr); 1982 ierr = PetscDrawSetTitle(draw,"Residual norm");CHKERRQ(ierr); 1983 x = (PetscReal) n; 1984 if (rnorm > 0.0) y = log10(rnorm); else y = -15.0; 1985 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 1986 if (n < 20 || !(n % 5)) { 1987 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 1988 } 1989 1990 ierr = PetscViewerDrawGetDrawLG(v,1,&lg);CHKERRQ(ierr); 1991 if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 1992 ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr); 1993 ierr = PetscDrawSetTitle(draw,"% elemts > .2*max elemt");CHKERRQ(ierr); 1994 ierr = SNESMonitorRange_Private(snes,n,&per);CHKERRQ(ierr); 1995 x = (PetscReal) n; 1996 y = 100.0*per; 1997 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 1998 if (n < 20 || !(n % 5)) { 1999 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 2000 } 2001 2002 ierr = PetscViewerDrawGetDrawLG(v,2,&lg);CHKERRQ(ierr); 2003 if (!n) {prev = rnorm;ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 2004 ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr); 2005 ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm");CHKERRQ(ierr); 2006 x = (PetscReal) n; 2007 y = (prev - rnorm)/prev; 2008 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 2009 if (n < 20 || !(n % 5)) { 2010 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 2011 } 2012 2013 ierr = PetscViewerDrawGetDrawLG(v,3,&lg);CHKERRQ(ierr); 2014 if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 2015 ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr); 2016 ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");CHKERRQ(ierr); 2017 x = (PetscReal) n; 2018 y = (prev - rnorm)/(prev*per); 2019 if (n > 2) { /*skip initial crazy value */ 2020 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 2021 } 2022 if (n < 20 || !(n % 5)) { 2023 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 2024 } 2025 prev = rnorm; 2026 PetscFunctionReturn(0); 2027 } 2028 2029 #undef __FUNCT__ 2030 #define __FUNCT__ "SNESMonitorLGRangeCreate" 2031 PetscErrorCode SNESMonitorLGRangeCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 2032 { 2033 PetscErrorCode ierr; 2034 2035 PetscFunctionBegin; 2036 ierr = KSPMonitorLGCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr); 2037 PetscFunctionReturn(0); 2038 } 2039 2040 #undef __FUNCT__ 2041 #define __FUNCT__ "SNESMonitorLGRangeDestroy" 2042 PetscErrorCode SNESMonitorLGRangeDestroy(PetscDrawLG *draw) 2043 { 2044 PetscErrorCode ierr; 2045 2046 PetscFunctionBegin; 2047 ierr = KSPMonitorLGDestroy(draw);CHKERRQ(ierr); 2048 PetscFunctionReturn(0); 2049 } 2050 2051 #undef __FUNCT__ 2052 #define __FUNCT__ "SNESMonitor" 2053 /* 2054 Runs the user provided monitor routines, if they exists. 2055 */ 2056 PetscErrorCode SNESMonitor(SNES snes,PetscInt iter,PetscReal rnorm) 2057 { 2058 PetscErrorCode ierr; 2059 PetscInt i,n = snes->numbermonitors; 2060 2061 PetscFunctionBegin; 2062 for (i=0; i<n; i++) { 2063 ierr = (*snes->monitor[i])(snes,iter,rnorm,snes->monitorcontext[i]);CHKERRQ(ierr); 2064 } 2065 PetscFunctionReturn(0); 2066 } 2067 2068 /* ------------ Routines to set performance monitoring options ----------- */ 2069 2070 #undef __FUNCT__ 2071 #define __FUNCT__ "SNESMonitorSet" 2072 /*@C 2073 SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every 2074 iteration of the nonlinear solver to display the iteration's 2075 progress. 2076 2077 Logically Collective on SNES 2078 2079 Input Parameters: 2080 + snes - the SNES context 2081 . func - monitoring routine 2082 . mctx - [optional] user-defined context for private data for the 2083 monitor routine (use PETSC_NULL if no context is desired) 2084 - monitordestroy - [optional] routine that frees monitor context 2085 (may be PETSC_NULL) 2086 2087 Calling sequence of func: 2088 $ int func(SNES snes,PetscInt its, PetscReal norm,void *mctx) 2089 2090 + snes - the SNES context 2091 . its - iteration number 2092 . norm - 2-norm function value (may be estimated) 2093 - mctx - [optional] monitoring context 2094 2095 Options Database Keys: 2096 + -snes_monitor - sets SNESMonitorDefault() 2097 . -snes_monitor_draw - sets line graph monitor, 2098 uses SNESMonitorLGCreate() 2099 - -snes_monitor_cancel - cancels all monitors that have 2100 been hardwired into a code by 2101 calls to SNESMonitorSet(), but 2102 does not cancel those set via 2103 the options database. 2104 2105 Notes: 2106 Several different monitoring routines may be set by calling 2107 SNESMonitorSet() multiple times; all will be called in the 2108 order in which they were set. 2109 2110 Fortran notes: Only a single monitor function can be set for each SNES object 2111 2112 Level: intermediate 2113 2114 .keywords: SNES, nonlinear, set, monitor 2115 2116 .seealso: SNESMonitorDefault(), SNESMonitorCancel() 2117 @*/ 2118 PetscErrorCode SNESMonitorSet(SNES snes,PetscErrorCode (*monitor)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**)) 2119 { 2120 PetscInt i; 2121 2122 PetscFunctionBegin; 2123 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2124 if (snes->numbermonitors >= MAXSNESMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 2125 for (i=0; i<snes->numbermonitors;i++) { 2126 if (monitor == snes->monitor[i] && monitordestroy == snes->monitordestroy[i] && mctx == snes->monitorcontext[i]) PetscFunctionReturn(0); 2127 2128 /* check if both default monitors that share common ASCII viewer */ 2129 if (monitor == snes->monitor[i] && monitor == SNESMonitorDefault) { 2130 if (mctx && snes->monitorcontext[i]) { 2131 PetscErrorCode ierr; 2132 PetscViewerASCIIMonitor viewer1 = (PetscViewerASCIIMonitor) mctx; 2133 PetscViewerASCIIMonitor viewer2 = (PetscViewerASCIIMonitor) snes->monitorcontext[i]; 2134 if (viewer1->viewer == viewer2->viewer) { 2135 ierr = (*monitordestroy)(&mctx);CHKERRQ(ierr); 2136 PetscFunctionReturn(0); 2137 } 2138 } 2139 } 2140 } 2141 snes->monitor[snes->numbermonitors] = monitor; 2142 snes->monitordestroy[snes->numbermonitors] = monitordestroy; 2143 snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 2144 PetscFunctionReturn(0); 2145 } 2146 2147 #undef __FUNCT__ 2148 #define __FUNCT__ "SNESMonitorCancel" 2149 /*@C 2150 SNESMonitorCancel - Clears all the monitor functions for a SNES object. 2151 2152 Logically Collective on SNES 2153 2154 Input Parameters: 2155 . snes - the SNES context 2156 2157 Options Database Key: 2158 . -snes_monitor_cancel - cancels all monitors that have been hardwired 2159 into a code by calls to SNESMonitorSet(), but does not cancel those 2160 set via the options database 2161 2162 Notes: 2163 There is no way to clear one specific monitor from a SNES object. 2164 2165 Level: intermediate 2166 2167 .keywords: SNES, nonlinear, set, monitor 2168 2169 .seealso: SNESMonitorDefault(), SNESMonitorSet() 2170 @*/ 2171 PetscErrorCode SNESMonitorCancel(SNES snes) 2172 { 2173 PetscErrorCode ierr; 2174 PetscInt i; 2175 2176 PetscFunctionBegin; 2177 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2178 for (i=0; i<snes->numbermonitors; i++) { 2179 if (snes->monitordestroy[i]) { 2180 ierr = (*snes->monitordestroy[i])(&snes->monitorcontext[i]);CHKERRQ(ierr); 2181 } 2182 } 2183 snes->numbermonitors = 0; 2184 PetscFunctionReturn(0); 2185 } 2186 2187 #undef __FUNCT__ 2188 #define __FUNCT__ "SNESSetConvergenceTest" 2189 /*@C 2190 SNESSetConvergenceTest - Sets the function that is to be used 2191 to test for convergence of the nonlinear iterative solution. 2192 2193 Logically Collective on SNES 2194 2195 Input Parameters: 2196 + snes - the SNES context 2197 . func - routine to test for convergence 2198 . cctx - [optional] context for private data for the convergence routine (may be PETSC_NULL) 2199 - destroy - [optional] destructor for the context (may be PETSC_NULL; PETSC_NULL_FUNCTION in Fortran) 2200 2201 Calling sequence of func: 2202 $ PetscErrorCode func (SNES snes,PetscInt it,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx) 2203 2204 + snes - the SNES context 2205 . it - current iteration (0 is the first and is before any Newton step) 2206 . cctx - [optional] convergence context 2207 . reason - reason for convergence/divergence 2208 . xnorm - 2-norm of current iterate 2209 . gnorm - 2-norm of current step 2210 - f - 2-norm of function 2211 2212 Level: advanced 2213 2214 .keywords: SNES, nonlinear, set, convergence, test 2215 2216 .seealso: SNESDefaultConverged(), SNESSkipConverged() 2217 @*/ 2218 PetscErrorCode SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*)) 2219 { 2220 PetscErrorCode ierr; 2221 2222 PetscFunctionBegin; 2223 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2224 if (!func) func = SNESSkipConverged; 2225 if (snes->ops->convergeddestroy) { 2226 ierr = (*snes->ops->convergeddestroy)(snes->cnvP);CHKERRQ(ierr); 2227 } 2228 snes->ops->converged = func; 2229 snes->ops->convergeddestroy = destroy; 2230 snes->cnvP = cctx; 2231 PetscFunctionReturn(0); 2232 } 2233 2234 #undef __FUNCT__ 2235 #define __FUNCT__ "SNESGetConvergedReason" 2236 /*@ 2237 SNESGetConvergedReason - Gets the reason the SNES iteration was stopped. 2238 2239 Not Collective 2240 2241 Input Parameter: 2242 . snes - the SNES context 2243 2244 Output Parameter: 2245 . reason - negative value indicates diverged, positive value converged, see SNESConvergedReason or the 2246 manual pages for the individual convergence tests for complete lists 2247 2248 Level: intermediate 2249 2250 Notes: Can only be called after the call the SNESSolve() is complete. 2251 2252 .keywords: SNES, nonlinear, set, convergence, test 2253 2254 .seealso: SNESSetConvergenceTest(), SNESConvergedReason 2255 @*/ 2256 PetscErrorCode SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason) 2257 { 2258 PetscFunctionBegin; 2259 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2260 PetscValidPointer(reason,2); 2261 *reason = snes->reason; 2262 PetscFunctionReturn(0); 2263 } 2264 2265 #undef __FUNCT__ 2266 #define __FUNCT__ "SNESSetConvergenceHistory" 2267 /*@ 2268 SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 2269 2270 Logically Collective on SNES 2271 2272 Input Parameters: 2273 + snes - iterative context obtained from SNESCreate() 2274 . a - array to hold history, this array will contain the function norms computed at each step 2275 . its - integer array holds the number of linear iterations for each solve. 2276 . na - size of a and its 2277 - reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero, 2278 else it continues storing new values for new nonlinear solves after the old ones 2279 2280 Notes: 2281 If 'a' and 'its' are PETSC_NULL then space is allocated for the history. If 'na' PETSC_DECIDE or PETSC_DEFAULT then a 2282 default array of length 10000 is allocated. 2283 2284 This routine is useful, e.g., when running a code for purposes 2285 of accurate performance monitoring, when no I/O should be done 2286 during the section of code that is being timed. 2287 2288 Level: intermediate 2289 2290 .keywords: SNES, set, convergence, history 2291 2292 .seealso: SNESGetConvergenceHistory() 2293 2294 @*/ 2295 PetscErrorCode SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscBool reset) 2296 { 2297 PetscErrorCode ierr; 2298 2299 PetscFunctionBegin; 2300 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2301 if (na) PetscValidScalarPointer(a,2); 2302 if (its) PetscValidIntPointer(its,3); 2303 if (na == PETSC_DECIDE || na == PETSC_DEFAULT || !a) { 2304 if (na == PETSC_DECIDE || na == PETSC_DEFAULT) na = 1000; 2305 ierr = PetscMalloc(na*sizeof(PetscReal),&a);CHKERRQ(ierr); 2306 ierr = PetscMalloc(na*sizeof(PetscInt),&its);CHKERRQ(ierr); 2307 snes->conv_malloc = PETSC_TRUE; 2308 } 2309 snes->conv_hist = a; 2310 snes->conv_hist_its = its; 2311 snes->conv_hist_max = na; 2312 snes->conv_hist_len = 0; 2313 snes->conv_hist_reset = reset; 2314 PetscFunctionReturn(0); 2315 } 2316 2317 #if defined(PETSC_HAVE_MATLAB_ENGINE) 2318 #include <engine.h> /* MATLAB include file */ 2319 #include <mex.h> /* MATLAB include file */ 2320 EXTERN_C_BEGIN 2321 #undef __FUNCT__ 2322 #define __FUNCT__ "SNESGetConvergenceHistoryMatlab" 2323 mxArray *SNESGetConvergenceHistoryMatlab(SNES snes) 2324 { 2325 mxArray *mat; 2326 PetscInt i; 2327 PetscReal *ar; 2328 2329 PetscFunctionBegin; 2330 mat = mxCreateDoubleMatrix(snes->conv_hist_len,1,mxREAL); 2331 ar = (PetscReal*) mxGetData(mat); 2332 for (i=0; i<snes->conv_hist_len; i++) { 2333 ar[i] = snes->conv_hist[i]; 2334 } 2335 PetscFunctionReturn(mat); 2336 } 2337 EXTERN_C_END 2338 #endif 2339 2340 2341 #undef __FUNCT__ 2342 #define __FUNCT__ "SNESGetConvergenceHistory" 2343 /*@C 2344 SNESGetConvergenceHistory - Gets the array used to hold the convergence history. 2345 2346 Not Collective 2347 2348 Input Parameter: 2349 . snes - iterative context obtained from SNESCreate() 2350 2351 Output Parameters: 2352 . a - array to hold history 2353 . its - integer array holds the number of linear iterations (or 2354 negative if not converged) for each solve. 2355 - na - size of a and its 2356 2357 Notes: 2358 The calling sequence for this routine in Fortran is 2359 $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr) 2360 2361 This routine is useful, e.g., when running a code for purposes 2362 of accurate performance monitoring, when no I/O should be done 2363 during the section of code that is being timed. 2364 2365 Level: intermediate 2366 2367 .keywords: SNES, get, convergence, history 2368 2369 .seealso: SNESSetConvergencHistory() 2370 2371 @*/ 2372 PetscErrorCode SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na) 2373 { 2374 PetscFunctionBegin; 2375 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2376 if (a) *a = snes->conv_hist; 2377 if (its) *its = snes->conv_hist_its; 2378 if (na) *na = snes->conv_hist_len; 2379 PetscFunctionReturn(0); 2380 } 2381 2382 #undef __FUNCT__ 2383 #define __FUNCT__ "SNESSetUpdate" 2384 /*@C 2385 SNESSetUpdate - Sets the general-purpose update function called 2386 at the beginning of every iteration of the nonlinear solve. Specifically 2387 it is called just before the Jacobian is "evaluated". 2388 2389 Logically Collective on SNES 2390 2391 Input Parameters: 2392 . snes - The nonlinear solver context 2393 . func - The function 2394 2395 Calling sequence of func: 2396 . func (SNES snes, PetscInt step); 2397 2398 . step - The current step of the iteration 2399 2400 Level: advanced 2401 2402 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() 2403 This is not used by most users. 2404 2405 .keywords: SNES, update 2406 2407 .seealso SNESDefaultUpdate(), SNESSetJacobian(), SNESSolve() 2408 @*/ 2409 PetscErrorCode SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt)) 2410 { 2411 PetscFunctionBegin; 2412 PetscValidHeaderSpecific(snes, SNES_CLASSID,1); 2413 snes->ops->update = func; 2414 PetscFunctionReturn(0); 2415 } 2416 2417 #undef __FUNCT__ 2418 #define __FUNCT__ "SNESDefaultUpdate" 2419 /*@ 2420 SNESDefaultUpdate - The default update function which does nothing. 2421 2422 Not collective 2423 2424 Input Parameters: 2425 . snes - The nonlinear solver context 2426 . step - The current step of the iteration 2427 2428 Level: intermediate 2429 2430 .keywords: SNES, update 2431 .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultShortolutionBC() 2432 @*/ 2433 PetscErrorCode SNESDefaultUpdate(SNES snes, PetscInt step) 2434 { 2435 PetscFunctionBegin; 2436 PetscFunctionReturn(0); 2437 } 2438 2439 #undef __FUNCT__ 2440 #define __FUNCT__ "SNESScaleStep_Private" 2441 /* 2442 SNESScaleStep_Private - Scales a step so that its length is less than the 2443 positive parameter delta. 2444 2445 Input Parameters: 2446 + snes - the SNES context 2447 . y - approximate solution of linear system 2448 . fnorm - 2-norm of current function 2449 - delta - trust region size 2450 2451 Output Parameters: 2452 + gpnorm - predicted function norm at the new point, assuming local 2453 linearization. The value is zero if the step lies within the trust 2454 region, and exceeds zero otherwise. 2455 - ynorm - 2-norm of the step 2456 2457 Note: 2458 For non-trust region methods such as SNESLS, the parameter delta 2459 is set to be the maximum allowable step size. 2460 2461 .keywords: SNES, nonlinear, scale, step 2462 */ 2463 PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm) 2464 { 2465 PetscReal nrm; 2466 PetscScalar cnorm; 2467 PetscErrorCode ierr; 2468 2469 PetscFunctionBegin; 2470 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2471 PetscValidHeaderSpecific(y,VEC_CLASSID,2); 2472 PetscCheckSameComm(snes,1,y,2); 2473 2474 ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr); 2475 if (nrm > *delta) { 2476 nrm = *delta/nrm; 2477 *gpnorm = (1.0 - nrm)*(*fnorm); 2478 cnorm = nrm; 2479 ierr = VecScale(y,cnorm);CHKERRQ(ierr); 2480 *ynorm = *delta; 2481 } else { 2482 *gpnorm = 0.0; 2483 *ynorm = nrm; 2484 } 2485 PetscFunctionReturn(0); 2486 } 2487 2488 #undef __FUNCT__ 2489 #define __FUNCT__ "SNESSolve" 2490 /*@C 2491 SNESSolve - Solves a nonlinear system F(x) = b. 2492 Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX(). 2493 2494 Collective on SNES 2495 2496 Input Parameters: 2497 + snes - the SNES context 2498 . b - the constant part of the equation, or PETSC_NULL to use zero. 2499 - x - the solution vector. 2500 2501 Notes: 2502 The user should initialize the vector,x, with the initial guess 2503 for the nonlinear solve prior to calling SNESSolve. In particular, 2504 to employ an initial guess of zero, the user should explicitly set 2505 this vector to zero by calling VecSet(). 2506 2507 Level: beginner 2508 2509 .keywords: SNES, nonlinear, solve 2510 2511 .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian() 2512 @*/ 2513 PetscErrorCode SNESSolve(SNES snes,Vec b,Vec x) 2514 { 2515 PetscErrorCode ierr; 2516 PetscBool flg; 2517 char filename[PETSC_MAX_PATH_LEN]; 2518 PetscViewer viewer; 2519 PetscInt grid; 2520 2521 PetscFunctionBegin; 2522 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2523 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 2524 PetscCheckSameComm(snes,1,x,3); 2525 if (b) PetscValidHeaderSpecific(b,VEC_CLASSID,2); 2526 if (b) PetscCheckSameComm(snes,1,b,2); 2527 2528 for (grid=0; grid<snes->gridsequence+1; grid++) { 2529 2530 /* set solution vector */ 2531 if (!grid) {ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);} 2532 ierr = VecDestroy(&snes->vec_sol);CHKERRQ(ierr); 2533 snes->vec_sol = x; 2534 /* set afine vector if provided */ 2535 if (b) { ierr = PetscObjectReference((PetscObject)b);CHKERRQ(ierr); } 2536 ierr = VecDestroy(&snes->vec_rhs);CHKERRQ(ierr); 2537 snes->vec_rhs = b; 2538 2539 ierr = SNESSetUp(snes);CHKERRQ(ierr); 2540 2541 if (!grid && snes->ops->computeinitialguess) { 2542 ierr = (*snes->ops->computeinitialguess)(snes,snes->vec_sol,snes->initialguessP);CHKERRQ(ierr); 2543 } 2544 2545 if (snes->conv_hist_reset) snes->conv_hist_len = 0; 2546 snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0; 2547 2548 ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 2549 ierr = (*snes->ops->solve)(snes);CHKERRQ(ierr); 2550 ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 2551 if (snes->domainerror){ 2552 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 2553 snes->domainerror = PETSC_FALSE; 2554 } 2555 if (!snes->reason) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason"); 2556 2557 ierr = PetscOptionsGetString(((PetscObject)snes)->prefix,"-snes_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 2558 if (flg && !PetscPreLoadingOn) { 2559 ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,filename,&viewer);CHKERRQ(ierr); 2560 ierr = SNESView(snes,viewer);CHKERRQ(ierr); 2561 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 2562 } 2563 2564 flg = PETSC_FALSE; 2565 ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_test_local_min",&flg,PETSC_NULL);CHKERRQ(ierr); 2566 if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); } 2567 if (snes->printreason) { 2568 if (snes->reason > 0) { 2569 ierr = PetscPrintf(((PetscObject)snes)->comm,"Nonlinear solve converged due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr); 2570 } else { 2571 ierr = PetscPrintf(((PetscObject)snes)->comm,"Nonlinear solve did not converge due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr); 2572 } 2573 } 2574 2575 if (snes->errorifnotconverged && snes->reason < 0) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_NOT_CONVERGED,"SNESSolve has not converged"); 2576 if (grid < snes->gridsequence) { 2577 DM fine; 2578 Vec xnew; 2579 Mat interp; 2580 2581 ierr = DMRefine(snes->dm,((PetscObject)snes)->comm,&fine);CHKERRQ(ierr); 2582 ierr = DMGetInterpolation(snes->dm,fine,&interp,PETSC_NULL);CHKERRQ(ierr); 2583 ierr = DMCreateGlobalVector(fine,&xnew);CHKERRQ(ierr); 2584 ierr = MatInterpolate(interp,x,xnew);CHKERRQ(ierr); 2585 ierr = MatDestroy(&interp);CHKERRQ(ierr); 2586 x = xnew; 2587 2588 ierr = SNESReset(snes);CHKERRQ(ierr); 2589 ierr = SNESSetDM(snes,fine);CHKERRQ(ierr); 2590 ierr = DMDestroy(&fine);CHKERRQ(ierr); 2591 } 2592 } 2593 PetscFunctionReturn(0); 2594 } 2595 2596 /* --------- Internal routines for SNES Package --------- */ 2597 2598 #undef __FUNCT__ 2599 #define __FUNCT__ "SNESSetType" 2600 /*@C 2601 SNESSetType - Sets the method for the nonlinear solver. 2602 2603 Collective on SNES 2604 2605 Input Parameters: 2606 + snes - the SNES context 2607 - type - a known method 2608 2609 Options Database Key: 2610 . -snes_type <type> - Sets the method; use -help for a list 2611 of available methods (for instance, ls or tr) 2612 2613 Notes: 2614 See "petsc/include/petscsnes.h" for available methods (for instance) 2615 + SNESLS - Newton's method with line search 2616 (systems of nonlinear equations) 2617 . SNESTR - Newton's method with trust region 2618 (systems of nonlinear equations) 2619 2620 Normally, it is best to use the SNESSetFromOptions() command and then 2621 set the SNES solver type from the options database rather than by using 2622 this routine. Using the options database provides the user with 2623 maximum flexibility in evaluating the many nonlinear solvers. 2624 The SNESSetType() routine is provided for those situations where it 2625 is necessary to set the nonlinear solver independently of the command 2626 line or options database. This might be the case, for example, when 2627 the choice of solver changes during the execution of the program, 2628 and the user's application is taking responsibility for choosing the 2629 appropriate method. 2630 2631 Level: intermediate 2632 2633 .keywords: SNES, set, type 2634 2635 .seealso: SNESType, SNESCreate() 2636 2637 @*/ 2638 PetscErrorCode SNESSetType(SNES snes,const SNESType type) 2639 { 2640 PetscErrorCode ierr,(*r)(SNES); 2641 PetscBool match; 2642 2643 PetscFunctionBegin; 2644 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2645 PetscValidCharPointer(type,2); 2646 2647 ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr); 2648 if (match) PetscFunctionReturn(0); 2649 2650 ierr = PetscFListFind(SNESList,((PetscObject)snes)->comm,type,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr); 2651 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type); 2652 /* Destroy the previous private SNES context */ 2653 if (snes->ops->destroy) { ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr); } 2654 /* Reinitialize function pointers in SNESOps structure */ 2655 snes->ops->setup = 0; 2656 snes->ops->solve = 0; 2657 snes->ops->view = 0; 2658 snes->ops->setfromoptions = 0; 2659 snes->ops->destroy = 0; 2660 /* Call the SNESCreate_XXX routine for this particular Nonlinear solver */ 2661 snes->setupcalled = PETSC_FALSE; 2662 ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr); 2663 ierr = (*r)(snes);CHKERRQ(ierr); 2664 #if defined(PETSC_HAVE_AMS) 2665 if (PetscAMSPublishAll) { 2666 ierr = PetscObjectAMSPublish((PetscObject)snes);CHKERRQ(ierr); 2667 } 2668 #endif 2669 PetscFunctionReturn(0); 2670 } 2671 2672 2673 /* --------------------------------------------------------------------- */ 2674 #undef __FUNCT__ 2675 #define __FUNCT__ "SNESRegisterDestroy" 2676 /*@ 2677 SNESRegisterDestroy - Frees the list of nonlinear solvers that were 2678 registered by SNESRegisterDynamic(). 2679 2680 Not Collective 2681 2682 Level: advanced 2683 2684 .keywords: SNES, nonlinear, register, destroy 2685 2686 .seealso: SNESRegisterAll(), SNESRegisterAll() 2687 @*/ 2688 PetscErrorCode SNESRegisterDestroy(void) 2689 { 2690 PetscErrorCode ierr; 2691 2692 PetscFunctionBegin; 2693 ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr); 2694 SNESRegisterAllCalled = PETSC_FALSE; 2695 PetscFunctionReturn(0); 2696 } 2697 2698 #undef __FUNCT__ 2699 #define __FUNCT__ "SNESGetType" 2700 /*@C 2701 SNESGetType - Gets the SNES method type and name (as a string). 2702 2703 Not Collective 2704 2705 Input Parameter: 2706 . snes - nonlinear solver context 2707 2708 Output Parameter: 2709 . type - SNES method (a character string) 2710 2711 Level: intermediate 2712 2713 .keywords: SNES, nonlinear, get, type, name 2714 @*/ 2715 PetscErrorCode SNESGetType(SNES snes,const SNESType *type) 2716 { 2717 PetscFunctionBegin; 2718 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2719 PetscValidPointer(type,2); 2720 *type = ((PetscObject)snes)->type_name; 2721 PetscFunctionReturn(0); 2722 } 2723 2724 #undef __FUNCT__ 2725 #define __FUNCT__ "SNESGetSolution" 2726 /*@ 2727 SNESGetSolution - Returns the vector where the approximate solution is 2728 stored. 2729 2730 Not Collective, but Vec is parallel if SNES is parallel 2731 2732 Input Parameter: 2733 . snes - the SNES context 2734 2735 Output Parameter: 2736 . x - the solution 2737 2738 Level: intermediate 2739 2740 .keywords: SNES, nonlinear, get, solution 2741 2742 .seealso: SNESGetSolutionUpdate(), SNESGetFunction() 2743 @*/ 2744 PetscErrorCode SNESGetSolution(SNES snes,Vec *x) 2745 { 2746 PetscFunctionBegin; 2747 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2748 PetscValidPointer(x,2); 2749 *x = snes->vec_sol; 2750 PetscFunctionReturn(0); 2751 } 2752 2753 #undef __FUNCT__ 2754 #define __FUNCT__ "SNESGetSolutionUpdate" 2755 /*@ 2756 SNESGetSolutionUpdate - Returns the vector where the solution update is 2757 stored. 2758 2759 Not Collective, but Vec is parallel if SNES is parallel 2760 2761 Input Parameter: 2762 . snes - the SNES context 2763 2764 Output Parameter: 2765 . x - the solution update 2766 2767 Level: advanced 2768 2769 .keywords: SNES, nonlinear, get, solution, update 2770 2771 .seealso: SNESGetSolution(), SNESGetFunction() 2772 @*/ 2773 PetscErrorCode SNESGetSolutionUpdate(SNES snes,Vec *x) 2774 { 2775 PetscFunctionBegin; 2776 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2777 PetscValidPointer(x,2); 2778 *x = snes->vec_sol_update; 2779 PetscFunctionReturn(0); 2780 } 2781 2782 #undef __FUNCT__ 2783 #define __FUNCT__ "SNESGetFunction" 2784 /*@C 2785 SNESGetFunction - Returns the vector where the function is stored. 2786 2787 Not Collective, but Vec is parallel if SNES is parallel 2788 2789 Input Parameter: 2790 . snes - the SNES context 2791 2792 Output Parameter: 2793 + r - the function (or PETSC_NULL) 2794 . func - the function (or PETSC_NULL) 2795 - ctx - the function context (or PETSC_NULL) 2796 2797 Level: advanced 2798 2799 .keywords: SNES, nonlinear, get, function 2800 2801 .seealso: SNESSetFunction(), SNESGetSolution() 2802 @*/ 2803 PetscErrorCode SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx) 2804 { 2805 PetscFunctionBegin; 2806 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2807 if (r) *r = snes->vec_func; 2808 if (func) *func = snes->ops->computefunction; 2809 if (ctx) *ctx = snes->funP; 2810 PetscFunctionReturn(0); 2811 } 2812 2813 #undef __FUNCT__ 2814 #define __FUNCT__ "SNESSetOptionsPrefix" 2815 /*@C 2816 SNESSetOptionsPrefix - Sets the prefix used for searching for all 2817 SNES options in the database. 2818 2819 Logically Collective on SNES 2820 2821 Input Parameter: 2822 + snes - the SNES context 2823 - prefix - the prefix to prepend to all option names 2824 2825 Notes: 2826 A hyphen (-) must NOT be given at the beginning of the prefix name. 2827 The first character of all runtime options is AUTOMATICALLY the hyphen. 2828 2829 Level: advanced 2830 2831 .keywords: SNES, set, options, prefix, database 2832 2833 .seealso: SNESSetFromOptions() 2834 @*/ 2835 PetscErrorCode SNESSetOptionsPrefix(SNES snes,const char prefix[]) 2836 { 2837 PetscErrorCode ierr; 2838 2839 PetscFunctionBegin; 2840 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2841 ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2842 if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);} 2843 ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 2844 PetscFunctionReturn(0); 2845 } 2846 2847 #undef __FUNCT__ 2848 #define __FUNCT__ "SNESAppendOptionsPrefix" 2849 /*@C 2850 SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 2851 SNES options in the database. 2852 2853 Logically Collective on SNES 2854 2855 Input Parameters: 2856 + snes - the SNES context 2857 - prefix - the prefix to prepend to all option names 2858 2859 Notes: 2860 A hyphen (-) must NOT be given at the beginning of the prefix name. 2861 The first character of all runtime options is AUTOMATICALLY the hyphen. 2862 2863 Level: advanced 2864 2865 .keywords: SNES, append, options, prefix, database 2866 2867 .seealso: SNESGetOptionsPrefix() 2868 @*/ 2869 PetscErrorCode SNESAppendOptionsPrefix(SNES snes,const char prefix[]) 2870 { 2871 PetscErrorCode ierr; 2872 2873 PetscFunctionBegin; 2874 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2875 ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2876 if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);} 2877 ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 2878 PetscFunctionReturn(0); 2879 } 2880 2881 #undef __FUNCT__ 2882 #define __FUNCT__ "SNESGetOptionsPrefix" 2883 /*@C 2884 SNESGetOptionsPrefix - Sets the prefix used for searching for all 2885 SNES options in the database. 2886 2887 Not Collective 2888 2889 Input Parameter: 2890 . snes - the SNES context 2891 2892 Output Parameter: 2893 . prefix - pointer to the prefix string used 2894 2895 Notes: On the fortran side, the user should pass in a string 'prefix' of 2896 sufficient length to hold the prefix. 2897 2898 Level: advanced 2899 2900 .keywords: SNES, get, options, prefix, database 2901 2902 .seealso: SNESAppendOptionsPrefix() 2903 @*/ 2904 PetscErrorCode SNESGetOptionsPrefix(SNES snes,const char *prefix[]) 2905 { 2906 PetscErrorCode ierr; 2907 2908 PetscFunctionBegin; 2909 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2910 ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2911 PetscFunctionReturn(0); 2912 } 2913 2914 2915 #undef __FUNCT__ 2916 #define __FUNCT__ "SNESRegister" 2917 /*@C 2918 SNESRegister - See SNESRegisterDynamic() 2919 2920 Level: advanced 2921 @*/ 2922 PetscErrorCode SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES)) 2923 { 2924 char fullname[PETSC_MAX_PATH_LEN]; 2925 PetscErrorCode ierr; 2926 2927 PetscFunctionBegin; 2928 ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 2929 ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 2930 PetscFunctionReturn(0); 2931 } 2932 2933 #undef __FUNCT__ 2934 #define __FUNCT__ "SNESTestLocalMin" 2935 PetscErrorCode SNESTestLocalMin(SNES snes) 2936 { 2937 PetscErrorCode ierr; 2938 PetscInt N,i,j; 2939 Vec u,uh,fh; 2940 PetscScalar value; 2941 PetscReal norm; 2942 2943 PetscFunctionBegin; 2944 ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr); 2945 ierr = VecDuplicate(u,&uh);CHKERRQ(ierr); 2946 ierr = VecDuplicate(u,&fh);CHKERRQ(ierr); 2947 2948 /* currently only works for sequential */ 2949 ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n"); 2950 ierr = VecGetSize(u,&N);CHKERRQ(ierr); 2951 for (i=0; i<N; i++) { 2952 ierr = VecCopy(u,uh);CHKERRQ(ierr); 2953 ierr = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr); 2954 for (j=-10; j<11; j++) { 2955 value = PetscSign(j)*exp(PetscAbs(j)-10.0); 2956 ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2957 ierr = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr); 2958 ierr = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr); 2959 ierr = PetscPrintf(PETSC_COMM_WORLD," j norm %D %18.16e\n",j,norm);CHKERRQ(ierr); 2960 value = -value; 2961 ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2962 } 2963 } 2964 ierr = VecDestroy(&uh);CHKERRQ(ierr); 2965 ierr = VecDestroy(&fh);CHKERRQ(ierr); 2966 PetscFunctionReturn(0); 2967 } 2968 2969 #undef __FUNCT__ 2970 #define __FUNCT__ "SNESKSPSetUseEW" 2971 /*@ 2972 SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for 2973 computing relative tolerance for linear solvers within an inexact 2974 Newton method. 2975 2976 Logically Collective on SNES 2977 2978 Input Parameters: 2979 + snes - SNES context 2980 - flag - PETSC_TRUE or PETSC_FALSE 2981 2982 Options Database: 2983 + -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence 2984 . -snes_ksp_ew_version ver - version of Eisenstat-Walker method 2985 . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0 2986 . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax 2987 . -snes_ksp_ew_gamma <gamma> - Sets gamma 2988 . -snes_ksp_ew_alpha <alpha> - Sets alpha 2989 . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2 2990 - -snes_ksp_ew_threshold <threshold> - Sets threshold 2991 2992 Notes: 2993 Currently, the default is to use a constant relative tolerance for 2994 the inner linear solvers. Alternatively, one can use the 2995 Eisenstat-Walker method, where the relative convergence tolerance 2996 is reset at each Newton iteration according progress of the nonlinear 2997 solver. 2998 2999 Level: advanced 3000 3001 Reference: 3002 S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an 3003 inexact Newton method", SISC 17 (1), pp.16-32, 1996. 3004 3005 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton 3006 3007 .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW() 3008 @*/ 3009 PetscErrorCode SNESKSPSetUseEW(SNES snes,PetscBool flag) 3010 { 3011 PetscFunctionBegin; 3012 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3013 PetscValidLogicalCollectiveBool(snes,flag,2); 3014 snes->ksp_ewconv = flag; 3015 PetscFunctionReturn(0); 3016 } 3017 3018 #undef __FUNCT__ 3019 #define __FUNCT__ "SNESKSPGetUseEW" 3020 /*@ 3021 SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method 3022 for computing relative tolerance for linear solvers within an 3023 inexact Newton method. 3024 3025 Not Collective 3026 3027 Input Parameter: 3028 . snes - SNES context 3029 3030 Output Parameter: 3031 . flag - PETSC_TRUE or PETSC_FALSE 3032 3033 Notes: 3034 Currently, the default is to use a constant relative tolerance for 3035 the inner linear solvers. Alternatively, one can use the 3036 Eisenstat-Walker method, where the relative convergence tolerance 3037 is reset at each Newton iteration according progress of the nonlinear 3038 solver. 3039 3040 Level: advanced 3041 3042 Reference: 3043 S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an 3044 inexact Newton method", SISC 17 (1), pp.16-32, 1996. 3045 3046 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton 3047 3048 .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW() 3049 @*/ 3050 PetscErrorCode SNESKSPGetUseEW(SNES snes, PetscBool *flag) 3051 { 3052 PetscFunctionBegin; 3053 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3054 PetscValidPointer(flag,2); 3055 *flag = snes->ksp_ewconv; 3056 PetscFunctionReturn(0); 3057 } 3058 3059 #undef __FUNCT__ 3060 #define __FUNCT__ "SNESKSPSetParametersEW" 3061 /*@ 3062 SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker 3063 convergence criteria for the linear solvers within an inexact 3064 Newton method. 3065 3066 Logically Collective on SNES 3067 3068 Input Parameters: 3069 + snes - SNES context 3070 . version - version 1, 2 (default is 2) or 3 3071 . rtol_0 - initial relative tolerance (0 <= rtol_0 < 1) 3072 . rtol_max - maximum relative tolerance (0 <= rtol_max < 1) 3073 . gamma - multiplicative factor for version 2 rtol computation 3074 (0 <= gamma2 <= 1) 3075 . alpha - power for version 2 rtol computation (1 < alpha <= 2) 3076 . alpha2 - power for safeguard 3077 - threshold - threshold for imposing safeguard (0 < threshold < 1) 3078 3079 Note: 3080 Version 3 was contributed by Luis Chacon, June 2006. 3081 3082 Use PETSC_DEFAULT to retain the default for any of the parameters. 3083 3084 Level: advanced 3085 3086 Reference: 3087 S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an 3088 inexact Newton method", Utah State University Math. Stat. Dept. Res. 3089 Report 6/94/75, June, 1994, to appear in SIAM J. Sci. Comput. 3090 3091 .keywords: SNES, KSP, Eisenstat, Walker, set, parameters 3092 3093 .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW() 3094 @*/ 3095 PetscErrorCode SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max, 3096 PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold) 3097 { 3098 SNESKSPEW *kctx; 3099 PetscFunctionBegin; 3100 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3101 kctx = (SNESKSPEW*)snes->kspconvctx; 3102 if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing"); 3103 PetscValidLogicalCollectiveInt(snes,version,2); 3104 PetscValidLogicalCollectiveReal(snes,rtol_0,3); 3105 PetscValidLogicalCollectiveReal(snes,rtol_max,4); 3106 PetscValidLogicalCollectiveReal(snes,gamma,5); 3107 PetscValidLogicalCollectiveReal(snes,alpha,6); 3108 PetscValidLogicalCollectiveReal(snes,alpha2,7); 3109 PetscValidLogicalCollectiveReal(snes,threshold,8); 3110 3111 if (version != PETSC_DEFAULT) kctx->version = version; 3112 if (rtol_0 != PETSC_DEFAULT) kctx->rtol_0 = rtol_0; 3113 if (rtol_max != PETSC_DEFAULT) kctx->rtol_max = rtol_max; 3114 if (gamma != PETSC_DEFAULT) kctx->gamma = gamma; 3115 if (alpha != PETSC_DEFAULT) kctx->alpha = alpha; 3116 if (alpha2 != PETSC_DEFAULT) kctx->alpha2 = alpha2; 3117 if (threshold != PETSC_DEFAULT) kctx->threshold = threshold; 3118 3119 if (kctx->version < 1 || kctx->version > 3) { 3120 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 and 3 are supported: %D",kctx->version); 3121 } 3122 if (kctx->rtol_0 < 0.0 || kctx->rtol_0 >= 1.0) { 3123 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_0 < 1.0: %G",kctx->rtol_0); 3124 } 3125 if (kctx->rtol_max < 0.0 || kctx->rtol_max >= 1.0) { 3126 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_max (%G) < 1.0\n",kctx->rtol_max); 3127 } 3128 if (kctx->gamma < 0.0 || kctx->gamma > 1.0) { 3129 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= gamma (%G) <= 1.0\n",kctx->gamma); 3130 } 3131 if (kctx->alpha <= 1.0 || kctx->alpha > 2.0) { 3132 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"1.0 < alpha (%G) <= 2.0\n",kctx->alpha); 3133 } 3134 if (kctx->threshold <= 0.0 || kctx->threshold >= 1.0) { 3135 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 < threshold (%G) < 1.0\n",kctx->threshold); 3136 } 3137 PetscFunctionReturn(0); 3138 } 3139 3140 #undef __FUNCT__ 3141 #define __FUNCT__ "SNESKSPGetParametersEW" 3142 /*@ 3143 SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker 3144 convergence criteria for the linear solvers within an inexact 3145 Newton method. 3146 3147 Not Collective 3148 3149 Input Parameters: 3150 snes - SNES context 3151 3152 Output Parameters: 3153 + version - version 1, 2 (default is 2) or 3 3154 . rtol_0 - initial relative tolerance (0 <= rtol_0 < 1) 3155 . rtol_max - maximum relative tolerance (0 <= rtol_max < 1) 3156 . gamma - multiplicative factor for version 2 rtol computation 3157 (0 <= gamma2 <= 1) 3158 . alpha - power for version 2 rtol computation (1 < alpha <= 2) 3159 . alpha2 - power for safeguard 3160 - threshold - threshold for imposing safeguard (0 < threshold < 1) 3161 3162 Level: advanced 3163 3164 .keywords: SNES, KSP, Eisenstat, Walker, get, parameters 3165 3166 .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW() 3167 @*/ 3168 PetscErrorCode SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max, 3169 PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold) 3170 { 3171 SNESKSPEW *kctx; 3172 PetscFunctionBegin; 3173 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3174 kctx = (SNESKSPEW*)snes->kspconvctx; 3175 if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing"); 3176 if(version) *version = kctx->version; 3177 if(rtol_0) *rtol_0 = kctx->rtol_0; 3178 if(rtol_max) *rtol_max = kctx->rtol_max; 3179 if(gamma) *gamma = kctx->gamma; 3180 if(alpha) *alpha = kctx->alpha; 3181 if(alpha2) *alpha2 = kctx->alpha2; 3182 if(threshold) *threshold = kctx->threshold; 3183 PetscFunctionReturn(0); 3184 } 3185 3186 #undef __FUNCT__ 3187 #define __FUNCT__ "SNESKSPEW_PreSolve" 3188 static PetscErrorCode SNESKSPEW_PreSolve(SNES snes, KSP ksp, Vec b, Vec x) 3189 { 3190 PetscErrorCode ierr; 3191 SNESKSPEW *kctx = (SNESKSPEW*)snes->kspconvctx; 3192 PetscReal rtol=PETSC_DEFAULT,stol; 3193 3194 PetscFunctionBegin; 3195 if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists"); 3196 if (!snes->iter) { /* first time in, so use the original user rtol */ 3197 rtol = kctx->rtol_0; 3198 } else { 3199 if (kctx->version == 1) { 3200 rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last; 3201 if (rtol < 0.0) rtol = -rtol; 3202 stol = pow(kctx->rtol_last,kctx->alpha2); 3203 if (stol > kctx->threshold) rtol = PetscMax(rtol,stol); 3204 } else if (kctx->version == 2) { 3205 rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha); 3206 stol = kctx->gamma * pow(kctx->rtol_last,kctx->alpha); 3207 if (stol > kctx->threshold) rtol = PetscMax(rtol,stol); 3208 } else if (kctx->version == 3) {/* contributed by Luis Chacon, June 2006. */ 3209 rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha); 3210 /* safeguard: avoid sharp decrease of rtol */ 3211 stol = kctx->gamma*pow(kctx->rtol_last,kctx->alpha); 3212 stol = PetscMax(rtol,stol); 3213 rtol = PetscMin(kctx->rtol_0,stol); 3214 /* safeguard: avoid oversolving */ 3215 stol = kctx->gamma*(snes->ttol)/snes->norm; 3216 stol = PetscMax(rtol,stol); 3217 rtol = PetscMin(kctx->rtol_0,stol); 3218 } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 or 3 are supported: %D",kctx->version); 3219 } 3220 /* safeguard: avoid rtol greater than one */ 3221 rtol = PetscMin(rtol,kctx->rtol_max); 3222 ierr = KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 3223 ierr = PetscInfo3(snes,"iter %D, Eisenstat-Walker (version %D) KSP rtol=%G\n",snes->iter,kctx->version,rtol);CHKERRQ(ierr); 3224 PetscFunctionReturn(0); 3225 } 3226 3227 #undef __FUNCT__ 3228 #define __FUNCT__ "SNESKSPEW_PostSolve" 3229 static PetscErrorCode SNESKSPEW_PostSolve(SNES snes, KSP ksp, Vec b, Vec x) 3230 { 3231 PetscErrorCode ierr; 3232 SNESKSPEW *kctx = (SNESKSPEW*)snes->kspconvctx; 3233 PCSide pcside; 3234 Vec lres; 3235 3236 PetscFunctionBegin; 3237 if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists"); 3238 ierr = KSPGetTolerances(ksp,&kctx->rtol_last,0,0,0);CHKERRQ(ierr); 3239 ierr = SNESGetFunctionNorm(snes,&kctx->norm_last);CHKERRQ(ierr); 3240 if (kctx->version == 1) { 3241 ierr = KSPGetPCSide(ksp,&pcside);CHKERRQ(ierr); 3242 if (pcside == PC_RIGHT) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */ 3243 /* KSP residual is true linear residual */ 3244 ierr = KSPGetResidualNorm(ksp,&kctx->lresid_last);CHKERRQ(ierr); 3245 } else { 3246 /* KSP residual is preconditioned residual */ 3247 /* compute true linear residual norm */ 3248 ierr = VecDuplicate(b,&lres);CHKERRQ(ierr); 3249 ierr = MatMult(snes->jacobian,x,lres);CHKERRQ(ierr); 3250 ierr = VecAYPX(lres,-1.0,b);CHKERRQ(ierr); 3251 ierr = VecNorm(lres,NORM_2,&kctx->lresid_last);CHKERRQ(ierr); 3252 ierr = VecDestroy(&lres);CHKERRQ(ierr); 3253 } 3254 } 3255 PetscFunctionReturn(0); 3256 } 3257 3258 #undef __FUNCT__ 3259 #define __FUNCT__ "SNES_KSPSolve" 3260 PetscErrorCode SNES_KSPSolve(SNES snes, KSP ksp, Vec b, Vec x) 3261 { 3262 PetscErrorCode ierr; 3263 3264 PetscFunctionBegin; 3265 if (snes->ksp_ewconv) { ierr = SNESKSPEW_PreSolve(snes,ksp,b,x);CHKERRQ(ierr); } 3266 ierr = KSPSolve(ksp,b,x);CHKERRQ(ierr); 3267 if (snes->ksp_ewconv) { ierr = SNESKSPEW_PostSolve(snes,ksp,b,x);CHKERRQ(ierr); } 3268 PetscFunctionReturn(0); 3269 } 3270 3271 #undef __FUNCT__ 3272 #define __FUNCT__ "SNESSetDM" 3273 /*@ 3274 SNESSetDM - Sets the DM that may be used by some preconditioners 3275 3276 Logically Collective on SNES 3277 3278 Input Parameters: 3279 + snes - the preconditioner context 3280 - dm - the dm 3281 3282 Level: intermediate 3283 3284 3285 .seealso: SNESGetDM(), KSPSetDM(), KSPGetDM() 3286 @*/ 3287 PetscErrorCode SNESSetDM(SNES snes,DM dm) 3288 { 3289 PetscErrorCode ierr; 3290 KSP ksp; 3291 3292 PetscFunctionBegin; 3293 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3294 if (dm) {ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);} 3295 ierr = DMDestroy(&snes->dm);CHKERRQ(ierr); 3296 snes->dm = dm; 3297 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 3298 ierr = KSPSetDM(ksp,dm);CHKERRQ(ierr); 3299 ierr = KSPSetDMActive(ksp,PETSC_FALSE);CHKERRQ(ierr); 3300 PetscFunctionReturn(0); 3301 } 3302 3303 #undef __FUNCT__ 3304 #define __FUNCT__ "SNESGetDM" 3305 /*@ 3306 SNESGetDM - Gets the DM that may be used by some preconditioners 3307 3308 Not Collective but DM obtained is parallel on SNES 3309 3310 Input Parameter: 3311 . snes - the preconditioner context 3312 3313 Output Parameter: 3314 . dm - the dm 3315 3316 Level: intermediate 3317 3318 3319 .seealso: SNESSetDM(), KSPSetDM(), KSPGetDM() 3320 @*/ 3321 PetscErrorCode SNESGetDM(SNES snes,DM *dm) 3322 { 3323 PetscFunctionBegin; 3324 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3325 *dm = snes->dm; 3326 PetscFunctionReturn(0); 3327 } 3328 3329 #if defined(PETSC_HAVE_MATLAB_ENGINE) 3330 #include <mex.h> 3331 3332 typedef struct {char *funcname; mxArray *ctx;} SNESMatlabContext; 3333 3334 #undef __FUNCT__ 3335 #define __FUNCT__ "SNESComputeFunction_Matlab" 3336 /* 3337 SNESComputeFunction_Matlab - Calls the function that has been set with 3338 SNESSetFunctionMatlab(). 3339 3340 Collective on SNES 3341 3342 Input Parameters: 3343 + snes - the SNES context 3344 - x - input vector 3345 3346 Output Parameter: 3347 . y - function vector, as set by SNESSetFunction() 3348 3349 Notes: 3350 SNESComputeFunction() is typically used within nonlinear solvers 3351 implementations, so most users would not generally call this routine 3352 themselves. 3353 3354 Level: developer 3355 3356 .keywords: SNES, nonlinear, compute, function 3357 3358 .seealso: SNESSetFunction(), SNESGetFunction() 3359 */ 3360 PetscErrorCode SNESComputeFunction_Matlab(SNES snes,Vec x,Vec y, void *ctx) 3361 { 3362 PetscErrorCode ierr; 3363 SNESMatlabContext *sctx = (SNESMatlabContext *)ctx; 3364 int nlhs = 1,nrhs = 5; 3365 mxArray *plhs[1],*prhs[5]; 3366 long long int lx = 0,ly = 0,ls = 0; 3367 3368 PetscFunctionBegin; 3369 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3370 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 3371 PetscValidHeaderSpecific(y,VEC_CLASSID,3); 3372 PetscCheckSameComm(snes,1,x,2); 3373 PetscCheckSameComm(snes,1,y,3); 3374 3375 /* call Matlab function in ctx with arguments x and y */ 3376 3377 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 3378 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 3379 ierr = PetscMemcpy(&ly,&y,sizeof(x));CHKERRQ(ierr); 3380 prhs[0] = mxCreateDoubleScalar((double)ls); 3381 prhs[1] = mxCreateDoubleScalar((double)lx); 3382 prhs[2] = mxCreateDoubleScalar((double)ly); 3383 prhs[3] = mxCreateString(sctx->funcname); 3384 prhs[4] = sctx->ctx; 3385 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeFunctionInternal");CHKERRQ(ierr); 3386 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 3387 mxDestroyArray(prhs[0]); 3388 mxDestroyArray(prhs[1]); 3389 mxDestroyArray(prhs[2]); 3390 mxDestroyArray(prhs[3]); 3391 mxDestroyArray(plhs[0]); 3392 PetscFunctionReturn(0); 3393 } 3394 3395 3396 #undef __FUNCT__ 3397 #define __FUNCT__ "SNESSetFunctionMatlab" 3398 /* 3399 SNESSetFunctionMatlab - Sets the function evaluation routine and function 3400 vector for use by the SNES routines in solving systems of nonlinear 3401 equations from MATLAB. Here the function is a string containing the name of a MATLAB function 3402 3403 Logically Collective on SNES 3404 3405 Input Parameters: 3406 + snes - the SNES context 3407 . r - vector to store function value 3408 - func - function evaluation routine 3409 3410 Calling sequence of func: 3411 $ func (SNES snes,Vec x,Vec f,void *ctx); 3412 3413 3414 Notes: 3415 The Newton-like methods typically solve linear systems of the form 3416 $ f'(x) x = -f(x), 3417 where f'(x) denotes the Jacobian matrix and f(x) is the function. 3418 3419 Level: beginner 3420 3421 .keywords: SNES, nonlinear, set, function 3422 3423 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction() 3424 */ 3425 PetscErrorCode SNESSetFunctionMatlab(SNES snes,Vec r,const char *func,mxArray *ctx) 3426 { 3427 PetscErrorCode ierr; 3428 SNESMatlabContext *sctx; 3429 3430 PetscFunctionBegin; 3431 /* currently sctx is memory bleed */ 3432 ierr = PetscMalloc(sizeof(SNESMatlabContext),&sctx);CHKERRQ(ierr); 3433 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 3434 /* 3435 This should work, but it doesn't 3436 sctx->ctx = ctx; 3437 mexMakeArrayPersistent(sctx->ctx); 3438 */ 3439 sctx->ctx = mxDuplicateArray(ctx); 3440 ierr = SNESSetFunction(snes,r,SNESComputeFunction_Matlab,sctx);CHKERRQ(ierr); 3441 PetscFunctionReturn(0); 3442 } 3443 3444 #undef __FUNCT__ 3445 #define __FUNCT__ "SNESComputeJacobian_Matlab" 3446 /* 3447 SNESComputeJacobian_Matlab - Calls the function that has been set with 3448 SNESSetJacobianMatlab(). 3449 3450 Collective on SNES 3451 3452 Input Parameters: 3453 + snes - the SNES context 3454 . x - input vector 3455 . A, B - the matrices 3456 - ctx - user context 3457 3458 Output Parameter: 3459 . flag - structure of the matrix 3460 3461 Level: developer 3462 3463 .keywords: SNES, nonlinear, compute, function 3464 3465 .seealso: SNESSetFunction(), SNESGetFunction() 3466 @*/ 3467 PetscErrorCode SNESComputeJacobian_Matlab(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag, void *ctx) 3468 { 3469 PetscErrorCode ierr; 3470 SNESMatlabContext *sctx = (SNESMatlabContext *)ctx; 3471 int nlhs = 2,nrhs = 6; 3472 mxArray *plhs[2],*prhs[6]; 3473 long long int lx = 0,lA = 0,ls = 0, lB = 0; 3474 3475 PetscFunctionBegin; 3476 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3477 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 3478 3479 /* call Matlab function in ctx with arguments x and y */ 3480 3481 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 3482 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 3483 ierr = PetscMemcpy(&lA,A,sizeof(x));CHKERRQ(ierr); 3484 ierr = PetscMemcpy(&lB,B,sizeof(x));CHKERRQ(ierr); 3485 prhs[0] = mxCreateDoubleScalar((double)ls); 3486 prhs[1] = mxCreateDoubleScalar((double)lx); 3487 prhs[2] = mxCreateDoubleScalar((double)lA); 3488 prhs[3] = mxCreateDoubleScalar((double)lB); 3489 prhs[4] = mxCreateString(sctx->funcname); 3490 prhs[5] = sctx->ctx; 3491 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeJacobianInternal");CHKERRQ(ierr); 3492 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 3493 *flag = (MatStructure) mxGetScalar(plhs[1]);CHKERRQ(ierr); 3494 mxDestroyArray(prhs[0]); 3495 mxDestroyArray(prhs[1]); 3496 mxDestroyArray(prhs[2]); 3497 mxDestroyArray(prhs[3]); 3498 mxDestroyArray(prhs[4]); 3499 mxDestroyArray(plhs[0]); 3500 mxDestroyArray(plhs[1]); 3501 PetscFunctionReturn(0); 3502 } 3503 3504 3505 #undef __FUNCT__ 3506 #define __FUNCT__ "SNESSetJacobianMatlab" 3507 /* 3508 SNESSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices 3509 vector for use by the SNES routines in solving systems of nonlinear 3510 equations from MATLAB. Here the function is a string containing the name of a MATLAB function 3511 3512 Logically Collective on SNES 3513 3514 Input Parameters: 3515 + snes - the SNES context 3516 . A,B - Jacobian matrices 3517 . func - function evaluation routine 3518 - ctx - user context 3519 3520 Calling sequence of func: 3521 $ flag = func (SNES snes,Vec x,Mat A,Mat B,void *ctx); 3522 3523 3524 Level: developer 3525 3526 .keywords: SNES, nonlinear, set, function 3527 3528 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction() 3529 */ 3530 PetscErrorCode SNESSetJacobianMatlab(SNES snes,Mat A,Mat B,const char *func,mxArray *ctx) 3531 { 3532 PetscErrorCode ierr; 3533 SNESMatlabContext *sctx; 3534 3535 PetscFunctionBegin; 3536 /* currently sctx is memory bleed */ 3537 ierr = PetscMalloc(sizeof(SNESMatlabContext),&sctx);CHKERRQ(ierr); 3538 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 3539 /* 3540 This should work, but it doesn't 3541 sctx->ctx = ctx; 3542 mexMakeArrayPersistent(sctx->ctx); 3543 */ 3544 sctx->ctx = mxDuplicateArray(ctx); 3545 ierr = SNESSetJacobian(snes,A,B,SNESComputeJacobian_Matlab,sctx);CHKERRQ(ierr); 3546 PetscFunctionReturn(0); 3547 } 3548 3549 #undef __FUNCT__ 3550 #define __FUNCT__ "SNESMonitor_Matlab" 3551 /* 3552 SNESMonitor_Matlab - Calls the function that has been set with SNESMonitorSetMatlab(). 3553 3554 Collective on SNES 3555 3556 .seealso: SNESSetFunction(), SNESGetFunction() 3557 @*/ 3558 PetscErrorCode SNESMonitor_Matlab(SNES snes,PetscInt it, PetscReal fnorm, void *ctx) 3559 { 3560 PetscErrorCode ierr; 3561 SNESMatlabContext *sctx = (SNESMatlabContext *)ctx; 3562 int nlhs = 1,nrhs = 6; 3563 mxArray *plhs[1],*prhs[6]; 3564 long long int lx = 0,ls = 0; 3565 Vec x=snes->vec_sol; 3566 3567 PetscFunctionBegin; 3568 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3569 3570 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 3571 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 3572 prhs[0] = mxCreateDoubleScalar((double)ls); 3573 prhs[1] = mxCreateDoubleScalar((double)it); 3574 prhs[2] = mxCreateDoubleScalar((double)fnorm); 3575 prhs[3] = mxCreateDoubleScalar((double)lx); 3576 prhs[4] = mxCreateString(sctx->funcname); 3577 prhs[5] = sctx->ctx; 3578 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESMonitorInternal");CHKERRQ(ierr); 3579 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 3580 mxDestroyArray(prhs[0]); 3581 mxDestroyArray(prhs[1]); 3582 mxDestroyArray(prhs[2]); 3583 mxDestroyArray(prhs[3]); 3584 mxDestroyArray(prhs[4]); 3585 mxDestroyArray(plhs[0]); 3586 PetscFunctionReturn(0); 3587 } 3588 3589 3590 #undef __FUNCT__ 3591 #define __FUNCT__ "SNESMonitorSetMatlab" 3592 /* 3593 SNESMonitorSetMatlab - Sets the monitor function from MATLAB 3594 3595 Level: developer 3596 3597 .keywords: SNES, nonlinear, set, function 3598 3599 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction() 3600 */ 3601 PetscErrorCode SNESMonitorSetMatlab(SNES snes,const char *func,mxArray *ctx) 3602 { 3603 PetscErrorCode ierr; 3604 SNESMatlabContext *sctx; 3605 3606 PetscFunctionBegin; 3607 /* currently sctx is memory bleed */ 3608 ierr = PetscMalloc(sizeof(SNESMatlabContext),&sctx);CHKERRQ(ierr); 3609 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 3610 /* 3611 This should work, but it doesn't 3612 sctx->ctx = ctx; 3613 mexMakeArrayPersistent(sctx->ctx); 3614 */ 3615 sctx->ctx = mxDuplicateArray(ctx); 3616 ierr = SNESMonitorSet(snes,SNESMonitor_Matlab,sctx,PETSC_NULL);CHKERRQ(ierr); 3617 PetscFunctionReturn(0); 3618 } 3619 3620 #endif 3621