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