xref: /petsc/src/mat/impls/baij/mpi/mpibaij.h (revision bebe2cf65d55febe21a5af8db2bd2e168caaa2e7)
1 
2 #if !defined(__MPIBAIJ_H)
3 #define __MPIBAIJ_H
4 #include <../src/mat/impls/baij/seq/baij.h>
5 #include <../src/mat/impls/aij/mpi/mpiaij.h>
6 #include <petscctable.h>
7 
8 #if defined(PETSC_USE_CTABLE)
9 #define PETSCTABLE PetscTable
10 #else
11 #define PETSCTABLE PetscInt*
12 #endif
13 
14 #define MPIBAIJHEADER \
15   PetscInt    *rangebs;                /* rmap->range/bs */                                    \
16   PetscInt    rstartbs,rendbs,cstartbs,cendbs;    /* map values / bs  */                       \
17   Mat         A,B;                     /* local submatrices: A (diag part), B (off-diag part) */ \
18   PetscMPIInt size;                    /* size of communicator */                              \
19   PetscMPIInt rank;                    /* rank of proc in communicator */                      \
20   PetscInt    bs2;                     /* block size, bs2 = bs*bs */                           \
21   PetscInt    Mbs,Nbs;                 /* number block rows/cols in matrix; M/bs, N/bs */      \
22   PetscInt    mbs,nbs;                 /* number block rows/cols on processor; m/bs, n/bs */   \
23                                                                                                \
24   /* The following variables are used for matrix assembly */                                   \
25                                                                                                \
26   PetscBool   donotstash;               /* if 1, off processor entries dropped */              \
27   MPI_Request *send_waits;              /* array of send requests */                           \
28   MPI_Request *recv_waits;              /* array of receive requests */                        \
29   PetscInt    nsends,nrecvs;           /* numbers of sends and receives */                     \
30   MatScalar   *svalues,*rvalues;       /* sending and receiving data */                        \
31   PetscInt    rmax;                     /* maximum message length */                           \
32   PETSCTABLE  colmap;                   /* local col number of off-diag col */                 \
33                                                                                                \
34   PetscInt *garray;                    /* work array */                                        \
35                                                                                                \
36   /* The following variable is used by blocked matrix assembly */                              \
37   MatScalar *barray;                    /* Block array of size bs2 */                          \
38                                                                                                \
39   /* The following variables are used for matrix-vector products */                            \
40                                                                                                \
41   Vec        lvec;                 /* local vector */                                          \
42   VecScatter Mvctx;                /* scatter context for vector */                            \
43   PetscBool  roworiented;          /* if true, row-oriented input, default true */             \
44                                                                                                \
45   /* The following variables are for MatGetRow() */                                            \
46                                                                                                \
47   PetscInt    *rowindices;         /* column indices for row */                                \
48   PetscScalar *rowvalues;          /* nonzero values in row */                                 \
49   PetscBool   getrowactive;        /* indicates MatGetRow(), not restored */                   \
50                                                                                                \
51   /* Some variables to make MatSetValues and others more efficient */                          \
52   PetscInt  rstart_bs,rend_bs;                                                                 \
53   PetscInt  cstart_bs,cend_bs;                                                                 \
54   PetscInt  *ht;                          /* Hash table to speed up matrix assembly */         \
55   MatScalar **hd;                         /* Hash table data */                                \
56   PetscInt  ht_size;                                                                           \
57   PetscInt  ht_total_ct,ht_insert_ct;     /* Hash table statistics */                          \
58   PetscBool ht_flag;                      /* Flag to indicate if hash tables are used */       \
59   double    ht_fact;                      /* Factor to determine the HT size */                \
60                                                                                                \
61   PetscInt  setvalueslen;       /* only used for single precision computations */              \
62   MatScalar *setvaluescopy;     /* area double precision values in MatSetValuesXXX() are copied*/ \
63                                 /* before calling MatSetValuesXXX_MPIBAIJ_MatScalar() */       \
64   PetscBool ijonly             /* used in  MatGetSubMatrices_MPIBAIJ_local() for getting ij structure only */
65 
66 typedef struct {
67   MPIBAIJHEADER;
68 } Mat_MPIBAIJ;
69 
70 PETSC_INTERN PetscErrorCode MatLoad_MPIBAIJ(Mat,PetscViewer);
71 PETSC_INTERN PetscErrorCode MatCreateColmap_MPIBAIJ_Private(Mat);
72 PETSC_INTERN PetscErrorCode MatGetSubMatrices_MPIBAIJ(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat*[]);
73 PETSC_INTERN PetscErrorCode MatGetSubMatrices_MPIBAIJ_local(Mat,PetscInt,const IS[],const IS[],MatReuse,PetscBool*,PetscBool*,Mat*);
74 PETSC_INTERN PetscErrorCode MatGetSubMatrix_MPIBAIJ_Private(Mat,IS,IS,PetscInt,MatReuse,Mat*);
75 PETSC_INTERN PetscErrorCode MatGetMultiProcBlock_MPIBAIJ(Mat,MPI_Comm,MatReuse,Mat*);
76 PETSC_INTERN PetscErrorCode MatIncreaseOverlap_MPIBAIJ(Mat,PetscInt,IS[],PetscInt);
77 PETSC_INTERN PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Once(Mat,PetscInt,IS*);
78 PETSC_INTERN PetscErrorCode MatFDColoringCreate_MPIXAIJ(Mat,ISColoring,MatFDColoring);
79 PETSC_INTERN PetscErrorCode MatFDColoringSetUp_MPIXAIJ(Mat,ISColoring,MatFDColoring);
80 PETSC_INTERN PetscErrorCode MatFDColoringApply_BAIJ(Mat,MatFDColoring,Vec,void*);
81 
82 PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_MPIBAIJ(Mat,const PetscInt *,Mat,const PetscInt*,PetscInt*);
83 #endif
84