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