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
 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
 21: -  PC_MG_KASKADE - like full multigrid except one never goes back to a coarser level
 22:                from a finer

 24: .seealso: PCMGSetType()

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

 31: #define PC_MG_V_CYCLE     1
 32: #define PC_MG_W_CYCLE     2

 34: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetType(PC,PCMGType);
 35: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetLevels(PC,PetscInt,MPI_Comm*);
 36: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetLevels(PC,PetscInt*);

 38: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothUp(PC,PetscInt);
 39: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothDown(PC,PetscInt);
 40: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetCycles(PC,PetscInt);
 41: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetCyclesOnLevel(PC,PetscInt,PetscInt);
 42: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetGalerkin(PC);
 43: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetGalerkin(PC,PetscTruth*);

 45: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetSmoother(PC,PetscInt,KSP*);
 46: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetSmootherDown(PC,PetscInt,KSP*);
 47: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetSmootherUp(PC,PetscInt,KSP*);
 48: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetCoarseSolve(PC,KSP*);


 51: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetRhs(PC,PetscInt,Vec);
 52: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetX(PC,PetscInt,Vec);
 53: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetR(PC,PetscInt,Vec);

 55: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetRestriction(PC,PetscInt,Mat);
 56: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetInterpolate(PC,PetscInt,Mat);
 57: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetResidual(PC,PetscInt,PetscErrorCode (*)(Mat,Vec,Vec,Vec),Mat);
 58: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCMGDefaultResidual(Mat,Vec,Vec,Vec);

 61: #endif