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