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