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