xref: /petsc/src/snes/impls/fas/fas.c (revision 34d65b3cd98b5a242667ff96254c2856b2e2b837)
1421d9b32SPeter Brune /* Defines the basic SNES object */
26038b46bSPeter Brune #include <../src/snes/impls/fas/fasimpls.h>    /*I  "petscsnesfas.h"  I*/
3421d9b32SPeter Brune 
4*34d65b3cSPeter Brune const char *const SNESFASTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","SNESFASType","SNES_FAS",0};
507144faaSPeter Brune 
6421d9b32SPeter Brune extern PetscErrorCode SNESDestroy_FAS(SNES snes);
7421d9b32SPeter Brune extern PetscErrorCode SNESSetUp_FAS(SNES snes);
8421d9b32SPeter Brune extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes);
9421d9b32SPeter Brune extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer);
10421d9b32SPeter Brune extern PetscErrorCode SNESSolve_FAS(SNES snes);
11421d9b32SPeter Brune extern PetscErrorCode SNESReset_FAS(SNES snes);
126273346dSPeter Brune extern PetscErrorCode SNESFASGalerkinDefaultFunction(SNES, Vec, Vec, void*);
13ab8d36c9SPeter Brune extern PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES, SNES*);
14421d9b32SPeter Brune 
151fbfccc6SJed Brown /*MC
161fbfccc6SJed Brown 
171fbfccc6SJed Brown SNESFAS - Full Approximation Scheme nonlinear multigrid solver.
181fbfccc6SJed Brown 
19d3bc2379SPeter Brune    The nonlinear problem is solved by correction using coarse versions
20d3bc2379SPeter Brune    of the nonlinear problem.  This problem is perturbed so that a projected
21d3bc2379SPeter Brune    solution of the fine problem elicits no correction from the coarse problem.
22d3bc2379SPeter Brune 
23d3bc2379SPeter Brune Options Database:
24d3bc2379SPeter Brune +   -snes_fas_levels -  The number of levels
25d3bc2379SPeter Brune .   -snes_fas_cycles<1> -  The number of cycles -- 1 for V, 2 for W
26d3bc2379SPeter Brune .   -snes_fas_type<additive, multiplicative>  -  Additive or multiplicative cycle
27d3bc2379SPeter Brune .   -snes_fas_galerkin<PETSC_FALSE> -  Form coarse problems by projection back upon the fine problem
28d3bc2379SPeter Brune .   -snes_fas_smoothup<1> -  The number of iterations of the post-smoother
29d3bc2379SPeter Brune .   -snes_fas_smoothdown<1> -  The number of iterations of the pre-smoother
30d3bc2379SPeter Brune .   -snes_fas_monitor -  Monitor progress of all of the levels
31d3bc2379SPeter Brune .   -fas_levels_snes_ -  SNES options for all smoothers
327d84e935SPeter Brune .   -fas_levels_cycle_snes_ -  SNES options for all cycles
33d3bc2379SPeter Brune .   -fas_levels_i_snes_ -  SNES options for the smoothers on level i
347d84e935SPeter Brune .   -fas_levels_i_cycle_snes_ - SNES options for the cycle on level i
35d3bc2379SPeter Brune -   -fas_coarse_snes_ -  SNES options for the coarsest smoother
36d3bc2379SPeter Brune 
37d3bc2379SPeter Brune Notes:
38d3bc2379SPeter Brune    The organization of the FAS solver is slightly different from the organization of PCMG
39d3bc2379SPeter Brune    As each level has smoother SNES instances(down and potentially up) and a cycle SNES instance.
40d3bc2379SPeter Brune    The cycle SNES instance may be used for monitoring convergence on a particular level.
411fbfccc6SJed Brown 
427d84e935SPeter Brune Level: beginner
431fbfccc6SJed Brown 
44d3bc2379SPeter Brune .seealso: PCMG, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
451fbfccc6SJed Brown M*/
46421d9b32SPeter Brune 
47421d9b32SPeter Brune #undef __FUNCT__
48421d9b32SPeter Brune #define __FUNCT__ "SNESCreate_FAS"
498cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_FAS(SNES snes)
50421d9b32SPeter Brune {
51421d9b32SPeter Brune   SNES_FAS       *fas;
52421d9b32SPeter Brune   PetscErrorCode ierr;
53421d9b32SPeter Brune 
54421d9b32SPeter Brune   PetscFunctionBegin;
55421d9b32SPeter Brune   snes->ops->destroy        = SNESDestroy_FAS;
56421d9b32SPeter Brune   snes->ops->setup          = SNESSetUp_FAS;
57421d9b32SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_FAS;
58421d9b32SPeter Brune   snes->ops->view           = SNESView_FAS;
59421d9b32SPeter Brune   snes->ops->solve          = SNESSolve_FAS;
60421d9b32SPeter Brune   snes->ops->reset          = SNESReset_FAS;
61421d9b32SPeter Brune 
62ed020824SBarry Smith   snes->usesksp = PETSC_FALSE;
63ed020824SBarry Smith   snes->usespc  = PETSC_FALSE;
64ed020824SBarry Smith 
6588976e71SPeter Brune   if (!snes->tolerancesset) {
660e444f03SPeter Brune     snes->max_funcs = 30000;
670e444f03SPeter Brune     snes->max_its   = 10000;
6888976e71SPeter Brune   }
690e444f03SPeter Brune 
70421d9b32SPeter Brune   ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr);
711aa26658SKarl Rupp 
72421d9b32SPeter Brune   snes->data                  = (void*) fas;
73421d9b32SPeter Brune   fas->level                  = 0;
74293a7e31SPeter Brune   fas->levels                 = 1;
75ee78dd50SPeter Brune   fas->n_cycles               = 1;
76ee78dd50SPeter Brune   fas->max_up_it              = 1;
77ee78dd50SPeter Brune   fas->max_down_it            = 1;
780298fd71SBarry Smith   fas->smoothu                = NULL;
790298fd71SBarry Smith   fas->smoothd                = NULL;
800298fd71SBarry Smith   fas->next                   = NULL;
810298fd71SBarry Smith   fas->previous               = NULL;
82ab8d36c9SPeter Brune   fas->fine                   = snes;
830298fd71SBarry Smith   fas->interpolate            = NULL;
840298fd71SBarry Smith   fas->restrct                = NULL;
850298fd71SBarry Smith   fas->inject                 = NULL;
860298fd71SBarry Smith   fas->monitor                = NULL;
87cc05f883SPeter Brune   fas->usedmfornumberoflevels = PETSC_FALSE;
88ddebd997SPeter Brune   fas->fastype                = SNES_FAS_MULTIPLICATIVE;
890dd27c6cSPeter Brune 
900298fd71SBarry Smith   fas->eventsmoothsetup    = 0;
910298fd71SBarry Smith   fas->eventsmoothsolve    = 0;
920298fd71SBarry Smith   fas->eventresidual       = 0;
930298fd71SBarry Smith   fas->eventinterprestrict = 0;
94efe1f98aSPeter Brune   PetscFunctionReturn(0);
95efe1f98aSPeter Brune }
96efe1f98aSPeter Brune 
97421d9b32SPeter Brune #undef __FUNCT__
98421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS"
99421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes)
100421d9b32SPeter Brune {
10177df8cc4SPeter Brune   PetscErrorCode ierr  = 0;
102421d9b32SPeter Brune   SNES_FAS       * fas = (SNES_FAS*)snes->data;
103421d9b32SPeter Brune 
104421d9b32SPeter Brune   PetscFunctionBegin;
105ab8d36c9SPeter Brune   ierr = SNESDestroy(&fas->smoothu);CHKERRQ(ierr);
106ab8d36c9SPeter Brune   ierr = SNESDestroy(&fas->smoothd);CHKERRQ(ierr);
1073dccd265SPeter Brune   ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
108bccf9bb3SJed Brown   ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
109bccf9bb3SJed Brown   ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
110bccf9bb3SJed Brown   ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
111742fe5e2SPeter Brune   if (fas->next) ierr = SNESReset(fas->next);CHKERRQ(ierr);
112421d9b32SPeter Brune   PetscFunctionReturn(0);
113421d9b32SPeter Brune }
114421d9b32SPeter Brune 
115421d9b32SPeter Brune #undef __FUNCT__
116421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS"
117421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes)
118421d9b32SPeter Brune {
119421d9b32SPeter Brune   SNES_FAS       * fas = (SNES_FAS*)snes->data;
120742fe5e2SPeter Brune   PetscErrorCode ierr  = 0;
121421d9b32SPeter Brune 
122421d9b32SPeter Brune   PetscFunctionBegin;
123421d9b32SPeter Brune   /* recursively resets and then destroys */
12479d9a41aSPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
1251aa26658SKarl Rupp   if (fas->next) {
1261aa26658SKarl Rupp     ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
1271aa26658SKarl Rupp   }
128421d9b32SPeter Brune   ierr = PetscFree(fas);CHKERRQ(ierr);
129421d9b32SPeter Brune   PetscFunctionReturn(0);
130421d9b32SPeter Brune }
131421d9b32SPeter Brune 
132421d9b32SPeter Brune #undef __FUNCT__
133421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS"
134421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes)
135421d9b32SPeter Brune {
13648bfdf8aSPeter Brune   SNES_FAS       *fas = (SNES_FAS*) snes->data;
137421d9b32SPeter Brune   PetscErrorCode ierr;
138efe1f98aSPeter Brune   VecScatter     injscatter;
139d1adcc6fSPeter Brune   PetscInt       dm_levels;
1403dccd265SPeter Brune   Vec            vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */
141ab8d36c9SPeter Brune   SNES           next;
142ab8d36c9SPeter Brune   PetscBool      isFine;
143f89ba88eSPeter Brune   SNESLineSearch linesearch;
144f89ba88eSPeter Brune   SNESLineSearch slinesearch;
145f89ba88eSPeter Brune   void           *lsprectx,*lspostctx;
1466b2b7091SBarry Smith   PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*);
1476b2b7091SBarry Smith   PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*);
148eff52c0eSPeter Brune 
1496b2b7091SBarry Smith   PetscFunctionBegin;
150ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
151ab8d36c9SPeter Brune   if (fas->usedmfornumberoflevels && isFine) {
152d1adcc6fSPeter Brune     ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr);
153d1adcc6fSPeter Brune     dm_levels++;
154cc05f883SPeter Brune     if (dm_levels > fas->levels) {
1552e8ce248SJed Brown       /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESFASSetLevels_FAS*/
1563dccd265SPeter Brune       vec_sol              = snes->vec_sol;
1573dccd265SPeter Brune       vec_func             = snes->vec_func;
1583dccd265SPeter Brune       vec_sol_update       = snes->vec_sol_update;
1593dccd265SPeter Brune       vec_rhs              = snes->vec_rhs;
1600298fd71SBarry Smith       snes->vec_sol        = NULL;
1610298fd71SBarry Smith       snes->vec_func       = NULL;
1620298fd71SBarry Smith       snes->vec_sol_update = NULL;
1630298fd71SBarry Smith       snes->vec_rhs        = NULL;
1643dccd265SPeter Brune 
1653dccd265SPeter Brune       /* reset the number of levels */
1660298fd71SBarry Smith       ierr = SNESFASSetLevels(snes,dm_levels,NULL);CHKERRQ(ierr);
167cc05f883SPeter Brune       ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
1683dccd265SPeter Brune 
1693dccd265SPeter Brune       snes->vec_sol        = vec_sol;
1703dccd265SPeter Brune       snes->vec_func       = vec_func;
1713dccd265SPeter Brune       snes->vec_rhs        = vec_rhs;
1723dccd265SPeter Brune       snes->vec_sol_update = vec_sol_update;
173d1adcc6fSPeter Brune     }
174d1adcc6fSPeter Brune   }
175ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
176ab8d36c9SPeter Brune   if (!isFine) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */
1773dccd265SPeter Brune 
178fa0ddf94SBarry Smith   ierr = SNESSetWorkVecs(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
179cc05f883SPeter Brune 
180ab8d36c9SPeter Brune   /* set up the smoothers if they haven't already been set up */
181ab8d36c9SPeter Brune   if (!fas->smoothd) {
182ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
183ab8d36c9SPeter Brune   }
184ab8d36c9SPeter Brune 
18579d9a41aSPeter Brune   if (snes->dm) {
186ab8d36c9SPeter Brune     /* set the smoother DMs properly */
187ab8d36c9SPeter Brune     if (fas->smoothu) ierr = SNESSetDM(fas->smoothu, snes->dm);CHKERRQ(ierr);
188ab8d36c9SPeter Brune     ierr = SNESSetDM(fas->smoothd, snes->dm);CHKERRQ(ierr);
18979d9a41aSPeter Brune     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
190ab8d36c9SPeter Brune     if (next) {
19179d9a41aSPeter Brune       /* for now -- assume the DM and the evaluation functions have been set externally */
192ab8d36c9SPeter Brune       if (!next->dm) {
193ce94432eSBarry Smith         ierr = DMCoarsen(snes->dm, PetscObjectComm((PetscObject)next), &next->dm);CHKERRQ(ierr);
194ab8d36c9SPeter Brune         ierr = SNESSetDM(next, next->dm);CHKERRQ(ierr);
19579d9a41aSPeter Brune       }
19679d9a41aSPeter Brune       /* set the interpolation and restriction from the DM */
19779d9a41aSPeter Brune       if (!fas->interpolate) {
198ab8d36c9SPeter Brune         ierr = DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
199bccf9bb3SJed Brown         if (!fas->restrct) {
200bccf9bb3SJed Brown           ierr         = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr);
20179d9a41aSPeter Brune           fas->restrct = fas->interpolate;
20279d9a41aSPeter Brune         }
203bccf9bb3SJed Brown       }
20479d9a41aSPeter Brune       /* set the injection from the DM */
20579d9a41aSPeter Brune       if (!fas->inject) {
206ab8d36c9SPeter Brune         ierr = DMCreateInjection(next->dm, snes->dm, &injscatter);CHKERRQ(ierr);
207ce94432eSBarry Smith         ierr = MatCreateScatter(PetscObjectComm((PetscObject)snes), injscatter, &fas->inject);CHKERRQ(ierr);
20879d9a41aSPeter Brune         ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr);
20979d9a41aSPeter Brune       }
21079d9a41aSPeter Brune     }
21179d9a41aSPeter Brune   }
21279d9a41aSPeter Brune   /*pass the smoother, function, and jacobian up to the next level if it's not user set already */
21379d9a41aSPeter Brune   if (fas->galerkin) {
2141aa26658SKarl Rupp     if (next) {
2150298fd71SBarry Smith       ierr = SNESSetFunction(next, NULL, SNESFASGalerkinDefaultFunction, next);CHKERRQ(ierr);
2161aa26658SKarl Rupp     }
2171aa26658SKarl Rupp     if (fas->smoothd && fas->level != fas->levels - 1) {
2180298fd71SBarry Smith       ierr = SNESSetFunction(fas->smoothd, NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr);
2191aa26658SKarl Rupp     }
2201aa26658SKarl Rupp     if (fas->smoothu && fas->level != fas->levels - 1) {
2210298fd71SBarry Smith       ierr = SNESSetFunction(fas->smoothu, NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr);
2221aa26658SKarl Rupp     }
22379d9a41aSPeter Brune   }
22479d9a41aSPeter Brune 
225534ebe21SPeter Brune   /* sets the down (pre) smoother's default norm and sets it from options */
226534ebe21SPeter Brune   if (fas->smoothd) {
227bc3f2f05SPeter Brune     if (fas->level == 0 && fas->levels != 1) {
228365a6726SPeter Brune       ierr = SNESSetNormSchedule(fas->smoothd, SNES_NORM_NONE);CHKERRQ(ierr);
229534ebe21SPeter Brune     } else {
230365a6726SPeter Brune       ierr = SNESSetNormSchedule(fas->smoothd, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr);
231534ebe21SPeter Brune     }
2327fce8c19SPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)fas->smoothd);CHKERRQ(ierr);
233534ebe21SPeter Brune     ierr = SNESSetFromOptions(fas->smoothd);CHKERRQ(ierr);
2347601faf0SJed Brown     ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr);
2357601faf0SJed Brown     ierr = SNESGetLineSearch(fas->smoothd,&slinesearch);CHKERRQ(ierr);
2366b2b7091SBarry Smith     ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr);
2376b2b7091SBarry Smith     ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr);
2386b2b7091SBarry Smith     ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr);
2396b2b7091SBarry Smith     ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr);
240f89ba88eSPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr);
2410dd27c6cSPeter Brune 
2420dd27c6cSPeter Brune     fas->smoothd->vec_sol        = snes->vec_sol;
2430dd27c6cSPeter Brune     ierr                         = PetscObjectReference((PetscObject)snes->vec_sol);CHKERRQ(ierr);
2440dd27c6cSPeter Brune     fas->smoothd->vec_sol_update = snes->vec_sol_update;
2450dd27c6cSPeter Brune     ierr                         = PetscObjectReference((PetscObject)snes->vec_sol_update);CHKERRQ(ierr);
2460dd27c6cSPeter Brune     fas->smoothd->vec_func       = snes->vec_func;
2470dd27c6cSPeter Brune     ierr                         = PetscObjectReference((PetscObject)snes->vec_func);CHKERRQ(ierr);
2480dd27c6cSPeter Brune 
2490dd27c6cSPeter Brune     if (fas->eventsmoothsetup) {ierr = PetscLogEventBegin(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
2500dd27c6cSPeter Brune     ierr = SNESSetUp(fas->smoothd);CHKERRQ(ierr);
2510dd27c6cSPeter Brune     if (fas->eventsmoothsetup) {ierr = PetscLogEventEnd(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
252534ebe21SPeter Brune   }
253534ebe21SPeter Brune 
254534ebe21SPeter Brune   /* sets the up (post) smoother's default norm and sets it from options */
255534ebe21SPeter Brune   if (fas->smoothu) {
256534ebe21SPeter Brune     if (fas->level != fas->levels - 1) {
257365a6726SPeter Brune       ierr = SNESSetNormSchedule(fas->smoothu, SNES_NORM_NONE);CHKERRQ(ierr);
258534ebe21SPeter Brune     } else {
259365a6726SPeter Brune       ierr = SNESSetNormSchedule(fas->smoothu, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr);
260534ebe21SPeter Brune     }
2617fce8c19SPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)fas->smoothu);CHKERRQ(ierr);
262534ebe21SPeter Brune     ierr = SNESSetFromOptions(fas->smoothu);CHKERRQ(ierr);
2637601faf0SJed Brown     ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr);
2647601faf0SJed Brown     ierr = SNESGetLineSearch(fas->smoothu,&slinesearch);CHKERRQ(ierr);
2656b2b7091SBarry Smith     ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr);
2666b2b7091SBarry Smith     ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr);
2676b2b7091SBarry Smith     ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr);
2686b2b7091SBarry Smith     ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr);
269f89ba88eSPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr);
2700dd27c6cSPeter Brune 
2710dd27c6cSPeter Brune     fas->smoothu->vec_sol        = snes->vec_sol;
2720dd27c6cSPeter Brune     ierr                         = PetscObjectReference((PetscObject)snes->vec_sol);CHKERRQ(ierr);
2730dd27c6cSPeter Brune     fas->smoothu->vec_sol_update = snes->vec_sol_update;
2740dd27c6cSPeter Brune     ierr                         = PetscObjectReference((PetscObject)snes->vec_sol_update);CHKERRQ(ierr);
2750dd27c6cSPeter Brune     fas->smoothu->vec_func       = snes->vec_func;
2760dd27c6cSPeter Brune     ierr                         = PetscObjectReference((PetscObject)snes->vec_func);CHKERRQ(ierr);
2770dd27c6cSPeter Brune 
2780dd27c6cSPeter Brune     if (fas->eventsmoothsetup) {ierr = PetscLogEventBegin(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
2790dd27c6cSPeter Brune     ierr = SNESSetUp(fas->smoothu);CHKERRQ(ierr);
2800dd27c6cSPeter Brune     if (fas->eventsmoothsetup) {ierr = PetscLogEventEnd(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
2810dd27c6cSPeter Brune 
282534ebe21SPeter Brune   }
283d06165b7SPeter Brune 
284ab8d36c9SPeter Brune   if (next) {
28579d9a41aSPeter Brune     /* gotta set up the solution vector for this to work */
286ab8d36c9SPeter Brune     if (!next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_sol);CHKERRQ(ierr);}
287ab8d36c9SPeter Brune     if (!next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_rhs);CHKERRQ(ierr);}
2887fce8c19SPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)next);CHKERRQ(ierr);
2897601faf0SJed Brown     ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr);
2907601faf0SJed Brown     ierr = SNESGetLineSearch(fas->next,&slinesearch);CHKERRQ(ierr);
2916b2b7091SBarry Smith     ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr);
2926b2b7091SBarry Smith     ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr);
2936b2b7091SBarry Smith     ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr);
2946b2b7091SBarry Smith     ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr);
295f89ba88eSPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr);
296ab8d36c9SPeter Brune     ierr = SNESSetUp(next);CHKERRQ(ierr);
29779d9a41aSPeter Brune   }
2986273346dSPeter Brune   /* setup FAS work vectors */
2996273346dSPeter Brune   if (fas->galerkin) {
3006273346dSPeter Brune     ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr);
3016273346dSPeter Brune     ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr);
3026273346dSPeter Brune   }
303421d9b32SPeter Brune   PetscFunctionReturn(0);
304421d9b32SPeter Brune }
305421d9b32SPeter Brune 
306421d9b32SPeter Brune #undef __FUNCT__
307421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS"
308421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
309421d9b32SPeter Brune {
310ee78dd50SPeter Brune   SNES_FAS       *fas   = (SNES_FAS*) snes->data;
311ee78dd50SPeter Brune   PetscInt       levels = 1;
3124d26bfa5SPeter Brune   PetscBool      flg    = PETSC_FALSE, upflg = PETSC_FALSE, downflg = PETSC_FALSE, monflg = PETSC_FALSE, galerkinflg = PETSC_FALSE;
313421d9b32SPeter Brune   PetscErrorCode ierr;
314ee78dd50SPeter Brune   char           monfilename[PETSC_MAX_PATH_LEN];
31507144faaSPeter Brune   SNESFASType    fastype;
316fde0ff24SPeter Brune   const char     *optionsprefix;
317f1c6b773SPeter Brune   SNESLineSearch linesearch;
31866585501SPeter Brune   PetscInt       m, n_up, n_down;
319ab8d36c9SPeter Brune   SNES           next;
320ab8d36c9SPeter Brune   PetscBool      isFine;
321421d9b32SPeter Brune 
322421d9b32SPeter Brune   PetscFunctionBegin;
323ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
324c90fad12SPeter Brune   ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr);
325ee78dd50SPeter Brune 
326ab8d36c9SPeter Brune   /* number of levels -- only process most options on the finest level */
327ab8d36c9SPeter Brune   if (isFine) {
328ee78dd50SPeter Brune     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
329c732cbdbSBarry Smith     if (!flg && snes->dm) {
330c732cbdbSBarry Smith       ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr);
331c732cbdbSBarry Smith       levels++;
332d1adcc6fSPeter Brune       fas->usedmfornumberoflevels = PETSC_TRUE;
333c732cbdbSBarry Smith     }
3340298fd71SBarry Smith     ierr    = SNESFASSetLevels(snes, levels, NULL);CHKERRQ(ierr);
33507144faaSPeter Brune     fastype = fas->fastype;
33607144faaSPeter Brune     ierr    = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr);
33707144faaSPeter Brune     if (flg) {
33807144faaSPeter Brune       ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr);
33907144faaSPeter Brune     }
340ee78dd50SPeter Brune 
341fde0ff24SPeter Brune     ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
342ab8d36c9SPeter Brune     ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&m,&flg);CHKERRQ(ierr);
343ab8d36c9SPeter Brune     if (flg) {
344ab8d36c9SPeter Brune       ierr = SNESFASSetCycles(snes, m);CHKERRQ(ierr);
345fde0ff24SPeter Brune     }
346fde0ff24SPeter Brune 
347ab8d36c9SPeter Brune     ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFASSetGalerkin",fas->galerkin,&galerkinflg,&flg);CHKERRQ(ierr);
348ab8d36c9SPeter Brune     if (flg) {
349ab8d36c9SPeter Brune       ierr = SNESFASSetGalerkin(snes, galerkinflg);CHKERRQ(ierr);
350ab8d36c9SPeter Brune     }
351ee78dd50SPeter Brune 
35266585501SPeter Brune     ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",fas->max_up_it,&n_up,&upflg);CHKERRQ(ierr);
353162d76ddSPeter Brune 
35466585501SPeter Brune     ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",fas->max_down_it,&n_down,&downflg);CHKERRQ(ierr);
355162d76ddSPeter Brune 
356c8c899caSPeter Brune     ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESFASSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr);
357c8c899caSPeter Brune     if (monflg) ierr = SNESFASSetMonitor(snes, PETSC_TRUE);CHKERRQ(ierr);
3580dd27c6cSPeter Brune 
3590dd27c6cSPeter Brune     flg    = PETSC_FALSE;
3600dd27c6cSPeter Brune     monflg = PETSC_TRUE;
3610dd27c6cSPeter Brune     ierr   = PetscOptionsBool("-snes_fas_log","Log times for each FAS level","SNESFASSetLog",monflg,&monflg,&flg);CHKERRQ(ierr);
3620dd27c6cSPeter Brune     if (flg) {ierr = SNESFASSetLog(snes,monflg);CHKERRQ(ierr);}
363ab8d36c9SPeter Brune   }
364ee78dd50SPeter Brune 
365421d9b32SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
3668cc86e31SPeter Brune   /* setup from the determined types if there is no pointwise procedure or smoother defined */
367162d76ddSPeter Brune   if (upflg) {
36866585501SPeter Brune     ierr = SNESFASSetNumberSmoothUp(snes,n_up);CHKERRQ(ierr);
369162d76ddSPeter Brune   }
370162d76ddSPeter Brune   if (downflg) {
37166585501SPeter Brune     ierr = SNESFASSetNumberSmoothDown(snes,n_down);CHKERRQ(ierr);
372162d76ddSPeter Brune   }
373eff52c0eSPeter Brune 
3749e764e56SPeter Brune   /* set up the default line search for coarse grid corrections */
3759e764e56SPeter Brune   if (fas->fastype == SNES_FAS_ADDITIVE) {
3769e764e56SPeter Brune     if (!snes->linesearch) {
3777601faf0SJed Brown       ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
3781a4f838cSPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr);
3799e764e56SPeter Brune     }
3809e764e56SPeter Brune   }
3819e764e56SPeter Brune 
382ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
383ee78dd50SPeter Brune   /* recursive option setting for the smoothers */
384ab8d36c9SPeter Brune   if (next) {ierr = SNESSetFromOptions(next);CHKERRQ(ierr);}
385421d9b32SPeter Brune   PetscFunctionReturn(0);
386421d9b32SPeter Brune }
387421d9b32SPeter Brune 
3889804daf3SBarry Smith #include <petscdraw.h>
389421d9b32SPeter Brune #undef __FUNCT__
390421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS"
391421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
392421d9b32SPeter Brune {
393421d9b32SPeter Brune   SNES_FAS       *fas = (SNES_FAS*) snes->data;
394656ede7eSPeter Brune   PetscBool      isFine,iascii,isdraw;
395ab8d36c9SPeter Brune   PetscInt       i;
396421d9b32SPeter Brune   PetscErrorCode ierr;
397ab8d36c9SPeter Brune   SNES           smoothu, smoothd, levelsnes;
398421d9b32SPeter Brune 
399421d9b32SPeter Brune   PetscFunctionBegin;
400ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
401ab8d36c9SPeter Brune   if (isFine) {
402251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
403656ede7eSPeter Brune     ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
404421d9b32SPeter Brune     if (iascii) {
405ab8d36c9SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "FAS: type is %s, levels=%D, cycles=%D\n",  SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr);
406ab8d36c9SPeter Brune       if (fas->galerkin) {
407ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
408421d9b32SPeter Brune       } else {
409ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
410421d9b32SPeter Brune       }
411ab8d36c9SPeter Brune       for (i=0; i<fas->levels; i++) {
412ab8d36c9SPeter Brune         ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr);
413ab8d36c9SPeter Brune         ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr);
414ab8d36c9SPeter Brune         ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr);
415ab8d36c9SPeter Brune         if (!i) {
416ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr);
417421d9b32SPeter Brune         } else {
418ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
419421d9b32SPeter Brune         }
420ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
421166b3ea4SJed Brown         if (smoothd) {
422ab8d36c9SPeter Brune           ierr = SNESView(smoothd,viewer);CHKERRQ(ierr);
423166b3ea4SJed Brown         } else {
424166b3ea4SJed Brown           ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr);
425166b3ea4SJed Brown         }
426ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
427ab8d36c9SPeter Brune         if (i && (smoothd == smoothu)) {
428ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
429ab8d36c9SPeter Brune         } else if (i) {
430ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
431ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
432166b3ea4SJed Brown           if (smoothu) {
433ab8d36c9SPeter Brune             ierr = SNESView(smoothu,viewer);CHKERRQ(ierr);
434166b3ea4SJed Brown           } else {
435166b3ea4SJed Brown             ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr);
436166b3ea4SJed Brown           }
437ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
438ab8d36c9SPeter Brune         }
439ab8d36c9SPeter Brune       }
440656ede7eSPeter Brune     } else if (isdraw) {
441656ede7eSPeter Brune       PetscDraw draw;
442b4375e8dSPeter Brune       PetscReal x,w,y,bottom,th,wth;
443656ede7eSPeter Brune       SNES_FAS  *curfas = fas;
444656ede7eSPeter Brune       ierr   = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
445656ede7eSPeter Brune       ierr   = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
446656ede7eSPeter Brune       ierr   = PetscDrawStringGetSize(draw,&wth,&th);CHKERRQ(ierr);
447656ede7eSPeter Brune       bottom = y - th;
448656ede7eSPeter Brune       while (curfas) {
449b4375e8dSPeter Brune         if (!curfas->smoothu) {
450656ede7eSPeter Brune           ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
451656ede7eSPeter Brune           if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr);
452656ede7eSPeter Brune           ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
453b4375e8dSPeter Brune         } else {
454b4375e8dSPeter Brune           w    = 0.5*PetscMin(1.0-x,x);
455b4375e8dSPeter Brune           ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
456b4375e8dSPeter Brune           if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr);
457b4375e8dSPeter Brune           ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
458b4375e8dSPeter Brune           ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
459b4375e8dSPeter Brune           if (curfas->smoothu) ierr = SNESView(curfas->smoothu,viewer);CHKERRQ(ierr);
460b4375e8dSPeter Brune           ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
461b4375e8dSPeter Brune         }
462656ede7eSPeter Brune         /* this is totally bogus but we have no way of knowing how low the previous one was draw to */
463656ede7eSPeter Brune         bottom -= 5*th;
4641aa26658SKarl Rupp         if (curfas->next) curfas = (SNES_FAS*)curfas->next->data;
4650298fd71SBarry Smith         else curfas = NULL;
466656ede7eSPeter Brune       }
467421d9b32SPeter Brune     }
468ab8d36c9SPeter Brune   }
469421d9b32SPeter Brune   PetscFunctionReturn(0);
470421d9b32SPeter Brune }
471421d9b32SPeter Brune 
472421d9b32SPeter Brune #undef __FUNCT__
47391f99d7cSPeter Brune #define __FUNCT__ "SNESFASDownSmooth_Private"
47439bd7f45SPeter Brune /*
47539bd7f45SPeter Brune Defines the action of the downsmoother
47639bd7f45SPeter Brune  */
47791f99d7cSPeter Brune PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
478b9c2fdf1SPeter Brune {
47939bd7f45SPeter Brune   PetscErrorCode      ierr = 0;
480742fe5e2SPeter Brune   SNESConvergedReason reason;
481ab8d36c9SPeter Brune   Vec                 FPC;
482ab8d36c9SPeter Brune   SNES                smoothd;
4830dd27c6cSPeter Brune   SNES_FAS            *fas = (SNES_FAS*) snes->data;
4846e111a19SKarl Rupp 
485421d9b32SPeter Brune   PetscFunctionBegin;
486ab8d36c9SPeter Brune   ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr);
487e4ed7901SPeter Brune   ierr = SNESSetInitialFunction(smoothd, F);CHKERRQ(ierr);
4880dd27c6cSPeter Brune   if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
489ab8d36c9SPeter Brune   ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr);
4900dd27c6cSPeter Brune   if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
491742fe5e2SPeter Brune   /* check convergence reason for the smoother */
492ab8d36c9SPeter Brune   ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr);
493e70c42e5SPeter Brune   if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN)) {
494742fe5e2SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
495742fe5e2SPeter Brune     PetscFunctionReturn(0);
496742fe5e2SPeter Brune   }
4970298fd71SBarry Smith   ierr = SNESGetFunction(smoothd, &FPC, NULL, NULL);CHKERRQ(ierr);
4984b32a720SPeter Brune   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
499b9c2fdf1SPeter Brune   ierr = SNESGetFunctionNorm(smoothd, fnorm);CHKERRQ(ierr);
50039bd7f45SPeter Brune   PetscFunctionReturn(0);
50139bd7f45SPeter Brune }
50239bd7f45SPeter Brune 
50339bd7f45SPeter Brune 
50439bd7f45SPeter Brune #undef __FUNCT__
50591f99d7cSPeter Brune #define __FUNCT__ "SNESFASUpSmooth_Private"
50639bd7f45SPeter Brune /*
50707144faaSPeter Brune Defines the action of the upsmoother
50839bd7f45SPeter Brune  */
5090adebc6cSBarry Smith PetscErrorCode SNESFASUpSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
5100adebc6cSBarry Smith {
51139bd7f45SPeter Brune   PetscErrorCode      ierr = 0;
51239bd7f45SPeter Brune   SNESConvergedReason reason;
513ab8d36c9SPeter Brune   Vec                 FPC;
514ab8d36c9SPeter Brune   SNES                smoothu;
5150dd27c6cSPeter Brune   SNES_FAS            *fas = (SNES_FAS*) snes->data;
516ab8d36c9SPeter Brune 
5176e111a19SKarl Rupp   PetscFunctionBegin;
518ab8d36c9SPeter Brune   ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr);
5190dd27c6cSPeter Brune   if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
520ab8d36c9SPeter Brune   ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr);
5210dd27c6cSPeter Brune   if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
52239bd7f45SPeter Brune   /* check convergence reason for the smoother */
523ab8d36c9SPeter Brune   ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr);
52439bd7f45SPeter Brune   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
52539bd7f45SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
52639bd7f45SPeter Brune     PetscFunctionReturn(0);
52739bd7f45SPeter Brune   }
5280298fd71SBarry Smith   ierr = SNESGetFunction(smoothu, &FPC, NULL, NULL);CHKERRQ(ierr);
5294b32a720SPeter Brune   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
530b9c2fdf1SPeter Brune   ierr = SNESGetFunctionNorm(smoothu, fnorm);CHKERRQ(ierr);
53139bd7f45SPeter Brune   PetscFunctionReturn(0);
53239bd7f45SPeter Brune }
53339bd7f45SPeter Brune 
53439bd7f45SPeter Brune #undef __FUNCT__
535938e4a01SJed Brown #define __FUNCT__ "SNESFASCreateCoarseVec"
536938e4a01SJed Brown /*@
537938e4a01SJed Brown    SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level
538938e4a01SJed Brown 
539938e4a01SJed Brown    Collective
540938e4a01SJed Brown 
541938e4a01SJed Brown    Input Arguments:
542938e4a01SJed Brown .  snes - SNESFAS
543938e4a01SJed Brown 
544938e4a01SJed Brown    Output Arguments:
545938e4a01SJed Brown .  Xcoarse - vector on level one coarser than snes
546938e4a01SJed Brown 
547938e4a01SJed Brown    Level: developer
548938e4a01SJed Brown 
549938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict()
550938e4a01SJed Brown @*/
551938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse)
552938e4a01SJed Brown {
553938e4a01SJed Brown   PetscErrorCode ierr;
554938e4a01SJed Brown   SNES_FAS       *fas = (SNES_FAS*)snes->data;
555938e4a01SJed Brown 
556938e4a01SJed Brown   PetscFunctionBegin;
5571aa26658SKarl Rupp   if (fas->rscale) {
5581aa26658SKarl Rupp     ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);
559f5af7f23SKarl Rupp   } else if (fas->inject) {
5600298fd71SBarry Smith     ierr = MatGetVecs(fas->inject,Xcoarse,NULL);CHKERRQ(ierr);
561ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr);
562938e4a01SJed Brown   PetscFunctionReturn(0);
563938e4a01SJed Brown }
564938e4a01SJed Brown 
565e9923e8dSJed Brown #undef __FUNCT__
566e9923e8dSJed Brown #define __FUNCT__ "SNESFASRestrict"
567e9923e8dSJed Brown /*@
568e9923e8dSJed Brown    SNESFASRestrict - restrict a Vec to the next coarser level
569e9923e8dSJed Brown 
570e9923e8dSJed Brown    Collective
571e9923e8dSJed Brown 
572e9923e8dSJed Brown    Input Arguments:
573e9923e8dSJed Brown +  fine - SNES from which to restrict
574e9923e8dSJed Brown -  Xfine - vector to restrict
575e9923e8dSJed Brown 
576e9923e8dSJed Brown    Output Arguments:
577e9923e8dSJed Brown .  Xcoarse - result of restriction
578e9923e8dSJed Brown 
579e9923e8dSJed Brown    Level: developer
580e9923e8dSJed Brown 
581e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection()
582e9923e8dSJed Brown @*/
583e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse)
584e9923e8dSJed Brown {
585e9923e8dSJed Brown   PetscErrorCode ierr;
586e9923e8dSJed Brown   SNES_FAS       *fas = (SNES_FAS*)fine->data;
587e9923e8dSJed Brown 
588e9923e8dSJed Brown   PetscFunctionBegin;
589e9923e8dSJed Brown   PetscValidHeaderSpecific(fine,SNES_CLASSID,1);
590e9923e8dSJed Brown   PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2);
591e9923e8dSJed Brown   PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3);
592e9923e8dSJed Brown   if (fas->inject) {
593e9923e8dSJed Brown     ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr);
594e9923e8dSJed Brown   } else {
595e9923e8dSJed Brown     ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr);
596e9923e8dSJed Brown     ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr);
597e9923e8dSJed Brown   }
598e9923e8dSJed Brown   PetscFunctionReturn(0);
599e9923e8dSJed Brown }
600e9923e8dSJed Brown 
601e9923e8dSJed Brown #undef __FUNCT__
6028c40d5fbSBarry Smith #define __FUNCT__ "SNESFASCoarseCorrection"
60339bd7f45SPeter Brune /*
60439bd7f45SPeter Brune 
60539bd7f45SPeter Brune Performs the FAS coarse correction as:
60639bd7f45SPeter Brune 
60739bd7f45SPeter Brune fine problem: F(x) = 0
60839bd7f45SPeter Brune coarse problem: F^c(x) = b^c
60939bd7f45SPeter Brune 
61039bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
61139bd7f45SPeter Brune 
61239bd7f45SPeter Brune  */
6130adebc6cSBarry Smith PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new)
6140adebc6cSBarry Smith {
61539bd7f45SPeter Brune   PetscErrorCode      ierr;
61639bd7f45SPeter Brune   Vec                 X_c, Xo_c, F_c, B_c;
61739bd7f45SPeter Brune   SNESConvergedReason reason;
618ab8d36c9SPeter Brune   SNES                next;
619ab8d36c9SPeter Brune   Mat                 restrct, interpolate;
6200dd27c6cSPeter Brune   SNES_FAS            *fasc;
6215fd66863SKarl Rupp 
62239bd7f45SPeter Brune   PetscFunctionBegin;
623ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
624ab8d36c9SPeter Brune   if (next) {
6250dd27c6cSPeter Brune     fasc = (SNES_FAS*)next->data;
6260dd27c6cSPeter Brune 
627ab8d36c9SPeter Brune     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
628ab8d36c9SPeter Brune     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
629ab8d36c9SPeter Brune 
630ab8d36c9SPeter Brune     X_c  = next->vec_sol;
631ab8d36c9SPeter Brune     Xo_c = next->work[0];
632ab8d36c9SPeter Brune     F_c  = next->vec_func;
633ab8d36c9SPeter Brune     B_c  = next->vec_rhs;
634efe1f98aSPeter Brune 
6350dd27c6cSPeter Brune     if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
636938e4a01SJed Brown     ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr);
637293a7e31SPeter Brune     /* restrict the defect */
638ab8d36c9SPeter Brune     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
6390dd27c6cSPeter Brune     if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
6400dd27c6cSPeter Brune 
6410dd27c6cSPeter Brune     if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);}
642ab8d36c9SPeter Brune     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
6430dd27c6cSPeter Brune     if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);}
6440dd27c6cSPeter Brune 
6450dd27c6cSPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */
646e4ed7901SPeter Brune     ierr = VecCopy(B_c, X_c);CHKERRQ(ierr);
647b9c2fdf1SPeter Brune     ierr = VecCopy(F_c, B_c);CHKERRQ(ierr);
648e4ed7901SPeter Brune     ierr = VecCopy(X_c, F_c);CHKERRQ(ierr);
649ee78dd50SPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
650ee78dd50SPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
651c90fad12SPeter Brune 
652c90fad12SPeter Brune     /* recurse to the next level */
653e4ed7901SPeter Brune     ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr);
654ab8d36c9SPeter Brune     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
655ab8d36c9SPeter Brune     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
656742fe5e2SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
657742fe5e2SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
658742fe5e2SPeter Brune       PetscFunctionReturn(0);
659742fe5e2SPeter Brune     }
660fa9694d7SPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
661fa9694d7SPeter Brune     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
6620dd27c6cSPeter Brune 
6630dd27c6cSPeter Brune     if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
664ab8d36c9SPeter Brune     ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr);
6650dd27c6cSPeter Brune     if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
666293a7e31SPeter Brune   }
66739bd7f45SPeter Brune   PetscFunctionReturn(0);
66839bd7f45SPeter Brune }
66939bd7f45SPeter Brune 
67039bd7f45SPeter Brune #undef __FUNCT__
6712cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Additive"
67239bd7f45SPeter Brune /*
67339bd7f45SPeter Brune 
67439bd7f45SPeter Brune The additive cycle looks like:
67539bd7f45SPeter Brune 
67607144faaSPeter Brune xhat = x
67707144faaSPeter Brune xhat = dS(x, b)
67807144faaSPeter Brune x = coarsecorrection(xhat, b_d)
67907144faaSPeter Brune x = x + nu*(xhat - x);
68039bd7f45SPeter Brune (optional) x = uS(x, b)
68139bd7f45SPeter Brune 
68239bd7f45SPeter Brune With the coarse RHS (defect correction) as below.
68339bd7f45SPeter Brune 
68439bd7f45SPeter Brune  */
6850adebc6cSBarry Smith PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X)
6860adebc6cSBarry Smith {
68707144faaSPeter Brune   Vec                 F, B, Xhat;
68822c1e704SPeter Brune   Vec                 X_c, Xo_c, F_c, B_c;
68939bd7f45SPeter Brune   PetscErrorCode      ierr;
69007144faaSPeter Brune   SNESConvergedReason reason;
69122c1e704SPeter Brune   PetscReal           xnorm, fnorm, ynorm;
69222c1e704SPeter Brune   PetscBool           lssuccess;
693ab8d36c9SPeter Brune   SNES                next;
694ab8d36c9SPeter Brune   Mat                 restrct, interpolate;
6950dd27c6cSPeter Brune   SNES_FAS            *fas = (SNES_FAS*)snes->data,*fasc;
6960dd27c6cSPeter Brune 
69739bd7f45SPeter Brune   PetscFunctionBegin;
698ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
69939bd7f45SPeter Brune   F    = snes->vec_func;
70039bd7f45SPeter Brune   B    = snes->vec_rhs;
701e7f468e7SPeter Brune   Xhat = snes->work[1];
70207144faaSPeter Brune   ierr = VecCopy(X, Xhat);CHKERRQ(ierr);
70307144faaSPeter Brune   /* recurse first */
704ab8d36c9SPeter Brune   if (next) {
7050dd27c6cSPeter Brune     fasc = (SNES_FAS*)next->data;
706ab8d36c9SPeter Brune     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
707ab8d36c9SPeter Brune     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
7080dd27c6cSPeter Brune     if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);}
70907144faaSPeter Brune     ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr);
7100dd27c6cSPeter Brune     if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);}
711c2a02606SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
712ab8d36c9SPeter Brune     X_c  = next->vec_sol;
713ab8d36c9SPeter Brune     Xo_c = next->work[0];
714ab8d36c9SPeter Brune     F_c  = next->vec_func;
715ab8d36c9SPeter Brune     B_c  = next->vec_rhs;
71639bd7f45SPeter Brune 
717938e4a01SJed Brown     ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr);
71807144faaSPeter Brune     /* restrict the defect */
719ab8d36c9SPeter Brune     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
72007144faaSPeter Brune 
72107144faaSPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
7220dd27c6cSPeter Brune     if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);}
723ab8d36c9SPeter Brune     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
7240dd27c6cSPeter Brune     if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);}
725e4ed7901SPeter Brune     ierr = VecCopy(B_c, X_c);CHKERRQ(ierr);
726b9c2fdf1SPeter Brune     ierr = VecCopy(F_c, B_c);CHKERRQ(ierr);
727e4ed7901SPeter Brune     ierr = VecCopy(X_c, F_c);CHKERRQ(ierr);
72807144faaSPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
72907144faaSPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
73007144faaSPeter Brune 
73107144faaSPeter Brune     /* recurse */
732e4ed7901SPeter Brune     ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr);
733ab8d36c9SPeter Brune     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
73407144faaSPeter Brune 
73507144faaSPeter Brune     /* smooth on this level */
73691f99d7cSPeter Brune     ierr = SNESFASDownSmooth_Private(snes, B, X, F, &fnorm);CHKERRQ(ierr);
73707144faaSPeter Brune 
738ab8d36c9SPeter Brune     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
73907144faaSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
74007144faaSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
74107144faaSPeter Brune       PetscFunctionReturn(0);
74207144faaSPeter Brune     }
74307144faaSPeter Brune 
74407144faaSPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
745c68acad4SPeter Brune     ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr);
746ab8d36c9SPeter Brune     ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr);
74707144faaSPeter Brune 
748ddebd997SPeter Brune     /* additive correction of the coarse direction*/
749f1c6b773SPeter Brune     ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr);
750f1c6b773SPeter Brune     ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssuccess);CHKERRQ(ierr);
7519e764e56SPeter Brune     if (!lssuccess) {
7529e764e56SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
7539e764e56SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
7549e764e56SPeter Brune         PetscFunctionReturn(0);
7559e764e56SPeter Brune       }
7569e764e56SPeter Brune     }
757b9c2fdf1SPeter Brune     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm);CHKERRQ(ierr);
75807144faaSPeter Brune   } else {
75991f99d7cSPeter Brune     ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
76007144faaSPeter Brune   }
76139bd7f45SPeter Brune   PetscFunctionReturn(0);
76239bd7f45SPeter Brune }
76339bd7f45SPeter Brune 
76439bd7f45SPeter Brune #undef __FUNCT__
7652cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Multiplicative"
76639bd7f45SPeter Brune /*
76739bd7f45SPeter Brune 
76839bd7f45SPeter Brune Defines the FAS cycle as:
76939bd7f45SPeter Brune 
77039bd7f45SPeter Brune fine problem: F(x) = 0
77139bd7f45SPeter Brune coarse problem: F^c(x) = b^c
77239bd7f45SPeter Brune 
77339bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
77439bd7f45SPeter Brune 
77539bd7f45SPeter Brune correction:
77639bd7f45SPeter Brune 
77739bd7f45SPeter Brune x = x + I(x^c - Rx)
77839bd7f45SPeter Brune 
77939bd7f45SPeter Brune  */
7800adebc6cSBarry Smith PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X)
7810adebc6cSBarry Smith {
78239bd7f45SPeter Brune 
78339bd7f45SPeter Brune   PetscErrorCode ierr;
78439bd7f45SPeter Brune   Vec            F,B;
785*34d65b3cSPeter Brune   SNES           next;
78639bd7f45SPeter Brune 
78739bd7f45SPeter Brune   PetscFunctionBegin;
78839bd7f45SPeter Brune   F = snes->vec_func;
78939bd7f45SPeter Brune   B = snes->vec_rhs;
79039bd7f45SPeter Brune   /* pre-smooth -- just update using the pre-smoother */
791*34d65b3cSPeter Brune   ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr);
79291f99d7cSPeter Brune   ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
793*34d65b3cSPeter Brune   if (next) {
7948c40d5fbSBarry Smith     ierr = SNESFASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr);
79591f99d7cSPeter Brune     ierr = SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
796fe6f9142SPeter Brune   }
797fa9694d7SPeter Brune   PetscFunctionReturn(0);
798421d9b32SPeter Brune }
799421d9b32SPeter Brune 
800421d9b32SPeter Brune #undef __FUNCT__
8018c94862eSPeter Brune #define __FUNCT__ "SNESFASCycleSetupPhase_Full"
8028c94862eSPeter Brune PetscErrorCode SNESFASCycleSetupPhase_Full(SNES snes)
8038c94862eSPeter Brune {
8048c94862eSPeter Brune   SNES           next;
8058c94862eSPeter Brune   SNES_FAS       *fas = (SNES_FAS*)snes->data;
8068c94862eSPeter Brune   PetscBool      isFine;
8078c94862eSPeter Brune   PetscErrorCode ierr;
8088c94862eSPeter Brune 
8098c94862eSPeter Brune   PetscFunctionBegin;
8108c94862eSPeter Brune   /* pre-smooth -- just update using the pre-smoother */
8118c94862eSPeter Brune   ierr = SNESFASCycleIsFine(snes,&isFine);CHKERRQ(ierr);
8128c94862eSPeter Brune   ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr);
8138c94862eSPeter Brune   fas->full_stage = 0;
8148c94862eSPeter Brune   if (next) {ierr = SNESFASCycleSetupPhase_Full(next);CHKERRQ(ierr);}
8158c94862eSPeter Brune   PetscFunctionReturn(0);
8168c94862eSPeter Brune }
8178c94862eSPeter Brune 
8188c94862eSPeter Brune #undef __FUNCT__
8198c94862eSPeter Brune #define __FUNCT__ "SNESFASCycle_Full"
8208c94862eSPeter Brune PetscErrorCode SNESFASCycle_Full(SNES snes, Vec X)
8218c94862eSPeter Brune {
8228c94862eSPeter Brune   PetscErrorCode ierr;
8238c94862eSPeter Brune   Vec            F,B;
8248c94862eSPeter Brune   SNES_FAS       *fas = (SNES_FAS*)snes->data;
8258c94862eSPeter Brune   PetscBool      isFine;
8268c94862eSPeter Brune   SNES           next;
8278c94862eSPeter Brune 
8288c94862eSPeter Brune   PetscFunctionBegin;
8298c94862eSPeter Brune   F = snes->vec_func;
8308c94862eSPeter Brune   B = snes->vec_rhs;
8318c94862eSPeter Brune   ierr = SNESFASCycleIsFine(snes,&isFine);CHKERRQ(ierr);
8328c94862eSPeter Brune   ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr);
8338c94862eSPeter Brune 
8348c94862eSPeter Brune   if (isFine) {
8358c94862eSPeter Brune     ierr = SNESFASCycleSetupPhase_Full(snes);CHKERRQ(ierr);
8368c94862eSPeter Brune   }
8378c94862eSPeter Brune 
8388c94862eSPeter Brune   if (fas->full_stage == 0) {
8398c94862eSPeter Brune     /* upsweep */
8408c94862eSPeter Brune     if (next) {
8418c94862eSPeter Brune       if (fas->level != 1) next->max_its += 1;
8428c94862eSPeter Brune       ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr);
8438c94862eSPeter Brune       ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
8448c94862eSPeter Brune       if (fas->level != 1) next->max_its -= 1;
8458c94862eSPeter Brune     } else {
8468c94862eSPeter Brune       ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
8478c94862eSPeter Brune     }
8488c94862eSPeter Brune     fas->full_stage = 1;
8498c94862eSPeter Brune   } else if (fas->full_stage == 1) {
8508c94862eSPeter Brune     if (snes->iter == 0) {ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);}
8518c94862eSPeter Brune     if (next) {
8528c94862eSPeter Brune       ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr);
8538c94862eSPeter Brune       ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
8548c94862eSPeter Brune     }
8558c94862eSPeter Brune   }
8568c94862eSPeter Brune   /* final v-cycle */
8578c94862eSPeter Brune   if (isFine) {
8588c94862eSPeter Brune     if (next) {
8598c94862eSPeter Brune       ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
8608c94862eSPeter Brune       ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr);
8618c94862eSPeter Brune       ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
8628c94862eSPeter Brune     }
8638c94862eSPeter Brune   }
8648c94862eSPeter Brune   PetscFunctionReturn(0);
8658c94862eSPeter Brune }
8668c94862eSPeter Brune 
8678c94862eSPeter Brune #undef __FUNCT__
868*34d65b3cSPeter Brune #define __FUNCT__ "SNESFASCycle_Kaskade"
869*34d65b3cSPeter Brune PetscErrorCode SNESFASCycle_Kaskade(SNES snes, Vec X)
870*34d65b3cSPeter Brune {
871*34d65b3cSPeter Brune 
872*34d65b3cSPeter Brune   PetscErrorCode ierr;
873*34d65b3cSPeter Brune   Vec            F,B;
874*34d65b3cSPeter Brune   SNES           next;
875*34d65b3cSPeter Brune 
876*34d65b3cSPeter Brune   PetscFunctionBegin;
877*34d65b3cSPeter Brune   F = snes->vec_func;
878*34d65b3cSPeter Brune   B = snes->vec_rhs;
879*34d65b3cSPeter Brune   ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr);
880*34d65b3cSPeter Brune   if (next) {
881*34d65b3cSPeter Brune     ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr);
882*34d65b3cSPeter Brune     ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
883*34d65b3cSPeter Brune   } else {
884*34d65b3cSPeter Brune     ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
885*34d65b3cSPeter Brune   }
886*34d65b3cSPeter Brune   PetscFunctionReturn(0);
887*34d65b3cSPeter Brune }
888*34d65b3cSPeter Brune 
889*34d65b3cSPeter Brune #undef __FUNCT__
890421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS"
891421d9b32SPeter Brune 
892421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes)
893421d9b32SPeter Brune {
894fa9694d7SPeter Brune   PetscErrorCode ierr;
895fe6f9142SPeter Brune   PetscInt       i, maxits;
896ddb5aff1SPeter Brune   Vec            X, F;
897fe6f9142SPeter Brune   PetscReal      fnorm;
898b17ce1afSJed Brown   SNES_FAS       *fas = (SNES_FAS*)snes->data,*ffas;
899b17ce1afSJed Brown   DM             dm;
900e70c42e5SPeter Brune   PetscBool      isFine;
901b17ce1afSJed Brown 
902421d9b32SPeter Brune   PetscFunctionBegin;
903fe6f9142SPeter Brune   maxits       = snes->max_its;      /* maximum number of iterations */
904fe6f9142SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
905fa9694d7SPeter Brune   X            = snes->vec_sol;
906f5a6d4f9SBarry Smith   F            = snes->vec_func;
907293a7e31SPeter Brune 
90818a66777SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
909293a7e31SPeter Brune   /*norm setup */
910ce8c27fbSBarry Smith   ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
911fe6f9142SPeter Brune   snes->iter = 0;
912fe6f9142SPeter Brune   snes->norm = 0.;
913ce8c27fbSBarry Smith   ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
914e4ed7901SPeter Brune   if (!snes->vec_func_init_set) {
9150dd27c6cSPeter Brune     if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);}
916fe6f9142SPeter Brune     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
9170dd27c6cSPeter Brune     if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);}
918fe6f9142SPeter Brune     if (snes->domainerror) {
919fe6f9142SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
920fe6f9142SPeter Brune       PetscFunctionReturn(0);
921fe6f9142SPeter Brune     }
9221aa26658SKarl Rupp   } else snes->vec_func_init_set = PETSC_FALSE;
923e4ed7901SPeter Brune 
924fe6f9142SPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
925189a9710SBarry Smith   if (PetscIsInfOrNanReal(fnorm)) {
926189a9710SBarry Smith     snes->reason = SNES_DIVERGED_FNORM_NAN;
927189a9710SBarry Smith     PetscFunctionReturn(0);
928189a9710SBarry Smith   }
929e4ed7901SPeter Brune 
9306c4a5ce2SBarry Smith   ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
931fe6f9142SPeter Brune   snes->norm = fnorm;
932ce8c27fbSBarry Smith   ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
933a71f0d7dSBarry Smith   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
934fe6f9142SPeter Brune   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
935fe6f9142SPeter Brune 
936fe6f9142SPeter Brune   /* test convergence */
937fe6f9142SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
938fe6f9142SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
939e4ed7901SPeter Brune 
940b17ce1afSJed Brown 
941b9c2fdf1SPeter Brune   if (isFine) {
942b9c2fdf1SPeter Brune     /* propagate scale-dependent data up the hierarchy */
943b17ce1afSJed Brown     ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
944b17ce1afSJed Brown     for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) {
945b17ce1afSJed Brown       DM dmcoarse;
946b17ce1afSJed Brown       ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr);
947b17ce1afSJed Brown       ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr);
948b17ce1afSJed Brown       dm   = dmcoarse;
949b17ce1afSJed Brown     }
950b9c2fdf1SPeter Brune   }
951b17ce1afSJed Brown 
952fe6f9142SPeter Brune   for (i = 0; i < maxits; i++) {
953fe6f9142SPeter Brune     /* Call general purpose update function */
954646217ecSPeter Brune 
955fe6f9142SPeter Brune     if (snes->ops->update) {
956fe6f9142SPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
957fe6f9142SPeter Brune     }
95807144faaSPeter Brune     if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
95991f99d7cSPeter Brune       ierr = SNESFASCycle_Multiplicative(snes, X);CHKERRQ(ierr);
9608c94862eSPeter Brune     } else if (fas->fastype == SNES_FAS_ADDITIVE) {
96191f99d7cSPeter Brune       ierr = SNESFASCycle_Additive(snes, X);CHKERRQ(ierr);
9628c94862eSPeter Brune     } else if (fas->fastype == SNES_FAS_FULL) {
9638c94862eSPeter Brune       ierr = SNESFASCycle_Full(snes, X);CHKERRQ(ierr);
964*34d65b3cSPeter Brune     } else if (fas->fastype ==SNES_FAS_KASKADE) {
965*34d65b3cSPeter Brune       ierr = SNESFASCycle_Kaskade(snes, X);CHKERRQ(ierr);
9668c94862eSPeter Brune     } else {
9678c94862eSPeter Brune       SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported FAS type");CHKERRQ(ierr);
96807144faaSPeter Brune     }
969742fe5e2SPeter Brune 
970742fe5e2SPeter Brune     /* check for FAS cycle divergence */
9711aa26658SKarl Rupp     if (snes->reason != SNES_CONVERGED_ITERATING) PetscFunctionReturn(0);
972b9c2fdf1SPeter Brune 
973c90fad12SPeter Brune     /* Monitor convergence */
974ce8c27fbSBarry Smith     ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
975c90fad12SPeter Brune     snes->iter = i+1;
976ce8c27fbSBarry Smith     ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
977a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
978c90fad12SPeter Brune     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
979c90fad12SPeter Brune     /* Test for convergence */
98066585501SPeter Brune     if (isFine) {
981b9c2fdf1SPeter Brune       ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,snes->norm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
982c90fad12SPeter Brune       if (snes->reason) break;
983fe6f9142SPeter Brune     }
98466585501SPeter Brune   }
985fe6f9142SPeter Brune   if (i == maxits) {
986fe6f9142SPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
987fe6f9142SPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
988fe6f9142SPeter Brune   }
989421d9b32SPeter Brune   PetscFunctionReturn(0);
990421d9b32SPeter Brune }
991