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