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