1 #define PETSCSNES_DLL 2 3 #include "include/private/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_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 PetscViewerASCIIMonitor 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 = PetscViewerASCIIMonitorCreate(snes->comm,monfilename,0,&monviewer);CHKERRQ(ierr); 246 ierr = SNESMonitorSet(snes,SNESMonitorDefault,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);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 = PetscViewerASCIIMonitorCreate(snes->comm,monfilename,0,&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 = PetscViewerASCIIMonitorCreate(snes->comm,monfilename,0,&monviewer);CHKERRQ(ierr); 258 ierr = SNESMonitorSet(snes,SNESMonitorDefaultShort,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);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 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 . r - vector to store function value 807 . func - function evaluation routine 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 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 1078 if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 1079 snes->jacobian = A; 1080 } 1081 if (B) { 1082 ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 1083 if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 1084 snes->jacobian_pre = B; 1085 } 1086 PetscFunctionReturn(0); 1087 } 1088 1089 #undef __FUNCT__ 1090 #define __FUNCT__ "SNESGetJacobian" 1091 /*@C 1092 SNESGetJacobian - Returns the Jacobian matrix and optionally the user 1093 provided context for evaluating the Jacobian. 1094 1095 Not Collective, but Mat object will be parallel if SNES object is 1096 1097 Input Parameter: 1098 . snes - the nonlinear solver context 1099 1100 Output Parameters: 1101 + A - location to stash Jacobian matrix (or PETSC_NULL) 1102 . B - location to stash preconditioner matrix (or PETSC_NULL) 1103 . func - location to put Jacobian function (or PETSC_NULL) 1104 - ctx - location to stash Jacobian ctx (or PETSC_NULL) 1105 1106 Level: advanced 1107 1108 .seealso: SNESSetJacobian(), SNESComputeJacobian() 1109 @*/ 1110 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetJacobian(SNES snes,Mat *A,Mat *B,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx) 1111 { 1112 PetscFunctionBegin; 1113 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1114 if (A) *A = snes->jacobian; 1115 if (B) *B = snes->jacobian_pre; 1116 if (func) *func = snes->ops->computejacobian; 1117 if (ctx) *ctx = snes->jacP; 1118 PetscFunctionReturn(0); 1119 } 1120 1121 /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 1122 EXTERN PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*); 1123 1124 #undef __FUNCT__ 1125 #define __FUNCT__ "SNESSetUp" 1126 /*@ 1127 SNESSetUp - Sets up the internal data structures for the later use 1128 of a nonlinear solver. 1129 1130 Collective on SNES 1131 1132 Input Parameters: 1133 . snes - the SNES context 1134 1135 Notes: 1136 For basic use of the SNES solvers the user need not explicitly call 1137 SNESSetUp(), since these actions will automatically occur during 1138 the call to SNESSolve(). However, if one wishes to control this 1139 phase separately, SNESSetUp() should be called after SNESCreate() 1140 and optional routines of the form SNESSetXXX(), but before SNESSolve(). 1141 1142 Level: advanced 1143 1144 .keywords: SNES, nonlinear, setup 1145 1146 .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 1147 @*/ 1148 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUp(SNES snes) 1149 { 1150 PetscErrorCode ierr; 1151 PetscTruth flg, iseqtr; 1152 1153 PetscFunctionBegin; 1154 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1155 if (snes->setupcalled) PetscFunctionReturn(0); 1156 1157 ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr); 1158 /* 1159 This version replaces the user provided Jacobian matrix with a 1160 matrix-free version but still employs the user-provided preconditioner matrix 1161 */ 1162 if (flg) { 1163 Mat J; 1164 ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1165 ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 1166 ierr = PetscInfo(snes,"Setting default matrix-free operator routines\n");CHKERRQ(ierr); 1167 ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 1168 ierr = MatDestroy(J);CHKERRQ(ierr); 1169 } 1170 1171 #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) && !defined(PETSC_USE_MAT_SINGLE) && !defined(PETSC_USE_LONG_DOUBLE) && !defined(PETSC_USE_INT) 1172 ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr); 1173 if (flg) { 1174 Mat J; 1175 ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1176 ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 1177 ierr = MatDestroy(J);CHKERRQ(ierr); 1178 } 1179 #endif 1180 1181 ierr = PetscOptionsHasName(snes->prefix,"-snes_mf",&flg);CHKERRQ(ierr); 1182 /* 1183 This version replaces both the user-provided Jacobian and the user- 1184 provided preconditioner matrix with the default matrix free version. 1185 */ 1186 if (flg) { 1187 Mat J; 1188 KSP ksp; 1189 PC pc; 1190 1191 ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1192 ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 1193 ierr = PetscInfo(snes,"Setting default matrix-free operator and preconditioner routines;\nThat is no preconditioner is being used.\n");CHKERRQ(ierr); 1194 ierr = SNESSetJacobian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr); 1195 ierr = MatDestroy(J);CHKERRQ(ierr); 1196 1197 /* force no preconditioner */ 1198 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1199 ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 1200 ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr); 1201 if (!flg) { 1202 ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 1203 } 1204 } 1205 1206 if (!snes->vec_func && !snes->afine) { 1207 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 1208 } 1209 if (!snes->ops->computefunction && !snes->afine) { 1210 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 1211 } 1212 if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option"); 1213 if (snes->vec_func == snes->vec_sol) { 1214 SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector"); 1215 } 1216 1217 /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 1218 ierr = PetscTypeCompare((PetscObject)snes,SNESTR,&iseqtr);CHKERRQ(ierr); 1219 if (snes->ksp_ewconv && !iseqtr) { 1220 KSP ksp; 1221 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1222 ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,snes);CHKERRQ(ierr); 1223 } 1224 1225 if (snes->ops->setup) {ierr = (*snes->ops->setup)(snes);CHKERRQ(ierr);} 1226 snes->setupcalled = PETSC_TRUE; 1227 PetscFunctionReturn(0); 1228 } 1229 1230 #undef __FUNCT__ 1231 #define __FUNCT__ "SNESDestroy" 1232 /*@ 1233 SNESDestroy - Destroys the nonlinear solver context that was created 1234 with SNESCreate(). 1235 1236 Collective on SNES 1237 1238 Input Parameter: 1239 . snes - the SNES context 1240 1241 Level: beginner 1242 1243 .keywords: SNES, nonlinear, destroy 1244 1245 .seealso: SNESCreate(), SNESSolve() 1246 @*/ 1247 PetscErrorCode PETSCSNES_DLLEXPORT SNESDestroy(SNES snes) 1248 { 1249 PetscErrorCode ierr; 1250 1251 PetscFunctionBegin; 1252 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1253 if (--snes->refct > 0) PetscFunctionReturn(0); 1254 1255 /* if memory was published with AMS then destroy it */ 1256 ierr = PetscObjectDepublish(snes);CHKERRQ(ierr); 1257 1258 if (snes->ops->destroy) {ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr);} 1259 ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr); 1260 if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 1261 if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 1262 if (snes->afine) {ierr = VecDestroy(snes->afine);CHKERRQ(ierr);} 1263 ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr); 1264 if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);} 1265 ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 1266 ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr); 1267 PetscFunctionReturn(0); 1268 } 1269 1270 /* ----------- Routines to set solver parameters ---------- */ 1271 1272 #undef __FUNCT__ 1273 #define __FUNCT__ "SNESSetTolerances" 1274 /*@ 1275 SNESSetTolerances - Sets various parameters used in convergence tests. 1276 1277 Collective on SNES 1278 1279 Input Parameters: 1280 + snes - the SNES context 1281 . abstol - absolute convergence tolerance 1282 . rtol - relative convergence tolerance 1283 . stol - convergence tolerance in terms of the norm 1284 of the change in the solution between steps 1285 . maxit - maximum number of iterations 1286 - maxf - maximum number of function evaluations 1287 1288 Options Database Keys: 1289 + -snes_atol <abstol> - Sets abstol 1290 . -snes_rtol <rtol> - Sets rtol 1291 . -snes_stol <stol> - Sets stol 1292 . -snes_max_it <maxit> - Sets maxit 1293 - -snes_max_funcs <maxf> - Sets maxf 1294 1295 Notes: 1296 The default maximum number of iterations is 50. 1297 The default maximum number of function evaluations is 1000. 1298 1299 Level: intermediate 1300 1301 .keywords: SNES, nonlinear, set, convergence, tolerances 1302 1303 .seealso: SNESSetTrustRegionTolerance() 1304 @*/ 1305 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf) 1306 { 1307 PetscFunctionBegin; 1308 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1309 if (abstol != PETSC_DEFAULT) snes->abstol = abstol; 1310 if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1311 if (stol != PETSC_DEFAULT) snes->xtol = stol; 1312 if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1313 if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 1314 PetscFunctionReturn(0); 1315 } 1316 1317 #undef __FUNCT__ 1318 #define __FUNCT__ "SNESGetTolerances" 1319 /*@ 1320 SNESGetTolerances - Gets various parameters used in convergence tests. 1321 1322 Not Collective 1323 1324 Input Parameters: 1325 + snes - the SNES context 1326 . abstol - absolute convergence tolerance 1327 . rtol - relative convergence tolerance 1328 . stol - convergence tolerance in terms of the norm 1329 of the change in the solution between steps 1330 . maxit - maximum number of iterations 1331 - maxf - maximum number of function evaluations 1332 1333 Notes: 1334 The user can specify PETSC_NULL for any parameter that is not needed. 1335 1336 Level: intermediate 1337 1338 .keywords: SNES, nonlinear, get, convergence, tolerances 1339 1340 .seealso: SNESSetTolerances() 1341 @*/ 1342 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetTolerances(SNES snes,PetscReal *abstol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf) 1343 { 1344 PetscFunctionBegin; 1345 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1346 if (abstol) *abstol = snes->abstol; 1347 if (rtol) *rtol = snes->rtol; 1348 if (stol) *stol = snes->xtol; 1349 if (maxit) *maxit = snes->max_its; 1350 if (maxf) *maxf = snes->max_funcs; 1351 PetscFunctionReturn(0); 1352 } 1353 1354 #undef __FUNCT__ 1355 #define __FUNCT__ "SNESSetTrustRegionTolerance" 1356 /*@ 1357 SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 1358 1359 Collective on SNES 1360 1361 Input Parameters: 1362 + snes - the SNES context 1363 - tol - tolerance 1364 1365 Options Database Key: 1366 . -snes_trtol <tol> - Sets tol 1367 1368 Level: intermediate 1369 1370 .keywords: SNES, nonlinear, set, trust region, tolerance 1371 1372 .seealso: SNESSetTolerances() 1373 @*/ 1374 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTrustRegionTolerance(SNES snes,PetscReal tol) 1375 { 1376 PetscFunctionBegin; 1377 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1378 snes->deltatol = tol; 1379 PetscFunctionReturn(0); 1380 } 1381 1382 /* 1383 Duplicate the lg monitors for SNES from KSP; for some reason with 1384 dynamic libraries things don't work under Sun4 if we just use 1385 macros instead of functions 1386 */ 1387 #undef __FUNCT__ 1388 #define __FUNCT__ "SNESMonitorLG" 1389 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLG(SNES snes,PetscInt it,PetscReal norm,void *ctx) 1390 { 1391 PetscErrorCode ierr; 1392 1393 PetscFunctionBegin; 1394 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1395 ierr = KSPMonitorLG((KSP)snes,it,norm,ctx);CHKERRQ(ierr); 1396 PetscFunctionReturn(0); 1397 } 1398 1399 #undef __FUNCT__ 1400 #define __FUNCT__ "SNESMonitorLGCreate" 1401 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1402 { 1403 PetscErrorCode ierr; 1404 1405 PetscFunctionBegin; 1406 ierr = KSPMonitorLGCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr); 1407 PetscFunctionReturn(0); 1408 } 1409 1410 #undef __FUNCT__ 1411 #define __FUNCT__ "SNESMonitorLGDestroy" 1412 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGDestroy(PetscDrawLG draw) 1413 { 1414 PetscErrorCode ierr; 1415 1416 PetscFunctionBegin; 1417 ierr = KSPMonitorLGDestroy(draw);CHKERRQ(ierr); 1418 PetscFunctionReturn(0); 1419 } 1420 1421 /* ------------ Routines to set performance monitoring options ----------- */ 1422 1423 #undef __FUNCT__ 1424 #define __FUNCT__ "SNESMonitorSet" 1425 /*@C 1426 SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every 1427 iteration of the nonlinear solver to display the iteration's 1428 progress. 1429 1430 Collective on SNES 1431 1432 Input Parameters: 1433 + snes - the SNES context 1434 . func - monitoring routine 1435 . mctx - [optional] user-defined context for private data for the 1436 monitor routine (use PETSC_NULL if no context is desired) 1437 - monitordestroy - [optional] routine that frees monitor context 1438 (may be PETSC_NULL) 1439 1440 Calling sequence of func: 1441 $ int func(SNES snes,PetscInt its, PetscReal norm,void *mctx) 1442 1443 + snes - the SNES context 1444 . its - iteration number 1445 . norm - 2-norm function value (may be estimated) 1446 - mctx - [optional] monitoring context 1447 1448 Options Database Keys: 1449 + -snes_monitor - sets SNESMonitorDefault() 1450 . -snes_monitor_draw - sets line graph monitor, 1451 uses SNESMonitorLGCreate() 1452 _ -snes_monitor_cancel - cancels all monitors that have 1453 been hardwired into a code by 1454 calls to SNESMonitorSet(), but 1455 does not cancel those set via 1456 the options database. 1457 1458 Notes: 1459 Several different monitoring routines may be set by calling 1460 SNESMonitorSet() multiple times; all will be called in the 1461 order in which they were set. 1462 1463 Level: intermediate 1464 1465 .keywords: SNES, nonlinear, set, monitor 1466 1467 .seealso: SNESMonitorDefault(), SNESMonitorCancel() 1468 @*/ 1469 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorSet(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void*)) 1470 { 1471 PetscFunctionBegin; 1472 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1473 if (snes->numbermonitors >= MAXSNESMONITORS) { 1474 SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1475 } 1476 snes->monitor[snes->numbermonitors] = func; 1477 snes->monitordestroy[snes->numbermonitors] = monitordestroy; 1478 snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 1479 PetscFunctionReturn(0); 1480 } 1481 1482 #undef __FUNCT__ 1483 #define __FUNCT__ "SNESMonitorCancel" 1484 /*@C 1485 SNESMonitorCancel - Clears all the monitor functions for a SNES object. 1486 1487 Collective on SNES 1488 1489 Input Parameters: 1490 . snes - the SNES context 1491 1492 Options Database Key: 1493 . -snes_monitor_cancel - cancels all monitors that have been hardwired 1494 into a code by calls to SNESMonitorSet(), but does not cancel those 1495 set via the options database 1496 1497 Notes: 1498 There is no way to clear one specific monitor from a SNES object. 1499 1500 Level: intermediate 1501 1502 .keywords: SNES, nonlinear, set, monitor 1503 1504 .seealso: SNESMonitorDefault(), SNESMonitorSet() 1505 @*/ 1506 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorCancel(SNES snes) 1507 { 1508 PetscErrorCode ierr; 1509 PetscInt i; 1510 1511 PetscFunctionBegin; 1512 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1513 for (i=0; i<snes->numbermonitors; i++) { 1514 if (snes->monitordestroy[i]) { 1515 ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr); 1516 } 1517 } 1518 snes->numbermonitors = 0; 1519 PetscFunctionReturn(0); 1520 } 1521 1522 #undef __FUNCT__ 1523 #define __FUNCT__ "SNESSetConvergenceTest" 1524 /*@C 1525 SNESSetConvergenceTest - Sets the function that is to be used 1526 to test for convergence of the nonlinear iterative solution. 1527 1528 Collective on SNES 1529 1530 Input Parameters: 1531 + snes - the SNES context 1532 . func - routine to test for convergence 1533 - cctx - [optional] context for private data for the convergence routine 1534 (may be PETSC_NULL) 1535 1536 Calling sequence of func: 1537 $ PetscErrorCode func (SNES snes,PetscInt it,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx) 1538 1539 + snes - the SNES context 1540 . it - current iteration (0 is the first and is before any Newton step) 1541 . cctx - [optional] convergence context 1542 . reason - reason for convergence/divergence 1543 . xnorm - 2-norm of current iterate 1544 . gnorm - 2-norm of current step 1545 - f - 2-norm of function 1546 1547 Level: advanced 1548 1549 .keywords: SNES, nonlinear, set, convergence, test 1550 1551 .seealso: SNESConverged_LS(), SNESConverged_TR() 1552 @*/ 1553 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx) 1554 { 1555 PetscFunctionBegin; 1556 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1557 (snes)->ops->converged = func; 1558 (snes)->cnvP = cctx; 1559 PetscFunctionReturn(0); 1560 } 1561 1562 #undef __FUNCT__ 1563 #define __FUNCT__ "SNESGetConvergedReason" 1564 /*@ 1565 SNESGetConvergedReason - Gets the reason the SNES iteration was stopped. 1566 1567 Not Collective 1568 1569 Input Parameter: 1570 . snes - the SNES context 1571 1572 Output Parameter: 1573 . reason - negative value indicates diverged, positive value converged, see petscsnes.h or the 1574 manual pages for the individual convergence tests for complete lists 1575 1576 Level: intermediate 1577 1578 Notes: Can only be called after the call the SNESSolve() is complete. 1579 1580 .keywords: SNES, nonlinear, set, convergence, test 1581 1582 .seealso: SNESSetConvergenceTest(), SNESConverged_LS(), SNESConverged_TR(), SNESConvergedReason 1583 @*/ 1584 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason) 1585 { 1586 PetscFunctionBegin; 1587 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1588 PetscValidPointer(reason,2); 1589 *reason = snes->reason; 1590 PetscFunctionReturn(0); 1591 } 1592 1593 #undef __FUNCT__ 1594 #define __FUNCT__ "SNESSetConvergenceHistory" 1595 /*@ 1596 SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1597 1598 Collective on SNES 1599 1600 Input Parameters: 1601 + snes - iterative context obtained from SNESCreate() 1602 . a - array to hold history 1603 . its - integer array holds the number of linear iterations for each solve. 1604 . na - size of a and its 1605 - reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero, 1606 else it continues storing new values for new nonlinear solves after the old ones 1607 1608 Notes: 1609 If set, this array will contain the function norms computed 1610 at each step. 1611 1612 This routine is useful, e.g., when running a code for purposes 1613 of accurate performance monitoring, when no I/O should be done 1614 during the section of code that is being timed. 1615 1616 Level: intermediate 1617 1618 .keywords: SNES, set, convergence, history 1619 1620 .seealso: SNESGetConvergenceHistory() 1621 1622 @*/ 1623 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt *its,PetscInt na,PetscTruth reset) 1624 { 1625 PetscFunctionBegin; 1626 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1627 if (na) PetscValidScalarPointer(a,2); 1628 snes->conv_hist = a; 1629 snes->conv_hist_its = its; 1630 snes->conv_hist_max = na; 1631 snes->conv_hist_reset = reset; 1632 PetscFunctionReturn(0); 1633 } 1634 1635 #undef __FUNCT__ 1636 #define __FUNCT__ "SNESGetConvergenceHistory" 1637 /*@C 1638 SNESGetConvergenceHistory - Gets the array used to hold the convergence history. 1639 1640 Collective on SNES 1641 1642 Input Parameter: 1643 . snes - iterative context obtained from SNESCreate() 1644 1645 Output Parameters: 1646 . a - array to hold history 1647 . its - integer array holds the number of linear iterations (or 1648 negative if not converged) for each solve. 1649 - na - size of a and its 1650 1651 Notes: 1652 The calling sequence for this routine in Fortran is 1653 $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr) 1654 1655 This routine is useful, e.g., when running a code for purposes 1656 of accurate performance monitoring, when no I/O should be done 1657 during the section of code that is being timed. 1658 1659 Level: intermediate 1660 1661 .keywords: SNES, get, convergence, history 1662 1663 .seealso: SNESSetConvergencHistory() 1664 1665 @*/ 1666 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na) 1667 { 1668 PetscFunctionBegin; 1669 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1670 if (a) *a = snes->conv_hist; 1671 if (its) *its = snes->conv_hist_its; 1672 if (na) *na = snes->conv_hist_len; 1673 PetscFunctionReturn(0); 1674 } 1675 1676 #undef __FUNCT__ 1677 #define __FUNCT__ "SNESSetUpdate" 1678 /*@C 1679 SNESSetUpdate - Sets the general-purpose update function called 1680 at the beginning o every iteration of the nonlinear solve. Specifically 1681 it is called just before the Jacobian is "evaluated". 1682 1683 Collective on SNES 1684 1685 Input Parameters: 1686 . snes - The nonlinear solver context 1687 . func - The function 1688 1689 Calling sequence of func: 1690 . func (SNES snes, PetscInt step); 1691 1692 . step - The current step of the iteration 1693 1694 Level: intermediate 1695 1696 .keywords: SNES, update 1697 1698 .seealso SNESDefaultUpdate(), SNESSetRhsBC(), SNESSetSolutionBC() 1699 @*/ 1700 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt)) 1701 { 1702 PetscFunctionBegin; 1703 PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 1704 snes->ops->update = func; 1705 PetscFunctionReturn(0); 1706 } 1707 1708 #undef __FUNCT__ 1709 #define __FUNCT__ "SNESDefaultUpdate" 1710 /*@ 1711 SNESDefaultUpdate - The default update function which does nothing. 1712 1713 Not collective 1714 1715 Input Parameters: 1716 . snes - The nonlinear solver context 1717 . step - The current step of the iteration 1718 1719 Level: intermediate 1720 1721 .keywords: SNES, update 1722 .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultShortolutionBC() 1723 @*/ 1724 PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultUpdate(SNES snes, PetscInt step) 1725 { 1726 PetscFunctionBegin; 1727 PetscFunctionReturn(0); 1728 } 1729 1730 #undef __FUNCT__ 1731 #define __FUNCT__ "SNESScaleStep_Private" 1732 /* 1733 SNESScaleStep_Private - Scales a step so that its length is less than the 1734 positive parameter delta. 1735 1736 Input Parameters: 1737 + snes - the SNES context 1738 . y - approximate solution of linear system 1739 . fnorm - 2-norm of current function 1740 - delta - trust region size 1741 1742 Output Parameters: 1743 + gpnorm - predicted function norm at the new point, assuming local 1744 linearization. The value is zero if the step lies within the trust 1745 region, and exceeds zero otherwise. 1746 - ynorm - 2-norm of the step 1747 1748 Note: 1749 For non-trust region methods such as SNESLS, the parameter delta 1750 is set to be the maximum allowable step size. 1751 1752 .keywords: SNES, nonlinear, scale, step 1753 */ 1754 PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm) 1755 { 1756 PetscReal nrm; 1757 PetscScalar cnorm; 1758 PetscErrorCode ierr; 1759 1760 PetscFunctionBegin; 1761 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1762 PetscValidHeaderSpecific(y,VEC_COOKIE,2); 1763 PetscCheckSameComm(snes,1,y,2); 1764 1765 ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr); 1766 if (nrm > *delta) { 1767 nrm = *delta/nrm; 1768 *gpnorm = (1.0 - nrm)*(*fnorm); 1769 cnorm = nrm; 1770 ierr = VecScale(y,cnorm);CHKERRQ(ierr); 1771 *ynorm = *delta; 1772 } else { 1773 *gpnorm = 0.0; 1774 *ynorm = nrm; 1775 } 1776 PetscFunctionReturn(0); 1777 } 1778 1779 #undef __FUNCT__ 1780 #define __FUNCT__ "SNESSolve" 1781 /*@C 1782 SNESSolve - Solves a nonlinear system F(x) = b. 1783 Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX(). 1784 1785 Collective on SNES 1786 1787 Input Parameters: 1788 + snes - the SNES context 1789 . b - the constant part of the equation, or PETSC_NULL to use zero. 1790 - x - the solution vector, or PETSC_NULL if it was set with SNESSetSolution() 1791 1792 Notes: 1793 The user should initialize the vector,x, with the initial guess 1794 for the nonlinear solve prior to calling SNESSolve. In particular, 1795 to employ an initial guess of zero, the user should explicitly set 1796 this vector to zero by calling VecSet(). 1797 1798 Level: beginner 1799 1800 .keywords: SNES, nonlinear, solve 1801 1802 .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetRhs(), SNESSetSolution() 1803 @*/ 1804 PetscErrorCode PETSCSNES_DLLEXPORT SNESSolve(SNES snes,Vec b,Vec x) 1805 { 1806 PetscErrorCode ierr; 1807 PetscTruth flg; 1808 char filename[PETSC_MAX_PATH_LEN]; 1809 PetscViewer viewer; 1810 1811 PetscFunctionBegin; 1812 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1813 if (!snes->ops->solve) SETERRQ(PETSC_ERR_ORDER,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()"); 1814 1815 if (b) { 1816 ierr = SNESSetRhs(snes, b); CHKERRQ(ierr); 1817 if (!snes->vec_func) { 1818 Vec r; 1819 1820 ierr = VecDuplicate(b, &r); CHKERRQ(ierr); 1821 ierr = SNESSetFunction(snes, r, PETSC_NULL, PETSC_NULL); CHKERRQ(ierr); 1822 } 1823 } 1824 if (x) { 1825 PetscValidHeaderSpecific(x,VEC_COOKIE,3); 1826 PetscCheckSameComm(snes,1,x,3); 1827 } else { 1828 ierr = SNESGetSolution(snes, &x); CHKERRQ(ierr); 1829 if (!x) { 1830 ierr = VecDuplicate(snes->vec_func_always, &x); CHKERRQ(ierr); 1831 } 1832 } 1833 snes->vec_sol = snes->vec_sol_always = x; 1834 if (!snes->setupcalled) { 1835 ierr = SNESSetUp(snes);CHKERRQ(ierr); 1836 } 1837 if (snes->conv_hist_reset) snes->conv_hist_len = 0; 1838 ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 1839 snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0; 1840 1841 ierr = PetscExceptionTry1((*(snes)->ops->solve)(snes),PETSC_ERR_ARG_DOMAIN); 1842 if (PetscExceptionValue(ierr)) { 1843 /* this means that a caller above me has also tryed this exception so I don't handle it here, pass it up */ 1844 PetscErrorCode pierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(pierr); 1845 } else if (PetscExceptionCaught(ierr,PETSC_ERR_ARG_DOMAIN)) { 1846 /* translate exception into SNES not converged reason */ 1847 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1848 ierr = 0; 1849 } 1850 CHKERRQ(ierr); 1851 1852 ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 1853 ierr = PetscOptionsGetString(snes->prefix,"-snes_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 1854 if (flg && !PetscPreLoadingOn) { 1855 ierr = PetscViewerASCIIOpen(snes->comm,filename,&viewer);CHKERRQ(ierr); 1856 ierr = SNESView(snes,viewer);CHKERRQ(ierr); 1857 ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr); 1858 } 1859 1860 ierr = PetscOptionsHasName(snes->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr); 1861 if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); } 1862 if (snes->printreason) { 1863 if (snes->reason > 0) { 1864 ierr = PetscPrintf(snes->comm,"Nonlinear solve converged due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr); 1865 } else { 1866 ierr = PetscPrintf(snes->comm,"Nonlinear solve did not converge due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr); 1867 } 1868 } 1869 1870 PetscFunctionReturn(0); 1871 } 1872 1873 /* --------- Internal routines for SNES Package --------- */ 1874 1875 #undef __FUNCT__ 1876 #define __FUNCT__ "SNESSetType" 1877 /*@C 1878 SNESSetType - Sets the method for the nonlinear solver. 1879 1880 Collective on SNES 1881 1882 Input Parameters: 1883 + snes - the SNES context 1884 - type - a known method 1885 1886 Options Database Key: 1887 . -snes_type <type> - Sets the method; use -help for a list 1888 of available methods (for instance, ls or tr) 1889 1890 Notes: 1891 See "petsc/include/petscsnes.h" for available methods (for instance) 1892 + SNESLS - Newton's method with line search 1893 (systems of nonlinear equations) 1894 . SNESTR - Newton's method with trust region 1895 (systems of nonlinear equations) 1896 1897 Normally, it is best to use the SNESSetFromOptions() command and then 1898 set the SNES solver type from the options database rather than by using 1899 this routine. Using the options database provides the user with 1900 maximum flexibility in evaluating the many nonlinear solvers. 1901 The SNESSetType() routine is provided for those situations where it 1902 is necessary to set the nonlinear solver independently of the command 1903 line or options database. This might be the case, for example, when 1904 the choice of solver changes during the execution of the program, 1905 and the user's application is taking responsibility for choosing the 1906 appropriate method. 1907 1908 Level: intermediate 1909 1910 .keywords: SNES, set, type 1911 1912 .seealso: SNESType, SNESCreate() 1913 1914 @*/ 1915 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetType(SNES snes,SNESType type) 1916 { 1917 PetscErrorCode ierr,(*r)(SNES); 1918 PetscTruth match; 1919 1920 PetscFunctionBegin; 1921 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1922 PetscValidCharPointer(type,2); 1923 1924 ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr); 1925 if (match) PetscFunctionReturn(0); 1926 1927 if (snes->setupcalled) { 1928 snes->setupcalled = PETSC_FALSE; 1929 ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr); 1930 snes->data = 0; 1931 } 1932 1933 /* Get the function pointers for the iterative method requested */ 1934 if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 1935 ierr = PetscFListFind(snes->comm,SNESList,type,(void (**)(void)) &r);CHKERRQ(ierr); 1936 if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type); 1937 ierr = PetscFree(snes->data);CHKERRQ(ierr); 1938 snes->data = 0; 1939 ierr = (*r)(snes);CHKERRQ(ierr); 1940 ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr); 1941 PetscFunctionReturn(0); 1942 } 1943 1944 1945 /* --------------------------------------------------------------------- */ 1946 #undef __FUNCT__ 1947 #define __FUNCT__ "SNESRegisterDestroy" 1948 /*@ 1949 SNESRegisterDestroy - Frees the list of nonlinear solvers that were 1950 registered by SNESRegisterDynamic(). 1951 1952 Not Collective 1953 1954 Level: advanced 1955 1956 .keywords: SNES, nonlinear, register, destroy 1957 1958 .seealso: SNESRegisterAll(), SNESRegisterAll() 1959 @*/ 1960 PetscErrorCode PETSCSNES_DLLEXPORT SNESRegisterDestroy(void) 1961 { 1962 PetscErrorCode ierr; 1963 1964 PetscFunctionBegin; 1965 if (SNESList) { 1966 ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr); 1967 SNESList = 0; 1968 } 1969 SNESRegisterAllCalled = PETSC_FALSE; 1970 PetscFunctionReturn(0); 1971 } 1972 1973 #undef __FUNCT__ 1974 #define __FUNCT__ "SNESGetType" 1975 /*@C 1976 SNESGetType - Gets the SNES method type and name (as a string). 1977 1978 Not Collective 1979 1980 Input Parameter: 1981 . snes - nonlinear solver context 1982 1983 Output Parameter: 1984 . type - SNES method (a character string) 1985 1986 Level: intermediate 1987 1988 .keywords: SNES, nonlinear, get, type, name 1989 @*/ 1990 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetType(SNES snes,SNESType *type) 1991 { 1992 PetscFunctionBegin; 1993 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1994 PetscValidPointer(type,2); 1995 *type = snes->type_name; 1996 PetscFunctionReturn(0); 1997 } 1998 1999 #undef __FUNCT__ 2000 #define __FUNCT__ "SNESGetSolution" 2001 /*@ 2002 SNESGetSolution - Returns the vector where the approximate solution is 2003 stored. 2004 2005 Not Collective, but Vec is parallel if SNES is parallel 2006 2007 Input Parameter: 2008 . snes - the SNES context 2009 2010 Output Parameter: 2011 . x - the solution 2012 2013 Level: intermediate 2014 2015 .keywords: SNES, nonlinear, get, solution 2016 2017 .seealso: SNESSetSolution(), SNESGetFunction(), SNESGetSolutionUpdate() 2018 @*/ 2019 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolution(SNES snes,Vec *x) 2020 { 2021 PetscFunctionBegin; 2022 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2023 PetscValidPointer(x,2); 2024 *x = snes->vec_sol_always; 2025 PetscFunctionReturn(0); 2026 } 2027 2028 #undef __FUNCT__ 2029 #define __FUNCT__ "SNESSetSolution" 2030 /*@ 2031 SNESSetSolution - Sets the vector where the approximate solution is stored. 2032 2033 Not Collective, but Vec is parallel if SNES is parallel 2034 2035 Input Parameters: 2036 + snes - the SNES context 2037 - x - the solution 2038 2039 Output Parameter: 2040 2041 Level: intermediate 2042 2043 Notes: this is not normally used, rather one simply calls SNESSolve() with 2044 the appropriate solution vector. 2045 2046 .keywords: SNES, nonlinear, set, solution 2047 2048 .seealso: SNESGetSolution(), SNESGetFunction(), SNESGetSolutionUpdate() 2049 @*/ 2050 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetSolution(SNES snes,Vec x) 2051 { 2052 PetscFunctionBegin; 2053 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2054 PetscValidHeaderSpecific(x,VEC_COOKIE,2); 2055 PetscCheckSameComm(snes,1,x,2); 2056 snes->vec_sol_always = x; 2057 PetscFunctionReturn(0); 2058 } 2059 2060 #undef __FUNCT__ 2061 #define __FUNCT__ "SNESGetSolutionUpdate" 2062 /*@ 2063 SNESGetSolutionUpdate - Returns the vector where the solution update is 2064 stored. 2065 2066 Not Collective, but Vec is parallel if SNES is parallel 2067 2068 Input Parameter: 2069 . snes - the SNES context 2070 2071 Output Parameter: 2072 . x - the solution update 2073 2074 Level: advanced 2075 2076 .keywords: SNES, nonlinear, get, solution, update 2077 2078 .seealso: SNESGetSolution(), SNESGetFunction 2079 @*/ 2080 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolutionUpdate(SNES snes,Vec *x) 2081 { 2082 PetscFunctionBegin; 2083 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2084 PetscValidPointer(x,2); 2085 *x = snes->vec_sol_update_always; 2086 PetscFunctionReturn(0); 2087 } 2088 2089 #undef __FUNCT__ 2090 #define __FUNCT__ "SNESGetFunction" 2091 /*@C 2092 SNESGetFunction - Returns the vector where the function is stored. 2093 2094 Not Collective, but Vec is parallel if SNES is parallel 2095 2096 Input Parameter: 2097 . snes - the SNES context 2098 2099 Output Parameter: 2100 + r - the function (or PETSC_NULL) 2101 . func - the function (or PETSC_NULL) 2102 - ctx - the function context (or PETSC_NULL) 2103 2104 Level: advanced 2105 2106 .keywords: SNES, nonlinear, get, function 2107 2108 .seealso: SNESSetFunction(), SNESGetSolution() 2109 @*/ 2110 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx) 2111 { 2112 PetscFunctionBegin; 2113 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2114 if (r) *r = snes->vec_func_always; 2115 if (func) *func = snes->ops->computefunction; 2116 if (ctx) *ctx = snes->funP; 2117 PetscFunctionReturn(0); 2118 } 2119 2120 #undef __FUNCT__ 2121 #define __FUNCT__ "SNESSetOptionsPrefix" 2122 /*@C 2123 SNESSetOptionsPrefix - Sets the prefix used for searching for all 2124 SNES options in the database. 2125 2126 Collective on SNES 2127 2128 Input Parameter: 2129 + snes - the SNES context 2130 - prefix - the prefix to prepend to all option names 2131 2132 Notes: 2133 A hyphen (-) must NOT be given at the beginning of the prefix name. 2134 The first character of all runtime options is AUTOMATICALLY the hyphen. 2135 2136 Level: advanced 2137 2138 .keywords: SNES, set, options, prefix, database 2139 2140 .seealso: SNESSetFromOptions() 2141 @*/ 2142 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetOptionsPrefix(SNES snes,const char prefix[]) 2143 { 2144 PetscErrorCode ierr; 2145 2146 PetscFunctionBegin; 2147 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2148 ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2149 ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 2150 PetscFunctionReturn(0); 2151 } 2152 2153 #undef __FUNCT__ 2154 #define __FUNCT__ "SNESAppendOptionsPrefix" 2155 /*@C 2156 SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 2157 SNES options in the database. 2158 2159 Collective on SNES 2160 2161 Input Parameters: 2162 + snes - the SNES context 2163 - prefix - the prefix to prepend to all option names 2164 2165 Notes: 2166 A hyphen (-) must NOT be given at the beginning of the prefix name. 2167 The first character of all runtime options is AUTOMATICALLY the hyphen. 2168 2169 Level: advanced 2170 2171 .keywords: SNES, append, options, prefix, database 2172 2173 .seealso: SNESGetOptionsPrefix() 2174 @*/ 2175 PetscErrorCode PETSCSNES_DLLEXPORT SNESAppendOptionsPrefix(SNES snes,const char prefix[]) 2176 { 2177 PetscErrorCode ierr; 2178 2179 PetscFunctionBegin; 2180 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2181 ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2182 ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 2183 PetscFunctionReturn(0); 2184 } 2185 2186 #undef __FUNCT__ 2187 #define __FUNCT__ "SNESGetOptionsPrefix" 2188 /*@C 2189 SNESGetOptionsPrefix - Sets the prefix used for searching for all 2190 SNES options in the database. 2191 2192 Not Collective 2193 2194 Input Parameter: 2195 . snes - the SNES context 2196 2197 Output Parameter: 2198 . prefix - pointer to the prefix string used 2199 2200 Notes: On the fortran side, the user should pass in a string 'prifix' of 2201 sufficient length to hold the prefix. 2202 2203 Level: advanced 2204 2205 .keywords: SNES, get, options, prefix, database 2206 2207 .seealso: SNESAppendOptionsPrefix() 2208 @*/ 2209 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetOptionsPrefix(SNES snes,const char *prefix[]) 2210 { 2211 PetscErrorCode ierr; 2212 2213 PetscFunctionBegin; 2214 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2215 ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2216 PetscFunctionReturn(0); 2217 } 2218 2219 2220 #undef __FUNCT__ 2221 #define __FUNCT__ "SNESRegister" 2222 /*@C 2223 SNESRegister - See SNESRegisterDynamic() 2224 2225 Level: advanced 2226 @*/ 2227 PetscErrorCode PETSCSNES_DLLEXPORT SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES)) 2228 { 2229 char fullname[PETSC_MAX_PATH_LEN]; 2230 PetscErrorCode ierr; 2231 2232 PetscFunctionBegin; 2233 ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 2234 ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 2235 PetscFunctionReturn(0); 2236 } 2237 2238 #undef __FUNCT__ 2239 #define __FUNCT__ "SNESTestLocalMin" 2240 PetscErrorCode PETSCSNES_DLLEXPORT SNESTestLocalMin(SNES snes) 2241 { 2242 PetscErrorCode ierr; 2243 PetscInt N,i,j; 2244 Vec u,uh,fh; 2245 PetscScalar value; 2246 PetscReal norm; 2247 2248 PetscFunctionBegin; 2249 ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr); 2250 ierr = VecDuplicate(u,&uh);CHKERRQ(ierr); 2251 ierr = VecDuplicate(u,&fh);CHKERRQ(ierr); 2252 2253 /* currently only works for sequential */ 2254 ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n"); 2255 ierr = VecGetSize(u,&N);CHKERRQ(ierr); 2256 for (i=0; i<N; i++) { 2257 ierr = VecCopy(u,uh);CHKERRQ(ierr); 2258 ierr = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr); 2259 for (j=-10; j<11; j++) { 2260 value = PetscSign(j)*exp(PetscAbs(j)-10.0); 2261 ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2262 ierr = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr); 2263 ierr = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr); 2264 ierr = PetscPrintf(PETSC_COMM_WORLD," j norm %D %18.16e\n",j,norm);CHKERRQ(ierr); 2265 value = -value; 2266 ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2267 } 2268 } 2269 ierr = VecDestroy(uh);CHKERRQ(ierr); 2270 ierr = VecDestroy(fh);CHKERRQ(ierr); 2271 PetscFunctionReturn(0); 2272 } 2273