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