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