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