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