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