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