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