1 2 #ifndef lint 3 static char vcid[] = "$Id: snes.c,v 1.11 1995/04/27 20:17:07 bsmith Exp bsmith $"; 4 #endif 5 6 #include "draw.h" 7 #include "snesimpl.h" 8 #include "options.h" 9 10 /*@ 11 SNESSetFromOptions - Sets various SLES parameters from user options. 12 13 Input Parameter: 14 . snes - the SNES context 15 16 .keywords: SNES, nonlinear, set, options, database 17 18 .seealso: SNESPrintHelp() 19 @*/ 20 int SNESSetFromOptions(SNES snes) 21 { 22 SNESMETHOD method; 23 double tmp; 24 SLES sles; 25 VALIDHEADER(snes,SNES_COOKIE); 26 if (SNESGetMethodFromOptions(snes,&method)) { 27 SNESSetMethod(snes,method); 28 } 29 if (OptionsHasName(0,0,"-help")) SNESPrintHelp(snes); 30 if (OptionsGetDouble(0,snes->prefix,"-snes_stol",&tmp)) { 31 SNESSetSolutionTolerance(snes,tmp); 32 } 33 if (OptionsGetDouble(0,snes->prefix,"-snes_ttol",&tmp)) { 34 SNESSetTruncationTolerance(snes,tmp); 35 } 36 if (OptionsGetDouble(0,snes->prefix,"-snes_atol",&tmp)) { 37 SNESSetAbsoluteTolerance(snes,tmp); 38 } 39 if (OptionsGetDouble(0,snes->prefix,"-snes_rtol",&tmp)) { 40 SNESSetRelativeTolerance(snes,tmp); 41 } 42 OptionsGetInt(0,snes->prefix,"-snes_max_it",&snes->max_its); 43 OptionsGetInt(0,snes->prefix,"-snes_max_resid",&snes->max_resids); 44 if (OptionsHasName(0,snes->prefix,"-snes_monitor")) { 45 SNESSetMonitor(snes,SNESDefaultMonitor,0); 46 } 47 if (OptionsHasName(0,snes->prefix,"-snes_xmonitor")){ 48 int ierr,mytid = 0; 49 DrawLGCtx lg; 50 MPI_Initialized(&mytid); 51 if (mytid) MPI_Comm_rank(snes->comm,&mytid); 52 if (!mytid) { 53 ierr = SNESLGMonitorCreate(0,0,0,0,300,300,&lg); CHKERR(ierr); 54 SNESSetMonitor(snes,SNESLGMonitor,(void *)lg); 55 } 56 } 57 SNESGetSLES(snes,&sles); 58 SLESSetFromOptions(sles); 59 if (!snes->SetFromOptions) return 0; 60 return (*snes->SetFromOptions)(snes); 61 } 62 63 /*@ 64 SNESPrintHelp - Prints all options for the SNES component. 65 66 Input Parameter: 67 . snes - the SNES context 68 69 .keywords: SNES, nonlinear, help 70 71 .seealso: SLESSetFromOptions() 72 @*/ 73 int SNESPrintHelp(SNES snes) 74 { 75 int rank; 76 char *prefix = "-"; 77 if (snes->prefix) prefix = snes->prefix; 78 VALIDHEADER(snes,SNES_COOKIE); 79 if (!snes->PrintHelp) return 0; 80 MPI_Comm_rank(snes->comm,&rank); if (rank) return 0; 81 fprintf(stderr,"SNES options ----------------------------\n"); 82 fprintf(stderr,"%ssnes_method [ls] \n",prefix); 83 fprintf(stderr,"%ssnes_stol tol (default %g)\n",prefix,snes->xtol); 84 fprintf(stderr,"%ssnes_atol tol (default %g)\n",prefix,snes->atol); 85 fprintf(stderr,"%ssnes_rtol tol (default %g)\n",prefix,snes->rtol); 86 fprintf(stderr,"%ssnes_ttol tol (default %g)\n",prefix,snes->trunctol); 87 fprintf(stderr,"%ssnes_max_it its (default %d)\n",prefix,snes->max_its); 88 fprintf(stderr,"%ssnes_monitor\n",prefix); 89 return (*snes->PrintHelp)(snes); 90 } 91 /*@ 92 SNESSetApplicationContext - Sets the optional user-defined context for 93 the nonlinear solvers. 94 95 Input Parameters: 96 . snes - the SNES context 97 . usrP - optional user context 98 99 .keywords: SNES, nonlinear, set, application, context 100 101 .seealso: SNESGetApplicationContext() 102 @*/ 103 int SNESSetApplicationContext(SNES snes,void *usrP) 104 { 105 VALIDHEADER(snes,SNES_COOKIE); 106 snes->user = usrP; 107 return 0; 108 } 109 /*@ 110 SNESGetApplicationContext - Gets the user-defined context for the 111 nonlinear solvers. 112 113 Input Parameter: 114 . snes - SNES context 115 116 Output Parameter: 117 . usrP - user context 118 119 .keywords: SNES, nonlinear, get, application, context 120 121 .seealso: SNESSetApplicationContext() 122 @*/ 123 int SNESGetApplicationContext( SNES snes, void **usrP ) 124 { 125 VALIDHEADER(snes,SNES_COOKIE); 126 *usrP = snes->user; 127 return 0; 128 } 129 130 /*@ 131 SNESGetSLES - Returns the SLES context for a SNES solver. 132 133 Input Parameter: 134 . snes - the SNES context 135 136 Output Parameter: 137 . sles - the SLES context 138 139 Notes: 140 The user can then directly manipulate the SLES context to set various 141 options, etc. Likewise, the user can then manipulate the KSP and PC 142 contexts as well. 143 144 .keywords: SNES, nonlinear, get, SLES, context 145 146 .seealso: SLESGetPC(), SLESGetKSP() 147 @*/ 148 int SNESGetSLES(SNES snes,SLES *sles) 149 { 150 VALIDHEADER(snes,SNES_COOKIE); 151 *sles = snes->sles; 152 return 0; 153 } 154 155 /* -----------------------------------------------------------*/ 156 157 /*@ 158 SNESCreate - Creates a nonlinear solver context. 159 160 Input Parameter: 161 . comm - MPI communicator 162 163 Output Parameter: 164 . outsnes - the new SNES context 165 166 .keywords: SNES, nonlinear, create, context 167 168 .seealso: SNESSetUp(), SNESSolve(), SNESDestroy() 169 @*/ 170 int SNESCreate(MPI_Comm comm, SNES *outsnes) 171 { 172 int ierr; 173 SNES snes; 174 *outsnes = 0; 175 PETSCHEADERCREATE(snes,_SNES,SNES_COOKIE,SNES_NLS,comm); 176 PLogObjectCreate(snes); 177 snes->max_its = 50; 178 snes->max_resids = 1000; 179 snes->max_funcs = 1000; 180 snes->norm = 0.0; 181 snes->rtol = 1.e-8; 182 snes->atol = 1.e-10; 183 snes->xtol = 1.e-8; 184 snes->trunctol = 1.e-12; 185 snes->nresids = 0; 186 snes->Monitor = 0; 187 ierr = SLESCreate(comm,&snes->sles); CHKERR(ierr); 188 PLogObjectParent(snes,snes->sles) 189 *outsnes = snes; 190 return 0; 191 } 192 193 /* --------------------------------------------------------------- */ 194 /*@C 195 SNESSetFunction - Sets the residual evaluation routine and residual 196 vector for use by the SNES routines. 197 198 Input Parameters: 199 . snes - the SNES context 200 . func - residual evaluation routine 201 . resid_neg - indicator whether func evaluates f or -f. 202 If resid_neg is nonzero, then func evaluates -f; otherwise, 203 func evaluates f. 204 . ctx - optional user-defined function context 205 . r - vector to store residual 206 207 Calling sequence of func: 208 . func (Vec x, Vec f, void *ctx); 209 210 . x - input vector 211 . f - residual vector or its negative 212 . ctx - optional user-defined context for private data for the 213 residual evaluation routine (may be null) 214 215 Notes: 216 The Newton-like methods typically solve linear systems of the form 217 $ f'(x) x = -f(x), 218 $ where f'(x) denotes the Jacobian matrix and f(x) is the residual. 219 By setting resid_neg = 1, the user can supply -f(x) directly. 220 221 .keywords: SNES, nonlinear, set, residual 222 223 .seealso: SNESGetFunction(), SNESSetJacobian(), SNESSetSolution() 224 @*/ 225 int SNESSetFunction( SNES snes, Vec r, int (*func)(Vec,Vec,void*), 226 void *ctx,int rneg) 227 { 228 VALIDHEADER(snes,SNES_COOKIE); 229 snes->ComputeResidual = func; 230 snes->rsign = rneg; 231 snes->vec_res = r; 232 snes->resP = ctx; 233 return 0; 234 } 235 236 int SNESComputeFunction(SNES snes,Vec x, Vec y) 237 { 238 int ierr; 239 Scalar mone = -1.0; 240 ierr = (*snes->ComputeResidual)(x,y,snes->resP); CHKERR(ierr); 241 snes->nresids++; 242 if (!snes->rsign) { 243 ierr = VecScale(&mone,y); CHKERR(ierr); 244 } 245 return 0; 246 } 247 248 /*@ 249 SNESSetJacobian - Sets the function to compute Jacobian as well as the 250 location to store it. 251 252 Input Parameters: 253 . snes - the SNES context 254 . A - Jacobian matrix 255 . func - Jacobian evaluation routine 256 . ctx - optional user-defined context for private data for the 257 Jacobian evaluation routine (may be null) 258 259 Calling sequence of func: 260 . func (Vec x, Mat *A, Mat *B, int *flag,void *ctx); 261 262 . x - input vector 263 . A - Jacobian matrix 264 . B - preconditioner matrix, usually the same as A 265 . flag - same as options to SLESSetOperators(). Usually 0 or 266 $ MAT_SAME_NONZERO_PATTERN 267 . ctx - optional user-defined Jacobian context 268 269 Notes: 270 The function func() takes a Mat * as an argument rather than a Mat. 271 This is to allow the Jacobian code to replace it with a new matrix 272 when appropriate, for instance, if the nonzero structure is changing. 273 274 .keywords: SNES, nonlinear, set, Jacobian, matrix 275 276 .seealso: SNESSetFunction(), SNESSetSolution() 277 @*/ 278 int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(Vec,Mat*,Mat*,int*,void*), 279 void *ctx) 280 { 281 VALIDHEADER(snes,SNES_COOKIE); 282 snes->ComputeJacobian = func; 283 snes->jacP = ctx; 284 snes->jacobian = A; 285 snes->jacobian_pre = B; 286 return 0; 287 } 288 289 /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 290 291 /*@ 292 SNESSetUp - Sets up the internal data structures for the later use 293 of a nonlinear solver. Call SNESSetUp() after calling SNESCreate() 294 and optional routines of the form SNESSetXXX(), but before calling 295 SNESSolve(). 296 297 Input Parameter: 298 . snes - the SNES context 299 300 .keywords: SNES, nonlinear, setup 301 302 .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 303 @*/ 304 int SNESSetUp(SNES snes) 305 { 306 VALIDHEADER(snes,SNES_COOKIE); 307 return (*(snes)->Setup)( snes ); 308 } 309 310 /*@ 311 SNESDestroy - Destroys the nonlinear solver context that was created 312 with SNESCreate(). 313 314 Input Parameter: 315 . snes - the SNES context 316 317 .keywords: SNES, nonlinear, destroy 318 319 .seealso: SNESCreate(), SNESSetUp(), SNESSolve() 320 @*/ 321 int SNESDestroy(SNES snes) 322 { 323 VALIDHEADER(snes,SNES_COOKIE); 324 return (*(snes)->destroy)( (PetscObject) snes ); 325 } 326 327 /* ----------- Routines to set solver parameters ---------- */ 328 329 /*@ 330 SNESSetMaxIterations - Sets the maximum number of global iterations to use. 331 332 Input Parameters: 333 . snes - the SNES context 334 . maxits - maximum number of iterations to use 335 336 Options Database Key: 337 $ -snes_max_it maxits 338 339 Note: 340 The default maximum number of iterations is 50. 341 342 .keywords: SNES, nonlinear, set, maximum, iterations 343 344 .seealso: SNESSetMaxResidualEvaluations() 345 @*/ 346 int SNESSetMaxIterations(SNES snes,int maxits) 347 { 348 VALIDHEADER(snes,SNES_COOKIE); 349 snes->max_its = maxits; 350 return 0; 351 } 352 353 /*@ 354 SNESSetMaxResidualEvaluations - Sets the maximum number of residual 355 evaluations to use. 356 357 Input Parameters: 358 . snes - the SNES context 359 . maxr - maximum number of residual evaluations 360 361 Options Database Key: 362 $ -snes_max_resid maxr 363 364 Note: 365 The default maximum number of residual evaluations is 1000. 366 367 .keywords: SNES, nonlinear, set, maximum, residual, evaluations 368 369 .seealso: SNESSetMaxIterations() 370 @*/ 371 int SNESSetMaxResidualEvaluations(SNES snes,int maxr) 372 { 373 VALIDHEADER(snes,SNES_COOKIE); 374 snes->max_resids = maxr; 375 return 0; 376 } 377 378 /*@ 379 SNESSetRelativeTolerance - Sets the relative convergence tolerance. 380 381 Input Parameters: 382 . snes - the SNES context 383 . rtol - tolerance 384 385 Options Database Key: 386 $ -snes_rtol tol 387 388 .keywords: SNES, nonlinear, set, relative, convergence, tolerance 389 390 .seealso: SNESSetAbsoluteTolerance(), SNESSetSolutionTolerance(), 391 SNESSetTruncationTolerance() 392 @*/ 393 int SNESSetRelativeTolerance(SNES snes,double rtol) 394 { 395 VALIDHEADER(snes,SNES_COOKIE); 396 snes->rtol = rtol; 397 return 0; 398 } 399 400 /*@ 401 SNESSetAbsoluteTolerance - Sets the absolute convergence tolerance. 402 403 Input Parameters: 404 . snes - the SNES context 405 . atol - tolerance 406 407 Options Database Key: 408 $ -snes_atol tol 409 410 Notes: 411 $ The following convergence monitoring routines use atol 412 413 .keywords: SNES, nonlinear, set, absolute, convergence, tolerance 414 415 .seealso: SNESSetRelativeTolerance(), SNESSetSolutionTolerance(), 416 SNESSetTruncationTolerance() 417 @*/ 418 int SNESSetAbsoluteTolerance(SNES snes,double atol) 419 { 420 VALIDHEADER(snes,SNES_COOKIE); 421 snes->atol = atol; 422 return 0; 423 } 424 425 /*@ 426 SNESSetTruncationTolerance - Sets the tolerance that may be used by the 427 step routines to control the accuracy of the step computation. 428 429 Input Parameters: 430 . snes - the SNES context 431 . tol - tolerance 432 433 Options Database Key: 434 $ -snes_ttol tol 435 436 Notes: 437 If the step computation involves an application of the inverse 438 Jacobian (or Hessian), this parameter may be used to control the 439 accuracy of that application. In particular, this tolerance is used 440 by SNESKSPDefaultConverged() and SNESKSPQuadraticConverged() to determine 441 the minimum convergence tolerance for the iterative linear solvers. 442 443 .keywords: SNES, nonlinear, set, truncation, tolerance 444 445 .seealso: SNESSetRelativeTolerance(), SNESSetSolutionTolerance(), 446 SNESSetAbsoluteTolerance() 447 @*/ 448 int SNESSetTruncationTolerance(SNES snes,double tol) 449 { 450 VALIDHEADER(snes,SNES_COOKIE); 451 snes->trunctol = tol; 452 return 0; 453 } 454 455 /*@ 456 SNESSetSolutionTolerance - Sets the convergence tolerance in terms of 457 the norm of the change in the solution between steps. 458 459 Input Parameters: 460 . snes - the SNES context 461 . tol - tolerance 462 463 Options Database Key: 464 $ -snes_stol tol 465 466 .keywords: SNES, nonlinear, set, solution, tolerance 467 468 .seealso: SNESSetTruncationTolerance(), SNESSetRelativeTolerance(), 469 SNESSetAbsoluteTolerance() 470 @*/ 471 int SNESSetSolutionTolerance( SNES snes, double tol ) 472 { 473 VALIDHEADER(snes,SNES_COOKIE); 474 snes->xtol = tol; 475 return 0; 476 } 477 478 /* ---------- Routines to set various aspects of nonlinear solver --------- */ 479 480 /*@ 481 SNESSetSolution - Sets the initial guess routine and solution vector 482 for use by the SNES routines. 483 484 Input Parameters: 485 . snes - the SNES context 486 . x - the solution vector 487 . func - optional routine to compute an initial guess (may be null) 488 . ctx - optional user-defined context for private data for the 489 initial guess routine (may be null) 490 491 Calling sequence of func: 492 int guess(Vec x, void *ctx) 493 494 . x - input vector 495 . ctx - optional user-defined initial guess context 496 497 Note: 498 If no initial guess routine is indicated, an initial guess of zero 499 will be used. 500 501 .keywords: SNES, nonlinear, set, solution, initial guess 502 503 .seealso: SNESGetSolution(), SNESSetJacobian(), SNESSetFunction() 504 @*/ 505 int SNESSetSolution(SNES snes,Vec x,int (*func)(Vec,void*),void *ctx) 506 { 507 VALIDHEADER(snes,SNES_COOKIE); 508 snes->vec_sol = x; 509 snes->ComputeInitialGuess = func; 510 snes->gusP = ctx; 511 return 0; 512 } 513 514 /* ------------ Routines to set performance monitoring options ----------- */ 515 516 /*@C 517 SNESSetMonitor - Sets the function that is to be used at every 518 iteration of the nonlinear solver to display the iteration's 519 progress. 520 521 Input Parameters: 522 . snes - the SNES context 523 . func - monitoring routine 524 . mctx - optional user-defined context for private data for the 525 monitor routine (may be null) 526 527 Calling sequence of func: 528 int func((SNES snes,int its, Vec x,Vec f,double fnorm,void *mctx) 529 530 . snes - the SNES context 531 . its - iteration number 532 . x - current iterate 533 . f - current residual (+/-). f is either the residual or its negative, 534 depending on the user's preference, as set with SNESSetFunction() 535 . fnorm - 2-norm residual value (may be estimated) 536 . mctx - optional monitoring context 537 538 .keywords: SNES, nonlinear, set, monitor 539 540 .seealso: SNESDefaultMonitor() 541 @*/ 542 int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*), 543 void *mctx ) 544 { 545 snes->Monitor = func; 546 snes->monP = (void*)mctx; 547 return 0; 548 } 549 550 /*@C 551 SNESSetConvergenceTest - Sets the function that is to be used 552 to test for convergence of the nonlinear iterative solution. 553 554 Input Parameters: 555 . snes - the SNES context 556 . func - routine to test for convergence 557 . cctx - optional context for private data for the convergence routine 558 (may be null) 559 560 Calling sequence of func: 561 int func (SNES snes,double xnorm,double pnorm,double fnorm, 562 void *cctx) 563 564 . snes - the SNES context 565 . xnorm - 2-norm of current iterate 566 . pnorm - 2-norm of current step 567 . fnorm - 2-norm of residual 568 . cctx - optional convergence context 569 570 .keywords: SNES, nonlinear, set, convergence, test 571 572 .seealso: SNESDefaultConverged() 573 @*/ 574 int SNESSetConvergenceTest(SNES nlP, 575 int (*func)(SNES,double,double,double,void*),void *cctx) 576 { 577 (nlP)->Converged = func; 578 (nlP)->cnvP = cctx; 579 return 0; 580 } 581 582 /*@ 583 SNESScaleStep - Scales a step so that its length is less than the 584 positive parameter delta. 585 586 Input Parameters: 587 . snes - the SNES context 588 . y - approximate solution of linear system 589 . fnorm - 2-norm of current residual 590 . delta - trust region size 591 592 Output Parameters: 593 . gpnorm - predicted residual norm at the new point, assuming local 594 linearization. The value is zero if the step lies within the trust 595 region, and exceeds zero otherwise. 596 . ynorm - 2-norm of the step 597 598 Note: 599 For non-trust region methods such as SNES_NLS, the parameter delta 600 is set to be the maximum allowable step size. 601 602 .keywords: SNES, nonlinear, scale, step 603 @*/ 604 int SNESScaleStep(SNES snes,Vec y,double *fnorm,double *delta, 605 double *gpnorm,double *ynorm) 606 { 607 double norm; 608 Scalar cnorm; 609 VecNorm(y, &norm ); 610 if (norm > *delta) { 611 norm = *delta/norm; 612 *gpnorm = (1.0 - norm)*(*fnorm); 613 cnorm = norm; 614 VecScale( &cnorm, y ); 615 *ynorm = *delta; 616 } else { 617 *gpnorm = 0.0; 618 *ynorm = norm; 619 } 620 return 0; 621 } 622 623 /*@ 624 SNESSolve - Solves a nonlinear system. Call SNESSolve after calling 625 SNESCreate(), optional routines of the form SNESSetXXX(), and SNESSetUp(). 626 627 Input Parameter: 628 . snes - the SNES context 629 630 Output Parameter: 631 its - number of iterations until termination 632 633 .keywords: SNES, nonlinear, solve 634 635 .seealso: SNESCreate(), SNESSetUp(), SNESDestroy() 636 @*/ 637 int SNESSolve(SNES snes,int *its) 638 { 639 int ierr; 640 PLogEventBegin(SNES_Solve,snes,0,0,0); 641 ierr = (*(snes)->Solver)( snes,its ); CHKERR(ierr); 642 PLogEventEnd(SNES_Solve,snes,0,0,0); 643 return 0; 644 } 645 646 /* --------- Internal routines for SNES Package --------- */ 647 648 /* 649 SNESComputeInitialGuess - Manages computation of initial approximation. 650 */ 651 int SNESComputeInitialGuess( SNES nlP,Vec x ) 652 { 653 int ierr; 654 Scalar zero = 0.0; 655 if (nlP->ComputeInitialGuess) { 656 ierr = (*nlP->ComputeInitialGuess)( x, nlP->gusP); CHKERR(ierr); 657 } 658 else VecSet(&zero, x ); 659 return 0; 660 } 661 662 /* ------------------------------------------------------------------ */ 663 664 665 #include "sys/nreg.h" 666 NRList *__NLList; 667 668 /*@ 669 SNESSetMethod - Sets the method for the nonlinear solver. 670 671 Input Parameters: 672 . snes - the SNES context 673 . method - choose from 674 675 Possible methods: 676 $ SNES_NLS - Newton's method with line search 677 $ SNES_NTR - Newton's method with trust region 678 @*/ 679 int SNESSetMethod( SNES snes, SNESMETHOD method) 680 { 681 int (*r)(SNES); 682 VALIDHEADER(snes,SNES_COOKIE); 683 /* Get the function pointers for the iterative method requested */ 684 if (!__NLList) {SNESRegisterAll();} 685 if (!__NLList) {SETERR(1,"Could not acquire list of SNES methods"); } 686 r = (int (*)(SNES))NRFindRoutine( __NLList, (int)method, (char *)0 ); 687 if (!r) {SETERR(1,"Unknown SNES method");} 688 return (*r)(snes); 689 } 690 691 /* --------------------------------------------------------------------- */ 692 /*@ 693 SNESRegister - Adds the method to the nonlinear solver package, given 694 a function pointer and a nonlinear solver name of the type SNESMETHOD. 695 696 Input Parameters: 697 . name - for instance SNES_NLS, SNES_NTR, ... 698 . sname - corresponding string for name 699 . create - routine to create method context 700 701 .keywords: SNES, nonlinear, register 702 703 .seealso: SNESRegisterAll(), SNESRegisterDestroy() 704 @*/ 705 int SNESRegister(int name, char *sname, int (*create)(SNES)) 706 { 707 int ierr; 708 if (!__NLList) {ierr = NRCreate(&__NLList); CHKERR(ierr);} 709 NRRegister( __NLList, name, sname, (int (*)(void*))create ); 710 return 0; 711 } 712 /* --------------------------------------------------------------------- */ 713 /*@ 714 SNESRegisterDestroy - Frees the list of nonlinear solvers that were 715 registered by SNESRegister(). 716 717 .keywords: SNES, nonlinear, register, destroy 718 719 .seealso: SNESRegisterAll(), SNESRegisterAll() 720 @*/ 721 int SNESRegisterDestroy() 722 { 723 if (__NLList) { 724 NRDestroy( __NLList ); 725 __NLList = 0; 726 } 727 return 0; 728 } 729 #include "options.h" 730 /*@C 731 SNESGetMethodFromOptions - Sets the selected method from the options 732 database. 733 734 Input parameters: 735 . ctx - the SNES context 736 737 Output Parameter: 738 . method - solver method 739 740 Returns: 741 Returns 1 if the method is found; 0 otherwise. 742 743 Options Database Key: 744 $ -snes_method method 745 746 .keywords: SNES, nonlinear, options, database, get, method 747 748 .seealso: SNESGetMethodName() 749 @*/ 750 int SNESGetMethodFromOptions(SNES ctx,SNESMETHOD *method) 751 { 752 char sbuf[50]; 753 if (OptionsGetString(0,ctx->prefix,"-snes_method", sbuf, 50 )) { 754 if (!__NLList) SNESRegisterAll(); 755 *method = (SNESMETHOD)NRFindID( __NLList, sbuf ); 756 return 1; 757 } 758 return 0; 759 } 760 761 /*@C 762 SNESGetMethodName - Gets the SNES method name (as a string) from 763 the method type. 764 765 Input Parameter: 766 . method - SNES method 767 768 Output Parameter: 769 . name - name of SNES method 770 771 .keywords: SNES, nonlinear, get, method, name 772 773 .seealso: SNESGetMethodFromOptions() 774 @*/ 775 int SNESGetMethodName(SNESMETHOD method,char **name) 776 { 777 if (!__NLList) SNESRegisterAll(); 778 *name = NRFindName( __NLList, (int) method ); 779 return 0; 780 } 781 782 #include <stdio.h> 783 /*@C 784 SNESPrintMethods - Prints the SNES methods available from the options 785 database. 786 787 Input Parameters: 788 . prefix - prefix (usually "-") 789 . name - the options database name (by default "snesmethod") 790 791 .keywords: SNES, nonlinear, print, methods, options, database 792 793 .seealso: SNESPrintHelp() 794 @*/ 795 int SNESPrintMethods(char* prefix,char *name) 796 { 797 FuncList *entry; 798 if (!__NLList) {SNESRegisterAll();} 799 entry = __NLList->head; 800 fprintf(stderr," %s%s (one of)",prefix,name); 801 while (entry) { 802 fprintf(stderr," %s",entry->name); 803 entry = entry->next; 804 } 805 fprintf(stderr,"\n"); 806 return 0; 807 } 808 809 /*@ 810 SNESGetSolution - Returns the vector where the approximate solution is 811 stored. 812 813 Input Parameter: 814 . snes - the SNES context 815 816 Output Parameter: 817 . x - the solution 818 819 .keywords: SNES, nonlinear, get, solution 820 821 .seealso: SNESSetSolution(), SNESGetFunction() 822 @*/ 823 int SNESGetSolution(SNES snes,Vec *x) 824 { 825 VALIDHEADER(snes,SNES_COOKIE); 826 *x = snes->vec_sol; 827 return 0; 828 } 829 830 /*@ 831 SNESGetFunction - Returns the vector where the residual is 832 stored. Actually usually returns the vector where the negative of 833 the residual is stored. 834 835 Input Parameter: 836 . snes - the SNES context 837 838 Output Parameter: 839 . r - the residual (or its negative) 840 841 .keywords: SNES, nonlinear, get residual 842 843 .seealso: SNESSetFunction(), SNESGetSolution() 844 @*/ 845 int SNESGetFunction(SNES snes,Vec *r) 846 { 847 VALIDHEADER(snes,SNES_COOKIE); 848 *r = snes->vec_res; 849 return 0; 850 } 851 852 853 854 855 856 857