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 PetscFunctionReturn(0); 1047 } 1048 1049 #undef __FUNCT__ 1050 #define __FUNCT__ "SNESSetFunction" 1051 /*@C 1052 SNESSetFunction - Sets the function evaluation routine and function 1053 vector for use by the SNES routines in solving systems of nonlinear 1054 equations. 1055 1056 Logically Collective on SNES 1057 1058 Input Parameters: 1059 + snes - the SNES context 1060 . r - vector to store function value 1061 . func - function evaluation routine 1062 - ctx - [optional] user-defined context for private data for the 1063 function evaluation routine (may be PETSC_NULL) 1064 1065 Calling sequence of func: 1066 $ func (SNES snes,Vec x,Vec f,void *ctx); 1067 1068 . f - function vector 1069 - ctx - optional user-defined function context 1070 1071 Notes: 1072 The Newton-like methods typically solve linear systems of the form 1073 $ f'(x) x = -f(x), 1074 where f'(x) denotes the Jacobian matrix and f(x) is the function. 1075 1076 Level: beginner 1077 1078 .keywords: SNES, nonlinear, set, function 1079 1080 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 1081 @*/ 1082 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx) 1083 { 1084 PetscErrorCode ierr; 1085 PetscFunctionBegin; 1086 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1087 if (r) PetscValidHeaderSpecific(r,VEC_CLASSID,2); 1088 if (r) PetscCheckSameComm(snes,1,r,2); 1089 if (!r && snes->dm) { 1090 ierr = DMCreateGlobalVector(snes->dm,&r);CHKERRQ(ierr); 1091 } else { 1092 ierr = PetscObjectReference((PetscObject)r);CHKERRQ(ierr); 1093 } 1094 if (snes->vec_func) { ierr = VecDestroy(snes->vec_func);CHKERRQ(ierr); } 1095 snes->ops->computefunction = func; 1096 snes->vec_func = r; 1097 snes->funP = ctx; 1098 PetscFunctionReturn(0); 1099 } 1100 1101 /* --------------------------------------------------------------- */ 1102 #undef __FUNCT__ 1103 #define __FUNCT__ "SNESGetRhs" 1104 /*@C 1105 SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set 1106 it assumes a zero right hand side. 1107 1108 Logically Collective on SNES 1109 1110 Input Parameter: 1111 . snes - the SNES context 1112 1113 Output Parameter: 1114 . rhs - the right hand side vector or PETSC_NULL if the right hand side vector is null 1115 1116 Level: intermediate 1117 1118 .keywords: SNES, nonlinear, get, function, right hand side 1119 1120 .seealso: SNESGetSolution(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction() 1121 @*/ 1122 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetRhs(SNES snes,Vec *rhs) 1123 { 1124 PetscFunctionBegin; 1125 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1126 PetscValidPointer(rhs,2); 1127 *rhs = snes->vec_rhs; 1128 PetscFunctionReturn(0); 1129 } 1130 1131 #undef __FUNCT__ 1132 #define __FUNCT__ "SNESComputeFunction" 1133 /*@ 1134 SNESComputeFunction - Calls the function that has been set with 1135 SNESSetFunction(). 1136 1137 Collective on SNES 1138 1139 Input Parameters: 1140 + snes - the SNES context 1141 - x - input vector 1142 1143 Output Parameter: 1144 . y - function vector, as set by SNESSetFunction() 1145 1146 Notes: 1147 SNESComputeFunction() is typically used within nonlinear solvers 1148 implementations, so most users would not generally call this routine 1149 themselves. 1150 1151 Level: developer 1152 1153 .keywords: SNES, nonlinear, compute, function 1154 1155 .seealso: SNESSetFunction(), SNESGetFunction() 1156 @*/ 1157 PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeFunction(SNES snes,Vec x,Vec y) 1158 { 1159 PetscErrorCode ierr; 1160 1161 PetscFunctionBegin; 1162 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1163 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 1164 PetscValidHeaderSpecific(y,VEC_CLASSID,3); 1165 PetscCheckSameComm(snes,1,x,2); 1166 PetscCheckSameComm(snes,1,y,3); 1167 1168 ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 1169 if (snes->ops->computefunction) { 1170 PetscStackPush("SNES user function"); 1171 ierr = (*snes->ops->computefunction)(snes,x,y,snes->funP); 1172 PetscStackPop; 1173 if (PetscExceptionValue(ierr)) { 1174 PetscErrorCode pierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(pierr); 1175 } 1176 CHKERRQ(ierr); 1177 } else if (snes->vec_rhs) { 1178 ierr = MatMult(snes->jacobian, x, y);CHKERRQ(ierr); 1179 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() before SNESComputeFunction(), likely called from SNESSolve()."); 1180 if (snes->vec_rhs) { 1181 ierr = VecAXPY(y,-1.0,snes->vec_rhs);CHKERRQ(ierr); 1182 } 1183 snes->nfuncs++; 1184 ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 1185 PetscFunctionReturn(0); 1186 } 1187 1188 #undef __FUNCT__ 1189 #define __FUNCT__ "SNESComputeJacobian" 1190 /*@ 1191 SNESComputeJacobian - Computes the Jacobian matrix that has been 1192 set with SNESSetJacobian(). 1193 1194 Collective on SNES and Mat 1195 1196 Input Parameters: 1197 + snes - the SNES context 1198 - x - input vector 1199 1200 Output Parameters: 1201 + A - Jacobian matrix 1202 . B - optional preconditioning matrix 1203 - flag - flag indicating matrix structure (one of, SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER) 1204 1205 Options Database Keys: 1206 + -snes_lag_preconditioner <lag> 1207 - -snes_lag_jacobian <lag> 1208 1209 Notes: 1210 Most users should not need to explicitly call this routine, as it 1211 is used internally within the nonlinear solvers. 1212 1213 See KSPSetOperators() for important information about setting the 1214 flag parameter. 1215 1216 Level: developer 1217 1218 .keywords: SNES, compute, Jacobian, matrix 1219 1220 .seealso: SNESSetJacobian(), KSPSetOperators(), MatStructure, SNESSetLagPreconditioner(), SNESSetLagJacobian() 1221 @*/ 1222 PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 1223 { 1224 PetscErrorCode ierr; 1225 PetscTruth flag; 1226 1227 PetscFunctionBegin; 1228 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1229 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 1230 PetscValidPointer(flg,5); 1231 PetscCheckSameComm(snes,1,X,2); 1232 if (!snes->ops->computejacobian) PetscFunctionReturn(0); 1233 1234 /* make sure that MatAssemblyBegin/End() is called on A matrix if it is matrix free */ 1235 1236 if (snes->lagjacobian == -2) { 1237 snes->lagjacobian = -1; 1238 ierr = PetscInfo(snes,"Recomputing Jacobian/preconditioner because lag is -2 (means compute Jacobian, but then never again) \n");CHKERRQ(ierr); 1239 } else if (snes->lagjacobian == -1) { 1240 *flg = SAME_PRECONDITIONER; 1241 ierr = PetscInfo(snes,"Reusing Jacobian/preconditioner because lag is -1\n");CHKERRQ(ierr); 1242 ierr = PetscTypeCompare((PetscObject)*A,MATMFFD,&flag);CHKERRQ(ierr); 1243 if (flag) { 1244 ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1245 ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1246 } 1247 PetscFunctionReturn(0); 1248 } else if (snes->lagjacobian > 1 && snes->iter % snes->lagjacobian) { 1249 *flg = SAME_PRECONDITIONER; 1250 ierr = PetscInfo2(snes,"Reusing Jacobian/preconditioner because lag is %D and SNES iteration is %D\n",snes->lagjacobian,snes->iter);CHKERRQ(ierr); 1251 ierr = PetscTypeCompare((PetscObject)*A,MATMFFD,&flag);CHKERRQ(ierr); 1252 if (flag) { 1253 ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1254 ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1255 } 1256 PetscFunctionReturn(0); 1257 } 1258 1259 *flg = DIFFERENT_NONZERO_PATTERN; 1260 ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 1261 PetscStackPush("SNES user Jacobian function"); 1262 ierr = (*snes->ops->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr); 1263 PetscStackPop; 1264 ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 1265 1266 if (snes->lagpreconditioner == -2) { 1267 ierr = PetscInfo(snes,"Rebuilding preconditioner exactly once since lag is -2\n");CHKERRQ(ierr); 1268 snes->lagpreconditioner = -1; 1269 } else if (snes->lagpreconditioner == -1) { 1270 *flg = SAME_PRECONDITIONER; 1271 ierr = PetscInfo(snes,"Reusing preconditioner because lag is -1\n");CHKERRQ(ierr); 1272 } else if (snes->lagpreconditioner > 1 && snes->iter % snes->lagpreconditioner) { 1273 *flg = SAME_PRECONDITIONER; 1274 ierr = PetscInfo2(snes,"Reusing preconditioner because lag is %D and SNES iteration is %D\n",snes->lagpreconditioner,snes->iter);CHKERRQ(ierr); 1275 } 1276 1277 /* make sure user returned a correct Jacobian and preconditioner */ 1278 /* PetscValidHeaderSpecific(*A,MAT_CLASSID,3); 1279 PetscValidHeaderSpecific(*B,MAT_CLASSID,4); */ 1280 PetscFunctionReturn(0); 1281 } 1282 1283 #undef __FUNCT__ 1284 #define __FUNCT__ "SNESSetJacobian" 1285 /*@C 1286 SNESSetJacobian - Sets the function to compute Jacobian as well as the 1287 location to store the matrix. 1288 1289 Logically Collective on SNES and Mat 1290 1291 Input Parameters: 1292 + snes - the SNES context 1293 . A - Jacobian matrix 1294 . B - preconditioner matrix (usually same as the Jacobian) 1295 . func - Jacobian evaluation routine 1296 - ctx - [optional] user-defined context for private data for the 1297 Jacobian evaluation routine (may be PETSC_NULL) 1298 1299 Calling sequence of func: 1300 $ func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 1301 1302 + x - input vector 1303 . A - Jacobian matrix 1304 . B - preconditioner matrix, usually the same as A 1305 . flag - flag indicating information about the preconditioner matrix 1306 structure (same as flag in KSPSetOperators()), one of SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER 1307 - ctx - [optional] user-defined Jacobian context 1308 1309 Notes: 1310 See KSPSetOperators() for important information about setting the flag 1311 output parameter in the routine func(). Be sure to read this information! 1312 1313 The routine func() takes Mat * as the matrix arguments rather than Mat. 1314 This allows the Jacobian evaluation routine to replace A and/or B with a 1315 completely new new matrix structure (not just different matrix elements) 1316 when appropriate, for instance, if the nonzero structure is changing 1317 throughout the global iterations. 1318 1319 If the A matrix and B matrix are different you must call MatAssemblyBegin/End() on 1320 each matrix. 1321 1322 Level: beginner 1323 1324 .keywords: SNES, nonlinear, set, Jacobian, matrix 1325 1326 .seealso: KSPSetOperators(), SNESSetFunction(), MatMFFDComputeJacobian(), SNESDefaultComputeJacobianColor(), MatStructure 1327 @*/ 1328 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 1329 { 1330 PetscErrorCode ierr; 1331 1332 PetscFunctionBegin; 1333 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1334 if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2); 1335 if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3); 1336 if (A) PetscCheckSameComm(snes,1,A,2); 1337 if (B) PetscCheckSameComm(snes,1,B,3); 1338 if (func) snes->ops->computejacobian = func; 1339 if (ctx) snes->jacP = ctx; 1340 if (A) { 1341 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 1342 if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 1343 snes->jacobian = A; 1344 } 1345 if (B) { 1346 ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 1347 if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 1348 snes->jacobian_pre = B; 1349 } 1350 PetscFunctionReturn(0); 1351 } 1352 1353 #undef __FUNCT__ 1354 #define __FUNCT__ "SNESGetJacobian" 1355 /*@C 1356 SNESGetJacobian - Returns the Jacobian matrix and optionally the user 1357 provided context for evaluating the Jacobian. 1358 1359 Not Collective, but Mat object will be parallel if SNES object is 1360 1361 Input Parameter: 1362 . snes - the nonlinear solver context 1363 1364 Output Parameters: 1365 + A - location to stash Jacobian matrix (or PETSC_NULL) 1366 . B - location to stash preconditioner matrix (or PETSC_NULL) 1367 . func - location to put Jacobian function (or PETSC_NULL) 1368 - ctx - location to stash Jacobian ctx (or PETSC_NULL) 1369 1370 Level: advanced 1371 1372 .seealso: SNESSetJacobian(), SNESComputeJacobian() 1373 @*/ 1374 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetJacobian(SNES snes,Mat *A,Mat *B,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx) 1375 { 1376 PetscFunctionBegin; 1377 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1378 if (A) *A = snes->jacobian; 1379 if (B) *B = snes->jacobian_pre; 1380 if (func) *func = snes->ops->computejacobian; 1381 if (ctx) *ctx = snes->jacP; 1382 PetscFunctionReturn(0); 1383 } 1384 1385 /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 1386 1387 #undef __FUNCT__ 1388 #define __FUNCT__ "SNESSetUp" 1389 /*@ 1390 SNESSetUp - Sets up the internal data structures for the later use 1391 of a nonlinear solver. 1392 1393 Collective on SNES 1394 1395 Input Parameters: 1396 . snes - the SNES context 1397 1398 Notes: 1399 For basic use of the SNES solvers the user need not explicitly call 1400 SNESSetUp(), since these actions will automatically occur during 1401 the call to SNESSolve(). However, if one wishes to control this 1402 phase separately, SNESSetUp() should be called after SNESCreate() 1403 and optional routines of the form SNESSetXXX(), but before SNESSolve(). 1404 1405 Level: advanced 1406 1407 .keywords: SNES, nonlinear, setup 1408 1409 .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 1410 @*/ 1411 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUp(SNES snes) 1412 { 1413 PetscErrorCode ierr; 1414 1415 PetscFunctionBegin; 1416 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1417 if (snes->setupcalled) PetscFunctionReturn(0); 1418 1419 if (!((PetscObject)snes)->type_name) { 1420 ierr = SNESSetType(snes,SNESLS);CHKERRQ(ierr); 1421 } 1422 1423 if (!snes->vec_func && !snes->vec_rhs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 1424 if (!snes->ops->computefunction && !snes->vec_rhs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 1425 if (snes->vec_func == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector"); 1426 if (!snes->ops->computejacobian && snes->dm) { 1427 Mat J; 1428 ISColoring coloring; 1429 MatFDColoring fd; 1430 1431 ierr = DMGetMatrix(snes->dm,MATAIJ,&J);CHKERRQ(ierr); 1432 ierr = DMGetColoring(snes->dm,IS_COLORING_GHOSTED,MATAIJ,&coloring);CHKERRQ(ierr); 1433 ierr = MatFDColoringCreate(J,coloring,&fd);CHKERRQ(ierr); 1434 ierr = MatFDColoringSetFunction(fd,(PetscErrorCode (*)(void))snes->ops->computefunction,snes->funP);CHKERRQ(ierr); 1435 ierr = SNESSetJacobian(snes,J,J,SNESDefaultComputeJacobianColor,fd);CHKERRQ(ierr); 1436 ierr = ISColoringDestroy(coloring);CHKERRQ(ierr); 1437 } 1438 if (!snes->ksp) {ierr = SNESGetKSP(snes, &snes->ksp);CHKERRQ(ierr);} 1439 1440 if (snes->ops->setup) { 1441 ierr = (*snes->ops->setup)(snes);CHKERRQ(ierr); 1442 } 1443 snes->setupcalled = PETSC_TRUE; 1444 PetscFunctionReturn(0); 1445 } 1446 1447 #undef __FUNCT__ 1448 #define __FUNCT__ "SNESDestroy" 1449 /*@ 1450 SNESDestroy - Destroys the nonlinear solver context that was created 1451 with SNESCreate(). 1452 1453 Collective on SNES 1454 1455 Input Parameter: 1456 . snes - the SNES context 1457 1458 Level: beginner 1459 1460 .keywords: SNES, nonlinear, destroy 1461 1462 .seealso: SNESCreate(), SNESSolve() 1463 @*/ 1464 PetscErrorCode PETSCSNES_DLLEXPORT SNESDestroy(SNES snes) 1465 { 1466 PetscErrorCode ierr; 1467 1468 PetscFunctionBegin; 1469 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1470 if (--((PetscObject)snes)->refct > 0) PetscFunctionReturn(0); 1471 1472 /* if memory was published with AMS then destroy it */ 1473 ierr = PetscObjectDepublish(snes);CHKERRQ(ierr); 1474 1475 if (snes->dm) {ierr = DMDestroy(snes->dm);CHKERRQ(ierr);} 1476 if (snes->ops->destroy) {ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr);} 1477 1478 if (snes->vec_rhs) {ierr = VecDestroy(snes->vec_rhs);CHKERRQ(ierr);} 1479 if (snes->vec_sol) {ierr = VecDestroy(snes->vec_sol);CHKERRQ(ierr);} 1480 if (snes->vec_func) {ierr = VecDestroy(snes->vec_func);CHKERRQ(ierr);} 1481 if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 1482 if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 1483 if (snes->ksp) {ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr);} 1484 ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr); 1485 if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);} 1486 ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 1487 if (snes->ops->convergeddestroy) {ierr = (*snes->ops->convergeddestroy)(snes->cnvP);CHKERRQ(ierr);} 1488 ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr); 1489 PetscFunctionReturn(0); 1490 } 1491 1492 /* ----------- Routines to set solver parameters ---------- */ 1493 1494 #undef __FUNCT__ 1495 #define __FUNCT__ "SNESSetLagPreconditioner" 1496 /*@ 1497 SNESSetLagPreconditioner - Determines when the preconditioner is rebuilt in the nonlinear solve. 1498 1499 Logically Collective on SNES 1500 1501 Input Parameters: 1502 + snes - the SNES context 1503 - lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time 1504 the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that 1505 1506 Options Database Keys: 1507 . -snes_lag_preconditioner <lag> 1508 1509 Notes: 1510 The default is 1 1511 The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 1512 If -1 is used before the very first nonlinear solve the preconditioner is still built because there is no previous preconditioner to use 1513 1514 Level: intermediate 1515 1516 .keywords: SNES, nonlinear, set, convergence, tolerances 1517 1518 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian() 1519 1520 @*/ 1521 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLagPreconditioner(SNES snes,PetscInt lag) 1522 { 1523 PetscFunctionBegin; 1524 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1525 if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater"); 1526 if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0"); 1527 PetscValidLogicalCollectiveInt(snes,lag,2); 1528 snes->lagpreconditioner = lag; 1529 PetscFunctionReturn(0); 1530 } 1531 1532 #undef __FUNCT__ 1533 #define __FUNCT__ "SNESGetLagPreconditioner" 1534 /*@ 1535 SNESGetLagPreconditioner - Indicates how often the preconditioner is rebuilt 1536 1537 Not Collective 1538 1539 Input Parameter: 1540 . snes - the SNES context 1541 1542 Output Parameter: 1543 . lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time 1544 the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that 1545 1546 Options Database Keys: 1547 . -snes_lag_preconditioner <lag> 1548 1549 Notes: 1550 The default is 1 1551 The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 1552 1553 Level: intermediate 1554 1555 .keywords: SNES, nonlinear, set, convergence, tolerances 1556 1557 .seealso: SNESSetTrustRegionTolerance(), SNESSetLagPreconditioner() 1558 1559 @*/ 1560 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLagPreconditioner(SNES snes,PetscInt *lag) 1561 { 1562 PetscFunctionBegin; 1563 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1564 *lag = snes->lagpreconditioner; 1565 PetscFunctionReturn(0); 1566 } 1567 1568 #undef __FUNCT__ 1569 #define __FUNCT__ "SNESSetLagJacobian" 1570 /*@ 1571 SNESSetLagJacobian - Determines when the Jacobian is rebuilt in the nonlinear solve. See SNESSetLagPreconditioner() for determining how 1572 often the preconditioner is rebuilt. 1573 1574 Logically Collective on SNES 1575 1576 Input Parameters: 1577 + snes - the SNES context 1578 - lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time 1579 the Jacobian is built etc. -2 means rebuild at next chance but then never again 1580 1581 Options Database Keys: 1582 . -snes_lag_jacobian <lag> 1583 1584 Notes: 1585 The default is 1 1586 The Jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 1587 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 1588 at the next Newton step but never again (unless it is reset to another value) 1589 1590 Level: intermediate 1591 1592 .keywords: SNES, nonlinear, set, convergence, tolerances 1593 1594 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagPreconditioner(), SNESGetLagJacobian() 1595 1596 @*/ 1597 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLagJacobian(SNES snes,PetscInt lag) 1598 { 1599 PetscFunctionBegin; 1600 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1601 if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater"); 1602 if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0"); 1603 PetscValidLogicalCollectiveInt(snes,lag,2); 1604 snes->lagjacobian = lag; 1605 PetscFunctionReturn(0); 1606 } 1607 1608 #undef __FUNCT__ 1609 #define __FUNCT__ "SNESGetLagJacobian" 1610 /*@ 1611 SNESGetLagJacobian - Indicates how often the Jacobian is rebuilt. See SNESGetLagPreconditioner() to determine when the preconditioner is rebuilt 1612 1613 Not Collective 1614 1615 Input Parameter: 1616 . snes - the SNES context 1617 1618 Output Parameter: 1619 . lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time 1620 the Jacobian is built etc. 1621 1622 Options Database Keys: 1623 . -snes_lag_jacobian <lag> 1624 1625 Notes: 1626 The default is 1 1627 The jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 1628 1629 Level: intermediate 1630 1631 .keywords: SNES, nonlinear, set, convergence, tolerances 1632 1633 .seealso: SNESSetTrustRegionTolerance(), SNESSetLagJacobian(), SNESSetLagPreconditioner(), SNESGetLagPreconditioner() 1634 1635 @*/ 1636 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLagJacobian(SNES snes,PetscInt *lag) 1637 { 1638 PetscFunctionBegin; 1639 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1640 *lag = snes->lagjacobian; 1641 PetscFunctionReturn(0); 1642 } 1643 1644 #undef __FUNCT__ 1645 #define __FUNCT__ "SNESSetTolerances" 1646 /*@ 1647 SNESSetTolerances - Sets various parameters used in convergence tests. 1648 1649 Logically Collective on SNES 1650 1651 Input Parameters: 1652 + snes - the SNES context 1653 . abstol - absolute convergence tolerance 1654 . rtol - relative convergence tolerance 1655 . stol - convergence tolerance in terms of the norm 1656 of the change in the solution between steps 1657 . maxit - maximum number of iterations 1658 - maxf - maximum number of function evaluations 1659 1660 Options Database Keys: 1661 + -snes_atol <abstol> - Sets abstol 1662 . -snes_rtol <rtol> - Sets rtol 1663 . -snes_stol <stol> - Sets stol 1664 . -snes_max_it <maxit> - Sets maxit 1665 - -snes_max_funcs <maxf> - Sets maxf 1666 1667 Notes: 1668 The default maximum number of iterations is 50. 1669 The default maximum number of function evaluations is 1000. 1670 1671 Level: intermediate 1672 1673 .keywords: SNES, nonlinear, set, convergence, tolerances 1674 1675 .seealso: SNESSetTrustRegionTolerance() 1676 @*/ 1677 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf) 1678 { 1679 PetscFunctionBegin; 1680 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1681 PetscValidLogicalCollectiveReal(snes,abstol,2); 1682 PetscValidLogicalCollectiveReal(snes,rtol,3); 1683 PetscValidLogicalCollectiveReal(snes,stol,4); 1684 PetscValidLogicalCollectiveInt(snes,maxit,5); 1685 PetscValidLogicalCollectiveInt(snes,maxf,6); 1686 1687 if (abstol != PETSC_DEFAULT) snes->abstol = abstol; 1688 if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1689 if (stol != PETSC_DEFAULT) snes->xtol = stol; 1690 if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1691 if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 1692 PetscFunctionReturn(0); 1693 } 1694 1695 #undef __FUNCT__ 1696 #define __FUNCT__ "SNESGetTolerances" 1697 /*@ 1698 SNESGetTolerances - Gets various parameters used in convergence tests. 1699 1700 Not Collective 1701 1702 Input Parameters: 1703 + snes - the SNES context 1704 . atol - absolute convergence tolerance 1705 . rtol - relative convergence tolerance 1706 . stol - convergence tolerance in terms of the norm 1707 of the change in the solution between steps 1708 . maxit - maximum number of iterations 1709 - maxf - maximum number of function evaluations 1710 1711 Notes: 1712 The user can specify PETSC_NULL for any parameter that is not needed. 1713 1714 Level: intermediate 1715 1716 .keywords: SNES, nonlinear, get, convergence, tolerances 1717 1718 .seealso: SNESSetTolerances() 1719 @*/ 1720 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf) 1721 { 1722 PetscFunctionBegin; 1723 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1724 if (atol) *atol = snes->abstol; 1725 if (rtol) *rtol = snes->rtol; 1726 if (stol) *stol = snes->xtol; 1727 if (maxit) *maxit = snes->max_its; 1728 if (maxf) *maxf = snes->max_funcs; 1729 PetscFunctionReturn(0); 1730 } 1731 1732 #undef __FUNCT__ 1733 #define __FUNCT__ "SNESSetTrustRegionTolerance" 1734 /*@ 1735 SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 1736 1737 Logically Collective on SNES 1738 1739 Input Parameters: 1740 + snes - the SNES context 1741 - tol - tolerance 1742 1743 Options Database Key: 1744 . -snes_trtol <tol> - Sets tol 1745 1746 Level: intermediate 1747 1748 .keywords: SNES, nonlinear, set, trust region, tolerance 1749 1750 .seealso: SNESSetTolerances() 1751 @*/ 1752 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTrustRegionTolerance(SNES snes,PetscReal tol) 1753 { 1754 PetscFunctionBegin; 1755 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1756 PetscValidLogicalCollectiveReal(snes,tol,2); 1757 snes->deltatol = tol; 1758 PetscFunctionReturn(0); 1759 } 1760 1761 /* 1762 Duplicate the lg monitors for SNES from KSP; for some reason with 1763 dynamic libraries things don't work under Sun4 if we just use 1764 macros instead of functions 1765 */ 1766 #undef __FUNCT__ 1767 #define __FUNCT__ "SNESMonitorLG" 1768 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLG(SNES snes,PetscInt it,PetscReal norm,void *ctx) 1769 { 1770 PetscErrorCode ierr; 1771 1772 PetscFunctionBegin; 1773 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1774 ierr = KSPMonitorLG((KSP)snes,it,norm,ctx);CHKERRQ(ierr); 1775 PetscFunctionReturn(0); 1776 } 1777 1778 #undef __FUNCT__ 1779 #define __FUNCT__ "SNESMonitorLGCreate" 1780 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1781 { 1782 PetscErrorCode ierr; 1783 1784 PetscFunctionBegin; 1785 ierr = KSPMonitorLGCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr); 1786 PetscFunctionReturn(0); 1787 } 1788 1789 #undef __FUNCT__ 1790 #define __FUNCT__ "SNESMonitorLGDestroy" 1791 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGDestroy(PetscDrawLG draw) 1792 { 1793 PetscErrorCode ierr; 1794 1795 PetscFunctionBegin; 1796 ierr = KSPMonitorLGDestroy(draw);CHKERRQ(ierr); 1797 PetscFunctionReturn(0); 1798 } 1799 1800 extern PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorRange_Private(SNES,PetscInt,PetscReal*); 1801 #undef __FUNCT__ 1802 #define __FUNCT__ "SNESMonitorLGRange" 1803 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGRange(SNES snes,PetscInt n,PetscReal rnorm,void *monctx) 1804 { 1805 PetscDrawLG lg; 1806 PetscErrorCode ierr; 1807 PetscReal x,y,per; 1808 PetscViewer v = (PetscViewer)monctx; 1809 static PetscReal prev; /* should be in the context */ 1810 PetscDraw draw; 1811 PetscFunctionBegin; 1812 if (!monctx) { 1813 MPI_Comm comm; 1814 1815 ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 1816 v = PETSC_VIEWER_DRAW_(comm); 1817 } 1818 ierr = PetscViewerDrawGetDrawLG(v,0,&lg);CHKERRQ(ierr); 1819 if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 1820 ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr); 1821 ierr = PetscDrawSetTitle(draw,"Residual norm");CHKERRQ(ierr); 1822 x = (PetscReal) n; 1823 if (rnorm > 0.0) y = log10(rnorm); else y = -15.0; 1824 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 1825 if (n < 20 || !(n % 5)) { 1826 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 1827 } 1828 1829 ierr = PetscViewerDrawGetDrawLG(v,1,&lg);CHKERRQ(ierr); 1830 if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 1831 ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr); 1832 ierr = PetscDrawSetTitle(draw,"% elemts > .2*max elemt");CHKERRQ(ierr); 1833 ierr = SNESMonitorRange_Private(snes,n,&per);CHKERRQ(ierr); 1834 x = (PetscReal) n; 1835 y = 100.0*per; 1836 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 1837 if (n < 20 || !(n % 5)) { 1838 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 1839 } 1840 1841 ierr = PetscViewerDrawGetDrawLG(v,2,&lg);CHKERRQ(ierr); 1842 if (!n) {prev = rnorm;ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 1843 ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr); 1844 ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm");CHKERRQ(ierr); 1845 x = (PetscReal) n; 1846 y = (prev - rnorm)/prev; 1847 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 1848 if (n < 20 || !(n % 5)) { 1849 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 1850 } 1851 1852 ierr = PetscViewerDrawGetDrawLG(v,3,&lg);CHKERRQ(ierr); 1853 if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 1854 ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr); 1855 ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");CHKERRQ(ierr); 1856 x = (PetscReal) n; 1857 y = (prev - rnorm)/(prev*per); 1858 if (n > 2) { /*skip initial crazy value */ 1859 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 1860 } 1861 if (n < 20 || !(n % 5)) { 1862 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 1863 } 1864 prev = rnorm; 1865 PetscFunctionReturn(0); 1866 } 1867 1868 #undef __FUNCT__ 1869 #define __FUNCT__ "SNESMonitorLGRangeCreate" 1870 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGRangeCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1871 { 1872 PetscErrorCode ierr; 1873 1874 PetscFunctionBegin; 1875 ierr = KSPMonitorLGCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr); 1876 PetscFunctionReturn(0); 1877 } 1878 1879 #undef __FUNCT__ 1880 #define __FUNCT__ "SNESMonitorLGRangeDestroy" 1881 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGRangeDestroy(PetscDrawLG draw) 1882 { 1883 PetscErrorCode ierr; 1884 1885 PetscFunctionBegin; 1886 ierr = KSPMonitorLGDestroy(draw);CHKERRQ(ierr); 1887 PetscFunctionReturn(0); 1888 } 1889 1890 /* ------------ Routines to set performance monitoring options ----------- */ 1891 1892 #undef __FUNCT__ 1893 #define __FUNCT__ "SNESMonitorSet" 1894 /*@C 1895 SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every 1896 iteration of the nonlinear solver to display the iteration's 1897 progress. 1898 1899 Logically Collective on SNES 1900 1901 Input Parameters: 1902 + snes - the SNES context 1903 . func - monitoring routine 1904 . mctx - [optional] user-defined context for private data for the 1905 monitor routine (use PETSC_NULL if no context is desired) 1906 - monitordestroy - [optional] routine that frees monitor context 1907 (may be PETSC_NULL) 1908 1909 Calling sequence of func: 1910 $ int func(SNES snes,PetscInt its, PetscReal norm,void *mctx) 1911 1912 + snes - the SNES context 1913 . its - iteration number 1914 . norm - 2-norm function value (may be estimated) 1915 - mctx - [optional] monitoring context 1916 1917 Options Database Keys: 1918 + -snes_monitor - sets SNESMonitorDefault() 1919 . -snes_monitor_draw - sets line graph monitor, 1920 uses SNESMonitorLGCreate() 1921 _ -snes_monitor_cancel - cancels all monitors that have 1922 been hardwired into a code by 1923 calls to SNESMonitorSet(), but 1924 does not cancel those set via 1925 the options database. 1926 1927 Notes: 1928 Several different monitoring routines may be set by calling 1929 SNESMonitorSet() multiple times; all will be called in the 1930 order in which they were set. 1931 1932 Fortran notes: Only a single monitor function can be set for each SNES object 1933 1934 Level: intermediate 1935 1936 .keywords: SNES, nonlinear, set, monitor 1937 1938 .seealso: SNESMonitorDefault(), SNESMonitorCancel() 1939 @*/ 1940 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorSet(SNES snes,PetscErrorCode (*monitor)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void*)) 1941 { 1942 PetscInt i; 1943 1944 PetscFunctionBegin; 1945 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1946 if (snes->numbermonitors >= MAXSNESMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1947 for (i=0; i<snes->numbermonitors;i++) { 1948 if (monitor == snes->monitor[i] && monitordestroy == snes->monitordestroy[i] && mctx == snes->monitorcontext[i]) PetscFunctionReturn(0); 1949 1950 /* check if both default monitors that share common ASCII viewer */ 1951 if (monitor == snes->monitor[i] && monitor == SNESMonitorDefault) { 1952 if (mctx && snes->monitorcontext[i]) { 1953 PetscErrorCode ierr; 1954 PetscViewerASCIIMonitor viewer1 = (PetscViewerASCIIMonitor) mctx; 1955 PetscViewerASCIIMonitor viewer2 = (PetscViewerASCIIMonitor) snes->monitorcontext[i]; 1956 if (viewer1->viewer == viewer2->viewer) { 1957 ierr = (*monitordestroy)(mctx);CHKERRQ(ierr); 1958 PetscFunctionReturn(0); 1959 } 1960 } 1961 } 1962 } 1963 snes->monitor[snes->numbermonitors] = monitor; 1964 snes->monitordestroy[snes->numbermonitors] = monitordestroy; 1965 snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 1966 PetscFunctionReturn(0); 1967 } 1968 1969 #undef __FUNCT__ 1970 #define __FUNCT__ "SNESMonitorCancel" 1971 /*@C 1972 SNESMonitorCancel - Clears all the monitor functions for a SNES object. 1973 1974 Logically Collective on SNES 1975 1976 Input Parameters: 1977 . snes - the SNES context 1978 1979 Options Database Key: 1980 . -snes_monitor_cancel - cancels all monitors that have been hardwired 1981 into a code by calls to SNESMonitorSet(), but does not cancel those 1982 set via the options database 1983 1984 Notes: 1985 There is no way to clear one specific monitor from a SNES object. 1986 1987 Level: intermediate 1988 1989 .keywords: SNES, nonlinear, set, monitor 1990 1991 .seealso: SNESMonitorDefault(), SNESMonitorSet() 1992 @*/ 1993 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorCancel(SNES snes) 1994 { 1995 PetscErrorCode ierr; 1996 PetscInt i; 1997 1998 PetscFunctionBegin; 1999 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2000 for (i=0; i<snes->numbermonitors; i++) { 2001 if (snes->monitordestroy[i]) { 2002 ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr); 2003 } 2004 } 2005 snes->numbermonitors = 0; 2006 PetscFunctionReturn(0); 2007 } 2008 2009 #undef __FUNCT__ 2010 #define __FUNCT__ "SNESSetConvergenceTest" 2011 /*@C 2012 SNESSetConvergenceTest - Sets the function that is to be used 2013 to test for convergence of the nonlinear iterative solution. 2014 2015 Logically Collective on SNES 2016 2017 Input Parameters: 2018 + snes - the SNES context 2019 . func - routine to test for convergence 2020 . cctx - [optional] context for private data for the convergence routine (may be PETSC_NULL) 2021 - destroy - [optional] destructor for the context (may be PETSC_NULL; PETSC_NULL_FUNCTION in Fortran) 2022 2023 Calling sequence of func: 2024 $ PetscErrorCode func (SNES snes,PetscInt it,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx) 2025 2026 + snes - the SNES context 2027 . it - current iteration (0 is the first and is before any Newton step) 2028 . cctx - [optional] convergence context 2029 . reason - reason for convergence/divergence 2030 . xnorm - 2-norm of current iterate 2031 . gnorm - 2-norm of current step 2032 - f - 2-norm of function 2033 2034 Level: advanced 2035 2036 .keywords: SNES, nonlinear, set, convergence, test 2037 2038 .seealso: SNESDefaultConverged(), SNESSkipConverged() 2039 @*/ 2040 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*)) 2041 { 2042 PetscErrorCode ierr; 2043 2044 PetscFunctionBegin; 2045 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2046 if (!func) func = SNESSkipConverged; 2047 if (snes->ops->convergeddestroy) { 2048 ierr = (*snes->ops->convergeddestroy)(snes->cnvP);CHKERRQ(ierr); 2049 } 2050 snes->ops->converged = func; 2051 snes->ops->convergeddestroy = destroy; 2052 snes->cnvP = cctx; 2053 PetscFunctionReturn(0); 2054 } 2055 2056 #undef __FUNCT__ 2057 #define __FUNCT__ "SNESGetConvergedReason" 2058 /*@ 2059 SNESGetConvergedReason - Gets the reason the SNES iteration was stopped. 2060 2061 Not Collective 2062 2063 Input Parameter: 2064 . snes - the SNES context 2065 2066 Output Parameter: 2067 . reason - negative value indicates diverged, positive value converged, see SNESConvergedReason or the 2068 manual pages for the individual convergence tests for complete lists 2069 2070 Level: intermediate 2071 2072 Notes: Can only be called after the call the SNESSolve() is complete. 2073 2074 .keywords: SNES, nonlinear, set, convergence, test 2075 2076 .seealso: SNESSetConvergenceTest(), SNESConvergedReason 2077 @*/ 2078 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason) 2079 { 2080 PetscFunctionBegin; 2081 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2082 PetscValidPointer(reason,2); 2083 *reason = snes->reason; 2084 PetscFunctionReturn(0); 2085 } 2086 2087 #undef __FUNCT__ 2088 #define __FUNCT__ "SNESSetConvergenceHistory" 2089 /*@ 2090 SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 2091 2092 Logically Collective on SNES 2093 2094 Input Parameters: 2095 + snes - iterative context obtained from SNESCreate() 2096 . a - array to hold history, this array will contain the function norms computed at each step 2097 . its - integer array holds the number of linear iterations for each solve. 2098 . na - size of a and its 2099 - reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero, 2100 else it continues storing new values for new nonlinear solves after the old ones 2101 2102 This routine is useful, e.g., when running a code for purposes 2103 of accurate performance monitoring, when no I/O should be done 2104 during the section of code that is being timed. 2105 2106 Level: intermediate 2107 2108 .keywords: SNES, set, convergence, history 2109 2110 .seealso: SNESGetConvergenceHistory() 2111 2112 @*/ 2113 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscTruth reset) 2114 { 2115 PetscFunctionBegin; 2116 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2117 if (na) PetscValidScalarPointer(a,2); 2118 if (its) PetscValidIntPointer(its,3); 2119 snes->conv_hist = a; 2120 snes->conv_hist_its = its; 2121 snes->conv_hist_max = na; 2122 snes->conv_hist_len = 0; 2123 snes->conv_hist_reset = reset; 2124 PetscFunctionReturn(0); 2125 } 2126 2127 #undef __FUNCT__ 2128 #define __FUNCT__ "SNESGetConvergenceHistory" 2129 /*@C 2130 SNESGetConvergenceHistory - Gets the array used to hold the convergence history. 2131 2132 Not Collective 2133 2134 Input Parameter: 2135 . snes - iterative context obtained from SNESCreate() 2136 2137 Output Parameters: 2138 . a - array to hold history 2139 . its - integer array holds the number of linear iterations (or 2140 negative if not converged) for each solve. 2141 - na - size of a and its 2142 2143 Notes: 2144 The calling sequence for this routine in Fortran is 2145 $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr) 2146 2147 This routine is useful, e.g., when running a code for purposes 2148 of accurate performance monitoring, when no I/O should be done 2149 during the section of code that is being timed. 2150 2151 Level: intermediate 2152 2153 .keywords: SNES, get, convergence, history 2154 2155 .seealso: SNESSetConvergencHistory() 2156 2157 @*/ 2158 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na) 2159 { 2160 PetscFunctionBegin; 2161 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2162 if (a) *a = snes->conv_hist; 2163 if (its) *its = snes->conv_hist_its; 2164 if (na) *na = snes->conv_hist_len; 2165 PetscFunctionReturn(0); 2166 } 2167 2168 #undef __FUNCT__ 2169 #define __FUNCT__ "SNESSetUpdate" 2170 /*@C 2171 SNESSetUpdate - Sets the general-purpose update function called 2172 at the beginning o every iteration of the nonlinear solve. Specifically 2173 it is called just before the Jacobian is "evaluated". 2174 2175 Logically Collective on SNES 2176 2177 Input Parameters: 2178 . snes - The nonlinear solver context 2179 . func - The function 2180 2181 Calling sequence of func: 2182 . func (SNES snes, PetscInt step); 2183 2184 . step - The current step of the iteration 2185 2186 Level: advanced 2187 2188 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() 2189 This is not used by most users. 2190 2191 .keywords: SNES, update 2192 2193 .seealso SNESDefaultUpdate(), SNESSetJacobian(), SNESSolve() 2194 @*/ 2195 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt)) 2196 { 2197 PetscFunctionBegin; 2198 PetscValidHeaderSpecific(snes, SNES_CLASSID,1); 2199 snes->ops->update = func; 2200 PetscFunctionReturn(0); 2201 } 2202 2203 #undef __FUNCT__ 2204 #define __FUNCT__ "SNESDefaultUpdate" 2205 /*@ 2206 SNESDefaultUpdate - The default update function which does nothing. 2207 2208 Not collective 2209 2210 Input Parameters: 2211 . snes - The nonlinear solver context 2212 . step - The current step of the iteration 2213 2214 Level: intermediate 2215 2216 .keywords: SNES, update 2217 .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultShortolutionBC() 2218 @*/ 2219 PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultUpdate(SNES snes, PetscInt step) 2220 { 2221 PetscFunctionBegin; 2222 PetscFunctionReturn(0); 2223 } 2224 2225 #undef __FUNCT__ 2226 #define __FUNCT__ "SNESScaleStep_Private" 2227 /* 2228 SNESScaleStep_Private - Scales a step so that its length is less than the 2229 positive parameter delta. 2230 2231 Input Parameters: 2232 + snes - the SNES context 2233 . y - approximate solution of linear system 2234 . fnorm - 2-norm of current function 2235 - delta - trust region size 2236 2237 Output Parameters: 2238 + gpnorm - predicted function norm at the new point, assuming local 2239 linearization. The value is zero if the step lies within the trust 2240 region, and exceeds zero otherwise. 2241 - ynorm - 2-norm of the step 2242 2243 Note: 2244 For non-trust region methods such as SNESLS, the parameter delta 2245 is set to be the maximum allowable step size. 2246 2247 .keywords: SNES, nonlinear, scale, step 2248 */ 2249 PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm) 2250 { 2251 PetscReal nrm; 2252 PetscScalar cnorm; 2253 PetscErrorCode ierr; 2254 2255 PetscFunctionBegin; 2256 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2257 PetscValidHeaderSpecific(y,VEC_CLASSID,2); 2258 PetscCheckSameComm(snes,1,y,2); 2259 2260 ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr); 2261 if (nrm > *delta) { 2262 nrm = *delta/nrm; 2263 *gpnorm = (1.0 - nrm)*(*fnorm); 2264 cnorm = nrm; 2265 ierr = VecScale(y,cnorm);CHKERRQ(ierr); 2266 *ynorm = *delta; 2267 } else { 2268 *gpnorm = 0.0; 2269 *ynorm = nrm; 2270 } 2271 PetscFunctionReturn(0); 2272 } 2273 2274 #undef __FUNCT__ 2275 #define __FUNCT__ "SNESSolve" 2276 /*@C 2277 SNESSolve - Solves a nonlinear system F(x) = b. 2278 Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX(). 2279 2280 Collective on SNES 2281 2282 Input Parameters: 2283 + snes - the SNES context 2284 . b - the constant part of the equation, or PETSC_NULL to use zero. 2285 - x - the solution vector. 2286 2287 Notes: 2288 The user should initialize the vector,x, with the initial guess 2289 for the nonlinear solve prior to calling SNESSolve. In particular, 2290 to employ an initial guess of zero, the user should explicitly set 2291 this vector to zero by calling VecSet(). 2292 2293 Level: beginner 2294 2295 .keywords: SNES, nonlinear, solve 2296 2297 .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian() 2298 @*/ 2299 PetscErrorCode PETSCSNES_DLLEXPORT SNESSolve(SNES snes,Vec b,Vec x) 2300 { 2301 PetscErrorCode ierr; 2302 PetscTruth flg; 2303 char filename[PETSC_MAX_PATH_LEN]; 2304 PetscViewer viewer; 2305 2306 PetscFunctionBegin; 2307 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2308 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 2309 PetscCheckSameComm(snes,1,x,3); 2310 if (b) PetscValidHeaderSpecific(b,VEC_CLASSID,2); 2311 if (b) PetscCheckSameComm(snes,1,b,2); 2312 2313 /* set solution vector */ 2314 ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 2315 if (snes->vec_sol) { ierr = VecDestroy(snes->vec_sol);CHKERRQ(ierr); } 2316 snes->vec_sol = x; 2317 /* set afine vector if provided */ 2318 if (b) { ierr = PetscObjectReference((PetscObject)b);CHKERRQ(ierr); } 2319 if (snes->vec_rhs) { ierr = VecDestroy(snes->vec_rhs);CHKERRQ(ierr); } 2320 snes->vec_rhs = b; 2321 2322 if (!snes->vec_func && snes->vec_rhs) { 2323 ierr = VecDuplicate(b, &snes->vec_func);CHKERRQ(ierr); 2324 } 2325 2326 ierr = SNESSetUp(snes);CHKERRQ(ierr); 2327 2328 if (snes->conv_hist_reset) snes->conv_hist_len = 0; 2329 snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0; 2330 2331 ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 2332 ierr = (*snes->ops->solve)(snes);CHKERRQ(ierr); 2333 ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 2334 if (snes->domainerror){ 2335 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 2336 snes->domainerror = PETSC_FALSE; 2337 } 2338 if (!snes->reason) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason"); 2339 2340 ierr = PetscOptionsGetString(((PetscObject)snes)->prefix,"-snes_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 2341 if (flg && !PetscPreLoadingOn) { 2342 ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,filename,&viewer);CHKERRQ(ierr); 2343 ierr = SNESView(snes,viewer);CHKERRQ(ierr); 2344 ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr); 2345 } 2346 2347 flg = PETSC_FALSE; 2348 ierr = PetscOptionsGetTruth(((PetscObject)snes)->prefix,"-snes_test_local_min",&flg,PETSC_NULL);CHKERRQ(ierr); 2349 if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); } 2350 if (snes->printreason) { 2351 if (snes->reason > 0) { 2352 ierr = PetscPrintf(((PetscObject)snes)->comm,"Nonlinear solve converged due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr); 2353 } else { 2354 ierr = PetscPrintf(((PetscObject)snes)->comm,"Nonlinear solve did not converge due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr); 2355 } 2356 } 2357 2358 if (snes->errorifnotconverged && snes->reason < 0) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_NOT_CONVERGED,"SNESSolve has not converged"); 2359 PetscFunctionReturn(0); 2360 } 2361 2362 /* --------- Internal routines for SNES Package --------- */ 2363 2364 #undef __FUNCT__ 2365 #define __FUNCT__ "SNESSetType" 2366 /*@C 2367 SNESSetType - Sets the method for the nonlinear solver. 2368 2369 Collective on SNES 2370 2371 Input Parameters: 2372 + snes - the SNES context 2373 - type - a known method 2374 2375 Options Database Key: 2376 . -snes_type <type> - Sets the method; use -help for a list 2377 of available methods (for instance, ls or tr) 2378 2379 Notes: 2380 See "petsc/include/petscsnes.h" for available methods (for instance) 2381 + SNESLS - Newton's method with line search 2382 (systems of nonlinear equations) 2383 . SNESTR - Newton's method with trust region 2384 (systems of nonlinear equations) 2385 2386 Normally, it is best to use the SNESSetFromOptions() command and then 2387 set the SNES solver type from the options database rather than by using 2388 this routine. Using the options database provides the user with 2389 maximum flexibility in evaluating the many nonlinear solvers. 2390 The SNESSetType() routine is provided for those situations where it 2391 is necessary to set the nonlinear solver independently of the command 2392 line or options database. This might be the case, for example, when 2393 the choice of solver changes during the execution of the program, 2394 and the user's application is taking responsibility for choosing the 2395 appropriate method. 2396 2397 Level: intermediate 2398 2399 .keywords: SNES, set, type 2400 2401 .seealso: SNESType, SNESCreate() 2402 2403 @*/ 2404 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetType(SNES snes,const SNESType type) 2405 { 2406 PetscErrorCode ierr,(*r)(SNES); 2407 PetscTruth match; 2408 2409 PetscFunctionBegin; 2410 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2411 PetscValidCharPointer(type,2); 2412 2413 ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr); 2414 if (match) PetscFunctionReturn(0); 2415 2416 ierr = PetscFListFind(SNESList,((PetscObject)snes)->comm,type,(void (**)(void)) &r);CHKERRQ(ierr); 2417 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type); 2418 /* Destroy the previous private SNES context */ 2419 if (snes->ops->destroy) { ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr); } 2420 /* Reinitialize function pointers in SNESOps structure */ 2421 snes->ops->setup = 0; 2422 snes->ops->solve = 0; 2423 snes->ops->view = 0; 2424 snes->ops->setfromoptions = 0; 2425 snes->ops->destroy = 0; 2426 /* Call the SNESCreate_XXX routine for this particular Nonlinear solver */ 2427 snes->setupcalled = PETSC_FALSE; 2428 ierr = (*r)(snes);CHKERRQ(ierr); 2429 ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr); 2430 #if defined(PETSC_HAVE_AMS) 2431 if (PetscAMSPublishAll) { 2432 ierr = PetscObjectAMSPublish((PetscObject)snes);CHKERRQ(ierr); 2433 } 2434 #endif 2435 PetscFunctionReturn(0); 2436 } 2437 2438 2439 /* --------------------------------------------------------------------- */ 2440 #undef __FUNCT__ 2441 #define __FUNCT__ "SNESRegisterDestroy" 2442 /*@ 2443 SNESRegisterDestroy - Frees the list of nonlinear solvers that were 2444 registered by SNESRegisterDynamic(). 2445 2446 Not Collective 2447 2448 Level: advanced 2449 2450 .keywords: SNES, nonlinear, register, destroy 2451 2452 .seealso: SNESRegisterAll(), SNESRegisterAll() 2453 @*/ 2454 PetscErrorCode PETSCSNES_DLLEXPORT SNESRegisterDestroy(void) 2455 { 2456 PetscErrorCode ierr; 2457 2458 PetscFunctionBegin; 2459 ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr); 2460 SNESRegisterAllCalled = PETSC_FALSE; 2461 PetscFunctionReturn(0); 2462 } 2463 2464 #undef __FUNCT__ 2465 #define __FUNCT__ "SNESGetType" 2466 /*@C 2467 SNESGetType - Gets the SNES method type and name (as a string). 2468 2469 Not Collective 2470 2471 Input Parameter: 2472 . snes - nonlinear solver context 2473 2474 Output Parameter: 2475 . type - SNES method (a character string) 2476 2477 Level: intermediate 2478 2479 .keywords: SNES, nonlinear, get, type, name 2480 @*/ 2481 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetType(SNES snes,const SNESType *type) 2482 { 2483 PetscFunctionBegin; 2484 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2485 PetscValidPointer(type,2); 2486 *type = ((PetscObject)snes)->type_name; 2487 PetscFunctionReturn(0); 2488 } 2489 2490 #undef __FUNCT__ 2491 #define __FUNCT__ "SNESGetSolution" 2492 /*@ 2493 SNESGetSolution - Returns the vector where the approximate solution is 2494 stored. 2495 2496 Not Collective, but Vec is parallel if SNES is parallel 2497 2498 Input Parameter: 2499 . snes - the SNES context 2500 2501 Output Parameter: 2502 . x - the solution 2503 2504 Level: intermediate 2505 2506 .keywords: SNES, nonlinear, get, solution 2507 2508 .seealso: SNESGetSolutionUpdate(), SNESGetFunction() 2509 @*/ 2510 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolution(SNES snes,Vec *x) 2511 { 2512 PetscFunctionBegin; 2513 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2514 PetscValidPointer(x,2); 2515 *x = snes->vec_sol; 2516 PetscFunctionReturn(0); 2517 } 2518 2519 #undef __FUNCT__ 2520 #define __FUNCT__ "SNESGetSolutionUpdate" 2521 /*@ 2522 SNESGetSolutionUpdate - Returns the vector where the solution update is 2523 stored. 2524 2525 Not Collective, but Vec is parallel if SNES is parallel 2526 2527 Input Parameter: 2528 . snes - the SNES context 2529 2530 Output Parameter: 2531 . x - the solution update 2532 2533 Level: advanced 2534 2535 .keywords: SNES, nonlinear, get, solution, update 2536 2537 .seealso: SNESGetSolution(), SNESGetFunction() 2538 @*/ 2539 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolutionUpdate(SNES snes,Vec *x) 2540 { 2541 PetscFunctionBegin; 2542 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2543 PetscValidPointer(x,2); 2544 *x = snes->vec_sol_update; 2545 PetscFunctionReturn(0); 2546 } 2547 2548 #undef __FUNCT__ 2549 #define __FUNCT__ "SNESGetFunction" 2550 /*@C 2551 SNESGetFunction - Returns the vector where the function is stored. 2552 2553 Not Collective, but Vec is parallel if SNES is parallel 2554 2555 Input Parameter: 2556 . snes - the SNES context 2557 2558 Output Parameter: 2559 + r - the function (or PETSC_NULL) 2560 . func - the function (or PETSC_NULL) 2561 - ctx - the function context (or PETSC_NULL) 2562 2563 Level: advanced 2564 2565 .keywords: SNES, nonlinear, get, function 2566 2567 .seealso: SNESSetFunction(), SNESGetSolution() 2568 @*/ 2569 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx) 2570 { 2571 PetscFunctionBegin; 2572 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2573 if (r) *r = snes->vec_func; 2574 if (func) *func = snes->ops->computefunction; 2575 if (ctx) *ctx = snes->funP; 2576 PetscFunctionReturn(0); 2577 } 2578 2579 #undef __FUNCT__ 2580 #define __FUNCT__ "SNESSetOptionsPrefix" 2581 /*@C 2582 SNESSetOptionsPrefix - Sets the prefix used for searching for all 2583 SNES options in the database. 2584 2585 Logically Collective on SNES 2586 2587 Input Parameter: 2588 + snes - the SNES context 2589 - prefix - the prefix to prepend to all option names 2590 2591 Notes: 2592 A hyphen (-) must NOT be given at the beginning of the prefix name. 2593 The first character of all runtime options is AUTOMATICALLY the hyphen. 2594 2595 Level: advanced 2596 2597 .keywords: SNES, set, options, prefix, database 2598 2599 .seealso: SNESSetFromOptions() 2600 @*/ 2601 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetOptionsPrefix(SNES snes,const char prefix[]) 2602 { 2603 PetscErrorCode ierr; 2604 2605 PetscFunctionBegin; 2606 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2607 ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2608 if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);} 2609 ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 2610 PetscFunctionReturn(0); 2611 } 2612 2613 #undef __FUNCT__ 2614 #define __FUNCT__ "SNESAppendOptionsPrefix" 2615 /*@C 2616 SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 2617 SNES options in the database. 2618 2619 Logically Collective on SNES 2620 2621 Input Parameters: 2622 + snes - the SNES context 2623 - prefix - the prefix to prepend to all option names 2624 2625 Notes: 2626 A hyphen (-) must NOT be given at the beginning of the prefix name. 2627 The first character of all runtime options is AUTOMATICALLY the hyphen. 2628 2629 Level: advanced 2630 2631 .keywords: SNES, append, options, prefix, database 2632 2633 .seealso: SNESGetOptionsPrefix() 2634 @*/ 2635 PetscErrorCode PETSCSNES_DLLEXPORT SNESAppendOptionsPrefix(SNES snes,const char prefix[]) 2636 { 2637 PetscErrorCode ierr; 2638 2639 PetscFunctionBegin; 2640 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2641 ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2642 if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);} 2643 ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 2644 PetscFunctionReturn(0); 2645 } 2646 2647 #undef __FUNCT__ 2648 #define __FUNCT__ "SNESGetOptionsPrefix" 2649 /*@C 2650 SNESGetOptionsPrefix - Sets the prefix used for searching for all 2651 SNES options in the database. 2652 2653 Not Collective 2654 2655 Input Parameter: 2656 . snes - the SNES context 2657 2658 Output Parameter: 2659 . prefix - pointer to the prefix string used 2660 2661 Notes: On the fortran side, the user should pass in a string 'prefix' of 2662 sufficient length to hold the prefix. 2663 2664 Level: advanced 2665 2666 .keywords: SNES, get, options, prefix, database 2667 2668 .seealso: SNESAppendOptionsPrefix() 2669 @*/ 2670 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetOptionsPrefix(SNES snes,const char *prefix[]) 2671 { 2672 PetscErrorCode ierr; 2673 2674 PetscFunctionBegin; 2675 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2676 ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2677 PetscFunctionReturn(0); 2678 } 2679 2680 2681 #undef __FUNCT__ 2682 #define __FUNCT__ "SNESRegister" 2683 /*@C 2684 SNESRegister - See SNESRegisterDynamic() 2685 2686 Level: advanced 2687 @*/ 2688 PetscErrorCode PETSCSNES_DLLEXPORT SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES)) 2689 { 2690 char fullname[PETSC_MAX_PATH_LEN]; 2691 PetscErrorCode ierr; 2692 2693 PetscFunctionBegin; 2694 ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 2695 ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 2696 PetscFunctionReturn(0); 2697 } 2698 2699 #undef __FUNCT__ 2700 #define __FUNCT__ "SNESTestLocalMin" 2701 PetscErrorCode PETSCSNES_DLLEXPORT SNESTestLocalMin(SNES snes) 2702 { 2703 PetscErrorCode ierr; 2704 PetscInt N,i,j; 2705 Vec u,uh,fh; 2706 PetscScalar value; 2707 PetscReal norm; 2708 2709 PetscFunctionBegin; 2710 ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr); 2711 ierr = VecDuplicate(u,&uh);CHKERRQ(ierr); 2712 ierr = VecDuplicate(u,&fh);CHKERRQ(ierr); 2713 2714 /* currently only works for sequential */ 2715 ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n"); 2716 ierr = VecGetSize(u,&N);CHKERRQ(ierr); 2717 for (i=0; i<N; i++) { 2718 ierr = VecCopy(u,uh);CHKERRQ(ierr); 2719 ierr = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr); 2720 for (j=-10; j<11; j++) { 2721 value = PetscSign(j)*exp(PetscAbs(j)-10.0); 2722 ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2723 ierr = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr); 2724 ierr = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr); 2725 ierr = PetscPrintf(PETSC_COMM_WORLD," j norm %D %18.16e\n",j,norm);CHKERRQ(ierr); 2726 value = -value; 2727 ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2728 } 2729 } 2730 ierr = VecDestroy(uh);CHKERRQ(ierr); 2731 ierr = VecDestroy(fh);CHKERRQ(ierr); 2732 PetscFunctionReturn(0); 2733 } 2734 2735 #undef __FUNCT__ 2736 #define __FUNCT__ "SNESKSPSetUseEW" 2737 /*@ 2738 SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for 2739 computing relative tolerance for linear solvers within an inexact 2740 Newton method. 2741 2742 Logically Collective on SNES 2743 2744 Input Parameters: 2745 + snes - SNES context 2746 - flag - PETSC_TRUE or PETSC_FALSE 2747 2748 Options Database: 2749 + -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence 2750 . -snes_ksp_ew_version ver - version of Eisenstat-Walker method 2751 . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0 2752 . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax 2753 . -snes_ksp_ew_gamma <gamma> - Sets gamma 2754 . -snes_ksp_ew_alpha <alpha> - Sets alpha 2755 . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2 2756 - -snes_ksp_ew_threshold <threshold> - Sets threshold 2757 2758 Notes: 2759 Currently, the default is to use a constant relative tolerance for 2760 the inner linear solvers. Alternatively, one can use the 2761 Eisenstat-Walker method, where the relative convergence tolerance 2762 is reset at each Newton iteration according progress of the nonlinear 2763 solver. 2764 2765 Level: advanced 2766 2767 Reference: 2768 S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an 2769 inexact Newton method", SISC 17 (1), pp.16-32, 1996. 2770 2771 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton 2772 2773 .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW() 2774 @*/ 2775 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPSetUseEW(SNES snes,PetscTruth flag) 2776 { 2777 PetscFunctionBegin; 2778 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2779 PetscValidLogicalCollectiveTruth(snes,flag,2); 2780 snes->ksp_ewconv = flag; 2781 PetscFunctionReturn(0); 2782 } 2783 2784 #undef __FUNCT__ 2785 #define __FUNCT__ "SNESKSPGetUseEW" 2786 /*@ 2787 SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method 2788 for computing relative tolerance for linear solvers within an 2789 inexact Newton method. 2790 2791 Not Collective 2792 2793 Input Parameter: 2794 . snes - SNES context 2795 2796 Output Parameter: 2797 . flag - PETSC_TRUE or PETSC_FALSE 2798 2799 Notes: 2800 Currently, the default is to use a constant relative tolerance for 2801 the inner linear solvers. Alternatively, one can use the 2802 Eisenstat-Walker method, where the relative convergence tolerance 2803 is reset at each Newton iteration according progress of the nonlinear 2804 solver. 2805 2806 Level: advanced 2807 2808 Reference: 2809 S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an 2810 inexact Newton method", SISC 17 (1), pp.16-32, 1996. 2811 2812 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton 2813 2814 .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW() 2815 @*/ 2816 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPGetUseEW(SNES snes, PetscTruth *flag) 2817 { 2818 PetscFunctionBegin; 2819 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2820 PetscValidPointer(flag,2); 2821 *flag = snes->ksp_ewconv; 2822 PetscFunctionReturn(0); 2823 } 2824 2825 #undef __FUNCT__ 2826 #define __FUNCT__ "SNESKSPSetParametersEW" 2827 /*@ 2828 SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker 2829 convergence criteria for the linear solvers within an inexact 2830 Newton method. 2831 2832 Logically Collective on SNES 2833 2834 Input Parameters: 2835 + snes - SNES context 2836 . version - version 1, 2 (default is 2) or 3 2837 . rtol_0 - initial relative tolerance (0 <= rtol_0 < 1) 2838 . rtol_max - maximum relative tolerance (0 <= rtol_max < 1) 2839 . gamma - multiplicative factor for version 2 rtol computation 2840 (0 <= gamma2 <= 1) 2841 . alpha - power for version 2 rtol computation (1 < alpha <= 2) 2842 . alpha2 - power for safeguard 2843 - threshold - threshold for imposing safeguard (0 < threshold < 1) 2844 2845 Note: 2846 Version 3 was contributed by Luis Chacon, June 2006. 2847 2848 Use PETSC_DEFAULT to retain the default for any of the parameters. 2849 2850 Level: advanced 2851 2852 Reference: 2853 S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an 2854 inexact Newton method", Utah State University Math. Stat. Dept. Res. 2855 Report 6/94/75, June, 1994, to appear in SIAM J. Sci. Comput. 2856 2857 .keywords: SNES, KSP, Eisenstat, Walker, set, parameters 2858 2859 .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW() 2860 @*/ 2861 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max, 2862 PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold) 2863 { 2864 SNESKSPEW *kctx; 2865 PetscFunctionBegin; 2866 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2867 kctx = (SNESKSPEW*)snes->kspconvctx; 2868 if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing"); 2869 PetscValidLogicalCollectiveInt(snes,version,2); 2870 PetscValidLogicalCollectiveReal(snes,rtol_0,3); 2871 PetscValidLogicalCollectiveReal(snes,rtol_max,4); 2872 PetscValidLogicalCollectiveReal(snes,gamma,5); 2873 PetscValidLogicalCollectiveReal(snes,alpha,6); 2874 PetscValidLogicalCollectiveReal(snes,alpha2,7); 2875 PetscValidLogicalCollectiveReal(snes,threshold,8); 2876 2877 if (version != PETSC_DEFAULT) kctx->version = version; 2878 if (rtol_0 != PETSC_DEFAULT) kctx->rtol_0 = rtol_0; 2879 if (rtol_max != PETSC_DEFAULT) kctx->rtol_max = rtol_max; 2880 if (gamma != PETSC_DEFAULT) kctx->gamma = gamma; 2881 if (alpha != PETSC_DEFAULT) kctx->alpha = alpha; 2882 if (alpha2 != PETSC_DEFAULT) kctx->alpha2 = alpha2; 2883 if (threshold != PETSC_DEFAULT) kctx->threshold = threshold; 2884 2885 if (kctx->version < 1 || kctx->version > 3) { 2886 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 and 3 are supported: %D",kctx->version); 2887 } 2888 if (kctx->rtol_0 < 0.0 || kctx->rtol_0 >= 1.0) { 2889 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_0 < 1.0: %G",kctx->rtol_0); 2890 } 2891 if (kctx->rtol_max < 0.0 || kctx->rtol_max >= 1.0) { 2892 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_max (%G) < 1.0\n",kctx->rtol_max); 2893 } 2894 if (kctx->gamma < 0.0 || kctx->gamma > 1.0) { 2895 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= gamma (%G) <= 1.0\n",kctx->gamma); 2896 } 2897 if (kctx->alpha <= 1.0 || kctx->alpha > 2.0) { 2898 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"1.0 < alpha (%G) <= 2.0\n",kctx->alpha); 2899 } 2900 if (kctx->threshold <= 0.0 || kctx->threshold >= 1.0) { 2901 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 < threshold (%G) < 1.0\n",kctx->threshold); 2902 } 2903 PetscFunctionReturn(0); 2904 } 2905 2906 #undef __FUNCT__ 2907 #define __FUNCT__ "SNESKSPGetParametersEW" 2908 /*@ 2909 SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker 2910 convergence criteria for the linear solvers within an inexact 2911 Newton method. 2912 2913 Not Collective 2914 2915 Input Parameters: 2916 snes - SNES context 2917 2918 Output Parameters: 2919 + version - version 1, 2 (default is 2) or 3 2920 . rtol_0 - initial relative tolerance (0 <= rtol_0 < 1) 2921 . rtol_max - maximum relative tolerance (0 <= rtol_max < 1) 2922 . gamma - multiplicative factor for version 2 rtol computation 2923 (0 <= gamma2 <= 1) 2924 . alpha - power for version 2 rtol computation (1 < alpha <= 2) 2925 . alpha2 - power for safeguard 2926 - threshold - threshold for imposing safeguard (0 < threshold < 1) 2927 2928 Level: advanced 2929 2930 .keywords: SNES, KSP, Eisenstat, Walker, get, parameters 2931 2932 .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW() 2933 @*/ 2934 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max, 2935 PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold) 2936 { 2937 SNESKSPEW *kctx; 2938 PetscFunctionBegin; 2939 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2940 kctx = (SNESKSPEW*)snes->kspconvctx; 2941 if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing"); 2942 if(version) *version = kctx->version; 2943 if(rtol_0) *rtol_0 = kctx->rtol_0; 2944 if(rtol_max) *rtol_max = kctx->rtol_max; 2945 if(gamma) *gamma = kctx->gamma; 2946 if(alpha) *alpha = kctx->alpha; 2947 if(alpha2) *alpha2 = kctx->alpha2; 2948 if(threshold) *threshold = kctx->threshold; 2949 PetscFunctionReturn(0); 2950 } 2951 2952 #undef __FUNCT__ 2953 #define __FUNCT__ "SNESKSPEW_PreSolve" 2954 static PetscErrorCode SNESKSPEW_PreSolve(SNES snes, KSP ksp, Vec b, Vec x) 2955 { 2956 PetscErrorCode ierr; 2957 SNESKSPEW *kctx = (SNESKSPEW*)snes->kspconvctx; 2958 PetscReal rtol=PETSC_DEFAULT,stol; 2959 2960 PetscFunctionBegin; 2961 if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists"); 2962 if (!snes->iter) { /* first time in, so use the original user rtol */ 2963 rtol = kctx->rtol_0; 2964 } else { 2965 if (kctx->version == 1) { 2966 rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last; 2967 if (rtol < 0.0) rtol = -rtol; 2968 stol = pow(kctx->rtol_last,kctx->alpha2); 2969 if (stol > kctx->threshold) rtol = PetscMax(rtol,stol); 2970 } else if (kctx->version == 2) { 2971 rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha); 2972 stol = kctx->gamma * pow(kctx->rtol_last,kctx->alpha); 2973 if (stol > kctx->threshold) rtol = PetscMax(rtol,stol); 2974 } else if (kctx->version == 3) {/* contributed by Luis Chacon, June 2006. */ 2975 rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha); 2976 /* safeguard: avoid sharp decrease of rtol */ 2977 stol = kctx->gamma*pow(kctx->rtol_last,kctx->alpha); 2978 stol = PetscMax(rtol,stol); 2979 rtol = PetscMin(kctx->rtol_0,stol); 2980 /* safeguard: avoid oversolving */ 2981 stol = kctx->gamma*(snes->ttol)/snes->norm; 2982 stol = PetscMax(rtol,stol); 2983 rtol = PetscMin(kctx->rtol_0,stol); 2984 } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 or 3 are supported: %D",kctx->version); 2985 } 2986 /* safeguard: avoid rtol greater than one */ 2987 rtol = PetscMin(rtol,kctx->rtol_max); 2988 ierr = KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 2989 ierr = PetscInfo3(snes,"iter %D, Eisenstat-Walker (version %D) KSP rtol=%G\n",snes->iter,kctx->version,rtol);CHKERRQ(ierr); 2990 PetscFunctionReturn(0); 2991 } 2992 2993 #undef __FUNCT__ 2994 #define __FUNCT__ "SNESKSPEW_PostSolve" 2995 static PetscErrorCode SNESKSPEW_PostSolve(SNES snes, KSP ksp, Vec b, Vec x) 2996 { 2997 PetscErrorCode ierr; 2998 SNESKSPEW *kctx = (SNESKSPEW*)snes->kspconvctx; 2999 PCSide pcside; 3000 Vec lres; 3001 3002 PetscFunctionBegin; 3003 if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists"); 3004 ierr = KSPGetTolerances(ksp,&kctx->rtol_last,0,0,0);CHKERRQ(ierr); 3005 ierr = SNESGetFunctionNorm(snes,&kctx->norm_last);CHKERRQ(ierr); 3006 if (kctx->version == 1) { 3007 ierr = KSPGetPCSide(ksp,&pcside);CHKERRQ(ierr); 3008 if (pcside == PC_RIGHT) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */ 3009 /* KSP residual is true linear residual */ 3010 ierr = KSPGetResidualNorm(ksp,&kctx->lresid_last);CHKERRQ(ierr); 3011 } else { 3012 /* KSP residual is preconditioned residual */ 3013 /* compute true linear residual norm */ 3014 ierr = VecDuplicate(b,&lres);CHKERRQ(ierr); 3015 ierr = MatMult(snes->jacobian,x,lres);CHKERRQ(ierr); 3016 ierr = VecAYPX(lres,-1.0,b);CHKERRQ(ierr); 3017 ierr = VecNorm(lres,NORM_2,&kctx->lresid_last);CHKERRQ(ierr); 3018 ierr = VecDestroy(lres);CHKERRQ(ierr); 3019 } 3020 } 3021 PetscFunctionReturn(0); 3022 } 3023 3024 #undef __FUNCT__ 3025 #define __FUNCT__ "SNES_KSPSolve" 3026 PetscErrorCode SNES_KSPSolve(SNES snes, KSP ksp, Vec b, Vec x) 3027 { 3028 PetscErrorCode ierr; 3029 3030 PetscFunctionBegin; 3031 if (snes->ksp_ewconv) { ierr = SNESKSPEW_PreSolve(snes,ksp,b,x);CHKERRQ(ierr); } 3032 ierr = KSPSolve(ksp,b,x);CHKERRQ(ierr); 3033 if (snes->ksp_ewconv) { ierr = SNESKSPEW_PostSolve(snes,ksp,b,x);CHKERRQ(ierr); } 3034 PetscFunctionReturn(0); 3035 } 3036 3037 #undef __FUNCT__ 3038 #define __FUNCT__ "SNESSetDM" 3039 /*@ 3040 SNESSetDM - Sets the DM that may be used by some preconditioners 3041 3042 Logically Collective on SNES 3043 3044 Input Parameters: 3045 + snes - the preconditioner context 3046 - dm - the dm 3047 3048 Level: intermediate 3049 3050 3051 .seealso: SNESGetDM(), KSPSetDM(), KSPGetDM() 3052 @*/ 3053 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetDM(SNES snes,DM dm) 3054 { 3055 PetscErrorCode ierr; 3056 KSP ksp; 3057 3058 PetscFunctionBegin; 3059 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3060 ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr); 3061 if (snes->dm) {ierr = DMDestroy(snes->dm);CHKERRQ(ierr);} 3062 snes->dm = dm; 3063 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 3064 ierr = KSPSetDM(ksp,dm);CHKERRQ(ierr); 3065 ierr = KSPSetDMActive(ksp,PETSC_FALSE);CHKERRQ(ierr); 3066 PetscFunctionReturn(0); 3067 } 3068 3069 #undef __FUNCT__ 3070 #define __FUNCT__ "SNESGetDM" 3071 /*@ 3072 SNESGetDM - Gets the DM that may be used by some preconditioners 3073 3074 Not Collective but DM obtained is parallel on SNES 3075 3076 Input Parameter: 3077 . snes - the preconditioner context 3078 3079 Output Parameter: 3080 . dm - the dm 3081 3082 Level: intermediate 3083 3084 3085 .seealso: SNESSetDM(), KSPSetDM(), KSPGetDM() 3086 @*/ 3087 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetDM(SNES snes,DM *dm) 3088 { 3089 PetscFunctionBegin; 3090 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3091 *dm = snes->dm; 3092 PetscFunctionReturn(0); 3093 } 3094