1 /*$Id: snes.c,v 1.221 2000/09/22 20:45:56 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","SNES");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 167 if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 168 ierr = OptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr); 169 if (flg) { 170 ierr = SNESSetType(snes,type);CHKERRQ(ierr); 171 } else if (!snes->type_name) { 172 ierr = SNESSetType(snes,deft);CHKERRQ(ierr); 173 } 174 175 ierr = OptionsDouble("-snes_stol","Stop if step length less then","SNESSetTolerances",snes->xtol,&snes->xtol,0);CHKERRQ(ierr); 176 ierr = OptionsDouble("-snes_atol","Stop if function norm less then","SNESSetTolerances",snes->atol,&snes->atol,0);CHKERRQ(ierr); 177 178 ierr = OptionsDouble("-snes_rtol","Stop if decrease in function norm less then","SNESSetTolerances",snes->rtol,&snes->rtol,0);CHKERRQ(ierr); 179 ierr = OptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,PETSC_NULL);CHKERRQ(ierr); 180 ierr = OptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr); 181 ierr = OptionsDouble("-snes_fmin","Minimization function tolerance","SNESSetMinimizationFunctionTolerance",snes->fmin,&snes->fmin,0);CHKERRQ(ierr); 182 183 ierr = OptionsName("-snes_ksp_ew_conv","Use Eisentat-Walker linear system convergence test","SNES_KSP_SetParametersEW",&snes->ksp_ewconv);CHKERRQ(ierr); 184 185 ierr = OptionsInt("-snes_ksp_ew_version","Version 1 or 2","SNES_KSP_SetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr); 186 ierr = OptionsDouble("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNES_KSP_SetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr); 187 ierr = OptionsDouble("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNES_KSP_SetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr); 188 ierr = OptionsDouble("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNES_KSP_SetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr); 189 ierr = OptionsDouble("-snes_ksp_ew_alpha","1 < alpha <= 2","SNES_KSP_SetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr); 190 ierr = OptionsDouble("-snes_ksp_ew_alpha2","alpha2","SNES_KSP_SetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr); 191 ierr = OptionsDouble("-snes_ksp_ew_threshold","0 < threshold < 1","SNES_KSP_SetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr); 192 193 ierr = OptionsName("-snes_no_convergence_test","Don't test for convergence","None",&flg);CHKERRQ(ierr); 194 if (flg) {snes->converged = 0;} 195 ierr = OptionsName("-snes_cancelmonitors","Remove all monitors","SNESClearMonitor",&flg);CHKERRQ(ierr); 196 if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);} 197 ierr = OptionsName("-snes_monitor","Monitor norm of function","SNESDefaultMonitor",&flg);CHKERRQ(ierr); 198 if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0,0);CHKERRQ(ierr);} 199 ierr = OptionsName("-snes_smonitor","Monitor norm of function (fewer digits)","SNESDefaultSMonitor",&flg);CHKERRQ(ierr); 200 if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0,0);CHKERRQ(ierr);} 201 ierr = OptionsName("-snes_vecmonitor","Plot solution at each iteration","SNESVecViewMonitor",&flg);CHKERRQ(ierr); 202 if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewMonitor,0,0);CHKERRQ(ierr);} 203 ierr = OptionsName("-snes_vecmonitor_update","Plot correction at each iteration","SNESVecViewUpdateMonitor",&flg);CHKERRQ(ierr); 204 if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewUpdateMonitor,0,0);CHKERRQ(ierr);} 205 ierr = OptionsName("-snes_xmonitor","Plot function norm at each iteration","SNESLGMonitor",&flg);CHKERRQ(ierr); 206 if (flg) {ierr = SNESSetMonitor(snes,SNESLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);} 207 208 ierr = OptionsName("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",&flg);CHKERRQ(ierr); 209 if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) { 210 ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr); 211 PLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Jacobian matrix\n"); 212 } else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 213 ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeHessian,snes->funP);CHKERRQ(ierr); 214 PLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Hessian matrix\n"); 215 } 216 217 if (snes->setfromoptions) { 218 ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr); 219 } 220 221 ierr = OptionsEnd();CHKERRQ(ierr); 222 223 ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr); 224 ierr = SLESSetFromOptions(sles);CHKERRQ(ierr); 225 226 PetscFunctionReturn(0); 227 } 228 229 230 #undef __FUNC__ 231 #define __FUNC__ /*<a name=""></a>*/"SNESSetApplicationContext" 232 /*@ 233 SNESSetApplicationContext - Sets the optional user-defined context for 234 the nonlinear solvers. 235 236 Collective on SNES 237 238 Input Parameters: 239 + snes - the SNES context 240 - usrP - optional user context 241 242 Level: intermediate 243 244 .keywords: SNES, nonlinear, set, application, context 245 246 .seealso: SNESGetApplicationContext() 247 @*/ 248 int SNESSetApplicationContext(SNES snes,void *usrP) 249 { 250 PetscFunctionBegin; 251 PetscValidHeaderSpecific(snes,SNES_COOKIE); 252 snes->user = usrP; 253 PetscFunctionReturn(0); 254 } 255 256 #undef __FUNC__ 257 #define __FUNC__ /*<a name=""></a>*/"SNESGetApplicationContext" 258 /*@C 259 SNESGetApplicationContext - Gets the user-defined context for the 260 nonlinear solvers. 261 262 Not Collective 263 264 Input Parameter: 265 . snes - SNES context 266 267 Output Parameter: 268 . usrP - user context 269 270 Level: intermediate 271 272 .keywords: SNES, nonlinear, get, application, context 273 274 .seealso: SNESSetApplicationContext() 275 @*/ 276 int SNESGetApplicationContext(SNES snes,void **usrP) 277 { 278 PetscFunctionBegin; 279 PetscValidHeaderSpecific(snes,SNES_COOKIE); 280 *usrP = snes->user; 281 PetscFunctionReturn(0); 282 } 283 284 #undef __FUNC__ 285 #define __FUNC__ /*<a name=""></a>*/"SNESGetIterationNumber" 286 /*@ 287 SNESGetIterationNumber - Gets the number of nonlinear iterations completed 288 at this time. 289 290 Not Collective 291 292 Input Parameter: 293 . snes - SNES context 294 295 Output Parameter: 296 . iter - iteration number 297 298 Notes: 299 For example, during the computation of iteration 2 this would return 1. 300 301 This is useful for using lagged Jacobians (where one does not recompute the 302 Jacobian at each SNES iteration). For example, the code 303 .vb 304 ierr = SNESGetIterationNumber(snes,&it); 305 if (!(it % 2)) { 306 [compute Jacobian here] 307 } 308 .ve 309 can be used in your ComputeJacobian() function to cause the Jacobian to be 310 recomputed every second SNES iteration. 311 312 Level: intermediate 313 314 .keywords: SNES, nonlinear, get, iteration, number 315 @*/ 316 int SNESGetIterationNumber(SNES snes,int* iter) 317 { 318 PetscFunctionBegin; 319 PetscValidHeaderSpecific(snes,SNES_COOKIE); 320 PetscValidIntPointer(iter); 321 *iter = snes->iter; 322 PetscFunctionReturn(0); 323 } 324 325 #undef __FUNC__ 326 #define __FUNC__ /*<a name=""></a>*/"SNESGetFunctionNorm" 327 /*@ 328 SNESGetFunctionNorm - Gets the norm of the current function that was set 329 with SNESSSetFunction(). 330 331 Collective on SNES 332 333 Input Parameter: 334 . snes - SNES context 335 336 Output Parameter: 337 . fnorm - 2-norm of function 338 339 Note: 340 SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only. 341 A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is 342 SNESGetGradientNorm(). 343 344 Level: intermediate 345 346 .keywords: SNES, nonlinear, get, function, norm 347 348 .seealso: SNESGetFunction() 349 @*/ 350 int SNESGetFunctionNorm(SNES snes,Scalar *fnorm) 351 { 352 PetscFunctionBegin; 353 PetscValidHeaderSpecific(snes,SNES_COOKIE); 354 PetscValidScalarPointer(fnorm); 355 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 356 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"For SNES_NONLINEAR_EQUATIONS only"); 357 } 358 *fnorm = snes->norm; 359 PetscFunctionReturn(0); 360 } 361 362 #undef __FUNC__ 363 #define __FUNC__ /*<a name=""></a>*/"SNESGetGradientNorm" 364 /*@ 365 SNESGetGradientNorm - Gets the norm of the current gradient that was set 366 with SNESSSetGradient(). 367 368 Collective on SNES 369 370 Input Parameter: 371 . snes - SNES context 372 373 Output Parameter: 374 . fnorm - 2-norm of gradient 375 376 Note: 377 SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION 378 methods only. A related routine for SNES_NONLINEAR_EQUATIONS methods 379 is SNESGetFunctionNorm(). 380 381 Level: intermediate 382 383 .keywords: SNES, nonlinear, get, gradient, norm 384 385 .seelso: SNESSetGradient() 386 @*/ 387 int SNESGetGradientNorm(SNES snes,Scalar *gnorm) 388 { 389 PetscFunctionBegin; 390 PetscValidHeaderSpecific(snes,SNES_COOKIE); 391 PetscValidScalarPointer(gnorm); 392 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 393 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 394 } 395 *gnorm = snes->norm; 396 PetscFunctionReturn(0); 397 } 398 399 #undef __FUNC__ 400 #define __FUNC__ /*<a name=""></a>*/"SNESGetNumberUnsuccessfulSteps" 401 /*@ 402 SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps 403 attempted by the nonlinear solver. 404 405 Not Collective 406 407 Input Parameter: 408 . snes - SNES context 409 410 Output Parameter: 411 . nfails - number of unsuccessful steps attempted 412 413 Notes: 414 This counter is reset to zero for each successive call to SNESSolve(). 415 416 Level: intermediate 417 418 .keywords: SNES, nonlinear, get, number, unsuccessful, steps 419 @*/ 420 int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails) 421 { 422 PetscFunctionBegin; 423 PetscValidHeaderSpecific(snes,SNES_COOKIE); 424 PetscValidIntPointer(nfails); 425 *nfails = snes->nfailures; 426 PetscFunctionReturn(0); 427 } 428 429 #undef __FUNC__ 430 #define __FUNC__ /*<a name=""></a>*/"SNESGetNumberLinearIterations" 431 /*@ 432 SNESGetNumberLinearIterations - Gets the total number of linear iterations 433 used by the nonlinear solver. 434 435 Not Collective 436 437 Input Parameter: 438 . snes - SNES context 439 440 Output Parameter: 441 . lits - number of linear iterations 442 443 Notes: 444 This counter is reset to zero for each successive call to SNESSolve(). 445 446 Level: intermediate 447 448 .keywords: SNES, nonlinear, get, number, linear, iterations 449 @*/ 450 int SNESGetNumberLinearIterations(SNES snes,int* lits) 451 { 452 PetscFunctionBegin; 453 PetscValidHeaderSpecific(snes,SNES_COOKIE); 454 PetscValidIntPointer(lits); 455 *lits = snes->linear_its; 456 PetscFunctionReturn(0); 457 } 458 459 #undef __FUNC__ 460 #define __FUNC__ /*<a name=""></a>*/"SNESGetSLES" 461 /*@C 462 SNESGetSLES - Returns the SLES context for a SNES solver. 463 464 Not Collective, but if SNES object is parallel, then SLES object is parallel 465 466 Input Parameter: 467 . snes - the SNES context 468 469 Output Parameter: 470 . sles - the SLES context 471 472 Notes: 473 The user can then directly manipulate the SLES context to set various 474 options, etc. Likewise, the user can then extract and manipulate the 475 KSP and PC contexts as well. 476 477 Level: beginner 478 479 .keywords: SNES, nonlinear, get, SLES, context 480 481 .seealso: SLESGetPC(), SLESGetKSP() 482 @*/ 483 int SNESGetSLES(SNES snes,SLES *sles) 484 { 485 PetscFunctionBegin; 486 PetscValidHeaderSpecific(snes,SNES_COOKIE); 487 *sles = snes->sles; 488 PetscFunctionReturn(0); 489 } 490 491 #undef __FUNC__ 492 #define __FUNC__ /*<a name=""></a>*/"SNESPublish_Petsc" 493 static int SNESPublish_Petsc(PetscObject obj) 494 { 495 #if defined(PETSC_HAVE_AMS) 496 SNES v = (SNES) obj; 497 int ierr; 498 #endif 499 500 PetscFunctionBegin; 501 502 #if defined(PETSC_HAVE_AMS) 503 /* if it is already published then return */ 504 if (v->amem >=0) PetscFunctionReturn(0); 505 506 ierr = PetscObjectPublishBaseBegin(obj);CHKERRQ(ierr); 507 ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Iteration",&v->iter,1,AMS_INT,AMS_READ, 508 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr); 509 ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Residual",&v->norm,1,AMS_DOUBLE,AMS_READ, 510 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr); 511 ierr = PetscObjectPublishBaseEnd(obj);CHKERRQ(ierr); 512 #endif 513 PetscFunctionReturn(0); 514 } 515 516 /* -----------------------------------------------------------*/ 517 #undef __FUNC__ 518 #define __FUNC__ /*<a name=""></a>*/"SNESCreate" 519 /*@C 520 SNESCreate - Creates a nonlinear solver context. 521 522 Collective on MPI_Comm 523 524 Input Parameters: 525 + comm - MPI communicator 526 - type - type of method, either 527 SNES_NONLINEAR_EQUATIONS (for systems of nonlinear equations) 528 or SNES_UNCONSTRAINED_MINIMIZATION (for unconstrained minimization) 529 530 Output Parameter: 531 . outsnes - the new SNES context 532 533 Options Database Keys: 534 + -snes_mf - Activates default matrix-free Jacobian-vector products, 535 and no preconditioning matrix 536 . -snes_mf_operator - Activates default matrix-free Jacobian-vector 537 products, and a user-provided preconditioning matrix 538 as set by SNESSetJacobian() 539 - -snes_fd - Uses (slow!) finite differences to compute Jacobian 540 541 Level: beginner 542 543 .keywords: SNES, nonlinear, create, context 544 545 .seealso: SNESSolve(), SNESDestroy() 546 @*/ 547 int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes) 548 { 549 int ierr; 550 SNES snes; 551 SNES_KSP_EW_ConvCtx *kctx; 552 553 PetscFunctionBegin; 554 *outsnes = 0; 555 if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS){ 556 SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"incorrect method type"); 557 } 558 PetscHeaderCreate(snes,_p_SNES,int,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView); 559 PLogObjectCreate(snes); 560 snes->bops->publish = SNESPublish_Petsc; 561 snes->max_its = 50; 562 snes->max_funcs = 10000; 563 snes->norm = 0.0; 564 if (type == SNES_UNCONSTRAINED_MINIMIZATION) { 565 snes->rtol = 1.e-8; 566 snes->ttol = 0.0; 567 snes->atol = 1.e-10; 568 } else { 569 snes->rtol = 1.e-8; 570 snes->ttol = 0.0; 571 snes->atol = 1.e-50; 572 } 573 snes->xtol = 1.e-8; 574 snes->trunctol = 1.e-12; /* no longer used */ 575 snes->nfuncs = 0; 576 snes->nfailures = 0; 577 snes->linear_its = 0; 578 snes->numbermonitors = 0; 579 snes->data = 0; 580 snes->view = 0; 581 snes->computeumfunction = 0; 582 snes->umfunP = 0; 583 snes->fc = 0; 584 snes->deltatol = 1.e-12; 585 snes->fmin = -1.e30; 586 snes->method_class = type; 587 snes->set_method_called = 0; 588 snes->setupcalled = 0; 589 snes->ksp_ewconv = PETSC_FALSE; 590 snes->vwork = 0; 591 snes->nwork = 0; 592 snes->conv_hist_len = 0; 593 snes->conv_hist_max = 0; 594 snes->conv_hist = PETSC_NULL; 595 snes->conv_hist_its = PETSC_NULL; 596 snes->conv_hist_reset = PETSC_TRUE; 597 snes->reason = SNES_CONVERGED_ITERATING; 598 599 /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 600 kctx = PetscNew(SNES_KSP_EW_ConvCtx);CHKPTRQ(kctx); 601 PLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx)); 602 snes->kspconvctx = (void*)kctx; 603 kctx->version = 2; 604 kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 605 this was too large for some test cases */ 606 kctx->rtol_last = 0; 607 kctx->rtol_max = .9; 608 kctx->gamma = 1.0; 609 kctx->alpha2 = .5*(1.0 + sqrt(5.0)); 610 kctx->alpha = kctx->alpha2; 611 kctx->threshold = .1; 612 kctx->lresid_last = 0; 613 kctx->norm_last = 0; 614 615 ierr = SLESCreate(comm,&snes->sles);CHKERRQ(ierr); 616 PLogObjectParent(snes,snes->sles) 617 618 *outsnes = snes; 619 ierr = PetscPublishAll(snes);CHKERRQ(ierr); 620 PetscFunctionReturn(0); 621 } 622 623 /* --------------------------------------------------------------- */ 624 #undef __FUNC__ 625 #define __FUNC__ /*<a name=""></a>*/"SNESSetFunction" 626 /*@C 627 SNESSetFunction - Sets the function evaluation routine and function 628 vector for use by the SNES routines in solving systems of nonlinear 629 equations. 630 631 Collective on SNES 632 633 Input Parameters: 634 + snes - the SNES context 635 . func - function evaluation routine 636 . r - vector to store function value 637 - ctx - [optional] user-defined context for private data for the 638 function evaluation routine (may be PETSC_NULL) 639 640 Calling sequence of func: 641 $ func (SNES snes,Vec x,Vec f,void *ctx); 642 643 . f - function vector 644 - ctx - optional user-defined function context 645 646 Notes: 647 The Newton-like methods typically solve linear systems of the form 648 $ f'(x) x = -f(x), 649 where f'(x) denotes the Jacobian matrix and f(x) is the function. 650 651 SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 652 Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 653 SNESSetMinimizationFunction() and SNESSetGradient(); 654 655 Level: beginner 656 657 .keywords: SNES, nonlinear, set, function 658 659 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 660 @*/ 661 int SNESSetFunction(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx) 662 { 663 PetscFunctionBegin; 664 PetscValidHeaderSpecific(snes,SNES_COOKIE); 665 PetscValidHeaderSpecific(r,VEC_COOKIE); 666 PetscCheckSameComm(snes,r); 667 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 668 SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only"); 669 } 670 671 snes->computefunction = func; 672 snes->vec_func = snes->vec_func_always = r; 673 snes->funP = ctx; 674 PetscFunctionReturn(0); 675 } 676 677 #undef __FUNC__ 678 #define __FUNC__ /*<a name=""></a>*/"SNESComputeFunction" 679 /*@ 680 SNESComputeFunction - Calls the function that has been set with 681 SNESSetFunction(). 682 683 Collective on SNES 684 685 Input Parameters: 686 + snes - the SNES context 687 - x - input vector 688 689 Output Parameter: 690 . y - function vector, as set by SNESSetFunction() 691 692 Notes: 693 SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 694 Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 695 SNESComputeMinimizationFunction() and SNESComputeGradient(); 696 697 SNESComputeFunction() is typically used within nonlinear solvers 698 implementations, so most users would not generally call this routine 699 themselves. 700 701 Level: developer 702 703 .keywords: SNES, nonlinear, compute, function 704 705 .seealso: SNESSetFunction(), SNESGetFunction() 706 @*/ 707 int SNESComputeFunction(SNES snes,Vec x,Vec y) 708 { 709 int ierr; 710 711 PetscFunctionBegin; 712 PetscValidHeaderSpecific(snes,SNES_COOKIE); 713 PetscValidHeaderSpecific(x,VEC_COOKIE); 714 PetscValidHeaderSpecific(y,VEC_COOKIE); 715 PetscCheckSameComm(snes,x); 716 PetscCheckSameComm(snes,y); 717 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 718 SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only"); 719 } 720 721 ierr = PLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 722 PetscStackPush("SNES user function"); 723 ierr = (*snes->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr); 724 PetscStackPop; 725 snes->nfuncs++; 726 ierr = PLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 727 PetscFunctionReturn(0); 728 } 729 730 #undef __FUNC__ 731 #define __FUNC__ /*<a name=""></a>*/"SNESSetMinimizationFunction" 732 /*@C 733 SNESSetMinimizationFunction - Sets the function evaluation routine for 734 unconstrained minimization. 735 736 Collective on SNES 737 738 Input Parameters: 739 + snes - the SNES context 740 . func - function evaluation routine 741 - ctx - [optional] user-defined context for private data for the 742 function evaluation routine (may be PETSC_NULL) 743 744 Calling sequence of func: 745 $ func (SNES snes,Vec x,PetscReal *f,void *ctx); 746 747 + x - input vector 748 . f - function 749 - ctx - [optional] user-defined function context 750 751 Level: beginner 752 753 Notes: 754 SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 755 methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 756 SNESSetFunction(). 757 758 .keywords: SNES, nonlinear, set, minimization, function 759 760 .seealso: SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(), 761 SNESSetHessian(), SNESSetGradient() 762 @*/ 763 int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,PetscReal*,void*),void *ctx) 764 { 765 PetscFunctionBegin; 766 PetscValidHeaderSpecific(snes,SNES_COOKIE); 767 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 768 SETERRQ(PETSC_ERR_ARG_WRONG,"Only for SNES_UNCONSTRAINED_MINIMIZATION"); 769 } 770 snes->computeumfunction = func; 771 snes->umfunP = ctx; 772 PetscFunctionReturn(0); 773 } 774 775 #undef __FUNC__ 776 #define __FUNC__ /*<a name=""></a>*/"SNESComputeMinimizationFunction" 777 /*@ 778 SNESComputeMinimizationFunction - Computes the function that has been 779 set with SNESSetMinimizationFunction(). 780 781 Collective on SNES 782 783 Input Parameters: 784 + snes - the SNES context 785 - x - input vector 786 787 Output Parameter: 788 . y - function value 789 790 Notes: 791 SNESComputeMinimizationFunction() is valid only for 792 SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 793 SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 794 795 SNESComputeMinimizationFunction() is typically used within minimization 796 implementations, so most users would not generally call this routine 797 themselves. 798 799 Level: developer 800 801 .keywords: SNES, nonlinear, compute, minimization, function 802 803 .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(), 804 SNESComputeGradient(), SNESComputeHessian() 805 @*/ 806 int SNESComputeMinimizationFunction(SNES snes,Vec x,PetscReal *y) 807 { 808 int ierr; 809 810 PetscFunctionBegin; 811 PetscValidHeaderSpecific(snes,SNES_COOKIE); 812 PetscValidHeaderSpecific(x,VEC_COOKIE); 813 PetscCheckSameComm(snes,x); 814 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 815 SETERRQ(PETSC_ERR_ARG_WRONG,"Only for SNES_UNCONSTRAINED_MINIMIZATION"); 816 } 817 818 ierr = PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0);CHKERRQ(ierr); 819 PetscStackPush("SNES user minimzation function"); 820 ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP);CHKERRQ(ierr); 821 PetscStackPop; 822 snes->nfuncs++; 823 ierr = PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0);CHKERRQ(ierr); 824 PetscFunctionReturn(0); 825 } 826 827 #undef __FUNC__ 828 #define __FUNC__ /*<a name=""></a>*/"SNESSetGradient" 829 /*@C 830 SNESSetGradient - Sets the gradient evaluation routine and gradient 831 vector for use by the SNES routines. 832 833 Collective on SNES 834 835 Input Parameters: 836 + snes - the SNES context 837 . func - function evaluation routine 838 . ctx - optional user-defined context for private data for the 839 gradient evaluation routine (may be PETSC_NULL) 840 - r - vector to store gradient value 841 842 Calling sequence of func: 843 $ func (SNES, Vec x, Vec g, void *ctx); 844 845 + x - input vector 846 . g - gradient vector 847 - ctx - optional user-defined gradient context 848 849 Notes: 850 SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 851 methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 852 SNESSetFunction(). 853 854 Level: beginner 855 856 .keywords: SNES, nonlinear, set, function 857 858 .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(), 859 SNESSetMinimizationFunction(), 860 @*/ 861 int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx) 862 { 863 PetscFunctionBegin; 864 PetscValidHeaderSpecific(snes,SNES_COOKIE); 865 PetscValidHeaderSpecific(r,VEC_COOKIE); 866 PetscCheckSameComm(snes,r); 867 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 868 SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 869 } 870 snes->computefunction = func; 871 snes->vec_func = snes->vec_func_always = r; 872 snes->funP = ctx; 873 PetscFunctionReturn(0); 874 } 875 876 #undef __FUNC__ 877 #define __FUNC__ /*<a name=""></a>*/"SNESComputeGradient" 878 /*@ 879 SNESComputeGradient - Computes the gradient that has been set with 880 SNESSetGradient(). 881 882 Collective on SNES 883 884 Input Parameters: 885 + snes - the SNES context 886 - x - input vector 887 888 Output Parameter: 889 . y - gradient vector 890 891 Notes: 892 SNESComputeGradient() is valid only for 893 SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 894 SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 895 896 SNESComputeGradient() is typically used within minimization 897 implementations, so most users would not generally call this routine 898 themselves. 899 900 Level: developer 901 902 .keywords: SNES, nonlinear, compute, gradient 903 904 .seealso: SNESSetGradient(), SNESGetGradient(), 905 SNESComputeMinimizationFunction(), SNESComputeHessian() 906 @*/ 907 int SNESComputeGradient(SNES snes,Vec x,Vec y) 908 { 909 int ierr; 910 911 PetscFunctionBegin; 912 PetscValidHeaderSpecific(snes,SNES_COOKIE); 913 PetscValidHeaderSpecific(x,VEC_COOKIE); 914 PetscValidHeaderSpecific(y,VEC_COOKIE); 915 PetscCheckSameComm(snes,x); 916 PetscCheckSameComm(snes,y); 917 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 918 SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 919 } 920 ierr = PLogEventBegin(SNES_GradientEval,snes,x,y,0);CHKERRQ(ierr); 921 PetscStackPush("SNES user gradient function"); 922 ierr = (*snes->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr); 923 PetscStackPop; 924 ierr = PLogEventEnd(SNES_GradientEval,snes,x,y,0);CHKERRQ(ierr); 925 PetscFunctionReturn(0); 926 } 927 928 #undef __FUNC__ 929 #define __FUNC__ /*<a name=""></a>*/"SNESComputeJacobian" 930 /*@ 931 SNESComputeJacobian - Computes the Jacobian matrix that has been 932 set with SNESSetJacobian(). 933 934 Collective on SNES and Mat 935 936 Input Parameters: 937 + snes - the SNES context 938 - x - input vector 939 940 Output Parameters: 941 + A - Jacobian matrix 942 . B - optional preconditioning matrix 943 - flag - flag indicating matrix structure 944 945 Notes: 946 Most users should not need to explicitly call this routine, as it 947 is used internally within the nonlinear solvers. 948 949 See SLESSetOperators() for important information about setting the 950 flag parameter. 951 952 SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS 953 methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION 954 methods is SNESComputeHessian(). 955 956 Level: developer 957 958 .keywords: SNES, compute, Jacobian, matrix 959 960 .seealso: SNESSetJacobian(), SLESSetOperators() 961 @*/ 962 int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 963 { 964 int ierr; 965 966 PetscFunctionBegin; 967 PetscValidHeaderSpecific(snes,SNES_COOKIE); 968 PetscValidHeaderSpecific(X,VEC_COOKIE); 969 PetscCheckSameComm(snes,X); 970 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 971 SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only"); 972 } 973 if (!snes->computejacobian) PetscFunctionReturn(0); 974 ierr = PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 975 *flg = DIFFERENT_NONZERO_PATTERN; 976 PetscStackPush("SNES user Jacobian function"); 977 ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr); 978 PetscStackPop; 979 ierr = PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 980 /* make sure user returned a correct Jacobian and preconditioner */ 981 PetscValidHeaderSpecific(*A,MAT_COOKIE); 982 PetscValidHeaderSpecific(*B,MAT_COOKIE); 983 PetscFunctionReturn(0); 984 } 985 986 #undef __FUNC__ 987 #define __FUNC__ /*<a name=""></a>*/"SNESComputeHessian" 988 /*@ 989 SNESComputeHessian - Computes the Hessian matrix that has been 990 set with SNESSetHessian(). 991 992 Collective on SNES and Mat 993 994 Input Parameters: 995 + snes - the SNES context 996 - x - input vector 997 998 Output Parameters: 999 + A - Hessian matrix 1000 . B - optional preconditioning matrix 1001 - flag - flag indicating matrix structure 1002 1003 Notes: 1004 Most users should not need to explicitly call this routine, as it 1005 is used internally within the nonlinear solvers. 1006 1007 See SLESSetOperators() for important information about setting the 1008 flag parameter. 1009 1010 SNESComputeHessian() is valid only for 1011 SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 1012 SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian(). 1013 1014 SNESComputeHessian() is typically used within minimization 1015 implementations, so most users would not generally call this routine 1016 themselves. 1017 1018 Level: developer 1019 1020 .keywords: SNES, compute, Hessian, matrix 1021 1022 .seealso: SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(), 1023 SNESComputeMinimizationFunction() 1024 @*/ 1025 int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag) 1026 { 1027 int ierr; 1028 1029 PetscFunctionBegin; 1030 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1031 PetscValidHeaderSpecific(x,VEC_COOKIE); 1032 PetscCheckSameComm(snes,x); 1033 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 1034 SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 1035 } 1036 if (!snes->computejacobian) PetscFunctionReturn(0); 1037 ierr = PLogEventBegin(SNES_HessianEval,snes,x,*A,*B);CHKERRQ(ierr); 1038 *flag = DIFFERENT_NONZERO_PATTERN; 1039 PetscStackPush("SNES user Hessian function"); 1040 ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP);CHKERRQ(ierr); 1041 PetscStackPop; 1042 ierr = PLogEventEnd(SNES_HessianEval,snes,x,*A,*B);CHKERRQ(ierr); 1043 /* make sure user returned a correct Jacobian and preconditioner */ 1044 PetscValidHeaderSpecific(*A,MAT_COOKIE); 1045 PetscValidHeaderSpecific(*B,MAT_COOKIE); 1046 PetscFunctionReturn(0); 1047 } 1048 1049 #undef __FUNC__ 1050 #define __FUNC__ /*<a name=""></a>*/"SNESSetJacobian" 1051 /*@C 1052 SNESSetJacobian - Sets the function to compute Jacobian as well as the 1053 location to store the matrix. 1054 1055 Collective on SNES and Mat 1056 1057 Input Parameters: 1058 + snes - the SNES context 1059 . A - Jacobian matrix 1060 . B - preconditioner matrix (usually same as the Jacobian) 1061 . func - Jacobian evaluation routine 1062 - ctx - [optional] user-defined context for private data for the 1063 Jacobian evaluation routine (may be PETSC_NULL) 1064 1065 Calling sequence of func: 1066 $ func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 1067 1068 + x - input vector 1069 . A - Jacobian matrix 1070 . B - preconditioner matrix, usually the same as A 1071 . flag - flag indicating information about the preconditioner matrix 1072 structure (same as flag in SLESSetOperators()) 1073 - ctx - [optional] user-defined Jacobian context 1074 1075 Notes: 1076 See SLESSetOperators() for important information about setting the flag 1077 output parameter in the routine func(). Be sure to read this information! 1078 1079 The routine func() takes Mat * as the matrix arguments rather than Mat. 1080 This allows the Jacobian evaluation routine to replace A and/or B with a 1081 completely new new matrix structure (not just different matrix elements) 1082 when appropriate, for instance, if the nonzero structure is changing 1083 throughout the global iterations. 1084 1085 Level: beginner 1086 1087 .keywords: SNES, nonlinear, set, Jacobian, matrix 1088 1089 .seealso: SLESSetOperators(), SNESSetFunction() 1090 @*/ 1091 int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 1092 { 1093 PetscFunctionBegin; 1094 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1095 if (A) PetscValidHeaderSpecific(A,MAT_COOKIE); 1096 if (B) PetscValidHeaderSpecific(B,MAT_COOKIE); 1097 if (A) PetscCheckSameComm(snes,A); 1098 if (B) PetscCheckSameComm(snes,B); 1099 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 1100 SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only"); 1101 } 1102 1103 snes->computejacobian = func; 1104 snes->jacP = ctx; 1105 snes->jacobian = A; 1106 snes->jacobian_pre = B; 1107 PetscFunctionReturn(0); 1108 } 1109 1110 #undef __FUNC__ 1111 #define __FUNC__ /*<a name=""></a>*/"SNESGetJacobian" 1112 /*@C 1113 SNESGetJacobian - Returns the Jacobian matrix and optionally the user 1114 provided context for evaluating the Jacobian. 1115 1116 Not Collective, but Mat object will be parallel if SNES object is 1117 1118 Input Parameter: 1119 . snes - the nonlinear solver context 1120 1121 Output Parameters: 1122 + A - location to stash Jacobian matrix (or PETSC_NULL) 1123 . B - location to stash preconditioner matrix (or PETSC_NULL) 1124 . ctx - location to stash Jacobian ctx (or PETSC_NULL) 1125 - func - location to put Jacobian function (or PETSC_NULL) 1126 1127 Level: advanced 1128 1129 .seealso: SNESSetJacobian(), SNESComputeJacobian() 1130 @*/ 1131 int SNESGetJacobian(SNES snes,Mat *A,Mat *B,void **ctx,int (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*)) 1132 { 1133 PetscFunctionBegin; 1134 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1135 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 1136 SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only"); 1137 } 1138 if (A) *A = snes->jacobian; 1139 if (B) *B = snes->jacobian_pre; 1140 if (ctx) *ctx = snes->jacP; 1141 if (func) *func = snes->computejacobian; 1142 PetscFunctionReturn(0); 1143 } 1144 1145 #undef __FUNC__ 1146 #define __FUNC__ /*<a name=""></a>*/"SNESSetHessian" 1147 /*@C 1148 SNESSetHessian - Sets the function to compute Hessian as well as the 1149 location to store the matrix. 1150 1151 Collective on SNES and Mat 1152 1153 Input Parameters: 1154 + snes - the SNES context 1155 . A - Hessian matrix 1156 . B - preconditioner matrix (usually same as the Hessian) 1157 . func - Jacobian evaluation routine 1158 - ctx - [optional] user-defined context for private data for the 1159 Hessian evaluation routine (may be PETSC_NULL) 1160 1161 Calling sequence of func: 1162 $ func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 1163 1164 + x - input vector 1165 . A - Hessian matrix 1166 . B - preconditioner matrix, usually the same as A 1167 . flag - flag indicating information about the preconditioner matrix 1168 structure (same as flag in SLESSetOperators()) 1169 - ctx - [optional] user-defined Hessian context 1170 1171 Notes: 1172 See SLESSetOperators() for important information about setting the flag 1173 output parameter in the routine func(). Be sure to read this information! 1174 1175 The function func() takes Mat * as the matrix arguments rather than Mat. 1176 This allows the Hessian evaluation routine to replace A and/or B with a 1177 completely new new matrix structure (not just different matrix elements) 1178 when appropriate, for instance, if the nonzero structure is changing 1179 throughout the global iterations. 1180 1181 Level: beginner 1182 1183 .keywords: SNES, nonlinear, set, Hessian, matrix 1184 1185 .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators() 1186 @*/ 1187 int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 1188 { 1189 PetscFunctionBegin; 1190 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1191 PetscValidHeaderSpecific(A,MAT_COOKIE); 1192 PetscValidHeaderSpecific(B,MAT_COOKIE); 1193 PetscCheckSameComm(snes,A); 1194 PetscCheckSameComm(snes,B); 1195 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 1196 SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 1197 } 1198 snes->computejacobian = func; 1199 snes->jacP = ctx; 1200 snes->jacobian = A; 1201 snes->jacobian_pre = B; 1202 PetscFunctionReturn(0); 1203 } 1204 1205 #undef __FUNC__ 1206 #define __FUNC__ /*<a name=""></a>*/"SNESGetHessian" 1207 /*@ 1208 SNESGetHessian - Returns the Hessian matrix and optionally the user 1209 provided context for evaluating the Hessian. 1210 1211 Not Collective, but Mat object is parallel if SNES object is parallel 1212 1213 Input Parameter: 1214 . snes - the nonlinear solver context 1215 1216 Output Parameters: 1217 + A - location to stash Hessian matrix (or PETSC_NULL) 1218 . B - location to stash preconditioner matrix (or PETSC_NULL) 1219 - ctx - location to stash Hessian ctx (or PETSC_NULL) 1220 1221 Level: advanced 1222 1223 .seealso: SNESSetHessian(), SNESComputeHessian() 1224 1225 .keywords: SNES, get, Hessian 1226 @*/ 1227 int SNESGetHessian(SNES snes,Mat *A,Mat *B,void **ctx) 1228 { 1229 PetscFunctionBegin; 1230 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1231 if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION){ 1232 SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 1233 } 1234 if (A) *A = snes->jacobian; 1235 if (B) *B = snes->jacobian_pre; 1236 if (ctx) *ctx = snes->jacP; 1237 PetscFunctionReturn(0); 1238 } 1239 1240 /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 1241 1242 #undef __FUNC__ 1243 #define __FUNC__ /*<a name=""></a>*/"SNESSetUp" 1244 /*@ 1245 SNESSetUp - Sets up the internal data structures for the later use 1246 of a nonlinear solver. 1247 1248 Collective on SNES 1249 1250 Input Parameters: 1251 + snes - the SNES context 1252 - x - the solution vector 1253 1254 Notes: 1255 For basic use of the SNES solvers the user need not explicitly call 1256 SNESSetUp(), since these actions will automatically occur during 1257 the call to SNESSolve(). However, if one wishes to control this 1258 phase separately, SNESSetUp() should be called after SNESCreate() 1259 and optional routines of the form SNESSetXXX(), but before SNESSolve(). 1260 1261 Level: advanced 1262 1263 .keywords: SNES, nonlinear, setup 1264 1265 .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 1266 @*/ 1267 int SNESSetUp(SNES snes,Vec x) 1268 { 1269 int ierr; 1270 PetscTruth flg; 1271 1272 PetscFunctionBegin; 1273 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1274 PetscValidHeaderSpecific(x,VEC_COOKIE); 1275 PetscCheckSameComm(snes,x); 1276 snes->vec_sol = snes->vec_sol_always = x; 1277 1278 ierr = OptionsHasName(snes->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr); 1279 /* 1280 This version replaces the user provided Jacobian matrix with a 1281 matrix-free version but still employs the user-provided preconditioner matrix 1282 */ 1283 if (flg) { 1284 Mat J; 1285 ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1286 PLogObjectParent(snes,J); 1287 snes->mfshell = J; 1288 snes->jacobian = J; 1289 if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 1290 PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n"); 1291 } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 1292 PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n"); 1293 } else { 1294 SETERRQ(PETSC_ERR_SUP,"Method class doesn't support matrix-free operator option"); 1295 } 1296 ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 1297 } 1298 ierr = OptionsHasName(snes->prefix,"-snes_mf",&flg);CHKERRQ(ierr); 1299 /* 1300 This version replaces both the user-provided Jacobian and the user- 1301 provided preconditioner matrix with the default matrix free version. 1302 */ 1303 if (flg) { 1304 Mat J; 1305 ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1306 PLogObjectParent(snes,J); 1307 snes->mfshell = J; 1308 if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 1309 ierr = SNESSetJacobian(snes,J,J,MatSNESMFFormJacobian,snes->funP);CHKERRQ(ierr); 1310 PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n"); 1311 } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 1312 ierr = SNESSetHessian(snes,J,J,MatSNESMFFormJacobian,snes->funP);CHKERRQ(ierr); 1313 PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n"); 1314 } else { 1315 SETERRQ(PETSC_ERR_SUP,"Method class doesn't support matrix-free option"); 1316 } 1317 ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 1318 } 1319 if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) { 1320 PetscTruth iseqtr; 1321 1322 if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 1323 if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 1324 if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option"); 1325 if (snes->vec_func == snes->vec_sol) { 1326 SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector"); 1327 } 1328 1329 /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 1330 ierr = PetscTypeCompare((PetscObject)snes,SNESEQTR,&iseqtr);CHKERRQ(ierr); 1331 if (snes->ksp_ewconv && !iseqtr) { 1332 SLES sles; KSP ksp; 1333 ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr); 1334 ierr = SLESGetKSP(sles,&ksp);CHKERRQ(ierr); 1335 ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,snes);CHKERRQ(ierr); 1336 } 1337 } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) { 1338 if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetGradient() first"); 1339 if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetGradient() first"); 1340 if (!snes->computeumfunction) { 1341 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetMinimizationFunction() first"); 1342 } 1343 if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetHessian()"); 1344 } else { 1345 SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Unknown method class"); 1346 } 1347 if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);} 1348 snes->setupcalled = 1; 1349 PetscFunctionReturn(0); 1350 } 1351 1352 #undef __FUNC__ 1353 #define __FUNC__ /*<a name=""></a>*/"SNESDestroy" 1354 /*@C 1355 SNESDestroy - Destroys the nonlinear solver context that was created 1356 with SNESCreate(). 1357 1358 Collective on SNES 1359 1360 Input Parameter: 1361 . snes - the SNES context 1362 1363 Level: beginner 1364 1365 .keywords: SNES, nonlinear, destroy 1366 1367 .seealso: SNESCreate(), SNESSolve() 1368 @*/ 1369 int SNESDestroy(SNES snes) 1370 { 1371 int i,ierr; 1372 1373 PetscFunctionBegin; 1374 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1375 if (--snes->refct > 0) PetscFunctionReturn(0); 1376 1377 /* if memory was published with AMS then destroy it */ 1378 ierr = PetscObjectDepublish(snes);CHKERRQ(ierr); 1379 1380 if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);} 1381 if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);} 1382 if (snes->mfshell) {ierr = MatDestroy(snes->mfshell);CHKERRQ(ierr);} 1383 ierr = SLESDestroy(snes->sles);CHKERRQ(ierr); 1384 if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);} 1385 for (i=0; i<snes->numbermonitors; i++) { 1386 if (snes->monitordestroy[i]) { 1387 ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr); 1388 } 1389 } 1390 PLogObjectDestroy((PetscObject)snes); 1391 PetscHeaderDestroy((PetscObject)snes); 1392 PetscFunctionReturn(0); 1393 } 1394 1395 /* ----------- Routines to set solver parameters ---------- */ 1396 1397 #undef __FUNC__ 1398 #define __FUNC__ /*<a name=""></a>*/"SNESSetTolerances" 1399 /*@ 1400 SNESSetTolerances - Sets various parameters used in convergence tests. 1401 1402 Collective on SNES 1403 1404 Input Parameters: 1405 + snes - the SNES context 1406 . atol - absolute convergence tolerance 1407 . rtol - relative convergence tolerance 1408 . stol - convergence tolerance in terms of the norm 1409 of the change in the solution between steps 1410 . maxit - maximum number of iterations 1411 - maxf - maximum number of function evaluations 1412 1413 Options Database Keys: 1414 + -snes_atol <atol> - Sets atol 1415 . -snes_rtol <rtol> - Sets rtol 1416 . -snes_stol <stol> - Sets stol 1417 . -snes_max_it <maxit> - Sets maxit 1418 - -snes_max_funcs <maxf> - Sets maxf 1419 1420 Notes: 1421 The default maximum number of iterations is 50. 1422 The default maximum number of function evaluations is 1000. 1423 1424 Level: intermediate 1425 1426 .keywords: SNES, nonlinear, set, convergence, tolerances 1427 1428 .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance() 1429 @*/ 1430 int SNESSetTolerances(SNES snes,PetscReal atol,PetscReal rtol,PetscReal stol,int maxit,int maxf) 1431 { 1432 PetscFunctionBegin; 1433 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1434 if (atol != PETSC_DEFAULT) snes->atol = atol; 1435 if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1436 if (stol != PETSC_DEFAULT) snes->xtol = stol; 1437 if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1438 if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 1439 PetscFunctionReturn(0); 1440 } 1441 1442 #undef __FUNC__ 1443 #define __FUNC__ /*<a name=""></a>*/"SNESGetTolerances" 1444 /*@ 1445 SNESGetTolerances - Gets various parameters used in convergence tests. 1446 1447 Not Collective 1448 1449 Input Parameters: 1450 + snes - the SNES context 1451 . atol - absolute convergence tolerance 1452 . rtol - relative convergence tolerance 1453 . stol - convergence tolerance in terms of the norm 1454 of the change in the solution between steps 1455 . maxit - maximum number of iterations 1456 - maxf - maximum number of function evaluations 1457 1458 Notes: 1459 The user can specify PETSC_NULL for any parameter that is not needed. 1460 1461 Level: intermediate 1462 1463 .keywords: SNES, nonlinear, get, convergence, tolerances 1464 1465 .seealso: SNESSetTolerances() 1466 @*/ 1467 int SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,int *maxit,int *maxf) 1468 { 1469 PetscFunctionBegin; 1470 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1471 if (atol) *atol = snes->atol; 1472 if (rtol) *rtol = snes->rtol; 1473 if (stol) *stol = snes->xtol; 1474 if (maxit) *maxit = snes->max_its; 1475 if (maxf) *maxf = snes->max_funcs; 1476 PetscFunctionReturn(0); 1477 } 1478 1479 #undef __FUNC__ 1480 #define __FUNC__ /*<a name=""></a>*/"SNESSetTrustRegionTolerance" 1481 /*@ 1482 SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 1483 1484 Collective on SNES 1485 1486 Input Parameters: 1487 + snes - the SNES context 1488 - tol - tolerance 1489 1490 Options Database Key: 1491 . -snes_trtol <tol> - Sets tol 1492 1493 Level: intermediate 1494 1495 .keywords: SNES, nonlinear, set, trust region, tolerance 1496 1497 .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance() 1498 @*/ 1499 int SNESSetTrustRegionTolerance(SNES snes,PetscReal tol) 1500 { 1501 PetscFunctionBegin; 1502 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1503 snes->deltatol = tol; 1504 PetscFunctionReturn(0); 1505 } 1506 1507 #undef __FUNC__ 1508 #define __FUNC__ /*<a name=""></a>*/"SNESSetMinimizationFunctionTolerance" 1509 /*@ 1510 SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance 1511 for unconstrained minimization solvers. 1512 1513 Collective on SNES 1514 1515 Input Parameters: 1516 + snes - the SNES context 1517 - ftol - minimum function tolerance 1518 1519 Options Database Key: 1520 . -snes_fmin <ftol> - Sets ftol 1521 1522 Note: 1523 SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION 1524 methods only. 1525 1526 Level: intermediate 1527 1528 .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance 1529 1530 .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance() 1531 @*/ 1532 int SNESSetMinimizationFunctionTolerance(SNES snes,PetscReal ftol) 1533 { 1534 PetscFunctionBegin; 1535 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1536 snes->fmin = ftol; 1537 PetscFunctionReturn(0); 1538 } 1539 /* 1540 Duplicate the lg monitors for SNES from KSP; for some reason with 1541 dynamic libraries things don't work under Sun4 if we just use 1542 macros instead of functions 1543 */ 1544 #undef __FUNC__ 1545 #define __FUNC__ /*<a name=""></a>*/"SNESLGMonitor" 1546 int SNESLGMonitor(SNES snes,int it,PetscReal norm,void *ctx) 1547 { 1548 int ierr; 1549 1550 PetscFunctionBegin; 1551 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1552 ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr); 1553 PetscFunctionReturn(0); 1554 } 1555 1556 #undef __FUNC__ 1557 #define __FUNC__ /*<a name=""></a>*/"SNESLGMonitorCreate" 1558 int SNESLGMonitorCreate(char *host,char *label,int x,int y,int m,int n,DrawLG *draw) 1559 { 1560 int ierr; 1561 1562 PetscFunctionBegin; 1563 ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr); 1564 PetscFunctionReturn(0); 1565 } 1566 1567 #undef __FUNC__ 1568 #define __FUNC__ /*<a name=""></a>*/"SNESLGMonitorDestroy" 1569 int SNESLGMonitorDestroy(DrawLG draw) 1570 { 1571 int ierr; 1572 1573 PetscFunctionBegin; 1574 ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr); 1575 PetscFunctionReturn(0); 1576 } 1577 1578 /* ------------ Routines to set performance monitoring options ----------- */ 1579 1580 #undef __FUNC__ 1581 #define __FUNC__ /*<a name=""></a>*/"SNESSetMonitor" 1582 /*@C 1583 SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every 1584 iteration of the nonlinear solver to display the iteration's 1585 progress. 1586 1587 Collective on SNES 1588 1589 Input Parameters: 1590 + snes - the SNES context 1591 . func - monitoring routine 1592 . mctx - [optional] user-defined context for private data for the 1593 monitor routine (use PETSC_NULL if no context is desitre) 1594 - monitordestroy - [optional] routine that frees monitor context 1595 (may be PETSC_NULL) 1596 1597 Calling sequence of func: 1598 $ int func(SNES snes,int its, PetscReal norm,void *mctx) 1599 1600 + snes - the SNES context 1601 . its - iteration number 1602 . norm - 2-norm function value (may be estimated) 1603 for SNES_NONLINEAR_EQUATIONS methods 1604 . norm - 2-norm gradient value (may be estimated) 1605 for SNES_UNCONSTRAINED_MINIMIZATION methods 1606 - mctx - [optional] monitoring context 1607 1608 Options Database Keys: 1609 + -snes_monitor - sets SNESDefaultMonitor() 1610 . -snes_xmonitor - sets line graph monitor, 1611 uses SNESLGMonitorCreate() 1612 _ -snes_cancelmonitors - cancels all monitors that have 1613 been hardwired into a code by 1614 calls to SNESSetMonitor(), but 1615 does not cancel those set via 1616 the options database. 1617 1618 Notes: 1619 Several different monitoring routines may be set by calling 1620 SNESSetMonitor() multiple times; all will be called in the 1621 order in which they were set. 1622 1623 Level: intermediate 1624 1625 .keywords: SNES, nonlinear, set, monitor 1626 1627 .seealso: SNESDefaultMonitor(), SNESClearMonitor() 1628 @*/ 1629 int SNESSetMonitor(SNES snes,int (*func)(SNES,int,PetscReal,void*),void *mctx,int (*monitordestroy)(void *)) 1630 { 1631 PetscFunctionBegin; 1632 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1633 if (snes->numbermonitors >= MAXSNESMONITORS) { 1634 SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1635 } 1636 1637 snes->monitor[snes->numbermonitors] = func; 1638 snes->monitordestroy[snes->numbermonitors] = monitordestroy; 1639 snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 1640 PetscFunctionReturn(0); 1641 } 1642 1643 #undef __FUNC__ 1644 #define __FUNC__ /*<a name=""></a>*/"SNESClearMonitor" 1645 /*@C 1646 SNESClearMonitor - Clears all the monitor functions for a SNES object. 1647 1648 Collective on SNES 1649 1650 Input Parameters: 1651 . snes - the SNES context 1652 1653 Options Database: 1654 . -snes_cancelmonitors - cancels all monitors that have been hardwired 1655 into a code by calls to SNESSetMonitor(), but does not cancel those 1656 set via the options database 1657 1658 Notes: 1659 There is no way to clear one specific monitor from a SNES object. 1660 1661 Level: intermediate 1662 1663 .keywords: SNES, nonlinear, set, monitor 1664 1665 .seealso: SNESDefaultMonitor(), SNESSetMonitor() 1666 @*/ 1667 int SNESClearMonitor(SNES snes) 1668 { 1669 PetscFunctionBegin; 1670 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1671 snes->numbermonitors = 0; 1672 PetscFunctionReturn(0); 1673 } 1674 1675 #undef __FUNC__ 1676 #define __FUNC__ /*<a name=""></a>*/"SNESSetConvergenceTest" 1677 /*@C 1678 SNESSetConvergenceTest - Sets the function that is to be used 1679 to test for convergence of the nonlinear iterative solution. 1680 1681 Collective on SNES 1682 1683 Input Parameters: 1684 + snes - the SNES context 1685 . func - routine to test for convergence 1686 - cctx - [optional] context for private data for the convergence routine 1687 (may be PETSC_NULL) 1688 1689 Calling sequence of func: 1690 $ int func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx) 1691 1692 + snes - the SNES context 1693 . cctx - [optional] convergence context 1694 . reason - reason for convergence/divergence 1695 . xnorm - 2-norm of current iterate 1696 . gnorm - 2-norm of current step (SNES_NONLINEAR_EQUATIONS methods) 1697 . f - 2-norm of function (SNES_NONLINEAR_EQUATIONS methods) 1698 . gnorm - 2-norm of current gradient (SNES_UNCONSTRAINED_MINIMIZATION methods) 1699 - f - function value (SNES_UNCONSTRAINED_MINIMIZATION methods) 1700 1701 Level: advanced 1702 1703 .keywords: SNES, nonlinear, set, convergence, test 1704 1705 .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(), 1706 SNESConverged_UM_LS(), SNESConverged_UM_TR() 1707 @*/ 1708 int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx) 1709 { 1710 PetscFunctionBegin; 1711 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1712 (snes)->converged = func; 1713 (snes)->cnvP = cctx; 1714 PetscFunctionReturn(0); 1715 } 1716 1717 #undef __FUNC__ 1718 #define __FUNC__ /*<a name=""></a>*/"SNESGetConvergedReason" 1719 /*@C 1720 SNESGetConvergedReason - Gets the reason the SNES iteration was stopped. 1721 1722 Not Collective 1723 1724 Input Parameter: 1725 . snes - the SNES context 1726 1727 Output Parameter: 1728 . reason - negative value indicates diverged, positive value converged, see petscsnes.h or the 1729 manual pages for the individual convergence tests for complete lists 1730 1731 Level: intermediate 1732 1733 Notes: Can only be called after the call the SNESSolve() is complete. 1734 1735 .keywords: SNES, nonlinear, set, convergence, test 1736 1737 .seealso: SNESSetConvergenceTest(), SNESConverged_EQ_LS(), SNESConverged_EQ_TR(), 1738 SNESConverged_UM_LS(), SNESConverged_UM_TR() 1739 @*/ 1740 int SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason) 1741 { 1742 PetscFunctionBegin; 1743 PetscValidHeaderSpecific(snes,SNES_COOKIE); 1744 *reason = snes->reason; 1745 PetscFunctionReturn(0); 1746 } 1747 1748 #undef __FUNC__ 1749 #define __FUNC__ /*<a name=""></a>*/"SNESSetConvergenceHistory" 1750 /*@ 1751 SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1752 1753 Collective on SNES 1754 1755 Input Parameters: 1756 + snes - iterative context obtained from SNESCreate() 1757 . a - array to hold history 1758 . its - integer array holds the number of linear iterations 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,"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 ierr = PLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 1926 snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0; 1927 ierr = (*(snes)->solve)(snes,its);CHKERRQ(ierr); 1928 ierr = PLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 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,"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,"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,"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,"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