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