xref: /petsc/src/snes/impls/fas/fas.c (revision 96a0c9949c19278363492ba7462c3ac2a88a3ed1)
1421d9b32SPeter Brune /* Defines the basic SNES object */
2a055c8caSBarry Smith #include <../src/snes/impls/fas/fasimpls.h>    /*I  "petscsnes.h"  I*/
3421d9b32SPeter Brune 
434d65b3cSPeter 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);
84416b707SBarry Smith extern PetscErrorCode SNESSetFromOptions_FAS(PetscOptionItems *PetscOptionsObject,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
26928e959bSPeter Brune .   -snes_fas_type<additive,multiplicative,full,kaskade>  -  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
31d39b9438SMatthew G. Knepley .   -snes_fas_full_downsweep<PETSC_FALSE> - call the downsmooth on the initial downsweep of full FAS
32d3bc2379SPeter Brune .   -fas_levels_snes_ -  SNES options for all smoothers
337d84e935SPeter Brune .   -fas_levels_cycle_snes_ -  SNES options for all cycles
34d3bc2379SPeter Brune .   -fas_levels_i_snes_ -  SNES options for the smoothers on level i
357d84e935SPeter Brune .   -fas_levels_i_cycle_snes_ - SNES options for the cycle on level i
36d3bc2379SPeter Brune -   -fas_coarse_snes_ -  SNES options for the coarsest smoother
37d3bc2379SPeter Brune 
38d3bc2379SPeter Brune Notes:
39d3bc2379SPeter Brune    The organization of the FAS solver is slightly different from the organization of PCMG
40d3bc2379SPeter Brune    As each level has smoother SNES instances(down and potentially up) and a cycle SNES instance.
41d3bc2379SPeter Brune    The cycle SNES instance may be used for monitoring convergence on a particular level.
421fbfccc6SJed Brown 
437d84e935SPeter Brune Level: beginner
441fbfccc6SJed Brown 
454f02bc6aSBarry Smith    References:
46*96a0c994SBarry Smith . 1. -  Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
474f02bc6aSBarry Smith    SIAM Review, 57(4), 2015
484f02bc6aSBarry Smith 
49d3bc2379SPeter Brune .seealso: PCMG, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
501fbfccc6SJed Brown M*/
51421d9b32SPeter Brune 
52421d9b32SPeter Brune #undef __FUNCT__
53421d9b32SPeter Brune #define __FUNCT__ "SNESCreate_FAS"
548cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_FAS(SNES snes)
55421d9b32SPeter Brune {
56421d9b32SPeter Brune   SNES_FAS       *fas;
57421d9b32SPeter Brune   PetscErrorCode ierr;
58421d9b32SPeter Brune 
59421d9b32SPeter Brune   PetscFunctionBegin;
60421d9b32SPeter Brune   snes->ops->destroy        = SNESDestroy_FAS;
61421d9b32SPeter Brune   snes->ops->setup          = SNESSetUp_FAS;
62421d9b32SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_FAS;
63421d9b32SPeter Brune   snes->ops->view           = SNESView_FAS;
64421d9b32SPeter Brune   snes->ops->solve          = SNESSolve_FAS;
65421d9b32SPeter Brune   snes->ops->reset          = SNESReset_FAS;
66421d9b32SPeter Brune 
67ed020824SBarry Smith   snes->usesksp = PETSC_FALSE;
68ed020824SBarry Smith   snes->usespc  = PETSC_FALSE;
69ed020824SBarry Smith 
7088976e71SPeter Brune   if (!snes->tolerancesset) {
710e444f03SPeter Brune     snes->max_funcs = 30000;
720e444f03SPeter Brune     snes->max_its   = 10000;
7388976e71SPeter Brune   }
740e444f03SPeter Brune 
75b00a9115SJed Brown   ierr = PetscNewLog(snes,&fas);CHKERRQ(ierr);
761aa26658SKarl Rupp 
77421d9b32SPeter Brune   snes->data                  = (void*) fas;
78421d9b32SPeter Brune   fas->level                  = 0;
79293a7e31SPeter Brune   fas->levels                 = 1;
80ee78dd50SPeter Brune   fas->n_cycles               = 1;
81ee78dd50SPeter Brune   fas->max_up_it              = 1;
82ee78dd50SPeter Brune   fas->max_down_it            = 1;
830298fd71SBarry Smith   fas->smoothu                = NULL;
840298fd71SBarry Smith   fas->smoothd                = NULL;
850298fd71SBarry Smith   fas->next                   = NULL;
860298fd71SBarry Smith   fas->previous               = NULL;
87ab8d36c9SPeter Brune   fas->fine                   = snes;
880298fd71SBarry Smith   fas->interpolate            = NULL;
890298fd71SBarry Smith   fas->restrct                = NULL;
900298fd71SBarry Smith   fas->inject                 = NULL;
910298fd71SBarry Smith   fas->monitor                = NULL;
92cc05f883SPeter Brune   fas->usedmfornumberoflevels = PETSC_FALSE;
93ddebd997SPeter Brune   fas->fastype                = SNES_FAS_MULTIPLICATIVE;
94928e959bSPeter Brune   fas->full_downsweep         = PETSC_FALSE;
950dd27c6cSPeter Brune 
960298fd71SBarry Smith   fas->eventsmoothsetup    = 0;
970298fd71SBarry Smith   fas->eventsmoothsolve    = 0;
980298fd71SBarry Smith   fas->eventresidual       = 0;
990298fd71SBarry Smith   fas->eventinterprestrict = 0;
100efe1f98aSPeter Brune   PetscFunctionReturn(0);
101efe1f98aSPeter Brune }
102efe1f98aSPeter Brune 
103421d9b32SPeter Brune #undef __FUNCT__
104421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS"
105421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes)
106421d9b32SPeter Brune {
10777df8cc4SPeter Brune   PetscErrorCode ierr  = 0;
108421d9b32SPeter Brune   SNES_FAS       * fas = (SNES_FAS*)snes->data;
109421d9b32SPeter Brune 
110421d9b32SPeter Brune   PetscFunctionBegin;
111ab8d36c9SPeter Brune   ierr = SNESDestroy(&fas->smoothu);CHKERRQ(ierr);
112ab8d36c9SPeter Brune   ierr = SNESDestroy(&fas->smoothd);CHKERRQ(ierr);
1133dccd265SPeter Brune   ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
114bccf9bb3SJed Brown   ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
115bccf9bb3SJed Brown   ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
116bccf9bb3SJed Brown   ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
1177cd33a7bSPeter Brune   if (fas->galerkin) {
1187cd33a7bSPeter Brune     ierr = VecDestroy(&fas->Xg);CHKERRQ(ierr);
1197cd33a7bSPeter Brune     ierr = VecDestroy(&fas->Fg);CHKERRQ(ierr);
1207cd33a7bSPeter Brune   }
121d477d801SPeter Brune   if (fas->next) {ierr = SNESReset(fas->next);CHKERRQ(ierr);}
122421d9b32SPeter Brune   PetscFunctionReturn(0);
123421d9b32SPeter Brune }
124421d9b32SPeter Brune 
125421d9b32SPeter Brune #undef __FUNCT__
126421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS"
127421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes)
128421d9b32SPeter Brune {
129421d9b32SPeter Brune   SNES_FAS       * fas = (SNES_FAS*)snes->data;
130742fe5e2SPeter Brune   PetscErrorCode ierr  = 0;
131421d9b32SPeter Brune 
132421d9b32SPeter Brune   PetscFunctionBegin;
133421d9b32SPeter Brune   /* recursively resets and then destroys */
13479d9a41aSPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
1351aa26658SKarl Rupp   if (fas->next) {
1361aa26658SKarl Rupp     ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
1371aa26658SKarl Rupp   }
138421d9b32SPeter Brune   ierr = PetscFree(fas);CHKERRQ(ierr);
139421d9b32SPeter Brune   PetscFunctionReturn(0);
140421d9b32SPeter Brune }
141421d9b32SPeter Brune 
142421d9b32SPeter Brune #undef __FUNCT__
143421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS"
144421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes)
145421d9b32SPeter Brune {
14648bfdf8aSPeter Brune   SNES_FAS       *fas = (SNES_FAS*) snes->data;
147421d9b32SPeter Brune   PetscErrorCode ierr;
148d1adcc6fSPeter Brune   PetscInt       dm_levels;
1493dccd265SPeter Brune   Vec            vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */
150ab8d36c9SPeter Brune   SNES           next;
151ab8d36c9SPeter Brune   PetscBool      isFine;
152f89ba88eSPeter Brune   SNESLineSearch linesearch;
153f89ba88eSPeter Brune   SNESLineSearch slinesearch;
154f89ba88eSPeter Brune   void           *lsprectx,*lspostctx;
1556b2b7091SBarry Smith   PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*);
1566b2b7091SBarry Smith   PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*);
157eff52c0eSPeter Brune 
1586b2b7091SBarry Smith   PetscFunctionBegin;
159ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
160ab8d36c9SPeter Brune   if (fas->usedmfornumberoflevels && isFine) {
161d1adcc6fSPeter Brune     ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr);
162d1adcc6fSPeter Brune     dm_levels++;
163cc05f883SPeter Brune     if (dm_levels > fas->levels) {
1642e8ce248SJed Brown       /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESFASSetLevels_FAS*/
1653dccd265SPeter Brune       vec_sol              = snes->vec_sol;
1663dccd265SPeter Brune       vec_func             = snes->vec_func;
1673dccd265SPeter Brune       vec_sol_update       = snes->vec_sol_update;
1683dccd265SPeter Brune       vec_rhs              = snes->vec_rhs;
1690298fd71SBarry Smith       snes->vec_sol        = NULL;
1700298fd71SBarry Smith       snes->vec_func       = NULL;
1710298fd71SBarry Smith       snes->vec_sol_update = NULL;
1720298fd71SBarry Smith       snes->vec_rhs        = NULL;
1733dccd265SPeter Brune 
1743dccd265SPeter Brune       /* reset the number of levels */
1750298fd71SBarry Smith       ierr = SNESFASSetLevels(snes,dm_levels,NULL);CHKERRQ(ierr);
176cc05f883SPeter Brune       ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
1773dccd265SPeter Brune 
1783dccd265SPeter Brune       snes->vec_sol        = vec_sol;
1793dccd265SPeter Brune       snes->vec_func       = vec_func;
1803dccd265SPeter Brune       snes->vec_rhs        = vec_rhs;
1813dccd265SPeter Brune       snes->vec_sol_update = vec_sol_update;
182d1adcc6fSPeter Brune     }
183d1adcc6fSPeter Brune   }
184ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
185ab8d36c9SPeter Brune   if (!isFine) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */
1863dccd265SPeter Brune 
187fa0ddf94SBarry Smith   ierr = SNESSetWorkVecs(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
188cc05f883SPeter Brune 
189ab8d36c9SPeter Brune   /* set up the smoothers if they haven't already been set up */
190ab8d36c9SPeter Brune   if (!fas->smoothd) {
191ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
192ab8d36c9SPeter Brune   }
193ab8d36c9SPeter Brune 
19479d9a41aSPeter Brune   if (snes->dm) {
195ab8d36c9SPeter Brune     /* set the smoother DMs properly */
196ab8d36c9SPeter Brune     if (fas->smoothu) ierr = SNESSetDM(fas->smoothu, snes->dm);CHKERRQ(ierr);
197ab8d36c9SPeter Brune     ierr = SNESSetDM(fas->smoothd, snes->dm);CHKERRQ(ierr);
19879d9a41aSPeter Brune     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
199ab8d36c9SPeter Brune     if (next) {
20079d9a41aSPeter Brune       /* for now -- assume the DM and the evaluation functions have been set externally */
201ab8d36c9SPeter Brune       if (!next->dm) {
202ce94432eSBarry Smith         ierr = DMCoarsen(snes->dm, PetscObjectComm((PetscObject)next), &next->dm);CHKERRQ(ierr);
203ab8d36c9SPeter Brune         ierr = SNESSetDM(next, next->dm);CHKERRQ(ierr);
20479d9a41aSPeter Brune       }
20579d9a41aSPeter Brune       /* set the interpolation and restriction from the DM */
20679d9a41aSPeter Brune       if (!fas->interpolate) {
207ab8d36c9SPeter Brune         ierr = DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
208bccf9bb3SJed Brown         if (!fas->restrct) {
209bccf9bb3SJed Brown           ierr         = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr);
21079d9a41aSPeter Brune           fas->restrct = fas->interpolate;
21179d9a41aSPeter Brune         }
212bccf9bb3SJed Brown       }
21379d9a41aSPeter Brune       /* set the injection from the DM */
21479d9a41aSPeter Brune       if (!fas->inject) {
2156dbf9973SLawrence Mitchell         ierr = DMCreateInjection(next->dm, snes->dm, &fas->inject);CHKERRQ(ierr);
21679d9a41aSPeter Brune       }
21779d9a41aSPeter Brune     }
21879d9a41aSPeter Brune   }
21979d9a41aSPeter Brune   /*pass the smoother, function, and jacobian up to the next level if it's not user set already */
22079d9a41aSPeter Brune   if (fas->galerkin) {
2211aa26658SKarl Rupp     if (next) {
2220298fd71SBarry Smith       ierr = SNESSetFunction(next, NULL, SNESFASGalerkinDefaultFunction, next);CHKERRQ(ierr);
2231aa26658SKarl Rupp     }
2241aa26658SKarl Rupp     if (fas->smoothd && fas->level != fas->levels - 1) {
2250298fd71SBarry Smith       ierr = SNESSetFunction(fas->smoothd, NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr);
2261aa26658SKarl Rupp     }
2271aa26658SKarl Rupp     if (fas->smoothu && fas->level != fas->levels - 1) {
2280298fd71SBarry Smith       ierr = SNESSetFunction(fas->smoothu, NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr);
2291aa26658SKarl Rupp     }
23079d9a41aSPeter Brune   }
23179d9a41aSPeter Brune 
232534ebe21SPeter Brune   /* sets the down (pre) smoother's default norm and sets it from options */
233534ebe21SPeter Brune   if (fas->smoothd) {
234bc3f2f05SPeter Brune     if (fas->level == 0 && fas->levels != 1) {
235365a6726SPeter Brune       ierr = SNESSetNormSchedule(fas->smoothd, SNES_NORM_NONE);CHKERRQ(ierr);
236534ebe21SPeter Brune     } else {
237365a6726SPeter Brune       ierr = SNESSetNormSchedule(fas->smoothd, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr);
238534ebe21SPeter Brune     }
2397fce8c19SPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)fas->smoothd);CHKERRQ(ierr);
240534ebe21SPeter Brune     ierr = SNESSetFromOptions(fas->smoothd);CHKERRQ(ierr);
2417601faf0SJed Brown     ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr);
2427601faf0SJed Brown     ierr = SNESGetLineSearch(fas->smoothd,&slinesearch);CHKERRQ(ierr);
2436b2b7091SBarry Smith     ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr);
2446b2b7091SBarry Smith     ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr);
2456b2b7091SBarry Smith     ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr);
2466b2b7091SBarry Smith     ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr);
247f89ba88eSPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr);
2480dd27c6cSPeter Brune 
2490dd27c6cSPeter Brune     fas->smoothd->vec_sol        = snes->vec_sol;
2500dd27c6cSPeter Brune     ierr                         = PetscObjectReference((PetscObject)snes->vec_sol);CHKERRQ(ierr);
2510dd27c6cSPeter Brune     fas->smoothd->vec_sol_update = snes->vec_sol_update;
2520dd27c6cSPeter Brune     ierr                         = PetscObjectReference((PetscObject)snes->vec_sol_update);CHKERRQ(ierr);
2530dd27c6cSPeter Brune     fas->smoothd->vec_func       = snes->vec_func;
2540dd27c6cSPeter Brune     ierr                         = PetscObjectReference((PetscObject)snes->vec_func);CHKERRQ(ierr);
2550dd27c6cSPeter Brune 
2560dd27c6cSPeter Brune     if (fas->eventsmoothsetup) {ierr = PetscLogEventBegin(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
2570dd27c6cSPeter Brune     ierr = SNESSetUp(fas->smoothd);CHKERRQ(ierr);
2580dd27c6cSPeter Brune     if (fas->eventsmoothsetup) {ierr = PetscLogEventEnd(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
259534ebe21SPeter Brune   }
260534ebe21SPeter Brune 
261534ebe21SPeter Brune   /* sets the up (post) smoother's default norm and sets it from options */
262534ebe21SPeter Brune   if (fas->smoothu) {
263534ebe21SPeter Brune     if (fas->level != fas->levels - 1) {
264365a6726SPeter Brune       ierr = SNESSetNormSchedule(fas->smoothu, SNES_NORM_NONE);CHKERRQ(ierr);
265534ebe21SPeter Brune     } else {
266365a6726SPeter Brune       ierr = SNESSetNormSchedule(fas->smoothu, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr);
267534ebe21SPeter Brune     }
2687fce8c19SPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)fas->smoothu);CHKERRQ(ierr);
269534ebe21SPeter Brune     ierr = SNESSetFromOptions(fas->smoothu);CHKERRQ(ierr);
2707601faf0SJed Brown     ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr);
2717601faf0SJed Brown     ierr = SNESGetLineSearch(fas->smoothu,&slinesearch);CHKERRQ(ierr);
2726b2b7091SBarry Smith     ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr);
2736b2b7091SBarry Smith     ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr);
2746b2b7091SBarry Smith     ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr);
2756b2b7091SBarry Smith     ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr);
276f89ba88eSPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr);
2770dd27c6cSPeter Brune 
2780dd27c6cSPeter Brune     fas->smoothu->vec_sol        = snes->vec_sol;
2790dd27c6cSPeter Brune     ierr                         = PetscObjectReference((PetscObject)snes->vec_sol);CHKERRQ(ierr);
2800dd27c6cSPeter Brune     fas->smoothu->vec_sol_update = snes->vec_sol_update;
2810dd27c6cSPeter Brune     ierr                         = PetscObjectReference((PetscObject)snes->vec_sol_update);CHKERRQ(ierr);
2820dd27c6cSPeter Brune     fas->smoothu->vec_func       = snes->vec_func;
2830dd27c6cSPeter Brune     ierr                         = PetscObjectReference((PetscObject)snes->vec_func);CHKERRQ(ierr);
2840dd27c6cSPeter Brune 
2850dd27c6cSPeter Brune     if (fas->eventsmoothsetup) {ierr = PetscLogEventBegin(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
2860dd27c6cSPeter Brune     ierr = SNESSetUp(fas->smoothu);CHKERRQ(ierr);
2870dd27c6cSPeter Brune     if (fas->eventsmoothsetup) {ierr = PetscLogEventEnd(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
2880dd27c6cSPeter Brune 
289534ebe21SPeter Brune   }
290d06165b7SPeter Brune 
291ab8d36c9SPeter Brune   if (next) {
29279d9a41aSPeter Brune     /* gotta set up the solution vector for this to work */
293ab8d36c9SPeter Brune     if (!next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_sol);CHKERRQ(ierr);}
294ab8d36c9SPeter Brune     if (!next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_rhs);CHKERRQ(ierr);}
2957fce8c19SPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)next);CHKERRQ(ierr);
2967601faf0SJed Brown     ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr);
2977601faf0SJed Brown     ierr = SNESGetLineSearch(fas->next,&slinesearch);CHKERRQ(ierr);
2986b2b7091SBarry Smith     ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr);
2996b2b7091SBarry Smith     ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr);
3006b2b7091SBarry Smith     ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr);
3016b2b7091SBarry Smith     ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr);
302f89ba88eSPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr);
303ab8d36c9SPeter Brune     ierr = SNESSetUp(next);CHKERRQ(ierr);
30479d9a41aSPeter Brune   }
3056273346dSPeter Brune   /* setup FAS work vectors */
3066273346dSPeter Brune   if (fas->galerkin) {
3076273346dSPeter Brune     ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr);
3086273346dSPeter Brune     ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr);
3096273346dSPeter Brune   }
310421d9b32SPeter Brune   PetscFunctionReturn(0);
311421d9b32SPeter Brune }
312421d9b32SPeter Brune 
313421d9b32SPeter Brune #undef __FUNCT__
314421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS"
3154416b707SBarry Smith PetscErrorCode SNESSetFromOptions_FAS(PetscOptionItems *PetscOptionsObject,SNES snes)
316421d9b32SPeter Brune {
317ee78dd50SPeter Brune   SNES_FAS       *fas   = (SNES_FAS*) snes->data;
318ee78dd50SPeter Brune   PetscInt       levels = 1;
31987f44e3fSPeter Brune   PetscBool      flg    = PETSC_FALSE, upflg = PETSC_FALSE, downflg = PETSC_FALSE, monflg = PETSC_FALSE, galerkinflg = PETSC_FALSE,continuationflg = PETSC_FALSE;
320421d9b32SPeter Brune   PetscErrorCode ierr;
321ee78dd50SPeter Brune   char           monfilename[PETSC_MAX_PATH_LEN];
32207144faaSPeter Brune   SNESFASType    fastype;
323fde0ff24SPeter Brune   const char     *optionsprefix;
324f1c6b773SPeter Brune   SNESLineSearch linesearch;
32566585501SPeter Brune   PetscInt       m, n_up, n_down;
326ab8d36c9SPeter Brune   SNES           next;
327ab8d36c9SPeter Brune   PetscBool      isFine;
328421d9b32SPeter Brune 
329421d9b32SPeter Brune   PetscFunctionBegin;
330ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
331e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SNESFAS Options-----------------------------------");CHKERRQ(ierr);
332ee78dd50SPeter Brune 
333ab8d36c9SPeter Brune   /* number of levels -- only process most options on the finest level */
334ab8d36c9SPeter Brune   if (isFine) {
335ee78dd50SPeter Brune     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
336c732cbdbSBarry Smith     if (!flg && snes->dm) {
337c732cbdbSBarry Smith       ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr);
338c732cbdbSBarry Smith       levels++;
339d1adcc6fSPeter Brune       fas->usedmfornumberoflevels = PETSC_TRUE;
340c732cbdbSBarry Smith     }
3410298fd71SBarry Smith     ierr    = SNESFASSetLevels(snes, levels, NULL);CHKERRQ(ierr);
34207144faaSPeter Brune     fastype = fas->fastype;
34307144faaSPeter Brune     ierr    = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr);
34407144faaSPeter Brune     if (flg) {
34507144faaSPeter Brune       ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr);
34607144faaSPeter Brune     }
347ee78dd50SPeter Brune 
348fde0ff24SPeter Brune     ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
349ab8d36c9SPeter Brune     ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&m,&flg);CHKERRQ(ierr);
350ab8d36c9SPeter Brune     if (flg) {
351ab8d36c9SPeter Brune       ierr = SNESFASSetCycles(snes, m);CHKERRQ(ierr);
352fde0ff24SPeter Brune     }
35387f44e3fSPeter Brune     ierr = PetscOptionsBool("-snes_fas_continuation","Corrected grid-sequence continuation","SNESFASSetContinuation",fas->continuation,&continuationflg,&flg);CHKERRQ(ierr);
35487f44e3fSPeter Brune     if (flg) {
35587f44e3fSPeter Brune       ierr = SNESFASSetContinuation(snes,continuationflg);CHKERRQ(ierr);
35687f44e3fSPeter Brune     }
357fde0ff24SPeter Brune 
358ab8d36c9SPeter Brune     ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFASSetGalerkin",fas->galerkin,&galerkinflg,&flg);CHKERRQ(ierr);
359ab8d36c9SPeter Brune     if (flg) {
360ab8d36c9SPeter Brune       ierr = SNESFASSetGalerkin(snes, galerkinflg);CHKERRQ(ierr);
361ab8d36c9SPeter Brune     }
362ee78dd50SPeter Brune 
363928e959bSPeter Brune     if (fas->fastype == SNES_FAS_FULL) {
364928e959bSPeter Brune       ierr   = PetscOptionsBool("-snes_fas_full_downsweep","Smooth on the initial upsweep for full FAS cycles","SNESFASFullSetDownSweep",fas->full_downsweep,&fas->full_downsweep,&flg);CHKERRQ(ierr);
365928e959bSPeter Brune       if (flg) {SNESFASFullSetDownSweep(snes,fas->full_downsweep);CHKERRQ(ierr);}
366928e959bSPeter Brune     }
367928e959bSPeter Brune 
36866585501SPeter Brune     ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",fas->max_up_it,&n_up,&upflg);CHKERRQ(ierr);
369162d76ddSPeter Brune 
37066585501SPeter Brune     ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",fas->max_down_it,&n_down,&downflg);CHKERRQ(ierr);
371162d76ddSPeter Brune 
372c8c899caSPeter Brune     ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESFASSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr);
373c8c899caSPeter Brune     if (monflg) ierr = SNESFASSetMonitor(snes, PETSC_TRUE);CHKERRQ(ierr);
3740dd27c6cSPeter Brune 
3750dd27c6cSPeter Brune     flg    = PETSC_FALSE;
3760dd27c6cSPeter Brune     monflg = PETSC_TRUE;
3770dd27c6cSPeter Brune     ierr   = PetscOptionsBool("-snes_fas_log","Log times for each FAS level","SNESFASSetLog",monflg,&monflg,&flg);CHKERRQ(ierr);
3780dd27c6cSPeter Brune     if (flg) {ierr = SNESFASSetLog(snes,monflg);CHKERRQ(ierr);}
379ab8d36c9SPeter Brune   }
380ee78dd50SPeter Brune 
381421d9b32SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
3828cc86e31SPeter Brune   /* setup from the determined types if there is no pointwise procedure or smoother defined */
383162d76ddSPeter Brune   if (upflg) {
38466585501SPeter Brune     ierr = SNESFASSetNumberSmoothUp(snes,n_up);CHKERRQ(ierr);
385162d76ddSPeter Brune   }
386162d76ddSPeter Brune   if (downflg) {
38766585501SPeter Brune     ierr = SNESFASSetNumberSmoothDown(snes,n_down);CHKERRQ(ierr);
388162d76ddSPeter Brune   }
389eff52c0eSPeter Brune 
3909e764e56SPeter Brune   /* set up the default line search for coarse grid corrections */
3919e764e56SPeter Brune   if (fas->fastype == SNES_FAS_ADDITIVE) {
3929e764e56SPeter Brune     if (!snes->linesearch) {
3937601faf0SJed Brown       ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
3941a4f838cSPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr);
3959e764e56SPeter Brune     }
3969e764e56SPeter Brune   }
3979e764e56SPeter Brune 
398ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
399ee78dd50SPeter Brune   /* recursive option setting for the smoothers */
400ab8d36c9SPeter Brune   if (next) {ierr = SNESSetFromOptions(next);CHKERRQ(ierr);}
401421d9b32SPeter Brune   PetscFunctionReturn(0);
402421d9b32SPeter Brune }
403421d9b32SPeter Brune 
4049804daf3SBarry Smith #include <petscdraw.h>
405421d9b32SPeter Brune #undef __FUNCT__
406421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS"
407421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
408421d9b32SPeter Brune {
409421d9b32SPeter Brune   SNES_FAS       *fas = (SNES_FAS*) snes->data;
410656ede7eSPeter Brune   PetscBool      isFine,iascii,isdraw;
411ab8d36c9SPeter Brune   PetscInt       i;
412421d9b32SPeter Brune   PetscErrorCode ierr;
413ab8d36c9SPeter Brune   SNES           smoothu, smoothd, levelsnes;
414421d9b32SPeter Brune 
415421d9b32SPeter Brune   PetscFunctionBegin;
416ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
417ab8d36c9SPeter Brune   if (isFine) {
418251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
419656ede7eSPeter Brune     ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
420421d9b32SPeter Brune     if (iascii) {
421ab8d36c9SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "FAS: type is %s, levels=%D, cycles=%D\n",  SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr);
422ab8d36c9SPeter Brune       if (fas->galerkin) {
423ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
424421d9b32SPeter Brune       } else {
425ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
426421d9b32SPeter Brune       }
427ab8d36c9SPeter Brune       for (i=0; i<fas->levels; i++) {
428ab8d36c9SPeter Brune         ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr);
429ab8d36c9SPeter Brune         ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr);
430ab8d36c9SPeter Brune         ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr);
431ab8d36c9SPeter Brune         if (!i) {
432ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr);
433421d9b32SPeter Brune         } else {
434ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
435421d9b32SPeter Brune         }
436ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
437166b3ea4SJed Brown         if (smoothd) {
438ab8d36c9SPeter Brune           ierr = SNESView(smoothd,viewer);CHKERRQ(ierr);
439166b3ea4SJed Brown         } else {
440166b3ea4SJed Brown           ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr);
441166b3ea4SJed Brown         }
442ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
443ab8d36c9SPeter Brune         if (i && (smoothd == smoothu)) {
444ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
445ab8d36c9SPeter Brune         } else if (i) {
446ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
447ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
448166b3ea4SJed Brown           if (smoothu) {
449ab8d36c9SPeter Brune             ierr = SNESView(smoothu,viewer);CHKERRQ(ierr);
450166b3ea4SJed Brown           } else {
451166b3ea4SJed Brown             ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr);
452166b3ea4SJed Brown           }
453ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
454ab8d36c9SPeter Brune         }
455ab8d36c9SPeter Brune       }
456656ede7eSPeter Brune     } else if (isdraw) {
457656ede7eSPeter Brune       PetscDraw draw;
458b4375e8dSPeter Brune       PetscReal x,w,y,bottom,th,wth;
459656ede7eSPeter Brune       SNES_FAS  *curfas = fas;
460656ede7eSPeter Brune       ierr   = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
461656ede7eSPeter Brune       ierr   = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
462656ede7eSPeter Brune       ierr   = PetscDrawStringGetSize(draw,&wth,&th);CHKERRQ(ierr);
463656ede7eSPeter Brune       bottom = y - th;
464656ede7eSPeter Brune       while (curfas) {
465b4375e8dSPeter Brune         if (!curfas->smoothu) {
466656ede7eSPeter Brune           ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
467656ede7eSPeter Brune           if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr);
468656ede7eSPeter Brune           ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
469b4375e8dSPeter Brune         } else {
470b4375e8dSPeter Brune           w    = 0.5*PetscMin(1.0-x,x);
471b4375e8dSPeter Brune           ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
472b4375e8dSPeter Brune           if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr);
473b4375e8dSPeter Brune           ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
474b4375e8dSPeter Brune           ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
475b4375e8dSPeter Brune           if (curfas->smoothu) ierr = SNESView(curfas->smoothu,viewer);CHKERRQ(ierr);
476b4375e8dSPeter Brune           ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
477b4375e8dSPeter Brune         }
478656ede7eSPeter Brune         /* this is totally bogus but we have no way of knowing how low the previous one was draw to */
479656ede7eSPeter Brune         bottom -= 5*th;
4801aa26658SKarl Rupp         if (curfas->next) curfas = (SNES_FAS*)curfas->next->data;
4810298fd71SBarry Smith         else curfas = NULL;
482656ede7eSPeter Brune       }
483421d9b32SPeter Brune     }
484ab8d36c9SPeter Brune   }
485421d9b32SPeter Brune   PetscFunctionReturn(0);
486421d9b32SPeter Brune }
487421d9b32SPeter Brune 
488421d9b32SPeter Brune #undef __FUNCT__
48991f99d7cSPeter Brune #define __FUNCT__ "SNESFASDownSmooth_Private"
49039bd7f45SPeter Brune /*
49139bd7f45SPeter Brune Defines the action of the downsmoother
49239bd7f45SPeter Brune  */
49391f99d7cSPeter Brune PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
494b9c2fdf1SPeter Brune {
49539bd7f45SPeter Brune   PetscErrorCode      ierr = 0;
496742fe5e2SPeter Brune   SNESConvergedReason reason;
497ab8d36c9SPeter Brune   Vec                 FPC;
498ab8d36c9SPeter Brune   SNES                smoothd;
4990dd27c6cSPeter Brune   SNES_FAS            *fas = (SNES_FAS*) snes->data;
5006e111a19SKarl Rupp 
501421d9b32SPeter Brune   PetscFunctionBegin;
502ab8d36c9SPeter Brune   ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr);
503e4ed7901SPeter Brune   ierr = SNESSetInitialFunction(smoothd, F);CHKERRQ(ierr);
5040dd27c6cSPeter Brune   if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
505ab8d36c9SPeter Brune   ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr);
5060dd27c6cSPeter Brune   if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
507742fe5e2SPeter Brune   /* check convergence reason for the smoother */
508ab8d36c9SPeter Brune   ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr);
5091ff805c3SPeter Brune   if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) {
510742fe5e2SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
511742fe5e2SPeter Brune     PetscFunctionReturn(0);
512742fe5e2SPeter Brune   }
5130298fd71SBarry Smith   ierr = SNESGetFunction(smoothd, &FPC, NULL, NULL);CHKERRQ(ierr);
5144b32a720SPeter Brune   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
51571dbe336SPeter Brune   if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);}
51639bd7f45SPeter Brune   PetscFunctionReturn(0);
51739bd7f45SPeter Brune }
51839bd7f45SPeter Brune 
51939bd7f45SPeter Brune 
52039bd7f45SPeter Brune #undef __FUNCT__
52191f99d7cSPeter Brune #define __FUNCT__ "SNESFASUpSmooth_Private"
52239bd7f45SPeter Brune /*
52307144faaSPeter Brune Defines the action of the upsmoother
52439bd7f45SPeter Brune  */
5250adebc6cSBarry Smith PetscErrorCode SNESFASUpSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
5260adebc6cSBarry Smith {
52739bd7f45SPeter Brune   PetscErrorCode      ierr = 0;
52839bd7f45SPeter Brune   SNESConvergedReason reason;
529ab8d36c9SPeter Brune   Vec                 FPC;
530ab8d36c9SPeter Brune   SNES                smoothu;
5310dd27c6cSPeter Brune   SNES_FAS            *fas = (SNES_FAS*) snes->data;
532ab8d36c9SPeter Brune 
5336e111a19SKarl Rupp   PetscFunctionBegin;
534ab8d36c9SPeter Brune   ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr);
5350dd27c6cSPeter Brune   if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
536ab8d36c9SPeter Brune   ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr);
5370dd27c6cSPeter Brune   if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
53839bd7f45SPeter Brune   /* check convergence reason for the smoother */
539ab8d36c9SPeter Brune   ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr);
5401ff805c3SPeter Brune   if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) {
54139bd7f45SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
54239bd7f45SPeter Brune     PetscFunctionReturn(0);
54339bd7f45SPeter Brune   }
5440298fd71SBarry Smith   ierr = SNESGetFunction(smoothu, &FPC, NULL, NULL);CHKERRQ(ierr);
5454b32a720SPeter Brune   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
54671dbe336SPeter Brune   if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);}
54739bd7f45SPeter Brune   PetscFunctionReturn(0);
54839bd7f45SPeter Brune }
54939bd7f45SPeter Brune 
55039bd7f45SPeter Brune #undef __FUNCT__
551938e4a01SJed Brown #define __FUNCT__ "SNESFASCreateCoarseVec"
552938e4a01SJed Brown /*@
553938e4a01SJed Brown    SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level
554938e4a01SJed Brown 
555938e4a01SJed Brown    Collective
556938e4a01SJed Brown 
557938e4a01SJed Brown    Input Arguments:
558938e4a01SJed Brown .  snes - SNESFAS
559938e4a01SJed Brown 
560938e4a01SJed Brown    Output Arguments:
561938e4a01SJed Brown .  Xcoarse - vector on level one coarser than snes
562938e4a01SJed Brown 
563938e4a01SJed Brown    Level: developer
564938e4a01SJed Brown 
565938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict()
566938e4a01SJed Brown @*/
567938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse)
568938e4a01SJed Brown {
569938e4a01SJed Brown   PetscErrorCode ierr;
570938e4a01SJed Brown   SNES_FAS       *fas = (SNES_FAS*)snes->data;
571938e4a01SJed Brown 
572938e4a01SJed Brown   PetscFunctionBegin;
5731aa26658SKarl Rupp   if (fas->rscale) {
5741aa26658SKarl Rupp     ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);
575f5af7f23SKarl Rupp   } else if (fas->inject) {
5762a7a6963SBarry Smith     ierr = MatCreateVecs(fas->inject,Xcoarse,NULL);CHKERRQ(ierr);
577ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr);
578938e4a01SJed Brown   PetscFunctionReturn(0);
579938e4a01SJed Brown }
580938e4a01SJed Brown 
581e9923e8dSJed Brown #undef __FUNCT__
582e9923e8dSJed Brown #define __FUNCT__ "SNESFASRestrict"
583e9923e8dSJed Brown /*@
584e9923e8dSJed Brown    SNESFASRestrict - restrict a Vec to the next coarser level
585e9923e8dSJed Brown 
586e9923e8dSJed Brown    Collective
587e9923e8dSJed Brown 
588e9923e8dSJed Brown    Input Arguments:
589e9923e8dSJed Brown +  fine - SNES from which to restrict
590e9923e8dSJed Brown -  Xfine - vector to restrict
591e9923e8dSJed Brown 
592e9923e8dSJed Brown    Output Arguments:
593e9923e8dSJed Brown .  Xcoarse - result of restriction
594e9923e8dSJed Brown 
595e9923e8dSJed Brown    Level: developer
596e9923e8dSJed Brown 
597e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection()
598e9923e8dSJed Brown @*/
599e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse)
600e9923e8dSJed Brown {
601e9923e8dSJed Brown   PetscErrorCode ierr;
602e9923e8dSJed Brown   SNES_FAS       *fas = (SNES_FAS*)fine->data;
603e9923e8dSJed Brown 
604e9923e8dSJed Brown   PetscFunctionBegin;
605e9923e8dSJed Brown   PetscValidHeaderSpecific(fine,SNES_CLASSID,1);
606e9923e8dSJed Brown   PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2);
607e9923e8dSJed Brown   PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3);
608e9923e8dSJed Brown   if (fas->inject) {
609e9923e8dSJed Brown     ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr);
610e9923e8dSJed Brown   } else {
611e9923e8dSJed Brown     ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr);
612e9923e8dSJed Brown     ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr);
613e9923e8dSJed Brown   }
614e9923e8dSJed Brown   PetscFunctionReturn(0);
615e9923e8dSJed Brown }
616e9923e8dSJed Brown 
617e9923e8dSJed Brown #undef __FUNCT__
6188c40d5fbSBarry Smith #define __FUNCT__ "SNESFASCoarseCorrection"
61939bd7f45SPeter Brune /*
62039bd7f45SPeter Brune 
62139bd7f45SPeter Brune Performs the FAS coarse correction as:
62239bd7f45SPeter Brune 
623b5527d98SMatthew G. Knepley fine problem:   F(x) = b
624b5527d98SMatthew G. Knepley coarse problem: F^c(x^c) = b^c
62539bd7f45SPeter Brune 
626b5527d98SMatthew G. Knepley b^c = F^c(Rx) - R(F(x) - b)
62739bd7f45SPeter Brune 
62839bd7f45SPeter Brune  */
6290adebc6cSBarry Smith PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new)
6300adebc6cSBarry Smith {
63139bd7f45SPeter Brune   PetscErrorCode      ierr;
63239bd7f45SPeter Brune   Vec                 X_c, Xo_c, F_c, B_c;
63339bd7f45SPeter Brune   SNESConvergedReason reason;
634ab8d36c9SPeter Brune   SNES                next;
635ab8d36c9SPeter Brune   Mat                 restrct, interpolate;
6360dd27c6cSPeter Brune   SNES_FAS            *fasc;
6375fd66863SKarl Rupp 
63839bd7f45SPeter Brune   PetscFunctionBegin;
639ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
640ab8d36c9SPeter Brune   if (next) {
6410dd27c6cSPeter Brune     fasc = (SNES_FAS*)next->data;
6420dd27c6cSPeter Brune 
643ab8d36c9SPeter Brune     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
644ab8d36c9SPeter Brune     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
645ab8d36c9SPeter Brune 
646ab8d36c9SPeter Brune     X_c  = next->vec_sol;
647ab8d36c9SPeter Brune     Xo_c = next->work[0];
648ab8d36c9SPeter Brune     F_c  = next->vec_func;
649ab8d36c9SPeter Brune     B_c  = next->vec_rhs;
650efe1f98aSPeter Brune 
6510dd27c6cSPeter Brune     if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
652938e4a01SJed Brown     ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr);
6535609cbf2SMatthew G. Knepley     /* restrict the defect: R(F(x) - b) */
654ab8d36c9SPeter Brune     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
6550dd27c6cSPeter Brune     if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
6560dd27c6cSPeter Brune 
6570dd27c6cSPeter Brune     if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);}
6585609cbf2SMatthew G. Knepley     /* F_c = F^c(Rx) - R(F(x) - b) since the second term was sitting in next->vec_rhs */
659ab8d36c9SPeter Brune     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
6600dd27c6cSPeter Brune     if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);}
6610dd27c6cSPeter Brune 
6620dd27c6cSPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */
663e4ed7901SPeter Brune     ierr = VecCopy(B_c, X_c);CHKERRQ(ierr);
664b9c2fdf1SPeter Brune     ierr = VecCopy(F_c, B_c);CHKERRQ(ierr);
665e4ed7901SPeter Brune     ierr = VecCopy(X_c, F_c);CHKERRQ(ierr);
666ee78dd50SPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
667ee78dd50SPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
668c90fad12SPeter Brune 
669c90fad12SPeter Brune     /* recurse to the next level */
670e4ed7901SPeter Brune     ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr);
671ab8d36c9SPeter Brune     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
672ab8d36c9SPeter Brune     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
673742fe5e2SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
674742fe5e2SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
675742fe5e2SPeter Brune       PetscFunctionReturn(0);
676742fe5e2SPeter Brune     }
677fa9694d7SPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
678fa9694d7SPeter Brune     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
6790dd27c6cSPeter Brune 
6800dd27c6cSPeter Brune     if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
681ab8d36c9SPeter Brune     ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr);
6820dd27c6cSPeter Brune     if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
683293a7e31SPeter Brune   }
68439bd7f45SPeter Brune   PetscFunctionReturn(0);
68539bd7f45SPeter Brune }
68639bd7f45SPeter Brune 
68739bd7f45SPeter Brune #undef __FUNCT__
6882cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Additive"
68939bd7f45SPeter Brune /*
69039bd7f45SPeter Brune 
69139bd7f45SPeter Brune The additive cycle looks like:
69239bd7f45SPeter Brune 
69307144faaSPeter Brune xhat = x
69407144faaSPeter Brune xhat = dS(x, b)
69507144faaSPeter Brune x = coarsecorrection(xhat, b_d)
69607144faaSPeter Brune x = x + nu*(xhat - x);
69739bd7f45SPeter Brune (optional) x = uS(x, b)
69839bd7f45SPeter Brune 
69939bd7f45SPeter Brune With the coarse RHS (defect correction) as below.
70039bd7f45SPeter Brune 
70139bd7f45SPeter Brune  */
7020adebc6cSBarry Smith PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X)
7030adebc6cSBarry Smith {
70407144faaSPeter Brune   Vec                  F, B, Xhat;
70522c1e704SPeter Brune   Vec                  X_c, Xo_c, F_c, B_c;
70639bd7f45SPeter Brune   PetscErrorCode       ierr;
70707144faaSPeter Brune   SNESConvergedReason  reason;
70822c1e704SPeter Brune   PetscReal            xnorm, fnorm, ynorm;
709422a814eSBarry Smith   SNESLineSearchReason lsresult;
710ab8d36c9SPeter Brune   SNES                 next;
711ab8d36c9SPeter Brune   Mat                  restrct, interpolate;
7120dd27c6cSPeter Brune   SNES_FAS             *fas = (SNES_FAS*)snes->data,*fasc;
7130dd27c6cSPeter Brune 
71439bd7f45SPeter Brune   PetscFunctionBegin;
715ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
71639bd7f45SPeter Brune   F    = snes->vec_func;
71739bd7f45SPeter Brune   B    = snes->vec_rhs;
718e7f468e7SPeter Brune   Xhat = snes->work[1];
71907144faaSPeter Brune   ierr = VecCopy(X, Xhat);CHKERRQ(ierr);
72007144faaSPeter Brune   /* recurse first */
721ab8d36c9SPeter Brune   if (next) {
7220dd27c6cSPeter Brune     fasc = (SNES_FAS*)next->data;
723ab8d36c9SPeter Brune     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
724ab8d36c9SPeter Brune     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
7250dd27c6cSPeter Brune     if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);}
72607144faaSPeter Brune     ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr);
7270dd27c6cSPeter Brune     if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);}
728c2a02606SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
729ab8d36c9SPeter Brune     X_c  = next->vec_sol;
730ab8d36c9SPeter Brune     Xo_c = next->work[0];
731ab8d36c9SPeter Brune     F_c  = next->vec_func;
732ab8d36c9SPeter Brune     B_c  = next->vec_rhs;
73339bd7f45SPeter Brune 
734938e4a01SJed Brown     ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr);
73507144faaSPeter Brune     /* restrict the defect */
736ab8d36c9SPeter Brune     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
73707144faaSPeter Brune 
73807144faaSPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
7390dd27c6cSPeter Brune     if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);}
740ab8d36c9SPeter Brune     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
7410dd27c6cSPeter Brune     if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);}
742e4ed7901SPeter Brune     ierr = VecCopy(B_c, X_c);CHKERRQ(ierr);
743b9c2fdf1SPeter Brune     ierr = VecCopy(F_c, B_c);CHKERRQ(ierr);
744e4ed7901SPeter Brune     ierr = VecCopy(X_c, F_c);CHKERRQ(ierr);
74507144faaSPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
74607144faaSPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
74707144faaSPeter Brune 
74807144faaSPeter Brune     /* recurse */
749e4ed7901SPeter Brune     ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr);
750ab8d36c9SPeter Brune     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
75107144faaSPeter Brune 
75207144faaSPeter Brune     /* smooth on this level */
75391f99d7cSPeter Brune     ierr = SNESFASDownSmooth_Private(snes, B, X, F, &fnorm);CHKERRQ(ierr);
75407144faaSPeter Brune 
755ab8d36c9SPeter Brune     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
75607144faaSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
75707144faaSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
75807144faaSPeter Brune       PetscFunctionReturn(0);
75907144faaSPeter Brune     }
76007144faaSPeter Brune 
76107144faaSPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
762c68acad4SPeter Brune     ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr);
763ab8d36c9SPeter Brune     ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr);
76407144faaSPeter Brune 
765ddebd997SPeter Brune     /* additive correction of the coarse direction*/
766f1c6b773SPeter Brune     ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr);
767422a814eSBarry Smith     ierr = SNESLineSearchGetReason(snes->linesearch, &lsresult);CHKERRQ(ierr);
768422a814eSBarry Smith     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm);CHKERRQ(ierr);
769422a814eSBarry Smith     if (lsresult) {
7709e764e56SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
7719e764e56SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
7729e764e56SPeter Brune         PetscFunctionReturn(0);
7739e764e56SPeter Brune       }
7749e764e56SPeter Brune     }
77507144faaSPeter Brune   } else {
77691f99d7cSPeter Brune     ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
77707144faaSPeter Brune   }
77839bd7f45SPeter Brune   PetscFunctionReturn(0);
77939bd7f45SPeter Brune }
78039bd7f45SPeter Brune 
78139bd7f45SPeter Brune #undef __FUNCT__
7822cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Multiplicative"
78339bd7f45SPeter Brune /*
78439bd7f45SPeter Brune 
78539bd7f45SPeter Brune Defines the FAS cycle as:
78639bd7f45SPeter Brune 
787b5527d98SMatthew G. Knepley fine problem: F(x) = b
78839bd7f45SPeter Brune coarse problem: F^c(x) = b^c
78939bd7f45SPeter Brune 
790b5527d98SMatthew G. Knepley b^c = F^c(Rx) - R(F(x) - b)
79139bd7f45SPeter Brune 
79239bd7f45SPeter Brune correction:
79339bd7f45SPeter Brune 
79439bd7f45SPeter Brune x = x + I(x^c - Rx)
79539bd7f45SPeter Brune 
79639bd7f45SPeter Brune  */
7970adebc6cSBarry Smith PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X)
7980adebc6cSBarry Smith {
79939bd7f45SPeter Brune 
80039bd7f45SPeter Brune   PetscErrorCode ierr;
80139bd7f45SPeter Brune   Vec            F,B;
80234d65b3cSPeter Brune   SNES           next;
80339bd7f45SPeter Brune 
80439bd7f45SPeter Brune   PetscFunctionBegin;
80539bd7f45SPeter Brune   F = snes->vec_func;
80639bd7f45SPeter Brune   B = snes->vec_rhs;
80739bd7f45SPeter Brune   /* pre-smooth -- just update using the pre-smoother */
80834d65b3cSPeter Brune   ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr);
80991f99d7cSPeter Brune   ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
81034d65b3cSPeter Brune   if (next) {
8118c40d5fbSBarry Smith     ierr = SNESFASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr);
81291f99d7cSPeter Brune     ierr = SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
813fe6f9142SPeter Brune   }
814fa9694d7SPeter Brune   PetscFunctionReturn(0);
815421d9b32SPeter Brune }
816421d9b32SPeter Brune 
817421d9b32SPeter Brune #undef __FUNCT__
8188c94862eSPeter Brune #define __FUNCT__ "SNESFASCycleSetupPhase_Full"
8198c94862eSPeter Brune PetscErrorCode SNESFASCycleSetupPhase_Full(SNES snes)
8208c94862eSPeter Brune {
8218c94862eSPeter Brune   SNES           next;
8228c94862eSPeter Brune   SNES_FAS       *fas = (SNES_FAS*)snes->data;
8238c94862eSPeter Brune   PetscBool      isFine;
8248c94862eSPeter Brune   PetscErrorCode ierr;
8258c94862eSPeter Brune 
8268c94862eSPeter Brune   PetscFunctionBegin;
8278c94862eSPeter Brune   /* pre-smooth -- just update using the pre-smoother */
8288c94862eSPeter Brune   ierr = SNESFASCycleIsFine(snes,&isFine);CHKERRQ(ierr);
8298c94862eSPeter Brune   ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr);
8308c94862eSPeter Brune   fas->full_stage = 0;
8318c94862eSPeter Brune   if (next) {ierr = SNESFASCycleSetupPhase_Full(next);CHKERRQ(ierr);}
8328c94862eSPeter Brune   PetscFunctionReturn(0);
8338c94862eSPeter Brune }
8348c94862eSPeter Brune 
8358c94862eSPeter Brune #undef __FUNCT__
8368c94862eSPeter Brune #define __FUNCT__ "SNESFASCycle_Full"
8378c94862eSPeter Brune PetscErrorCode SNESFASCycle_Full(SNES snes, Vec X)
8388c94862eSPeter Brune {
8398c94862eSPeter Brune   PetscErrorCode ierr;
8408c94862eSPeter Brune   Vec            F,B;
8418c94862eSPeter Brune   SNES_FAS       *fas = (SNES_FAS*)snes->data;
8428c94862eSPeter Brune   PetscBool      isFine;
8438c94862eSPeter Brune   SNES           next;
8448c94862eSPeter Brune 
8458c94862eSPeter Brune   PetscFunctionBegin;
8468c94862eSPeter Brune   F = snes->vec_func;
8478c94862eSPeter Brune   B = snes->vec_rhs;
8488c94862eSPeter Brune   ierr = SNESFASCycleIsFine(snes,&isFine);CHKERRQ(ierr);
8498c94862eSPeter Brune   ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr);
8508c94862eSPeter Brune 
8518c94862eSPeter Brune   if (isFine) {
8528c94862eSPeter Brune     ierr = SNESFASCycleSetupPhase_Full(snes);CHKERRQ(ierr);
8538c94862eSPeter Brune   }
8548c94862eSPeter Brune 
8558c94862eSPeter Brune   if (fas->full_stage == 0) {
856928e959bSPeter Brune     /* downsweep */
8578c94862eSPeter Brune     if (next) {
8588c94862eSPeter Brune       if (fas->level != 1) next->max_its += 1;
859928e959bSPeter Brune       if (fas->full_downsweep||isFine) {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       if (fas->level != 1) next->max_its -= 1;
8638c94862eSPeter Brune     } else {
864a3a80b83SMatthew G. Knepley       /* The smoother on the coarse level is the coarse solver */
8658c94862eSPeter Brune       ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
8668c94862eSPeter Brune     }
8678c94862eSPeter Brune     fas->full_stage = 1;
8688c94862eSPeter Brune   } else if (fas->full_stage == 1) {
8698c94862eSPeter Brune     if (snes->iter == 0) {ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);}
8708c94862eSPeter Brune     if (next) {
8718c94862eSPeter Brune       ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr);
8728c94862eSPeter Brune       ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
8738c94862eSPeter Brune     }
8748c94862eSPeter Brune   }
8758c94862eSPeter Brune   /* final v-cycle */
8768c94862eSPeter Brune   if (isFine) {
8778c94862eSPeter Brune     if (next) {
8788c94862eSPeter Brune       ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr);
8798c94862eSPeter Brune       ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
8808c94862eSPeter Brune     }
8818c94862eSPeter Brune   }
8828c94862eSPeter Brune   PetscFunctionReturn(0);
8838c94862eSPeter Brune }
8848c94862eSPeter Brune 
8858c94862eSPeter Brune #undef __FUNCT__
88634d65b3cSPeter Brune #define __FUNCT__ "SNESFASCycle_Kaskade"
88734d65b3cSPeter Brune PetscErrorCode SNESFASCycle_Kaskade(SNES snes, Vec X)
88834d65b3cSPeter Brune {
88934d65b3cSPeter Brune   PetscErrorCode ierr;
89034d65b3cSPeter Brune   Vec            F,B;
89134d65b3cSPeter Brune   SNES           next;
89234d65b3cSPeter Brune 
89334d65b3cSPeter Brune   PetscFunctionBegin;
89434d65b3cSPeter Brune   F = snes->vec_func;
89534d65b3cSPeter Brune   B = snes->vec_rhs;
89634d65b3cSPeter Brune   ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr);
89734d65b3cSPeter Brune   if (next) {
89834d65b3cSPeter Brune     ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr);
89934d65b3cSPeter Brune     ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
90034d65b3cSPeter Brune   } else {
90134d65b3cSPeter Brune     ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
90234d65b3cSPeter Brune   }
90334d65b3cSPeter Brune   PetscFunctionReturn(0);
90434d65b3cSPeter Brune }
90534d65b3cSPeter Brune 
906fffbeea8SBarry Smith PetscBool SNEScite = PETSC_FALSE;
907fffbeea8SBarry Smith const char SNESCitation[] = "@techreport{pbmkbsxt2012,\n"
908fffbeea8SBarry Smith                             "  title = {Composing Scalable Nonlinear Algebraic Solvers},\n"
909fffbeea8SBarry Smith                             "  author = {Peter Brune and Mathew Knepley and Barry Smith and Xuemin Tu},\n"
910fffbeea8SBarry Smith                             "  year = 2013,\n"
911fffbeea8SBarry Smith                             "  type = Preprint,\n"
912fffbeea8SBarry Smith                             "  number = {ANL/MCS-P2010-0112},\n"
913fffbeea8SBarry Smith                             "  institution = {Argonne National Laboratory}\n}\n";
914fffbeea8SBarry Smith 
91534d65b3cSPeter Brune #undef __FUNCT__
916421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS"
917421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes)
918421d9b32SPeter Brune {
919fa9694d7SPeter Brune   PetscErrorCode ierr;
920fe6f9142SPeter Brune   PetscInt       i, maxits;
921ddb5aff1SPeter Brune   Vec            X, F;
922fe6f9142SPeter Brune   PetscReal      fnorm;
923b17ce1afSJed Brown   SNES_FAS       *fas = (SNES_FAS*)snes->data,*ffas;
924b17ce1afSJed Brown   DM             dm;
925e70c42e5SPeter Brune   PetscBool      isFine;
926b17ce1afSJed Brown 
927421d9b32SPeter Brune   PetscFunctionBegin;
928c579b300SPatrick Farrell 
9296c4ed002SBarry Smith   if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
930c579b300SPatrick Farrell 
931fffbeea8SBarry Smith   ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr);
932fe6f9142SPeter Brune   maxits       = snes->max_its;      /* maximum number of iterations */
933fe6f9142SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
934fa9694d7SPeter Brune   X            = snes->vec_sol;
935f5a6d4f9SBarry Smith   F            = snes->vec_func;
936293a7e31SPeter Brune 
93718a66777SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
938293a7e31SPeter Brune   /*norm setup */
939e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
940fe6f9142SPeter Brune   snes->iter = 0;
941fe6f9142SPeter Brune   snes->norm = 0.;
942e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
943e4ed7901SPeter Brune   if (!snes->vec_func_init_set) {
9440dd27c6cSPeter Brune     if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);}
945fe6f9142SPeter Brune     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
9460dd27c6cSPeter Brune     if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);}
9471aa26658SKarl Rupp   } else snes->vec_func_init_set = PETSC_FALSE;
948e4ed7901SPeter Brune 
949fe6f9142SPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
950422a814eSBarry Smith   SNESCheckFunctionNorm(snes,fnorm);
951e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
952fe6f9142SPeter Brune   snes->norm = fnorm;
953e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
954a71f0d7dSBarry Smith   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
955fe6f9142SPeter Brune   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
956fe6f9142SPeter Brune 
957fe6f9142SPeter Brune   /* test convergence */
958fe6f9142SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
959fe6f9142SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
960e4ed7901SPeter Brune 
961b17ce1afSJed Brown 
962b9c2fdf1SPeter Brune   if (isFine) {
963b9c2fdf1SPeter Brune     /* propagate scale-dependent data up the hierarchy */
964b17ce1afSJed Brown     ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
965b17ce1afSJed Brown     for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) {
966b17ce1afSJed Brown       DM dmcoarse;
967b17ce1afSJed Brown       ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr);
968b17ce1afSJed Brown       ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr);
969b17ce1afSJed Brown       dm   = dmcoarse;
970b17ce1afSJed Brown     }
971b9c2fdf1SPeter Brune   }
972b17ce1afSJed Brown 
973fe6f9142SPeter Brune   for (i = 0; i < maxits; i++) {
974fe6f9142SPeter Brune     /* Call general purpose update function */
975646217ecSPeter Brune 
976fe6f9142SPeter Brune     if (snes->ops->update) {
977fe6f9142SPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
978fe6f9142SPeter Brune     }
97907144faaSPeter Brune     if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
98091f99d7cSPeter Brune       ierr = SNESFASCycle_Multiplicative(snes, X);CHKERRQ(ierr);
9818c94862eSPeter Brune     } else if (fas->fastype == SNES_FAS_ADDITIVE) {
98291f99d7cSPeter Brune       ierr = SNESFASCycle_Additive(snes, X);CHKERRQ(ierr);
9838c94862eSPeter Brune     } else if (fas->fastype == SNES_FAS_FULL) {
9848c94862eSPeter Brune       ierr = SNESFASCycle_Full(snes, X);CHKERRQ(ierr);
98534d65b3cSPeter Brune     } else if (fas->fastype ==SNES_FAS_KASKADE) {
98634d65b3cSPeter Brune       ierr = SNESFASCycle_Kaskade(snes, X);CHKERRQ(ierr);
9876c4ed002SBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported FAS type");
988742fe5e2SPeter Brune 
989742fe5e2SPeter Brune     /* check for FAS cycle divergence */
9901aa26658SKarl Rupp     if (snes->reason != SNES_CONVERGED_ITERATING) PetscFunctionReturn(0);
991b9c2fdf1SPeter Brune 
992c90fad12SPeter Brune     /* Monitor convergence */
993e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
994c90fad12SPeter Brune     snes->iter = i+1;
995e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
996a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
997c90fad12SPeter Brune     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
998c90fad12SPeter Brune     /* Test for convergence */
99966585501SPeter Brune     if (isFine) {
1000b9c2fdf1SPeter Brune       ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,snes->norm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1001c90fad12SPeter Brune       if (snes->reason) break;
1002fe6f9142SPeter Brune     }
100366585501SPeter Brune   }
1004fe6f9142SPeter Brune   if (i == maxits) {
1005fe6f9142SPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
1006fe6f9142SPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
1007fe6f9142SPeter Brune   }
1008421d9b32SPeter Brune   PetscFunctionReturn(0);
1009421d9b32SPeter Brune }
1010