xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision 66bfb163c73c2c7e13b37dc600d5cce69ff4f6af)
1 
2 /*
3   Defines matrix-matrix product routines for pairs of MPIAIJ matrices
4           C = A * B
5 */
6 #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/
7 #include <../src/mat/utils/freespace.h>
8 #include <../src/mat/impls/aij/mpi/mpiaij.h>
9 #include <petscbt.h>
10 #include <../src/mat/impls/dense/mpi/mpidense.h>
11 
12 #undef __FUNCT__
13 #define __FUNCT__ "MatMatMult_MPIAIJ_MPIAIJ"
14 PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C)
15 {
16   PetscErrorCode ierr;
17 
18   PetscFunctionBegin;
19   if (scall == MAT_INITIAL_MATRIX){
20     ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);/* numeric product is computed as well */
21   } else if (scall == MAT_REUSE_MATRIX){
22     ierr = MatMatMultNumeric_MPIAIJ_MPIAIJ(A,B,*C);CHKERRQ(ierr);
23   } else {
24     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall);
25   }
26   PetscFunctionReturn(0);
27 }
28 
29 #undef __FUNCT__
30 #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatMultMPI"
31 PetscErrorCode PetscContainerDestroy_Mat_MatMatMultMPI(void *ptr)
32 {
33   PetscErrorCode       ierr;
34   Mat_MatMatMultMPI    *mult=(Mat_MatMatMultMPI*)ptr;
35 
36   PetscFunctionBegin;
37   ierr = PetscFree2(mult->startsj,mult->startsj_r);CHKERRQ(ierr);
38   ierr = PetscFree(mult->bufa);CHKERRQ(ierr);
39   ierr = ISDestroy(&mult->isrowa);CHKERRQ(ierr);
40   ierr = ISDestroy(&mult->isrowb);CHKERRQ(ierr);
41   ierr = ISDestroy(&mult->iscolb);CHKERRQ(ierr);
42   ierr = MatDestroy(&mult->C_seq);CHKERRQ(ierr);
43   ierr = MatDestroy(&mult->A_loc);CHKERRQ(ierr);
44   ierr = MatDestroy(&mult->B_seq);CHKERRQ(ierr);
45   ierr = MatDestroy(&mult->B_loc);CHKERRQ(ierr);
46   ierr = MatDestroy(&mult->B_oth);CHKERRQ(ierr);
47   ierr = PetscFree(mult->abi);CHKERRQ(ierr);
48   ierr = PetscFree(mult->abj);CHKERRQ(ierr);
49   ierr = PetscFree(mult);CHKERRQ(ierr);
50   PetscFunctionReturn(0);
51 }
52 
53 extern PetscErrorCode MatDestroy_AIJ(Mat);
54 #undef __FUNCT__
55 #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult"
56 PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
57 {
58   PetscErrorCode     ierr;
59   PetscContainer     container;
60   Mat_MatMatMultMPI  *mult=PETSC_NULL;
61 
62   PetscFunctionBegin;
63   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
64   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
65   ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
66   A->ops->destroy   = mult->destroy;
67   A->ops->duplicate = mult->duplicate;
68   if (A->ops->destroy) {
69     ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
70   }
71   ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatMultMPI",0);CHKERRQ(ierr);
72   PetscFunctionReturn(0);
73 }
74 
75 #undef __FUNCT__
76 #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult"
77 PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M)
78 {
79   PetscErrorCode     ierr;
80   Mat_MatMatMultMPI  *mult;
81   PetscContainer     container;
82 
83   PetscFunctionBegin;
84   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
85   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
86   ierr  = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
87   /* Note: the container is not duplicated, because it requires deep copying of
88      several large data sets (see PetscContainerDestroy_Mat_MatMatMultMPI()).
89      These data sets are only used for repeated calling of MatMatMultNumeric().
90      *M is unlikely being used in this way. Thus we create *M with pure mpiaij format */
91   ierr = (*mult->duplicate)(A,op,M);CHKERRQ(ierr);
92   (*M)->ops->destroy   = mult->destroy;   /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's container! */
93   (*M)->ops->duplicate = mult->duplicate; /* = MatDuplicate_MPIAIJ */
94   PetscFunctionReturn(0);
95 }
96 
97 #undef __FUNCT__
98 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
99 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
100 {
101   PetscErrorCode     ierr;
102   PetscInt           start,end;
103   Mat_MatMatMultMPI  *mult;
104   PetscContainer     container;
105   Mat                *aloc;
106 
107   PetscFunctionBegin;
108   if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend){
109     SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, (%D, %D) != (%D,%D)",A->cmap->rstart,A->cmap->rend,B->rmap->rstart,B->rmap->rend);
110   }
111   ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr);
112 
113   /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */
114   ierr = MatGetBrowsOfAcols(A,B,MAT_INITIAL_MATRIX,&mult->isrowb,&mult->iscolb,&mult->B_seq);CHKERRQ(ierr);
115 
116   /*  create a seq matrix A_seq = submatrix of A by taking all local rows of A */
117   start = A->rmap->rstart; end = A->rmap->rend;
118   ierr = ISCreateStride(PETSC_COMM_SELF,end-start,start,1,&mult->isrowa);CHKERRQ(ierr);
119   ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_INITIAL_MATRIX,&aloc);CHKERRQ(ierr);
120   mult->A_loc = aloc[0];
121   ierr = PetscFree(aloc);CHKERRQ(ierr);
122 
123   /* compute C_seq = A_seq * B_seq */
124   ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,MAT_INITIAL_MATRIX,fill,&mult->C_seq);CHKERRQ(ierr);
125 
126   /* create mpi matrix C by concatinating C_seq */
127   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); /* prevent C_seq being destroyed by MatMerge() */
128   ierr = MatMerge(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,MAT_INITIAL_MATRIX,C);CHKERRQ(ierr);
129 
130   /* attach the supporting struct to C for reuse of symbolic C */
131   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
132   ierr = PetscContainerSetPointer(container,mult);CHKERRQ(ierr);
133   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatMultMPI);CHKERRQ(ierr);
134   ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr);
135   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
136   mult->destroy   = (*C)->ops->destroy;
137   mult->duplicate = (*C)->ops->duplicate;
138   (*C)->ops->destroy   = MatDestroy_MPIAIJ_MatMatMult;
139   (*C)->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult;
140 
141   PetscFunctionReturn(0);
142 }
143 
144 /* This routine is called ONLY in the case of reusing previously computed symbolic C */
145 #undef __FUNCT__
146 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ"
147 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat B,Mat C)
148 {
149   PetscErrorCode     ierr;
150   Mat                *seq;
151   Mat_MatMatMultMPI  *mult;
152   PetscContainer     container;
153 
154   PetscFunctionBegin;
155   ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
156   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
157   ierr  = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
158   seq = &mult->B_seq;
159   ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr);
160   mult->B_seq = *seq;
161 
162   seq = &mult->A_loc;
163   ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr);
164   mult->A_loc = *seq;
165 
166   ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,MAT_REUSE_MATRIX,0.0,&mult->C_seq);CHKERRQ(ierr);
167 
168   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr);
169   ierr = MatMerge(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr);
170   PetscFunctionReturn(0);
171 }
172 
173 #undef __FUNCT__
174 #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense"
175 PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
176 {
177   PetscErrorCode ierr;
178 
179   PetscFunctionBegin;
180   if (scall == MAT_INITIAL_MATRIX){
181     ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr);
182   }
183   ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr);
184   PetscFunctionReturn(0);
185 }
186 
187 typedef struct {
188   Mat         workB;
189   PetscScalar *rvalues,*svalues;
190   MPI_Request *rwaits,*swaits;
191 } MPIAIJ_MPIDense;
192 
193 #undef __FUNCT__
194 #define __FUNCT__ "MPIAIJ_MPIDenseDestroy"
195 PetscErrorCode MPIAIJ_MPIDenseDestroy(void *ctx)
196 {
197   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx;
198   PetscErrorCode  ierr;
199 
200   PetscFunctionBegin;
201   ierr = MatDestroy(&contents->workB);CHKERRQ(ierr);
202   ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr);
203   ierr = PetscFree(contents);CHKERRQ(ierr);
204   PetscFunctionReturn(0);
205 }
206 
207 #undef __FUNCT__
208 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense"
209 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
210 {
211   PetscErrorCode         ierr;
212   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
213   PetscInt               nz = aij->B->cmap->n;
214   PetscContainer         container;
215   MPIAIJ_MPIDense        *contents;
216   VecScatter             ctx = aij->Mvctx;
217   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
218   VecScatter_MPI_General *to   = ( VecScatter_MPI_General*) ctx->todata;
219   PetscInt               m=A->rmap->n,n=B->cmap->n;
220 
221   PetscFunctionBegin;
222   ierr = MatCreate(((PetscObject)B)->comm,C);CHKERRQ(ierr);
223   ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr);
224   ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr);
225   ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
226   ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
227 
228   ierr = PetscNew(MPIAIJ_MPIDense,&contents);CHKERRQ(ierr);
229   /* Create work matrix used to store off processor rows of B needed for local product */
230   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,PETSC_NULL,&contents->workB);CHKERRQ(ierr);
231   /* Create work arrays needed */
232   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],PetscScalar,&contents->rvalues,
233                       B->cmap->N*to->starts[to->n],PetscScalar,&contents->svalues,
234                       from->n,MPI_Request,&contents->rwaits,
235                       to->n,MPI_Request,&contents->swaits);CHKERRQ(ierr);
236 
237   ierr = PetscContainerCreate(((PetscObject)A)->comm,&container);CHKERRQ(ierr);
238   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
239   ierr = PetscContainerSetUserDestroy(container,MPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
240   ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr);
241   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
242   PetscFunctionReturn(0);
243 }
244 
245 #undef __FUNCT__
246 #define __FUNCT__ "MatMPIDenseScatter"
247 /*
248     Performs an efficient scatter on the rows of B needed by this process; this is
249     a modification of the VecScatterBegin_() routines.
250 */
251 PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB)
252 {
253   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)A->data;
254   PetscErrorCode         ierr;
255   PetscScalar            *b,*w,*svalues,*rvalues;
256   VecScatter             ctx = aij->Mvctx;
257   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
258   VecScatter_MPI_General *to   = ( VecScatter_MPI_General*) ctx->todata;
259   PetscInt               i,j,k;
260   PetscInt               *sindices,*sstarts,*rindices,*rstarts;
261   PetscMPIInt            *sprocs,*rprocs,nrecvs;
262   MPI_Request            *swaits,*rwaits;
263   MPI_Comm               comm = ((PetscObject)A)->comm;
264   PetscMPIInt            tag = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n;
265   MPI_Status             status;
266   MPIAIJ_MPIDense        *contents;
267   PetscContainer         container;
268   Mat                    workB;
269 
270   PetscFunctionBegin;
271   ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr);
272   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
273   ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr);
274 
275   workB = *outworkB = contents->workB;
276   if (nrows != workB->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number of rows of workB %D not equal to columns of aij->B %D",nrows,workB->cmap->n);
277   sindices  = to->indices;
278   sstarts   = to->starts;
279   sprocs    = to->procs;
280   swaits    = contents->swaits;
281   svalues   = contents->svalues;
282 
283   rindices  = from->indices;
284   rstarts   = from->starts;
285   rprocs    = from->procs;
286   rwaits    = contents->rwaits;
287   rvalues   = contents->rvalues;
288 
289   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
290   ierr = MatGetArray(workB,&w);CHKERRQ(ierr);
291 
292   for (i=0; i<from->n; i++) {
293     ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr);
294   }
295 
296   for (i=0; i<to->n; i++) {
297     /* pack a message at a time */
298     CHKMEMQ;
299     for (j=0; j<sstarts[i+1]-sstarts[i]; j++){
300       for (k=0; k<ncols; k++) {
301         svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k];
302       }
303     }
304     CHKMEMQ;
305     ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr);
306   }
307 
308   nrecvs = from->n;
309   while (nrecvs) {
310     ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr);
311     nrecvs--;
312     /* unpack a message at a time */
313     CHKMEMQ;
314     for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++){
315       for (k=0; k<ncols; k++) {
316         w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k];
317       }
318     }
319     CHKMEMQ;
320   }
321   if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);}
322 
323   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
324   ierr = MatRestoreArray(workB,&w);CHKERRQ(ierr);
325   ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
326   ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
327   PetscFunctionReturn(0);
328 }
329 extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat);
330 
331 #undef __FUNCT__
332 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense"
333 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
334 {
335   PetscErrorCode       ierr;
336   Mat_MPIAIJ           *aij = (Mat_MPIAIJ*)A->data;
337   Mat_MPIDense         *bdense = (Mat_MPIDense*)B->data;
338   Mat_MPIDense         *cdense = (Mat_MPIDense*)C->data;
339   Mat                  workB;
340 
341   PetscFunctionBegin;
342 
343   /* diagonal block of A times all local rows of B*/
344   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr);
345 
346   /* get off processor parts of B needed to complete the product */
347   ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr);
348 
349   /* off-diagonal block of A times nonlocal rows of B */
350   ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr);
351   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
352   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
353   PetscFunctionReturn(0);
354 }
355 
356