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