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