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