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