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