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