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