1 #ifdef PETSC_RCS_HEADER 2 static char vcid[] = "$Id: snes.c,v 1.131 1997/10/19 03:29:25 bsmith Exp bsmith $"; 3 #endif 4 5 #include "src/snes/snesimpl.h" /*I "snes.h" I*/ 6 #include "src/sys/nreg.h" 7 #include "pinclude/pviewer.h" 8 #include <math.h> 9 10 extern int SNESGetTypeFromOptions_Private(SNES,SNESType*,int*); 11 extern int SNESPrintTypes_Private(MPI_Comm,char*,char*); 12 13 int SNESRegisterAllCalled = 0; 14 15 #undef __FUNC__ 16 #define __FUNC__ "SNESView" 17 /*@ 18 SNESView - Prints the SNES data structure. 19 20 Input Parameters: 21 . SNES - the SNES context 22 . viewer - visualization context 23 24 Options Database Key: 25 $ -snes_view : calls SNESView() at end of SNESSolve() 26 27 Notes: 28 The available visualization contexts include 29 $ VIEWER_STDOUT_SELF - standard output (default) 30 $ VIEWER_STDOUT_WORLD - synchronized standard 31 $ output where only the first processor opens 32 $ the file. All other processors send their 33 $ data to the first processor to print. 34 35 The user can open alternative vistualization contexts with 36 $ ViewerFileOpenASCII() - output to a specified file 37 38 .keywords: SNES, view 39 40 .seealso: ViewerFileOpenASCII() 41 @*/ 42 int SNESView(SNES snes,Viewer viewer) 43 { 44 SNES_KSP_EW_ConvCtx *kctx; 45 FILE *fd; 46 int ierr; 47 SLES sles; 48 char *method; 49 ViewerType vtype; 50 51 PetscFunctionBegin; 52 PetscValidHeaderSpecific(snes,SNES_COOKIE); 53 if (viewer) {PetscValidHeader(viewer);} 54 else { viewer = VIEWER_STDOUT_SELF; } 55 56 ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 57 if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 58 ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 59 PetscFPrintf(snes->comm,fd,"SNES Object:\n"); 60 SNESGetType(snes,PETSC_NULL,&method); 61 PetscFPrintf(snes->comm,fd," method: %s\n",method); 62 if (snes->view) (*snes->view)((PetscObject)snes,viewer); 63 PetscFPrintf(snes->comm,fd, 64 " maximum iterations=%d, maximum function evaluations=%d\n", 65 snes->max_its,snes->max_funcs); 66 PetscFPrintf(snes->comm,fd, 67 " tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n", 68 snes->rtol, snes->atol, snes->trunctol, snes->xtol); 69 PetscFPrintf(snes->comm,fd, 70 " total number of linear solver iterations=%d\n",snes->linear_its); 71 PetscFPrintf(snes->comm,fd, 72 " total number of function evaluations=%d\n",snes->nfuncs); 73 if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) 74 PetscFPrintf(snes->comm,fd," min function tolerance=%g\n",snes->fmin); 75 if (snes->ksp_ewconv) { 76 kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 77 if (kctx) { 78 PetscFPrintf(snes->comm,fd, 79 " Eisenstat-Walker computation of KSP relative tolerance (version %d)\n", 80 kctx->version); 81 PetscFPrintf(snes->comm,fd, 82 " rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0, 83 kctx->rtol_max,kctx->threshold); 84 PetscFPrintf(snes->comm,fd," gamma=%g, alpha=%g, alpha2=%g\n", 85 kctx->gamma,kctx->alpha,kctx->alpha2); 86 } 87 } 88 } else if (vtype == STRING_VIEWER) { 89 SNESGetType(snes,PETSC_NULL,&method); 90 ViewerStringSPrintf(viewer," %-3.3s",method); 91 } 92 SNESGetSLES(snes,&sles); 93 ierr = SLESView(sles,viewer); CHKERRQ(ierr); 94 PetscFunctionReturn(0); 95 } 96 97 /* 98 We retain a list of functions that also take SNES command 99 line options. These are called at the end SNESSetFromOptions() 100 */ 101 #define MAXSETFROMOPTIONS 5 102 static int numberofsetfromoptions; 103 static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES); 104 105 #undef __FUNC__ 106 #define __FUNC__ "SNESAddOptionsChecker" 107 /*@ 108 SNESAddOptionsChecker - Adds an additional function to check for SNES options. 109 110 Input Parameter: 111 . snescheck - function that checks for options 112 113 .seealso: SNESSetFromOptions() 114 @*/ 115 int SNESAddOptionsChecker(int (*snescheck)(SNES) ) 116 { 117 PetscFunctionBegin; 118 if (numberofsetfromoptions >= MAXSETFROMOPTIONS) { 119 SETERRQ(1,0,"Too many options checkers, only 5 allowed"); 120 } 121 122 othersetfromoptions[numberofsetfromoptions++] = snescheck; 123 PetscFunctionReturn(0); 124 } 125 126 #undef __FUNC__ 127 #define __FUNC__ "SNESSetFromOptions" 128 /*@ 129 SNESSetFromOptions - Sets various SNES and SLES parameters from user options. 130 131 Input Parameter: 132 . snes - the SNES context 133 134 Notes: 135 To see all options, run your program with the -help option or consult 136 the users manual. 137 138 .keywords: SNES, nonlinear, set, options, database 139 140 .seealso: SNESPrintHelp(), SNESSetOptionsPrefix() 141 @*/ 142 int SNESSetFromOptions(SNES snes) 143 { 144 SNESType method; 145 double tmp; 146 SLES sles; 147 int ierr, flg,i,loc[4],nmax = 4; 148 int version = PETSC_DEFAULT; 149 double rtol_0 = PETSC_DEFAULT; 150 double rtol_max = PETSC_DEFAULT; 151 double gamma2 = PETSC_DEFAULT; 152 double alpha = PETSC_DEFAULT; 153 double alpha2 = PETSC_DEFAULT; 154 double threshold = PETSC_DEFAULT; 155 156 PetscFunctionBegin; 157 loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300; 158 159 PetscValidHeaderSpecific(snes,SNES_COOKIE); 160 if (snes->setup_called) SETERRQ(1,0,"Must call prior to SNESSetUp"); 161 ierr = SNESGetTypeFromOptions_Private(snes,&method,&flg); CHKERRQ(ierr); 162 if (flg) { 163 ierr = SNESSetType(snes,method); CHKERRQ(ierr); 164 } 165 else if (!snes->set_method_called) { 166 if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 167 ierr = SNESSetType(snes,SNES_EQ_LS); CHKERRQ(ierr); 168 } else { 169 ierr = SNESSetType(snes,SNES_UM_TR); CHKERRQ(ierr); 170 } 171 } 172 173 ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 174 if (flg) { ierr = SNESPrintHelp(snes); CHKERRQ(ierr);} 175 ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg); CHKERRQ(ierr); 176 if (flg) { 177 ierr = SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 178 } 179 ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg); CHKERRQ(ierr); 180 if (flg) { 181 ierr = SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 182 } 183 ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg); CHKERRQ(ierr); 184 if (flg) { 185 ierr = SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 186 } 187 ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, &flg); CHKERRQ(ierr); 188 ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, &flg);CHKERRQ(ierr); 189 ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg); CHKERRQ(ierr); 190 if (flg) { ierr = SNESSetTrustRegionTolerance(snes,tmp); CHKERRQ(ierr); } 191 ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg); CHKERRQ(ierr); 192 if (flg) { ierr = SNESSetMinimizationFunctionTolerance(snes,tmp); CHKERRQ(ierr);} 193 ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg); CHKERRQ(ierr); 194 if (flg) { snes->ksp_ewconv = 1; } 195 ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,&flg); CHKERRQ(ierr); 196 ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,&flg); CHKERRQ(ierr); 197 ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,&flg);CHKERRQ(ierr); 198 ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,&flg); CHKERRQ(ierr); 199 ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,&flg); CHKERRQ(ierr); 200 ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,&flg); CHKERRQ(ierr); 201 ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,&flg);CHKERRQ(ierr); 202 203 ierr = OptionsHasName(snes->prefix,"-snes_no_convergence_test",&flg); CHKERRQ(ierr); 204 if (flg) {snes->converged = 0;} 205 206 ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha, 207 alpha2,threshold); CHKERRQ(ierr); 208 ierr = OptionsHasName(snes->prefix,"-snes_cancelmonitors",&flg); CHKERRQ(ierr); 209 if (flg) {ierr = SNESSetMonitor(snes,0,0);CHKERRQ(ierr);} 210 ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg); CHKERRQ(ierr); 211 if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0);CHKERRQ(ierr);} 212 ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg); CHKERRQ(ierr); 213 if (flg) { 214 ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0); CHKERRQ(ierr); 215 } 216 ierr = OptionsGetIntArray(snes->prefix,"-snes_xmonitor",loc,&nmax,&flg);CHKERRQ(ierr); 217 if (flg) { 218 int rank = 0; 219 DrawLG lg; 220 MPI_Initialized(&rank); 221 if (rank) MPI_Comm_rank(snes->comm,&rank); 222 if (!rank) { 223 ierr = SNESLGMonitorCreate(0,0,loc[0],loc[1],loc[2],loc[3],&lg); CHKERRQ(ierr); 224 ierr = SNESSetMonitor(snes,SNESLGMonitor,(void *)lg); CHKERRQ(ierr); 225 snes->xmonitor = lg; 226 PLogObjectParent(snes,lg); 227 } 228 } 229 ierr = OptionsHasName(snes->prefix,"-snes_fd", &flg); CHKERRQ(ierr); 230 if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) { 231 ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre, 232 SNESDefaultComputeJacobian,snes->funP); CHKERRQ(ierr); 233 PLogInfo(snes, 234 "SNESSetFromOptions: Setting default finite difference Jacobian matrix\n"); 235 } 236 else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 237 ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre, 238 SNESDefaultComputeHessian,snes->funP); CHKERRQ(ierr); 239 PLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Hessian matrix\n"); 240 } 241 242 for ( i=0; i<numberofsetfromoptions; i++ ) { 243 ierr = (*othersetfromoptions[i])(snes); CHKERRQ(ierr); 244 } 245 246 ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 247 ierr = SLESSetFromOptions(sles); CHKERRQ(ierr); 248 if (snes->setfromoptions) { 249 ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr); 250 } 251 PetscFunctionReturn(0); 252 } 253 254 #undef __FUNC__ 255 #define __FUNC__ "SNESPrintHelp" 256 /*@ 257 SNESPrintHelp - Prints all options for the SNES component. 258 259 Input Parameter: 260 . snes - the SNES context 261 262 Options Database Keys: 263 $ -help, -h 264 265 .keywords: SNES, nonlinear, help 266 267 .seealso: SNESSetFromOptions() 268 @*/ 269 int SNESPrintHelp(SNES snes) 270 { 271 char p[64]; 272 SNES_KSP_EW_ConvCtx *kctx; 273 int ierr; 274 275 PetscFunctionBegin; 276 PetscValidHeaderSpecific(snes,SNES_COOKIE); 277 278 PetscStrcpy(p,"-"); 279 if (snes->prefix) PetscStrcat(p, snes->prefix); 280 281 kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 282 283 PetscPrintf(snes->comm,"SNES options ------------------------------------------------\n"); 284 SNESPrintTypes_Private(snes->comm,p,"snes_type"); 285 PetscPrintf(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p); 286 PetscPrintf(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its); 287 PetscPrintf(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs); 288 PetscPrintf(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol); 289 PetscPrintf(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol); 290 PetscPrintf(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol); 291 PetscPrintf(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol); 292 PetscPrintf(snes->comm," SNES Monitoring Options: Choose any of the following\n"); 293 PetscPrintf(snes->comm," %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p); 294 PetscPrintf(snes->comm," %ssnes_monitor: use default SNES convergence monitor, prints\n\ 295 residual norm at each iteration.\n",p); 296 PetscPrintf(snes->comm," %ssnes_smonitor: same as the above, but prints fewer digits of the\n\ 297 residual norm for small residual norms. This is useful to conceal\n\ 298 meaningless digits that may be different on different machines.\n",p); 299 PetscPrintf(snes->comm," %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p); 300 if (snes->type == SNES_NONLINEAR_EQUATIONS) { 301 PetscPrintf(snes->comm, 302 " Options for solving systems of nonlinear equations only:\n"); 303 PetscPrintf(snes->comm," %ssnes_fd: use finite differences for Jacobian\n",p); 304 PetscPrintf(snes->comm," %ssnes_mf: use matrix-free Jacobian\n",p); 305 PetscPrintf(snes->comm," %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p); 306 PetscPrintf(snes->comm," %ssnes_no_convergence_test: Do not test for convergence, always run to SNES max its\n",p); 307 PetscPrintf(snes->comm," %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p); 308 PetscPrintf(snes->comm, 309 " %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version); 310 PetscPrintf(snes->comm, 311 " %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0); 312 PetscPrintf(snes->comm, 313 " %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max); 314 PetscPrintf(snes->comm, 315 " %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma); 316 PetscPrintf(snes->comm, 317 " %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha); 318 PetscPrintf(snes->comm, 319 " %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2); 320 PetscPrintf(snes->comm, 321 " %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold); 322 } else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) { 323 PetscPrintf(snes->comm, 324 " Options for solving unconstrained minimization problems only:\n"); 325 PetscPrintf(snes->comm," %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin); 326 PetscPrintf(snes->comm," %ssnes_fd: use finite differences for Hessian\n",p); 327 PetscPrintf(snes->comm," %ssnes_mf: use matrix-free Hessian\n",p); 328 } 329 PetscPrintf(snes->comm," Run program with -help %ssnes_type <method> for help on ",p); 330 PetscPrintf(snes->comm,"a particular method\n"); 331 if (snes->printhelp) { 332 ierr = (*snes->printhelp)(snes,p);CHKERRQ(ierr); 333 } 334 PetscFunctionReturn(0); 335 } 336 337 #undef __FUNC__ 338 #define __FUNC__ "SNESSetApplicationContext" 339 /*@ 340 SNESSetApplicationContext - Sets the optional user-defined context for 341 the nonlinear solvers. 342 343 Input Parameters: 344 . snes - the SNES context 345 . usrP - optional user context 346 347 .keywords: SNES, nonlinear, set, application, context 348 349 .seealso: SNESGetApplicationContext() 350 @*/ 351 int SNESSetApplicationContext(SNES snes,void *usrP) 352 { 353 PetscFunctionBegin; 354 PetscValidHeaderSpecific(snes,SNES_COOKIE); 355 snes->user = usrP; 356 PetscFunctionReturn(0); 357 } 358 359 #undef __FUNC__ 360 #define __FUNC__ "SNESGetApplicationContext" 361 /*@C 362 SNESGetApplicationContext - Gets the user-defined context for the 363 nonlinear solvers. 364 365 Input Parameter: 366 . snes - SNES context 367 368 Output Parameter: 369 . usrP - user context 370 371 .keywords: SNES, nonlinear, get, application, context 372 373 .seealso: SNESSetApplicationContext() 374 @*/ 375 int SNESGetApplicationContext( SNES snes, void **usrP ) 376 { 377 PetscFunctionBegin; 378 PetscValidHeaderSpecific(snes,SNES_COOKIE); 379 *usrP = snes->user; 380 PetscFunctionReturn(0); 381 } 382 383 #undef __FUNC__ 384 #define __FUNC__ "SNESGetIterationNumber" 385 /*@ 386 SNESGetIterationNumber - Gets the current iteration number of the 387 nonlinear solver. 388 389 Input Parameter: 390 . snes - SNES context 391 392 Output Parameter: 393 . iter - iteration number 394 395 .keywords: SNES, nonlinear, get, iteration, number 396 @*/ 397 int SNESGetIterationNumber(SNES snes,int* iter) 398 { 399 PetscFunctionBegin; 400 PetscValidHeaderSpecific(snes,SNES_COOKIE); 401 PetscValidIntPointer(iter); 402 *iter = snes->iter; 403 PetscFunctionReturn(0); 404 } 405 406 #undef __FUNC__ 407 #define __FUNC__ "SNESGetFunctionNorm" 408 /*@ 409 SNESGetFunctionNorm - Gets the norm of the current function that was set 410 with SNESSSetFunction(). 411 412 Input Parameter: 413 . snes - SNES context 414 415 Output Parameter: 416 . fnorm - 2-norm of function 417 418 Note: 419 SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only. 420 A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is 421 SNESGetGradientNorm(). 422 423 .keywords: SNES, nonlinear, get, function, norm 424 425 .seealso: SNESSetFunction() 426 @*/ 427 int SNESGetFunctionNorm(SNES snes,Scalar *fnorm) 428 { 429 PetscFunctionBegin; 430 PetscValidHeaderSpecific(snes,SNES_COOKIE); 431 PetscValidScalarPointer(fnorm); 432 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 433 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_NONLINEAR_EQUATIONS only"); 434 } 435 *fnorm = snes->norm; 436 PetscFunctionReturn(0); 437 } 438 439 #undef __FUNC__ 440 #define __FUNC__ "SNESGetGradientNorm" 441 /*@ 442 SNESGetGradientNorm - Gets the norm of the current gradient that was set 443 with SNESSSetGradient(). 444 445 Input Parameter: 446 . snes - SNES context 447 448 Output Parameter: 449 . fnorm - 2-norm of gradient 450 451 Note: 452 SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION 453 methods only. A related routine for SNES_NONLINEAR_EQUATIONS methods 454 is SNESGetFunctionNorm(). 455 456 .keywords: SNES, nonlinear, get, gradient, norm 457 458 .seelso: SNESSetGradient() 459 @*/ 460 int SNESGetGradientNorm(SNES snes,Scalar *gnorm) 461 { 462 PetscFunctionBegin; 463 PetscValidHeaderSpecific(snes,SNES_COOKIE); 464 PetscValidScalarPointer(gnorm); 465 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 466 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 467 } 468 *gnorm = snes->norm; 469 PetscFunctionReturn(0); 470 } 471 472 #undef __FUNC__ 473 #define __FUNC__ "SNESGetNumberUnsuccessfulSteps" 474 /*@ 475 SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps 476 attempted by the nonlinear solver. 477 478 Input Parameter: 479 . snes - SNES context 480 481 Output Parameter: 482 . nfails - number of unsuccessful steps attempted 483 484 Notes: 485 This counter is reset to zero for each successive call to SNESSolve(). 486 487 .keywords: SNES, nonlinear, get, number, unsuccessful, steps 488 @*/ 489 int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails) 490 { 491 PetscFunctionBegin; 492 PetscValidHeaderSpecific(snes,SNES_COOKIE); 493 PetscValidIntPointer(nfails); 494 *nfails = snes->nfailures; 495 PetscFunctionReturn(0); 496 } 497 498 #undef __FUNC__ 499 #define __FUNC__ "SNESGetNumberLinearIterations" 500 /*@ 501 SNESGetNumberLinearIterations - Gets the total number of linear iterations 502 used by the nonlinear solver. 503 504 Input Parameter: 505 . snes - SNES context 506 507 Output Parameter: 508 . lits - number of linear iterations 509 510 Notes: 511 This counter is reset to zero for each successive call to SNESSolve(). 512 513 .keywords: SNES, nonlinear, get, number, linear, iterations 514 @*/ 515 int SNESGetNumberLinearIterations(SNES snes,int* lits) 516 { 517 PetscFunctionBegin; 518 PetscValidHeaderSpecific(snes,SNES_COOKIE); 519 PetscValidIntPointer(lits); 520 *lits = snes->linear_its; 521 PetscFunctionReturn(0); 522 } 523 524 #undef __FUNC__ 525 #define __FUNC__ "SNESGetSLES" 526 /*@C 527 SNESGetSLES - Returns the SLES context for a SNES solver. 528 529 Input Parameter: 530 . snes - the SNES context 531 532 Output Parameter: 533 . sles - the SLES context 534 535 Notes: 536 The user can then directly manipulate the SLES context to set various 537 options, etc. Likewise, the user can then extract and manipulate the 538 KSP and PC contexts as well. 539 540 .keywords: SNES, nonlinear, get, SLES, context 541 542 .seealso: SLESGetPC(), SLESGetKSP() 543 @*/ 544 int SNESGetSLES(SNES snes,SLES *sles) 545 { 546 PetscFunctionBegin; 547 PetscValidHeaderSpecific(snes,SNES_COOKIE); 548 *sles = snes->sles; 549 PetscFunctionReturn(0); 550 } 551 /* -----------------------------------------------------------*/ 552 #undef __FUNC__ 553 #define __FUNC__ "SNESCreate" 554 /*@C 555 SNESCreate - Creates a nonlinear solver context. 556 557 Input Parameter: 558 . comm - MPI communicator 559 . type - type of method, one of 560 $ SNES_NONLINEAR_EQUATIONS 561 $ (for systems of nonlinear equations) 562 $ SNES_UNCONSTRAINED_MINIMIZATION 563 $ (for unconstrained minimization) 564 565 Output Parameter: 566 . outsnes - the new SNES context 567 568 Options Database Key: 569 $ -snes_mf - use default matrix-free Jacobian-vector products, 570 $ and no preconditioning matrix 571 $ -snes_mf_operator - use default matrix-free Jacobian-vector 572 $ products, and a user-provided preconditioning matrix 573 $ as set by SNESSetJacobian() 574 $ -snes_fd - use (slow!) finite differences to compute Jacobian 575 576 .keywords: SNES, nonlinear, create, context 577 578 .seealso: SNESSolve(), SNESDestroy() 579 @*/ 580 int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes) 581 { 582 int ierr; 583 SNES snes; 584 SNES_KSP_EW_ConvCtx *kctx; 585 586 PetscFunctionBegin; 587 *outsnes = 0; 588 if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS){ 589 SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"incorrect method type"); 590 } 591 PetscHeaderCreate(snes,_p_SNES,SNES_COOKIE,SNES_UNKNOWN_METHOD,comm,SNESDestroy,SNESView); 592 PLogObjectCreate(snes); 593 snes->max_its = 50; 594 snes->max_funcs = 10000; 595 snes->norm = 0.0; 596 if (type == SNES_UNCONSTRAINED_MINIMIZATION) { 597 snes->rtol = 1.e-8; 598 snes->ttol = 0.0; 599 snes->atol = 1.e-10; 600 } else { 601 snes->rtol = 1.e-8; 602 snes->ttol = 0.0; 603 snes->atol = 1.e-50; 604 } 605 snes->xtol = 1.e-8; 606 snes->trunctol = 1.e-12; /* no longer used */ 607 snes->nfuncs = 0; 608 snes->nfailures = 0; 609 snes->linear_its = 0; 610 snes->numbermonitors = 0; 611 snes->data = 0; 612 snes->view = 0; 613 snes->computeumfunction = 0; 614 snes->umfunP = 0; 615 snes->fc = 0; 616 snes->deltatol = 1.e-12; 617 snes->fmin = -1.e30; 618 snes->method_class = type; 619 snes->set_method_called = 0; 620 snes->setup_called = 0; 621 snes->ksp_ewconv = 0; 622 snes->type = -1; 623 snes->xmonitor = 0; 624 snes->vwork = 0; 625 snes->nwork = 0; 626 snes->conv_hist_size = 0; 627 snes->conv_act_size = 0; 628 snes->conv_hist = 0; 629 630 /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 631 kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx); 632 PLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx)); 633 snes->kspconvctx = (void*)kctx; 634 kctx->version = 2; 635 kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 636 this was too large for some test cases */ 637 kctx->rtol_last = 0; 638 kctx->rtol_max = .9; 639 kctx->gamma = 1.0; 640 kctx->alpha2 = .5*(1.0 + sqrt(5.0)); 641 kctx->alpha = kctx->alpha2; 642 kctx->threshold = .1; 643 kctx->lresid_last = 0; 644 kctx->norm_last = 0; 645 646 ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr); 647 PLogObjectParent(snes,snes->sles) 648 649 *outsnes = snes; 650 PetscFunctionReturn(0); 651 } 652 653 /* --------------------------------------------------------------- */ 654 #undef __FUNC__ 655 #define __FUNC__ "SNESSetFunction" 656 /*@C 657 SNESSetFunction - Sets the function evaluation routine and function 658 vector for use by the SNES routines in solving systems of nonlinear 659 equations. 660 661 Input Parameters: 662 . snes - the SNES context 663 . func - function evaluation routine 664 . r - vector to store function value 665 . ctx - [optional] user-defined context for private data for the 666 function evaluation routine (may be PETSC_NULL) 667 668 Calling sequence of func: 669 . func (SNES snes,Vec x,Vec f,void *ctx); 670 671 . x - input vector 672 . f - function vector 673 . ctx - optional user-defined function context 674 675 Notes: 676 The Newton-like methods typically solve linear systems of the form 677 $ f'(x) x = -f(x), 678 $ where f'(x) denotes the Jacobian matrix and f(x) is the function. 679 680 SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 681 Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 682 SNESSetMinimizationFunction() and SNESSetGradient(); 683 684 .keywords: SNES, nonlinear, set, function 685 686 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 687 @*/ 688 int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx) 689 { 690 PetscFunctionBegin; 691 PetscValidHeaderSpecific(snes,SNES_COOKIE); 692 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 693 snes->computefunction = func; 694 snes->vec_func = snes->vec_func_always = r; 695 snes->funP = ctx; 696 PetscFunctionReturn(0); 697 } 698 699 #undef __FUNC__ 700 #define __FUNC__ "SNESComputeFunction" 701 /*@ 702 SNESComputeFunction - Computes the function that has been set with 703 SNESSetFunction(). 704 705 Input Parameters: 706 . snes - the SNES context 707 . x - input vector 708 709 Output Parameter: 710 . y - function vector, as set by SNESSetFunction() 711 712 Notes: 713 SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 714 Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 715 SNESComputeMinimizationFunction() and SNESComputeGradient(); 716 717 .keywords: SNES, nonlinear, compute, function 718 719 .seealso: SNESSetFunction(), SNESGetFunction() 720 @*/ 721 int SNESComputeFunction(SNES snes,Vec x, Vec y) 722 { 723 int ierr; 724 725 PetscFunctionBegin; 726 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 727 SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 728 } 729 PLogEventBegin(SNES_FunctionEval,snes,x,y,0); 730 PetscStackPush("SNES user function"); 731 ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 732 PetscStackPop; 733 snes->nfuncs++; 734 PLogEventEnd(SNES_FunctionEval,snes,x,y,0); 735 PetscFunctionReturn(0); 736 } 737 738 #undef __FUNC__ 739 #define __FUNC__ "SNESSetMinimizationFunction" 740 /*@C 741 SNESSetMinimizationFunction - Sets the function evaluation routine for 742 unconstrained minimization. 743 744 Input Parameters: 745 . snes - the SNES context 746 . func - function evaluation routine 747 . ctx - [optional] user-defined context for private data for the 748 function evaluation routine (may be PETSC_NULL) 749 750 Calling sequence of func: 751 . func (SNES snes,Vec x,double *f,void *ctx); 752 753 . x - input vector 754 . f - function 755 . ctx - [optional] user-defined function context 756 757 Notes: 758 SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 759 methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 760 SNESSetFunction(). 761 762 .keywords: SNES, nonlinear, set, minimization, function 763 764 .seealso: SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(), 765 SNESSetHessian(), SNESSetGradient() 766 @*/ 767 int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*), 768 void *ctx) 769 { 770 PetscFunctionBegin; 771 PetscValidHeaderSpecific(snes,SNES_COOKIE); 772 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION"); 773 snes->computeumfunction = func; 774 snes->umfunP = ctx; 775 PetscFunctionReturn(0); 776 } 777 778 #undef __FUNC__ 779 #define __FUNC__ "SNESComputeMinimizationFunction" 780 /*@ 781 SNESComputeMinimizationFunction - Computes the function that has been 782 set with SNESSetMinimizationFunction(). 783 784 Input Parameters: 785 . snes - the SNES context 786 . x - input vector 787 788 Output Parameter: 789 . y - function value 790 791 Notes: 792 SNESComputeMinimizationFunction() is valid only for 793 SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 794 SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 795 796 .keywords: SNES, nonlinear, compute, minimization, function 797 798 .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(), 799 SNESComputeGradient(), SNESComputeHessian() 800 @*/ 801 int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y) 802 { 803 int ierr; 804 805 PetscFunctionBegin; 806 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION"); 807 PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0); 808 PetscStackPush("SNES user minimzation function"); 809 ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr); 810 PetscStackPop; 811 snes->nfuncs++; 812 PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0); 813 PetscFunctionReturn(0); 814 } 815 816 #undef __FUNC__ 817 #define __FUNC__ "SNESSetGradient" 818 /*@C 819 SNESSetGradient - Sets the gradient evaluation routine and gradient 820 vector for use by the SNES routines. 821 822 Input Parameters: 823 . snes - the SNES context 824 . func - function evaluation routine 825 . ctx - optional user-defined context for private data for the 826 gradient evaluation routine (may be PETSC_NULL) 827 . r - vector to store gradient value 828 829 Calling sequence of func: 830 . func (SNES, Vec x, Vec g, void *ctx); 831 832 . x - input vector 833 . g - gradient vector 834 . ctx - optional user-defined gradient context 835 836 Notes: 837 SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 838 methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 839 SNESSetFunction(). 840 841 .keywords: SNES, nonlinear, set, function 842 843 .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(), 844 SNESSetMinimizationFunction(), 845 @*/ 846 int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx) 847 { 848 PetscFunctionBegin; 849 PetscValidHeaderSpecific(snes,SNES_COOKIE); 850 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 851 snes->computefunction = func; 852 snes->vec_func = snes->vec_func_always = r; 853 snes->funP = ctx; 854 PetscFunctionReturn(0); 855 } 856 857 #undef __FUNC__ 858 #define __FUNC__ "SNESComputeGradient" 859 /*@ 860 SNESComputeGradient - Computes the gradient that has been set with 861 SNESSetGradient(). 862 863 Input Parameters: 864 . snes - the SNES context 865 . x - input vector 866 867 Output Parameter: 868 . y - gradient vector 869 870 Notes: 871 SNESComputeGradient() is valid only for 872 SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 873 SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 874 875 .keywords: SNES, nonlinear, compute, gradient 876 877 .seealso: SNESSetGradient(), SNESGetGradient(), 878 SNESComputeMinimizationFunction(), SNESComputeHessian() 879 @*/ 880 int SNESComputeGradient(SNES snes,Vec x, Vec y) 881 { 882 int ierr; 883 884 PetscFunctionBegin; 885 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 886 SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 887 } 888 PLogEventBegin(SNES_GradientEval,snes,x,y,0); 889 PetscStackPush("SNES user gradient function"); 890 ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 891 PetscStackPop; 892 PLogEventEnd(SNES_GradientEval,snes,x,y,0); 893 PetscFunctionReturn(0); 894 } 895 896 #undef __FUNC__ 897 #define __FUNC__ "SNESComputeJacobian" 898 /*@ 899 SNESComputeJacobian - Computes the Jacobian matrix that has been 900 set with SNESSetJacobian(). 901 902 Input Parameters: 903 . snes - the SNES context 904 . x - input vector 905 906 Output Parameters: 907 . A - Jacobian matrix 908 . B - optional preconditioning matrix 909 . flag - flag indicating matrix structure 910 911 Notes: 912 Most users should not need to explicitly call this routine, as it 913 is used internally within the nonlinear solvers. 914 915 See SLESSetOperators() for important information about setting the 916 flag parameter. 917 918 SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS 919 methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION 920 methods is SNESComputeHessian(). 921 922 .keywords: SNES, compute, Jacobian, matrix 923 924 .seealso: SNESSetJacobian(), SLESSetOperators() 925 @*/ 926 int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 927 { 928 int ierr; 929 930 PetscFunctionBegin; 931 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 932 SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 933 } 934 if (!snes->computejacobian) PetscFunctionReturn(0); 935 PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B); 936 *flg = DIFFERENT_NONZERO_PATTERN; 937 PetscStackPush("SNES user Jacobian function"); 938 ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr); 939 PetscStackPop; 940 PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B); 941 /* make sure user returned a correct Jacobian and preconditioner */ 942 PetscValidHeaderSpecific(*A,MAT_COOKIE); 943 PetscValidHeaderSpecific(*B,MAT_COOKIE); 944 PetscFunctionReturn(0); 945 } 946 947 #undef __FUNC__ 948 #define __FUNC__ "SNESComputeHessian" 949 /*@ 950 SNESComputeHessian - Computes the Hessian matrix that has been 951 set with SNESSetHessian(). 952 953 Input Parameters: 954 . snes - the SNES context 955 . x - input vector 956 957 Output Parameters: 958 . A - Hessian matrix 959 . B - optional preconditioning matrix 960 . flag - flag indicating matrix structure 961 962 Notes: 963 Most users should not need to explicitly call this routine, as it 964 is used internally within the nonlinear solvers. 965 966 See SLESSetOperators() for important information about setting the 967 flag parameter. 968 969 SNESComputeHessian() is valid only for 970 SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 971 SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian(). 972 973 .keywords: SNES, compute, Hessian, matrix 974 975 .seealso: SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(), 976 SNESComputeMinimizationFunction() 977 @*/ 978 int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag) 979 { 980 int ierr; 981 982 PetscFunctionBegin; 983 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 984 SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 985 } 986 if (!snes->computejacobian) PetscFunctionReturn(0); 987 PLogEventBegin(SNES_HessianEval,snes,x,*A,*B); 988 *flag = DIFFERENT_NONZERO_PATTERN; 989 PetscStackPush("SNES user Hessian function"); 990 ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr); 991 PetscStackPop; 992 PLogEventEnd(SNES_HessianEval,snes,x,*A,*B); 993 /* make sure user returned a correct Jacobian and preconditioner */ 994 PetscValidHeaderSpecific(*A,MAT_COOKIE); 995 PetscValidHeaderSpecific(*B,MAT_COOKIE); 996 PetscFunctionReturn(0); 997 } 998 999 #undef __FUNC__ 1000 #define __FUNC__ "SNESSetJacobian" 1001 /*@C 1002 SNESSetJacobian - Sets the function to compute Jacobian as well as the 1003 location to store the matrix. 1004 1005 Input Parameters: 1006 . snes - the SNES context 1007 . A - Jacobian matrix 1008 . B - preconditioner matrix (usually same as the Jacobian) 1009 . func - Jacobian evaluation routine 1010 . ctx - [optional] user-defined context for private data for the 1011 Jacobian evaluation routine (may be PETSC_NULL) 1012 1013 Calling sequence of func: 1014 . func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 1015 1016 . x - input vector 1017 . A - Jacobian matrix 1018 . B - preconditioner matrix, usually the same as A 1019 . flag - flag indicating information about the preconditioner matrix 1020 structure (same as flag in SLESSetOperators()) 1021 . ctx - [optional] user-defined Jacobian context 1022 1023 Notes: 1024 See SLESSetOperators() for important information about setting the flag 1025 output parameter in the routine func(). Be sure to read this information! 1026 1027 The routine func() takes Mat * as the matrix arguments rather than Mat. 1028 This allows the Jacobian evaluation routine to replace A and/or B with a 1029 completely new new matrix structure (not just different matrix elements) 1030 when appropriate, for instance, if the nonzero structure is changing 1031 throughout the global iterations. 1032 1033 .keywords: SNES, nonlinear, set, Jacobian, matrix 1034 1035 .seealso: SLESSetOperators(), SNESSetFunction() 1036 @*/ 1037 int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 1038 MatStructure*,void*),void *ctx) 1039 { 1040 PetscFunctionBegin; 1041 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1042 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 1043 snes->computejacobian = func; 1044 snes->jacP = ctx; 1045 snes->jacobian = A; 1046 snes->jacobian_pre = B; 1047 PetscFunctionReturn(0); 1048 } 1049 1050 #undef __FUNC__ 1051 #define __FUNC__ "SNESGetJacobian" 1052 /*@ 1053 SNESGetJacobian - Returns the Jacobian matrix and optionally the user 1054 provided context for evaluating the Jacobian. 1055 1056 Input Parameter: 1057 . snes - the nonlinear solver context 1058 1059 Output Parameters: 1060 . A - location to stash Jacobian matrix (or PETSC_NULL) 1061 . B - location to stash preconditioner matrix (or PETSC_NULL) 1062 . ctx - location to stash Jacobian ctx (or PETSC_NULL) 1063 1064 .seealso: SNESSetJacobian(), SNESComputeJacobian() 1065 @*/ 1066 int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx) 1067 { 1068 PetscFunctionBegin; 1069 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 1070 SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 1071 } 1072 if (A) *A = snes->jacobian; 1073 if (B) *B = snes->jacobian_pre; 1074 if (ctx) *ctx = snes->jacP; 1075 PetscFunctionReturn(0); 1076 } 1077 1078 #undef __FUNC__ 1079 #define __FUNC__ "SNESSetHessian" 1080 /*@C 1081 SNESSetHessian - Sets the function to compute Hessian as well as the 1082 location to store the matrix. 1083 1084 Input Parameters: 1085 . snes - the SNES context 1086 . A - Hessian matrix 1087 . B - preconditioner matrix (usually same as the Hessian) 1088 . func - Jacobian evaluation routine 1089 . ctx - [optional] user-defined context for private data for the 1090 Hessian evaluation routine (may be PETSC_NULL) 1091 1092 Calling sequence of func: 1093 . func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 1094 1095 . x - input vector 1096 . A - Hessian matrix 1097 . B - preconditioner matrix, usually the same as A 1098 . flag - flag indicating information about the preconditioner matrix 1099 structure (same as flag in SLESSetOperators()) 1100 . ctx - [optional] user-defined Hessian context 1101 1102 Notes: 1103 See SLESSetOperators() for important information about setting the flag 1104 output parameter in the routine func(). Be sure to read this information! 1105 1106 The function func() takes Mat * as the matrix arguments rather than Mat. 1107 This allows the Hessian evaluation routine to replace A and/or B with a 1108 completely new new matrix structure (not just different matrix elements) 1109 when appropriate, for instance, if the nonzero structure is changing 1110 throughout the global iterations. 1111 1112 .keywords: SNES, nonlinear, set, Hessian, matrix 1113 1114 .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators() 1115 @*/ 1116 int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 1117 MatStructure*,void*),void *ctx) 1118 { 1119 PetscFunctionBegin; 1120 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1121 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 1122 SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 1123 } 1124 snes->computejacobian = func; 1125 snes->jacP = ctx; 1126 snes->jacobian = A; 1127 snes->jacobian_pre = B; 1128 PetscFunctionReturn(0); 1129 } 1130 1131 #undef __FUNC__ 1132 #define __FUNC__ "SNESGetHessian" 1133 /*@ 1134 SNESGetHessian - Returns the Hessian matrix and optionally the user 1135 provided context for evaluating the Hessian. 1136 1137 Input Parameter: 1138 . snes - the nonlinear solver context 1139 1140 Output Parameters: 1141 . A - location to stash Hessian matrix (or PETSC_NULL) 1142 . B - location to stash preconditioner matrix (or PETSC_NULL) 1143 . ctx - location to stash Hessian ctx (or PETSC_NULL) 1144 1145 .seealso: SNESSetHessian(), SNESComputeHessian() 1146 @*/ 1147 int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx) 1148 { 1149 PetscFunctionBegin; 1150 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION){ 1151 SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 1152 } 1153 if (A) *A = snes->jacobian; 1154 if (B) *B = snes->jacobian_pre; 1155 if (ctx) *ctx = snes->jacP; 1156 PetscFunctionReturn(0); 1157 } 1158 1159 /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 1160 1161 #undef __FUNC__ 1162 #define __FUNC__ "SNESSetUp" 1163 /*@ 1164 SNESSetUp - Sets up the internal data structures for the later use 1165 of a nonlinear solver. 1166 1167 Input Parameter: 1168 . snes - the SNES context 1169 . x - the solution vector 1170 1171 Notes: 1172 For basic use of the SNES solvers the user need not explicitly call 1173 SNESSetUp(), since these actions will automatically occur during 1174 the call to SNESSolve(). However, if one wishes to control this 1175 phase separately, SNESSetUp() should be called after SNESCreate() 1176 and optional routines of the form SNESSetXXX(), but before SNESSolve(). 1177 1178 .keywords: SNES, nonlinear, setup 1179 1180 .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 1181 @*/ 1182 int SNESSetUp(SNES snes,Vec x) 1183 { 1184 int ierr, flg; 1185 1186 PetscFunctionBegin; 1187 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1188 PetscValidHeaderSpecific(x,VEC_COOKIE); 1189 snes->vec_sol = snes->vec_sol_always = x; 1190 1191 ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg); CHKERRQ(ierr); 1192 /* 1193 This version replaces the user provided Jacobian matrix with a 1194 matrix-free version but still employs the user-provided preconditioner matrix 1195 */ 1196 if (flg) { 1197 Mat J; 1198 ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1199 PLogObjectParent(snes,J); 1200 snes->mfshell = J; 1201 if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 1202 snes->jacobian = J; 1203 PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n"); 1204 } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 1205 snes->jacobian = J; 1206 PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n"); 1207 } else { 1208 SETERRQ(1,0,"Method class doesn't support matrix-free operator option"); 1209 } 1210 } 1211 ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg); CHKERRQ(ierr); 1212 /* 1213 This version replaces both the user-provided Jacobian and the user- 1214 provided preconditioner matrix with the default matrix free version. 1215 */ 1216 if (flg) { 1217 Mat J; 1218 ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1219 PLogObjectParent(snes,J); 1220 snes->mfshell = J; 1221 if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 1222 ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 1223 PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n"); 1224 } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 1225 ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 1226 PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n"); 1227 } else { 1228 SETERRQ(1,0,"Method class doesn't support matrix-free option"); 1229 } 1230 } 1231 if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) { 1232 if (!snes->vec_func) SETERRQ(1,0,"Must call SNESSetFunction() first"); 1233 if (!snes->computefunction) SETERRQ(1,0,"Must call SNESSetFunction() first"); 1234 if (!snes->jacobian) SETERRQ(1,0,"Must call SNESSetJacobian() first"); 1235 if (snes->vec_func == snes->vec_sol) SETERRQ(1,0,"Solution vector cannot be function vector"); 1236 1237 /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 1238 if (snes->ksp_ewconv && snes->type != SNES_EQ_TR) { 1239 SLES sles; KSP ksp; 1240 ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 1241 ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr); 1242 ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private, 1243 (void *)snes); CHKERRQ(ierr); 1244 } 1245 } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) { 1246 if (!snes->vec_func) SETERRQ(1,0,"Must call SNESSetGradient() first"); 1247 if (!snes->computefunction) SETERRQ(1,0,"Must call SNESSetGradient() first"); 1248 if (!snes->computeumfunction) SETERRQ(1,0,"Must call SNESSetMinimizationFunction() first"); 1249 if (!snes->jacobian) SETERRQ(1,0,"Must call SNESSetHessian() first"); 1250 } else { 1251 SETERRQ(1,0,"Unknown method class"); 1252 } 1253 if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);} 1254 snes->setup_called = 1; 1255 PetscFunctionReturn(0); 1256 } 1257 1258 #undef __FUNC__ 1259 #define __FUNC__ "SNESDestroy" 1260 /*@C 1261 SNESDestroy - Destroys the nonlinear solver context that was created 1262 with SNESCreate(). 1263 1264 Input Parameter: 1265 . snes - the SNES context 1266 1267 .keywords: SNES, nonlinear, destroy 1268 1269 .seealso: SNESCreate(), SNESSolve() 1270 @*/ 1271 int SNESDestroy(SNES snes) 1272 { 1273 int ierr; 1274 1275 PetscFunctionBegin; 1276 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1277 if (--snes->refct > 0) PetscFunctionReturn(0); 1278 1279 if (snes->destroy) {ierr = (*(snes)->destroy)((PetscObject)snes); CHKERRQ(ierr);} 1280 if (snes->kspconvctx) PetscFree(snes->kspconvctx); 1281 if (snes->mfshell) MatDestroy(snes->mfshell); 1282 ierr = SLESDestroy(snes->sles); CHKERRQ(ierr); 1283 if (snes->xmonitor) SNESLGMonitorDestroy(snes->xmonitor); 1284 if (snes->vwork) VecDestroyVecs(snes->vwork,snes->nvwork); 1285 PLogObjectDestroy((PetscObject)snes); 1286 PetscHeaderDestroy((PetscObject)snes); 1287 PetscFunctionReturn(0); 1288 } 1289 1290 /* ----------- Routines to set solver parameters ---------- */ 1291 1292 #undef __FUNC__ 1293 #define __FUNC__ "SNESSetTolerances" 1294 /*@ 1295 SNESSetTolerances - Sets various parameters used in convergence tests. 1296 1297 Input Parameters: 1298 . snes - the SNES context 1299 . atol - absolute convergence tolerance 1300 . rtol - relative convergence tolerance 1301 . stol - convergence tolerance in terms of the norm 1302 of the change in the solution between steps 1303 . maxit - maximum number of iterations 1304 . maxf - maximum number of function evaluations 1305 1306 Options Database Keys: 1307 $ -snes_atol <atol> 1308 $ -snes_rtol <rtol> 1309 $ -snes_stol <stol> 1310 $ -snes_max_it <maxit> 1311 $ -snes_max_funcs <maxf> 1312 1313 Notes: 1314 The default maximum number of iterations is 50. 1315 The default maximum number of function evaluations is 1000. 1316 1317 .keywords: SNES, nonlinear, set, convergence, tolerances 1318 1319 .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance() 1320 @*/ 1321 int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf) 1322 { 1323 PetscFunctionBegin; 1324 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1325 if (atol != PETSC_DEFAULT) snes->atol = atol; 1326 if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1327 if (stol != PETSC_DEFAULT) snes->xtol = stol; 1328 if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1329 if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 1330 PetscFunctionReturn(0); 1331 } 1332 1333 #undef __FUNC__ 1334 #define __FUNC__ "SNESGetTolerances" 1335 /*@ 1336 SNESGetTolerances - Gets various parameters used in convergence tests. 1337 1338 Input Parameters: 1339 . snes - the SNES context 1340 . atol - absolute convergence tolerance 1341 . rtol - relative convergence tolerance 1342 . stol - convergence tolerance in terms of the norm 1343 of the change in the solution between steps 1344 . maxit - maximum number of iterations 1345 . maxf - maximum number of function evaluations 1346 1347 Notes: 1348 The user can specify PETSC_NULL for any parameter that is not needed. 1349 1350 .keywords: SNES, nonlinear, get, convergence, tolerances 1351 1352 .seealso: SNESSetTolerances() 1353 @*/ 1354 int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf) 1355 { 1356 PetscFunctionBegin; 1357 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1358 if (atol) *atol = snes->atol; 1359 if (rtol) *rtol = snes->rtol; 1360 if (stol) *stol = snes->xtol; 1361 if (maxit) *maxit = snes->max_its; 1362 if (maxf) *maxf = snes->max_funcs; 1363 PetscFunctionReturn(0); 1364 } 1365 1366 #undef __FUNC__ 1367 #define __FUNC__ "SNESSetTrustRegionTolerance" 1368 /*@ 1369 SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 1370 1371 Input Parameters: 1372 . snes - the SNES context 1373 . tol - tolerance 1374 1375 Options Database Key: 1376 $ -snes_trtol <tol> 1377 1378 .keywords: SNES, nonlinear, set, trust region, tolerance 1379 1380 .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance() 1381 @*/ 1382 int SNESSetTrustRegionTolerance(SNES snes,double tol) 1383 { 1384 PetscFunctionBegin; 1385 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1386 snes->deltatol = tol; 1387 PetscFunctionReturn(0); 1388 } 1389 1390 #undef __FUNC__ 1391 #define __FUNC__ "SNESSetMinimizationFunctionTolerance" 1392 /*@ 1393 SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance 1394 for unconstrained minimization solvers. 1395 1396 Input Parameters: 1397 . snes - the SNES context 1398 . ftol - minimum function tolerance 1399 1400 Options Database Key: 1401 $ -snes_fmin <ftol> 1402 1403 Note: 1404 SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION 1405 methods only. 1406 1407 .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance 1408 1409 .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance() 1410 @*/ 1411 int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol) 1412 { 1413 PetscFunctionBegin; 1414 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1415 snes->fmin = ftol; 1416 PetscFunctionReturn(0); 1417 } 1418 1419 /* ------------ Routines to set performance monitoring options ----------- */ 1420 1421 #undef __FUNC__ 1422 #define __FUNC__ "SNESSetMonitor" 1423 /*@C 1424 SNESSetMonitor - Sets the function that is to be used at every 1425 iteration of the nonlinear solver to display the iteration's 1426 progress. 1427 1428 Input Parameters: 1429 . snes - the SNES context 1430 . func - monitoring routine 1431 . mctx - [optional] user-defined context for private data for the 1432 monitor routine (may be PETSC_NULL) 1433 1434 Calling sequence of func: 1435 int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx) 1436 1437 $ snes - the SNES context 1438 $ its - iteration number 1439 $ mctx - [optional] monitoring context 1440 $ 1441 $ SNES_NONLINEAR_EQUATIONS methods: 1442 $ norm - 2-norm function value (may be estimated) 1443 $ 1444 $ SNES_UNCONSTRAINED_MINIMIZATION methods: 1445 $ norm - 2-norm gradient value (may be estimated) 1446 1447 Options Database Keys: 1448 $ -snes_monitor : sets SNESDefaultMonitor() 1449 $ -snes_xmonitor : sets line graph monitor, 1450 $ uses SNESLGMonitorCreate() 1451 $ -snes_cancelmonitors : cancels all monitors that have 1452 $ been hardwired into a code by 1453 $ calls to SNESSetMonitor(), but 1454 $ does not cancel those set via 1455 $ the options database. 1456 1457 1458 Notes: 1459 Several different monitoring routines may be set by calling 1460 SNESSetMonitor() multiple times; all will be called in the 1461 order in which they were set. 1462 1463 .keywords: SNES, nonlinear, set, monitor 1464 1465 .seealso: SNESDefaultMonitor() 1466 @*/ 1467 int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx ) 1468 { 1469 PetscFunctionBegin; 1470 if (!func) { 1471 snes->numbermonitors = 0; 1472 PetscFunctionReturn(0); 1473 } 1474 if (snes->numbermonitors >= MAXSNESMONITORS) { 1475 SETERRQ(1,0,"Too many monitors set"); 1476 } 1477 1478 snes->monitor[snes->numbermonitors] = func; 1479 snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 1480 PetscFunctionReturn(0); 1481 } 1482 1483 #undef __FUNC__ 1484 #define __FUNC__ "SNESSetConvergenceTest" 1485 /*@C 1486 SNESSetConvergenceTest - Sets the function that is to be used 1487 to test for convergence of the nonlinear iterative solution. 1488 1489 Input Parameters: 1490 . snes - the SNES context 1491 . func - routine to test for convergence 1492 . cctx - [optional] context for private data for the convergence routine 1493 (may be PETSC_NULL) 1494 1495 Calling sequence of func: 1496 int func (SNES snes,double xnorm,double gnorm, 1497 double f,void *cctx) 1498 1499 $ snes - the SNES context 1500 $ cctx - [optional] convergence context 1501 $ xnorm - 2-norm of current iterate 1502 $ 1503 $ SNES_NONLINEAR_EQUATIONS methods: 1504 $ gnorm - 2-norm of current step 1505 $ f - 2-norm of function 1506 $ 1507 $ SNES_UNCONSTRAINED_MINIMIZATION methods: 1508 $ gnorm - 2-norm of current gradient 1509 $ f - function value 1510 1511 .keywords: SNES, nonlinear, set, convergence, test 1512 1513 .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(), 1514 SNESConverged_UM_LS(), SNESConverged_UM_TR() 1515 @*/ 1516 int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx) 1517 { 1518 PetscFunctionBegin; 1519 (snes)->converged = func; 1520 (snes)->cnvP = cctx; 1521 PetscFunctionReturn(0); 1522 } 1523 1524 #undef __FUNC__ 1525 #define __FUNC__ "SNESSetConvergenceHistory" 1526 /*@ 1527 SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1528 1529 Input Parameters: 1530 . snes - iterative context obtained from SNESCreate() 1531 . a - array to hold history 1532 . na - size of a 1533 1534 Notes: 1535 If set, this array will contain the function norms (for 1536 SNES_NONLINEAR_EQUATIONS methods) or gradient norms 1537 (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed 1538 at each step. 1539 1540 This routine is useful, e.g., when running a code for purposes 1541 of accurate performance monitoring, when no I/O should be done 1542 during the section of code that is being timed. 1543 1544 .keywords: SNES, set, convergence, history 1545 @*/ 1546 int SNESSetConvergenceHistory(SNES snes, double *a, int na) 1547 { 1548 PetscFunctionBegin; 1549 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1550 if (na) PetscValidScalarPointer(a); 1551 snes->conv_hist = a; 1552 snes->conv_hist_size = na; 1553 PetscFunctionReturn(0); 1554 } 1555 1556 #undef __FUNC__ 1557 #define __FUNC__ "SNESScaleStep_Private" 1558 /* 1559 SNESScaleStep_Private - Scales a step so that its length is less than the 1560 positive parameter delta. 1561 1562 Input Parameters: 1563 . snes - the SNES context 1564 . y - approximate solution of linear system 1565 . fnorm - 2-norm of current function 1566 . delta - trust region size 1567 1568 Output Parameters: 1569 . gpnorm - predicted function norm at the new point, assuming local 1570 linearization. The value is zero if the step lies within the trust 1571 region, and exceeds zero otherwise. 1572 . ynorm - 2-norm of the step 1573 1574 Note: 1575 For non-trust region methods such as SNES_EQ_LS, the parameter delta 1576 is set to be the maximum allowable step size. 1577 1578 .keywords: SNES, nonlinear, scale, step 1579 */ 1580 int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta, 1581 double *gpnorm,double *ynorm) 1582 { 1583 double norm; 1584 Scalar cnorm; 1585 int ierr; 1586 1587 PetscFunctionBegin; 1588 ierr = VecNorm(y,NORM_2, &norm );CHKERRQ(ierr); 1589 if (norm > *delta) { 1590 norm = *delta/norm; 1591 *gpnorm = (1.0 - norm)*(*fnorm); 1592 cnorm = norm; 1593 VecScale( &cnorm, y ); 1594 *ynorm = *delta; 1595 } else { 1596 *gpnorm = 0.0; 1597 *ynorm = norm; 1598 } 1599 PetscFunctionReturn(0); 1600 } 1601 1602 #undef __FUNC__ 1603 #define __FUNC__ "SNESSolve" 1604 /*@ 1605 SNESSolve - Solves a nonlinear system. Call SNESSolve after calling 1606 SNESCreate() and optional routines of the form SNESSetXXX(). 1607 1608 Input Parameter: 1609 . snes - the SNES context 1610 . x - the solution vector 1611 1612 Output Parameter: 1613 its - number of iterations until termination 1614 1615 Note: 1616 The user should initialize the vector, x, with the initial guess 1617 for the nonlinear solve prior to calling SNESSolve. In particular, 1618 to employ an initial guess of zero, the user should explicitly set 1619 this vector to zero by calling VecSet(). 1620 1621 .keywords: SNES, nonlinear, solve 1622 1623 .seealso: SNESCreate(), SNESDestroy() 1624 @*/ 1625 int SNESSolve(SNES snes,Vec x,int *its) 1626 { 1627 int ierr, flg; 1628 1629 PetscFunctionBegin; 1630 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1631 PetscValidIntPointer(its); 1632 if (!snes->setup_called) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);} 1633 else {snes->vec_sol = snes->vec_sol_always = x;} 1634 PLogEventBegin(SNES_Solve,snes,0,0,0); 1635 snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0; 1636 ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr); 1637 PLogEventEnd(SNES_Solve,snes,0,0,0); 1638 ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr); 1639 if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); } 1640 PetscFunctionReturn(0); 1641 } 1642 1643 /* --------- Internal routines for SNES Package --------- */ 1644 static NRList *__SNESList = 0; 1645 1646 #undef __FUNC__ 1647 #define __FUNC__ "SNESSetType" 1648 /*@ 1649 SNESSetType - Sets the method for the nonlinear solver. 1650 1651 Input Parameters: 1652 . snes - the SNES context 1653 . method - a known method 1654 1655 Options Database Command: 1656 $ -snes_type <method> 1657 $ Use -help for a list of available methods 1658 $ (for instance, ls or tr) 1659 1660 Notes: 1661 See "petsc/include/snes.h" for available methods (for instance) 1662 $ Systems of nonlinear equations: 1663 $ SNES_EQ_LS - Newton's method with line search 1664 $ SNES_EQ_TR - Newton's method with trust region 1665 $ Unconstrained minimization: 1666 $ SNES_UM_TR - Newton's method with trust region 1667 $ SNES_UM_LS - Newton's method with line search 1668 1669 Normally, it is best to use the SNESSetFromOptions() command and then 1670 set the SNES solver type from the options database rather than by using 1671 this routine. Using the options database provides the user with 1672 maximum flexibility in evaluating the many nonlinear solvers. 1673 The SNESSetType() routine is provided for those situations where it 1674 is necessary to set the nonlinear solver independently of the command 1675 line or options database. This might be the case, for example, when 1676 the choice of solver changes during the execution of the program, 1677 and the user's application is taking responsibility for choosing the 1678 appropriate method. In other words, this routine is for the advanced user. 1679 1680 .keywords: SNES, set, method 1681 @*/ 1682 int SNESSetType(SNES snes,SNESType method) 1683 { 1684 int ierr; 1685 int (*r)(SNES); 1686 1687 PetscFunctionBegin; 1688 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1689 if (snes->type == (int) method) PetscFunctionReturn(0); 1690 1691 /* Get the function pointers for the iterative method requested */ 1692 if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(); CHKERRQ(ierr);} 1693 if (!__SNESList) {SETERRQ(1,0,"Could not get methods");} 1694 r = (int (*)(SNES))NRFindRoutine( __SNESList, (int)method, (char *)0 ); 1695 if (!r) {SETERRQ(1,0,"Unknown method");} 1696 if (snes->data) PetscFree(snes->data); 1697 snes->set_method_called = 1; 1698 ierr = (*r)(snes);CHKERRQ(ierr); 1699 PetscFunctionReturn(0); 1700 } 1701 1702 /* --------------------------------------------------------------------- */ 1703 #undef __FUNC__ 1704 #define __FUNC__ "SNESRegister" 1705 /*@C 1706 SNESRegister - Adds the method to the nonlinear solver package, given 1707 a function pointer and a nonlinear solver name of the type SNESType. 1708 1709 Input Parameters: 1710 . name - either a predefined name such as SNES_EQ_LS, or SNES_NEW 1711 to indicate a new user-defined solver 1712 . sname - corresponding string for name 1713 . create - routine to create method context 1714 1715 Output Parameter: 1716 . oname - type associated with this new solver 1717 1718 Notes: 1719 Multiple user-defined nonlinear solvers can be added by calling 1720 SNESRegister() with the input parameter "name" set to be SNES_NEW; 1721 each call will return a unique solver type in the output 1722 parameter "oname". 1723 1724 .keywords: SNES, nonlinear, register 1725 1726 .seealso: SNESRegisterAll(), SNESRegisterDestroy() 1727 @*/ 1728 int SNESRegister(SNESType name,SNESType *oname, char *sname, int (*create)(SNES)) 1729 { 1730 int ierr; 1731 static int numberregistered = 0; 1732 1733 PetscFunctionBegin; 1734 if (name == SNES_NEW) name = (SNESType) ((int) SNES_NEW + numberregistered++); 1735 1736 if (oname) *oname = name; 1737 if (!__SNESList) {ierr = NRCreate(&__SNESList); CHKERRQ(ierr);} 1738 NRRegister( __SNESList, (int) name, sname, (int (*)(void*))create ); 1739 PetscFunctionReturn(0); 1740 } 1741 1742 /* --------------------------------------------------------------------- */ 1743 #undef __FUNC__ 1744 #define __FUNC__ "SNESRegisterDestroy" 1745 /*@C 1746 SNESRegisterDestroy - Frees the list of nonlinear solvers that were 1747 registered by SNESRegister(). 1748 1749 .keywords: SNES, nonlinear, register, destroy 1750 1751 .seealso: SNESRegisterAll(), SNESRegisterAll() 1752 @*/ 1753 int SNESRegisterDestroy() 1754 { 1755 PetscFunctionBegin; 1756 if (__SNESList) { 1757 NRDestroy( __SNESList ); 1758 __SNESList = 0; 1759 } 1760 SNESRegisterAllCalled = 0; 1761 PetscFunctionReturn(0); 1762 } 1763 1764 #undef __FUNC__ 1765 #define __FUNC__ "SNESGetTypeFromOptions_Private" 1766 /* 1767 SNESGetTypeFromOptions_Private - Sets the selected method from the 1768 options database. 1769 1770 Input Parameter: 1771 . ctx - the SNES context 1772 1773 Output Parameter: 1774 . method - solver method 1775 . flg - indicate if method found 1776 1777 Options Database Key: 1778 $ -snes_type method 1779 */ 1780 int SNESGetTypeFromOptions_Private(SNES ctx,SNESType *method,int *flg) 1781 { 1782 int ierr; 1783 char sbuf[50]; 1784 1785 PetscFunctionBegin; 1786 ierr = OptionsGetString(ctx->prefix,"-snes_type",sbuf,50,flg);CHKERRQ(ierr); 1787 if (*flg) { 1788 if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);} 1789 *method = (SNESType)NRFindID( __SNESList, sbuf ); 1790 if (*method == (SNESType) -1) SETERRQ(1,1,"Invalid SNES type"); 1791 } 1792 PetscFunctionReturn(0); 1793 } 1794 1795 #undef __FUNC__ 1796 #define __FUNC__ "SNESGetType" 1797 /*@C 1798 SNESGetType - Gets the SNES method type and name (as a string). 1799 1800 Input Parameter: 1801 . snes - nonlinear solver context 1802 1803 Output Parameter: 1804 . method - SNES method (or use PETSC_NULL) 1805 . name - name of SNES method (or use PETSC_NULL) 1806 1807 .keywords: SNES, nonlinear, get, method, name 1808 @*/ 1809 int SNESGetType(SNES snes, SNESType *method,char **name) 1810 { 1811 int ierr; 1812 1813 PetscFunctionBegin; 1814 if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);} 1815 if (method) *method = (SNESType) snes->type; 1816 if (name) *name = NRFindName( __SNESList, (int) snes->type ); 1817 PetscFunctionReturn(0); 1818 } 1819 1820 #undef __FUNC__ 1821 #define __FUNC__ "SNESPrintTypes_Private" 1822 /* 1823 SNESPrintTypes_Private - Prints the SNES methods available from the 1824 options database. 1825 1826 Input Parameters: 1827 . comm - communicator (usually PETSC_COMM_WORLD) 1828 . prefix - prefix (usually "-") 1829 . name - the options database name (by default "snes_type") 1830 */ 1831 int SNESPrintTypes_Private(MPI_Comm comm,char* prefix,char *name) 1832 { 1833 FuncList *entry; 1834 1835 PetscFunctionBegin; 1836 if (!__SNESList) {SNESRegisterAll();} 1837 entry = __SNESList->head; 1838 PetscPrintf(comm," %s%s (one of)",prefix,name); 1839 while (entry) { 1840 PetscPrintf(comm," %s",entry->name); 1841 entry = entry->next; 1842 } 1843 PetscPrintf(comm,"\n"); 1844 PetscFunctionReturn(0); 1845 } 1846 1847 #undef __FUNC__ 1848 #define __FUNC__ "SNESGetSolution" 1849 /*@C 1850 SNESGetSolution - Returns the vector where the approximate solution is 1851 stored. 1852 1853 Input Parameter: 1854 . snes - the SNES context 1855 1856 Output Parameter: 1857 . x - the solution 1858 1859 .keywords: SNES, nonlinear, get, solution 1860 1861 .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate() 1862 @*/ 1863 int SNESGetSolution(SNES snes,Vec *x) 1864 { 1865 PetscFunctionBegin; 1866 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1867 *x = snes->vec_sol_always; 1868 PetscFunctionReturn(0); 1869 } 1870 1871 #undef __FUNC__ 1872 #define __FUNC__ "SNESGetSolutionUpdate" 1873 /*@C 1874 SNESGetSolutionUpdate - Returns the vector where the solution update is 1875 stored. 1876 1877 Input Parameter: 1878 . snes - the SNES context 1879 1880 Output Parameter: 1881 . x - the solution update 1882 1883 Notes: 1884 This vector is implementation dependent. 1885 1886 .keywords: SNES, nonlinear, get, solution, update 1887 1888 .seealso: SNESGetSolution(), SNESGetFunction 1889 @*/ 1890 int SNESGetSolutionUpdate(SNES snes,Vec *x) 1891 { 1892 PetscFunctionBegin; 1893 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1894 *x = snes->vec_sol_update_always; 1895 PetscFunctionReturn(0); 1896 } 1897 1898 #undef __FUNC__ 1899 #define __FUNC__ "SNESGetFunction" 1900 /*@C 1901 SNESGetFunction - Returns the vector where the function is stored. 1902 1903 Input Parameter: 1904 . snes - the SNES context 1905 1906 Output Parameter: 1907 . r - the function 1908 1909 Notes: 1910 SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only 1911 Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 1912 SNESGetMinimizationFunction() and SNESGetGradient(); 1913 1914 .keywords: SNES, nonlinear, get, function 1915 1916 .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(), 1917 SNESGetGradient() 1918 @*/ 1919 int SNESGetFunction(SNES snes,Vec *r) 1920 { 1921 PetscFunctionBegin; 1922 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1923 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 1924 *r = snes->vec_func_always; 1925 PetscFunctionReturn(0); 1926 } 1927 1928 #undef __FUNC__ 1929 #define __FUNC__ "SNESGetGradient" 1930 /*@C 1931 SNESGetGradient - Returns the vector where the gradient is stored. 1932 1933 Input Parameter: 1934 . snes - the SNES context 1935 1936 Output Parameter: 1937 . r - the gradient 1938 1939 Notes: 1940 SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods 1941 only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 1942 SNESGetFunction(). 1943 1944 .keywords: SNES, nonlinear, get, gradient 1945 1946 .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction() 1947 @*/ 1948 int SNESGetGradient(SNES snes,Vec *r) 1949 { 1950 PetscFunctionBegin; 1951 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1952 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 1953 SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 1954 } 1955 *r = snes->vec_func_always; 1956 PetscFunctionReturn(0); 1957 } 1958 1959 #undef __FUNC__ 1960 #define __FUNC__ "SNESGetMinimizationFunction" 1961 /*@ 1962 SNESGetMinimizationFunction - Returns the scalar function value for 1963 unconstrained minimization problems. 1964 1965 Input Parameter: 1966 . snes - the SNES context 1967 1968 Output Parameter: 1969 . r - the function 1970 1971 Notes: 1972 SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 1973 methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 1974 SNESGetFunction(). 1975 1976 .keywords: SNES, nonlinear, get, function 1977 1978 .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction() 1979 @*/ 1980 int SNESGetMinimizationFunction(SNES snes,double *r) 1981 { 1982 PetscFunctionBegin; 1983 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1984 PetscValidScalarPointer(r); 1985 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 1986 SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 1987 } 1988 *r = snes->fc; 1989 PetscFunctionReturn(0); 1990 } 1991 1992 #undef __FUNC__ 1993 #define __FUNC__ "SNESSetOptionsPrefix" 1994 /*@C 1995 SNESSetOptionsPrefix - Sets the prefix used for searching for all 1996 SNES options in the database. 1997 1998 Input Parameter: 1999 . snes - the SNES context 2000 . prefix - the prefix to prepend to all option names 2001 2002 Notes: 2003 A hyphen (-) must NOT be given at the beginning of the prefix name. 2004 The first character of all runtime options is AUTOMATICALLY the 2005 hyphen. 2006 2007 .keywords: SNES, set, options, prefix, database 2008 2009 .seealso: SNESSetFromOptions() 2010 @*/ 2011 int SNESSetOptionsPrefix(SNES snes,char *prefix) 2012 { 2013 int ierr; 2014 2015 PetscFunctionBegin; 2016 PetscValidHeaderSpecific(snes,SNES_COOKIE); 2017 ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 2018 ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 2019 PetscFunctionReturn(0); 2020 } 2021 2022 #undef __FUNC__ 2023 #define __FUNC__ "SNESAppendOptionsPrefix" 2024 /*@C 2025 SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 2026 SNES options in the database. 2027 2028 Input Parameter: 2029 . snes - the SNES context 2030 . prefix - the prefix to prepend to all option names 2031 2032 Notes: 2033 A hyphen (-) must NOT be given at the beginning of the prefix name. 2034 The first character of all runtime options is AUTOMATICALLY the 2035 hyphen. 2036 2037 .keywords: SNES, append, options, prefix, database 2038 2039 .seealso: SNESGetOptionsPrefix() 2040 @*/ 2041 int SNESAppendOptionsPrefix(SNES snes,char *prefix) 2042 { 2043 int ierr; 2044 2045 PetscFunctionBegin; 2046 PetscValidHeaderSpecific(snes,SNES_COOKIE); 2047 ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 2048 ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 2049 PetscFunctionReturn(0); 2050 } 2051 2052 #undef __FUNC__ 2053 #define __FUNC__ "SNESGetOptionsPrefix" 2054 /*@ 2055 SNESGetOptionsPrefix - Sets the prefix used for searching for all 2056 SNES options in the database. 2057 2058 Input Parameter: 2059 . snes - the SNES context 2060 2061 Output Parameter: 2062 . prefix - pointer to the prefix string used 2063 2064 .keywords: SNES, get, options, prefix, database 2065 2066 .seealso: SNESAppendOptionsPrefix() 2067 @*/ 2068 int SNESGetOptionsPrefix(SNES snes,char **prefix) 2069 { 2070 int ierr; 2071 2072 PetscFunctionBegin; 2073 PetscValidHeaderSpecific(snes,SNES_COOKIE); 2074 ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 2075 PetscFunctionReturn(0); 2076 } 2077 2078 2079 2080 2081 2082 2083