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