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