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