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