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