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