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