Actual source code: pcmgimpl.h

petsc-3.8.4 2018-03-24
Report Typos and Errors
  1: /*
  2:       Data structure used for Multigrid preconditioner.
  3: */
  6:  #include <petsc/private/pcimpl.h>
  7:  #include <petscksp.h>

  9: /*
 10:      Each level has its own copy of this data.
 11:      Level (0) is always the coarsest level and Level (levels-1) is the finest.
 12: */
 13: typedef struct {
 14:   PetscInt cycles;                             /* Type of cycle to run: 1 V 2 W */
 15:   PetscInt level;                              /* level = 0 coarsest level */
 16:   PetscInt levels;                             /* number of active levels used */
 17:   Vec      b;                                  /* Right hand side */
 18:   Vec      x;                                  /* Solution */
 19:   Vec      r;                                  /* Residual */

 21:   PetscErrorCode (*residual)(Mat,Vec,Vec,Vec);

 23:   Mat           A;                             /* matrix used in forming residual*/
 24:   KSP           smoothd;                       /* pre smoother */
 25:   KSP           smoothu;                       /* post smoother */
 26:   Mat           interpolate;
 27:   Mat           restrct;                       /* restrict is a reserved word in C99 and on Cray */
 28:   Vec           rscale;                        /* scaling of restriction matrix */
 29:   PetscLogEvent eventsmoothsetup;              /* if logging times for each level */
 30:   PetscLogEvent eventsmoothsolve;
 31:   PetscLogEvent eventresidual;
 32:   PetscLogEvent eventinterprestrict;
 33: } PC_MG_Levels;

 35: /*
 36:     This data structure is shared by all the levels.
 37: */
 38: typedef struct {
 39:   PCMGType         am;                        /* Multiplicative, additive or full */
 40:   PetscInt         cyclesperpcapply;          /* Number of cycles to use in each PCApply(), multiplicative only*/
 41:   PetscInt         maxlevels;                 /* total number of levels allocated */
 42:   PCMGGalerkinType galerkin;                  /* use Galerkin process to compute coarser matrices */
 43:   PetscBool        usedmfornumberoflevels;    /* sets the number of levels by getting this information out of the DM */

 45:   PetscInt     nlevels;
 46:   PC_MG_Levels **levels;
 47:   PetscInt     default_smoothu;               /* number of smooths per level if not over-ridden */
 48:   PetscInt     default_smoothd;               /*  with calls to KSPSetTolerances() */
 49:   PetscReal    rtol,abstol,dtol,ttol;         /* tolerances for when running with PCApplyRichardson_MG */

 51:   void          *innerctx;                    /* optional data for preconditioner, like PCEXOTIC that inherits off of PCMG */
 52:   PetscLogStage stageApply;
 53:   PetscErrorCode (*view)(PC,PetscViewer);     /* GAMG and other objects that use PCMG can set their own viewer here */
 54: } PC_MG;

 56: PETSC_INTERN PetscErrorCode PCSetUp_MG(PC);
 57: PETSC_INTERN PetscErrorCode PCDestroy_MG(PC);
 58: PETSC_INTERN PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC);
 59: PETSC_INTERN PetscErrorCode PCView_MG(PC,PetscViewer);
 60: PETSC_DEPRECATED("Use PCMGResidualDefault()") PETSC_STATIC_INLINE PetscErrorCode PCMGResidual_Default(Mat A,Vec b,Vec x,Vec r) {
 61:   return PCMGResidualDefault(A,b,x,r);
 62: }

 64: #endif