1 #ifndef lint 2 static char vcid[] = "$Id: snes.c,v 1.84 1996/08/26 02:35:13 bsmith 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 for Jacobian and preconditioner 449 $ -snes_mf_operator - use default matrix free for Jacobian but not for preconditioner 450 $ -snes_fd - use slow finite differences to compute Jacobian 451 452 .keywords: SNES, nonlinear, create, context 453 454 .seealso: SNESSolve(), SNESDestroy() 455 @*/ 456 int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes) 457 { 458 int ierr; 459 SNES snes; 460 SNES_KSP_EW_ConvCtx *kctx; 461 462 *outsnes = 0; 463 if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS) 464 SETERRQ(1,"SNESCreate:incorrect method type"); 465 PetscHeaderCreate(snes,_SNES,SNES_COOKIE,SNES_UNKNOWN_METHOD,comm); 466 PLogObjectCreate(snes); 467 snes->max_its = 50; 468 snes->max_funcs = 1000; 469 snes->norm = 0.0; 470 if (type == SNES_UNCONSTRAINED_MINIMIZATION) { 471 snes->rtol = 1.e-8; 472 snes->ttol = 0.0; 473 snes->atol = 1.e-10; 474 } 475 else { 476 snes->rtol = 1.e-8; 477 snes->ttol = 0.0; 478 snes->atol = 1.e-50; 479 } 480 snes->xtol = 1.e-8; 481 snes->trunctol = 1.e-12; /* no longer used */ 482 snes->nfuncs = 0; 483 snes->nfailures = 0; 484 snes->monitor = 0; 485 snes->data = 0; 486 snes->view = 0; 487 snes->computeumfunction = 0; 488 snes->umfunP = 0; 489 snes->fc = 0; 490 snes->deltatol = 1.e-12; 491 snes->fmin = -1.e30; 492 snes->method_class = type; 493 snes->set_method_called = 0; 494 snes->setup_called = 0; 495 snes->ksp_ewconv = 0; 496 snes->type = -1; 497 snes->xmonitor = 0; 498 snes->vwork = 0; 499 snes->nwork = 0; 500 501 /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 502 kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx); 503 snes->kspconvctx = (void*)kctx; 504 kctx->version = 2; 505 kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 506 this was too large for some test cases */ 507 kctx->rtol_last = 0; 508 kctx->rtol_max = .9; 509 kctx->gamma = 1.0; 510 kctx->alpha2 = .5*(1.0 + sqrt(5.0)); 511 kctx->alpha = kctx->alpha2; 512 kctx->threshold = .1; 513 kctx->lresid_last = 0; 514 kctx->norm_last = 0; 515 516 ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr); 517 PLogObjectParent(snes,snes->sles) 518 519 *outsnes = snes; 520 return 0; 521 } 522 523 /* --------------------------------------------------------------- */ 524 /*@C 525 SNESSetFunction - Sets the function evaluation routine and function 526 vector for use by the SNES routines in solving systems of nonlinear 527 equations. 528 529 Input Parameters: 530 . snes - the SNES context 531 . func - function evaluation routine 532 . ctx - optional user-defined function context 533 . r - vector to store function value 534 535 Calling sequence of func: 536 . func (SNES, Vec x, Vec f, void *ctx); 537 538 . x - input vector 539 . f - vector function 540 . ctx - optional user-defined context for private data for the 541 function evaluation routine (may be null) 542 543 Notes: 544 The Newton-like methods typically solve linear systems of the form 545 $ f'(x) x = -f(x), 546 $ where f'(x) denotes the Jacobian matrix and f(x) is the function. 547 548 SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 549 Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 550 SNESSetMinimizationFunction() and SNESSetGradient(); 551 552 .keywords: SNES, nonlinear, set, function 553 554 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 555 @*/ 556 int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx) 557 { 558 PetscValidHeaderSpecific(snes,SNES_COOKIE); 559 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1, 560 "SNESSetFunction:For SNES_NONLINEAR_EQUATIONS only"); 561 snes->computefunction = func; 562 snes->vec_func = snes->vec_func_always = r; 563 snes->funP = ctx; 564 return 0; 565 } 566 567 /*@ 568 SNESComputeFunction - Computes the function that has been set with 569 SNESSetFunction(). 570 571 Input Parameters: 572 . snes - the SNES context 573 . x - input vector 574 575 Output Parameter: 576 . y - function vector, as set by SNESSetFunction() 577 578 Notes: 579 SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 580 Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 581 SNESComputeMinimizationFunction() and SNESComputeGradient(); 582 583 .keywords: SNES, nonlinear, compute, function 584 585 .seealso: SNESSetFunction(), SNESGetFunction() 586 @*/ 587 int SNESComputeFunction(SNES snes,Vec x, Vec y) 588 { 589 int ierr; 590 591 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 592 SETERRQ(1,"SNESComputeFunction: For SNES_NONLINEAR_EQUATIONS only"); 593 PLogEventBegin(SNES_FunctionEval,snes,x,y,0); 594 ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 595 PLogEventEnd(SNES_FunctionEval,snes,x,y,0); 596 return 0; 597 } 598 599 /*@C 600 SNESSetMinimizationFunction - Sets the function evaluation routine for 601 unconstrained minimization. 602 603 Input Parameters: 604 . snes - the SNES context 605 . func - function evaluation routine 606 . ctx - optional user-defined function context 607 608 Calling sequence of func: 609 . func (SNES snes,Vec x,double *f,void *ctx); 610 611 . x - input vector 612 . f - function 613 . ctx - optional user-defined context for private data for the 614 function evaluation routine (may be null) 615 616 Notes: 617 SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 618 methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 619 SNESSetFunction(). 620 621 .keywords: SNES, nonlinear, set, minimization, function 622 623 .seealso: SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(), 624 SNESSetHessian(), SNESSetGradient() 625 @*/ 626 int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*), 627 void *ctx) 628 { 629 PetscValidHeaderSpecific(snes,SNES_COOKIE); 630 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1, 631 "SNESSetMinimizationFunction:Only for SNES_UNCONSTRAINED_MINIMIZATION"); 632 snes->computeumfunction = func; 633 snes->umfunP = ctx; 634 return 0; 635 } 636 637 /*@ 638 SNESComputeMinimizationFunction - Computes the function that has been 639 set with SNESSetMinimizationFunction(). 640 641 Input Parameters: 642 . snes - the SNES context 643 . x - input vector 644 645 Output Parameter: 646 . y - function value 647 648 Notes: 649 SNESComputeMinimizationFunction() is valid only for 650 SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 651 SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 652 653 .keywords: SNES, nonlinear, compute, minimization, function 654 655 .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(), 656 SNESComputeGradient(), SNESComputeHessian() 657 @*/ 658 int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y) 659 { 660 int ierr; 661 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1, 662 "SNESComputeMinimizationFunction:Only for SNES_UNCONSTRAINED_MINIMIZATION"); 663 PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0); 664 ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr); 665 PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0); 666 return 0; 667 } 668 669 /*@C 670 SNESSetGradient - Sets the gradient evaluation routine and gradient 671 vector for use by the SNES routines. 672 673 Input Parameters: 674 . snes - the SNES context 675 . func - function evaluation routine 676 . ctx - optional user-defined function context 677 . r - vector to store gradient value 678 679 Calling sequence of func: 680 . func (SNES, Vec x, Vec g, void *ctx); 681 682 . x - input vector 683 . g - gradient vector 684 . ctx - optional user-defined context for private data for the 685 function evaluation routine (may be null) 686 687 Notes: 688 SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 689 methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 690 SNESSetFunction(). 691 692 .keywords: SNES, nonlinear, set, function 693 694 .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(), 695 SNESSetMinimizationFunction(), 696 @*/ 697 int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx) 698 { 699 PetscValidHeaderSpecific(snes,SNES_COOKIE); 700 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1, 701 "SNESSetGradient:For SNES_UNCONSTRAINED_MINIMIZATION only"); 702 snes->computefunction = func; 703 snes->vec_func = snes->vec_func_always = r; 704 snes->funP = ctx; 705 return 0; 706 } 707 708 /*@ 709 SNESComputeGradient - Computes the gradient that has been set with 710 SNESSetGradient(). 711 712 Input Parameters: 713 . snes - the SNES context 714 . x - input vector 715 716 Output Parameter: 717 . y - gradient vector 718 719 Notes: 720 SNESComputeGradient() is valid only for 721 SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 722 SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 723 724 .keywords: SNES, nonlinear, compute, gradient 725 726 .seealso: SNESSetGradient(), SNESGetGradient(), 727 SNESComputeMinimizationFunction(), SNESComputeHessian() 728 @*/ 729 int SNESComputeGradient(SNES snes,Vec x, Vec y) 730 { 731 int ierr; 732 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 733 SETERRQ(1,"SNESComputeGradient:For SNES_UNCONSTRAINED_MINIMIZATION only"); 734 PLogEventBegin(SNES_GradientEval,snes,x,y,0); 735 ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 736 PLogEventEnd(SNES_GradientEval,snes,x,y,0); 737 return 0; 738 } 739 740 /*@ 741 SNESComputeJacobian - Computes the Jacobian matrix that has been 742 set with SNESSetJacobian(). 743 744 Input Parameters: 745 . snes - the SNES context 746 . x - input vector 747 748 Output Parameters: 749 . A - Jacobian matrix 750 . B - optional preconditioning matrix 751 . flag - flag indicating matrix structure 752 753 Notes: 754 Most users should not need to explicitly call this routine, as it 755 is used internally within the nonlinear solvers. 756 757 See SLESSetOperators() for information about setting the flag parameter. 758 759 SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS 760 methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION 761 methods is SNESComputeJacobian(). 762 763 .keywords: SNES, compute, Jacobian, matrix 764 765 .seealso: SNESSetJacobian(), SLESSetOperators() 766 @*/ 767 int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 768 { 769 int ierr; 770 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 771 SETERRQ(1,"SNESComputeJacobian: For SNES_NONLINEAR_EQUATIONS only"); 772 if (!snes->computejacobian) return 0; 773 PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B); 774 *flg = DIFFERENT_NONZERO_PATTERN; 775 ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr); 776 PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B); 777 /* make sure user returned a correct Jacobian and preconditioner */ 778 PetscValidHeaderSpecific(*A,MAT_COOKIE); 779 PetscValidHeaderSpecific(*B,MAT_COOKIE); 780 return 0; 781 } 782 783 /*@ 784 SNESComputeHessian - Computes the Hessian matrix that has been 785 set with SNESSetHessian(). 786 787 Input Parameters: 788 . snes - the SNES context 789 . x - input vector 790 791 Output Parameters: 792 . A - Hessian matrix 793 . B - optional preconditioning matrix 794 . flag - flag indicating matrix structure 795 796 Notes: 797 Most users should not need to explicitly call this routine, as it 798 is used internally within the nonlinear solvers. 799 800 See SLESSetOperators() for information about setting the flag parameter. 801 802 SNESComputeHessian() is valid only for 803 SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 804 SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian(). 805 806 .keywords: SNES, compute, Hessian, matrix 807 808 .seealso: SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(), 809 SNESComputeMinimizationFunction() 810 @*/ 811 int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag) 812 { 813 int ierr; 814 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 815 SETERRQ(1,"SNESComputeHessian:For SNES_UNCONSTRAINED_MINIMIZATION only"); 816 if (!snes->computejacobian) return 0; 817 PLogEventBegin(SNES_HessianEval,snes,x,*A,*B); 818 *flag = DIFFERENT_NONZERO_PATTERN; 819 ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr); 820 PLogEventEnd(SNES_HessianEval,snes,x,*A,*B); 821 /* make sure user returned a correct Jacobian and preconditioner */ 822 PetscValidHeaderSpecific(*A,MAT_COOKIE); 823 PetscValidHeaderSpecific(*B,MAT_COOKIE); 824 return 0; 825 } 826 827 /*@C 828 SNESSetJacobian - Sets the function to compute Jacobian as well as the 829 location to store it. 830 831 Input Parameters: 832 . snes - the SNES context 833 . A - Jacobian matrix 834 . B - preconditioner matrix (usually same as the Jacobian) 835 . func - Jacobian evaluation routine 836 . ctx - optional user-defined context for private data for the 837 Jacobian evaluation routine (may be null) 838 839 Calling sequence of func: 840 . func (SNES,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 841 842 . x - input vector 843 . A - Jacobian matrix 844 . B - preconditioner matrix, usually the same as A 845 . flag - flag indicating information about the preconditioner matrix 846 structure (same as flag in SLESSetOperators()) 847 . ctx - optional user-defined Jacobian context 848 849 Notes: 850 See SLESSetOperators() for information about setting the flag input 851 parameter in the routine func(). Be sure to read this information! 852 853 The routine func() takes Mat * as the matrix arguments rather than Mat. 854 This allows the Jacobian evaluation routine to replace A and/or B with a 855 completely new new matrix structure (not just different matrix elements) 856 when appropriate, for instance, if the nonzero structure is changing 857 throughout the global iterations. 858 859 .keywords: SNES, nonlinear, set, Jacobian, matrix 860 861 .seealso: SLESSetOperators(), SNESSetFunction() 862 @*/ 863 int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 864 MatStructure*,void*),void *ctx) 865 { 866 PetscValidHeaderSpecific(snes,SNES_COOKIE); 867 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 868 SETERRQ(1,"SNESSetJacobian:For SNES_NONLINEAR_EQUATIONS only"); 869 snes->computejacobian = func; 870 snes->jacP = ctx; 871 snes->jacobian = A; 872 snes->jacobian_pre = B; 873 return 0; 874 } 875 876 /*@ 877 SNESGetJacobian - Returns the Jacobian matrix and optionally the user 878 provided context for evaluating the Jacobian. 879 880 Input Parameter: 881 . snes - the nonlinear solver context 882 883 Output Parameters: 884 . A - location to stash Jacobian matrix (or PETSC_NULL) 885 . B - location to stash preconditioner matrix (or PETSC_NULL) 886 . ctx - location to stash Jacobian ctx (or PETSC_NULL) 887 888 .seealso: SNESSetJacobian(), SNESComputeJacobian() 889 @*/ 890 int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx) 891 { 892 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 893 SETERRQ(1,"SNESSetJacobian:For SNES_NONLINEAR_EQUATIONS only"); 894 if (A) *A = snes->jacobian; 895 if (B) *B = snes->jacobian_pre; 896 if (ctx) *ctx = snes->jacP; 897 return 0; 898 } 899 900 /*@C 901 SNESSetHessian - Sets the function to compute Hessian as well as the 902 location to store it. 903 904 Input Parameters: 905 . snes - the SNES context 906 . A - Hessian matrix 907 . B - preconditioner matrix (usually same as the Hessian) 908 . func - Jacobian evaluation routine 909 . ctx - optional user-defined context for private data for the 910 Hessian evaluation routine (may be null) 911 912 Calling sequence of func: 913 . func (SNES,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 914 915 . x - input vector 916 . A - Hessian matrix 917 . B - preconditioner matrix, usually the same as A 918 . flag - flag indicating information about the preconditioner matrix 919 structure (same as flag in SLESSetOperators()) 920 . ctx - optional user-defined Hessian context 921 922 Notes: 923 See SLESSetOperators() for information about setting the flag input 924 parameter in the routine func(). Be sure to read this information! 925 926 The function func() takes Mat * as the matrix arguments rather than Mat. 927 This allows the Hessian evaluation routine to replace A and/or B with a 928 completely new new matrix structure (not just different matrix elements) 929 when appropriate, for instance, if the nonzero structure is changing 930 throughout the global iterations. 931 932 .keywords: SNES, nonlinear, set, Hessian, matrix 933 934 .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators() 935 @*/ 936 int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 937 MatStructure*,void*),void *ctx) 938 { 939 PetscValidHeaderSpecific(snes,SNES_COOKIE); 940 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 941 SETERRQ(1,"SNESSetHessian:For SNES_UNCONSTRAINED_MINIMIZATION only"); 942 snes->computejacobian = func; 943 snes->jacP = ctx; 944 snes->jacobian = A; 945 snes->jacobian_pre = B; 946 return 0; 947 } 948 949 /*@ 950 SNESGetHessian - Returns the Hessian matrix and optionally the user 951 provided context for evaluating the Hessian. 952 953 Input Parameter: 954 . snes - the nonlinear solver context 955 956 Output Parameters: 957 . A - location to stash Hessian matrix (or PETSC_NULL) 958 . B - location to stash preconditioner matrix (or PETSC_NULL) 959 . ctx - location to stash Hessian ctx (or PETSC_NULL) 960 961 .seealso: SNESSetHessian(), SNESComputeHessian() 962 @*/ 963 int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx) 964 { 965 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 966 SETERRQ(1,"SNESSetHessian:For SNES_UNCONSTRAINED_MINIMIZATION only"); 967 if (A) *A = snes->jacobian; 968 if (B) *B = snes->jacobian_pre; 969 if (ctx) *ctx = snes->jacP; 970 return 0; 971 } 972 973 /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 974 975 /*@ 976 SNESSetUp - Sets up the internal data structures for the later use 977 of a nonlinear solver. 978 979 Input Parameter: 980 . snes - the SNES context 981 . x - the solution vector 982 983 Notes: 984 For basic use of the SNES solvers the user need not explicitly call 985 SNESSetUp(), since these actions will automatically occur during 986 the call to SNESSolve(). However, if one wishes to control this 987 phase separately, SNESSetUp() should be called after SNESCreate() 988 and optional routines of the form SNESSetXXX(), but before SNESSolve(). 989 990 .keywords: SNES, nonlinear, setup 991 992 .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 993 @*/ 994 int SNESSetUp(SNES snes,Vec x) 995 { 996 int ierr, flg; 997 PetscValidHeaderSpecific(snes,SNES_COOKIE); 998 PetscValidHeaderSpecific(x,VEC_COOKIE); 999 snes->vec_sol = snes->vec_sol_always = x; 1000 1001 ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg); CHKERRQ(ierr); 1002 /* 1003 This version replaces the user provided Jacobian matrix with a 1004 matrix free version but still uses the user provided preconditioner matrix 1005 */ 1006 if (flg) { 1007 Mat J; 1008 ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1009 PLogObjectParent(snes,J); 1010 snes->mfshell = J; 1011 if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 1012 snes->jacobian = J; 1013 PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n"); 1014 } 1015 else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 1016 snes->jacobian = J; 1017 PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n"); 1018 } else SETERRQ(1,"SNESSetUp:Method class doesn't support matrix-free operator option"); 1019 } 1020 ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg); CHKERRQ(ierr); 1021 /* 1022 This version replaces both the user provided Jacobian and the user 1023 provided preconditioner matrix with the default matrix free version. 1024 */ 1025 if (flg) { 1026 Mat J; 1027 ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1028 PLogObjectParent(snes,J); 1029 snes->mfshell = J; 1030 if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 1031 ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 1032 PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n"); 1033 } 1034 else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 1035 ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 1036 PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n"); 1037 } else SETERRQ(1,"SNESSetUp:Method class doesn't support matrix-free option"); 1038 } 1039 if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) { 1040 if (!snes->vec_func) SETERRQ(1,"SNESSetUp:Must call SNESSetFunction() first"); 1041 if (!snes->computefunction) SETERRQ(1,"SNESSetUp:Must call SNESSetFunction() first"); 1042 if (!snes->jacobian) SETERRQ(1,"SNESSetUp:Must call SNESSetJacobian() first"); 1043 if (snes->vec_func == snes->vec_sol) SETERRQ(1,"SNESSetUp:Solution vector cannot be function vector"); 1044 1045 /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 1046 if (snes->ksp_ewconv && snes->type != SNES_EQ_TR) { 1047 SLES sles; KSP ksp; 1048 ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 1049 ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr); 1050 ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private, 1051 (void *)snes); CHKERRQ(ierr); 1052 } 1053 } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) { 1054 if (!snes->vec_func) SETERRQ(1,"SNESSetUp:Must call SNESSetGradient() first"); 1055 if (!snes->computefunction) SETERRQ(1,"SNESSetUp:Must call SNESSetGradient() first"); 1056 if (!snes->computeumfunction) 1057 SETERRQ(1,"SNESSetUp:Must call SNESSetMinimizationFunction() first"); 1058 if (!snes->jacobian) SETERRQ(1,"SNESSetUp:Must call SNESSetHessian() first"); 1059 } else SETERRQ(1,"SNESSetUp:Unknown method class"); 1060 if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);} 1061 snes->setup_called = 1; 1062 return 0; 1063 } 1064 1065 /*@C 1066 SNESDestroy - Destroys the nonlinear solver context that was created 1067 with SNESCreate(). 1068 1069 Input Parameter: 1070 . snes - the SNES context 1071 1072 .keywords: SNES, nonlinear, destroy 1073 1074 .seealso: SNESCreate(), SNESSolve() 1075 @*/ 1076 int SNESDestroy(SNES snes) 1077 { 1078 int ierr; 1079 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1080 ierr = (*(snes)->destroy)((PetscObject)snes); CHKERRQ(ierr); 1081 if (snes->kspconvctx) PetscFree(snes->kspconvctx); 1082 if (snes->mfshell) MatDestroy(snes->mfshell); 1083 ierr = SLESDestroy(snes->sles); CHKERRQ(ierr); 1084 if (snes->xmonitor) SNESLGMonitorDestroy(snes->xmonitor); 1085 if (snes->vwork) VecDestroyVecs(snes->vwork,snes->nvwork); 1086 PLogObjectDestroy((PetscObject)snes); 1087 PetscHeaderDestroy((PetscObject)snes); 1088 return 0; 1089 } 1090 1091 /* ----------- Routines to set solver parameters ---------- */ 1092 1093 /*@ 1094 SNESSetTolerances - Sets various parameters used in convergence tests. 1095 1096 Input Parameters: 1097 . snes - the SNES context 1098 . atol - absolute convergence tolerance 1099 . rtol - relative convergence tolerance 1100 . stol - convergence tolerance in terms of the norm 1101 of the change in the solution between steps 1102 . maxit - maximum number of iterations 1103 . maxf - maximum number of function evaluations 1104 1105 Options Database Keys: 1106 $ -snes_atol <atol> 1107 $ -snes_rtol <rtol> 1108 $ -snes_stol <stol> 1109 $ -snes_max_it <maxit> 1110 $ -snes_max_funcs <maxf> 1111 1112 Notes: 1113 The default maximum number of iterations is 50. 1114 The default maximum number of function evaluations is 1000. 1115 1116 .keywords: SNES, nonlinear, set, convergence, tolerances 1117 1118 .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance() 1119 @*/ 1120 int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf) 1121 { 1122 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1123 if (atol != PETSC_DEFAULT) snes->atol = atol; 1124 if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1125 if (stol != PETSC_DEFAULT) snes->xtol = stol; 1126 if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1127 if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 1128 return 0; 1129 } 1130 1131 /*@ 1132 SNESGetTolerances - Gets various parameters used in convergence tests. 1133 1134 Input Parameters: 1135 . snes - the SNES context 1136 . atol - absolute convergence tolerance 1137 . rtol - relative convergence tolerance 1138 . stol - convergence tolerance in terms of the norm 1139 of the change in the solution between steps 1140 . maxit - maximum number of iterations 1141 . maxf - maximum number of function evaluations 1142 1143 Notes: 1144 The user can specify PETSC_NULL for any parameter that is not needed. 1145 1146 .keywords: SNES, nonlinear, get, convergence, tolerances 1147 1148 .seealso: SNESSetTolerances() 1149 @*/ 1150 int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf) 1151 { 1152 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1153 if (atol) *atol = snes->atol; 1154 if (rtol) *rtol = snes->rtol; 1155 if (stol) *stol = snes->xtol; 1156 if (maxit) *maxit = snes->max_its; 1157 if (maxf) *maxf = snes->max_funcs; 1158 return 0; 1159 } 1160 1161 /*@ 1162 SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 1163 1164 Input Parameters: 1165 . snes - the SNES context 1166 . tol - tolerance 1167 1168 Options Database Key: 1169 $ -snes_trtol <tol> 1170 1171 .keywords: SNES, nonlinear, set, trust region, tolerance 1172 1173 .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance() 1174 @*/ 1175 int SNESSetTrustRegionTolerance(SNES snes,double tol) 1176 { 1177 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1178 snes->deltatol = tol; 1179 return 0; 1180 } 1181 1182 /*@ 1183 SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance 1184 for unconstrained minimization solvers. 1185 1186 Input Parameters: 1187 . snes - the SNES context 1188 . ftol - minimum function tolerance 1189 1190 Options Database Key: 1191 $ -snes_fmin <ftol> 1192 1193 Note: 1194 SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION 1195 methods only. 1196 1197 .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance 1198 1199 .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance() 1200 @*/ 1201 int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol) 1202 { 1203 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1204 snes->fmin = ftol; 1205 return 0; 1206 } 1207 1208 /* ------------ Routines to set performance monitoring options ----------- */ 1209 1210 /*@C 1211 SNESSetMonitor - Sets the function that is to be used at every 1212 iteration of the nonlinear solver to display the iteration's 1213 progress. 1214 1215 Input Parameters: 1216 . snes - the SNES context 1217 . func - monitoring routine 1218 . mctx - optional user-defined context for private data for the 1219 monitor routine (may be null) 1220 1221 Calling sequence of func: 1222 int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx) 1223 1224 $ snes - the SNES context 1225 $ its - iteration number 1226 $ mctx - optional monitoring context 1227 $ 1228 $ SNES_NONLINEAR_EQUATIONS methods: 1229 $ norm - 2-norm function value (may be estimated) 1230 $ 1231 $ SNES_UNCONSTRAINED_MINIMIZATION methods: 1232 $ norm - 2-norm gradient value (may be estimated) 1233 1234 .keywords: SNES, nonlinear, set, monitor 1235 1236 .seealso: SNESDefaultMonitor() 1237 @*/ 1238 int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx ) 1239 { 1240 snes->monitor = func; 1241 snes->monP = mctx; 1242 return 0; 1243 } 1244 1245 /*@C 1246 SNESSetConvergenceTest - Sets the function that is to be used 1247 to test for convergence of the nonlinear iterative solution. 1248 1249 Input Parameters: 1250 . snes - the SNES context 1251 . func - routine to test for convergence 1252 . cctx - optional context for private data for the convergence routine 1253 (may be null) 1254 1255 Calling sequence of func: 1256 int func (SNES snes,double xnorm,double gnorm, 1257 double f,void *cctx) 1258 1259 $ snes - the SNES context 1260 $ cctx - optional convergence context 1261 $ xnorm - 2-norm of current iterate 1262 $ 1263 $ SNES_NONLINEAR_EQUATIONS methods: 1264 $ gnorm - 2-norm of current step 1265 $ f - 2-norm of function 1266 $ 1267 $ SNES_UNCONSTRAINED_MINIMIZATION methods: 1268 $ gnorm - 2-norm of current gradient 1269 $ f - function value 1270 1271 .keywords: SNES, nonlinear, set, convergence, test 1272 1273 .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(), 1274 SNESConverged_UM_LS(), SNESConverged_UM_TR() 1275 @*/ 1276 int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx) 1277 { 1278 (snes)->converged = func; 1279 (snes)->cnvP = cctx; 1280 return 0; 1281 } 1282 1283 /* 1284 SNESScaleStep_Private - Scales a step so that its length is less than the 1285 positive parameter delta. 1286 1287 Input Parameters: 1288 . snes - the SNES context 1289 . y - approximate solution of linear system 1290 . fnorm - 2-norm of current function 1291 . delta - trust region size 1292 1293 Output Parameters: 1294 . gpnorm - predicted function norm at the new point, assuming local 1295 linearization. The value is zero if the step lies within the trust 1296 region, and exceeds zero otherwise. 1297 . ynorm - 2-norm of the step 1298 1299 Note: 1300 For non-trust region methods such as SNES_EQ_LS, the parameter delta 1301 is set to be the maximum allowable step size. 1302 1303 .keywords: SNES, nonlinear, scale, step 1304 */ 1305 int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta, 1306 double *gpnorm,double *ynorm) 1307 { 1308 double norm; 1309 Scalar cnorm; 1310 VecNorm(y,NORM_2, &norm ); 1311 if (norm > *delta) { 1312 norm = *delta/norm; 1313 *gpnorm = (1.0 - norm)*(*fnorm); 1314 cnorm = norm; 1315 VecScale( &cnorm, y ); 1316 *ynorm = *delta; 1317 } else { 1318 *gpnorm = 0.0; 1319 *ynorm = norm; 1320 } 1321 return 0; 1322 } 1323 1324 /*@ 1325 SNESSolve - Solves a nonlinear system. Call SNESSolve after calling 1326 SNESCreate() and optional routines of the form SNESSetXXX(). 1327 1328 Input Parameter: 1329 . snes - the SNES context 1330 . x - the solution vector 1331 1332 Output Parameter: 1333 its - number of iterations until termination 1334 1335 Note: 1336 The user should initialize the vector, x, with the initial guess 1337 for the nonlinear solve prior to calling SNESSolve. In particular, 1338 to employ an initial guess of zero, the user should explicitly set 1339 this vector to zero by calling VecSet(). 1340 1341 .keywords: SNES, nonlinear, solve 1342 1343 .seealso: SNESCreate(), SNESDestroy() 1344 @*/ 1345 int SNESSolve(SNES snes,Vec x,int *its) 1346 { 1347 int ierr, flg; 1348 1349 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1350 PetscValidIntPointer(its); 1351 if (!snes->setup_called) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);} 1352 else {snes->vec_sol = snes->vec_sol_always = x;} 1353 PLogEventBegin(SNES_Solve,snes,0,0,0); 1354 ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr); 1355 PLogEventEnd(SNES_Solve,snes,0,0,0); 1356 ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr); 1357 if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); } 1358 return 0; 1359 } 1360 1361 /* --------- Internal routines for SNES Package --------- */ 1362 static NRList *__SNESList = 0; 1363 1364 /*@ 1365 SNESSetType - Sets the method for the nonlinear solver. 1366 1367 Input Parameters: 1368 . snes - the SNES context 1369 . method - a known method 1370 1371 Notes: 1372 See "petsc/include/snes.h" for available methods (for instance) 1373 $ Systems of nonlinear equations: 1374 $ SNES_EQ_LS - Newton's method with line search 1375 $ SNES_EQ_TR - Newton's method with trust region 1376 $ Unconstrained minimization: 1377 $ SNES_UM_TR - Newton's method with trust region 1378 $ SNES_UM_LS - Newton's method with line search 1379 1380 Options Database Command: 1381 $ -snes_type <method> 1382 $ Use -help for a list of available methods 1383 $ (for instance, ls or tr) 1384 1385 .keysords: SNES, set, method 1386 @*/ 1387 int SNESSetType(SNES snes,SNESType method) 1388 { 1389 int (*r)(SNES); 1390 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1391 if (snes->type == (int) method) return 0; 1392 1393 /* Get the function pointers for the iterative method requested */ 1394 if (!__SNESList) {SNESRegisterAll();} 1395 if (!__SNESList) {SETERRQ(1,"SNESSetType:Could not get methods");} 1396 r = (int (*)(SNES))NRFindRoutine( __SNESList, (int)method, (char *)0 ); 1397 if (!r) {SETERRQ(1,"SNESSetType:Unknown method");} 1398 if (snes->data) PetscFree(snes->data); 1399 snes->set_method_called = 1; 1400 return (*r)(snes); 1401 } 1402 1403 /* --------------------------------------------------------------------- */ 1404 /*@C 1405 SNESRegister - Adds the method to the nonlinear solver package, given 1406 a function pointer and a nonlinear solver name of the type SNESType. 1407 1408 Input Parameters: 1409 . name - for instance SNES_EQ_LS, SNES_EQ_TR, ... 1410 . sname - corresponding string for name 1411 . create - routine to create method context 1412 1413 .keywords: SNES, nonlinear, register 1414 1415 .seealso: SNESRegisterAll(), SNESRegisterDestroy() 1416 @*/ 1417 int SNESRegister(int name, char *sname, int (*create)(SNES)) 1418 { 1419 int ierr; 1420 if (!__SNESList) {ierr = NRCreate(&__SNESList); CHKERRQ(ierr);} 1421 NRRegister( __SNESList, name, sname, (int (*)(void*))create ); 1422 return 0; 1423 } 1424 /* --------------------------------------------------------------------- */ 1425 /*@C 1426 SNESRegisterDestroy - Frees the list of nonlinear solvers that were 1427 registered by SNESRegister(). 1428 1429 .keywords: SNES, nonlinear, register, destroy 1430 1431 .seealso: SNESRegisterAll(), SNESRegisterAll() 1432 @*/ 1433 int SNESRegisterDestroy() 1434 { 1435 if (__SNESList) { 1436 NRDestroy( __SNESList ); 1437 __SNESList = 0; 1438 } 1439 return 0; 1440 } 1441 1442 /* 1443 SNESGetTypeFromOptions_Private - Sets the selected method from the 1444 options database. 1445 1446 Input Parameter: 1447 . ctx - the SNES context 1448 1449 Output Parameter: 1450 . method - solver method 1451 1452 Returns: 1453 Returns 1 if the method is found; 0 otherwise. 1454 1455 Options Database Key: 1456 $ -snes_type method 1457 */ 1458 int SNESGetTypeFromOptions_Private(SNES ctx,SNESType *method,int *flg) 1459 { 1460 int ierr; 1461 char sbuf[50]; 1462 1463 ierr = OptionsGetString(ctx->prefix,"-snes_type",sbuf,50,flg);CHKERRQ(ierr); 1464 if (*flg) { 1465 if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);} 1466 *method = (SNESType)NRFindID( __SNESList, sbuf ); 1467 } 1468 return 0; 1469 } 1470 1471 /*@C 1472 SNESGetType - Gets the SNES method type and name (as a string). 1473 1474 Input Parameter: 1475 . snes - nonlinear solver context 1476 1477 Output Parameter: 1478 . method - SNES method (or use PETSC_NULL) 1479 . name - name of SNES method (or use PETSC_NULL) 1480 1481 .keywords: SNES, nonlinear, get, method, name 1482 @*/ 1483 int SNESGetType(SNES snes, SNESType *method,char **name) 1484 { 1485 int ierr; 1486 if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);} 1487 if (method) *method = (SNESType) snes->type; 1488 if (name) *name = NRFindName( __SNESList, (int) snes->type ); 1489 return 0; 1490 } 1491 1492 #include <stdio.h> 1493 /* 1494 SNESPrintTypes_Private - Prints the SNES methods available from the 1495 options database. 1496 1497 Input Parameters: 1498 . comm - communicator (usually MPI_COMM_WORLD) 1499 . prefix - prefix (usually "-") 1500 . name - the options database name (by default "snes_type") 1501 */ 1502 int SNESPrintTypes_Private(MPI_Comm comm,char* prefix,char *name) 1503 { 1504 FuncList *entry; 1505 if (!__SNESList) {SNESRegisterAll();} 1506 entry = __SNESList->head; 1507 PetscPrintf(comm," %s%s (one of)",prefix,name); 1508 while (entry) { 1509 PetscPrintf(comm," %s",entry->name); 1510 entry = entry->next; 1511 } 1512 PetscPrintf(comm,"\n"); 1513 return 0; 1514 } 1515 1516 /*@C 1517 SNESGetSolution - Returns the vector where the approximate solution is 1518 stored. 1519 1520 Input Parameter: 1521 . snes - the SNES context 1522 1523 Output Parameter: 1524 . x - the solution 1525 1526 .keywords: SNES, nonlinear, get, solution 1527 1528 .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate() 1529 @*/ 1530 int SNESGetSolution(SNES snes,Vec *x) 1531 { 1532 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1533 *x = snes->vec_sol_always; 1534 return 0; 1535 } 1536 1537 /*@C 1538 SNESGetSolutionUpdate - Returns the vector where the solution update is 1539 stored. 1540 1541 Input Parameter: 1542 . snes - the SNES context 1543 1544 Output Parameter: 1545 . x - the solution update 1546 1547 Notes: 1548 This vector is implementation dependent. 1549 1550 .keywords: SNES, nonlinear, get, solution, update 1551 1552 .seealso: SNESGetSolution(), SNESGetFunction 1553 @*/ 1554 int SNESGetSolutionUpdate(SNES snes,Vec *x) 1555 { 1556 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1557 *x = snes->vec_sol_update_always; 1558 return 0; 1559 } 1560 1561 /*@C 1562 SNESGetFunction - Returns the vector where the function is stored. 1563 1564 Input Parameter: 1565 . snes - the SNES context 1566 1567 Output Parameter: 1568 . r - the function 1569 1570 Notes: 1571 SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only 1572 Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 1573 SNESGetMinimizationFunction() and SNESGetGradient(); 1574 1575 .keywords: SNES, nonlinear, get, function 1576 1577 .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(), 1578 SNESGetGradient() 1579 @*/ 1580 int SNESGetFunction(SNES snes,Vec *r) 1581 { 1582 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1583 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1, 1584 "SNESGetFunction:For SNES_NONLINEAR_EQUATIONS only"); 1585 *r = snes->vec_func_always; 1586 return 0; 1587 } 1588 1589 /*@C 1590 SNESGetGradient - Returns the vector where the gradient is stored. 1591 1592 Input Parameter: 1593 . snes - the SNES context 1594 1595 Output Parameter: 1596 . r - the gradient 1597 1598 Notes: 1599 SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods 1600 only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 1601 SNESGetFunction(). 1602 1603 .keywords: SNES, nonlinear, get, gradient 1604 1605 .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction() 1606 @*/ 1607 int SNESGetGradient(SNES snes,Vec *r) 1608 { 1609 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1610 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1, 1611 "SNESGetGradient:For SNES_UNCONSTRAINED_MINIMIZATION only"); 1612 *r = snes->vec_func_always; 1613 return 0; 1614 } 1615 1616 /*@ 1617 SNESGetMinimizationFunction - Returns the scalar function value for 1618 unconstrained minimization problems. 1619 1620 Input Parameter: 1621 . snes - the SNES context 1622 1623 Output Parameter: 1624 . r - the function 1625 1626 Notes: 1627 SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 1628 methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 1629 SNESGetFunction(). 1630 1631 .keywords: SNES, nonlinear, get, function 1632 1633 .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction() 1634 @*/ 1635 int SNESGetMinimizationFunction(SNES snes,double *r) 1636 { 1637 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1638 PetscValidScalarPointer(r); 1639 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1, 1640 "SNESGetMinimizationFunction:For SNES_UNCONSTRAINED_MINIMIZATION only"); 1641 *r = snes->fc; 1642 return 0; 1643 } 1644 1645 /*@C 1646 SNESSetOptionsPrefix - Sets the prefix used for searching for all 1647 SNES options in the database. 1648 1649 Input Parameter: 1650 . snes - the SNES context 1651 . prefix - the prefix to prepend to all option names 1652 1653 .keywords: SNES, set, options, prefix, database 1654 1655 .seealso: SNESSetFromOptions() 1656 @*/ 1657 int SNESSetOptionsPrefix(SNES snes,char *prefix) 1658 { 1659 int ierr; 1660 1661 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1662 ierr = PetscObjectSetPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 1663 ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 1664 return 0; 1665 } 1666 1667 /*@C 1668 SNESAppendOptionsPrefix - Append to the prefix used for searching for all 1669 SNES options in the database. 1670 1671 Input Parameter: 1672 . snes - the SNES context 1673 . prefix - the prefix to prepend to all option names 1674 1675 .keywords: SNES, append, options, prefix, database 1676 1677 .seealso: SNESGetOptionsPrefix() 1678 @*/ 1679 int SNESAppendOptionsPrefix(SNES snes,char *prefix) 1680 { 1681 int ierr; 1682 1683 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1684 ierr = PetscObjectAppendPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 1685 ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 1686 return 0; 1687 } 1688 1689 /*@ 1690 SNESGetOptionsPrefix - Sets the prefix used for searching for all 1691 SNES options in the database. 1692 1693 Input Parameter: 1694 . snes - the SNES context 1695 1696 Output Parameter: 1697 . prefix - pointer to the prefix string used 1698 1699 .keywords: SNES, get, options, prefix, database 1700 1701 .seealso: SNESAppendOptionsPrefix() 1702 @*/ 1703 int SNESGetOptionsPrefix(SNES snes,char **prefix) 1704 { 1705 int ierr; 1706 1707 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1708 ierr = PetscObjectGetPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 1709 return 0; 1710 } 1711 1712 1713 1714 1715 1716