Actual source code: fasimpls.h

  1: #ifndef PETSC_SNES_FASIMPLS_H
  2: #define PETSC_SNES_FASIMPLS_H

  4: #include <petsc/private/snesimpl.h>
  5: #include <petsc/private/linesearchimpl.h>
  6: #include <petscdm.h>

  8: typedef struct {
  9:   /* flags for knowing the global place of this FAS object */
 10:   PetscInt level;  /* level = 0 coarsest level */
 11:   PetscInt levels; /* if level + 1 = levels; we're the last turtle */

 13:   /* smoothing objects */
 14:   SNES smoothu; /* the SNES for presmoothing */
 15:   SNES smoothd; /* the SNES for postsmoothing */

 17:   /* coarse grid correction objects */
 18:   SNES next;        /* the SNES instance for the next coarser level in the hierarchy */
 19:   SNES fine;        /* the finest SNES instance; used as a reference for prefixes */
 20:   SNES previous;    /* the SNES instance for the next finer level in the hierarchy */
 21:   Mat  interpolate; /* interpolation */
 22:   Mat  inject;      /* injection operator (unscaled) */
 23:   Mat  restrct;     /* restriction operator */
 24:   Vec  rscale;      /* the pointwise scaling of the restriction operator */

 26:   /* method parameters */
 27:   PetscInt    n_cycles;               /* number of cycles on this level */
 28:   SNESFASType fastype;                /* FAS type */
 29:   PetscInt    max_up_it;              /* number of pre-smooths */
 30:   PetscInt    max_down_it;            /* number of post-smooth cycles */
 31:   PetscBool   usedmfornumberoflevels; /* uses a DM to generate a number of the levels */
 32:   PetscBool   full_downsweep;         /* smooth on the initial full downsweep */
 33:   PetscBool   full_total;             /* use total residual restriction and total solution interpolation on the initial downsweep and upsweep */
 34:   PetscBool   continuation;           /* sets the setup to default to continuation */
 35:   PetscInt    full_stage;             /* stage of the full cycle -- 0 is the upswing, 1 is the downsweep and final V-cycle */

 37:   /* Galerkin FAS state */
 38:   PetscBool galerkin; /* use Galerkin formation of the coarse problem */
 39:   Vec       Xg;       /* Galerkin solution projection */
 40:   Vec       Fg;       /* Galerkin function projection */

 42:   /* if logging times for each level */
 43:   PetscLogEvent eventsmoothsetup;    /* level setup */
 44:   PetscLogEvent eventsmoothsolve;    /* level smoother solves */
 45:   PetscLogEvent eventresidual;       /* level residual evaluation */
 46:   PetscLogEvent eventinterprestrict; /* level interpolation and restriction */
 47: } SNES_FAS;

 49: PETSC_INTERN PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES, SNES *);

 51: #endif // PETSC_SNES_FASIMPLS_H