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