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