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