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