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