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 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 = NULL; 79 fas->smoothd = NULL; 80 fas->next = NULL; 81 fas->previous = NULL; 82 fas->fine = snes; 83 fas->interpolate = NULL; 84 fas->restrct = NULL; 85 fas->inject = NULL; 86 fas->monitor = NULL; 87 fas->usedmfornumberoflevels = PETSC_FALSE; 88 fas->fastype = SNES_FAS_MULTIPLICATIVE; 89 90 fas->eventsmoothsetup = 0; 91 fas->eventsmoothsolve = 0; 92 fas->eventresidual = 0; 93 fas->eventinterprestrict = 0; 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 = NULL; 161 snes->vec_func = NULL; 162 snes->vec_sol_update = NULL; 163 snes->vec_rhs = NULL; 164 165 /* reset the number of levels */ 166 ierr = SNESFASSetLevels(snes,dm_levels,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, PetscObjectComm((PetscObject)next), &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(PetscObjectComm((PetscObject)snes), 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, NULL, SNESFASGalerkinDefaultFunction, next);CHKERRQ(ierr); 220 } 221 if (fas->smoothd && fas->level != fas->levels - 1) { 222 ierr = SNESSetFunction(fas->smoothd, NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr); 223 } 224 if (fas->smoothu && fas->level != fas->levels - 1) { 225 ierr = SNESSetFunction(fas->smoothu, 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, 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 #include <petscdraw.h> 393 #undef __FUNCT__ 394 #define __FUNCT__ "SNESView_FAS" 395 PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 396 { 397 SNES_FAS *fas = (SNES_FAS*) snes->data; 398 PetscBool isFine,iascii,isdraw; 399 PetscInt i; 400 PetscErrorCode ierr; 401 SNES smoothu, smoothd, levelsnes; 402 403 PetscFunctionBegin; 404 ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 405 if (isFine) { 406 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 407 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 408 if (iascii) { 409 ierr = PetscViewerASCIIPrintf(viewer, "FAS: type is %s, levels=%D, cycles=%D\n", SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr); 410 if (fas->galerkin) { 411 ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr); 412 } else { 413 ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr); 414 } 415 for (i=0; i<fas->levels; i++) { 416 ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr); 417 ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr); 418 ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr); 419 if (!i) { 420 ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr); 421 } else { 422 ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 423 } 424 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 425 ierr = SNESView(smoothd,viewer);CHKERRQ(ierr); 426 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 427 if (i && (smoothd == smoothu)) { 428 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 429 } else if (i) { 430 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 431 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 432 ierr = SNESView(smoothu,viewer);CHKERRQ(ierr); 433 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 434 } 435 } 436 } else if (isdraw) { 437 PetscDraw draw; 438 PetscReal x,w,y,bottom,th,wth; 439 SNES_FAS *curfas = fas; 440 ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 441 ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 442 ierr = PetscDrawStringGetSize(draw,&wth,&th);CHKERRQ(ierr); 443 bottom = y - th; 444 while (curfas) { 445 if (!curfas->smoothu) { 446 ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr); 447 if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr); 448 ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 449 } else { 450 w = 0.5*PetscMin(1.0-x,x); 451 ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr); 452 if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr); 453 ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 454 ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr); 455 if (curfas->smoothu) ierr = SNESView(curfas->smoothu,viewer);CHKERRQ(ierr); 456 ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 457 } 458 /* this is totally bogus but we have no way of knowing how low the previous one was draw to */ 459 bottom -= 5*th; 460 if (curfas->next) curfas = (SNES_FAS*)curfas->next->data; 461 else curfas = NULL; 462 } 463 } 464 } 465 PetscFunctionReturn(0); 466 } 467 468 #undef __FUNCT__ 469 #define __FUNCT__ "SNESFASDownSmooth_Private" 470 /* 471 Defines the action of the downsmoother 472 */ 473 PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) 474 { 475 PetscErrorCode ierr = 0; 476 SNESConvergedReason reason; 477 Vec FPC; 478 SNES smoothd; 479 SNES_FAS *fas = (SNES_FAS*) snes->data; 480 481 PetscFunctionBegin; 482 ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr); 483 ierr = SNESSetInitialFunction(smoothd, F);CHKERRQ(ierr); 484 ierr = SNESSetInitialFunctionNorm(smoothd, *fnorm);CHKERRQ(ierr); 485 if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 486 ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr); 487 if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 488 /* check convergence reason for the smoother */ 489 ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr); 490 if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN)) { 491 snes->reason = SNES_DIVERGED_INNER; 492 PetscFunctionReturn(0); 493 } 494 ierr = SNESGetFunction(smoothd, &FPC, NULL, NULL);CHKERRQ(ierr); 495 ierr = VecCopy(FPC, F);CHKERRQ(ierr); 496 ierr = SNESGetFunctionNorm(smoothd, fnorm);CHKERRQ(ierr); 497 PetscFunctionReturn(0); 498 } 499 500 501 #undef __FUNCT__ 502 #define __FUNCT__ "SNESFASUpSmooth_Private" 503 /* 504 Defines the action of the upsmoother 505 */ 506 PetscErrorCode SNESFASUpSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) 507 { 508 PetscErrorCode ierr = 0; 509 SNESConvergedReason reason; 510 Vec FPC; 511 SNES smoothu; 512 SNES_FAS *fas = (SNES_FAS*) snes->data; 513 514 PetscFunctionBegin; 515 ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr); 516 if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 517 ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr); 518 if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 519 /* check convergence reason for the smoother */ 520 ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr); 521 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 522 snes->reason = SNES_DIVERGED_INNER; 523 PetscFunctionReturn(0); 524 } 525 ierr = SNESGetFunction(smoothu, &FPC, NULL, NULL);CHKERRQ(ierr); 526 ierr = VecCopy(FPC, F);CHKERRQ(ierr); 527 ierr = SNESGetFunctionNorm(smoothu, fnorm);CHKERRQ(ierr); 528 PetscFunctionReturn(0); 529 } 530 531 #undef __FUNCT__ 532 #define __FUNCT__ "SNESFASCreateCoarseVec" 533 /*@ 534 SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level 535 536 Collective 537 538 Input Arguments: 539 . snes - SNESFAS 540 541 Output Arguments: 542 . Xcoarse - vector on level one coarser than snes 543 544 Level: developer 545 546 .seealso: SNESFASSetRestriction(), SNESFASRestrict() 547 @*/ 548 PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse) 549 { 550 PetscErrorCode ierr; 551 SNES_FAS *fas = (SNES_FAS*)snes->data; 552 553 PetscFunctionBegin; 554 if (fas->rscale) { 555 ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr); 556 } else if (fas->inject) { 557 ierr = MatGetVecs(fas->inject,Xcoarse,NULL);CHKERRQ(ierr); 558 } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr); 559 PetscFunctionReturn(0); 560 } 561 562 #undef __FUNCT__ 563 #define __FUNCT__ "SNESFASRestrict" 564 /*@ 565 SNESFASRestrict - restrict a Vec to the next coarser level 566 567 Collective 568 569 Input Arguments: 570 + fine - SNES from which to restrict 571 - Xfine - vector to restrict 572 573 Output Arguments: 574 . Xcoarse - result of restriction 575 576 Level: developer 577 578 .seealso: SNESFASSetRestriction(), SNESFASSetInjection() 579 @*/ 580 PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse) 581 { 582 PetscErrorCode ierr; 583 SNES_FAS *fas = (SNES_FAS*)fine->data; 584 585 PetscFunctionBegin; 586 PetscValidHeaderSpecific(fine,SNES_CLASSID,1); 587 PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2); 588 PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3); 589 if (fas->inject) { 590 ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr); 591 } else { 592 ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr); 593 ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr); 594 } 595 PetscFunctionReturn(0); 596 } 597 598 #undef __FUNCT__ 599 #define __FUNCT__ "SNESFASCoarseCorrection" 600 /* 601 602 Performs the FAS coarse correction as: 603 604 fine problem: F(x) = 0 605 coarse problem: F^c(x) = b^c 606 607 b^c = F^c(I^c_fx^f - I^c_fF(x)) 608 609 */ 610 PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) 611 { 612 PetscErrorCode ierr; 613 Vec X_c, Xo_c, F_c, B_c; 614 SNESConvergedReason reason; 615 SNES next; 616 Mat restrct, interpolate; 617 SNES_FAS *fasc; 618 619 PetscFunctionBegin; 620 ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 621 if (next) { 622 fasc = (SNES_FAS*)next->data; 623 624 ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 625 ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 626 627 X_c = next->vec_sol; 628 Xo_c = next->work[0]; 629 F_c = next->vec_func; 630 B_c = next->vec_rhs; 631 632 if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 633 ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr); 634 /* restrict the defect */ 635 ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr); 636 if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 637 638 if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 639 ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr); 640 if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 641 642 /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */ 643 ierr = VecCopy(B_c, X_c);CHKERRQ(ierr); 644 ierr = VecCopy(F_c, B_c);CHKERRQ(ierr); 645 ierr = VecCopy(X_c, F_c);CHKERRQ(ierr); 646 /* set initial guess of the coarse problem to the projected fine solution */ 647 ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 648 649 /* recurse to the next level */ 650 ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr); 651 ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 652 ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 653 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 654 snes->reason = SNES_DIVERGED_INNER; 655 PetscFunctionReturn(0); 656 } 657 /* correct as x <- x + I(x^c - Rx)*/ 658 ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 659 660 if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 661 ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr); 662 if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 663 } 664 PetscFunctionReturn(0); 665 } 666 667 #undef __FUNCT__ 668 #define __FUNCT__ "SNESFASCycle_Additive" 669 /* 670 671 The additive cycle looks like: 672 673 xhat = x 674 xhat = dS(x, b) 675 x = coarsecorrection(xhat, b_d) 676 x = x + nu*(xhat - x); 677 (optional) x = uS(x, b) 678 679 With the coarse RHS (defect correction) as below. 680 681 */ 682 PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X) 683 { 684 Vec F, B, Xhat; 685 Vec X_c, Xo_c, F_c, B_c; 686 PetscErrorCode ierr; 687 SNESConvergedReason reason; 688 PetscReal xnorm, fnorm, ynorm; 689 PetscBool lssuccess; 690 SNES next; 691 Mat restrct, interpolate; 692 SNES_FAS *fas = (SNES_FAS*)snes->data,*fasc; 693 694 PetscFunctionBegin; 695 ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 696 F = snes->vec_func; 697 B = snes->vec_rhs; 698 Xhat = snes->work[1]; 699 ierr = VecCopy(X, Xhat);CHKERRQ(ierr); 700 /* recurse first */ 701 if (next) { 702 fasc = (SNES_FAS*)next->data; 703 ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 704 ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 705 if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 706 ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr); 707 if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 708 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 709 X_c = next->vec_sol; 710 Xo_c = next->work[0]; 711 F_c = next->vec_func; 712 B_c = next->vec_rhs; 713 714 ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr); 715 /* restrict the defect */ 716 ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr); 717 718 /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 719 if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 720 ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr); 721 if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 722 ierr = VecCopy(B_c, X_c);CHKERRQ(ierr); 723 ierr = VecCopy(F_c, B_c);CHKERRQ(ierr); 724 ierr = VecCopy(X_c, F_c);CHKERRQ(ierr); 725 /* set initial guess of the coarse problem to the projected fine solution */ 726 ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 727 728 /* recurse */ 729 ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr); 730 ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 731 732 /* smooth on this level */ 733 ierr = SNESFASDownSmooth_Private(snes, B, X, F, &fnorm);CHKERRQ(ierr); 734 735 ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 736 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 737 snes->reason = SNES_DIVERGED_INNER; 738 PetscFunctionReturn(0); 739 } 740 741 /* correct as x <- x + I(x^c - Rx)*/ 742 ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr); 743 ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr); 744 745 /* additive correction of the coarse direction*/ 746 ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr); 747 ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssuccess);CHKERRQ(ierr); 748 if (!lssuccess) { 749 if (++snes->numFailures >= snes->maxFailures) { 750 snes->reason = SNES_DIVERGED_LINE_SEARCH; 751 PetscFunctionReturn(0); 752 } 753 } 754 ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm);CHKERRQ(ierr); 755 } else { 756 ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 757 } 758 PetscFunctionReturn(0); 759 } 760 761 #undef __FUNCT__ 762 #define __FUNCT__ "SNESFASCycle_Multiplicative" 763 /* 764 765 Defines the FAS cycle as: 766 767 fine problem: F(x) = 0 768 coarse problem: F^c(x) = b^c 769 770 b^c = F^c(I^c_fx^f - I^c_fF(x)) 771 772 correction: 773 774 x = x + I(x^c - Rx) 775 776 */ 777 PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X) 778 { 779 780 PetscErrorCode ierr; 781 Vec F,B; 782 SNES_FAS *fas = (SNES_FAS*)snes->data; 783 784 PetscFunctionBegin; 785 F = snes->vec_func; 786 B = snes->vec_rhs; 787 /* pre-smooth -- just update using the pre-smoother */ 788 ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 789 if (fas->level != 0) { 790 ierr = SNESFASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr); 791 ierr = SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 792 } 793 PetscFunctionReturn(0); 794 } 795 796 #undef __FUNCT__ 797 #define __FUNCT__ "SNESSolve_FAS" 798 799 PetscErrorCode SNESSolve_FAS(SNES snes) 800 { 801 PetscErrorCode ierr; 802 PetscInt i, maxits; 803 Vec X, F; 804 PetscReal fnorm; 805 SNES_FAS *fas = (SNES_FAS*)snes->data,*ffas; 806 DM dm; 807 PetscBool isFine; 808 809 PetscFunctionBegin; 810 maxits = snes->max_its; /* maximum number of iterations */ 811 snes->reason = SNES_CONVERGED_ITERATING; 812 X = snes->vec_sol; 813 F = snes->vec_func; 814 815 ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 816 /*norm setup */ 817 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 818 snes->iter = 0; 819 snes->norm = 0.; 820 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 821 if (!snes->vec_func_init_set) { 822 if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 823 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 824 if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 825 if (snes->domainerror) { 826 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 827 PetscFunctionReturn(0); 828 } 829 } else snes->vec_func_init_set = PETSC_FALSE; 830 831 if (!snes->norm_init_set) { 832 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 833 if (PetscIsInfOrNanReal(fnorm)) { 834 snes->reason = SNES_DIVERGED_FNORM_NAN; 835 PetscFunctionReturn(0); 836 } 837 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 838 } else { 839 fnorm = snes->norm_init; 840 snes->norm_init_set = PETSC_FALSE; 841 } 842 843 snes->norm = fnorm; 844 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 845 ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 846 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 847 848 /* set parameter for default relative tolerance convergence test */ 849 snes->ttol = fnorm*snes->rtol; 850 /* test convergence */ 851 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 852 if (snes->reason) PetscFunctionReturn(0); 853 854 855 if (isFine) { 856 /* propagate scale-dependent data up the hierarchy */ 857 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 858 for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) { 859 DM dmcoarse; 860 ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr); 861 ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr); 862 dm = dmcoarse; 863 } 864 } 865 866 for (i = 0; i < maxits; i++) { 867 /* Call general purpose update function */ 868 869 if (snes->ops->update) { 870 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 871 } 872 if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 873 ierr = SNESFASCycle_Multiplicative(snes, X);CHKERRQ(ierr); 874 } else { 875 ierr = SNESFASCycle_Additive(snes, X);CHKERRQ(ierr); 876 } 877 878 /* check for FAS cycle divergence */ 879 if (snes->reason != SNES_CONVERGED_ITERATING) PetscFunctionReturn(0); 880 881 /* Monitor convergence */ 882 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 883 snes->iter = i+1; 884 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 885 ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 886 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 887 /* Test for convergence */ 888 if (isFine) { 889 ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,snes->norm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 890 if (snes->reason) break; 891 } 892 } 893 if (i == maxits) { 894 ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 895 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 896 } 897 PetscFunctionReturn(0); 898 } 899