1 #define PETSCSNES_DLL 2 3 #include "private/snesimpl.h" /*I "petscsnes.h" I*/ 4 5 PetscTruth SNESRegisterAllCalled = PETSC_FALSE; 6 PetscFList SNESList = PETSC_NULL; 7 8 /* Logging support */ 9 PetscClassId PETSCSNES_DLLEXPORT SNES_CLASSID; 10 PetscLogEvent SNES_Solve, SNES_LineSearch, SNES_FunctionEval, SNES_JacobianEval; 11 12 #undef __FUNCT__ 13 #define __FUNCT__ "SNESSetErrorIfNotConverged" 14 /*@ 15 SNESSetErrorIfNotConverged - Causes SNESSolve() to generate an error if the solver has not converged. 16 17 Logically Collective on SNES 18 19 Input Parameters: 20 + snes - iterative context obtained from SNESCreate() 21 - flg - PETSC_TRUE indicates you want the error generated 22 23 Options database keys: 24 . -snes_error_if_not_converged : this takes an optional truth value (0/1/no/yes/true/false) 25 26 Level: intermediate 27 28 Notes: 29 Normally PETSc continues if a linear solver fails to converge, you can call SNESGetConvergedReason() after a SNESSolve() 30 to determine if it has converged. 31 32 .keywords: SNES, set, initial guess, nonzero 33 34 .seealso: SNESGetErrorIfNotConverged(), KSPGetErrorIfNotConverged(), KSPSetErrorIFNotConverged() 35 @*/ 36 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetErrorIfNotConverged(SNES snes,PetscTruth flg) 37 { 38 PetscFunctionBegin; 39 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 40 PetscValidLogicalCollectiveTruth(snes,flg,2); 41 snes->errorifnotconverged = flg; 42 PetscFunctionReturn(0); 43 } 44 45 #undef __FUNCT__ 46 #define __FUNCT__ "SNESGetErrorIfNotConverged" 47 /*@ 48 SNESGetErrorIfNotConverged - Will SNESSolve() generate an error if the solver does not converge? 49 50 Not Collective 51 52 Input Parameter: 53 . snes - iterative context obtained from SNESCreate() 54 55 Output Parameter: 56 . flag - PETSC_TRUE if it will generate an error, else PETSC_FALSE 57 58 Level: intermediate 59 60 .keywords: SNES, set, initial guess, nonzero 61 62 .seealso: SNESSetErrorIfNotConverged(), KSPGetErrorIfNotConverged(), KSPSetErrorIFNotConverged() 63 @*/ 64 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetErrorIfNotConverged(SNES snes,PetscTruth *flag) 65 { 66 PetscFunctionBegin; 67 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 68 PetscValidPointer(flag,2); 69 *flag = snes->errorifnotconverged; 70 PetscFunctionReturn(0); 71 } 72 73 #undef __FUNCT__ 74 #define __FUNCT__ "SNESSetFunctionDomainError" 75 /*@ 76 SNESSetFunctionDomainError - tells SNES that the input vector to your FormFunction is not 77 in the functions domain. For example, negative pressure. 78 79 Logically Collective on SNES 80 81 Input Parameters: 82 . SNES - the SNES context 83 84 Level: advanced 85 86 .keywords: SNES, view 87 88 .seealso: SNESCreate(), SNESSetFunction() 89 @*/ 90 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFunctionDomainError(SNES snes) 91 { 92 PetscFunctionBegin; 93 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 94 snes->domainerror = PETSC_TRUE; 95 PetscFunctionReturn(0); 96 } 97 98 #undef __FUNCT__ 99 #define __FUNCT__ "SNESView" 100 /*@C 101 SNESView - Prints the SNES data structure. 102 103 Collective on SNES 104 105 Input Parameters: 106 + SNES - the SNES context 107 - viewer - visualization context 108 109 Options Database Key: 110 . -snes_view - Calls SNESView() at end of SNESSolve() 111 112 Notes: 113 The available visualization contexts include 114 + PETSC_VIEWER_STDOUT_SELF - standard output (default) 115 - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 116 output where only the first processor opens 117 the file. All other processors send their 118 data to the first processor to print. 119 120 The user can open an alternative visualization context with 121 PetscViewerASCIIOpen() - output to a specified file. 122 123 Level: beginner 124 125 .keywords: SNES, view 126 127 .seealso: PetscViewerASCIIOpen() 128 @*/ 129 PetscErrorCode PETSCSNES_DLLEXPORT SNESView(SNES snes,PetscViewer viewer) 130 { 131 SNESKSPEW *kctx; 132 PetscErrorCode ierr; 133 KSP ksp; 134 const SNESType type; 135 PetscTruth iascii,isstring; 136 137 PetscFunctionBegin; 138 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 139 if (!viewer) { 140 ierr = PetscViewerASCIIGetStdout(((PetscObject)snes)->comm,&viewer);CHKERRQ(ierr); 141 } 142 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 143 PetscCheckSameComm(snes,1,viewer,2); 144 145 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 146 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 147 if (iascii) { 148 if (((PetscObject)snes)->prefix) { 149 ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:(%s)\n",((PetscObject)snes)->prefix);CHKERRQ(ierr); 150 } else { 151 ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:\n");CHKERRQ(ierr); 152 } 153 ierr = SNESGetType(snes,&type);CHKERRQ(ierr); 154 if (type) { 155 ierr = PetscViewerASCIIPrintf(viewer," type: %s\n",type);CHKERRQ(ierr); 156 } else { 157 ierr = PetscViewerASCIIPrintf(viewer," type: not set yet\n");CHKERRQ(ierr); 158 } 159 if (snes->ops->view) { 160 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 161 ierr = (*snes->ops->view)(snes,viewer);CHKERRQ(ierr); 162 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 163 } 164 ierr = PetscViewerASCIIPrintf(viewer," maximum iterations=%D, maximum function evaluations=%D\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr); 165 ierr = PetscViewerASCIIPrintf(viewer," tolerances: relative=%G, absolute=%G, solution=%G\n", 166 snes->rtol,snes->abstol,snes->xtol);CHKERRQ(ierr); 167 ierr = PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",snes->linear_its);CHKERRQ(ierr); 168 ierr = PetscViewerASCIIPrintf(viewer," total number of function evaluations=%D\n",snes->nfuncs);CHKERRQ(ierr); 169 if (snes->ksp_ewconv) { 170 kctx = (SNESKSPEW *)snes->kspconvctx; 171 if (kctx) { 172 ierr = PetscViewerASCIIPrintf(viewer," Eisenstat-Walker computation of KSP relative tolerance (version %D)\n",kctx->version);CHKERRQ(ierr); 173 ierr = PetscViewerASCIIPrintf(viewer," rtol_0=%G, rtol_max=%G, threshold=%G\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr); 174 ierr = PetscViewerASCIIPrintf(viewer," gamma=%G, alpha=%G, alpha2=%G\n",kctx->gamma,kctx->alpha,kctx->alpha2);CHKERRQ(ierr); 175 } 176 } 177 if (snes->lagpreconditioner == -1) { 178 ierr = PetscViewerASCIIPrintf(viewer," Preconditioned is never rebuilt\n");CHKERRQ(ierr); 179 } else if (snes->lagpreconditioner > 1) { 180 ierr = PetscViewerASCIIPrintf(viewer," Preconditioned is rebuilt every %D new Jacobians\n",snes->lagpreconditioner);CHKERRQ(ierr); 181 } 182 if (snes->lagjacobian == -1) { 183 ierr = PetscViewerASCIIPrintf(viewer," Jacobian is never rebuilt\n");CHKERRQ(ierr); 184 } else if (snes->lagjacobian > 1) { 185 ierr = PetscViewerASCIIPrintf(viewer," Jacobian is rebuilt every %D Newton iterations\n",snes->lagjacobian);CHKERRQ(ierr); 186 } 187 } else if (isstring) { 188 ierr = SNESGetType(snes,&type);CHKERRQ(ierr); 189 ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr); 190 } 191 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 192 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 193 ierr = KSPView(ksp,viewer);CHKERRQ(ierr); 194 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 195 PetscFunctionReturn(0); 196 } 197 198 /* 199 We retain a list of functions that also take SNES command 200 line options. These are called at the end SNESSetFromOptions() 201 */ 202 #define MAXSETFROMOPTIONS 5 203 static PetscInt numberofsetfromoptions; 204 static PetscErrorCode (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES); 205 206 #undef __FUNCT__ 207 #define __FUNCT__ "SNESAddOptionsChecker" 208 /*@C 209 SNESAddOptionsChecker - Adds an additional function to check for SNES options. 210 211 Not Collective 212 213 Input Parameter: 214 . snescheck - function that checks for options 215 216 Level: developer 217 218 .seealso: SNESSetFromOptions() 219 @*/ 220 PetscErrorCode PETSCSNES_DLLEXPORT SNESAddOptionsChecker(PetscErrorCode (*snescheck)(SNES)) 221 { 222 PetscFunctionBegin; 223 if (numberofsetfromoptions >= MAXSETFROMOPTIONS) { 224 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %D allowed", MAXSETFROMOPTIONS); 225 } 226 othersetfromoptions[numberofsetfromoptions++] = snescheck; 227 PetscFunctionReturn(0); 228 } 229 230 EXTERN PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*); 231 232 #undef __FUNCT__ 233 #define __FUNCT__ "SNESSetUpMatrixFree_Private" 234 static PetscErrorCode SNESSetUpMatrixFree_Private(SNES snes, PetscTruth hasOperator, PetscInt version) 235 { 236 Mat J; 237 KSP ksp; 238 PC pc; 239 PetscTruth match; 240 PetscErrorCode ierr; 241 242 PetscFunctionBegin; 243 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 244 245 if(!snes->vec_func && (snes->jacobian || snes->jacobian_pre)) { 246 Mat A = snes->jacobian, B = snes->jacobian_pre; 247 ierr = MatGetVecs(A ? A : B, PETSC_NULL,&snes->vec_func);CHKERRQ(ierr); 248 } 249 250 if (version == 1) { 251 ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr); 252 ierr = MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);CHKERRQ(ierr); 253 ierr = MatMFFDSetFromOptions(J);CHKERRQ(ierr); 254 } else if (version == 2) { 255 if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"SNESSetFunction() must be called first"); 256 #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SCALAR_SINGLE) && !defined(PETSC_USE_SCALAR_MAT_SINGLE) && !defined(PETSC_USE_SCALAR_LONG_DOUBLE) && !defined(PETSC_USE_SCALAR_INT) 257 ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_func,&J);CHKERRQ(ierr); 258 #else 259 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "matrix-free operator rutines (version 2)"); 260 #endif 261 } else { 262 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "matrix-free operator rutines, only version 1 and 2"); 263 } 264 265 ierr = PetscInfo1(snes,"Setting default matrix-free operator routines (version %D)\n", version);CHKERRQ(ierr); 266 if (hasOperator) { 267 /* This version replaces the user provided Jacobian matrix with a 268 matrix-free version but still employs the user-provided preconditioner matrix. */ 269 ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 270 } else { 271 /* This version replaces both the user-provided Jacobian and the user- 272 provided preconditioner matrix with the default matrix free version. */ 273 ierr = SNESSetJacobian(snes,J,J,MatMFFDComputeJacobian,snes->funP);CHKERRQ(ierr); 274 /* Force no preconditioner */ 275 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 276 ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 277 ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&match);CHKERRQ(ierr); 278 if (!match) { 279 ierr = PetscInfo(snes,"Setting default matrix-free preconditioner routines\nThat is no preconditioner is being used\n");CHKERRQ(ierr); 280 ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 281 } 282 } 283 ierr = MatDestroy(J);CHKERRQ(ierr); 284 285 PetscFunctionReturn(0); 286 } 287 288 #undef __FUNCT__ 289 #define __FUNCT__ "SNESSetFromOptions" 290 /*@ 291 SNESSetFromOptions - Sets various SNES and KSP parameters from user options. 292 293 Collective on SNES 294 295 Input Parameter: 296 . snes - the SNES context 297 298 Options Database Keys: 299 + -snes_type <type> - ls, tr, umls, umtr, test 300 . -snes_stol - convergence tolerance in terms of the norm 301 of the change in the solution between steps 302 . -snes_atol <abstol> - absolute tolerance of residual norm 303 . -snes_rtol <rtol> - relative decrease in tolerance norm from initial 304 . -snes_max_it <max_it> - maximum number of iterations 305 . -snes_max_funcs <max_funcs> - maximum number of function evaluations 306 . -snes_max_fail <max_fail> - maximum number of failures 307 . -snes_max_linear_solve_fail - number of linear solver failures before SNESSolve() stops 308 . -snes_lag_preconditioner <lag> - how often preconditioner is rebuilt (use -1 to never rebuild) 309 . -snes_lag_jacobian <lag> - how often Jacobian is rebuilt (use -1 to never rebuild) 310 . -snes_trtol <trtol> - trust region tolerance 311 . -snes_no_convergence_test - skip convergence test in nonlinear 312 solver; hence iterations will continue until max_it 313 or some other criterion is reached. Saves expense 314 of convergence test 315 . -snes_monitor <optional filename> - prints residual norm at each iteration. if no 316 filename given prints to stdout 317 . -snes_monitor_solution - plots solution at each iteration 318 . -snes_monitor_residual - plots residual (not its norm) at each iteration 319 . -snes_monitor_solution_update - plots update to solution at each iteration 320 . -snes_monitor_draw - plots residual norm at each iteration 321 . -snes_fd - use finite differences to compute Jacobian; very slow, only for testing 322 . -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration 323 - -snes_converged_reason - print the reason for convergence/divergence after each solve 324 325 Options Database for Eisenstat-Walker method: 326 + -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence 327 . -snes_ksp_ew_version ver - version of Eisenstat-Walker method 328 . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0 329 . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax 330 . -snes_ksp_ew_gamma <gamma> - Sets gamma 331 . -snes_ksp_ew_alpha <alpha> - Sets alpha 332 . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2 333 - -snes_ksp_ew_threshold <threshold> - Sets threshold 334 335 Notes: 336 To see all options, run your program with the -help option or consult 337 the <A href="../../docs/manual.pdf#nameddest=ch_snes">SNES chapter of the users manual</A>. 338 339 Level: beginner 340 341 .keywords: SNES, nonlinear, set, options, database 342 343 .seealso: SNESSetOptionsPrefix() 344 @*/ 345 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFromOptions(SNES snes) 346 { 347 PetscTruth flg,mf,mf_operator; 348 PetscInt i,indx,lag,mf_version; 349 MatStructure matflag; 350 const char *deft = SNESLS; 351 const char *convtests[] = {"default","skip"}; 352 SNESKSPEW *kctx = NULL; 353 char type[256], monfilename[PETSC_MAX_PATH_LEN]; 354 PetscViewerASCIIMonitor monviewer; 355 PetscErrorCode ierr; 356 357 PetscFunctionBegin; 358 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 359 360 if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 361 ierr = PetscOptionsBegin(((PetscObject)snes)->comm,((PetscObject)snes)->prefix,"Nonlinear solver (SNES) options","SNES");CHKERRQ(ierr); 362 if (((PetscObject)snes)->type_name) { deft = ((PetscObject)snes)->type_name; } 363 ierr = PetscOptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr); 364 if (flg) { 365 ierr = SNESSetType(snes,type);CHKERRQ(ierr); 366 } else if (!((PetscObject)snes)->type_name) { 367 ierr = SNESSetType(snes,deft);CHKERRQ(ierr); 368 } 369 /* not used here, but called so will go into help messaage */ 370 ierr = PetscOptionsName("-snes_view","Print detailed information on solver used","SNESView",0);CHKERRQ(ierr); 371 372 ierr = PetscOptionsReal("-snes_stol","Stop if step length less than","SNESSetTolerances",snes->xtol,&snes->xtol,0);CHKERRQ(ierr); 373 ierr = PetscOptionsReal("-snes_atol","Stop if function norm less than","SNESSetTolerances",snes->abstol,&snes->abstol,0);CHKERRQ(ierr); 374 375 ierr = PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less than","SNESSetTolerances",snes->rtol,&snes->rtol,0);CHKERRQ(ierr); 376 ierr = PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,PETSC_NULL);CHKERRQ(ierr); 377 ierr = PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr); 378 ierr = PetscOptionsInt("-snes_max_fail","Maximum failures","SNESSetTolerances",snes->maxFailures,&snes->maxFailures,PETSC_NULL);CHKERRQ(ierr); 379 ierr = PetscOptionsInt("-snes_max_linear_solve_fail","Maximum failures in linear solves allowed","SNESSetMaxLinearSolveFailures",snes->maxLinearSolveFailures,&snes->maxLinearSolveFailures,PETSC_NULL);CHKERRQ(ierr); 380 ierr = PetscOptionsTruth("-snes_error_if_not_converged","Generate error if solver does not converge","SNESSetErrorIfNotConverged",snes->errorifnotconverged,&snes->errorifnotconverged,PETSC_NULL);CHKERRQ(ierr); 381 382 ierr = PetscOptionsInt("-snes_lag_preconditioner","How often to rebuild preconditioner","SNESSetLagPreconditioner",snes->lagpreconditioner,&lag,&flg);CHKERRQ(ierr); 383 if (flg) { 384 ierr = SNESSetLagPreconditioner(snes,lag);CHKERRQ(ierr); 385 } 386 ierr = PetscOptionsInt("-snes_lag_jacobian","How often to rebuild Jacobian","SNESSetLagJacobian",snes->lagjacobian,&lag,&flg);CHKERRQ(ierr); 387 if (flg) { 388 ierr = SNESSetLagJacobian(snes,lag);CHKERRQ(ierr); 389 } 390 391 ierr = PetscOptionsEList("-snes_convergence_test","Convergence test","SNESSetConvergenceTest",convtests,2,"default",&indx,&flg);CHKERRQ(ierr); 392 if (flg) { 393 switch (indx) { 394 case 0: ierr = SNESSetConvergenceTest(snes,SNESDefaultConverged,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); break; 395 case 1: ierr = SNESSetConvergenceTest(snes,SNESSkipConverged,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); break; 396 } 397 } 398 399 ierr = PetscOptionsTruth("-snes_converged_reason","Print reason for converged or diverged","SNESSolve",snes->printreason,&snes->printreason,PETSC_NULL);CHKERRQ(ierr); 400 401 kctx = (SNESKSPEW *)snes->kspconvctx; 402 403 ierr = PetscOptionsTruth("-snes_ksp_ew","Use Eisentat-Walker linear system convergence test","SNESKSPSetUseEW",snes->ksp_ewconv,&snes->ksp_ewconv,PETSC_NULL);CHKERRQ(ierr); 404 405 ierr = PetscOptionsInt("-snes_ksp_ew_version","Version 1, 2 or 3","SNESKSPSetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr); 406 ierr = PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNESKSPSetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr); 407 ierr = PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNESKSPSetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr); 408 ierr = PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNESKSPSetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr); 409 ierr = PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNESKSPSetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr); 410 ierr = PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNESKSPSetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr); 411 ierr = PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNESKSPSetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr); 412 413 flg = PETSC_FALSE; 414 ierr = PetscOptionsTruth("-snes_monitor_cancel","Remove all monitors","SNESMonitorCancel",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 415 if (flg) {ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);} 416 417 ierr = PetscOptionsString("-snes_monitor","Monitor norm of function","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 418 if (flg) { 419 ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,monfilename,((PetscObject)snes)->tablevel,&monviewer);CHKERRQ(ierr); 420 ierr = SNESMonitorSet(snes,SNESMonitorDefault,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 421 } 422 423 ierr = PetscOptionsString("-snes_monitor_range","Monitor range of elements of function","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 424 if (flg) { 425 ierr = SNESMonitorSet(snes,SNESMonitorRange,0,0);CHKERRQ(ierr); 426 } 427 428 ierr = PetscOptionsString("-snes_ratiomonitor","Monitor ratios of norms of function","SNESMonitorSetRatio","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 429 if (flg) { 430 ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,monfilename,((PetscObject)snes)->tablevel,&monviewer);CHKERRQ(ierr); 431 ierr = SNESMonitorSetRatio(snes,monviewer);CHKERRQ(ierr); 432 } 433 434 ierr = PetscOptionsString("-snes_monitor_short","Monitor norm of function (fewer digits)","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 435 if (flg) { 436 ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,monfilename,((PetscObject)snes)->tablevel,&monviewer);CHKERRQ(ierr); 437 ierr = SNESMonitorSet(snes,SNESMonitorDefaultShort,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 438 } 439 440 flg = PETSC_FALSE; 441 ierr = PetscOptionsTruth("-snes_monitor_solution","Plot solution at each iteration","SNESMonitorSolution",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 442 if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorSolution,0,0);CHKERRQ(ierr);} 443 flg = PETSC_FALSE; 444 ierr = PetscOptionsTruth("-snes_monitor_solution_update","Plot correction at each iteration","SNESMonitorSolutionUpdate",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 445 if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorSolutionUpdate,0,0);CHKERRQ(ierr);} 446 flg = PETSC_FALSE; 447 ierr = PetscOptionsTruth("-snes_monitor_residual","Plot residual at each iteration","SNESMonitorResidual",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 448 if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorResidual,0,0);CHKERRQ(ierr);} 449 flg = PETSC_FALSE; 450 ierr = PetscOptionsTruth("-snes_monitor_draw","Plot function norm at each iteration","SNESMonitorLG",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 451 if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);} 452 flg = PETSC_FALSE; 453 ierr = PetscOptionsTruth("-snes_monitor_range_draw","Plot function range at each iteration","SNESMonitorLG",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 454 if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorLGRange,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);} 455 456 flg = PETSC_FALSE; 457 ierr = PetscOptionsTruth("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 458 if (flg) { 459 ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr); 460 ierr = PetscInfo(snes,"Setting default finite difference Jacobian matrix\n");CHKERRQ(ierr); 461 } 462 463 mf = mf_operator = PETSC_FALSE; 464 flg = PETSC_FALSE; 465 ierr = PetscOptionsTruth("-snes_mf_operator","Use a Matrix-Free Jacobian with user-provided preconditioner matrix","MatCreateSNESMF",PETSC_FALSE,&mf_operator,&flg);CHKERRQ(ierr); 466 if (flg && mf_operator) mf = PETSC_TRUE; 467 flg = PETSC_FALSE; 468 ierr = PetscOptionsTruth("-snes_mf","Use a Matrix-Free Jacobian with no preconditioner matrix","MatCreateSNESMF",PETSC_FALSE,&mf,&flg);CHKERRQ(ierr); 469 if (!flg && mf_operator) mf = PETSC_TRUE; 470 mf_version = 1; 471 ierr = PetscOptionsInt("-snes_mf_version","Matrix-Free routines version 1 or 2","None",mf_version,&mf_version,0);CHKERRQ(ierr); 472 473 for(i = 0; i < numberofsetfromoptions; i++) { 474 ierr = (*othersetfromoptions[i])(snes);CHKERRQ(ierr); 475 } 476 477 if (snes->ops->setfromoptions) { 478 ierr = (*snes->ops->setfromoptions)(snes);CHKERRQ(ierr); 479 } 480 481 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 482 ierr = PetscObjectProcessOptionsHandlers((PetscObject)snes);CHKERRQ(ierr); 483 ierr = PetscOptionsEnd();CHKERRQ(ierr); 484 485 if (mf) { ierr = SNESSetUpMatrixFree_Private(snes, mf_operator, mf_version);CHKERRQ(ierr); } 486 487 if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);} 488 ierr = KSPGetOperators(snes->ksp,PETSC_NULL,PETSC_NULL,&matflag); 489 ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,matflag);CHKERRQ(ierr); 490 ierr = KSPSetFromOptions(snes->ksp);CHKERRQ(ierr); 491 492 PetscFunctionReturn(0); 493 } 494 495 496 #undef __FUNCT__ 497 #define __FUNCT__ "SNESSetApplicationContext" 498 /*@ 499 SNESSetApplicationContext - Sets the optional user-defined context for 500 the nonlinear solvers. 501 502 Logically Collective on SNES 503 504 Input Parameters: 505 + snes - the SNES context 506 - usrP - optional user context 507 508 Level: intermediate 509 510 .keywords: SNES, nonlinear, set, application, context 511 512 .seealso: SNESGetApplicationContext() 513 @*/ 514 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetApplicationContext(SNES snes,void *usrP) 515 { 516 PetscFunctionBegin; 517 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 518 snes->user = usrP; 519 PetscFunctionReturn(0); 520 } 521 522 #undef __FUNCT__ 523 #define __FUNCT__ "SNESGetApplicationContext" 524 /*@C 525 SNESGetApplicationContext - Gets the user-defined context for the 526 nonlinear solvers. 527 528 Not Collective 529 530 Input Parameter: 531 . snes - SNES context 532 533 Output Parameter: 534 . usrP - user context 535 536 Level: intermediate 537 538 .keywords: SNES, nonlinear, get, application, context 539 540 .seealso: SNESSetApplicationContext() 541 @*/ 542 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetApplicationContext(SNES snes,void **usrP) 543 { 544 PetscFunctionBegin; 545 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 546 *usrP = snes->user; 547 PetscFunctionReturn(0); 548 } 549 550 #undef __FUNCT__ 551 #define __FUNCT__ "SNESGetIterationNumber" 552 /*@ 553 SNESGetIterationNumber - Gets the number of nonlinear iterations completed 554 at this time. 555 556 Not Collective 557 558 Input Parameter: 559 . snes - SNES context 560 561 Output Parameter: 562 . iter - iteration number 563 564 Notes: 565 For example, during the computation of iteration 2 this would return 1. 566 567 This is useful for using lagged Jacobians (where one does not recompute the 568 Jacobian at each SNES iteration). For example, the code 569 .vb 570 ierr = SNESGetIterationNumber(snes,&it); 571 if (!(it % 2)) { 572 [compute Jacobian here] 573 } 574 .ve 575 can be used in your ComputeJacobian() function to cause the Jacobian to be 576 recomputed every second SNES iteration. 577 578 Level: intermediate 579 580 .keywords: SNES, nonlinear, get, iteration, number, 581 582 .seealso: SNESGetFunctionNorm(), SNESGetLinearSolveIterations() 583 @*/ 584 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetIterationNumber(SNES snes,PetscInt* iter) 585 { 586 PetscFunctionBegin; 587 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 588 PetscValidIntPointer(iter,2); 589 *iter = snes->iter; 590 PetscFunctionReturn(0); 591 } 592 593 #undef __FUNCT__ 594 #define __FUNCT__ "SNESGetFunctionNorm" 595 /*@ 596 SNESGetFunctionNorm - Gets the norm of the current function that was set 597 with SNESSSetFunction(). 598 599 Collective on SNES 600 601 Input Parameter: 602 . snes - SNES context 603 604 Output Parameter: 605 . fnorm - 2-norm of function 606 607 Level: intermediate 608 609 .keywords: SNES, nonlinear, get, function, norm 610 611 .seealso: SNESGetFunction(), SNESGetIterationNumber(), SNESGetLinearSolveIterations() 612 @*/ 613 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunctionNorm(SNES snes,PetscReal *fnorm) 614 { 615 PetscFunctionBegin; 616 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 617 PetscValidScalarPointer(fnorm,2); 618 *fnorm = snes->norm; 619 PetscFunctionReturn(0); 620 } 621 622 #undef __FUNCT__ 623 #define __FUNCT__ "SNESGetNonlinearStepFailures" 624 /*@ 625 SNESGetNonlinearStepFailures - Gets the number of unsuccessful steps 626 attempted by the nonlinear solver. 627 628 Not Collective 629 630 Input Parameter: 631 . snes - SNES context 632 633 Output Parameter: 634 . nfails - number of unsuccessful steps attempted 635 636 Notes: 637 This counter is reset to zero for each successive call to SNESSolve(). 638 639 Level: intermediate 640 641 .keywords: SNES, nonlinear, get, number, unsuccessful, steps 642 643 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), 644 SNESSetMaxNonlinearStepFailures(), SNESGetMaxNonlinearStepFailures() 645 @*/ 646 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetNonlinearStepFailures(SNES snes,PetscInt* nfails) 647 { 648 PetscFunctionBegin; 649 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 650 PetscValidIntPointer(nfails,2); 651 *nfails = snes->numFailures; 652 PetscFunctionReturn(0); 653 } 654 655 #undef __FUNCT__ 656 #define __FUNCT__ "SNESSetMaxNonlinearStepFailures" 657 /*@ 658 SNESSetMaxNonlinearStepFailures - Sets the maximum number of unsuccessful steps 659 attempted by the nonlinear solver before it gives up. 660 661 Not Collective 662 663 Input Parameters: 664 + snes - SNES context 665 - maxFails - maximum of unsuccessful steps 666 667 Level: intermediate 668 669 .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps 670 671 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), 672 SNESGetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures() 673 @*/ 674 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMaxNonlinearStepFailures(SNES snes, PetscInt maxFails) 675 { 676 PetscFunctionBegin; 677 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 678 snes->maxFailures = maxFails; 679 PetscFunctionReturn(0); 680 } 681 682 #undef __FUNCT__ 683 #define __FUNCT__ "SNESGetMaxNonlinearStepFailures" 684 /*@ 685 SNESGetMaxNonlinearStepFailures - Gets the maximum number of unsuccessful steps 686 attempted by the nonlinear solver before it gives up. 687 688 Not Collective 689 690 Input Parameter: 691 . snes - SNES context 692 693 Output Parameter: 694 . maxFails - maximum of unsuccessful steps 695 696 Level: intermediate 697 698 .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps 699 700 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), 701 SNESSetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures() 702 703 @*/ 704 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetMaxNonlinearStepFailures(SNES snes, PetscInt *maxFails) 705 { 706 PetscFunctionBegin; 707 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 708 PetscValidIntPointer(maxFails,2); 709 *maxFails = snes->maxFailures; 710 PetscFunctionReturn(0); 711 } 712 713 #undef __FUNCT__ 714 #define __FUNCT__ "SNESGetNumberFunctionEvals" 715 /*@ 716 SNESGetNumberFunctionEvals - Gets the number of user provided function evaluations 717 done by SNES. 718 719 Not Collective 720 721 Input Parameter: 722 . snes - SNES context 723 724 Output Parameter: 725 . nfuncs - number of evaluations 726 727 Level: intermediate 728 729 .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps 730 731 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures() 732 @*/ 733 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetNumberFunctionEvals(SNES snes, PetscInt *nfuncs) 734 { 735 PetscFunctionBegin; 736 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 737 PetscValidIntPointer(nfuncs,2); 738 *nfuncs = snes->nfuncs; 739 PetscFunctionReturn(0); 740 } 741 742 #undef __FUNCT__ 743 #define __FUNCT__ "SNESGetLinearSolveFailures" 744 /*@ 745 SNESGetLinearSolveFailures - Gets the number of failed (non-converged) 746 linear solvers. 747 748 Not Collective 749 750 Input Parameter: 751 . snes - SNES context 752 753 Output Parameter: 754 . nfails - number of failed solves 755 756 Notes: 757 This counter is reset to zero for each successive call to SNESSolve(). 758 759 Level: intermediate 760 761 .keywords: SNES, nonlinear, get, number, unsuccessful, steps 762 763 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures() 764 @*/ 765 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLinearSolveFailures(SNES snes,PetscInt* nfails) 766 { 767 PetscFunctionBegin; 768 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 769 PetscValidIntPointer(nfails,2); 770 *nfails = snes->numLinearSolveFailures; 771 PetscFunctionReturn(0); 772 } 773 774 #undef __FUNCT__ 775 #define __FUNCT__ "SNESSetMaxLinearSolveFailures" 776 /*@ 777 SNESSetMaxLinearSolveFailures - the number of failed linear solve attempts 778 allowed before SNES returns with a diverged reason of SNES_DIVERGED_LINEAR_SOLVE 779 780 Logically Collective on SNES 781 782 Input Parameters: 783 + snes - SNES context 784 - maxFails - maximum allowed linear solve failures 785 786 Level: intermediate 787 788 Notes: By default this is 0; that is SNES returns on the first failed linear solve 789 790 .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps 791 792 .seealso: SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations() 793 @*/ 794 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMaxLinearSolveFailures(SNES snes, PetscInt maxFails) 795 { 796 PetscFunctionBegin; 797 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 798 PetscValidLogicalCollectiveInt(snes,maxFails,2); 799 snes->maxLinearSolveFailures = maxFails; 800 PetscFunctionReturn(0); 801 } 802 803 #undef __FUNCT__ 804 #define __FUNCT__ "SNESGetMaxLinearSolveFailures" 805 /*@ 806 SNESGetMaxLinearSolveFailures - gets the maximum number of linear solve failures that 807 are allowed before SNES terminates 808 809 Not Collective 810 811 Input Parameter: 812 . snes - SNES context 813 814 Output Parameter: 815 . maxFails - maximum of unsuccessful solves allowed 816 817 Level: intermediate 818 819 Notes: By default this is 1; that is SNES returns on the first failed linear solve 820 821 .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps 822 823 .seealso: SNESGetLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), 824 @*/ 825 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetMaxLinearSolveFailures(SNES snes, PetscInt *maxFails) 826 { 827 PetscFunctionBegin; 828 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 829 PetscValidIntPointer(maxFails,2); 830 *maxFails = snes->maxLinearSolveFailures; 831 PetscFunctionReturn(0); 832 } 833 834 #undef __FUNCT__ 835 #define __FUNCT__ "SNESGetLinearSolveIterations" 836 /*@ 837 SNESGetLinearSolveIterations - Gets the total number of linear iterations 838 used by the nonlinear solver. 839 840 Not Collective 841 842 Input Parameter: 843 . snes - SNES context 844 845 Output Parameter: 846 . lits - number of linear iterations 847 848 Notes: 849 This counter is reset to zero for each successive call to SNESSolve(). 850 851 Level: intermediate 852 853 .keywords: SNES, nonlinear, get, number, linear, iterations 854 855 .seealso: SNESGetIterationNumber(), SNESGetFunctionNorm(), SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures() 856 @*/ 857 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLinearSolveIterations(SNES snes,PetscInt* lits) 858 { 859 PetscFunctionBegin; 860 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 861 PetscValidIntPointer(lits,2); 862 *lits = snes->linear_its; 863 PetscFunctionReturn(0); 864 } 865 866 #undef __FUNCT__ 867 #define __FUNCT__ "SNESGetKSP" 868 /*@ 869 SNESGetKSP - Returns the KSP context for a SNES solver. 870 871 Not Collective, but if SNES object is parallel, then KSP object is parallel 872 873 Input Parameter: 874 . snes - the SNES context 875 876 Output Parameter: 877 . ksp - the KSP context 878 879 Notes: 880 The user can then directly manipulate the KSP context to set various 881 options, etc. Likewise, the user can then extract and manipulate the 882 PC contexts as well. 883 884 Level: beginner 885 886 .keywords: SNES, nonlinear, get, KSP, context 887 888 .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP() 889 @*/ 890 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetKSP(SNES snes,KSP *ksp) 891 { 892 PetscErrorCode ierr; 893 894 PetscFunctionBegin; 895 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 896 PetscValidPointer(ksp,2); 897 898 if (!snes->ksp) { 899 ierr = KSPCreate(((PetscObject)snes)->comm,&snes->ksp);CHKERRQ(ierr); 900 ierr = PetscObjectIncrementTabLevel((PetscObject)snes->ksp,(PetscObject)snes,1);CHKERRQ(ierr); 901 ierr = PetscLogObjectParent(snes,snes->ksp);CHKERRQ(ierr); 902 } 903 *ksp = snes->ksp; 904 PetscFunctionReturn(0); 905 } 906 907 #undef __FUNCT__ 908 #define __FUNCT__ "SNESSetKSP" 909 /*@ 910 SNESSetKSP - Sets a KSP context for the SNES object to use 911 912 Not Collective, but the SNES and KSP objects must live on the same MPI_Comm 913 914 Input Parameters: 915 + snes - the SNES context 916 - ksp - the KSP context 917 918 Notes: 919 The SNES object already has its KSP object, you can obtain with SNESGetKSP() 920 so this routine is rarely needed. 921 922 The KSP object that is already in the SNES object has its reference count 923 decreased by one. 924 925 Level: developer 926 927 .keywords: SNES, nonlinear, get, KSP, context 928 929 .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP() 930 @*/ 931 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetKSP(SNES snes,KSP ksp) 932 { 933 PetscErrorCode ierr; 934 935 PetscFunctionBegin; 936 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 937 PetscValidHeaderSpecific(ksp,KSP_CLASSID,2); 938 PetscCheckSameComm(snes,1,ksp,2); 939 ierr = PetscObjectReference((PetscObject)ksp);CHKERRQ(ierr); 940 if (snes->ksp) {ierr = PetscObjectDereference((PetscObject)snes->ksp);CHKERRQ(ierr);} 941 snes->ksp = ksp; 942 PetscFunctionReturn(0); 943 } 944 945 #if 0 946 #undef __FUNCT__ 947 #define __FUNCT__ "SNESPublish_Petsc" 948 static PetscErrorCode SNESPublish_Petsc(PetscObject obj) 949 { 950 PetscFunctionBegin; 951 PetscFunctionReturn(0); 952 } 953 #endif 954 955 /* -----------------------------------------------------------*/ 956 #undef __FUNCT__ 957 #define __FUNCT__ "SNESCreate" 958 /*@ 959 SNESCreate - Creates a nonlinear solver context. 960 961 Collective on MPI_Comm 962 963 Input Parameters: 964 . comm - MPI communicator 965 966 Output Parameter: 967 . outsnes - the new SNES context 968 969 Options Database Keys: 970 + -snes_mf - Activates default matrix-free Jacobian-vector products, 971 and no preconditioning matrix 972 . -snes_mf_operator - Activates default matrix-free Jacobian-vector 973 products, and a user-provided preconditioning matrix 974 as set by SNESSetJacobian() 975 - -snes_fd - Uses (slow!) finite differences to compute Jacobian 976 977 Level: beginner 978 979 .keywords: SNES, nonlinear, create, context 980 981 .seealso: SNESSolve(), SNESDestroy(), SNES, SNESSetLagPreconditioner() 982 983 @*/ 984 PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate(MPI_Comm comm,SNES *outsnes) 985 { 986 PetscErrorCode ierr; 987 SNES snes; 988 SNESKSPEW *kctx; 989 990 PetscFunctionBegin; 991 PetscValidPointer(outsnes,2); 992 *outsnes = PETSC_NULL; 993 #ifndef PETSC_USE_DYNAMIC_LIBRARIES 994 ierr = SNESInitializePackage(PETSC_NULL);CHKERRQ(ierr); 995 #endif 996 997 ierr = PetscHeaderCreate(snes,_p_SNES,struct _SNESOps,SNES_CLASSID,0,"SNES",comm,SNESDestroy,SNESView);CHKERRQ(ierr); 998 999 snes->ops->converged = SNESDefaultConverged; 1000 snes->max_its = 50; 1001 snes->max_funcs = 10000; 1002 snes->norm = 0.0; 1003 snes->rtol = 1.e-8; 1004 snes->ttol = 0.0; 1005 snes->abstol = 1.e-50; 1006 snes->xtol = 1.e-8; 1007 snes->deltatol = 1.e-12; 1008 snes->nfuncs = 0; 1009 snes->numFailures = 0; 1010 snes->maxFailures = 1; 1011 snes->linear_its = 0; 1012 snes->lagjacobian = 1; 1013 snes->lagpreconditioner = 1; 1014 snes->numbermonitors = 0; 1015 snes->data = 0; 1016 snes->setupcalled = PETSC_FALSE; 1017 snes->ksp_ewconv = PETSC_FALSE; 1018 snes->vwork = 0; 1019 snes->nwork = 0; 1020 snes->conv_hist_len = 0; 1021 snes->conv_hist_max = 0; 1022 snes->conv_hist = PETSC_NULL; 1023 snes->conv_hist_its = PETSC_NULL; 1024 snes->conv_hist_reset = PETSC_TRUE; 1025 snes->reason = SNES_CONVERGED_ITERATING; 1026 1027 snes->numLinearSolveFailures = 0; 1028 snes->maxLinearSolveFailures = 1; 1029 1030 /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 1031 ierr = PetscNewLog(snes,SNESKSPEW,&kctx);CHKERRQ(ierr); 1032 snes->kspconvctx = (void*)kctx; 1033 kctx->version = 2; 1034 kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 1035 this was too large for some test cases */ 1036 kctx->rtol_last = 0.0; 1037 kctx->rtol_max = .9; 1038 kctx->gamma = 1.0; 1039 kctx->alpha = .5*(1.0 + sqrt(5.0)); 1040 kctx->alpha2 = kctx->alpha; 1041 kctx->threshold = .1; 1042 kctx->lresid_last = 0.0; 1043 kctx->norm_last = 0.0; 1044 1045 *outsnes = snes; 1046 ierr = PetscPublishAll(snes);CHKERRQ(ierr); 1047 PetscFunctionReturn(0); 1048 } 1049 1050 #undef __FUNCT__ 1051 #define __FUNCT__ "SNESSetFunction" 1052 /*@C 1053 SNESSetFunction - Sets the function evaluation routine and function 1054 vector for use by the SNES routines in solving systems of nonlinear 1055 equations. 1056 1057 Logically Collective on SNES 1058 1059 Input Parameters: 1060 + snes - the SNES context 1061 . r - vector to store function value 1062 . func - function evaluation routine 1063 - ctx - [optional] user-defined context for private data for the 1064 function evaluation routine (may be PETSC_NULL) 1065 1066 Calling sequence of func: 1067 $ func (SNES snes,Vec x,Vec f,void *ctx); 1068 1069 . f - function vector 1070 - ctx - optional user-defined function context 1071 1072 Notes: 1073 The Newton-like methods typically solve linear systems of the form 1074 $ f'(x) x = -f(x), 1075 where f'(x) denotes the Jacobian matrix and f(x) is the function. 1076 1077 Level: beginner 1078 1079 .keywords: SNES, nonlinear, set, function 1080 1081 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 1082 @*/ 1083 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx) 1084 { 1085 PetscErrorCode ierr; 1086 PetscFunctionBegin; 1087 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1088 if (r) PetscValidHeaderSpecific(r,VEC_CLASSID,2); 1089 if (r) PetscCheckSameComm(snes,1,r,2); 1090 if (!r && snes->dm) { 1091 ierr = DMCreateGlobalVector(snes->dm,&r);CHKERRQ(ierr); 1092 } else { 1093 ierr = PetscObjectReference((PetscObject)r);CHKERRQ(ierr); 1094 } 1095 if (snes->vec_func) { ierr = VecDestroy(snes->vec_func);CHKERRQ(ierr); } 1096 snes->ops->computefunction = func; 1097 snes->vec_func = r; 1098 snes->funP = ctx; 1099 PetscFunctionReturn(0); 1100 } 1101 1102 /* --------------------------------------------------------------- */ 1103 #undef __FUNCT__ 1104 #define __FUNCT__ "SNESGetRhs" 1105 /*@C 1106 SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set 1107 it assumes a zero right hand side. 1108 1109 Logically Collective on SNES 1110 1111 Input Parameter: 1112 . snes - the SNES context 1113 1114 Output Parameter: 1115 . rhs - the right hand side vector or PETSC_NULL if the right hand side vector is null 1116 1117 Level: intermediate 1118 1119 .keywords: SNES, nonlinear, get, function, right hand side 1120 1121 .seealso: SNESGetSolution(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction() 1122 @*/ 1123 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetRhs(SNES snes,Vec *rhs) 1124 { 1125 PetscFunctionBegin; 1126 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1127 PetscValidPointer(rhs,2); 1128 *rhs = snes->vec_rhs; 1129 PetscFunctionReturn(0); 1130 } 1131 1132 #undef __FUNCT__ 1133 #define __FUNCT__ "SNESComputeFunction" 1134 /*@ 1135 SNESComputeFunction - Calls the function that has been set with 1136 SNESSetFunction(). 1137 1138 Collective on SNES 1139 1140 Input Parameters: 1141 + snes - the SNES context 1142 - x - input vector 1143 1144 Output Parameter: 1145 . y - function vector, as set by SNESSetFunction() 1146 1147 Notes: 1148 SNESComputeFunction() is typically used within nonlinear solvers 1149 implementations, so most users would not generally call this routine 1150 themselves. 1151 1152 Level: developer 1153 1154 .keywords: SNES, nonlinear, compute, function 1155 1156 .seealso: SNESSetFunction(), SNESGetFunction() 1157 @*/ 1158 PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeFunction(SNES snes,Vec x,Vec y) 1159 { 1160 PetscErrorCode ierr; 1161 1162 PetscFunctionBegin; 1163 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1164 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 1165 PetscValidHeaderSpecific(y,VEC_CLASSID,3); 1166 PetscCheckSameComm(snes,1,x,2); 1167 PetscCheckSameComm(snes,1,y,3); 1168 1169 ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 1170 if (snes->ops->computefunction) { 1171 PetscStackPush("SNES user function"); 1172 ierr = (*snes->ops->computefunction)(snes,x,y,snes->funP); 1173 PetscStackPop; 1174 if (PetscExceptionValue(ierr)) { 1175 PetscErrorCode pierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(pierr); 1176 } 1177 CHKERRQ(ierr); 1178 } else if (snes->vec_rhs) { 1179 ierr = MatMult(snes->jacobian, x, y);CHKERRQ(ierr); 1180 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() before SNESComputeFunction(), likely called from SNESSolve()."); 1181 if (snes->vec_rhs) { 1182 ierr = VecAXPY(y,-1.0,snes->vec_rhs);CHKERRQ(ierr); 1183 } 1184 snes->nfuncs++; 1185 ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 1186 PetscFunctionReturn(0); 1187 } 1188 1189 #undef __FUNCT__ 1190 #define __FUNCT__ "SNESComputeJacobian" 1191 /*@ 1192 SNESComputeJacobian - Computes the Jacobian matrix that has been 1193 set with SNESSetJacobian(). 1194 1195 Collective on SNES and Mat 1196 1197 Input Parameters: 1198 + snes - the SNES context 1199 - x - input vector 1200 1201 Output Parameters: 1202 + A - Jacobian matrix 1203 . B - optional preconditioning matrix 1204 - flag - flag indicating matrix structure (one of, SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER) 1205 1206 Options Database Keys: 1207 + -snes_lag_preconditioner <lag> 1208 - -snes_lag_jacobian <lag> 1209 1210 Notes: 1211 Most users should not need to explicitly call this routine, as it 1212 is used internally within the nonlinear solvers. 1213 1214 See KSPSetOperators() for important information about setting the 1215 flag parameter. 1216 1217 Level: developer 1218 1219 .keywords: SNES, compute, Jacobian, matrix 1220 1221 .seealso: SNESSetJacobian(), KSPSetOperators(), MatStructure, SNESSetLagPreconditioner(), SNESSetLagJacobian() 1222 @*/ 1223 PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 1224 { 1225 PetscErrorCode ierr; 1226 PetscTruth flag; 1227 1228 PetscFunctionBegin; 1229 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1230 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 1231 PetscValidPointer(flg,5); 1232 PetscCheckSameComm(snes,1,X,2); 1233 if (!snes->ops->computejacobian) PetscFunctionReturn(0); 1234 1235 /* make sure that MatAssemblyBegin/End() is called on A matrix if it is matrix free */ 1236 1237 if (snes->lagjacobian == -2) { 1238 snes->lagjacobian = -1; 1239 ierr = PetscInfo(snes,"Recomputing Jacobian/preconditioner because lag is -2 (means compute Jacobian, but then never again) \n");CHKERRQ(ierr); 1240 } else if (snes->lagjacobian == -1) { 1241 *flg = SAME_PRECONDITIONER; 1242 ierr = PetscInfo(snes,"Reusing Jacobian/preconditioner because lag is -1\n");CHKERRQ(ierr); 1243 ierr = PetscTypeCompare((PetscObject)*A,MATMFFD,&flag);CHKERRQ(ierr); 1244 if (flag) { 1245 ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1246 ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1247 } 1248 PetscFunctionReturn(0); 1249 } else if (snes->lagjacobian > 1 && snes->iter % snes->lagjacobian) { 1250 *flg = SAME_PRECONDITIONER; 1251 ierr = PetscInfo2(snes,"Reusing Jacobian/preconditioner because lag is %D and SNES iteration is %D\n",snes->lagjacobian,snes->iter);CHKERRQ(ierr); 1252 ierr = PetscTypeCompare((PetscObject)*A,MATMFFD,&flag);CHKERRQ(ierr); 1253 if (flag) { 1254 ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1255 ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1256 } 1257 PetscFunctionReturn(0); 1258 } 1259 1260 *flg = DIFFERENT_NONZERO_PATTERN; 1261 ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 1262 PetscStackPush("SNES user Jacobian function"); 1263 ierr = (*snes->ops->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr); 1264 PetscStackPop; 1265 ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 1266 1267 if (snes->lagpreconditioner == -2) { 1268 ierr = PetscInfo(snes,"Rebuilding preconditioner exactly once since lag is -2\n");CHKERRQ(ierr); 1269 snes->lagpreconditioner = -1; 1270 } else if (snes->lagpreconditioner == -1) { 1271 *flg = SAME_PRECONDITIONER; 1272 ierr = PetscInfo(snes,"Reusing preconditioner because lag is -1\n");CHKERRQ(ierr); 1273 } else if (snes->lagpreconditioner > 1 && snes->iter % snes->lagpreconditioner) { 1274 *flg = SAME_PRECONDITIONER; 1275 ierr = PetscInfo2(snes,"Reusing preconditioner because lag is %D and SNES iteration is %D\n",snes->lagpreconditioner,snes->iter);CHKERRQ(ierr); 1276 } 1277 1278 /* make sure user returned a correct Jacobian and preconditioner */ 1279 /* PetscValidHeaderSpecific(*A,MAT_CLASSID,3); 1280 PetscValidHeaderSpecific(*B,MAT_CLASSID,4); */ 1281 PetscFunctionReturn(0); 1282 } 1283 1284 #undef __FUNCT__ 1285 #define __FUNCT__ "SNESSetJacobian" 1286 /*@C 1287 SNESSetJacobian - Sets the function to compute Jacobian as well as the 1288 location to store the matrix. 1289 1290 Logically Collective on SNES and Mat 1291 1292 Input Parameters: 1293 + snes - the SNES context 1294 . A - Jacobian matrix 1295 . B - preconditioner matrix (usually same as the Jacobian) 1296 . func - Jacobian evaluation routine 1297 - ctx - [optional] user-defined context for private data for the 1298 Jacobian evaluation routine (may be PETSC_NULL) 1299 1300 Calling sequence of func: 1301 $ func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 1302 1303 + x - input vector 1304 . A - Jacobian matrix 1305 . B - preconditioner matrix, usually the same as A 1306 . flag - flag indicating information about the preconditioner matrix 1307 structure (same as flag in KSPSetOperators()), one of SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER 1308 - ctx - [optional] user-defined Jacobian context 1309 1310 Notes: 1311 See KSPSetOperators() for important information about setting the flag 1312 output parameter in the routine func(). Be sure to read this information! 1313 1314 The routine func() takes Mat * as the matrix arguments rather than Mat. 1315 This allows the Jacobian evaluation routine to replace A and/or B with a 1316 completely new new matrix structure (not just different matrix elements) 1317 when appropriate, for instance, if the nonzero structure is changing 1318 throughout the global iterations. 1319 1320 If the A matrix and B matrix are different you must call MatAssemblyBegin/End() on 1321 each matrix. 1322 1323 Level: beginner 1324 1325 .keywords: SNES, nonlinear, set, Jacobian, matrix 1326 1327 .seealso: KSPSetOperators(), SNESSetFunction(), MatMFFDComputeJacobian(), SNESDefaultComputeJacobianColor(), MatStructure 1328 @*/ 1329 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 1330 { 1331 PetscErrorCode ierr; 1332 1333 PetscFunctionBegin; 1334 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1335 if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2); 1336 if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3); 1337 if (A) PetscCheckSameComm(snes,1,A,2); 1338 if (B) PetscCheckSameComm(snes,1,B,3); 1339 if (func) snes->ops->computejacobian = func; 1340 if (ctx) snes->jacP = ctx; 1341 if (A) { 1342 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 1343 if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 1344 snes->jacobian = A; 1345 } 1346 if (B) { 1347 ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 1348 if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 1349 snes->jacobian_pre = B; 1350 } 1351 PetscFunctionReturn(0); 1352 } 1353 1354 #undef __FUNCT__ 1355 #define __FUNCT__ "SNESGetJacobian" 1356 /*@C 1357 SNESGetJacobian - Returns the Jacobian matrix and optionally the user 1358 provided context for evaluating the Jacobian. 1359 1360 Not Collective, but Mat object will be parallel if SNES object is 1361 1362 Input Parameter: 1363 . snes - the nonlinear solver context 1364 1365 Output Parameters: 1366 + A - location to stash Jacobian matrix (or PETSC_NULL) 1367 . B - location to stash preconditioner matrix (or PETSC_NULL) 1368 . func - location to put Jacobian function (or PETSC_NULL) 1369 - ctx - location to stash Jacobian ctx (or PETSC_NULL) 1370 1371 Level: advanced 1372 1373 .seealso: SNESSetJacobian(), SNESComputeJacobian() 1374 @*/ 1375 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetJacobian(SNES snes,Mat *A,Mat *B,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx) 1376 { 1377 PetscFunctionBegin; 1378 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1379 if (A) *A = snes->jacobian; 1380 if (B) *B = snes->jacobian_pre; 1381 if (func) *func = snes->ops->computejacobian; 1382 if (ctx) *ctx = snes->jacP; 1383 PetscFunctionReturn(0); 1384 } 1385 1386 /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 1387 1388 #undef __FUNCT__ 1389 #define __FUNCT__ "SNESSetUp" 1390 /*@ 1391 SNESSetUp - Sets up the internal data structures for the later use 1392 of a nonlinear solver. 1393 1394 Collective on SNES 1395 1396 Input Parameters: 1397 . snes - the SNES context 1398 1399 Notes: 1400 For basic use of the SNES solvers the user need not explicitly call 1401 SNESSetUp(), since these actions will automatically occur during 1402 the call to SNESSolve(). However, if one wishes to control this 1403 phase separately, SNESSetUp() should be called after SNESCreate() 1404 and optional routines of the form SNESSetXXX(), but before SNESSolve(). 1405 1406 Level: advanced 1407 1408 .keywords: SNES, nonlinear, setup 1409 1410 .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 1411 @*/ 1412 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUp(SNES snes) 1413 { 1414 PetscErrorCode ierr; 1415 1416 PetscFunctionBegin; 1417 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1418 if (snes->setupcalled) PetscFunctionReturn(0); 1419 1420 if (!((PetscObject)snes)->type_name) { 1421 ierr = SNESSetType(snes,SNESLS);CHKERRQ(ierr); 1422 } 1423 1424 if (!snes->vec_func && !snes->vec_rhs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 1425 if (!snes->ops->computefunction && !snes->vec_rhs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 1426 if (snes->vec_func == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector"); 1427 if (!snes->ops->computejacobian && snes->dm) { 1428 Mat J; 1429 ISColoring coloring; 1430 MatFDColoring fd; 1431 1432 ierr = DMGetMatrix(snes->dm,MATAIJ,&J);CHKERRQ(ierr); 1433 ierr = DMGetColoring(snes->dm,IS_COLORING_GHOSTED,MATAIJ,&coloring);CHKERRQ(ierr); 1434 ierr = MatFDColoringCreate(J,coloring,&fd);CHKERRQ(ierr); 1435 ierr = MatFDColoringSetFunction(fd,(PetscErrorCode (*)(void))snes->ops->computefunction,snes->funP);CHKERRQ(ierr); 1436 ierr = SNESSetJacobian(snes,J,J,SNESDefaultComputeJacobianColor,fd);CHKERRQ(ierr); 1437 ierr = ISColoringDestroy(coloring);CHKERRQ(ierr); 1438 } 1439 if (!snes->ksp) {ierr = SNESGetKSP(snes, &snes->ksp);CHKERRQ(ierr);} 1440 1441 if (snes->ops->setup) { 1442 ierr = (*snes->ops->setup)(snes);CHKERRQ(ierr); 1443 } 1444 snes->setupcalled = PETSC_TRUE; 1445 PetscFunctionReturn(0); 1446 } 1447 1448 #undef __FUNCT__ 1449 #define __FUNCT__ "SNESDestroy" 1450 /*@ 1451 SNESDestroy - Destroys the nonlinear solver context that was created 1452 with SNESCreate(). 1453 1454 Collective on SNES 1455 1456 Input Parameter: 1457 . snes - the SNES context 1458 1459 Level: beginner 1460 1461 .keywords: SNES, nonlinear, destroy 1462 1463 .seealso: SNESCreate(), SNESSolve() 1464 @*/ 1465 PetscErrorCode PETSCSNES_DLLEXPORT SNESDestroy(SNES snes) 1466 { 1467 PetscErrorCode ierr; 1468 1469 PetscFunctionBegin; 1470 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1471 if (--((PetscObject)snes)->refct > 0) PetscFunctionReturn(0); 1472 1473 /* if memory was published with AMS then destroy it */ 1474 ierr = PetscObjectDepublish(snes);CHKERRQ(ierr); 1475 1476 if (snes->dm) {ierr = DMDestroy(snes->dm);CHKERRQ(ierr);} 1477 if (snes->ops->destroy) {ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr);} 1478 1479 if (snes->vec_rhs) {ierr = VecDestroy(snes->vec_rhs);CHKERRQ(ierr);} 1480 if (snes->vec_sol) {ierr = VecDestroy(snes->vec_sol);CHKERRQ(ierr);} 1481 if (snes->vec_func) {ierr = VecDestroy(snes->vec_func);CHKERRQ(ierr);} 1482 if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 1483 if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 1484 if (snes->ksp) {ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr);} 1485 ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr); 1486 if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);} 1487 ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 1488 if (snes->ops->convergeddestroy) {ierr = (*snes->ops->convergeddestroy)(snes->cnvP);CHKERRQ(ierr);} 1489 ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr); 1490 PetscFunctionReturn(0); 1491 } 1492 1493 /* ----------- Routines to set solver parameters ---------- */ 1494 1495 #undef __FUNCT__ 1496 #define __FUNCT__ "SNESSetLagPreconditioner" 1497 /*@ 1498 SNESSetLagPreconditioner - Determines when the preconditioner is rebuilt in the nonlinear solve. 1499 1500 Logically Collective on SNES 1501 1502 Input Parameters: 1503 + snes - the SNES context 1504 - lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time 1505 the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that 1506 1507 Options Database Keys: 1508 . -snes_lag_preconditioner <lag> 1509 1510 Notes: 1511 The default is 1 1512 The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 1513 If -1 is used before the very first nonlinear solve the preconditioner is still built because there is no previous preconditioner to use 1514 1515 Level: intermediate 1516 1517 .keywords: SNES, nonlinear, set, convergence, tolerances 1518 1519 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian() 1520 1521 @*/ 1522 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLagPreconditioner(SNES snes,PetscInt lag) 1523 { 1524 PetscFunctionBegin; 1525 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1526 if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater"); 1527 if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0"); 1528 PetscValidLogicalCollectiveInt(snes,lag,2); 1529 snes->lagpreconditioner = lag; 1530 PetscFunctionReturn(0); 1531 } 1532 1533 #undef __FUNCT__ 1534 #define __FUNCT__ "SNESGetLagPreconditioner" 1535 /*@ 1536 SNESGetLagPreconditioner - Indicates how often the preconditioner is rebuilt 1537 1538 Not Collective 1539 1540 Input Parameter: 1541 . snes - the SNES context 1542 1543 Output Parameter: 1544 . lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time 1545 the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that 1546 1547 Options Database Keys: 1548 . -snes_lag_preconditioner <lag> 1549 1550 Notes: 1551 The default is 1 1552 The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 1553 1554 Level: intermediate 1555 1556 .keywords: SNES, nonlinear, set, convergence, tolerances 1557 1558 .seealso: SNESSetTrustRegionTolerance(), SNESSetLagPreconditioner() 1559 1560 @*/ 1561 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLagPreconditioner(SNES snes,PetscInt *lag) 1562 { 1563 PetscFunctionBegin; 1564 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1565 *lag = snes->lagpreconditioner; 1566 PetscFunctionReturn(0); 1567 } 1568 1569 #undef __FUNCT__ 1570 #define __FUNCT__ "SNESSetLagJacobian" 1571 /*@ 1572 SNESSetLagJacobian - Determines when the Jacobian is rebuilt in the nonlinear solve. See SNESSetLagPreconditioner() for determining how 1573 often the preconditioner is rebuilt. 1574 1575 Logically Collective on SNES 1576 1577 Input Parameters: 1578 + snes - the SNES context 1579 - lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time 1580 the Jacobian is built etc. -2 means rebuild at next chance but then never again 1581 1582 Options Database Keys: 1583 . -snes_lag_jacobian <lag> 1584 1585 Notes: 1586 The default is 1 1587 The Jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 1588 If -1 is used before the very first nonlinear solve the CODE WILL FAIL! because no Jacobian is used, use -2 to indicate you want it recomputed 1589 at the next Newton step but never again (unless it is reset to another value) 1590 1591 Level: intermediate 1592 1593 .keywords: SNES, nonlinear, set, convergence, tolerances 1594 1595 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagPreconditioner(), SNESGetLagJacobian() 1596 1597 @*/ 1598 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLagJacobian(SNES snes,PetscInt lag) 1599 { 1600 PetscFunctionBegin; 1601 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1602 if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater"); 1603 if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0"); 1604 PetscValidLogicalCollectiveInt(snes,lag,2); 1605 snes->lagjacobian = lag; 1606 PetscFunctionReturn(0); 1607 } 1608 1609 #undef __FUNCT__ 1610 #define __FUNCT__ "SNESGetLagJacobian" 1611 /*@ 1612 SNESGetLagJacobian - Indicates how often the Jacobian is rebuilt. See SNESGetLagPreconditioner() to determine when the preconditioner is rebuilt 1613 1614 Not Collective 1615 1616 Input Parameter: 1617 . snes - the SNES context 1618 1619 Output Parameter: 1620 . lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time 1621 the Jacobian is built etc. 1622 1623 Options Database Keys: 1624 . -snes_lag_jacobian <lag> 1625 1626 Notes: 1627 The default is 1 1628 The jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 1629 1630 Level: intermediate 1631 1632 .keywords: SNES, nonlinear, set, convergence, tolerances 1633 1634 .seealso: SNESSetTrustRegionTolerance(), SNESSetLagJacobian(), SNESSetLagPreconditioner(), SNESGetLagPreconditioner() 1635 1636 @*/ 1637 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLagJacobian(SNES snes,PetscInt *lag) 1638 { 1639 PetscFunctionBegin; 1640 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1641 *lag = snes->lagjacobian; 1642 PetscFunctionReturn(0); 1643 } 1644 1645 #undef __FUNCT__ 1646 #define __FUNCT__ "SNESSetTolerances" 1647 /*@ 1648 SNESSetTolerances - Sets various parameters used in convergence tests. 1649 1650 Logically Collective on SNES 1651 1652 Input Parameters: 1653 + snes - the SNES context 1654 . abstol - absolute convergence tolerance 1655 . rtol - relative convergence tolerance 1656 . stol - convergence tolerance in terms of the norm 1657 of the change in the solution between steps 1658 . maxit - maximum number of iterations 1659 - maxf - maximum number of function evaluations 1660 1661 Options Database Keys: 1662 + -snes_atol <abstol> - Sets abstol 1663 . -snes_rtol <rtol> - Sets rtol 1664 . -snes_stol <stol> - Sets stol 1665 . -snes_max_it <maxit> - Sets maxit 1666 - -snes_max_funcs <maxf> - Sets maxf 1667 1668 Notes: 1669 The default maximum number of iterations is 50. 1670 The default maximum number of function evaluations is 1000. 1671 1672 Level: intermediate 1673 1674 .keywords: SNES, nonlinear, set, convergence, tolerances 1675 1676 .seealso: SNESSetTrustRegionTolerance() 1677 @*/ 1678 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf) 1679 { 1680 PetscFunctionBegin; 1681 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1682 PetscValidLogicalCollectiveReal(snes,abstol,2); 1683 PetscValidLogicalCollectiveReal(snes,rtol,3); 1684 PetscValidLogicalCollectiveReal(snes,stol,4); 1685 PetscValidLogicalCollectiveInt(snes,maxit,5); 1686 PetscValidLogicalCollectiveInt(snes,maxf,6); 1687 1688 if (abstol != PETSC_DEFAULT) snes->abstol = abstol; 1689 if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1690 if (stol != PETSC_DEFAULT) snes->xtol = stol; 1691 if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1692 if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 1693 PetscFunctionReturn(0); 1694 } 1695 1696 #undef __FUNCT__ 1697 #define __FUNCT__ "SNESGetTolerances" 1698 /*@ 1699 SNESGetTolerances - Gets various parameters used in convergence tests. 1700 1701 Not Collective 1702 1703 Input Parameters: 1704 + snes - the SNES context 1705 . atol - absolute convergence tolerance 1706 . rtol - relative convergence tolerance 1707 . stol - convergence tolerance in terms of the norm 1708 of the change in the solution between steps 1709 . maxit - maximum number of iterations 1710 - maxf - maximum number of function evaluations 1711 1712 Notes: 1713 The user can specify PETSC_NULL for any parameter that is not needed. 1714 1715 Level: intermediate 1716 1717 .keywords: SNES, nonlinear, get, convergence, tolerances 1718 1719 .seealso: SNESSetTolerances() 1720 @*/ 1721 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf) 1722 { 1723 PetscFunctionBegin; 1724 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1725 if (atol) *atol = snes->abstol; 1726 if (rtol) *rtol = snes->rtol; 1727 if (stol) *stol = snes->xtol; 1728 if (maxit) *maxit = snes->max_its; 1729 if (maxf) *maxf = snes->max_funcs; 1730 PetscFunctionReturn(0); 1731 } 1732 1733 #undef __FUNCT__ 1734 #define __FUNCT__ "SNESSetTrustRegionTolerance" 1735 /*@ 1736 SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 1737 1738 Logically Collective on SNES 1739 1740 Input Parameters: 1741 + snes - the SNES context 1742 - tol - tolerance 1743 1744 Options Database Key: 1745 . -snes_trtol <tol> - Sets tol 1746 1747 Level: intermediate 1748 1749 .keywords: SNES, nonlinear, set, trust region, tolerance 1750 1751 .seealso: SNESSetTolerances() 1752 @*/ 1753 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTrustRegionTolerance(SNES snes,PetscReal tol) 1754 { 1755 PetscFunctionBegin; 1756 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1757 PetscValidLogicalCollectiveReal(snes,tol,2); 1758 snes->deltatol = tol; 1759 PetscFunctionReturn(0); 1760 } 1761 1762 /* 1763 Duplicate the lg monitors for SNES from KSP; for some reason with 1764 dynamic libraries things don't work under Sun4 if we just use 1765 macros instead of functions 1766 */ 1767 #undef __FUNCT__ 1768 #define __FUNCT__ "SNESMonitorLG" 1769 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLG(SNES snes,PetscInt it,PetscReal norm,void *ctx) 1770 { 1771 PetscErrorCode ierr; 1772 1773 PetscFunctionBegin; 1774 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1775 ierr = KSPMonitorLG((KSP)snes,it,norm,ctx);CHKERRQ(ierr); 1776 PetscFunctionReturn(0); 1777 } 1778 1779 #undef __FUNCT__ 1780 #define __FUNCT__ "SNESMonitorLGCreate" 1781 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1782 { 1783 PetscErrorCode ierr; 1784 1785 PetscFunctionBegin; 1786 ierr = KSPMonitorLGCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr); 1787 PetscFunctionReturn(0); 1788 } 1789 1790 #undef __FUNCT__ 1791 #define __FUNCT__ "SNESMonitorLGDestroy" 1792 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGDestroy(PetscDrawLG draw) 1793 { 1794 PetscErrorCode ierr; 1795 1796 PetscFunctionBegin; 1797 ierr = KSPMonitorLGDestroy(draw);CHKERRQ(ierr); 1798 PetscFunctionReturn(0); 1799 } 1800 1801 extern PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorRange_Private(SNES,PetscInt,PetscReal*); 1802 #undef __FUNCT__ 1803 #define __FUNCT__ "SNESMonitorLGRange" 1804 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGRange(SNES snes,PetscInt n,PetscReal rnorm,void *monctx) 1805 { 1806 PetscDrawLG lg; 1807 PetscErrorCode ierr; 1808 PetscReal x,y,per; 1809 PetscViewer v = (PetscViewer)monctx; 1810 static PetscReal prev; /* should be in the context */ 1811 PetscDraw draw; 1812 PetscFunctionBegin; 1813 if (!monctx) { 1814 MPI_Comm comm; 1815 1816 ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 1817 v = PETSC_VIEWER_DRAW_(comm); 1818 } 1819 ierr = PetscViewerDrawGetDrawLG(v,0,&lg);CHKERRQ(ierr); 1820 if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 1821 ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr); 1822 ierr = PetscDrawSetTitle(draw,"Residual norm");CHKERRQ(ierr); 1823 x = (PetscReal) n; 1824 if (rnorm > 0.0) y = log10(rnorm); else y = -15.0; 1825 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 1826 if (n < 20 || !(n % 5)) { 1827 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 1828 } 1829 1830 ierr = PetscViewerDrawGetDrawLG(v,1,&lg);CHKERRQ(ierr); 1831 if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 1832 ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr); 1833 ierr = PetscDrawSetTitle(draw,"% elemts > .2*max elemt");CHKERRQ(ierr); 1834 ierr = SNESMonitorRange_Private(snes,n,&per);CHKERRQ(ierr); 1835 x = (PetscReal) n; 1836 y = 100.0*per; 1837 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 1838 if (n < 20 || !(n % 5)) { 1839 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 1840 } 1841 1842 ierr = PetscViewerDrawGetDrawLG(v,2,&lg);CHKERRQ(ierr); 1843 if (!n) {prev = rnorm;ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 1844 ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr); 1845 ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm");CHKERRQ(ierr); 1846 x = (PetscReal) n; 1847 y = (prev - rnorm)/prev; 1848 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 1849 if (n < 20 || !(n % 5)) { 1850 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 1851 } 1852 1853 ierr = PetscViewerDrawGetDrawLG(v,3,&lg);CHKERRQ(ierr); 1854 if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 1855 ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr); 1856 ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");CHKERRQ(ierr); 1857 x = (PetscReal) n; 1858 y = (prev - rnorm)/(prev*per); 1859 if (n > 2) { /*skip initial crazy value */ 1860 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 1861 } 1862 if (n < 20 || !(n % 5)) { 1863 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 1864 } 1865 prev = rnorm; 1866 PetscFunctionReturn(0); 1867 } 1868 1869 #undef __FUNCT__ 1870 #define __FUNCT__ "SNESMonitorLGRangeCreate" 1871 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGRangeCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1872 { 1873 PetscErrorCode ierr; 1874 1875 PetscFunctionBegin; 1876 ierr = KSPMonitorLGCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr); 1877 PetscFunctionReturn(0); 1878 } 1879 1880 #undef __FUNCT__ 1881 #define __FUNCT__ "SNESMonitorLGRangeDestroy" 1882 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGRangeDestroy(PetscDrawLG draw) 1883 { 1884 PetscErrorCode ierr; 1885 1886 PetscFunctionBegin; 1887 ierr = KSPMonitorLGDestroy(draw);CHKERRQ(ierr); 1888 PetscFunctionReturn(0); 1889 } 1890 1891 /* ------------ Routines to set performance monitoring options ----------- */ 1892 1893 #undef __FUNCT__ 1894 #define __FUNCT__ "SNESMonitorSet" 1895 /*@C 1896 SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every 1897 iteration of the nonlinear solver to display the iteration's 1898 progress. 1899 1900 Logically Collective on SNES 1901 1902 Input Parameters: 1903 + snes - the SNES context 1904 . func - monitoring routine 1905 . mctx - [optional] user-defined context for private data for the 1906 monitor routine (use PETSC_NULL if no context is desired) 1907 - monitordestroy - [optional] routine that frees monitor context 1908 (may be PETSC_NULL) 1909 1910 Calling sequence of func: 1911 $ int func(SNES snes,PetscInt its, PetscReal norm,void *mctx) 1912 1913 + snes - the SNES context 1914 . its - iteration number 1915 . norm - 2-norm function value (may be estimated) 1916 - mctx - [optional] monitoring context 1917 1918 Options Database Keys: 1919 + -snes_monitor - sets SNESMonitorDefault() 1920 . -snes_monitor_draw - sets line graph monitor, 1921 uses SNESMonitorLGCreate() 1922 _ -snes_monitor_cancel - cancels all monitors that have 1923 been hardwired into a code by 1924 calls to SNESMonitorSet(), but 1925 does not cancel those set via 1926 the options database. 1927 1928 Notes: 1929 Several different monitoring routines may be set by calling 1930 SNESMonitorSet() multiple times; all will be called in the 1931 order in which they were set. 1932 1933 Fortran notes: Only a single monitor function can be set for each SNES object 1934 1935 Level: intermediate 1936 1937 .keywords: SNES, nonlinear, set, monitor 1938 1939 .seealso: SNESMonitorDefault(), SNESMonitorCancel() 1940 @*/ 1941 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorSet(SNES snes,PetscErrorCode (*monitor)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void*)) 1942 { 1943 PetscInt i; 1944 1945 PetscFunctionBegin; 1946 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1947 if (snes->numbermonitors >= MAXSNESMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1948 for (i=0; i<snes->numbermonitors;i++) { 1949 if (monitor == snes->monitor[i] && monitordestroy == snes->monitordestroy[i] && mctx == snes->monitorcontext[i]) PetscFunctionReturn(0); 1950 1951 /* check if both default monitors that share common ASCII viewer */ 1952 if (monitor == snes->monitor[i] && monitor == SNESMonitorDefault) { 1953 if (mctx && snes->monitorcontext[i]) { 1954 PetscErrorCode ierr; 1955 PetscViewerASCIIMonitor viewer1 = (PetscViewerASCIIMonitor) mctx; 1956 PetscViewerASCIIMonitor viewer2 = (PetscViewerASCIIMonitor) snes->monitorcontext[i]; 1957 if (viewer1->viewer == viewer2->viewer) { 1958 ierr = (*monitordestroy)(mctx);CHKERRQ(ierr); 1959 PetscFunctionReturn(0); 1960 } 1961 } 1962 } 1963 } 1964 snes->monitor[snes->numbermonitors] = monitor; 1965 snes->monitordestroy[snes->numbermonitors] = monitordestroy; 1966 snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 1967 PetscFunctionReturn(0); 1968 } 1969 1970 #undef __FUNCT__ 1971 #define __FUNCT__ "SNESMonitorCancel" 1972 /*@C 1973 SNESMonitorCancel - Clears all the monitor functions for a SNES object. 1974 1975 Logically Collective on SNES 1976 1977 Input Parameters: 1978 . snes - the SNES context 1979 1980 Options Database Key: 1981 . -snes_monitor_cancel - cancels all monitors that have been hardwired 1982 into a code by calls to SNESMonitorSet(), but does not cancel those 1983 set via the options database 1984 1985 Notes: 1986 There is no way to clear one specific monitor from a SNES object. 1987 1988 Level: intermediate 1989 1990 .keywords: SNES, nonlinear, set, monitor 1991 1992 .seealso: SNESMonitorDefault(), SNESMonitorSet() 1993 @*/ 1994 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorCancel(SNES snes) 1995 { 1996 PetscErrorCode ierr; 1997 PetscInt i; 1998 1999 PetscFunctionBegin; 2000 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2001 for (i=0; i<snes->numbermonitors; i++) { 2002 if (snes->monitordestroy[i]) { 2003 ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr); 2004 } 2005 } 2006 snes->numbermonitors = 0; 2007 PetscFunctionReturn(0); 2008 } 2009 2010 #undef __FUNCT__ 2011 #define __FUNCT__ "SNESSetConvergenceTest" 2012 /*@C 2013 SNESSetConvergenceTest - Sets the function that is to be used 2014 to test for convergence of the nonlinear iterative solution. 2015 2016 Logically Collective on SNES 2017 2018 Input Parameters: 2019 + snes - the SNES context 2020 . func - routine to test for convergence 2021 . cctx - [optional] context for private data for the convergence routine (may be PETSC_NULL) 2022 - destroy - [optional] destructor for the context (may be PETSC_NULL; PETSC_NULL_FUNCTION in Fortran) 2023 2024 Calling sequence of func: 2025 $ PetscErrorCode func (SNES snes,PetscInt it,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx) 2026 2027 + snes - the SNES context 2028 . it - current iteration (0 is the first and is before any Newton step) 2029 . cctx - [optional] convergence context 2030 . reason - reason for convergence/divergence 2031 . xnorm - 2-norm of current iterate 2032 . gnorm - 2-norm of current step 2033 - f - 2-norm of function 2034 2035 Level: advanced 2036 2037 .keywords: SNES, nonlinear, set, convergence, test 2038 2039 .seealso: SNESDefaultConverged(), SNESSkipConverged() 2040 @*/ 2041 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*)) 2042 { 2043 PetscErrorCode ierr; 2044 2045 PetscFunctionBegin; 2046 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2047 if (!func) func = SNESSkipConverged; 2048 if (snes->ops->convergeddestroy) { 2049 ierr = (*snes->ops->convergeddestroy)(snes->cnvP);CHKERRQ(ierr); 2050 } 2051 snes->ops->converged = func; 2052 snes->ops->convergeddestroy = destroy; 2053 snes->cnvP = cctx; 2054 PetscFunctionReturn(0); 2055 } 2056 2057 #undef __FUNCT__ 2058 #define __FUNCT__ "SNESGetConvergedReason" 2059 /*@ 2060 SNESGetConvergedReason - Gets the reason the SNES iteration was stopped. 2061 2062 Not Collective 2063 2064 Input Parameter: 2065 . snes - the SNES context 2066 2067 Output Parameter: 2068 . reason - negative value indicates diverged, positive value converged, see SNESConvergedReason or the 2069 manual pages for the individual convergence tests for complete lists 2070 2071 Level: intermediate 2072 2073 Notes: Can only be called after the call the SNESSolve() is complete. 2074 2075 .keywords: SNES, nonlinear, set, convergence, test 2076 2077 .seealso: SNESSetConvergenceTest(), SNESConvergedReason 2078 @*/ 2079 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason) 2080 { 2081 PetscFunctionBegin; 2082 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2083 PetscValidPointer(reason,2); 2084 *reason = snes->reason; 2085 PetscFunctionReturn(0); 2086 } 2087 2088 #undef __FUNCT__ 2089 #define __FUNCT__ "SNESSetConvergenceHistory" 2090 /*@ 2091 SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 2092 2093 Logically Collective on SNES 2094 2095 Input Parameters: 2096 + snes - iterative context obtained from SNESCreate() 2097 . a - array to hold history, this array will contain the function norms computed at each step 2098 . its - integer array holds the number of linear iterations for each solve. 2099 . na - size of a and its 2100 - reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero, 2101 else it continues storing new values for new nonlinear solves after the old ones 2102 2103 This routine is useful, e.g., when running a code for purposes 2104 of accurate performance monitoring, when no I/O should be done 2105 during the section of code that is being timed. 2106 2107 Level: intermediate 2108 2109 .keywords: SNES, set, convergence, history 2110 2111 .seealso: SNESGetConvergenceHistory() 2112 2113 @*/ 2114 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscTruth reset) 2115 { 2116 PetscFunctionBegin; 2117 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2118 if (na) PetscValidScalarPointer(a,2); 2119 if (its) PetscValidIntPointer(its,3); 2120 snes->conv_hist = a; 2121 snes->conv_hist_its = its; 2122 snes->conv_hist_max = na; 2123 snes->conv_hist_len = 0; 2124 snes->conv_hist_reset = reset; 2125 PetscFunctionReturn(0); 2126 } 2127 2128 #undef __FUNCT__ 2129 #define __FUNCT__ "SNESGetConvergenceHistory" 2130 /*@C 2131 SNESGetConvergenceHistory - Gets the array used to hold the convergence history. 2132 2133 Not Collective 2134 2135 Input Parameter: 2136 . snes - iterative context obtained from SNESCreate() 2137 2138 Output Parameters: 2139 . a - array to hold history 2140 . its - integer array holds the number of linear iterations (or 2141 negative if not converged) for each solve. 2142 - na - size of a and its 2143 2144 Notes: 2145 The calling sequence for this routine in Fortran is 2146 $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr) 2147 2148 This routine is useful, e.g., when running a code for purposes 2149 of accurate performance monitoring, when no I/O should be done 2150 during the section of code that is being timed. 2151 2152 Level: intermediate 2153 2154 .keywords: SNES, get, convergence, history 2155 2156 .seealso: SNESSetConvergencHistory() 2157 2158 @*/ 2159 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na) 2160 { 2161 PetscFunctionBegin; 2162 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2163 if (a) *a = snes->conv_hist; 2164 if (its) *its = snes->conv_hist_its; 2165 if (na) *na = snes->conv_hist_len; 2166 PetscFunctionReturn(0); 2167 } 2168 2169 #undef __FUNCT__ 2170 #define __FUNCT__ "SNESSetUpdate" 2171 /*@C 2172 SNESSetUpdate - Sets the general-purpose update function called 2173 at the beginning o every iteration of the nonlinear solve. Specifically 2174 it is called just before the Jacobian is "evaluated". 2175 2176 Logically Collective on SNES 2177 2178 Input Parameters: 2179 . snes - The nonlinear solver context 2180 . func - The function 2181 2182 Calling sequence of func: 2183 . func (SNES snes, PetscInt step); 2184 2185 . step - The current step of the iteration 2186 2187 Level: advanced 2188 2189 Note: This is NOT what one uses to update the ghost points before a function evaluation, that should be done at the beginning of your FormFunction() 2190 This is not used by most users. 2191 2192 .keywords: SNES, update 2193 2194 .seealso SNESDefaultUpdate(), SNESSetJacobian(), SNESSolve() 2195 @*/ 2196 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt)) 2197 { 2198 PetscFunctionBegin; 2199 PetscValidHeaderSpecific(snes, SNES_CLASSID,1); 2200 snes->ops->update = func; 2201 PetscFunctionReturn(0); 2202 } 2203 2204 #undef __FUNCT__ 2205 #define __FUNCT__ "SNESDefaultUpdate" 2206 /*@ 2207 SNESDefaultUpdate - The default update function which does nothing. 2208 2209 Not collective 2210 2211 Input Parameters: 2212 . snes - The nonlinear solver context 2213 . step - The current step of the iteration 2214 2215 Level: intermediate 2216 2217 .keywords: SNES, update 2218 .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultShortolutionBC() 2219 @*/ 2220 PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultUpdate(SNES snes, PetscInt step) 2221 { 2222 PetscFunctionBegin; 2223 PetscFunctionReturn(0); 2224 } 2225 2226 #undef __FUNCT__ 2227 #define __FUNCT__ "SNESScaleStep_Private" 2228 /* 2229 SNESScaleStep_Private - Scales a step so that its length is less than the 2230 positive parameter delta. 2231 2232 Input Parameters: 2233 + snes - the SNES context 2234 . y - approximate solution of linear system 2235 . fnorm - 2-norm of current function 2236 - delta - trust region size 2237 2238 Output Parameters: 2239 + gpnorm - predicted function norm at the new point, assuming local 2240 linearization. The value is zero if the step lies within the trust 2241 region, and exceeds zero otherwise. 2242 - ynorm - 2-norm of the step 2243 2244 Note: 2245 For non-trust region methods such as SNESLS, the parameter delta 2246 is set to be the maximum allowable step size. 2247 2248 .keywords: SNES, nonlinear, scale, step 2249 */ 2250 PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm) 2251 { 2252 PetscReal nrm; 2253 PetscScalar cnorm; 2254 PetscErrorCode ierr; 2255 2256 PetscFunctionBegin; 2257 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2258 PetscValidHeaderSpecific(y,VEC_CLASSID,2); 2259 PetscCheckSameComm(snes,1,y,2); 2260 2261 ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr); 2262 if (nrm > *delta) { 2263 nrm = *delta/nrm; 2264 *gpnorm = (1.0 - nrm)*(*fnorm); 2265 cnorm = nrm; 2266 ierr = VecScale(y,cnorm);CHKERRQ(ierr); 2267 *ynorm = *delta; 2268 } else { 2269 *gpnorm = 0.0; 2270 *ynorm = nrm; 2271 } 2272 PetscFunctionReturn(0); 2273 } 2274 2275 #undef __FUNCT__ 2276 #define __FUNCT__ "SNESSolve" 2277 /*@C 2278 SNESSolve - Solves a nonlinear system F(x) = b. 2279 Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX(). 2280 2281 Collective on SNES 2282 2283 Input Parameters: 2284 + snes - the SNES context 2285 . b - the constant part of the equation, or PETSC_NULL to use zero. 2286 - x - the solution vector. 2287 2288 Notes: 2289 The user should initialize the vector,x, with the initial guess 2290 for the nonlinear solve prior to calling SNESSolve. In particular, 2291 to employ an initial guess of zero, the user should explicitly set 2292 this vector to zero by calling VecSet(). 2293 2294 Level: beginner 2295 2296 .keywords: SNES, nonlinear, solve 2297 2298 .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian() 2299 @*/ 2300 PetscErrorCode PETSCSNES_DLLEXPORT SNESSolve(SNES snes,Vec b,Vec x) 2301 { 2302 PetscErrorCode ierr; 2303 PetscTruth flg; 2304 char filename[PETSC_MAX_PATH_LEN]; 2305 PetscViewer viewer; 2306 2307 PetscFunctionBegin; 2308 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2309 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 2310 PetscCheckSameComm(snes,1,x,3); 2311 if (b) PetscValidHeaderSpecific(b,VEC_CLASSID,2); 2312 if (b) PetscCheckSameComm(snes,1,b,2); 2313 2314 /* set solution vector */ 2315 ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 2316 if (snes->vec_sol) { ierr = VecDestroy(snes->vec_sol);CHKERRQ(ierr); } 2317 snes->vec_sol = x; 2318 /* set afine vector if provided */ 2319 if (b) { ierr = PetscObjectReference((PetscObject)b);CHKERRQ(ierr); } 2320 if (snes->vec_rhs) { ierr = VecDestroy(snes->vec_rhs);CHKERRQ(ierr); } 2321 snes->vec_rhs = b; 2322 2323 if (!snes->vec_func && snes->vec_rhs) { 2324 ierr = VecDuplicate(b, &snes->vec_func);CHKERRQ(ierr); 2325 } 2326 2327 ierr = SNESSetUp(snes);CHKERRQ(ierr); 2328 2329 if (snes->conv_hist_reset) snes->conv_hist_len = 0; 2330 snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0; 2331 2332 ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 2333 ierr = (*snes->ops->solve)(snes);CHKERRQ(ierr); 2334 ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 2335 if (snes->domainerror){ 2336 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 2337 snes->domainerror = PETSC_FALSE; 2338 } 2339 if (!snes->reason) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason"); 2340 2341 ierr = PetscOptionsGetString(((PetscObject)snes)->prefix,"-snes_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 2342 if (flg && !PetscPreLoadingOn) { 2343 ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,filename,&viewer);CHKERRQ(ierr); 2344 ierr = SNESView(snes,viewer);CHKERRQ(ierr); 2345 ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr); 2346 } 2347 2348 flg = PETSC_FALSE; 2349 ierr = PetscOptionsGetTruth(((PetscObject)snes)->prefix,"-snes_test_local_min",&flg,PETSC_NULL);CHKERRQ(ierr); 2350 if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); } 2351 if (snes->printreason) { 2352 if (snes->reason > 0) { 2353 ierr = PetscPrintf(((PetscObject)snes)->comm,"Nonlinear solve converged due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr); 2354 } else { 2355 ierr = PetscPrintf(((PetscObject)snes)->comm,"Nonlinear solve did not converge due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr); 2356 } 2357 } 2358 2359 if (snes->errorifnotconverged && snes->reason < 0) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_NOT_CONVERGED,"SNESSolve has not converged"); 2360 PetscFunctionReturn(0); 2361 } 2362 2363 /* --------- Internal routines for SNES Package --------- */ 2364 2365 #undef __FUNCT__ 2366 #define __FUNCT__ "SNESSetType" 2367 /*@C 2368 SNESSetType - Sets the method for the nonlinear solver. 2369 2370 Collective on SNES 2371 2372 Input Parameters: 2373 + snes - the SNES context 2374 - type - a known method 2375 2376 Options Database Key: 2377 . -snes_type <type> - Sets the method; use -help for a list 2378 of available methods (for instance, ls or tr) 2379 2380 Notes: 2381 See "petsc/include/petscsnes.h" for available methods (for instance) 2382 + SNESLS - Newton's method with line search 2383 (systems of nonlinear equations) 2384 . SNESTR - Newton's method with trust region 2385 (systems of nonlinear equations) 2386 2387 Normally, it is best to use the SNESSetFromOptions() command and then 2388 set the SNES solver type from the options database rather than by using 2389 this routine. Using the options database provides the user with 2390 maximum flexibility in evaluating the many nonlinear solvers. 2391 The SNESSetType() routine is provided for those situations where it 2392 is necessary to set the nonlinear solver independently of the command 2393 line or options database. This might be the case, for example, when 2394 the choice of solver changes during the execution of the program, 2395 and the user's application is taking responsibility for choosing the 2396 appropriate method. 2397 2398 Level: intermediate 2399 2400 .keywords: SNES, set, type 2401 2402 .seealso: SNESType, SNESCreate() 2403 2404 @*/ 2405 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetType(SNES snes,const SNESType type) 2406 { 2407 PetscErrorCode ierr,(*r)(SNES); 2408 PetscTruth match; 2409 2410 PetscFunctionBegin; 2411 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2412 PetscValidCharPointer(type,2); 2413 2414 ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr); 2415 if (match) PetscFunctionReturn(0); 2416 2417 ierr = PetscFListFind(SNESList,((PetscObject)snes)->comm,type,(void (**)(void)) &r);CHKERRQ(ierr); 2418 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type); 2419 /* Destroy the previous private SNES context */ 2420 if (snes->ops->destroy) { ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr); } 2421 /* Reinitialize function pointers in SNESOps structure */ 2422 snes->ops->setup = 0; 2423 snes->ops->solve = 0; 2424 snes->ops->view = 0; 2425 snes->ops->setfromoptions = 0; 2426 snes->ops->destroy = 0; 2427 /* Call the SNESCreate_XXX routine for this particular Nonlinear solver */ 2428 snes->setupcalled = PETSC_FALSE; 2429 ierr = (*r)(snes);CHKERRQ(ierr); 2430 ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr); 2431 PetscFunctionReturn(0); 2432 } 2433 2434 2435 /* --------------------------------------------------------------------- */ 2436 #undef __FUNCT__ 2437 #define __FUNCT__ "SNESRegisterDestroy" 2438 /*@ 2439 SNESRegisterDestroy - Frees the list of nonlinear solvers that were 2440 registered by SNESRegisterDynamic(). 2441 2442 Not Collective 2443 2444 Level: advanced 2445 2446 .keywords: SNES, nonlinear, register, destroy 2447 2448 .seealso: SNESRegisterAll(), SNESRegisterAll() 2449 @*/ 2450 PetscErrorCode PETSCSNES_DLLEXPORT SNESRegisterDestroy(void) 2451 { 2452 PetscErrorCode ierr; 2453 2454 PetscFunctionBegin; 2455 ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr); 2456 SNESRegisterAllCalled = PETSC_FALSE; 2457 PetscFunctionReturn(0); 2458 } 2459 2460 #undef __FUNCT__ 2461 #define __FUNCT__ "SNESGetType" 2462 /*@C 2463 SNESGetType - Gets the SNES method type and name (as a string). 2464 2465 Not Collective 2466 2467 Input Parameter: 2468 . snes - nonlinear solver context 2469 2470 Output Parameter: 2471 . type - SNES method (a character string) 2472 2473 Level: intermediate 2474 2475 .keywords: SNES, nonlinear, get, type, name 2476 @*/ 2477 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetType(SNES snes,const SNESType *type) 2478 { 2479 PetscFunctionBegin; 2480 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2481 PetscValidPointer(type,2); 2482 *type = ((PetscObject)snes)->type_name; 2483 PetscFunctionReturn(0); 2484 } 2485 2486 #undef __FUNCT__ 2487 #define __FUNCT__ "SNESGetSolution" 2488 /*@ 2489 SNESGetSolution - Returns the vector where the approximate solution is 2490 stored. 2491 2492 Not Collective, but Vec is parallel if SNES is parallel 2493 2494 Input Parameter: 2495 . snes - the SNES context 2496 2497 Output Parameter: 2498 . x - the solution 2499 2500 Level: intermediate 2501 2502 .keywords: SNES, nonlinear, get, solution 2503 2504 .seealso: SNESGetSolutionUpdate(), SNESGetFunction() 2505 @*/ 2506 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolution(SNES snes,Vec *x) 2507 { 2508 PetscFunctionBegin; 2509 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2510 PetscValidPointer(x,2); 2511 *x = snes->vec_sol; 2512 PetscFunctionReturn(0); 2513 } 2514 2515 #undef __FUNCT__ 2516 #define __FUNCT__ "SNESGetSolutionUpdate" 2517 /*@ 2518 SNESGetSolutionUpdate - Returns the vector where the solution update is 2519 stored. 2520 2521 Not Collective, but Vec is parallel if SNES is parallel 2522 2523 Input Parameter: 2524 . snes - the SNES context 2525 2526 Output Parameter: 2527 . x - the solution update 2528 2529 Level: advanced 2530 2531 .keywords: SNES, nonlinear, get, solution, update 2532 2533 .seealso: SNESGetSolution(), SNESGetFunction() 2534 @*/ 2535 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolutionUpdate(SNES snes,Vec *x) 2536 { 2537 PetscFunctionBegin; 2538 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2539 PetscValidPointer(x,2); 2540 *x = snes->vec_sol_update; 2541 PetscFunctionReturn(0); 2542 } 2543 2544 #undef __FUNCT__ 2545 #define __FUNCT__ "SNESGetFunction" 2546 /*@C 2547 SNESGetFunction - Returns the vector where the function is stored. 2548 2549 Not Collective, but Vec is parallel if SNES is parallel 2550 2551 Input Parameter: 2552 . snes - the SNES context 2553 2554 Output Parameter: 2555 + r - the function (or PETSC_NULL) 2556 . func - the function (or PETSC_NULL) 2557 - ctx - the function context (or PETSC_NULL) 2558 2559 Level: advanced 2560 2561 .keywords: SNES, nonlinear, get, function 2562 2563 .seealso: SNESSetFunction(), SNESGetSolution() 2564 @*/ 2565 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx) 2566 { 2567 PetscFunctionBegin; 2568 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2569 if (r) *r = snes->vec_func; 2570 if (func) *func = snes->ops->computefunction; 2571 if (ctx) *ctx = snes->funP; 2572 PetscFunctionReturn(0); 2573 } 2574 2575 #undef __FUNCT__ 2576 #define __FUNCT__ "SNESSetOptionsPrefix" 2577 /*@C 2578 SNESSetOptionsPrefix - Sets the prefix used for searching for all 2579 SNES options in the database. 2580 2581 Logically Collective on SNES 2582 2583 Input Parameter: 2584 + snes - the SNES context 2585 - prefix - the prefix to prepend to all option names 2586 2587 Notes: 2588 A hyphen (-) must NOT be given at the beginning of the prefix name. 2589 The first character of all runtime options is AUTOMATICALLY the hyphen. 2590 2591 Level: advanced 2592 2593 .keywords: SNES, set, options, prefix, database 2594 2595 .seealso: SNESSetFromOptions() 2596 @*/ 2597 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetOptionsPrefix(SNES snes,const char prefix[]) 2598 { 2599 PetscErrorCode ierr; 2600 2601 PetscFunctionBegin; 2602 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2603 ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2604 if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);} 2605 ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 2606 PetscFunctionReturn(0); 2607 } 2608 2609 #undef __FUNCT__ 2610 #define __FUNCT__ "SNESAppendOptionsPrefix" 2611 /*@C 2612 SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 2613 SNES options in the database. 2614 2615 Logically Collective on SNES 2616 2617 Input Parameters: 2618 + snes - the SNES context 2619 - prefix - the prefix to prepend to all option names 2620 2621 Notes: 2622 A hyphen (-) must NOT be given at the beginning of the prefix name. 2623 The first character of all runtime options is AUTOMATICALLY the hyphen. 2624 2625 Level: advanced 2626 2627 .keywords: SNES, append, options, prefix, database 2628 2629 .seealso: SNESGetOptionsPrefix() 2630 @*/ 2631 PetscErrorCode PETSCSNES_DLLEXPORT SNESAppendOptionsPrefix(SNES snes,const char prefix[]) 2632 { 2633 PetscErrorCode ierr; 2634 2635 PetscFunctionBegin; 2636 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2637 ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2638 if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);} 2639 ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 2640 PetscFunctionReturn(0); 2641 } 2642 2643 #undef __FUNCT__ 2644 #define __FUNCT__ "SNESGetOptionsPrefix" 2645 /*@C 2646 SNESGetOptionsPrefix - Sets the prefix used for searching for all 2647 SNES options in the database. 2648 2649 Not Collective 2650 2651 Input Parameter: 2652 . snes - the SNES context 2653 2654 Output Parameter: 2655 . prefix - pointer to the prefix string used 2656 2657 Notes: On the fortran side, the user should pass in a string 'prefix' of 2658 sufficient length to hold the prefix. 2659 2660 Level: advanced 2661 2662 .keywords: SNES, get, options, prefix, database 2663 2664 .seealso: SNESAppendOptionsPrefix() 2665 @*/ 2666 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetOptionsPrefix(SNES snes,const char *prefix[]) 2667 { 2668 PetscErrorCode ierr; 2669 2670 PetscFunctionBegin; 2671 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2672 ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2673 PetscFunctionReturn(0); 2674 } 2675 2676 2677 #undef __FUNCT__ 2678 #define __FUNCT__ "SNESRegister" 2679 /*@C 2680 SNESRegister - See SNESRegisterDynamic() 2681 2682 Level: advanced 2683 @*/ 2684 PetscErrorCode PETSCSNES_DLLEXPORT SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES)) 2685 { 2686 char fullname[PETSC_MAX_PATH_LEN]; 2687 PetscErrorCode ierr; 2688 2689 PetscFunctionBegin; 2690 ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 2691 ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 2692 PetscFunctionReturn(0); 2693 } 2694 2695 #undef __FUNCT__ 2696 #define __FUNCT__ "SNESTestLocalMin" 2697 PetscErrorCode PETSCSNES_DLLEXPORT SNESTestLocalMin(SNES snes) 2698 { 2699 PetscErrorCode ierr; 2700 PetscInt N,i,j; 2701 Vec u,uh,fh; 2702 PetscScalar value; 2703 PetscReal norm; 2704 2705 PetscFunctionBegin; 2706 ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr); 2707 ierr = VecDuplicate(u,&uh);CHKERRQ(ierr); 2708 ierr = VecDuplicate(u,&fh);CHKERRQ(ierr); 2709 2710 /* currently only works for sequential */ 2711 ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n"); 2712 ierr = VecGetSize(u,&N);CHKERRQ(ierr); 2713 for (i=0; i<N; i++) { 2714 ierr = VecCopy(u,uh);CHKERRQ(ierr); 2715 ierr = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr); 2716 for (j=-10; j<11; j++) { 2717 value = PetscSign(j)*exp(PetscAbs(j)-10.0); 2718 ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2719 ierr = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr); 2720 ierr = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr); 2721 ierr = PetscPrintf(PETSC_COMM_WORLD," j norm %D %18.16e\n",j,norm);CHKERRQ(ierr); 2722 value = -value; 2723 ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2724 } 2725 } 2726 ierr = VecDestroy(uh);CHKERRQ(ierr); 2727 ierr = VecDestroy(fh);CHKERRQ(ierr); 2728 PetscFunctionReturn(0); 2729 } 2730 2731 #undef __FUNCT__ 2732 #define __FUNCT__ "SNESKSPSetUseEW" 2733 /*@ 2734 SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for 2735 computing relative tolerance for linear solvers within an inexact 2736 Newton method. 2737 2738 Logically Collective on SNES 2739 2740 Input Parameters: 2741 + snes - SNES context 2742 - flag - PETSC_TRUE or PETSC_FALSE 2743 2744 Options Database: 2745 + -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence 2746 . -snes_ksp_ew_version ver - version of Eisenstat-Walker method 2747 . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0 2748 . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax 2749 . -snes_ksp_ew_gamma <gamma> - Sets gamma 2750 . -snes_ksp_ew_alpha <alpha> - Sets alpha 2751 . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2 2752 - -snes_ksp_ew_threshold <threshold> - Sets threshold 2753 2754 Notes: 2755 Currently, the default is to use a constant relative tolerance for 2756 the inner linear solvers. Alternatively, one can use the 2757 Eisenstat-Walker method, where the relative convergence tolerance 2758 is reset at each Newton iteration according progress of the nonlinear 2759 solver. 2760 2761 Level: advanced 2762 2763 Reference: 2764 S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an 2765 inexact Newton method", SISC 17 (1), pp.16-32, 1996. 2766 2767 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton 2768 2769 .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW() 2770 @*/ 2771 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPSetUseEW(SNES snes,PetscTruth flag) 2772 { 2773 PetscFunctionBegin; 2774 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2775 PetscValidLogicalCollectiveTruth(snes,flag,2); 2776 snes->ksp_ewconv = flag; 2777 PetscFunctionReturn(0); 2778 } 2779 2780 #undef __FUNCT__ 2781 #define __FUNCT__ "SNESKSPGetUseEW" 2782 /*@ 2783 SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method 2784 for computing relative tolerance for linear solvers within an 2785 inexact Newton method. 2786 2787 Not Collective 2788 2789 Input Parameter: 2790 . snes - SNES context 2791 2792 Output Parameter: 2793 . flag - PETSC_TRUE or PETSC_FALSE 2794 2795 Notes: 2796 Currently, the default is to use a constant relative tolerance for 2797 the inner linear solvers. Alternatively, one can use the 2798 Eisenstat-Walker method, where the relative convergence tolerance 2799 is reset at each Newton iteration according progress of the nonlinear 2800 solver. 2801 2802 Level: advanced 2803 2804 Reference: 2805 S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an 2806 inexact Newton method", SISC 17 (1), pp.16-32, 1996. 2807 2808 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton 2809 2810 .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW() 2811 @*/ 2812 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPGetUseEW(SNES snes, PetscTruth *flag) 2813 { 2814 PetscFunctionBegin; 2815 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2816 PetscValidPointer(flag,2); 2817 *flag = snes->ksp_ewconv; 2818 PetscFunctionReturn(0); 2819 } 2820 2821 #undef __FUNCT__ 2822 #define __FUNCT__ "SNESKSPSetParametersEW" 2823 /*@ 2824 SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker 2825 convergence criteria for the linear solvers within an inexact 2826 Newton method. 2827 2828 Logically Collective on SNES 2829 2830 Input Parameters: 2831 + snes - SNES context 2832 . version - version 1, 2 (default is 2) or 3 2833 . rtol_0 - initial relative tolerance (0 <= rtol_0 < 1) 2834 . rtol_max - maximum relative tolerance (0 <= rtol_max < 1) 2835 . gamma - multiplicative factor for version 2 rtol computation 2836 (0 <= gamma2 <= 1) 2837 . alpha - power for version 2 rtol computation (1 < alpha <= 2) 2838 . alpha2 - power for safeguard 2839 - threshold - threshold for imposing safeguard (0 < threshold < 1) 2840 2841 Note: 2842 Version 3 was contributed by Luis Chacon, June 2006. 2843 2844 Use PETSC_DEFAULT to retain the default for any of the parameters. 2845 2846 Level: advanced 2847 2848 Reference: 2849 S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an 2850 inexact Newton method", Utah State University Math. Stat. Dept. Res. 2851 Report 6/94/75, June, 1994, to appear in SIAM J. Sci. Comput. 2852 2853 .keywords: SNES, KSP, Eisenstat, Walker, set, parameters 2854 2855 .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW() 2856 @*/ 2857 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max, 2858 PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold) 2859 { 2860 SNESKSPEW *kctx; 2861 PetscFunctionBegin; 2862 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2863 kctx = (SNESKSPEW*)snes->kspconvctx; 2864 if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing"); 2865 PetscValidLogicalCollectiveInt(snes,version,2); 2866 PetscValidLogicalCollectiveReal(snes,rtol_0,3); 2867 PetscValidLogicalCollectiveReal(snes,rtol_max,4); 2868 PetscValidLogicalCollectiveReal(snes,gamma,5); 2869 PetscValidLogicalCollectiveReal(snes,alpha,6); 2870 PetscValidLogicalCollectiveReal(snes,alpha2,7); 2871 PetscValidLogicalCollectiveReal(snes,threshold,8); 2872 2873 if (version != PETSC_DEFAULT) kctx->version = version; 2874 if (rtol_0 != PETSC_DEFAULT) kctx->rtol_0 = rtol_0; 2875 if (rtol_max != PETSC_DEFAULT) kctx->rtol_max = rtol_max; 2876 if (gamma != PETSC_DEFAULT) kctx->gamma = gamma; 2877 if (alpha != PETSC_DEFAULT) kctx->alpha = alpha; 2878 if (alpha2 != PETSC_DEFAULT) kctx->alpha2 = alpha2; 2879 if (threshold != PETSC_DEFAULT) kctx->threshold = threshold; 2880 2881 if (kctx->version < 1 || kctx->version > 3) { 2882 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 and 3 are supported: %D",kctx->version); 2883 } 2884 if (kctx->rtol_0 < 0.0 || kctx->rtol_0 >= 1.0) { 2885 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_0 < 1.0: %G",kctx->rtol_0); 2886 } 2887 if (kctx->rtol_max < 0.0 || kctx->rtol_max >= 1.0) { 2888 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_max (%G) < 1.0\n",kctx->rtol_max); 2889 } 2890 if (kctx->gamma < 0.0 || kctx->gamma > 1.0) { 2891 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= gamma (%G) <= 1.0\n",kctx->gamma); 2892 } 2893 if (kctx->alpha <= 1.0 || kctx->alpha > 2.0) { 2894 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"1.0 < alpha (%G) <= 2.0\n",kctx->alpha); 2895 } 2896 if (kctx->threshold <= 0.0 || kctx->threshold >= 1.0) { 2897 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 < threshold (%G) < 1.0\n",kctx->threshold); 2898 } 2899 PetscFunctionReturn(0); 2900 } 2901 2902 #undef __FUNCT__ 2903 #define __FUNCT__ "SNESKSPGetParametersEW" 2904 /*@ 2905 SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker 2906 convergence criteria for the linear solvers within an inexact 2907 Newton method. 2908 2909 Not Collective 2910 2911 Input Parameters: 2912 snes - SNES context 2913 2914 Output Parameters: 2915 + version - version 1, 2 (default is 2) or 3 2916 . rtol_0 - initial relative tolerance (0 <= rtol_0 < 1) 2917 . rtol_max - maximum relative tolerance (0 <= rtol_max < 1) 2918 . gamma - multiplicative factor for version 2 rtol computation 2919 (0 <= gamma2 <= 1) 2920 . alpha - power for version 2 rtol computation (1 < alpha <= 2) 2921 . alpha2 - power for safeguard 2922 - threshold - threshold for imposing safeguard (0 < threshold < 1) 2923 2924 Level: advanced 2925 2926 .keywords: SNES, KSP, Eisenstat, Walker, get, parameters 2927 2928 .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW() 2929 @*/ 2930 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max, 2931 PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold) 2932 { 2933 SNESKSPEW *kctx; 2934 PetscFunctionBegin; 2935 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2936 kctx = (SNESKSPEW*)snes->kspconvctx; 2937 if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing"); 2938 if(version) *version = kctx->version; 2939 if(rtol_0) *rtol_0 = kctx->rtol_0; 2940 if(rtol_max) *rtol_max = kctx->rtol_max; 2941 if(gamma) *gamma = kctx->gamma; 2942 if(alpha) *alpha = kctx->alpha; 2943 if(alpha2) *alpha2 = kctx->alpha2; 2944 if(threshold) *threshold = kctx->threshold; 2945 PetscFunctionReturn(0); 2946 } 2947 2948 #undef __FUNCT__ 2949 #define __FUNCT__ "SNESKSPEW_PreSolve" 2950 static PetscErrorCode SNESKSPEW_PreSolve(SNES snes, KSP ksp, Vec b, Vec x) 2951 { 2952 PetscErrorCode ierr; 2953 SNESKSPEW *kctx = (SNESKSPEW*)snes->kspconvctx; 2954 PetscReal rtol=PETSC_DEFAULT,stol; 2955 2956 PetscFunctionBegin; 2957 if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists"); 2958 if (!snes->iter) { /* first time in, so use the original user rtol */ 2959 rtol = kctx->rtol_0; 2960 } else { 2961 if (kctx->version == 1) { 2962 rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last; 2963 if (rtol < 0.0) rtol = -rtol; 2964 stol = pow(kctx->rtol_last,kctx->alpha2); 2965 if (stol > kctx->threshold) rtol = PetscMax(rtol,stol); 2966 } else if (kctx->version == 2) { 2967 rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha); 2968 stol = kctx->gamma * pow(kctx->rtol_last,kctx->alpha); 2969 if (stol > kctx->threshold) rtol = PetscMax(rtol,stol); 2970 } else if (kctx->version == 3) {/* contributed by Luis Chacon, June 2006. */ 2971 rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha); 2972 /* safeguard: avoid sharp decrease of rtol */ 2973 stol = kctx->gamma*pow(kctx->rtol_last,kctx->alpha); 2974 stol = PetscMax(rtol,stol); 2975 rtol = PetscMin(kctx->rtol_0,stol); 2976 /* safeguard: avoid oversolving */ 2977 stol = kctx->gamma*(snes->ttol)/snes->norm; 2978 stol = PetscMax(rtol,stol); 2979 rtol = PetscMin(kctx->rtol_0,stol); 2980 } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 or 3 are supported: %D",kctx->version); 2981 } 2982 /* safeguard: avoid rtol greater than one */ 2983 rtol = PetscMin(rtol,kctx->rtol_max); 2984 ierr = KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 2985 ierr = PetscInfo3(snes,"iter %D, Eisenstat-Walker (version %D) KSP rtol=%G\n",snes->iter,kctx->version,rtol);CHKERRQ(ierr); 2986 PetscFunctionReturn(0); 2987 } 2988 2989 #undef __FUNCT__ 2990 #define __FUNCT__ "SNESKSPEW_PostSolve" 2991 static PetscErrorCode SNESKSPEW_PostSolve(SNES snes, KSP ksp, Vec b, Vec x) 2992 { 2993 PetscErrorCode ierr; 2994 SNESKSPEW *kctx = (SNESKSPEW*)snes->kspconvctx; 2995 PCSide pcside; 2996 Vec lres; 2997 2998 PetscFunctionBegin; 2999 if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists"); 3000 ierr = KSPGetTolerances(ksp,&kctx->rtol_last,0,0,0);CHKERRQ(ierr); 3001 ierr = SNESGetFunctionNorm(snes,&kctx->norm_last);CHKERRQ(ierr); 3002 if (kctx->version == 1) { 3003 ierr = KSPGetPCSide(ksp,&pcside);CHKERRQ(ierr); 3004 if (pcside == PC_RIGHT) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */ 3005 /* KSP residual is true linear residual */ 3006 ierr = KSPGetResidualNorm(ksp,&kctx->lresid_last);CHKERRQ(ierr); 3007 } else { 3008 /* KSP residual is preconditioned residual */ 3009 /* compute true linear residual norm */ 3010 ierr = VecDuplicate(b,&lres);CHKERRQ(ierr); 3011 ierr = MatMult(snes->jacobian,x,lres);CHKERRQ(ierr); 3012 ierr = VecAYPX(lres,-1.0,b);CHKERRQ(ierr); 3013 ierr = VecNorm(lres,NORM_2,&kctx->lresid_last);CHKERRQ(ierr); 3014 ierr = VecDestroy(lres);CHKERRQ(ierr); 3015 } 3016 } 3017 PetscFunctionReturn(0); 3018 } 3019 3020 #undef __FUNCT__ 3021 #define __FUNCT__ "SNES_KSPSolve" 3022 PetscErrorCode SNES_KSPSolve(SNES snes, KSP ksp, Vec b, Vec x) 3023 { 3024 PetscErrorCode ierr; 3025 3026 PetscFunctionBegin; 3027 if (snes->ksp_ewconv) { ierr = SNESKSPEW_PreSolve(snes,ksp,b,x);CHKERRQ(ierr); } 3028 ierr = KSPSolve(ksp,b,x);CHKERRQ(ierr); 3029 if (snes->ksp_ewconv) { ierr = SNESKSPEW_PostSolve(snes,ksp,b,x);CHKERRQ(ierr); } 3030 PetscFunctionReturn(0); 3031 } 3032 3033 #undef __FUNCT__ 3034 #define __FUNCT__ "SNESSetDM" 3035 /*@ 3036 SNESSetDM - Sets the DM that may be used by some preconditioners 3037 3038 Logically Collective on SNES 3039 3040 Input Parameters: 3041 + snes - the preconditioner context 3042 - dm - the dm 3043 3044 Level: intermediate 3045 3046 3047 .seealso: SNESGetDM(), KSPSetDM(), KSPGetDM() 3048 @*/ 3049 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetDM(SNES snes,DM dm) 3050 { 3051 PetscErrorCode ierr; 3052 KSP ksp; 3053 3054 PetscFunctionBegin; 3055 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3056 ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr); 3057 if (snes->dm) {ierr = DMDestroy(snes->dm);CHKERRQ(ierr);} 3058 snes->dm = dm; 3059 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 3060 ierr = KSPSetDM(ksp,dm);CHKERRQ(ierr); 3061 ierr = KSPSetDMActive(ksp,PETSC_FALSE);CHKERRQ(ierr); 3062 PetscFunctionReturn(0); 3063 } 3064 3065 #undef __FUNCT__ 3066 #define __FUNCT__ "SNESGetDM" 3067 /*@ 3068 SNESGetDM - Gets the DM that may be used by some preconditioners 3069 3070 Not Collective but DM obtained is parallel on SNES 3071 3072 Input Parameter: 3073 . snes - the preconditioner context 3074 3075 Output Parameter: 3076 . dm - the dm 3077 3078 Level: intermediate 3079 3080 3081 .seealso: SNESSetDM(), KSPSetDM(), KSPGetDM() 3082 @*/ 3083 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetDM(SNES snes,DM *dm) 3084 { 3085 PetscFunctionBegin; 3086 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3087 *dm = snes->dm; 3088 PetscFunctionReturn(0); 3089 } 3090