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