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