xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision fcb1c9af1a39adbb75e1011fa083d9d197867fec)
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 
106   PetscFunctionBegin;
107   if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend){
108     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);
109   }
110   ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr);
111 
112   /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */
113   ierr = MatGetBrowsOfAcols(A,B,MAT_INITIAL_MATRIX,&mult->isrowb,&mult->iscolb,&mult->brstart,&mult->B_seq);CHKERRQ(ierr);
114 
115   /*  create a seq matrix A_seq = submatrix of A by taking all local rows of A */
116   start = A->rmap->rstart; end = A->rmap->rend;
117   ierr = ISCreateStride(PETSC_COMM_SELF,end-start,start,1,&mult->isrowa);CHKERRQ(ierr);
118   ierr = MatMPIAIJGetLocalMatCondensed(A,MAT_INITIAL_MATRIX,&mult->isrowa,&mult->isrowb,&mult->A_loc);CHKERRQ(ierr);
119 
120   /* compute C_seq = A_seq * B_seq */
121   ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,MAT_INITIAL_MATRIX,fill,&mult->C_seq);CHKERRQ(ierr);
122 
123   /* create mpi matrix C by concatinating C_seq */
124   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); /* prevent C_seq being destroyed by MatMerge() */
125   ierr = MatMerge(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,MAT_INITIAL_MATRIX,C);CHKERRQ(ierr);
126 
127   /* attach the supporting struct to C for reuse of symbolic C */
128   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
129   ierr = PetscContainerSetPointer(container,mult);CHKERRQ(ierr);
130   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatMultMPI);CHKERRQ(ierr);
131   ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr);
132   mult->destroy   = (*C)->ops->destroy;
133   mult->duplicate = (*C)->ops->duplicate;
134   (*C)->ops->destroy   = MatDestroy_MPIAIJ_MatMatMult;
135   (*C)->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult;
136   PetscFunctionReturn(0);
137 }
138 
139 /* This routine is called ONLY in the case of reusing previously computed symbolic C */
140 #undef __FUNCT__
141 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ"
142 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat B,Mat C)
143 {
144   PetscErrorCode     ierr;
145   Mat                *seq;
146   Mat_MatMatMultMPI  *mult;
147   PetscContainer     container;
148 
149   PetscFunctionBegin;
150   ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
151   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
152   ierr  = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
153 
154   seq = &mult->B_seq;
155   ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr);
156   mult->B_seq = *seq;
157 
158   seq = &mult->A_loc;
159   ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr);
160   mult->A_loc = *seq;
161 
162   ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,MAT_REUSE_MATRIX,0.0,&mult->C_seq);CHKERRQ(ierr);
163 
164   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr);
165   ierr = MatMerge(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr);
166   PetscFunctionReturn(0);
167 }
168 
169 #undef __FUNCT__
170 #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense"
171 PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
172 {
173   PetscErrorCode ierr;
174 
175   PetscFunctionBegin;
176   if (scall == MAT_INITIAL_MATRIX){
177     ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr);
178   }
179   ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr);
180   PetscFunctionReturn(0);
181 }
182 
183 typedef struct {
184   Mat         workB;
185   PetscScalar *rvalues,*svalues;
186   MPI_Request *rwaits,*swaits;
187 } MPIAIJ_MPIDense;
188 
189 #undef __FUNCT__
190 #define __FUNCT__ "MPIAIJ_MPIDenseDestroy"
191 PetscErrorCode MPIAIJ_MPIDenseDestroy(void *ctx)
192 {
193   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx;
194   PetscErrorCode  ierr;
195 
196   PetscFunctionBegin;
197   ierr = MatDestroy(&contents->workB);CHKERRQ(ierr);
198   ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr);
199   ierr = PetscFree(contents);CHKERRQ(ierr);
200   PetscFunctionReturn(0);
201 }
202 
203 #undef __FUNCT__
204 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense"
205 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
206 {
207   PetscErrorCode         ierr;
208   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
209   PetscInt               nz = aij->B->cmap->n;
210   PetscContainer         container;
211   MPIAIJ_MPIDense        *contents;
212   VecScatter             ctx = aij->Mvctx;
213   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
214   VecScatter_MPI_General *to   = ( VecScatter_MPI_General*) ctx->todata;
215   PetscInt               m=A->rmap->n,n=B->cmap->n;
216 
217   PetscFunctionBegin;
218   ierr = MatCreate(((PetscObject)B)->comm,C);CHKERRQ(ierr);
219   ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr);
220   ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr);
221   ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
222   ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
223 
224   ierr = PetscNew(MPIAIJ_MPIDense,&contents);CHKERRQ(ierr);
225   /* Create work matrix used to store off processor rows of B needed for local product */
226   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,PETSC_NULL,&contents->workB);CHKERRQ(ierr);
227   /* Create work arrays needed */
228   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],PetscScalar,&contents->rvalues,
229                       B->cmap->N*to->starts[to->n],PetscScalar,&contents->svalues,
230                       from->n,MPI_Request,&contents->rwaits,
231                       to->n,MPI_Request,&contents->swaits);CHKERRQ(ierr);
232 
233   ierr = PetscContainerCreate(((PetscObject)A)->comm,&container);CHKERRQ(ierr);
234   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
235   ierr = PetscContainerSetUserDestroy(container,MPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
236   ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr);
237   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
238   PetscFunctionReturn(0);
239 }
240 
241 #undef __FUNCT__
242 #define __FUNCT__ "MatMPIDenseScatter"
243 /*
244     Performs an efficient scatter on the rows of B needed by this process; this is
245     a modification of the VecScatterBegin_() routines.
246 */
247 PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB)
248 {
249   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)A->data;
250   PetscErrorCode         ierr;
251   PetscScalar            *b,*w,*svalues,*rvalues;
252   VecScatter             ctx = aij->Mvctx;
253   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
254   VecScatter_MPI_General *to   = ( VecScatter_MPI_General*) ctx->todata;
255   PetscInt               i,j,k;
256   PetscInt               *sindices,*sstarts,*rindices,*rstarts;
257   PetscMPIInt            *sprocs,*rprocs,nrecvs;
258   MPI_Request            *swaits,*rwaits;
259   MPI_Comm               comm = ((PetscObject)A)->comm;
260   PetscMPIInt            tag = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n;
261   MPI_Status             status;
262   MPIAIJ_MPIDense        *contents;
263   PetscContainer         container;
264   Mat                    workB;
265 
266   PetscFunctionBegin;
267   ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr);
268   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
269   ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr);
270 
271   workB = *outworkB = contents->workB;
272   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);
273   sindices  = to->indices;
274   sstarts   = to->starts;
275   sprocs    = to->procs;
276   swaits    = contents->swaits;
277   svalues   = contents->svalues;
278 
279   rindices  = from->indices;
280   rstarts   = from->starts;
281   rprocs    = from->procs;
282   rwaits    = contents->rwaits;
283   rvalues   = contents->rvalues;
284 
285   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
286   ierr = MatGetArray(workB,&w);CHKERRQ(ierr);
287 
288   for (i=0; i<from->n; i++) {
289     ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr);
290   }
291 
292   for (i=0; i<to->n; i++) {
293     /* pack a message at a time */
294     CHKMEMQ;
295     for (j=0; j<sstarts[i+1]-sstarts[i]; j++){
296       for (k=0; k<ncols; k++) {
297         svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k];
298       }
299     }
300     CHKMEMQ;
301     ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr);
302   }
303 
304   nrecvs = from->n;
305   while (nrecvs) {
306     ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr);
307     nrecvs--;
308     /* unpack a message at a time */
309     CHKMEMQ;
310     for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++){
311       for (k=0; k<ncols; k++) {
312         w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k];
313       }
314     }
315     CHKMEMQ;
316   }
317   if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);}
318 
319   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
320   ierr = MatRestoreArray(workB,&w);CHKERRQ(ierr);
321   ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
322   ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
323   PetscFunctionReturn(0);
324 }
325 extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat);
326 
327 #undef __FUNCT__
328 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense"
329 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
330 {
331   PetscErrorCode       ierr;
332   Mat_MPIAIJ           *aij = (Mat_MPIAIJ*)A->data;
333   Mat_MPIDense         *bdense = (Mat_MPIDense*)B->data;
334   Mat_MPIDense         *cdense = (Mat_MPIDense*)C->data;
335   Mat                  workB;
336 
337   PetscFunctionBegin;
338 
339   /* diagonal block of A times all local rows of B*/
340   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr);
341 
342   /* get off processor parts of B needed to complete the product */
343   ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr);
344 
345   /* off-diagonal block of A times nonlocal rows of B */
346   ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr);
347   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
348   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
349   PetscFunctionReturn(0);
350 }
351 
352