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