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