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