Actual source code: petscmg.h

  1: /* $Id: petscmg.h,v 1.21 2001/04/04 21:07:36 bsmith Exp $ */
  2: /*
  3:       Structure used for Multigrid preconditioners 
  4: */
 7:  #include petscsles.h

  9: /*E
 10:     MGType - Determines the type of multigrid method that is run.

 12:    Level: beginner

 14:    Values:
 15: +  MGMULTIPLICATIVE (default) - traditional V or W cycle as determined by MGSetCycles()
 16: .  MGADDITIVE - the additive multigrid preconditioner where all levels are
 17:                 smoothed before updating the residual
 18: .  MGFULL - same as multiplicative except one also performs grid sequencing, 
 19:             that is starts on the coarsest grid, performs a cycle, interpolates
 20:             to the next, performs a cycle etc
 21: -  MGKASKADE - like full multigrid except one never goes back to a coarser level
 22:                from a finer

 24: .seealso: MGSetType()

 26: E*/
 27: typedef enum { MGMULTIPLICATIVE,MGADDITIVE,MGFULL,MGKASKADE } MGType;
 28: #define MGCASCADE MGKASKADE;

 30: #define MG_V_CYCLE     1
 31: #define MG_W_CYCLE     2

 33: EXTERN int MGSetType(PC,MGType);
 34: EXTERN int MGCheck(PC);
 35: EXTERN int MGSetLevels(PC,int,MPI_Comm*);
 36: EXTERN int MGGetLevels(PC,int*);

 38: EXTERN int MGSetNumberSmoothUp(PC,int);
 39: EXTERN int MGSetNumberSmoothDown(PC,int);
 40: EXTERN int MGSetCycles(PC,int);
 41: EXTERN int MGSetCyclesOnLevel(PC,int,int);

 43: EXTERN int MGGetSmoother(PC,int,SLES*);
 44: EXTERN int MGGetSmootherDown(PC,int,SLES*);
 45: EXTERN int MGGetSmootherUp(PC,int,SLES*);
 46: EXTERN int MGGetCoarseSolve(PC,SLES*);

 48: EXTERN int MGSetRhs(PC,int,Vec);
 49: EXTERN int MGSetX(PC,int,Vec);
 50: EXTERN int MGSetR(PC,int,Vec);

 52: EXTERN int MGSetRestriction(PC,int,Mat);
 53: EXTERN int MGSetInterpolate(PC,int,Mat);
 54: EXTERN int MGSetResidual(PC,int,int (*)(Mat,Vec,Vec,Vec),Mat);
 55: EXTERN int MGDefaultResidual(Mat,Vec,Vec,Vec);


 58: #endif