Actual source code: petscmg.h

  1: /*
  2:       Structure used for Multigrid preconditioners 
  3: */
 6:  #include petscksp.h

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

 12:    Level: beginner

 14:    Values:
 15: +  PC_MG_MULTIPLICATIVE (default) - traditional V or W cycle as determined by PCMGSetCycles()
 16: .  PC_MG_ADDITIVE - the additive multigrid preconditioner where all levels are
 17:                 smoothed before updating the residual. This only uses the 
 18:                 down smoother, in the preconditioner the upper smoother is ignored
 19: .  PC_MG_FULL - same as multiplicative except one also performs grid sequencing, 
 20:             that is starts on the coarsest grid, performs a cycle, interpolates
 21:             to the next, performs a cycle etc
 22: -  PC_MG_KASKADE - like full multigrid except one never goes back to a coarser level
 23:                from a finer

 25: .seealso: PCMGSetType()

 27: E*/
 28: typedef enum { PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE } PCMGType;
 30: #define PC_MG_CASCADE PC_MG_KASKADE;

 32: /*E
 33:     PCMGCyleType - Use V-cycle or W-cycle

 35:    Level: beginner

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

 41: .seealso: PCMGSetCycleType()

 43: E*/
 44: typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType;

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

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

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


 66: EXTERN PetscErrorCode  PCMGSetRhs(PC,PetscInt,Vec);
 67: EXTERN PetscErrorCode  PCMGSetX(PC,PetscInt,Vec);
 68: EXTERN PetscErrorCode  PCMGSetR(PC,PetscInt,Vec);

 70: EXTERN PetscErrorCode  PCMGSetRestriction(PC,PetscInt,Mat);
 71: EXTERN PetscErrorCode  PCMGSetInterpolate(PC,PetscInt,Mat);
 72: EXTERN PetscErrorCode  PCMGSetResidual(PC,PetscInt,PetscErrorCode (*)(Mat,Vec,Vec,Vec),Mat);
 73: EXTERN PetscErrorCode  PCMGDefaultResidual(Mat,Vec,Vec,Vec);

 76: #endif