1 2 #include <private/snesimpl.h> /*I "petscsnes.h" I*/ 3 4 PetscBool SNESRegisterAllCalled = PETSC_FALSE; 5 PetscFList SNESList = PETSC_NULL; 6 7 /* Logging support */ 8 PetscClassId SNES_CLASSID; 9 PetscLogEvent SNES_Solve, SNES_LineSearch, SNES_FunctionEval, SNES_JacobianEval, SNES_GSEval; 10 11 #undef __FUNCT__ 12 #define __FUNCT__ "SNESDMComputeJacobian" 13 /* 14 Translates from a SNES call to a DM call in computing a Jacobian 15 */ 16 PetscErrorCode SNESDMComputeJacobian(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr) 17 { 18 PetscErrorCode ierr; 19 DM dm; 20 21 PetscFunctionBegin; 22 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 23 ierr = DMComputeJacobian(dm,X,*J,*B,flag);CHKERRQ(ierr); 24 PetscFunctionReturn(0); 25 } 26 27 #undef __FUNCT__ 28 #define __FUNCT__ "SNESSetErrorIfNotConverged" 29 /*@ 30 SNESSetErrorIfNotConverged - Causes SNESSolve() to generate an error if the solver has not converged. 31 32 Logically Collective on SNES 33 34 Input Parameters: 35 + snes - iterative context obtained from SNESCreate() 36 - flg - PETSC_TRUE indicates you want the error generated 37 38 Options database keys: 39 . -snes_error_if_not_converged : this takes an optional truth value (0/1/no/yes/true/false) 40 41 Level: intermediate 42 43 Notes: 44 Normally PETSc continues if a linear solver fails to converge, you can call SNESGetConvergedReason() after a SNESSolve() 45 to determine if it has converged. 46 47 .keywords: SNES, set, initial guess, nonzero 48 49 .seealso: SNESGetErrorIfNotConverged(), KSPGetErrorIfNotConverged(), KSPSetErrorIFNotConverged() 50 @*/ 51 PetscErrorCode SNESSetErrorIfNotConverged(SNES snes,PetscBool flg) 52 { 53 PetscFunctionBegin; 54 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 55 PetscValidLogicalCollectiveBool(snes,flg,2); 56 snes->errorifnotconverged = flg; 57 58 PetscFunctionReturn(0); 59 } 60 61 #undef __FUNCT__ 62 #define __FUNCT__ "SNESGetErrorIfNotConverged" 63 /*@ 64 SNESGetErrorIfNotConverged - Will SNESSolve() generate an error if the solver does not converge? 65 66 Not Collective 67 68 Input Parameter: 69 . snes - iterative context obtained from SNESCreate() 70 71 Output Parameter: 72 . flag - PETSC_TRUE if it will generate an error, else PETSC_FALSE 73 74 Level: intermediate 75 76 .keywords: SNES, set, initial guess, nonzero 77 78 .seealso: SNESSetErrorIfNotConverged(), KSPGetErrorIfNotConverged(), KSPSetErrorIFNotConverged() 79 @*/ 80 PetscErrorCode SNESGetErrorIfNotConverged(SNES snes,PetscBool *flag) 81 { 82 PetscFunctionBegin; 83 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 84 PetscValidPointer(flag,2); 85 *flag = snes->errorifnotconverged; 86 PetscFunctionReturn(0); 87 } 88 89 #undef __FUNCT__ 90 #define __FUNCT__ "SNESSetFunctionDomainError" 91 /*@ 92 SNESSetFunctionDomainError - tells SNES that the input vector to your FormFunction is not 93 in the functions domain. For example, negative pressure. 94 95 Logically Collective on SNES 96 97 Input Parameters: 98 . SNES - the SNES context 99 100 Level: advanced 101 102 .keywords: SNES, view 103 104 .seealso: SNESCreate(), SNESSetFunction() 105 @*/ 106 PetscErrorCode SNESSetFunctionDomainError(SNES snes) 107 { 108 PetscFunctionBegin; 109 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 110 snes->domainerror = PETSC_TRUE; 111 PetscFunctionReturn(0); 112 } 113 114 #undef __FUNCT__ 115 #define __FUNCT__ "SNESView" 116 /*@C 117 SNESView - Prints the SNES data structure. 118 119 Collective on SNES 120 121 Input Parameters: 122 + SNES - the SNES context 123 - viewer - visualization context 124 125 Options Database Key: 126 . -snes_view - Calls SNESView() at end of SNESSolve() 127 128 Notes: 129 The available visualization contexts include 130 + PETSC_VIEWER_STDOUT_SELF - standard output (default) 131 - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 132 output where only the first processor opens 133 the file. All other processors send their 134 data to the first processor to print. 135 136 The user can open an alternative visualization context with 137 PetscViewerASCIIOpen() - output to a specified file. 138 139 Level: beginner 140 141 .keywords: SNES, view 142 143 .seealso: PetscViewerASCIIOpen() 144 @*/ 145 PetscErrorCode SNESView(SNES snes,PetscViewer viewer) 146 { 147 SNESKSPEW *kctx; 148 PetscErrorCode ierr; 149 KSP ksp; 150 PetscBool iascii,isstring; 151 152 PetscFunctionBegin; 153 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 154 if (!viewer) { 155 ierr = PetscViewerASCIIGetStdout(((PetscObject)snes)->comm,&viewer);CHKERRQ(ierr); 156 } 157 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 158 PetscCheckSameComm(snes,1,viewer,2); 159 160 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 161 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 162 if (iascii) { 163 ierr = PetscObjectPrintClassNamePrefixType((PetscObject)snes,viewer,"SNES Object");CHKERRQ(ierr); 164 if (snes->ops->view) { 165 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 166 ierr = (*snes->ops->view)(snes,viewer);CHKERRQ(ierr); 167 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 168 } 169 ierr = PetscViewerASCIIPrintf(viewer," maximum iterations=%D, maximum function evaluations=%D\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr); 170 ierr = PetscViewerASCIIPrintf(viewer," tolerances: relative=%G, absolute=%G, solution=%G\n", 171 snes->rtol,snes->abstol,snes->xtol);CHKERRQ(ierr); 172 ierr = PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",snes->linear_its);CHKERRQ(ierr); 173 ierr = PetscViewerASCIIPrintf(viewer," total number of function evaluations=%D\n",snes->nfuncs);CHKERRQ(ierr); 174 if (snes->ksp_ewconv) { 175 kctx = (SNESKSPEW *)snes->kspconvctx; 176 if (kctx) { 177 ierr = PetscViewerASCIIPrintf(viewer," Eisenstat-Walker computation of KSP relative tolerance (version %D)\n",kctx->version);CHKERRQ(ierr); 178 ierr = PetscViewerASCIIPrintf(viewer," rtol_0=%G, rtol_max=%G, threshold=%G\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr); 179 ierr = PetscViewerASCIIPrintf(viewer," gamma=%G, alpha=%G, alpha2=%G\n",kctx->gamma,kctx->alpha,kctx->alpha2);CHKERRQ(ierr); 180 } 181 } 182 if (snes->lagpreconditioner == -1) { 183 ierr = PetscViewerASCIIPrintf(viewer," Preconditioned is never rebuilt\n");CHKERRQ(ierr); 184 } else if (snes->lagpreconditioner > 1) { 185 ierr = PetscViewerASCIIPrintf(viewer," Preconditioned is rebuilt every %D new Jacobians\n",snes->lagpreconditioner);CHKERRQ(ierr); 186 } 187 if (snes->lagjacobian == -1) { 188 ierr = PetscViewerASCIIPrintf(viewer," Jacobian is never rebuilt\n");CHKERRQ(ierr); 189 } else if (snes->lagjacobian > 1) { 190 ierr = PetscViewerASCIIPrintf(viewer," Jacobian is rebuilt every %D SNES iterations\n",snes->lagjacobian);CHKERRQ(ierr); 191 } 192 } else if (isstring) { 193 const char *type; 194 ierr = SNESGetType(snes,&type);CHKERRQ(ierr); 195 ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr); 196 } 197 if (snes->pc && snes->usespc) { 198 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 199 ierr = SNESView(snes->pc, viewer);CHKERRQ(ierr); 200 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 201 } 202 if (snes->usesksp) { 203 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 204 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 205 ierr = KSPView(ksp,viewer);CHKERRQ(ierr); 206 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 207 } 208 PetscFunctionReturn(0); 209 } 210 211 /* 212 We retain a list of functions that also take SNES command 213 line options. These are called at the end SNESSetFromOptions() 214 */ 215 #define MAXSETFROMOPTIONS 5 216 static PetscInt numberofsetfromoptions; 217 static PetscErrorCode (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES); 218 219 #undef __FUNCT__ 220 #define __FUNCT__ "SNESAddOptionsChecker" 221 /*@C 222 SNESAddOptionsChecker - Adds an additional function to check for SNES options. 223 224 Not Collective 225 226 Input Parameter: 227 . snescheck - function that checks for options 228 229 Level: developer 230 231 .seealso: SNESSetFromOptions() 232 @*/ 233 PetscErrorCode SNESAddOptionsChecker(PetscErrorCode (*snescheck)(SNES)) 234 { 235 PetscFunctionBegin; 236 if (numberofsetfromoptions >= MAXSETFROMOPTIONS) { 237 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %D allowed", MAXSETFROMOPTIONS); 238 } 239 othersetfromoptions[numberofsetfromoptions++] = snescheck; 240 PetscFunctionReturn(0); 241 } 242 243 extern PetscErrorCode SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*); 244 245 #undef __FUNCT__ 246 #define __FUNCT__ "SNESSetUpMatrixFree_Private" 247 static PetscErrorCode SNESSetUpMatrixFree_Private(SNES snes, PetscBool hasOperator, PetscInt version) 248 { 249 Mat J; 250 KSP ksp; 251 PC pc; 252 PetscBool match; 253 PetscErrorCode ierr; 254 255 PetscFunctionBegin; 256 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 257 258 if(!snes->vec_func && (snes->jacobian || snes->jacobian_pre)) { 259 Mat A = snes->jacobian, B = snes->jacobian_pre; 260 ierr = MatGetVecs(A ? A : B, PETSC_NULL,&snes->vec_func);CHKERRQ(ierr); 261 } 262 263 if (version == 1) { 264 ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr); 265 ierr = MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);CHKERRQ(ierr); 266 ierr = MatSetFromOptions(J);CHKERRQ(ierr); 267 } else if (version == 2) { 268 if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"SNESSetFunction() must be called first"); 269 #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128) 270 ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_func,&J);CHKERRQ(ierr); 271 #else 272 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "matrix-free operator rutines (version 2)"); 273 #endif 274 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "matrix-free operator rutines, only version 1 and 2"); 275 276 ierr = PetscInfo1(snes,"Setting default matrix-free operator routines (version %D)\n", version);CHKERRQ(ierr); 277 if (hasOperator) { 278 /* This version replaces the user provided Jacobian matrix with a 279 matrix-free version but still employs the user-provided preconditioner matrix. */ 280 ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 281 } else { 282 /* This version replaces both the user-provided Jacobian and the user- 283 provided preconditioner matrix with the default matrix free version. */ 284 ierr = SNESSetJacobian(snes,J,J,MatMFFDComputeJacobian,snes->funP);CHKERRQ(ierr); 285 /* Force no preconditioner */ 286 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 287 ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 288 ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&match);CHKERRQ(ierr); 289 if (!match) { 290 ierr = PetscInfo(snes,"Setting default matrix-free preconditioner routines\nThat is no preconditioner is being used\n");CHKERRQ(ierr); 291 ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 292 } 293 } 294 ierr = MatDestroy(&J);CHKERRQ(ierr); 295 PetscFunctionReturn(0); 296 } 297 298 #undef __FUNCT__ 299 #define __FUNCT__ "SNESSetFromOptions" 300 /*@ 301 SNESSetFromOptions - Sets various SNES and KSP parameters from user options. 302 303 Collective on SNES 304 305 Input Parameter: 306 . snes - the SNES context 307 308 Options Database Keys: 309 + -snes_type <type> - ls, tr, ngmres, ncg, richardson, qn, vi, fas 310 . -snes_stol - convergence tolerance in terms of the norm 311 of the change in the solution between steps 312 . -snes_atol <abstol> - absolute tolerance of residual norm 313 . -snes_rtol <rtol> - relative decrease in tolerance norm from initial 314 . -snes_max_it <max_it> - maximum number of iterations 315 . -snes_max_funcs <max_funcs> - maximum number of function evaluations 316 . -snes_max_fail <max_fail> - maximum number of line search failures allowed before stopping, default is none 317 . -snes_max_linear_solve_fail - number of linear solver failures before SNESSolve() stops 318 . -snes_lag_preconditioner <lag> - how often preconditioner is rebuilt (use -1 to never rebuild) 319 . -snes_lag_jacobian <lag> - how often Jacobian is rebuilt (use -1 to never rebuild) 320 . -snes_trtol <trtol> - trust region tolerance 321 . -snes_no_convergence_test - skip convergence test in nonlinear 322 solver; hence iterations will continue until max_it 323 or some other criterion is reached. Saves expense 324 of convergence test 325 . -snes_monitor <optional filename> - prints residual norm at each iteration. if no 326 filename given prints to stdout 327 . -snes_monitor_solution - plots solution at each iteration 328 . -snes_monitor_residual - plots residual (not its norm) at each iteration 329 . -snes_monitor_solution_update - plots update to solution at each iteration 330 . -snes_monitor_draw - plots residual norm at each iteration 331 . -snes_fd - use finite differences to compute Jacobian; very slow, only for testing 332 . -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration 333 - -snes_converged_reason - print the reason for convergence/divergence after each solve 334 335 Options Database for Eisenstat-Walker method: 336 + -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence 337 . -snes_ksp_ew_version ver - version of Eisenstat-Walker method 338 . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0 339 . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax 340 . -snes_ksp_ew_gamma <gamma> - Sets gamma 341 . -snes_ksp_ew_alpha <alpha> - Sets alpha 342 . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2 343 - -snes_ksp_ew_threshold <threshold> - Sets threshold 344 345 Notes: 346 To see all options, run your program with the -help option or consult 347 the <A href="../../docs/manual.pdf#nameddest=ch_snes">SNES chapter of the users manual</A>. 348 349 Level: beginner 350 351 .keywords: SNES, nonlinear, set, options, database 352 353 .seealso: SNESSetOptionsPrefix() 354 @*/ 355 PetscErrorCode SNESSetFromOptions(SNES snes) 356 { 357 PetscBool flg,set,mf,mf_operator,pcset; 358 PetscInt i,indx,lag,mf_version,grids; 359 MatStructure matflag; 360 const char *deft = SNESLS; 361 const char *convtests[] = {"default","skip"}; 362 SNESKSPEW *kctx = NULL; 363 char type[256], monfilename[PETSC_MAX_PATH_LEN]; 364 const char *optionsprefix; 365 PetscViewer monviewer; 366 PetscErrorCode ierr; 367 368 PetscFunctionBegin; 369 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 370 371 if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 372 ierr = PetscObjectOptionsBegin((PetscObject)snes);CHKERRQ(ierr); 373 if (((PetscObject)snes)->type_name) { deft = ((PetscObject)snes)->type_name; } 374 ierr = PetscOptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr); 375 if (flg) { 376 ierr = SNESSetType(snes,type);CHKERRQ(ierr); 377 } else if (!((PetscObject)snes)->type_name) { 378 ierr = SNESSetType(snes,deft);CHKERRQ(ierr); 379 } 380 /* not used here, but called so will go into help messaage */ 381 ierr = PetscOptionsName("-snes_view","Print detailed information on solver used","SNESView",0);CHKERRQ(ierr); 382 383 ierr = PetscOptionsReal("-snes_stol","Stop if step length less than","SNESSetTolerances",snes->xtol,&snes->xtol,0);CHKERRQ(ierr); 384 ierr = PetscOptionsReal("-snes_atol","Stop if function norm less than","SNESSetTolerances",snes->abstol,&snes->abstol,0);CHKERRQ(ierr); 385 386 ierr = PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less than","SNESSetTolerances",snes->rtol,&snes->rtol,0);CHKERRQ(ierr); 387 ierr = PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,PETSC_NULL);CHKERRQ(ierr); 388 ierr = PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr); 389 ierr = PetscOptionsInt("-snes_max_fail","Maximum failures","SNESSetTolerances",snes->maxFailures,&snes->maxFailures,PETSC_NULL);CHKERRQ(ierr); 390 ierr = PetscOptionsInt("-snes_max_linear_solve_fail","Maximum failures in linear solves allowed","SNESSetMaxLinearSolveFailures",snes->maxLinearSolveFailures,&snes->maxLinearSolveFailures,PETSC_NULL);CHKERRQ(ierr); 391 ierr = PetscOptionsBool("-snes_error_if_not_converged","Generate error if solver does not converge","SNESSetErrorIfNotConverged",snes->errorifnotconverged,&snes->errorifnotconverged,PETSC_NULL);CHKERRQ(ierr); 392 393 ierr = PetscOptionsInt("-snes_lag_preconditioner","How often to rebuild preconditioner","SNESSetLagPreconditioner",snes->lagpreconditioner,&lag,&flg);CHKERRQ(ierr); 394 if (flg) { 395 ierr = SNESSetLagPreconditioner(snes,lag);CHKERRQ(ierr); 396 } 397 ierr = PetscOptionsInt("-snes_lag_jacobian","How often to rebuild Jacobian","SNESSetLagJacobian",snes->lagjacobian,&lag,&flg);CHKERRQ(ierr); 398 if (flg) { 399 ierr = SNESSetLagJacobian(snes,lag);CHKERRQ(ierr); 400 } 401 ierr = PetscOptionsInt("-snes_grid_sequence","Use grid sequencing to generate initial guess","SNESSetGridSequence",snes->gridsequence,&grids,&flg);CHKERRQ(ierr); 402 if (flg) { 403 ierr = SNESSetGridSequence(snes,grids);CHKERRQ(ierr); 404 } 405 406 ierr = PetscOptionsEList("-snes_convergence_test","Convergence test","SNESSetConvergenceTest",convtests,2,"default",&indx,&flg);CHKERRQ(ierr); 407 if (flg) { 408 switch (indx) { 409 case 0: ierr = SNESSetConvergenceTest(snes,SNESDefaultConverged,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); break; 410 case 1: ierr = SNESSetConvergenceTest(snes,SNESSkipConverged,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); break; 411 } 412 } 413 414 ierr = PetscOptionsBool("-snes_converged_reason","Print reason for converged or diverged","SNESSolve",snes->printreason,&snes->printreason,PETSC_NULL);CHKERRQ(ierr); 415 416 kctx = (SNESKSPEW *)snes->kspconvctx; 417 418 ierr = PetscOptionsBool("-snes_ksp_ew","Use Eisentat-Walker linear system convergence test","SNESKSPSetUseEW",snes->ksp_ewconv,&snes->ksp_ewconv,PETSC_NULL);CHKERRQ(ierr); 419 420 ierr = PetscOptionsInt("-snes_ksp_ew_version","Version 1, 2 or 3","SNESKSPSetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr); 421 ierr = PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNESKSPSetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr); 422 ierr = PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNESKSPSetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr); 423 ierr = PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNESKSPSetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr); 424 ierr = PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNESKSPSetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr); 425 ierr = PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNESKSPSetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr); 426 ierr = PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNESKSPSetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr); 427 428 flg = PETSC_FALSE; 429 ierr = PetscOptionsBool("-snes_monitor_cancel","Remove all monitors","SNESMonitorCancel",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 430 if (flg) {ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);} 431 432 ierr = PetscOptionsString("-snes_monitor","Monitor norm of function","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 433 if (flg) { 434 ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,monfilename,&monviewer);CHKERRQ(ierr); 435 ierr = SNESMonitorSet(snes,SNESMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 436 } 437 438 ierr = PetscOptionsString("-snes_monitor_range","Monitor range of elements of function","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 439 if (flg) { 440 ierr = SNESMonitorSet(snes,SNESMonitorRange,0,0);CHKERRQ(ierr); 441 } 442 443 ierr = PetscOptionsString("-snes_ratiomonitor","Monitor ratios of norms of function","SNESMonitorSetRatio","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 444 if (flg) { 445 ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,monfilename,&monviewer);CHKERRQ(ierr); 446 ierr = SNESMonitorSetRatio(snes,monviewer);CHKERRQ(ierr); 447 } 448 449 ierr = PetscOptionsString("-snes_monitor_short","Monitor norm of function (fewer digits)","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 450 if (flg) { 451 ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,monfilename,&monviewer);CHKERRQ(ierr); 452 ierr = SNESMonitorSet(snes,SNESMonitorDefaultShort,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 453 } 454 455 ierr = PetscOptionsString("-snes_monitor_python","Use Python function","SNESMonitorSet",0,monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 456 if (flg) {ierr = PetscPythonMonitorSet((PetscObject)snes,monfilename);CHKERRQ(ierr);} 457 458 flg = PETSC_FALSE; 459 ierr = PetscOptionsBool("-snes_monitor_solution","Plot solution at each iteration","SNESMonitorSolution",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 460 if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorSolution,0,0);CHKERRQ(ierr);} 461 flg = PETSC_FALSE; 462 ierr = PetscOptionsBool("-snes_monitor_solution_update","Plot correction at each iteration","SNESMonitorSolutionUpdate",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 463 if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorSolutionUpdate,0,0);CHKERRQ(ierr);} 464 flg = PETSC_FALSE; 465 ierr = PetscOptionsBool("-snes_monitor_residual","Plot residual at each iteration","SNESMonitorResidual",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 466 if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorResidual,0,0);CHKERRQ(ierr);} 467 flg = PETSC_FALSE; 468 ierr = PetscOptionsBool("-snes_monitor_draw","Plot function norm at each iteration","SNESMonitorLG",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 469 if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);} 470 flg = PETSC_FALSE; 471 ierr = PetscOptionsBool("-snes_monitor_range_draw","Plot function range at each iteration","SNESMonitorLG",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 472 if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorLGRange,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);} 473 474 flg = PETSC_FALSE; 475 ierr = PetscOptionsBool("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 476 if (flg) { 477 ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr); 478 ierr = PetscInfo(snes,"Setting default finite difference Jacobian matrix\n");CHKERRQ(ierr); 479 } 480 481 mf = mf_operator = PETSC_FALSE; 482 flg = PETSC_FALSE; 483 ierr = PetscOptionsBool("-snes_mf_operator","Use a Matrix-Free Jacobian with user-provided preconditioner matrix","MatCreateSNESMF",PETSC_FALSE,&mf_operator,&flg);CHKERRQ(ierr); 484 if (flg && mf_operator) { 485 snes->mf_operator = PETSC_TRUE; 486 mf = PETSC_TRUE; 487 } 488 flg = PETSC_FALSE; 489 ierr = PetscOptionsBool("-snes_mf","Use a Matrix-Free Jacobian with no preconditioner matrix","MatCreateSNESMF",PETSC_FALSE,&mf,&flg);CHKERRQ(ierr); 490 if (!flg && mf_operator) mf = PETSC_TRUE; 491 mf_version = 1; 492 ierr = PetscOptionsInt("-snes_mf_version","Matrix-Free routines version 1 or 2","None",mf_version,&mf_version,0);CHKERRQ(ierr); 493 494 495 /* GS Options */ 496 ierr = PetscOptionsInt("-snes_gs_sweeps","Number of sweeps of GS to apply","SNESComputeGS",snes->gssweeps,&snes->gssweeps,PETSC_NULL);CHKERRQ(ierr); 497 498 /* line search options */ 499 ierr = PetscOptionsReal("-snes_ls_alpha","Constant function norm must decrease by","None",snes->ls_alpha,&snes->ls_alpha,0);CHKERRQ(ierr); 500 ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",snes->maxstep,&snes->maxstep,0);CHKERRQ(ierr); 501 ierr = PetscOptionsReal("-snes_ls_steptol","Minimum lambda allowed","None",snes->steptol,&snes->steptol,0);CHKERRQ(ierr); 502 ierr = PetscOptionsReal("-snes_ls_damping","Damping parameter","SNES",snes->damping,&snes->damping,&flg);CHKERRQ(ierr); 503 ierr = PetscOptionsInt("-snes_ls_it" ,"Line search iterations","SNES",snes->ls_its,&snes->ls_its,&flg);CHKERRQ(ierr); 504 ierr = PetscOptionsBool("-snes_ls_monitor","Print progress of line searches","SNESLineSearchSetMonitor",snes->ls_monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 505 if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);} 506 ierr = PetscOptionsEnum("-snes_ls","Line search used","SNESLineSearchSet",SNESLineSearchTypes,(PetscEnum)snes->ls_type,(PetscEnum*)&indx,&flg);CHKERRQ(ierr); 507 if (flg) { 508 ierr = SNESLineSearchSetType(snes,(SNESLineSearchType)indx);CHKERRQ(ierr); 509 } 510 flg = snes->ops->precheckstep == SNESLineSearchPreCheckPicard ? PETSC_TRUE : PETSC_FALSE; 511 ierr = PetscOptionsBool("-snes_ls_precheck_picard","Use a correction that sometimes improves convergence of Picard iteration","SNESLineSearchPreCheckPicard",flg,&flg,&set);CHKERRQ(ierr); 512 if (set) { 513 if (flg) { 514 snes->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */ 515 ierr = PetscOptionsReal("-snes_ls_precheck_picard_angle","Maximum angle at which to activate the correction","none",snes->precheck_picard_angle,&snes->precheck_picard_angle,PETSC_NULL);CHKERRQ(ierr); 516 ierr = SNESLineSearchSetPreCheck(snes,SNESLineSearchPreCheckPicard,&snes->precheck_picard_angle);CHKERRQ(ierr); 517 } else { 518 ierr = SNESLineSearchSetPreCheck(snes,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 519 } 520 } 521 522 for(i = 0; i < numberofsetfromoptions; i++) { 523 ierr = (*othersetfromoptions[i])(snes);CHKERRQ(ierr); 524 } 525 526 if (snes->ops->setfromoptions) { 527 ierr = (*snes->ops->setfromoptions)(snes);CHKERRQ(ierr); 528 } 529 530 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 531 ierr = PetscObjectProcessOptionsHandlers((PetscObject)snes);CHKERRQ(ierr); 532 ierr = PetscOptionsEnd();CHKERRQ(ierr); 533 534 if (mf) { ierr = SNESSetUpMatrixFree_Private(snes, mf_operator, mf_version);CHKERRQ(ierr); } 535 536 if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);} 537 ierr = KSPGetOperators(snes->ksp,PETSC_NULL,PETSC_NULL,&matflag); 538 ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,matflag);CHKERRQ(ierr); 539 ierr = KSPSetFromOptions(snes->ksp);CHKERRQ(ierr); 540 541 /* if someone has set the SNES PC type, create it. */ 542 ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 543 ierr = PetscOptionsHasName(optionsprefix, "-npc_snes_type", &pcset);CHKERRQ(ierr); 544 if (pcset && (!snes->pc)) { 545 ierr = SNESGetPC(snes, &snes->pc);CHKERRQ(ierr); 546 } 547 if (snes->pc) { 548 ierr = SNESSetOptionsPrefix(snes->pc, optionsprefix);CHKERRQ(ierr); 549 ierr = SNESAppendOptionsPrefix(snes->pc, "npc_");CHKERRQ(ierr); 550 ierr = SNESSetDM(snes->pc, snes->dm);CHKERRQ(ierr); 551 ierr = SNESSetGS(snes->pc, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 552 /* Should we make a duplicate vector and matrix? Leave the DM to make it? */ 553 ierr = SNESSetFunction(snes->pc, snes->vec_func, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 554 ierr = SNESSetJacobian(snes->pc, snes->jacobian, snes->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 555 ierr = SNESSetFromOptions(snes->pc);CHKERRQ(ierr); 556 } 557 PetscFunctionReturn(0); 558 } 559 560 #undef __FUNCT__ 561 #define __FUNCT__ "SNESSetComputeApplicationContext" 562 /*@ 563 SNESSetComputeApplicationContext - Sets an optional function to compute a user-defined context for 564 the nonlinear solvers. 565 566 Logically Collective on SNES 567 568 Input Parameters: 569 + snes - the SNES context 570 . compute - function to compute the context 571 - destroy - function to destroy the context 572 573 Level: intermediate 574 575 .keywords: SNES, nonlinear, set, application, context 576 577 .seealso: SNESGetApplicationContext(), SNESSetComputeApplicationContext(), SNESGetApplicationContext() 578 @*/ 579 PetscErrorCode SNESSetComputeApplicationContext(SNES snes,PetscErrorCode (*compute)(SNES,void**),PetscErrorCode (*destroy)(void**)) 580 { 581 PetscFunctionBegin; 582 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 583 snes->ops->usercompute = compute; 584 snes->ops->userdestroy = destroy; 585 PetscFunctionReturn(0); 586 } 587 588 #undef __FUNCT__ 589 #define __FUNCT__ "SNESSetApplicationContext" 590 /*@ 591 SNESSetApplicationContext - Sets the optional user-defined context for 592 the nonlinear solvers. 593 594 Logically Collective on SNES 595 596 Input Parameters: 597 + snes - the SNES context 598 - usrP - optional user context 599 600 Level: intermediate 601 602 .keywords: SNES, nonlinear, set, application, context 603 604 .seealso: SNESGetApplicationContext(), SNESSetApplicationContext() 605 @*/ 606 PetscErrorCode SNESSetApplicationContext(SNES snes,void *usrP) 607 { 608 PetscErrorCode ierr; 609 KSP ksp; 610 611 PetscFunctionBegin; 612 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 613 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 614 ierr = KSPSetApplicationContext(ksp,usrP);CHKERRQ(ierr); 615 snes->user = usrP; 616 PetscFunctionReturn(0); 617 } 618 619 #undef __FUNCT__ 620 #define __FUNCT__ "SNESGetApplicationContext" 621 /*@ 622 SNESGetApplicationContext - Gets the user-defined context for the 623 nonlinear solvers. 624 625 Not Collective 626 627 Input Parameter: 628 . snes - SNES context 629 630 Output Parameter: 631 . usrP - user context 632 633 Level: intermediate 634 635 .keywords: SNES, nonlinear, get, application, context 636 637 .seealso: SNESSetApplicationContext() 638 @*/ 639 PetscErrorCode SNESGetApplicationContext(SNES snes,void *usrP) 640 { 641 PetscFunctionBegin; 642 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 643 *(void**)usrP = snes->user; 644 PetscFunctionReturn(0); 645 } 646 647 #undef __FUNCT__ 648 #define __FUNCT__ "SNESGetIterationNumber" 649 /*@ 650 SNESGetIterationNumber - Gets the number of nonlinear iterations completed 651 at this time. 652 653 Not Collective 654 655 Input Parameter: 656 . snes - SNES context 657 658 Output Parameter: 659 . iter - iteration number 660 661 Notes: 662 For example, during the computation of iteration 2 this would return 1. 663 664 This is useful for using lagged Jacobians (where one does not recompute the 665 Jacobian at each SNES iteration). For example, the code 666 .vb 667 ierr = SNESGetIterationNumber(snes,&it); 668 if (!(it % 2)) { 669 [compute Jacobian here] 670 } 671 .ve 672 can be used in your ComputeJacobian() function to cause the Jacobian to be 673 recomputed every second SNES iteration. 674 675 Level: intermediate 676 677 .keywords: SNES, nonlinear, get, iteration, number, 678 679 .seealso: SNESGetFunctionNorm(), SNESGetLinearSolveIterations() 680 @*/ 681 PetscErrorCode SNESGetIterationNumber(SNES snes,PetscInt* iter) 682 { 683 PetscFunctionBegin; 684 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 685 PetscValidIntPointer(iter,2); 686 *iter = snes->iter; 687 PetscFunctionReturn(0); 688 } 689 690 #undef __FUNCT__ 691 #define __FUNCT__ "SNESGetFunctionNorm" 692 /*@ 693 SNESGetFunctionNorm - Gets the norm of the current function that was set 694 with SNESSSetFunction(). 695 696 Collective on SNES 697 698 Input Parameter: 699 . snes - SNES context 700 701 Output Parameter: 702 . fnorm - 2-norm of function 703 704 Level: intermediate 705 706 .keywords: SNES, nonlinear, get, function, norm 707 708 .seealso: SNESGetFunction(), SNESGetIterationNumber(), SNESGetLinearSolveIterations() 709 @*/ 710 PetscErrorCode SNESGetFunctionNorm(SNES snes,PetscReal *fnorm) 711 { 712 PetscFunctionBegin; 713 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 714 PetscValidScalarPointer(fnorm,2); 715 *fnorm = snes->norm; 716 PetscFunctionReturn(0); 717 } 718 719 #undef __FUNCT__ 720 #define __FUNCT__ "SNESGetNonlinearStepFailures" 721 /*@ 722 SNESGetNonlinearStepFailures - Gets the number of unsuccessful steps 723 attempted by the nonlinear solver. 724 725 Not Collective 726 727 Input Parameter: 728 . snes - SNES context 729 730 Output Parameter: 731 . nfails - number of unsuccessful steps attempted 732 733 Notes: 734 This counter is reset to zero for each successive call to SNESSolve(). 735 736 Level: intermediate 737 738 .keywords: SNES, nonlinear, get, number, unsuccessful, steps 739 740 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), 741 SNESSetMaxNonlinearStepFailures(), SNESGetMaxNonlinearStepFailures() 742 @*/ 743 PetscErrorCode SNESGetNonlinearStepFailures(SNES snes,PetscInt* nfails) 744 { 745 PetscFunctionBegin; 746 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 747 PetscValidIntPointer(nfails,2); 748 *nfails = snes->numFailures; 749 PetscFunctionReturn(0); 750 } 751 752 #undef __FUNCT__ 753 #define __FUNCT__ "SNESSetMaxNonlinearStepFailures" 754 /*@ 755 SNESSetMaxNonlinearStepFailures - Sets the maximum number of unsuccessful steps 756 attempted by the nonlinear solver before it gives up. 757 758 Not Collective 759 760 Input Parameters: 761 + snes - SNES context 762 - maxFails - maximum of unsuccessful steps 763 764 Level: intermediate 765 766 .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps 767 768 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), 769 SNESGetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures() 770 @*/ 771 PetscErrorCode SNESSetMaxNonlinearStepFailures(SNES snes, PetscInt maxFails) 772 { 773 PetscFunctionBegin; 774 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 775 snes->maxFailures = maxFails; 776 PetscFunctionReturn(0); 777 } 778 779 #undef __FUNCT__ 780 #define __FUNCT__ "SNESGetMaxNonlinearStepFailures" 781 /*@ 782 SNESGetMaxNonlinearStepFailures - Gets the maximum number of unsuccessful steps 783 attempted by the nonlinear solver before it gives up. 784 785 Not Collective 786 787 Input Parameter: 788 . snes - SNES context 789 790 Output Parameter: 791 . maxFails - maximum of unsuccessful steps 792 793 Level: intermediate 794 795 .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps 796 797 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), 798 SNESSetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures() 799 800 @*/ 801 PetscErrorCode SNESGetMaxNonlinearStepFailures(SNES snes, PetscInt *maxFails) 802 { 803 PetscFunctionBegin; 804 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 805 PetscValidIntPointer(maxFails,2); 806 *maxFails = snes->maxFailures; 807 PetscFunctionReturn(0); 808 } 809 810 #undef __FUNCT__ 811 #define __FUNCT__ "SNESGetNumberFunctionEvals" 812 /*@ 813 SNESGetNumberFunctionEvals - Gets the number of user provided function evaluations 814 done by SNES. 815 816 Not Collective 817 818 Input Parameter: 819 . snes - SNES context 820 821 Output Parameter: 822 . nfuncs - number of evaluations 823 824 Level: intermediate 825 826 .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps 827 828 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures() 829 @*/ 830 PetscErrorCode SNESGetNumberFunctionEvals(SNES snes, PetscInt *nfuncs) 831 { 832 PetscFunctionBegin; 833 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 834 PetscValidIntPointer(nfuncs,2); 835 *nfuncs = snes->nfuncs; 836 PetscFunctionReturn(0); 837 } 838 839 #undef __FUNCT__ 840 #define __FUNCT__ "SNESGetLinearSolveFailures" 841 /*@ 842 SNESGetLinearSolveFailures - Gets the number of failed (non-converged) 843 linear solvers. 844 845 Not Collective 846 847 Input Parameter: 848 . snes - SNES context 849 850 Output Parameter: 851 . nfails - number of failed solves 852 853 Notes: 854 This counter is reset to zero for each successive call to SNESSolve(). 855 856 Level: intermediate 857 858 .keywords: SNES, nonlinear, get, number, unsuccessful, steps 859 860 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures() 861 @*/ 862 PetscErrorCode SNESGetLinearSolveFailures(SNES snes,PetscInt* nfails) 863 { 864 PetscFunctionBegin; 865 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 866 PetscValidIntPointer(nfails,2); 867 *nfails = snes->numLinearSolveFailures; 868 PetscFunctionReturn(0); 869 } 870 871 #undef __FUNCT__ 872 #define __FUNCT__ "SNESSetMaxLinearSolveFailures" 873 /*@ 874 SNESSetMaxLinearSolveFailures - the number of failed linear solve attempts 875 allowed before SNES returns with a diverged reason of SNES_DIVERGED_LINEAR_SOLVE 876 877 Logically Collective on SNES 878 879 Input Parameters: 880 + snes - SNES context 881 - maxFails - maximum allowed linear solve failures 882 883 Level: intermediate 884 885 Notes: By default this is 0; that is SNES returns on the first failed linear solve 886 887 .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps 888 889 .seealso: SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations() 890 @*/ 891 PetscErrorCode SNESSetMaxLinearSolveFailures(SNES snes, PetscInt maxFails) 892 { 893 PetscFunctionBegin; 894 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 895 PetscValidLogicalCollectiveInt(snes,maxFails,2); 896 snes->maxLinearSolveFailures = maxFails; 897 PetscFunctionReturn(0); 898 } 899 900 #undef __FUNCT__ 901 #define __FUNCT__ "SNESGetMaxLinearSolveFailures" 902 /*@ 903 SNESGetMaxLinearSolveFailures - gets the maximum number of linear solve failures that 904 are allowed before SNES terminates 905 906 Not Collective 907 908 Input Parameter: 909 . snes - SNES context 910 911 Output Parameter: 912 . maxFails - maximum of unsuccessful solves allowed 913 914 Level: intermediate 915 916 Notes: By default this is 1; that is SNES returns on the first failed linear solve 917 918 .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps 919 920 .seealso: SNESGetLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), 921 @*/ 922 PetscErrorCode SNESGetMaxLinearSolveFailures(SNES snes, PetscInt *maxFails) 923 { 924 PetscFunctionBegin; 925 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 926 PetscValidIntPointer(maxFails,2); 927 *maxFails = snes->maxLinearSolveFailures; 928 PetscFunctionReturn(0); 929 } 930 931 #undef __FUNCT__ 932 #define __FUNCT__ "SNESGetLinearSolveIterations" 933 /*@ 934 SNESGetLinearSolveIterations - Gets the total number of linear iterations 935 used by the nonlinear solver. 936 937 Not Collective 938 939 Input Parameter: 940 . snes - SNES context 941 942 Output Parameter: 943 . lits - number of linear iterations 944 945 Notes: 946 This counter is reset to zero for each successive call to SNESSolve(). 947 948 Level: intermediate 949 950 .keywords: SNES, nonlinear, get, number, linear, iterations 951 952 .seealso: SNESGetIterationNumber(), SNESGetFunctionNorm(), SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures() 953 @*/ 954 PetscErrorCode SNESGetLinearSolveIterations(SNES snes,PetscInt* lits) 955 { 956 PetscFunctionBegin; 957 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 958 PetscValidIntPointer(lits,2); 959 *lits = snes->linear_its; 960 PetscFunctionReturn(0); 961 } 962 963 #undef __FUNCT__ 964 #define __FUNCT__ "SNESGetKSP" 965 /*@ 966 SNESGetKSP - Returns the KSP context for a SNES solver. 967 968 Not Collective, but if SNES object is parallel, then KSP object is parallel 969 970 Input Parameter: 971 . snes - the SNES context 972 973 Output Parameter: 974 . ksp - the KSP context 975 976 Notes: 977 The user can then directly manipulate the KSP context to set various 978 options, etc. Likewise, the user can then extract and manipulate the 979 PC contexts as well. 980 981 Level: beginner 982 983 .keywords: SNES, nonlinear, get, KSP, context 984 985 .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP() 986 @*/ 987 PetscErrorCode SNESGetKSP(SNES snes,KSP *ksp) 988 { 989 PetscErrorCode ierr; 990 991 PetscFunctionBegin; 992 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 993 PetscValidPointer(ksp,2); 994 995 if (!snes->ksp) { 996 ierr = KSPCreate(((PetscObject)snes)->comm,&snes->ksp);CHKERRQ(ierr); 997 ierr = PetscObjectIncrementTabLevel((PetscObject)snes->ksp,(PetscObject)snes,1);CHKERRQ(ierr); 998 ierr = PetscLogObjectParent(snes,snes->ksp);CHKERRQ(ierr); 999 } 1000 *ksp = snes->ksp; 1001 PetscFunctionReturn(0); 1002 } 1003 1004 #undef __FUNCT__ 1005 #define __FUNCT__ "SNESSetKSP" 1006 /*@ 1007 SNESSetKSP - Sets a KSP context for the SNES object to use 1008 1009 Not Collective, but the SNES and KSP objects must live on the same MPI_Comm 1010 1011 Input Parameters: 1012 + snes - the SNES context 1013 - ksp - the KSP context 1014 1015 Notes: 1016 The SNES object already has its KSP object, you can obtain with SNESGetKSP() 1017 so this routine is rarely needed. 1018 1019 The KSP object that is already in the SNES object has its reference count 1020 decreased by one. 1021 1022 Level: developer 1023 1024 .keywords: SNES, nonlinear, get, KSP, context 1025 1026 .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP() 1027 @*/ 1028 PetscErrorCode SNESSetKSP(SNES snes,KSP ksp) 1029 { 1030 PetscErrorCode ierr; 1031 1032 PetscFunctionBegin; 1033 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1034 PetscValidHeaderSpecific(ksp,KSP_CLASSID,2); 1035 PetscCheckSameComm(snes,1,ksp,2); 1036 ierr = PetscObjectReference((PetscObject)ksp);CHKERRQ(ierr); 1037 if (snes->ksp) {ierr = PetscObjectDereference((PetscObject)snes->ksp);CHKERRQ(ierr);} 1038 snes->ksp = ksp; 1039 PetscFunctionReturn(0); 1040 } 1041 1042 #if 0 1043 #undef __FUNCT__ 1044 #define __FUNCT__ "SNESPublish_Petsc" 1045 static PetscErrorCode SNESPublish_Petsc(PetscObject obj) 1046 { 1047 PetscFunctionBegin; 1048 PetscFunctionReturn(0); 1049 } 1050 #endif 1051 1052 /* -----------------------------------------------------------*/ 1053 #undef __FUNCT__ 1054 #define __FUNCT__ "SNESCreate" 1055 /*@ 1056 SNESCreate - Creates a nonlinear solver context. 1057 1058 Collective on MPI_Comm 1059 1060 Input Parameters: 1061 . comm - MPI communicator 1062 1063 Output Parameter: 1064 . outsnes - the new SNES context 1065 1066 Options Database Keys: 1067 + -snes_mf - Activates default matrix-free Jacobian-vector products, 1068 and no preconditioning matrix 1069 . -snes_mf_operator - Activates default matrix-free Jacobian-vector 1070 products, and a user-provided preconditioning matrix 1071 as set by SNESSetJacobian() 1072 - -snes_fd - Uses (slow!) finite differences to compute Jacobian 1073 1074 Level: beginner 1075 1076 .keywords: SNES, nonlinear, create, context 1077 1078 .seealso: SNESSolve(), SNESDestroy(), SNES, SNESSetLagPreconditioner() 1079 1080 @*/ 1081 PetscErrorCode SNESCreate(MPI_Comm comm,SNES *outsnes) 1082 { 1083 PetscErrorCode ierr; 1084 SNES snes; 1085 SNESKSPEW *kctx; 1086 1087 PetscFunctionBegin; 1088 PetscValidPointer(outsnes,2); 1089 *outsnes = PETSC_NULL; 1090 #ifndef PETSC_USE_DYNAMIC_LIBRARIES 1091 ierr = SNESInitializePackage(PETSC_NULL);CHKERRQ(ierr); 1092 #endif 1093 1094 ierr = PetscHeaderCreate(snes,_p_SNES,struct _SNESOps,SNES_CLASSID,0,"SNES","Nonlinear solver","SNES",comm,SNESDestroy,SNESView);CHKERRQ(ierr); 1095 1096 snes->ops->converged = SNESDefaultConverged; 1097 snes->usesksp = PETSC_TRUE; 1098 snes->max_its = 50; 1099 snes->max_funcs = 10000; 1100 snes->norm = 0.0; 1101 snes->rtol = 1.e-8; 1102 snes->ttol = 0.0; 1103 snes->abstol = 1.e-50; 1104 snes->xtol = 1.e-8; 1105 snes->deltatol = 1.e-12; 1106 snes->nfuncs = 0; 1107 snes->numFailures = 0; 1108 snes->maxFailures = 1; 1109 snes->linear_its = 0; 1110 snes->lagjacobian = 1; 1111 snes->lagpreconditioner = 1; 1112 snes->numbermonitors = 0; 1113 snes->data = 0; 1114 snes->setupcalled = PETSC_FALSE; 1115 snes->ksp_ewconv = PETSC_FALSE; 1116 snes->nwork = 0; 1117 snes->work = 0; 1118 snes->nvwork = 0; 1119 snes->vwork = 0; 1120 snes->conv_hist_len = 0; 1121 snes->conv_hist_max = 0; 1122 snes->conv_hist = PETSC_NULL; 1123 snes->conv_hist_its = PETSC_NULL; 1124 snes->conv_hist_reset = PETSC_TRUE; 1125 snes->reason = SNES_CONVERGED_ITERATING; 1126 snes->gssweeps = 1; 1127 1128 /* initialize the line search options */ 1129 snes->ls_type = SNES_LS_BASIC; 1130 snes->ls_its = 1; 1131 snes->damping = 1.0; 1132 snes->maxstep = 1e8; 1133 snes->steptol = 1e-12; 1134 snes->ls_alpha = 1e-4; 1135 snes->ls_monitor = PETSC_NULL; 1136 1137 snes->ops->linesearch = PETSC_NULL; 1138 snes->precheck = PETSC_NULL; 1139 snes->ops->precheckstep = PETSC_NULL; 1140 snes->postcheck = PETSC_NULL; 1141 snes->ops->postcheckstep= PETSC_NULL; 1142 1143 snes->numLinearSolveFailures = 0; 1144 snes->maxLinearSolveFailures = 1; 1145 1146 /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 1147 ierr = PetscNewLog(snes,SNESKSPEW,&kctx);CHKERRQ(ierr); 1148 snes->kspconvctx = (void*)kctx; 1149 kctx->version = 2; 1150 kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 1151 this was too large for some test cases */ 1152 kctx->rtol_last = 0.0; 1153 kctx->rtol_max = .9; 1154 kctx->gamma = 1.0; 1155 kctx->alpha = .5*(1.0 + PetscSqrtReal(5.0)); 1156 kctx->alpha2 = kctx->alpha; 1157 kctx->threshold = .1; 1158 kctx->lresid_last = 0.0; 1159 kctx->norm_last = 0.0; 1160 1161 *outsnes = snes; 1162 PetscFunctionReturn(0); 1163 } 1164 1165 #undef __FUNCT__ 1166 #define __FUNCT__ "SNESSetFunction" 1167 /*@C 1168 SNESSetFunction - Sets the function evaluation routine and function 1169 vector for use by the SNES routines in solving systems of nonlinear 1170 equations. 1171 1172 Logically Collective on SNES 1173 1174 Input Parameters: 1175 + snes - the SNES context 1176 . r - vector to store function value 1177 . func - function evaluation routine 1178 - ctx - [optional] user-defined context for private data for the 1179 function evaluation routine (may be PETSC_NULL) 1180 1181 Calling sequence of func: 1182 $ func (SNES snes,Vec x,Vec f,void *ctx); 1183 1184 . f - function vector 1185 - ctx - optional user-defined function context 1186 1187 Notes: 1188 The Newton-like methods typically solve linear systems of the form 1189 $ f'(x) x = -f(x), 1190 where f'(x) denotes the Jacobian matrix and f(x) is the function. 1191 1192 Level: beginner 1193 1194 .keywords: SNES, nonlinear, set, function 1195 1196 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetPicard() 1197 @*/ 1198 PetscErrorCode SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx) 1199 { 1200 PetscErrorCode ierr; 1201 PetscFunctionBegin; 1202 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1203 if (r) { 1204 PetscValidHeaderSpecific(r,VEC_CLASSID,2); 1205 PetscCheckSameComm(snes,1,r,2); 1206 ierr = PetscObjectReference((PetscObject)r);CHKERRQ(ierr); 1207 ierr = VecDestroy(&snes->vec_func);CHKERRQ(ierr); 1208 snes->vec_func = r; 1209 } else if (!snes->vec_func && snes->dm) { 1210 ierr = DMCreateGlobalVector(snes->dm,&snes->vec_func);CHKERRQ(ierr); 1211 } 1212 if (func) snes->ops->computefunction = func; 1213 if (ctx) snes->funP = ctx; 1214 PetscFunctionReturn(0); 1215 } 1216 1217 1218 #undef __FUNCT__ 1219 #define __FUNCT__ "SNESSetGS" 1220 /*@C 1221 SNESSetGS - Sets the user nonlinear Gauss-Seidel routine for 1222 use with composed nonlinear solvers. 1223 1224 Input Parameters: 1225 + snes - the SNES context 1226 . gsfunc - function evaluation routine 1227 - ctx - [optional] user-defined context for private data for the 1228 smoother evaluation routine (may be PETSC_NULL) 1229 1230 Calling sequence of func: 1231 $ func (SNES snes,Vec x,Vec b,void *ctx); 1232 1233 + X - solution vector 1234 . B - RHS vector 1235 - ctx - optional user-defined Gauss-Seidel context 1236 1237 Notes: 1238 The GS routines are used by the composed nonlinear solver to generate 1239 a problem appropriate update to the solution, particularly FAS. 1240 1241 Level: intermediate 1242 1243 .keywords: SNES, nonlinear, set, Gauss-Seidel 1244 1245 .seealso: SNESGetFunction(), SNESComputeGS() 1246 @*/ 1247 PetscErrorCode SNESSetGS(SNES snes, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx) { 1248 PetscFunctionBegin; 1249 if (gsfunc) snes->ops->computegs = gsfunc; 1250 if (ctx) snes->gsP = ctx; 1251 PetscFunctionReturn(0); 1252 } 1253 1254 #undef __FUNCT__ 1255 #define __FUNCT__ "SNESSetGSSweeps" 1256 /*@ 1257 SNESSetGSSweeps - Sets the number of sweeps of GS to use. 1258 1259 Input Parameters: 1260 + snes - the SNES context 1261 - sweeps - the number of sweeps of GS to perform. 1262 1263 Level: intermediate 1264 1265 .keywords: SNES, nonlinear, set, Gauss-Siedel 1266 1267 .seealso: SNESSetGS(), SNESGetGS(), SNESSetPC(), SNESGetGSSweeps() 1268 @*/ 1269 1270 PetscErrorCode SNESSetGSSweeps(SNES snes, PetscInt sweeps) { 1271 PetscFunctionBegin; 1272 snes->gssweeps = sweeps; 1273 PetscFunctionReturn(0); 1274 } 1275 1276 1277 #undef __FUNCT__ 1278 #define __FUNCT__ "SNESGetGSSweeps" 1279 /*@ 1280 SNESGetGSSweeps - Gets the number of sweeps GS will use. 1281 1282 Input Parameters: 1283 . snes - the SNES context 1284 1285 Output Parameters: 1286 . sweeps - the number of sweeps of GS to perform. 1287 1288 Level: intermediate 1289 1290 .keywords: SNES, nonlinear, set, Gauss-Siedel 1291 1292 .seealso: SNESSetGS(), SNESGetGS(), SNESSetPC(), SNESSetGSSweeps() 1293 @*/ 1294 PetscErrorCode SNESGetGSSweeps(SNES snes, PetscInt * sweeps) { 1295 PetscFunctionBegin; 1296 *sweeps = snes->gssweeps; 1297 PetscFunctionReturn(0); 1298 } 1299 1300 #undef __FUNCT__ 1301 #define __FUNCT__ "SNESPicardComputeFunction" 1302 PetscErrorCode SNESPicardComputeFunction(SNES snes,Vec x,Vec f,void *ctx) 1303 { 1304 PetscErrorCode ierr; 1305 PetscFunctionBegin; 1306 /* A(x)*x - b(x) */ 1307 ierr = (*snes->ops->computepjacobian)(snes,x,&snes->jacobian,&snes->jacobian_pre,&snes->matstruct,snes->jacP);CHKERRQ(ierr); 1308 ierr = (*snes->ops->computepfunction)(snes,x,f,snes->funP);CHKERRQ(ierr); 1309 ierr = VecView(x,PETSC_VIEWER_BINARY_WORLD);CHKERRQ(ierr); 1310 ierr = VecView(f,PETSC_VIEWER_BINARY_WORLD);CHKERRQ(ierr); 1311 ierr = VecScale(f,-1.0);CHKERRQ(ierr); 1312 ierr = MatMultAdd(snes->jacobian_pre,x,f,f);CHKERRQ(ierr); 1313 PetscFunctionReturn(0); 1314 } 1315 1316 #undef __FUNCT__ 1317 #define __FUNCT__ "SNESPicardComputeJacobian" 1318 PetscErrorCode SNESPicardComputeJacobian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx) 1319 { 1320 PetscFunctionBegin; 1321 *flag = snes->matstruct; 1322 PetscFunctionReturn(0); 1323 } 1324 1325 #undef __FUNCT__ 1326 #define __FUNCT__ "SNESSetPicard" 1327 /*@C 1328 SNESSetPicard - Use SNES to solve the semilinear-system A(x) x = b(x) via a Picard type iteration (Picard linearization) 1329 1330 Logically Collective on SNES 1331 1332 Input Parameters: 1333 + snes - the SNES context 1334 . r - vector to store function value 1335 . func - function evaluation routine 1336 . jmat - normally the same as mat but you can pass another matrix for which you compute the Jacobian of A(x) x - b(x) (see jmat below) 1337 . mat - matrix to store A 1338 . mfunc - function to compute matrix value 1339 - ctx - [optional] user-defined context for private data for the 1340 function evaluation routine (may be PETSC_NULL) 1341 1342 Calling sequence of func: 1343 $ func (SNES snes,Vec x,Vec f,void *ctx); 1344 1345 + f - function vector 1346 - ctx - optional user-defined function context 1347 1348 Calling sequence of mfunc: 1349 $ mfunc (SNES snes,Vec x,Mat *jmat,Mat *mat,int *flag,void *ctx); 1350 1351 + x - input vector 1352 . jmat - Form Jacobian matrix of A(x) x - b(x) if available, not there is really no reason to use it in this way since then you can just use SNESSetJacobian(), 1353 normally just pass mat in this location 1354 . mat - form A(x) matrix 1355 . flag - flag indicating information about the preconditioner matrix 1356 structure (same as flag in KSPSetOperators()), one of SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER 1357 - ctx - [optional] user-defined Jacobian context 1358 1359 Notes: 1360 One can call SNESSetPicard() or SNESSetFunction() (and possibly SNESSetJacobian()) but cannot call both 1361 1362 $ Solves the equation A(x) x = b(x) via the defect correction algorithm A(x^{n}) (x^{n+1} - x^{n}) = b(x^{n}) - A(x^{n})x^{n} 1363 $ Note that when an exact solver is used this corresponds to the "classic" Picard A(x^{n}) x^{n+1} = b(x^{n}) iteration. 1364 1365 Run with -snes_mf_operator to solve the system with Newton's method using A(x^{n}) to construct the preconditioner. 1366 1367 We implement the defect correction form of the Picard iteration because it converges much more generally when inexact linear solvers are used then 1368 the direct Picard iteration A(x^n) x^{n+1} = b(x^n) 1369 1370 There is some controversity over the definition of a Picard iteration for nonlinear systems but almost everyone agrees that it involves a linear solve and some 1371 believe it is the iteration A(x^{n}) x^{n+1} = b(x^{n}) hence we use the name Picard. If anyone has an authoritative reference that defines the Picard iteration 1372 different please contact us at petsc-dev@mcs.anl.gov and we'll have an entirely new argument :-). 1373 1374 Level: beginner 1375 1376 .keywords: SNES, nonlinear, set, function 1377 1378 .seealso: SNESGetFunction(), SNESSetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESGetPicard(), SNESLineSearchPreCheckPicard() 1379 @*/ 1380 PetscErrorCode SNESSetPicard(SNES snes,Vec r,PetscErrorCode (*func)(SNES,Vec,Vec,void*),Mat jmat, Mat mat, PetscErrorCode (*mfunc)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 1381 { 1382 PetscErrorCode ierr; 1383 PetscFunctionBegin; 1384 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1385 snes->ops->computepfunction = func; 1386 snes->ops->computepjacobian = mfunc; 1387 ierr = SNESSetFunction(snes,r,SNESPicardComputeFunction,ctx);CHKERRQ(ierr); 1388 ierr = SNESSetJacobian(snes,jmat,mat,SNESPicardComputeJacobian,ctx);CHKERRQ(ierr); 1389 PetscFunctionReturn(0); 1390 } 1391 1392 #undef __FUNCT__ 1393 #define __FUNCT__ "SNESSetComputeInitialGuess" 1394 /*@C 1395 SNESSetComputeInitialGuess - Sets a routine used to compute an initial guess for the problem 1396 1397 Logically Collective on SNES 1398 1399 Input Parameters: 1400 + snes - the SNES context 1401 . func - function evaluation routine 1402 - ctx - [optional] user-defined context for private data for the 1403 function evaluation routine (may be PETSC_NULL) 1404 1405 Calling sequence of func: 1406 $ func (SNES snes,Vec x,void *ctx); 1407 1408 . f - function vector 1409 - ctx - optional user-defined function context 1410 1411 Level: intermediate 1412 1413 .keywords: SNES, nonlinear, set, function 1414 1415 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 1416 @*/ 1417 PetscErrorCode SNESSetComputeInitialGuess(SNES snes,PetscErrorCode (*func)(SNES,Vec,void*),void *ctx) 1418 { 1419 PetscFunctionBegin; 1420 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1421 if (func) snes->ops->computeinitialguess = func; 1422 if (ctx) snes->initialguessP = ctx; 1423 PetscFunctionReturn(0); 1424 } 1425 1426 /* --------------------------------------------------------------- */ 1427 #undef __FUNCT__ 1428 #define __FUNCT__ "SNESGetRhs" 1429 /*@C 1430 SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set 1431 it assumes a zero right hand side. 1432 1433 Logically Collective on SNES 1434 1435 Input Parameter: 1436 . snes - the SNES context 1437 1438 Output Parameter: 1439 . rhs - the right hand side vector or PETSC_NULL if the right hand side vector is null 1440 1441 Level: intermediate 1442 1443 .keywords: SNES, nonlinear, get, function, right hand side 1444 1445 .seealso: SNESGetSolution(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction() 1446 @*/ 1447 PetscErrorCode SNESGetRhs(SNES snes,Vec *rhs) 1448 { 1449 PetscFunctionBegin; 1450 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1451 PetscValidPointer(rhs,2); 1452 *rhs = snes->vec_rhs; 1453 PetscFunctionReturn(0); 1454 } 1455 1456 #undef __FUNCT__ 1457 #define __FUNCT__ "SNESComputeFunction" 1458 /*@ 1459 SNESComputeFunction - Calls the function that has been set with 1460 SNESSetFunction(). 1461 1462 Collective on SNES 1463 1464 Input Parameters: 1465 + snes - the SNES context 1466 - x - input vector 1467 1468 Output Parameter: 1469 . y - function vector, as set by SNESSetFunction() 1470 1471 Notes: 1472 SNESComputeFunction() is typically used within nonlinear solvers 1473 implementations, so most users would not generally call this routine 1474 themselves. 1475 1476 Level: developer 1477 1478 .keywords: SNES, nonlinear, compute, function 1479 1480 .seealso: SNESSetFunction(), SNESGetFunction() 1481 @*/ 1482 PetscErrorCode SNESComputeFunction(SNES snes,Vec x,Vec y) 1483 { 1484 PetscErrorCode ierr; 1485 1486 PetscFunctionBegin; 1487 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1488 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 1489 PetscValidHeaderSpecific(y,VEC_CLASSID,3); 1490 PetscCheckSameComm(snes,1,x,2); 1491 PetscCheckSameComm(snes,1,y,3); 1492 VecValidValues(x,2,PETSC_TRUE); 1493 1494 ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 1495 if (snes->ops->computefunction) { 1496 PetscStackPush("SNES user function"); 1497 ierr = (*snes->ops->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr); 1498 PetscStackPop; 1499 } else if (snes->dm) { 1500 ierr = DMComputeFunction(snes->dm,x,y);CHKERRQ(ierr); 1501 } else if (snes->vec_rhs) { 1502 ierr = MatMult(snes->jacobian, x, y);CHKERRQ(ierr); 1503 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() or SNESSetDM() before SNESComputeFunction(), likely called from SNESSolve()."); 1504 if (snes->vec_rhs) { 1505 ierr = VecAXPY(y,-1.0,snes->vec_rhs);CHKERRQ(ierr); 1506 } 1507 snes->nfuncs++; 1508 ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 1509 VecValidValues(y,3,PETSC_FALSE); 1510 PetscFunctionReturn(0); 1511 } 1512 1513 #undef __FUNCT__ 1514 #define __FUNCT__ "SNESComputeGS" 1515 /*@ 1516 SNESComputeGS - Calls the Gauss-Seidel function that has been set with 1517 SNESSetGS(). 1518 1519 Collective on SNES 1520 1521 Input Parameters: 1522 + snes - the SNES context 1523 . x - input vector 1524 - b - rhs vector 1525 1526 Output Parameter: 1527 . x - new solution vector 1528 1529 Notes: 1530 SNESComputeGS() is typically used within composed nonlinear solver 1531 implementations, so most users would not generally call this routine 1532 themselves. 1533 1534 Level: developer 1535 1536 .keywords: SNES, nonlinear, compute, function 1537 1538 .seealso: SNESSetGS(), SNESComputeFunction() 1539 @*/ 1540 PetscErrorCode SNESComputeGS(SNES snes,Vec b,Vec x) 1541 { 1542 PetscErrorCode ierr; 1543 PetscInt i; 1544 1545 PetscFunctionBegin; 1546 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1547 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 1548 if (b) PetscValidHeaderSpecific(b,VEC_CLASSID,3); 1549 PetscCheckSameComm(snes,1,x,2); 1550 if(b) PetscCheckSameComm(snes,1,b,3); 1551 if (b) VecValidValues(b,2,PETSC_TRUE); 1552 ierr = PetscLogEventBegin(SNES_GSEval,snes,x,b,0);CHKERRQ(ierr); 1553 if (snes->ops->computegs) { 1554 for (i = 0; i < snes->gssweeps; i++) { 1555 PetscStackPush("SNES user GS"); 1556 ierr = (*snes->ops->computegs)(snes,x,b,snes->gsP);CHKERRQ(ierr); 1557 PetscStackPop; 1558 } 1559 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetGS() before SNESComputeGS(), likely called from SNESSolve()."); 1560 ierr = PetscLogEventEnd(SNES_GSEval,snes,x,b,0);CHKERRQ(ierr); 1561 VecValidValues(x,3,PETSC_FALSE); 1562 PetscFunctionReturn(0); 1563 } 1564 1565 1566 #undef __FUNCT__ 1567 #define __FUNCT__ "SNESComputeJacobian" 1568 /*@ 1569 SNESComputeJacobian - Computes the Jacobian matrix that has been 1570 set with SNESSetJacobian(). 1571 1572 Collective on SNES and Mat 1573 1574 Input Parameters: 1575 + snes - the SNES context 1576 - x - input vector 1577 1578 Output Parameters: 1579 + A - Jacobian matrix 1580 . B - optional preconditioning matrix 1581 - flag - flag indicating matrix structure (one of, SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER) 1582 1583 Options Database Keys: 1584 + -snes_lag_preconditioner <lag> 1585 . -snes_lag_jacobian <lag> 1586 . -snes_compare_explicit - Compare the computed Jacobian to the finite difference Jacobian and output the differences 1587 . -snes_compare_explicit_draw - Compare the computed Jacobian to the finite difference Jacobian and draw the result 1588 . -snes_compare_explicit_contour - Compare the computed Jacobian to the finite difference Jacobian and draw a contour plot with the result 1589 . -snes_compare_operator - Make the comparison options above use the operator instead of the preconditioning matrix 1590 . -snes_compare_coloring - Compute the finite differece Jacobian using coloring and display norms of difference 1591 . -snes_compare_coloring_display - Compute the finite differece Jacobian using coloring and display verbose differences 1592 . -snes_compare_coloring_threshold - Display only those matrix entries that differ by more than a given threshold 1593 . -snes_compare_coloring_threshold_atol - Absolute tolerance for difference in matrix entries to be displayed by -snes_compare_coloring_threshold 1594 . -snes_compare_coloring_threshold_rtol - Relative tolerance for difference in matrix entries to be displayed by -snes_compare_coloring_threshold 1595 . -snes_compare_coloring_draw - Compute the finite differece Jacobian using coloring and draw differences 1596 - -snes_compare_coloring_draw_contour - Compute the finite differece Jacobian using coloring and show contours of matrices and differences 1597 1598 1599 Notes: 1600 Most users should not need to explicitly call this routine, as it 1601 is used internally within the nonlinear solvers. 1602 1603 See KSPSetOperators() for important information about setting the 1604 flag parameter. 1605 1606 Level: developer 1607 1608 .keywords: SNES, compute, Jacobian, matrix 1609 1610 .seealso: SNESSetJacobian(), KSPSetOperators(), MatStructure, SNESSetLagPreconditioner(), SNESSetLagJacobian() 1611 @*/ 1612 PetscErrorCode SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 1613 { 1614 PetscErrorCode ierr; 1615 PetscBool flag; 1616 1617 PetscFunctionBegin; 1618 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1619 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 1620 PetscValidPointer(flg,5); 1621 PetscCheckSameComm(snes,1,X,2); 1622 VecValidValues(X,2,PETSC_TRUE); 1623 if (!snes->ops->computejacobian) PetscFunctionReturn(0); 1624 1625 /* make sure that MatAssemblyBegin/End() is called on A matrix if it is matrix free */ 1626 1627 if (snes->lagjacobian == -2) { 1628 snes->lagjacobian = -1; 1629 ierr = PetscInfo(snes,"Recomputing Jacobian/preconditioner because lag is -2 (means compute Jacobian, but then never again) \n");CHKERRQ(ierr); 1630 } else if (snes->lagjacobian == -1) { 1631 *flg = SAME_PRECONDITIONER; 1632 ierr = PetscInfo(snes,"Reusing Jacobian/preconditioner because lag is -1\n");CHKERRQ(ierr); 1633 ierr = PetscTypeCompare((PetscObject)*A,MATMFFD,&flag);CHKERRQ(ierr); 1634 if (flag) { 1635 ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1636 ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1637 } 1638 PetscFunctionReturn(0); 1639 } else if (snes->lagjacobian > 1 && snes->iter % snes->lagjacobian) { 1640 *flg = SAME_PRECONDITIONER; 1641 ierr = PetscInfo2(snes,"Reusing Jacobian/preconditioner because lag is %D and SNES iteration is %D\n",snes->lagjacobian,snes->iter);CHKERRQ(ierr); 1642 ierr = PetscTypeCompare((PetscObject)*A,MATMFFD,&flag);CHKERRQ(ierr); 1643 if (flag) { 1644 ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1645 ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1646 } 1647 PetscFunctionReturn(0); 1648 } 1649 1650 *flg = DIFFERENT_NONZERO_PATTERN; 1651 ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 1652 PetscStackPush("SNES user Jacobian function"); 1653 ierr = (*snes->ops->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr); 1654 PetscStackPop; 1655 ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 1656 1657 if (snes->lagpreconditioner == -2) { 1658 ierr = PetscInfo(snes,"Rebuilding preconditioner exactly once since lag is -2\n");CHKERRQ(ierr); 1659 snes->lagpreconditioner = -1; 1660 } else if (snes->lagpreconditioner == -1) { 1661 *flg = SAME_PRECONDITIONER; 1662 ierr = PetscInfo(snes,"Reusing preconditioner because lag is -1\n");CHKERRQ(ierr); 1663 } else if (snes->lagpreconditioner > 1 && snes->iter % snes->lagpreconditioner) { 1664 *flg = SAME_PRECONDITIONER; 1665 ierr = PetscInfo2(snes,"Reusing preconditioner because lag is %D and SNES iteration is %D\n",snes->lagpreconditioner,snes->iter);CHKERRQ(ierr); 1666 } 1667 1668 /* make sure user returned a correct Jacobian and preconditioner */ 1669 /* PetscValidHeaderSpecific(*A,MAT_CLASSID,3); 1670 PetscValidHeaderSpecific(*B,MAT_CLASSID,4); */ 1671 { 1672 PetscBool flag = PETSC_FALSE,flag_draw = PETSC_FALSE,flag_contour = PETSC_FALSE,flag_operator = PETSC_FALSE; 1673 ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_explicit",&flag,PETSC_NULL);CHKERRQ(ierr); 1674 ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_explicit_draw",&flag_draw,PETSC_NULL);CHKERRQ(ierr); 1675 ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_explicit_draw_contour",&flag_contour,PETSC_NULL);CHKERRQ(ierr); 1676 ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_operator",&flag_operator,PETSC_NULL);CHKERRQ(ierr); 1677 if (flag || flag_draw || flag_contour) { 1678 Mat Bexp_mine = PETSC_NULL,Bexp,FDexp; 1679 MatStructure mstruct; 1680 PetscViewer vdraw,vstdout; 1681 PetscBool flg; 1682 if (flag_operator) { 1683 ierr = MatComputeExplicitOperator(*A,&Bexp_mine);CHKERRQ(ierr); 1684 Bexp = Bexp_mine; 1685 } else { 1686 /* See if the preconditioning matrix can be viewed and added directly */ 1687 ierr = PetscTypeCompareAny((PetscObject)*B,&flg,MATSEQAIJ,MATMPIAIJ,MATSEQDENSE,MATMPIDENSE,MATSEQBAIJ,MATMPIBAIJ,MATSEQSBAIJ,MATMPIBAIJ,"");CHKERRQ(ierr); 1688 if (flg) Bexp = *B; 1689 else { 1690 /* If the "preconditioning" matrix is itself MATSHELL or some other type without direct support */ 1691 ierr = MatComputeExplicitOperator(*B,&Bexp_mine);CHKERRQ(ierr); 1692 Bexp = Bexp_mine; 1693 } 1694 } 1695 ierr = MatConvert(Bexp,MATSAME,MAT_INITIAL_MATRIX,&FDexp);CHKERRQ(ierr); 1696 ierr = SNESDefaultComputeJacobian(snes,X,&FDexp,&FDexp,&mstruct,NULL);CHKERRQ(ierr); 1697 ierr = PetscViewerASCIIGetStdout(((PetscObject)snes)->comm,&vstdout);CHKERRQ(ierr); 1698 if (flag_draw || flag_contour) { 1699 ierr = PetscViewerDrawOpen(((PetscObject)snes)->comm,0,"Explicit Jacobians",PETSC_DECIDE,PETSC_DECIDE,300,300,&vdraw);CHKERRQ(ierr); 1700 if (flag_contour) {ierr = PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);CHKERRQ(ierr);} 1701 } else vdraw = PETSC_NULL; 1702 ierr = PetscViewerASCIIPrintf(vstdout,"Explicit %s\n",flag_operator?"Jacobian":"preconditioning Jacobian");CHKERRQ(ierr); 1703 if (flag) {ierr = MatView(Bexp,vstdout);CHKERRQ(ierr);} 1704 if (vdraw) {ierr = MatView(Bexp,vdraw);CHKERRQ(ierr);} 1705 ierr = PetscViewerASCIIPrintf(vstdout,"Finite difference Jacobian\n");CHKERRQ(ierr); 1706 if (flag) {ierr = MatView(FDexp,vstdout);CHKERRQ(ierr);} 1707 if (vdraw) {ierr = MatView(FDexp,vdraw);CHKERRQ(ierr);} 1708 ierr = MatAYPX(FDexp,-1.0,Bexp,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 1709 ierr = PetscViewerASCIIPrintf(vstdout,"User-provided matrix minus finite difference Jacobian\n");CHKERRQ(ierr); 1710 if (flag) {ierr = MatView(FDexp,vstdout);CHKERRQ(ierr);} 1711 if (vdraw) { /* Always use contour for the difference */ 1712 ierr = PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);CHKERRQ(ierr); 1713 ierr = MatView(FDexp,vdraw);CHKERRQ(ierr); 1714 ierr = PetscViewerPopFormat(vdraw);CHKERRQ(ierr); 1715 } 1716 if (flag_contour) {ierr = PetscViewerPopFormat(vdraw);CHKERRQ(ierr);} 1717 ierr = PetscViewerDestroy(&vdraw);CHKERRQ(ierr); 1718 ierr = MatDestroy(&Bexp_mine);CHKERRQ(ierr); 1719 ierr = MatDestroy(&FDexp);CHKERRQ(ierr); 1720 } 1721 } 1722 { 1723 PetscBool flag = PETSC_FALSE,flag_display = PETSC_FALSE,flag_draw = PETSC_FALSE,flag_contour = PETSC_FALSE,flag_threshold = PETSC_FALSE; 1724 PetscReal threshold_atol = PETSC_SQRT_MACHINE_EPSILON,threshold_rtol = 10*PETSC_SQRT_MACHINE_EPSILON; 1725 ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_coloring",&flag,PETSC_NULL);CHKERRQ(ierr); 1726 ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_coloring_display",&flag_display,PETSC_NULL);CHKERRQ(ierr); 1727 ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_coloring_draw",&flag_draw,PETSC_NULL);CHKERRQ(ierr); 1728 ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_coloring_draw_contour",&flag_contour,PETSC_NULL);CHKERRQ(ierr); 1729 ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold",&flag_threshold,PETSC_NULL);CHKERRQ(ierr); 1730 ierr = PetscOptionsGetReal(((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold_rtol",&threshold_rtol,PETSC_NULL);CHKERRQ(ierr); 1731 ierr = PetscOptionsGetReal(((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold_atol",&threshold_atol,PETSC_NULL);CHKERRQ(ierr); 1732 if (flag || flag_display || flag_draw || flag_contour || flag_threshold) { 1733 Mat Bfd; 1734 MatStructure mstruct; 1735 PetscViewer vdraw,vstdout; 1736 ISColoring iscoloring; 1737 MatFDColoring matfdcoloring; 1738 PetscErrorCode (*func)(SNES,Vec,Vec,void*); 1739 void *funcctx; 1740 PetscReal norm1,norm2,normmax; 1741 1742 ierr = MatDuplicate(*B,MAT_DO_NOT_COPY_VALUES,&Bfd);CHKERRQ(ierr); 1743 ierr = MatGetColoring(Bfd,MATCOLORINGSL,&iscoloring);CHKERRQ(ierr); 1744 ierr = MatFDColoringCreate(Bfd,iscoloring,&matfdcoloring);CHKERRQ(ierr); 1745 ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 1746 1747 /* This method of getting the function is currently unreliable since it doesn't work for DM local functions. */ 1748 ierr = SNESGetFunction(snes,PETSC_NULL,&func,&funcctx);CHKERRQ(ierr); 1749 ierr = MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode(*)(void))func,funcctx);CHKERRQ(ierr); 1750 ierr = PetscObjectSetOptionsPrefix((PetscObject)matfdcoloring,((PetscObject)snes)->prefix);CHKERRQ(ierr); 1751 ierr = PetscObjectAppendOptionsPrefix((PetscObject)matfdcoloring,"coloring_");CHKERRQ(ierr); 1752 ierr = MatFDColoringSetFromOptions(matfdcoloring);CHKERRQ(ierr); 1753 ierr = MatFDColoringApply(Bfd,matfdcoloring,X,&mstruct,snes);CHKERRQ(ierr); 1754 ierr = MatFDColoringDestroy(&matfdcoloring);CHKERRQ(ierr); 1755 1756 ierr = PetscViewerASCIIGetStdout(((PetscObject)snes)->comm,&vstdout);CHKERRQ(ierr); 1757 if (flag_draw || flag_contour) { 1758 ierr = PetscViewerDrawOpen(((PetscObject)snes)->comm,0,"Colored Jacobians",PETSC_DECIDE,PETSC_DECIDE,300,300,&vdraw);CHKERRQ(ierr); 1759 if (flag_contour) {ierr = PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);CHKERRQ(ierr);} 1760 } else vdraw = PETSC_NULL; 1761 ierr = PetscViewerASCIIPrintf(vstdout,"Explicit preconditioning Jacobian\n");CHKERRQ(ierr); 1762 if (flag_display) {ierr = MatView(*B,vstdout);CHKERRQ(ierr);} 1763 if (vdraw) {ierr = MatView(*B,vdraw);CHKERRQ(ierr);} 1764 ierr = PetscViewerASCIIPrintf(vstdout,"Colored Finite difference Jacobian\n");CHKERRQ(ierr); 1765 if (flag_display) {ierr = MatView(Bfd,vstdout);CHKERRQ(ierr);} 1766 if (vdraw) {ierr = MatView(Bfd,vdraw);CHKERRQ(ierr);} 1767 ierr = MatAYPX(Bfd,-1.0,*B,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 1768 ierr = MatNorm(Bfd,NORM_1,&norm1);CHKERRQ(ierr); 1769 ierr = MatNorm(Bfd,NORM_FROBENIUS,&norm2);CHKERRQ(ierr); 1770 ierr = MatNorm(Bfd,NORM_MAX,&normmax);CHKERRQ(ierr); 1771 ierr = PetscViewerASCIIPrintf(vstdout,"User-provided matrix minus finite difference Jacobian, norm1=%G normFrob=%G normmax=%G\n",norm1,norm2,normmax);CHKERRQ(ierr); 1772 if (flag_display) {ierr = MatView(Bfd,vstdout);CHKERRQ(ierr);} 1773 if (vdraw) { /* Always use contour for the difference */ 1774 ierr = PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);CHKERRQ(ierr); 1775 ierr = MatView(Bfd,vdraw);CHKERRQ(ierr); 1776 ierr = PetscViewerPopFormat(vdraw);CHKERRQ(ierr); 1777 } 1778 if (flag_contour) {ierr = PetscViewerPopFormat(vdraw);CHKERRQ(ierr);} 1779 1780 if (flag_threshold) { 1781 PetscInt bs,rstart,rend,i; 1782 ierr = MatGetBlockSize(*B,&bs);CHKERRQ(ierr); 1783 ierr = MatGetOwnershipRange(*B,&rstart,&rend);CHKERRQ(ierr); 1784 for (i=rstart; i<rend; i++) { 1785 const PetscScalar *ba,*ca; 1786 const PetscInt *bj,*cj; 1787 PetscInt bn,cn,j,maxentrycol = -1,maxdiffcol = -1,maxrdiffcol = -1; 1788 PetscReal maxentry = 0,maxdiff = 0,maxrdiff = 0; 1789 ierr = MatGetRow(*B,i,&bn,&bj,&ba);CHKERRQ(ierr); 1790 ierr = MatGetRow(Bfd,i,&cn,&cj,&ca);CHKERRQ(ierr); 1791 if (bn != cn) SETERRQ(((PetscObject)*A)->comm,PETSC_ERR_PLIB,"Unexpected different nonzero pattern in -snes_compare_coloring_threshold"); 1792 for (j=0; j<bn; j++) { 1793 PetscReal rdiff = PetscAbsScalar(ca[j]) / (threshold_atol + threshold_rtol*PetscAbsScalar(ba[j])); 1794 if (PetscAbsScalar(ba[j]) > PetscAbs(maxentry)) { 1795 maxentrycol = bj[j]; 1796 maxentry = PetscRealPart(ba[j]); 1797 } 1798 if (PetscAbsScalar(ca[j]) > PetscAbs(maxdiff)) { 1799 maxdiffcol = bj[j]; 1800 maxdiff = PetscRealPart(ca[j]); 1801 } 1802 if (rdiff > maxrdiff) { 1803 maxrdiffcol = bj[j]; 1804 maxrdiff = rdiff; 1805 } 1806 } 1807 if (maxrdiff > 1) { 1808 ierr = PetscViewerASCIIPrintf(vstdout,"row %D (maxentry=%G at %D, maxdiff=%G at %D, maxrdiff=%G at %D):",i,maxentry,maxentrycol,maxdiff,maxdiffcol,maxrdiff,maxrdiffcol);CHKERRQ(ierr); 1809 for (j=0; j<bn; j++) { 1810 PetscReal rdiff; 1811 rdiff = PetscAbsScalar(ca[j]) / (threshold_atol + threshold_rtol*PetscAbsScalar(ba[j])); 1812 if (rdiff > 1) { 1813 ierr = PetscViewerASCIIPrintf(vstdout," (%D,%G:%G)",bj[j],PetscRealPart(ba[j]),PetscRealPart(ca[j]));CHKERRQ(ierr); 1814 } 1815 } 1816 ierr = PetscViewerASCIIPrintf(vstdout,"\n",i,maxentry,maxdiff,maxrdiff);CHKERRQ(ierr); 1817 } 1818 ierr = MatRestoreRow(*B,i,&bn,&bj,&ba);CHKERRQ(ierr); 1819 ierr = MatRestoreRow(Bfd,i,&cn,&cj,&ca);CHKERRQ(ierr); 1820 } 1821 } 1822 ierr = PetscViewerDestroy(&vdraw);CHKERRQ(ierr); 1823 ierr = MatDestroy(&Bfd);CHKERRQ(ierr); 1824 } 1825 } 1826 PetscFunctionReturn(0); 1827 } 1828 1829 #undef __FUNCT__ 1830 #define __FUNCT__ "SNESSetJacobian" 1831 /*@C 1832 SNESSetJacobian - Sets the function to compute Jacobian as well as the 1833 location to store the matrix. 1834 1835 Logically Collective on SNES and Mat 1836 1837 Input Parameters: 1838 + snes - the SNES context 1839 . A - Jacobian matrix 1840 . B - preconditioner matrix (usually same as the Jacobian) 1841 . func - Jacobian evaluation routine (if PETSC_NULL then SNES retains any previously set value) 1842 - ctx - [optional] user-defined context for private data for the 1843 Jacobian evaluation routine (may be PETSC_NULL) (if PETSC_NULL then SNES retains any previously set value) 1844 1845 Calling sequence of func: 1846 $ func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 1847 1848 + x - input vector 1849 . A - Jacobian matrix 1850 . B - preconditioner matrix, usually the same as A 1851 . flag - flag indicating information about the preconditioner matrix 1852 structure (same as flag in KSPSetOperators()), one of SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER 1853 - ctx - [optional] user-defined Jacobian context 1854 1855 Notes: 1856 See KSPSetOperators() for important information about setting the flag 1857 output parameter in the routine func(). Be sure to read this information! 1858 1859 The routine func() takes Mat * as the matrix arguments rather than Mat. 1860 This allows the Jacobian evaluation routine to replace A and/or B with a 1861 completely new new matrix structure (not just different matrix elements) 1862 when appropriate, for instance, if the nonzero structure is changing 1863 throughout the global iterations. 1864 1865 If the A matrix and B matrix are different you must call MatAssemblyBegin/End() on 1866 each matrix. 1867 1868 If using SNESDefaultComputeJacobianColor() to assemble a Jacobian, the ctx argument 1869 must be a MatFDColoring. 1870 1871 Other defect-correction schemes can be used by computing a different matrix in place of the Jacobian. One common 1872 example is to use the "Picard linearization" which only differentiates through the highest order parts of each term. 1873 1874 Level: beginner 1875 1876 .keywords: SNES, nonlinear, set, Jacobian, matrix 1877 1878 .seealso: KSPSetOperators(), SNESSetFunction(), MatMFFDComputeJacobian(), SNESDefaultComputeJacobianColor(), MatStructure 1879 @*/ 1880 PetscErrorCode SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 1881 { 1882 PetscErrorCode ierr; 1883 1884 PetscFunctionBegin; 1885 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1886 if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2); 1887 if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3); 1888 if (A) PetscCheckSameComm(snes,1,A,2); 1889 if (B) PetscCheckSameComm(snes,1,B,3); 1890 if (func) snes->ops->computejacobian = func; 1891 if (ctx) snes->jacP = ctx; 1892 if (A) { 1893 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 1894 ierr = MatDestroy(&snes->jacobian);CHKERRQ(ierr); 1895 snes->jacobian = A; 1896 } 1897 if (B) { 1898 ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 1899 ierr = MatDestroy(&snes->jacobian_pre);CHKERRQ(ierr); 1900 snes->jacobian_pre = B; 1901 } 1902 PetscFunctionReturn(0); 1903 } 1904 1905 #undef __FUNCT__ 1906 #define __FUNCT__ "SNESGetJacobian" 1907 /*@C 1908 SNESGetJacobian - Returns the Jacobian matrix and optionally the user 1909 provided context for evaluating the Jacobian. 1910 1911 Not Collective, but Mat object will be parallel if SNES object is 1912 1913 Input Parameter: 1914 . snes - the nonlinear solver context 1915 1916 Output Parameters: 1917 + A - location to stash Jacobian matrix (or PETSC_NULL) 1918 . B - location to stash preconditioner matrix (or PETSC_NULL) 1919 . func - location to put Jacobian function (or PETSC_NULL) 1920 - ctx - location to stash Jacobian ctx (or PETSC_NULL) 1921 1922 Level: advanced 1923 1924 .seealso: SNESSetJacobian(), SNESComputeJacobian() 1925 @*/ 1926 PetscErrorCode SNESGetJacobian(SNES snes,Mat *A,Mat *B,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx) 1927 { 1928 PetscFunctionBegin; 1929 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1930 if (A) *A = snes->jacobian; 1931 if (B) *B = snes->jacobian_pre; 1932 if (func) *func = snes->ops->computejacobian; 1933 if (ctx) *ctx = snes->jacP; 1934 PetscFunctionReturn(0); 1935 } 1936 1937 /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 1938 1939 #undef __FUNCT__ 1940 #define __FUNCT__ "SNESSetUp" 1941 /*@ 1942 SNESSetUp - Sets up the internal data structures for the later use 1943 of a nonlinear solver. 1944 1945 Collective on SNES 1946 1947 Input Parameters: 1948 . snes - the SNES context 1949 1950 Notes: 1951 For basic use of the SNES solvers the user need not explicitly call 1952 SNESSetUp(), since these actions will automatically occur during 1953 the call to SNESSolve(). However, if one wishes to control this 1954 phase separately, SNESSetUp() should be called after SNESCreate() 1955 and optional routines of the form SNESSetXXX(), but before SNESSolve(). 1956 1957 Level: advanced 1958 1959 .keywords: SNES, nonlinear, setup 1960 1961 .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 1962 @*/ 1963 PetscErrorCode SNESSetUp(SNES snes) 1964 { 1965 PetscErrorCode ierr; 1966 1967 PetscFunctionBegin; 1968 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1969 if (snes->setupcalled) PetscFunctionReturn(0); 1970 1971 if (!((PetscObject)snes)->type_name) { 1972 ierr = SNESSetType(snes,SNESLS);CHKERRQ(ierr); 1973 } 1974 1975 ierr = SNESGetFunction(snes,&snes->vec_func,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 1976 if (snes->vec_func == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector"); 1977 if (snes->vec_rhs == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be right hand side vector"); 1978 1979 if (!snes->vec_sol_update /* && snes->vec_sol */) { 1980 ierr = VecDuplicate(snes->vec_sol,&snes->vec_sol_update);CHKERRQ(ierr); 1981 ierr = PetscLogObjectParent(snes,snes->vec_sol_update);CHKERRQ(ierr); 1982 } 1983 1984 if (!snes->ops->computejacobian && snes->dm) { 1985 Mat J,B; 1986 ierr = DMCreateMatrix(snes->dm,MATAIJ,&B);CHKERRQ(ierr); 1987 if (snes->mf_operator) { 1988 ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr); 1989 ierr = MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);CHKERRQ(ierr); 1990 ierr = MatSetFromOptions(J);CHKERRQ(ierr); 1991 } else { 1992 J = B; 1993 ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr); 1994 } 1995 ierr = SNESSetJacobian(snes,J,B,SNESDMComputeJacobian,PETSC_NULL);CHKERRQ(ierr); 1996 ierr = MatDestroy(&J);CHKERRQ(ierr); 1997 ierr = MatDestroy(&B);CHKERRQ(ierr); 1998 } else if (!snes->jacobian && snes->ops->computejacobian == MatMFFDComputeJacobian) { 1999 Mat J; 2000 ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr); 2001 ierr = MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);CHKERRQ(ierr); 2002 ierr = MatSetFromOptions(J);CHKERRQ(ierr); 2003 ierr = SNESSetJacobian(snes,J,J,MatMFFDComputeJacobian,snes->funP);CHKERRQ(ierr); 2004 ierr = MatDestroy(&J);CHKERRQ(ierr); 2005 } else if (snes->dm && snes->mf_operator && !snes->jacobian_pre && !snes->jacobian) { 2006 Mat J,B; 2007 ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr); 2008 ierr = MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);CHKERRQ(ierr); 2009 ierr = MatSetFromOptions(J);CHKERRQ(ierr); 2010 ierr = DMCreateMatrix(snes->dm,MATAIJ,&B);CHKERRQ(ierr); 2011 ierr = SNESSetJacobian(snes,J,B,SNESDMComputeJacobian,snes->funP);CHKERRQ(ierr); 2012 ierr = MatDestroy(&J);CHKERRQ(ierr); 2013 ierr = MatDestroy(&B);CHKERRQ(ierr); 2014 } else if (snes->dm && !snes->jacobian_pre){ 2015 Mat J; 2016 ierr = DMCreateMatrix(snes->dm,MATAIJ,&J);CHKERRQ(ierr); 2017 ierr = SNESSetJacobian(snes,J,J,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 2018 ierr = MatDestroy(&J);CHKERRQ(ierr); 2019 } 2020 if (!snes->ops->computefunction && !snes->dm) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must provide a residual function with SNESSetFunction() or call SNESSetDM()"); 2021 if (!snes->vec_func) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must provide a vector when calling SNESSetFunction() or call SNESSetDM()"); 2022 2023 if (!snes->ksp) {ierr = SNESGetKSP(snes, &snes->ksp);CHKERRQ(ierr);} 2024 2025 if (snes->ops->usercompute && !snes->user) { 2026 ierr = (*snes->ops->usercompute)(snes,(void**)&snes->user);CHKERRQ(ierr); 2027 } 2028 2029 if (snes->ops->setup) { 2030 ierr = (*snes->ops->setup)(snes);CHKERRQ(ierr); 2031 } 2032 2033 snes->setupcalled = PETSC_TRUE; 2034 PetscFunctionReturn(0); 2035 } 2036 2037 #undef __FUNCT__ 2038 #define __FUNCT__ "SNESReset" 2039 /*@ 2040 SNESReset - Resets a SNES context to the snessetupcalled = 0 state and removes any allocated Vecs and Mats 2041 2042 Collective on SNES 2043 2044 Input Parameter: 2045 . snes - iterative context obtained from SNESCreate() 2046 2047 Level: intermediate 2048 2049 Notes: Also calls the application context destroy routine set with SNESSetComputeApplicationContext() 2050 2051 .keywords: SNES, destroy 2052 2053 .seealso: SNESCreate(), SNESSetUp(), SNESSolve() 2054 @*/ 2055 PetscErrorCode SNESReset(SNES snes) 2056 { 2057 PetscErrorCode ierr; 2058 2059 PetscFunctionBegin; 2060 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2061 if (snes->ops->userdestroy && snes->user) { 2062 ierr = (*snes->ops->userdestroy)((void**)&snes->user);CHKERRQ(ierr); 2063 snes->user = PETSC_NULL; 2064 } 2065 if (snes->pc) { 2066 ierr = SNESReset(snes->pc);CHKERRQ(ierr); 2067 } 2068 2069 if (snes->ops->reset) { 2070 ierr = (*snes->ops->reset)(snes);CHKERRQ(ierr); 2071 } 2072 if (snes->ksp) {ierr = KSPReset(snes->ksp);CHKERRQ(ierr);} 2073 ierr = VecDestroy(&snes->vec_rhs);CHKERRQ(ierr); 2074 ierr = VecDestroy(&snes->vec_sol);CHKERRQ(ierr); 2075 ierr = VecDestroy(&snes->vec_sol_update);CHKERRQ(ierr); 2076 ierr = VecDestroy(&snes->vec_func);CHKERRQ(ierr); 2077 ierr = MatDestroy(&snes->jacobian);CHKERRQ(ierr); 2078 ierr = MatDestroy(&snes->jacobian_pre);CHKERRQ(ierr); 2079 if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);} 2080 if (snes->vwork) {ierr = VecDestroyVecs(snes->nvwork,&snes->vwork);CHKERRQ(ierr);} 2081 snes->nwork = snes->nvwork = 0; 2082 snes->setupcalled = PETSC_FALSE; 2083 PetscFunctionReturn(0); 2084 } 2085 2086 #undef __FUNCT__ 2087 #define __FUNCT__ "SNESDestroy" 2088 /*@ 2089 SNESDestroy - Destroys the nonlinear solver context that was created 2090 with SNESCreate(). 2091 2092 Collective on SNES 2093 2094 Input Parameter: 2095 . snes - the SNES context 2096 2097 Level: beginner 2098 2099 .keywords: SNES, nonlinear, destroy 2100 2101 .seealso: SNESCreate(), SNESSolve() 2102 @*/ 2103 PetscErrorCode SNESDestroy(SNES *snes) 2104 { 2105 PetscErrorCode ierr; 2106 2107 PetscFunctionBegin; 2108 if (!*snes) PetscFunctionReturn(0); 2109 PetscValidHeaderSpecific((*snes),SNES_CLASSID,1); 2110 if (--((PetscObject)(*snes))->refct > 0) {*snes = 0; PetscFunctionReturn(0);} 2111 2112 ierr = SNESReset((*snes));CHKERRQ(ierr); 2113 ierr = SNESDestroy(&(*snes)->pc);CHKERRQ(ierr); 2114 2115 /* if memory was published with AMS then destroy it */ 2116 ierr = PetscObjectDepublish((*snes));CHKERRQ(ierr); 2117 if ((*snes)->ops->destroy) {ierr = (*((*snes))->ops->destroy)((*snes));CHKERRQ(ierr);} 2118 2119 ierr = DMDestroy(&(*snes)->dm);CHKERRQ(ierr); 2120 ierr = KSPDestroy(&(*snes)->ksp);CHKERRQ(ierr); 2121 2122 ierr = PetscFree((*snes)->kspconvctx);CHKERRQ(ierr); 2123 if ((*snes)->ops->convergeddestroy) { 2124 ierr = (*(*snes)->ops->convergeddestroy)((*snes)->cnvP);CHKERRQ(ierr); 2125 } 2126 if ((*snes)->conv_malloc) { 2127 ierr = PetscFree((*snes)->conv_hist);CHKERRQ(ierr); 2128 ierr = PetscFree((*snes)->conv_hist_its);CHKERRQ(ierr); 2129 } 2130 ierr = PetscViewerDestroy(&(*snes)->ls_monitor);CHKERRQ(ierr); 2131 ierr = SNESMonitorCancel((*snes));CHKERRQ(ierr); 2132 ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr); 2133 PetscFunctionReturn(0); 2134 } 2135 2136 /* ----------- Routines to set solver parameters ---------- */ 2137 2138 #undef __FUNCT__ 2139 #define __FUNCT__ "SNESSetLagPreconditioner" 2140 /*@ 2141 SNESSetLagPreconditioner - Determines when the preconditioner is rebuilt in the nonlinear solve. 2142 2143 Logically Collective on SNES 2144 2145 Input Parameters: 2146 + snes - the SNES context 2147 - lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time 2148 the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that 2149 2150 Options Database Keys: 2151 . -snes_lag_preconditioner <lag> 2152 2153 Notes: 2154 The default is 1 2155 The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 2156 If -1 is used before the very first nonlinear solve the preconditioner is still built because there is no previous preconditioner to use 2157 2158 Level: intermediate 2159 2160 .keywords: SNES, nonlinear, set, convergence, tolerances 2161 2162 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian() 2163 2164 @*/ 2165 PetscErrorCode SNESSetLagPreconditioner(SNES snes,PetscInt lag) 2166 { 2167 PetscFunctionBegin; 2168 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2169 if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater"); 2170 if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0"); 2171 PetscValidLogicalCollectiveInt(snes,lag,2); 2172 snes->lagpreconditioner = lag; 2173 PetscFunctionReturn(0); 2174 } 2175 2176 #undef __FUNCT__ 2177 #define __FUNCT__ "SNESSetGridSequence" 2178 /*@ 2179 SNESSetGridSequence - sets the number of steps of grid sequencing that SNES does 2180 2181 Logically Collective on SNES 2182 2183 Input Parameters: 2184 + snes - the SNES context 2185 - steps - the number of refinements to do, defaults to 0 2186 2187 Options Database Keys: 2188 . -snes_grid_sequence <steps> 2189 2190 Level: intermediate 2191 2192 Notes: 2193 Use SNESGetSolution() to extract the fine grid solution after grid sequencing. 2194 2195 .keywords: SNES, nonlinear, set, convergence, tolerances 2196 2197 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian() 2198 2199 @*/ 2200 PetscErrorCode SNESSetGridSequence(SNES snes,PetscInt steps) 2201 { 2202 PetscFunctionBegin; 2203 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2204 PetscValidLogicalCollectiveInt(snes,steps,2); 2205 snes->gridsequence = steps; 2206 PetscFunctionReturn(0); 2207 } 2208 2209 #undef __FUNCT__ 2210 #define __FUNCT__ "SNESGetLagPreconditioner" 2211 /*@ 2212 SNESGetLagPreconditioner - Indicates how often the preconditioner is rebuilt 2213 2214 Not Collective 2215 2216 Input Parameter: 2217 . snes - the SNES context 2218 2219 Output Parameter: 2220 . lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time 2221 the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that 2222 2223 Options Database Keys: 2224 . -snes_lag_preconditioner <lag> 2225 2226 Notes: 2227 The default is 1 2228 The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 2229 2230 Level: intermediate 2231 2232 .keywords: SNES, nonlinear, set, convergence, tolerances 2233 2234 .seealso: SNESSetTrustRegionTolerance(), SNESSetLagPreconditioner() 2235 2236 @*/ 2237 PetscErrorCode SNESGetLagPreconditioner(SNES snes,PetscInt *lag) 2238 { 2239 PetscFunctionBegin; 2240 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2241 *lag = snes->lagpreconditioner; 2242 PetscFunctionReturn(0); 2243 } 2244 2245 #undef __FUNCT__ 2246 #define __FUNCT__ "SNESSetLagJacobian" 2247 /*@ 2248 SNESSetLagJacobian - Determines when the Jacobian is rebuilt in the nonlinear solve. See SNESSetLagPreconditioner() for determining how 2249 often the preconditioner is rebuilt. 2250 2251 Logically Collective on SNES 2252 2253 Input Parameters: 2254 + snes - the SNES context 2255 - lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time 2256 the Jacobian is built etc. -2 means rebuild at next chance but then never again 2257 2258 Options Database Keys: 2259 . -snes_lag_jacobian <lag> 2260 2261 Notes: 2262 The default is 1 2263 The Jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 2264 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 2265 at the next Newton step but never again (unless it is reset to another value) 2266 2267 Level: intermediate 2268 2269 .keywords: SNES, nonlinear, set, convergence, tolerances 2270 2271 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagPreconditioner(), SNESGetLagJacobian() 2272 2273 @*/ 2274 PetscErrorCode SNESSetLagJacobian(SNES snes,PetscInt lag) 2275 { 2276 PetscFunctionBegin; 2277 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2278 if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater"); 2279 if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0"); 2280 PetscValidLogicalCollectiveInt(snes,lag,2); 2281 snes->lagjacobian = lag; 2282 PetscFunctionReturn(0); 2283 } 2284 2285 #undef __FUNCT__ 2286 #define __FUNCT__ "SNESGetLagJacobian" 2287 /*@ 2288 SNESGetLagJacobian - Indicates how often the Jacobian is rebuilt. See SNESGetLagPreconditioner() to determine when the preconditioner is rebuilt 2289 2290 Not Collective 2291 2292 Input Parameter: 2293 . snes - the SNES context 2294 2295 Output Parameter: 2296 . lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time 2297 the Jacobian is built etc. 2298 2299 Options Database Keys: 2300 . -snes_lag_jacobian <lag> 2301 2302 Notes: 2303 The default is 1 2304 The jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 2305 2306 Level: intermediate 2307 2308 .keywords: SNES, nonlinear, set, convergence, tolerances 2309 2310 .seealso: SNESSetTrustRegionTolerance(), SNESSetLagJacobian(), SNESSetLagPreconditioner(), SNESGetLagPreconditioner() 2311 2312 @*/ 2313 PetscErrorCode SNESGetLagJacobian(SNES snes,PetscInt *lag) 2314 { 2315 PetscFunctionBegin; 2316 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2317 *lag = snes->lagjacobian; 2318 PetscFunctionReturn(0); 2319 } 2320 2321 #undef __FUNCT__ 2322 #define __FUNCT__ "SNESSetTolerances" 2323 /*@ 2324 SNESSetTolerances - Sets various parameters used in convergence tests. 2325 2326 Logically Collective on SNES 2327 2328 Input Parameters: 2329 + snes - the SNES context 2330 . abstol - absolute convergence tolerance 2331 . rtol - relative convergence tolerance 2332 . stol - convergence tolerance in terms of the norm 2333 of the change in the solution between steps 2334 . maxit - maximum number of iterations 2335 - maxf - maximum number of function evaluations 2336 2337 Options Database Keys: 2338 + -snes_atol <abstol> - Sets abstol 2339 . -snes_rtol <rtol> - Sets rtol 2340 . -snes_stol <stol> - Sets stol 2341 . -snes_max_it <maxit> - Sets maxit 2342 - -snes_max_funcs <maxf> - Sets maxf 2343 2344 Notes: 2345 The default maximum number of iterations is 50. 2346 The default maximum number of function evaluations is 1000. 2347 2348 Level: intermediate 2349 2350 .keywords: SNES, nonlinear, set, convergence, tolerances 2351 2352 .seealso: SNESSetTrustRegionTolerance() 2353 @*/ 2354 PetscErrorCode SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf) 2355 { 2356 PetscFunctionBegin; 2357 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2358 PetscValidLogicalCollectiveReal(snes,abstol,2); 2359 PetscValidLogicalCollectiveReal(snes,rtol,3); 2360 PetscValidLogicalCollectiveReal(snes,stol,4); 2361 PetscValidLogicalCollectiveInt(snes,maxit,5); 2362 PetscValidLogicalCollectiveInt(snes,maxf,6); 2363 2364 if (abstol != PETSC_DEFAULT) { 2365 if (abstol < 0.0) SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %G must be non-negative",abstol); 2366 snes->abstol = abstol; 2367 } 2368 if (rtol != PETSC_DEFAULT) { 2369 if (rtol < 0.0 || 1.0 <= rtol) SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Relative tolerance %G must be non-negative and less than 1.0",rtol); 2370 snes->rtol = rtol; 2371 } 2372 if (stol != PETSC_DEFAULT) { 2373 if (stol < 0.0) SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Step tolerance %G must be non-negative",stol); 2374 snes->xtol = stol; 2375 } 2376 if (maxit != PETSC_DEFAULT) { 2377 if (maxit < 0) SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxit); 2378 snes->max_its = maxit; 2379 } 2380 if (maxf != PETSC_DEFAULT) { 2381 if (maxf < 0) SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of function evaluations %D must be non-negative",maxf); 2382 snes->max_funcs = maxf; 2383 } 2384 PetscFunctionReturn(0); 2385 } 2386 2387 #undef __FUNCT__ 2388 #define __FUNCT__ "SNESGetTolerances" 2389 /*@ 2390 SNESGetTolerances - Gets various parameters used in convergence tests. 2391 2392 Not Collective 2393 2394 Input Parameters: 2395 + snes - the SNES context 2396 . atol - absolute convergence tolerance 2397 . rtol - relative convergence tolerance 2398 . stol - convergence tolerance in terms of the norm 2399 of the change in the solution between steps 2400 . maxit - maximum number of iterations 2401 - maxf - maximum number of function evaluations 2402 2403 Notes: 2404 The user can specify PETSC_NULL for any parameter that is not needed. 2405 2406 Level: intermediate 2407 2408 .keywords: SNES, nonlinear, get, convergence, tolerances 2409 2410 .seealso: SNESSetTolerances() 2411 @*/ 2412 PetscErrorCode SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf) 2413 { 2414 PetscFunctionBegin; 2415 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2416 if (atol) *atol = snes->abstol; 2417 if (rtol) *rtol = snes->rtol; 2418 if (stol) *stol = snes->xtol; 2419 if (maxit) *maxit = snes->max_its; 2420 if (maxf) *maxf = snes->max_funcs; 2421 PetscFunctionReturn(0); 2422 } 2423 2424 #undef __FUNCT__ 2425 #define __FUNCT__ "SNESSetTrustRegionTolerance" 2426 /*@ 2427 SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 2428 2429 Logically Collective on SNES 2430 2431 Input Parameters: 2432 + snes - the SNES context 2433 - tol - tolerance 2434 2435 Options Database Key: 2436 . -snes_trtol <tol> - Sets tol 2437 2438 Level: intermediate 2439 2440 .keywords: SNES, nonlinear, set, trust region, tolerance 2441 2442 .seealso: SNESSetTolerances() 2443 @*/ 2444 PetscErrorCode SNESSetTrustRegionTolerance(SNES snes,PetscReal tol) 2445 { 2446 PetscFunctionBegin; 2447 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2448 PetscValidLogicalCollectiveReal(snes,tol,2); 2449 snes->deltatol = tol; 2450 PetscFunctionReturn(0); 2451 } 2452 2453 /* 2454 Duplicate the lg monitors for SNES from KSP; for some reason with 2455 dynamic libraries things don't work under Sun4 if we just use 2456 macros instead of functions 2457 */ 2458 #undef __FUNCT__ 2459 #define __FUNCT__ "SNESMonitorLG" 2460 PetscErrorCode SNESMonitorLG(SNES snes,PetscInt it,PetscReal norm,void *ctx) 2461 { 2462 PetscErrorCode ierr; 2463 2464 PetscFunctionBegin; 2465 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2466 ierr = KSPMonitorLG((KSP)snes,it,norm,ctx);CHKERRQ(ierr); 2467 PetscFunctionReturn(0); 2468 } 2469 2470 #undef __FUNCT__ 2471 #define __FUNCT__ "SNESMonitorLGCreate" 2472 PetscErrorCode SNESMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 2473 { 2474 PetscErrorCode ierr; 2475 2476 PetscFunctionBegin; 2477 ierr = KSPMonitorLGCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr); 2478 PetscFunctionReturn(0); 2479 } 2480 2481 #undef __FUNCT__ 2482 #define __FUNCT__ "SNESMonitorLGDestroy" 2483 PetscErrorCode SNESMonitorLGDestroy(PetscDrawLG *draw) 2484 { 2485 PetscErrorCode ierr; 2486 2487 PetscFunctionBegin; 2488 ierr = KSPMonitorLGDestroy(draw);CHKERRQ(ierr); 2489 PetscFunctionReturn(0); 2490 } 2491 2492 extern PetscErrorCode SNESMonitorRange_Private(SNES,PetscInt,PetscReal*); 2493 #undef __FUNCT__ 2494 #define __FUNCT__ "SNESMonitorLGRange" 2495 PetscErrorCode SNESMonitorLGRange(SNES snes,PetscInt n,PetscReal rnorm,void *monctx) 2496 { 2497 PetscDrawLG lg; 2498 PetscErrorCode ierr; 2499 PetscReal x,y,per; 2500 PetscViewer v = (PetscViewer)monctx; 2501 static PetscReal prev; /* should be in the context */ 2502 PetscDraw draw; 2503 PetscFunctionBegin; 2504 if (!monctx) { 2505 MPI_Comm comm; 2506 2507 ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 2508 v = PETSC_VIEWER_DRAW_(comm); 2509 } 2510 ierr = PetscViewerDrawGetDrawLG(v,0,&lg);CHKERRQ(ierr); 2511 if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 2512 ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr); 2513 ierr = PetscDrawSetTitle(draw,"Residual norm");CHKERRQ(ierr); 2514 x = (PetscReal) n; 2515 if (rnorm > 0.0) y = log10(rnorm); else y = -15.0; 2516 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 2517 if (n < 20 || !(n % 5)) { 2518 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 2519 } 2520 2521 ierr = PetscViewerDrawGetDrawLG(v,1,&lg);CHKERRQ(ierr); 2522 if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 2523 ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr); 2524 ierr = PetscDrawSetTitle(draw,"% elemts > .2*max elemt");CHKERRQ(ierr); 2525 ierr = SNESMonitorRange_Private(snes,n,&per);CHKERRQ(ierr); 2526 x = (PetscReal) n; 2527 y = 100.0*per; 2528 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 2529 if (n < 20 || !(n % 5)) { 2530 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 2531 } 2532 2533 ierr = PetscViewerDrawGetDrawLG(v,2,&lg);CHKERRQ(ierr); 2534 if (!n) {prev = rnorm;ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 2535 ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr); 2536 ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm");CHKERRQ(ierr); 2537 x = (PetscReal) n; 2538 y = (prev - rnorm)/prev; 2539 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 2540 if (n < 20 || !(n % 5)) { 2541 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 2542 } 2543 2544 ierr = PetscViewerDrawGetDrawLG(v,3,&lg);CHKERRQ(ierr); 2545 if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 2546 ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr); 2547 ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");CHKERRQ(ierr); 2548 x = (PetscReal) n; 2549 y = (prev - rnorm)/(prev*per); 2550 if (n > 2) { /*skip initial crazy value */ 2551 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 2552 } 2553 if (n < 20 || !(n % 5)) { 2554 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 2555 } 2556 prev = rnorm; 2557 PetscFunctionReturn(0); 2558 } 2559 2560 #undef __FUNCT__ 2561 #define __FUNCT__ "SNESMonitorLGRangeCreate" 2562 PetscErrorCode SNESMonitorLGRangeCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 2563 { 2564 PetscErrorCode ierr; 2565 2566 PetscFunctionBegin; 2567 ierr = KSPMonitorLGCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr); 2568 PetscFunctionReturn(0); 2569 } 2570 2571 #undef __FUNCT__ 2572 #define __FUNCT__ "SNESMonitorLGRangeDestroy" 2573 PetscErrorCode SNESMonitorLGRangeDestroy(PetscDrawLG *draw) 2574 { 2575 PetscErrorCode ierr; 2576 2577 PetscFunctionBegin; 2578 ierr = KSPMonitorLGDestroy(draw);CHKERRQ(ierr); 2579 PetscFunctionReturn(0); 2580 } 2581 2582 #undef __FUNCT__ 2583 #define __FUNCT__ "SNESMonitor" 2584 /*@ 2585 SNESMonitor - runs the user provided monitor routines, if they exist 2586 2587 Collective on SNES 2588 2589 Input Parameters: 2590 + snes - nonlinear solver context obtained from SNESCreate() 2591 . iter - iteration number 2592 - rnorm - relative norm of the residual 2593 2594 Notes: 2595 This routine is called by the SNES implementations. 2596 It does not typically need to be called by the user. 2597 2598 Level: developer 2599 2600 .seealso: SNESMonitorSet() 2601 @*/ 2602 PetscErrorCode SNESMonitor(SNES snes,PetscInt iter,PetscReal rnorm) 2603 { 2604 PetscErrorCode ierr; 2605 PetscInt i,n = snes->numbermonitors; 2606 2607 PetscFunctionBegin; 2608 for (i=0; i<n; i++) { 2609 ierr = (*snes->monitor[i])(snes,iter,rnorm,snes->monitorcontext[i]);CHKERRQ(ierr); 2610 } 2611 PetscFunctionReturn(0); 2612 } 2613 2614 /* ------------ Routines to set performance monitoring options ----------- */ 2615 2616 #undef __FUNCT__ 2617 #define __FUNCT__ "SNESMonitorSet" 2618 /*@C 2619 SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every 2620 iteration of the nonlinear solver to display the iteration's 2621 progress. 2622 2623 Logically Collective on SNES 2624 2625 Input Parameters: 2626 + snes - the SNES context 2627 . func - monitoring routine 2628 . mctx - [optional] user-defined context for private data for the 2629 monitor routine (use PETSC_NULL if no context is desired) 2630 - monitordestroy - [optional] routine that frees monitor context 2631 (may be PETSC_NULL) 2632 2633 Calling sequence of func: 2634 $ int func(SNES snes,PetscInt its, PetscReal norm,void *mctx) 2635 2636 + snes - the SNES context 2637 . its - iteration number 2638 . norm - 2-norm function value (may be estimated) 2639 - mctx - [optional] monitoring context 2640 2641 Options Database Keys: 2642 + -snes_monitor - sets SNESMonitorDefault() 2643 . -snes_monitor_draw - sets line graph monitor, 2644 uses SNESMonitorLGCreate() 2645 - -snes_monitor_cancel - cancels all monitors that have 2646 been hardwired into a code by 2647 calls to SNESMonitorSet(), but 2648 does not cancel those set via 2649 the options database. 2650 2651 Notes: 2652 Several different monitoring routines may be set by calling 2653 SNESMonitorSet() multiple times; all will be called in the 2654 order in which they were set. 2655 2656 Fortran notes: Only a single monitor function can be set for each SNES object 2657 2658 Level: intermediate 2659 2660 .keywords: SNES, nonlinear, set, monitor 2661 2662 .seealso: SNESMonitorDefault(), SNESMonitorCancel() 2663 @*/ 2664 PetscErrorCode SNESMonitorSet(SNES snes,PetscErrorCode (*monitor)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**)) 2665 { 2666 PetscInt i; 2667 PetscErrorCode ierr; 2668 2669 PetscFunctionBegin; 2670 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2671 if (snes->numbermonitors >= MAXSNESMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 2672 for (i=0; i<snes->numbermonitors;i++) { 2673 if (monitor == snes->monitor[i] && monitordestroy == snes->monitordestroy[i] && mctx == snes->monitorcontext[i]) { 2674 if (monitordestroy) { 2675 ierr = (*monitordestroy)(&mctx);CHKERRQ(ierr); 2676 } 2677 PetscFunctionReturn(0); 2678 } 2679 } 2680 snes->monitor[snes->numbermonitors] = monitor; 2681 snes->monitordestroy[snes->numbermonitors] = monitordestroy; 2682 snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 2683 PetscFunctionReturn(0); 2684 } 2685 2686 #undef __FUNCT__ 2687 #define __FUNCT__ "SNESMonitorCancel" 2688 /*@C 2689 SNESMonitorCancel - Clears all the monitor functions for a SNES object. 2690 2691 Logically Collective on SNES 2692 2693 Input Parameters: 2694 . snes - the SNES context 2695 2696 Options Database Key: 2697 . -snes_monitor_cancel - cancels all monitors that have been hardwired 2698 into a code by calls to SNESMonitorSet(), but does not cancel those 2699 set via the options database 2700 2701 Notes: 2702 There is no way to clear one specific monitor from a SNES object. 2703 2704 Level: intermediate 2705 2706 .keywords: SNES, nonlinear, set, monitor 2707 2708 .seealso: SNESMonitorDefault(), SNESMonitorSet() 2709 @*/ 2710 PetscErrorCode SNESMonitorCancel(SNES snes) 2711 { 2712 PetscErrorCode ierr; 2713 PetscInt i; 2714 2715 PetscFunctionBegin; 2716 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2717 for (i=0; i<snes->numbermonitors; i++) { 2718 if (snes->monitordestroy[i]) { 2719 ierr = (*snes->monitordestroy[i])(&snes->monitorcontext[i]);CHKERRQ(ierr); 2720 } 2721 } 2722 snes->numbermonitors = 0; 2723 PetscFunctionReturn(0); 2724 } 2725 2726 #undef __FUNCT__ 2727 #define __FUNCT__ "SNESSetConvergenceTest" 2728 /*@C 2729 SNESSetConvergenceTest - Sets the function that is to be used 2730 to test for convergence of the nonlinear iterative solution. 2731 2732 Logically Collective on SNES 2733 2734 Input Parameters: 2735 + snes - the SNES context 2736 . func - routine to test for convergence 2737 . cctx - [optional] context for private data for the convergence routine (may be PETSC_NULL) 2738 - destroy - [optional] destructor for the context (may be PETSC_NULL; PETSC_NULL_FUNCTION in Fortran) 2739 2740 Calling sequence of func: 2741 $ PetscErrorCode func (SNES snes,PetscInt it,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx) 2742 2743 + snes - the SNES context 2744 . it - current iteration (0 is the first and is before any Newton step) 2745 . cctx - [optional] convergence context 2746 . reason - reason for convergence/divergence 2747 . xnorm - 2-norm of current iterate 2748 . gnorm - 2-norm of current step 2749 - f - 2-norm of function 2750 2751 Level: advanced 2752 2753 .keywords: SNES, nonlinear, set, convergence, test 2754 2755 .seealso: SNESDefaultConverged(), SNESSkipConverged() 2756 @*/ 2757 PetscErrorCode SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*)) 2758 { 2759 PetscErrorCode ierr; 2760 2761 PetscFunctionBegin; 2762 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2763 if (!func) func = SNESSkipConverged; 2764 if (snes->ops->convergeddestroy) { 2765 ierr = (*snes->ops->convergeddestroy)(snes->cnvP);CHKERRQ(ierr); 2766 } 2767 snes->ops->converged = func; 2768 snes->ops->convergeddestroy = destroy; 2769 snes->cnvP = cctx; 2770 PetscFunctionReturn(0); 2771 } 2772 2773 #undef __FUNCT__ 2774 #define __FUNCT__ "SNESGetConvergedReason" 2775 /*@ 2776 SNESGetConvergedReason - Gets the reason the SNES iteration was stopped. 2777 2778 Not Collective 2779 2780 Input Parameter: 2781 . snes - the SNES context 2782 2783 Output Parameter: 2784 . reason - negative value indicates diverged, positive value converged, see SNESConvergedReason or the 2785 manual pages for the individual convergence tests for complete lists 2786 2787 Level: intermediate 2788 2789 Notes: Can only be called after the call the SNESSolve() is complete. 2790 2791 .keywords: SNES, nonlinear, set, convergence, test 2792 2793 .seealso: SNESSetConvergenceTest(), SNESConvergedReason 2794 @*/ 2795 PetscErrorCode SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason) 2796 { 2797 PetscFunctionBegin; 2798 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2799 PetscValidPointer(reason,2); 2800 *reason = snes->reason; 2801 PetscFunctionReturn(0); 2802 } 2803 2804 #undef __FUNCT__ 2805 #define __FUNCT__ "SNESSetConvergenceHistory" 2806 /*@ 2807 SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 2808 2809 Logically Collective on SNES 2810 2811 Input Parameters: 2812 + snes - iterative context obtained from SNESCreate() 2813 . a - array to hold history, this array will contain the function norms computed at each step 2814 . its - integer array holds the number of linear iterations for each solve. 2815 . na - size of a and its 2816 - reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero, 2817 else it continues storing new values for new nonlinear solves after the old ones 2818 2819 Notes: 2820 If 'a' and 'its' are PETSC_NULL then space is allocated for the history. If 'na' PETSC_DECIDE or PETSC_DEFAULT then a 2821 default array of length 10000 is allocated. 2822 2823 This routine is useful, e.g., when running a code for purposes 2824 of accurate performance monitoring, when no I/O should be done 2825 during the section of code that is being timed. 2826 2827 Level: intermediate 2828 2829 .keywords: SNES, set, convergence, history 2830 2831 .seealso: SNESGetConvergenceHistory() 2832 2833 @*/ 2834 PetscErrorCode SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscBool reset) 2835 { 2836 PetscErrorCode ierr; 2837 2838 PetscFunctionBegin; 2839 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2840 if (na) PetscValidScalarPointer(a,2); 2841 if (its) PetscValidIntPointer(its,3); 2842 if (na == PETSC_DECIDE || na == PETSC_DEFAULT || !a) { 2843 if (na == PETSC_DECIDE || na == PETSC_DEFAULT) na = 1000; 2844 ierr = PetscMalloc(na*sizeof(PetscReal),&a);CHKERRQ(ierr); 2845 ierr = PetscMalloc(na*sizeof(PetscInt),&its);CHKERRQ(ierr); 2846 snes->conv_malloc = PETSC_TRUE; 2847 } 2848 snes->conv_hist = a; 2849 snes->conv_hist_its = its; 2850 snes->conv_hist_max = na; 2851 snes->conv_hist_len = 0; 2852 snes->conv_hist_reset = reset; 2853 PetscFunctionReturn(0); 2854 } 2855 2856 #if defined(PETSC_HAVE_MATLAB_ENGINE) 2857 #include <engine.h> /* MATLAB include file */ 2858 #include <mex.h> /* MATLAB include file */ 2859 EXTERN_C_BEGIN 2860 #undef __FUNCT__ 2861 #define __FUNCT__ "SNESGetConvergenceHistoryMatlab" 2862 mxArray *SNESGetConvergenceHistoryMatlab(SNES snes) 2863 { 2864 mxArray *mat; 2865 PetscInt i; 2866 PetscReal *ar; 2867 2868 PetscFunctionBegin; 2869 mat = mxCreateDoubleMatrix(snes->conv_hist_len,1,mxREAL); 2870 ar = (PetscReal*) mxGetData(mat); 2871 for (i=0; i<snes->conv_hist_len; i++) { 2872 ar[i] = snes->conv_hist[i]; 2873 } 2874 PetscFunctionReturn(mat); 2875 } 2876 EXTERN_C_END 2877 #endif 2878 2879 2880 #undef __FUNCT__ 2881 #define __FUNCT__ "SNESGetConvergenceHistory" 2882 /*@C 2883 SNESGetConvergenceHistory - Gets the array used to hold the convergence history. 2884 2885 Not Collective 2886 2887 Input Parameter: 2888 . snes - iterative context obtained from SNESCreate() 2889 2890 Output Parameters: 2891 . a - array to hold history 2892 . its - integer array holds the number of linear iterations (or 2893 negative if not converged) for each solve. 2894 - na - size of a and its 2895 2896 Notes: 2897 The calling sequence for this routine in Fortran is 2898 $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr) 2899 2900 This routine is useful, e.g., when running a code for purposes 2901 of accurate performance monitoring, when no I/O should be done 2902 during the section of code that is being timed. 2903 2904 Level: intermediate 2905 2906 .keywords: SNES, get, convergence, history 2907 2908 .seealso: SNESSetConvergencHistory() 2909 2910 @*/ 2911 PetscErrorCode SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na) 2912 { 2913 PetscFunctionBegin; 2914 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2915 if (a) *a = snes->conv_hist; 2916 if (its) *its = snes->conv_hist_its; 2917 if (na) *na = snes->conv_hist_len; 2918 PetscFunctionReturn(0); 2919 } 2920 2921 #undef __FUNCT__ 2922 #define __FUNCT__ "SNESSetUpdate" 2923 /*@C 2924 SNESSetUpdate - Sets the general-purpose update function called 2925 at the beginning of every iteration of the nonlinear solve. Specifically 2926 it is called just before the Jacobian is "evaluated". 2927 2928 Logically Collective on SNES 2929 2930 Input Parameters: 2931 . snes - The nonlinear solver context 2932 . func - The function 2933 2934 Calling sequence of func: 2935 . func (SNES snes, PetscInt step); 2936 2937 . step - The current step of the iteration 2938 2939 Level: advanced 2940 2941 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() 2942 This is not used by most users. 2943 2944 .keywords: SNES, update 2945 2946 .seealso SNESDefaultUpdate(), SNESSetJacobian(), SNESSolve() 2947 @*/ 2948 PetscErrorCode SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt)) 2949 { 2950 PetscFunctionBegin; 2951 PetscValidHeaderSpecific(snes, SNES_CLASSID,1); 2952 snes->ops->update = func; 2953 PetscFunctionReturn(0); 2954 } 2955 2956 #undef __FUNCT__ 2957 #define __FUNCT__ "SNESDefaultUpdate" 2958 /*@ 2959 SNESDefaultUpdate - The default update function which does nothing. 2960 2961 Not collective 2962 2963 Input Parameters: 2964 . snes - The nonlinear solver context 2965 . step - The current step of the iteration 2966 2967 Level: intermediate 2968 2969 .keywords: SNES, update 2970 .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultShortolutionBC() 2971 @*/ 2972 PetscErrorCode SNESDefaultUpdate(SNES snes, PetscInt step) 2973 { 2974 PetscFunctionBegin; 2975 PetscFunctionReturn(0); 2976 } 2977 2978 #undef __FUNCT__ 2979 #define __FUNCT__ "SNESScaleStep_Private" 2980 /* 2981 SNESScaleStep_Private - Scales a step so that its length is less than the 2982 positive parameter delta. 2983 2984 Input Parameters: 2985 + snes - the SNES context 2986 . y - approximate solution of linear system 2987 . fnorm - 2-norm of current function 2988 - delta - trust region size 2989 2990 Output Parameters: 2991 + gpnorm - predicted function norm at the new point, assuming local 2992 linearization. The value is zero if the step lies within the trust 2993 region, and exceeds zero otherwise. 2994 - ynorm - 2-norm of the step 2995 2996 Note: 2997 For non-trust region methods such as SNESLS, the parameter delta 2998 is set to be the maximum allowable step size. 2999 3000 .keywords: SNES, nonlinear, scale, step 3001 */ 3002 PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm) 3003 { 3004 PetscReal nrm; 3005 PetscScalar cnorm; 3006 PetscErrorCode ierr; 3007 3008 PetscFunctionBegin; 3009 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3010 PetscValidHeaderSpecific(y,VEC_CLASSID,2); 3011 PetscCheckSameComm(snes,1,y,2); 3012 3013 ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr); 3014 if (nrm > *delta) { 3015 nrm = *delta/nrm; 3016 *gpnorm = (1.0 - nrm)*(*fnorm); 3017 cnorm = nrm; 3018 ierr = VecScale(y,cnorm);CHKERRQ(ierr); 3019 *ynorm = *delta; 3020 } else { 3021 *gpnorm = 0.0; 3022 *ynorm = nrm; 3023 } 3024 PetscFunctionReturn(0); 3025 } 3026 3027 #undef __FUNCT__ 3028 #define __FUNCT__ "SNESSolve" 3029 /*@C 3030 SNESSolve - Solves a nonlinear system F(x) = b. 3031 Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX(). 3032 3033 Collective on SNES 3034 3035 Input Parameters: 3036 + snes - the SNES context 3037 . b - the constant part of the equation F(x) = b, or PETSC_NULL to use zero. 3038 - x - the solution vector. 3039 3040 Notes: 3041 The user should initialize the vector,x, with the initial guess 3042 for the nonlinear solve prior to calling SNESSolve. In particular, 3043 to employ an initial guess of zero, the user should explicitly set 3044 this vector to zero by calling VecSet(). 3045 3046 Level: beginner 3047 3048 .keywords: SNES, nonlinear, solve 3049 3050 .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetGridSequence(), SNESGetSolution() 3051 @*/ 3052 PetscErrorCode SNESSolve(SNES snes,Vec b,Vec x) 3053 { 3054 PetscErrorCode ierr; 3055 PetscBool flg; 3056 char filename[PETSC_MAX_PATH_LEN]; 3057 PetscViewer viewer; 3058 PetscInt grid; 3059 Vec xcreated = PETSC_NULL; 3060 3061 PetscFunctionBegin; 3062 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3063 if (x) PetscValidHeaderSpecific(x,VEC_CLASSID,3); 3064 if (x) PetscCheckSameComm(snes,1,x,3); 3065 if (b) PetscValidHeaderSpecific(b,VEC_CLASSID,2); 3066 if (b) PetscCheckSameComm(snes,1,b,2); 3067 3068 if (!x && snes->dm) { 3069 ierr = DMCreateGlobalVector(snes->dm,&xcreated);CHKERRQ(ierr); 3070 x = xcreated; 3071 } 3072 3073 for (grid=0; grid<snes->gridsequence; grid++) {ierr = PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm));CHKERRQ(ierr);} 3074 for (grid=0; grid<snes->gridsequence+1; grid++) { 3075 3076 /* set solution vector */ 3077 if (!grid) {ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);} 3078 ierr = VecDestroy(&snes->vec_sol);CHKERRQ(ierr); 3079 snes->vec_sol = x; 3080 /* set afine vector if provided */ 3081 if (b) { ierr = PetscObjectReference((PetscObject)b);CHKERRQ(ierr); } 3082 ierr = VecDestroy(&snes->vec_rhs);CHKERRQ(ierr); 3083 snes->vec_rhs = b; 3084 3085 ierr = SNESSetUp(snes);CHKERRQ(ierr); 3086 3087 if (!grid) { 3088 if (snes->ops->computeinitialguess) { 3089 ierr = (*snes->ops->computeinitialguess)(snes,snes->vec_sol,snes->initialguessP);CHKERRQ(ierr); 3090 } else if (snes->dm) { 3091 PetscBool ig; 3092 ierr = DMHasInitialGuess(snes->dm,&ig);CHKERRQ(ierr); 3093 if (ig) { 3094 ierr = DMComputeInitialGuess(snes->dm,snes->vec_sol);CHKERRQ(ierr); 3095 } 3096 } 3097 } 3098 3099 if (snes->conv_hist_reset) snes->conv_hist_len = 0; 3100 snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0; 3101 3102 ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 3103 ierr = (*snes->ops->solve)(snes);CHKERRQ(ierr); 3104 ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 3105 if (snes->domainerror){ 3106 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 3107 snes->domainerror = PETSC_FALSE; 3108 } 3109 if (!snes->reason) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason"); 3110 3111 ierr = PetscOptionsGetString(((PetscObject)snes)->prefix,"-snes_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 3112 if (flg && !PetscPreLoadingOn) { 3113 ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,filename,&viewer);CHKERRQ(ierr); 3114 ierr = SNESView(snes,viewer);CHKERRQ(ierr); 3115 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 3116 } 3117 3118 flg = PETSC_FALSE; 3119 ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_test_local_min",&flg,PETSC_NULL);CHKERRQ(ierr); 3120 if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); } 3121 if (snes->printreason) { 3122 ierr = PetscViewerASCIIAddTab(PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm),((PetscObject)snes)->tablevel);CHKERRQ(ierr); 3123 if (snes->reason > 0) { 3124 ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm),"Nonlinear solve converged due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr); 3125 } else { 3126 ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm),"Nonlinear solve did not converge due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr); 3127 } 3128 ierr = PetscViewerASCIISubtractTab(PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm),((PetscObject)snes)->tablevel);CHKERRQ(ierr); 3129 } 3130 3131 flg = PETSC_FALSE; 3132 ierr = PetscOptionsGetString(((PetscObject)snes)->prefix,"-snes_view_solution_vtk",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 3133 if (flg) { 3134 PetscViewer viewer; 3135 ierr = PetscViewerCreate(((PetscObject)snes)->comm,&viewer);CHKERRQ(ierr); 3136 ierr = PetscViewerSetType(viewer,PETSCVIEWERVTK);CHKERRQ(ierr); 3137 ierr = PetscViewerFileSetName(viewer,filename);CHKERRQ(ierr); 3138 ierr = VecView(snes->vec_sol,viewer);CHKERRQ(ierr); 3139 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 3140 } 3141 3142 if (snes->errorifnotconverged && snes->reason < 0) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_NOT_CONVERGED,"SNESSolve has not converged"); 3143 if (grid < snes->gridsequence) { 3144 DM fine; 3145 Vec xnew; 3146 Mat interp; 3147 3148 ierr = DMRefine(snes->dm,((PetscObject)snes)->comm,&fine);CHKERRQ(ierr); 3149 ierr = DMCreateInterpolation(snes->dm,fine,&interp,PETSC_NULL);CHKERRQ(ierr); 3150 ierr = DMCreateGlobalVector(fine,&xnew);CHKERRQ(ierr); 3151 ierr = MatInterpolate(interp,x,xnew);CHKERRQ(ierr); 3152 ierr = MatDestroy(&interp);CHKERRQ(ierr); 3153 x = xnew; 3154 3155 ierr = SNESReset(snes);CHKERRQ(ierr); 3156 ierr = SNESSetDM(snes,fine);CHKERRQ(ierr); 3157 ierr = DMDestroy(&fine);CHKERRQ(ierr); 3158 ierr = PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm));CHKERRQ(ierr); 3159 } 3160 } 3161 ierr = VecDestroy(&xcreated);CHKERRQ(ierr); 3162 PetscFunctionReturn(0); 3163 } 3164 3165 /* --------- Internal routines for SNES Package --------- */ 3166 3167 #undef __FUNCT__ 3168 #define __FUNCT__ "SNESSetType" 3169 /*@C 3170 SNESSetType - Sets the method for the nonlinear solver. 3171 3172 Collective on SNES 3173 3174 Input Parameters: 3175 + snes - the SNES context 3176 - type - a known method 3177 3178 Options Database Key: 3179 . -snes_type <type> - Sets the method; use -help for a list 3180 of available methods (for instance, ls or tr) 3181 3182 Notes: 3183 See "petsc/include/petscsnes.h" for available methods (for instance) 3184 + SNESLS - Newton's method with line search 3185 (systems of nonlinear equations) 3186 . SNESTR - Newton's method with trust region 3187 (systems of nonlinear equations) 3188 3189 Normally, it is best to use the SNESSetFromOptions() command and then 3190 set the SNES solver type from the options database rather than by using 3191 this routine. Using the options database provides the user with 3192 maximum flexibility in evaluating the many nonlinear solvers. 3193 The SNESSetType() routine is provided for those situations where it 3194 is necessary to set the nonlinear solver independently of the command 3195 line or options database. This might be the case, for example, when 3196 the choice of solver changes during the execution of the program, 3197 and the user's application is taking responsibility for choosing the 3198 appropriate method. 3199 3200 Level: intermediate 3201 3202 .keywords: SNES, set, type 3203 3204 .seealso: SNESType, SNESCreate() 3205 3206 @*/ 3207 PetscErrorCode SNESSetType(SNES snes,const SNESType type) 3208 { 3209 PetscErrorCode ierr,(*r)(SNES); 3210 PetscBool match; 3211 3212 PetscFunctionBegin; 3213 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3214 PetscValidCharPointer(type,2); 3215 3216 ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr); 3217 if (match) PetscFunctionReturn(0); 3218 3219 ierr = PetscFListFind(SNESList,((PetscObject)snes)->comm,type,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr); 3220 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type); 3221 /* Destroy the previous private SNES context */ 3222 if (snes->ops->destroy) { ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr); } 3223 /* Reinitialize function pointers in SNESOps structure */ 3224 snes->ops->setup = 0; 3225 snes->ops->solve = 0; 3226 snes->ops->view = 0; 3227 snes->ops->setfromoptions = 0; 3228 snes->ops->destroy = 0; 3229 /* Call the SNESCreate_XXX routine for this particular Nonlinear solver */ 3230 snes->setupcalled = PETSC_FALSE; 3231 ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr); 3232 ierr = (*r)(snes);CHKERRQ(ierr); 3233 #if defined(PETSC_HAVE_AMS) 3234 if (PetscAMSPublishAll) { 3235 ierr = PetscObjectAMSPublish((PetscObject)snes);CHKERRQ(ierr); 3236 } 3237 #endif 3238 PetscFunctionReturn(0); 3239 } 3240 3241 3242 /* --------------------------------------------------------------------- */ 3243 #undef __FUNCT__ 3244 #define __FUNCT__ "SNESRegisterDestroy" 3245 /*@ 3246 SNESRegisterDestroy - Frees the list of nonlinear solvers that were 3247 registered by SNESRegisterDynamic(). 3248 3249 Not Collective 3250 3251 Level: advanced 3252 3253 .keywords: SNES, nonlinear, register, destroy 3254 3255 .seealso: SNESRegisterAll(), SNESRegisterAll() 3256 @*/ 3257 PetscErrorCode SNESRegisterDestroy(void) 3258 { 3259 PetscErrorCode ierr; 3260 3261 PetscFunctionBegin; 3262 ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr); 3263 SNESRegisterAllCalled = PETSC_FALSE; 3264 PetscFunctionReturn(0); 3265 } 3266 3267 #undef __FUNCT__ 3268 #define __FUNCT__ "SNESGetType" 3269 /*@C 3270 SNESGetType - Gets the SNES method type and name (as a string). 3271 3272 Not Collective 3273 3274 Input Parameter: 3275 . snes - nonlinear solver context 3276 3277 Output Parameter: 3278 . type - SNES method (a character string) 3279 3280 Level: intermediate 3281 3282 .keywords: SNES, nonlinear, get, type, name 3283 @*/ 3284 PetscErrorCode SNESGetType(SNES snes,const SNESType *type) 3285 { 3286 PetscFunctionBegin; 3287 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3288 PetscValidPointer(type,2); 3289 *type = ((PetscObject)snes)->type_name; 3290 PetscFunctionReturn(0); 3291 } 3292 3293 #undef __FUNCT__ 3294 #define __FUNCT__ "SNESGetSolution" 3295 /*@ 3296 SNESGetSolution - Returns the vector where the approximate solution is 3297 stored. This is the fine grid solution when using SNESSetGridSequence(). 3298 3299 Not Collective, but Vec is parallel if SNES is parallel 3300 3301 Input Parameter: 3302 . snes - the SNES context 3303 3304 Output Parameter: 3305 . x - the solution 3306 3307 Level: intermediate 3308 3309 .keywords: SNES, nonlinear, get, solution 3310 3311 .seealso: SNESGetSolutionUpdate(), SNESGetFunction() 3312 @*/ 3313 PetscErrorCode SNESGetSolution(SNES snes,Vec *x) 3314 { 3315 PetscFunctionBegin; 3316 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3317 PetscValidPointer(x,2); 3318 *x = snes->vec_sol; 3319 PetscFunctionReturn(0); 3320 } 3321 3322 #undef __FUNCT__ 3323 #define __FUNCT__ "SNESGetSolutionUpdate" 3324 /*@ 3325 SNESGetSolutionUpdate - Returns the vector where the solution update is 3326 stored. 3327 3328 Not Collective, but Vec is parallel if SNES is parallel 3329 3330 Input Parameter: 3331 . snes - the SNES context 3332 3333 Output Parameter: 3334 . x - the solution update 3335 3336 Level: advanced 3337 3338 .keywords: SNES, nonlinear, get, solution, update 3339 3340 .seealso: SNESGetSolution(), SNESGetFunction() 3341 @*/ 3342 PetscErrorCode SNESGetSolutionUpdate(SNES snes,Vec *x) 3343 { 3344 PetscFunctionBegin; 3345 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3346 PetscValidPointer(x,2); 3347 *x = snes->vec_sol_update; 3348 PetscFunctionReturn(0); 3349 } 3350 3351 #undef __FUNCT__ 3352 #define __FUNCT__ "SNESGetFunction" 3353 /*@C 3354 SNESGetFunction - Returns the vector where the function is stored. 3355 3356 Not Collective, but Vec is parallel if SNES is parallel. Collective if Vec is requested, but has not been created yet. 3357 3358 Input Parameter: 3359 . snes - the SNES context 3360 3361 Output Parameter: 3362 + r - the function (or PETSC_NULL) 3363 . func - the function (or PETSC_NULL) 3364 - ctx - the function context (or PETSC_NULL) 3365 3366 Level: advanced 3367 3368 .keywords: SNES, nonlinear, get, function 3369 3370 .seealso: SNESSetFunction(), SNESGetSolution() 3371 @*/ 3372 PetscErrorCode SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx) 3373 { 3374 PetscErrorCode ierr; 3375 3376 PetscFunctionBegin; 3377 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3378 if (r) { 3379 if (!snes->vec_func) { 3380 if (snes->vec_rhs) { 3381 ierr = VecDuplicate(snes->vec_rhs,&snes->vec_func);CHKERRQ(ierr); 3382 } else if (snes->vec_sol) { 3383 ierr = VecDuplicate(snes->vec_sol,&snes->vec_func);CHKERRQ(ierr); 3384 } else if (snes->dm) { 3385 ierr = DMCreateGlobalVector(snes->dm,&snes->vec_func);CHKERRQ(ierr); 3386 } 3387 } 3388 *r = snes->vec_func; 3389 } 3390 if (func) *func = snes->ops->computefunction; 3391 if (ctx) *ctx = snes->funP; 3392 PetscFunctionReturn(0); 3393 } 3394 3395 /*@C 3396 SNESGetGS - Returns the GS function and context. 3397 3398 Input Parameter: 3399 . snes - the SNES context 3400 3401 Output Parameter: 3402 + gsfunc - the function (or PETSC_NULL) 3403 - ctx - the function context (or PETSC_NULL) 3404 3405 Level: advanced 3406 3407 .keywords: SNES, nonlinear, get, function 3408 3409 .seealso: SNESSetGS(), SNESGetFunction() 3410 @*/ 3411 3412 #undef __FUNCT__ 3413 #define __FUNCT__ "SNESGetGS" 3414 PetscErrorCode SNESGetGS (SNES snes, PetscErrorCode(**func)(SNES, Vec, Vec, void*), void ** ctx) 3415 { 3416 PetscFunctionBegin; 3417 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3418 if (func) *func = snes->ops->computegs; 3419 if (ctx) *ctx = snes->funP; 3420 PetscFunctionReturn(0); 3421 } 3422 3423 #undef __FUNCT__ 3424 #define __FUNCT__ "SNESSetOptionsPrefix" 3425 /*@C 3426 SNESSetOptionsPrefix - Sets the prefix used for searching for all 3427 SNES options in the database. 3428 3429 Logically Collective on SNES 3430 3431 Input Parameter: 3432 + snes - the SNES context 3433 - prefix - the prefix to prepend to all option names 3434 3435 Notes: 3436 A hyphen (-) must NOT be given at the beginning of the prefix name. 3437 The first character of all runtime options is AUTOMATICALLY the hyphen. 3438 3439 Level: advanced 3440 3441 .keywords: SNES, set, options, prefix, database 3442 3443 .seealso: SNESSetFromOptions() 3444 @*/ 3445 PetscErrorCode SNESSetOptionsPrefix(SNES snes,const char prefix[]) 3446 { 3447 PetscErrorCode ierr; 3448 3449 PetscFunctionBegin; 3450 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3451 ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 3452 if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);} 3453 ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 3454 PetscFunctionReturn(0); 3455 } 3456 3457 #undef __FUNCT__ 3458 #define __FUNCT__ "SNESAppendOptionsPrefix" 3459 /*@C 3460 SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 3461 SNES options in the database. 3462 3463 Logically Collective on SNES 3464 3465 Input Parameters: 3466 + snes - the SNES context 3467 - prefix - the prefix to prepend to all option names 3468 3469 Notes: 3470 A hyphen (-) must NOT be given at the beginning of the prefix name. 3471 The first character of all runtime options is AUTOMATICALLY the hyphen. 3472 3473 Level: advanced 3474 3475 .keywords: SNES, append, options, prefix, database 3476 3477 .seealso: SNESGetOptionsPrefix() 3478 @*/ 3479 PetscErrorCode SNESAppendOptionsPrefix(SNES snes,const char prefix[]) 3480 { 3481 PetscErrorCode ierr; 3482 3483 PetscFunctionBegin; 3484 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3485 ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 3486 if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);} 3487 ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 3488 PetscFunctionReturn(0); 3489 } 3490 3491 #undef __FUNCT__ 3492 #define __FUNCT__ "SNESGetOptionsPrefix" 3493 /*@C 3494 SNESGetOptionsPrefix - Sets the prefix used for searching for all 3495 SNES options in the database. 3496 3497 Not Collective 3498 3499 Input Parameter: 3500 . snes - the SNES context 3501 3502 Output Parameter: 3503 . prefix - pointer to the prefix string used 3504 3505 Notes: On the fortran side, the user should pass in a string 'prefix' of 3506 sufficient length to hold the prefix. 3507 3508 Level: advanced 3509 3510 .keywords: SNES, get, options, prefix, database 3511 3512 .seealso: SNESAppendOptionsPrefix() 3513 @*/ 3514 PetscErrorCode SNESGetOptionsPrefix(SNES snes,const char *prefix[]) 3515 { 3516 PetscErrorCode ierr; 3517 3518 PetscFunctionBegin; 3519 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3520 ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 3521 PetscFunctionReturn(0); 3522 } 3523 3524 3525 #undef __FUNCT__ 3526 #define __FUNCT__ "SNESRegister" 3527 /*@C 3528 SNESRegister - See SNESRegisterDynamic() 3529 3530 Level: advanced 3531 @*/ 3532 PetscErrorCode SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES)) 3533 { 3534 char fullname[PETSC_MAX_PATH_LEN]; 3535 PetscErrorCode ierr; 3536 3537 PetscFunctionBegin; 3538 ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 3539 ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 3540 PetscFunctionReturn(0); 3541 } 3542 3543 #undef __FUNCT__ 3544 #define __FUNCT__ "SNESTestLocalMin" 3545 PetscErrorCode SNESTestLocalMin(SNES snes) 3546 { 3547 PetscErrorCode ierr; 3548 PetscInt N,i,j; 3549 Vec u,uh,fh; 3550 PetscScalar value; 3551 PetscReal norm; 3552 3553 PetscFunctionBegin; 3554 ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr); 3555 ierr = VecDuplicate(u,&uh);CHKERRQ(ierr); 3556 ierr = VecDuplicate(u,&fh);CHKERRQ(ierr); 3557 3558 /* currently only works for sequential */ 3559 ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n"); 3560 ierr = VecGetSize(u,&N);CHKERRQ(ierr); 3561 for (i=0; i<N; i++) { 3562 ierr = VecCopy(u,uh);CHKERRQ(ierr); 3563 ierr = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr); 3564 for (j=-10; j<11; j++) { 3565 value = PetscSign(j)*exp(PetscAbs(j)-10.0); 3566 ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 3567 ierr = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr); 3568 ierr = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr); 3569 ierr = PetscPrintf(PETSC_COMM_WORLD," j norm %D %18.16e\n",j,norm);CHKERRQ(ierr); 3570 value = -value; 3571 ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 3572 } 3573 } 3574 ierr = VecDestroy(&uh);CHKERRQ(ierr); 3575 ierr = VecDestroy(&fh);CHKERRQ(ierr); 3576 PetscFunctionReturn(0); 3577 } 3578 3579 #undef __FUNCT__ 3580 #define __FUNCT__ "SNESKSPSetUseEW" 3581 /*@ 3582 SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for 3583 computing relative tolerance for linear solvers within an inexact 3584 Newton method. 3585 3586 Logically Collective on SNES 3587 3588 Input Parameters: 3589 + snes - SNES context 3590 - flag - PETSC_TRUE or PETSC_FALSE 3591 3592 Options Database: 3593 + -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence 3594 . -snes_ksp_ew_version ver - version of Eisenstat-Walker method 3595 . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0 3596 . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax 3597 . -snes_ksp_ew_gamma <gamma> - Sets gamma 3598 . -snes_ksp_ew_alpha <alpha> - Sets alpha 3599 . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2 3600 - -snes_ksp_ew_threshold <threshold> - Sets threshold 3601 3602 Notes: 3603 Currently, the default is to use a constant relative tolerance for 3604 the inner linear solvers. Alternatively, one can use the 3605 Eisenstat-Walker method, where the relative convergence tolerance 3606 is reset at each Newton iteration according progress of the nonlinear 3607 solver. 3608 3609 Level: advanced 3610 3611 Reference: 3612 S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an 3613 inexact Newton method", SISC 17 (1), pp.16-32, 1996. 3614 3615 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton 3616 3617 .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW() 3618 @*/ 3619 PetscErrorCode SNESKSPSetUseEW(SNES snes,PetscBool flag) 3620 { 3621 PetscFunctionBegin; 3622 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3623 PetscValidLogicalCollectiveBool(snes,flag,2); 3624 snes->ksp_ewconv = flag; 3625 PetscFunctionReturn(0); 3626 } 3627 3628 #undef __FUNCT__ 3629 #define __FUNCT__ "SNESKSPGetUseEW" 3630 /*@ 3631 SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method 3632 for computing relative tolerance for linear solvers within an 3633 inexact Newton method. 3634 3635 Not Collective 3636 3637 Input Parameter: 3638 . snes - SNES context 3639 3640 Output Parameter: 3641 . flag - PETSC_TRUE or PETSC_FALSE 3642 3643 Notes: 3644 Currently, the default is to use a constant relative tolerance for 3645 the inner linear solvers. Alternatively, one can use the 3646 Eisenstat-Walker method, where the relative convergence tolerance 3647 is reset at each Newton iteration according progress of the nonlinear 3648 solver. 3649 3650 Level: advanced 3651 3652 Reference: 3653 S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an 3654 inexact Newton method", SISC 17 (1), pp.16-32, 1996. 3655 3656 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton 3657 3658 .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW() 3659 @*/ 3660 PetscErrorCode SNESKSPGetUseEW(SNES snes, PetscBool *flag) 3661 { 3662 PetscFunctionBegin; 3663 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3664 PetscValidPointer(flag,2); 3665 *flag = snes->ksp_ewconv; 3666 PetscFunctionReturn(0); 3667 } 3668 3669 #undef __FUNCT__ 3670 #define __FUNCT__ "SNESKSPSetParametersEW" 3671 /*@ 3672 SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker 3673 convergence criteria for the linear solvers within an inexact 3674 Newton method. 3675 3676 Logically Collective on SNES 3677 3678 Input Parameters: 3679 + snes - SNES context 3680 . version - version 1, 2 (default is 2) or 3 3681 . rtol_0 - initial relative tolerance (0 <= rtol_0 < 1) 3682 . rtol_max - maximum relative tolerance (0 <= rtol_max < 1) 3683 . gamma - multiplicative factor for version 2 rtol computation 3684 (0 <= gamma2 <= 1) 3685 . alpha - power for version 2 rtol computation (1 < alpha <= 2) 3686 . alpha2 - power for safeguard 3687 - threshold - threshold for imposing safeguard (0 < threshold < 1) 3688 3689 Note: 3690 Version 3 was contributed by Luis Chacon, June 2006. 3691 3692 Use PETSC_DEFAULT to retain the default for any of the parameters. 3693 3694 Level: advanced 3695 3696 Reference: 3697 S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an 3698 inexact Newton method", Utah State University Math. Stat. Dept. Res. 3699 Report 6/94/75, June, 1994, to appear in SIAM J. Sci. Comput. 3700 3701 .keywords: SNES, KSP, Eisenstat, Walker, set, parameters 3702 3703 .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW() 3704 @*/ 3705 PetscErrorCode SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max, 3706 PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold) 3707 { 3708 SNESKSPEW *kctx; 3709 PetscFunctionBegin; 3710 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3711 kctx = (SNESKSPEW*)snes->kspconvctx; 3712 if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing"); 3713 PetscValidLogicalCollectiveInt(snes,version,2); 3714 PetscValidLogicalCollectiveReal(snes,rtol_0,3); 3715 PetscValidLogicalCollectiveReal(snes,rtol_max,4); 3716 PetscValidLogicalCollectiveReal(snes,gamma,5); 3717 PetscValidLogicalCollectiveReal(snes,alpha,6); 3718 PetscValidLogicalCollectiveReal(snes,alpha2,7); 3719 PetscValidLogicalCollectiveReal(snes,threshold,8); 3720 3721 if (version != PETSC_DEFAULT) kctx->version = version; 3722 if (rtol_0 != PETSC_DEFAULT) kctx->rtol_0 = rtol_0; 3723 if (rtol_max != PETSC_DEFAULT) kctx->rtol_max = rtol_max; 3724 if (gamma != PETSC_DEFAULT) kctx->gamma = gamma; 3725 if (alpha != PETSC_DEFAULT) kctx->alpha = alpha; 3726 if (alpha2 != PETSC_DEFAULT) kctx->alpha2 = alpha2; 3727 if (threshold != PETSC_DEFAULT) kctx->threshold = threshold; 3728 3729 if (kctx->version < 1 || kctx->version > 3) { 3730 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 and 3 are supported: %D",kctx->version); 3731 } 3732 if (kctx->rtol_0 < 0.0 || kctx->rtol_0 >= 1.0) { 3733 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_0 < 1.0: %G",kctx->rtol_0); 3734 } 3735 if (kctx->rtol_max < 0.0 || kctx->rtol_max >= 1.0) { 3736 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_max (%G) < 1.0\n",kctx->rtol_max); 3737 } 3738 if (kctx->gamma < 0.0 || kctx->gamma > 1.0) { 3739 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= gamma (%G) <= 1.0\n",kctx->gamma); 3740 } 3741 if (kctx->alpha <= 1.0 || kctx->alpha > 2.0) { 3742 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"1.0 < alpha (%G) <= 2.0\n",kctx->alpha); 3743 } 3744 if (kctx->threshold <= 0.0 || kctx->threshold >= 1.0) { 3745 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 < threshold (%G) < 1.0\n",kctx->threshold); 3746 } 3747 PetscFunctionReturn(0); 3748 } 3749 3750 #undef __FUNCT__ 3751 #define __FUNCT__ "SNESKSPGetParametersEW" 3752 /*@ 3753 SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker 3754 convergence criteria for the linear solvers within an inexact 3755 Newton method. 3756 3757 Not Collective 3758 3759 Input Parameters: 3760 snes - SNES context 3761 3762 Output Parameters: 3763 + version - version 1, 2 (default is 2) or 3 3764 . rtol_0 - initial relative tolerance (0 <= rtol_0 < 1) 3765 . rtol_max - maximum relative tolerance (0 <= rtol_max < 1) 3766 . gamma - multiplicative factor for version 2 rtol computation 3767 (0 <= gamma2 <= 1) 3768 . alpha - power for version 2 rtol computation (1 < alpha <= 2) 3769 . alpha2 - power for safeguard 3770 - threshold - threshold for imposing safeguard (0 < threshold < 1) 3771 3772 Level: advanced 3773 3774 .keywords: SNES, KSP, Eisenstat, Walker, get, parameters 3775 3776 .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW() 3777 @*/ 3778 PetscErrorCode SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max, 3779 PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold) 3780 { 3781 SNESKSPEW *kctx; 3782 PetscFunctionBegin; 3783 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3784 kctx = (SNESKSPEW*)snes->kspconvctx; 3785 if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing"); 3786 if(version) *version = kctx->version; 3787 if(rtol_0) *rtol_0 = kctx->rtol_0; 3788 if(rtol_max) *rtol_max = kctx->rtol_max; 3789 if(gamma) *gamma = kctx->gamma; 3790 if(alpha) *alpha = kctx->alpha; 3791 if(alpha2) *alpha2 = kctx->alpha2; 3792 if(threshold) *threshold = kctx->threshold; 3793 PetscFunctionReturn(0); 3794 } 3795 3796 #undef __FUNCT__ 3797 #define __FUNCT__ "SNESKSPEW_PreSolve" 3798 static PetscErrorCode SNESKSPEW_PreSolve(SNES snes, KSP ksp, Vec b, Vec x) 3799 { 3800 PetscErrorCode ierr; 3801 SNESKSPEW *kctx = (SNESKSPEW*)snes->kspconvctx; 3802 PetscReal rtol=PETSC_DEFAULT,stol; 3803 3804 PetscFunctionBegin; 3805 if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists"); 3806 if (!snes->iter) { /* first time in, so use the original user rtol */ 3807 rtol = kctx->rtol_0; 3808 } else { 3809 if (kctx->version == 1) { 3810 rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last; 3811 if (rtol < 0.0) rtol = -rtol; 3812 stol = pow(kctx->rtol_last,kctx->alpha2); 3813 if (stol > kctx->threshold) rtol = PetscMax(rtol,stol); 3814 } else if (kctx->version == 2) { 3815 rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha); 3816 stol = kctx->gamma * pow(kctx->rtol_last,kctx->alpha); 3817 if (stol > kctx->threshold) rtol = PetscMax(rtol,stol); 3818 } else if (kctx->version == 3) {/* contributed by Luis Chacon, June 2006. */ 3819 rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha); 3820 /* safeguard: avoid sharp decrease of rtol */ 3821 stol = kctx->gamma*pow(kctx->rtol_last,kctx->alpha); 3822 stol = PetscMax(rtol,stol); 3823 rtol = PetscMin(kctx->rtol_0,stol); 3824 /* safeguard: avoid oversolving */ 3825 stol = kctx->gamma*(snes->ttol)/snes->norm; 3826 stol = PetscMax(rtol,stol); 3827 rtol = PetscMin(kctx->rtol_0,stol); 3828 } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 or 3 are supported: %D",kctx->version); 3829 } 3830 /* safeguard: avoid rtol greater than one */ 3831 rtol = PetscMin(rtol,kctx->rtol_max); 3832 ierr = KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 3833 ierr = PetscInfo3(snes,"iter %D, Eisenstat-Walker (version %D) KSP rtol=%G\n",snes->iter,kctx->version,rtol);CHKERRQ(ierr); 3834 PetscFunctionReturn(0); 3835 } 3836 3837 #undef __FUNCT__ 3838 #define __FUNCT__ "SNESKSPEW_PostSolve" 3839 static PetscErrorCode SNESKSPEW_PostSolve(SNES snes, KSP ksp, Vec b, Vec x) 3840 { 3841 PetscErrorCode ierr; 3842 SNESKSPEW *kctx = (SNESKSPEW*)snes->kspconvctx; 3843 PCSide pcside; 3844 Vec lres; 3845 3846 PetscFunctionBegin; 3847 if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists"); 3848 ierr = KSPGetTolerances(ksp,&kctx->rtol_last,0,0,0);CHKERRQ(ierr); 3849 ierr = SNESGetFunctionNorm(snes,&kctx->norm_last);CHKERRQ(ierr); 3850 if (kctx->version == 1) { 3851 ierr = KSPGetPCSide(ksp,&pcside);CHKERRQ(ierr); 3852 if (pcside == PC_RIGHT) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */ 3853 /* KSP residual is true linear residual */ 3854 ierr = KSPGetResidualNorm(ksp,&kctx->lresid_last);CHKERRQ(ierr); 3855 } else { 3856 /* KSP residual is preconditioned residual */ 3857 /* compute true linear residual norm */ 3858 ierr = VecDuplicate(b,&lres);CHKERRQ(ierr); 3859 ierr = MatMult(snes->jacobian,x,lres);CHKERRQ(ierr); 3860 ierr = VecAYPX(lres,-1.0,b);CHKERRQ(ierr); 3861 ierr = VecNorm(lres,NORM_2,&kctx->lresid_last);CHKERRQ(ierr); 3862 ierr = VecDestroy(&lres);CHKERRQ(ierr); 3863 } 3864 } 3865 PetscFunctionReturn(0); 3866 } 3867 3868 #undef __FUNCT__ 3869 #define __FUNCT__ "SNES_KSPSolve" 3870 PetscErrorCode SNES_KSPSolve(SNES snes, KSP ksp, Vec b, Vec x) 3871 { 3872 PetscErrorCode ierr; 3873 3874 PetscFunctionBegin; 3875 if (snes->ksp_ewconv) { ierr = SNESKSPEW_PreSolve(snes,ksp,b,x);CHKERRQ(ierr); } 3876 ierr = KSPSolve(ksp,b,x);CHKERRQ(ierr); 3877 if (snes->ksp_ewconv) { ierr = SNESKSPEW_PostSolve(snes,ksp,b,x);CHKERRQ(ierr); } 3878 PetscFunctionReturn(0); 3879 } 3880 3881 #undef __FUNCT__ 3882 #define __FUNCT__ "SNESSetDM" 3883 /*@ 3884 SNESSetDM - Sets the DM that may be used by some preconditioners 3885 3886 Logically Collective on SNES 3887 3888 Input Parameters: 3889 + snes - the preconditioner context 3890 - dm - the dm 3891 3892 Level: intermediate 3893 3894 3895 .seealso: SNESGetDM(), KSPSetDM(), KSPGetDM() 3896 @*/ 3897 PetscErrorCode SNESSetDM(SNES snes,DM dm) 3898 { 3899 PetscErrorCode ierr; 3900 KSP ksp; 3901 3902 PetscFunctionBegin; 3903 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3904 if (dm) {ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);} 3905 ierr = DMDestroy(&snes->dm);CHKERRQ(ierr); 3906 snes->dm = dm; 3907 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 3908 ierr = KSPSetDM(ksp,dm);CHKERRQ(ierr); 3909 ierr = KSPSetDMActive(ksp,PETSC_FALSE);CHKERRQ(ierr); 3910 if (snes->pc) { 3911 ierr = SNESSetDM(snes->pc, snes->dm);CHKERRQ(ierr); 3912 } 3913 PetscFunctionReturn(0); 3914 } 3915 3916 #undef __FUNCT__ 3917 #define __FUNCT__ "SNESGetDM" 3918 /*@ 3919 SNESGetDM - Gets the DM that may be used by some preconditioners 3920 3921 Not Collective but DM obtained is parallel on SNES 3922 3923 Input Parameter: 3924 . snes - the preconditioner context 3925 3926 Output Parameter: 3927 . dm - the dm 3928 3929 Level: intermediate 3930 3931 3932 .seealso: SNESSetDM(), KSPSetDM(), KSPGetDM() 3933 @*/ 3934 PetscErrorCode SNESGetDM(SNES snes,DM *dm) 3935 { 3936 PetscFunctionBegin; 3937 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 3938 *dm = snes->dm; 3939 PetscFunctionReturn(0); 3940 } 3941 3942 #undef __FUNCT__ 3943 #define __FUNCT__ "SNESSetPC" 3944 /*@ 3945 SNESSetPC - Sets the nonlinear preconditioner to be used. 3946 3947 Collective on SNES 3948 3949 Input Parameters: 3950 + snes - iterative context obtained from SNESCreate() 3951 - pc - the preconditioner object 3952 3953 Notes: 3954 Use SNESGetPC() to retrieve the preconditioner context (for example, 3955 to configure it using the API). 3956 3957 Level: developer 3958 3959 .keywords: SNES, set, precondition 3960 .seealso: SNESGetPC() 3961 @*/ 3962 PetscErrorCode SNESSetPC(SNES snes, SNES pc) 3963 { 3964 PetscErrorCode ierr; 3965 3966 PetscFunctionBegin; 3967 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 3968 PetscValidHeaderSpecific(pc, SNES_CLASSID, 2); 3969 PetscCheckSameComm(snes, 1, pc, 2); 3970 ierr = PetscObjectReference((PetscObject) pc);CHKERRQ(ierr); 3971 ierr = SNESDestroy(&snes->pc);CHKERRQ(ierr); 3972 snes->pc = pc; 3973 ierr = PetscLogObjectParent(snes, snes->pc);CHKERRQ(ierr); 3974 PetscFunctionReturn(0); 3975 } 3976 3977 #undef __FUNCT__ 3978 #define __FUNCT__ "SNESGetPC" 3979 /*@ 3980 SNESGetPC - Returns a pointer to the nonlinear preconditioning context set with SNESSetPC(). 3981 3982 Not Collective 3983 3984 Input Parameter: 3985 . snes - iterative context obtained from SNESCreate() 3986 3987 Output Parameter: 3988 . pc - preconditioner context 3989 3990 Level: developer 3991 3992 .keywords: SNES, get, preconditioner 3993 .seealso: SNESSetPC() 3994 @*/ 3995 PetscErrorCode SNESGetPC(SNES snes, SNES *pc) 3996 { 3997 PetscErrorCode ierr; 3998 3999 PetscFunctionBegin; 4000 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 4001 PetscValidPointer(pc, 2); 4002 if (!snes->pc) { 4003 ierr = SNESCreate(((PetscObject) snes)->comm, &snes->pc);CHKERRQ(ierr); 4004 ierr = PetscObjectIncrementTabLevel((PetscObject) snes->pc, (PetscObject) snes, 1);CHKERRQ(ierr); 4005 ierr = PetscLogObjectParent(snes, snes->pc);CHKERRQ(ierr); 4006 } 4007 *pc = snes->pc; 4008 PetscFunctionReturn(0); 4009 } 4010 4011 #if defined(PETSC_HAVE_MATLAB_ENGINE) 4012 #include <mex.h> 4013 4014 typedef struct {char *funcname; mxArray *ctx;} SNESMatlabContext; 4015 4016 #undef __FUNCT__ 4017 #define __FUNCT__ "SNESComputeFunction_Matlab" 4018 /* 4019 SNESComputeFunction_Matlab - Calls the function that has been set with 4020 SNESSetFunctionMatlab(). 4021 4022 Collective on SNES 4023 4024 Input Parameters: 4025 + snes - the SNES context 4026 - x - input vector 4027 4028 Output Parameter: 4029 . y - function vector, as set by SNESSetFunction() 4030 4031 Notes: 4032 SNESComputeFunction() is typically used within nonlinear solvers 4033 implementations, so most users would not generally call this routine 4034 themselves. 4035 4036 Level: developer 4037 4038 .keywords: SNES, nonlinear, compute, function 4039 4040 .seealso: SNESSetFunction(), SNESGetFunction() 4041 */ 4042 PetscErrorCode SNESComputeFunction_Matlab(SNES snes,Vec x,Vec y, void *ctx) 4043 { 4044 PetscErrorCode ierr; 4045 SNESMatlabContext *sctx = (SNESMatlabContext *)ctx; 4046 int nlhs = 1,nrhs = 5; 4047 mxArray *plhs[1],*prhs[5]; 4048 long long int lx = 0,ly = 0,ls = 0; 4049 4050 PetscFunctionBegin; 4051 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 4052 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 4053 PetscValidHeaderSpecific(y,VEC_CLASSID,3); 4054 PetscCheckSameComm(snes,1,x,2); 4055 PetscCheckSameComm(snes,1,y,3); 4056 4057 /* call Matlab function in ctx with arguments x and y */ 4058 4059 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 4060 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 4061 ierr = PetscMemcpy(&ly,&y,sizeof(x));CHKERRQ(ierr); 4062 prhs[0] = mxCreateDoubleScalar((double)ls); 4063 prhs[1] = mxCreateDoubleScalar((double)lx); 4064 prhs[2] = mxCreateDoubleScalar((double)ly); 4065 prhs[3] = mxCreateString(sctx->funcname); 4066 prhs[4] = sctx->ctx; 4067 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeFunctionInternal");CHKERRQ(ierr); 4068 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 4069 mxDestroyArray(prhs[0]); 4070 mxDestroyArray(prhs[1]); 4071 mxDestroyArray(prhs[2]); 4072 mxDestroyArray(prhs[3]); 4073 mxDestroyArray(plhs[0]); 4074 PetscFunctionReturn(0); 4075 } 4076 4077 4078 #undef __FUNCT__ 4079 #define __FUNCT__ "SNESSetFunctionMatlab" 4080 /* 4081 SNESSetFunctionMatlab - Sets the function evaluation routine and function 4082 vector for use by the SNES routines in solving systems of nonlinear 4083 equations from MATLAB. Here the function is a string containing the name of a MATLAB function 4084 4085 Logically Collective on SNES 4086 4087 Input Parameters: 4088 + snes - the SNES context 4089 . r - vector to store function value 4090 - func - function evaluation routine 4091 4092 Calling sequence of func: 4093 $ func (SNES snes,Vec x,Vec f,void *ctx); 4094 4095 4096 Notes: 4097 The Newton-like methods typically solve linear systems of the form 4098 $ f'(x) x = -f(x), 4099 where f'(x) denotes the Jacobian matrix and f(x) is the function. 4100 4101 Level: beginner 4102 4103 .keywords: SNES, nonlinear, set, function 4104 4105 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction() 4106 */ 4107 PetscErrorCode SNESSetFunctionMatlab(SNES snes,Vec r,const char *func,mxArray *ctx) 4108 { 4109 PetscErrorCode ierr; 4110 SNESMatlabContext *sctx; 4111 4112 PetscFunctionBegin; 4113 /* currently sctx is memory bleed */ 4114 ierr = PetscMalloc(sizeof(SNESMatlabContext),&sctx);CHKERRQ(ierr); 4115 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 4116 /* 4117 This should work, but it doesn't 4118 sctx->ctx = ctx; 4119 mexMakeArrayPersistent(sctx->ctx); 4120 */ 4121 sctx->ctx = mxDuplicateArray(ctx); 4122 ierr = SNESSetFunction(snes,r,SNESComputeFunction_Matlab,sctx);CHKERRQ(ierr); 4123 PetscFunctionReturn(0); 4124 } 4125 4126 #undef __FUNCT__ 4127 #define __FUNCT__ "SNESComputeJacobian_Matlab" 4128 /* 4129 SNESComputeJacobian_Matlab - Calls the function that has been set with 4130 SNESSetJacobianMatlab(). 4131 4132 Collective on SNES 4133 4134 Input Parameters: 4135 + snes - the SNES context 4136 . x - input vector 4137 . A, B - the matrices 4138 - ctx - user context 4139 4140 Output Parameter: 4141 . flag - structure of the matrix 4142 4143 Level: developer 4144 4145 .keywords: SNES, nonlinear, compute, function 4146 4147 .seealso: SNESSetFunction(), SNESGetFunction() 4148 @*/ 4149 PetscErrorCode SNESComputeJacobian_Matlab(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag, void *ctx) 4150 { 4151 PetscErrorCode ierr; 4152 SNESMatlabContext *sctx = (SNESMatlabContext *)ctx; 4153 int nlhs = 2,nrhs = 6; 4154 mxArray *plhs[2],*prhs[6]; 4155 long long int lx = 0,lA = 0,ls = 0, lB = 0; 4156 4157 PetscFunctionBegin; 4158 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 4159 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 4160 4161 /* call Matlab function in ctx with arguments x and y */ 4162 4163 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 4164 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 4165 ierr = PetscMemcpy(&lA,A,sizeof(x));CHKERRQ(ierr); 4166 ierr = PetscMemcpy(&lB,B,sizeof(x));CHKERRQ(ierr); 4167 prhs[0] = mxCreateDoubleScalar((double)ls); 4168 prhs[1] = mxCreateDoubleScalar((double)lx); 4169 prhs[2] = mxCreateDoubleScalar((double)lA); 4170 prhs[3] = mxCreateDoubleScalar((double)lB); 4171 prhs[4] = mxCreateString(sctx->funcname); 4172 prhs[5] = sctx->ctx; 4173 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeJacobianInternal");CHKERRQ(ierr); 4174 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 4175 *flag = (MatStructure) mxGetScalar(plhs[1]);CHKERRQ(ierr); 4176 mxDestroyArray(prhs[0]); 4177 mxDestroyArray(prhs[1]); 4178 mxDestroyArray(prhs[2]); 4179 mxDestroyArray(prhs[3]); 4180 mxDestroyArray(prhs[4]); 4181 mxDestroyArray(plhs[0]); 4182 mxDestroyArray(plhs[1]); 4183 PetscFunctionReturn(0); 4184 } 4185 4186 4187 #undef __FUNCT__ 4188 #define __FUNCT__ "SNESSetJacobianMatlab" 4189 /* 4190 SNESSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices 4191 vector for use by the SNES routines in solving systems of nonlinear 4192 equations from MATLAB. Here the function is a string containing the name of a MATLAB function 4193 4194 Logically Collective on SNES 4195 4196 Input Parameters: 4197 + snes - the SNES context 4198 . A,B - Jacobian matrices 4199 . func - function evaluation routine 4200 - ctx - user context 4201 4202 Calling sequence of func: 4203 $ flag = func (SNES snes,Vec x,Mat A,Mat B,void *ctx); 4204 4205 4206 Level: developer 4207 4208 .keywords: SNES, nonlinear, set, function 4209 4210 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction() 4211 */ 4212 PetscErrorCode SNESSetJacobianMatlab(SNES snes,Mat A,Mat B,const char *func,mxArray *ctx) 4213 { 4214 PetscErrorCode ierr; 4215 SNESMatlabContext *sctx; 4216 4217 PetscFunctionBegin; 4218 /* currently sctx is memory bleed */ 4219 ierr = PetscMalloc(sizeof(SNESMatlabContext),&sctx);CHKERRQ(ierr); 4220 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 4221 /* 4222 This should work, but it doesn't 4223 sctx->ctx = ctx; 4224 mexMakeArrayPersistent(sctx->ctx); 4225 */ 4226 sctx->ctx = mxDuplicateArray(ctx); 4227 ierr = SNESSetJacobian(snes,A,B,SNESComputeJacobian_Matlab,sctx);CHKERRQ(ierr); 4228 PetscFunctionReturn(0); 4229 } 4230 4231 #undef __FUNCT__ 4232 #define __FUNCT__ "SNESMonitor_Matlab" 4233 /* 4234 SNESMonitor_Matlab - Calls the function that has been set with SNESMonitorSetMatlab(). 4235 4236 Collective on SNES 4237 4238 .seealso: SNESSetFunction(), SNESGetFunction() 4239 @*/ 4240 PetscErrorCode SNESMonitor_Matlab(SNES snes,PetscInt it, PetscReal fnorm, void *ctx) 4241 { 4242 PetscErrorCode ierr; 4243 SNESMatlabContext *sctx = (SNESMatlabContext *)ctx; 4244 int nlhs = 1,nrhs = 6; 4245 mxArray *plhs[1],*prhs[6]; 4246 long long int lx = 0,ls = 0; 4247 Vec x=snes->vec_sol; 4248 4249 PetscFunctionBegin; 4250 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 4251 4252 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 4253 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 4254 prhs[0] = mxCreateDoubleScalar((double)ls); 4255 prhs[1] = mxCreateDoubleScalar((double)it); 4256 prhs[2] = mxCreateDoubleScalar((double)fnorm); 4257 prhs[3] = mxCreateDoubleScalar((double)lx); 4258 prhs[4] = mxCreateString(sctx->funcname); 4259 prhs[5] = sctx->ctx; 4260 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESMonitorInternal");CHKERRQ(ierr); 4261 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 4262 mxDestroyArray(prhs[0]); 4263 mxDestroyArray(prhs[1]); 4264 mxDestroyArray(prhs[2]); 4265 mxDestroyArray(prhs[3]); 4266 mxDestroyArray(prhs[4]); 4267 mxDestroyArray(plhs[0]); 4268 PetscFunctionReturn(0); 4269 } 4270 4271 4272 #undef __FUNCT__ 4273 #define __FUNCT__ "SNESMonitorSetMatlab" 4274 /* 4275 SNESMonitorSetMatlab - Sets the monitor function from MATLAB 4276 4277 Level: developer 4278 4279 .keywords: SNES, nonlinear, set, function 4280 4281 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction() 4282 */ 4283 PetscErrorCode SNESMonitorSetMatlab(SNES snes,const char *func,mxArray *ctx) 4284 { 4285 PetscErrorCode ierr; 4286 SNESMatlabContext *sctx; 4287 4288 PetscFunctionBegin; 4289 /* currently sctx is memory bleed */ 4290 ierr = PetscMalloc(sizeof(SNESMatlabContext),&sctx);CHKERRQ(ierr); 4291 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 4292 /* 4293 This should work, but it doesn't 4294 sctx->ctx = ctx; 4295 mexMakeArrayPersistent(sctx->ctx); 4296 */ 4297 sctx->ctx = mxDuplicateArray(ctx); 4298 ierr = SNESMonitorSet(snes,SNESMonitor_Matlab,sctx,PETSC_NULL);CHKERRQ(ierr); 4299 PetscFunctionReturn(0); 4300 } 4301 4302 #endif 4303