1 #ifndef lint 2 static char vcid[] = "$Id: snes.c,v 1.86 1996/09/08 22:52:35 curfman Exp curfman $"; 3 #endif 4 5 #include "draw.h" /*I "draw.h" I*/ 6 #include "src/snes/snesimpl.h" /*I "snes.h" I*/ 7 #include "src/sys/nreg.h" 8 #include "pinclude/pviewer.h" 9 #include <math.h> 10 11 extern int SNESGetTypeFromOptions_Private(SNES,SNESType*,int*); 12 extern int SNESPrintTypes_Private(MPI_Comm,char*,char*); 13 14 /*@ 15 SNESView - Prints the SNES data structure. 16 17 Input Parameters: 18 . SNES - the SNES context 19 . viewer - visualization context 20 21 Options Database Key: 22 $ -snes_view : calls SNESView() at end of SNESSolve() 23 24 Notes: 25 The available visualization contexts include 26 $ VIEWER_STDOUT_SELF - standard output (default) 27 $ VIEWER_STDOUT_WORLD - synchronized standard 28 $ output where only the first processor opens 29 $ the file. All other processors send their 30 $ data to the first processor to print. 31 32 The user can open alternative vistualization contexts with 33 $ ViewerFileOpenASCII() - output to a specified file 34 35 .keywords: SNES, view 36 37 .seealso: ViewerFileOpenASCII() 38 @*/ 39 int SNESView(SNES snes,Viewer viewer) 40 { 41 SNES_KSP_EW_ConvCtx *kctx; 42 FILE *fd; 43 int ierr; 44 SLES sles; 45 char *method; 46 ViewerType vtype; 47 48 PetscValidHeaderSpecific(snes,SNES_COOKIE); 49 if (viewer) {PetscValidHeader(viewer);} 50 else { viewer = VIEWER_STDOUT_SELF; } 51 52 ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 53 if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 54 ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 55 PetscFPrintf(snes->comm,fd,"SNES Object:\n"); 56 SNESGetType(snes,PETSC_NULL,&method); 57 PetscFPrintf(snes->comm,fd," method: %s\n",method); 58 if (snes->view) (*snes->view)((PetscObject)snes,viewer); 59 PetscFPrintf(snes->comm,fd, 60 " maximum iterations=%d, maximum function evaluations=%d\n", 61 snes->max_its,snes->max_funcs); 62 PetscFPrintf(snes->comm,fd, 63 " tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n", 64 snes->rtol, snes->atol, snes->trunctol, snes->xtol); 65 if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) 66 PetscFPrintf(snes->comm,fd," min function tolerance=%g\n",snes->fmin); 67 if (snes->ksp_ewconv) { 68 kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 69 if (kctx) { 70 PetscFPrintf(snes->comm,fd, 71 " Eisenstat-Walker computation of KSP relative tolerance (version %d)\n", 72 kctx->version); 73 PetscFPrintf(snes->comm,fd, 74 " rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0, 75 kctx->rtol_max,kctx->threshold); 76 PetscFPrintf(snes->comm,fd," gamma=%g, alpha=%g, alpha2=%g\n", 77 kctx->gamma,kctx->alpha,kctx->alpha2); 78 } 79 } 80 } else if (vtype == STRING_VIEWER) { 81 SNESGetType(snes,PETSC_NULL,&method); 82 ViewerStringSPrintf(viewer," %-3.3s",method); 83 } 84 SNESGetSLES(snes,&sles); 85 ierr = SLESView(sles,viewer); CHKERRQ(ierr); 86 return 0; 87 } 88 89 /*@ 90 SNESSetFromOptions - Sets various SNES and SLES parameters from user options. 91 92 Input Parameter: 93 . snes - the SNES context 94 95 .keywords: SNES, nonlinear, set, options, database 96 97 .seealso: SNESPrintHelp(), SNESSetOptionsPrefix() 98 @*/ 99 int SNESSetFromOptions(SNES snes) 100 { 101 SNESType method; 102 double tmp; 103 SLES sles; 104 int ierr, flg, mset = 0; 105 int version = PETSC_DEFAULT; 106 double rtol_0 = PETSC_DEFAULT; 107 double rtol_max = PETSC_DEFAULT; 108 double gamma2 = PETSC_DEFAULT; 109 double alpha = PETSC_DEFAULT; 110 double alpha2 = PETSC_DEFAULT; 111 double threshold = PETSC_DEFAULT; 112 113 PetscValidHeaderSpecific(snes,SNES_COOKIE); 114 if (snes->setup_called) SETERRQ(1,"SNESSetFromOptions:Must call prior to SNESSetUp"); 115 ierr = SNESGetTypeFromOptions_Private(snes,&method,&flg); CHKERRQ(ierr); 116 if (flg) { 117 ierr = SNESSetType(snes,method); CHKERRQ(ierr); 118 } 119 else if (!snes->set_method_called) { 120 if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 121 ierr = SNESSetType(snes,SNES_EQ_LS); CHKERRQ(ierr); 122 } 123 else { 124 ierr = SNESSetType(snes,SNES_UM_TR); CHKERRQ(ierr); 125 } 126 } 127 128 ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 129 if (flg) { SNESPrintHelp(snes); } 130 ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg); CHKERRQ(ierr); 131 if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);} 132 ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg); CHKERRQ(ierr); 133 if (flg) {SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);} 134 ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg); CHKERRQ(ierr); 135 if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT); } 136 ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, &flg); CHKERRQ(ierr); 137 ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, &flg);CHKERRQ(ierr); 138 ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg); CHKERRQ(ierr); 139 if (flg) { SNESSetTrustRegionTolerance(snes,tmp); } 140 ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg); CHKERRQ(ierr); 141 if (flg) { SNESSetMinimizationFunctionTolerance(snes,tmp); } 142 ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg); CHKERRQ(ierr); 143 if (flg) { snes->ksp_ewconv = 1; } 144 ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,&flg); CHKERRQ(ierr); 145 ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,&flg); CHKERRQ(ierr); 146 ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,&flg);CHKERRQ(ierr); 147 ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,&flg); CHKERRQ(ierr); 148 ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,&flg); CHKERRQ(ierr); 149 ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,&flg); CHKERRQ(ierr); 150 ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,&flg);CHKERRQ(ierr); 151 ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha, 152 alpha2,threshold); CHKERRQ(ierr); 153 ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg); CHKERRQ(ierr); 154 if (flg) {SNESSetMonitor(snes,SNESDefaultMonitor,0); mset = 1;} 155 ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg); CHKERRQ(ierr); 156 if (flg) { 157 if (mset) SETERRQ(1,"SNESSetFromOptions:Monitor for SNES is already set"); 158 SNESSetMonitor(snes,SNESDefaultSMonitor,0); mset = 1; 159 } 160 ierr = OptionsHasName(snes->prefix,"-snes_xmonitor",&flg); CHKERRQ(ierr); 161 if (flg) { 162 int rank = 0; 163 DrawLG lg; 164 if (mset) SETERRQ(1,"SNESSetFromOptions:Monitor for SNES is already set"); 165 MPI_Initialized(&rank); 166 if (rank) MPI_Comm_rank(snes->comm,&rank); 167 if (!rank) { 168 ierr = SNESLGMonitorCreate(0,0,0,0,300,300,&lg); CHKERRQ(ierr); 169 ierr = SNESSetMonitor(snes,SNESLGMonitor,(void *)lg); CHKERRQ(ierr); 170 snes->xmonitor = lg; 171 PLogObjectParent(snes,lg); 172 } 173 mset = 1; 174 } 175 ierr = OptionsHasName(snes->prefix,"-snes_fd", &flg); CHKERRQ(ierr); 176 if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) { 177 ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre, 178 SNESDefaultComputeJacobian,snes->funP); CHKERRQ(ierr); 179 PLogInfo(snes, 180 "SNESSetFromOptions: Setting default finite difference Jacobian matrix\n"); 181 } 182 else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 183 ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre, 184 SNESDefaultComputeHessian,snes->funP); CHKERRQ(ierr); 185 PLogInfo(snes, 186 "SNESSetFromOptions: Setting default finite difference Hessian matrix\n"); 187 } 188 ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 189 ierr = SLESSetFromOptions(sles); CHKERRQ(ierr); 190 if (!snes->setfromoptions) return 0; 191 return (*snes->setfromoptions)(snes); 192 } 193 194 /*@ 195 SNESPrintHelp - Prints all options for the SNES component. 196 197 Input Parameter: 198 . snes - the SNES context 199 200 Options Database Keys: 201 $ -help, -h 202 203 .keywords: SNES, nonlinear, help 204 205 .seealso: SNESSetFromOptions() 206 @*/ 207 int SNESPrintHelp(SNES snes) 208 { 209 char p[64]; 210 SNES_KSP_EW_ConvCtx *kctx; 211 212 PetscValidHeaderSpecific(snes,SNES_COOKIE); 213 214 PetscStrcpy(p,"-"); 215 if (snes->prefix) PetscStrcat(p, snes->prefix); 216 217 kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 218 219 PetscPrintf(snes->comm,"SNES options ------------------------------------------------\n"); 220 SNESPrintTypes_Private(snes->comm,p,"snes_type"); 221 PetscPrintf(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p); 222 PetscPrintf(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its); 223 PetscPrintf(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs); 224 PetscPrintf(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol); 225 PetscPrintf(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol); 226 PetscPrintf(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol); 227 PetscPrintf(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol); 228 PetscPrintf(snes->comm," SNES Monitoring Options: Choose 1 of the following\n"); 229 PetscPrintf(snes->comm," %ssnes_monitor: use default SNES convergence monitor\n",p); 230 PetscPrintf(snes->comm," %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p); 231 if (snes->type == SNES_NONLINEAR_EQUATIONS) { 232 PetscPrintf(snes->comm, 233 " Options for solving systems of nonlinear equations only:\n"); 234 PetscPrintf(snes->comm," %ssnes_fd: use finite differences for Jacobian\n",p); 235 PetscPrintf(snes->comm," %ssnes_mf: use matrix-free Jacobian\n",p); 236 PetscPrintf(snes->comm," %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p); 237 PetscPrintf(snes->comm," %ssnes_mf_err: relative error in function evaluation for matrix-free Jacobian\n",p); 238 PetscPrintf(snes->comm," %ssnes_mf_umin: minimum iterate parameter for matrix-free Jacobian\n",p); 239 PetscPrintf(snes->comm," %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p); 240 PetscPrintf(snes->comm, 241 " %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version); 242 PetscPrintf(snes->comm, 243 " %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0); 244 PetscPrintf(snes->comm, 245 " %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max); 246 PetscPrintf(snes->comm, 247 " %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma); 248 PetscPrintf(snes->comm, 249 " %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha); 250 PetscPrintf(snes->comm, 251 " %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2); 252 PetscPrintf(snes->comm, 253 " %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold); 254 } 255 else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) { 256 PetscPrintf(snes->comm, 257 " Options for solving unconstrained minimization problems only:\n"); 258 PetscPrintf(snes->comm," %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin); 259 PetscPrintf(snes->comm," %ssnes_fd: use finite differences for Hessian\n",p); 260 PetscPrintf(snes->comm," %ssnes_mf: use matrix-free Hessian\n",p); 261 PetscPrintf(snes->comm," %ssnes_mf_err: relative error in gradient evaluation for matrix-free Hessian\n",p); 262 PetscPrintf(snes->comm," %ssnes_mf_umin: minimum iterate parameter for matrix-free Hessian\n",p); 263 } 264 PetscPrintf(snes->comm," Run program with -help %ssnes_type <method> for help on ",p); 265 PetscPrintf(snes->comm,"a particular method\n"); 266 if (snes->printhelp) (*snes->printhelp)(snes,p); 267 return 0; 268 } 269 /*@ 270 SNESSetApplicationContext - Sets the optional user-defined context for 271 the nonlinear solvers. 272 273 Input Parameters: 274 . snes - the SNES context 275 . usrP - optional user context 276 277 .keywords: SNES, nonlinear, set, application, context 278 279 .seealso: SNESGetApplicationContext() 280 @*/ 281 int SNESSetApplicationContext(SNES snes,void *usrP) 282 { 283 PetscValidHeaderSpecific(snes,SNES_COOKIE); 284 snes->user = usrP; 285 return 0; 286 } 287 288 /*@C 289 SNESGetApplicationContext - Gets the user-defined context for the 290 nonlinear solvers. 291 292 Input Parameter: 293 . snes - SNES context 294 295 Output Parameter: 296 . usrP - user context 297 298 .keywords: SNES, nonlinear, get, application, context 299 300 .seealso: SNESSetApplicationContext() 301 @*/ 302 int SNESGetApplicationContext( SNES snes, void **usrP ) 303 { 304 PetscValidHeaderSpecific(snes,SNES_COOKIE); 305 *usrP = snes->user; 306 return 0; 307 } 308 309 /*@ 310 SNESGetIterationNumber - Gets the current iteration number of the 311 nonlinear solver. 312 313 Input Parameter: 314 . snes - SNES context 315 316 Output Parameter: 317 . iter - iteration number 318 319 .keywords: SNES, nonlinear, get, iteration, number 320 @*/ 321 int SNESGetIterationNumber(SNES snes,int* iter) 322 { 323 PetscValidHeaderSpecific(snes,SNES_COOKIE); 324 PetscValidIntPointer(iter); 325 *iter = snes->iter; 326 return 0; 327 } 328 329 /*@ 330 SNESGetFunctionNorm - Gets the norm of the current function that was set 331 with SNESSSetFunction(). 332 333 Input Parameter: 334 . snes - SNES context 335 336 Output Parameter: 337 . fnorm - 2-norm of function 338 339 Note: 340 SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only. 341 A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is 342 SNESGetGradientNorm(). 343 344 .keywords: SNES, nonlinear, get, function, norm 345 346 .seealso: SNESSetFunction() 347 @*/ 348 int SNESGetFunctionNorm(SNES snes,Scalar *fnorm) 349 { 350 PetscValidHeaderSpecific(snes,SNES_COOKIE); 351 PetscValidScalarPointer(fnorm); 352 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 353 SETERRQ(1,"SNESGetFunctionNorm:For SNES_NONLINEAR_EQUATIONS only"); 354 } 355 *fnorm = snes->norm; 356 return 0; 357 } 358 359 /*@ 360 SNESGetGradientNorm - Gets the norm of the current gradient that was set 361 with SNESSSetGradient(). 362 363 Input Parameter: 364 . snes - SNES context 365 366 Output Parameter: 367 . fnorm - 2-norm of gradient 368 369 Note: 370 SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION 371 methods only. A related routine for SNES_NONLINEAR_EQUATIONS methods 372 is SNESGetFunctionNorm(). 373 374 .keywords: SNES, nonlinear, get, gradient, norm 375 376 .seelso: SNESSetGradient() 377 @*/ 378 int SNESGetGradientNorm(SNES snes,Scalar *gnorm) 379 { 380 PetscValidHeaderSpecific(snes,SNES_COOKIE); 381 PetscValidScalarPointer(gnorm); 382 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 383 SETERRQ(1,"SNESGetGradientNorm:For SNES_UNCONSTRAINED_MINIMIZATION only"); 384 } 385 *gnorm = snes->norm; 386 return 0; 387 } 388 389 /*@ 390 SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps 391 attempted by the nonlinear solver. 392 393 Input Parameter: 394 . snes - SNES context 395 396 Output Parameter: 397 . nfails - number of unsuccessful steps attempted 398 399 .keywords: SNES, nonlinear, get, number, unsuccessful, steps 400 @*/ 401 int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails) 402 { 403 PetscValidHeaderSpecific(snes,SNES_COOKIE); 404 PetscValidIntPointer(nfails); 405 *nfails = snes->nfailures; 406 return 0; 407 } 408 /*@C 409 SNESGetSLES - Returns the SLES context for a SNES solver. 410 411 Input Parameter: 412 . snes - the SNES context 413 414 Output Parameter: 415 . sles - the SLES context 416 417 Notes: 418 The user can then directly manipulate the SLES context to set various 419 options, etc. Likewise, the user can then extract and manipulate the 420 KSP and PC contexts as well. 421 422 .keywords: SNES, nonlinear, get, SLES, context 423 424 .seealso: SLESGetPC(), SLESGetKSP() 425 @*/ 426 int SNESGetSLES(SNES snes,SLES *sles) 427 { 428 PetscValidHeaderSpecific(snes,SNES_COOKIE); 429 *sles = snes->sles; 430 return 0; 431 } 432 /* -----------------------------------------------------------*/ 433 /*@C 434 SNESCreate - Creates a nonlinear solver context. 435 436 Input Parameter: 437 . comm - MPI communicator 438 . type - type of method, one of 439 $ SNES_NONLINEAR_EQUATIONS 440 $ (for systems of nonlinear equations) 441 $ SNES_UNCONSTRAINED_MINIMIZATION 442 $ (for unconstrained minimization) 443 444 Output Parameter: 445 . outsnes - the new SNES context 446 447 Options Database Key: 448 $ -snes_mf - use default matrix-free Jacobian-vector products, 449 $ and no preconditioning matrix 450 $ -snes_mf_operator - use default matrix-free Jacobian-vector 451 $ products, and a user-provided preconditioning matrix 452 $ as set by SNESSetJacobian() 453 $ -snes_fd - use (slow!) finite differences to compute Jacobian 454 455 .keywords: SNES, nonlinear, create, context 456 457 .seealso: SNESSolve(), SNESDestroy() 458 @*/ 459 int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes) 460 { 461 int ierr; 462 SNES snes; 463 SNES_KSP_EW_ConvCtx *kctx; 464 465 *outsnes = 0; 466 if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS) 467 SETERRQ(1,"SNESCreate:incorrect method type"); 468 PetscHeaderCreate(snes,_SNES,SNES_COOKIE,SNES_UNKNOWN_METHOD,comm); 469 PLogObjectCreate(snes); 470 snes->max_its = 50; 471 snes->max_funcs = 1000; 472 snes->norm = 0.0; 473 if (type == SNES_UNCONSTRAINED_MINIMIZATION) { 474 snes->rtol = 1.e-8; 475 snes->ttol = 0.0; 476 snes->atol = 1.e-10; 477 } 478 else { 479 snes->rtol = 1.e-8; 480 snes->ttol = 0.0; 481 snes->atol = 1.e-50; 482 } 483 snes->xtol = 1.e-8; 484 snes->trunctol = 1.e-12; /* no longer used */ 485 snes->nfuncs = 0; 486 snes->nfailures = 0; 487 snes->monitor = 0; 488 snes->data = 0; 489 snes->view = 0; 490 snes->computeumfunction = 0; 491 snes->umfunP = 0; 492 snes->fc = 0; 493 snes->deltatol = 1.e-12; 494 snes->fmin = -1.e30; 495 snes->method_class = type; 496 snes->set_method_called = 0; 497 snes->setup_called = 0; 498 snes->ksp_ewconv = 0; 499 snes->type = -1; 500 snes->xmonitor = 0; 501 snes->vwork = 0; 502 snes->nwork = 0; 503 504 /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 505 kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx); 506 snes->kspconvctx = (void*)kctx; 507 kctx->version = 2; 508 kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 509 this was too large for some test cases */ 510 kctx->rtol_last = 0; 511 kctx->rtol_max = .9; 512 kctx->gamma = 1.0; 513 kctx->alpha2 = .5*(1.0 + sqrt(5.0)); 514 kctx->alpha = kctx->alpha2; 515 kctx->threshold = .1; 516 kctx->lresid_last = 0; 517 kctx->norm_last = 0; 518 519 ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr); 520 PLogObjectParent(snes,snes->sles) 521 522 *outsnes = snes; 523 return 0; 524 } 525 526 /* --------------------------------------------------------------- */ 527 /*@C 528 SNESSetFunction - Sets the function evaluation routine and function 529 vector for use by the SNES routines in solving systems of nonlinear 530 equations. 531 532 Input Parameters: 533 . snes - the SNES context 534 . func - function evaluation routine 535 . ctx - optional user-defined function context 536 . r - vector to store function value 537 538 Calling sequence of func: 539 . func (SNES, Vec x, Vec f, void *ctx); 540 541 . x - input vector 542 . f - vector function 543 . ctx - optional user-defined context for private data for the 544 function evaluation routine (may be null) 545 546 Notes: 547 The Newton-like methods typically solve linear systems of the form 548 $ f'(x) x = -f(x), 549 $ where f'(x) denotes the Jacobian matrix and f(x) is the function. 550 551 SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 552 Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 553 SNESSetMinimizationFunction() and SNESSetGradient(); 554 555 .keywords: SNES, nonlinear, set, function 556 557 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 558 @*/ 559 int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx) 560 { 561 PetscValidHeaderSpecific(snes,SNES_COOKIE); 562 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1, 563 "SNESSetFunction:For SNES_NONLINEAR_EQUATIONS only"); 564 snes->computefunction = func; 565 snes->vec_func = snes->vec_func_always = r; 566 snes->funP = ctx; 567 return 0; 568 } 569 570 /*@ 571 SNESComputeFunction - Computes the function that has been set with 572 SNESSetFunction(). 573 574 Input Parameters: 575 . snes - the SNES context 576 . x - input vector 577 578 Output Parameter: 579 . y - function vector, as set by SNESSetFunction() 580 581 Notes: 582 SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 583 Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 584 SNESComputeMinimizationFunction() and SNESComputeGradient(); 585 586 .keywords: SNES, nonlinear, compute, function 587 588 .seealso: SNESSetFunction(), SNESGetFunction() 589 @*/ 590 int SNESComputeFunction(SNES snes,Vec x, Vec y) 591 { 592 int ierr; 593 594 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 595 SETERRQ(1,"SNESComputeFunction: For SNES_NONLINEAR_EQUATIONS only"); 596 PLogEventBegin(SNES_FunctionEval,snes,x,y,0); 597 ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 598 PLogEventEnd(SNES_FunctionEval,snes,x,y,0); 599 return 0; 600 } 601 602 /*@C 603 SNESSetMinimizationFunction - Sets the function evaluation routine for 604 unconstrained minimization. 605 606 Input Parameters: 607 . snes - the SNES context 608 . func - function evaluation routine 609 . ctx - optional user-defined function context 610 611 Calling sequence of func: 612 . func (SNES snes,Vec x,double *f,void *ctx); 613 614 . x - input vector 615 . f - function 616 . ctx - optional user-defined context for private data for the 617 function evaluation routine (may be null) 618 619 Notes: 620 SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 621 methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 622 SNESSetFunction(). 623 624 .keywords: SNES, nonlinear, set, minimization, function 625 626 .seealso: SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(), 627 SNESSetHessian(), SNESSetGradient() 628 @*/ 629 int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*), 630 void *ctx) 631 { 632 PetscValidHeaderSpecific(snes,SNES_COOKIE); 633 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1, 634 "SNESSetMinimizationFunction:Only for SNES_UNCONSTRAINED_MINIMIZATION"); 635 snes->computeumfunction = func; 636 snes->umfunP = ctx; 637 return 0; 638 } 639 640 /*@ 641 SNESComputeMinimizationFunction - Computes the function that has been 642 set with SNESSetMinimizationFunction(). 643 644 Input Parameters: 645 . snes - the SNES context 646 . x - input vector 647 648 Output Parameter: 649 . y - function value 650 651 Notes: 652 SNESComputeMinimizationFunction() is valid only for 653 SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 654 SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 655 656 .keywords: SNES, nonlinear, compute, minimization, function 657 658 .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(), 659 SNESComputeGradient(), SNESComputeHessian() 660 @*/ 661 int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y) 662 { 663 int ierr; 664 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1, 665 "SNESComputeMinimizationFunction:Only for SNES_UNCONSTRAINED_MINIMIZATION"); 666 PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0); 667 ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr); 668 PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0); 669 return 0; 670 } 671 672 /*@C 673 SNESSetGradient - Sets the gradient evaluation routine and gradient 674 vector for use by the SNES routines. 675 676 Input Parameters: 677 . snes - the SNES context 678 . func - function evaluation routine 679 . ctx - optional user-defined function context 680 . r - vector to store gradient value 681 682 Calling sequence of func: 683 . func (SNES, Vec x, Vec g, void *ctx); 684 685 . x - input vector 686 . g - gradient vector 687 . ctx - optional user-defined context for private data for the 688 function evaluation routine (may be null) 689 690 Notes: 691 SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 692 methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 693 SNESSetFunction(). 694 695 .keywords: SNES, nonlinear, set, function 696 697 .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(), 698 SNESSetMinimizationFunction(), 699 @*/ 700 int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx) 701 { 702 PetscValidHeaderSpecific(snes,SNES_COOKIE); 703 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1, 704 "SNESSetGradient:For SNES_UNCONSTRAINED_MINIMIZATION only"); 705 snes->computefunction = func; 706 snes->vec_func = snes->vec_func_always = r; 707 snes->funP = ctx; 708 return 0; 709 } 710 711 /*@ 712 SNESComputeGradient - Computes the gradient that has been set with 713 SNESSetGradient(). 714 715 Input Parameters: 716 . snes - the SNES context 717 . x - input vector 718 719 Output Parameter: 720 . y - gradient vector 721 722 Notes: 723 SNESComputeGradient() is valid only for 724 SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 725 SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 726 727 .keywords: SNES, nonlinear, compute, gradient 728 729 .seealso: SNESSetGradient(), SNESGetGradient(), 730 SNESComputeMinimizationFunction(), SNESComputeHessian() 731 @*/ 732 int SNESComputeGradient(SNES snes,Vec x, Vec y) 733 { 734 int ierr; 735 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 736 SETERRQ(1,"SNESComputeGradient:For SNES_UNCONSTRAINED_MINIMIZATION only"); 737 PLogEventBegin(SNES_GradientEval,snes,x,y,0); 738 ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 739 PLogEventEnd(SNES_GradientEval,snes,x,y,0); 740 return 0; 741 } 742 743 /*@ 744 SNESComputeJacobian - Computes the Jacobian matrix that has been 745 set with SNESSetJacobian(). 746 747 Input Parameters: 748 . snes - the SNES context 749 . x - input vector 750 751 Output Parameters: 752 . A - Jacobian matrix 753 . B - optional preconditioning matrix 754 . flag - flag indicating matrix structure 755 756 Notes: 757 Most users should not need to explicitly call this routine, as it 758 is used internally within the nonlinear solvers. 759 760 See SLESSetOperators() for information about setting the flag parameter. 761 762 SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS 763 methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION 764 methods is SNESComputeJacobian(). 765 766 .keywords: SNES, compute, Jacobian, matrix 767 768 .seealso: SNESSetJacobian(), SLESSetOperators() 769 @*/ 770 int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 771 { 772 int ierr; 773 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 774 SETERRQ(1,"SNESComputeJacobian: For SNES_NONLINEAR_EQUATIONS only"); 775 if (!snes->computejacobian) return 0; 776 PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B); 777 *flg = DIFFERENT_NONZERO_PATTERN; 778 ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr); 779 PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B); 780 /* make sure user returned a correct Jacobian and preconditioner */ 781 PetscValidHeaderSpecific(*A,MAT_COOKIE); 782 PetscValidHeaderSpecific(*B,MAT_COOKIE); 783 return 0; 784 } 785 786 /*@ 787 SNESComputeHessian - Computes the Hessian matrix that has been 788 set with SNESSetHessian(). 789 790 Input Parameters: 791 . snes - the SNES context 792 . x - input vector 793 794 Output Parameters: 795 . A - Hessian matrix 796 . B - optional preconditioning matrix 797 . flag - flag indicating matrix structure 798 799 Notes: 800 Most users should not need to explicitly call this routine, as it 801 is used internally within the nonlinear solvers. 802 803 See SLESSetOperators() for information about setting the flag parameter. 804 805 SNESComputeHessian() is valid only for 806 SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 807 SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian(). 808 809 .keywords: SNES, compute, Hessian, matrix 810 811 .seealso: SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(), 812 SNESComputeMinimizationFunction() 813 @*/ 814 int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag) 815 { 816 int ierr; 817 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 818 SETERRQ(1,"SNESComputeHessian:For SNES_UNCONSTRAINED_MINIMIZATION only"); 819 if (!snes->computejacobian) return 0; 820 PLogEventBegin(SNES_HessianEval,snes,x,*A,*B); 821 *flag = DIFFERENT_NONZERO_PATTERN; 822 ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr); 823 PLogEventEnd(SNES_HessianEval,snes,x,*A,*B); 824 /* make sure user returned a correct Jacobian and preconditioner */ 825 PetscValidHeaderSpecific(*A,MAT_COOKIE); 826 PetscValidHeaderSpecific(*B,MAT_COOKIE); 827 return 0; 828 } 829 830 /*@C 831 SNESSetJacobian - Sets the function to compute Jacobian as well as the 832 location to store it. 833 834 Input Parameters: 835 . snes - the SNES context 836 . A - Jacobian matrix 837 . B - preconditioner matrix (usually same as the Jacobian) 838 . func - Jacobian evaluation routine 839 . ctx - optional user-defined context for private data for the 840 Jacobian evaluation routine (may be null) 841 842 Calling sequence of func: 843 . func (SNES,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 844 845 . x - input vector 846 . A - Jacobian matrix 847 . B - preconditioner matrix, usually the same as A 848 . flag - flag indicating information about the preconditioner matrix 849 structure (same as flag in SLESSetOperators()) 850 . ctx - optional user-defined Jacobian context 851 852 Notes: 853 See SLESSetOperators() for information about setting the flag input 854 parameter in the routine func(). Be sure to read this information! 855 856 The routine func() takes Mat * as the matrix arguments rather than Mat. 857 This allows the Jacobian evaluation routine to replace A and/or B with a 858 completely new new matrix structure (not just different matrix elements) 859 when appropriate, for instance, if the nonzero structure is changing 860 throughout the global iterations. 861 862 .keywords: SNES, nonlinear, set, Jacobian, matrix 863 864 .seealso: SLESSetOperators(), SNESSetFunction() 865 @*/ 866 int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 867 MatStructure*,void*),void *ctx) 868 { 869 PetscValidHeaderSpecific(snes,SNES_COOKIE); 870 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 871 SETERRQ(1,"SNESSetJacobian:For SNES_NONLINEAR_EQUATIONS only"); 872 snes->computejacobian = func; 873 snes->jacP = ctx; 874 snes->jacobian = A; 875 snes->jacobian_pre = B; 876 return 0; 877 } 878 879 /*@ 880 SNESGetJacobian - Returns the Jacobian matrix and optionally the user 881 provided context for evaluating the Jacobian. 882 883 Input Parameter: 884 . snes - the nonlinear solver context 885 886 Output Parameters: 887 . A - location to stash Jacobian matrix (or PETSC_NULL) 888 . B - location to stash preconditioner matrix (or PETSC_NULL) 889 . ctx - location to stash Jacobian ctx (or PETSC_NULL) 890 891 .seealso: SNESSetJacobian(), SNESComputeJacobian() 892 @*/ 893 int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx) 894 { 895 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 896 SETERRQ(1,"SNESSetJacobian:For SNES_NONLINEAR_EQUATIONS only"); 897 if (A) *A = snes->jacobian; 898 if (B) *B = snes->jacobian_pre; 899 if (ctx) *ctx = snes->jacP; 900 return 0; 901 } 902 903 /*@C 904 SNESSetHessian - Sets the function to compute Hessian as well as the 905 location to store it. 906 907 Input Parameters: 908 . snes - the SNES context 909 . A - Hessian matrix 910 . B - preconditioner matrix (usually same as the Hessian) 911 . func - Jacobian evaluation routine 912 . ctx - optional user-defined context for private data for the 913 Hessian evaluation routine (may be null) 914 915 Calling sequence of func: 916 . func (SNES,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 917 918 . x - input vector 919 . A - Hessian matrix 920 . B - preconditioner matrix, usually the same as A 921 . flag - flag indicating information about the preconditioner matrix 922 structure (same as flag in SLESSetOperators()) 923 . ctx - optional user-defined Hessian context 924 925 Notes: 926 See SLESSetOperators() for information about setting the flag input 927 parameter in the routine func(). Be sure to read this information! 928 929 The function func() takes Mat * as the matrix arguments rather than Mat. 930 This allows the Hessian evaluation routine to replace A and/or B with a 931 completely new new matrix structure (not just different matrix elements) 932 when appropriate, for instance, if the nonzero structure is changing 933 throughout the global iterations. 934 935 .keywords: SNES, nonlinear, set, Hessian, matrix 936 937 .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators() 938 @*/ 939 int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 940 MatStructure*,void*),void *ctx) 941 { 942 PetscValidHeaderSpecific(snes,SNES_COOKIE); 943 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 944 SETERRQ(1,"SNESSetHessian:For SNES_UNCONSTRAINED_MINIMIZATION only"); 945 snes->computejacobian = func; 946 snes->jacP = ctx; 947 snes->jacobian = A; 948 snes->jacobian_pre = B; 949 return 0; 950 } 951 952 /*@ 953 SNESGetHessian - Returns the Hessian matrix and optionally the user 954 provided context for evaluating the Hessian. 955 956 Input Parameter: 957 . snes - the nonlinear solver context 958 959 Output Parameters: 960 . A - location to stash Hessian matrix (or PETSC_NULL) 961 . B - location to stash preconditioner matrix (or PETSC_NULL) 962 . ctx - location to stash Hessian ctx (or PETSC_NULL) 963 964 .seealso: SNESSetHessian(), SNESComputeHessian() 965 @*/ 966 int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx) 967 { 968 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 969 SETERRQ(1,"SNESSetHessian:For SNES_UNCONSTRAINED_MINIMIZATION only"); 970 if (A) *A = snes->jacobian; 971 if (B) *B = snes->jacobian_pre; 972 if (ctx) *ctx = snes->jacP; 973 return 0; 974 } 975 976 /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 977 978 /*@ 979 SNESSetUp - Sets up the internal data structures for the later use 980 of a nonlinear solver. 981 982 Input Parameter: 983 . snes - the SNES context 984 . x - the solution vector 985 986 Notes: 987 For basic use of the SNES solvers the user need not explicitly call 988 SNESSetUp(), since these actions will automatically occur during 989 the call to SNESSolve(). However, if one wishes to control this 990 phase separately, SNESSetUp() should be called after SNESCreate() 991 and optional routines of the form SNESSetXXX(), but before SNESSolve(). 992 993 .keywords: SNES, nonlinear, setup 994 995 .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 996 @*/ 997 int SNESSetUp(SNES snes,Vec x) 998 { 999 int ierr, flg; 1000 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1001 PetscValidHeaderSpecific(x,VEC_COOKIE); 1002 snes->vec_sol = snes->vec_sol_always = x; 1003 1004 ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg); CHKERRQ(ierr); 1005 /* 1006 This version replaces the user provided Jacobian matrix with a 1007 matrix-free version but still employs the user-provided preconditioner matrix 1008 */ 1009 if (flg) { 1010 Mat J; 1011 ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1012 PLogObjectParent(snes,J); 1013 snes->mfshell = J; 1014 if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 1015 snes->jacobian = J; 1016 PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n"); 1017 } 1018 else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 1019 snes->jacobian = J; 1020 PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n"); 1021 } else SETERRQ(1,"SNESSetUp:Method class doesn't support matrix-free operator option"); 1022 } 1023 ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg); CHKERRQ(ierr); 1024 /* 1025 This version replaces both the user-provided Jacobian and the user- 1026 provided preconditioner matrix with the default matrix free version. 1027 */ 1028 if (flg) { 1029 Mat J; 1030 ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1031 PLogObjectParent(snes,J); 1032 snes->mfshell = J; 1033 if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 1034 ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 1035 PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n"); 1036 } 1037 else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 1038 ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 1039 PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n"); 1040 } else SETERRQ(1,"SNESSetUp:Method class doesn't support matrix-free option"); 1041 } 1042 if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) { 1043 if (!snes->vec_func) SETERRQ(1,"SNESSetUp:Must call SNESSetFunction() first"); 1044 if (!snes->computefunction) SETERRQ(1,"SNESSetUp:Must call SNESSetFunction() first"); 1045 if (!snes->jacobian) SETERRQ(1,"SNESSetUp:Must call SNESSetJacobian() first"); 1046 if (snes->vec_func == snes->vec_sol) SETERRQ(1,"SNESSetUp:Solution vector cannot be function vector"); 1047 1048 /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 1049 if (snes->ksp_ewconv && snes->type != SNES_EQ_TR) { 1050 SLES sles; KSP ksp; 1051 ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 1052 ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr); 1053 ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private, 1054 (void *)snes); CHKERRQ(ierr); 1055 } 1056 } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) { 1057 if (!snes->vec_func) SETERRQ(1,"SNESSetUp:Must call SNESSetGradient() first"); 1058 if (!snes->computefunction) SETERRQ(1,"SNESSetUp:Must call SNESSetGradient() first"); 1059 if (!snes->computeumfunction) 1060 SETERRQ(1,"SNESSetUp:Must call SNESSetMinimizationFunction() first"); 1061 if (!snes->jacobian) SETERRQ(1,"SNESSetUp:Must call SNESSetHessian() first"); 1062 } else SETERRQ(1,"SNESSetUp:Unknown method class"); 1063 if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);} 1064 snes->setup_called = 1; 1065 return 0; 1066 } 1067 1068 /*@C 1069 SNESDestroy - Destroys the nonlinear solver context that was created 1070 with SNESCreate(). 1071 1072 Input Parameter: 1073 . snes - the SNES context 1074 1075 .keywords: SNES, nonlinear, destroy 1076 1077 .seealso: SNESCreate(), SNESSolve() 1078 @*/ 1079 int SNESDestroy(SNES snes) 1080 { 1081 int ierr; 1082 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1083 ierr = (*(snes)->destroy)((PetscObject)snes); CHKERRQ(ierr); 1084 if (snes->kspconvctx) PetscFree(snes->kspconvctx); 1085 if (snes->mfshell) MatDestroy(snes->mfshell); 1086 ierr = SLESDestroy(snes->sles); CHKERRQ(ierr); 1087 if (snes->xmonitor) SNESLGMonitorDestroy(snes->xmonitor); 1088 if (snes->vwork) VecDestroyVecs(snes->vwork,snes->nvwork); 1089 PLogObjectDestroy((PetscObject)snes); 1090 PetscHeaderDestroy((PetscObject)snes); 1091 return 0; 1092 } 1093 1094 /* ----------- Routines to set solver parameters ---------- */ 1095 1096 /*@ 1097 SNESSetTolerances - Sets various parameters used in convergence tests. 1098 1099 Input Parameters: 1100 . snes - the SNES context 1101 . atol - absolute convergence tolerance 1102 . rtol - relative convergence tolerance 1103 . stol - convergence tolerance in terms of the norm 1104 of the change in the solution between steps 1105 . maxit - maximum number of iterations 1106 . maxf - maximum number of function evaluations 1107 1108 Options Database Keys: 1109 $ -snes_atol <atol> 1110 $ -snes_rtol <rtol> 1111 $ -snes_stol <stol> 1112 $ -snes_max_it <maxit> 1113 $ -snes_max_funcs <maxf> 1114 1115 Notes: 1116 The default maximum number of iterations is 50. 1117 The default maximum number of function evaluations is 1000. 1118 1119 .keywords: SNES, nonlinear, set, convergence, tolerances 1120 1121 .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance() 1122 @*/ 1123 int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf) 1124 { 1125 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1126 if (atol != PETSC_DEFAULT) snes->atol = atol; 1127 if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1128 if (stol != PETSC_DEFAULT) snes->xtol = stol; 1129 if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1130 if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 1131 return 0; 1132 } 1133 1134 /*@ 1135 SNESGetTolerances - Gets various parameters used in convergence tests. 1136 1137 Input Parameters: 1138 . snes - the SNES context 1139 . atol - absolute convergence tolerance 1140 . rtol - relative convergence tolerance 1141 . stol - convergence tolerance in terms of the norm 1142 of the change in the solution between steps 1143 . maxit - maximum number of iterations 1144 . maxf - maximum number of function evaluations 1145 1146 Notes: 1147 The user can specify PETSC_NULL for any parameter that is not needed. 1148 1149 .keywords: SNES, nonlinear, get, convergence, tolerances 1150 1151 .seealso: SNESSetTolerances() 1152 @*/ 1153 int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf) 1154 { 1155 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1156 if (atol) *atol = snes->atol; 1157 if (rtol) *rtol = snes->rtol; 1158 if (stol) *stol = snes->xtol; 1159 if (maxit) *maxit = snes->max_its; 1160 if (maxf) *maxf = snes->max_funcs; 1161 return 0; 1162 } 1163 1164 /*@ 1165 SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 1166 1167 Input Parameters: 1168 . snes - the SNES context 1169 . tol - tolerance 1170 1171 Options Database Key: 1172 $ -snes_trtol <tol> 1173 1174 .keywords: SNES, nonlinear, set, trust region, tolerance 1175 1176 .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance() 1177 @*/ 1178 int SNESSetTrustRegionTolerance(SNES snes,double tol) 1179 { 1180 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1181 snes->deltatol = tol; 1182 return 0; 1183 } 1184 1185 /*@ 1186 SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance 1187 for unconstrained minimization solvers. 1188 1189 Input Parameters: 1190 . snes - the SNES context 1191 . ftol - minimum function tolerance 1192 1193 Options Database Key: 1194 $ -snes_fmin <ftol> 1195 1196 Note: 1197 SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION 1198 methods only. 1199 1200 .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance 1201 1202 .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance() 1203 @*/ 1204 int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol) 1205 { 1206 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1207 snes->fmin = ftol; 1208 return 0; 1209 } 1210 1211 /* ------------ Routines to set performance monitoring options ----------- */ 1212 1213 /*@C 1214 SNESSetMonitor - Sets the function that is to be used at every 1215 iteration of the nonlinear solver to display the iteration's 1216 progress. 1217 1218 Input Parameters: 1219 . snes - the SNES context 1220 . func - monitoring routine 1221 . mctx - optional user-defined context for private data for the 1222 monitor routine (may be null) 1223 1224 Calling sequence of func: 1225 int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx) 1226 1227 $ snes - the SNES context 1228 $ its - iteration number 1229 $ mctx - optional monitoring context 1230 $ 1231 $ SNES_NONLINEAR_EQUATIONS methods: 1232 $ norm - 2-norm function value (may be estimated) 1233 $ 1234 $ SNES_UNCONSTRAINED_MINIMIZATION methods: 1235 $ norm - 2-norm gradient value (may be estimated) 1236 1237 .keywords: SNES, nonlinear, set, monitor 1238 1239 .seealso: SNESDefaultMonitor() 1240 @*/ 1241 int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx ) 1242 { 1243 snes->monitor = func; 1244 snes->monP = mctx; 1245 return 0; 1246 } 1247 1248 /*@C 1249 SNESSetConvergenceTest - Sets the function that is to be used 1250 to test for convergence of the nonlinear iterative solution. 1251 1252 Input Parameters: 1253 . snes - the SNES context 1254 . func - routine to test for convergence 1255 . cctx - optional context for private data for the convergence routine 1256 (may be null) 1257 1258 Calling sequence of func: 1259 int func (SNES snes,double xnorm,double gnorm, 1260 double f,void *cctx) 1261 1262 $ snes - the SNES context 1263 $ cctx - optional convergence context 1264 $ xnorm - 2-norm of current iterate 1265 $ 1266 $ SNES_NONLINEAR_EQUATIONS methods: 1267 $ gnorm - 2-norm of current step 1268 $ f - 2-norm of function 1269 $ 1270 $ SNES_UNCONSTRAINED_MINIMIZATION methods: 1271 $ gnorm - 2-norm of current gradient 1272 $ f - function value 1273 1274 .keywords: SNES, nonlinear, set, convergence, test 1275 1276 .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(), 1277 SNESConverged_UM_LS(), SNESConverged_UM_TR() 1278 @*/ 1279 int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx) 1280 { 1281 (snes)->converged = func; 1282 (snes)->cnvP = cctx; 1283 return 0; 1284 } 1285 1286 /* 1287 SNESScaleStep_Private - Scales a step so that its length is less than the 1288 positive parameter delta. 1289 1290 Input Parameters: 1291 . snes - the SNES context 1292 . y - approximate solution of linear system 1293 . fnorm - 2-norm of current function 1294 . delta - trust region size 1295 1296 Output Parameters: 1297 . gpnorm - predicted function norm at the new point, assuming local 1298 linearization. The value is zero if the step lies within the trust 1299 region, and exceeds zero otherwise. 1300 . ynorm - 2-norm of the step 1301 1302 Note: 1303 For non-trust region methods such as SNES_EQ_LS, the parameter delta 1304 is set to be the maximum allowable step size. 1305 1306 .keywords: SNES, nonlinear, scale, step 1307 */ 1308 int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta, 1309 double *gpnorm,double *ynorm) 1310 { 1311 double norm; 1312 Scalar cnorm; 1313 VecNorm(y,NORM_2, &norm ); 1314 if (norm > *delta) { 1315 norm = *delta/norm; 1316 *gpnorm = (1.0 - norm)*(*fnorm); 1317 cnorm = norm; 1318 VecScale( &cnorm, y ); 1319 *ynorm = *delta; 1320 } else { 1321 *gpnorm = 0.0; 1322 *ynorm = norm; 1323 } 1324 return 0; 1325 } 1326 1327 /*@ 1328 SNESSolve - Solves a nonlinear system. Call SNESSolve after calling 1329 SNESCreate() and optional routines of the form SNESSetXXX(). 1330 1331 Input Parameter: 1332 . snes - the SNES context 1333 . x - the solution vector 1334 1335 Output Parameter: 1336 its - number of iterations until termination 1337 1338 Note: 1339 The user should initialize the vector, x, with the initial guess 1340 for the nonlinear solve prior to calling SNESSolve. In particular, 1341 to employ an initial guess of zero, the user should explicitly set 1342 this vector to zero by calling VecSet(). 1343 1344 .keywords: SNES, nonlinear, solve 1345 1346 .seealso: SNESCreate(), SNESDestroy() 1347 @*/ 1348 int SNESSolve(SNES snes,Vec x,int *its) 1349 { 1350 int ierr, flg; 1351 1352 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1353 PetscValidIntPointer(its); 1354 if (!snes->setup_called) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);} 1355 else {snes->vec_sol = snes->vec_sol_always = x;} 1356 PLogEventBegin(SNES_Solve,snes,0,0,0); 1357 ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr); 1358 PLogEventEnd(SNES_Solve,snes,0,0,0); 1359 ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr); 1360 if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); } 1361 return 0; 1362 } 1363 1364 /* --------- Internal routines for SNES Package --------- */ 1365 static NRList *__SNESList = 0; 1366 1367 /*@ 1368 SNESSetType - Sets the method for the nonlinear solver. 1369 1370 Input Parameters: 1371 . snes - the SNES context 1372 . method - a known method 1373 1374 Notes: 1375 See "petsc/include/snes.h" for available methods (for instance) 1376 $ Systems of nonlinear equations: 1377 $ SNES_EQ_LS - Newton's method with line search 1378 $ SNES_EQ_TR - Newton's method with trust region 1379 $ Unconstrained minimization: 1380 $ SNES_UM_TR - Newton's method with trust region 1381 $ SNES_UM_LS - Newton's method with line search 1382 1383 Options Database Command: 1384 $ -snes_type <method> 1385 $ Use -help for a list of available methods 1386 $ (for instance, ls or tr) 1387 1388 .keysords: SNES, set, method 1389 @*/ 1390 int SNESSetType(SNES snes,SNESType method) 1391 { 1392 int (*r)(SNES); 1393 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1394 if (snes->type == (int) method) return 0; 1395 1396 /* Get the function pointers for the iterative method requested */ 1397 if (!__SNESList) {SNESRegisterAll();} 1398 if (!__SNESList) {SETERRQ(1,"SNESSetType:Could not get methods");} 1399 r = (int (*)(SNES))NRFindRoutine( __SNESList, (int)method, (char *)0 ); 1400 if (!r) {SETERRQ(1,"SNESSetType:Unknown method");} 1401 if (snes->data) PetscFree(snes->data); 1402 snes->set_method_called = 1; 1403 return (*r)(snes); 1404 } 1405 1406 /* --------------------------------------------------------------------- */ 1407 /*@C 1408 SNESRegister - Adds the method to the nonlinear solver package, given 1409 a function pointer and a nonlinear solver name of the type SNESType. 1410 1411 Input Parameters: 1412 . name - for instance SNES_EQ_LS, SNES_EQ_TR, ... 1413 . sname - corresponding string for name 1414 . create - routine to create method context 1415 1416 .keywords: SNES, nonlinear, register 1417 1418 .seealso: SNESRegisterAll(), SNESRegisterDestroy() 1419 @*/ 1420 int SNESRegister(int name, char *sname, int (*create)(SNES)) 1421 { 1422 int ierr; 1423 if (!__SNESList) {ierr = NRCreate(&__SNESList); CHKERRQ(ierr);} 1424 NRRegister( __SNESList, name, sname, (int (*)(void*))create ); 1425 return 0; 1426 } 1427 /* --------------------------------------------------------------------- */ 1428 /*@C 1429 SNESRegisterDestroy - Frees the list of nonlinear solvers that were 1430 registered by SNESRegister(). 1431 1432 .keywords: SNES, nonlinear, register, destroy 1433 1434 .seealso: SNESRegisterAll(), SNESRegisterAll() 1435 @*/ 1436 int SNESRegisterDestroy() 1437 { 1438 if (__SNESList) { 1439 NRDestroy( __SNESList ); 1440 __SNESList = 0; 1441 } 1442 return 0; 1443 } 1444 1445 /* 1446 SNESGetTypeFromOptions_Private - Sets the selected method from the 1447 options database. 1448 1449 Input Parameter: 1450 . ctx - the SNES context 1451 1452 Output Parameter: 1453 . method - solver method 1454 1455 Returns: 1456 Returns 1 if the method is found; 0 otherwise. 1457 1458 Options Database Key: 1459 $ -snes_type method 1460 */ 1461 int SNESGetTypeFromOptions_Private(SNES ctx,SNESType *method,int *flg) 1462 { 1463 int ierr; 1464 char sbuf[50]; 1465 1466 ierr = OptionsGetString(ctx->prefix,"-snes_type",sbuf,50,flg);CHKERRQ(ierr); 1467 if (*flg) { 1468 if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);} 1469 *method = (SNESType)NRFindID( __SNESList, sbuf ); 1470 } 1471 return 0; 1472 } 1473 1474 /*@C 1475 SNESGetType - Gets the SNES method type and name (as a string). 1476 1477 Input Parameter: 1478 . snes - nonlinear solver context 1479 1480 Output Parameter: 1481 . method - SNES method (or use PETSC_NULL) 1482 . name - name of SNES method (or use PETSC_NULL) 1483 1484 .keywords: SNES, nonlinear, get, method, name 1485 @*/ 1486 int SNESGetType(SNES snes, SNESType *method,char **name) 1487 { 1488 int ierr; 1489 if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);} 1490 if (method) *method = (SNESType) snes->type; 1491 if (name) *name = NRFindName( __SNESList, (int) snes->type ); 1492 return 0; 1493 } 1494 1495 #include <stdio.h> 1496 /* 1497 SNESPrintTypes_Private - Prints the SNES methods available from the 1498 options database. 1499 1500 Input Parameters: 1501 . comm - communicator (usually MPI_COMM_WORLD) 1502 . prefix - prefix (usually "-") 1503 . name - the options database name (by default "snes_type") 1504 */ 1505 int SNESPrintTypes_Private(MPI_Comm comm,char* prefix,char *name) 1506 { 1507 FuncList *entry; 1508 if (!__SNESList) {SNESRegisterAll();} 1509 entry = __SNESList->head; 1510 PetscPrintf(comm," %s%s (one of)",prefix,name); 1511 while (entry) { 1512 PetscPrintf(comm," %s",entry->name); 1513 entry = entry->next; 1514 } 1515 PetscPrintf(comm,"\n"); 1516 return 0; 1517 } 1518 1519 /*@C 1520 SNESGetSolution - Returns the vector where the approximate solution is 1521 stored. 1522 1523 Input Parameter: 1524 . snes - the SNES context 1525 1526 Output Parameter: 1527 . x - the solution 1528 1529 .keywords: SNES, nonlinear, get, solution 1530 1531 .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate() 1532 @*/ 1533 int SNESGetSolution(SNES snes,Vec *x) 1534 { 1535 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1536 *x = snes->vec_sol_always; 1537 return 0; 1538 } 1539 1540 /*@C 1541 SNESGetSolutionUpdate - Returns the vector where the solution update is 1542 stored. 1543 1544 Input Parameter: 1545 . snes - the SNES context 1546 1547 Output Parameter: 1548 . x - the solution update 1549 1550 Notes: 1551 This vector is implementation dependent. 1552 1553 .keywords: SNES, nonlinear, get, solution, update 1554 1555 .seealso: SNESGetSolution(), SNESGetFunction 1556 @*/ 1557 int SNESGetSolutionUpdate(SNES snes,Vec *x) 1558 { 1559 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1560 *x = snes->vec_sol_update_always; 1561 return 0; 1562 } 1563 1564 /*@C 1565 SNESGetFunction - Returns the vector where the function is stored. 1566 1567 Input Parameter: 1568 . snes - the SNES context 1569 1570 Output Parameter: 1571 . r - the function 1572 1573 Notes: 1574 SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only 1575 Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 1576 SNESGetMinimizationFunction() and SNESGetGradient(); 1577 1578 .keywords: SNES, nonlinear, get, function 1579 1580 .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(), 1581 SNESGetGradient() 1582 @*/ 1583 int SNESGetFunction(SNES snes,Vec *r) 1584 { 1585 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1586 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1, 1587 "SNESGetFunction:For SNES_NONLINEAR_EQUATIONS only"); 1588 *r = snes->vec_func_always; 1589 return 0; 1590 } 1591 1592 /*@C 1593 SNESGetGradient - Returns the vector where the gradient is stored. 1594 1595 Input Parameter: 1596 . snes - the SNES context 1597 1598 Output Parameter: 1599 . r - the gradient 1600 1601 Notes: 1602 SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods 1603 only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 1604 SNESGetFunction(). 1605 1606 .keywords: SNES, nonlinear, get, gradient 1607 1608 .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction() 1609 @*/ 1610 int SNESGetGradient(SNES snes,Vec *r) 1611 { 1612 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1613 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1, 1614 "SNESGetGradient:For SNES_UNCONSTRAINED_MINIMIZATION only"); 1615 *r = snes->vec_func_always; 1616 return 0; 1617 } 1618 1619 /*@ 1620 SNESGetMinimizationFunction - Returns the scalar function value for 1621 unconstrained minimization problems. 1622 1623 Input Parameter: 1624 . snes - the SNES context 1625 1626 Output Parameter: 1627 . r - the function 1628 1629 Notes: 1630 SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 1631 methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 1632 SNESGetFunction(). 1633 1634 .keywords: SNES, nonlinear, get, function 1635 1636 .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction() 1637 @*/ 1638 int SNESGetMinimizationFunction(SNES snes,double *r) 1639 { 1640 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1641 PetscValidScalarPointer(r); 1642 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1, 1643 "SNESGetMinimizationFunction:For SNES_UNCONSTRAINED_MINIMIZATION only"); 1644 *r = snes->fc; 1645 return 0; 1646 } 1647 1648 /*@C 1649 SNESSetOptionsPrefix - Sets the prefix used for searching for all 1650 SNES options in the database. 1651 1652 Input Parameter: 1653 . snes - the SNES context 1654 . prefix - the prefix to prepend to all option names 1655 1656 .keywords: SNES, set, options, prefix, database 1657 1658 .seealso: SNESSetFromOptions() 1659 @*/ 1660 int SNESSetOptionsPrefix(SNES snes,char *prefix) 1661 { 1662 int ierr; 1663 1664 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1665 ierr = PetscObjectSetPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 1666 ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 1667 return 0; 1668 } 1669 1670 /*@C 1671 SNESAppendOptionsPrefix - Append to the prefix used for searching for all 1672 SNES options in the database. 1673 1674 Input Parameter: 1675 . snes - the SNES context 1676 . prefix - the prefix to prepend to all option names 1677 1678 .keywords: SNES, append, options, prefix, database 1679 1680 .seealso: SNESGetOptionsPrefix() 1681 @*/ 1682 int SNESAppendOptionsPrefix(SNES snes,char *prefix) 1683 { 1684 int ierr; 1685 1686 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1687 ierr = PetscObjectAppendPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 1688 ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 1689 return 0; 1690 } 1691 1692 /*@ 1693 SNESGetOptionsPrefix - Sets the prefix used for searching for all 1694 SNES options in the database. 1695 1696 Input Parameter: 1697 . snes - the SNES context 1698 1699 Output Parameter: 1700 . prefix - pointer to the prefix string used 1701 1702 .keywords: SNES, get, options, prefix, database 1703 1704 .seealso: SNESAppendOptionsPrefix() 1705 @*/ 1706 int SNESGetOptionsPrefix(SNES snes,char **prefix) 1707 { 1708 int ierr; 1709 1710 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1711 ierr = PetscObjectGetPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 1712 return 0; 1713 } 1714 1715 1716 1717 1718 1719