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