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