xref: /petsc/src/mat/impls/aij/mpi/mpiaij.h (revision b41ce5d507ea9a58bfa83cf403107a702e77a67d)
1 
2 #if !defined(__MPIAIJ_H)
3 #define __MPIAIJ_H
4 
5 #include <../src/mat/impls/aij/seq/aij.h>
6 
7 typedef struct { /* used by MatCreateMPIAIJSumSeqAIJ for reusing the merged matrix */
8   PetscLayout rowmap;
9   PetscInt    **buf_ri,**buf_rj;
10   PetscMPIInt *len_s,*len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
11   PetscMPIInt nsend,nrecv;
12   PetscInt    *bi,*bj;    /* i and j array of the local portion of mpi C (matrix product) - rename to ci, cj! */
13   PetscInt    *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
14   PetscErrorCode (*destroy)(Mat);
15   PetscErrorCode (*duplicate)(Mat,MatDuplicateOption,Mat*);
16 } Mat_Merge_SeqsToMPI;
17 
18 typedef struct { /* used by MatPtAP_MPIAIJ_MPIAIJ() and MatMatMult_MPIAIJ_MPIAIJ() */
19   PetscInt    *startsj_s,*startsj_r;    /* used by MatGetBrowsOfAoCols_MPIAIJ */
20   PetscScalar *bufa;                    /* used by MatGetBrowsOfAoCols_MPIAIJ */
21   Mat         P_loc,P_oth;     /* partial B_seq -- intend to replace B_seq */
22   PetscInt    *api,*apj;       /* symbolic i and j arrays of the local product A_loc*B_seq */
23   PetscScalar *apv;
24   MatReuse    reuse;           /* flag to skip MatGetBrowsOfAoCols_MPIAIJ() and MatMPIAIJGetLocalMat() in 1st call of MatPtAPNumeric_MPIAIJ_MPIAIJ() */
25   PetscScalar *apa;            /* tmp array for store a row of A*P used in MatMatMult() */
26   Mat         A_loc;           /* used by MatTransposeMatMult(), contains api and apj */
27   Mat         Pt;              /* used by MatTransposeMatMult(), Pt = P^T */
28   PetscBool   scalable;        /* flag determines scalable or non-scalable implementation */
29   Mat         Rd,Ro,AP_loc,C_loc,C_oth;
30   PetscInt    algType;         /* implementation algorithm */
31 
32   Mat_Merge_SeqsToMPI *merge;
33   PetscErrorCode (*destroy)(Mat);
34   PetscErrorCode (*duplicate)(Mat,MatDuplicateOption,Mat*);
35   PetscErrorCode (*view)(Mat,PetscViewer);
36 } Mat_PtAPMPI;
37 
38 typedef struct {
39   Mat A,B;                             /* local submatrices: A (diag part),
40                                            B (off-diag part) */
41   PetscMPIInt size;                     /* size of communicator */
42   PetscMPIInt rank;                     /* rank of proc in communicator */
43 
44   /* The following variables are used for matrix assembly */
45   PetscBool   donotstash;               /* PETSC_TRUE if off processor entries dropped */
46   MPI_Request *send_waits;              /* array of send requests */
47   MPI_Request *recv_waits;              /* array of receive requests */
48   PetscInt    nsends,nrecvs;           /* numbers of sends and receives */
49   PetscScalar *svalues,*rvalues;       /* sending and receiving data */
50   PetscInt    rmax;                     /* maximum message length */
51 #if defined(PETSC_USE_CTABLE)
52   PetscTable colmap;
53 #else
54   PetscInt *colmap;                     /* local col number of off-diag col */
55 #endif
56   PetscInt *garray;                     /* global index of all off-processor columns */
57 
58   /* The following variables are used for matrix-vector products */
59   Vec        lvec;                 /* local vector */
60   Vec        diag;
61   VecScatter Mvctx;                /* scatter context for vector */
62   PetscBool  roworiented;          /* if true, row-oriented input, default true */
63 
64   /* The following variables are for MatGetRow() */
65   PetscInt    *rowindices;         /* column indices for row */
66   PetscScalar *rowvalues;          /* nonzero values in row */
67   PetscBool   getrowactive;        /* indicates MatGetRow(), not restored */
68 
69   /* Used by MatDistribute_MPIAIJ() to allow reuse of previous matrix allocation  and nonzero pattern */
70   PetscInt *ld;                    /* number of entries per row left of diagona block */
71 
72   /* Used by MatMatMult() and MatPtAP() */
73   Mat_PtAPMPI *ptap;
74 
75   /* used by MatMatMatMult() */
76   Mat_MatMatMatMult *matmatmatmult;
77 
78   /* Used by MPICUSP and MPICUSPARSE classes */
79   void * spptr;
80 
81 } Mat_MPIAIJ;
82 
83 PETSC_EXTERN PetscErrorCode MatCreate_MPIAIJ(Mat);
84 
85 PETSC_INTERN PetscErrorCode MatAssemblyEnd_MPIAIJ(Mat,MatAssemblyType);
86 
87 PETSC_INTERN PetscErrorCode MatSetUpMultiply_MPIAIJ(Mat);
88 PETSC_INTERN PetscErrorCode MatDisAssemble_MPIAIJ(Mat);
89 PETSC_INTERN PetscErrorCode MatDuplicate_MPIAIJ(Mat,MatDuplicateOption,Mat*);
90 PETSC_INTERN PetscErrorCode MatIncreaseOverlap_MPIAIJ(Mat,PetscInt,IS [],PetscInt);
91 PETSC_INTERN PetscErrorCode MatIncreaseOverlap_MPIAIJ_Scalable(Mat,PetscInt,IS [],PetscInt);
92 PETSC_INTERN PetscErrorCode MatFDColoringCreate_MPIXAIJ(Mat,ISColoring,MatFDColoring);
93 PETSC_INTERN PetscErrorCode MatFDColoringSetUp_MPIXAIJ(Mat,ISColoring,MatFDColoring);
94 PETSC_INTERN PetscErrorCode MatCreateSubMatrices_MPIAIJ (Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
95 PETSC_INTERN PetscErrorCode MatCreateSubMatricesMPI_MPIAIJ (Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
96 PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_All(Mat,MatCreateSubMatrixOption,MatReuse,Mat *[]);
97 
98 
99 PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ(Mat,IS,IS,MatReuse,Mat*);
100 PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_nonscalable(Mat,IS,IS,PetscInt,MatReuse,Mat*);
101 PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_SameRowDist(Mat,IS,IS,MatReuse,Mat*);
102 PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_SameRowColDist(Mat,IS,IS,MatReuse,Mat*);
103 PETSC_INTERN PetscErrorCode MatGetMultiProcBlock_MPIAIJ(Mat,MPI_Comm,MatReuse,Mat*);
104 
105 PETSC_INTERN PetscErrorCode MatLoad_MPIAIJ(Mat,PetscViewer);
106 PETSC_INTERN PetscErrorCode MatCreateColmap_MPIAIJ_Private(Mat);
107 PETSC_INTERN PetscErrorCode MatMatMult_MPIDense_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
108 PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
109 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,PetscReal,Mat*);
110 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*);
111 PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat);
112 PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,Mat);
113 
114 PETSC_INTERN PetscErrorCode MatMatMatMult_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
115 PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,PetscReal,Mat*);
116 PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,Mat);
117 
118 PETSC_INTERN PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
119 PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*);
120 PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat);
121 
122 PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(Mat,Mat,PetscReal,Mat*);
123 PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap(Mat,Mat,Mat);
124 
125 PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat);
126 PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ(Mat);
127 
128 PETSC_INTERN PetscErrorCode MatRARt_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
129 
130 PETSC_INTERN PetscErrorCode MatGetBrowsOfAoCols_MPIAIJ(Mat,Mat,MatReuse,PetscInt**,PetscInt**,MatScalar**,Mat*);
131 PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar [],InsertMode);
132 PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat);
133 PETSC_INTERN PetscErrorCode PetscContainerDestroy_Mat_MatMatMultMPI(void*);
134 PETSC_INTERN PetscErrorCode MatSetOption_MPIAIJ(Mat,MatOption,PetscBool);
135 
136 PETSC_INTERN PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
137 PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,PetscReal,Mat*);
138 PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*);
139 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat);
140 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,Mat);
141 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat,Mat,Mat);
142 PETSC_INTERN PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*);
143 PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat*);
144 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIDense(Mat,Mat,Mat);
145 PETSC_INTERN PetscErrorCode MatGetSeqNonzeroStructure_MPIAIJ(Mat,Mat*);
146 
147 PETSC_INTERN PetscErrorCode MatSetFromOptions_MPIAIJ(PetscOptionItems*,Mat);
148 PETSC_INTERN PetscErrorCode MatMPIAIJSetPreallocation_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
149 
150 #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128) && !defined(PETSC_USE_REAL___FP16)
151 PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_MPIAIJ_TFS(Mat,IS,IS,const MatFactorInfo*,Mat*);
152 #endif
153 PETSC_INTERN PetscErrorCode MatSolve_MPIAIJ(Mat,Vec,Vec);
154 PETSC_INTERN PetscErrorCode MatILUFactor_MPIAIJ(Mat,IS,IS,const MatFactorInfo*);
155 
156 PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_MPIX_private(PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,PetscInt*);
157 
158 extern PetscErrorCode MatGetDiagonalBlock_MPIAIJ(Mat,Mat*);
159 extern PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat,Vec);
160 
161 PETSC_INTERN PetscErrorCode MatGetSeqMats_MPIAIJ(Mat,Mat*,Mat*);
162 PETSC_INTERN PetscErrorCode MatSetSeqMats_MPIAIJ(Mat,IS,IS,IS,MatStructure,Mat,Mat);
163 
164 /* compute apa = A[i,:]*P = Ad[i,:]*P_loc + Ao*[i,:]*P_oth using sparse axpy */
165 #define AProw_scalable(i,ad,ao,p_loc,p_oth,api,apj,apa) \
166 {\
167   PetscInt    _anz,_pnz,_j,_k,*_ai,*_aj,_row,*_pi,*_pj,_nextp,*_apJ;      \
168   PetscScalar *_aa,_valtmp,*_pa;                             \
169   _apJ = apj + api[i];\
170   /* diagonal portion of A */\
171   _ai  = ad->i;\
172   _anz = _ai[i+1] - _ai[i];\
173   _aj  = ad->j + _ai[i];\
174   _aa  = ad->a + _ai[i];\
175   for (_j=0; _j<_anz; _j++) {\
176     _row = _aj[_j]; \
177     _pi  = p_loc->i;                                 \
178     _pnz = _pi[_row+1] - _pi[_row];         \
179     _pj  = p_loc->j + _pi[_row];                 \
180     _pa  = p_loc->a + _pi[_row];                 \
181     /* perform sparse axpy */                    \
182     _valtmp = _aa[_j];                           \
183     _nextp  = 0; \
184     for (_k=0; _nextp<_pnz; _k++) {                    \
185       if (_apJ[_k] == _pj[_nextp]) { /* column of AP == column of P */   \
186         apa[_k] += _valtmp*_pa[_nextp++];                                \
187       } \
188     }                                           \
189     (void)PetscLogFlops(2.0*_pnz);              \
190   }                                             \
191   /* off-diagonal portion of A */               \
192   _ai  = ao->i;\
193   _anz = _ai[i+1] - _ai[i];                     \
194   _aj  = ao->j + _ai[i];                         \
195   _aa  = ao->a + _ai[i];                         \
196   for (_j=0; _j<_anz; _j++) {                      \
197     _row = _aj[_j];    \
198     _pi  = p_oth->i;                         \
199     _pnz = _pi[_row+1] - _pi[_row];          \
200     _pj  = p_oth->j + _pi[_row];                  \
201     _pa  = p_oth->a + _pi[_row];                  \
202     /* perform sparse axpy */                     \
203     _valtmp = _aa[_j];                             \
204     _nextp  = 0; \
205     for (_k=0; _nextp<_pnz; _k++) {                     \
206       if (_apJ[_k] == _pj[_nextp]) { /* column of AP == column of P */\
207         apa[_k] += _valtmp*_pa[_nextp++];                       \
208       }                                                     \
209     }                                            \
210     (void)PetscLogFlops(2.0*_pnz);               \
211   } \
212 }
213 
214 #define AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa) \
215 {\
216   PetscInt    _anz,_pnz,_j,_k,*_ai,*_aj,_row,*_pi,*_pj;      \
217   PetscScalar *_aa,_valtmp,*_pa;                             \
218   /* diagonal portion of A */\
219   _ai  = ad->i;\
220   _anz = _ai[i+1] - _ai[i];\
221   _aj  = ad->j + _ai[i];\
222   _aa  = ad->a + _ai[i];\
223   for (_j=0; _j<_anz; _j++) {\
224     _row = _aj[_j]; \
225     _pi  = p_loc->i;                                 \
226     _pnz = _pi[_row+1] - _pi[_row];         \
227     _pj  = p_loc->j + _pi[_row];                 \
228     _pa  = p_loc->a + _pi[_row];                 \
229     /* perform dense axpy */                    \
230     _valtmp = _aa[_j];                           \
231     for (_k=0; _k<_pnz; _k++) {                    \
232       apa[_pj[_k]] += _valtmp*_pa[_k];               \
233     }                                           \
234     (void)PetscLogFlops(2.0*_pnz);              \
235   }                                             \
236   /* off-diagonal portion of A */               \
237   _ai  = ao->i;\
238   _anz = _ai[i+1] - _ai[i];                     \
239   _aj  = ao->j + _ai[i];                         \
240   _aa  = ao->a + _ai[i];                         \
241   for (_j=0; _j<_anz; _j++) {                      \
242     _row = _aj[_j];    \
243     _pi  = p_oth->i;                         \
244     _pnz = _pi[_row+1] - _pi[_row];          \
245     _pj  = p_oth->j + _pi[_row];                  \
246     _pa  = p_oth->a + _pi[_row];                  \
247     /* perform dense axpy */                     \
248     _valtmp = _aa[_j];                             \
249     for (_k=0; _k<_pnz; _k++) {                     \
250       apa[_pj[_k]] += _valtmp*_pa[_k];                \
251     }                                            \
252     (void)PetscLogFlops(2.0*_pnz);               \
253   } \
254 }
255 
256 #endif
257