Actual source code: pcmgimpl.h
petsc-3.11.4 2019-09-28
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: Mat inject; /* Used for moving state if provided. */
29: Vec rscale; /* scaling of restriction matrix */
30: PetscLogEvent eventsmoothsetup; /* if logging times for each level */
31: PetscLogEvent eventsmoothsolve;
32: PetscLogEvent eventresidual;
33: PetscLogEvent eventinterprestrict;
34: } PC_MG_Levels;
36: /*
37: This data structure is shared by all the levels.
38: */
39: typedef struct {
40: PCMGType am; /* Multiplicative, additive or full */
41: PetscInt cyclesperpcapply; /* Number of cycles to use in each PCApply(), multiplicative only*/
42: PetscInt maxlevels; /* total number of levels allocated */
43: PCMGGalerkinType galerkin; /* use Galerkin process to compute coarser matrices */
44: PetscBool usedmfornumberoflevels; /* sets the number of levels by getting this information out of the DM */
46: PetscInt nlevels;
47: PC_MG_Levels **levels;
48: PetscInt default_smoothu; /* number of smooths per level if not over-ridden */
49: PetscInt default_smoothd; /* with calls to KSPSetTolerances() */
50: PetscReal rtol,abstol,dtol,ttol; /* tolerances for when running with PCApplyRichardson_MG */
52: void *innerctx; /* optional data for preconditioner, like PCEXOTIC that inherits off of PCMG */
53: PetscLogStage stageApply;
54: PetscErrorCode (*view)(PC,PetscViewer); /* GAMG and other objects that use PCMG can set their own viewer here */
55: } PC_MG;
57: PETSC_INTERN PetscErrorCode PCSetUp_MG(PC);
58: PETSC_INTERN PetscErrorCode PCDestroy_MG(PC);
59: PETSC_INTERN PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC);
60: PETSC_INTERN PetscErrorCode PCView_MG(PC,PetscViewer);
61: PETSC_INTERN PetscErrorCode PCMGGetLevels_MG(PC,PetscInt *);
62: PETSC_INTERN PetscErrorCode PCMGSetLevels_MG(PC,PetscInt,MPI_Comm *);
63: PETSC_DEPRECATED("Use PCMGResidualDefault()") PETSC_STATIC_INLINE PetscErrorCode PCMGResidual_Default(Mat A,Vec b,Vec x,Vec r) {
64: return PCMGResidualDefault(A,b,x,r);
65: }
67: #endif