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