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