1 2 #ifdef PETSC_RCS_HEADER 3 static char vcid[] = "$Id: snes.c,v 1.144 1998/04/13 17:55:33 bsmith Exp curfman $"; 4 #endif 5 6 #include "src/snes/snesimpl.h" /*I "snes.h" I*/ 7 #include "src/sys/nreg.h" 8 #include "pinclude/pviewer.h" 9 #include <math.h> 10 11 int SNESRegisterAllCalled = 0; 12 DLList SNESList = 0; 13 14 15 #undef __FUNC__ 16 #define __FUNC__ "SNESView" 17 /*@ 18 SNESView - Prints the SNES data structure. 19 20 Input Parameters: 21 . SNES - the SNES context 22 . viewer - visualization context 23 24 Collective on SNES, unless Viewer is VIEWER_STDOUT_SELF 25 26 Options Database Key: 27 . -snes_view - Calls SNESView() at end of SNESSolve() 28 29 Notes: 30 The available visualization contexts include 31 . VIEWER_STDOUT_SELF - standard output (default) 32 . VIEWER_STDOUT_WORLD - synchronized standard 33 output where only the first processor opens 34 the file. All other processors send their 35 data to the first processor to print. 36 37 The user can open alternative vistualization contexts with 38 $ ViewerFileOpenASCII() - output to a specified file 39 40 .keywords: SNES, view 41 42 .seealso: ViewerFileOpenASCII() 43 @*/ 44 int SNESView(SNES snes,Viewer viewer) 45 { 46 SNES_KSP_EW_ConvCtx *kctx; 47 FILE *fd; 48 int ierr; 49 SLES sles; 50 char *method; 51 ViewerType vtype; 52 53 PetscFunctionBegin; 54 PetscValidHeaderSpecific(snes,SNES_COOKIE); 55 if (viewer) {PetscValidHeader(viewer);} 56 else { viewer = VIEWER_STDOUT_SELF; } 57 58 ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 59 if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 60 ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 61 PetscFPrintf(snes->comm,fd,"SNES Object:\n"); 62 SNESGetType(snes,&method); 63 PetscFPrintf(snes->comm,fd," method: %s\n",method); 64 if (snes->view) (*snes->view)(snes,viewer); 65 PetscFPrintf(snes->comm,fd, 66 " maximum iterations=%d, maximum function evaluations=%d\n", 67 snes->max_its,snes->max_funcs); 68 PetscFPrintf(snes->comm,fd, 69 " tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n", 70 snes->rtol, snes->atol, snes->trunctol, snes->xtol); 71 PetscFPrintf(snes->comm,fd, 72 " total number of linear solver iterations=%d\n",snes->linear_its); 73 PetscFPrintf(snes->comm,fd, 74 " total number of function evaluations=%d\n",snes->nfuncs); 75 if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) 76 PetscFPrintf(snes->comm,fd," min function tolerance=%g\n",snes->fmin); 77 if (snes->ksp_ewconv) { 78 kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 79 if (kctx) { 80 PetscFPrintf(snes->comm,fd, 81 " Eisenstat-Walker computation of KSP relative tolerance (version %d)\n", 82 kctx->version); 83 PetscFPrintf(snes->comm,fd, 84 " rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0, 85 kctx->rtol_max,kctx->threshold); 86 PetscFPrintf(snes->comm,fd," gamma=%g, alpha=%g, alpha2=%g\n", 87 kctx->gamma,kctx->alpha,kctx->alpha2); 88 } 89 } 90 } else if (vtype == STRING_VIEWER) { 91 SNESGetType(snes,&method); 92 ViewerStringSPrintf(viewer," %-3.3s",method); 93 } 94 SNESGetSLES(snes,&sles); 95 ierr = SLESView(sles,viewer); CHKERRQ(ierr); 96 PetscFunctionReturn(0); 97 } 98 99 /* 100 We retain a list of functions that also take SNES command 101 line options. These are called at the end SNESSetFromOptions() 102 */ 103 #define MAXSETFROMOPTIONS 5 104 static int numberofsetfromoptions; 105 static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES); 106 107 #undef __FUNC__ 108 #define __FUNC__ "SNESAddOptionsChecker" 109 /*@ 110 SNESAddOptionsChecker - Adds an additional function to check for SNES options. 111 112 Input Parameter: 113 . snescheck - function that checks for options 114 115 Not Collective 116 117 .seealso: SNESSetFromOptions() 118 @*/ 119 int SNESAddOptionsChecker(int (*snescheck)(SNES) ) 120 { 121 PetscFunctionBegin; 122 if (numberofsetfromoptions >= MAXSETFROMOPTIONS) { 123 SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many options checkers, only 5 allowed"); 124 } 125 126 othersetfromoptions[numberofsetfromoptions++] = snescheck; 127 PetscFunctionReturn(0); 128 } 129 130 #undef __FUNC__ 131 #define __FUNC__ "SNESSetFromOptions" 132 /*@ 133 SNESSetFromOptions - Sets various SNES and SLES parameters from user options. 134 135 Input Parameter: 136 . snes - the SNES context 137 138 Collective on SNES 139 Notes: 140 To see all options, run your program with the -help option or consult 141 the users manual. 142 143 .keywords: SNES, nonlinear, set, options, database 144 145 .seealso: SNESPrintHelp(), SNESSetOptionsPrefix() 146 @*/ 147 int SNESSetFromOptions(SNES snes) 148 { 149 char method[256]; 150 double tmp; 151 SLES sles; 152 int ierr, flg,i,loc[4],nmax = 4; 153 int version = PETSC_DEFAULT; 154 double rtol_0 = PETSC_DEFAULT; 155 double rtol_max = PETSC_DEFAULT; 156 double gamma2 = PETSC_DEFAULT; 157 double alpha = PETSC_DEFAULT; 158 double alpha2 = PETSC_DEFAULT; 159 double threshold = PETSC_DEFAULT; 160 161 PetscFunctionBegin; 162 loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300; 163 164 PetscValidHeaderSpecific(snes,SNES_COOKIE); 165 if (snes->setupcalled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call prior to SNESSetUp"); 166 167 if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 168 ierr = OptionsGetString(snes->prefix,"-snes_type",method,256,&flg); 169 if (flg) { 170 ierr = SNESSetType(snes,(SNESType) method); CHKERRQ(ierr); 171 } 172 ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg); CHKERRQ(ierr); 173 if (flg) { 174 ierr = SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 175 } 176 ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg); CHKERRQ(ierr); 177 if (flg) { 178 ierr = SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 179 } 180 ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg); CHKERRQ(ierr); 181 if (flg) { 182 ierr = SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 183 } 184 ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, &flg); CHKERRQ(ierr); 185 ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, &flg);CHKERRQ(ierr); 186 ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg); CHKERRQ(ierr); 187 if (flg) { ierr = SNESSetTrustRegionTolerance(snes,tmp); CHKERRQ(ierr); } 188 ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg); CHKERRQ(ierr); 189 if (flg) { ierr = SNESSetMinimizationFunctionTolerance(snes,tmp); CHKERRQ(ierr);} 190 ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg); CHKERRQ(ierr); 191 if (flg) { snes->ksp_ewconv = 1; } 192 ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,&flg); CHKERRQ(ierr); 193 ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,&flg); CHKERRQ(ierr); 194 ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,&flg);CHKERRQ(ierr); 195 ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,&flg); CHKERRQ(ierr); 196 ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,&flg); CHKERRQ(ierr); 197 ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,&flg); CHKERRQ(ierr); 198 ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,&flg);CHKERRQ(ierr); 199 200 ierr = OptionsHasName(snes->prefix,"-snes_no_convergence_test",&flg); CHKERRQ(ierr); 201 if (flg) {snes->converged = 0;} 202 203 ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha, 204 alpha2,threshold); CHKERRQ(ierr); 205 ierr = OptionsHasName(snes->prefix,"-snes_cancelmonitors",&flg); CHKERRQ(ierr); 206 if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);} 207 ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg); CHKERRQ(ierr); 208 if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0);CHKERRQ(ierr);} 209 ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg); CHKERRQ(ierr); 210 if (flg) { 211 ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0); CHKERRQ(ierr); 212 } 213 ierr = OptionsGetIntArray(snes->prefix,"-snes_xmonitor",loc,&nmax,&flg);CHKERRQ(ierr); 214 if (flg) { 215 int rank = 0; 216 DrawLG lg; 217 MPI_Initialized(&rank); 218 if (rank) MPI_Comm_rank(snes->comm,&rank); 219 if (!rank) { 220 ierr = SNESLGMonitorCreate(0,0,loc[0],loc[1],loc[2],loc[3],&lg); CHKERRQ(ierr); 221 ierr = SNESSetMonitor(snes,SNESLGMonitor,(void *)lg); CHKERRQ(ierr); 222 snes->xmonitor = lg; 223 PLogObjectParent(snes,lg); 224 } 225 } 226 ierr = OptionsHasName(snes->prefix,"-snes_fd", &flg); CHKERRQ(ierr); 227 if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) { 228 ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre, 229 SNESDefaultComputeJacobian,snes->funP); CHKERRQ(ierr); 230 PLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Jacobian matrix\n"); 231 } else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 232 ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre, 233 SNESDefaultComputeHessian,snes->funP); CHKERRQ(ierr); 234 PLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Hessian matrix\n"); 235 } 236 237 for ( i=0; i<numberofsetfromoptions; i++ ) { 238 ierr = (*othersetfromoptions[i])(snes); CHKERRQ(ierr); 239 } 240 241 ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 242 ierr = SLESSetFromOptions(sles); CHKERRQ(ierr); 243 /* 244 Since the private setfromoptions requires the type to all ready have 245 been set we make sure a type is set by this time 246 */ 247 if (!snes->type_name) { 248 if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 249 ierr = SNESSetType(snes,SNES_EQ_LS); CHKERRQ(ierr); 250 } else { 251 ierr = SNESSetType(snes,SNES_UM_TR); CHKERRQ(ierr); 252 } 253 } 254 255 ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 256 if (flg) { ierr = SNESPrintHelp(snes); CHKERRQ(ierr);} 257 258 if (snes->setfromoptions) { 259 ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr); 260 } 261 PetscFunctionReturn(0); 262 } 263 264 265 #undef __FUNC__ 266 #define __FUNC__ "SNESSetApplicationContext" 267 /*@ 268 SNESSetApplicationContext - Sets the optional user-defined context for 269 the nonlinear solvers. 270 271 Input Parameters: 272 . snes - the SNES context 273 . usrP - optional user context 274 275 Collective on SNES 276 277 .keywords: SNES, nonlinear, set, application, context 278 279 .seealso: SNESGetApplicationContext() 280 @*/ 281 int SNESSetApplicationContext(SNES snes,void *usrP) 282 { 283 PetscFunctionBegin; 284 PetscValidHeaderSpecific(snes,SNES_COOKIE); 285 snes->user = usrP; 286 PetscFunctionReturn(0); 287 } 288 289 #undef __FUNC__ 290 #define __FUNC__ "SNESGetApplicationContext" 291 /*@C 292 SNESGetApplicationContext - Gets the user-defined context for the 293 nonlinear solvers. 294 295 Input Parameter: 296 . snes - SNES context 297 298 Output Parameter: 299 . usrP - user context 300 301 Not Collective 302 303 .keywords: SNES, nonlinear, get, application, context 304 305 .seealso: SNESSetApplicationContext() 306 @*/ 307 int SNESGetApplicationContext( SNES snes, void **usrP ) 308 { 309 PetscFunctionBegin; 310 PetscValidHeaderSpecific(snes,SNES_COOKIE); 311 *usrP = snes->user; 312 PetscFunctionReturn(0); 313 } 314 315 #undef __FUNC__ 316 #define __FUNC__ "SNESGetIterationNumber" 317 /*@ 318 SNESGetIterationNumber - Gets the current iteration number of the 319 nonlinear solver. 320 321 Input Parameter: 322 . snes - SNES context 323 324 Output Parameter: 325 . iter - iteration number 326 327 Not Collective 328 329 .keywords: SNES, nonlinear, get, iteration, number 330 @*/ 331 int SNESGetIterationNumber(SNES snes,int* iter) 332 { 333 PetscFunctionBegin; 334 PetscValidHeaderSpecific(snes,SNES_COOKIE); 335 PetscValidIntPointer(iter); 336 *iter = snes->iter; 337 PetscFunctionReturn(0); 338 } 339 340 #undef __FUNC__ 341 #define __FUNC__ "SNESGetFunctionNorm" 342 /*@ 343 SNESGetFunctionNorm - Gets the norm of the current function that was set 344 with SNESSSetFunction(). 345 346 Input Parameter: 347 . snes - SNES context 348 349 Output Parameter: 350 . fnorm - 2-norm of function 351 352 Collective on SNES 353 354 Note: 355 SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only. 356 A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is 357 SNESGetGradientNorm(). 358 359 .keywords: SNES, nonlinear, get, function, norm 360 361 .seealso: SNESSetFunction() 362 @*/ 363 int SNESGetFunctionNorm(SNES snes,Scalar *fnorm) 364 { 365 PetscFunctionBegin; 366 PetscValidHeaderSpecific(snes,SNES_COOKIE); 367 PetscValidScalarPointer(fnorm); 368 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 369 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_NONLINEAR_EQUATIONS only"); 370 } 371 *fnorm = snes->norm; 372 PetscFunctionReturn(0); 373 } 374 375 #undef __FUNC__ 376 #define __FUNC__ "SNESGetGradientNorm" 377 /*@ 378 SNESGetGradientNorm - Gets the norm of the current gradient that was set 379 with SNESSSetGradient(). 380 381 Input Parameter: 382 . snes - SNES context 383 384 Output Parameter: 385 . fnorm - 2-norm of gradient 386 387 Collective on SNES 388 389 Note: 390 SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION 391 methods only. A related routine for SNES_NONLINEAR_EQUATIONS methods 392 is SNESGetFunctionNorm(). 393 394 .keywords: SNES, nonlinear, get, gradient, norm 395 396 .seelso: SNESSetGradient() 397 @*/ 398 int SNESGetGradientNorm(SNES snes,Scalar *gnorm) 399 { 400 PetscFunctionBegin; 401 PetscValidHeaderSpecific(snes,SNES_COOKIE); 402 PetscValidScalarPointer(gnorm); 403 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 404 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 405 } 406 *gnorm = snes->norm; 407 PetscFunctionReturn(0); 408 } 409 410 #undef __FUNC__ 411 #define __FUNC__ "SNESGetNumberUnsuccessfulSteps" 412 /*@ 413 SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps 414 attempted by the nonlinear solver. 415 416 Input Parameter: 417 . snes - SNES context 418 419 Output Parameter: 420 . nfails - number of unsuccessful steps attempted 421 422 Not Collective 423 424 Notes: 425 This counter is reset to zero for each successive call to SNESSolve(). 426 427 .keywords: SNES, nonlinear, get, number, unsuccessful, steps 428 @*/ 429 int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails) 430 { 431 PetscFunctionBegin; 432 PetscValidHeaderSpecific(snes,SNES_COOKIE); 433 PetscValidIntPointer(nfails); 434 *nfails = snes->nfailures; 435 PetscFunctionReturn(0); 436 } 437 438 #undef __FUNC__ 439 #define __FUNC__ "SNESGetNumberLinearIterations" 440 /*@ 441 SNESGetNumberLinearIterations - Gets the total number of linear iterations 442 used by the nonlinear solver. 443 444 Input Parameter: 445 . snes - SNES context 446 447 Output Parameter: 448 . lits - number of linear iterations 449 450 Not Collective 451 452 Notes: 453 This counter is reset to zero for each successive call to SNESSolve(). 454 455 .keywords: SNES, nonlinear, get, number, linear, iterations 456 @*/ 457 int SNESGetNumberLinearIterations(SNES snes,int* lits) 458 { 459 PetscFunctionBegin; 460 PetscValidHeaderSpecific(snes,SNES_COOKIE); 461 PetscValidIntPointer(lits); 462 *lits = snes->linear_its; 463 PetscFunctionReturn(0); 464 } 465 466 #undef __FUNC__ 467 #define __FUNC__ "SNESGetSLES" 468 /*@C 469 SNESGetSLES - Returns the SLES context for a SNES solver. 470 471 Input Parameter: 472 . snes - the SNES context 473 474 Output Parameter: 475 . sles - the SLES context 476 477 Not Collective, but if SNES object is parallel, then SLES object is parallel 478 479 Notes: 480 The user can then directly manipulate the SLES context to set various 481 options, etc. Likewise, the user can then extract and manipulate the 482 KSP and PC contexts as well. 483 484 .keywords: SNES, nonlinear, get, SLES, context 485 486 .seealso: SLESGetPC(), SLESGetKSP() 487 @*/ 488 int SNESGetSLES(SNES snes,SLES *sles) 489 { 490 PetscFunctionBegin; 491 PetscValidHeaderSpecific(snes,SNES_COOKIE); 492 *sles = snes->sles; 493 PetscFunctionReturn(0); 494 } 495 496 /* -----------------------------------------------------------*/ 497 #undef __FUNC__ 498 #define __FUNC__ "SNESCreate" 499 /*@C 500 SNESCreate - Creates a nonlinear solver context. 501 502 Input Parameter: 503 . comm - MPI communicator 504 . type - type of method, one of 505 $ SNES_NONLINEAR_EQUATIONS 506 $ (for systems of nonlinear equations) 507 $ SNES_UNCONSTRAINED_MINIMIZATION 508 $ (for unconstrained minimization) 509 510 Output Parameter: 511 . outsnes - the new SNES context 512 513 Collective on MPI_Comm 514 515 Options Database Key: 516 $ -snes_mf - use default matrix-free Jacobian-vector products, 517 $ and no preconditioning matrix 518 $ -snes_mf_operator - use default matrix-free Jacobian-vector 519 $ products, and a user-provided preconditioning matrix 520 $ as set by SNESSetJacobian() 521 $ -snes_fd - use (slow!) finite differences to compute Jacobian 522 523 .keywords: SNES, nonlinear, create, context 524 525 .seealso: SNESSolve(), SNESDestroy() 526 @*/ 527 int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes) 528 { 529 int ierr; 530 SNES snes; 531 SNES_KSP_EW_ConvCtx *kctx; 532 533 PetscFunctionBegin; 534 *outsnes = 0; 535 if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS){ 536 SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"incorrect method type"); 537 } 538 PetscHeaderCreate(snes,_p_SNES,int,SNES_COOKIE,0,comm,SNESDestroy,SNESView); 539 PLogObjectCreate(snes); 540 snes->max_its = 50; 541 snes->max_funcs = 10000; 542 snes->norm = 0.0; 543 if (type == SNES_UNCONSTRAINED_MINIMIZATION) { 544 snes->rtol = 1.e-8; 545 snes->ttol = 0.0; 546 snes->atol = 1.e-10; 547 } else { 548 snes->rtol = 1.e-8; 549 snes->ttol = 0.0; 550 snes->atol = 1.e-50; 551 } 552 snes->xtol = 1.e-8; 553 snes->trunctol = 1.e-12; /* no longer used */ 554 snes->nfuncs = 0; 555 snes->nfailures = 0; 556 snes->linear_its = 0; 557 snes->numbermonitors = 0; 558 snes->data = 0; 559 snes->view = 0; 560 snes->computeumfunction = 0; 561 snes->umfunP = 0; 562 snes->fc = 0; 563 snes->deltatol = 1.e-12; 564 snes->fmin = -1.e30; 565 snes->method_class = type; 566 snes->set_method_called = 0; 567 snes->setupcalled = 0; 568 snes->ksp_ewconv = 0; 569 snes->xmonitor = 0; 570 snes->vwork = 0; 571 snes->nwork = 0; 572 snes->conv_hist_size = 0; 573 snes->conv_act_size = 0; 574 snes->conv_hist = 0; 575 576 /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 577 kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx); 578 PLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx)); 579 snes->kspconvctx = (void*)kctx; 580 kctx->version = 2; 581 kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 582 this was too large for some test cases */ 583 kctx->rtol_last = 0; 584 kctx->rtol_max = .9; 585 kctx->gamma = 1.0; 586 kctx->alpha2 = .5*(1.0 + sqrt(5.0)); 587 kctx->alpha = kctx->alpha2; 588 kctx->threshold = .1; 589 kctx->lresid_last = 0; 590 kctx->norm_last = 0; 591 592 ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr); 593 PLogObjectParent(snes,snes->sles) 594 595 *outsnes = snes; 596 PetscFunctionReturn(0); 597 } 598 599 /* --------------------------------------------------------------- */ 600 #undef __FUNC__ 601 #define __FUNC__ "SNESSetFunction" 602 /*@C 603 SNESSetFunction - Sets the function evaluation routine and function 604 vector for use by the SNES routines in solving systems of nonlinear 605 equations. 606 607 Input Parameters: 608 . snes - the SNES context 609 . func - function evaluation routine 610 . r - vector to store function value 611 . ctx - [optional] user-defined context for private data for the 612 function evaluation routine (may be PETSC_NULL) 613 614 Collective on SNES 615 616 Calling sequence of func: 617 . func (SNES snes,Vec x,Vec f,void *ctx); 618 619 . x - input vector 620 . f - function vector 621 . ctx - optional user-defined function context 622 623 Notes: 624 The Newton-like methods typically solve linear systems of the form 625 $ f'(x) x = -f(x), 626 $ where f'(x) denotes the Jacobian matrix and f(x) is the function. 627 628 SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 629 Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 630 SNESSetMinimizationFunction() and SNESSetGradient(); 631 632 .keywords: SNES, nonlinear, set, function 633 634 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 635 @*/ 636 int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx) 637 { 638 PetscFunctionBegin; 639 PetscValidHeaderSpecific(snes,SNES_COOKIE); 640 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 641 SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 642 } 643 snes->computefunction = func; 644 snes->vec_func = snes->vec_func_always = r; 645 snes->funP = ctx; 646 PetscFunctionReturn(0); 647 } 648 649 #undef __FUNC__ 650 #define __FUNC__ "SNESComputeFunction" 651 /*@ 652 SNESComputeFunction - Computes the function that has been set with 653 SNESSetFunction(). 654 655 Input Parameters: 656 . snes - the SNES context 657 . x - input vector 658 659 Output Parameter: 660 . y - function vector, as set by SNESSetFunction() 661 662 Collective on SNES 663 664 Notes: 665 SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 666 Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 667 SNESComputeMinimizationFunction() and SNESComputeGradient(); 668 669 .keywords: SNES, nonlinear, compute, function 670 671 .seealso: SNESSetFunction(), SNESGetFunction() 672 @*/ 673 int SNESComputeFunction(SNES snes,Vec x, Vec y) 674 { 675 int ierr; 676 677 PetscFunctionBegin; 678 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 679 SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 680 } 681 PLogEventBegin(SNES_FunctionEval,snes,x,y,0); 682 PetscStackPush("SNES user function"); 683 ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 684 PetscStackPop; 685 snes->nfuncs++; 686 PLogEventEnd(SNES_FunctionEval,snes,x,y,0); 687 PetscFunctionReturn(0); 688 } 689 690 #undef __FUNC__ 691 #define __FUNC__ "SNESSetMinimizationFunction" 692 /*@C 693 SNESSetMinimizationFunction - Sets the function evaluation routine for 694 unconstrained minimization. 695 696 Input Parameters: 697 . snes - the SNES context 698 . func - function evaluation routine 699 . ctx - [optional] user-defined context for private data for the 700 function evaluation routine (may be PETSC_NULL) 701 702 Collective on SNES 703 704 Calling sequence of func: 705 . func (SNES snes,Vec x,double *f,void *ctx); 706 707 . x - input vector 708 . f - function 709 . ctx - [optional] user-defined function context 710 711 Notes: 712 SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 713 methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 714 SNESSetFunction(). 715 716 .keywords: SNES, nonlinear, set, minimization, function 717 718 .seealso: SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(), 719 SNESSetHessian(), SNESSetGradient() 720 @*/ 721 int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*), 722 void *ctx) 723 { 724 PetscFunctionBegin; 725 PetscValidHeaderSpecific(snes,SNES_COOKIE); 726 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 727 SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION"); 728 } 729 snes->computeumfunction = func; 730 snes->umfunP = ctx; 731 PetscFunctionReturn(0); 732 } 733 734 #undef __FUNC__ 735 #define __FUNC__ "SNESComputeMinimizationFunction" 736 /*@ 737 SNESComputeMinimizationFunction - Computes the function that has been 738 set with SNESSetMinimizationFunction(). 739 740 Input Parameters: 741 . snes - the SNES context 742 . x - input vector 743 744 Output Parameter: 745 . y - function value 746 747 Collective on SNES 748 749 Notes: 750 SNESComputeMinimizationFunction() is valid only for 751 SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 752 SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 753 754 .keywords: SNES, nonlinear, compute, minimization, function 755 756 .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(), 757 SNESComputeGradient(), SNESComputeHessian() 758 @*/ 759 int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y) 760 { 761 int ierr; 762 763 PetscFunctionBegin; 764 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 765 SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION"); 766 } 767 PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0); 768 PetscStackPush("SNES user minimzation function"); 769 ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr); 770 PetscStackPop; 771 snes->nfuncs++; 772 PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0); 773 PetscFunctionReturn(0); 774 } 775 776 #undef __FUNC__ 777 #define __FUNC__ "SNESSetGradient" 778 /*@C 779 SNESSetGradient - Sets the gradient evaluation routine and gradient 780 vector for use by the SNES routines. 781 782 Input Parameters: 783 . snes - the SNES context 784 . func - function evaluation routine 785 . ctx - optional user-defined context for private data for the 786 gradient evaluation routine (may be PETSC_NULL) 787 . r - vector to store gradient value 788 789 Collective on SNES 790 791 Calling sequence of func: 792 . func (SNES, Vec x, Vec g, void *ctx); 793 794 . x - input vector 795 . g - gradient vector 796 . ctx - optional user-defined gradient context 797 798 Notes: 799 SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 800 methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 801 SNESSetFunction(). 802 803 .keywords: SNES, nonlinear, set, function 804 805 .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(), 806 SNESSetMinimizationFunction(), 807 @*/ 808 int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx) 809 { 810 PetscFunctionBegin; 811 PetscValidHeaderSpecific(snes,SNES_COOKIE); 812 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 813 SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 814 } 815 snes->computefunction = func; 816 snes->vec_func = snes->vec_func_always = r; 817 snes->funP = ctx; 818 PetscFunctionReturn(0); 819 } 820 821 #undef __FUNC__ 822 #define __FUNC__ "SNESComputeGradient" 823 /*@ 824 SNESComputeGradient - Computes the gradient that has been set with 825 SNESSetGradient(). 826 827 Input Parameters: 828 . snes - the SNES context 829 . x - input vector 830 831 Output Parameter: 832 . y - gradient vector 833 834 Collective on SNES 835 836 Notes: 837 SNESComputeGradient() is valid only for 838 SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 839 SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 840 841 .keywords: SNES, nonlinear, compute, gradient 842 843 .seealso: SNESSetGradient(), SNESGetGradient(), 844 SNESComputeMinimizationFunction(), SNESComputeHessian() 845 @*/ 846 int SNESComputeGradient(SNES snes,Vec x, Vec y) 847 { 848 int ierr; 849 850 PetscFunctionBegin; 851 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 852 SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 853 } 854 PLogEventBegin(SNES_GradientEval,snes,x,y,0); 855 PetscStackPush("SNES user gradient function"); 856 ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 857 PetscStackPop; 858 PLogEventEnd(SNES_GradientEval,snes,x,y,0); 859 PetscFunctionReturn(0); 860 } 861 862 #undef __FUNC__ 863 #define __FUNC__ "SNESComputeJacobian" 864 /*@ 865 SNESComputeJacobian - Computes the Jacobian matrix that has been 866 set with SNESSetJacobian(). 867 868 Input Parameters: 869 . snes - the SNES context 870 . x - input vector 871 872 Output Parameters: 873 . A - Jacobian matrix 874 . B - optional preconditioning matrix 875 . flag - flag indicating matrix structure 876 877 Collective on SNES and Mat 878 879 Notes: 880 Most users should not need to explicitly call this routine, as it 881 is used internally within the nonlinear solvers. 882 883 See SLESSetOperators() for important information about setting the 884 flag parameter. 885 886 SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS 887 methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION 888 methods is SNESComputeHessian(). 889 890 .keywords: SNES, compute, Jacobian, matrix 891 892 .seealso: SNESSetJacobian(), SLESSetOperators() 893 @*/ 894 int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 895 { 896 int ierr; 897 898 PetscFunctionBegin; 899 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 900 SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 901 } 902 if (!snes->computejacobian) PetscFunctionReturn(0); 903 PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B); 904 *flg = DIFFERENT_NONZERO_PATTERN; 905 PetscStackPush("SNES user Jacobian function"); 906 ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr); 907 PetscStackPop; 908 PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B); 909 /* make sure user returned a correct Jacobian and preconditioner */ 910 PetscValidHeaderSpecific(*A,MAT_COOKIE); 911 PetscValidHeaderSpecific(*B,MAT_COOKIE); 912 PetscFunctionReturn(0); 913 } 914 915 #undef __FUNC__ 916 #define __FUNC__ "SNESComputeHessian" 917 /*@ 918 SNESComputeHessian - Computes the Hessian matrix that has been 919 set with SNESSetHessian(). 920 921 Input Parameters: 922 . snes - the SNES context 923 . x - input vector 924 925 Output Parameters: 926 . A - Hessian matrix 927 . B - optional preconditioning matrix 928 . flag - flag indicating matrix structure 929 930 Collective on SNES and Mat 931 932 Notes: 933 Most users should not need to explicitly call this routine, as it 934 is used internally within the nonlinear solvers. 935 936 See SLESSetOperators() for important information about setting the 937 flag parameter. 938 939 SNESComputeHessian() is valid only for 940 SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 941 SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian(). 942 943 .keywords: SNES, compute, Hessian, matrix 944 945 .seealso: SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(), 946 SNESComputeMinimizationFunction() 947 @*/ 948 int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag) 949 { 950 int ierr; 951 952 PetscFunctionBegin; 953 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 954 SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 955 } 956 if (!snes->computejacobian) PetscFunctionReturn(0); 957 PLogEventBegin(SNES_HessianEval,snes,x,*A,*B); 958 *flag = DIFFERENT_NONZERO_PATTERN; 959 PetscStackPush("SNES user Hessian function"); 960 ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr); 961 PetscStackPop; 962 PLogEventEnd(SNES_HessianEval,snes,x,*A,*B); 963 /* make sure user returned a correct Jacobian and preconditioner */ 964 PetscValidHeaderSpecific(*A,MAT_COOKIE); 965 PetscValidHeaderSpecific(*B,MAT_COOKIE); 966 PetscFunctionReturn(0); 967 } 968 969 #undef __FUNC__ 970 #define __FUNC__ "SNESSetJacobian" 971 /*@C 972 SNESSetJacobian - Sets the function to compute Jacobian as well as the 973 location to store the matrix. 974 975 Input Parameters: 976 . snes - the SNES context 977 . A - Jacobian matrix 978 . B - preconditioner matrix (usually same as the Jacobian) 979 . func - Jacobian evaluation routine 980 . ctx - [optional] user-defined context for private data for the 981 Jacobian evaluation routine (may be PETSC_NULL) 982 983 Collective on SNES and Mat 984 985 Calling sequence of func: 986 . func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 987 988 . x - input vector 989 . A - Jacobian matrix 990 . B - preconditioner matrix, usually the same as A 991 . flag - flag indicating information about the preconditioner matrix 992 structure (same as flag in SLESSetOperators()) 993 . ctx - [optional] user-defined Jacobian context 994 995 Notes: 996 See SLESSetOperators() for important information about setting the flag 997 output parameter in the routine func(). Be sure to read this information! 998 999 The routine func() takes Mat * as the matrix arguments rather than Mat. 1000 This allows the Jacobian evaluation routine to replace A and/or B with a 1001 completely new new matrix structure (not just different matrix elements) 1002 when appropriate, for instance, if the nonzero structure is changing 1003 throughout the global iterations. 1004 1005 .keywords: SNES, nonlinear, set, Jacobian, matrix 1006 1007 .seealso: SLESSetOperators(), SNESSetFunction() 1008 @*/ 1009 int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 1010 MatStructure*,void*),void *ctx) 1011 { 1012 PetscFunctionBegin; 1013 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1014 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 1015 SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 1016 } 1017 snes->computejacobian = func; 1018 snes->jacP = ctx; 1019 snes->jacobian = A; 1020 snes->jacobian_pre = B; 1021 PetscFunctionReturn(0); 1022 } 1023 1024 #undef __FUNC__ 1025 #define __FUNC__ "SNESGetJacobian" 1026 /*@ 1027 SNESGetJacobian - Returns the Jacobian matrix and optionally the user 1028 provided context for evaluating the Jacobian. 1029 1030 Input Parameter: 1031 . snes - the nonlinear solver context 1032 1033 Output Parameters: 1034 . A - location to stash Jacobian matrix (or PETSC_NULL) 1035 . B - location to stash preconditioner matrix (or PETSC_NULL) 1036 . ctx - location to stash Jacobian ctx (or PETSC_NULL) 1037 1038 Not Collective, but Mat object will be parallel if SNES object is 1039 1040 .seealso: SNESSetJacobian(), SNESComputeJacobian() 1041 @*/ 1042 int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx) 1043 { 1044 PetscFunctionBegin; 1045 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 1046 SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 1047 } 1048 if (A) *A = snes->jacobian; 1049 if (B) *B = snes->jacobian_pre; 1050 if (ctx) *ctx = snes->jacP; 1051 PetscFunctionReturn(0); 1052 } 1053 1054 #undef __FUNC__ 1055 #define __FUNC__ "SNESSetHessian" 1056 /*@C 1057 SNESSetHessian - Sets the function to compute Hessian as well as the 1058 location to store the matrix. 1059 1060 Input Parameters: 1061 . snes - the SNES context 1062 . A - Hessian matrix 1063 . B - preconditioner matrix (usually same as the Hessian) 1064 . func - Jacobian evaluation routine 1065 . ctx - [optional] user-defined context for private data for the 1066 Hessian evaluation routine (may be PETSC_NULL) 1067 1068 Collective on SNES and Mat 1069 1070 Calling sequence of func: 1071 . func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 1072 1073 . x - input vector 1074 . A - Hessian matrix 1075 . B - preconditioner matrix, usually the same as A 1076 . flag - flag indicating information about the preconditioner matrix 1077 structure (same as flag in SLESSetOperators()) 1078 . ctx - [optional] user-defined Hessian context 1079 1080 Notes: 1081 See SLESSetOperators() for important information about setting the flag 1082 output parameter in the routine func(). Be sure to read this information! 1083 1084 The function func() takes Mat * as the matrix arguments rather than Mat. 1085 This allows the Hessian evaluation routine to replace A and/or B with a 1086 completely new new matrix structure (not just different matrix elements) 1087 when appropriate, for instance, if the nonzero structure is changing 1088 throughout the global iterations. 1089 1090 .keywords: SNES, nonlinear, set, Hessian, matrix 1091 1092 .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators() 1093 @*/ 1094 int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 1095 MatStructure*,void*),void *ctx) 1096 { 1097 PetscFunctionBegin; 1098 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1099 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 1100 SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 1101 } 1102 snes->computejacobian = func; 1103 snes->jacP = ctx; 1104 snes->jacobian = A; 1105 snes->jacobian_pre = B; 1106 PetscFunctionReturn(0); 1107 } 1108 1109 #undef __FUNC__ 1110 #define __FUNC__ "SNESGetHessian" 1111 /*@ 1112 SNESGetHessian - Returns the Hessian matrix and optionally the user 1113 provided context for evaluating the Hessian. 1114 1115 Input Parameter: 1116 . snes - the nonlinear solver context 1117 1118 Output Parameters: 1119 . A - location to stash Hessian matrix (or PETSC_NULL) 1120 . B - location to stash preconditioner matrix (or PETSC_NULL) 1121 . ctx - location to stash Hessian ctx (or PETSC_NULL) 1122 1123 Not Collective, but Mat object is parallel if SNES object is parallel 1124 1125 .seealso: SNESSetHessian(), SNESComputeHessian() 1126 @*/ 1127 int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx) 1128 { 1129 PetscFunctionBegin; 1130 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION){ 1131 SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 1132 } 1133 if (A) *A = snes->jacobian; 1134 if (B) *B = snes->jacobian_pre; 1135 if (ctx) *ctx = snes->jacP; 1136 PetscFunctionReturn(0); 1137 } 1138 1139 /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 1140 1141 #undef __FUNC__ 1142 #define __FUNC__ "SNESSetUp" 1143 /*@ 1144 SNESSetUp - Sets up the internal data structures for the later use 1145 of a nonlinear solver. 1146 1147 Input Parameter: 1148 . snes - the SNES context 1149 . x - the solution vector 1150 1151 Collective on SNES 1152 1153 Notes: 1154 For basic use of the SNES solvers the user need not explicitly call 1155 SNESSetUp(), since these actions will automatically occur during 1156 the call to SNESSolve(). However, if one wishes to control this 1157 phase separately, SNESSetUp() should be called after SNESCreate() 1158 and optional routines of the form SNESSetXXX(), but before SNESSolve(). 1159 1160 .keywords: SNES, nonlinear, setup 1161 1162 .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 1163 @*/ 1164 int SNESSetUp(SNES snes,Vec x) 1165 { 1166 int ierr, flg; 1167 1168 PetscFunctionBegin; 1169 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1170 PetscValidHeaderSpecific(x,VEC_COOKIE); 1171 snes->vec_sol = snes->vec_sol_always = x; 1172 1173 ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg); CHKERRQ(ierr); 1174 /* 1175 This version replaces the user provided Jacobian matrix with a 1176 matrix-free version but still employs the user-provided preconditioner matrix 1177 */ 1178 if (flg) { 1179 Mat J; 1180 ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1181 PLogObjectParent(snes,J); 1182 snes->mfshell = J; 1183 if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 1184 snes->jacobian = J; 1185 PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n"); 1186 } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 1187 snes->jacobian = J; 1188 PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n"); 1189 } else { 1190 SETERRQ(PETSC_ERR_SUP,0,"Method class doesn't support matrix-free operator option"); 1191 } 1192 } 1193 ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg); CHKERRQ(ierr); 1194 /* 1195 This version replaces both the user-provided Jacobian and the user- 1196 provided preconditioner matrix with the default matrix free version. 1197 */ 1198 if (flg) { 1199 Mat J; 1200 ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1201 PLogObjectParent(snes,J); 1202 snes->mfshell = J; 1203 if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 1204 ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 1205 PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n"); 1206 } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 1207 ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 1208 PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n"); 1209 } else { 1210 SETERRQ(PETSC_ERR_SUP,0,"Method class doesn't support matrix-free option"); 1211 } 1212 } 1213 if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) { 1214 if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetFunction() first"); 1215 if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetFunction() first"); 1216 if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetJacobian() first"); 1217 if (snes->vec_func == snes->vec_sol) { 1218 SETERRQ(PETSC_ERR_ARG_IDN,0,"Solution vector cannot be function vector"); 1219 } 1220 1221 /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 1222 if (snes->ksp_ewconv && PetscStrcmp(snes->type_name,SNES_EQ_TR)) { 1223 SLES sles; KSP ksp; 1224 ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 1225 ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr); 1226 ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private, 1227 (void *)snes); CHKERRQ(ierr); 1228 } 1229 } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) { 1230 if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetGradient() first"); 1231 if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetGradient() first"); 1232 if (!snes->computeumfunction) { 1233 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetMinimizationFunction() first"); 1234 } 1235 if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetHessian() first"); 1236 } else { 1237 SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown method class"); 1238 } 1239 if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);} 1240 snes->setupcalled = 1; 1241 PetscFunctionReturn(0); 1242 } 1243 1244 #undef __FUNC__ 1245 #define __FUNC__ "SNESDestroy" 1246 /*@C 1247 SNESDestroy - Destroys the nonlinear solver context that was created 1248 with SNESCreate(). 1249 1250 Input Parameter: 1251 . snes - the SNES context 1252 1253 Collective on SNES 1254 1255 .keywords: SNES, nonlinear, destroy 1256 1257 .seealso: SNESCreate(), SNESSolve() 1258 @*/ 1259 int SNESDestroy(SNES snes) 1260 { 1261 int ierr; 1262 1263 PetscFunctionBegin; 1264 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1265 if (--snes->refct > 0) PetscFunctionReturn(0); 1266 1267 if (snes->destroy) {ierr = (*(snes)->destroy)(snes); CHKERRQ(ierr);} 1268 if (snes->kspconvctx) PetscFree(snes->kspconvctx); 1269 if (snes->mfshell) {ierr = MatDestroy(snes->mfshell);CHKERRQ(ierr);} 1270 ierr = SLESDestroy(snes->sles); CHKERRQ(ierr); 1271 if (snes->xmonitor) {ierr = SNESLGMonitorDestroy(snes->xmonitor);CHKERRQ(ierr);} 1272 if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);} 1273 PLogObjectDestroy((PetscObject)snes); 1274 PetscHeaderDestroy((PetscObject)snes); 1275 PetscFunctionReturn(0); 1276 } 1277 1278 /* ----------- Routines to set solver parameters ---------- */ 1279 1280 #undef __FUNC__ 1281 #define __FUNC__ "SNESSetTolerances" 1282 /*@ 1283 SNESSetTolerances - Sets various parameters used in convergence tests. 1284 1285 Input Parameters: 1286 . snes - the SNES context 1287 . atol - absolute convergence tolerance 1288 . rtol - relative convergence tolerance 1289 . stol - convergence tolerance in terms of the norm 1290 of the change in the solution between steps 1291 . maxit - maximum number of iterations 1292 . maxf - maximum number of function evaluations 1293 1294 Collective on SNES 1295 1296 Options Database Keys: 1297 $ -snes_atol <atol> 1298 $ -snes_rtol <rtol> 1299 $ -snes_stol <stol> 1300 $ -snes_max_it <maxit> 1301 $ -snes_max_funcs <maxf> 1302 1303 Notes: 1304 The default maximum number of iterations is 50. 1305 The default maximum number of function evaluations is 1000. 1306 1307 .keywords: SNES, nonlinear, set, convergence, tolerances 1308 1309 .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance() 1310 @*/ 1311 int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf) 1312 { 1313 PetscFunctionBegin; 1314 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1315 if (atol != PETSC_DEFAULT) snes->atol = atol; 1316 if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1317 if (stol != PETSC_DEFAULT) snes->xtol = stol; 1318 if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1319 if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 1320 PetscFunctionReturn(0); 1321 } 1322 1323 #undef __FUNC__ 1324 #define __FUNC__ "SNESGetTolerances" 1325 /*@ 1326 SNESGetTolerances - Gets various parameters used in convergence tests. 1327 1328 Input Parameters: 1329 . snes - the SNES context 1330 . atol - absolute convergence tolerance 1331 . rtol - relative convergence tolerance 1332 . stol - convergence tolerance in terms of the norm 1333 of the change in the solution between steps 1334 . maxit - maximum number of iterations 1335 . maxf - maximum number of function evaluations 1336 1337 Not Collective 1338 1339 Notes: 1340 The user can specify PETSC_NULL for any parameter that is not needed. 1341 1342 .keywords: SNES, nonlinear, get, convergence, tolerances 1343 1344 .seealso: SNESSetTolerances() 1345 @*/ 1346 int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf) 1347 { 1348 PetscFunctionBegin; 1349 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1350 if (atol) *atol = snes->atol; 1351 if (rtol) *rtol = snes->rtol; 1352 if (stol) *stol = snes->xtol; 1353 if (maxit) *maxit = snes->max_its; 1354 if (maxf) *maxf = snes->max_funcs; 1355 PetscFunctionReturn(0); 1356 } 1357 1358 #undef __FUNC__ 1359 #define __FUNC__ "SNESSetTrustRegionTolerance" 1360 /*@ 1361 SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 1362 1363 Input Parameters: 1364 . snes - the SNES context 1365 . tol - tolerance 1366 1367 Collective on SNES 1368 1369 Options Database Key: 1370 $ -snes_trtol <tol> 1371 1372 .keywords: SNES, nonlinear, set, trust region, tolerance 1373 1374 .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance() 1375 @*/ 1376 int SNESSetTrustRegionTolerance(SNES snes,double tol) 1377 { 1378 PetscFunctionBegin; 1379 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1380 snes->deltatol = tol; 1381 PetscFunctionReturn(0); 1382 } 1383 1384 #undef __FUNC__ 1385 #define __FUNC__ "SNESSetMinimizationFunctionTolerance" 1386 /*@ 1387 SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance 1388 for unconstrained minimization solvers. 1389 1390 Input Parameters: 1391 . snes - the SNES context 1392 . ftol - minimum function tolerance 1393 1394 Collective on SNES 1395 1396 Options Database Key: 1397 $ -snes_fmin <ftol> 1398 1399 Note: 1400 SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION 1401 methods only. 1402 1403 .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance 1404 1405 .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance() 1406 @*/ 1407 int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol) 1408 { 1409 PetscFunctionBegin; 1410 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1411 snes->fmin = ftol; 1412 PetscFunctionReturn(0); 1413 } 1414 1415 /* ------------ Routines to set performance monitoring options ----------- */ 1416 1417 #undef __FUNC__ 1418 #define __FUNC__ "SNESSetMonitor" 1419 /*@C 1420 SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every 1421 iteration of the nonlinear solver to display the iteration's 1422 progress. 1423 1424 Input Parameters: 1425 . snes - the SNES context 1426 . func - monitoring routine 1427 . mctx - [optional] user-defined context for private data for the 1428 monitor routine (may be PETSC_NULL) 1429 1430 Collective on SNES 1431 1432 Calling sequence of func: 1433 int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx) 1434 1435 $ snes - the SNES context 1436 $ its - iteration number 1437 $ mctx - [optional] monitoring context 1438 $ 1439 $ SNES_NONLINEAR_EQUATIONS methods: 1440 $ norm - 2-norm function value (may be estimated) 1441 $ 1442 $ SNES_UNCONSTRAINED_MINIMIZATION methods: 1443 $ norm - 2-norm gradient value (may be estimated) 1444 1445 Options Database Keys: 1446 $ -snes_monitor : sets SNESDefaultMonitor() 1447 $ -snes_xmonitor : sets line graph monitor, 1448 $ uses SNESLGMonitorCreate() 1449 $ -snes_cancelmonitors : cancels all monitors that have 1450 $ been hardwired into a code by 1451 $ calls to SNESSetMonitor(), but 1452 $ does not cancel those set via 1453 $ the options database. 1454 1455 1456 Notes: 1457 Several different monitoring routines may be set by calling 1458 SNESSetMonitor() multiple times; all will be called in the 1459 order in which they were set. 1460 1461 .keywords: SNES, nonlinear, set, monitor 1462 1463 .seealso: SNESDefaultMonitor(), SNESClearMonitor() 1464 @*/ 1465 int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx ) 1466 { 1467 PetscFunctionBegin; 1468 if (snes->numbermonitors >= MAXSNESMONITORS) { 1469 SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many monitors set"); 1470 } 1471 1472 snes->monitor[snes->numbermonitors] = func; 1473 snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 1474 PetscFunctionReturn(0); 1475 } 1476 1477 #undef __FUNC__ 1478 #define __FUNC__ "SNESClearMonitor" 1479 /*@C 1480 SNESClearMonitor - Clears all the monitor functions for a SNES object. 1481 1482 Input Parameters: 1483 . snes - the SNES context 1484 1485 Collective on SNES 1486 1487 Options Database: 1488 $ -snes_cancelmonitors : cancels all monitors that have 1489 $ been hardwired into a code by 1490 $ calls to SNESSetMonitor(), but 1491 $ does not cancel those set via 1492 $ the options database. 1493 1494 Notes: 1495 There is no way to clear one specific monitor from a SNES object. 1496 1497 .keywords: SNES, nonlinear, set, monitor 1498 1499 .seealso: SNESDefaultMonitor(), SNESSetMonitor() 1500 @*/ 1501 int SNESClearMonitor( SNES snes ) 1502 { 1503 PetscFunctionBegin; 1504 snes->numbermonitors = 0; 1505 PetscFunctionReturn(0); 1506 } 1507 1508 #undef __FUNC__ 1509 #define __FUNC__ "SNESSetConvergenceTest" 1510 /*@C 1511 SNESSetConvergenceTest - Sets the function that is to be used 1512 to test for convergence of the nonlinear iterative solution. 1513 1514 Input Parameters: 1515 . snes - the SNES context 1516 . func - routine to test for convergence 1517 . cctx - [optional] context for private data for the convergence routine 1518 (may be PETSC_NULL) 1519 1520 Collective on SNES 1521 1522 Calling sequence of func: 1523 int func (SNES snes,double xnorm,double gnorm, 1524 double f,void *cctx) 1525 1526 $ snes - the SNES context 1527 $ cctx - [optional] convergence context 1528 $ xnorm - 2-norm of current iterate 1529 $ 1530 $ SNES_NONLINEAR_EQUATIONS methods: 1531 $ gnorm - 2-norm of current step 1532 $ f - 2-norm of function 1533 $ 1534 $ SNES_UNCONSTRAINED_MINIMIZATION methods: 1535 $ gnorm - 2-norm of current gradient 1536 $ f - function value 1537 1538 .keywords: SNES, nonlinear, set, convergence, test 1539 1540 .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(), 1541 SNESConverged_UM_LS(), SNESConverged_UM_TR() 1542 @*/ 1543 int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx) 1544 { 1545 PetscFunctionBegin; 1546 (snes)->converged = func; 1547 (snes)->cnvP = cctx; 1548 PetscFunctionReturn(0); 1549 } 1550 1551 #undef __FUNC__ 1552 #define __FUNC__ "SNESSetConvergenceHistory" 1553 /*@ 1554 SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1555 1556 Input Parameters: 1557 . snes - iterative context obtained from SNESCreate() 1558 . a - array to hold history 1559 . na - size of a 1560 1561 Collective on SNES 1562 1563 Notes: 1564 If set, this array will contain the function norms (for 1565 SNES_NONLINEAR_EQUATIONS methods) or gradient norms 1566 (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed 1567 at each step. 1568 1569 This routine is useful, e.g., when running a code for purposes 1570 of accurate performance monitoring, when no I/O should be done 1571 during the section of code that is being timed. 1572 1573 .keywords: SNES, set, convergence, history 1574 @*/ 1575 int SNESSetConvergenceHistory(SNES snes, double *a, int na) 1576 { 1577 PetscFunctionBegin; 1578 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1579 if (na) PetscValidScalarPointer(a); 1580 snes->conv_hist = a; 1581 snes->conv_hist_size = na; 1582 PetscFunctionReturn(0); 1583 } 1584 1585 #undef __FUNC__ 1586 #define __FUNC__ "SNESScaleStep_Private" 1587 /* 1588 SNESScaleStep_Private - Scales a step so that its length is less than the 1589 positive parameter delta. 1590 1591 Input Parameters: 1592 . snes - the SNES context 1593 . y - approximate solution of linear system 1594 . fnorm - 2-norm of current function 1595 . delta - trust region size 1596 1597 Output Parameters: 1598 . gpnorm - predicted function norm at the new point, assuming local 1599 linearization. The value is zero if the step lies within the trust 1600 region, and exceeds zero otherwise. 1601 . ynorm - 2-norm of the step 1602 1603 Note: 1604 For non-trust region methods such as SNES_EQ_LS, the parameter delta 1605 is set to be the maximum allowable step size. 1606 1607 .keywords: SNES, nonlinear, scale, step 1608 */ 1609 int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta, 1610 double *gpnorm,double *ynorm) 1611 { 1612 double norm; 1613 Scalar cnorm; 1614 int ierr; 1615 1616 PetscFunctionBegin; 1617 ierr = VecNorm(y,NORM_2, &norm );CHKERRQ(ierr); 1618 if (norm > *delta) { 1619 norm = *delta/norm; 1620 *gpnorm = (1.0 - norm)*(*fnorm); 1621 cnorm = norm; 1622 VecScale( &cnorm, y ); 1623 *ynorm = *delta; 1624 } else { 1625 *gpnorm = 0.0; 1626 *ynorm = norm; 1627 } 1628 PetscFunctionReturn(0); 1629 } 1630 1631 #undef __FUNC__ 1632 #define __FUNC__ "SNESSolve" 1633 /*@ 1634 SNESSolve - Solves a nonlinear system. Call SNESSolve after calling 1635 SNESCreate() and optional routines of the form SNESSetXXX(). 1636 1637 Input Parameter: 1638 . snes - the SNES context 1639 . x - the solution vector 1640 1641 Output Parameter: 1642 its - number of iterations until termination 1643 1644 Collective on SNES 1645 1646 Note: 1647 The user should initialize the vector, x, with the initial guess 1648 for the nonlinear solve prior to calling SNESSolve. In particular, 1649 to employ an initial guess of zero, the user should explicitly set 1650 this vector to zero by calling VecSet(). 1651 1652 .keywords: SNES, nonlinear, solve 1653 1654 .seealso: SNESCreate(), SNESDestroy() 1655 @*/ 1656 int SNESSolve(SNES snes,Vec x,int *its) 1657 { 1658 int ierr, flg; 1659 1660 PetscFunctionBegin; 1661 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1662 PetscValidIntPointer(its); 1663 if (!snes->setupcalled) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);} 1664 else {snes->vec_sol = snes->vec_sol_always = x;} 1665 PLogEventBegin(SNES_Solve,snes,0,0,0); 1666 snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0; 1667 ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr); 1668 PLogEventEnd(SNES_Solve,snes,0,0,0); 1669 ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr); 1670 if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); } 1671 PetscFunctionReturn(0); 1672 } 1673 1674 /* --------- Internal routines for SNES Package --------- */ 1675 1676 #undef __FUNC__ 1677 #define __FUNC__ "SNESSetType" 1678 /*@C 1679 SNESSetType - Sets the method for the nonlinear solver. 1680 1681 Input Parameters: 1682 . snes - the SNES context 1683 . method - a known method 1684 1685 Collective on SNES 1686 1687 Options Database Command: 1688 $ -snes_type <method> 1689 $ Use -help for a list of available methods 1690 $ (for instance, ls or tr) 1691 1692 Notes: 1693 See "petsc/include/snes.h" for available methods (for instance) 1694 $ Systems of nonlinear equations: 1695 $ SNES_EQ_LS - Newton's method with line search 1696 $ SNES_EQ_TR - Newton's method with trust region 1697 $ Unconstrained minimization: 1698 $ SNES_UM_TR - Newton's method with trust region 1699 $ SNES_UM_LS - Newton's method with line search 1700 1701 Normally, it is best to use the SNESSetFromOptions() command and then 1702 set the SNES solver type from the options database rather than by using 1703 this routine. Using the options database provides the user with 1704 maximum flexibility in evaluating the many nonlinear solvers. 1705 The SNESSetType() routine is provided for those situations where it 1706 is necessary to set the nonlinear solver independently of the command 1707 line or options database. This might be the case, for example, when 1708 the choice of solver changes during the execution of the program, 1709 and the user's application is taking responsibility for choosing the 1710 appropriate method. In other words, this routine is for the advanced user. 1711 1712 .keywords: SNES, set, method 1713 @*/ 1714 int SNESSetType(SNES snes,SNESType method) 1715 { 1716 int ierr; 1717 int (*r)(SNES); 1718 1719 PetscFunctionBegin; 1720 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1721 1722 if (!PetscStrcmp(snes->type_name,method)) PetscFunctionReturn(0); 1723 1724 if (snes->setupcalled) { 1725 ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr); 1726 snes->data = 0; 1727 } 1728 1729 /* Get the function pointers for the iterative method requested */ 1730 if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL); CHKERRQ(ierr);} 1731 1732 ierr = DLRegisterFind(snes->comm, SNESList, method,(int (**)(void *)) &r );CHKERRQ(ierr); 1733 1734 if (snes->data) PetscFree(snes->data); 1735 snes->data = 0; 1736 ierr = (*r)(snes); CHKERRQ(ierr); 1737 1738 if (snes->type_name) PetscFree(snes->type_name); 1739 snes->type_name = (char *) PetscMalloc((PetscStrlen(method)+1)*sizeof(char));CHKPTRQ(snes->type_name); 1740 PetscStrcpy(snes->type_name,method); 1741 snes->set_method_called = 1; 1742 1743 PetscFunctionReturn(0); 1744 } 1745 1746 1747 /* --------------------------------------------------------------------- */ 1748 #undef __FUNC__ 1749 #define __FUNC__ "SNESRegisterDestroy" 1750 /*@C 1751 SNESRegisterDestroy - Frees the list of nonlinear solvers that were 1752 registered by SNESRegister(). 1753 1754 Not Collective 1755 1756 .keywords: SNES, nonlinear, register, destroy 1757 1758 .seealso: SNESRegisterAll(), SNESRegisterAll() 1759 @*/ 1760 int SNESRegisterDestroy(void) 1761 { 1762 int ierr; 1763 1764 PetscFunctionBegin; 1765 if (SNESList) { 1766 ierr = DLRegisterDestroy( SNESList );CHKERRQ(ierr); 1767 SNESList = 0; 1768 } 1769 SNESRegisterAllCalled = 0; 1770 PetscFunctionReturn(0); 1771 } 1772 1773 #undef __FUNC__ 1774 #define __FUNC__ "SNESGetType" 1775 /*@C 1776 SNESGetType - Gets the SNES method type and name (as a string). 1777 1778 Input Parameter: 1779 . snes - nonlinear solver context 1780 1781 Output Parameter: 1782 . method - SNES method (a charactor string) 1783 1784 Not Collective 1785 1786 .keywords: SNES, nonlinear, get, method, name 1787 @*/ 1788 int SNESGetType(SNES snes, SNESType *method) 1789 { 1790 PetscFunctionBegin; 1791 *method = snes->type_name; 1792 PetscFunctionReturn(0); 1793 } 1794 1795 #undef __FUNC__ 1796 #define __FUNC__ "SNESGetSolution" 1797 /*@C 1798 SNESGetSolution - Returns the vector where the approximate solution is 1799 stored. 1800 1801 Input Parameter: 1802 . snes - the SNES context 1803 1804 Output Parameter: 1805 . x - the solution 1806 1807 Not Collective, but Vec is parallel if SNES is parallel 1808 1809 .keywords: SNES, nonlinear, get, solution 1810 1811 .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate() 1812 @*/ 1813 int SNESGetSolution(SNES snes,Vec *x) 1814 { 1815 PetscFunctionBegin; 1816 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1817 *x = snes->vec_sol_always; 1818 PetscFunctionReturn(0); 1819 } 1820 1821 #undef __FUNC__ 1822 #define __FUNC__ "SNESGetSolutionUpdate" 1823 /*@C 1824 SNESGetSolutionUpdate - Returns the vector where the solution update is 1825 stored. 1826 1827 Input Parameter: 1828 . snes - the SNES context 1829 1830 Output Parameter: 1831 . x - the solution update 1832 1833 Not Collective, but Vec is parallel if SNES is parallel 1834 1835 .keywords: SNES, nonlinear, get, solution, update 1836 1837 .seealso: SNESGetSolution(), SNESGetFunction 1838 @*/ 1839 int SNESGetSolutionUpdate(SNES snes,Vec *x) 1840 { 1841 PetscFunctionBegin; 1842 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1843 *x = snes->vec_sol_update_always; 1844 PetscFunctionReturn(0); 1845 } 1846 1847 #undef __FUNC__ 1848 #define __FUNC__ "SNESGetFunction" 1849 /*@C 1850 SNESGetFunction - Returns the vector where the function is stored. 1851 1852 Input Parameter: 1853 . snes - the SNES context 1854 1855 Output Parameter: 1856 . r - the function 1857 1858 Not Collective, but Vec is parallel if SNES is parallel 1859 1860 Notes: 1861 SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only 1862 Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 1863 SNESGetMinimizationFunction() and SNESGetGradient(); 1864 1865 .keywords: SNES, nonlinear, get, function 1866 1867 .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(), 1868 SNESGetGradient() 1869 @*/ 1870 int SNESGetFunction(SNES snes,Vec *r) 1871 { 1872 PetscFunctionBegin; 1873 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1874 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 1875 SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 1876 } 1877 *r = snes->vec_func_always; 1878 PetscFunctionReturn(0); 1879 } 1880 1881 #undef __FUNC__ 1882 #define __FUNC__ "SNESGetGradient" 1883 /*@C 1884 SNESGetGradient - Returns the vector where the gradient is stored. 1885 1886 Input Parameter: 1887 . snes - the SNES context 1888 1889 Output Parameter: 1890 . r - the gradient 1891 1892 Not Collective, but Vec is parallel if SNES is parallel 1893 1894 Notes: 1895 SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods 1896 only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 1897 SNESGetFunction(). 1898 1899 .keywords: SNES, nonlinear, get, gradient 1900 1901 .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction() 1902 @*/ 1903 int SNESGetGradient(SNES snes,Vec *r) 1904 { 1905 PetscFunctionBegin; 1906 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1907 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 1908 SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 1909 } 1910 *r = snes->vec_func_always; 1911 PetscFunctionReturn(0); 1912 } 1913 1914 #undef __FUNC__ 1915 #define __FUNC__ "SNESGetMinimizationFunction" 1916 /*@ 1917 SNESGetMinimizationFunction - Returns the scalar function value for 1918 unconstrained minimization problems. 1919 1920 Input Parameter: 1921 . snes - the SNES context 1922 1923 Output Parameter: 1924 . r - the function 1925 1926 Not Collective, but Vec is parallel if SNES is parallel 1927 1928 Notes: 1929 SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 1930 methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 1931 SNESGetFunction(). 1932 1933 .keywords: SNES, nonlinear, get, function 1934 1935 .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction() 1936 @*/ 1937 int SNESGetMinimizationFunction(SNES snes,double *r) 1938 { 1939 PetscFunctionBegin; 1940 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1941 PetscValidScalarPointer(r); 1942 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 1943 SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 1944 } 1945 *r = snes->fc; 1946 PetscFunctionReturn(0); 1947 } 1948 1949 #undef __FUNC__ 1950 #define __FUNC__ "SNESSetOptionsPrefix" 1951 /*@C 1952 SNESSetOptionsPrefix - Sets the prefix used for searching for all 1953 SNES options in the database. 1954 1955 Input Parameter: 1956 . snes - the SNES context 1957 . prefix - the prefix to prepend to all option names 1958 1959 Collective on SNES 1960 1961 Notes: 1962 A hyphen (-) must NOT be given at the beginning of the prefix name. 1963 The first character of all runtime options is AUTOMATICALLY the 1964 hyphen. 1965 1966 .keywords: SNES, set, options, prefix, database 1967 1968 .seealso: SNESSetFromOptions() 1969 @*/ 1970 int SNESSetOptionsPrefix(SNES snes,char *prefix) 1971 { 1972 int ierr; 1973 1974 PetscFunctionBegin; 1975 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1976 ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 1977 ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 1978 PetscFunctionReturn(0); 1979 } 1980 1981 #undef __FUNC__ 1982 #define __FUNC__ "SNESAppendOptionsPrefix" 1983 /*@C 1984 SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 1985 SNES options in the database. 1986 1987 Input Parameter: 1988 . snes - the SNES context 1989 . prefix - the prefix to prepend to all option names 1990 1991 Collective on SNES 1992 1993 Notes: 1994 A hyphen (-) must NOT be given at the beginning of the prefix name. 1995 The first character of all runtime options is AUTOMATICALLY the 1996 hyphen. 1997 1998 .keywords: SNES, append, options, prefix, database 1999 2000 .seealso: SNESGetOptionsPrefix() 2001 @*/ 2002 int SNESAppendOptionsPrefix(SNES snes,char *prefix) 2003 { 2004 int ierr; 2005 2006 PetscFunctionBegin; 2007 PetscValidHeaderSpecific(snes,SNES_COOKIE); 2008 ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 2009 ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 2010 PetscFunctionReturn(0); 2011 } 2012 2013 #undef __FUNC__ 2014 #define __FUNC__ "SNESGetOptionsPrefix" 2015 /*@ 2016 SNESGetOptionsPrefix - Sets the prefix used for searching for all 2017 SNES options in the database. 2018 2019 Input Parameter: 2020 . snes - the SNES context 2021 2022 Output Parameter: 2023 . prefix - pointer to the prefix string used 2024 2025 Not Collective 2026 2027 .keywords: SNES, get, options, prefix, database 2028 2029 .seealso: SNESAppendOptionsPrefix() 2030 @*/ 2031 int SNESGetOptionsPrefix(SNES snes,char **prefix) 2032 { 2033 int ierr; 2034 2035 PetscFunctionBegin; 2036 PetscValidHeaderSpecific(snes,SNES_COOKIE); 2037 ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 2038 PetscFunctionReturn(0); 2039 } 2040 2041 #undef __FUNC__ 2042 #define __FUNC__ "SNESPrintHelp" 2043 /*@ 2044 SNESPrintHelp - Prints all options for the SNES component. 2045 2046 Input Parameter: 2047 . snes - the SNES context 2048 2049 Collective on SNES 2050 2051 Options Database Keys: 2052 $ -help, -h 2053 2054 .keywords: SNES, nonlinear, help 2055 2056 .seealso: SNESSetFromOptions() 2057 @*/ 2058 int SNESPrintHelp(SNES snes) 2059 { 2060 char p[64]; 2061 SNES_KSP_EW_ConvCtx *kctx; 2062 int ierr; 2063 2064 PetscFunctionBegin; 2065 PetscValidHeaderSpecific(snes,SNES_COOKIE); 2066 2067 PetscStrcpy(p,"-"); 2068 if (snes->prefix) PetscStrcat(p, snes->prefix); 2069 2070 kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 2071 2072 (*PetscHelpPrintf)(snes->comm,"SNES options ------------------------------------------------\n"); 2073 ierr = DLRegisterPrintTypes(snes->comm,stdout,snes->prefix,"snes_type",SNESList);CHKERRQ(ierr); 2074 (*PetscHelpPrintf)(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p); 2075 (*PetscHelpPrintf)(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its); 2076 (*PetscHelpPrintf)(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs); 2077 (*PetscHelpPrintf)(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol); 2078 (*PetscHelpPrintf)(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol); 2079 (*PetscHelpPrintf)(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol); 2080 (*PetscHelpPrintf)(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol); 2081 (*PetscHelpPrintf)(snes->comm," SNES Monitoring Options: Choose any of the following\n"); 2082 (*PetscHelpPrintf)(snes->comm," %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p); 2083 (*PetscHelpPrintf)(snes->comm," %ssnes_monitor: use default SNES convergence monitor, prints\n\ 2084 residual norm at each iteration.\n",p); 2085 (*PetscHelpPrintf)(snes->comm," %ssnes_smonitor: same as the above, but prints fewer digits of the\n\ 2086 residual norm for small residual norms. This is useful to conceal\n\ 2087 meaningless digits that may be different on different machines.\n",p); 2088 (*PetscHelpPrintf)(snes->comm," %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p); 2089 if (snes->type == SNES_NONLINEAR_EQUATIONS) { 2090 (*PetscHelpPrintf)(snes->comm, 2091 " Options for solving systems of nonlinear equations only:\n"); 2092 (*PetscHelpPrintf)(snes->comm," %ssnes_fd: use finite differences for Jacobian\n",p); 2093 (*PetscHelpPrintf)(snes->comm," %ssnes_mf: use matrix-free Jacobian\n",p); 2094 (*PetscHelpPrintf)(snes->comm," %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p); 2095 (*PetscHelpPrintf)(snes->comm," %ssnes_no_convergence_test: Do not test for convergence, always run to SNES max its\n",p); 2096 (*PetscHelpPrintf)(snes->comm," %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p); 2097 (*PetscHelpPrintf)(snes->comm, 2098 " %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version); 2099 (*PetscHelpPrintf)(snes->comm, 2100 " %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0); 2101 (*PetscHelpPrintf)(snes->comm, 2102 " %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max); 2103 (*PetscHelpPrintf)(snes->comm, 2104 " %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma); 2105 (*PetscHelpPrintf)(snes->comm, 2106 " %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha); 2107 (*PetscHelpPrintf)(snes->comm, 2108 " %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2); 2109 (*PetscHelpPrintf)(snes->comm, 2110 " %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold); 2111 } else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) { 2112 (*PetscHelpPrintf)(snes->comm, 2113 " Options for solving unconstrained minimization problems only:\n"); 2114 (*PetscHelpPrintf)(snes->comm," %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin); 2115 (*PetscHelpPrintf)(snes->comm," %ssnes_fd: use finite differences for Hessian\n",p); 2116 (*PetscHelpPrintf)(snes->comm," %ssnes_mf: use matrix-free Hessian\n",p); 2117 } 2118 (*PetscHelpPrintf)(snes->comm," Run program with -help %ssnes_type <method> for help on ",p); 2119 (*PetscHelpPrintf)(snes->comm,"a particular method\n"); 2120 if (snes->printhelp) { 2121 ierr = (*snes->printhelp)(snes,p);CHKERRQ(ierr); 2122 } 2123 PetscFunctionReturn(0); 2124 } 2125 2126 2127 2128 2129 2130