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