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