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