xref: /petsc/src/ksp/pc/impls/bjacobi/bjacobi.h (revision 226f8a8a5081bc6ad7227cd631662400f0d6e2a0)
1a4963045SJacob Faibussowitsch #pragma once
24b9ad928SBarry Smith /*
34b9ad928SBarry Smith     Private data for block Jacobi and block Gauss-Seidel preconditioner.
44b9ad928SBarry Smith */
5af0996ceSBarry Smith #include <petsc/private/pcimpl.h>
64b9ad928SBarry Smith 
74b9ad928SBarry Smith /*
84b9ad928SBarry Smith        This data is general for all implementations
94b9ad928SBarry Smith */
104b9ad928SBarry Smith typedef struct {
115a7e1895SHong Zhang   PetscInt     n;              /* number of global blocks */
125a7e1895SHong Zhang   PetscInt     n_local;        /* number of blocks in this subcommunicator or in this process */
1313f74950SBarry Smith   PetscInt     first_local;    /* number of first block on processor */
14*7addb90fSBarry Smith   PetscBool    use_true_local; /* use block from true matrix, not the matrix used to construct the preconditioner for local MatMult() */
152bbb50dfSHong Zhang   KSP         *ksp;            /* KSP contexts for blocks or for subcommunicator */
164b9ad928SBarry Smith   void        *data;           /* implementation-specific data */
1713f74950SBarry Smith   PetscInt    *l_lens;         /* lens of each block */
1813f74950SBarry Smith   PetscInt    *g_lens;
195a7e1895SHong Zhang   PetscSubcomm psubcomm; /* for multiple processors per block */
204b9ad928SBarry Smith } PC_BJacobi;
214b9ad928SBarry Smith 
224b9ad928SBarry Smith /*
234b9ad928SBarry Smith        This data is specific for certain implementations
244b9ad928SBarry Smith */
254b9ad928SBarry Smith 
264b9ad928SBarry Smith /*  This is for multiple blocks per processor */
274b9ad928SBarry Smith typedef struct {
284b9ad928SBarry Smith   Vec      *x, *y;      /* work vectors for solves on each block */
2913f74950SBarry Smith   PetscInt *starts;     /* starting point of each block */
304b9ad928SBarry Smith   Mat      *mat, *pmat; /* submatrices for each block */
314b9ad928SBarry Smith   IS       *is;         /* for gathering the submatrices */
324b9ad928SBarry Smith } PC_BJacobi_Multiblock;
334b9ad928SBarry Smith 
344b9ad928SBarry Smith /*  This is for a single block per processor */
354b9ad928SBarry Smith typedef struct {
364b9ad928SBarry Smith   Vec x, y;
374b9ad928SBarry Smith } PC_BJacobi_Singleblock;
384b9ad928SBarry Smith 
395a7e1895SHong Zhang /*  This is for multiple processors per block */
405a7e1895SHong Zhang typedef struct {
415a7e1895SHong Zhang   PC           pc;         /* preconditioner used on each subcommunicator */
42ce94432eSBarry Smith   Vec          xsub, ysub; /* vectors of a subcommunicator to hold parallel vectors of PetscObjectComm((PetscObject)pc) */
43*7addb90fSBarry Smith   Mat          submats;    /* the matrices belong to a subcommunicator */
44d6037b41SHong Zhang   PetscSubcomm psubcomm;
455a7e1895SHong Zhang } PC_BJacobi_Multiproc;
46