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