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