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