xref: /petsc/src/ksp/pc/impls/bjacobi/bjacobi.h (revision 5a7e1895ba1221176dc63d0bbec6b671c7ff4efe)
14b9ad928SBarry Smith 
24b9ad928SBarry Smith #if !defined(__BJACOBI_H)
34b9ad928SBarry Smith #define __BJACOBI_H
44b9ad928SBarry Smith /*
54b9ad928SBarry Smith     Private data for block Jacobi and block Gauss-Seidel preconditioner.
64b9ad928SBarry Smith */
74b9ad928SBarry Smith #include "petscksp.h"
86356e834SBarry Smith #include "private/pcimpl.h"
94b9ad928SBarry Smith 
104b9ad928SBarry Smith /*
114b9ad928SBarry Smith        This data is general for all implementations
124b9ad928SBarry Smith */
134b9ad928SBarry Smith typedef struct {
14*5a7e1895SHong Zhang   PetscInt   n;                 /* number of global blocks */
15*5a7e1895SHong Zhang   PetscInt   n_local;           /* number of blocks in this subcommunicator or in this process */
1613f74950SBarry Smith   PetscInt   first_local;       /* number of first block on processor */
174b9ad928SBarry Smith   PetscTruth use_true_local;    /* use block from true matrix, not preconditioner matrix for local MatMult() */
184b9ad928SBarry Smith   KSP        *ksp;              /* KSP contexts for blocks */
194b9ad928SBarry Smith   void       *data;             /* implementation-specific data */
204b9ad928SBarry Smith   PetscTruth same_local_solves; /* flag indicating whether all local solvers are same (used for PCView()) */
2113f74950SBarry Smith   PetscInt   *l_lens;           /* lens of each block */
2213f74950SBarry Smith   PetscInt   *g_lens;
234b9ad928SBarry Smith   Mat        tp_mat,tp_pmat;    /* diagonal block of matrix for this processor */
24*5a7e1895SHong Zhang   PetscSubcomm psubcomm;        /* for multiple processors per block */
254b9ad928SBarry Smith } PC_BJacobi;
264b9ad928SBarry Smith 
274b9ad928SBarry Smith /*
284b9ad928SBarry Smith        This data is specific for certain implementations
294b9ad928SBarry Smith */
304b9ad928SBarry Smith 
314b9ad928SBarry Smith /*  This is for multiple blocks per processor */
324b9ad928SBarry Smith typedef struct {
334b9ad928SBarry Smith   Vec              *x,*y;             /* work vectors for solves on each block */
3413f74950SBarry Smith   PetscInt         *starts;           /* starting point of each block */
354b9ad928SBarry Smith   Mat              *mat,*pmat;        /* submatrices for each block */
364b9ad928SBarry Smith   IS               *is;               /* for gathering the submatrices */
374b9ad928SBarry Smith } PC_BJacobi_Multiblock;
384b9ad928SBarry Smith 
394b9ad928SBarry Smith /*  This is for a single block per processor */
404b9ad928SBarry Smith typedef struct {
414b9ad928SBarry Smith   Vec  x,y;
424b9ad928SBarry Smith } PC_BJacobi_Singleblock;
434b9ad928SBarry Smith 
44*5a7e1895SHong Zhang /*  This is for multiple processors per block */
45*5a7e1895SHong Zhang typedef struct {
46*5a7e1895SHong Zhang   KSP          ksp;                /* ksp used on each subcommunicator */
47*5a7e1895SHong Zhang   PC           pc;                 /* preconditioner used on each subcommunicator */
48*5a7e1895SHong Zhang   Vec          xsub,ysub;          /* vectors of a subcommunicator to hold parallel vectors of ((PetscObject)pc)->comm */
49*5a7e1895SHong Zhang   Mat          submats;            /* matrix and optional preconditioner matrix belong to a subcommunicator */
50*5a7e1895SHong Zhang   PetscSubcomm psubcomm;           /* only psubcomm->comm is needed! */
51*5a7e1895SHong Zhang   PetscInt     nsubcomm;           /* total number of blocks, which is same as number of subcommunicators */
52*5a7e1895SHong Zhang } PC_BJacobi_Multiproc;
534b9ad928SBarry Smith #endif
544b9ad928SBarry Smith 
554b9ad928SBarry Smith 
56