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