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