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