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