1 /* Defines the basic SNES object */ 2 #include <../src/snes/impls/fas/fasimpls.h> /*I "petscsnesfas.h" I*/ 3 4 const char *const SNESFASTypes[] = {"MULTIPLICATIVE","ADDITIVE","SNESFASType","SNES_FAS",0}; 5 6 extern PetscErrorCode SNESDestroy_FAS(SNES snes); 7 extern PetscErrorCode SNESSetUp_FAS(SNES snes); 8 extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes); 9 extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer); 10 extern PetscErrorCode SNESSolve_FAS(SNES snes); 11 extern PetscErrorCode SNESReset_FAS(SNES snes); 12 extern PetscErrorCode SNESFASGalerkinDefaultFunction(SNES, Vec, Vec, void *); 13 extern PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES, SNES *); 14 15 /*MC 16 17 SNESFAS - Full Approximation Scheme nonlinear multigrid solver. 18 19 The nonlinear problem is solved by correction using coarse versions 20 of the nonlinear problem. This problem is perturbed so that a projected 21 solution of the fine problem elicits no correction from the coarse problem. 22 23 Options Database: 24 + -snes_fas_levels - The number of levels 25 . -snes_fas_cycles<1> - The number of cycles -- 1 for V, 2 for W 26 . -snes_fas_type<additive, multiplicative> - Additive or multiplicative cycle 27 . -snes_fas_galerkin<PETSC_FALSE> - Form coarse problems by projection back upon the fine problem 28 . -snes_fas_smoothup<1> - The number of iterations of the post-smoother 29 . -snes_fas_smoothdown<1> - The number of iterations of the pre-smoother 30 . -snes_fas_monitor - Monitor progress of all of the levels 31 . -fas_levels_snes_ - SNES options for all smoothers 32 . -fas_levels_cycle_snes_ - SNES options for all cycles 33 . -fas_levels_i_snes_ - SNES options for the smoothers on level i 34 . -fas_levels_i_cycle_snes_ - SNES options for the cycle on level i 35 - -fas_coarse_snes_ - SNES options for the coarsest smoother 36 37 Notes: 38 The organization of the FAS solver is slightly different from the organization of PCMG 39 As each level has smoother SNES instances(down and potentially up) and a cycle SNES instance. 40 The cycle SNES instance may be used for monitoring convergence on a particular level. 41 42 Level: beginner 43 44 .seealso: PCMG, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 45 M*/ 46 47 #undef __FUNCT__ 48 #define __FUNCT__ "SNESCreate_FAS" 49 PETSC_EXTERN_C PetscErrorCode SNESCreate_FAS(SNES snes) 50 { 51 SNES_FAS * fas; 52 PetscErrorCode ierr; 53 54 PetscFunctionBegin; 55 snes->ops->destroy = SNESDestroy_FAS; 56 snes->ops->setup = SNESSetUp_FAS; 57 snes->ops->setfromoptions = SNESSetFromOptions_FAS; 58 snes->ops->view = SNESView_FAS; 59 snes->ops->solve = SNESSolve_FAS; 60 snes->ops->reset = SNESReset_FAS; 61 62 snes->usesksp = PETSC_FALSE; 63 snes->usespc = PETSC_FALSE; 64 65 if (!snes->tolerancesset) { 66 snes->max_funcs = 30000; 67 snes->max_its = 10000; 68 } 69 70 ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr); 71 snes->data = (void*) fas; 72 fas->level = 0; 73 fas->levels = 1; 74 fas->n_cycles = 1; 75 fas->max_up_it = 1; 76 fas->max_down_it = 1; 77 fas->smoothu = PETSC_NULL; 78 fas->smoothd = PETSC_NULL; 79 fas->next = PETSC_NULL; 80 fas->previous = PETSC_NULL; 81 fas->fine = snes; 82 fas->interpolate = PETSC_NULL; 83 fas->restrct = PETSC_NULL; 84 fas->inject = PETSC_NULL; 85 fas->monitor = PETSC_NULL; 86 fas->usedmfornumberoflevels = PETSC_FALSE; 87 fas->fastype = SNES_FAS_MULTIPLICATIVE; 88 89 fas->eventsmoothsetup = PETSC_FALSE; 90 fas->eventsmoothsolve = PETSC_FALSE; 91 fas->eventresidual = PETSC_FALSE; 92 fas->eventinterprestrict = PETSC_FALSE; 93 PetscFunctionReturn(0); 94 } 95 96 #undef __FUNCT__ 97 #define __FUNCT__ "SNESReset_FAS" 98 PetscErrorCode SNESReset_FAS(SNES snes) 99 { 100 PetscErrorCode ierr = 0; 101 SNES_FAS * fas = (SNES_FAS *)snes->data; 102 103 PetscFunctionBegin; 104 ierr = SNESDestroy(&fas->smoothu);CHKERRQ(ierr); 105 ierr = SNESDestroy(&fas->smoothd);CHKERRQ(ierr); 106 ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 107 ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 108 ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 109 ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 110 if (fas->next) ierr = SNESReset(fas->next);CHKERRQ(ierr); 111 PetscFunctionReturn(0); 112 } 113 114 #undef __FUNCT__ 115 #define __FUNCT__ "SNESDestroy_FAS" 116 PetscErrorCode SNESDestroy_FAS(SNES snes) 117 { 118 SNES_FAS * fas = (SNES_FAS *)snes->data; 119 PetscErrorCode ierr = 0; 120 121 PetscFunctionBegin; 122 /* recursively resets and then destroys */ 123 ierr = SNESReset(snes);CHKERRQ(ierr); 124 if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 125 ierr = PetscFree(fas);CHKERRQ(ierr); 126 PetscFunctionReturn(0); 127 } 128 129 #undef __FUNCT__ 130 #define __FUNCT__ "SNESSetUp_FAS" 131 PetscErrorCode SNESSetUp_FAS(SNES snes) 132 { 133 SNES_FAS *fas = (SNES_FAS *) snes->data; 134 PetscErrorCode ierr; 135 VecScatter injscatter; 136 PetscInt dm_levels; 137 Vec vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */ 138 SNES next; 139 PetscBool isFine; 140 SNESLineSearch linesearch; 141 SNESLineSearch slinesearch; 142 void *lsprectx,*lspostctx; 143 PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*); 144 PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool *,PetscBool *,void*); 145 146 PetscFunctionBegin; 147 ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 148 if (fas->usedmfornumberoflevels && isFine) { 149 ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr); 150 dm_levels++; 151 if (dm_levels > fas->levels) { 152 /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESFASSetLevels_FAS*/ 153 vec_sol = snes->vec_sol; 154 vec_func = snes->vec_func; 155 vec_sol_update = snes->vec_sol_update; 156 vec_rhs = snes->vec_rhs; 157 snes->vec_sol = PETSC_NULL; 158 snes->vec_func = PETSC_NULL; 159 snes->vec_sol_update = PETSC_NULL; 160 snes->vec_rhs = PETSC_NULL; 161 162 /* reset the number of levels */ 163 ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr); 164 ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 165 166 snes->vec_sol = vec_sol; 167 snes->vec_func = vec_func; 168 snes->vec_rhs = vec_rhs; 169 snes->vec_sol_update = vec_sol_update; 170 } 171 } 172 ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 173 if (!isFine) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */ 174 175 if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 176 ierr = SNESDefaultGetWork(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 177 } else { 178 ierr = SNESDefaultGetWork(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 179 } 180 181 /* set up the smoothers if they haven't already been set up */ 182 if (!fas->smoothd) { 183 ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); 184 } 185 186 if (snes->dm) { 187 /* set the smoother DMs properly */ 188 if (fas->smoothu) ierr = SNESSetDM(fas->smoothu, snes->dm);CHKERRQ(ierr); 189 ierr = SNESSetDM(fas->smoothd, snes->dm);CHKERRQ(ierr); 190 /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 191 if (next) { 192 /* for now -- assume the DM and the evaluation functions have been set externally */ 193 if (!next->dm) { 194 ierr = DMCoarsen(snes->dm, ((PetscObject)next)->comm, &next->dm);CHKERRQ(ierr); 195 ierr = SNESSetDM(next, next->dm);CHKERRQ(ierr); 196 } 197 /* set the interpolation and restriction from the DM */ 198 if (!fas->interpolate) { 199 ierr = DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 200 if (!fas->restrct) { 201 ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr); 202 fas->restrct = fas->interpolate; 203 } 204 } 205 /* set the injection from the DM */ 206 if (!fas->inject) { 207 ierr = DMCreateInjection(next->dm, snes->dm, &injscatter);CHKERRQ(ierr); 208 ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr); 209 ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr); 210 } 211 } 212 } 213 /*pass the smoother, function, and jacobian up to the next level if it's not user set already */ 214 if (fas->galerkin) { 215 if (next) 216 ierr = SNESSetFunction(next, PETSC_NULL, SNESFASGalerkinDefaultFunction, next);CHKERRQ(ierr); 217 if (fas->smoothd && fas->level != fas->levels - 1) ierr = SNESSetFunction(fas->smoothd, PETSC_NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr); 218 if (fas->smoothu && fas->level != fas->levels - 1) ierr = SNESSetFunction(fas->smoothu, PETSC_NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr); 219 } 220 221 /* sets the down (pre) smoother's default norm and sets it from options */ 222 if (fas->smoothd) { 223 if (fas->level == 0 && fas->levels != 1) { 224 ierr = SNESSetNormType(fas->smoothd, SNES_NORM_NONE);CHKERRQ(ierr); 225 } else { 226 ierr = SNESSetNormType(fas->smoothd, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr); 227 } 228 ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)fas->smoothd);CHKERRQ(ierr); 229 ierr = SNESSetFromOptions(fas->smoothd);CHKERRQ(ierr); 230 ierr = SNESGetSNESLineSearch(snes,&linesearch);CHKERRQ(ierr); 231 ierr = SNESGetSNESLineSearch(fas->smoothd,&slinesearch);CHKERRQ(ierr); 232 ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr); 233 ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr); 234 ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr); 235 ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr); 236 ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr); 237 238 fas->smoothd->vec_sol = snes->vec_sol; 239 ierr = PetscObjectReference((PetscObject)snes->vec_sol);CHKERRQ(ierr); 240 fas->smoothd->vec_sol_update = snes->vec_sol_update; 241 ierr = PetscObjectReference((PetscObject)snes->vec_sol_update);CHKERRQ(ierr); 242 fas->smoothd->vec_func = snes->vec_func; 243 ierr = PetscObjectReference((PetscObject)snes->vec_func);CHKERRQ(ierr); 244 245 if (fas->eventsmoothsetup) {ierr = PetscLogEventBegin(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 246 ierr = SNESSetUp(fas->smoothd);CHKERRQ(ierr); 247 if (fas->eventsmoothsetup) {ierr = PetscLogEventEnd(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 248 } 249 250 /* sets the up (post) smoother's default norm and sets it from options */ 251 if (fas->smoothu) { 252 if (fas->level != fas->levels - 1) { 253 ierr = SNESSetNormType(fas->smoothu, SNES_NORM_NONE);CHKERRQ(ierr); 254 } else { 255 ierr = SNESSetNormType(fas->smoothu, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr); 256 } 257 ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)fas->smoothu);CHKERRQ(ierr); 258 ierr = SNESSetFromOptions(fas->smoothu);CHKERRQ(ierr); 259 ierr = SNESGetSNESLineSearch(snes,&linesearch);CHKERRQ(ierr); 260 ierr = SNESGetSNESLineSearch(fas->smoothu,&slinesearch);CHKERRQ(ierr); 261 ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr); 262 ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr); 263 ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr); 264 ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr); 265 ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr); 266 267 fas->smoothu->vec_sol = snes->vec_sol; 268 ierr = PetscObjectReference((PetscObject)snes->vec_sol);CHKERRQ(ierr); 269 fas->smoothu->vec_sol_update = snes->vec_sol_update; 270 ierr = PetscObjectReference((PetscObject)snes->vec_sol_update);CHKERRQ(ierr); 271 fas->smoothu->vec_func = snes->vec_func; 272 ierr = PetscObjectReference((PetscObject)snes->vec_func);CHKERRQ(ierr); 273 274 if (fas->eventsmoothsetup) {ierr = PetscLogEventBegin(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 275 ierr = SNESSetUp(fas->smoothu);CHKERRQ(ierr); 276 if (fas->eventsmoothsetup) {ierr = PetscLogEventEnd(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 277 278 } 279 280 if (next) { 281 /* gotta set up the solution vector for this to work */ 282 if (!next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_sol);CHKERRQ(ierr);} 283 if (!next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_rhs);CHKERRQ(ierr);} 284 ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)next);CHKERRQ(ierr); 285 ierr = SNESGetSNESLineSearch(snes,&linesearch);CHKERRQ(ierr); 286 ierr = SNESGetSNESLineSearch(fas->next,&slinesearch);CHKERRQ(ierr); 287 ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr); 288 ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr); 289 ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr); 290 ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr); 291 ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr); 292 ierr = SNESSetUp(next);CHKERRQ(ierr); 293 } 294 /* setup FAS work vectors */ 295 if (fas->galerkin) { 296 ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr); 297 ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr); 298 } 299 PetscFunctionReturn(0); 300 } 301 302 #undef __FUNCT__ 303 #define __FUNCT__ "SNESSetFromOptions_FAS" 304 PetscErrorCode SNESSetFromOptions_FAS(SNES snes) 305 { 306 SNES_FAS *fas = (SNES_FAS *) snes->data; 307 PetscInt levels = 1; 308 PetscBool flg = PETSC_FALSE, upflg = PETSC_FALSE, downflg = PETSC_FALSE, monflg = PETSC_FALSE, galerkinflg = PETSC_FALSE; 309 PetscErrorCode ierr; 310 char monfilename[PETSC_MAX_PATH_LEN]; 311 SNESFASType fastype; 312 const char *optionsprefix; 313 SNESLineSearch linesearch; 314 PetscInt m, n_up, n_down; 315 SNES next; 316 PetscBool isFine; 317 318 PetscFunctionBegin; 319 ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 320 ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr); 321 322 /* number of levels -- only process most options on the finest level */ 323 if (isFine) { 324 ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 325 if (!flg && snes->dm) { 326 ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 327 levels++; 328 fas->usedmfornumberoflevels = PETSC_TRUE; 329 } 330 ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr); 331 fastype = fas->fastype; 332 ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr); 333 if (flg) { 334 ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr); 335 } 336 337 ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 338 ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&m,&flg);CHKERRQ(ierr); 339 if (flg) { 340 ierr = SNESFASSetCycles(snes, m);CHKERRQ(ierr); 341 } 342 343 ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFASSetGalerkin",fas->galerkin,&galerkinflg,&flg);CHKERRQ(ierr); 344 if (flg) { 345 ierr = SNESFASSetGalerkin(snes, galerkinflg);CHKERRQ(ierr); 346 } 347 348 ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",fas->max_up_it,&n_up,&upflg);CHKERRQ(ierr); 349 350 ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",fas->max_down_it,&n_down,&downflg);CHKERRQ(ierr); 351 352 ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESFASSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 353 if (monflg) ierr = SNESFASSetMonitor(snes, PETSC_TRUE);CHKERRQ(ierr); 354 355 flg = PETSC_FALSE; 356 monflg = PETSC_TRUE; 357 ierr = PetscOptionsBool("-snes_fas_log","Log times for each FAS level","SNESFASSetLog",monflg,&monflg,&flg);CHKERRQ(ierr); 358 if (flg) {ierr = SNESFASSetLog(snes,monflg);CHKERRQ(ierr);} 359 } 360 361 ierr = PetscOptionsTail();CHKERRQ(ierr); 362 /* setup from the determined types if there is no pointwise procedure or smoother defined */ 363 if (upflg) { 364 ierr = SNESFASSetNumberSmoothUp(snes,n_up);CHKERRQ(ierr); 365 } 366 if (downflg) { 367 ierr = SNESFASSetNumberSmoothDown(snes,n_down);CHKERRQ(ierr); 368 } 369 370 /* set up the default line search for coarse grid corrections */ 371 if (fas->fastype == SNES_FAS_ADDITIVE) { 372 if (!snes->linesearch) { 373 ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr); 374 ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr); 375 } 376 } 377 378 ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 379 /* recursive option setting for the smoothers */ 380 if (next) {ierr = SNESSetFromOptions(next);CHKERRQ(ierr);} 381 PetscFunctionReturn(0); 382 } 383 384 #undef __FUNCT__ 385 #define __FUNCT__ "SNESView_FAS" 386 PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 387 { 388 SNES_FAS *fas = (SNES_FAS *) snes->data; 389 PetscBool isFine,iascii,isdraw; 390 PetscInt i; 391 PetscErrorCode ierr; 392 SNES smoothu, smoothd, levelsnes; 393 394 PetscFunctionBegin; 395 ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 396 if (isFine) { 397 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 398 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 399 if (iascii) { 400 ierr = PetscViewerASCIIPrintf(viewer, "FAS: type is %s, levels=%D, cycles=%D\n", SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr); 401 if (fas->galerkin) { 402 ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr); 403 } else { 404 ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr); 405 } 406 for (i=0; i<fas->levels; i++) { 407 ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr); 408 ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr); 409 ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr); 410 if (!i) { 411 ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr); 412 } else { 413 ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 414 } 415 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 416 ierr = SNESView(smoothd,viewer);CHKERRQ(ierr); 417 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 418 if (i && (smoothd == smoothu)) { 419 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 420 } else if (i) { 421 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 422 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 423 ierr = SNESView(smoothu,viewer);CHKERRQ(ierr); 424 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 425 } 426 } 427 } else if (isdraw) { 428 PetscDraw draw; 429 PetscReal x,w,y,bottom,th,wth; 430 SNES_FAS *curfas = fas; 431 ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 432 ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 433 ierr = PetscDrawStringGetSize(draw,&wth,&th);CHKERRQ(ierr); 434 bottom = y - th; 435 while (curfas) { 436 if (!curfas->smoothu) { 437 ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr); 438 if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr); 439 ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 440 } else { 441 w = 0.5*PetscMin(1.0-x,x); 442 ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr); 443 if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr); 444 ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 445 ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr); 446 if (curfas->smoothu) ierr = SNESView(curfas->smoothu,viewer);CHKERRQ(ierr); 447 ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 448 } 449 /* this is totally bogus but we have no way of knowing how low the previous one was draw to */ 450 bottom -= 5*th; 451 if (curfas->next) { 452 curfas = (SNES_FAS*)curfas->next->data; 453 } else { 454 curfas = PETSC_NULL; 455 } 456 } 457 } 458 } 459 PetscFunctionReturn(0); 460 } 461 462 #undef __FUNCT__ 463 #define __FUNCT__ "SNESFASDownSmooth_Private" 464 /* 465 Defines the action of the downsmoother 466 */ 467 PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) 468 { 469 PetscErrorCode ierr = 0; 470 SNESConvergedReason reason; 471 Vec FPC; 472 SNES smoothd; 473 SNES_FAS *fas = (SNES_FAS *) snes->data; 474 475 PetscFunctionBegin; 476 ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr); 477 ierr = SNESSetInitialFunction(smoothd, F);CHKERRQ(ierr); 478 ierr = SNESSetInitialFunctionNorm(smoothd, *fnorm);CHKERRQ(ierr); 479 if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 480 ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr); 481 if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 482 /* check convergence reason for the smoother */ 483 ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr); 484 if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN)) { 485 snes->reason = SNES_DIVERGED_INNER; 486 PetscFunctionReturn(0); 487 } 488 ierr = SNESGetFunction(smoothd, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 489 ierr = VecCopy(FPC, F);CHKERRQ(ierr); 490 ierr = SNESGetFunctionNorm(smoothd, fnorm);CHKERRQ(ierr); 491 PetscFunctionReturn(0); 492 } 493 494 495 #undef __FUNCT__ 496 #define __FUNCT__ "SNESFASUpSmooth_Private" 497 /* 498 Defines the action of the upsmoother 499 */ 500 PetscErrorCode SNESFASUpSmooth_Private (SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) 501 { 502 PetscErrorCode ierr = 0; 503 SNESConvergedReason reason; 504 Vec FPC; 505 SNES smoothu; 506 SNES_FAS *fas = (SNES_FAS *) snes->data; 507 508 PetscFunctionBegin; 509 ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr); 510 if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 511 ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr); 512 if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 513 /* check convergence reason for the smoother */ 514 ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr); 515 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 516 snes->reason = SNES_DIVERGED_INNER; 517 PetscFunctionReturn(0); 518 } 519 ierr = SNESGetFunction(smoothu, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 520 ierr = VecCopy(FPC, F);CHKERRQ(ierr); 521 ierr = SNESGetFunctionNorm(smoothu, fnorm);CHKERRQ(ierr); 522 PetscFunctionReturn(0); 523 } 524 525 #undef __FUNCT__ 526 #define __FUNCT__ "SNESFASCreateCoarseVec" 527 /*@ 528 SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level 529 530 Collective 531 532 Input Arguments: 533 . snes - SNESFAS 534 535 Output Arguments: 536 . Xcoarse - vector on level one coarser than snes 537 538 Level: developer 539 540 .seealso: SNESFASSetRestriction(), SNESFASRestrict() 541 @*/ 542 PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse) 543 { 544 PetscErrorCode ierr; 545 SNES_FAS *fas = (SNES_FAS*)snes->data; 546 547 PetscFunctionBegin; 548 if (fas->rscale) {ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);} 549 else if (fas->inject) {ierr = MatGetVecs(fas->inject,Xcoarse,PETSC_NULL);CHKERRQ(ierr);} 550 else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr); 551 PetscFunctionReturn(0); 552 } 553 554 #undef __FUNCT__ 555 #define __FUNCT__ "SNESFASRestrict" 556 /*@ 557 SNESFASRestrict - restrict a Vec to the next coarser level 558 559 Collective 560 561 Input Arguments: 562 + fine - SNES from which to restrict 563 - Xfine - vector to restrict 564 565 Output Arguments: 566 . Xcoarse - result of restriction 567 568 Level: developer 569 570 .seealso: SNESFASSetRestriction(), SNESFASSetInjection() 571 @*/ 572 PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse) 573 { 574 PetscErrorCode ierr; 575 SNES_FAS *fas = (SNES_FAS*)fine->data; 576 577 PetscFunctionBegin; 578 PetscValidHeaderSpecific(fine,SNES_CLASSID,1); 579 PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2); 580 PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3); 581 if (fas->inject) { 582 ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr); 583 } else { 584 ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr); 585 ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr); 586 } 587 PetscFunctionReturn(0); 588 } 589 590 #undef __FUNCT__ 591 #define __FUNCT__ "SNESFASCoarseCorrection" 592 /* 593 594 Performs the FAS coarse correction as: 595 596 fine problem: F(x) = 0 597 coarse problem: F^c(x) = b^c 598 599 b^c = F^c(I^c_fx^f - I^c_fF(x)) 600 601 */ 602 PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) 603 { 604 PetscErrorCode ierr; 605 Vec X_c, Xo_c, F_c, B_c; 606 SNESConvergedReason reason; 607 SNES next; 608 Mat restrct, interpolate; 609 SNES_FAS *fasc; 610 611 PetscFunctionBegin; 612 ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 613 if (next) { 614 fasc = (SNES_FAS *)next->data; 615 616 ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 617 ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 618 619 X_c = next->vec_sol; 620 Xo_c = next->work[0]; 621 F_c = next->vec_func; 622 B_c = next->vec_rhs; 623 624 if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 625 ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr); 626 /* restrict the defect */ 627 ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr); 628 if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 629 630 if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 631 ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr); 632 if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 633 634 /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */ 635 ierr = VecCopy(B_c, X_c);CHKERRQ(ierr); 636 ierr = VecCopy(F_c, B_c);CHKERRQ(ierr); 637 ierr = VecCopy(X_c, F_c);CHKERRQ(ierr); 638 /* set initial guess of the coarse problem to the projected fine solution */ 639 ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 640 641 /* recurse to the next level */ 642 ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr); 643 ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 644 ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 645 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 646 snes->reason = SNES_DIVERGED_INNER; 647 PetscFunctionReturn(0); 648 } 649 /* correct as x <- x + I(x^c - Rx)*/ 650 ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 651 652 if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 653 ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr); 654 if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 655 } 656 PetscFunctionReturn(0); 657 } 658 659 #undef __FUNCT__ 660 #define __FUNCT__ "SNESFASCycle_Additive" 661 /* 662 663 The additive cycle looks like: 664 665 xhat = x 666 xhat = dS(x, b) 667 x = coarsecorrection(xhat, b_d) 668 x = x + nu*(xhat - x); 669 (optional) x = uS(x, b) 670 671 With the coarse RHS (defect correction) as below. 672 673 */ 674 PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X) 675 { 676 Vec F, B, Xhat; 677 Vec X_c, Xo_c, F_c, B_c; 678 PetscErrorCode ierr; 679 SNESConvergedReason reason; 680 PetscReal xnorm, fnorm, ynorm; 681 PetscBool lssuccess; 682 SNES next; 683 Mat restrct, interpolate; 684 SNES_FAS *fas = (SNES_FAS *)snes->data,*fasc; 685 686 PetscFunctionBegin; 687 ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 688 F = snes->vec_func; 689 B = snes->vec_rhs; 690 Xhat = snes->work[1]; 691 ierr = VecCopy(X, Xhat);CHKERRQ(ierr); 692 /* recurse first */ 693 if (next) { 694 fasc = (SNES_FAS*)next->data; 695 ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 696 ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 697 if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 698 ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr); 699 if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 700 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 701 X_c = next->vec_sol; 702 Xo_c = next->work[0]; 703 F_c = next->vec_func; 704 B_c = next->vec_rhs; 705 706 ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr); 707 /* restrict the defect */ 708 ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr); 709 710 /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 711 if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 712 ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr); 713 if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 714 ierr = VecCopy(B_c, X_c);CHKERRQ(ierr); 715 ierr = VecCopy(F_c, B_c);CHKERRQ(ierr); 716 ierr = VecCopy(X_c, F_c);CHKERRQ(ierr); 717 /* set initial guess of the coarse problem to the projected fine solution */ 718 ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 719 720 /* recurse */ 721 ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr); 722 ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 723 724 /* smooth on this level */ 725 ierr = SNESFASDownSmooth_Private(snes, B, X, F, &fnorm);CHKERRQ(ierr); 726 727 ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 728 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 729 snes->reason = SNES_DIVERGED_INNER; 730 PetscFunctionReturn(0); 731 } 732 733 /* correct as x <- x + I(x^c - Rx)*/ 734 ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr); 735 ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr); 736 737 /* additive correction of the coarse direction*/ 738 ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr); 739 ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssuccess);CHKERRQ(ierr); 740 if (!lssuccess) { 741 if (++snes->numFailures >= snes->maxFailures) { 742 snes->reason = SNES_DIVERGED_LINE_SEARCH; 743 PetscFunctionReturn(0); 744 } 745 } 746 ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm);CHKERRQ(ierr); 747 } else { 748 ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 749 } 750 PetscFunctionReturn(0); 751 } 752 753 #undef __FUNCT__ 754 #define __FUNCT__ "SNESFASCycle_Multiplicative" 755 /* 756 757 Defines the FAS cycle as: 758 759 fine problem: F(x) = 0 760 coarse problem: F^c(x) = b^c 761 762 b^c = F^c(I^c_fx^f - I^c_fF(x)) 763 764 correction: 765 766 x = x + I(x^c - Rx) 767 768 */ 769 PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X) 770 { 771 772 PetscErrorCode ierr; 773 Vec F,B; 774 SNES_FAS *fas = (SNES_FAS *)snes->data; 775 776 PetscFunctionBegin; 777 F = snes->vec_func; 778 B = snes->vec_rhs; 779 /* pre-smooth -- just update using the pre-smoother */ 780 ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 781 if (fas->level != 0) { 782 ierr = SNESFASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr); 783 ierr = SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 784 } 785 PetscFunctionReturn(0); 786 } 787 788 #undef __FUNCT__ 789 #define __FUNCT__ "SNESSolve_FAS" 790 791 PetscErrorCode SNESSolve_FAS(SNES snes) 792 { 793 PetscErrorCode ierr; 794 PetscInt i, maxits; 795 Vec X, F; 796 PetscReal fnorm; 797 SNES_FAS *fas = (SNES_FAS *)snes->data,*ffas; 798 DM dm; 799 PetscBool isFine; 800 801 PetscFunctionBegin; 802 maxits = snes->max_its; /* maximum number of iterations */ 803 snes->reason = SNES_CONVERGED_ITERATING; 804 X = snes->vec_sol; 805 F = snes->vec_func; 806 807 ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 808 /*norm setup */ 809 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 810 snes->iter = 0; 811 snes->norm = 0.; 812 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 813 if (!snes->vec_func_init_set) { 814 if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 815 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 816 if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 817 if (snes->domainerror) { 818 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 819 PetscFunctionReturn(0); 820 } 821 } else { 822 snes->vec_func_init_set = PETSC_FALSE; 823 } 824 825 if (!snes->norm_init_set) { 826 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 827 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 828 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 829 } else { 830 fnorm = snes->norm_init; 831 snes->norm_init_set = PETSC_FALSE; 832 } 833 834 snes->norm = fnorm; 835 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 836 SNESLogConvHistory(snes,fnorm,0); 837 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 838 839 /* set parameter for default relative tolerance convergence test */ 840 snes->ttol = fnorm*snes->rtol; 841 /* test convergence */ 842 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 843 if (snes->reason) PetscFunctionReturn(0); 844 845 846 if (isFine) { 847 /* propagate scale-dependent data up the hierarchy */ 848 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 849 for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) { 850 DM dmcoarse; 851 ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr); 852 ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr); 853 dm = dmcoarse; 854 } 855 } 856 857 for (i = 0; i < maxits; i++) { 858 /* Call general purpose update function */ 859 860 if (snes->ops->update) { 861 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 862 } 863 if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 864 ierr = SNESFASCycle_Multiplicative(snes, X);CHKERRQ(ierr); 865 } else { 866 ierr = SNESFASCycle_Additive(snes, X);CHKERRQ(ierr); 867 } 868 869 /* check for FAS cycle divergence */ 870 if (snes->reason != SNES_CONVERGED_ITERATING) { 871 PetscFunctionReturn(0); 872 } 873 874 /* Monitor convergence */ 875 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 876 snes->iter = i+1; 877 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 878 SNESLogConvHistory(snes,snes->norm,0); 879 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 880 /* Test for convergence */ 881 if (isFine) { 882 ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,snes->norm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 883 if (snes->reason) break; 884 } 885 } 886 if (i == maxits) { 887 ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 888 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 889 } 890 PetscFunctionReturn(0); 891 } 892