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