xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision bed7c2b961e7f796158d2bfe4e13e61fa2401fe6)
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 = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
21     ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);/* numeric product is computed as well in old version! */
22     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
23   }
24   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
25   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
26   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
27   PetscFunctionReturn(0);
28 }
29 
30 #undef __FUNCT__
31 #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatMultMPI"
32 PetscErrorCode PetscContainerDestroy_Mat_MatMatMultMPI(void *ptr)
33 {
34   PetscErrorCode       ierr;
35   Mat_MatMatMultMPI    *mult=(Mat_MatMatMultMPI*)ptr;
36 
37   PetscFunctionBegin;
38   ierr = ISDestroy(&mult->isrowa);CHKERRQ(ierr);
39   ierr = ISDestroy(&mult->isrowb);CHKERRQ(ierr);
40   ierr = ISDestroy(&mult->iscolb);CHKERRQ(ierr);
41   ierr = MatDestroy(&mult->C_seq);CHKERRQ(ierr);
42   ierr = MatDestroy(&mult->A_loc);CHKERRQ(ierr);
43   ierr = MatDestroy(&mult->B_seq);CHKERRQ(ierr);
44   ierr = PetscFree(mult);CHKERRQ(ierr);
45   PetscFunctionReturn(0);
46 }
47 
48 #undef __FUNCT__
49 #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult_new"
50 PetscErrorCode MatDestroy_MPIAIJ_MatMatMult_new(Mat A)
51 {
52   PetscErrorCode ierr;
53   Mat_MPIAIJ     *a=(Mat_MPIAIJ*)A->data;
54   Mat_PtAPMPI    *ptap=a->ptap;
55 
56   PetscFunctionBegin;
57   ierr = PetscFree2(ptap->startsj,ptap->startsj_r);CHKERRQ(ierr);
58   ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
59   ierr = MatDestroy(&ptap->B_loc);CHKERRQ(ierr);
60   ierr = MatDestroy(&ptap->B_oth);CHKERRQ(ierr);
61   ierr = PetscFree(ptap->abi);CHKERRQ(ierr);
62   ierr = PetscFree(ptap->abj);CHKERRQ(ierr);
63   ierr = PetscFree(ptap->apa);CHKERRQ(ierr);
64   ierr = ptap->destroy(A);CHKERRQ(ierr);
65   ierr = PetscFree(ptap);CHKERRQ(ierr);
66   PetscFunctionReturn(0);
67 }
68 
69 #undef __FUNCT__
70 #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult"
71 PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
72 {
73   PetscErrorCode     ierr;
74   PetscContainer     container;
75   Mat_MatMatMultMPI  *mult=PETSC_NULL;
76 
77   PetscFunctionBegin;
78   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
79   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
80   ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
81   A->ops->destroy   = mult->destroy;
82   A->ops->duplicate = mult->duplicate;
83   if (A->ops->destroy) {
84     ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
85   }
86   ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatMultMPI",0);CHKERRQ(ierr);
87   PetscFunctionReturn(0);
88 }
89 
90 #undef __FUNCT__
91 #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult"
92 PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M)
93 {
94   PetscErrorCode     ierr;
95   Mat_MatMatMultMPI  *mult;
96   PetscContainer     container;
97 
98   PetscFunctionBegin;
99   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
100   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
101   ierr  = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
102   /* Note: the container is not duplicated, because it requires deep copying of
103      several large data sets (see PetscContainerDestroy_Mat_MatMatMultMPI()).
104      These data sets are only used for repeated calling of MatMatMultNumeric().
105      *M is unlikely being used in this way. Thus we create *M with pure mpiaij format */
106   ierr = (*mult->duplicate)(A,op,M);CHKERRQ(ierr);
107   (*M)->ops->destroy   = mult->destroy;   /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's container! */
108   (*M)->ops->duplicate = mult->duplicate; /* = MatDuplicate_MPIAIJ */
109   PetscFunctionReturn(0);
110 }
111 
112 #undef __FUNCT__
113 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_new"
114 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_new(Mat A,Mat P,Mat C)
115 {
116   PetscErrorCode     ierr;
117   Mat_MPIAIJ         *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
118   Mat_SeqAIJ         *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
119   Mat_SeqAIJ         *cd=(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
120   PetscInt           *adi=ad->i,*adj,*aoi=ao->i,*aoj;
121   PetscScalar        *ada,*aoa,*cda=cd->a,*coa=co->a;
122   Mat_SeqAIJ         *p_loc,*p_oth;
123   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
124   PetscScalar        *pa_loc,*pa_oth,*pa,*apa,valtmp,*ca;
125   PetscInt           cm=C->rmap->n,anz,pnz;
126   Mat_PtAPMPI        *ptap=c->ptap;
127   PetscInt           *api,*apj,*apJ,cnz,i,j,k,row;
128   MPI_Comm           comm=((PetscObject)A)->comm;
129   PetscMPIInt        rank;
130   PetscInt           rstart=C->rmap->rstart,cstart=C->cmap->rstart,cend=C->cmap->rend;
131   PetscInt           cdnz,conz,k0,k1;
132 
133   PetscFunctionBegin;
134   /* 1) get P_oth = ptap->B_oth  and P_loc = ptap->B_loc */
135   /*-----------------------------------------------------*/
136   if (ptap->reuse == MAT_INITIAL_MATRIX){
137     ptap->reuse = MAT_REUSE_MATRIX;
138   } else { /* update numerical values of P_oth and P_loc */
139     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->B_oth);CHKERRQ(ierr);
140     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->B_loc);CHKERRQ(ierr);
141   }
142 
143   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
144   /*----------------------------------------------------------*/
145   /* get data from symbolic products */
146   p_loc = (Mat_SeqAIJ*)(ptap->B_loc)->data;
147   p_oth = (Mat_SeqAIJ*)(ptap->B_oth)->data;
148   pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a;
149   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
150 
151   /* get apa for storing dense row A[i,:]*P */
152   apa = ptap->apa;
153 
154   for (i=0; i<cm; i++) {
155     /* diagonal portion of A */
156     anz = adi[i+1] - adi[i];
157     adj = ad->j + adi[i];
158     ada = ad->a + adi[i];
159     for (j=0; j<anz; j++) {
160       row = adj[j];
161       pnz = pi_loc[row+1] - pi_loc[row];
162       pj  = pj_loc + pi_loc[row];
163       pa  = pa_loc + pi_loc[row];
164 
165       /* perform dense axpy */
166       valtmp = ada[j];
167       for (k=0; k<pnz; k++){
168         apa[pj[k]] += valtmp*pa[k];
169       }
170       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
171     }
172 
173     /* off-diagonal portion of A */
174     anz = aoi[i+1] - aoi[i];
175     aoj = ao->j + aoi[i];
176     aoa = ao->a + aoi[i];
177     for (j=0; j<anz; j++) {
178       row = aoj[j];
179       pnz = pi_oth[row+1] - pi_oth[row];
180       pj  = pj_oth + pi_oth[row];
181       pa  = pa_oth + pi_oth[row];
182 
183       /* perform dense axpy */
184       valtmp = aoa[j];
185       for (k=0; k<pnz; k++){
186         apa[pj[k]] += valtmp*pa[k];
187       }
188       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
189     }
190 
191     /* set values in C */
192     row = rstart + i;
193     api = ptap->abi;
194     apj = ptap->abj;
195     apJ = apj + api[i];
196     cnz = api[i+1] - api[i];
197     cdnz = cd->i[i+1] - cd->i[i];
198     conz = co->i[i+1] - co->i[i];
199 
200     /* 1st off-diagoanl part of C */
201     ca = coa + co->i[i];
202     k  = 0;
203     for (k0=0; k0<conz; k0++){
204       if (apJ[k] >= cstart) break;
205       ca[k0]      = apa[apJ[k]];
206       apa[apJ[k]] = 0.0;
207       k++;
208     }
209 
210     /* diagonal part of C */
211     ca = cda + cd->i[i];
212     for (k1=0; k1<cdnz; k1++){
213       ca[k1]      = apa[apJ[k]];
214       apa[apJ[k]] = 0.0;
215       k++;
216     }
217 
218     /* 2nd off-diagoanl part of C */
219     ca = coa + co->i[i];
220     for (; k0<conz; k0++){
221       ca[k0]      = apa[apJ[k]];
222       apa[apJ[k]] = 0.0;
223       k++;
224     }
225   }
226   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
227   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
228   PetscFunctionReturn(0);
229 }
230 
231 #undef __FUNCT__
232 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_new"
233 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_new(Mat A,Mat P,PetscReal fill,Mat *C)
234 {
235   PetscErrorCode       ierr;
236   MPI_Comm             comm=((PetscObject)A)->comm;
237   PetscMPIInt          size,rank;
238   Mat                  B_mpi;
239   Mat_PtAPMPI          *ptap;
240   PetscFreeSpaceList   free_space=PETSC_NULL,current_space=PETSC_NULL;
241   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
242   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
243   Mat_SeqAIJ           *p_loc,*p_oth;
244   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ,*dnz,*onz;
245   PetscInt             *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,nnz,rstart=A->rmap->rstart;
246   PetscInt             nlnk,*lnk,*owners_co,*coi,*coj,i,k,pnz,row;
247   PetscInt             am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
248   PetscBT              lnkbt;
249   PetscInt             nzi,*bi,*bj;
250   PetscInt             *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j;
251   PetscScalar          *apa;
252 
253   PetscFunctionBegin;
254   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
255   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
256 
257   /* create struct Mat_PtAPMPI and attached it to C later */
258   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
259   ptap->abnz_max = 0;
260   ptap->reuse    = MAT_INITIAL_MATRIX;
261 
262   /* malloc apa to store dense row A[i,:]*P */
263   ierr = PetscMalloc((P->cmap->N)*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
264   ierr = PetscMemzero(apa,P->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
265   ptap->apa = apa;
266 
267 
268   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
269   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->B_oth);CHKERRQ(ierr);
270   /* get P_loc by taking all local rows of P */
271   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->B_loc);CHKERRQ(ierr);
272 
273   p_loc = (Mat_SeqAIJ*)(ptap->B_loc)->data;
274   p_oth = (Mat_SeqAIJ*)(ptap->B_oth)->data;
275   pi_loc = p_loc->i; pj_loc = p_loc->j;
276   pi_oth = p_oth->i; pj_oth = p_oth->j;
277 
278   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
279   /*-------------------------------------------------------------------*/
280   ierr  = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr);
281   ptap->abi = api;
282   api[0]    = 0;
283 
284   /* create and initialize a linked list */
285   nlnk = pN+1;
286   ierr = PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
287 
288   /* Initial FreeSpace size is fill*nnz(A) */
289   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am])),&free_space);CHKERRQ(ierr);
290   current_space = free_space;
291 
292   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
293   for (i=0;i<am;i++) {
294     apnz = 0;
295     /* diagonal portion of A */
296     nzi = adi[i+1] - adi[i];
297     for (j=0; j<nzi; j++){
298       row = *adj++;
299       pnz = pi_loc[row+1] - pi_loc[row];
300       Jptr  = pj_loc + pi_loc[row];
301       /* add non-zero cols of P into the sorted linked list lnk */
302       ierr = PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
303       apnz += nlnk;
304     }
305     /* off-diagonal portion of A */
306     nzi = aoi[i+1] - aoi[i];
307     for (j=0; j<nzi; j++){
308       row = *aoj++;
309       pnz = pi_oth[row+1] - pi_oth[row];
310       Jptr  = pj_oth + pi_oth[row];
311       ierr = PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
312       apnz += nlnk;
313     }
314 
315     api[i+1] = api[i] + apnz;
316     if (ptap->abnz_max < apnz) ptap->abnz_max = apnz;
317 
318     /* if free space is not available, double the total space in the list */
319     if (current_space->local_remaining<apnz) {
320       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
321       nspacedouble++;
322     }
323 
324     /* Copy data into free space, then initialize lnk */
325     ierr = PetscLLClean(pN,pN,apnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
326     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
327     current_space->array           += apnz;
328     current_space->local_used      += apnz;
329     current_space->local_remaining -= apnz;
330   }
331 
332   /* Allocate space for apj, initialize apj, and */
333   /* destroy list of free space and other temporary array(s) */
334   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ptap->abj);CHKERRQ(ierr);
335   apj = ptap->abj;
336   ierr = PetscFreeSpaceContiguous(&free_space,ptap->abj);CHKERRQ(ierr);
337   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
338 
339   /* create and assemble symbolic parallel matrix B_mpi */
340   /*----------------------------------------------------*/
341   ierr = MatCreate(comm,&B_mpi);CHKERRQ(ierr);
342   ierr = MatSetSizes(B_mpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
343   ierr = MatSetType(B_mpi,MATMPIAIJ);CHKERRQ(ierr);
344   ierr = MatMPIAIJSetPreallocation(B_mpi,0,dnz,0,onz);CHKERRQ(ierr);
345   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
346   ierr = MatSetBlockSize(B_mpi,1);CHKERRQ(ierr);
347   for (i=0; i<am; i++){
348     row = i + rstart;
349     apnz = api[i+1] - api[i];
350     ierr = MatSetValues(B_mpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
351     apj += apnz;
352   }
353   ierr = MatAssemblyBegin(B_mpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
354   ierr = MatAssemblyEnd(B_mpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
355   //ierr = MatView(B_mpi,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
356 
357   ptap->destroy              = B_mpi->ops->destroy;
358   B_mpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_new;
359   B_mpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult_new;
360 
361   /* attach the supporting struct to B_mpi for reuse */
362   c = (Mat_MPIAIJ*)B_mpi->data;
363   c->ptap  = ptap;
364 
365   *C = B_mpi;
366   PetscFunctionReturn(0);
367 }
368 
369 #undef __FUNCT__
370 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
371 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
372 {
373   PetscErrorCode     ierr;
374   Mat_MatMatMultMPI  *mult;
375   PetscContainer     container;
376   Mat                AB,*seq;
377   Mat_MPIAIJ         *a=(Mat_MPIAIJ*)A->data;
378   PetscInt           *idx,i,start,ncols,nzA,nzB,*cmap,imark;
379   PetscBool          matmatmult_new=PETSC_FALSE;
380 
381   PetscFunctionBegin;
382   if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend){
383     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);
384   }
385   ierr = PetscOptionsGetBool(PETSC_NULL,"-matmatmult_new",&matmatmult_new,PETSC_NULL);CHKERRQ(ierr);
386   if (matmatmult_new){
387     ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_new(A,B,fill,C);;CHKERRQ(ierr);
388     PetscFunctionReturn(0);
389   }
390 
391   ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr);
392 
393   /* get isrowb: nonzero col of A */
394   start = A->cmap->rstart;
395   cmap  = a->garray;
396   nzA   = a->A->cmap->n;
397   nzB   = a->B->cmap->n;
398   ierr  = PetscMalloc((nzA+nzB)*sizeof(PetscInt), &idx);CHKERRQ(ierr);
399   ncols = 0;
400   for (i=0; i<nzB; i++) {  /* row < local row index */
401     if (cmap[i] < start) idx[ncols++] = cmap[i];
402     else break;
403   }
404   imark = i;
405   for (i=0; i<nzA; i++) idx[ncols++] = start + i;  /* local rows */
406   for (i=imark; i<nzB; i++) idx[ncols++] = cmap[i]; /* row > local row index */
407   ierr = ISCreateGeneral(PETSC_COMM_SELF,ncols,idx,PETSC_OWN_POINTER,&mult->isrowb);CHKERRQ(ierr);
408   ierr = ISCreateStride(PETSC_COMM_SELF,B->cmap->N,0,1,&mult->iscolb);CHKERRQ(ierr);
409 
410   /*  get isrowa: all local rows of A */
411   ierr = ISCreateStride(PETSC_COMM_SELF,A->rmap->n,A->rmap->rstart,1,&mult->isrowa);CHKERRQ(ierr);
412 
413   /* Below should go to MatMatMultNumeric_MPIAIJ_MPIAIJ() - How to generate C there? */
414   /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */
415   ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_INITIAL_MATRIX,&seq);CHKERRQ(ierr);
416   mult->B_seq = *seq;
417   ierr = PetscFree(seq);CHKERRQ(ierr);
418 
419   /*  create a seq matrix A_seq = submatrix of A by taking all local rows of A */
420   ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_INITIAL_MATRIX,&seq);CHKERRQ(ierr);
421   mult->A_loc = *seq;
422   ierr = PetscFree(seq);CHKERRQ(ierr);
423 
424   /* compute C_seq = A_seq * B_seq */
425   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,fill,&mult->C_seq);CHKERRQ(ierr);
426   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(mult->A_loc,mult->B_seq,mult->C_seq);CHKERRQ(ierr);
427 
428   /* create mpi matrix C by concatinating C_seq */
429   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); /* prevent C_seq being destroyed by MatMerge() */
430   ierr = MatMergeSymbolic(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,&AB);CHKERRQ(ierr);
431   ierr = MatMergeNumeric(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,AB);CHKERRQ(ierr);
432 
433   /* attach the supporting struct to C for reuse of symbolic C */
434   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
435   ierr = PetscContainerSetPointer(container,mult);CHKERRQ(ierr);
436   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatMultMPI);CHKERRQ(ierr);
437   ierr = PetscObjectCompose((PetscObject)AB,"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr);
438   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
439   mult->destroy      = AB->ops->destroy;
440   mult->duplicate    = AB->ops->duplicate;
441   AB->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ;
442   AB->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
443   AB->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult;
444   AB->ops->matmult        = MatMatMult_MPIAIJ_MPIAIJ;
445   *C                      = AB;
446   PetscFunctionReturn(0);
447 }
448 
449 /* This routine is called ONLY in the case of reusing previously computed symbolic C */
450 #undef __FUNCT__
451 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ"
452 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat B,Mat C)
453 {
454   PetscErrorCode     ierr;
455   Mat                *seq;
456   Mat_MatMatMultMPI  *mult;
457   PetscContainer     container;
458 
459   PetscFunctionBegin;
460   ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
461   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
462   ierr  = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
463   seq = &mult->B_seq;
464   ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr);
465   mult->B_seq = *seq;
466 
467   seq = &mult->A_loc;
468   ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr);
469   mult->A_loc = *seq;
470 
471   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_SparseAxpy(mult->A_loc,mult->B_seq,mult->C_seq);CHKERRQ(ierr);
472 
473   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr);
474   ierr = MatMerge(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr);
475   PetscFunctionReturn(0);
476 }
477 
478 #undef __FUNCT__
479 #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense"
480 PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
481 {
482   PetscErrorCode ierr;
483 
484   PetscFunctionBegin;
485   if (scall == MAT_INITIAL_MATRIX){
486     ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr);
487   }
488   ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr);
489   PetscFunctionReturn(0);
490 }
491 
492 typedef struct {
493   Mat         workB;
494   PetscScalar *rvalues,*svalues;
495   MPI_Request *rwaits,*swaits;
496 } MPIAIJ_MPIDense;
497 
498 #undef __FUNCT__
499 #define __FUNCT__ "MPIAIJ_MPIDenseDestroy"
500 PetscErrorCode MPIAIJ_MPIDenseDestroy(void *ctx)
501 {
502   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx;
503   PetscErrorCode  ierr;
504 
505   PetscFunctionBegin;
506   ierr = MatDestroy(&contents->workB);CHKERRQ(ierr);
507   ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr);
508   ierr = PetscFree(contents);CHKERRQ(ierr);
509   PetscFunctionReturn(0);
510 }
511 
512 #undef __FUNCT__
513 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense"
514 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
515 {
516   PetscErrorCode         ierr;
517   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
518   PetscInt               nz = aij->B->cmap->n;
519   PetscContainer         container;
520   MPIAIJ_MPIDense        *contents;
521   VecScatter             ctx = aij->Mvctx;
522   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
523   VecScatter_MPI_General *to   = ( VecScatter_MPI_General*) ctx->todata;
524   PetscInt               m=A->rmap->n,n=B->cmap->n;
525 
526   PetscFunctionBegin;
527   ierr = MatCreate(((PetscObject)B)->comm,C);CHKERRQ(ierr);
528   ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr);
529   ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr);
530   ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
531   ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
532   (*C)->ops->matmult = MatMatMult_MPIAIJ_MPIDense;
533 
534   ierr = PetscNew(MPIAIJ_MPIDense,&contents);CHKERRQ(ierr);
535   /* Create work matrix used to store off processor rows of B needed for local product */
536   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,PETSC_NULL,&contents->workB);CHKERRQ(ierr);
537   /* Create work arrays needed */
538   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],PetscScalar,&contents->rvalues,
539                       B->cmap->N*to->starts[to->n],PetscScalar,&contents->svalues,
540                       from->n,MPI_Request,&contents->rwaits,
541                       to->n,MPI_Request,&contents->swaits);CHKERRQ(ierr);
542 
543   ierr = PetscContainerCreate(((PetscObject)A)->comm,&container);CHKERRQ(ierr);
544   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
545   ierr = PetscContainerSetUserDestroy(container,MPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
546   ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr);
547   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
548   PetscFunctionReturn(0);
549 }
550 
551 #undef __FUNCT__
552 #define __FUNCT__ "MatMPIDenseScatter"
553 /*
554     Performs an efficient scatter on the rows of B needed by this process; this is
555     a modification of the VecScatterBegin_() routines.
556 */
557 PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB)
558 {
559   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)A->data;
560   PetscErrorCode         ierr;
561   PetscScalar            *b,*w,*svalues,*rvalues;
562   VecScatter             ctx = aij->Mvctx;
563   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
564   VecScatter_MPI_General *to   = ( VecScatter_MPI_General*) ctx->todata;
565   PetscInt               i,j,k;
566   PetscInt               *sindices,*sstarts,*rindices,*rstarts;
567   PetscMPIInt            *sprocs,*rprocs,nrecvs;
568   MPI_Request            *swaits,*rwaits;
569   MPI_Comm               comm = ((PetscObject)A)->comm;
570   PetscMPIInt            tag = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n;
571   MPI_Status             status;
572   MPIAIJ_MPIDense        *contents;
573   PetscContainer         container;
574   Mat                    workB;
575 
576   PetscFunctionBegin;
577   ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr);
578   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
579   ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr);
580 
581   workB = *outworkB = contents->workB;
582   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);
583   sindices  = to->indices;
584   sstarts   = to->starts;
585   sprocs    = to->procs;
586   swaits    = contents->swaits;
587   svalues   = contents->svalues;
588 
589   rindices  = from->indices;
590   rstarts   = from->starts;
591   rprocs    = from->procs;
592   rwaits    = contents->rwaits;
593   rvalues   = contents->rvalues;
594 
595   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
596   ierr = MatGetArray(workB,&w);CHKERRQ(ierr);
597 
598   for (i=0; i<from->n; i++) {
599     ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr);
600   }
601 
602   for (i=0; i<to->n; i++) {
603     /* pack a message at a time */
604     CHKMEMQ;
605     for (j=0; j<sstarts[i+1]-sstarts[i]; j++){
606       for (k=0; k<ncols; k++) {
607         svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k];
608       }
609     }
610     CHKMEMQ;
611     ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr);
612   }
613 
614   nrecvs = from->n;
615   while (nrecvs) {
616     ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr);
617     nrecvs--;
618     /* unpack a message at a time */
619     CHKMEMQ;
620     for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++){
621       for (k=0; k<ncols; k++) {
622         w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k];
623       }
624     }
625     CHKMEMQ;
626   }
627   if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);}
628 
629   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
630   ierr = MatRestoreArray(workB,&w);CHKERRQ(ierr);
631   ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
632   ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
633   PetscFunctionReturn(0);
634 }
635 extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat);
636 
637 #undef __FUNCT__
638 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense"
639 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
640 {
641   PetscErrorCode       ierr;
642   Mat_MPIAIJ           *aij = (Mat_MPIAIJ*)A->data;
643   Mat_MPIDense         *bdense = (Mat_MPIDense*)B->data;
644   Mat_MPIDense         *cdense = (Mat_MPIDense*)C->data;
645   Mat                  workB;
646 
647   PetscFunctionBegin;
648 
649   /* diagonal block of A times all local rows of B*/
650   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr);
651 
652   /* get off processor parts of B needed to complete the product */
653   ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr);
654 
655   /* off-diagonal block of A times nonlocal rows of B */
656   ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr);
657   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
658   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
659   PetscFunctionReturn(0);
660 }
661 
662