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