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