Actual source code: petscpcmg.h

petsc-3.3-p7 2013-05-11
  1: /*
  2:       Structure used for Multigrid preconditioners 
  3: */
 6:  #include petscksp.h

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

 11:    Level: beginner

 13:    Values:
 14: +  PC_MG_MULTIPLICATIVE (default) - traditional V or W cycle as determined by PCMGSetCycles()
 15: .  PC_MG_ADDITIVE - the additive multigrid preconditioner where all levels are
 16:                 smoothed before updating the residual. This only uses the 
 17:                 down smoother, in the preconditioner the upper smoother is ignored
 18: .  PC_MG_FULL - 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. This is much like the F-cycle presented in "Multigrid" by Trottenberg, Oosterlee, Schuller page 49, but that
 21:             algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm
 22:             calls the V-cycle only on the coarser level and has a post-smoother instead.
 23: -  PC_MG_KASKADE - like full multigrid except one never goes back to a coarser level
 24:                from a finer

 26: .seealso: PCMGSetType()

 28: E*/
 29: typedef enum { PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE } PCMGType;
 30: PETSC_EXTERN const char *PCMGTypes[];
 31: #define PC_MG_CASCADE PC_MG_KASKADE;

 33: /*E
 34:     PCMGCycleType - Use V-cycle or W-cycle

 36:    Level: beginner

 38:    Values:
 39: +  PC_MG_V_CYCLE
 40: -  PC_MG_W_CYCLE

 42: .seealso: PCMGSetCycleType()

 44: E*/
 45: typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType;
 46: PETSC_EXTERN const char *PCMGCycleTypes[];

 48: PETSC_EXTERN PetscErrorCode PCMGSetType(PC,PCMGType);
 49: PETSC_EXTERN PetscErrorCode PCMGSetLevels(PC,PetscInt,MPI_Comm*);
 50: PETSC_EXTERN PetscErrorCode PCMGGetLevels(PC,PetscInt*);

 52: PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothUp(PC,PetscInt);
 53: PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothDown(PC,PetscInt);
 54: PETSC_EXTERN PetscErrorCode PCMGSetCycleType(PC,PCMGCycleType);
 55: PETSC_EXTERN PetscErrorCode PCMGSetCycleTypeOnLevel(PC,PetscInt,PCMGCycleType);
 56: PETSC_EXTERN PetscErrorCode PCMGSetCyclesOnLevel(PC,PetscInt,PetscInt);
 57: PETSC_EXTERN PetscErrorCode PCMGMultiplicativeSetCycles(PC,PetscInt);
 58: PETSC_EXTERN PetscErrorCode PCMGSetGalerkin(PC,PetscBool);
 59: PETSC_EXTERN PetscErrorCode PCMGGetGalerkin(PC,PetscBool*);

 61: PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*);
 62: PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*);
 63: PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*);
 64: PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*);


 67: PETSC_EXTERN PetscErrorCode PCMGSetRhs(PC,PetscInt,Vec);
 68: PETSC_EXTERN PetscErrorCode PCMGSetX(PC,PetscInt,Vec);
 69: PETSC_EXTERN PetscErrorCode PCMGSetR(PC,PetscInt,Vec);

 71: PETSC_EXTERN PetscErrorCode PCMGSetRestriction(PC,PetscInt,Mat);
 72: PETSC_EXTERN PetscErrorCode PCMGGetRestriction(PC,PetscInt,Mat*);
 73: PETSC_EXTERN PetscErrorCode PCMGSetInterpolation(PC,PetscInt,Mat);
 74: PETSC_EXTERN PetscErrorCode PCMGGetInterpolation(PC,PetscInt,Mat*);
 75: PETSC_EXTERN PetscErrorCode PCMGSetRScale(PC,PetscInt,Vec);
 76: PETSC_EXTERN PetscErrorCode PCMGGetRScale(PC,PetscInt,Vec*);
 77: PETSC_EXTERN PetscErrorCode PCMGSetResidual(PC,PetscInt,PetscErrorCode (*)(Mat,Vec,Vec,Vec),Mat);
 78: PETSC_EXTERN PetscErrorCode PCMGDefaultResidual(Mat,Vec,Vec,Vec);

 80: /*E
 81:     PCExoticType - Face based or wirebasket based coarse grid space

 83:    Level: beginner

 85: .seealso: PCExoticSetType(), PCEXOTIC
 86: E*/
 87: typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType;
 88: PETSC_EXTERN const char *PCExoticTypes[];
 89: PETSC_EXTERN PetscErrorCode PCExoticSetType(PC,PCExoticType);


 92: #endif