xref: /petsc/src/snes/impls/fas/fasimpls.h (revision 6dfbafefc36c5de14c2c975d80905386b64e67ea)
1519f805aSKarl Rupp #if !defined(_SNES_FASIMPLS)
2421d9b32SPeter Brune #define _SNES_FASIMPLS
3421d9b32SPeter Brune 
4af0996ceSBarry Smith #include <petsc/private/snesimpl.h>
5af0996ceSBarry Smith #include <petsc/private/linesearchimpl.h>
6539c167fSBarry Smith #include <petscdm.h>
7421d9b32SPeter Brune 
8421d9b32SPeter Brune typedef struct {
9421d9b32SPeter Brune 
10421d9b32SPeter Brune   /* flags for knowing the global place of this FAS object */
11421d9b32SPeter Brune   PetscInt level;                              /* level = 0 coarsest level */
12421d9b32SPeter Brune   PetscInt levels;                             /* if level + 1 = levels; we're the last turtle */
13421d9b32SPeter Brune 
14421d9b32SPeter Brune   /* smoothing objects */
15ab8d36c9SPeter Brune   SNES smoothu;                                /* the SNES for presmoothing */
16ab8d36c9SPeter Brune   SNES smoothd;                                /* the SNES for postsmoothing */
1722c1e704SPeter Brune 
18421d9b32SPeter Brune   /* coarse grid correction objects */
196273346dSPeter Brune   SNES next;                                   /* the SNES instance for the next coarser level in the hierarchy */
20ab8d36c9SPeter Brune   SNES fine;                                   /* the finest SNES instance; used as a reference for prefixes */
216273346dSPeter Brune   SNES previous;                               /* the SNES instance for the next finer level in the hierarchy */
22421d9b32SPeter Brune   Mat  interpolate;                            /* interpolation */
23efe1f98aSPeter Brune   Mat  inject;                                 /* injection operator (unscaled) */
24421d9b32SPeter Brune   Mat  restrct;                                /* restriction operator */
25421d9b32SPeter Brune   Vec  rscale;                                 /* the pointwise scaling of the restriction operator */
26421d9b32SPeter Brune 
27ee78dd50SPeter Brune   /* method parameters */
28ee78dd50SPeter Brune   PetscInt    n_cycles;                        /* number of cycles on this level */
2907144faaSPeter Brune   SNESFASType fastype;                         /* FAS type */
30ee78dd50SPeter Brune   PetscInt    max_up_it;                       /* number of pre-smooths */
31ee78dd50SPeter Brune   PetscInt    max_down_it;                     /* number of post-smooth cycles */
32d1adcc6fSPeter Brune   PetscBool   usedmfornumberoflevels;          /* uses a DM to generate a number of the levels */
33928e959bSPeter Brune   PetscBool   full_downsweep;                  /* smooth on the initial full downsweep */
34*6dfbafefSToby Isaac   PetscBool   full_total;                      /* use total residual restriction and total solution interpolation on the initial downsweep and upsweep */
3587f44e3fSPeter Brune   PetscBool   continuation;                    /* sets the setup to default to continuation */
368c94862eSPeter Brune   PetscInt    full_stage;                      /* stage of the full cycle -- 0 is the upswing, 1 is the downsweep and final V-cycle */
376273346dSPeter Brune 
386273346dSPeter Brune   /* Galerkin FAS state */
396273346dSPeter Brune   PetscBool galerkin;                          /* use Galerkin formation of the coarse problem */
406273346dSPeter Brune   Vec       Xg;                                /* Galerkin solution projection */
416273346dSPeter Brune   Vec       Fg;                                /* Galerkin function projection */
426273346dSPeter Brune 
430dd27c6cSPeter Brune   /* if logging times for each level */
440dd27c6cSPeter Brune   PetscLogEvent eventsmoothsetup;              /* level setup */
450dd27c6cSPeter Brune   PetscLogEvent eventsmoothsolve;              /* level smoother solves */
460dd27c6cSPeter Brune   PetscLogEvent eventresidual;                 /* level residual evaluation */
470dd27c6cSPeter Brune   PetscLogEvent eventinterprestrict;           /* level interpolation and restriction */
48421d9b32SPeter Brune } SNES_FAS;
49421d9b32SPeter Brune 
501943db53SBarry Smith PETSC_INTERN PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES,SNES*);
511943db53SBarry Smith 
52421d9b32SPeter Brune #endif
53