xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision 2ec85ea456d4da44b4dd09fd939e615b8c7fb52b)
1 
2 /*
3   Defines projective product routines where A is a MPIAIJ matrix
4           C = P^T * A * P
5 */
6 
7 #include <../src/mat/impls/aij/seq/aij.h>   /*I "petscmat.h" I*/
8 #include <../src/mat/utils/freespace.h>
9 #include <../src/mat/impls/aij/mpi/mpiaij.h>
10 #include <petscbt.h>
11 #include <petsctime.h>
12 
13 /* #define PTAP_PROFILE */
14 
15 extern PetscErrorCode MatDestroy_MPIAIJ(Mat);
16 PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A)
17 {
18   PetscErrorCode ierr;
19   Mat_MPIAIJ     *a=(Mat_MPIAIJ*)A->data;
20   Mat_PtAPMPI    *ptap=a->ptap;
21 
22   PetscFunctionBegin;
23   if (ptap) {
24     Mat_Merge_SeqsToMPI *merge=ptap->merge;
25     ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
26     ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
27     ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
28     ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
29     ierr = MatDestroy(&ptap->A_loc);CHKERRQ(ierr); /* used by MatTransposeMatMult() */
30     ierr = MatDestroy(&ptap->Rd);CHKERRQ(ierr);
31     ierr = MatDestroy(&ptap->Ro);CHKERRQ(ierr);
32     if (ptap->AP_loc) { /* used by alg_rap */
33       Mat_SeqAIJ *ap = (Mat_SeqAIJ*)(ptap->AP_loc)->data;
34       ierr = PetscFree(ap->i);CHKERRQ(ierr);
35       ierr = PetscFree2(ap->j,ap->a);CHKERRQ(ierr);
36       ierr = MatDestroy(&ptap->AP_loc);CHKERRQ(ierr);
37     } else { /* used by alg_ptap */
38       ierr = PetscFree(ptap->api);CHKERRQ(ierr);
39       ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
40     }
41     ierr = MatDestroy(&ptap->C_loc);CHKERRQ(ierr);
42     ierr = MatDestroy(&ptap->C_oth);CHKERRQ(ierr);
43     if (ptap->apa) {ierr = PetscFree(ptap->apa);CHKERRQ(ierr);}
44 
45     if (merge) { /* used by alg_ptap */
46       ierr = PetscFree(merge->id_r);CHKERRQ(ierr);
47       ierr = PetscFree(merge->len_s);CHKERRQ(ierr);
48       ierr = PetscFree(merge->len_r);CHKERRQ(ierr);
49       ierr = PetscFree(merge->bi);CHKERRQ(ierr);
50       ierr = PetscFree(merge->bj);CHKERRQ(ierr);
51       ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr);
52       ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr);
53       ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr);
54       ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr);
55       ierr = PetscFree(merge->coi);CHKERRQ(ierr);
56       ierr = PetscFree(merge->coj);CHKERRQ(ierr);
57       ierr = PetscFree(merge->owners_co);CHKERRQ(ierr);
58       ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr);
59       ierr = PetscFree(ptap->merge);CHKERRQ(ierr);
60     }
61 
62     ierr = ptap->destroy(A);CHKERRQ(ierr);
63     ierr = PetscFree(ptap);CHKERRQ(ierr);
64   }
65   PetscFunctionReturn(0);
66 }
67 
68 PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M)
69 {
70   PetscErrorCode ierr;
71   Mat_MPIAIJ     *a     = (Mat_MPIAIJ*)A->data;
72   Mat_PtAPMPI    *ptap  = a->ptap;
73 
74   PetscFunctionBegin;
75   ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr);
76   (*M)->ops->destroy   = ptap->destroy;
77   (*M)->ops->duplicate = ptap->duplicate;
78   PetscFunctionReturn(0);
79 }
80 
81 PETSC_INTERN PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
82 {
83   PetscErrorCode ierr;
84   PetscBool      rap=PETSC_TRUE; /* do R=P^T locally, then C=R*A*P */
85   MPI_Comm       comm;
86 
87   PetscFunctionBegin;
88   /* check if matrix local sizes are compatible */
89   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
90   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, Arow (%D, %D) != Prow (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend);
91   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, Acol (%D, %D) != Prow (%D,%D)",A->cmap->rstart,A->cmap->rend,P->rmap->rstart,P->rmap->rend);
92 
93   ierr = PetscOptionsGetBool(NULL,NULL,"-matrap",&rap,NULL);CHKERRQ(ierr);
94   if (scall == MAT_INITIAL_MATRIX) {
95     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
96     if (rap) { /* do R=P^T locally, then C=R*A*P */
97       ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);
98     } else {       /* do P^T*A*P */
99       ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(A,P,fill,C);CHKERRQ(ierr);
100     }
101     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
102   }
103   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
104   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
105   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
106   PetscFunctionReturn(0);
107 }
108 
109 #if defined(PETSC_HAVE_HYPRE)
110 PETSC_INTERN PetscErrorCode MatPtAPSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat*);
111 #endif
112 
113 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
114 {
115   PetscErrorCode      ierr;
116   Mat_PtAPMPI         *ptap;
117   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
118   MPI_Comm            comm;
119   PetscMPIInt         size,rank;
120   Mat                 Cmpi;
121   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
122   PetscInt            am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
123   PetscInt            *lnk,i,k,pnz,row,nsend;
124   PetscBT             lnkbt;
125   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
126   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
127   PetscInt            len,proc,*dnz,*onz,*owners,nzi,nspacedouble;
128   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
129   MPI_Request         *swaits,*rwaits;
130   MPI_Status          *sstatus,rstatus;
131   PetscLayout         rowmap;
132   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
133   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
134   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0,Crmax,*aj,*ai,*pi;
135   Mat_SeqAIJ          *p_loc,*p_oth,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*c_loc,*c_oth;
136   PetscScalar         *apv;
137   PetscTable          ta;
138 #if defined(PETSC_HAVE_HYPRE)
139   const char          *algTypes[3] = {"scalable","nonscalable","hypre"};
140   PetscInt            nalg = 3;
141 #else
142   const char          *algTypes[2] = {"scalable","nonscalable"};
143   PetscInt            nalg = 2;
144 #endif
145   PetscInt            alg = 1; /* set default algorithm */
146 #if defined(PETSC_USE_INFO)
147   PetscReal           apfill;
148 #endif
149 #if defined(PTAP_PROFILE)
150   PetscLogDouble      t0,t1,t11,t12,t2,t3,t4;
151 #endif
152 
153   PetscFunctionBegin;
154   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
155   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
156   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
157 
158   /* create symbolic parallel matrix Cmpi */
159   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
160   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
161 
162   ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
163   ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[1],&alg,NULL);CHKERRQ(ierr);
164   ierr = PetscOptionsEnd();CHKERRQ(ierr);
165 
166   if (alg == 1) {
167     /* Do dense axpy in MatPtAPNumeric_MPIAIJ_MPIAIJ() */
168     Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ;
169   } else if (alg == 0) {
170     Cmpi->ops->ptapnumeric = 0; /* not done yet */
171     SETERRQ(comm,PETSC_ERR_ARG_SIZ,"Not done yet");
172 #if defined(PETSC_HAVE_HYPRE)
173   } else if (alg == 2) {
174     /* Use boomerAMGBuildCoarseOperator */
175     ierr = MatDestroy(&Cmpi);CHKERRQ(ierr);
176     ierr = MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);CHKERRQ(ierr);
177     PetscFunctionReturn(0);
178 #endif
179   } else {
180     SETERRQ1(comm,PETSC_ERR_PLIB,"Unknown algorithm number %D\n",alg);
181   }
182 
183 #if defined(PTAP_PROFILE)
184   ierr = PetscTime(&t0);CHKERRQ(ierr);
185 #endif
186 
187   /* create struct Mat_PtAPMPI and attached it to C later */
188   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
189   ptap->reuse = MAT_INITIAL_MATRIX;
190 
191   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
192   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
193   /* get P_loc by taking all local rows of P */
194   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
195 
196   /* (0) compute Rd = Pd^T, Ro = Po^T  */
197   /* --------------------------------- */
198   ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
199   ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
200 #if defined(PTAP_PROFILE)
201   ierr = PetscTime(&t11);CHKERRQ(ierr);
202 #endif
203 
204   /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */
205   /* ----------------------------------------------------------------- */
206   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
207   p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
208 
209   /* create and initialize a linked list */
210   ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */
211   MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta);
212   MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta);
213   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */
214   /* printf("[%d] est %d, Crmax %d; pN %d\n",rank,5*(p_loc->rmax+p_oth->rmax + (PetscInt)(1.e-2*pN)),Crmax,pN); */
215 
216   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
217 
218   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
219   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr);
220   current_space = free_space;
221   nspacedouble  = 0;
222 
223   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
224   api[0] = 0;
225   for (i=0; i<am; i++) {
226     /* diagonal portion: Ad[i,:]*P */
227     ai = ad->i; pi = p_loc->i;
228     nzi = ai[i+1] - ai[i];
229     aj  = ad->j + ai[i];
230     for (j=0; j<nzi; j++) {
231       row  = aj[j];
232       pnz  = pi[row+1] - pi[row];
233       Jptr = p_loc->j + pi[row];
234       /* add non-zero cols of P into the sorted linked list lnk */
235       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
236     }
237     /* off-diagonal portion: Ao[i,:]*P */
238     ai = ao->i; pi = p_oth->i;
239     nzi = ai[i+1] - ai[i];
240     aj  = ao->j + ai[i];
241     for (j=0; j<nzi; j++) {
242       row  = aj[j];
243       pnz  = pi[row+1] - pi[row];
244       Jptr = p_oth->j + pi[row];
245       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
246     }
247     apnz     = lnk[0];
248     api[i+1] = api[i] + apnz;
249     if (ap_rmax < apnz) ap_rmax = apnz;
250 
251     /* if free space is not available, double the total space in the list */
252     if (current_space->local_remaining<apnz) {
253       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
254       nspacedouble++;
255     }
256 
257     /* Copy data into free space, then initialize lnk */
258     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
259 
260     current_space->array           += apnz;
261     current_space->local_used      += apnz;
262     current_space->local_remaining -= apnz;
263   }
264   /* Allocate space for apj and apv, initialize apj, and */
265   /* destroy list of free space and other temporary array(s) */
266   ierr   = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr);
267   ierr   = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
268   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
269 
270   /* Create AP_loc for reuse */
271   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr);
272 
273 #if defined(PETSC_USE_INFO)
274   apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1);
275   ptap->AP_loc->info.mallocs           = nspacedouble;
276   ptap->AP_loc->info.fill_ratio_given  = fill;
277   ptap->AP_loc->info.fill_ratio_needed = apfill;
278 
279   if (api[am]) {
280     ierr = PetscInfo3(ptap->AP_loc,"AP_loc reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)apfill);CHKERRQ(ierr);
281     ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr);
282   } else {
283     ierr = PetscInfo(ptap->AP_loc,"AP_loc is empty \n");CHKERRQ(ierr);
284   }
285 #endif
286 
287 #if defined(PTAP_PROFILE)
288   ierr = PetscTime(&t12);CHKERRQ(ierr);
289 #endif
290 
291   /* (2-1) compute symbolic Co = Ro*AP_loc  */
292   /* ------------------------------------ */
293   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
294 #if defined(PTAP_PROFILE)
295   ierr = PetscTime(&t1);CHKERRQ(ierr);
296 #endif
297 
298   /* (3) send coj of C_oth to other processors  */
299   /* ------------------------------------------ */
300   /* determine row ownership */
301   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
302   rowmap->n  = pn;
303   rowmap->bs = 1;
304   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
305   owners = rowmap->range;
306 
307   /* determine the number of messages to send, their lengths */
308   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
309   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
310   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
311 
312   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
313   coi   = c_oth->i; coj = c_oth->j;
314   con   = ptap->C_oth->rmap->n;
315   proc  = 0;
316   for (i=0; i<con; i++) {
317     while (prmap[i] >= owners[proc+1]) proc++;
318     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
319     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
320   }
321 
322   len          = 0; /* max length of buf_si[], see (4) */
323   owners_co[0] = 0;
324   nsend        = 0;
325   for (proc=0; proc<size; proc++) {
326     owners_co[proc+1] = owners_co[proc] + len_si[proc];
327     if (len_s[proc]) {
328       nsend++;
329       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
330       len         += len_si[proc];
331     }
332   }
333 
334   /* determine the number and length of messages to receive for coi and coj  */
335   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
336   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
337 
338   /* post the Irecv and Isend of coj */
339   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
340   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
341   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
342   for (proc=0, k=0; proc<size; proc++) {
343     if (!len_s[proc]) continue;
344     i    = owners_co[proc];
345     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
346     k++;
347   }
348 
349   /* (2-2) compute symbolic C_loc = Rd*AP_loc */
350   /* ---------------------------------------- */
351   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr);
352   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
353 
354   /* receives coj are complete */
355   for (i=0; i<nrecv; i++) {
356     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
357   }
358   ierr = PetscFree(rwaits);CHKERRQ(ierr);
359   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
360 
361   /* add received column indices into ta to update Crmax */
362   for (k=0; k<nrecv; k++) {/* k-th received message */
363     Jptr = buf_rj[k];
364     for (j=0; j<len_r[k]; j++) {
365       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
366     }
367   }
368   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
369   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
370 
371   /* (4) send and recv coi */
372   /*-----------------------*/
373   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
374   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
375   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
376   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
377   for (proc=0,k=0; proc<size; proc++) {
378     if (!len_s[proc]) continue;
379     /* form outgoing message for i-structure:
380          buf_si[0]:                 nrows to be sent
381                [1:nrows]:           row index (global)
382                [nrows+1:2*nrows+1]: i-structure index
383     */
384     /*-------------------------------------------*/
385     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
386     buf_si_i    = buf_si + nrows+1;
387     buf_si[0]   = nrows;
388     buf_si_i[0] = 0;
389     nrows       = 0;
390     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
391       nzi = coi[i+1] - coi[i];
392       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
393       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
394       nrows++;
395     }
396     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
397     k++;
398     buf_si += len_si[proc];
399   }
400   for (i=0; i<nrecv; i++) {
401     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
402   }
403   ierr = PetscFree(rwaits);CHKERRQ(ierr);
404   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
405 
406   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
407   ierr = PetscFree(len_ri);CHKERRQ(ierr);
408   ierr = PetscFree(swaits);CHKERRQ(ierr);
409   ierr = PetscFree(buf_s);CHKERRQ(ierr);
410 #if defined(PTAP_PROFILE)
411   ierr = PetscTime(&t2);CHKERRQ(ierr);
412 #endif
413   /* (5) compute the local portion of Cmpi      */
414   /* ------------------------------------------ */
415   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
416   ierr          = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr);
417   current_space = free_space;
418 
419   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
420   for (k=0; k<nrecv; k++) {
421     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
422     nrows       = *buf_ri_k[k];
423     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
424     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
425   }
426 
427   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
428   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
429   for (i=0; i<pn; i++) {
430     /* add C_loc into Cmpi */
431     nzi  = c_loc->i[i+1] - c_loc->i[i];
432     Jptr = c_loc->j + c_loc->i[i];
433     ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
434 
435     /* add received col data into lnk */
436     for (k=0; k<nrecv; k++) { /* k-th received message */
437       if (i == *nextrow[k]) { /* i-th row */
438         nzi  = *(nextci[k]+1) - *nextci[k];
439         Jptr = buf_rj[k] + *nextci[k];
440         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
441         nextrow[k]++; nextci[k]++;
442       }
443     }
444     nzi = lnk[0];
445 
446     /* copy data into free space, then initialize lnk */
447     ierr = PetscLLCondensedClean(pN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
448     ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr);
449   }
450   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
451   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
452   ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr);
453 #if defined(PTAP_PROFILE)
454   ierr = PetscTime(&t3);CHKERRQ(ierr);
455 #endif
456 
457   /* local sizes and preallocation */
458   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
459   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
460   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
461   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
462 
463   /* members in merge */
464   ierr = PetscFree(id_r);CHKERRQ(ierr);
465   ierr = PetscFree(len_r);CHKERRQ(ierr);
466   ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr);
467   ierr = PetscFree(buf_ri);CHKERRQ(ierr);
468   ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr);
469   ierr = PetscFree(buf_rj);CHKERRQ(ierr);
470   ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr);
471 
472   /* attach the supporting struct to Cmpi for reuse */
473   c = (Mat_MPIAIJ*)Cmpi->data;
474   c->ptap         = ptap;
475   ptap->duplicate = Cmpi->ops->duplicate;
476   ptap->destroy   = Cmpi->ops->destroy;
477 
478   if (alg == 1) {
479     ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
480   }
481 
482   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
483   Cmpi->assembled        = PETSC_FALSE;
484   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
485   Cmpi->ops->duplicate   = MatDuplicate_MPIAIJ_MatPtAP;
486   *C                     = Cmpi;
487 
488 #if defined(PTAP_PROFILE)
489   ierr = PetscTime(&t4);CHKERRQ(ierr);
490   if (rank == 1) {
491     printf("PtAPSym: %g + %g + %g + %g + %g + %g = %g\n",t11-t0,t1-t11,t12-t11,t2-t2,t3-t2,t4-t3,t4-t0);
492   }
493 #endif
494   PetscFunctionReturn(0);
495 }
496 
497 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
498 {
499   PetscErrorCode    ierr;
500   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
501   Mat_SeqAIJ        *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
502   Mat_SeqAIJ        *ap,*p_loc,*p_oth,*c_seq;
503   Mat_PtAPMPI       *ptap = c->ptap;
504   Mat               AP_loc,C_loc,C_oth;
505   PetscInt          i,rstart,rend,cm,ncols,row;
506   PetscInt          *api,*apj,am = A->rmap->n,j,col,apnz;
507   PetscScalar       *apa;
508   const PetscInt    *cols;
509   const PetscScalar *vals;
510 #if defined(PTAP_PROFILE)
511   PetscMPIInt       rank;
512   MPI_Comm          comm;
513   PetscLogDouble    t0,t1,t2,t3,t4,eR,eAP,eCseq,eCmpi;
514 #endif
515 
516   PetscFunctionBegin;
517   ierr = MatZeroEntries(C);CHKERRQ(ierr);
518 
519   /* 1) get R = Pd^T,Ro = Po^T */
520 #if defined(PTAP_PROFILE)
521   ierr = PetscTime(&t0);CHKERRQ(ierr);
522 #endif
523   if (ptap->reuse == MAT_REUSE_MATRIX) {
524     ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
525     ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
526   }
527 #if defined(PTAP_PROFILE)
528   ierr = PetscTime(&t1);CHKERRQ(ierr);
529   eR = t1 - t0;
530 #endif
531 
532   /* 2) get AP_loc */
533   AP_loc = ptap->AP_loc;
534   ap = (Mat_SeqAIJ*)AP_loc->data;
535 
536   /* 2-1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
537   /*-----------------------------------------------------*/
538   if (ptap->reuse == MAT_REUSE_MATRIX) {
539     /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */
540     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
541     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
542   }
543 
544   /* 2-2) compute numeric A_loc*P - dominating part */
545   /* ---------------------------------------------- */
546   /* get data from symbolic products */
547   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
548   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
549   apa   = ptap->apa;
550   api   = ap->i;
551   apj   = ap->j;
552   for (i=0; i<am; i++) {
553     /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */
554     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
555     apnz = api[i+1] - api[i];
556     for (j=0; j<apnz; j++) {
557       col = apj[j+api[i]];
558       ap->a[j+ap->i[i]] = apa[col];
559       apa[col] = 0.0;
560     }
561     ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
562   }
563 #if defined(PTAP_PROFILE)
564   ierr = PetscTime(&t2);CHKERRQ(ierr);
565   eAP = t2 - t1;
566 #endif
567 
568   /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */
569   ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr);
570   ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr);
571   C_loc = ptap->C_loc;
572   C_oth = ptap->C_oth;
573 #if defined(PTAP_PROFILE)
574   ierr = PetscTime(&t3);CHKERRQ(ierr);
575   eCseq = t3 - t2;
576 #endif
577 
578   /* add C_loc and Co to to C */
579   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
580 
581   /* C_loc -> C */
582   cm    = C_loc->rmap->N;
583   c_seq = (Mat_SeqAIJ*)C_loc->data;
584   cols = c_seq->j;
585   vals = c_seq->a;
586   for (i=0; i<cm; i++) {
587     ncols = c_seq->i[i+1] - c_seq->i[i];
588     row = rstart + i;
589     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
590     cols += ncols; vals += ncols;
591   }
592 
593   /* Co -> C, off-processor part */
594   cm = C_oth->rmap->N;
595   c_seq = (Mat_SeqAIJ*)C_oth->data;
596   cols = c_seq->j;
597   vals = c_seq->a;
598   for (i=0; i<cm; i++) {
599     ncols = c_seq->i[i+1] - c_seq->i[i];
600     row = p->garray[i];
601     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
602     cols += ncols; vals += ncols;
603   }
604   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
605   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
606 #if defined(PTAP_PROFILE)
607   ierr = PetscTime(&t4);CHKERRQ(ierr);
608   eCmpi = t4 - t3;
609 
610   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
611   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
612   if (rank==1) {
613     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);
614   }
615 #endif
616   ptap->reuse = MAT_REUSE_MATRIX;
617   PetscFunctionReturn(0);
618 }
619 
620 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(Mat A,Mat P,PetscReal fill,Mat *C)
621 {
622   PetscErrorCode      ierr;
623   Mat                 Cmpi;
624   Mat_PtAPMPI         *ptap;
625   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
626   Mat_MPIAIJ          *a        =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
627   Mat_SeqAIJ          *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
628   Mat_SeqAIJ          *p_loc,*p_oth;
629   PetscInt            *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ;
630   PetscInt            *adi=ad->i,*aj,*aoi=ao->i,nnz;
631   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
632   PetscInt            am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n;
633   PetscBT             lnkbt;
634   MPI_Comm            comm;
635   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0;
636   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
637   PetscInt            len,proc,*dnz,*onz,*owners;
638   PetscInt            nzi,*pti,*ptj;
639   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
640   MPI_Request         *swaits,*rwaits;
641   MPI_Status          *sstatus,rstatus;
642   Mat_Merge_SeqsToMPI *merge;
643   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0;
644   PetscReal           afill=1.0,afill_tmp;
645 
646   PetscFunctionBegin;
647   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
648   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
649   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
650 
651   /* create struct Mat_PtAPMPI and attached it to C later */
652   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
653   ierr        = PetscNew(&merge);CHKERRQ(ierr);
654   ptap->merge = merge;
655   ptap->reuse = MAT_INITIAL_MATRIX;
656 
657   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
658   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
659 
660   /* get P_loc by taking all local rows of P */
661   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
662 
663   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
664   p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
665   pi_loc = p_loc->i; pj_loc = p_loc->j;
666   pi_oth = p_oth->i; pj_oth = p_oth->j;
667 
668   /* (1) compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth (api,apj) */
669   /*--------------------------------------------------------------------------*/
670   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
671   api[0] = 0;
672 
673   /* create and initialize a linked list */
674   ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr);
675 
676   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) -OOM for ex56, np=8k on Intrepid! */
677   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr);
678   current_space = free_space;
679 
680   for (i=0; i<am; i++) {
681     /* diagonal portion of A */
682     nzi = adi[i+1] - adi[i];
683     aj  = ad->j + adi[i];
684     for (j=0; j<nzi; j++) {
685       row  = aj[j];
686       pnz  = pi_loc[row+1] - pi_loc[row];
687       Jptr = pj_loc + pi_loc[row];
688       /* add non-zero cols of P into the sorted linked list lnk */
689       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
690     }
691     /* off-diagonal portion of A */
692     nzi = aoi[i+1] - aoi[i];
693     aj  = ao->j + aoi[i];
694     for (j=0; j<nzi; j++) {
695       row  = aj[j];
696       pnz  = pi_oth[row+1] - pi_oth[row];
697       Jptr = pj_oth + pi_oth[row];
698       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
699     }
700     apnz     = lnk[0];
701     api[i+1] = api[i] + apnz;
702     if (ap_rmax < apnz) ap_rmax = apnz;
703 
704     /* if free space is not available, double the total space in the list */
705     if (current_space->local_remaining<apnz) {
706       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
707       nspacedouble++;
708     }
709 
710     /* Copy data into free space, then initialize lnk */
711     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
712 
713     current_space->array           += apnz;
714     current_space->local_used      += apnz;
715     current_space->local_remaining -= apnz;
716   }
717 
718   /* Allocate space for apj, initialize apj, and */
719   /* destroy list of free space and other temporary array(s) */
720   ierr      = PetscMalloc1(api[am]+1,&apj);CHKERRQ(ierr);
721   ierr      = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
722   afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1);
723   if (afill_tmp > afill) afill = afill_tmp;
724 
725   /* (2) determine symbolic Co=(p->B)^T*AP - send to others (coi,coj)*/
726   /*-----------------------------------------------------------------*/
727   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
728 
729   /* then, compute symbolic Co = (p->B)^T*AP */
730   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
731                          >= (num of nonzero rows of C_seq) - pn */
732   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
733   coi[0] = 0;
734 
735   /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */
736   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(poti[pon],api[am]));
737   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
738   current_space = free_space;
739 
740   for (i=0; i<pon; i++) {
741     pnz = poti[i+1] - poti[i];
742     ptJ = potj + poti[i];
743     for (j=0; j<pnz; j++) {
744       row  = ptJ[j]; /* row of AP == col of Pot */
745       apnz = api[row+1] - api[row];
746       Jptr = apj + api[row];
747       /* add non-zero cols of AP into the sorted linked list lnk */
748       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
749     }
750     nnz = lnk[0];
751 
752     /* If free space is not available, double the total space in the list */
753     if (current_space->local_remaining<nnz) {
754       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
755       nspacedouble++;
756     }
757 
758     /* Copy data into free space, and zero out denserows */
759     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
760 
761     current_space->array           += nnz;
762     current_space->local_used      += nnz;
763     current_space->local_remaining -= nnz;
764 
765     coi[i+1] = coi[i] + nnz;
766   }
767 
768   ierr      = PetscMalloc1(coi[pon],&coj);CHKERRQ(ierr);
769   ierr      = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
770   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1);
771   if (afill_tmp > afill) afill = afill_tmp;
772   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
773 
774   /* (3) send j-array (coj) of Co to other processors */
775   /*--------------------------------------------------*/
776   ierr = PetscCalloc1(size,&merge->len_s);CHKERRQ(ierr);
777   len_s        = merge->len_s;
778   merge->nsend = 0;
779 
780 
781   /* determine row ownership */
782   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
783   merge->rowmap->n  = pn;
784   merge->rowmap->bs = 1;
785 
786   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
787   owners = merge->rowmap->range;
788 
789   /* determine the number of messages to send, their lengths */
790   ierr = PetscMalloc2(size,&len_si,size,&sstatus);CHKERRQ(ierr);
791   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
792   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
793 
794   proc = 0;
795   for (i=0; i<pon; i++) {
796     while (prmap[i] >= owners[proc+1]) proc++;
797     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
798     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
799   }
800 
801   len          = 0; /* max length of buf_si[], see (4) */
802   owners_co[0] = 0;
803   for (proc=0; proc<size; proc++) {
804     owners_co[proc+1] = owners_co[proc] + len_si[proc];
805     if (len_s[proc]) {
806       merge->nsend++;
807       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
808       len         += len_si[proc];
809     }
810   }
811 
812   /* determine the number and length of messages to receive for coi and coj  */
813   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
814   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
815 
816   /* post the Irecv and Isend of coj */
817   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
818   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
819   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
820   for (proc=0, k=0; proc<size; proc++) {
821     if (!len_s[proc]) continue;
822     i    = owners_co[proc];
823     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
824     k++;
825   }
826 
827   /* receives and sends of coj are complete */
828   for (i=0; i<merge->nrecv; i++) {
829     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
830   }
831   ierr = PetscFree(rwaits);CHKERRQ(ierr);
832   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
833 
834   /* (4) send and recv coi */
835   /*-----------------------*/
836   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
837   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
838   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
839   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
840   for (proc=0,k=0; proc<size; proc++) {
841     if (!len_s[proc]) continue;
842     /* form outgoing message for i-structure:
843          buf_si[0]:                 nrows to be sent
844                [1:nrows]:           row index (global)
845                [nrows+1:2*nrows+1]: i-structure index
846     */
847     /*-------------------------------------------*/
848     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
849     buf_si_i    = buf_si + nrows+1;
850     buf_si[0]   = nrows;
851     buf_si_i[0] = 0;
852     nrows       = 0;
853     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
854       nzi = coi[i+1] - coi[i];
855       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
856       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
857       nrows++;
858     }
859     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
860     k++;
861     buf_si += len_si[proc];
862   }
863   i = merge->nrecv;
864   while (i--) {
865     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
866   }
867   ierr = PetscFree(rwaits);CHKERRQ(ierr);
868   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
869 
870   ierr = PetscFree2(len_si,sstatus);CHKERRQ(ierr);
871   ierr = PetscFree(len_ri);CHKERRQ(ierr);
872   ierr = PetscFree(swaits);CHKERRQ(ierr);
873   ierr = PetscFree(buf_s);CHKERRQ(ierr);
874 
875   /* (5) compute the local portion of C (mpi mat) */
876   /*----------------------------------------------*/
877   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
878 
879   /* allocate pti array and free space for accumulating nonzero column info */
880   ierr   = PetscMalloc1(pn+1,&pti);CHKERRQ(ierr);
881   pti[0] = 0;
882 
883   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
884   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(pi_loc[pm],api[am]));
885   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
886   current_space = free_space;
887 
888   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
889   for (k=0; k<merge->nrecv; k++) {
890     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
891     nrows       = *buf_ri_k[k];
892     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
893     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
894   }
895   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
896 
897   for (i=0; i<pn; i++) {
898     /* add pdt[i,:]*AP into lnk */
899     pnz = pdti[i+1] - pdti[i];
900     ptJ = pdtj + pdti[i];
901     for (j=0; j<pnz; j++) {
902       row  = ptJ[j];  /* row of AP == col of Pt */
903       apnz = api[row+1] - api[row];
904       Jptr = apj + api[row];
905       /* add non-zero cols of AP into the sorted linked list lnk */
906       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
907     }
908 
909     /* add received col data into lnk */
910     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
911       if (i == *nextrow[k]) { /* i-th row */
912         nzi  = *(nextci[k]+1) - *nextci[k];
913         Jptr = buf_rj[k] + *nextci[k];
914         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
915         nextrow[k]++; nextci[k]++;
916       }
917     }
918     nnz = lnk[0];
919 
920     /* if free space is not available, make more free space */
921     if (current_space->local_remaining<nnz) {
922       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
923       nspacedouble++;
924     }
925     /* copy data into free space, then initialize lnk */
926     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
927     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
928 
929     current_space->array           += nnz;
930     current_space->local_used      += nnz;
931     current_space->local_remaining -= nnz;
932 
933     pti[i+1] = pti[i] + nnz;
934   }
935   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
936   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
937 
938   ierr      = PetscMalloc1(pti[pn]+1,&ptj);CHKERRQ(ierr);
939   ierr      = PetscFreeSpaceContiguous(&free_space,ptj);CHKERRQ(ierr);
940   afill_tmp = (PetscReal)pti[pn]/(pi_loc[pm] + api[am]+1);
941   if (afill_tmp > afill) afill = afill_tmp;
942   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
943 
944   /* (6) create symbolic parallel matrix Cmpi */
945   /*------------------------------------------*/
946   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
947   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
948   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
949   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
950   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
951   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
952 
953   merge->bi        = pti;      /* Cseq->i */
954   merge->bj        = ptj;      /* Cseq->j */
955   merge->coi       = coi;      /* Co->i   */
956   merge->coj       = coj;      /* Co->j   */
957   merge->buf_ri    = buf_ri;
958   merge->buf_rj    = buf_rj;
959   merge->owners_co = owners_co;
960   merge->destroy   = Cmpi->ops->destroy;
961 
962   /* attach the supporting struct to Cmpi for reuse */
963   c           = (Mat_MPIAIJ*)Cmpi->data;
964   c->ptap     = ptap;
965   ptap->api   = api;
966   ptap->apj   = apj;
967   ptap->duplicate = Cmpi->ops->duplicate;
968   ptap->destroy   = Cmpi->ops->destroy;
969 
970   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
971   Cmpi->assembled        = PETSC_FALSE;
972   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
973   Cmpi->ops->duplicate   = MatDuplicate_MPIAIJ_MatPtAP;
974   Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap;
975   *C                     = Cmpi;
976 
977   /* flag 'scalable' determines which implementations to be used:
978        0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa;
979        1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */
980   /* set default scalable */
981   ptap->scalable = PETSC_FALSE; /* PETSC_TRUE; */
982 
983   ierr = PetscOptionsGetBool(((PetscObject)Cmpi)->options,((PetscObject)Cmpi)->prefix,"-matptap_scalable",&ptap->scalable,NULL);CHKERRQ(ierr);
984   if (!ptap->scalable) {  /* Do dense axpy */
985     ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
986   } else {
987     ierr = PetscCalloc1(ap_rmax+1,&ptap->apa);CHKERRQ(ierr);
988   }
989 
990 #if defined(PETSC_USE_INFO)
991   if (pti[pn] != 0) {
992     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
993     ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
994   } else {
995     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
996   }
997 #endif
998   PetscFunctionReturn(0);
999 }
1000 
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; 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