Actual source code: pcmgimpl.h

  1: /*
  2:       Data structure used for Multigrid preconditioner.
  3: */
  4: #pragma once

  6: #include <petsc/private/pcimpl.h>
  7: #include <petscksp.h>
  8: #include <petscpctypes.h>
  9: #define PETSC_MG_MAXLEVELS 10
 10: /*
 11:      Each level has its own copy of this data.
 12:      Level (0) is always the coarsest level and Level (levels-1) is the finest.
 13: */
 14: typedef struct {
 15:   PetscInt cycles; /* Type of cycle to run: 1 V 2 W */
 16:   PetscInt level;  /* level = 0 coarsest level */
 17:   PetscInt levels; /* number of active levels used */
 18:   Vec      b;      /* Right hand side */
 19:   Vec      x;      /* Solution */
 20:   Vec      r;      /* Residual */
 21:   Mat      B;
 22:   Mat      X;
 23:   Mat      R;
 24:   Mat      coarseSpace; /* A vector space which should be accurately captured by the next coarser mesh,
 25:                                                   and thus accurately interpolated. The columns of this dense matrix
 26:                                                   correspond to the same function discretized in
 27:                                                   the sequence of spaces. */

 29:   PetscErrorCode (*residual)(Mat, Vec, Vec, Vec);
 30:   PetscErrorCode (*residualtranspose)(Mat, Vec, Vec, Vec);
 31:   PetscErrorCode (*matresidual)(Mat, Mat, Mat, Mat);
 32:   PetscErrorCode (*matresidualtranspose)(Mat, Mat, Mat, Mat);

 34:   Mat           A;       /* matrix used in forming residual*/
 35:   KSP           smoothd; /* pre smoother */
 36:   KSP           smoothu; /* post smoother */
 37:   KSP           cr;      /* post compatible relaxation (cr) */
 38:   Vec           crx;     /* cr solution */
 39:   Vec           crb;     /* cr rhs */
 40:   Mat           interpolate;
 41:   Mat           restrct;          /* restrict is a reserved word in C99 and on Cray */
 42:   Mat           inject;           /* Used for moving state if provided. */
 43:   Vec           rscale;           /* scaling of restriction matrix */
 44:   PetscLogEvent eventsmoothsetup; /* if logging times for each level */
 45:   PetscLogEvent eventsmoothsolve;
 46:   PetscLogEvent eventresidual;
 47:   PetscLogEvent eventinterprestrict;
 48: } PC_MG_Levels;

 50: /*
 51:     This data structure is shared by all the levels.
 52: */
 53: typedef struct {
 54:   PCMGType         am;                     /* Multiplicative, additive or full */
 55:   PetscInt         cyclesperpcapply;       /* Number of cycles to use in each PCApply(), multiplicative only*/
 56:   PetscInt         maxlevels;              /* total number of levels allocated */
 57:   PCMGGalerkinType galerkin;               /* use Galerkin process to compute coarser matrices */
 58:   PetscBool        usedmfornumberoflevels; /* sets the number of levels by getting this information out of the DM */

 60:   PetscBool           adaptInterpolation; /* flag to adapt the interpolator based upon the coarseSpace */
 61:   PCMGCoarseSpaceType coarseSpaceType;    /* Type of coarse space: polynomials, harmonics, eigenvectors, ... */
 62:   PetscInt            Nc;                 /* The number of vectors in coarseSpace */
 63:   PetscInt            eigenvalue;         /* Key for storing the eigenvalue as a scalar in the eigenvector Vec */
 64:   PetscBool           mespMonitor;        /* flag to monitor the multilevel eigensolver */

 66:   PetscBool compatibleRelaxation; /* flag to monitor the coarse space quality using an auxiliary solve with compatible relaxation */

 68:   PetscInt       nlevels;
 69:   PC_MG_Levels **levels;
 70:   PetscInt       default_smoothu;          /* number of smooths per level if not over-ridden */
 71:   PetscInt       default_smoothd;          /*  with calls to KSPSetTolerances() */
 72:   PetscReal      rtol, abstol, dtol, ttol; /* tolerances for when running with PCApplyRichardson_MG */

 74:   void         *innerctx; /* optional data for preconditioner, like PCEXOTIC that inherits off of PCMG */
 75:   PetscLogStage stageApply;
 76:   PetscErrorCode (*view)(PC, PetscViewer); /* GAMG and other objects that use PCMG can set their own viewer here */
 77:   PetscReal min_eigen_DinvA[PETSC_MG_MAXLEVELS];
 78:   PetscReal max_eigen_DinvA[PETSC_MG_MAXLEVELS];
 79: } PC_MG;

 81: PETSC_INTERN PetscErrorCode PCSetUp_MG(PC);
 82: PETSC_INTERN PetscErrorCode PCDestroy_MG(PC);
 83: PETSC_INTERN PetscErrorCode PCSetFromOptions_MG(PC, PetscOptionItems *PetscOptionsObject);
 84: PETSC_INTERN PetscErrorCode PCView_MG(PC, PetscViewer);
 85: PETSC_INTERN PetscErrorCode PCMGGetLevels_MG(PC, PetscInt *);
 86: PETSC_INTERN PetscErrorCode PCMGSetLevels_MG(PC, PetscInt, MPI_Comm *);
 87: PETSC_DEPRECATED_FUNCTION(3, 5, 0, "PCMGResidualDefault()", ) static inline PetscErrorCode PCMGResidual_Default(Mat A, Vec b, Vec x, Vec r)
 88: {
 89:   return PCMGResidualDefault(A, b, x, r);
 90: }

 92: PETSC_INTERN PetscErrorCode DMSetBasisFunction_Internal(PetscInt, PetscBool, PetscInt, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *));
 93: PETSC_INTERN PetscErrorCode PCMGComputeCoarseSpace_Internal(PC, PetscInt, PCMGCoarseSpaceType, PetscInt, Mat, Mat *);
 94: PETSC_INTERN PetscErrorCode PCMGAdaptInterpolator_Internal(PC, PetscInt, KSP, KSP, Mat, Mat);
 95: PETSC_INTERN PetscErrorCode PCMGRecomputeLevelOperators_Internal(PC, PetscInt);
 96: PETSC_INTERN PetscErrorCode PCMGACycle_Private(PC, PC_MG_Levels **, PetscBool, PetscBool);
 97: PETSC_INTERN PetscErrorCode PCMGFCycle_Private(PC, PC_MG_Levels **, PetscBool, PetscBool);
 98: PETSC_INTERN PetscErrorCode PCMGKCycle_Private(PC, PC_MG_Levels **, PetscBool, PetscBool);
 99: PETSC_INTERN PetscErrorCode PCMGMCycle_Private(PC, PC_MG_Levels **, PetscBool, PetscBool, PCRichardsonConvergedReason *);

101: PETSC_INTERN PetscErrorCode PCMGGDSWCreateCoarseSpace_Private(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *);