xref: /petsc/src/mat/impls/aij/mpi/mpimatmatmult.c (revision dcffc833068d8d86d19a4b82060221dc7997d85b)
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 #define DEBUG_MATTrMatMult
16  */
17 
18 #undef __FUNCT__
19 #define __FUNCT__ "MatMatMult_MPIAIJ_MPIAIJ"
20 PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C)
21 {
22   PetscErrorCode ierr;
23 
24   PetscFunctionBegin;
25   if (scall == MAT_INITIAL_MATRIX){
26     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
27     ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);
28     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
29   }
30   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
31   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
32   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
33   PetscFunctionReturn(0);
34 }
35 
36 #undef __FUNCT__
37 #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatMultMPI"
38 PetscErrorCode PetscContainerDestroy_Mat_MatMatMultMPI(void *ptr)
39 {
40   PetscErrorCode       ierr;
41   Mat_MatMatMultMPI    *mult=(Mat_MatMatMultMPI*)ptr;
42 
43   PetscFunctionBegin;
44   ierr = ISDestroy(&mult->isrowa);CHKERRQ(ierr);
45   ierr = ISDestroy(&mult->isrowb);CHKERRQ(ierr);
46   ierr = ISDestroy(&mult->iscolb);CHKERRQ(ierr);
47   ierr = MatDestroy(&mult->C_seq);CHKERRQ(ierr);
48   ierr = MatDestroy(&mult->A_loc);CHKERRQ(ierr);
49   ierr = MatDestroy(&mult->B_seq);CHKERRQ(ierr);
50   ierr = PetscFree(mult);CHKERRQ(ierr);
51   PetscFunctionReturn(0);
52 }
53 
54 #undef __FUNCT__
55 #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult"
56 PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
57 {
58   PetscErrorCode ierr;
59   Mat_MPIAIJ     *a=(Mat_MPIAIJ*)A->data;
60   Mat_PtAPMPI    *ptap=a->ptap;
61 
62   PetscFunctionBegin;
63   ierr = PetscFree2(ptap->startsj,ptap->startsj_r);CHKERRQ(ierr);
64   ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
65   ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
66   ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
67   ierr = PetscFree(ptap->api);CHKERRQ(ierr);
68   ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
69   ierr = PetscFree(ptap->apa);CHKERRQ(ierr);
70   ierr = ptap->destroy(A);CHKERRQ(ierr);
71   ierr = PetscFree(ptap);CHKERRQ(ierr);
72   PetscFunctionReturn(0);
73 }
74 
75 #undef __FUNCT__
76 #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult_32"
77 PetscErrorCode MatDestroy_MPIAIJ_MatMatMult_32(Mat A)
78 {
79   PetscErrorCode     ierr;
80   PetscContainer     container;
81   Mat_MatMatMultMPI  *mult=PETSC_NULL;
82 
83   PetscFunctionBegin;
84   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
85   if (!container) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_PLIB,"Container does not exist");
86   ierr = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
87   A->ops->destroy   = mult->destroy;
88   A->ops->duplicate = mult->duplicate;
89   if (A->ops->destroy) {
90     ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
91   }
92   ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatMultMPI",0);CHKERRQ(ierr);
93   PetscFunctionReturn(0);
94 }
95 
96 #undef __FUNCT__
97 #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult_32"
98 PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult_32(Mat A, MatDuplicateOption op, Mat *M)
99 {
100   PetscErrorCode     ierr;
101   Mat_MatMatMultMPI  *mult;
102   PetscContainer     container;
103 
104   PetscFunctionBegin;
105   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
106   if (!container) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_PLIB,"Container does not exist");
107   ierr  = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
108   /* Note: the container is not duplicated, because it requires deep copying of
109      several large data sets (see PetscContainerDestroy_Mat_MatMatMultMPI()).
110      These data sets are only used for repeated calling of MatMatMultNumeric().
111      *M is unlikely being used in this way. Thus we create *M with pure mpiaij format */
112   ierr = (*mult->duplicate)(A,op,M);CHKERRQ(ierr);
113   (*M)->ops->destroy   = mult->destroy;   /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's container! */
114   (*M)->ops->duplicate = mult->duplicate; /* = MatDuplicate_MPIAIJ */
115   PetscFunctionReturn(0);
116 }
117 
118 #undef __FUNCT__
119 #define __FUNCT__ "MatDuplicate_MPIAIJ_MatMatMult"
120 PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M)
121 {
122   PetscErrorCode     ierr;
123   Mat_MPIAIJ         *a=(Mat_MPIAIJ*)A->data;
124   Mat_PtAPMPI        *ptap=a->ptap;
125 
126   PetscFunctionBegin;
127   ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr);
128   (*M)->ops->destroy   = ptap->destroy;   /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's special structure! */
129   (*M)->ops->duplicate = ptap->duplicate; /* = MatDuplicate_MPIAIJ */
130   PetscFunctionReturn(0);
131 }
132 
133 #undef __FUNCT__
134 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ"
135 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
136 {
137   PetscErrorCode     ierr;
138   Mat_MPIAIJ         *a=(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
139   Mat_SeqAIJ         *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
140   Mat_SeqAIJ         *cd=(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
141   PetscInt           *adi=ad->i,*adj,*aoi=ao->i,*aoj;
142   PetscScalar        *ada,*aoa,*cda=cd->a,*coa=co->a;
143   Mat_SeqAIJ         *p_loc,*p_oth;
144   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
145   PetscScalar        *pa_loc,*pa_oth,*pa,*apa,valtmp,*ca;
146   PetscInt           cm=C->rmap->n,anz,pnz;
147   Mat_PtAPMPI        *ptap=c->ptap;
148   PetscInt           *api,*apj,*apJ,i,j,k,row;
149   PetscInt           cstart=C->cmap->rstart;
150   PetscInt           cdnz,conz,k0,k1;
151 #if defined(DEBUG_MATMATMULT)
152   PetscMPIInt        rank=a->rank;
153 #endif
154 
155   PetscFunctionBegin;
156 #if defined(DEBUG_MATMATMULT)
157   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultNumeric_MPIAIJ_MPIAIJ()...\n",rank);
158 #endif
159 
160   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
161   /*-----------------------------------------------------*/
162   /* update numerical values of P_oth and P_loc */
163   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
164   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
165 #if defined(DEBUG_MATMATMULT)
166   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] got P_oth and P_loc...\n",rank);
167 #endif
168 
169   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
170   /*----------------------------------------------------------*/
171   /* get data from symbolic products */
172   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
173   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
174   pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a;
175   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
176 
177   /* get apa for storing dense row A[i,:]*P */
178   apa = ptap->apa;
179 
180   api = ptap->api;
181   apj = ptap->apj;
182   for (i=0; i<cm; i++) {
183     /* diagonal portion of A */
184     anz = adi[i+1] - adi[i];
185     adj = ad->j + adi[i];
186     ada = ad->a + adi[i];
187     for (j=0; j<anz; j++) {
188       row = adj[j];
189       pnz = pi_loc[row+1] - pi_loc[row];
190       pj  = pj_loc + pi_loc[row];
191       pa  = pa_loc + pi_loc[row];
192 
193       /* perform dense axpy */
194       valtmp = ada[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     /* off-diagonal portion of A */
202     anz = aoi[i+1] - aoi[i];
203     aoj = ao->j + aoi[i];
204     aoa = ao->a + aoi[i];
205     for (j=0; j<anz; j++) {
206       row = aoj[j];
207       pnz = pi_oth[row+1] - pi_oth[row];
208       pj  = pj_oth + pi_oth[row];
209       pa  = pa_oth + pi_oth[row];
210 
211       /* perform dense axpy */
212       valtmp = aoa[j];
213       for (k=0; k<pnz; k++){
214         apa[pj[k]] += valtmp*pa[k];
215       }
216       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
217     }
218 
219     /* set values in C */
220     apJ = apj + api[i];
221     cdnz = cd->i[i+1] - cd->i[i];
222     conz = co->i[i+1] - co->i[i];
223 
224     /* 1st off-diagoanl part of C */
225     ca = coa + co->i[i];
226     k  = 0;
227     for (k0=0; k0<conz; k0++){
228       if (apJ[k] >= cstart) break;
229       ca[k0]      = apa[apJ[k]];
230       apa[apJ[k]] = 0.0;
231       k++;
232     }
233 
234     /* diagonal part of C */
235     ca = cda + cd->i[i];
236     for (k1=0; k1<cdnz; k1++){
237       ca[k1]      = apa[apJ[k]];
238       apa[apJ[k]] = 0.0;
239       k++;
240     }
241 
242     /* 2nd off-diagoanl part of C */
243     ca = coa + co->i[i];
244     for (; k0<conz; k0++){
245       ca[k0]      = apa[apJ[k]];
246       apa[apJ[k]] = 0.0;
247       k++;
248     }
249   }
250   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
251   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
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             *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            scalable=PETSC_FALSE;
274   PetscInt             nlnk_max,armax,prmax;
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(comm,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 
284   ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
285 
286     ierr = PetscOptionsBool("-matmatmult_32","Use a scalable but slower C=A*B","",scalable,&scalable,PETSC_NULL);CHKERRQ(ierr);
287     if (scalable){
288       ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable_32(A,P,fill,C);;CHKERRQ(ierr);
289       PetscFunctionReturn(0);
290     }
291     ierr = PetscOptionsBool("-matmatmult_scalable","Use a scalable but slower C=A*B","",scalable,&scalable,PETSC_NULL);CHKERRQ(ierr);
292     if (scalable){
293       ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(A,P,fill,C);;CHKERRQ(ierr);
294       PetscFunctionReturn(0);
295     }
296   ierr = PetscOptionsEnd();CHKERRQ(ierr);
297 
298 #if defined(DEBUG_MATMATMULT)
299   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultSymbolic_MPIAIJ_MPIAIJ()...\n",rank);
300 #endif
301 
302   /* create struct Mat_PtAPMPI and attached it to C later */
303   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
304 
305   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
306   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
307 #if defined(DEBUG_MATMATMULT)
308   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_oth is done...\n",rank);
309 #endif
310   /* get P_loc by taking all local rows of P */
311   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
312 
313   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
314   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
315   pi_loc = p_loc->i; pj_loc = p_loc->j;
316   pi_oth = p_oth->i; pj_oth = p_oth->j;
317 
318 #if defined(DEBUG_MATMATMULT)
319   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_loc is done, Annz %D * P_locnnz %D = %D...\n",rank,ad->i[am]+ao->i[am],p_loc->rmax,(ad->i[am]+ao->i[am])*p_loc->rmax);
320 #endif
321 
322   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
323   /*-------------------------------------------------------------------*/
324   ierr  = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr);
325   ptap->api = api;
326   api[0]    = 0;
327 
328   /* create and initialize a linked list */
329   armax = ad->rmax+ao->rmax;
330   prmax = PetscMax(p_loc->rmax,p_oth->rmax);
331   nlnk_max = armax*prmax;
332   if (!nlnk_max || nlnk_max > pN) nlnk_max = pN;
333 #if defined(DEBUG_MATMATMULT)
334   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] pN %d; nlnk_max %d; Armax %d+%d=%d; Prmax Max(%d,%d)=%d\n",rank,pN,nlnk_max,ad->rmax,ao->rmax,armax,p_loc->rmax,p_oth->rmax,prmax);
335 #endif
336   ierr = PetscLLCondensedCreate(nlnk_max,pN,&lnk,&lnkbt);CHKERRQ(ierr);
337 
338   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
339   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
340   current_space = free_space;
341 
342   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
343   for (i=0; i<am; i++) {
344     apnz = 0;
345     /* diagonal portion of A */
346     nzi = adi[i+1] - adi[i];
347     for (j=0; j<nzi; j++){
348       row = *adj++;
349       pnz = pi_loc[row+1] - pi_loc[row];
350       Jptr  = pj_loc + pi_loc[row];
351       /* add non-zero cols of P into the sorted linked list lnk */
352       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
353     }
354     /* off-diagonal portion of A */
355     nzi = aoi[i+1] - aoi[i];
356     for (j=0; j<nzi; j++){
357       row = *aoj++;
358       pnz = pi_oth[row+1] - pi_oth[row];
359       Jptr  = pj_oth + pi_oth[row];
360       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
361     }
362 
363     apnz     = lnk[0];
364     api[i+1] = api[i] + apnz;
365 
366     /* if free space is not available, double the total space in the list */
367     if (current_space->local_remaining<apnz) {
368       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
369       nspacedouble++;
370     }
371 
372     /* Copy data into free space, then initialize lnk */
373     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
374     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
375     current_space->array           += apnz;
376     current_space->local_used      += apnz;
377     current_space->local_remaining -= apnz;
378   }
379 
380   /* Allocate space for apj, initialize apj, and */
381   /* destroy list of free space and other temporary array(s) */
382   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ptap->apj);CHKERRQ(ierr);
383   apj  = ptap->apj;
384   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
385   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
386 #if defined(DEBUG_MATMATMULT)
387   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] AP is done...\n",rank);
388 #endif
389 
390   /* malloc apa to store dense row A[i,:]*P */
391   ierr = PetscMalloc(pN*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
392   ierr = PetscMemzero(apa,pN*sizeof(PetscScalar));CHKERRQ(ierr);
393   ptap->apa = apa;
394 #if defined(DEBUG_MATMATMULT)
395   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Malloc apa pN %D is done...\n",rank,pN);
396 #endif
397 
398   /* create and assemble symbolic parallel matrix Cmpi */
399   /*----------------------------------------------------*/
400   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
401   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
402   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
403   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
404   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
405   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
406   for (i=0; i<am; i++){
407     row  = i + rstart;
408     apnz = api[i+1] - api[i];
409     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
410     apj += apnz;
411   }
412   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
413   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
414 
415   ptap->destroy             = Cmpi->ops->destroy;
416   ptap->duplicate           = Cmpi->ops->duplicate;
417   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ;
418   Cmpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
419   Cmpi->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult;
420 
421   /* attach the supporting struct to Cmpi for reuse */
422   c = (Mat_MPIAIJ*)Cmpi->data;
423   c->ptap  = ptap;
424 
425   *C = Cmpi;
426 
427   /* set MatInfo */
428   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]) + 1.e-5;
429   if (afill < 1.0) afill = 1.0;
430   Cmpi->info.mallocs           = nspacedouble;
431   Cmpi->info.fill_ratio_given  = fill;
432   Cmpi->info.fill_ratio_needed = afill;
433 
434 #if defined(PETSC_USE_INFO)
435   if (api[am]) {
436     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
437     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
438   } else {
439     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
440   }
441 #endif
442   PetscFunctionReturn(0);
443 }
444 
445 /* implementation used in PETSc-3.2 */
446 /* This routine is called ONLY in the case of reusing previously computed symbolic C */
447 #undef __FUNCT__
448 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable_32"
449 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable_32(Mat A,Mat B,Mat C)
450 {
451   PetscErrorCode     ierr;
452   Mat                *seq;
453   Mat_MatMatMultMPI  *mult;
454   PetscContainer     container;
455 
456   PetscFunctionBegin;
457   ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
458   if (!container) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_PLIB,"Container does not exist");
459   ierr  = PetscContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
460 
461   if (mult->skipNumeric){ /* first numeric product is done during symbolic product */
462     mult->skipNumeric = PETSC_FALSE;
463     PetscFunctionReturn(0);
464   }
465 #if defined(DEBUG_MATMATMULT)
466   PetscMPIInt rank;
467   MPI_Comm    comm = ((PetscObject)C)->comm;
468   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
469   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
470   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable_32()...\n",rank);
471 #endif
472 
473   seq = &mult->B_seq;
474   ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr);
475   mult->B_seq = *seq;
476 
477   seq = &mult->A_loc;
478   ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&seq);CHKERRQ(ierr);
479   mult->A_loc = *seq;
480 
481   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(mult->A_loc,mult->B_seq,mult->C_seq);CHKERRQ(ierr);
482 
483   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr);
484   ierr = MatMerge(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr);
485   PetscFunctionReturn(0);
486 }
487 
488 /* numeric product is computed as well */
489 #undef __FUNCT__
490 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable_32"
491 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable_32(Mat A,Mat B,PetscReal fill,Mat *C)
492 {
493   PetscErrorCode     ierr;
494   Mat_MatMatMultMPI  *mult;
495   PetscContainer     container;
496   Mat                AB,*seq;
497   Mat_MPIAIJ         *a=(Mat_MPIAIJ*)A->data;
498   PetscInt           *idx,i,start,ncols,nzA,nzB,*cmap,imark;
499 #if defined(DEBUG_MATMATMULT)
500   MPI_Comm             comm = ((PetscObject)A)->comm;
501   PetscMPIInt          rank=a->rank;
502 #endif
503 
504   PetscFunctionBegin;
505 #if defined(DEBUG_MATMATMULT)
506   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable_32()...\n",rank);
507 #endif
508   if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend){
509     SETERRQ4(((PetscObject)A)->comm,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);
510   }
511 
512   ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr);
513 
514   /* get isrowb: nonzero col of A */
515   start = A->cmap->rstart;
516   cmap  = a->garray;
517   nzA   = a->A->cmap->n;
518   nzB   = a->B->cmap->n;
519   ierr  = PetscMalloc((nzA+nzB)*sizeof(PetscInt), &idx);CHKERRQ(ierr);
520   ncols = 0;
521   for (i=0; i<nzB; i++) {  /* row < local row index */
522     if (cmap[i] < start) idx[ncols++] = cmap[i];
523     else break;
524   }
525   imark = i;
526   for (i=0; i<nzA; i++) idx[ncols++] = start + i;  /* local rows */
527   for (i=imark; i<nzB; i++) idx[ncols++] = cmap[i]; /* row > local row index */
528   ierr = ISCreateGeneral(PETSC_COMM_SELF,ncols,idx,PETSC_OWN_POINTER,&mult->isrowb);CHKERRQ(ierr);
529   ierr = ISCreateStride(PETSC_COMM_SELF,B->cmap->N,0,1,&mult->iscolb);CHKERRQ(ierr);
530 
531   /*  get isrowa: all local rows of A */
532   ierr = ISCreateStride(PETSC_COMM_SELF,A->rmap->n,A->rmap->rstart,1,&mult->isrowa);CHKERRQ(ierr);
533 
534   /* Below should go to MatMatMultNumeric_MPIAIJ_MPIAIJ() - How to generate C there? */
535   /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */
536   ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_INITIAL_MATRIX,&seq);CHKERRQ(ierr);
537   mult->B_seq = *seq;
538   ierr = PetscFree(seq);CHKERRQ(ierr);
539 #if defined(DEBUG_MATMATMULT)
540   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] B_seq is done...\n",rank);
541 #endif
542 
543   /*  create a seq matrix A_seq = submatrix of A by taking all local rows of A */
544   ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_INITIAL_MATRIX,&seq);CHKERRQ(ierr);
545   mult->A_loc = *seq;
546   ierr = PetscFree(seq);CHKERRQ(ierr);
547 #if defined(DEBUG_MATMATMULT)
548   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] A_loc is done...\n",rank);
549 #endif
550 
551   /* compute C_seq = A_seq * B_seq */
552   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(mult->A_loc,mult->B_seq,fill,&mult->C_seq);CHKERRQ(ierr);
553 #if defined(DEBUG_MATMATMULT)
554   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
555   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] C_seq Symbolic is done...\n",rank);
556 #endif
557   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(mult->A_loc,mult->B_seq,mult->C_seq);CHKERRQ(ierr);
558 #if defined(DEBUG_MATMATMULT)
559   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] C_seq Numeric is done...\n",rank);
560 #endif
561 
562   /* create mpi matrix C by concatinating C_seq */
563   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); /* prevent C_seq being destroyed by MatMerge() */
564   ierr = MatMergeSymbolic(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,&AB);CHKERRQ(ierr);
565   ierr = MatMergeNumeric(((PetscObject)A)->comm,mult->C_seq,B->cmap->n,AB);CHKERRQ(ierr);
566 #if defined(DEBUG_MATMATMULT)
567   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Merge is done...\n",rank);
568 #endif
569 
570   /* attach the supporting struct to C for reuse of symbolic C */
571   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
572   ierr = PetscContainerSetPointer(container,mult);CHKERRQ(ierr);
573   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatMultMPI);CHKERRQ(ierr);
574   ierr = PetscObjectCompose((PetscObject)AB,"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr);
575   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
576   mult->skipNumeric  = PETSC_TRUE; /* a numeric product is done here */
577   mult->destroy      = AB->ops->destroy;
578   mult->duplicate    = AB->ops->duplicate;
579   AB->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable_32;
580   AB->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult_32;
581   AB->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult_32;
582   AB->ops->matmult        = MatMatMult_MPIAIJ_MPIAIJ;
583   *C                      = AB;
584   PetscFunctionReturn(0);
585 }
586 
587 #undef __FUNCT__
588 #define __FUNCT__ "MatMatMult_MPIAIJ_MPIDense"
589 PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
590 {
591   PetscErrorCode ierr;
592 
593   PetscFunctionBegin;
594   if (scall == MAT_INITIAL_MATRIX){
595     ierr = MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr);
596   }
597   ierr = MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr);
598   PetscFunctionReturn(0);
599 }
600 
601 typedef struct {
602   Mat         workB;
603   PetscScalar *rvalues,*svalues;
604   MPI_Request *rwaits,*swaits;
605 } MPIAIJ_MPIDense;
606 
607 #undef __FUNCT__
608 #define __FUNCT__ "MPIAIJ_MPIDenseDestroy"
609 PetscErrorCode MPIAIJ_MPIDenseDestroy(void *ctx)
610 {
611   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx;
612   PetscErrorCode  ierr;
613 
614   PetscFunctionBegin;
615   ierr = MatDestroy(&contents->workB);CHKERRQ(ierr);
616   ierr = PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);CHKERRQ(ierr);
617   ierr = PetscFree(contents);CHKERRQ(ierr);
618   PetscFunctionReturn(0);
619 }
620 
621 #undef __FUNCT__
622 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIDense"
623 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
624 {
625   PetscErrorCode         ierr;
626   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
627   PetscInt               nz = aij->B->cmap->n;
628   PetscContainer         container;
629   MPIAIJ_MPIDense        *contents;
630   VecScatter             ctx = aij->Mvctx;
631   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
632   VecScatter_MPI_General *to   = ( VecScatter_MPI_General*) ctx->todata;
633   PetscInt               m=A->rmap->n,n=B->cmap->n;
634 
635   PetscFunctionBegin;
636   ierr = MatCreate(((PetscObject)B)->comm,C);CHKERRQ(ierr);
637   ierr = MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);CHKERRQ(ierr);
638   ierr = MatSetType(*C,MATMPIDENSE);CHKERRQ(ierr);
639   ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
640   ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
641   (*C)->ops->matmult = MatMatMult_MPIAIJ_MPIDense;
642 
643   ierr = PetscNew(MPIAIJ_MPIDense,&contents);CHKERRQ(ierr);
644   /* Create work matrix used to store off processor rows of B needed for local product */
645   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,PETSC_NULL,&contents->workB);CHKERRQ(ierr);
646   /* Create work arrays needed */
647   ierr = PetscMalloc4(B->cmap->N*from->starts[from->n],PetscScalar,&contents->rvalues,
648                       B->cmap->N*to->starts[to->n],PetscScalar,&contents->svalues,
649                       from->n,MPI_Request,&contents->rwaits,
650                       to->n,MPI_Request,&contents->swaits);CHKERRQ(ierr);
651 
652   ierr = PetscContainerCreate(((PetscObject)A)->comm,&container);CHKERRQ(ierr);
653   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
654   ierr = PetscContainerSetUserDestroy(container,MPIAIJ_MPIDenseDestroy);CHKERRQ(ierr);
655   ierr = PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);CHKERRQ(ierr);
656   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
657   PetscFunctionReturn(0);
658 }
659 
660 #undef __FUNCT__
661 #define __FUNCT__ "MatMPIDenseScatter"
662 /*
663     Performs an efficient scatter on the rows of B needed by this process; this is
664     a modification of the VecScatterBegin_() routines.
665 */
666 PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB)
667 {
668   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)A->data;
669   PetscErrorCode         ierr;
670   PetscScalar            *b,*w,*svalues,*rvalues;
671   VecScatter             ctx = aij->Mvctx;
672   VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
673   VecScatter_MPI_General *to   = ( VecScatter_MPI_General*) ctx->todata;
674   PetscInt               i,j,k;
675   PetscInt               *sindices,*sstarts,*rindices,*rstarts;
676   PetscMPIInt            *sprocs,*rprocs,nrecvs;
677   MPI_Request            *swaits,*rwaits;
678   MPI_Comm               comm = ((PetscObject)A)->comm;
679   PetscMPIInt            tag = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n;
680   MPI_Status             status;
681   MPIAIJ_MPIDense        *contents;
682   PetscContainer         container;
683   Mat                    workB;
684 
685   PetscFunctionBegin;
686   ierr = PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);CHKERRQ(ierr);
687   if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist");
688   ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr);
689 
690   workB = *outworkB = contents->workB;
691   if (nrows != workB->rmap->n) SETERRQ2(comm,PETSC_ERR_PLIB,"Number of rows of workB %D not equal to columns of aij->B %D",nrows,workB->cmap->n);
692   sindices  = to->indices;
693   sstarts   = to->starts;
694   sprocs    = to->procs;
695   swaits    = contents->swaits;
696   svalues   = contents->svalues;
697 
698   rindices  = from->indices;
699   rstarts   = from->starts;
700   rprocs    = from->procs;
701   rwaits    = contents->rwaits;
702   rvalues   = contents->rvalues;
703 
704   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
705   ierr = MatGetArray(workB,&w);CHKERRQ(ierr);
706 
707   for (i=0; i<from->n; i++) {
708     ierr = MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);CHKERRQ(ierr);
709   }
710 
711   for (i=0; i<to->n; i++) {
712     /* pack a message at a time */
713     CHKMEMQ;
714     for (j=0; j<sstarts[i+1]-sstarts[i]; j++){
715       for (k=0; k<ncols; k++) {
716         svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k];
717       }
718     }
719     CHKMEMQ;
720     ierr = MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);CHKERRQ(ierr);
721   }
722 
723   nrecvs = from->n;
724   while (nrecvs) {
725     ierr = MPI_Waitany(from->n,rwaits,&imdex,&status);CHKERRQ(ierr);
726     nrecvs--;
727     /* unpack a message at a time */
728     CHKMEMQ;
729     for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++){
730       for (k=0; k<ncols; k++) {
731         w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k];
732       }
733     }
734     CHKMEMQ;
735   }
736   if (to->n) {ierr = MPI_Waitall(to->n,swaits,to->sstatus);CHKERRQ(ierr);}
737 
738   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
739   ierr = MatRestoreArray(workB,&w);CHKERRQ(ierr);
740   ierr = MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
741   ierr = MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
742   PetscFunctionReturn(0);
743 }
744 extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat);
745 
746 #undef __FUNCT__
747 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIDense"
748 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
749 {
750   PetscErrorCode       ierr;
751   Mat_MPIAIJ           *aij = (Mat_MPIAIJ*)A->data;
752   Mat_MPIDense         *bdense = (Mat_MPIDense*)B->data;
753   Mat_MPIDense         *cdense = (Mat_MPIDense*)C->data;
754   Mat                  workB;
755 
756   PetscFunctionBegin;
757 
758   /* diagonal block of A times all local rows of B*/
759   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);CHKERRQ(ierr);
760 
761   /* get off processor parts of B needed to complete the product */
762   ierr = MatMPIDenseScatter(A,B,C,&workB);CHKERRQ(ierr);
763 
764   /* off-diagonal block of A times nonlocal rows of B */
765   ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);CHKERRQ(ierr);
766   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
767   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
768   PetscFunctionReturn(0);
769 }
770 
771 #undef __FUNCT__
772 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable"
773 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable(Mat A,Mat P,Mat C)
774 {
775   PetscErrorCode     ierr;
776   Mat_MPIAIJ         *a=(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
777   Mat_SeqAIJ         *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
778   Mat_SeqAIJ         *cd=(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
779   PetscInt           *adi=ad->i,*adj,*aoi=ao->i,*aoj;
780   PetscScalar        *ada,*aoa,*cda=cd->a,*coa=co->a;
781   Mat_SeqAIJ         *p_loc,*p_oth;
782   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
783   PetscScalar        *pa_loc,*pa_oth,*pa,valtmp,*ca;
784   PetscInt           cm=C->rmap->n,anz,pnz;
785   Mat_PtAPMPI        *ptap=c->ptap;
786   PetscScalar        *apa_sparse=ptap->apa;
787   PetscInt           *api,*apj,*apJ,i,j,k,row;
788   PetscInt           cstart=C->cmap->rstart;
789   PetscInt           cdnz,conz,k0,k1,nextp;
790 #if defined(DEBUG_MATMATMULT)
791   PetscMPIInt        rank=a->rank;
792 #endif
793 
794   PetscFunctionBegin;
795 #if defined(DEBUG_MATMATMULT)
796   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable()...\n",rank);
797 #endif
798 
799   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
800   /*-----------------------------------------------------*/
801   /* update numerical values of P_oth and P_loc */
802   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
803   ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
804 #if defined(DEBUG_MATMATMULT)
805   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] got P_oth and P_loc...\n",rank);
806 #endif
807 
808   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
809   /*----------------------------------------------------------*/
810   /* get data from symbolic products */
811   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
812   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
813   pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a;
814   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
815 
816   api = ptap->api;
817   apj = ptap->apj;
818   for (i=0; i<cm; i++) {
819     apJ = apj + api[i];
820 
821     /* diagonal portion of A */
822     anz = adi[i+1] - adi[i];
823     adj = ad->j + adi[i];
824     ada = ad->a + adi[i];
825     for (j=0; j<anz; j++) {
826       row = adj[j];
827       pnz = pi_loc[row+1] - pi_loc[row];
828       pj  = pj_loc + pi_loc[row];
829       pa  = pa_loc + pi_loc[row];
830       /* perform sparse axpy */
831       valtmp = ada[j];
832       nextp  = 0;
833       for (k=0; nextp<pnz; k++) {
834         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
835           apa_sparse[k] += valtmp*pa[nextp++];
836         }
837       }
838       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
839     }
840 
841     /* off-diagonal portion of A */
842     anz = aoi[i+1] - aoi[i];
843     aoj = ao->j + aoi[i];
844     aoa = ao->a + aoi[i];
845     for (j=0; j<anz; j++) {
846       row = aoj[j];
847       pnz = pi_oth[row+1] - pi_oth[row];
848       pj  = pj_oth + pi_oth[row];
849       pa  = pa_oth + pi_oth[row];
850       /* perform sparse axpy */
851       valtmp = aoa[j];
852       nextp  = 0;
853       for (k=0; nextp<pnz; k++) {
854         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
855           apa_sparse[k] += valtmp*pa[nextp++];
856         }
857       }
858       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
859     }
860 
861     /* set values in C */
862     cdnz = cd->i[i+1] - cd->i[i];
863     conz = co->i[i+1] - co->i[i];
864 
865     /* 1st off-diagoanl part of C */
866     ca = coa + co->i[i];
867     k  = 0;
868     for (k0=0; k0<conz; k0++){
869       if (apJ[k] >= cstart) break;
870       ca[k0]      = apa_sparse[k];
871       apa_sparse[k] = 0.0;
872       k++;
873     }
874 
875     /* diagonal part of C */
876     ca = cda + cd->i[i];
877     for (k1=0; k1<cdnz; k1++){
878       ca[k1]      = apa_sparse[k];
879       apa_sparse[k] = 0.0;
880       k++;
881     }
882 
883     /* 2nd off-diagoanl part of C */
884     ca = coa + co->i[i];
885     for (; k0<conz; k0++){
886       ca[k0]      = apa_sparse[k];
887       apa_sparse[k] = 0.0;
888       k++;
889     }
890   }
891   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
892   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
893   PetscFunctionReturn(0);
894 }
895 
896 /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ(), except using LLCondensed to avoid O(BN) memory requirement */
897 #undef __FUNCT__
898 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable"
899 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable(Mat A,Mat P,PetscReal fill,Mat *C)
900 {
901   PetscErrorCode       ierr;
902   MPI_Comm             comm=((PetscObject)A)->comm;
903   Mat                  Cmpi;
904   Mat_PtAPMPI          *ptap;
905   PetscFreeSpaceList   free_space=PETSC_NULL,current_space=PETSC_NULL;
906   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*c;
907   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
908   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
909   PetscInt             *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
910   PetscInt             i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=0;
911   PetscInt             am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
912   PetscInt             nlnk_max,armax,prmax;
913   PetscReal            afill;
914   PetscScalar          *apa;
915 #if defined(DEBUG_MATMATMULT)
916   PetscMPIInt          rank=a->rank;
917 #endif
918 
919   PetscFunctionBegin;
920 #if defined(DEBUG_MATMATMULT)
921   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
922   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call MatMatMultSymbolic_MPIAIJ_MPIAIJ_Scalable()...\n",rank);
923 #endif
924 
925   /* create struct Mat_PtAPMPI and attached it to C later */
926   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
927 
928   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
929   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
930 #if defined(DEBUG_MATMATMULT)
931   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_oth is done...\n",rank);
932 #endif
933   /* get P_loc by taking all local rows of P */
934   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
935 
936   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
937   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
938   pi_loc = p_loc->i; pj_loc = p_loc->j;
939   pi_oth = p_oth->i; pj_oth = p_oth->j;
940 
941 #if defined(DEBUG_MATMATMULT)
942   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_loc is done, Annz %D * P_locnnz %D = %D...\n",rank,ad->i[am]+ao->i[am],p_loc->rmax,(ad->i[am]+ao->i[am])*p_loc->rmax);
943 #endif
944 
945   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
946   /*-------------------------------------------------------------------*/
947   ierr  = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr);
948   ptap->api = api;
949   api[0]    = 0;
950 
951   /* create and initialize a linked list */
952   armax = ad->rmax+ao->rmax;
953   prmax = PetscMax(p_loc->rmax,p_oth->rmax);
954   nlnk_max = armax*prmax;
955   if (!nlnk_max || nlnk_max > pN) nlnk_max = pN;
956 #if defined(DEBUG_MATMATMULT)
957   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] pN %d; nlnk_max %d; Armax %d+%d=%d; Prmax Max(%d,%d)=%d\n",rank,pN,nlnk_max,ad->rmax,ao->rmax,armax,p_loc->rmax,p_oth->rmax,prmax);
958 #endif
959   ierr = PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);CHKERRQ(ierr);
960 
961   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
962   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
963   current_space = free_space;
964 
965   ierr = MatPreallocateInitialize(comm,am,pn,dnz,onz);CHKERRQ(ierr);
966   for (i=0; i<am; i++) {
967     apnz = 0;
968     /* diagonal portion of A */
969     nzi = adi[i+1] - adi[i];
970     for (j=0; j<nzi; j++){
971       row = *adj++;
972       pnz = pi_loc[row+1] - pi_loc[row];
973       Jptr  = pj_loc + pi_loc[row];
974       /* add non-zero cols of P into the sorted linked list lnk */
975       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
976     }
977     /* off-diagonal portion of A */
978     nzi = aoi[i+1] - aoi[i];
979     for (j=0; j<nzi; j++){
980       row = *aoj++;
981       pnz = pi_oth[row+1] - pi_oth[row];
982       Jptr  = pj_oth + pi_oth[row];
983       ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr);
984     }
985 
986     apnz     = *lnk;
987     api[i+1] = api[i] + apnz;
988     if (apnz > apnz_max) apnz_max = apnz;
989 
990     /* if free space is not available, double the total space in the list */
991     if (current_space->local_remaining<apnz) {
992       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
993       nspacedouble++;
994     }
995 
996     /* Copy data into free space, then initialize lnk */
997     ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr);
998     ierr = MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);CHKERRQ(ierr);
999     current_space->array           += apnz;
1000     current_space->local_used      += apnz;
1001     current_space->local_remaining -= apnz;
1002   }
1003 
1004   /* Allocate space for apj, initialize apj, and */
1005   /* destroy list of free space and other temporary array(s) */
1006   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ptap->apj);CHKERRQ(ierr);
1007   apj  = ptap->apj;
1008   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
1009   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
1010 #if defined(DEBUG_MATMATMULT)
1011   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] AP is done..., apnz_max %d\n",rank,apnz_max);
1012 #endif
1013 
1014   /* create and assemble symbolic parallel matrix Cmpi */
1015   /*----------------------------------------------------*/
1016   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1017   ierr = MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1018   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1019   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1020   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1021   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1022 
1023   /* malloc apa for assembly Cmpi */
1024   ierr = PetscMalloc(apnz_max*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
1025   ierr = PetscMemzero(apa,apnz_max*sizeof(PetscScalar));CHKERRQ(ierr);
1026   ptap->apa = apa;
1027   for (i=0; i<am; i++){
1028     row  = i + rstart;
1029     apnz = api[i+1] - api[i];
1030     ierr = MatSetValues(Cmpi,1,&row,apnz,apj,apa,INSERT_VALUES);CHKERRQ(ierr);
1031     apj += apnz;
1032   }
1033   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1034   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1035 
1036   ptap->destroy             = Cmpi->ops->destroy;
1037   ptap->duplicate           = Cmpi->ops->duplicate;
1038   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_Scalable;
1039   Cmpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
1040   Cmpi->ops->duplicate      = MatDuplicate_MPIAIJ_MatMatMult;
1041 
1042   /* attach the supporting struct to Cmpi for reuse */
1043   c = (Mat_MPIAIJ*)Cmpi->data;
1044   c->ptap  = ptap;
1045 
1046   *C = Cmpi;
1047 
1048   /* set MatInfo */
1049   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]) + 1.e-5;
1050   if (afill < 1.0) afill = 1.0;
1051   Cmpi->info.mallocs           = nspacedouble;
1052   Cmpi->info.fill_ratio_given  = fill;
1053   Cmpi->info.fill_ratio_needed = afill;
1054 
1055 #if defined(PETSC_USE_INFO)
1056   if (api[am]) {
1057     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
1058     ierr = PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
1059   } else {
1060     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1061   }
1062 #endif
1063   PetscFunctionReturn(0);
1064 }
1065 
1066 /*-------------------------------------------------------------------------*/
1067 #undef __FUNCT__
1068 #define __FUNCT__ "MatTransposeMatMult_MPIAIJ_MPIAIJ"
1069 PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C)
1070 {
1071   PetscErrorCode ierr;
1072 
1073   PetscFunctionBegin;
1074   if (scall == MAT_INITIAL_MATRIX){
1075     ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);CHKERRQ(ierr);
1076   }
1077   ierr = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(P,A,*C);CHKERRQ(ierr);
1078   PetscFunctionReturn(0);
1079 }
1080 
1081 #undef __FUNCT__
1082 #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ"
1083 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C)
1084 {
1085   PetscErrorCode       ierr;
1086   Mat_Merge_SeqsToMPI  *merge;
1087   Mat_MPIAIJ           *p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1088   Mat_SeqAIJ           *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1089   Mat_PtAPMPI          *ptap;
1090   PetscInt             *adj,*apJ;
1091   PetscInt             i,j,k,anz,pnz,row,*cj;
1092   MatScalar            *ada,*apa,*ca,valtmp;
1093   PetscInt             am=A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1094   MPI_Comm             comm=((PetscObject)C)->comm;
1095   PetscMPIInt          size,rank,taga,*len_s;
1096   PetscInt             *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1097   PetscInt             **buf_ri,**buf_rj;
1098   PetscInt             cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */
1099   MPI_Request          *s_waits,*r_waits;
1100   MPI_Status           *status;
1101   MatScalar            **abuf_r,*ba_i,*pA,*coa,*ba;
1102   PetscInt             *api,*apj,*coi,*coj;
1103   PetscInt             *poJ=po->j,*pdJ=pd->j;
1104   PetscInt             sparse_axpy;
1105   Mat                  A_loc;
1106   Mat_SeqAIJ           *a_loc;
1107 
1108   PetscFunctionBegin;
1109   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1110   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1111   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
1112 
1113   ptap  = c->ptap;
1114   merge = ptap->merge;
1115 
1116   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
1117   /*--------------------------------------------------------------*/
1118   /* get data from symbolic products */
1119   coi = merge->coi; coj = merge->coj;
1120   ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr);
1121   ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr);
1122 
1123   bi     = merge->bi; bj = merge->bj;
1124   owners = merge->rowmap->range;
1125   ierr   = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
1126   ierr   = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr);
1127 
1128   /* get A_loc by taking all local rows of A */
1129   A_loc = ptap->A_loc;
1130   ierr = MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);CHKERRQ(ierr);
1131   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1132   api   = a_loc->i;
1133   apj   = a_loc->j;
1134 
1135   ierr = PetscMalloc((A->cmap->N)*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
1136   ierr = PetscMemzero(apa,A->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
1137 
1138     for (i=0; i<am; i++) {
1139       /* 2-a) put A[i,:] to dense array apa */
1140       anz = api[i+1] - api[i];
1141       adj = apj + api[i];
1142       ada = a_loc->a + api[i];
1143       for (j=0; j<anz; j++){
1144         apa[adj[j]] = ada[j];
1145       }
1146 
1147       /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1148       /*--------------------------------------------------------------*/
1149       /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1150       pnz = po->i[i+1] - po->i[i];
1151       poJ = po->j + po->i[i];
1152       pA  = po->a + po->i[i];
1153       for (j=0; j<pnz; j++){
1154         row = poJ[j];
1155         cnz = coi[row+1] - coi[row];
1156         cj  = coj + coi[row];
1157         ca  = coa + coi[row];
1158         /* perform dense axpy */
1159         valtmp = pA[j];
1160         for (k=0; k<cnz; k++) {
1161           ca[k] += valtmp*apa[cj[k]];
1162         }
1163         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1164       }
1165 
1166       /* put the value into Cd (diagonal part) */
1167       pnz = pd->i[i+1] - pd->i[i];
1168       pdJ = pd->j + pd->i[i];
1169       pA  = pd->a + pd->i[i];
1170       for (j=0; j<pnz; j++){
1171         row = pdJ[j];
1172         cnz = bi[row+1] - bi[row];
1173         cj  = bj + bi[row];
1174         ca  = ba + bi[row];
1175         /* perform dense axpy */
1176         valtmp = pA[j];
1177         for (k=0; k<cnz; k++) {
1178           ca[k] += valtmp*apa[cj[k]];
1179         }
1180         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1181       }
1182 
1183       /* zero the current row of A*P */
1184       apJ = apj + api[i];
1185       for (k=0; k<anz; k++) apa[apJ[k]] = 0.0;
1186     }
1187 
1188   /* 3) send and recv matrix values coa */
1189   /*------------------------------------*/
1190   buf_ri = merge->buf_ri;
1191   buf_rj = merge->buf_rj;
1192   len_s  = merge->len_s;
1193   ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1194   ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1195 
1196   ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr);
1197   for (proc=0,k=0; proc<size; proc++){
1198     if (!len_s[proc]) continue;
1199     i = merge->owners_co[proc];
1200     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1201     k++;
1202   }
1203   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1204   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1205 
1206   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1207   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1208   ierr = PetscFree(coa);CHKERRQ(ierr);
1209 
1210   /* 4) insert local Cseq and received values into Cmpi */
1211   /*----------------------------------------------------*/
1212   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
1213   for (k=0; k<merge->nrecv; k++){
1214     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure  -- memory corruption ???*/
1215     nrows       = *(buf_ri_k[k]);
1216     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1217     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
1218   }
1219 
1220   for (i=0; i<cm; i++) {
1221     row = owners[rank] + i; /* global row index of C_seq */
1222     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1223     ba_i = ba + bi[i];
1224     bnz  = bi[i+1] - bi[i];
1225     /* add received vals into ba */
1226     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
1227       /* i-th row */
1228       if (i == *nextrow[k]) {
1229         cnz = *(nextci[k]+1) - *nextci[k];
1230         cj  = buf_rj[k] + *(nextci[k]);
1231         ca  = abuf_r[k] + *(nextci[k]);
1232         nextcj = 0;
1233         for (j=0; nextcj<cnz; j++){
1234           if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */
1235             ba_i[j] += ca[nextcj++];
1236           }
1237         }
1238         nextrow[k]++; nextci[k]++;
1239         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1240       }
1241     }
1242     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1243   }
1244   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1245   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1246 
1247   ierr = PetscFree(ba);CHKERRQ(ierr);
1248   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1249   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1250   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1251   ierr = PetscFree(apa);CHKERRQ(ierr);
1252   PetscFunctionReturn(0);
1253 }
1254 
1255 /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1256 #undef __FUNCT__
1257 #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ"
1258 PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C)
1259 {
1260   PetscErrorCode       ierr;
1261   Mat                  Cmpi;
1262   Mat_PtAPMPI          *ptap;
1263   PetscFreeSpaceList   free_space=PETSC_NULL,current_space=PETSC_NULL;
1264   Mat_MPIAIJ           *p=(Mat_MPIAIJ*)P->data,*c;
1265   PetscInt             *pdti,*pdtj,*poti,*potj,*ptJ;
1266   PetscInt             nnz;
1267   PetscInt             nlnk,*lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1268   PetscInt             am=A->rmap->n,pn=P->cmap->n;
1269   PetscBT              lnkbt;
1270   MPI_Comm             comm=((PetscObject)A)->comm;
1271   PetscMPIInt          size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1272   PetscInt             **buf_rj,**buf_ri,**buf_ri_k;
1273   PetscInt             len,proc,*dnz,*onz,*owners;
1274   PetscInt             nzi,*bi,*bj;
1275   PetscInt             nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1276   MPI_Request          *swaits,*rwaits;
1277   MPI_Status           *sstatus,rstatus;
1278   Mat_Merge_SeqsToMPI  *merge;
1279   PetscInt             *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0;
1280   PetscReal            afill=1.0,afill_tmp;
1281   PetscInt             rstart = P->cmap->rstart,rmax,aN=A->cmap->N;
1282   PetscScalar          *vals;
1283   Mat                  A_loc;
1284   Mat_SeqAIJ           *a_loc;
1285 
1286   PetscFunctionBegin;
1287   /* check if matrix local sizes are compatible */
1288   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend){
1289     SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A (%D, %D) != P (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend);
1290   }
1291 
1292   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1293   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1294 #if defined(DEBUG_MATTrMatMult)
1295   ierr = PetscSynchronizedPrintf(comm,"[%d] TransposeMatMultSymbolic P: %d %d, %d %d; A %d %d, %d %d\n",rank,P->rmap->N,pN,P->rmap->n,P->cmap->n,A->rmap->N,aN,A->rmap->n,A->cmap->n);
1296   ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr);
1297 #endif
1298 
1299   /* create struct Mat_PtAPMPI and attached it to C later */
1300   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
1301 
1302   /* get A_loc by taking all local rows of A */
1303   ierr = MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);CHKERRQ(ierr);
1304   ptap->A_loc = A_loc;
1305   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1306   api   = a_loc->i;
1307   apj   = a_loc->j;
1308 
1309   /* determine symbolic Co=(p->B)^T*A - send to others */
1310   /*----------------------------------------------------*/
1311   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
1312 
1313   /* then, compute symbolic Co = (p->B)^T*A */
1314   pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1315                          >= (num of nonzero rows of C_seq) - pn */
1316   ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr);
1317   coi[0] = 0;
1318 
1319   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1320   nnz           = fill*(poti[pon] + api[am]);
1321   ierr          = PetscFreeSpaceGet(nnz,&free_space);
1322   current_space = free_space;
1323 #if defined(DEBUG_MATTrMatMult)
1324   ierr = PetscSynchronizedPrintf(comm, "[%d] nnz = fill %g *(%d + %d)\n",rank,fill,poti[pon],api[am]);CHKERRQ(ierr);
1325   ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr);
1326 #endif
1327   /* create and initialize a linked list */
1328   nlnk = aN+1;
1329   ierr = PetscLLCreate(aN,aN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1330 
1331   for (i=0; i<pon; i++) {
1332     nnz = 0;
1333     pnz = poti[i+1] - poti[i];
1334     ptJ = potj + poti[i];
1335     for (j=0; j<pnz; j++){
1336       row  = ptJ[j]; /* row of A_loc == col of Pot */
1337       apnz = api[row+1] - api[row];
1338       Jptr = apj + api[row];
1339       /* add non-zero cols of AP into the sorted linked list lnk */
1340       ierr = PetscLLAddSorted(apnz,Jptr,aN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1341       nnz += nlnk;
1342     }
1343 
1344     /* If free space is not available, double the total space in the list */
1345     if (current_space->local_remaining<nnz) {
1346       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1347       nspacedouble++;
1348     }
1349 
1350     /* Copy data into free space, and zero out denserows */
1351     ierr = PetscLLClean(aN,aN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
1352     current_space->array           += nnz;
1353     current_space->local_used      += nnz;
1354     current_space->local_remaining -= nnz;
1355     coi[i+1] = coi[i] + nnz;
1356   }
1357 
1358   ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr);
1359   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
1360   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]);
1361   if (afill_tmp > afill) afill = afill_tmp;
1362 
1363   /* send j-array (coj) of Co to other processors */
1364   /*----------------------------------------------*/
1365   /* determine row ownership */
1366   ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr);
1367   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
1368   merge->rowmap->n = pn;
1369   merge->rowmap->bs = 1;
1370   ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
1371   owners = merge->rowmap->range;
1372 
1373   /* determine the number of messages to send, their lengths */
1374   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr);
1375   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1376   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr);
1377   len_s = merge->len_s;
1378   merge->nsend = 0;
1379 
1380   ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr);
1381   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1382 
1383   proc = 0;
1384   for (i=0; i<pon; i++){
1385     while (prmap[i] >= owners[proc+1]) proc++;
1386     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1387     len_s[proc] += coi[i+1] - coi[i];
1388   }
1389 
1390   len   = 0;  /* max length of buf_si[] */
1391   owners_co[0] = 0;
1392   for (proc=0; proc<size; proc++){
1393     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1394     if (len_si[proc]){
1395       merge->nsend++;
1396       len_si[proc] = 2*(len_si[proc] + 1);
1397       len += len_si[proc];
1398     }
1399   }
1400 
1401   /* determine the number and length of messages to receive for coi and coj  */
1402   ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
1403   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
1404 
1405   /* post the Irecv and Isend of coj */
1406   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1407   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1408   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr);
1409   for (proc=0, k=0; proc<size; proc++){
1410     if (!len_s[proc]) continue;
1411     i = owners_co[proc];
1412     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
1413     k++;
1414   }
1415 
1416   /* receives and sends of coj are complete */
1417   ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr);
1418   for (i=0; i<merge->nrecv; i++){
1419     ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr);
1420   }
1421   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1422   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1423 
1424   /* send and recv coi */
1425   /*-------------------*/
1426   ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1427   ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1428   ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr);
1429   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1430   for (proc=0,k=0; proc<size; proc++){
1431     if (!len_s[proc]) continue;
1432     /* form outgoing message for i-structure:
1433          buf_si[0]:                 nrows to be sent
1434                [1:nrows]:           row index (global)
1435                [nrows+1:2*nrows+1]: i-structure index
1436     */
1437     /*-------------------------------------------*/
1438     nrows = len_si[proc]/2 - 1;
1439     buf_si_i    = buf_si + nrows+1;
1440     buf_si[0]   = nrows;
1441     buf_si_i[0] = 0;
1442     nrows = 0;
1443     for (i=owners_co[proc]; i<owners_co[proc+1]; i++){
1444       nzi = coi[i+1] - coi[i];
1445       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
1446       buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */
1447       nrows++;
1448     }
1449     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
1450     k++;
1451     buf_si += len_si[proc];
1452   }
1453   i = merge->nrecv;
1454   while (i--) {
1455     ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr);
1456   }
1457   ierr = PetscFree(rwaits);CHKERRQ(ierr);
1458   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1459   ierr = PetscFree(len_si);CHKERRQ(ierr);
1460   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1461   ierr = PetscFree(swaits);CHKERRQ(ierr);
1462   ierr = PetscFree(sstatus);CHKERRQ(ierr);
1463   ierr = PetscFree(buf_s);CHKERRQ(ierr);
1464 
1465   /* compute the local portion of C (mpi mat) */
1466   /*------------------------------------------*/
1467   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
1468 
1469   /* allocate bi array and free space for accumulating nonzero column info */
1470   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
1471   bi[0] = 0;
1472 
1473   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1474   nnz           = fill*(pdti[pn] + poti[pon] + api[am]);
1475   ierr          = PetscFreeSpaceGet(nnz,&free_space);
1476   current_space = free_space;
1477 
1478 #if defined(DEBUG_MATTrMatMult)
1479   ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] nnz=%d + %d + %d\n",rank,pdti[pn],poti[pon],api[am]);
1480 #endif
1481 
1482   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
1483   for (k=0; k<merge->nrecv; k++){
1484     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1485     nrows       = *buf_ri_k[k];
1486     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1487     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
1488   }
1489 
1490   ierr = MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);CHKERRQ(ierr);
1491   rmax = 0;
1492   for (i=0; i<pn; i++) {
1493     /* add pdt[i,:]*AP into lnk */
1494     nnz = 0;
1495     pnz = pdti[i+1] - pdti[i];
1496     ptJ = pdtj + pdti[i];
1497     for (j=0; j<pnz; j++){
1498       row  = ptJ[j];  /* row of AP == col of Pt */
1499       apnz = api[row+1] - api[row];
1500       Jptr = apj + api[row];
1501       /* add non-zero cols of AP into the sorted linked list lnk */
1502       ierr = PetscLLAddSorted(apnz,Jptr,aN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1503       nnz += nlnk;
1504     }
1505     /* add received col data into lnk */
1506     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
1507       if (i == *nextrow[k]) { /* i-th row */
1508         nzi = *(nextci[k]+1) - *nextci[k];
1509         Jptr  = buf_rj[k] + *nextci[k];
1510         ierr = PetscLLAddSorted(nzi,Jptr,aN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1511         nnz += nlnk;
1512         nextrow[k]++; nextci[k]++;
1513       }
1514     }
1515 
1516     /* if free space is not available, make more free space */
1517     if (current_space->local_remaining<nnz) {
1518       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
1519       nspacedouble++;
1520     }
1521     /* copy data into free space, then initialize lnk */
1522     ierr = PetscLLClean(aN,aN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
1523     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
1524     current_space->array           += nnz;
1525     current_space->local_used      += nnz;
1526     current_space->local_remaining -= nnz;
1527     bi[i+1] = bi[i] + nnz;
1528     if (nnz > rmax) rmax = nnz;
1529   }
1530   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1531 
1532   ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1533   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
1534   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + api[am]);
1535   if (afill_tmp > afill) afill = afill_tmp;
1536   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1537   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
1538   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
1539 
1540   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
1541   /*----------------------------------------------------------------------------------*/
1542   ierr = PetscMalloc((rmax+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
1543   ierr = PetscMemzero(vals,rmax*sizeof(PetscScalar));CHKERRQ(ierr);
1544 
1545   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
1546   ierr = MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1547   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
1548   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
1549   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1550   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
1551   for (i=0; i<pn; i++){
1552     row = i + rstart;
1553     nnz = bi[i+1] - bi[i];
1554     Jptr = bj + bi[i];
1555     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
1556   }
1557   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1558   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1559   ierr = PetscFree(vals);CHKERRQ(ierr);
1560 
1561   merge->bi            = bi;
1562   merge->bj            = bj;
1563   merge->coi           = coi;
1564   merge->coj           = coj;
1565   merge->buf_ri        = buf_ri;
1566   merge->buf_rj        = buf_rj;
1567   merge->owners_co     = owners_co;
1568   merge->destroy       = Cmpi->ops->destroy;
1569   merge->duplicate     = Cmpi->ops->duplicate;
1570 
1571   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_PtAP;
1572 
1573   /* attach the supporting struct to Cmpi for reuse */
1574   c = (Mat_MPIAIJ*)Cmpi->data;
1575   c->ptap        = ptap;
1576   ptap->api      = PETSC_NULL;
1577   ptap->apj      = PETSC_NULL;
1578   ptap->merge    = merge;
1579   ptap->apnz_max = ap_rmax;
1580 
1581   *C = Cmpi;
1582 #if defined(PETSC_USE_INFO)
1583   if (bi[pn] != 0) {
1584     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
1585     ierr = PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
1586   } else {
1587     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1588   }
1589 #endif
1590   PetscFunctionReturn(0);
1591 }
1592