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