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