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