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