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