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