xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision 78d30b94ba0e23c9ea405755060bbf738523281e)
1 
2 /*
3   Defines projective product routines where A is a MPIAIJ matrix
4           C = P^T * A * P
5 */
6 
7 #include <../src/mat/impls/aij/seq/aij.h>   /*I "petscmat.h" I*/
8 #include <../src/mat/utils/freespace.h>
9 #include <../src/mat/impls/aij/mpi/mpiaij.h>
10 #include <petscbt.h>
11 #include <petsctime.h>
12 
13 /* #define PTAP_PROFILE */
14 
15 extern PetscErrorCode MatDestroy_MPIAIJ(Mat);
16 #undef __FUNCT__
17 #define __FUNCT__ "MatDestroy_MPIAIJ_PtAP"
18 PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A)
19 {
20   PetscErrorCode ierr;
21   Mat_MPIAIJ     *a=(Mat_MPIAIJ*)A->data;
22   Mat_PtAPMPI    *ptap=a->ptap;
23 
24   PetscFunctionBegin;
25   if (ptap) {
26     Mat_Merge_SeqsToMPI *merge=ptap->merge;
27     ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
28     ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
29     ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
30     ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
31     ierr = MatDestroy(&ptap->A_loc);CHKERRQ(ierr); /* used by MatTransposeMatMult() */
32     ierr = MatDestroy(&ptap->Rd);CHKERRQ(ierr);
33     ierr = MatDestroy(&ptap->Ro);CHKERRQ(ierr);
34     if (ptap->AP_loc) { /* used by alg_rap */
35       Mat_SeqAIJ *ap = (Mat_SeqAIJ*)(ptap->AP_loc)->data;
36       ierr = PetscFree(ap->i);CHKERRQ(ierr);
37       ierr = PetscFree2(ap->j,ap->a);CHKERRQ(ierr);
38       ierr = MatDestroy(&ptap->AP_loc);CHKERRQ(ierr);
39     } else { /* used by alg_ptap */
40       ierr = PetscFree(ptap->api);CHKERRQ(ierr);
41       ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
42     }
43     ierr = MatDestroy(&ptap->C_loc);CHKERRQ(ierr);
44     ierr = MatDestroy(&ptap->C_oth);CHKERRQ(ierr);
45     if (ptap->apa) {ierr = PetscFree(ptap->apa);CHKERRQ(ierr);}
46     if (merge) { /* used by alg_ptap */
47       ierr = PetscFree(merge->id_r);CHKERRQ(ierr);
48       ierr = PetscFree(merge->len_s);CHKERRQ(ierr);
49       ierr = PetscFree(merge->len_r);CHKERRQ(ierr);
50       ierr = PetscFree(merge->bi);CHKERRQ(ierr);
51       ierr = PetscFree(merge->bj);CHKERRQ(ierr);
52       ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr);
53       ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr);
54       ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr);
55       ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr);
56       ierr = PetscFree(merge->coi);CHKERRQ(ierr);
57       ierr = PetscFree(merge->coj);CHKERRQ(ierr);
58       ierr = PetscFree(merge->owners_co);CHKERRQ(ierr);
59       ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr);
60       ierr = merge->destroy(A);CHKERRQ(ierr);
61       ierr = PetscFree(ptap->merge);CHKERRQ(ierr);
62     } else {
63       ierr = MatDestroy_MPIAIJ(A);CHKERRQ(ierr);
64     }
65     ierr = PetscFree(ptap);CHKERRQ(ierr);
66   }
67   PetscFunctionReturn(0);
68 }
69 
70 #undef __FUNCT__
71 #define __FUNCT__ "MatDuplicate_MPIAIJ_MatPtAP"
72 PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M)
73 {
74   PetscErrorCode ierr;
75   Mat_MPIAIJ     *a     = (Mat_MPIAIJ*)A->data;
76   Mat_PtAPMPI    *ptap  = a->ptap;
77 
78   PetscFunctionBegin;
79   ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr);
80   (*M)->ops->destroy   = ptap->destroy;
81   (*M)->ops->duplicate = ptap->duplicate;
82   PetscFunctionReturn(0);
83 }
84 
85 #undef __FUNCT__
86 #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ"
87 PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
88 {
89   PetscErrorCode ierr;
90   PetscBool      rap=PETSC_TRUE; /* do R=P^T locally, then C=R*A*P */
91   MPI_Comm       comm;
92 
93   PetscFunctionBegin;
94   /* check if matrix local sizes are compatible */
95   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
96   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
97     SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, Arow (%D, %D) != Prow (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend);
98   }
99   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) {
100     SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, Acol (%D, %D) != Prow (%D,%D)",A->cmap->rstart,A->cmap->rend,P->rmap->rstart,P->rmap->rend);
101   }
102 
103   ierr = PetscOptionsGetBool(NULL,"-matrap",&rap,NULL);CHKERRQ(ierr);
104   if (scall == MAT_INITIAL_MATRIX) {
105     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
106     if (rap) { /* do R=P^T locally, then C=R*A*P */
107       ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);
108     } else {       /* do P^T*A*P */
109       ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(A,P,fill,C);CHKERRQ(ierr);
110     }
111     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
112   }
113   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
114   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
115   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
116   PetscFunctionReturn(0);
117 }
118 
119 /* Does not requires array of size pN, thus scalable */
120 #undef __FUNCT__
121 #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ"
122 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
123 {
124   PetscErrorCode      ierr;
125   Mat_PtAPMPI         *ptap;
126   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
127   MPI_Comm            comm;
128   PetscMPIInt         size,rank;
129   Mat                 Cmpi;
130   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
131   PetscInt            am=A->rmap->n,pm=P->rmap->n;
132   PetscInt            *lnk,i,k,pnz,row,nnz;
133   PetscInt            pN=P->cmap->N,pn=P->cmap->n;
134   PetscBT             lnkbt;
135   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
136   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
137   PetscInt            len,proc,*dnz,*onz,*owners;
138   PetscInt            nzi,nspacedouble;
139   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
140   MPI_Request         *swaits,*rwaits;
141   MPI_Status          *sstatus,rstatus;
142   PetscLayout         rowmap;
143   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
144   PetscInt            nsend;
145   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
146   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0,Crmax;
147   PetscInt            rmax,*aj,*ai,*pi;
148   Mat_SeqAIJ          *p_loc,*p_oth,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*c_loc,*c_oth;
149   PetscScalar         *apv;
150   PetscTable          ta;
151   const char          *algTypes[2] = {"scalable","nonscalable"};
152   PetscInt            alg=1; /* set default algorithm */
153 #if defined(PETSC_USE_INFO)
154   PetscReal           apfill;
155 #endif
156 #if defined(PTAP_PROFILE)
157   PetscLogDouble      t0,t1,t11,t12,t2,t3,t4;
158 #endif
159 
160   PetscFunctionBegin;
161   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
162   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
163   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
164 #if defined(PTAP_PROFILE)
165   ierr = PetscTime(&t0);CHKERRQ(ierr);
166 #endif
167 
168   /* create struct Mat_PtAPMPI and attached it to C later */
169   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
170   ptap->reuse = MAT_INITIAL_MATRIX;
171 
172   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
173   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
174   /* get P_loc by taking all local rows of P */
175   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
176 
177   /* (0) compute Rd = Pd^T, Ro = Po^T  */
178   /* --------------------------------- */
179   ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
180   ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
181 #if defined(PTAP_PROFILE)
182   ierr = PetscTime(&t11);CHKERRQ(ierr);
183 #endif
184 
185   /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */
186   /* ----------------------------------------------------------------- */
187   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
188   p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
189 
190   /* create and initialize a linked list */
191   Crmax = p_loc->rmax+p_oth->rmax;
192   ierr = PetscTableCreate(2*Crmax,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */
193 
194   for (row=0; row<ptap->P_loc->rmap->N; row++) {
195     nzi = p_loc->i[row+1] - p_loc->i[row];
196     for (j=0; j<nzi; j++) {
197       Jptr = j + p_loc->j + p_loc->i[row];
198       ierr = PetscTableAdd(ta,*Jptr+1,1,INSERT_VALUES);CHKERRQ(ierr); /* key must be >0 */
199     }
200   }
201   for (row=0; row<ptap->P_oth->rmap->N; row++) {
202     nzi = p_oth->i[row+1] - p_oth->i[row];
203     for (j=0; j<nzi; j++) {
204       Jptr = j + p_oth->j + p_oth->i[row];
205       ierr = PetscTableAdd(ta,*Jptr+1,1,INSERT_VALUES);CHKERRQ(ierr); /* key must be >0 */
206     }
207   }
208   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */
209 
210   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
211 
212   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
213   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ad->i[am]+ao->i[am]+p_loc->i[pm])),&free_space);CHKERRQ(ierr);
214   current_space = free_space;
215   nspacedouble  = 0;
216 
217   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
218   api[0] = 0;
219   for (i=0; i<am; i++) {
220     /* diagonal portion: Ad[i,:]*P */
221     ai = ad->i; pi = p_loc->i;
222     nzi = ai[i+1] - ai[i];
223     aj  = ad->j + ai[i];
224     for (j=0; j<nzi; j++) {
225       row  = aj[j];
226       pnz  = pi[row+1] - pi[row];
227       Jptr = p_loc->j + pi[row];
228       /* add non-zero cols of P into the sorted linked list lnk */
229       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
230     }
231     /* off-diagonal portion: Ao[i,:]*P */
232     ai = ao->i; pi = p_oth->i;
233     nzi = ai[i+1] - ai[i];
234     aj  = ao->j + ai[i];
235     for (j=0; j<nzi; j++) {
236       row  = aj[j];
237       pnz  = pi[row+1] - pi[row];
238       Jptr = p_oth->j + pi[row];
239       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
240     }
241     apnz     = lnk[0];
242     api[i+1] = api[i] + apnz;
243     if (ap_rmax < apnz) ap_rmax = apnz;
244 
245     /* if free space is not available, double the total space in the list */
246     if (current_space->local_remaining<apnz) {
247       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
248       nspacedouble++;
249     }
250 
251     /* Copy data into free space, then initialize lnk */
252     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
253 
254     current_space->array           += apnz;
255     current_space->local_used      += apnz;
256     current_space->local_remaining -= apnz;
257   }
258   /* Allocate space for apj and apv, initialize apj, and */
259   /* destroy list of free space and other temporary array(s) */
260   ierr   = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr);
261   ierr   = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
262   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
263 
264   /* Create AP_loc for reuse */
265   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr);
266 
267 #if defined(PETSC_USE_INFO)
268   apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1);
269   ptap->AP_loc->info.mallocs           = nspacedouble;
270   ptap->AP_loc->info.fill_ratio_given  = fill;
271   ptap->AP_loc->info.fill_ratio_needed = apfill;
272 
273   if (api[am]) {
274     ierr = PetscInfo3(ptap->AP_loc,"AP_loc reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)apfill);CHKERRQ(ierr);
275     ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr);
276   } else {
277     ierr = PetscInfo(ptap->AP_loc,"AP_loc is empty \n");CHKERRQ(ierr);
278   }
279 #endif
280 
281 #if defined(PTAP_PROFILE)
282   ierr = PetscTime(&t12);CHKERRQ(ierr);
283 #endif
284 
285   /* (2-1) compute symbolic Co = Ro*AP_loc  */
286   /* ------------------------------------ */
287   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
288   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
289 #if defined(PTAP_PROFILE)
290   ierr = PetscTime(&t1);CHKERRQ(ierr);
291 #endif
292 
293   /* (3) send coj of C_oth to other processors  */
294   /* ------------------------------------------ */
295   /* determine row ownership */
296   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
297   rowmap->n  = pn;
298   rowmap->bs = 1;
299   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
300   owners = rowmap->range;
301 
302   /* determine the number of messages to send, their lengths */
303   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
304   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
305   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
306 
307   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
308   coi   = c_oth->i; coj = c_oth->j;
309   con   = ptap->C_oth->rmap->n;
310   proc  = 0;
311   for (i=0; i<con; i++) {
312     while (prmap[i] >= owners[proc+1]) proc++;
313     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
314     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
315   }
316 
317   len          = 0; /* max length of buf_si[], see (4) */
318   owners_co[0] = 0;
319   nsend        = 0;
320   for (proc=0; proc<size; proc++) {
321     owners_co[proc+1] = owners_co[proc] + len_si[proc];
322     if (len_s[proc]) {
323       nsend++;
324       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
325       len         += len_si[proc];
326     }
327   }
328 
329   /* determine the number and length of messages to receive for coi and coj  */
330   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
331   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
332 
333   /* post the Irecv and Isend of coj */
334   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
335   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
336   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
337   for (proc=0, k=0; proc<size; proc++) {
338     if (!len_s[proc]) continue;
339     i    = owners_co[proc];
340     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
341     k++;
342   }
343 
344   /* (2-2) compute symbolic C_loc = Rd*AP_loc */
345   /* ---------------------------------------- */
346   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr);
347   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
348 
349   /* receives coj are complete */
350   for (i=0; i<nrecv; i++) {
351     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
352   }
353   ierr = PetscFree(rwaits);CHKERRQ(ierr);
354   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
355 
356   /* add received column indices into ta to update Crmax */
357   for (k=0; k<nrecv; k++) {/* k-th received message */
358     Jptr = buf_rj[k];
359     for (j=0; j<len_r[k]; j++) {
360       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
361     }
362   }
363   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
364   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
365 
366   /* (4) send and recv coi */
367   /*-----------------------*/
368   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
369   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
370   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
371   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
372   for (proc=0,k=0; proc<size; proc++) {
373     if (!len_s[proc]) continue;
374     /* form outgoing message for i-structure:
375          buf_si[0]:                 nrows to be sent
376                [1:nrows]:           row index (global)
377                [nrows+1:2*nrows+1]: i-structure index
378     */
379     /*-------------------------------------------*/
380     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
381     buf_si_i    = buf_si + nrows+1;
382     buf_si[0]   = nrows;
383     buf_si_i[0] = 0;
384     nrows       = 0;
385     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
386       nzi = coi[i+1] - coi[i];
387       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
388       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
389       nrows++;
390     }
391     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
392     k++;
393     buf_si += len_si[proc];
394   }
395   for (i=0; i<nrecv; i++) {
396     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
397   }
398   ierr = PetscFree(rwaits);CHKERRQ(ierr);
399   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
400 
401   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
402   ierr = PetscFree(len_ri);CHKERRQ(ierr);
403   ierr = PetscFree(swaits);CHKERRQ(ierr);
404   ierr = PetscFree(buf_s);CHKERRQ(ierr);
405 #if defined(PTAP_PROFILE)
406   ierr = PetscTime(&t2);CHKERRQ(ierr);
407 #endif
408   /* (5) compute the local portion of Cmpi      */
409   /* ------------------------------------------ */
410   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
411   ierr          = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr);
412   current_space = free_space;
413 
414   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
415   for (k=0; k<nrecv; k++) {
416     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
417     nrows       = *buf_ri_k[k];
418     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
419     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
420   }
421 
422   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
423   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
424   rmax = 0;
425   for (i=0; i<pn; i++) {
426     /* add C_loc into Cmpi */
427     nzi  = c_loc->i[i+1] - c_loc->i[i];
428     Jptr = c_loc->j + c_loc->i[i];
429     ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
430 
431     /* add received col data into lnk */
432     for (k=0; k<nrecv; k++) { /* k-th received message */
433       if (i == *nextrow[k]) { /* i-th row */
434         nzi  = *(nextci[k]+1) - *nextci[k];
435         Jptr = buf_rj[k] + *nextci[k];
436         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
437         nextrow[k]++; nextci[k]++;
438       }
439     }
440     nnz = lnk[0];
441 
442     /* copy data into free space, then initialize lnk */
443     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
444     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
445     if (nnz > rmax) rmax = nnz;
446   }
447   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
448   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
449   ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr);
450 #if defined(PTAP_PROFILE)
451   ierr = PetscTime(&t3);CHKERRQ(ierr);
452 #endif
453 
454   /* (6) create symbolic parallel matrix Cmpi */
455   /*------------------------------------------*/
456   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
457   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
458   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
459   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
460   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
461   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
462 
463   /* members in merge */
464   ierr = PetscFree(id_r);CHKERRQ(ierr);
465   ierr = PetscFree(len_r);CHKERRQ(ierr);
466   ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr);
467   ierr = PetscFree(buf_ri);CHKERRQ(ierr);
468   ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr);
469   ierr = PetscFree(buf_rj);CHKERRQ(ierr);
470   ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr);
471 
472   /* attach the supporting struct to Cmpi for reuse */
473   c = (Mat_MPIAIJ*)Cmpi->data;
474   c->ptap         = ptap;
475   ptap->rmax      = ap_rmax;
476   ptap->duplicate = Cmpi->ops->duplicate;
477   ptap->destroy   = Cmpi->ops->destroy;
478 
479   ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
480   ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,2,algTypes[1],&alg,NULL);CHKERRQ(ierr);
481   ierr = PetscOptionsEnd();CHKERRQ(ierr);
482 
483   if (alg == 1) {
484     /* Do dense axpy in MatPtAPNumeric_MPIAIJ_MPIAIJ() */
485     ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
486     Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ;
487   } else {
488     SETERRQ(comm,PETSC_ERR_ARG_SIZ,"Not done yet");
489   }
490 
491 
492   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
493   Cmpi->assembled        = PETSC_FALSE;
494   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
495   Cmpi->ops->duplicate   = MatDuplicate_MPIAIJ_MatPtAP;
496   *C                     = Cmpi;
497 
498 #if defined(PTAP_PROFILE)
499   ierr = PetscTime(&t4);CHKERRQ(ierr);
500   if (rank == 1) {
501     printf("PtAPSym: %g + %g + %g + %g + %g + %g = %g\n",t11-t0,t1-t11,t12-t11,t2-t2,t3-t2,t4-t3,t4-t0);
502   }
503 #endif
504   PetscFunctionReturn(0);
505 }
506 
507 #undef __FUNCT__
508 #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ"
509 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
510 {
511   PetscErrorCode    ierr;
512   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
513   Mat_SeqAIJ        *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
514   Mat_SeqAIJ        *ap,*p_loc,*p_oth,*c_seq;
515   Mat_PtAPMPI       *ptap = c->ptap;
516   Mat               AP_loc,C_loc,C_oth;
517   PetscInt          i,rstart,rend,cm,ncols,row;
518   PetscInt          *api,*apj,am = A->rmap->n,j,col,apnz;
519   PetscScalar       *apa;
520   const PetscInt    *cols;
521   const PetscScalar *vals;
522 #if defined(PTAP_PROFILE)
523   PetscMPIInt       rank;
524   MPI_Comm          comm;
525   PetscLogDouble    t0,t1,t2,t3,t4,eR,eAP,eCseq,eCmpi;
526 #endif
527 
528   PetscFunctionBegin;
529   ierr = MatZeroEntries(C);CHKERRQ(ierr);
530 
531   /* 1) get R = Pd^T,Ro = Po^T */
532 #if defined(PTAP_PROFILE)
533   ierr = PetscTime(&t0);CHKERRQ(ierr);
534 #endif
535   if (ptap->reuse == MAT_REUSE_MATRIX) {
536     ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
537     ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
538   }
539 #if defined(PTAP_PROFILE)
540   ierr = PetscTime(&t1);CHKERRQ(ierr);
541   eR = t1 - t0;
542 #endif
543 
544   /* 2) get AP_loc */
545   AP_loc = ptap->AP_loc;
546   ap = (Mat_SeqAIJ*)AP_loc->data;
547 
548   /* 2-1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
549   /*-----------------------------------------------------*/
550   if (ptap->reuse == MAT_REUSE_MATRIX) {
551     /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */
552     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
553     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
554   }
555 
556   /* 2-2) compute numeric A_loc*P - dominating part */
557   /* ---------------------------------------------- */
558   /* get data from symbolic products */
559   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
560   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
561   apa   = ptap->apa;
562   api   = ap->i;
563   apj   = ap->j;
564   for (i=0; i<am; i++) {
565     /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */
566     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
567     apnz = api[i+1] - api[i];
568     for (j=0; j<apnz; j++) {
569       col = apj[j+api[i]];
570       ap->a[j+ap->i[i]] = apa[col];
571       apa[col] = 0.0;
572     }
573     ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
574   }
575 #if defined(PTAP_PROFILE)
576   ierr = PetscTime(&t2);CHKERRQ(ierr);
577   eAP = t2 - t1;
578 #endif
579 
580   /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */
581   ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr);
582   ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr);
583   C_loc = ptap->C_loc;
584   C_oth = ptap->C_oth;
585 #if defined(PTAP_PROFILE)
586   ierr = PetscTime(&t3);CHKERRQ(ierr);
587   eCseq = t3 - t2;
588 #endif
589 
590   /* add C_loc and Co to to C */
591   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
592 
593   /* C_loc -> C */
594   cm    = C_loc->rmap->N;
595   c_seq = (Mat_SeqAIJ*)C_loc->data;
596   cols = c_seq->j;
597   vals = c_seq->a;
598   for (i=0; i<cm; i++) {
599     ncols = c_seq->i[i+1] - c_seq->i[i];
600     row = rstart + i;
601     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
602     cols += ncols; vals += ncols;
603   }
604 
605   /* Co -> C, off-processor part */
606   cm = C_oth->rmap->N;
607   c_seq = (Mat_SeqAIJ*)C_oth->data;
608   cols = c_seq->j;
609   vals = c_seq->a;
610   for (i=0; i<cm; i++) {
611     ncols = c_seq->i[i+1] - c_seq->i[i];
612     row = p->garray[i];
613     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
614     cols += ncols; vals += ncols;
615   }
616   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
617   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
618 #if defined(PTAP_PROFILE)
619   ierr = PetscTime(&t4);CHKERRQ(ierr);
620   eCmpi = t4 - t3;
621 
622   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
623   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
624   if (rank==1) {
625     ierr = PetscPrintf(MPI_COMM_SELF," R %g, AP %g, Cseq %g, Cmpi %g = %g\n", eR,eAP,eCseq,eCmpi,eR+eAP+eCseq+eCmpi);CHKERRQ(ierr);
626   }
627 #endif
628   ptap->reuse = MAT_REUSE_MATRIX;
629   PetscFunctionReturn(0);
630 }
631 
632 #undef __FUNCT__
633 #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap"
634 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(Mat A,Mat P,PetscReal fill,Mat *C)
635 {
636   PetscErrorCode      ierr;
637   Mat                 Cmpi;
638   Mat_PtAPMPI         *ptap;
639   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
640   Mat_MPIAIJ          *a        =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
641   Mat_SeqAIJ          *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
642   Mat_SeqAIJ          *p_loc,*p_oth;
643   PetscInt            *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ;
644   PetscInt            *adi=ad->i,*aj,*aoi=ao->i,nnz;
645   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
646   PetscInt            am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n;
647   PetscBT             lnkbt;
648   MPI_Comm            comm;
649   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0;
650   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
651   PetscInt            len,proc,*dnz,*onz,*owners;
652   PetscInt            nzi,*pti,*ptj;
653   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
654   MPI_Request         *swaits,*rwaits;
655   MPI_Status          *sstatus,rstatus;
656   Mat_Merge_SeqsToMPI *merge;
657   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0;
658   PetscReal           afill=1.0,afill_tmp;
659   PetscInt            rmax;
660 
661   PetscFunctionBegin;
662   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
663   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
664   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
665 
666   /* create struct Mat_PtAPMPI and attached it to C later */
667   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
668   ierr        = PetscNew(&merge);CHKERRQ(ierr);
669   ptap->merge = merge;
670   ptap->reuse = MAT_INITIAL_MATRIX;
671 
672   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
673   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
674 
675   /* get P_loc by taking all local rows of P */
676   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
677 
678   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
679   p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
680   pi_loc = p_loc->i; pj_loc = p_loc->j;
681   pi_oth = p_oth->i; pj_oth = p_oth->j;
682 
683   /* (1) compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth (api,apj) */
684   /*--------------------------------------------------------------------------*/
685   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
686   api[0] = 0;
687 
688   /* create and initialize a linked list */
689   ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr);
690 
691   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) -OOM for ex56, np=8k on Intrepid! */
692   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
693   current_space = free_space;
694 
695   for (i=0; i<am; i++) {
696     /* diagonal portion of A */
697     nzi = adi[i+1] - adi[i];
698     aj  = ad->j + adi[i];
699     for (j=0; j<nzi; j++) {
700       row  = aj[j];
701       pnz  = pi_loc[row+1] - pi_loc[row];
702       Jptr = pj_loc + pi_loc[row];
703       /* add non-zero cols of P into the sorted linked list lnk */
704       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
705     }
706     /* off-diagonal portion of A */
707     nzi = aoi[i+1] - aoi[i];
708     aj  = ao->j + aoi[i];
709     for (j=0; j<nzi; j++) {
710       row  = aj[j];
711       pnz  = pi_oth[row+1] - pi_oth[row];
712       Jptr = pj_oth + pi_oth[row];
713       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
714     }
715     apnz     = lnk[0];
716     api[i+1] = api[i] + apnz;
717     if (ap_rmax < apnz) ap_rmax = apnz;
718 
719     /* if free space is not available, double the total space in the list */
720     if (current_space->local_remaining<apnz) {
721       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
722       nspacedouble++;
723     }
724 
725     /* Copy data into free space, then initialize lnk */
726     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
727 
728     current_space->array           += apnz;
729     current_space->local_used      += apnz;
730     current_space->local_remaining -= apnz;
731   }
732 
733   /* Allocate space for apj, initialize apj, and */
734   /* destroy list of free space and other temporary array(s) */
735   ierr      = PetscMalloc1(api[am]+1,&apj);CHKERRQ(ierr);
736   ierr      = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
737   afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1);
738   if (afill_tmp > afill) afill = afill_tmp;
739 
740   /* (2) determine symbolic Co=(p->B)^T*AP - send to others (coi,coj)*/
741   /*-----------------------------------------------------------------*/
742   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
743 
744   /* then, compute symbolic Co = (p->B)^T*AP */
745   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
746                          >= (num of nonzero rows of C_seq) - pn */
747   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
748   coi[0] = 0;
749 
750   /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */
751   nnz           = fill*(poti[pon] + api[am]);
752   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
753   current_space = free_space;
754 
755   for (i=0; i<pon; i++) {
756     pnz = poti[i+1] - poti[i];
757     ptJ = potj + poti[i];
758     for (j=0; j<pnz; j++) {
759       row  = ptJ[j]; /* row of AP == col of Pot */
760       apnz = api[row+1] - api[row];
761       Jptr = apj + api[row];
762       /* add non-zero cols of AP into the sorted linked list lnk */
763       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
764     }
765     nnz = lnk[0];
766 
767     /* If free space is not available, double the total space in the list */
768     if (current_space->local_remaining<nnz) {
769       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
770       nspacedouble++;
771     }
772 
773     /* Copy data into free space, and zero out denserows */
774     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
775 
776     current_space->array           += nnz;
777     current_space->local_used      += nnz;
778     current_space->local_remaining -= nnz;
779 
780     coi[i+1] = coi[i] + nnz;
781   }
782 
783   ierr      = PetscMalloc1(coi[pon],&coj);CHKERRQ(ierr);
784   ierr      = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
785   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1);
786   if (afill_tmp > afill) afill = afill_tmp;
787   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
788 
789   /* (3) send j-array (coj) of Co to other processors */
790   /*--------------------------------------------------*/
791   ierr = PetscCalloc1(size,&merge->len_s);CHKERRQ(ierr);
792   len_s        = merge->len_s;
793   merge->nsend = 0;
794 
795 
796   /* determine row ownership */
797   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
798   merge->rowmap->n  = pn;
799   merge->rowmap->bs = 1;
800 
801   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
802   owners = merge->rowmap->range;
803 
804   /* determine the number of messages to send, their lengths */
805   ierr = PetscMalloc2(size,&len_si,size,&sstatus);CHKERRQ(ierr);
806   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
807   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
808 
809   proc = 0;
810   for (i=0; i<pon; i++) {
811     while (prmap[i] >= owners[proc+1]) proc++;
812     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
813     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
814   }
815 
816   len          = 0; /* max length of buf_si[], see (4) */
817   owners_co[0] = 0;
818   for (proc=0; proc<size; proc++) {
819     owners_co[proc+1] = owners_co[proc] + len_si[proc];
820     if (len_s[proc]) {
821       merge->nsend++;
822       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
823       len         += len_si[proc];
824     }
825   }
826 
827   /* determine the number and length of messages to receive for coi and coj  */
828   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
829   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
830 
831   /* post the Irecv and Isend of coj */
832   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
833   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
834   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
835   for (proc=0, k=0; proc<size; proc++) {
836     if (!len_s[proc]) continue;
837     i    = owners_co[proc];
838     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
839     k++;
840   }
841 
842   /* receives and sends of coj are complete */
843   for (i=0; i<merge->nrecv; i++) {
844     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
845   }
846   ierr = PetscFree(rwaits);CHKERRQ(ierr);
847   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
848 
849   /* (4) send and recv coi */
850   /*-----------------------*/
851   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
852   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
853   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
854   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
855   for (proc=0,k=0; proc<size; proc++) {
856     if (!len_s[proc]) continue;
857     /* form outgoing message for i-structure:
858          buf_si[0]:                 nrows to be sent
859                [1:nrows]:           row index (global)
860                [nrows+1:2*nrows+1]: i-structure index
861     */
862     /*-------------------------------------------*/
863     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
864     buf_si_i    = buf_si + nrows+1;
865     buf_si[0]   = nrows;
866     buf_si_i[0] = 0;
867     nrows       = 0;
868     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
869       nzi = coi[i+1] - coi[i];
870       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
871       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
872       nrows++;
873     }
874     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
875     k++;
876     buf_si += len_si[proc];
877   }
878   i = merge->nrecv;
879   while (i--) {
880     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
881   }
882   ierr = PetscFree(rwaits);CHKERRQ(ierr);
883   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
884 
885   ierr = PetscFree2(len_si,sstatus);CHKERRQ(ierr);
886   ierr = PetscFree(len_ri);CHKERRQ(ierr);
887   ierr = PetscFree(swaits);CHKERRQ(ierr);
888   ierr = PetscFree(buf_s);CHKERRQ(ierr);
889 
890   /* (5) compute the local portion of C (mpi mat) */
891   /*----------------------------------------------*/
892   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
893 
894   /* allocate pti array and free space for accumulating nonzero column info */
895   ierr   = PetscMalloc1(pn+1,&pti);CHKERRQ(ierr);
896   pti[0] = 0;
897 
898   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
899   nnz           = fill*(pi_loc[pm] + api[am]);
900   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
901   current_space = free_space;
902 
903   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
904   for (k=0; k<merge->nrecv; k++) {
905     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
906     nrows       = *buf_ri_k[k];
907     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
908     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
909   }
910   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
911   rmax = 0;
912   for (i=0; i<pn; i++) {
913     /* add pdt[i,:]*AP into lnk */
914     pnz = pdti[i+1] - pdti[i];
915     ptJ = pdtj + pdti[i];
916     for (j=0; j<pnz; j++) {
917       row  = ptJ[j];  /* row of AP == col of Pt */
918       apnz = api[row+1] - api[row];
919       Jptr = apj + api[row];
920       /* add non-zero cols of AP into the sorted linked list lnk */
921       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
922     }
923 
924     /* add received col data into lnk */
925     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
926       if (i == *nextrow[k]) { /* i-th row */
927         nzi  = *(nextci[k]+1) - *nextci[k];
928         Jptr = buf_rj[k] + *nextci[k];
929         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
930         nextrow[k]++; nextci[k]++;
931       }
932     }
933     nnz = lnk[0];
934 
935     /* if free space is not available, make more free space */
936     if (current_space->local_remaining<nnz) {
937       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
938       nspacedouble++;
939     }
940     /* copy data into free space, then initialize lnk */
941     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
942     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
943 
944     current_space->array           += nnz;
945     current_space->local_used      += nnz;
946     current_space->local_remaining -= nnz;
947 
948     pti[i+1] = pti[i] + nnz;
949     if (nnz > rmax) rmax = nnz;
950   }
951   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
952   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
953 
954   ierr      = PetscMalloc1(pti[pn]+1,&ptj);CHKERRQ(ierr);
955   ierr      = PetscFreeSpaceContiguous(&free_space,ptj);CHKERRQ(ierr);
956   afill_tmp = (PetscReal)pti[pn]/(pi_loc[pm] + api[am]+1);
957   if (afill_tmp > afill) afill = afill_tmp;
958   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
959 
960   /* (6) create symbolic parallel matrix Cmpi */
961   /*------------------------------------------*/
962   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
963   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
964   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
965   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
966   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
967   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
968 
969   merge->bi        = pti;      /* Cseq->i */
970   merge->bj        = ptj;      /* Cseq->j */
971   merge->coi       = coi;      /* Co->i   */
972   merge->coj       = coj;      /* Co->j   */
973   merge->buf_ri    = buf_ri;
974   merge->buf_rj    = buf_rj;
975   merge->owners_co = owners_co;
976   merge->destroy   = Cmpi->ops->destroy;
977 
978   /* attach the supporting struct to Cmpi for reuse */
979   c           = (Mat_MPIAIJ*)Cmpi->data;
980   c->ptap     = ptap;
981   ptap->api   = api;
982   ptap->apj   = apj;
983   ptap->rmax  = ap_rmax;
984   ptap->duplicate = Cmpi->ops->duplicate;
985   ptap->destroy   = Cmpi->ops->destroy;
986 
987   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
988   Cmpi->assembled        = PETSC_FALSE;
989   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
990   Cmpi->ops->duplicate   = MatDuplicate_MPIAIJ_MatPtAP;
991   Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap;
992   *C                     = Cmpi;
993 
994   /* flag 'scalable' determines which implementations to be used:
995        0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa;
996        1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */
997   /* set default scalable */
998   ptap->scalable = PETSC_FALSE; /* PETSC_TRUE; */
999 
1000   ierr = PetscOptionsGetBool(((PetscObject)Cmpi)->prefix,"-matptap_scalable",&ptap->scalable,NULL);CHKERRQ(ierr);
1001   if (!ptap->scalable) {  /* Do dense axpy */
1002     ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
1003   } else {
1004     ierr = PetscCalloc1(ap_rmax+1,&ptap->apa);CHKERRQ(ierr);
1005   }
1006 
1007 #if defined(PETSC_USE_INFO)
1008   if (pti[pn] != 0) {
1009     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
1010     ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1011   } else {
1012     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1013   }
1014 #endif
1015   PetscFunctionReturn(0);
1016 }
1017 
1018 #undef __FUNCT__
1019 #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap"
1020 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap(Mat A,Mat P,Mat C)
1021 {
1022   PetscErrorCode      ierr;
1023   Mat_MPIAIJ          *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1024   Mat_SeqAIJ          *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
1025   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1026   Mat_SeqAIJ          *p_loc,*p_oth;
1027   Mat_PtAPMPI         *ptap;
1028   Mat_Merge_SeqsToMPI *merge;
1029   PetscInt            *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp;
1030   PetscInt            *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj;
1031   PetscInt            i,j,k,anz,pnz,apnz,nextap,row,*cj;
1032   MatScalar           *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp;
1033   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1034   MPI_Comm            comm;
1035   PetscMPIInt         size,rank,taga,*len_s;
1036   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1037   PetscInt            **buf_ri,**buf_rj;
1038   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
1039   MPI_Request         *s_waits,*r_waits;
1040   MPI_Status          *status;
1041   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
1042   PetscInt            *api,*apj,*coi,*coj;
1043   PetscInt            *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend;
1044   PetscBool           scalable;
1045 #if defined(PTAP_PROFILE)
1046   PetscLogDouble t0,t1,t2,eP,t3,t4,et2_AP=0.0,ePtAP=0.0,t2_0,t2_1,t2_2;
1047 #endif
1048 
1049   PetscFunctionBegin;
1050   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
1051 #if defined(PTAP_PROFILE)
1052   ierr = PetscTime(&t0);CHKERRQ(ierr);
1053 #endif
1054   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1055   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1056 
1057   ptap = c->ptap;
1058   if (!ptap) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_INCOMP,"MatPtAP() has not been called to create matrix C yet, cannot use MAT_REUSE_MATRIX");
1059   merge    = ptap->merge;
1060   apa      = ptap->apa;
1061   scalable = ptap->scalable;
1062 
1063   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
1064   /*-----------------------------------------------------*/
1065   if (ptap->reuse == MAT_INITIAL_MATRIX) {
1066     /* P_oth and P_loc are obtained in MatPtASymbolic(), skip calling MatGetBrowsOfAoCols() and MatMPIAIJGetLocalMat() */
1067     ptap->reuse = MAT_REUSE_MATRIX;
1068   } else { /* update numerical values of P_oth and P_loc */
1069     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
1070     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
1071   }
1072 #if defined(PTAP_PROFILE)
1073   ierr = PetscTime(&t1);CHKERRQ(ierr);
1074   eP = t1-t0;
1075 #endif
1076   /*
1077   printf("[%d] Ad: %d, %d; Ao: %d, %d; P_loc: %d, %d; P_oth %d, %d;\n",rank,
1078          a->A->rmap->N,a->A->cmap->N,a->B->rmap->N,a->B->cmap->N,
1079          ptap->P_loc->rmap->N,ptap->P_loc->cmap->N,
1080          ptap->P_oth->rmap->N,ptap->P_oth->cmap->N);
1081    */
1082 
1083   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
1084   /*--------------------------------------------------------------*/
1085   /* get data from symbolic products */
1086   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
1087   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
1088   pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a;
1089   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
1090 
1091   coi  = merge->coi; coj = merge->coj;
1092   ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
1093 
1094   bi     = merge->bi; bj = merge->bj;
1095   owners = merge->rowmap->range;
1096   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);  /* ba: Cseq->a */
1097 
1098   api = ptap->api; apj = ptap->apj;
1099 
1100   if (!scalable) { /* Do dense axpy on apa (length of pN, stores A[i,:]*P) - nonscalable, but faster (could take 1/3 scalable time) */
1101     ierr = PetscInfo(C,"Using non-scalable dense axpy\n");CHKERRQ(ierr);
1102 #if defined(PTAP_PROFILE)
1103     ierr = PetscTime(&t1);CHKERRQ(ierr);
1104 #endif
1105     for (i=0; i<am; i++) {
1106 #if defined(PTAP_PROFILE)
1107       ierr = PetscTime(&t2_0);CHKERRQ(ierr);
1108 #endif
1109       /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */
1110       /*------------------------------------------------------------*/
1111       apJ = apj + api[i];
1112 
1113       /* diagonal portion of A */
1114       anz = adi[i+1] - adi[i];
1115       adj = ad->j + adi[i];
1116       ada = ad->a + adi[i];
1117       for (j=0; j<anz; j++) {
1118         row = adj[j];
1119         pnz = pi_loc[row+1] - pi_loc[row];
1120         pj  = pj_loc + pi_loc[row];
1121         pa  = pa_loc + pi_loc[row];
1122 
1123         /* perform dense axpy */
1124         valtmp = ada[j];
1125         for (k=0; k<pnz; k++) {
1126           apa[pj[k]] += valtmp*pa[k];
1127         }
1128         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
1129       }
1130 
1131       /* off-diagonal portion of A */
1132       anz = aoi[i+1] - aoi[i];
1133       aoj = ao->j + aoi[i];
1134       aoa = ao->a + aoi[i];
1135       for (j=0; j<anz; j++) {
1136         row = aoj[j];
1137         pnz = pi_oth[row+1] - pi_oth[row];
1138         pj  = pj_oth + pi_oth[row];
1139         pa  = pa_oth + pi_oth[row];
1140 
1141         /* perform dense axpy */
1142         valtmp = aoa[j];
1143         for (k=0; k<pnz; k++) {
1144           apa[pj[k]] += valtmp*pa[k];
1145         }
1146         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
1147       }
1148 #if defined(PTAP_PROFILE)
1149       ierr    = PetscTime(&t2_1);CHKERRQ(ierr);
1150       et2_AP += t2_1 - t2_0;
1151 #endif
1152 
1153       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
1154       /*--------------------------------------------------------------*/
1155       apnz = api[i+1] - api[i];
1156       /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */
1157       pnz = po->i[i+1] - po->i[i];
1158       poJ = po->j + po->i[i];
1159       pA  = po->a + po->i[i];
1160       for (j=0; j<pnz; j++) {
1161         row = poJ[j];
1162         cnz = coi[row+1] - coi[row];
1163         cj  = coj + coi[row];
1164         ca  = coa + coi[row];
1165         /* perform dense axpy */
1166         valtmp = pA[j];
1167         for (k=0; k<cnz; k++) {
1168           ca[k] += valtmp*apa[cj[k]];
1169         }
1170         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1171       }
1172       /* put the value into Cd (diagonal part) */
1173       pnz = pd->i[i+1] - pd->i[i];
1174       pdJ = pd->j + pd->i[i];
1175       pA  = pd->a + pd->i[i];
1176       for (j=0; j<pnz; j++) {
1177         row = pdJ[j];
1178         cnz = bi[row+1] - bi[row];
1179         cj  = bj + bi[row];
1180         ca  = ba + bi[row];
1181         /* perform dense axpy */
1182         valtmp = pA[j];
1183         for (k=0; k<cnz; k++) {
1184           ca[k] += valtmp*apa[cj[k]];
1185         }
1186         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1187       }
1188 
1189       /* zero the current row of A*P */
1190       for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0;
1191 #if defined(PTAP_PROFILE)
1192       ierr      = PetscTime(&t2_2);CHKERRQ(ierr);
1193       ePtAP += t2_2 - t2_1;
1194 #endif
1195     }
1196   } else { /* Do sparse axpy on apa (length of ap_rmax, stores A[i,:]*P) - scalable, but slower */
1197     ierr = PetscInfo(C,"Using scalable sparse axpy\n");CHKERRQ(ierr);
1198     /*-----------------------------------------------------------------------------------------*/
1199     pA=pa_loc;
1200     for (i=0; i<am; i++) {
1201 #if defined(PTAP_PROFILE)
1202       ierr = PetscTime(&t2_0);CHKERRQ(ierr);
1203 #endif
1204       /* form i-th sparse row of A*P */
1205       apnz = api[i+1] - api[i];
1206       apJ  = apj + api[i];
1207       /* diagonal portion of A */
1208       anz = adi[i+1] - adi[i];
1209       adj = ad->j + adi[i];
1210       ada = ad->a + adi[i];
1211       for (j=0; j<anz; j++) {
1212         row    = adj[j];
1213         pnz    = pi_loc[row+1] - pi_loc[row];
1214         pj     = pj_loc + pi_loc[row];
1215         pa     = pa_loc + pi_loc[row];
1216         valtmp = ada[j];
1217         nextp  = 0;
1218         for (k=0; nextp<pnz; k++) {
1219           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
1220             apa[k] += valtmp*pa[nextp++];
1221           }
1222         }
1223         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
1224       }
1225       /* off-diagonal portion of A */
1226       anz = aoi[i+1] - aoi[i];
1227       aoj = ao->j + aoi[i];
1228       aoa = ao->a + aoi[i];
1229       for (j=0; j<anz; j++) {
1230         row    = aoj[j];
1231         pnz    = pi_oth[row+1] - pi_oth[row];
1232         pj     = pj_oth + pi_oth[row];
1233         pa     = pa_oth + pi_oth[row];
1234         valtmp = aoa[j];
1235         nextp  = 0;
1236         for (k=0; nextp<pnz; k++) {
1237           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
1238             apa[k] += valtmp*pa[nextp++];
1239           }
1240         }
1241         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
1242       }
1243 #if defined(PTAP_PROFILE)
1244       ierr    = PetscTime(&t2_1);CHKERRQ(ierr);
1245       et2_AP += t2_1 - t2_0;
1246 #endif
1247 
1248       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
1249       /*--------------------------------------------------------------*/
1250       pnz = pi_loc[i+1] - pi_loc[i];
1251       pJ  = pj_loc + pi_loc[i];
1252       for (j=0; j<pnz; j++) {
1253         nextap = 0;
1254         row    = pJ[j]; /* global index */
1255         if (row < pcstart || row >=pcend) { /* put the value into Co */
1256           row = *poJ;
1257           cj  = coj + coi[row];
1258           ca  = coa + coi[row]; poJ++;
1259         } else {                            /* put the value into Cd */
1260           row = *pdJ;
1261           cj  = bj + bi[row];
1262           ca  = ba + bi[row]; pdJ++;
1263         }
1264         valtmp = pA[j];
1265         for (k=0; nextap<apnz; k++) {
1266           if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++];
1267         }
1268         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
1269       }
1270       pA += pnz;
1271       /* zero the current row info for A*P */
1272       ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr);
1273 #if defined(PTAP_PROFILE)
1274       ierr      = PetscTime(&t2_2);CHKERRQ(ierr);
1275       ePtAP += t2_2 - t2_1;
1276 #endif
1277     }
1278   }
1279 #if defined(PTAP_PROFILE)
1280   ierr = PetscTime(&t2);CHKERRQ(ierr);
1281 #endif
1282 
1283   /* 3) send and recv matrix values coa */
1284   /*------------------------------------*/
1285   buf_ri = merge->buf_ri;
1286   buf_rj = merge->buf_rj;
1287   len_s  = merge->len_s;
1288   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
1289   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
1290 
1291   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
1292   for (proc=0,k=0; proc<size; proc++) {
1293     if (!len_s[proc]) continue;
1294     i    = merge->owners_co[proc];
1295     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
1296     k++;
1297   }
1298   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
1299   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
1300 
1301   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
1302   ierr = PetscFree(r_waits);CHKERRQ(ierr);
1303   ierr = PetscFree(coa);CHKERRQ(ierr);
1304 #if defined(PTAP_PROFILE)
1305   ierr = PetscTime(&t3);CHKERRQ(ierr);
1306 #endif
1307 
1308   /* 4) insert local Cseq and received values into Cmpi */
1309   /*------------------------------------------------------*/
1310   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
1311   for (k=0; k<merge->nrecv; k++) {
1312     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1313     nrows       = *(buf_ri_k[k]);
1314     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1315     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1316   }
1317 
1318   for (i=0; i<cm; i++) {
1319     row  = owners[rank] + i; /* global row index of C_seq */
1320     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1321     ba_i = ba + bi[i];
1322     bnz  = bi[i+1] - bi[i];
1323     /* add received vals into ba */
1324     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1325       /* i-th row */
1326       if (i == *nextrow[k]) {
1327         cnz    = *(nextci[k]+1) - *nextci[k];
1328         cj     = buf_rj[k] + *(nextci[k]);
1329         ca     = abuf_r[k] + *(nextci[k]);
1330         nextcj = 0;
1331         for (j=0; nextcj<cnz; j++) {
1332           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1333             ba_i[j] += ca[nextcj++];
1334           }
1335         }
1336         nextrow[k]++; nextci[k]++;
1337         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1338       }
1339     }
1340     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
1341   }
1342   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1343   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1344 
1345   ierr = PetscFree(ba);CHKERRQ(ierr);
1346   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
1347   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
1348   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
1349 #if defined(PTAP_PROFILE)
1350   ierr = PetscTime(&t4);CHKERRQ(ierr);
1351   if (rank==1) {
1352     ierr = PetscPrintf(MPI_COMM_SELF,"  [%d] PtAPNum %g/P + %g/PtAP( %g/A*P + %g/Pt*AP ) + %g/comm + %g/Cloc = %g\n\n",rank,eP,t2-t1,et2_AP,ePtAP,t3-t2,t4-t3,t4-t0);CHKERRQ(ierr);
1353   }
1354 #endif
1355   PetscFunctionReturn(0);
1356 }
1357