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