1 #ifndef lint 2 static char vcid[] = "$Id: snes.c,v 1.80 1996/08/12 03:42:55 bsmith Exp curfman $"; 3 #endif 4 5 #include "draw.h" /*I "draw.h" I*/ 6 #include "src/snes/snesimpl.h" /*I "snes.h" I*/ 7 #include "src/sys/nreg.h" 8 #include "pinclude/pviewer.h" 9 #include <math.h> 10 11 extern int SNESGetTypeFromOptions_Private(SNES,SNESType*,int*); 12 extern int SNESPrintTypes_Private(MPI_Comm,char*,char*); 13 14 /*@ 15 SNESView - Prints the SNES data structure. 16 17 Input Parameters: 18 . SNES - the SNES context 19 . viewer - visualization context 20 21 Options Database Key: 22 $ -snes_view : calls SNESView() at end of SNESSolve() 23 24 Notes: 25 The available visualization contexts include 26 $ VIEWER_STDOUT_SELF - standard output (default) 27 $ VIEWER_STDOUT_WORLD - synchronized standard 28 $ output where only the first processor opens 29 $ the file. All other processors send their 30 $ data to the first processor to print. 31 32 The user can open alternative vistualization contexts with 33 $ ViewerFileOpenASCII() - output to a specified file 34 35 .keywords: SNES, view 36 37 .seealso: ViewerFileOpenASCII() 38 @*/ 39 int SNESView(SNES snes,Viewer viewer) 40 { 41 SNES_KSP_EW_ConvCtx *kctx; 42 FILE *fd; 43 int ierr; 44 SLES sles; 45 char *method; 46 ViewerType vtype; 47 48 PetscValidHeaderSpecific(snes,SNES_COOKIE); 49 if (viewer) {PetscValidHeader(viewer);} 50 else { viewer = VIEWER_STDOUT_SELF; } 51 52 ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 53 if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 54 ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 55 PetscFPrintf(snes->comm,fd,"SNES Object:\n"); 56 SNESGetType(snes,PETSC_NULL,&method); 57 PetscFPrintf(snes->comm,fd," method: %s\n",method); 58 if (snes->view) (*snes->view)((PetscObject)snes,viewer); 59 PetscFPrintf(snes->comm,fd, 60 " maximum iterations=%d, maximum function evaluations=%d\n", 61 snes->max_its,snes->max_funcs); 62 PetscFPrintf(snes->comm,fd, 63 " tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n", 64 snes->rtol, snes->atol, snes->trunctol, snes->xtol); 65 if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) 66 PetscFPrintf(snes->comm,fd," min function tolerance=%g\n",snes->fmin); 67 if (snes->ksp_ewconv) { 68 kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 69 if (kctx) { 70 PetscFPrintf(snes->comm,fd, 71 " Eisenstat-Walker computation of KSP relative tolerance (version %d)\n", 72 kctx->version); 73 PetscFPrintf(snes->comm,fd, 74 " rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0, 75 kctx->rtol_max,kctx->threshold); 76 PetscFPrintf(snes->comm,fd," gamma=%g, alpha=%g, alpha2=%g\n", 77 kctx->gamma,kctx->alpha,kctx->alpha2); 78 } 79 } 80 } else if (vtype == STRING_VIEWER) { 81 SNESGetType(snes,PETSC_NULL,&method); 82 ViewerStringSPrintf(viewer," %-3.3s",method); 83 } 84 SNESGetSLES(snes,&sles); 85 ierr = SLESView(sles,viewer); CHKERRQ(ierr); 86 return 0; 87 } 88 89 /*@ 90 SNESSetFromOptions - Sets various SNES and SLES parameters from user options. 91 92 Input Parameter: 93 . snes - the SNES context 94 95 .keywords: SNES, nonlinear, set, options, database 96 97 .seealso: SNESPrintHelp(), SNESSetOptionsPrefix() 98 @*/ 99 int SNESSetFromOptions(SNES snes) 100 { 101 SNESType method; 102 double tmp; 103 SLES sles; 104 int ierr, flg, mset = 0; 105 int version = PETSC_DEFAULT; 106 double rtol_0 = PETSC_DEFAULT; 107 double rtol_max = PETSC_DEFAULT; 108 double gamma2 = PETSC_DEFAULT; 109 double alpha = PETSC_DEFAULT; 110 double alpha2 = PETSC_DEFAULT; 111 double threshold = PETSC_DEFAULT; 112 113 PetscValidHeaderSpecific(snes,SNES_COOKIE); 114 if (snes->setup_called) SETERRQ(1,"SNESSetFromOptions:Must call prior to SNESSetUp"); 115 ierr = SNESGetTypeFromOptions_Private(snes,&method,&flg); CHKERRQ(ierr); 116 if (flg) { 117 ierr = SNESSetType(snes,method); CHKERRQ(ierr); 118 } 119 else if (!snes->set_method_called) { 120 if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 121 ierr = SNESSetType(snes,SNES_EQ_LS); CHKERRQ(ierr); 122 } 123 else { 124 ierr = SNESSetType(snes,SNES_UM_TR); CHKERRQ(ierr); 125 } 126 } 127 128 ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 129 if (flg) { SNESPrintHelp(snes); } 130 ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg); CHKERRQ(ierr); 131 if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);} 132 ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg); CHKERRQ(ierr); 133 if (flg) {SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);} 134 ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg); CHKERRQ(ierr); 135 if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT); } 136 ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, &flg); CHKERRQ(ierr); 137 ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, &flg);CHKERRQ(ierr); 138 ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg); CHKERRQ(ierr); 139 if (flg) { SNESSetTrustRegionTolerance(snes,tmp); } 140 ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg); CHKERRQ(ierr); 141 if (flg) { SNESSetMinimizationFunctionTolerance(snes,tmp); } 142 ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg); CHKERRQ(ierr); 143 if (flg) { snes->ksp_ewconv = 1; } 144 ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,&flg); CHKERRQ(ierr); 145 ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,&flg); CHKERRQ(ierr); 146 ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,&flg);CHKERRQ(ierr); 147 ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,&flg); CHKERRQ(ierr); 148 ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,&flg); CHKERRQ(ierr); 149 ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,&flg); CHKERRQ(ierr); 150 ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,&flg);CHKERRQ(ierr); 151 ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha, 152 alpha2,threshold); CHKERRQ(ierr); 153 ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg); CHKERRQ(ierr); 154 if (flg) {SNESSetMonitor(snes,SNESDefaultMonitor,0); mset = 1;} 155 ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg); CHKERRQ(ierr); 156 if (flg) { 157 if (mset) SETERRQ(1,"SNESSetFromOptions:Monitor for SNES is already set"); 158 SNESSetMonitor(snes,SNESDefaultSMonitor,0); mset = 1; 159 } 160 ierr = OptionsHasName(snes->prefix,"-snes_xmonitor",&flg); CHKERRQ(ierr); 161 if (flg) { 162 int rank = 0; 163 DrawLG lg; 164 if (mset) SETERRQ(1,"SNESSetFromOptions:Monitor for SNES is already set"); 165 MPI_Initialized(&rank); 166 if (rank) MPI_Comm_rank(snes->comm,&rank); 167 if (!rank) { 168 ierr = SNESLGMonitorCreate(0,0,0,0,300,300,&lg); CHKERRQ(ierr); 169 ierr = SNESSetMonitor(snes,SNESLGMonitor,(void *)lg); CHKERRQ(ierr); 170 snes->xmonitor = lg; 171 PLogObjectParent(snes,lg); 172 } 173 mset = 1; 174 } 175 ierr = OptionsHasName(snes->prefix,"-snes_fd", &flg); CHKERRQ(ierr); 176 if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) { 177 ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre, 178 SNESDefaultComputeJacobian,snes->funP); CHKERRQ(ierr); 179 PLogInfo(snes, 180 "SNESSetFromOptions: Setting default finite difference Jacobian matrix\n"); 181 } 182 else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 183 ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre, 184 SNESDefaultComputeHessian,snes->funP); CHKERRQ(ierr); 185 PLogInfo(snes, 186 "SNESSetFromOptions: Setting default finite difference Hessian matrix\n"); 187 } 188 ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 189 ierr = SLESSetFromOptions(sles); CHKERRQ(ierr); 190 if (!snes->setfromoptions) return 0; 191 return (*snes->setfromoptions)(snes); 192 } 193 194 /*@ 195 SNESPrintHelp - Prints all options for the SNES component. 196 197 Input Parameter: 198 . snes - the SNES context 199 200 Options Database Keys: 201 $ -help, -h 202 203 .keywords: SNES, nonlinear, help 204 205 .seealso: SNESSetFromOptions() 206 @*/ 207 int SNESPrintHelp(SNES snes) 208 { 209 char p[64]; 210 SNES_KSP_EW_ConvCtx *kctx; 211 212 PetscValidHeaderSpecific(snes,SNES_COOKIE); 213 214 PetscStrcpy(p,"-"); 215 if (snes->prefix) PetscStrcat(p, snes->prefix); 216 217 kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 218 219 PetscPrintf(snes->comm,"SNES options ------------------------------------------------\n"); 220 SNESPrintTypes_Private(snes->comm,p,"snes_type"); 221 PetscPrintf(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p); 222 PetscPrintf(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its); 223 PetscPrintf(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs); 224 PetscPrintf(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol); 225 PetscPrintf(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol); 226 PetscPrintf(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol); 227 PetscPrintf(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol); 228 PetscPrintf(snes->comm," SNES Monitoring Options: Choose 1 of the following\n"); 229 PetscPrintf(snes->comm," %ssnes_monitor: use default SNES convergence monitor\n",p); 230 PetscPrintf(snes->comm," %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p); 231 if (snes->type == SNES_NONLINEAR_EQUATIONS) { 232 PetscPrintf(snes->comm, 233 " Options for solving systems of nonlinear equations only:\n"); 234 PetscPrintf(snes->comm," %ssnes_fd: use finite differences for Jacobian\n",p); 235 PetscPrintf(snes->comm," %ssnes_mf: use matrix-free Jacobian\n",p); 236 PetscPrintf(snes->comm," %ssnes_mf_err: relative error in function evaluation for matrix-free Jacobian\n",p); 237 PetscPrintf(snes->comm," %ssnes_mf_umin: minimum iterate parameter for matrix-free Jacobian\n",p); 238 PetscPrintf(snes->comm," %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p); 239 PetscPrintf(snes->comm, 240 " %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version); 241 PetscPrintf(snes->comm, 242 " %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0); 243 PetscPrintf(snes->comm, 244 " %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max); 245 PetscPrintf(snes->comm, 246 " %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma); 247 PetscPrintf(snes->comm, 248 " %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha); 249 PetscPrintf(snes->comm, 250 " %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2); 251 PetscPrintf(snes->comm, 252 " %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold); 253 } 254 else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) { 255 PetscPrintf(snes->comm, 256 " Options for solving unconstrained minimization problems only:\n"); 257 PetscPrintf(snes->comm," %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin); 258 PetscPrintf(snes->comm," %ssnes_fd: use finite differences for Hessian\n",p); 259 PetscPrintf(snes->comm," %ssnes_mf: use matrix-free Hessian\n",p); 260 PetscPrintf(snes->comm," %ssnes_mf_err: relative error in gradient evaluation for matrix-free Hessian\n",p); 261 PetscPrintf(snes->comm," %ssnes_mf_umin: minimum iterate parameter for matrix-free Hessian\n",p); 262 } 263 PetscPrintf(snes->comm," Run program with -help %ssnes_type <method> for help on ",p); 264 PetscPrintf(snes->comm,"a particular method\n"); 265 if (snes->printhelp) (*snes->printhelp)(snes,p); 266 return 0; 267 } 268 /*@ 269 SNESSetApplicationContext - Sets the optional user-defined context for 270 the nonlinear solvers. 271 272 Input Parameters: 273 . snes - the SNES context 274 . usrP - optional user context 275 276 .keywords: SNES, nonlinear, set, application, context 277 278 .seealso: SNESGetApplicationContext() 279 @*/ 280 int SNESSetApplicationContext(SNES snes,void *usrP) 281 { 282 PetscValidHeaderSpecific(snes,SNES_COOKIE); 283 snes->user = usrP; 284 return 0; 285 } 286 287 /*@C 288 SNESGetApplicationContext - Gets the user-defined context for the 289 nonlinear solvers. 290 291 Input Parameter: 292 . snes - SNES context 293 294 Output Parameter: 295 . usrP - user context 296 297 .keywords: SNES, nonlinear, get, application, context 298 299 .seealso: SNESSetApplicationContext() 300 @*/ 301 int SNESGetApplicationContext( SNES snes, void **usrP ) 302 { 303 PetscValidHeaderSpecific(snes,SNES_COOKIE); 304 *usrP = snes->user; 305 return 0; 306 } 307 308 /*@ 309 SNESGetIterationNumber - Gets the current iteration number of the 310 nonlinear solver. 311 312 Input Parameter: 313 . snes - SNES context 314 315 Output Parameter: 316 . iter - iteration number 317 318 .keywords: SNES, nonlinear, get, iteration, number 319 @*/ 320 int SNESGetIterationNumber(SNES snes,int* iter) 321 { 322 PetscValidHeaderSpecific(snes,SNES_COOKIE); 323 PetscValidIntPointer(iter); 324 *iter = snes->iter; 325 return 0; 326 } 327 328 /*@ 329 SNESGetFunctionNorm - Gets the norm of the current function that was set 330 with SNESSSetFunction(). 331 332 Input Parameter: 333 . snes - SNES context 334 335 Output Parameter: 336 . fnorm - 2-norm of function 337 338 Note: 339 SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only. 340 A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is 341 SNESGetGradientNorm(). 342 343 .keywords: SNES, nonlinear, get, function, norm 344 345 .seealso: SNESSetFunction() 346 @*/ 347 int SNESGetFunctionNorm(SNES snes,Scalar *fnorm) 348 { 349 PetscValidHeaderSpecific(snes,SNES_COOKIE); 350 PetscValidScalarPointer(fnorm); 351 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 352 SETERRQ(1,"SNESGetFunctionNorm:For SNES_NONLINEAR_EQUATIONS only"); 353 } 354 *fnorm = snes->norm; 355 return 0; 356 } 357 358 /*@ 359 SNESGetGradientNorm - Gets the norm of the current gradient that was set 360 with SNESSSetGradient(). 361 362 Input Parameter: 363 . snes - SNES context 364 365 Output Parameter: 366 . fnorm - 2-norm of gradient 367 368 Note: 369 SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION 370 methods only. A related routine for SNES_NONLINEAR_EQUATIONS methods 371 is SNESGetFunctionNorm(). 372 373 .keywords: SNES, nonlinear, get, gradient, norm 374 375 .seelso: SNESSetGradient() 376 @*/ 377 int SNESGetGradientNorm(SNES snes,Scalar *gnorm) 378 { 379 PetscValidHeaderSpecific(snes,SNES_COOKIE); 380 PetscValidScalarPointer(gnorm); 381 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 382 SETERRQ(1,"SNESGetGradientNorm:For SNES_UNCONSTRAINED_MINIMIZATION only"); 383 } 384 *gnorm = snes->norm; 385 return 0; 386 } 387 388 /*@ 389 SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps 390 attempted by the nonlinear solver. 391 392 Input Parameter: 393 . snes - SNES context 394 395 Output Parameter: 396 . nfails - number of unsuccessful steps attempted 397 398 .keywords: SNES, nonlinear, get, number, unsuccessful, steps 399 @*/ 400 int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails) 401 { 402 PetscValidHeaderSpecific(snes,SNES_COOKIE); 403 PetscValidIntPointer(nfails); 404 *nfails = snes->nfailures; 405 return 0; 406 } 407 /*@C 408 SNESGetSLES - Returns the SLES context for a SNES solver. 409 410 Input Parameter: 411 . snes - the SNES context 412 413 Output Parameter: 414 . sles - the SLES context 415 416 Notes: 417 The user can then directly manipulate the SLES context to set various 418 options, etc. Likewise, the user can then extract and manipulate the 419 KSP and PC contexts as well. 420 421 .keywords: SNES, nonlinear, get, SLES, context 422 423 .seealso: SLESGetPC(), SLESGetKSP() 424 @*/ 425 int SNESGetSLES(SNES snes,SLES *sles) 426 { 427 PetscValidHeaderSpecific(snes,SNES_COOKIE); 428 *sles = snes->sles; 429 return 0; 430 } 431 /* -----------------------------------------------------------*/ 432 /*@C 433 SNESCreate - Creates a nonlinear solver context. 434 435 Input Parameter: 436 . comm - MPI communicator 437 . type - type of method, one of 438 $ SNES_NONLINEAR_EQUATIONS 439 $ (for systems of nonlinear equations) 440 $ SNES_UNCONSTRAINED_MINIMIZATION 441 $ (for unconstrained minimization) 442 443 Output Parameter: 444 . outsnes - the new SNES context 445 446 .keywords: SNES, nonlinear, create, context 447 448 .seealso: SNESSolve(), SNESDestroy() 449 @*/ 450 int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes) 451 { 452 int ierr; 453 SNES snes; 454 SNES_KSP_EW_ConvCtx *kctx; 455 456 *outsnes = 0; 457 if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS) 458 SETERRQ(1,"SNESCreate:incorrect method type"); 459 PetscHeaderCreate(snes,_SNES,SNES_COOKIE,SNES_UNKNOWN_METHOD,comm); 460 PLogObjectCreate(snes); 461 snes->max_its = 50; 462 snes->max_funcs = 1000; 463 snes->norm = 0.0; 464 if (type == SNES_UNCONSTRAINED_MINIMIZATION) { 465 snes->rtol = 1.e-8; 466 snes->ttol = 0.0; 467 snes->atol = 1.e-10; 468 } 469 else { 470 snes->rtol = 1.e-8; 471 snes->ttol = 0.0; 472 snes->atol = 1.e-50; 473 } 474 snes->xtol = 1.e-8; 475 snes->trunctol = 1.e-12; /* no longer used */ 476 snes->nfuncs = 0; 477 snes->nfailures = 0; 478 snes->monitor = 0; 479 snes->data = 0; 480 snes->view = 0; 481 snes->computeumfunction = 0; 482 snes->umfunP = 0; 483 snes->fc = 0; 484 snes->deltatol = 1.e-12; 485 snes->fmin = -1.e30; 486 snes->method_class = type; 487 snes->set_method_called = 0; 488 snes->setup_called = 0; 489 snes->ksp_ewconv = 0; 490 snes->type = -1; 491 snes->xmonitor = 0; 492 snes->vwork = 0; 493 snes->nwork = 0; 494 495 /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 496 kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx); 497 snes->kspconvctx = (void*)kctx; 498 kctx->version = 2; 499 kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 500 this was too large for some test cases */ 501 kctx->rtol_last = 0; 502 kctx->rtol_max = .9; 503 kctx->gamma = 1.0; 504 kctx->alpha2 = .5*(1.0 + sqrt(5.0)); 505 kctx->alpha = kctx->alpha2; 506 kctx->threshold = .1; 507 kctx->lresid_last = 0; 508 kctx->norm_last = 0; 509 510 ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr); 511 PLogObjectParent(snes,snes->sles) 512 513 *outsnes = snes; 514 return 0; 515 } 516 517 /* --------------------------------------------------------------- */ 518 /*@C 519 SNESSetFunction - Sets the function evaluation routine and function 520 vector for use by the SNES routines in solving systems of nonlinear 521 equations. 522 523 Input Parameters: 524 . snes - the SNES context 525 . func - function evaluation routine 526 . ctx - optional user-defined function context 527 . r - vector to store function value 528 529 Calling sequence of func: 530 . func (SNES, Vec x, Vec f, void *ctx); 531 532 . x - input vector 533 . f - vector function 534 . ctx - optional user-defined context for private data for the 535 function evaluation routine (may be null) 536 537 Notes: 538 The Newton-like methods typically solve linear systems of the form 539 $ f'(x) x = -f(x), 540 $ where f'(x) denotes the Jacobian matrix and f(x) is the function. 541 542 SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 543 Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 544 SNESSetMinimizationFunction() and SNESSetGradient(); 545 546 .keywords: SNES, nonlinear, set, function 547 548 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 549 @*/ 550 int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx) 551 { 552 PetscValidHeaderSpecific(snes,SNES_COOKIE); 553 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1, 554 "SNESSetFunction:For SNES_NONLINEAR_EQUATIONS only"); 555 snes->computefunction = func; 556 snes->vec_func = snes->vec_func_always = r; 557 snes->funP = ctx; 558 return 0; 559 } 560 561 /*@ 562 SNESComputeFunction - Computes the function that has been set with 563 SNESSetFunction(). 564 565 Input Parameters: 566 . snes - the SNES context 567 . x - input vector 568 569 Output Parameter: 570 . y - function vector, as set by SNESSetFunction() 571 572 n Notes: 573 SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 574 Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 575 SNESComputeMinimizationFunction() and SNESComputeGradient(); 576 577 .keywords: SNES, nonlinear, compute, function 578 579 .seealso: SNESSetFunction(), SNESGetFunction() 580 @*/ 581 int SNESComputeFunction(SNES snes,Vec x, Vec y) 582 { 583 int ierr; 584 585 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 586 SETERRQ(1,"SNESComputeFunction: For SNES_NONLINEAR_EQUATIONS only"); 587 PLogEventBegin(SNES_FunctionEval,snes,x,y,0); 588 ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 589 PLogEventEnd(SNES_FunctionEval,snes,x,y,0); 590 return 0; 591 } 592 593 /*@C 594 SNESSetMinimizationFunction - Sets the function evaluation routine for 595 unconstrained minimization. 596 597 Input Parameters: 598 . snes - the SNES context 599 . func - function evaluation routine 600 . ctx - optional user-defined function context 601 602 Calling sequence of func: 603 . func (SNES snes,Vec x,double *f,void *ctx); 604 605 . x - input vector 606 . f - function 607 . ctx - optional user-defined context for private data for the 608 function evaluation routine (may be null) 609 610 Notes: 611 SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 612 methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 613 SNESSetFunction(). 614 615 .keywords: SNES, nonlinear, set, minimization, function 616 617 .seealso: SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(), 618 SNESSetHessian(), SNESSetGradient() 619 @*/ 620 int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*), 621 void *ctx) 622 { 623 PetscValidHeaderSpecific(snes,SNES_COOKIE); 624 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1, 625 "SNESSetMinimizationFunction:Only for SNES_UNCONSTRAINED_MINIMIZATION"); 626 snes->computeumfunction = func; 627 snes->umfunP = ctx; 628 return 0; 629 } 630 631 /*@ 632 SNESComputeMinimizationFunction - Computes the function that has been 633 set with SNESSetMinimizationFunction(). 634 635 Input Parameters: 636 . snes - the SNES context 637 . x - input vector 638 639 Output Parameter: 640 . y - function value 641 642 Notes: 643 SNESComputeMinimizationFunction() is valid only for 644 SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 645 SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 646 647 .keywords: SNES, nonlinear, compute, minimization, function 648 649 .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(), 650 SNESComputeGradient(), SNESComputeHessian() 651 @*/ 652 int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y) 653 { 654 int ierr; 655 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1, 656 "SNESComputeMinimizationFunction:Only for SNES_UNCONSTRAINED_MINIMIZATION"); 657 PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0); 658 ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr); 659 PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0); 660 return 0; 661 } 662 663 /*@C 664 SNESSetGradient - Sets the gradient evaluation routine and gradient 665 vector for use by the SNES routines. 666 667 Input Parameters: 668 . snes - the SNES context 669 . func - function evaluation routine 670 . ctx - optional user-defined function context 671 . r - vector to store gradient value 672 673 Calling sequence of func: 674 . func (SNES, Vec x, Vec g, void *ctx); 675 676 . x - input vector 677 . g - gradient vector 678 . ctx - optional user-defined context for private data for the 679 function evaluation routine (may be null) 680 681 Notes: 682 SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 683 methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 684 SNESSetFunction(). 685 686 .keywords: SNES, nonlinear, set, function 687 688 .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(), 689 SNESSetMinimizationFunction(), 690 @*/ 691 int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx) 692 { 693 PetscValidHeaderSpecific(snes,SNES_COOKIE); 694 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1, 695 "SNESSetGradient:For SNES_UNCONSTRAINED_MINIMIZATION only"); 696 snes->computefunction = func; 697 snes->vec_func = snes->vec_func_always = r; 698 snes->funP = ctx; 699 return 0; 700 } 701 702 /*@ 703 SNESComputeGradient - Computes the gradient that has been set with 704 SNESSetGradient(). 705 706 Input Parameters: 707 . snes - the SNES context 708 . x - input vector 709 710 Output Parameter: 711 . y - gradient vector 712 713 Notes: 714 SNESComputeGradient() is valid only for 715 SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 716 SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 717 718 .keywords: SNES, nonlinear, compute, gradient 719 720 .seealso: SNESSetGradient(), SNESGetGradient(), 721 SNESComputeMinimizationFunction(), SNESComputeHessian() 722 @*/ 723 int SNESComputeGradient(SNES snes,Vec x, Vec y) 724 { 725 int ierr; 726 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 727 SETERRQ(1,"SNESComputeGradient:For SNES_UNCONSTRAINED_MINIMIZATION only"); 728 PLogEventBegin(SNES_GradientEval,snes,x,y,0); 729 ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 730 PLogEventEnd(SNES_GradientEval,snes,x,y,0); 731 return 0; 732 } 733 734 /*@ 735 SNESComputeJacobian - Computes the Jacobian matrix that has been 736 set with SNESSetJacobian(). 737 738 Input Parameters: 739 . snes - the SNES context 740 . x - input vector 741 742 Output Parameters: 743 . A - Jacobian matrix 744 . B - optional preconditioning matrix 745 . flag - flag indicating matrix structure 746 747 Notes: 748 Most users should not need to explicitly call this routine, as it 749 is used internally within the nonlinear solvers. 750 751 See SLESSetOperators() for information about setting the flag parameter. 752 753 SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS 754 methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION 755 methods is SNESComputeJacobian(). 756 757 .keywords: SNES, compute, Jacobian, matrix 758 759 .seealso: SNESSetJacobian(), SLESSetOperators() 760 @*/ 761 int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 762 { 763 int ierr; 764 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 765 SETERRQ(1,"SNESComputeJacobian: For SNES_NONLINEAR_EQUATIONS only"); 766 if (!snes->computejacobian) return 0; 767 PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B); 768 *flg = DIFFERENT_NONZERO_PATTERN; 769 ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr); 770 PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B); 771 /* make sure user returned a correct Jacobian and preconditioner */ 772 PetscValidHeaderSpecific(*A,MAT_COOKIE); 773 PetscValidHeaderSpecific(*B,MAT_COOKIE); 774 return 0; 775 } 776 777 /*@ 778 SNESComputeHessian - Computes the Hessian matrix that has been 779 set with SNESSetHessian(). 780 781 Input Parameters: 782 . snes - the SNES context 783 . x - input vector 784 785 Output Parameters: 786 . A - Hessian matrix 787 . B - optional preconditioning matrix 788 . flag - flag indicating matrix structure 789 790 Notes: 791 Most users should not need to explicitly call this routine, as it 792 is used internally within the nonlinear solvers. 793 794 See SLESSetOperators() for information about setting the flag parameter. 795 796 SNESComputeHessian() is valid only for 797 SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 798 SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian(). 799 800 .keywords: SNES, compute, Hessian, matrix 801 802 .seealso: SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(), 803 SNESComputeMinimizationFunction() 804 @*/ 805 int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag) 806 { 807 int ierr; 808 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 809 SETERRQ(1,"SNESComputeHessian:For SNES_UNCONSTRAINED_MINIMIZATION only"); 810 if (!snes->computejacobian) return 0; 811 PLogEventBegin(SNES_HessianEval,snes,x,*A,*B); 812 *flag = DIFFERENT_NONZERO_PATTERN; 813 ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr); 814 PLogEventEnd(SNES_HessianEval,snes,x,*A,*B); 815 /* make sure user returned a correct Jacobian and preconditioner */ 816 PetscValidHeaderSpecific(*A,MAT_COOKIE); 817 PetscValidHeaderSpecific(*B,MAT_COOKIE); 818 return 0; 819 } 820 821 /*@C 822 SNESSetJacobian - Sets the function to compute Jacobian as well as the 823 location to store it. 824 825 Input Parameters: 826 . snes - the SNES context 827 . A - Jacobian matrix 828 . B - preconditioner matrix (usually same as the Jacobian) 829 . func - Jacobian evaluation routine 830 . ctx - optional user-defined context for private data for the 831 Jacobian evaluation routine (may be null) 832 833 Calling sequence of func: 834 . func (SNES,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 835 836 . x - input vector 837 . A - Jacobian matrix 838 . B - preconditioner matrix, usually the same as A 839 . flag - flag indicating information about the preconditioner matrix 840 structure (same as flag in SLESSetOperators()) 841 . ctx - optional user-defined Jacobian context 842 843 Notes: 844 See SLESSetOperators() for information about setting the flag input 845 parameter in the routine func(). Be sure to read this information! 846 847 The routine func() takes Mat * as the matrix arguments rather than Mat. 848 This allows the Jacobian evaluation routine to replace A and/or B with a 849 completely new new matrix structure (not just different matrix elements) 850 when appropriate, for instance, if the nonzero structure is changing 851 throughout the global iterations. 852 853 .keywords: SNES, nonlinear, set, Jacobian, matrix 854 855 .seealso: SLESSetOperators(), SNESSetFunction() 856 @*/ 857 int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 858 MatStructure*,void*),void *ctx) 859 { 860 PetscValidHeaderSpecific(snes,SNES_COOKIE); 861 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 862 SETERRQ(1,"SNESSetJacobian:For SNES_NONLINEAR_EQUATIONS only"); 863 snes->computejacobian = func; 864 snes->jacP = ctx; 865 snes->jacobian = A; 866 snes->jacobian_pre = B; 867 return 0; 868 } 869 870 /*@ 871 SNESGetJacobian - Returns the Jacobian matrix and optionally the user 872 provided context for evaluating the Jacobian. 873 874 Input Parameter: 875 . snes - the nonlinear solver context 876 877 Output Parameters: 878 . A - location to stash Jacobian matrix (or PETSC_NULL) 879 . B - location to stash preconditioner matrix (or PETSC_NULL) 880 . ctx - location to stash Jacobian ctx (or PETSC_NULL) 881 882 .seealso: SNESSetJacobian(), SNESComputeJacobian() 883 @*/ 884 int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx) 885 { 886 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 887 SETERRQ(1,"SNESSetJacobian:For SNES_NONLINEAR_EQUATIONS only"); 888 if (A) *A = snes->jacobian; 889 if (B) *B = snes->jacobian_pre; 890 if (ctx) *ctx = snes->jacP; 891 return 0; 892 } 893 894 /*@C 895 SNESSetHessian - Sets the function to compute Hessian as well as the 896 location to store it. 897 898 Input Parameters: 899 . snes - the SNES context 900 . A - Hessian matrix 901 . B - preconditioner matrix (usually same as the Hessian) 902 . func - Jacobian evaluation routine 903 . ctx - optional user-defined context for private data for the 904 Hessian evaluation routine (may be null) 905 906 Calling sequence of func: 907 . func (SNES,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 908 909 . x - input vector 910 . A - Hessian matrix 911 . B - preconditioner matrix, usually the same as A 912 . flag - flag indicating information about the preconditioner matrix 913 structure (same as flag in SLESSetOperators()) 914 . ctx - optional user-defined Hessian context 915 916 Notes: 917 See SLESSetOperators() for information about setting the flag input 918 parameter in the routine func(). Be sure to read this information! 919 920 The function func() takes Mat * as the matrix arguments rather than Mat. 921 This allows the Hessian evaluation routine to replace A and/or B with a 922 completely new new matrix structure (not just different matrix elements) 923 when appropriate, for instance, if the nonzero structure is changing 924 throughout the global iterations. 925 926 .keywords: SNES, nonlinear, set, Hessian, matrix 927 928 .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators() 929 @*/ 930 int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 931 MatStructure*,void*),void *ctx) 932 { 933 PetscValidHeaderSpecific(snes,SNES_COOKIE); 934 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 935 SETERRQ(1,"SNESSetHessian:For SNES_UNCONSTRAINED_MINIMIZATION only"); 936 snes->computejacobian = func; 937 snes->jacP = ctx; 938 snes->jacobian = A; 939 snes->jacobian_pre = B; 940 return 0; 941 } 942 943 /*@ 944 SNESGetHessian - Returns the Hessian matrix and optionally the user 945 provided context for evaluating the Hessian. 946 947 Input Parameter: 948 . snes - the nonlinear solver context 949 950 Output Parameters: 951 . A - location to stash Hessian matrix (or PETSC_NULL) 952 . B - location to stash preconditioner matrix (or PETSC_NULL) 953 . ctx - location to stash Hessian ctx (or PETSC_NULL) 954 955 .seealso: SNESSetHessian(), SNESComputeHessian() 956 @*/ 957 int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx) 958 { 959 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 960 SETERRQ(1,"SNESSetHessian:For SNES_UNCONSTRAINED_MINIMIZATION only"); 961 if (A) *A = snes->jacobian; 962 if (B) *B = snes->jacobian_pre; 963 if (ctx) *ctx = snes->jacP; 964 return 0; 965 } 966 967 /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 968 969 /*@ 970 SNESSetUp - Sets up the internal data structures for the later use 971 of a nonlinear solver. 972 973 Input Parameter: 974 . snes - the SNES context 975 . x - the solution vector 976 977 Notes: 978 For basic use of the SNES solvers the user need not explicitly call 979 SNESSetUp(), since these actions will automatically occur during 980 the call to SNESSolve(). However, if one wishes to control this 981 phase separately, SNESSetUp() should be called after SNESCreate() 982 and optional routines of the form SNESSetXXX(), but before SNESSolve(). 983 984 .keywords: SNES, nonlinear, setup 985 986 .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 987 @*/ 988 int SNESSetUp(SNES snes,Vec x) 989 { 990 int ierr, flg; 991 PetscValidHeaderSpecific(snes,SNES_COOKIE); 992 PetscValidHeaderSpecific(x,VEC_COOKIE); 993 snes->vec_sol = snes->vec_sol_always = x; 994 995 ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg); CHKERRQ(ierr); 996 if (flg) { 997 Mat J; 998 ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr); 999 PLogObjectParent(snes,J); 1000 snes->mfshell = J; 1001 if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 1002 ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 1003 PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n"); 1004 } 1005 else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 1006 ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 1007 PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n"); 1008 } else SETERRQ(1,"SNESSetUp:Method class doesn't support matrix-free option"); 1009 } 1010 if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) { 1011 if (!snes->vec_func) SETERRQ(1,"SNESSetUp:Must call SNESSetFunction() first"); 1012 if (!snes->computefunction) SETERRQ(1,"SNESSetUp:Must call SNESSetFunction() first"); 1013 if (!snes->jacobian) SETERRQ(1,"SNESSetUp:Must call SNESSetJacobian() first"); 1014 if (snes->vec_func == snes->vec_sol) SETERRQ(1,"SNESSetUp:Solution vector cannot be function vector"); 1015 1016 /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 1017 if (snes->ksp_ewconv && snes->type != SNES_EQ_TR) { 1018 SLES sles; KSP ksp; 1019 ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 1020 ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr); 1021 ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private, 1022 (void *)snes); CHKERRQ(ierr); 1023 } 1024 } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) { 1025 if (!snes->vec_func) SETERRQ(1,"SNESSetUp:Must call SNESSetGradient() first"); 1026 if (!snes->computefunction) SETERRQ(1,"SNESSetUp:Must call SNESSetGradient() first"); 1027 if (!snes->computeumfunction) 1028 SETERRQ(1,"SNESSetUp:Must call SNESSetMinimizationFunction() first"); 1029 if (!snes->jacobian) SETERRQ(1,"SNESSetUp:Must call SNESSetHessian() first"); 1030 } else SETERRQ(1,"SNESSetUp:Unknown method class"); 1031 if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);} 1032 snes->setup_called = 1; 1033 return 0; 1034 } 1035 1036 /*@C 1037 SNESDestroy - Destroys the nonlinear solver context that was created 1038 with SNESCreate(). 1039 1040 Input Parameter: 1041 . snes - the SNES context 1042 1043 .keywords: SNES, nonlinear, destroy 1044 1045 .seealso: SNESCreate(), SNESSolve() 1046 @*/ 1047 int SNESDestroy(SNES snes) 1048 { 1049 int ierr; 1050 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1051 ierr = (*(snes)->destroy)((PetscObject)snes); CHKERRQ(ierr); 1052 if (snes->kspconvctx) PetscFree(snes->kspconvctx); 1053 if (snes->mfshell) MatDestroy(snes->mfshell); 1054 ierr = SLESDestroy(snes->sles); CHKERRQ(ierr); 1055 if (snes->xmonitor) SNESLGMonitorDestroy(snes->xmonitor); 1056 if (snes->vwork) VecDestroyVecs(snes->vwork,snes->nvwork); 1057 PLogObjectDestroy((PetscObject)snes); 1058 PetscHeaderDestroy((PetscObject)snes); 1059 return 0; 1060 } 1061 1062 /* ----------- Routines to set solver parameters ---------- */ 1063 1064 /*@ 1065 SNESSetTolerances - Sets various parameters used in convergence tests. 1066 1067 Input Parameters: 1068 . snes - the SNES context 1069 . atol - tolerance 1070 1071 Options Database Key: 1072 $ -snes_atol <atol> - absolute convergence tolerance 1073 $ -snes_rtol <rtol> - relative convergence tolerance 1074 $ -snes_stol <stol> - convergence tolerance in terms of 1075 $ the norm of the change in the solution between steps 1076 $ -snes_max_it <maxit> - maximum number of iterations 1077 $ -snes_max_funcs <maxf> - maximum number of function evaluations 1078 1079 Notes: 1080 The default maximum number of iterations is 50. 1081 The default maximum number of function evaluations is 1000. 1082 1083 .keywords: SNES, nonlinear, set, absolute, convergence, tolerance 1084 1085 .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance() 1086 @*/ 1087 int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf) 1088 { 1089 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1090 if (atol != PETSC_DEFAULT) snes->atol = atol; 1091 if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1092 if (stol != PETSC_DEFAULT) snes->xtol = stol; 1093 if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1094 if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 1095 return 0; 1096 } 1097 1098 /*@ 1099 SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 1100 1101 Input Parameters: 1102 . snes - the SNES context 1103 . tol - tolerance 1104 1105 Options Database Key: 1106 $ -snes_trtol <tol> 1107 1108 .keywords: SNES, nonlinear, set, trust region, tolerance 1109 1110 .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance() 1111 @*/ 1112 int SNESSetTrustRegionTolerance(SNES snes,double tol) 1113 { 1114 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1115 snes->deltatol = tol; 1116 return 0; 1117 } 1118 1119 /*@ 1120 SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance 1121 for unconstrained minimization solvers. 1122 1123 Input Parameters: 1124 . snes - the SNES context 1125 . ftol - minimum function tolerance 1126 1127 Options Database Key: 1128 $ -snes_fmin <ftol> 1129 1130 Note: 1131 SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION 1132 methods only. 1133 1134 .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance 1135 1136 .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance() 1137 @*/ 1138 int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol) 1139 { 1140 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1141 snes->fmin = ftol; 1142 return 0; 1143 } 1144 1145 /* ------------ Routines to set performance monitoring options ----------- */ 1146 1147 /*@C 1148 SNESSetMonitor - Sets the function that is to be used at every 1149 iteration of the nonlinear solver to display the iteration's 1150 progress. 1151 1152 Input Parameters: 1153 . snes - the SNES context 1154 . func - monitoring routine 1155 . mctx - optional user-defined context for private data for the 1156 monitor routine (may be null) 1157 1158 Calling sequence of func: 1159 int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx) 1160 1161 $ snes - the SNES context 1162 $ its - iteration number 1163 $ mctx - optional monitoring context 1164 $ 1165 $ SNES_NONLINEAR_EQUATIONS methods: 1166 $ norm - 2-norm function value (may be estimated) 1167 $ 1168 $ SNES_UNCONSTRAINED_MINIMIZATION methods: 1169 $ norm - 2-norm gradient value (may be estimated) 1170 1171 .keywords: SNES, nonlinear, set, monitor 1172 1173 .seealso: SNESDefaultMonitor() 1174 @*/ 1175 int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx ) 1176 { 1177 snes->monitor = func; 1178 snes->monP = mctx; 1179 return 0; 1180 } 1181 1182 /*@C 1183 SNESSetConvergenceTest - Sets the function that is to be used 1184 to test for convergence of the nonlinear iterative solution. 1185 1186 Input Parameters: 1187 . snes - the SNES context 1188 . func - routine to test for convergence 1189 . cctx - optional context for private data for the convergence routine 1190 (may be null) 1191 1192 Calling sequence of func: 1193 int func (SNES snes,double xnorm,double gnorm, 1194 double f,void *cctx) 1195 1196 $ snes - the SNES context 1197 $ cctx - optional convergence context 1198 $ xnorm - 2-norm of current iterate 1199 $ 1200 $ SNES_NONLINEAR_EQUATIONS methods: 1201 $ gnorm - 2-norm of current step 1202 $ f - 2-norm of function 1203 $ 1204 $ SNES_UNCONSTRAINED_MINIMIZATION methods: 1205 $ gnorm - 2-norm of current gradient 1206 $ f - function value 1207 1208 .keywords: SNES, nonlinear, set, convergence, test 1209 1210 .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(), 1211 SNESConverged_UM_LS(), SNESConverged_UM_TR() 1212 @*/ 1213 int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx) 1214 { 1215 (snes)->converged = func; 1216 (snes)->cnvP = cctx; 1217 return 0; 1218 } 1219 1220 /* 1221 SNESScaleStep_Private - Scales a step so that its length is less than the 1222 positive parameter delta. 1223 1224 Input Parameters: 1225 . snes - the SNES context 1226 . y - approximate solution of linear system 1227 . fnorm - 2-norm of current function 1228 . delta - trust region size 1229 1230 Output Parameters: 1231 . gpnorm - predicted function norm at the new point, assuming local 1232 linearization. The value is zero if the step lies within the trust 1233 region, and exceeds zero otherwise. 1234 . ynorm - 2-norm of the step 1235 1236 Note: 1237 For non-trust region methods such as SNES_EQ_LS, the parameter delta 1238 is set to be the maximum allowable step size. 1239 1240 .keywords: SNES, nonlinear, scale, step 1241 */ 1242 int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta, 1243 double *gpnorm,double *ynorm) 1244 { 1245 double norm; 1246 Scalar cnorm; 1247 VecNorm(y,NORM_2, &norm ); 1248 if (norm > *delta) { 1249 norm = *delta/norm; 1250 *gpnorm = (1.0 - norm)*(*fnorm); 1251 cnorm = norm; 1252 VecScale( &cnorm, y ); 1253 *ynorm = *delta; 1254 } else { 1255 *gpnorm = 0.0; 1256 *ynorm = norm; 1257 } 1258 return 0; 1259 } 1260 1261 /*@ 1262 SNESSolve - Solves a nonlinear system. Call SNESSolve after calling 1263 SNESCreate() and optional routines of the form SNESSetXXX(). 1264 1265 Input Parameter: 1266 . snes - the SNES context 1267 . x - the solution vector 1268 1269 Output Parameter: 1270 its - number of iterations until termination 1271 1272 Note: 1273 The user should initialize the vector, x, with the initial guess 1274 for the nonlinear solve prior to calling SNESSolve. In particular, 1275 to employ an initial guess of zero, the user should explicitly set 1276 this vector to zero by calling VecSet(). 1277 1278 .keywords: SNES, nonlinear, solve 1279 1280 .seealso: SNESCreate(), SNESDestroy() 1281 @*/ 1282 int SNESSolve(SNES snes,Vec x,int *its) 1283 { 1284 int ierr, flg; 1285 1286 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1287 PetscValidIntPointer(its); 1288 if (!snes->setup_called) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);} 1289 else {snes->vec_sol = snes->vec_sol_always = x;} 1290 PLogEventBegin(SNES_Solve,snes,0,0,0); 1291 ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr); 1292 PLogEventEnd(SNES_Solve,snes,0,0,0); 1293 ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr); 1294 if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); } 1295 return 0; 1296 } 1297 1298 /* --------- Internal routines for SNES Package --------- */ 1299 static NRList *__SNESList = 0; 1300 1301 /*@ 1302 SNESSetType - Sets the method for the nonlinear solver. 1303 1304 Input Parameters: 1305 . snes - the SNES context 1306 . method - a known method 1307 1308 Notes: 1309 See "petsc/include/snes.h" for available methods (for instance) 1310 $ Systems of nonlinear equations: 1311 $ SNES_EQ_LS - Newton's method with line search 1312 $ SNES_EQ_TR - Newton's method with trust region 1313 $ Unconstrained minimization: 1314 $ SNES_UM_TR - Newton's method with trust region 1315 $ SNES_UM_LS - Newton's method with line search 1316 1317 Options Database Command: 1318 $ -snes_type <method> 1319 $ Use -help for a list of available methods 1320 $ (for instance, ls or tr) 1321 1322 .keysords: SNES, set, method 1323 @*/ 1324 int SNESSetType(SNES snes,SNESType method) 1325 { 1326 int (*r)(SNES); 1327 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1328 if (snes->type == (int) method) return 0; 1329 1330 /* Get the function pointers for the iterative method requested */ 1331 if (!__SNESList) {SNESRegisterAll();} 1332 if (!__SNESList) {SETERRQ(1,"SNESSetType:Could not get methods");} 1333 r = (int (*)(SNES))NRFindRoutine( __SNESList, (int)method, (char *)0 ); 1334 if (!r) {SETERRQ(1,"SNESSetType:Unknown method");} 1335 if (snes->data) PetscFree(snes->data); 1336 snes->set_method_called = 1; 1337 return (*r)(snes); 1338 } 1339 1340 /* --------------------------------------------------------------------- */ 1341 /*@C 1342 SNESRegister - Adds the method to the nonlinear solver package, given 1343 a function pointer and a nonlinear solver name of the type SNESType. 1344 1345 Input Parameters: 1346 . name - for instance SNES_EQ_LS, SNES_EQ_TR, ... 1347 . sname - corresponding string for name 1348 . create - routine to create method context 1349 1350 .keywords: SNES, nonlinear, register 1351 1352 .seealso: SNESRegisterAll(), SNESRegisterDestroy() 1353 @*/ 1354 int SNESRegister(int name, char *sname, int (*create)(SNES)) 1355 { 1356 int ierr; 1357 if (!__SNESList) {ierr = NRCreate(&__SNESList); CHKERRQ(ierr);} 1358 NRRegister( __SNESList, name, sname, (int (*)(void*))create ); 1359 return 0; 1360 } 1361 /* --------------------------------------------------------------------- */ 1362 /*@C 1363 SNESRegisterDestroy - Frees the list of nonlinear solvers that were 1364 registered by SNESRegister(). 1365 1366 .keywords: SNES, nonlinear, register, destroy 1367 1368 .seealso: SNESRegisterAll(), SNESRegisterAll() 1369 @*/ 1370 int SNESRegisterDestroy() 1371 { 1372 if (__SNESList) { 1373 NRDestroy( __SNESList ); 1374 __SNESList = 0; 1375 } 1376 return 0; 1377 } 1378 1379 /* 1380 SNESGetTypeFromOptions_Private - Sets the selected method from the 1381 options database. 1382 1383 Input Parameter: 1384 . ctx - the SNES context 1385 1386 Output Parameter: 1387 . method - solver method 1388 1389 Returns: 1390 Returns 1 if the method is found; 0 otherwise. 1391 1392 Options Database Key: 1393 $ -snes_type method 1394 */ 1395 int SNESGetTypeFromOptions_Private(SNES ctx,SNESType *method,int *flg) 1396 { 1397 int ierr; 1398 char sbuf[50]; 1399 1400 ierr = OptionsGetString(ctx->prefix,"-snes_type",sbuf,50,flg);CHKERRQ(ierr); 1401 if (*flg) { 1402 if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);} 1403 *method = (SNESType)NRFindID( __SNESList, sbuf ); 1404 } 1405 return 0; 1406 } 1407 1408 /*@C 1409 SNESGetType - Gets the SNES method type and name (as a string). 1410 1411 Input Parameter: 1412 . snes - nonlinear solver context 1413 1414 Output Parameter: 1415 . method - SNES method (or use PETSC_NULL) 1416 . name - name of SNES method (or use PETSC_NULL) 1417 1418 .keywords: SNES, nonlinear, get, method, name 1419 @*/ 1420 int SNESGetType(SNES snes, SNESType *method,char **name) 1421 { 1422 int ierr; 1423 if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);} 1424 if (method) *method = (SNESType) snes->type; 1425 if (name) *name = NRFindName( __SNESList, (int) snes->type ); 1426 return 0; 1427 } 1428 1429 #include <stdio.h> 1430 /* 1431 SNESPrintTypes_Private - Prints the SNES methods available from the 1432 options database. 1433 1434 Input Parameters: 1435 . comm - communicator (usually MPI_COMM_WORLD) 1436 . prefix - prefix (usually "-") 1437 . name - the options database name (by default "snes_type") 1438 */ 1439 int SNESPrintTypes_Private(MPI_Comm comm,char* prefix,char *name) 1440 { 1441 FuncList *entry; 1442 if (!__SNESList) {SNESRegisterAll();} 1443 entry = __SNESList->head; 1444 PetscPrintf(comm," %s%s (one of)",prefix,name); 1445 while (entry) { 1446 PetscPrintf(comm," %s",entry->name); 1447 entry = entry->next; 1448 } 1449 PetscPrintf(comm,"\n"); 1450 return 0; 1451 } 1452 1453 /*@C 1454 SNESGetSolution - Returns the vector where the approximate solution is 1455 stored. 1456 1457 Input Parameter: 1458 . snes - the SNES context 1459 1460 Output Parameter: 1461 . x - the solution 1462 1463 .keywords: SNES, nonlinear, get, solution 1464 1465 .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate() 1466 @*/ 1467 int SNESGetSolution(SNES snes,Vec *x) 1468 { 1469 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1470 *x = snes->vec_sol_always; 1471 return 0; 1472 } 1473 1474 /*@C 1475 SNESGetSolutionUpdate - Returns the vector where the solution update is 1476 stored. 1477 1478 Input Parameter: 1479 . snes - the SNES context 1480 1481 Output Parameter: 1482 . x - the solution update 1483 1484 Notes: 1485 This vector is implementation dependent. 1486 1487 .keywords: SNES, nonlinear, get, solution, update 1488 1489 .seealso: SNESGetSolution(), SNESGetFunction 1490 @*/ 1491 int SNESGetSolutionUpdate(SNES snes,Vec *x) 1492 { 1493 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1494 *x = snes->vec_sol_update_always; 1495 return 0; 1496 } 1497 1498 /*@C 1499 SNESGetFunction - Returns the vector where the function is stored. 1500 1501 Input Parameter: 1502 . snes - the SNES context 1503 1504 Output Parameter: 1505 . r - the function 1506 1507 Notes: 1508 SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only 1509 Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 1510 SNESGetMinimizationFunction() and SNESGetGradient(); 1511 1512 .keywords: SNES, nonlinear, get, function 1513 1514 .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(), 1515 SNESGetGradient() 1516 @*/ 1517 int SNESGetFunction(SNES snes,Vec *r) 1518 { 1519 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1520 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1, 1521 "SNESGetFunction:For SNES_NONLINEAR_EQUATIONS only"); 1522 *r = snes->vec_func_always; 1523 return 0; 1524 } 1525 1526 /*@C 1527 SNESGetGradient - Returns the vector where the gradient is stored. 1528 1529 Input Parameter: 1530 . snes - the SNES context 1531 1532 Output Parameter: 1533 . r - the gradient 1534 1535 Notes: 1536 SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods 1537 only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 1538 SNESGetFunction(). 1539 1540 .keywords: SNES, nonlinear, get, gradient 1541 1542 .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction() 1543 @*/ 1544 int SNESGetGradient(SNES snes,Vec *r) 1545 { 1546 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1547 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1, 1548 "SNESGetGradient:For SNES_UNCONSTRAINED_MINIMIZATION only"); 1549 *r = snes->vec_func_always; 1550 return 0; 1551 } 1552 1553 /*@ 1554 SNESGetMinimizationFunction - Returns the scalar function value for 1555 unconstrained minimization problems. 1556 1557 Input Parameter: 1558 . snes - the SNES context 1559 1560 Output Parameter: 1561 . r - the function 1562 1563 Notes: 1564 SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 1565 methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 1566 SNESGetFunction(). 1567 1568 .keywords: SNES, nonlinear, get, function 1569 1570 .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction() 1571 @*/ 1572 int SNESGetMinimizationFunction(SNES snes,double *r) 1573 { 1574 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1575 PetscValidScalarPointer(r); 1576 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1, 1577 "SNESGetMinimizationFunction:For SNES_UNCONSTRAINED_MINIMIZATION only"); 1578 *r = snes->fc; 1579 return 0; 1580 } 1581 1582 /*@C 1583 SNESSetOptionsPrefix - Sets the prefix used for searching for all 1584 SNES options in the database. 1585 1586 Input Parameter: 1587 . snes - the SNES context 1588 . prefix - the prefix to prepend to all option names 1589 1590 .keywords: SNES, set, options, prefix, database 1591 1592 .seealso: SNESSetFromOptions() 1593 @*/ 1594 int SNESSetOptionsPrefix(SNES snes,char *prefix) 1595 { 1596 int ierr; 1597 1598 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1599 ierr = PetscObjectSetPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 1600 ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 1601 return 0; 1602 } 1603 1604 /*@C 1605 SNESAppendOptionsPrefix - Append to the prefix used for searching for all 1606 SNES options in the database. 1607 1608 Input Parameter: 1609 . snes - the SNES context 1610 . prefix - the prefix to prepend to all option names 1611 1612 .keywords: SNES, append, options, prefix, database 1613 1614 .seealso: SNESGetOptionsPrefix() 1615 @*/ 1616 int SNESAppendOptionsPrefix(SNES snes,char *prefix) 1617 { 1618 int ierr; 1619 1620 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1621 ierr = PetscObjectAppendPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 1622 ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 1623 return 0; 1624 } 1625 1626 /*@ 1627 SNESGetOptionsPrefix - Sets the prefix used for searching for all 1628 SNES options in the database. 1629 1630 Input Parameter: 1631 . snes - the SNES context 1632 1633 Output Parameter: 1634 . prefix - pointer to the prefix string used 1635 1636 .keywords: SNES, get, options, prefix, database 1637 1638 .seealso: SNESAppendOptionsPrefix() 1639 @*/ 1640 int SNESGetOptionsPrefix(SNES snes,char **prefix) 1641 { 1642 int ierr; 1643 1644 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1645 ierr = PetscObjectGetPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 1646 return 0; 1647 } 1648 1649 1650 1651 1652 1653