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 /*@ 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 /*@ 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 ierr = VecAXPY(y,-1.0,snes->afine);CHKERRQ(ierr); 813 } 814 snes->nfuncs++; 815 ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 816 PetscFunctionReturn(0); 817 } 818 819 #undef __FUNCT__ 820 #define __FUNCT__ "SNESComputeJacobian" 821 /*@ 822 SNESComputeJacobian - Computes the Jacobian matrix that has been 823 set with SNESSetJacobian(). 824 825 Collective on SNES and Mat 826 827 Input Parameters: 828 + snes - the SNES context 829 - x - input vector 830 831 Output Parameters: 832 + A - Jacobian matrix 833 . B - optional preconditioning matrix 834 - flag - flag indicating matrix structure 835 836 Notes: 837 Most users should not need to explicitly call this routine, as it 838 is used internally within the nonlinear solvers. 839 840 See KSPSetOperators() for important information about setting the 841 flag parameter. 842 843 Level: developer 844 845 .keywords: SNES, compute, Jacobian, matrix 846 847 .seealso: SNESSetJacobian(), KSPSetOperators() 848 @*/ 849 PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 850 { 851 PetscErrorCode ierr; 852 853 PetscFunctionBegin; 854 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 855 PetscValidHeaderSpecific(X,VEC_COOKIE,2); 856 PetscValidPointer(flg,5); 857 PetscCheckSameComm(snes,1,X,2); 858 if (!snes->computejacobian) PetscFunctionReturn(0); 859 ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 860 *flg = DIFFERENT_NONZERO_PATTERN; 861 PetscStackPush("SNES user Jacobian function"); 862 CHKMEMQ; 863 ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr); 864 CHKMEMQ; 865 PetscStackPop; 866 ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 867 /* make sure user returned a correct Jacobian and preconditioner */ 868 PetscValidHeaderSpecific(*A,MAT_COOKIE,3); 869 PetscValidHeaderSpecific(*B,MAT_COOKIE,4); 870 PetscFunctionReturn(0); 871 } 872 873 #undef __FUNCT__ 874 #define __FUNCT__ "SNESSetJacobian" 875 /*@C 876 SNESSetJacobian - Sets the function to compute Jacobian as well as the 877 location to store the matrix. 878 879 Collective on SNES and Mat 880 881 Input Parameters: 882 + snes - the SNES context 883 . A - Jacobian matrix 884 . B - preconditioner matrix (usually same as the Jacobian) 885 . func - Jacobian evaluation routine 886 - ctx - [optional] user-defined context for private data for the 887 Jacobian evaluation routine (may be PETSC_NULL) 888 889 Calling sequence of func: 890 $ func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 891 892 + x - input vector 893 . A - Jacobian matrix 894 . B - preconditioner matrix, usually the same as A 895 . flag - flag indicating information about the preconditioner matrix 896 structure (same as flag in KSPSetOperators()) 897 - ctx - [optional] user-defined Jacobian context 898 899 Notes: 900 See KSPSetOperators() for important information about setting the flag 901 output parameter in the routine func(). Be sure to read this information! 902 903 The routine func() takes Mat * as the matrix arguments rather than Mat. 904 This allows the Jacobian evaluation routine to replace A and/or B with a 905 completely new new matrix structure (not just different matrix elements) 906 when appropriate, for instance, if the nonzero structure is changing 907 throughout the global iterations. 908 909 Level: beginner 910 911 .keywords: SNES, nonlinear, set, Jacobian, matrix 912 913 .seealso: KSPSetOperators(), SNESSetFunction(), MatSNESMFComputeJacobian(), SNESDefaultComputeJacobianColor() 914 @*/ 915 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 916 { 917 PetscErrorCode ierr; 918 919 PetscFunctionBegin; 920 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 921 if (A) PetscValidHeaderSpecific(A,MAT_COOKIE,2); 922 if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3); 923 if (A) PetscCheckSameComm(snes,1,A,2); 924 if (B) PetscCheckSameComm(snes,1,B,2); 925 if (func) snes->computejacobian = func; 926 if (ctx) snes->jacP = ctx; 927 if (A) { 928 if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 929 snes->jacobian = A; 930 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 931 } 932 if (B) { 933 if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 934 snes->jacobian_pre = B; 935 ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 936 } 937 ierr = KSPSetOperators(snes->ksp,A,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 938 PetscFunctionReturn(0); 939 } 940 941 #undef __FUNCT__ 942 #define __FUNCT__ "SNESGetJacobian" 943 /*@C 944 SNESGetJacobian - Returns the Jacobian matrix and optionally the user 945 provided context for evaluating the Jacobian. 946 947 Not Collective, but Mat object will be parallel if SNES object is 948 949 Input Parameter: 950 . snes - the nonlinear solver context 951 952 Output Parameters: 953 + A - location to stash Jacobian matrix (or PETSC_NULL) 954 . B - location to stash preconditioner matrix (or PETSC_NULL) 955 . func - location to put Jacobian function (or PETSC_NULL) 956 - ctx - location to stash Jacobian ctx (or PETSC_NULL) 957 958 Level: advanced 959 960 .seealso: SNESSetJacobian(), SNESComputeJacobian() 961 @*/ 962 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetJacobian(SNES snes,Mat *A,Mat *B,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx) 963 { 964 PetscFunctionBegin; 965 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 966 if (A) *A = snes->jacobian; 967 if (B) *B = snes->jacobian_pre; 968 if (func) *func = snes->computejacobian; 969 if (ctx) *ctx = snes->jacP; 970 PetscFunctionReturn(0); 971 } 972 973 /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 974 EXTERN PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*); 975 976 #undef __FUNCT__ 977 #define __FUNCT__ "SNESSetUp" 978 /*@ 979 SNESSetUp - Sets up the internal data structures for the later use 980 of a nonlinear solver. 981 982 Collective on SNES 983 984 Input Parameters: 985 . snes - the SNES context 986 987 Notes: 988 For basic use of the SNES solvers the user need not explicitly call 989 SNESSetUp(), since these actions will automatically occur during 990 the call to SNESSolve(). However, if one wishes to control this 991 phase separately, SNESSetUp() should be called after SNESCreate() 992 and optional routines of the form SNESSetXXX(), but before SNESSolve(). 993 994 Level: advanced 995 996 .keywords: SNES, nonlinear, setup 997 998 .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 999 @*/ 1000 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUp(SNES snes) 1001 { 1002 PetscErrorCode ierr; 1003 PetscTruth flg, iseqtr; 1004 1005 PetscFunctionBegin; 1006 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1007 if (snes->setupcalled) PetscFunctionReturn(0); 1008 1009 ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr); 1010 /* 1011 This version replaces the user provided Jacobian matrix with a 1012 matrix-free version but still employs the user-provided preconditioner matrix 1013 */ 1014 if (flg) { 1015 Mat J; 1016 ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1017 ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 1018 ierr = PetscLogInfo((snes,"SNESSetUp: Setting default matrix-free operator routines\n"));CHKERRQ(ierr); 1019 ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 1020 ierr = MatDestroy(J);CHKERRQ(ierr); 1021 } 1022 1023 #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) && !defined(PETSC_USE_MAT_SINGLE) 1024 ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr); 1025 if (flg) { 1026 Mat J; 1027 ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1028 ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 1029 ierr = MatDestroy(J);CHKERRQ(ierr); 1030 } 1031 #endif 1032 1033 ierr = PetscOptionsHasName(snes->prefix,"-snes_mf",&flg);CHKERRQ(ierr); 1034 /* 1035 This version replaces both the user-provided Jacobian and the user- 1036 provided preconditioner matrix with the default matrix free version. 1037 */ 1038 if (flg) { 1039 Mat J; 1040 KSP ksp; 1041 PC pc; 1042 1043 ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1044 ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 1045 ierr = PetscLogInfo((snes,"SNESSetUp: Setting default matrix-free operator and preconditioner routines\n"));CHKERRQ(ierr); 1046 ierr = SNESSetJacobian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr); 1047 ierr = MatDestroy(J);CHKERRQ(ierr); 1048 1049 /* force no preconditioner */ 1050 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1051 ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 1052 ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr); 1053 if (!flg) { 1054 ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 1055 } 1056 } 1057 1058 if (!snes->vec_func && !snes->afine) { 1059 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 1060 } 1061 if (!snes->computefunction && !snes->afine) { 1062 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 1063 } 1064 if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option"); 1065 if (snes->vec_func == snes->vec_sol) { 1066 SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector"); 1067 } 1068 1069 /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 1070 ierr = PetscTypeCompare((PetscObject)snes,SNESTR,&iseqtr);CHKERRQ(ierr); 1071 if (snes->ksp_ewconv && !iseqtr) { 1072 KSP ksp; 1073 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1074 ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,snes);CHKERRQ(ierr); 1075 } 1076 1077 if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);} 1078 snes->setupcalled = PETSC_TRUE; 1079 PetscFunctionReturn(0); 1080 } 1081 1082 #undef __FUNCT__ 1083 #define __FUNCT__ "SNESDestroy" 1084 /*@ 1085 SNESDestroy - Destroys the nonlinear solver context that was created 1086 with SNESCreate(). 1087 1088 Collective on SNES 1089 1090 Input Parameter: 1091 . snes - the SNES context 1092 1093 Level: beginner 1094 1095 .keywords: SNES, nonlinear, destroy 1096 1097 .seealso: SNESCreate(), SNESSolve() 1098 @*/ 1099 PetscErrorCode PETSCSNES_DLLEXPORT SNESDestroy(SNES snes) 1100 { 1101 PetscInt i; 1102 PetscErrorCode ierr; 1103 1104 PetscFunctionBegin; 1105 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1106 if (--snes->refct > 0) PetscFunctionReturn(0); 1107 1108 /* if memory was published with AMS then destroy it */ 1109 ierr = PetscObjectDepublish(snes);CHKERRQ(ierr); 1110 1111 if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);} 1112 if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);} 1113 if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 1114 if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 1115 if (snes->afine) {ierr = VecDestroy(snes->afine);CHKERRQ(ierr);} 1116 ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr); 1117 if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);} 1118 for (i=0; i<snes->numbermonitors; i++) { 1119 if (snes->monitordestroy[i]) { 1120 ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr); 1121 } 1122 } 1123 ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr); 1124 PetscFunctionReturn(0); 1125 } 1126 1127 /* ----------- Routines to set solver parameters ---------- */ 1128 1129 #undef __FUNCT__ 1130 #define __FUNCT__ "SNESSetTolerances" 1131 /*@ 1132 SNESSetTolerances - Sets various parameters used in convergence tests. 1133 1134 Collective on SNES 1135 1136 Input Parameters: 1137 + snes - the SNES context 1138 . abstol - absolute convergence tolerance 1139 . rtol - relative convergence tolerance 1140 . stol - convergence tolerance in terms of the norm 1141 of the change in the solution between steps 1142 . maxit - maximum number of iterations 1143 - maxf - maximum number of function evaluations 1144 1145 Options Database Keys: 1146 + -snes_atol <abstol> - Sets abstol 1147 . -snes_rtol <rtol> - Sets rtol 1148 . -snes_stol <stol> - Sets stol 1149 . -snes_max_it <maxit> - Sets maxit 1150 - -snes_max_funcs <maxf> - Sets maxf 1151 1152 Notes: 1153 The default maximum number of iterations is 50. 1154 The default maximum number of function evaluations is 1000. 1155 1156 Level: intermediate 1157 1158 .keywords: SNES, nonlinear, set, convergence, tolerances 1159 1160 .seealso: SNESSetTrustRegionTolerance() 1161 @*/ 1162 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf) 1163 { 1164 PetscFunctionBegin; 1165 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1166 if (abstol != PETSC_DEFAULT) snes->abstol = abstol; 1167 if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1168 if (stol != PETSC_DEFAULT) snes->xtol = stol; 1169 if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1170 if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 1171 PetscFunctionReturn(0); 1172 } 1173 1174 #undef __FUNCT__ 1175 #define __FUNCT__ "SNESGetTolerances" 1176 /*@ 1177 SNESGetTolerances - Gets various parameters used in convergence tests. 1178 1179 Not Collective 1180 1181 Input Parameters: 1182 + snes - the SNES context 1183 . abstol - absolute convergence tolerance 1184 . rtol - relative convergence tolerance 1185 . stol - convergence tolerance in terms of the norm 1186 of the change in the solution between steps 1187 . maxit - maximum number of iterations 1188 - maxf - maximum number of function evaluations 1189 1190 Notes: 1191 The user can specify PETSC_NULL for any parameter that is not needed. 1192 1193 Level: intermediate 1194 1195 .keywords: SNES, nonlinear, get, convergence, tolerances 1196 1197 .seealso: SNESSetTolerances() 1198 @*/ 1199 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetTolerances(SNES snes,PetscReal *abstol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf) 1200 { 1201 PetscFunctionBegin; 1202 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1203 if (abstol) *abstol = snes->abstol; 1204 if (rtol) *rtol = snes->rtol; 1205 if (stol) *stol = snes->xtol; 1206 if (maxit) *maxit = snes->max_its; 1207 if (maxf) *maxf = snes->max_funcs; 1208 PetscFunctionReturn(0); 1209 } 1210 1211 #undef __FUNCT__ 1212 #define __FUNCT__ "SNESSetTrustRegionTolerance" 1213 /*@ 1214 SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 1215 1216 Collective on SNES 1217 1218 Input Parameters: 1219 + snes - the SNES context 1220 - tol - tolerance 1221 1222 Options Database Key: 1223 . -snes_trtol <tol> - Sets tol 1224 1225 Level: intermediate 1226 1227 .keywords: SNES, nonlinear, set, trust region, tolerance 1228 1229 .seealso: SNESSetTolerances() 1230 @*/ 1231 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTrustRegionTolerance(SNES snes,PetscReal tol) 1232 { 1233 PetscFunctionBegin; 1234 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1235 snes->deltatol = tol; 1236 PetscFunctionReturn(0); 1237 } 1238 1239 /* 1240 Duplicate the lg monitors for SNES from KSP; for some reason with 1241 dynamic libraries things don't work under Sun4 if we just use 1242 macros instead of functions 1243 */ 1244 #undef __FUNCT__ 1245 #define __FUNCT__ "SNESLGMonitor" 1246 PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitor(SNES snes,PetscInt it,PetscReal norm,void *ctx) 1247 { 1248 PetscErrorCode ierr; 1249 1250 PetscFunctionBegin; 1251 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1252 ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr); 1253 PetscFunctionReturn(0); 1254 } 1255 1256 #undef __FUNCT__ 1257 #define __FUNCT__ "SNESLGMonitorCreate" 1258 PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1259 { 1260 PetscErrorCode ierr; 1261 1262 PetscFunctionBegin; 1263 ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr); 1264 PetscFunctionReturn(0); 1265 } 1266 1267 #undef __FUNCT__ 1268 #define __FUNCT__ "SNESLGMonitorDestroy" 1269 PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitorDestroy(PetscDrawLG draw) 1270 { 1271 PetscErrorCode ierr; 1272 1273 PetscFunctionBegin; 1274 ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr); 1275 PetscFunctionReturn(0); 1276 } 1277 1278 /* ------------ Routines to set performance monitoring options ----------- */ 1279 1280 #undef __FUNCT__ 1281 #define __FUNCT__ "SNESSetMonitor" 1282 /*@C 1283 SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every 1284 iteration of the nonlinear solver to display the iteration's 1285 progress. 1286 1287 Collective on SNES 1288 1289 Input Parameters: 1290 + snes - the SNES context 1291 . func - monitoring routine 1292 . mctx - [optional] user-defined context for private data for the 1293 monitor routine (use PETSC_NULL if no context is desitre) 1294 - monitordestroy - [optional] routine that frees monitor context 1295 (may be PETSC_NULL) 1296 1297 Calling sequence of func: 1298 $ int func(SNES snes,PetscInt its, PetscReal norm,void *mctx) 1299 1300 + snes - the SNES context 1301 . its - iteration number 1302 . norm - 2-norm function value (may be estimated) 1303 - mctx - [optional] monitoring context 1304 1305 Options Database Keys: 1306 + -snes_monitor - sets SNESDefaultMonitor() 1307 . -snes_xmonitor - sets line graph monitor, 1308 uses SNESLGMonitorCreate() 1309 _ -snes_cancelmonitors - cancels all monitors that have 1310 been hardwired into a code by 1311 calls to SNESSetMonitor(), but 1312 does not cancel those set via 1313 the options database. 1314 1315 Notes: 1316 Several different monitoring routines may be set by calling 1317 SNESSetMonitor() multiple times; all will be called in the 1318 order in which they were set. 1319 1320 Level: intermediate 1321 1322 .keywords: SNES, nonlinear, set, monitor 1323 1324 .seealso: SNESDefaultMonitor(), SNESClearMonitor() 1325 @*/ 1326 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMonitor(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void*)) 1327 { 1328 PetscFunctionBegin; 1329 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1330 if (snes->numbermonitors >= MAXSNESMONITORS) { 1331 SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1332 } 1333 snes->monitor[snes->numbermonitors] = func; 1334 snes->monitordestroy[snes->numbermonitors] = monitordestroy; 1335 snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 1336 PetscFunctionReturn(0); 1337 } 1338 1339 #undef __FUNCT__ 1340 #define __FUNCT__ "SNESClearMonitor" 1341 /*@C 1342 SNESClearMonitor - Clears all the monitor functions for a SNES object. 1343 1344 Collective on SNES 1345 1346 Input Parameters: 1347 . snes - the SNES context 1348 1349 Options Database Key: 1350 . -snes_cancelmonitors - cancels all monitors that have been hardwired 1351 into a code by calls to SNESSetMonitor(), but does not cancel those 1352 set via the options database 1353 1354 Notes: 1355 There is no way to clear one specific monitor from a SNES object. 1356 1357 Level: intermediate 1358 1359 .keywords: SNES, nonlinear, set, monitor 1360 1361 .seealso: SNESDefaultMonitor(), SNESSetMonitor() 1362 @*/ 1363 PetscErrorCode PETSCSNES_DLLEXPORT SNESClearMonitor(SNES snes) 1364 { 1365 PetscFunctionBegin; 1366 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1367 snes->numbermonitors = 0; 1368 PetscFunctionReturn(0); 1369 } 1370 1371 #undef __FUNCT__ 1372 #define __FUNCT__ "SNESSetConvergenceTest" 1373 /*@C 1374 SNESSetConvergenceTest - Sets the function that is to be used 1375 to test for convergence of the nonlinear iterative solution. 1376 1377 Collective on SNES 1378 1379 Input Parameters: 1380 + snes - the SNES context 1381 . func - routine to test for convergence 1382 - cctx - [optional] context for private data for the convergence routine 1383 (may be PETSC_NULL) 1384 1385 Calling sequence of func: 1386 $ PetscErrorCode func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx) 1387 1388 + snes - the SNES context 1389 . cctx - [optional] convergence context 1390 . reason - reason for convergence/divergence 1391 . xnorm - 2-norm of current iterate 1392 . gnorm - 2-norm of current step 1393 - f - 2-norm of function 1394 1395 Level: advanced 1396 1397 .keywords: SNES, nonlinear, set, convergence, test 1398 1399 .seealso: SNESConverged_LS(), SNESConverged_TR() 1400 @*/ 1401 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx) 1402 { 1403 PetscFunctionBegin; 1404 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1405 (snes)->converged = func; 1406 (snes)->cnvP = cctx; 1407 PetscFunctionReturn(0); 1408 } 1409 1410 #undef __FUNCT__ 1411 #define __FUNCT__ "SNESGetConvergedReason" 1412 /*@ 1413 SNESGetConvergedReason - Gets the reason the SNES iteration was stopped. 1414 1415 Not Collective 1416 1417 Input Parameter: 1418 . snes - the SNES context 1419 1420 Output Parameter: 1421 . reason - negative value indicates diverged, positive value converged, see petscsnes.h or the 1422 manual pages for the individual convergence tests for complete lists 1423 1424 Level: intermediate 1425 1426 Notes: Can only be called after the call the SNESSolve() is complete. 1427 1428 .keywords: SNES, nonlinear, set, convergence, test 1429 1430 .seealso: SNESSetConvergenceTest(), SNESConverged_LS(), SNESConverged_TR(), SNESConvergedReason 1431 @*/ 1432 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason) 1433 { 1434 PetscFunctionBegin; 1435 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1436 PetscValidPointer(reason,2); 1437 *reason = snes->reason; 1438 PetscFunctionReturn(0); 1439 } 1440 1441 #undef __FUNCT__ 1442 #define __FUNCT__ "SNESSetConvergenceHistory" 1443 /*@ 1444 SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1445 1446 Collective on SNES 1447 1448 Input Parameters: 1449 + snes - iterative context obtained from SNESCreate() 1450 . a - array to hold history 1451 . its - integer array holds the number of linear iterations for each solve. 1452 . na - size of a and its 1453 - reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero, 1454 else it continues storing new values for new nonlinear solves after the old ones 1455 1456 Notes: 1457 If set, this array will contain the function norms computed 1458 at each step. 1459 1460 This routine is useful, e.g., when running a code for purposes 1461 of accurate performance monitoring, when no I/O should be done 1462 during the section of code that is being timed. 1463 1464 Level: intermediate 1465 1466 .keywords: SNES, set, convergence, history 1467 1468 .seealso: SNESGetConvergenceHistory() 1469 1470 @*/ 1471 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt *its,PetscInt na,PetscTruth reset) 1472 { 1473 PetscFunctionBegin; 1474 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1475 if (na) PetscValidScalarPointer(a,2); 1476 snes->conv_hist = a; 1477 snes->conv_hist_its = its; 1478 snes->conv_hist_max = na; 1479 snes->conv_hist_reset = reset; 1480 PetscFunctionReturn(0); 1481 } 1482 1483 #undef __FUNCT__ 1484 #define __FUNCT__ "SNESGetConvergenceHistory" 1485 /*@C 1486 SNESGetConvergenceHistory - Gets the array used to hold the convergence history. 1487 1488 Collective on SNES 1489 1490 Input Parameter: 1491 . snes - iterative context obtained from SNESCreate() 1492 1493 Output Parameters: 1494 . a - array to hold history 1495 . its - integer array holds the number of linear iterations (or 1496 negative if not converged) for each solve. 1497 - na - size of a and its 1498 1499 Notes: 1500 The calling sequence for this routine in Fortran is 1501 $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr) 1502 1503 This routine is useful, e.g., when running a code for purposes 1504 of accurate performance monitoring, when no I/O should be done 1505 during the section of code that is being timed. 1506 1507 Level: intermediate 1508 1509 .keywords: SNES, get, convergence, history 1510 1511 .seealso: SNESSetConvergencHistory() 1512 1513 @*/ 1514 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na) 1515 { 1516 PetscFunctionBegin; 1517 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1518 if (a) *a = snes->conv_hist; 1519 if (its) *its = snes->conv_hist_its; 1520 if (na) *na = snes->conv_hist_len; 1521 PetscFunctionReturn(0); 1522 } 1523 1524 #undef __FUNCT__ 1525 #define __FUNCT__ "SNESSetUpdate" 1526 /*@C 1527 SNESSetUpdate - Sets the general-purpose update function called 1528 at the beginning of every step of the iteration. 1529 1530 Collective on SNES 1531 1532 Input Parameters: 1533 . snes - The nonlinear solver context 1534 . func - The function 1535 1536 Calling sequence of func: 1537 . func (SNES snes, PetscInt step); 1538 1539 . step - The current step of the iteration 1540 1541 Level: intermediate 1542 1543 .keywords: SNES, update 1544 1545 .seealso SNESDefaultUpdate(), SNESSetRhsBC(), SNESSetSolutionBC() 1546 @*/ 1547 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt)) 1548 { 1549 PetscFunctionBegin; 1550 PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 1551 snes->update = func; 1552 PetscFunctionReturn(0); 1553 } 1554 1555 #undef __FUNCT__ 1556 #define __FUNCT__ "SNESDefaultUpdate" 1557 /*@ 1558 SNESDefaultUpdate - The default update function which does nothing. 1559 1560 Not collective 1561 1562 Input Parameters: 1563 . snes - The nonlinear solver context 1564 . step - The current step of the iteration 1565 1566 Level: intermediate 1567 1568 .keywords: SNES, update 1569 .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultSolutionBC() 1570 @*/ 1571 PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultUpdate(SNES snes, PetscInt step) 1572 { 1573 PetscFunctionBegin; 1574 PetscFunctionReturn(0); 1575 } 1576 1577 #undef __FUNCT__ 1578 #define __FUNCT__ "SNESScaleStep_Private" 1579 /* 1580 SNESScaleStep_Private - Scales a step so that its length is less than the 1581 positive parameter delta. 1582 1583 Input Parameters: 1584 + snes - the SNES context 1585 . y - approximate solution of linear system 1586 . fnorm - 2-norm of current function 1587 - delta - trust region size 1588 1589 Output Parameters: 1590 + gpnorm - predicted function norm at the new point, assuming local 1591 linearization. The value is zero if the step lies within the trust 1592 region, and exceeds zero otherwise. 1593 - ynorm - 2-norm of the step 1594 1595 Note: 1596 For non-trust region methods such as SNESLS, the parameter delta 1597 is set to be the maximum allowable step size. 1598 1599 .keywords: SNES, nonlinear, scale, step 1600 */ 1601 PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm) 1602 { 1603 PetscReal nrm; 1604 PetscScalar cnorm; 1605 PetscErrorCode ierr; 1606 1607 PetscFunctionBegin; 1608 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1609 PetscValidHeaderSpecific(y,VEC_COOKIE,2); 1610 PetscCheckSameComm(snes,1,y,2); 1611 1612 ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr); 1613 if (nrm > *delta) { 1614 nrm = *delta/nrm; 1615 *gpnorm = (1.0 - nrm)*(*fnorm); 1616 cnorm = nrm; 1617 ierr = VecScale(y,cnorm);CHKERRQ(ierr); 1618 *ynorm = *delta; 1619 } else { 1620 *gpnorm = 0.0; 1621 *ynorm = nrm; 1622 } 1623 PetscFunctionReturn(0); 1624 } 1625 1626 #undef __FUNCT__ 1627 #define __FUNCT__ "SNESSolve" 1628 /*@ 1629 SNESSolve - Solves a nonlinear system F(x) = b. 1630 Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX(). 1631 1632 Collective on SNES 1633 1634 Input Parameters: 1635 + snes - the SNES context 1636 . b - the constant part of the equation, or PETSC_NULL to use zero. 1637 - x - the solution vector, or PETSC_NULL if it was set with SNESSetSolution() 1638 1639 Notes: 1640 The user should initialize the vector,x, with the initial guess 1641 for the nonlinear solve prior to calling SNESSolve. In particular, 1642 to employ an initial guess of zero, the user should explicitly set 1643 this vector to zero by calling VecSet(). 1644 1645 Level: beginner 1646 1647 .keywords: SNES, nonlinear, solve 1648 1649 .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetRhs(), SNESSetSolution() 1650 @*/ 1651 PetscErrorCode PETSCSNES_DLLEXPORT SNESSolve(SNES snes,Vec b,Vec x) 1652 { 1653 PetscErrorCode ierr; 1654 PetscTruth flg; 1655 1656 PetscFunctionBegin; 1657 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1658 if (!snes->solve) SETERRQ(PETSC_ERR_ORDER,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()"); 1659 1660 if (b) { 1661 ierr = SNESSetRhs(snes, b); CHKERRQ(ierr); 1662 if (!snes->vec_func) { 1663 Vec r; 1664 1665 ierr = VecDuplicate(b, &r); CHKERRQ(ierr); 1666 ierr = SNESSetFunction(snes, r, PETSC_NULL, PETSC_NULL); CHKERRQ(ierr); 1667 } 1668 } 1669 if (x) { 1670 PetscValidHeaderSpecific(x,VEC_COOKIE,3); 1671 PetscCheckSameComm(snes,1,x,3); 1672 } else { 1673 ierr = SNESGetSolution(snes, &x); CHKERRQ(ierr); 1674 if (!x) { 1675 ierr = VecDuplicate(snes->vec_func_always, &x); CHKERRQ(ierr); 1676 } 1677 } 1678 snes->vec_sol = snes->vec_sol_always = x; 1679 if (!snes->setupcalled) { 1680 ierr = SNESSetUp(snes);CHKERRQ(ierr); 1681 } 1682 if (snes->conv_hist_reset) snes->conv_hist_len = 0; 1683 ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 1684 snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0; 1685 1686 ierr = PetscExceptionTry1((*(snes)->solve)(snes),PETSC_ERR_ARG_DOMAIN); 1687 if (PetscExceptionValue(ierr)) { 1688 /* this means that a caller above me has also tryed this exception so I don't handle it here, pass it up */ 1689 PetscErrorCode pierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(pierr); 1690 } else if (PetscExceptionCaught(ierr,PETSC_ERR_ARG_DOMAIN)) { 1691 /* translate exception into SNES not converged reason */ 1692 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1693 ierr = 0; 1694 } 1695 CHKERRQ(ierr); 1696 1697 ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 1698 ierr = PetscOptionsHasName(snes->prefix,"-snes_view",&flg);CHKERRQ(ierr); 1699 if (flg && !PetscPreLoadingOn) { ierr = SNESView(snes,PETSC_VIEWER_STDOUT_(snes->comm));CHKERRQ(ierr); } 1700 ierr = PetscOptionsHasName(snes->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr); 1701 if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); } 1702 if (snes->printreason) { 1703 if (snes->reason > 0) { 1704 ierr = PetscPrintf(snes->comm,"Nonlinear solve converged due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr); 1705 } else { 1706 ierr = PetscPrintf(snes->comm,"Nonlinear solve did not converge due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr); 1707 } 1708 } 1709 1710 PetscFunctionReturn(0); 1711 } 1712 1713 /* --------- Internal routines for SNES Package --------- */ 1714 1715 #undef __FUNCT__ 1716 #define __FUNCT__ "SNESSetType" 1717 /*@C 1718 SNESSetType - Sets the method for the nonlinear solver. 1719 1720 Collective on SNES 1721 1722 Input Parameters: 1723 + snes - the SNES context 1724 - type - a known method 1725 1726 Options Database Key: 1727 . -snes_type <type> - Sets the method; use -help for a list 1728 of available methods (for instance, ls or tr) 1729 1730 Notes: 1731 See "petsc/include/petscsnes.h" for available methods (for instance) 1732 + SNESLS - Newton's method with line search 1733 (systems of nonlinear equations) 1734 . SNESTR - Newton's method with trust region 1735 (systems of nonlinear equations) 1736 1737 Normally, it is best to use the SNESSetFromOptions() command and then 1738 set the SNES solver type from the options database rather than by using 1739 this routine. Using the options database provides the user with 1740 maximum flexibility in evaluating the many nonlinear solvers. 1741 The SNESSetType() routine is provided for those situations where it 1742 is necessary to set the nonlinear solver independently of the command 1743 line or options database. This might be the case, for example, when 1744 the choice of solver changes during the execution of the program, 1745 and the user's application is taking responsibility for choosing the 1746 appropriate method. 1747 1748 Level: intermediate 1749 1750 .keywords: SNES, set, type 1751 1752 .seealso: SNESType, SNESCreate() 1753 1754 @*/ 1755 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetType(SNES snes,SNESType type) 1756 { 1757 PetscErrorCode ierr,(*r)(SNES); 1758 PetscTruth match; 1759 1760 PetscFunctionBegin; 1761 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1762 PetscValidCharPointer(type,2); 1763 1764 ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr); 1765 if (match) PetscFunctionReturn(0); 1766 1767 if (snes->setupcalled) { 1768 snes->setupcalled = PETSC_FALSE; 1769 ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr); 1770 snes->data = 0; 1771 } 1772 1773 /* Get the function pointers for the iterative method requested */ 1774 if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 1775 ierr = PetscFListFind(snes->comm,SNESList,type,(void (**)(void)) &r);CHKERRQ(ierr); 1776 if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type); 1777 if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);} 1778 snes->data = 0; 1779 ierr = (*r)(snes);CHKERRQ(ierr); 1780 ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr); 1781 PetscFunctionReturn(0); 1782 } 1783 1784 1785 /* --------------------------------------------------------------------- */ 1786 #undef __FUNCT__ 1787 #define __FUNCT__ "SNESRegisterDestroy" 1788 /*@ 1789 SNESRegisterDestroy - Frees the list of nonlinear solvers that were 1790 registered by SNESRegisterDynamic(). 1791 1792 Not Collective 1793 1794 Level: advanced 1795 1796 .keywords: SNES, nonlinear, register, destroy 1797 1798 .seealso: SNESRegisterAll(), SNESRegisterAll() 1799 @*/ 1800 PetscErrorCode PETSCSNES_DLLEXPORT SNESRegisterDestroy(void) 1801 { 1802 PetscErrorCode ierr; 1803 1804 PetscFunctionBegin; 1805 if (SNESList) { 1806 ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr); 1807 SNESList = 0; 1808 } 1809 SNESRegisterAllCalled = PETSC_FALSE; 1810 PetscFunctionReturn(0); 1811 } 1812 1813 #undef __FUNCT__ 1814 #define __FUNCT__ "SNESGetType" 1815 /*@C 1816 SNESGetType - Gets the SNES method type and name (as a string). 1817 1818 Not Collective 1819 1820 Input Parameter: 1821 . snes - nonlinear solver context 1822 1823 Output Parameter: 1824 . type - SNES method (a character string) 1825 1826 Level: intermediate 1827 1828 .keywords: SNES, nonlinear, get, type, name 1829 @*/ 1830 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetType(SNES snes,SNESType *type) 1831 { 1832 PetscFunctionBegin; 1833 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1834 PetscValidPointer(type,2); 1835 *type = snes->type_name; 1836 PetscFunctionReturn(0); 1837 } 1838 1839 #undef __FUNCT__ 1840 #define __FUNCT__ "SNESGetSolution" 1841 /*@ 1842 SNESGetSolution - Returns the vector where the approximate solution is 1843 stored. 1844 1845 Not Collective, but Vec is parallel if SNES is parallel 1846 1847 Input Parameter: 1848 . snes - the SNES context 1849 1850 Output Parameter: 1851 . x - the solution 1852 1853 Level: intermediate 1854 1855 .keywords: SNES, nonlinear, get, solution 1856 1857 .seealso: SNESSetSolution(), SNESGetFunction(), SNESGetSolutionUpdate() 1858 @*/ 1859 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolution(SNES snes,Vec *x) 1860 { 1861 PetscFunctionBegin; 1862 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1863 PetscValidPointer(x,2); 1864 *x = snes->vec_sol_always; 1865 PetscFunctionReturn(0); 1866 } 1867 1868 #undef __FUNCT__ 1869 #define __FUNCT__ "SNESSetSolution" 1870 /*@C 1871 SNESSetSolution - Sets the vector where the approximate solution is stored. 1872 1873 Not Collective, but Vec is parallel if SNES is parallel 1874 1875 Input Parameters: 1876 + snes - the SNES context 1877 - x - the solution 1878 1879 Output Parameter: 1880 1881 Level: intermediate 1882 1883 Notes: this is not normally used, rather one simply calls SNESSolve() with 1884 the appropriate solution vector. 1885 1886 .keywords: SNES, nonlinear, set, solution 1887 1888 .seealso: SNESGetSolution(), SNESGetFunction(), SNESGetSolutionUpdate() 1889 @*/ 1890 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetSolution(SNES snes,Vec x) 1891 { 1892 PetscFunctionBegin; 1893 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1894 PetscValidHeaderSpecific(x,VEC_COOKIE,2); 1895 PetscCheckSameComm(snes,1,x,2); 1896 snes->vec_sol_always = x; 1897 PetscFunctionReturn(0); 1898 } 1899 1900 #undef __FUNCT__ 1901 #define __FUNCT__ "SNESGetSolutionUpdate" 1902 /*@ 1903 SNESGetSolutionUpdate - Returns the vector where the solution update is 1904 stored. 1905 1906 Not Collective, but Vec is parallel if SNES is parallel 1907 1908 Input Parameter: 1909 . snes - the SNES context 1910 1911 Output Parameter: 1912 . x - the solution update 1913 1914 Level: advanced 1915 1916 .keywords: SNES, nonlinear, get, solution, update 1917 1918 .seealso: SNESGetSolution(), SNESGetFunction 1919 @*/ 1920 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolutionUpdate(SNES snes,Vec *x) 1921 { 1922 PetscFunctionBegin; 1923 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1924 PetscValidPointer(x,2); 1925 *x = snes->vec_sol_update_always; 1926 PetscFunctionReturn(0); 1927 } 1928 1929 #undef __FUNCT__ 1930 #define __FUNCT__ "SNESGetFunction" 1931 /*@C 1932 SNESGetFunction - Returns the vector where the function is stored. 1933 1934 Not Collective, but Vec is parallel if SNES is parallel 1935 1936 Input Parameter: 1937 . snes - the SNES context 1938 1939 Output Parameter: 1940 + r - the function (or PETSC_NULL) 1941 . func - the function (or PETSC_NULL) 1942 - ctx - the function context (or PETSC_NULL) 1943 1944 Level: advanced 1945 1946 .keywords: SNES, nonlinear, get, function 1947 1948 .seealso: SNESSetFunction(), SNESGetSolution() 1949 @*/ 1950 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx) 1951 { 1952 PetscFunctionBegin; 1953 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1954 if (r) *r = snes->vec_func_always; 1955 if (func) *func = snes->computefunction; 1956 if (ctx) *ctx = snes->funP; 1957 PetscFunctionReturn(0); 1958 } 1959 1960 #undef __FUNCT__ 1961 #define __FUNCT__ "SNESSetOptionsPrefix" 1962 /*@C 1963 SNESSetOptionsPrefix - Sets the prefix used for searching for all 1964 SNES options in the database. 1965 1966 Collective on SNES 1967 1968 Input Parameter: 1969 + snes - the SNES context 1970 - prefix - the prefix to prepend to all option names 1971 1972 Notes: 1973 A hyphen (-) must NOT be given at the beginning of the prefix name. 1974 The first character of all runtime options is AUTOMATICALLY the hyphen. 1975 1976 Level: advanced 1977 1978 .keywords: SNES, set, options, prefix, database 1979 1980 .seealso: SNESSetFromOptions() 1981 @*/ 1982 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetOptionsPrefix(SNES snes,const char prefix[]) 1983 { 1984 PetscErrorCode ierr; 1985 1986 PetscFunctionBegin; 1987 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1988 ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 1989 ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 1990 PetscFunctionReturn(0); 1991 } 1992 1993 #undef __FUNCT__ 1994 #define __FUNCT__ "SNESAppendOptionsPrefix" 1995 /*@C 1996 SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 1997 SNES options in the database. 1998 1999 Collective on SNES 2000 2001 Input Parameters: 2002 + snes - the SNES context 2003 - prefix - the prefix to prepend to all option names 2004 2005 Notes: 2006 A hyphen (-) must NOT be given at the beginning of the prefix name. 2007 The first character of all runtime options is AUTOMATICALLY the hyphen. 2008 2009 Level: advanced 2010 2011 .keywords: SNES, append, options, prefix, database 2012 2013 .seealso: SNESGetOptionsPrefix() 2014 @*/ 2015 PetscErrorCode PETSCSNES_DLLEXPORT SNESAppendOptionsPrefix(SNES snes,const char prefix[]) 2016 { 2017 PetscErrorCode ierr; 2018 2019 PetscFunctionBegin; 2020 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2021 ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2022 ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 2023 PetscFunctionReturn(0); 2024 } 2025 2026 #undef __FUNCT__ 2027 #define __FUNCT__ "SNESGetOptionsPrefix" 2028 /*@C 2029 SNESGetOptionsPrefix - Sets the prefix used for searching for all 2030 SNES options in the database. 2031 2032 Not Collective 2033 2034 Input Parameter: 2035 . snes - the SNES context 2036 2037 Output Parameter: 2038 . prefix - pointer to the prefix string used 2039 2040 Notes: On the fortran side, the user should pass in a string 'prifix' of 2041 sufficient length to hold the prefix. 2042 2043 Level: advanced 2044 2045 .keywords: SNES, get, options, prefix, database 2046 2047 .seealso: SNESAppendOptionsPrefix() 2048 @*/ 2049 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetOptionsPrefix(SNES snes,const char *prefix[]) 2050 { 2051 PetscErrorCode ierr; 2052 2053 PetscFunctionBegin; 2054 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2055 ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2056 PetscFunctionReturn(0); 2057 } 2058 2059 2060 #undef __FUNCT__ 2061 #define __FUNCT__ "SNESRegister" 2062 /*@C 2063 SNESRegister - See SNESRegisterDynamic() 2064 2065 Level: advanced 2066 @*/ 2067 PetscErrorCode PETSCSNES_DLLEXPORT SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES)) 2068 { 2069 char fullname[PETSC_MAX_PATH_LEN]; 2070 PetscErrorCode ierr; 2071 2072 PetscFunctionBegin; 2073 ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 2074 ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 2075 PetscFunctionReturn(0); 2076 } 2077 2078 #undef __FUNCT__ 2079 #define __FUNCT__ "SNESTestLocalMin" 2080 PetscErrorCode PETSCSNES_DLLEXPORT SNESTestLocalMin(SNES snes) 2081 { 2082 PetscErrorCode ierr; 2083 PetscInt N,i,j; 2084 Vec u,uh,fh; 2085 PetscScalar value; 2086 PetscReal norm; 2087 2088 PetscFunctionBegin; 2089 ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr); 2090 ierr = VecDuplicate(u,&uh);CHKERRQ(ierr); 2091 ierr = VecDuplicate(u,&fh);CHKERRQ(ierr); 2092 2093 /* currently only works for sequential */ 2094 ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n"); 2095 ierr = VecGetSize(u,&N);CHKERRQ(ierr); 2096 for (i=0; i<N; i++) { 2097 ierr = VecCopy(u,uh);CHKERRQ(ierr); 2098 ierr = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr); 2099 for (j=-10; j<11; j++) { 2100 value = PetscSign(j)*exp(PetscAbs(j)-10.0); 2101 ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2102 ierr = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr); 2103 ierr = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr); 2104 ierr = PetscPrintf(PETSC_COMM_WORLD," j norm %D %18.16e\n",j,norm);CHKERRQ(ierr); 2105 value = -value; 2106 ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2107 } 2108 } 2109 ierr = VecDestroy(uh);CHKERRQ(ierr); 2110 ierr = VecDestroy(fh);CHKERRQ(ierr); 2111 PetscFunctionReturn(0); 2112 } 2113