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