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