xref: /petsc/src/snes/impls/fas/fas.c (revision 5fd668637986a8d8518383a9159eebc368e1d5b4)
1421d9b32SPeter Brune /* Defines the basic SNES object */
26038b46bSPeter Brune #include <../src/snes/impls/fas/fasimpls.h>    /*I  "petscsnesfas.h"  I*/
3421d9b32SPeter Brune 
46a6fc655SJed Brown const char *const SNESFASTypes[] = {"MULTIPLICATIVE","ADDITIVE","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"
491fbfccc6SJed Brown PETSC_EXTERN_C 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);
71421d9b32SPeter Brune   snes->data                  = (void*) fas;
72421d9b32SPeter Brune   fas->level                  = 0;
73293a7e31SPeter Brune   fas->levels                 = 1;
74ee78dd50SPeter Brune   fas->n_cycles               = 1;
75ee78dd50SPeter Brune   fas->max_up_it              = 1;
76ee78dd50SPeter Brune   fas->max_down_it            = 1;
77ab8d36c9SPeter Brune   fas->smoothu                = PETSC_NULL;
78ab8d36c9SPeter Brune   fas->smoothd                = PETSC_NULL;
79421d9b32SPeter Brune   fas->next                   = PETSC_NULL;
806273346dSPeter Brune   fas->previous               = PETSC_NULL;
81ab8d36c9SPeter Brune   fas->fine                   = snes;
82421d9b32SPeter Brune   fas->interpolate            = PETSC_NULL;
83421d9b32SPeter Brune   fas->restrct                = PETSC_NULL;
84efe1f98aSPeter Brune   fas->inject                 = PETSC_NULL;
85cc05f883SPeter Brune   fas->monitor                = PETSC_NULL;
86cc05f883SPeter Brune   fas->usedmfornumberoflevels = PETSC_FALSE;
87ddebd997SPeter Brune   fas->fastype                = SNES_FAS_MULTIPLICATIVE;
880dd27c6cSPeter Brune 
890dd27c6cSPeter Brune   fas->eventsmoothsetup       = PETSC_FALSE;
900dd27c6cSPeter Brune   fas->eventsmoothsolve       = PETSC_FALSE;
910dd27c6cSPeter Brune   fas->eventresidual          = PETSC_FALSE;
920dd27c6cSPeter Brune   fas->eventinterprestrict    = PETSC_FALSE;
930dd27c6cSPeter Brune 
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);
11222c1e704SPeter Brune 
113421d9b32SPeter Brune   PetscFunctionReturn(0);
114421d9b32SPeter Brune }
115421d9b32SPeter Brune 
116421d9b32SPeter Brune #undef __FUNCT__
117421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS"
118421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes)
119421d9b32SPeter Brune {
120421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
121742fe5e2SPeter Brune   PetscErrorCode ierr = 0;
122421d9b32SPeter Brune 
123421d9b32SPeter Brune   PetscFunctionBegin;
124421d9b32SPeter Brune   /* recursively resets and then destroys */
12579d9a41aSPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
126421d9b32SPeter Brune   if (fas->next)         ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
127421d9b32SPeter Brune   ierr = PetscFree(fas);CHKERRQ(ierr);
128ee1fd11aSPeter Brune 
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;
1603dccd265SPeter Brune       snes->vec_sol = PETSC_NULL;
1613dccd265SPeter Brune       snes->vec_func = PETSC_NULL;
1623dccd265SPeter Brune       snes->vec_sol_update = PETSC_NULL;
1633dccd265SPeter Brune       snes->vec_rhs = PETSC_NULL;
1643dccd265SPeter Brune 
1653dccd265SPeter Brune       /* reset the number of levels */
166d1adcc6fSPeter Brune       ierr = SNESFASSetLevels(snes,dm_levels,PETSC_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 
17807144faaSPeter Brune   if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
179e4ed7901SPeter Brune     ierr = SNESDefaultGetWork(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
18007144faaSPeter Brune   } else {
181e7f468e7SPeter Brune     ierr = SNESDefaultGetWork(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
18207144faaSPeter Brune   }
183cc05f883SPeter Brune 
184ab8d36c9SPeter Brune   /* set up the smoothers if they haven't already been set up */
185ab8d36c9SPeter Brune   if (!fas->smoothd) {
186ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
187ab8d36c9SPeter Brune   }
188ab8d36c9SPeter Brune 
18979d9a41aSPeter Brune   if (snes->dm) {
190ab8d36c9SPeter Brune     /* set the smoother DMs properly */
191ab8d36c9SPeter Brune     if (fas->smoothu) ierr = SNESSetDM(fas->smoothu, snes->dm);CHKERRQ(ierr);
192ab8d36c9SPeter Brune     ierr = SNESSetDM(fas->smoothd, snes->dm);CHKERRQ(ierr);
19379d9a41aSPeter Brune     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
194ab8d36c9SPeter Brune     if (next) {
19579d9a41aSPeter Brune       /* for now -- assume the DM and the evaluation functions have been set externally */
196ab8d36c9SPeter Brune       if (!next->dm) {
197ab8d36c9SPeter Brune         ierr = DMCoarsen(snes->dm, ((PetscObject)next)->comm, &next->dm);CHKERRQ(ierr);
198ab8d36c9SPeter Brune         ierr = SNESSetDM(next, next->dm);CHKERRQ(ierr);
19979d9a41aSPeter Brune       }
20079d9a41aSPeter Brune       /* set the interpolation and restriction from the DM */
20179d9a41aSPeter Brune       if (!fas->interpolate) {
202ab8d36c9SPeter Brune         ierr = DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
203bccf9bb3SJed Brown         if (!fas->restrct) {
204bccf9bb3SJed Brown           ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr);
20579d9a41aSPeter Brune           fas->restrct = fas->interpolate;
20679d9a41aSPeter Brune         }
207bccf9bb3SJed Brown       }
20879d9a41aSPeter Brune       /* set the injection from the DM */
20979d9a41aSPeter Brune       if (!fas->inject) {
210ab8d36c9SPeter Brune         ierr = DMCreateInjection(next->dm, snes->dm, &injscatter);CHKERRQ(ierr);
21179d9a41aSPeter Brune         ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr);
21279d9a41aSPeter Brune         ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr);
21379d9a41aSPeter Brune       }
21479d9a41aSPeter Brune     }
21579d9a41aSPeter Brune   }
21679d9a41aSPeter Brune   /*pass the smoother, function, and jacobian up to the next level if it's not user set already */
21779d9a41aSPeter Brune   if (fas->galerkin) {
218ab8d36c9SPeter Brune     if (next)
219ab8d36c9SPeter Brune       ierr = SNESSetFunction(next, PETSC_NULL, SNESFASGalerkinDefaultFunction, next);CHKERRQ(ierr);
220ab8d36c9SPeter Brune     if (fas->smoothd && fas->level != fas->levels - 1) ierr = SNESSetFunction(fas->smoothd, PETSC_NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr);
221ab8d36c9SPeter Brune     if (fas->smoothu && fas->level != fas->levels - 1) ierr = SNESSetFunction(fas->smoothu, PETSC_NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr);
22279d9a41aSPeter Brune   }
22379d9a41aSPeter Brune 
224534ebe21SPeter Brune   /* sets the down (pre) smoother's default norm and sets it from options */
225534ebe21SPeter Brune   if (fas->smoothd) {
226bc3f2f05SPeter Brune     if (fas->level == 0 && fas->levels != 1) {
227534ebe21SPeter Brune       ierr = SNESSetNormType(fas->smoothd, SNES_NORM_NONE);CHKERRQ(ierr);
228534ebe21SPeter Brune      } else {
229534ebe21SPeter Brune       ierr = SNESSetNormType(fas->smoothd, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr);
230534ebe21SPeter Brune     }
2317fce8c19SPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)fas->smoothd);CHKERRQ(ierr);
232534ebe21SPeter Brune     ierr = SNESSetFromOptions(fas->smoothd);CHKERRQ(ierr);
233f89ba88eSPeter Brune     ierr = SNESGetSNESLineSearch(snes,&linesearch);CHKERRQ(ierr);
234f89ba88eSPeter Brune     ierr = SNESGetSNESLineSearch(fas->smoothd,&slinesearch);CHKERRQ(ierr);
2356b2b7091SBarry Smith     ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr);
2366b2b7091SBarry Smith     ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr);
2376b2b7091SBarry Smith     ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr);
2386b2b7091SBarry Smith     ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr);
239f89ba88eSPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr);
2400dd27c6cSPeter Brune 
2410dd27c6cSPeter Brune     fas->smoothd->vec_sol = snes->vec_sol;
2420dd27c6cSPeter Brune     ierr = PetscObjectReference((PetscObject)snes->vec_sol);CHKERRQ(ierr);
2430dd27c6cSPeter Brune     fas->smoothd->vec_sol_update = snes->vec_sol_update;
2440dd27c6cSPeter Brune     ierr = PetscObjectReference((PetscObject)snes->vec_sol_update);CHKERRQ(ierr);
2450dd27c6cSPeter Brune     fas->smoothd->vec_func = snes->vec_func;
2460dd27c6cSPeter Brune     ierr = PetscObjectReference((PetscObject)snes->vec_func);CHKERRQ(ierr);
2470dd27c6cSPeter Brune 
2480dd27c6cSPeter Brune     if (fas->eventsmoothsetup) {ierr = PetscLogEventBegin(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
2490dd27c6cSPeter Brune     ierr = SNESSetUp(fas->smoothd);CHKERRQ(ierr);
2500dd27c6cSPeter Brune     if (fas->eventsmoothsetup) {ierr = PetscLogEventEnd(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
251534ebe21SPeter Brune   }
252534ebe21SPeter Brune 
253534ebe21SPeter Brune   /* sets the up (post) smoother's default norm and sets it from options */
254534ebe21SPeter Brune   if (fas->smoothu) {
255534ebe21SPeter Brune     if (fas->level != fas->levels - 1) {
256534ebe21SPeter Brune       ierr = SNESSetNormType(fas->smoothu, SNES_NORM_NONE);CHKERRQ(ierr);
257534ebe21SPeter Brune     } else {
258534ebe21SPeter Brune       ierr = SNESSetNormType(fas->smoothu, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr);
259534ebe21SPeter Brune     }
2607fce8c19SPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)fas->smoothu);CHKERRQ(ierr);
261534ebe21SPeter Brune     ierr = SNESSetFromOptions(fas->smoothu);CHKERRQ(ierr);
262f89ba88eSPeter Brune     ierr = SNESGetSNESLineSearch(snes,&linesearch);CHKERRQ(ierr);
263f89ba88eSPeter Brune     ierr = SNESGetSNESLineSearch(fas->smoothu,&slinesearch);CHKERRQ(ierr);
2646b2b7091SBarry Smith     ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr);
2656b2b7091SBarry Smith     ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr);
2666b2b7091SBarry Smith     ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr);
2676b2b7091SBarry Smith     ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr);
268f89ba88eSPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr);
2690dd27c6cSPeter Brune 
2700dd27c6cSPeter Brune     fas->smoothu->vec_sol = snes->vec_sol;
2710dd27c6cSPeter Brune     ierr = PetscObjectReference((PetscObject)snes->vec_sol);CHKERRQ(ierr);
2720dd27c6cSPeter Brune     fas->smoothu->vec_sol_update = snes->vec_sol_update;
2730dd27c6cSPeter Brune     ierr = PetscObjectReference((PetscObject)snes->vec_sol_update);CHKERRQ(ierr);
2740dd27c6cSPeter Brune     fas->smoothu->vec_func = snes->vec_func;
2750dd27c6cSPeter Brune     ierr = PetscObjectReference((PetscObject)snes->vec_func);CHKERRQ(ierr);
2760dd27c6cSPeter Brune 
2770dd27c6cSPeter Brune     if (fas->eventsmoothsetup) {ierr = PetscLogEventBegin(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
2780dd27c6cSPeter Brune     ierr = SNESSetUp(fas->smoothu);CHKERRQ(ierr);
2790dd27c6cSPeter Brune     if (fas->eventsmoothsetup) {ierr = PetscLogEventEnd(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
2800dd27c6cSPeter Brune 
281534ebe21SPeter Brune   }
282d06165b7SPeter Brune 
283ab8d36c9SPeter Brune   if (next) {
28479d9a41aSPeter Brune     /* gotta set up the solution vector for this to work */
285ab8d36c9SPeter Brune     if (!next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_sol);CHKERRQ(ierr);}
286ab8d36c9SPeter Brune     if (!next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_rhs);CHKERRQ(ierr);}
2877fce8c19SPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)next);CHKERRQ(ierr);
288f89ba88eSPeter Brune     ierr = SNESGetSNESLineSearch(snes,&linesearch);CHKERRQ(ierr);
289f89ba88eSPeter Brune     ierr = SNESGetSNESLineSearch(fas->next,&slinesearch);CHKERRQ(ierr);
2906b2b7091SBarry Smith     ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr);
2916b2b7091SBarry Smith     ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr);
2926b2b7091SBarry Smith     ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr);
2936b2b7091SBarry Smith     ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr);
294f89ba88eSPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr);
295ab8d36c9SPeter Brune     ierr = SNESSetUp(next);CHKERRQ(ierr);
29679d9a41aSPeter Brune   }
2976273346dSPeter Brune   /* setup FAS work vectors */
2986273346dSPeter Brune   if (fas->galerkin) {
2996273346dSPeter Brune     ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr);
3006273346dSPeter Brune     ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr);
3016273346dSPeter Brune   }
302421d9b32SPeter Brune   PetscFunctionReturn(0);
303421d9b32SPeter Brune }
304421d9b32SPeter Brune 
305421d9b32SPeter Brune #undef __FUNCT__
306421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS"
307421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
308421d9b32SPeter Brune {
309ee78dd50SPeter Brune   SNES_FAS       *fas = (SNES_FAS *) snes->data;
310ee78dd50SPeter Brune   PetscInt       levels = 1;
3114d26bfa5SPeter Brune   PetscBool      flg = PETSC_FALSE, upflg = PETSC_FALSE, downflg = PETSC_FALSE, monflg = PETSC_FALSE, galerkinflg = PETSC_FALSE;
312421d9b32SPeter Brune   PetscErrorCode ierr;
313ee78dd50SPeter Brune   char           monfilename[PETSC_MAX_PATH_LEN];
31407144faaSPeter Brune   SNESFASType    fastype;
315fde0ff24SPeter Brune   const char     *optionsprefix;
316f1c6b773SPeter Brune   SNESLineSearch linesearch;
31766585501SPeter Brune   PetscInt       m, n_up, n_down;
318ab8d36c9SPeter Brune   SNES           next;
319ab8d36c9SPeter Brune   PetscBool      isFine;
320421d9b32SPeter Brune 
321421d9b32SPeter Brune   PetscFunctionBegin;
322ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
323c90fad12SPeter Brune   ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr);
324ee78dd50SPeter Brune 
325ab8d36c9SPeter Brune   /* number of levels -- only process most options on the finest level */
326ab8d36c9SPeter Brune   if (isFine) {
327ee78dd50SPeter Brune     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
328c732cbdbSBarry Smith     if (!flg && snes->dm) {
329c732cbdbSBarry Smith       ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr);
330c732cbdbSBarry Smith       levels++;
331d1adcc6fSPeter Brune       fas->usedmfornumberoflevels = PETSC_TRUE;
332c732cbdbSBarry Smith     }
333ee78dd50SPeter Brune     ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr);
33407144faaSPeter Brune     fastype = fas->fastype;
33507144faaSPeter Brune     ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr);
33607144faaSPeter Brune     if (flg) {
33707144faaSPeter Brune       ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr);
33807144faaSPeter Brune     }
339ee78dd50SPeter Brune 
340fde0ff24SPeter Brune     ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
341ab8d36c9SPeter Brune     ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&m,&flg);CHKERRQ(ierr);
342ab8d36c9SPeter Brune     if (flg) {
343ab8d36c9SPeter Brune       ierr = SNESFASSetCycles(snes, m);CHKERRQ(ierr);
344fde0ff24SPeter Brune     }
345fde0ff24SPeter Brune 
346ab8d36c9SPeter Brune     ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFASSetGalerkin",fas->galerkin,&galerkinflg,&flg);CHKERRQ(ierr);
347ab8d36c9SPeter Brune     if (flg) {
348ab8d36c9SPeter Brune       ierr = SNESFASSetGalerkin(snes, galerkinflg);CHKERRQ(ierr);
349ab8d36c9SPeter Brune     }
350ee78dd50SPeter Brune 
35166585501SPeter Brune     ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",fas->max_up_it,&n_up,&upflg);CHKERRQ(ierr);
352162d76ddSPeter Brune 
35366585501SPeter Brune     ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",fas->max_down_it,&n_down,&downflg);CHKERRQ(ierr);
354162d76ddSPeter Brune 
355c8c899caSPeter Brune     ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESFASSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr);
356c8c899caSPeter Brune     if (monflg) ierr = SNESFASSetMonitor(snes, PETSC_TRUE);CHKERRQ(ierr);
3570dd27c6cSPeter Brune 
3580dd27c6cSPeter Brune     flg = PETSC_FALSE;
3590dd27c6cSPeter Brune     monflg = PETSC_TRUE;
3600dd27c6cSPeter Brune     ierr = PetscOptionsBool("-snes_fas_log","Log times for each FAS level","SNESFASSetLog",monflg,&monflg,&flg);CHKERRQ(ierr);
3610dd27c6cSPeter Brune     if (flg){ierr = SNESFASSetLog(snes,monflg);CHKERRQ(ierr);}
362ab8d36c9SPeter Brune   }
363ee78dd50SPeter Brune 
364421d9b32SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
3658cc86e31SPeter Brune   /* setup from the determined types if there is no pointwise procedure or smoother defined */
366162d76ddSPeter Brune   if (upflg) {
36766585501SPeter Brune     ierr = SNESFASSetNumberSmoothUp(snes,n_up);CHKERRQ(ierr);
368162d76ddSPeter Brune   }
369162d76ddSPeter Brune   if (downflg) {
37066585501SPeter Brune     ierr = SNESFASSetNumberSmoothDown(snes,n_down);CHKERRQ(ierr);
371162d76ddSPeter Brune   }
372eff52c0eSPeter Brune 
3739e764e56SPeter Brune   /* set up the default line search for coarse grid corrections */
3749e764e56SPeter Brune   if (fas->fastype == SNES_FAS_ADDITIVE) {
3759e764e56SPeter Brune     if (!snes->linesearch) {
376f1c6b773SPeter Brune       ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr);
3771a4f838cSPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr);
3789e764e56SPeter Brune     }
3799e764e56SPeter Brune   }
3809e764e56SPeter Brune 
381ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
382ee78dd50SPeter Brune   /* recursive option setting for the smoothers */
383ab8d36c9SPeter Brune   if (next) {ierr = SNESSetFromOptions(next);CHKERRQ(ierr);}
384421d9b32SPeter Brune   PetscFunctionReturn(0);
385421d9b32SPeter Brune }
386421d9b32SPeter Brune 
387421d9b32SPeter Brune #undef __FUNCT__
388421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS"
389421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
390421d9b32SPeter Brune {
391421d9b32SPeter Brune   SNES_FAS       *fas = (SNES_FAS *) snes->data;
392656ede7eSPeter Brune   PetscBool      isFine,iascii,isdraw;
393ab8d36c9SPeter Brune   PetscInt       i;
394421d9b32SPeter Brune   PetscErrorCode ierr;
395ab8d36c9SPeter Brune   SNES           smoothu, smoothd, levelsnes;
396421d9b32SPeter Brune 
397421d9b32SPeter Brune   PetscFunctionBegin;
398ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
399ab8d36c9SPeter Brune   if (isFine) {
400251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
401656ede7eSPeter Brune     ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
402421d9b32SPeter Brune     if (iascii) {
403ab8d36c9SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "FAS: type is %s, levels=%D, cycles=%D\n",  SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr);
404ab8d36c9SPeter Brune       if (fas->galerkin) {
405ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
406421d9b32SPeter Brune       } else {
407ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
408421d9b32SPeter Brune       }
409ab8d36c9SPeter Brune       for (i=0; i<fas->levels; i++) {
410ab8d36c9SPeter Brune         ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr);
411ab8d36c9SPeter Brune         ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr);
412ab8d36c9SPeter Brune         ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr);
413ab8d36c9SPeter Brune         if (!i) {
414ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr);
415421d9b32SPeter Brune         } else {
416ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
417421d9b32SPeter Brune         }
418ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
419ab8d36c9SPeter Brune         ierr = SNESView(smoothd,viewer);CHKERRQ(ierr);
420ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
421ab8d36c9SPeter Brune         if (i && (smoothd == smoothu)) {
422ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
423ab8d36c9SPeter Brune         } else if (i) {
424ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
425ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
426ab8d36c9SPeter Brune           ierr = SNESView(smoothu,viewer);CHKERRQ(ierr);
427ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
428ab8d36c9SPeter Brune         }
429ab8d36c9SPeter Brune       }
430656ede7eSPeter Brune     } else if (isdraw) {
431656ede7eSPeter Brune       PetscDraw draw;
432b4375e8dSPeter Brune       PetscReal x,w,y,bottom,th,wth;
433656ede7eSPeter Brune       SNES_FAS *curfas = fas;
434656ede7eSPeter Brune       ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
435656ede7eSPeter Brune       ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
436656ede7eSPeter Brune       ierr = PetscDrawStringGetSize(draw,&wth,&th);CHKERRQ(ierr);
437656ede7eSPeter Brune       bottom = y - th;
438656ede7eSPeter Brune       while (curfas) {
439b4375e8dSPeter Brune         if (!curfas->smoothu) {
440656ede7eSPeter Brune           ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
441656ede7eSPeter Brune           if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr);
442656ede7eSPeter Brune           ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
443b4375e8dSPeter Brune         } else {
444b4375e8dSPeter Brune           w = 0.5*PetscMin(1.0-x,x);
445b4375e8dSPeter Brune           ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
446b4375e8dSPeter Brune           if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr);
447b4375e8dSPeter Brune           ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
448b4375e8dSPeter Brune           ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
449b4375e8dSPeter Brune           if (curfas->smoothu) ierr = SNESView(curfas->smoothu,viewer);CHKERRQ(ierr);
450b4375e8dSPeter Brune           ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
451b4375e8dSPeter Brune         }
452656ede7eSPeter Brune         /* this is totally bogus but we have no way of knowing how low the previous one was draw to */
453656ede7eSPeter Brune         bottom -= 5*th;
454656ede7eSPeter Brune         if (curfas->next) {
455656ede7eSPeter Brune           curfas = (SNES_FAS*)curfas->next->data;
456656ede7eSPeter Brune         } else {
457656ede7eSPeter Brune           curfas = PETSC_NULL;
458656ede7eSPeter Brune         }
459656ede7eSPeter Brune       }
460421d9b32SPeter Brune     }
461ab8d36c9SPeter Brune   }
462421d9b32SPeter Brune   PetscFunctionReturn(0);
463421d9b32SPeter Brune }
464421d9b32SPeter Brune 
465421d9b32SPeter Brune #undef __FUNCT__
46691f99d7cSPeter Brune #define __FUNCT__ "SNESFASDownSmooth_Private"
46739bd7f45SPeter Brune /*
46839bd7f45SPeter Brune Defines the action of the downsmoother
46939bd7f45SPeter Brune  */
47091f99d7cSPeter Brune PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
471b9c2fdf1SPeter Brune {
47239bd7f45SPeter Brune   PetscErrorCode      ierr = 0;
473742fe5e2SPeter Brune   SNESConvergedReason reason;
474ab8d36c9SPeter Brune   Vec                 FPC;
475ab8d36c9SPeter Brune   SNES                smoothd;
4760dd27c6cSPeter Brune   SNES_FAS            *fas = (SNES_FAS *) snes->data;
4776e111a19SKarl Rupp 
478421d9b32SPeter Brune   PetscFunctionBegin;
479ab8d36c9SPeter Brune   ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr);
480e4ed7901SPeter Brune   ierr = SNESSetInitialFunction(smoothd, F);CHKERRQ(ierr);
481e4ed7901SPeter Brune   ierr = SNESSetInitialFunctionNorm(smoothd, *fnorm);CHKERRQ(ierr);
4820dd27c6cSPeter Brune   if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
483ab8d36c9SPeter Brune   ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr);
4840dd27c6cSPeter Brune   if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
485742fe5e2SPeter Brune   /* check convergence reason for the smoother */
486ab8d36c9SPeter Brune   ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr);
487e70c42e5SPeter Brune   if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN)) {
488742fe5e2SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
489742fe5e2SPeter Brune     PetscFunctionReturn(0);
490742fe5e2SPeter Brune   }
491ab8d36c9SPeter Brune   ierr = SNESGetFunction(smoothd, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
4924b32a720SPeter Brune   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
493b9c2fdf1SPeter Brune   ierr = SNESGetFunctionNorm(smoothd, fnorm);CHKERRQ(ierr);
49439bd7f45SPeter Brune   PetscFunctionReturn(0);
49539bd7f45SPeter Brune }
49639bd7f45SPeter Brune 
49739bd7f45SPeter Brune 
49839bd7f45SPeter Brune #undef __FUNCT__
49991f99d7cSPeter Brune #define __FUNCT__ "SNESFASUpSmooth_Private"
50039bd7f45SPeter Brune /*
50107144faaSPeter Brune Defines the action of the upsmoother
50239bd7f45SPeter Brune  */
5030adebc6cSBarry Smith PetscErrorCode SNESFASUpSmooth_Private (SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
5040adebc6cSBarry Smith {
50539bd7f45SPeter Brune   PetscErrorCode      ierr = 0;
50639bd7f45SPeter Brune   SNESConvergedReason reason;
507ab8d36c9SPeter Brune   Vec                 FPC;
508ab8d36c9SPeter Brune   SNES                smoothu;
5090dd27c6cSPeter Brune   SNES_FAS            *fas = (SNES_FAS *) snes->data;
510ab8d36c9SPeter Brune 
5116e111a19SKarl Rupp   PetscFunctionBegin;
512ab8d36c9SPeter Brune   ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr);
5130dd27c6cSPeter Brune   if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
514ab8d36c9SPeter Brune   ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr);
5150dd27c6cSPeter Brune   if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
51639bd7f45SPeter Brune   /* check convergence reason for the smoother */
517ab8d36c9SPeter Brune   ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr);
51839bd7f45SPeter Brune   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
51939bd7f45SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
52039bd7f45SPeter Brune     PetscFunctionReturn(0);
52139bd7f45SPeter Brune   }
522ab8d36c9SPeter Brune   ierr = SNESGetFunction(smoothu, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
5234b32a720SPeter Brune   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
524b9c2fdf1SPeter Brune   ierr = SNESGetFunctionNorm(smoothu, fnorm);CHKERRQ(ierr);
52539bd7f45SPeter Brune   PetscFunctionReturn(0);
52639bd7f45SPeter Brune }
52739bd7f45SPeter Brune 
52839bd7f45SPeter Brune #undef __FUNCT__
529938e4a01SJed Brown #define __FUNCT__ "SNESFASCreateCoarseVec"
530938e4a01SJed Brown /*@
531938e4a01SJed Brown    SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level
532938e4a01SJed Brown 
533938e4a01SJed Brown    Collective
534938e4a01SJed Brown 
535938e4a01SJed Brown    Input Arguments:
536938e4a01SJed Brown .  snes - SNESFAS
537938e4a01SJed Brown 
538938e4a01SJed Brown    Output Arguments:
539938e4a01SJed Brown .  Xcoarse - vector on level one coarser than snes
540938e4a01SJed Brown 
541938e4a01SJed Brown    Level: developer
542938e4a01SJed Brown 
543938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict()
544938e4a01SJed Brown @*/
545938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse)
546938e4a01SJed Brown {
547938e4a01SJed Brown   PetscErrorCode ierr;
548938e4a01SJed Brown   SNES_FAS       *fas = (SNES_FAS*)snes->data;
549938e4a01SJed Brown 
550938e4a01SJed Brown   PetscFunctionBegin;
551938e4a01SJed Brown   if (fas->rscale) {ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);}
552938e4a01SJed Brown   else if (fas->inject) {ierr = MatGetVecs(fas->inject,Xcoarse,PETSC_NULL);CHKERRQ(ierr);}
553938e4a01SJed Brown   else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr);
554938e4a01SJed Brown   PetscFunctionReturn(0);
555938e4a01SJed Brown }
556938e4a01SJed Brown 
557e9923e8dSJed Brown #undef __FUNCT__
558e9923e8dSJed Brown #define __FUNCT__ "SNESFASRestrict"
559e9923e8dSJed Brown /*@
560e9923e8dSJed Brown    SNESFASRestrict - restrict a Vec to the next coarser level
561e9923e8dSJed Brown 
562e9923e8dSJed Brown    Collective
563e9923e8dSJed Brown 
564e9923e8dSJed Brown    Input Arguments:
565e9923e8dSJed Brown +  fine - SNES from which to restrict
566e9923e8dSJed Brown -  Xfine - vector to restrict
567e9923e8dSJed Brown 
568e9923e8dSJed Brown    Output Arguments:
569e9923e8dSJed Brown .  Xcoarse - result of restriction
570e9923e8dSJed Brown 
571e9923e8dSJed Brown    Level: developer
572e9923e8dSJed Brown 
573e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection()
574e9923e8dSJed Brown @*/
575e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse)
576e9923e8dSJed Brown {
577e9923e8dSJed Brown   PetscErrorCode ierr;
578e9923e8dSJed Brown   SNES_FAS       *fas = (SNES_FAS*)fine->data;
579e9923e8dSJed Brown 
580e9923e8dSJed Brown   PetscFunctionBegin;
581e9923e8dSJed Brown   PetscValidHeaderSpecific(fine,SNES_CLASSID,1);
582e9923e8dSJed Brown   PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2);
583e9923e8dSJed Brown   PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3);
584e9923e8dSJed Brown   if (fas->inject) {
585e9923e8dSJed Brown     ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr);
586e9923e8dSJed Brown   } else {
587e9923e8dSJed Brown     ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr);
588e9923e8dSJed Brown     ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr);
589e9923e8dSJed Brown   }
590e9923e8dSJed Brown   PetscFunctionReturn(0);
591e9923e8dSJed Brown }
592e9923e8dSJed Brown 
593e9923e8dSJed Brown #undef __FUNCT__
5948c40d5fbSBarry Smith #define __FUNCT__ "SNESFASCoarseCorrection"
59539bd7f45SPeter Brune /*
59639bd7f45SPeter Brune 
59739bd7f45SPeter Brune Performs the FAS coarse correction as:
59839bd7f45SPeter Brune 
59939bd7f45SPeter Brune fine problem: F(x) = 0
60039bd7f45SPeter Brune coarse problem: F^c(x) = b^c
60139bd7f45SPeter Brune 
60239bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
60339bd7f45SPeter Brune 
60439bd7f45SPeter Brune  */
6050adebc6cSBarry Smith PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new)
6060adebc6cSBarry Smith {
60739bd7f45SPeter Brune   PetscErrorCode      ierr;
60839bd7f45SPeter Brune   Vec                 X_c, Xo_c, F_c, B_c;
60939bd7f45SPeter Brune   SNESConvergedReason reason;
610ab8d36c9SPeter Brune   SNES                next;
611ab8d36c9SPeter Brune   Mat                 restrct, interpolate;
6120dd27c6cSPeter Brune   SNES_FAS            *fasc;
613*5fd66863SKarl Rupp 
61439bd7f45SPeter Brune   PetscFunctionBegin;
615ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
616ab8d36c9SPeter Brune   if (next) {
6170dd27c6cSPeter Brune     fasc = (SNES_FAS *)next->data;
6180dd27c6cSPeter Brune 
619ab8d36c9SPeter Brune     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
620ab8d36c9SPeter Brune     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
621ab8d36c9SPeter Brune 
622ab8d36c9SPeter Brune     X_c  = next->vec_sol;
623ab8d36c9SPeter Brune     Xo_c = next->work[0];
624ab8d36c9SPeter Brune     F_c  = next->vec_func;
625ab8d36c9SPeter Brune     B_c  = next->vec_rhs;
626efe1f98aSPeter Brune 
6270dd27c6cSPeter Brune     if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
628938e4a01SJed Brown     ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr);
629293a7e31SPeter Brune     /* restrict the defect */
630ab8d36c9SPeter Brune     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
6310dd27c6cSPeter Brune     if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
6320dd27c6cSPeter Brune 
6330dd27c6cSPeter Brune     if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);}
634ab8d36c9SPeter Brune     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
6350dd27c6cSPeter Brune     if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);}
6360dd27c6cSPeter Brune 
6370dd27c6cSPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */
638e4ed7901SPeter Brune     ierr = VecCopy(B_c, X_c);CHKERRQ(ierr);
639b9c2fdf1SPeter Brune     ierr = VecCopy(F_c, B_c);CHKERRQ(ierr);
640e4ed7901SPeter Brune     ierr = VecCopy(X_c, F_c);CHKERRQ(ierr);
641ee78dd50SPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
642ee78dd50SPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
643c90fad12SPeter Brune 
644c90fad12SPeter Brune     /* recurse to the next level */
645e4ed7901SPeter Brune     ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr);
646ab8d36c9SPeter Brune     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
647ab8d36c9SPeter Brune     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
648742fe5e2SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
649742fe5e2SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
650742fe5e2SPeter Brune       PetscFunctionReturn(0);
651742fe5e2SPeter Brune     }
652fa9694d7SPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
653fa9694d7SPeter Brune     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
6540dd27c6cSPeter Brune 
6550dd27c6cSPeter Brune     if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
656ab8d36c9SPeter Brune     ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr);
6570dd27c6cSPeter Brune     if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
658293a7e31SPeter Brune   }
65939bd7f45SPeter Brune   PetscFunctionReturn(0);
66039bd7f45SPeter Brune }
66139bd7f45SPeter Brune 
66239bd7f45SPeter Brune #undef __FUNCT__
6632cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Additive"
66439bd7f45SPeter Brune /*
66539bd7f45SPeter Brune 
66639bd7f45SPeter Brune The additive cycle looks like:
66739bd7f45SPeter Brune 
66807144faaSPeter Brune xhat = x
66907144faaSPeter Brune xhat = dS(x, b)
67007144faaSPeter Brune x = coarsecorrection(xhat, b_d)
67107144faaSPeter Brune x = x + nu*(xhat - x);
67239bd7f45SPeter Brune (optional) x = uS(x, b)
67339bd7f45SPeter Brune 
67439bd7f45SPeter Brune With the coarse RHS (defect correction) as below.
67539bd7f45SPeter Brune 
67639bd7f45SPeter Brune  */
6770adebc6cSBarry Smith PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X)
6780adebc6cSBarry Smith {
67907144faaSPeter Brune   Vec                 F, B, Xhat;
68022c1e704SPeter Brune   Vec                 X_c, Xo_c, F_c, B_c;
68139bd7f45SPeter Brune   PetscErrorCode      ierr;
68207144faaSPeter Brune   SNESConvergedReason reason;
68322c1e704SPeter Brune   PetscReal           xnorm, fnorm, ynorm;
68422c1e704SPeter Brune   PetscBool           lssuccess;
685ab8d36c9SPeter Brune   SNES                next;
686ab8d36c9SPeter Brune   Mat                 restrct, interpolate;
6870dd27c6cSPeter Brune   SNES_FAS            *fas = (SNES_FAS *)snes->data,*fasc;
6880dd27c6cSPeter Brune 
68939bd7f45SPeter Brune   PetscFunctionBegin;
690ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
69139bd7f45SPeter Brune   F = snes->vec_func;
69239bd7f45SPeter Brune   B = snes->vec_rhs;
693e7f468e7SPeter Brune   Xhat = snes->work[1];
69407144faaSPeter Brune   ierr = VecCopy(X, Xhat);CHKERRQ(ierr);
69507144faaSPeter Brune   /* recurse first */
696ab8d36c9SPeter Brune   if (next) {
6970dd27c6cSPeter Brune     fasc = (SNES_FAS*)next->data;
698ab8d36c9SPeter Brune     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
699ab8d36c9SPeter Brune     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
7000dd27c6cSPeter Brune     if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);}
70107144faaSPeter Brune     ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr);
7020dd27c6cSPeter Brune     if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);}
703c2a02606SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
704ab8d36c9SPeter Brune     X_c  = next->vec_sol;
705ab8d36c9SPeter Brune     Xo_c = next->work[0];
706ab8d36c9SPeter Brune     F_c  = next->vec_func;
707ab8d36c9SPeter Brune     B_c  = next->vec_rhs;
70839bd7f45SPeter Brune 
709938e4a01SJed Brown     ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr);
71007144faaSPeter Brune     /* restrict the defect */
711ab8d36c9SPeter Brune     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
71207144faaSPeter Brune 
71307144faaSPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
7140dd27c6cSPeter Brune     if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);}
715ab8d36c9SPeter Brune     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
7160dd27c6cSPeter Brune     if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);}
717e4ed7901SPeter Brune     ierr = VecCopy(B_c, X_c);CHKERRQ(ierr);
718b9c2fdf1SPeter Brune     ierr = VecCopy(F_c, B_c);CHKERRQ(ierr);
719e4ed7901SPeter Brune     ierr = VecCopy(X_c, F_c);CHKERRQ(ierr);
72007144faaSPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
72107144faaSPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
72207144faaSPeter Brune 
72307144faaSPeter Brune     /* recurse */
724e4ed7901SPeter Brune     ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr);
725ab8d36c9SPeter Brune     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
72607144faaSPeter Brune 
72707144faaSPeter Brune     /* smooth on this level */
72891f99d7cSPeter Brune     ierr = SNESFASDownSmooth_Private(snes, B, X, F, &fnorm);CHKERRQ(ierr);
72907144faaSPeter Brune 
730ab8d36c9SPeter Brune     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
73107144faaSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
73207144faaSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
73307144faaSPeter Brune       PetscFunctionReturn(0);
73407144faaSPeter Brune     }
73507144faaSPeter Brune 
73607144faaSPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
737c68acad4SPeter Brune     ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr);
738ab8d36c9SPeter Brune     ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr);
73907144faaSPeter Brune 
740ddebd997SPeter Brune     /* additive correction of the coarse direction*/
741f1c6b773SPeter Brune     ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr);
742f1c6b773SPeter Brune     ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssuccess);CHKERRQ(ierr);
7439e764e56SPeter Brune     if (!lssuccess) {
7449e764e56SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
7459e764e56SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
7469e764e56SPeter Brune         PetscFunctionReturn(0);
7479e764e56SPeter Brune       }
7489e764e56SPeter Brune     }
749b9c2fdf1SPeter Brune     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm);CHKERRQ(ierr);
75007144faaSPeter Brune   } else {
75191f99d7cSPeter Brune     ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
75207144faaSPeter Brune   }
75339bd7f45SPeter Brune   PetscFunctionReturn(0);
75439bd7f45SPeter Brune }
75539bd7f45SPeter Brune 
75639bd7f45SPeter Brune #undef __FUNCT__
7572cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Multiplicative"
75839bd7f45SPeter Brune /*
75939bd7f45SPeter Brune 
76039bd7f45SPeter Brune Defines the FAS cycle as:
76139bd7f45SPeter Brune 
76239bd7f45SPeter Brune fine problem: F(x) = 0
76339bd7f45SPeter Brune coarse problem: F^c(x) = b^c
76439bd7f45SPeter Brune 
76539bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
76639bd7f45SPeter Brune 
76739bd7f45SPeter Brune correction:
76839bd7f45SPeter Brune 
76939bd7f45SPeter Brune x = x + I(x^c - Rx)
77039bd7f45SPeter Brune 
77139bd7f45SPeter Brune  */
7720adebc6cSBarry Smith PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X)
7730adebc6cSBarry Smith {
77439bd7f45SPeter Brune 
77539bd7f45SPeter Brune   PetscErrorCode      ierr;
77639bd7f45SPeter Brune   Vec                 F,B;
77739bd7f45SPeter Brune   SNES_FAS            *fas = (SNES_FAS *)snes->data;
77839bd7f45SPeter Brune 
77939bd7f45SPeter Brune   PetscFunctionBegin;
78039bd7f45SPeter Brune   F = snes->vec_func;
78139bd7f45SPeter Brune   B = snes->vec_rhs;
78239bd7f45SPeter Brune   /* pre-smooth -- just update using the pre-smoother */
78391f99d7cSPeter Brune   ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
784c90fad12SPeter Brune   if (fas->level != 0) {
7858c40d5fbSBarry Smith     ierr = SNESFASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr);
78691f99d7cSPeter Brune     ierr = SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
787fe6f9142SPeter Brune   }
788fa9694d7SPeter Brune   PetscFunctionReturn(0);
789421d9b32SPeter Brune }
790421d9b32SPeter Brune 
791421d9b32SPeter Brune #undef __FUNCT__
792421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS"
793421d9b32SPeter Brune 
794421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes)
795421d9b32SPeter Brune {
796fa9694d7SPeter Brune   PetscErrorCode ierr;
797fe6f9142SPeter Brune   PetscInt       i, maxits;
798ddb5aff1SPeter Brune   Vec            X, F;
799fe6f9142SPeter Brune   PetscReal      fnorm;
800b17ce1afSJed Brown   SNES_FAS       *fas = (SNES_FAS *)snes->data,*ffas;
801b17ce1afSJed Brown   DM             dm;
802e70c42e5SPeter Brune   PetscBool      isFine;
803b17ce1afSJed Brown 
804421d9b32SPeter Brune   PetscFunctionBegin;
805fe6f9142SPeter Brune   maxits = snes->max_its;            /* maximum number of iterations */
806fe6f9142SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
807fa9694d7SPeter Brune   X = snes->vec_sol;
808f5a6d4f9SBarry Smith   F = snes->vec_func;
809293a7e31SPeter Brune 
81018a66777SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
811293a7e31SPeter Brune   /*norm setup */
812fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
813fe6f9142SPeter Brune   snes->iter = 0;
814fe6f9142SPeter Brune   snes->norm = 0.;
815fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
816e4ed7901SPeter Brune   if (!snes->vec_func_init_set) {
8170dd27c6cSPeter Brune     if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);}
818fe6f9142SPeter Brune     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
8190dd27c6cSPeter Brune     if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);}
820fe6f9142SPeter Brune     if (snes->domainerror) {
821fe6f9142SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
822fe6f9142SPeter Brune       PetscFunctionReturn(0);
823fe6f9142SPeter Brune     }
824e4ed7901SPeter Brune   } else {
825e4ed7901SPeter Brune     snes->vec_func_init_set = PETSC_FALSE;
826e4ed7901SPeter Brune   }
827e4ed7901SPeter Brune 
828e4ed7901SPeter Brune   if (!snes->norm_init_set) {
829fe6f9142SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
830fe6f9142SPeter Brune     if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
831fe6f9142SPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
832e4ed7901SPeter Brune   } else {
833e4ed7901SPeter Brune     fnorm = snes->norm_init;
834e4ed7901SPeter Brune     snes->norm_init_set = PETSC_FALSE;
835e4ed7901SPeter Brune   }
836e4ed7901SPeter Brune 
837fe6f9142SPeter Brune   snes->norm = fnorm;
838fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
839fe6f9142SPeter Brune   SNESLogConvHistory(snes,fnorm,0);
840fe6f9142SPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
841fe6f9142SPeter Brune 
842fe6f9142SPeter Brune   /* set parameter for default relative tolerance convergence test */
843fe6f9142SPeter Brune   snes->ttol = fnorm*snes->rtol;
844fe6f9142SPeter Brune   /* test convergence */
845fe6f9142SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
846fe6f9142SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
847e4ed7901SPeter Brune 
848b17ce1afSJed Brown 
849b9c2fdf1SPeter Brune   if (isFine) {
850b9c2fdf1SPeter Brune     /* propagate scale-dependent data up the hierarchy */
851b17ce1afSJed Brown     ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
852b17ce1afSJed Brown     for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) {
853b17ce1afSJed Brown       DM dmcoarse;
854b17ce1afSJed Brown       ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr);
855b17ce1afSJed Brown       ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr);
856b17ce1afSJed Brown       dm = dmcoarse;
857b17ce1afSJed Brown     }
858b9c2fdf1SPeter Brune   }
859b17ce1afSJed Brown 
860fe6f9142SPeter Brune   for (i = 0; i < maxits; i++) {
861fe6f9142SPeter Brune     /* Call general purpose update function */
862646217ecSPeter Brune 
863fe6f9142SPeter Brune     if (snes->ops->update) {
864fe6f9142SPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
865fe6f9142SPeter Brune     }
86607144faaSPeter Brune     if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
86791f99d7cSPeter Brune       ierr = SNESFASCycle_Multiplicative(snes, X);CHKERRQ(ierr);
86807144faaSPeter Brune     } else {
86991f99d7cSPeter Brune       ierr = SNESFASCycle_Additive(snes, X);CHKERRQ(ierr);
87007144faaSPeter Brune     }
871742fe5e2SPeter Brune 
872742fe5e2SPeter Brune     /* check for FAS cycle divergence */
873742fe5e2SPeter Brune     if (snes->reason != SNES_CONVERGED_ITERATING) {
874742fe5e2SPeter Brune       PetscFunctionReturn(0);
875742fe5e2SPeter Brune     }
876b9c2fdf1SPeter Brune 
877c90fad12SPeter Brune     /* Monitor convergence */
878c90fad12SPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
879c90fad12SPeter Brune     snes->iter = i+1;
880c90fad12SPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
881c90fad12SPeter Brune     SNESLogConvHistory(snes,snes->norm,0);
882c90fad12SPeter Brune     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
883c90fad12SPeter Brune     /* Test for convergence */
88466585501SPeter Brune     if (isFine) {
885b9c2fdf1SPeter Brune       ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,snes->norm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
886c90fad12SPeter Brune       if (snes->reason) break;
887fe6f9142SPeter Brune     }
88866585501SPeter Brune   }
889fe6f9142SPeter Brune   if (i == maxits) {
890fe6f9142SPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
891fe6f9142SPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
892fe6f9142SPeter Brune   }
893421d9b32SPeter Brune   PetscFunctionReturn(0);
894421d9b32SPeter Brune }
895