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