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