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