1 /*$Id: snes.c,v 1.233 2001/08/06 21:17:07 bsmith Exp balay $*/ 2 3 #include "src/snes/snesimpl.h" /*I "petscsnes.h" I*/ 4 5 PetscTruth SNESRegisterAllCalled = PETSC_FALSE; 6 PetscFList SNESList = 0; 7 8 #undef __FUNCT__ 9 #define __FUNCT__ "SNESGetProblemType" 10 /*@C 11 SNESGetProblemType -Indicates if SNES is solving a nonlinear system or a minimization 12 13 Not Collective 14 15 Input Parameter: 16 . SNES - the SNES context 17 18 Output Parameter: 19 . type - SNES_NONLINEAR_EQUATIONS (for systems of nonlinear equations) 20 or SNES_UNCONSTRAINED_MINIMIZATION (for unconstrained minimization) 21 22 Level: intermediate 23 24 .keywords: SNES, problem type 25 26 .seealso: SNESCreate() 27 @*/ 28 int SNESGetProblemType(SNES snes,SNESProblemType *type) 29 { 30 PetscFunctionBegin; 31 PetscValidHeaderSpecific(snes,SNES_COOKIE); 32 *type = snes->method_class; 33 PetscFunctionReturn(0); 34 } 35 36 #undef __FUNCT__ 37 #define __FUNCT__ "SNESView" 38 /*@C 39 SNESView - Prints the SNES data structure. 40 41 Collective on SNES 42 43 Input Parameters: 44 + SNES - the SNES context 45 - viewer - visualization context 46 47 Options Database Key: 48 . -snes_view - Calls SNESView() at end of SNESSolve() 49 50 Notes: 51 The available visualization contexts include 52 + PETSC_VIEWER_STDOUT_SELF - standard output (default) 53 - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 54 output where only the first processor opens 55 the file. All other processors send their 56 data to the first processor to print. 57 58 The user can open an alternative visualization context with 59 PetscViewerASCIIOpen() - output to a specified file. 60 61 Level: beginner 62 63 .keywords: SNES, view 64 65 .seealso: PetscViewerASCIIOpen() 66 @*/ 67 int SNESView(SNES snes,PetscViewer viewer) 68 { 69 SNES_KSP_EW_ConvCtx *kctx; 70 int ierr; 71 SLES sles; 72 char *type; 73 PetscTruth isascii,isstring; 74 75 PetscFunctionBegin; 76 PetscValidHeaderSpecific(snes,SNES_COOKIE); 77 if (!viewer) viewer = PETSC_VIEWER_STDOUT_(snes->comm); 78 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE); 79 PetscCheckSameComm(snes,viewer); 80 81 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 82 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 83 if (isascii) { 84 if (snes->prefix) { 85 ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:(%s)\n",snes->prefix);CHKERRQ(ierr); 86 } else { 87 ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:\n");CHKERRQ(ierr); 88 } 89 ierr = SNESGetType(snes,&type);CHKERRQ(ierr); 90 if (type) { 91 ierr = PetscViewerASCIIPrintf(viewer," type: %s\n",type);CHKERRQ(ierr); 92 } else { 93 ierr = PetscViewerASCIIPrintf(viewer," type: not set yet\n");CHKERRQ(ierr); 94 } 95 if (snes->view) { 96 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 97 ierr = (*snes->view)(snes,viewer);CHKERRQ(ierr); 98 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 99 } 100 ierr = PetscViewerASCIIPrintf(viewer," maximum iterations=%d, maximum function evaluations=%d\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr); 101 ierr = PetscViewerASCIIPrintf(viewer," tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n", 102 snes->rtol,snes->atol,snes->trunctol,snes->xtol);CHKERRQ(ierr); 103 ierr = PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%d\n",snes->linear_its);CHKERRQ(ierr); 104 ierr = PetscViewerASCIIPrintf(viewer," total number of function evaluations=%d\n",snes->nfuncs);CHKERRQ(ierr); 105 if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 106 ierr = PetscViewerASCIIPrintf(viewer," min function tolerance=%g\n",snes->fmin);CHKERRQ(ierr); 107 } 108 if (snes->ksp_ewconv) { 109 kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 110 if (kctx) { 111 ierr = PetscViewerASCIIPrintf(viewer," Eisenstat-Walker computation of KSP relative tolerance (version %d)\n",kctx->version);CHKERRQ(ierr); 112 ierr = PetscViewerASCIIPrintf(viewer," rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr); 113 ierr = PetscViewerASCIIPrintf(viewer," gamma=%g, alpha=%g, alpha2=%g\n",kctx->gamma,kctx->alpha,kctx->alpha2);CHKERRQ(ierr); 114 } 115 } 116 } else if (isstring) { 117 ierr = SNESGetType(snes,&type);CHKERRQ(ierr); 118 ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr); 119 } 120 ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr); 121 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 122 ierr = SLESView(sles,viewer);CHKERRQ(ierr); 123 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 124 PetscFunctionReturn(0); 125 } 126 127 #undef __FUNCT__ 128 #define __FUNCT__ "SNESSetFromOptions" 129 /*@ 130 SNESSetFromOptions - Sets various SNES and SLES parameters from user options. 131 132 Collective on SNES 133 134 Input Parameter: 135 . snes - the SNES context 136 137 Options Database Keys: 138 + -snes_type <type> - ls, tr, umls, umtr, test 139 . -snes_stol - convergence tolerance in terms of the norm 140 of the change in the solution between steps 141 . -snes_atol <atol> - absolute tolerance of residual norm 142 . -snes_rtol <rtol> - relative decrease in tolerance norm from initial 143 . -snes_max_it <max_it> - maximum number of iterations 144 . -snes_max_funcs <max_funcs> - maximum number of function evaluations 145 . -snes_trtol <trtol> - trust region tolerance 146 . -snes_no_convergence_test - skip convergence test in nonlinear or minimization 147 solver; hence iterations will continue until max_it 148 or some other criterion is reached. Saves expense 149 of convergence test 150 . -snes_monitor - prints residual norm at each iteration 151 . -snes_vecmonitor - plots solution at each iteration 152 . -snes_vecmonitor_update - plots update to solution at each iteration 153 . -snes_xmonitor - plots residual norm at each iteration 154 . -snes_fd - use finite differences to compute Jacobian; very slow, only for testing 155 - -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration 156 157 Options Database for Eisenstat-Walker method: 158 + -snes_ksp_eq_conv - use Eisenstat-Walker method for determining linear system convergence 159 . -snes_ksp_eq_version ver - version of Eisenstat-Walker method 160 . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0 161 . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax 162 . -snes_ksp_ew_gamma <gamma> - Sets gamma 163 . -snes_ksp_ew_alpha <alpha> - Sets alpha 164 . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2 165 - -snes_ksp_ew_threshold <threshold> - Sets threshold 166 167 Notes: 168 To see all options, run your program with the -help option or consult 169 the users manual. 170 171 Level: beginner 172 173 .keywords: SNES, nonlinear, set, options, database 174 175 .seealso: SNESSetOptionsPrefix() 176 @*/ 177 int SNESSetFromOptions(SNES snes) 178 { 179 SLES sles; 180 SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 181 PetscTruth flg; 182 int ierr; 183 char *deft,type[256]; 184 185 PetscFunctionBegin; 186 PetscValidHeaderSpecific(snes,SNES_COOKIE); 187 188 ierr = PetscOptionsBegin(snes->comm,snes->prefix,"Nonlinear solver (SNES) options","SNES");CHKERRQ(ierr); 189 if (snes->type_name) { 190 deft = snes->type_name; 191 } else { 192 if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 193 deft = SNESEQLS; 194 } else { 195 deft = SNESUMTR; 196 } 197 } 198 199 if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 200 ierr = PetscOptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr); 201 if (flg) { 202 ierr = SNESSetType(snes,type);CHKERRQ(ierr); 203 } else if (!snes->type_name) { 204 ierr = SNESSetType(snes,deft);CHKERRQ(ierr); 205 } 206 207 ierr = PetscOptionsReal("-snes_stol","Stop if step length less then","SNESSetTolerances",snes->xtol,&snes->xtol,0);CHKERRQ(ierr); 208 ierr = PetscOptionsReal("-snes_atol","Stop if function norm less then","SNESSetTolerances",snes->atol,&snes->atol,0);CHKERRQ(ierr); 209 210 ierr = PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less then","SNESSetTolerances",snes->rtol,&snes->rtol,0);CHKERRQ(ierr); 211 ierr = PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,PETSC_NULL);CHKERRQ(ierr); 212 ierr = PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr); 213 ierr = PetscOptionsReal("-snes_fmin","Minimization function tolerance","SNESSetMinimizationFunctionTolerance",snes->fmin,&snes->fmin,0);CHKERRQ(ierr); 214 215 ierr = PetscOptionsName("-snes_ksp_ew_conv","Use Eisentat-Walker linear system convergence test","SNES_KSP_SetParametersEW",&snes->ksp_ewconv);CHKERRQ(ierr); 216 217 ierr = PetscOptionsInt("-snes_ksp_ew_version","Version 1 or 2","SNES_KSP_SetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr); 218 ierr = PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNES_KSP_SetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr); 219 ierr = PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNES_KSP_SetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr); 220 ierr = PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNES_KSP_SetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr); 221 ierr = PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNES_KSP_SetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr); 222 ierr = PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNES_KSP_SetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr); 223 ierr = PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNES_KSP_SetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr); 224 225 ierr = PetscOptionsName("-snes_no_convergence_test","Don't test for convergence","None",&flg);CHKERRQ(ierr); 226 if (flg) {snes->converged = 0;} 227 ierr = PetscOptionsName("-snes_cancelmonitors","Remove all monitors","SNESClearMonitor",&flg);CHKERRQ(ierr); 228 if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);} 229 ierr = PetscOptionsName("-snes_monitor","Monitor norm of function","SNESDefaultMonitor",&flg);CHKERRQ(ierr); 230 if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0,0);CHKERRQ(ierr);} 231 ierr = PetscOptionsName("-snes_ratiomonitor","Monitor norm of function","SNESSetRatioMonitor",&flg);CHKERRQ(ierr); 232 if (flg) {ierr = SNESSetRatioMonitor(snes);CHKERRQ(ierr);} 233 ierr = PetscOptionsName("-snes_smonitor","Monitor norm of function (fewer digits)","SNESDefaultSMonitor",&flg);CHKERRQ(ierr); 234 if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0,0);CHKERRQ(ierr);} 235 ierr = PetscOptionsName("-snes_vecmonitor","Plot solution at each iteration","SNESVecViewMonitor",&flg);CHKERRQ(ierr); 236 if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewMonitor,0,0);CHKERRQ(ierr);} 237 ierr = PetscOptionsName("-snes_vecmonitor_update","Plot correction at each iteration","SNESVecViewUpdateMonitor",&flg);CHKERRQ(ierr); 238 if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewUpdateMonitor,0,0);CHKERRQ(ierr);} 239 ierr = PetscOptionsName("-snes_xmonitor","Plot function norm at each iteration","SNESLGMonitor",&flg);CHKERRQ(ierr); 240 if (flg) {ierr = SNESSetMonitor(snes,SNESLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);} 241 242 ierr = PetscOptionsName("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",&flg);CHKERRQ(ierr); 243 if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) { 244 ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr); 245 PetscLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Jacobian matrix\n"); 246 } else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 247 ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeHessian,snes->funP);CHKERRQ(ierr); 248 PetscLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Hessian matrix\n"); 249 } 250 251 if (snes->setfromoptions) { 252 ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr); 253 } 254 255 ierr = PetscOptionsEnd();CHKERRQ(ierr); 256 257 ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr); 258 ierr = SLESSetFromOptions(sles);CHKERRQ(ierr); 259 260 PetscFunctionReturn(0); 261 } 262 263 264 #undef __FUNCT__ 265 #define __FUNCT__ "SNESSetApplicationContext" 266 /*@ 267 SNESSetApplicationContext - Sets the optional user-defined context for 268 the nonlinear solvers. 269 270 Collective on SNES 271 272 Input Parameters: 273 + snes - the SNES context 274 - usrP - optional user context 275 276 Level: intermediate 277 278 .keywords: SNES, nonlinear, set, application, context 279 280 .seealso: SNESGetApplicationContext() 281 @*/ 282 int SNESSetApplicationContext(SNES snes,void *usrP) 283 { 284 PetscFunctionBegin; 285 PetscValidHeaderSpecific(snes,SNES_COOKIE); 286 snes->user = usrP; 287 PetscFunctionReturn(0); 288 } 289 290 #undef __FUNCT__ 291 #define __FUNCT__ "SNESGetApplicationContext" 292 /*@C 293 SNESGetApplicationContext - Gets the user-defined context for the 294 nonlinear solvers. 295 296 Not Collective 297 298 Input Parameter: 299 . snes - SNES context 300 301 Output Parameter: 302 . usrP - user context 303 304 Level: intermediate 305 306 .keywords: SNES, nonlinear, get, application, context 307 308 .seealso: SNESSetApplicationContext() 309 @*/ 310 int SNESGetApplicationContext(SNES snes,void **usrP) 311 { 312 PetscFunctionBegin; 313 PetscValidHeaderSpecific(snes,SNES_COOKIE); 314 *usrP = snes->user; 315 PetscFunctionReturn(0); 316 } 317 318 #undef __FUNCT__ 319 #define __FUNCT__ "SNESGetIterationNumber" 320 /*@ 321 SNESGetIterationNumber - Gets the number of nonlinear iterations completed 322 at this time. 323 324 Not Collective 325 326 Input Parameter: 327 . snes - SNES context 328 329 Output Parameter: 330 . iter - iteration number 331 332 Notes: 333 For example, during the computation of iteration 2 this would return 1. 334 335 This is useful for using lagged Jacobians (where one does not recompute the 336 Jacobian at each SNES iteration). For example, the code 337 .vb 338 ierr = SNESGetIterationNumber(snes,&it); 339 if (!(it % 2)) { 340 [compute Jacobian here] 341 } 342 .ve 343 can be used in your ComputeJacobian() function to cause the Jacobian to be 344 recomputed every second SNES iteration. 345 346 Level: intermediate 347 348 .keywords: SNES, nonlinear, get, iteration, number 349 @*/ 350 int SNESGetIterationNumber(SNES snes,int* iter) 351 { 352 PetscFunctionBegin; 353 PetscValidHeaderSpecific(snes,SNES_COOKIE); 354 PetscValidIntPointer(iter); 355 *iter = snes->iter; 356 PetscFunctionReturn(0); 357 } 358 359 #undef __FUNCT__ 360 #define __FUNCT__ "SNESGetFunctionNorm" 361 /*@ 362 SNESGetFunctionNorm - Gets the norm of the current function that was set 363 with SNESSSetFunction(). 364 365 Collective on SNES 366 367 Input Parameter: 368 . snes - SNES context 369 370 Output Parameter: 371 . fnorm - 2-norm of function 372 373 Note: 374 SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only. 375 A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is 376 SNESGetGradientNorm(). 377 378 Level: intermediate 379 380 .keywords: SNES, nonlinear, get, function, norm 381 382 .seealso: SNESGetFunction() 383 @*/ 384 int SNESGetFunctionNorm(SNES snes,PetscScalar *fnorm) 385 { 386 PetscFunctionBegin; 387 PetscValidHeaderSpecific(snes,SNES_COOKIE); 388 PetscValidScalarPointer(fnorm); 389 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 390 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"For SNES_NONLINEAR_EQUATIONS only"); 391 } 392 *fnorm = snes->norm; 393 PetscFunctionReturn(0); 394 } 395 396 #undef __FUNCT__ 397 #define __FUNCT__ "SNESGetGradientNorm" 398 /*@ 399 SNESGetGradientNorm - Gets the norm of the current gradient that was set 400 with SNESSSetGradient(). 401 402 Collective on SNES 403 404 Input Parameter: 405 . snes - SNES context 406 407 Output Parameter: 408 . fnorm - 2-norm of gradient 409 410 Note: 411 SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION 412 methods only. A related routine for SNES_NONLINEAR_EQUATIONS methods 413 is SNESGetFunctionNorm(). 414 415 Level: intermediate 416 417 .keywords: SNES, nonlinear, get, gradient, norm 418 419 .seelso: SNESSetGradient() 420 @*/ 421 int SNESGetGradientNorm(SNES snes,PetscScalar *gnorm) 422 { 423 PetscFunctionBegin; 424 PetscValidHeaderSpecific(snes,SNES_COOKIE); 425 PetscValidScalarPointer(gnorm); 426 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 427 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 428 } 429 *gnorm = snes->norm; 430 PetscFunctionReturn(0); 431 } 432 433 #undef __FUNCT__ 434 #define __FUNCT__ "SNESGetNumberUnsuccessfulSteps" 435 /*@ 436 SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps 437 attempted by the nonlinear solver. 438 439 Not Collective 440 441 Input Parameter: 442 . snes - SNES context 443 444 Output Parameter: 445 . nfails - number of unsuccessful steps attempted 446 447 Notes: 448 This counter is reset to zero for each successive call to SNESSolve(). 449 450 Level: intermediate 451 452 .keywords: SNES, nonlinear, get, number, unsuccessful, steps 453 @*/ 454 int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails) 455 { 456 PetscFunctionBegin; 457 PetscValidHeaderSpecific(snes,SNES_COOKIE); 458 PetscValidIntPointer(nfails); 459 *nfails = snes->nfailures; 460 PetscFunctionReturn(0); 461 } 462 463 #undef __FUNCT__ 464 #define __FUNCT__ "SNESGetNumberLinearIterations" 465 /*@ 466 SNESGetNumberLinearIterations - Gets the total number of linear iterations 467 used by the nonlinear solver. 468 469 Not Collective 470 471 Input Parameter: 472 . snes - SNES context 473 474 Output Parameter: 475 . lits - number of linear iterations 476 477 Notes: 478 This counter is reset to zero for each successive call to SNESSolve(). 479 480 Level: intermediate 481 482 .keywords: SNES, nonlinear, get, number, linear, iterations 483 @*/ 484 int SNESGetNumberLinearIterations(SNES snes,int* lits) 485 { 486 PetscFunctionBegin; 487 PetscValidHeaderSpecific(snes,SNES_COOKIE); 488 PetscValidIntPointer(lits); 489 *lits = snes->linear_its; 490 PetscFunctionReturn(0); 491 } 492 493 #undef __FUNCT__ 494 #define __FUNCT__ "SNESGetSLES" 495 /*@C 496 SNESGetSLES - Returns the SLES context for a SNES solver. 497 498 Not Collective, but if SNES object is parallel, then SLES object is parallel 499 500 Input Parameter: 501 . snes - the SNES context 502 503 Output Parameter: 504 . sles - the SLES context 505 506 Notes: 507 The user can then directly manipulate the SLES context to set various 508 options, etc. Likewise, the user can then extract and manipulate the 509 KSP and PC contexts as well. 510 511 Level: beginner 512 513 .keywords: SNES, nonlinear, get, SLES, context 514 515 .seealso: SLESGetPC(), SLESGetKSP() 516 @*/ 517 int SNESGetSLES(SNES snes,SLES *sles) 518 { 519 PetscFunctionBegin; 520 PetscValidHeaderSpecific(snes,SNES_COOKIE); 521 *sles = snes->sles; 522 PetscFunctionReturn(0); 523 } 524 525 #undef __FUNCT__ 526 #define __FUNCT__ "SNESPublish_Petsc" 527 static int SNESPublish_Petsc(PetscObject obj) 528 { 529 #if defined(PETSC_HAVE_AMS) 530 SNES v = (SNES) obj; 531 int ierr; 532 #endif 533 534 PetscFunctionBegin; 535 536 #if defined(PETSC_HAVE_AMS) 537 /* if it is already published then return */ 538 if (v->amem >=0) PetscFunctionReturn(0); 539 540 ierr = PetscObjectPublishBaseBegin(obj);CHKERRQ(ierr); 541 ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Iteration",&v->iter,1,AMS_INT,AMS_READ, 542 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr); 543 ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Residual",&v->norm,1,AMS_DOUBLE,AMS_READ, 544 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr); 545 ierr = PetscObjectPublishBaseEnd(obj);CHKERRQ(ierr); 546 #endif 547 PetscFunctionReturn(0); 548 } 549 550 /* -----------------------------------------------------------*/ 551 #undef __FUNCT__ 552 #define __FUNCT__ "SNESCreate" 553 /*@C 554 SNESCreate - Creates a nonlinear solver context. 555 556 Collective on MPI_Comm 557 558 Input Parameters: 559 + comm - MPI communicator 560 - type - type of method, either 561 SNES_NONLINEAR_EQUATIONS (for systems of nonlinear equations) 562 or SNES_UNCONSTRAINED_MINIMIZATION (for unconstrained minimization) 563 564 Output Parameter: 565 . outsnes - the new SNES context 566 567 Options Database Keys: 568 + -snes_mf - Activates default matrix-free Jacobian-vector products, 569 and no preconditioning matrix 570 . -snes_mf_operator - Activates default matrix-free Jacobian-vector 571 products, and a user-provided preconditioning matrix 572 as set by SNESSetJacobian() 573 - -snes_fd - Uses (slow!) finite differences to compute Jacobian 574 575 Level: beginner 576 577 .keywords: SNES, nonlinear, create, context 578 579 .seealso: SNESSolve(), SNESDestroy(), SNESProblemType, SNES 580 @*/ 581 int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes) 582 { 583 int ierr; 584 SNES snes; 585 SNES_KSP_EW_ConvCtx *kctx; 586 587 PetscFunctionBegin; 588 *outsnes = 0; 589 if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS){ 590 SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"incorrect method type"); 591 } 592 PetscHeaderCreate(snes,_p_SNES,int,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView); 593 PetscLogObjectCreate(snes); 594 snes->bops->publish = SNESPublish_Petsc; 595 snes->max_its = 50; 596 snes->max_funcs = 10000; 597 snes->norm = 0.0; 598 if (type == SNES_UNCONSTRAINED_MINIMIZATION) { 599 snes->rtol = 1.e-8; 600 snes->ttol = 0.0; 601 snes->atol = 1.e-10; 602 } else { 603 snes->rtol = 1.e-8; 604 snes->ttol = 0.0; 605 snes->atol = 1.e-50; 606 } 607 snes->xtol = 1.e-8; 608 snes->trunctol = 1.e-12; /* no longer used */ 609 snes->nfuncs = 0; 610 snes->nfailures = 0; 611 snes->linear_its = 0; 612 snes->numbermonitors = 0; 613 snes->data = 0; 614 snes->view = 0; 615 snes->computeumfunction = 0; 616 snes->umfunP = 0; 617 snes->fc = 0; 618 snes->deltatol = 1.e-12; 619 snes->fmin = -1.e30; 620 snes->method_class = type; 621 snes->set_method_called = 0; 622 snes->setupcalled = 0; 623 snes->ksp_ewconv = PETSC_FALSE; 624 snes->vwork = 0; 625 snes->nwork = 0; 626 snes->conv_hist_len = 0; 627 snes->conv_hist_max = 0; 628 snes->conv_hist = PETSC_NULL; 629 snes->conv_hist_its = PETSC_NULL; 630 snes->conv_hist_reset = PETSC_TRUE; 631 snes->reason = SNES_CONVERGED_ITERATING; 632 633 /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 634 ierr = PetscNew(SNES_KSP_EW_ConvCtx,&kctx);CHKERRQ(ierr); 635 PetscLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx)); 636 snes->kspconvctx = (void*)kctx; 637 kctx->version = 2; 638 kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 639 this was too large for some test cases */ 640 kctx->rtol_last = 0; 641 kctx->rtol_max = .9; 642 kctx->gamma = 1.0; 643 kctx->alpha2 = .5*(1.0 + sqrt(5.0)); 644 kctx->alpha = kctx->alpha2; 645 kctx->threshold = .1; 646 kctx->lresid_last = 0; 647 kctx->norm_last = 0; 648 649 ierr = SLESCreate(comm,&snes->sles);CHKERRQ(ierr); 650 PetscLogObjectParent(snes,snes->sles) 651 652 *outsnes = snes; 653 ierr = PetscPublishAll(snes);CHKERRQ(ierr); 654 PetscFunctionReturn(0); 655 } 656 657 /* --------------------------------------------------------------- */ 658 #undef __FUNCT__ 659 #define __FUNCT__ "SNESSetFunction" 660 /*@C 661 SNESSetFunction - Sets the function evaluation routine and function 662 vector for use by the SNES routines in solving systems of nonlinear 663 equations. 664 665 Collective on SNES 666 667 Input Parameters: 668 + snes - the SNES context 669 . func - function evaluation routine 670 . r - vector to store function value 671 - ctx - [optional] user-defined context for private data for the 672 function evaluation routine (may be PETSC_NULL) 673 674 Calling sequence of func: 675 $ func (SNES snes,Vec x,Vec f,void *ctx); 676 677 . f - function vector 678 - ctx - optional user-defined function context 679 680 Notes: 681 The Newton-like methods typically solve linear systems of the form 682 $ f'(x) x = -f(x), 683 where f'(x) denotes the Jacobian matrix and f(x) is the function. 684 685 SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 686 Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 687 SNESSetMinimizationFunction() and SNESSetGradient(); 688 689 Level: beginner 690 691 .keywords: SNES, nonlinear, set, function 692 693 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 694 @*/ 695 int SNESSetFunction(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx) 696 { 697 PetscFunctionBegin; 698 PetscValidHeaderSpecific(snes,SNES_COOKIE); 699 PetscValidHeaderSpecific(r,VEC_COOKIE); 700 PetscCheckSameComm(snes,r); 701 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 702 SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only"); 703 } 704 705 snes->computefunction = func; 706 snes->vec_func = snes->vec_func_always = r; 707 snes->funP = ctx; 708 PetscFunctionReturn(0); 709 } 710 711 #undef __FUNCT__ 712 #define __FUNCT__ "SNESComputeFunction" 713 /*@ 714 SNESComputeFunction - Calls the function that has been set with 715 SNESSetFunction(). 716 717 Collective on SNES 718 719 Input Parameters: 720 + snes - the SNES context 721 - x - input vector 722 723 Output Parameter: 724 . y - function vector, as set by SNESSetFunction() 725 726 Notes: 727 SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 728 Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 729 SNESComputeMinimizationFunction() and SNESComputeGradient(); 730 731 SNESComputeFunction() is typically used within nonlinear solvers 732 implementations, so most users would not generally call this routine 733 themselves. 734 735 Level: developer 736 737 .keywords: SNES, nonlinear, compute, function 738 739 .seealso: SNESSetFunction(), SNESGetFunction() 740 @*/ 741 int SNESComputeFunction(SNES snes,Vec x,Vec y) 742 { 743 int ierr; 744 745 PetscFunctionBegin; 746 PetscValidHeaderSpecific(snes,SNES_COOKIE); 747 PetscValidHeaderSpecific(x,VEC_COOKIE); 748 PetscValidHeaderSpecific(y,VEC_COOKIE); 749 PetscCheckSameComm(snes,x); 750 PetscCheckSameComm(snes,y); 751 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 752 SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only"); 753 } 754 755 ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 756 PetscStackPush("SNES user function"); 757 ierr = (*snes->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr); 758 PetscStackPop; 759 snes->nfuncs++; 760 ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 761 PetscFunctionReturn(0); 762 } 763 764 #undef __FUNCT__ 765 #define __FUNCT__ "SNESSetMinimizationFunction" 766 /*@C 767 SNESSetMinimizationFunction - Sets the function evaluation routine for 768 unconstrained minimization. 769 770 Collective on SNES 771 772 Input Parameters: 773 + snes - the SNES context 774 . func - function evaluation routine 775 - ctx - [optional] user-defined context for private data for the 776 function evaluation routine (may be PETSC_NULL) 777 778 Calling sequence of func: 779 $ func (SNES snes,Vec x,PetscReal *f,void *ctx); 780 781 + x - input vector 782 . f - function 783 - ctx - [optional] user-defined function context 784 785 Level: beginner 786 787 Notes: 788 SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 789 methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 790 SNESSetFunction(). 791 792 .keywords: SNES, nonlinear, set, minimization, function 793 794 .seealso: SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(), 795 SNESSetHessian(), SNESSetGradient() 796 @*/ 797 int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,PetscReal*,void*),void *ctx) 798 { 799 PetscFunctionBegin; 800 PetscValidHeaderSpecific(snes,SNES_COOKIE); 801 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 802 SETERRQ(PETSC_ERR_ARG_WRONG,"Only for SNES_UNCONSTRAINED_MINIMIZATION"); 803 } 804 snes->computeumfunction = func; 805 snes->umfunP = ctx; 806 PetscFunctionReturn(0); 807 } 808 809 #undef __FUNCT__ 810 #define __FUNCT__ "SNESComputeMinimizationFunction" 811 /*@ 812 SNESComputeMinimizationFunction - Computes the function that has been 813 set with SNESSetMinimizationFunction(). 814 815 Collective on SNES 816 817 Input Parameters: 818 + snes - the SNES context 819 - x - input vector 820 821 Output Parameter: 822 . y - function value 823 824 Notes: 825 SNESComputeMinimizationFunction() is valid only for 826 SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 827 SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 828 829 SNESComputeMinimizationFunction() is typically used within minimization 830 implementations, so most users would not generally call this routine 831 themselves. 832 833 Level: developer 834 835 .keywords: SNES, nonlinear, compute, minimization, function 836 837 .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(), 838 SNESComputeGradient(), SNESComputeHessian() 839 @*/ 840 int SNESComputeMinimizationFunction(SNES snes,Vec x,PetscReal *y) 841 { 842 int ierr; 843 844 PetscFunctionBegin; 845 PetscValidHeaderSpecific(snes,SNES_COOKIE); 846 PetscValidHeaderSpecific(x,VEC_COOKIE); 847 PetscCheckSameComm(snes,x); 848 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 849 SETERRQ(PETSC_ERR_ARG_WRONG,"Only for SNES_UNCONSTRAINED_MINIMIZATION"); 850 } 851 852 ierr = PetscLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0);CHKERRQ(ierr); 853 PetscStackPush("SNES user minimzation function"); 854 ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP);CHKERRQ(ierr); 855 PetscStackPop; 856 snes->nfuncs++; 857 ierr = PetscLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0);CHKERRQ(ierr); 858 PetscFunctionReturn(0); 859 } 860 861 #undef __FUNCT__ 862 #define __FUNCT__ "SNESSetGradient" 863 /*@C 864 SNESSetGradient - Sets the gradient evaluation routine and gradient 865 vector for use by the SNES routines. 866 867 Collective on SNES 868 869 Input Parameters: 870 + snes - the SNES context 871 . func - function evaluation routine 872 . ctx - optional user-defined context for private data for the 873 gradient evaluation routine (may be PETSC_NULL) 874 - r - vector to store gradient value 875 876 Calling sequence of func: 877 $ func (SNES, Vec x, Vec g, void *ctx); 878 879 + x - input vector 880 . g - gradient vector 881 - ctx - optional user-defined gradient context 882 883 Notes: 884 SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 885 methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 886 SNESSetFunction(). 887 888 Level: beginner 889 890 .keywords: SNES, nonlinear, set, function 891 892 .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(), 893 SNESSetMinimizationFunction(), 894 @*/ 895 int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx) 896 { 897 PetscFunctionBegin; 898 PetscValidHeaderSpecific(snes,SNES_COOKIE); 899 PetscValidHeaderSpecific(r,VEC_COOKIE); 900 PetscCheckSameComm(snes,r); 901 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 902 SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 903 } 904 snes->computefunction = func; 905 snes->vec_func = snes->vec_func_always = r; 906 snes->funP = ctx; 907 PetscFunctionReturn(0); 908 } 909 910 #undef __FUNCT__ 911 #define __FUNCT__ "SNESComputeGradient" 912 /*@ 913 SNESComputeGradient - Computes the gradient that has been set with 914 SNESSetGradient(). 915 916 Collective on SNES 917 918 Input Parameters: 919 + snes - the SNES context 920 - x - input vector 921 922 Output Parameter: 923 . y - gradient vector 924 925 Notes: 926 SNESComputeGradient() is valid only for 927 SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 928 SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 929 930 SNESComputeGradient() is typically used within minimization 931 implementations, so most users would not generally call this routine 932 themselves. 933 934 Level: developer 935 936 .keywords: SNES, nonlinear, compute, gradient 937 938 .seealso: SNESSetGradient(), SNESGetGradient(), 939 SNESComputeMinimizationFunction(), SNESComputeHessian() 940 @*/ 941 int SNESComputeGradient(SNES snes,Vec x,Vec y) 942 { 943 int ierr; 944 945 PetscFunctionBegin; 946 PetscValidHeaderSpecific(snes,SNES_COOKIE); 947 PetscValidHeaderSpecific(x,VEC_COOKIE); 948 PetscValidHeaderSpecific(y,VEC_COOKIE); 949 PetscCheckSameComm(snes,x); 950 PetscCheckSameComm(snes,y); 951 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 952 SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 953 } 954 ierr = PetscLogEventBegin(SNES_GradientEval,snes,x,y,0);CHKERRQ(ierr); 955 PetscStackPush("SNES user gradient function"); 956 ierr = (*snes->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr); 957 PetscStackPop; 958 ierr = PetscLogEventEnd(SNES_GradientEval,snes,x,y,0);CHKERRQ(ierr); 959 PetscFunctionReturn(0); 960 } 961 962 #undef __FUNCT__ 963 #define __FUNCT__ "SNESComputeJacobian" 964 /*@ 965 SNESComputeJacobian - Computes the Jacobian matrix that has been 966 set with SNESSetJacobian(). 967 968 Collective on SNES and Mat 969 970 Input Parameters: 971 + snes - the SNES context 972 - x - input vector 973 974 Output Parameters: 975 + A - Jacobian matrix 976 . B - optional preconditioning matrix 977 - flag - flag indicating matrix structure 978 979 Notes: 980 Most users should not need to explicitly call this routine, as it 981 is used internally within the nonlinear solvers. 982 983 See SLESSetOperators() for important information about setting the 984 flag parameter. 985 986 SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS 987 methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION 988 methods is SNESComputeHessian(). 989 990 Level: developer 991 992 .keywords: SNES, compute, Jacobian, matrix 993 994 .seealso: SNESSetJacobian(), SLESSetOperators() 995 @*/ 996 int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 997 { 998 int ierr; 999 1000 PetscFunctionBegin; 1001 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1002 PetscValidHeaderSpecific(X,VEC_COOKIE); 1003 PetscCheckSameComm(snes,X); 1004 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 1005 SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only"); 1006 } 1007 if (!snes->computejacobian) PetscFunctionReturn(0); 1008 ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 1009 *flg = DIFFERENT_NONZERO_PATTERN; 1010 PetscStackPush("SNES user Jacobian function"); 1011 ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr); 1012 PetscStackPop; 1013 ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 1014 /* make sure user returned a correct Jacobian and preconditioner */ 1015 PetscValidHeaderSpecific(*A,MAT_COOKIE); 1016 PetscValidHeaderSpecific(*B,MAT_COOKIE); 1017 PetscFunctionReturn(0); 1018 } 1019 1020 #undef __FUNCT__ 1021 #define __FUNCT__ "SNESComputeHessian" 1022 /*@ 1023 SNESComputeHessian - Computes the Hessian matrix that has been 1024 set with SNESSetHessian(). 1025 1026 Collective on SNES and Mat 1027 1028 Input Parameters: 1029 + snes - the SNES context 1030 - x - input vector 1031 1032 Output Parameters: 1033 + A - Hessian matrix 1034 . B - optional preconditioning matrix 1035 - flag - flag indicating matrix structure 1036 1037 Notes: 1038 Most users should not need to explicitly call this routine, as it 1039 is used internally within the nonlinear solvers. 1040 1041 See SLESSetOperators() for important information about setting the 1042 flag parameter. 1043 1044 SNESComputeHessian() is valid only for 1045 SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 1046 SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian(). 1047 1048 SNESComputeHessian() is typically used within minimization 1049 implementations, so most users would not generally call this routine 1050 themselves. 1051 1052 Level: developer 1053 1054 .keywords: SNES, compute, Hessian, matrix 1055 1056 .seealso: SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(), 1057 SNESComputeMinimizationFunction() 1058 @*/ 1059 int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag) 1060 { 1061 int ierr; 1062 1063 PetscFunctionBegin; 1064 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1065 PetscValidHeaderSpecific(x,VEC_COOKIE); 1066 PetscCheckSameComm(snes,x); 1067 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 1068 SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 1069 } 1070 if (!snes->computejacobian) PetscFunctionReturn(0); 1071 ierr = PetscLogEventBegin(SNES_HessianEval,snes,x,*A,*B);CHKERRQ(ierr); 1072 *flag = DIFFERENT_NONZERO_PATTERN; 1073 PetscStackPush("SNES user Hessian function"); 1074 ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP);CHKERRQ(ierr); 1075 PetscStackPop; 1076 ierr = PetscLogEventEnd(SNES_HessianEval,snes,x,*A,*B);CHKERRQ(ierr); 1077 /* make sure user returned a correct Jacobian and preconditioner */ 1078 PetscValidHeaderSpecific(*A,MAT_COOKIE); 1079 PetscValidHeaderSpecific(*B,MAT_COOKIE); 1080 PetscFunctionReturn(0); 1081 } 1082 1083 #undef __FUNCT__ 1084 #define __FUNCT__ "SNESSetJacobian" 1085 /*@C 1086 SNESSetJacobian - Sets the function to compute Jacobian as well as the 1087 location to store the matrix. 1088 1089 Collective on SNES and Mat 1090 1091 Input Parameters: 1092 + snes - the SNES context 1093 . A - Jacobian matrix 1094 . B - preconditioner matrix (usually same as the Jacobian) 1095 . func - Jacobian evaluation routine 1096 - ctx - [optional] user-defined context for private data for the 1097 Jacobian evaluation routine (may be PETSC_NULL) 1098 1099 Calling sequence of func: 1100 $ func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 1101 1102 + x - input vector 1103 . A - Jacobian matrix 1104 . B - preconditioner matrix, usually the same as A 1105 . flag - flag indicating information about the preconditioner matrix 1106 structure (same as flag in SLESSetOperators()) 1107 - ctx - [optional] user-defined Jacobian context 1108 1109 Notes: 1110 See SLESSetOperators() for important information about setting the flag 1111 output parameter in the routine func(). Be sure to read this information! 1112 1113 The routine func() takes Mat * as the matrix arguments rather than Mat. 1114 This allows the Jacobian evaluation routine to replace A and/or B with a 1115 completely new new matrix structure (not just different matrix elements) 1116 when appropriate, for instance, if the nonzero structure is changing 1117 throughout the global iterations. 1118 1119 Level: beginner 1120 1121 .keywords: SNES, nonlinear, set, Jacobian, matrix 1122 1123 .seealso: SLESSetOperators(), SNESSetFunction() 1124 @*/ 1125 int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 1126 { 1127 int ierr; 1128 1129 PetscFunctionBegin; 1130 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1131 if (A) PetscValidHeaderSpecific(A,MAT_COOKIE); 1132 if (B) PetscValidHeaderSpecific(B,MAT_COOKIE); 1133 if (A) PetscCheckSameComm(snes,A); 1134 if (B) PetscCheckSameComm(snes,B); 1135 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 1136 SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only"); 1137 } 1138 1139 if (func) snes->computejacobian = func; 1140 if (ctx) snes->jacP = ctx; 1141 if (A) { 1142 if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 1143 snes->jacobian = A; 1144 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 1145 } 1146 if (B) { 1147 if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 1148 snes->jacobian_pre = B; 1149 ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 1150 } 1151 PetscFunctionReturn(0); 1152 } 1153 1154 #undef __FUNCT__ 1155 #define __FUNCT__ "SNESGetJacobian" 1156 /*@C 1157 SNESGetJacobian - Returns the Jacobian matrix and optionally the user 1158 provided context for evaluating the Jacobian. 1159 1160 Not Collective, but Mat object will be parallel if SNES object is 1161 1162 Input Parameter: 1163 . snes - the nonlinear solver context 1164 1165 Output Parameters: 1166 + A - location to stash Jacobian matrix (or PETSC_NULL) 1167 . B - location to stash preconditioner matrix (or PETSC_NULL) 1168 . ctx - location to stash Jacobian ctx (or PETSC_NULL) 1169 - func - location to put Jacobian function (or PETSC_NULL) 1170 1171 Level: advanced 1172 1173 .seealso: SNESSetJacobian(), SNESComputeJacobian() 1174 @*/ 1175 int SNESGetJacobian(SNES snes,Mat *A,Mat *B,void **ctx,int (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*)) 1176 { 1177 PetscFunctionBegin; 1178 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1179 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 1180 SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only"); 1181 } 1182 if (A) *A = snes->jacobian; 1183 if (B) *B = snes->jacobian_pre; 1184 if (ctx) *ctx = snes->jacP; 1185 if (func) *func = snes->computejacobian; 1186 PetscFunctionReturn(0); 1187 } 1188 1189 #undef __FUNCT__ 1190 #define __FUNCT__ "SNESSetHessian" 1191 /*@C 1192 SNESSetHessian - Sets the function to compute Hessian as well as the 1193 location to store the matrix. 1194 1195 Collective on SNES and Mat 1196 1197 Input Parameters: 1198 + snes - the SNES context 1199 . A - Hessian matrix 1200 . B - preconditioner matrix (usually same as the Hessian) 1201 . func - Jacobian evaluation routine 1202 - ctx - [optional] user-defined context for private data for the 1203 Hessian evaluation routine (may be PETSC_NULL) 1204 1205 Calling sequence of func: 1206 $ func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 1207 1208 + x - input vector 1209 . A - Hessian matrix 1210 . B - preconditioner matrix, usually the same as A 1211 . flag - flag indicating information about the preconditioner matrix 1212 structure (same as flag in SLESSetOperators()) 1213 - ctx - [optional] user-defined Hessian context 1214 1215 Notes: 1216 See SLESSetOperators() for important information about setting the flag 1217 output parameter in the routine func(). Be sure to read this information! 1218 1219 The function func() takes Mat * as the matrix arguments rather than Mat. 1220 This allows the Hessian evaluation routine to replace A and/or B with a 1221 completely new new matrix structure (not just different matrix elements) 1222 when appropriate, for instance, if the nonzero structure is changing 1223 throughout the global iterations. 1224 1225 Level: beginner 1226 1227 .keywords: SNES, nonlinear, set, Hessian, matrix 1228 1229 .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators() 1230 @*/ 1231 int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 1232 { 1233 int ierr; 1234 1235 PetscFunctionBegin; 1236 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1237 PetscValidHeaderSpecific(A,MAT_COOKIE); 1238 PetscValidHeaderSpecific(B,MAT_COOKIE); 1239 PetscCheckSameComm(snes,A); 1240 PetscCheckSameComm(snes,B); 1241 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 1242 SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 1243 } 1244 if (func) snes->computejacobian = func; 1245 if (ctx) snes->jacP = ctx; 1246 if (A) { 1247 if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 1248 snes->jacobian = A; 1249 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 1250 } 1251 if (B) { 1252 if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 1253 snes->jacobian_pre = B; 1254 ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 1255 } 1256 PetscFunctionReturn(0); 1257 } 1258 1259 #undef __FUNCT__ 1260 #define __FUNCT__ "SNESGetHessian" 1261 /*@ 1262 SNESGetHessian - Returns the Hessian matrix and optionally the user 1263 provided context for evaluating the Hessian. 1264 1265 Not Collective, but Mat object is parallel if SNES object is parallel 1266 1267 Input Parameter: 1268 . snes - the nonlinear solver context 1269 1270 Output Parameters: 1271 + A - location to stash Hessian matrix (or PETSC_NULL) 1272 . B - location to stash preconditioner matrix (or PETSC_NULL) 1273 - ctx - location to stash Hessian ctx (or PETSC_NULL) 1274 1275 Level: advanced 1276 1277 .seealso: SNESSetHessian(), SNESComputeHessian() 1278 1279 .keywords: SNES, get, Hessian 1280 @*/ 1281 int SNESGetHessian(SNES snes,Mat *A,Mat *B,void **ctx) 1282 { 1283 PetscFunctionBegin; 1284 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1285 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION){ 1286 SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 1287 } 1288 if (A) *A = snes->jacobian; 1289 if (B) *B = snes->jacobian_pre; 1290 if (ctx) *ctx = snes->jacP; 1291 PetscFunctionReturn(0); 1292 } 1293 1294 /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 1295 1296 #undef __FUNCT__ 1297 #define __FUNCT__ "SNESSetUp" 1298 /*@ 1299 SNESSetUp - Sets up the internal data structures for the later use 1300 of a nonlinear solver. 1301 1302 Collective on SNES 1303 1304 Input Parameters: 1305 + snes - the SNES context 1306 - x - the solution vector 1307 1308 Notes: 1309 For basic use of the SNES solvers the user need not explicitly call 1310 SNESSetUp(), since these actions will automatically occur during 1311 the call to SNESSolve(). However, if one wishes to control this 1312 phase separately, SNESSetUp() should be called after SNESCreate() 1313 and optional routines of the form SNESSetXXX(), but before SNESSolve(). 1314 1315 Level: advanced 1316 1317 .keywords: SNES, nonlinear, setup 1318 1319 .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 1320 @*/ 1321 int SNESSetUp(SNES snes,Vec x) 1322 { 1323 int ierr; 1324 PetscTruth flg; 1325 1326 PetscFunctionBegin; 1327 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1328 PetscValidHeaderSpecific(x,VEC_COOKIE); 1329 PetscCheckSameComm(snes,x); 1330 snes->vec_sol = snes->vec_sol_always = x; 1331 1332 ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr); 1333 /* 1334 This version replaces the user provided Jacobian matrix with a 1335 matrix-free version but still employs the user-provided preconditioner matrix 1336 */ 1337 if (flg) { 1338 Mat J; 1339 ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1340 ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 1341 PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator routines\n"); 1342 ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 1343 ierr = MatDestroy(J);CHKERRQ(ierr); 1344 } 1345 ierr = PetscOptionsHasName(snes->prefix,"-snes_mf",&flg);CHKERRQ(ierr); 1346 /* 1347 This version replaces both the user-provided Jacobian and the user- 1348 provided preconditioner matrix with the default matrix free version. 1349 */ 1350 if (flg) { 1351 Mat J; 1352 SLES sles; 1353 PC pc; 1354 1355 ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1356 ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 1357 PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator and preconditioner routines\n"); 1358 if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 1359 ierr = SNESSetJacobian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr); 1360 } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 1361 ierr = SNESSetHessian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr); 1362 } else { 1363 SETERRQ(PETSC_ERR_SUP,"Method class doesn't support matrix-free option"); 1364 } 1365 ierr = MatDestroy(J);CHKERRQ(ierr); 1366 1367 /* force no preconditioner */ 1368 ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr); 1369 ierr = SLESGetPC(sles,&pc);CHKERRQ(ierr); 1370 ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 1371 } 1372 1373 if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) { 1374 PetscTruth iseqtr; 1375 1376 if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 1377 if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 1378 if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option"); 1379 if (snes->vec_func == snes->vec_sol) { 1380 SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector"); 1381 } 1382 1383 /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 1384 ierr = PetscTypeCompare((PetscObject)snes,SNESEQTR,&iseqtr);CHKERRQ(ierr); 1385 if (snes->ksp_ewconv && !iseqtr) { 1386 SLES sles; KSP ksp; 1387 ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr); 1388 ierr = SLESGetKSP(sles,&ksp);CHKERRQ(ierr); 1389 ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,snes);CHKERRQ(ierr); 1390 } 1391 } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) { 1392 if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetGradient() first"); 1393 if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetGradient() first"); 1394 if (!snes->computeumfunction) { 1395 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetMinimizationFunction() first"); 1396 } 1397 if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetHessian()"); 1398 } else { 1399 SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Unknown method class"); 1400 } 1401 if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);} 1402 snes->setupcalled = 1; 1403 PetscFunctionReturn(0); 1404 } 1405 1406 #undef __FUNCT__ 1407 #define __FUNCT__ "SNESDestroy" 1408 /*@C 1409 SNESDestroy - Destroys the nonlinear solver context that was created 1410 with SNESCreate(). 1411 1412 Collective on SNES 1413 1414 Input Parameter: 1415 . snes - the SNES context 1416 1417 Level: beginner 1418 1419 .keywords: SNES, nonlinear, destroy 1420 1421 .seealso: SNESCreate(), SNESSolve() 1422 @*/ 1423 int SNESDestroy(SNES snes) 1424 { 1425 int i,ierr; 1426 1427 PetscFunctionBegin; 1428 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1429 if (--snes->refct > 0) PetscFunctionReturn(0); 1430 1431 /* if memory was published with AMS then destroy it */ 1432 ierr = PetscObjectDepublish(snes);CHKERRQ(ierr); 1433 1434 if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);} 1435 if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);} 1436 if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 1437 if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 1438 ierr = SLESDestroy(snes->sles);CHKERRQ(ierr); 1439 if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);} 1440 for (i=0; i<snes->numbermonitors; i++) { 1441 if (snes->monitordestroy[i]) { 1442 ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr); 1443 } 1444 } 1445 PetscLogObjectDestroy((PetscObject)snes); 1446 PetscHeaderDestroy((PetscObject)snes); 1447 PetscFunctionReturn(0); 1448 } 1449 1450 /* ----------- Routines to set solver parameters ---------- */ 1451 1452 #undef __FUNCT__ 1453 #define __FUNCT__ "SNESSetTolerances" 1454 /*@ 1455 SNESSetTolerances - Sets various parameters used in convergence tests. 1456 1457 Collective on SNES 1458 1459 Input Parameters: 1460 + snes - the SNES context 1461 . atol - absolute convergence tolerance 1462 . rtol - relative convergence tolerance 1463 . stol - convergence tolerance in terms of the norm 1464 of the change in the solution between steps 1465 . maxit - maximum number of iterations 1466 - maxf - maximum number of function evaluations 1467 1468 Options Database Keys: 1469 + -snes_atol <atol> - Sets atol 1470 . -snes_rtol <rtol> - Sets rtol 1471 . -snes_stol <stol> - Sets stol 1472 . -snes_max_it <maxit> - Sets maxit 1473 - -snes_max_funcs <maxf> - Sets maxf 1474 1475 Notes: 1476 The default maximum number of iterations is 50. 1477 The default maximum number of function evaluations is 1000. 1478 1479 Level: intermediate 1480 1481 .keywords: SNES, nonlinear, set, convergence, tolerances 1482 1483 .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance() 1484 @*/ 1485 int SNESSetTolerances(SNES snes,PetscReal atol,PetscReal rtol,PetscReal stol,int maxit,int maxf) 1486 { 1487 PetscFunctionBegin; 1488 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1489 if (atol != PETSC_DEFAULT) snes->atol = atol; 1490 if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1491 if (stol != PETSC_DEFAULT) snes->xtol = stol; 1492 if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1493 if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 1494 PetscFunctionReturn(0); 1495 } 1496 1497 #undef __FUNCT__ 1498 #define __FUNCT__ "SNESGetTolerances" 1499 /*@ 1500 SNESGetTolerances - Gets various parameters used in convergence tests. 1501 1502 Not Collective 1503 1504 Input Parameters: 1505 + snes - the SNES context 1506 . atol - absolute convergence tolerance 1507 . rtol - relative convergence tolerance 1508 . stol - convergence tolerance in terms of the norm 1509 of the change in the solution between steps 1510 . maxit - maximum number of iterations 1511 - maxf - maximum number of function evaluations 1512 1513 Notes: 1514 The user can specify PETSC_NULL for any parameter that is not needed. 1515 1516 Level: intermediate 1517 1518 .keywords: SNES, nonlinear, get, convergence, tolerances 1519 1520 .seealso: SNESSetTolerances() 1521 @*/ 1522 int SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,int *maxit,int *maxf) 1523 { 1524 PetscFunctionBegin; 1525 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1526 if (atol) *atol = snes->atol; 1527 if (rtol) *rtol = snes->rtol; 1528 if (stol) *stol = snes->xtol; 1529 if (maxit) *maxit = snes->max_its; 1530 if (maxf) *maxf = snes->max_funcs; 1531 PetscFunctionReturn(0); 1532 } 1533 1534 #undef __FUNCT__ 1535 #define __FUNCT__ "SNESSetTrustRegionTolerance" 1536 /*@ 1537 SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 1538 1539 Collective on SNES 1540 1541 Input Parameters: 1542 + snes - the SNES context 1543 - tol - tolerance 1544 1545 Options Database Key: 1546 . -snes_trtol <tol> - Sets tol 1547 1548 Level: intermediate 1549 1550 .keywords: SNES, nonlinear, set, trust region, tolerance 1551 1552 .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance() 1553 @*/ 1554 int SNESSetTrustRegionTolerance(SNES snes,PetscReal tol) 1555 { 1556 PetscFunctionBegin; 1557 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1558 snes->deltatol = tol; 1559 PetscFunctionReturn(0); 1560 } 1561 1562 #undef __FUNCT__ 1563 #define __FUNCT__ "SNESSetMinimizationFunctionTolerance" 1564 /*@ 1565 SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance 1566 for unconstrained minimization solvers. 1567 1568 Collective on SNES 1569 1570 Input Parameters: 1571 + snes - the SNES context 1572 - ftol - minimum function tolerance 1573 1574 Options Database Key: 1575 . -snes_fmin <ftol> - Sets ftol 1576 1577 Note: 1578 SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION 1579 methods only. 1580 1581 Level: intermediate 1582 1583 .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance 1584 1585 .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance() 1586 @*/ 1587 int SNESSetMinimizationFunctionTolerance(SNES snes,PetscReal ftol) 1588 { 1589 PetscFunctionBegin; 1590 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1591 snes->fmin = ftol; 1592 PetscFunctionReturn(0); 1593 } 1594 /* 1595 Duplicate the lg monitors for SNES from KSP; for some reason with 1596 dynamic libraries things don't work under Sun4 if we just use 1597 macros instead of functions 1598 */ 1599 #undef __FUNCT__ 1600 #define __FUNCT__ "SNESLGMonitor" 1601 int SNESLGMonitor(SNES snes,int it,PetscReal norm,void *ctx) 1602 { 1603 int ierr; 1604 1605 PetscFunctionBegin; 1606 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1607 ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr); 1608 PetscFunctionReturn(0); 1609 } 1610 1611 #undef __FUNCT__ 1612 #define __FUNCT__ "SNESLGMonitorCreate" 1613 int SNESLGMonitorCreate(char *host,char *label,int x,int y,int m,int n,PetscDrawLG *draw) 1614 { 1615 int ierr; 1616 1617 PetscFunctionBegin; 1618 ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr); 1619 PetscFunctionReturn(0); 1620 } 1621 1622 #undef __FUNCT__ 1623 #define __FUNCT__ "SNESLGMonitorDestroy" 1624 int SNESLGMonitorDestroy(PetscDrawLG draw) 1625 { 1626 int ierr; 1627 1628 PetscFunctionBegin; 1629 ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr); 1630 PetscFunctionReturn(0); 1631 } 1632 1633 /* ------------ Routines to set performance monitoring options ----------- */ 1634 1635 #undef __FUNCT__ 1636 #define __FUNCT__ "SNESSetMonitor" 1637 /*@C 1638 SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every 1639 iteration of the nonlinear solver to display the iteration's 1640 progress. 1641 1642 Collective on SNES 1643 1644 Input Parameters: 1645 + snes - the SNES context 1646 . func - monitoring routine 1647 . mctx - [optional] user-defined context for private data for the 1648 monitor routine (use PETSC_NULL if no context is desitre) 1649 - monitordestroy - [optional] routine that frees monitor context 1650 (may be PETSC_NULL) 1651 1652 Calling sequence of func: 1653 $ int func(SNES snes,int its, PetscReal norm,void *mctx) 1654 1655 + snes - the SNES context 1656 . its - iteration number 1657 . norm - 2-norm function value (may be estimated) 1658 for SNES_NONLINEAR_EQUATIONS methods 1659 . norm - 2-norm gradient value (may be estimated) 1660 for SNES_UNCONSTRAINED_MINIMIZATION methods 1661 - mctx - [optional] monitoring context 1662 1663 Options Database Keys: 1664 + -snes_monitor - sets SNESDefaultMonitor() 1665 . -snes_xmonitor - sets line graph monitor, 1666 uses SNESLGMonitorCreate() 1667 _ -snes_cancelmonitors - cancels all monitors that have 1668 been hardwired into a code by 1669 calls to SNESSetMonitor(), but 1670 does not cancel those set via 1671 the options database. 1672 1673 Notes: 1674 Several different monitoring routines may be set by calling 1675 SNESSetMonitor() multiple times; all will be called in the 1676 order in which they were set. 1677 1678 Level: intermediate 1679 1680 .keywords: SNES, nonlinear, set, monitor 1681 1682 .seealso: SNESDefaultMonitor(), SNESClearMonitor() 1683 @*/ 1684 int SNESSetMonitor(SNES snes,int (*func)(SNES,int,PetscReal,void*),void *mctx,int (*monitordestroy)(void *)) 1685 { 1686 PetscFunctionBegin; 1687 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1688 if (snes->numbermonitors >= MAXSNESMONITORS) { 1689 SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1690 } 1691 1692 snes->monitor[snes->numbermonitors] = func; 1693 snes->monitordestroy[snes->numbermonitors] = monitordestroy; 1694 snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 1695 PetscFunctionReturn(0); 1696 } 1697 1698 #undef __FUNCT__ 1699 #define __FUNCT__ "SNESClearMonitor" 1700 /*@C 1701 SNESClearMonitor - Clears all the monitor functions for a SNES object. 1702 1703 Collective on SNES 1704 1705 Input Parameters: 1706 . snes - the SNES context 1707 1708 Options Database: 1709 . -snes_cancelmonitors - cancels all monitors that have been hardwired 1710 into a code by calls to SNESSetMonitor(), but does not cancel those 1711 set via the options database 1712 1713 Notes: 1714 There is no way to clear one specific monitor from a SNES object. 1715 1716 Level: intermediate 1717 1718 .keywords: SNES, nonlinear, set, monitor 1719 1720 .seealso: SNESDefaultMonitor(), SNESSetMonitor() 1721 @*/ 1722 int SNESClearMonitor(SNES snes) 1723 { 1724 PetscFunctionBegin; 1725 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1726 snes->numbermonitors = 0; 1727 PetscFunctionReturn(0); 1728 } 1729 1730 #undef __FUNCT__ 1731 #define __FUNCT__ "SNESSetConvergenceTest" 1732 /*@C 1733 SNESSetConvergenceTest - Sets the function that is to be used 1734 to test for convergence of the nonlinear iterative solution. 1735 1736 Collective on SNES 1737 1738 Input Parameters: 1739 + snes - the SNES context 1740 . func - routine to test for convergence 1741 - cctx - [optional] context for private data for the convergence routine 1742 (may be PETSC_NULL) 1743 1744 Calling sequence of func: 1745 $ int func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx) 1746 1747 + snes - the SNES context 1748 . cctx - [optional] convergence context 1749 . reason - reason for convergence/divergence 1750 . xnorm - 2-norm of current iterate 1751 . gnorm - 2-norm of current step (SNES_NONLINEAR_EQUATIONS methods) 1752 . f - 2-norm of function (SNES_NONLINEAR_EQUATIONS methods) 1753 . gnorm - 2-norm of current gradient (SNES_UNCONSTRAINED_MINIMIZATION methods) 1754 - f - function value (SNES_UNCONSTRAINED_MINIMIZATION methods) 1755 1756 Level: advanced 1757 1758 .keywords: SNES, nonlinear, set, convergence, test 1759 1760 .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(), 1761 SNESConverged_UM_LS(), SNESConverged_UM_TR() 1762 @*/ 1763 int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx) 1764 { 1765 PetscFunctionBegin; 1766 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1767 (snes)->converged = func; 1768 (snes)->cnvP = cctx; 1769 PetscFunctionReturn(0); 1770 } 1771 1772 #undef __FUNCT__ 1773 #define __FUNCT__ "SNESGetConvergedReason" 1774 /*@C 1775 SNESGetConvergedReason - Gets the reason the SNES iteration was stopped. 1776 1777 Not Collective 1778 1779 Input Parameter: 1780 . snes - the SNES context 1781 1782 Output Parameter: 1783 . reason - negative value indicates diverged, positive value converged, see petscsnes.h or the 1784 manual pages for the individual convergence tests for complete lists 1785 1786 Level: intermediate 1787 1788 Notes: Can only be called after the call the SNESSolve() is complete. 1789 1790 .keywords: SNES, nonlinear, set, convergence, test 1791 1792 .seealso: SNESSetConvergenceTest(), SNESConverged_EQ_LS(), SNESConverged_EQ_TR(), 1793 SNESConverged_UM_LS(), SNESConverged_UM_TR(), SNESConvergedReason 1794 @*/ 1795 int SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason) 1796 { 1797 PetscFunctionBegin; 1798 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1799 *reason = snes->reason; 1800 PetscFunctionReturn(0); 1801 } 1802 1803 #undef __FUNCT__ 1804 #define __FUNCT__ "SNESSetConvergenceHistory" 1805 /*@ 1806 SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1807 1808 Collective on SNES 1809 1810 Input Parameters: 1811 + snes - iterative context obtained from SNESCreate() 1812 . a - array to hold history 1813 . its - integer array holds the number of linear iterations for each solve. 1814 . na - size of a and its 1815 - reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero, 1816 else it continues storing new values for new nonlinear solves after the old ones 1817 1818 Notes: 1819 If set, this array will contain the function norms (for 1820 SNES_NONLINEAR_EQUATIONS methods) or gradient norms 1821 (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed 1822 at each step. 1823 1824 This routine is useful, e.g., when running a code for purposes 1825 of accurate performance monitoring, when no I/O should be done 1826 during the section of code that is being timed. 1827 1828 Level: intermediate 1829 1830 .keywords: SNES, set, convergence, history 1831 1832 .seealso: SNESGetConvergenceHistory() 1833 1834 @*/ 1835 int SNESSetConvergenceHistory(SNES snes,PetscReal *a,int *its,int na,PetscTruth reset) 1836 { 1837 PetscFunctionBegin; 1838 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1839 if (na) PetscValidScalarPointer(a); 1840 snes->conv_hist = a; 1841 snes->conv_hist_its = its; 1842 snes->conv_hist_max = na; 1843 snes->conv_hist_reset = reset; 1844 PetscFunctionReturn(0); 1845 } 1846 1847 #undef __FUNCT__ 1848 #define __FUNCT__ "SNESGetConvergenceHistory" 1849 /*@C 1850 SNESGetConvergenceHistory - Gets the array used to hold the convergence history. 1851 1852 Collective on SNES 1853 1854 Input Parameter: 1855 . snes - iterative context obtained from SNESCreate() 1856 1857 Output Parameters: 1858 . a - array to hold history 1859 . its - integer array holds the number of linear iterations (or 1860 negative if not converged) for each solve. 1861 - na - size of a and its 1862 1863 Notes: 1864 The calling sequence for this routine in Fortran is 1865 $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr) 1866 1867 This routine is useful, e.g., when running a code for purposes 1868 of accurate performance monitoring, when no I/O should be done 1869 during the section of code that is being timed. 1870 1871 Level: intermediate 1872 1873 .keywords: SNES, get, convergence, history 1874 1875 .seealso: SNESSetConvergencHistory() 1876 1877 @*/ 1878 int SNESGetConvergenceHistory(SNES snes,PetscReal **a,int **its,int *na) 1879 { 1880 PetscFunctionBegin; 1881 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1882 if (a) *a = snes->conv_hist; 1883 if (its) *its = snes->conv_hist_its; 1884 if (na) *na = snes->conv_hist_len; 1885 PetscFunctionReturn(0); 1886 } 1887 1888 #undef __FUNCT__ 1889 #define __FUNCT__ "SNESScaleStep_Private" 1890 /* 1891 SNESScaleStep_Private - Scales a step so that its length is less than the 1892 positive parameter delta. 1893 1894 Input Parameters: 1895 + snes - the SNES context 1896 . y - approximate solution of linear system 1897 . fnorm - 2-norm of current function 1898 - delta - trust region size 1899 1900 Output Parameters: 1901 + gpnorm - predicted function norm at the new point, assuming local 1902 linearization. The value is zero if the step lies within the trust 1903 region, and exceeds zero otherwise. 1904 - ynorm - 2-norm of the step 1905 1906 Note: 1907 For non-trust region methods such as SNESEQLS, the parameter delta 1908 is set to be the maximum allowable step size. 1909 1910 .keywords: SNES, nonlinear, scale, step 1911 */ 1912 int SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta, 1913 PetscReal *gpnorm,PetscReal *ynorm) 1914 { 1915 PetscReal norm; 1916 PetscScalar cnorm; 1917 int ierr; 1918 1919 PetscFunctionBegin; 1920 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1921 PetscValidHeaderSpecific(y,VEC_COOKIE); 1922 PetscCheckSameComm(snes,y); 1923 1924 ierr = VecNorm(y,NORM_2,&norm);CHKERRQ(ierr); 1925 if (norm > *delta) { 1926 norm = *delta/norm; 1927 *gpnorm = (1.0 - norm)*(*fnorm); 1928 cnorm = norm; 1929 ierr = VecScale(&cnorm,y);CHKERRQ(ierr); 1930 *ynorm = *delta; 1931 } else { 1932 *gpnorm = 0.0; 1933 *ynorm = norm; 1934 } 1935 PetscFunctionReturn(0); 1936 } 1937 1938 #undef __FUNCT__ 1939 #define __FUNCT__ "SNESSolve" 1940 /*@ 1941 SNESSolve - Solves a nonlinear system. Call SNESSolve after calling 1942 SNESCreate() and optional routines of the form SNESSetXXX(). 1943 1944 Collective on SNES 1945 1946 Input Parameters: 1947 + snes - the SNES context 1948 - x - the solution vector 1949 1950 Output Parameter: 1951 . its - number of iterations until termination 1952 1953 Notes: 1954 The user should initialize the vector,x, with the initial guess 1955 for the nonlinear solve prior to calling SNESSolve. In particular, 1956 to employ an initial guess of zero, the user should explicitly set 1957 this vector to zero by calling VecSet(). 1958 1959 Level: beginner 1960 1961 .keywords: SNES, nonlinear, solve 1962 1963 .seealso: SNESCreate(), SNESDestroy() 1964 @*/ 1965 int SNESSolve(SNES snes,Vec x,int *its) 1966 { 1967 int ierr; 1968 PetscTruth flg; 1969 1970 PetscFunctionBegin; 1971 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1972 PetscValidHeaderSpecific(x,VEC_COOKIE); 1973 PetscCheckSameComm(snes,x); 1974 PetscValidIntPointer(its); 1975 if (!snes->solve) SETERRQ(1,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()"); 1976 1977 if (!snes->setupcalled) {ierr = SNESSetUp(snes,x);CHKERRQ(ierr);} 1978 else {snes->vec_sol = snes->vec_sol_always = x;} 1979 if (snes->conv_hist_reset == PETSC_TRUE) snes->conv_hist_len = 0; 1980 ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 1981 snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0; 1982 ierr = (*(snes)->solve)(snes,its);CHKERRQ(ierr); 1983 ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 1984 ierr = PetscOptionsHasName(snes->prefix,"-snes_view",&flg);CHKERRQ(ierr); 1985 if (flg && !PetscPreLoadingOn) { ierr = SNESView(snes,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); } 1986 PetscFunctionReturn(0); 1987 } 1988 1989 /* --------- Internal routines for SNES Package --------- */ 1990 1991 #undef __FUNCT__ 1992 #define __FUNCT__ "SNESSetType" 1993 /*@C 1994 SNESSetType - Sets the method for the nonlinear solver. 1995 1996 Collective on SNES 1997 1998 Input Parameters: 1999 + snes - the SNES context 2000 - type - a known method 2001 2002 Options Database Key: 2003 . -snes_type <type> - Sets the method; use -help for a list 2004 of available methods (for instance, ls or tr) 2005 2006 Notes: 2007 See "petsc/include/petscsnes.h" for available methods (for instance) 2008 + SNESEQLS - Newton's method with line search 2009 (systems of nonlinear equations) 2010 . SNESEQTR - Newton's method with trust region 2011 (systems of nonlinear equations) 2012 . SNESUMTR - Newton's method with trust region 2013 (unconstrained minimization) 2014 - SNESUMLS - Newton's method with line search 2015 (unconstrained minimization) 2016 2017 Normally, it is best to use the SNESSetFromOptions() command and then 2018 set the SNES solver type from the options database rather than by using 2019 this routine. Using the options database provides the user with 2020 maximum flexibility in evaluating the many nonlinear solvers. 2021 The SNESSetType() routine is provided for those situations where it 2022 is necessary to set the nonlinear solver independently of the command 2023 line or options database. This might be the case, for example, when 2024 the choice of solver changes during the execution of the program, 2025 and the user's application is taking responsibility for choosing the 2026 appropriate method. 2027 2028 Level: intermediate 2029 2030 .keywords: SNES, set, type 2031 2032 .seealso: SNESType, SNESCreate() 2033 2034 @*/ 2035 int SNESSetType(SNES snes,SNESType type) 2036 { 2037 int ierr,(*r)(SNES); 2038 PetscTruth match; 2039 2040 PetscFunctionBegin; 2041 PetscValidHeaderSpecific(snes,SNES_COOKIE); 2042 PetscValidCharPointer(type); 2043 2044 ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr); 2045 if (match) PetscFunctionReturn(0); 2046 2047 if (snes->setupcalled) { 2048 ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr); 2049 snes->data = 0; 2050 } 2051 2052 /* Get the function pointers for the iterative method requested */ 2053 if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 2054 2055 ierr = PetscFListFind(snes->comm,SNESList,type,(void (**)(void)) &r);CHKERRQ(ierr); 2056 2057 if (!r) SETERRQ1(1,"Unable to find requested SNES type %s",type); 2058 2059 if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);} 2060 snes->data = 0; 2061 ierr = (*r)(snes);CHKERRQ(ierr); 2062 2063 ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr); 2064 snes->set_method_called = 1; 2065 2066 PetscFunctionReturn(0); 2067 } 2068 2069 2070 /* --------------------------------------------------------------------- */ 2071 #undef __FUNCT__ 2072 #define __FUNCT__ "SNESRegisterDestroy" 2073 /*@C 2074 SNESRegisterDestroy - Frees the list of nonlinear solvers that were 2075 registered by SNESRegisterDynamic(). 2076 2077 Not Collective 2078 2079 Level: advanced 2080 2081 .keywords: SNES, nonlinear, register, destroy 2082 2083 .seealso: SNESRegisterAll(), SNESRegisterAll() 2084 @*/ 2085 int SNESRegisterDestroy(void) 2086 { 2087 int ierr; 2088 2089 PetscFunctionBegin; 2090 if (SNESList) { 2091 ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr); 2092 SNESList = 0; 2093 } 2094 SNESRegisterAllCalled = PETSC_FALSE; 2095 PetscFunctionReturn(0); 2096 } 2097 2098 #undef __FUNCT__ 2099 #define __FUNCT__ "SNESGetType" 2100 /*@C 2101 SNESGetType - Gets the SNES method type and name (as a string). 2102 2103 Not Collective 2104 2105 Input Parameter: 2106 . snes - nonlinear solver context 2107 2108 Output Parameter: 2109 . type - SNES method (a character string) 2110 2111 Level: intermediate 2112 2113 .keywords: SNES, nonlinear, get, type, name 2114 @*/ 2115 int SNESGetType(SNES snes,SNESType *type) 2116 { 2117 PetscFunctionBegin; 2118 PetscValidHeaderSpecific(snes,SNES_COOKIE); 2119 *type = snes->type_name; 2120 PetscFunctionReturn(0); 2121 } 2122 2123 #undef __FUNCT__ 2124 #define __FUNCT__ "SNESGetSolution" 2125 /*@C 2126 SNESGetSolution - Returns the vector where the approximate solution is 2127 stored. 2128 2129 Not Collective, but Vec is parallel if SNES is parallel 2130 2131 Input Parameter: 2132 . snes - the SNES context 2133 2134 Output Parameter: 2135 . x - the solution 2136 2137 Level: advanced 2138 2139 .keywords: SNES, nonlinear, get, solution 2140 2141 .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate() 2142 @*/ 2143 int SNESGetSolution(SNES snes,Vec *x) 2144 { 2145 PetscFunctionBegin; 2146 PetscValidHeaderSpecific(snes,SNES_COOKIE); 2147 *x = snes->vec_sol_always; 2148 PetscFunctionReturn(0); 2149 } 2150 2151 #undef __FUNCT__ 2152 #define __FUNCT__ "SNESGetSolutionUpdate" 2153 /*@C 2154 SNESGetSolutionUpdate - Returns the vector where the solution update is 2155 stored. 2156 2157 Not Collective, but Vec is parallel if SNES is parallel 2158 2159 Input Parameter: 2160 . snes - the SNES context 2161 2162 Output Parameter: 2163 . x - the solution update 2164 2165 Level: advanced 2166 2167 .keywords: SNES, nonlinear, get, solution, update 2168 2169 .seealso: SNESGetSolution(), SNESGetFunction 2170 @*/ 2171 int SNESGetSolutionUpdate(SNES snes,Vec *x) 2172 { 2173 PetscFunctionBegin; 2174 PetscValidHeaderSpecific(snes,SNES_COOKIE); 2175 *x = snes->vec_sol_update_always; 2176 PetscFunctionReturn(0); 2177 } 2178 2179 #undef __FUNCT__ 2180 #define __FUNCT__ "SNESGetFunction" 2181 /*@C 2182 SNESGetFunction - Returns the vector where the function is stored. 2183 2184 Not Collective, but Vec is parallel if SNES is parallel 2185 2186 Input Parameter: 2187 . snes - the SNES context 2188 2189 Output Parameter: 2190 + r - the function (or PETSC_NULL) 2191 . ctx - the function context (or PETSC_NULL) 2192 - func - the function (or PETSC_NULL) 2193 2194 Notes: 2195 SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only 2196 Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 2197 SNESGetMinimizationFunction() and SNESGetGradient(); 2198 2199 Level: advanced 2200 2201 .keywords: SNES, nonlinear, get, function 2202 2203 .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(), 2204 SNESGetGradient() 2205 2206 @*/ 2207 int SNESGetFunction(SNES snes,Vec *r,void **ctx,int (**func)(SNES,Vec,Vec,void*)) 2208 { 2209 PetscFunctionBegin; 2210 PetscValidHeaderSpecific(snes,SNES_COOKIE); 2211 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 2212 SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only"); 2213 } 2214 if (r) *r = snes->vec_func_always; 2215 if (ctx) *ctx = snes->funP; 2216 if (func) *func = snes->computefunction; 2217 PetscFunctionReturn(0); 2218 } 2219 2220 #undef __FUNCT__ 2221 #define __FUNCT__ "SNESGetGradient" 2222 /*@C 2223 SNESGetGradient - Returns the vector where the gradient is stored. 2224 2225 Not Collective, but Vec is parallel if SNES is parallel 2226 2227 Input Parameter: 2228 . snes - the SNES context 2229 2230 Output Parameter: 2231 + r - the gradient (or PETSC_NULL) 2232 - ctx - the gradient context (or PETSC_NULL) 2233 2234 Notes: 2235 SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods 2236 only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 2237 SNESGetFunction(). 2238 2239 Level: advanced 2240 2241 .keywords: SNES, nonlinear, get, gradient 2242 2243 .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction(), 2244 SNESSetGradient(), SNESSetFunction() 2245 2246 @*/ 2247 int SNESGetGradient(SNES snes,Vec *r,void **ctx) 2248 { 2249 PetscFunctionBegin; 2250 PetscValidHeaderSpecific(snes,SNES_COOKIE); 2251 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 2252 SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 2253 } 2254 if (r) *r = snes->vec_func_always; 2255 if (ctx) *ctx = snes->funP; 2256 PetscFunctionReturn(0); 2257 } 2258 2259 #undef __FUNCT__ 2260 #define __FUNCT__ "SNESGetMinimizationFunction" 2261 /*@C 2262 SNESGetMinimizationFunction - Returns the scalar function value for 2263 unconstrained minimization problems. 2264 2265 Not Collective 2266 2267 Input Parameter: 2268 . snes - the SNES context 2269 2270 Output Parameter: 2271 + r - the function (or PETSC_NULL) 2272 - ctx - the function context (or PETSC_NULL) 2273 2274 Notes: 2275 SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 2276 methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 2277 SNESGetFunction(). 2278 2279 Level: advanced 2280 2281 .keywords: SNES, nonlinear, get, function 2282 2283 .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction(), SNESSetFunction() 2284 2285 @*/ 2286 int SNESGetMinimizationFunction(SNES snes,PetscReal *r,void **ctx) 2287 { 2288 PetscFunctionBegin; 2289 PetscValidHeaderSpecific(snes,SNES_COOKIE); 2290 PetscValidScalarPointer(r); 2291 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 2292 SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 2293 } 2294 if (r) *r = snes->fc; 2295 if (ctx) *ctx = snes->umfunP; 2296 PetscFunctionReturn(0); 2297 } 2298 2299 #undef __FUNCT__ 2300 #define __FUNCT__ "SNESSetOptionsPrefix" 2301 /*@C 2302 SNESSetOptionsPrefix - Sets the prefix used for searching for all 2303 SNES options in the database. 2304 2305 Collective on SNES 2306 2307 Input Parameter: 2308 + snes - the SNES context 2309 - prefix - the prefix to prepend to all option names 2310 2311 Notes: 2312 A hyphen (-) must NOT be given at the beginning of the prefix name. 2313 The first character of all runtime options is AUTOMATICALLY the hyphen. 2314 2315 Level: advanced 2316 2317 .keywords: SNES, set, options, prefix, database 2318 2319 .seealso: SNESSetFromOptions() 2320 @*/ 2321 int SNESSetOptionsPrefix(SNES snes,char *prefix) 2322 { 2323 int ierr; 2324 2325 PetscFunctionBegin; 2326 PetscValidHeaderSpecific(snes,SNES_COOKIE); 2327 ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2328 ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 2329 PetscFunctionReturn(0); 2330 } 2331 2332 #undef __FUNCT__ 2333 #define __FUNCT__ "SNESAppendOptionsPrefix" 2334 /*@C 2335 SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 2336 SNES options in the database. 2337 2338 Collective on SNES 2339 2340 Input Parameters: 2341 + snes - the SNES context 2342 - prefix - the prefix to prepend to all option names 2343 2344 Notes: 2345 A hyphen (-) must NOT be given at the beginning of the prefix name. 2346 The first character of all runtime options is AUTOMATICALLY the hyphen. 2347 2348 Level: advanced 2349 2350 .keywords: SNES, append, options, prefix, database 2351 2352 .seealso: SNESGetOptionsPrefix() 2353 @*/ 2354 int SNESAppendOptionsPrefix(SNES snes,char *prefix) 2355 { 2356 int ierr; 2357 2358 PetscFunctionBegin; 2359 PetscValidHeaderSpecific(snes,SNES_COOKIE); 2360 ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2361 ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 2362 PetscFunctionReturn(0); 2363 } 2364 2365 #undef __FUNCT__ 2366 #define __FUNCT__ "SNESGetOptionsPrefix" 2367 /*@C 2368 SNESGetOptionsPrefix - Sets the prefix used for searching for all 2369 SNES options in the database. 2370 2371 Not Collective 2372 2373 Input Parameter: 2374 . snes - the SNES context 2375 2376 Output Parameter: 2377 . prefix - pointer to the prefix string used 2378 2379 Notes: On the fortran side, the user should pass in a string 'prifix' of 2380 sufficient length to hold the prefix. 2381 2382 Level: advanced 2383 2384 .keywords: SNES, get, options, prefix, database 2385 2386 .seealso: SNESAppendOptionsPrefix() 2387 @*/ 2388 int SNESGetOptionsPrefix(SNES snes,char **prefix) 2389 { 2390 int ierr; 2391 2392 PetscFunctionBegin; 2393 PetscValidHeaderSpecific(snes,SNES_COOKIE); 2394 ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2395 PetscFunctionReturn(0); 2396 } 2397 2398 /*MC 2399 SNESRegisterDynamic - Adds a method to the nonlinear solver package. 2400 2401 Synopsis: 2402 int SNESRegisterDynamic(char *name_solver,char *path,char *name_create,int (*routine_create)(SNES)) 2403 2404 Not collective 2405 2406 Input Parameters: 2407 + name_solver - name of a new user-defined solver 2408 . path - path (either absolute or relative) the library containing this solver 2409 . name_create - name of routine to create method context 2410 - routine_create - routine to create method context 2411 2412 Notes: 2413 SNESRegisterDynamic() may be called multiple times to add several user-defined solvers. 2414 2415 If dynamic libraries are used, then the fourth input argument (routine_create) 2416 is ignored. 2417 2418 Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, ${BOPT}, 2419 and others of the form ${any_environmental_variable} occuring in pathname will be 2420 replaced with appropriate values. 2421 2422 Sample usage: 2423 .vb 2424 SNESRegisterDynamic("my_solver",/home/username/my_lib/lib/libg/solaris/mylib.a, 2425 "MySolverCreate",MySolverCreate); 2426 .ve 2427 2428 Then, your solver can be chosen with the procedural interface via 2429 $ SNESSetType(snes,"my_solver") 2430 or at runtime via the option 2431 $ -snes_type my_solver 2432 2433 Level: advanced 2434 2435 .keywords: SNES, nonlinear, register 2436 2437 .seealso: SNESRegisterAll(), SNESRegisterDestroy() 2438 M*/ 2439 2440 #undef __FUNCT__ 2441 #define __FUNCT__ "SNESRegister" 2442 int SNESRegister(char *sname,char *path,char *name,int (*function)(SNES)) 2443 { 2444 char fullname[256]; 2445 int ierr; 2446 2447 PetscFunctionBegin; 2448 ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 2449 ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)())function);CHKERRQ(ierr); 2450 PetscFunctionReturn(0); 2451 } 2452