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