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