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