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