xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision a2f3521de2a9d7869700f4c17e26c23fcfeaa6f6)
182412ba7SHong Zhang 
242c57489SHong Zhang /*
342c57489SHong Zhang   Defines projective product routines where A is a MPIAIJ matrix
442c57489SHong Zhang           C = P^T * A * P
542c57489SHong Zhang */
642c57489SHong Zhang 
7c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>   /*I "petscmat.h" I*/
8c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
9c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h>
10c6db04a5SJed Brown #include <petscbt.h>
1142c57489SHong Zhang 
12fe615159SHong Zhang /*
13fe615159SHong Zhang #define DEBUG_MATPTAP
14fe615159SHong Zhang  */
15fe615159SHong Zhang 
1609573ac7SBarry Smith extern PetscErrorCode MatDestroy_MPIAIJ(Mat);
17cc6cb787SHong Zhang #undef __FUNCT__
18f8487c73SHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_PtAP"
19f8487c73SHong Zhang PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A)
20cc6cb787SHong Zhang {
21cc6cb787SHong Zhang   PetscErrorCode       ierr;
22f8487c73SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data;
23f8487c73SHong Zhang   Mat_PtAPMPI          *ptap=a->ptap;
24cc6cb787SHong Zhang 
25cc6cb787SHong Zhang   PetscFunctionBegin;
26f8487c73SHong Zhang   if (ptap){
27c8b0795cSMark F. Adams     Mat_Merge_SeqsToMPI  *merge=ptap->merge;
28b7f45c76SHong Zhang     ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
29f8487c73SHong Zhang     ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
30a1a4e84aSHong Zhang     ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
31a1a4e84aSHong Zhang     ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
32c5af039cSHong Zhang     ierr = MatDestroy(&ptap->A_loc);CHKERRQ(ierr); /* used by MatTransposeMatMult() */
33c5af039cSHong Zhang     if (ptap->api){ierr = PetscFree(ptap->api);CHKERRQ(ierr);}
34c5af039cSHong Zhang     if (ptap->apj){ierr = PetscFree(ptap->apj);CHKERRQ(ierr);}
35c8b0795cSMark F. Adams     if (merge) {
36cc6cb787SHong Zhang       ierr = PetscFree(merge->id_r);CHKERRQ(ierr);
37cc6cb787SHong Zhang       ierr = PetscFree(merge->len_s);CHKERRQ(ierr);
38cc6cb787SHong Zhang       ierr = PetscFree(merge->len_r);CHKERRQ(ierr);
39cc6cb787SHong Zhang       ierr = PetscFree(merge->bi);CHKERRQ(ierr);
40cc6cb787SHong Zhang       ierr = PetscFree(merge->bj);CHKERRQ(ierr);
41c05d87d6SBarry Smith       ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr);
42cc6cb787SHong Zhang       ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr);
43c05d87d6SBarry Smith       ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr);
44cc6cb787SHong Zhang       ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr);
4505b42c5fSBarry Smith       ierr = PetscFree(merge->coi);CHKERRQ(ierr);
4605b42c5fSBarry Smith       ierr = PetscFree(merge->coj);CHKERRQ(ierr);
4705b42c5fSBarry Smith       ierr = PetscFree(merge->owners_co);CHKERRQ(ierr);
486bf464f9SBarry Smith       ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr);
49dce485f0SBarry Smith       ierr = merge->destroy(A);CHKERRQ(ierr);
50f8487c73SHong Zhang       ierr = PetscFree(ptap->merge);CHKERRQ(ierr);
51bf0cc555SLisandro Dalcin     }
52f8487c73SHong Zhang     ierr = PetscFree(ptap);CHKERRQ(ierr);
53c8b0795cSMark F. Adams   }
54c8b0795cSMark F. Adams 
55cc6cb787SHong Zhang   PetscFunctionReturn(0);
56cc6cb787SHong Zhang }
57cc6cb787SHong Zhang 
58cc6cb787SHong Zhang #undef __FUNCT__
59cc6cb787SHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatPtAP"
60f0c0a3a6SBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M)
61f0c0a3a6SBarry Smith {
62cc6cb787SHong Zhang   PetscErrorCode       ierr;
6311a6856fSHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data;
6411a6856fSHong Zhang   Mat_PtAPMPI          *ptap = a->ptap;
6511a6856fSHong Zhang   Mat_Merge_SeqsToMPI  *merge = ptap->merge;
66cc6cb787SHong Zhang 
67cc6cb787SHong Zhang   PetscFunctionBegin;
68dce485f0SBarry Smith   ierr = (*merge->duplicate)(A,op,M);CHKERRQ(ierr);
6911a6856fSHong Zhang   (*M)->ops->destroy   = merge->destroy;
7011a6856fSHong Zhang   (*M)->ops->duplicate = merge->duplicate;
71cc6cb787SHong Zhang   PetscFunctionReturn(0);
72cc6cb787SHong Zhang }
73cc6cb787SHong Zhang 
7442c57489SHong Zhang #undef __FUNCT__
757a7894deSKris Buschelman #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ"
767a7894deSKris Buschelman PetscErrorCode MatPtAPSymbolic_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
7742c57489SHong Zhang {
7842c57489SHong Zhang   PetscErrorCode ierr;
7942c57489SHong Zhang 
8042c57489SHong Zhang   PetscFunctionBegin;
8165e19b50SBarry Smith   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);
827a7894deSKris Buschelman   ierr = (*P->ops->ptapsymbolic_mpiaij)(A,P,fill,C);CHKERRQ(ierr);
837a7894deSKris Buschelman   PetscFunctionReturn(0);
847a7894deSKris Buschelman }
857a7894deSKris Buschelman 
867a7894deSKris Buschelman #undef __FUNCT__
877a7894deSKris Buschelman #define __FUNCT__ "MatPtAPNumeric_MPIAIJ"
887a7894deSKris Buschelman PetscErrorCode MatPtAPNumeric_MPIAIJ(Mat A,Mat P,Mat C)
897a7894deSKris Buschelman {
907a7894deSKris Buschelman   PetscErrorCode ierr;
917a7894deSKris Buschelman 
927a7894deSKris Buschelman   PetscFunctionBegin;
9365e19b50SBarry Smith   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);
947a7894deSKris Buschelman   ierr = (*P->ops->ptapnumeric_mpiaij)(A,P,C);CHKERRQ(ierr);
9542c57489SHong Zhang   PetscFunctionReturn(0);
9642c57489SHong Zhang }
9742c57489SHong Zhang 
9842c57489SHong Zhang #undef __FUNCT__
9942c57489SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ"
10042c57489SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
10142c57489SHong Zhang {
10242c57489SHong Zhang   PetscErrorCode       ierr;
10377584fe6SHong Zhang   Mat                  Cmpi;
104f8487c73SHong Zhang   Mat_PtAPMPI          *ptap;
105a1a86e44SBarry Smith   PetscFreeSpaceList   free_space=PETSC_NULL,current_space=PETSC_NULL;
106f8487c73SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
10742c57489SHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
10842c57489SHong Zhang   Mat_SeqAIJ           *p_loc,*p_oth;
109ded4f561SHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ;
11077584fe6SHong Zhang   PetscInt             *adi=ad->i,*aj,*aoi=ao->i,nnz;
11113f74950SBarry Smith   PetscInt             nlnk,*lnk,*owners_co,*coi,*coj,i,k,pnz,row;
112d16ca5f0SHong Zhang   PetscInt             am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n;
11342c57489SHong Zhang   PetscBT              lnkbt;
1147adad957SLisandro Dalcin   MPI_Comm             comm=((PetscObject)A)->comm;
115529bc5dcSHong Zhang   PetscMPIInt          size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
11642c57489SHong Zhang   PetscInt             **buf_rj,**buf_ri,**buf_ri_k;
11742c57489SHong Zhang   PetscInt             len,proc,*dnz,*onz,*owners;
118ded4f561SHong Zhang   PetscInt             nzi,*bi,*bj;
11942c57489SHong Zhang   PetscInt             nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
120ded4f561SHong Zhang   MPI_Request          *swaits,*rwaits;
121ded4f561SHong Zhang   MPI_Status           *sstatus,rstatus;
12242c57489SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
12377584fe6SHong Zhang   PetscInt             *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0;
124d16ca5f0SHong Zhang   PetscReal            afill=1.0,afill_tmp;
12508cb4532SHong Zhang   PetscInt             rstart = P->cmap->rstart,rmax;
12608cb4532SHong Zhang   PetscScalar          *vals;
12742c57489SHong Zhang 
12842c57489SHong Zhang   PetscFunctionBegin;
129c5af039cSHong Zhang   /* check if matrix local sizes are compatible */
130c5af039cSHong Zhang   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend){
131c5af039cSHong Zhang     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);
132c5af039cSHong Zhang   }
133c5af039cSHong Zhang   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend){
134c5af039cSHong Zhang     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);
135c5af039cSHong Zhang   }
136c5af039cSHong Zhang 
13783445d95SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
13883445d95SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
13983445d95SHong Zhang 
140f8487c73SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
141f8487c73SHong Zhang   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
142f8487c73SHong Zhang   ptap->reuse    = MAT_INITIAL_MATRIX;
1436c6e00e0SHong Zhang 
14477584fe6SHong Zhang 
1456c6e00e0SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
146b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
147fe615159SHong Zhang #if defined(DEBUG_MATPTAP)
148fe615159SHong Zhang    if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_oth is done...\n",rank);
149fe615159SHong Zhang #endif
150fe615159SHong Zhang 
1516c6e00e0SHong Zhang   /* get P_loc by taking all local rows of P */
152a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
153fe615159SHong Zhang #if defined(DEBUG_MATPTAP)
154fe615159SHong Zhang   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
155fe615159SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_loc is done...\n",rank);
156fe615159SHong Zhang #endif
1576c6e00e0SHong Zhang 
158a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
159a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
1606c6e00e0SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
1616c6e00e0SHong Zhang   pi_oth = p_oth->i; pj_oth = p_oth->j;
16242c57489SHong Zhang 
163fd0ff01cSHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
164fd0ff01cSHong Zhang   /*-------------------------------------------------------------------*/
16577584fe6SHong Zhang   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&api);CHKERRQ(ierr);
16682412ba7SHong Zhang   api[0] = 0;
16783445d95SHong Zhang 
16883445d95SHong Zhang   /* create and initialize a linked list */
16983445d95SHong Zhang   nlnk = pN+1;
170fe615159SHong Zhang #if defined(DEBUG_MATPTAP)
171fe615159SHong Zhang   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
172fe615159SHong Zhang   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]);
173fe615159SHong Zhang #endif
17483445d95SHong Zhang   ierr = PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
175fe615159SHong Zhang #if defined(DEBUG_MATPTAP)
176fe615159SHong Zhang   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
177fe615159SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call PetscLLCreate() is done\n",rank);
178fe615159SHong Zhang #endif
17983445d95SHong Zhang 
180d16ca5f0SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
181d16ca5f0SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
182f4a743e1SHong Zhang   current_space = free_space;
183f4a743e1SHong Zhang 
184f4a743e1SHong Zhang   for (i=0; i<am; i++) {
18582412ba7SHong Zhang     apnz = 0;
186f4a743e1SHong Zhang     /* diagonal portion of A */
187ded4f561SHong Zhang     nzi = adi[i+1] - adi[i];
18877584fe6SHong Zhang     aj  = ad->j + adi[i];
189ded4f561SHong Zhang     for (j=0; j<nzi; j++){
19077584fe6SHong Zhang       row  = aj[j];
19182412ba7SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
192ded4f561SHong Zhang       Jptr = pj_loc + pi_loc[row];
19383445d95SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
194dadf0e6bSHong Zhang       ierr = PetscLLAddSorted(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
19582412ba7SHong Zhang       apnz += nlnk;
196f4a743e1SHong Zhang     }
197f4a743e1SHong Zhang     /* off-diagonal portion of A */
198ded4f561SHong Zhang     nzi = aoi[i+1] - aoi[i];
19977584fe6SHong Zhang     aj  = ao->j + aoi[i];
200ded4f561SHong Zhang     for (j=0; j<nzi; j++){
20177584fe6SHong Zhang       row = aj[j];
20282412ba7SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
203ded4f561SHong Zhang       Jptr  = pj_oth + pi_oth[row];
204dadf0e6bSHong Zhang       ierr = PetscLLAddSorted(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
20582412ba7SHong Zhang       apnz += nlnk;
206f4a743e1SHong Zhang     }
207f4a743e1SHong Zhang 
20882412ba7SHong Zhang     api[i+1] = api[i] + apnz;
20977584fe6SHong Zhang     if (ap_rmax < apnz) ap_rmax = apnz;
210f4a743e1SHong Zhang 
21183445d95SHong Zhang     /* if free space is not available, double the total space in the list */
21282412ba7SHong Zhang     if (current_space->local_remaining<apnz) {
2132ba1acd5SHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
214f2b054eeSHong Zhang       nspacedouble++;
215f4a743e1SHong Zhang     }
216f4a743e1SHong Zhang 
217f4a743e1SHong Zhang     /* Copy data into free space, then initialize lnk */
21882412ba7SHong Zhang     ierr = PetscLLClean(pN,pN,apnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
21982412ba7SHong Zhang     current_space->array           += apnz;
22082412ba7SHong Zhang     current_space->local_used      += apnz;
22182412ba7SHong Zhang     current_space->local_remaining -= apnz;
222f4a743e1SHong Zhang   }
22382412ba7SHong Zhang   /* Allocate space for apj, initialize apj, and */
224f4a743e1SHong Zhang   /* destroy list of free space and other temporary array(s) */
22577584fe6SHong Zhang   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&apj);CHKERRQ(ierr);
22677584fe6SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
227d16ca5f0SHong Zhang   afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]);
228d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
229f4a743e1SHong Zhang 
230ded4f561SHong Zhang   /* determine symbolic Co=(p->B)^T*AP - send to others */
231ded4f561SHong Zhang   /*----------------------------------------------------*/
232de0260b3SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
23342c57489SHong Zhang 
234ded4f561SHong Zhang   /* then, compute symbolic Co = (p->B)^T*AP */
235d0f46423SBarry Smith   pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors
23683445d95SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
237de0260b3SHong Zhang   ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr);
238de0260b3SHong Zhang   coi[0] = 0;
23942c57489SHong Zhang 
240d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */
241d16ca5f0SHong Zhang   nnz           = fill*(poti[pon] + api[am]);
242a1a86e44SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);
24342c57489SHong Zhang   current_space = free_space;
24442c57489SHong Zhang 
245de0260b3SHong Zhang   for (i=0; i<pon; i++) {
246ded4f561SHong Zhang     nnz = 0;
24782412ba7SHong Zhang     pnz = poti[i+1] - poti[i];
24877584fe6SHong Zhang     ptJ = potj + poti[i];
24977584fe6SHong Zhang     for (j=0; j<pnz; j++){
25077584fe6SHong Zhang       row  = ptJ[j]; /* row of AP == col of Pot */
25182412ba7SHong Zhang       apnz = api[row+1] - api[row];
252ded4f561SHong Zhang       Jptr = apj + api[row];
25383445d95SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
254dadf0e6bSHong Zhang       ierr = PetscLLAddSorted(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
255ded4f561SHong Zhang       nnz += nlnk;
25642c57489SHong Zhang     }
25742c57489SHong Zhang 
25883445d95SHong Zhang     /* If free space is not available, double the total space in the list */
259ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
2604238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
261d16ca5f0SHong Zhang       nspacedouble++;
26242c57489SHong Zhang     }
26342c57489SHong Zhang 
26442c57489SHong Zhang     /* Copy data into free space, and zero out denserows */
265ded4f561SHong Zhang     ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
266ded4f561SHong Zhang     current_space->array           += nnz;
267ded4f561SHong Zhang     current_space->local_used      += nnz;
268ded4f561SHong Zhang     current_space->local_remaining -= nnz;
269ded4f561SHong Zhang     coi[i+1] = coi[i] + nnz;
27042c57489SHong Zhang   }
271de0260b3SHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr);
272a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
273d16ca5f0SHong Zhang   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]);
274d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
275de0260b3SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
27642c57489SHong Zhang 
277fe615159SHong Zhang #if defined(DEBUG_MATPTAP)
278fe615159SHong Zhang   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
279fe615159SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] local Co is done\n",rank);
280fe615159SHong Zhang #endif
281fe615159SHong Zhang 
282ded4f561SHong Zhang   /* send j-array (coj) of Co to other processors */
283ded4f561SHong Zhang   /*----------------------------------------------*/
28442c57489SHong Zhang   /* determine row ownership */
28583445d95SHong Zhang   ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr);
28626283091SBarry Smith   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
2877a2fc3feSBarry Smith   merge->rowmap->n = pn;
2887a2fc3feSBarry Smith   merge->rowmap->bs = 1;
28926283091SBarry Smith   ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
2907a2fc3feSBarry Smith   owners = merge->rowmap->range;
29142c57489SHong Zhang 
29242c57489SHong Zhang   /* determine the number of messages to send, their lengths */
29342c57489SHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr);
29483445d95SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
29542c57489SHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr);
29642c57489SHong Zhang   len_s = merge->len_s;
29742c57489SHong Zhang   merge->nsend = 0;
29883445d95SHong Zhang 
299de0260b3SHong Zhang   ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr);
300de0260b3SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
301de0260b3SHong Zhang 
30283445d95SHong Zhang   proc = 0;
303de0260b3SHong Zhang   for (i=0; i<pon; i++){
304de0260b3SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
305de0260b3SHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
306de0260b3SHong Zhang     len_s[proc] += coi[i+1] - coi[i];
30783445d95SHong Zhang   }
308de0260b3SHong Zhang 
309de0260b3SHong Zhang   len   = 0;  /* max length of buf_si[] */
310de0260b3SHong Zhang   owners_co[0] = 0;
31142c57489SHong Zhang   for (proc=0; proc<size; proc++){
312de0260b3SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
31383445d95SHong Zhang     if (len_si[proc]){
31442c57489SHong Zhang       merge->nsend++;
31583445d95SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
31642c57489SHong Zhang       len += len_si[proc];
31742c57489SHong Zhang     }
31842c57489SHong Zhang   }
31942c57489SHong Zhang 
320ded4f561SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
32142c57489SHong Zhang   ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
32242c57489SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
32342c57489SHong Zhang 
324ded4f561SHong Zhang   /* post the Irecv and Isend of coj */
325529bc5dcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
326529bc5dcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
327ded4f561SHong Zhang   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr);
32842c57489SHong Zhang   for (proc=0, k=0; proc<size; proc++){
32942c57489SHong Zhang     if (!len_s[proc]) continue;
330de0260b3SHong Zhang     i = owners_co[proc];
331529bc5dcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
33242c57489SHong Zhang     k++;
33342c57489SHong Zhang   }
33442c57489SHong Zhang 
335ded4f561SHong Zhang   /* receives and sends of coj are complete */
336ded4f561SHong Zhang   ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr);
33777584fe6SHong Zhang   for (i=0; i<merge->nrecv; i++){
338c280ed6aSJed Brown     PetscMPIInt icompleted;
339c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
340ded4f561SHong Zhang   }
341ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
3420c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
34382412ba7SHong Zhang 
344ded4f561SHong Zhang   /* send and recv coi */
345ded4f561SHong Zhang   /*-------------------*/
346529bc5dcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
347529bc5dcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
34842c57489SHong Zhang   ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr);
34942c57489SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
35042c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++){
35142c57489SHong Zhang     if (!len_s[proc]) continue;
35242c57489SHong Zhang     /* form outgoing message for i-structure:
35342c57489SHong Zhang          buf_si[0]:                 nrows to be sent
35442c57489SHong Zhang                [1:nrows]:           row index (global)
35542c57489SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
35642c57489SHong Zhang     */
35742c57489SHong Zhang     /*-------------------------------------------*/
35842c57489SHong Zhang     nrows = len_si[proc]/2 - 1;
35942c57489SHong Zhang     buf_si_i    = buf_si + nrows+1;
36042c57489SHong Zhang     buf_si[0]   = nrows;
36142c57489SHong Zhang     buf_si_i[0] = 0;
36242c57489SHong Zhang     nrows = 0;
363de0260b3SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++){
364ded4f561SHong Zhang       nzi = coi[i+1] - coi[i];
365ded4f561SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
366de0260b3SHong Zhang       buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */
36742c57489SHong Zhang       nrows++;
36842c57489SHong Zhang     }
369529bc5dcSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
37042c57489SHong Zhang     k++;
37142c57489SHong Zhang     buf_si += len_si[proc];
37242c57489SHong Zhang   }
373ded4f561SHong Zhang   i = merge->nrecv;
374ded4f561SHong Zhang   while (i--) {
375c280ed6aSJed Brown     PetscMPIInt icompleted;
376c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
377ded4f561SHong Zhang   }
378ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
3790c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
380f2b054eeSHong Zhang   /*
381ae15b995SBarry Smith   ierr = PetscInfo2(A,"nsend: %d, nrecv: %d\n",merge->nsend,merge->nrecv);CHKERRQ(ierr);
38242c57489SHong Zhang   for (i=0; i<merge->nrecv; i++){
383ae15b995SBarry Smith     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);
38442c57489SHong Zhang   }
385f2b054eeSHong Zhang   */
38642c57489SHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
38742c57489SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
388ded4f561SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
389ded4f561SHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
39042c57489SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
39142c57489SHong Zhang 
392fe615159SHong Zhang #if defined(DEBUG_MATPTAP)
393fe615159SHong Zhang   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
394fe615159SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Co is done\n",rank);
395fe615159SHong Zhang #endif
396fe615159SHong Zhang 
397de0260b3SHong Zhang   /* compute the local portion of C (mpi mat) */
398de0260b3SHong Zhang   /*------------------------------------------*/
399ded4f561SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
400ded4f561SHong Zhang 
40142c57489SHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
40242c57489SHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
40342c57489SHong Zhang   bi[0] = 0;
40442c57489SHong Zhang 
405d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
406d16ca5f0SHong Zhang   nnz           = fill*(pi_loc[pm] + api[am]);
407a1a86e44SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);
40842c57489SHong Zhang   current_space = free_space;
40942c57489SHong Zhang 
410687f1162SBarry Smith   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
41142c57489SHong Zhang   for (k=0; k<merge->nrecv; k++){
41242c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
41342c57489SHong Zhang     nrows       = *buf_ri_k[k];
41442c57489SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
41542c57489SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
41642c57489SHong Zhang   }
41742c57489SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
41808cb4532SHong Zhang   rmax = 0;
41942c57489SHong Zhang   for (i=0; i<pn; i++) {
420ded4f561SHong Zhang     /* add pdt[i,:]*AP into lnk */
421ded4f561SHong Zhang     nnz = 0;
422ded4f561SHong Zhang     pnz = pdti[i+1] - pdti[i];
42377584fe6SHong Zhang     ptJ = pdtj + pdti[i];
42477584fe6SHong Zhang     for (j=0; j<pnz; j++){
42577584fe6SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
426ded4f561SHong Zhang       apnz = api[row+1] - api[row];
427ded4f561SHong Zhang       Jptr = apj + api[row];
428ded4f561SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
429dadf0e6bSHong Zhang       ierr = PetscLLAddSorted(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
430ded4f561SHong Zhang       nnz += nlnk;
431ded4f561SHong Zhang     }
43242c57489SHong Zhang     /* add received col data into lnk */
43342c57489SHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
43442c57489SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
435ded4f561SHong Zhang         nzi = *(nextci[k]+1) - *nextci[k];
436ded4f561SHong Zhang         Jptr  = buf_rj[k] + *nextci[k];
437dadf0e6bSHong Zhang         ierr = PetscLLAddSorted(nzi,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
438ded4f561SHong Zhang         nnz += nlnk;
43942c57489SHong Zhang         nextrow[k]++; nextci[k]++;
44042c57489SHong Zhang       }
44142c57489SHong Zhang     }
44242c57489SHong Zhang 
44342c57489SHong Zhang     /* if free space is not available, make more free space */
444ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
4454238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
446d16ca5f0SHong Zhang       nspacedouble++;
44742c57489SHong Zhang     }
44842c57489SHong Zhang     /* copy data into free space, then initialize lnk */
449ded4f561SHong Zhang     ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
450ded4f561SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
451ded4f561SHong Zhang     current_space->array           += nnz;
452ded4f561SHong Zhang     current_space->local_used      += nnz;
453ded4f561SHong Zhang     current_space->local_remaining -= nnz;
454ded4f561SHong Zhang     bi[i+1] = bi[i] + nnz;
45508cb4532SHong Zhang     if (nnz > rmax) rmax = nnz;
45642c57489SHong Zhang   }
457ded4f561SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
4580572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
45942c57489SHong Zhang 
46042c57489SHong Zhang   ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
461a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
462d16ca5f0SHong Zhang   afill_tmp = (PetscReal)bi[pn]/(pi_loc[pm] + api[am]);
463d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
46442c57489SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
46542c57489SHong Zhang 
46608cb4532SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part - check runex55_SA ?  */
46708cb4532SHong Zhang   /*------------------------------------------------------------------------------------------------------*/
46808cb4532SHong Zhang   ierr = PetscMalloc((rmax+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
46908cb4532SHong Zhang   ierr = PetscMemzero(vals,rmax*sizeof(PetscScalar));CHKERRQ(ierr);
47008cb4532SHong Zhang 
47177584fe6SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
47277584fe6SHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
473*a2f3521dSMark F. Adams   ierr = MatSetBlockSizes(Cmpi,P->cmap->bs,P->cmap->bs);CHKERRQ(ierr);
47477584fe6SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
47577584fe6SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
47642c57489SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
477*a2f3521dSMark F. Adams 
47808cb4532SHong Zhang   for (i=0; i<pn; i++){
47908cb4532SHong Zhang     row = i + rstart;
48008cb4532SHong Zhang     nnz = bi[i+1] - bi[i];
48108cb4532SHong Zhang     Jptr = bj + bi[i];
48208cb4532SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
48308cb4532SHong Zhang   }
48408cb4532SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
48508cb4532SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
48608cb4532SHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
48742c57489SHong Zhang 
48842c57489SHong Zhang   merge->bi            = bi;
48942c57489SHong Zhang   merge->bj            = bj;
490de0260b3SHong Zhang   merge->coi           = coi;
491de0260b3SHong Zhang   merge->coj           = coj;
49242c57489SHong Zhang   merge->buf_ri        = buf_ri;
49342c57489SHong Zhang   merge->buf_rj        = buf_rj;
494de0260b3SHong Zhang   merge->owners_co     = owners_co;
49577584fe6SHong Zhang   merge->destroy       = Cmpi->ops->destroy;
49677584fe6SHong Zhang   merge->duplicate     = Cmpi->ops->duplicate;
497cc6cb787SHong Zhang 
49877584fe6SHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
49908cb4532SHong Zhang   /* Cmpi->assembled      = PETSC_FALSE; */
50077584fe6SHong Zhang   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_PtAP;
50177584fe6SHong Zhang   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP;
50242c57489SHong Zhang 
50377584fe6SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
50477584fe6SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
505f8487c73SHong Zhang   c->ptap        = ptap;
50677584fe6SHong Zhang   ptap->api      = api;
50777584fe6SHong Zhang   ptap->apj      = apj;
508f8487c73SHong Zhang   ptap->merge    = merge;
509d6ab1dc0SHong Zhang   ptap->rmax     = ap_rmax;
510f8487c73SHong Zhang 
51177584fe6SHong Zhang   *C = Cmpi;
512f2b054eeSHong Zhang #if defined(PETSC_USE_INFO)
513f2b054eeSHong Zhang   if (bi[pn] != 0) {
51477584fe6SHong Zhang     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
51577584fe6SHong Zhang     ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
516f2b054eeSHong Zhang   } else {
51777584fe6SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
518f2b054eeSHong Zhang   }
519f2b054eeSHong Zhang #endif
52042c57489SHong Zhang   PetscFunctionReturn(0);
52142c57489SHong Zhang }
52242c57489SHong Zhang 
52342c57489SHong Zhang #undef __FUNCT__
52442c57489SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ"
52542c57489SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
52642c57489SHong Zhang {
52742c57489SHong Zhang   PetscErrorCode       ierr;
52842c57489SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
529f8487c73SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
53042c57489SHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
531de0260b3SHong Zhang   Mat_SeqAIJ           *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
53242c57489SHong Zhang   Mat_SeqAIJ           *p_loc,*p_oth;
533f8487c73SHong Zhang   Mat_PtAPMPI          *ptap;
5341d633602SHong Zhang   PetscInt             *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp;
53582412ba7SHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj;
53682412ba7SHong Zhang   PetscInt             i,j,k,anz,pnz,apnz,nextap,row,*cj;
537e900a757SHong Zhang   MatScalar            *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp;
538d0f46423SBarry Smith   PetscInt             am=A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
5397adad957SLisandro Dalcin   MPI_Comm             comm=((PetscObject)C)->comm;
54042c57489SHong Zhang   PetscMPIInt          size,rank,taga,*len_s;
54182412ba7SHong Zhang   PetscInt             *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
54242c57489SHong Zhang   PetscInt             **buf_ri,**buf_rj;
54350f5bf86SHong Zhang   PetscInt             cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */
54442c57489SHong Zhang   MPI_Request          *s_waits,*r_waits;
54542c57489SHong Zhang   MPI_Status           *status;
54682412ba7SHong Zhang   MatScalar            **abuf_r,*ba_i,*pA,*coa,*ba;
54782412ba7SHong Zhang   PetscInt             *api,*apj,*coi,*coj;
548d0f46423SBarry Smith   PetscInt             *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend;
549a9d06231SHong Zhang   PetscInt             sparse_axpy;
55042c57489SHong Zhang 
55142c57489SHong Zhang   PetscFunctionBegin;
55242c57489SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
55342c57489SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
55442c57489SHong Zhang 
555f8487c73SHong Zhang   ptap  = c->ptap;
5561c4805f8SJed Brown   if (!ptap) SETERRQ(((PetscObject)C)->comm,PETSC_ERR_ARG_INCOMP,"MatPtAP() has not been called to create matrix C yet, cannot use MAT_REUSE_MATRIX");
557f8487c73SHong Zhang   merge = ptap->merge;
5586c6e00e0SHong Zhang 
559a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
560f9473cf7SHong Zhang   /*--------------------------------------------------*/
561f8487c73SHong Zhang   if (ptap->reuse == MAT_INITIAL_MATRIX){
562f8487c73SHong Zhang     ptap->reuse = MAT_REUSE_MATRIX;
5636c6e00e0SHong Zhang   } else { /* update numerical values of P_oth and P_loc */
564b7f45c76SHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
565a1a4e84aSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
5666c6e00e0SHong Zhang   }
567f8487c73SHong Zhang 
568f9473cf7SHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
569f9473cf7SHong Zhang   /*--------------------------------------------------------------*/
57042c57489SHong Zhang   /* get data from symbolic products */
571a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
572a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
573a9d06231SHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a;
57442c57489SHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
57542c57489SHong Zhang 
576de0260b3SHong Zhang   coi = merge->coi; coj = merge->coj;
577de0260b3SHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr);
578de0260b3SHong Zhang   ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr);
579de0260b3SHong Zhang 
580de0260b3SHong Zhang   bi     = merge->bi; bj = merge->bj;
5817a2fc3feSBarry Smith   owners = merge->rowmap->range;
582de0260b3SHong Zhang   ierr   = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
583de0260b3SHong Zhang   ierr   = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr);
58442c57489SHong Zhang 
585a1a4e84aSHong Zhang   api = ptap->api; apj = ptap->apj;
5860496f32cSHong Zhang   /* flag 'sparse_axpy' determines which implementations to be used:
5870496f32cSHong Zhang        0: do dense axpy in MatPtAPNumeric() - fastest, but requires storage of a dense array apa; (default)
5880496f32cSHong Zhang        1: do one sparse axpy - uses same memory as sparse_axpy=0 and might execute less flops
5890496f32cSHong Zhang           (apnz vs. cnz in the outerproduct), slower than case '0' when cnz is not too large than apnz;
5900496f32cSHong Zhang        2: do two sparse axpy in MatPtAPNumeric() - slowest, uses a sparse array apa */
591a9d06231SHong Zhang   /* set default sparse_axpy */
592e900a757SHong Zhang   sparse_axpy = 0;
5930496f32cSHong Zhang   ierr = PetscOptionsGetInt(PETSC_NULL,"-matptap_sparseaxpy",&sparse_axpy,PETSC_NULL);CHKERRQ(ierr);
594a9d06231SHong Zhang   if (sparse_axpy == 0){  /* Do not perform sparse axpy */
595a9d06231SHong Zhang     /*--------------------------------------------------*/
596a9d06231SHong Zhang     /* malloc apa to store dense row A[i,:]*P */
597a9d06231SHong Zhang     ierr = PetscMalloc((P->cmap->N)*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
598a9d06231SHong Zhang     ierr = PetscMemzero(apa,P->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
599a9d06231SHong Zhang 
600a9d06231SHong Zhang     for (i=0; i<am; i++) {
601a9d06231SHong Zhang       /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */
602a9d06231SHong Zhang       /*------------------------------------------------------------*/
603a9d06231SHong Zhang       apJ = apj + api[i];
604a9d06231SHong Zhang 
605a9d06231SHong Zhang       /* diagonal portion of A */
606a9d06231SHong Zhang       anz = adi[i+1] - adi[i];
607a9d06231SHong Zhang       adj = ad->j + adi[i];
608a9d06231SHong Zhang       ada = ad->a + adi[i];
609a9d06231SHong Zhang       for (j=0; j<anz; j++) {
610a9d06231SHong Zhang         row = adj[j];
611a9d06231SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
612a9d06231SHong Zhang         pj  = pj_loc + pi_loc[row];
613a9d06231SHong Zhang         pa  = pa_loc + pi_loc[row];
614a9d06231SHong Zhang 
615a9d06231SHong Zhang         /* perform dense axpy */
616e900a757SHong Zhang         valtmp = ada[j];
617a9d06231SHong Zhang         for (k=0; k<pnz; k++){
618e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
619a9d06231SHong Zhang         }
620a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
621a9d06231SHong Zhang       }
622a9d06231SHong Zhang 
623a9d06231SHong Zhang       /* off-diagonal portion of A */
624a9d06231SHong Zhang       anz = aoi[i+1] - aoi[i];
625a9d06231SHong Zhang       aoj = ao->j + aoi[i];
626a9d06231SHong Zhang       aoa = ao->a + aoi[i];
627a9d06231SHong Zhang       for (j=0; j<anz; j++) {
628a9d06231SHong Zhang         row = aoj[j];
629a9d06231SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
630a9d06231SHong Zhang         pj  = pj_oth + pi_oth[row];
631a9d06231SHong Zhang         pa  = pa_oth + pi_oth[row];
632a9d06231SHong Zhang 
633a9d06231SHong Zhang         /* perform dense axpy */
634e900a757SHong Zhang         valtmp = aoa[j];
635a9d06231SHong Zhang         for (k=0; k<pnz; k++){
636e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
637a9d06231SHong Zhang         }
638a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
639a9d06231SHong Zhang       }
640a9d06231SHong Zhang 
641a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
642a9d06231SHong Zhang       /*--------------------------------------------------------------*/
643a9d06231SHong Zhang       apnz = api[i+1] - api[i];
644a9d06231SHong Zhang       /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */
645a9d06231SHong Zhang       pnz = po->i[i+1] - po->i[i];
646e900a757SHong Zhang       poJ = po->j + po->i[i];
647a9d06231SHong Zhang       pA  = po->a + po->i[i];
648a9d06231SHong Zhang       for (j=0; j<pnz; j++){
649e900a757SHong Zhang         row = poJ[j];
650e900a757SHong Zhang         cnz = coi[row+1] - coi[row];
651e900a757SHong Zhang         cj  = coj + coi[row];
652e900a757SHong Zhang         ca  = coa + coi[row];
653a9d06231SHong Zhang         /* perform dense axpy */
654e900a757SHong Zhang         valtmp = pA[j];
655a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
656e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
657a9d06231SHong Zhang         }
658a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
659a9d06231SHong Zhang       }
660a9d06231SHong Zhang 
661a9d06231SHong Zhang       /* put the value into Cd (diagonal part) */
662a9d06231SHong Zhang       pnz = pd->i[i+1] - pd->i[i];
663e900a757SHong Zhang       pdJ = pd->j + pd->i[i];
664a9d06231SHong Zhang       pA  = pd->a + pd->i[i];
665a9d06231SHong Zhang       for (j=0; j<pnz; j++){
666e900a757SHong Zhang         row = pdJ[j];
667e900a757SHong Zhang         cnz = bi[row+1] - bi[row];
668e900a757SHong Zhang         cj  = bj + bi[row];
669e900a757SHong Zhang         ca  = ba + bi[row];
670a9d06231SHong Zhang         /* perform dense axpy */
671e900a757SHong Zhang         valtmp = pA[j];
672a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
673e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
674a9d06231SHong Zhang         }
675a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
676a9d06231SHong Zhang       }
677a9d06231SHong Zhang 
678a9d06231SHong Zhang       /* zero the current row of A*P */
679a9d06231SHong Zhang       for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0;
680a9d06231SHong Zhang     }
681a9d06231SHong Zhang   } else if (sparse_axpy == 1){  /* Perform one sparse axpy */
682a9d06231SHong Zhang     /*------------------------------------------------------*/
6831d633602SHong Zhang     /* malloc apa to store dense row A[i,:]*P */
6841d633602SHong Zhang     ierr = PetscMalloc((P->cmap->N)*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
6851d633602SHong Zhang     ierr = PetscMemzero(apa,P->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
6861d633602SHong Zhang 
6871d633602SHong Zhang     for (i=0; i<am; i++) {
688f9473cf7SHong Zhang       /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */
689f9473cf7SHong Zhang       /*------------------------------------------------------------*/
6901d633602SHong Zhang       apJ = apj + api[i];
6911d633602SHong Zhang 
6921d633602SHong Zhang       /* diagonal portion of A */
6931d633602SHong Zhang       anz = adi[i+1] - adi[i];
6941d633602SHong Zhang       adj = ad->j + adi[i];
6951d633602SHong Zhang       ada = ad->a + adi[i];
6961d633602SHong Zhang       for (j=0; j<anz; j++) {
6971d633602SHong Zhang         row = adj[j];
6981d633602SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
6991d633602SHong Zhang         pj  = pj_loc + pi_loc[row];
7001d633602SHong Zhang         pa  = pa_loc + pi_loc[row];
7011d633602SHong Zhang 
7021d633602SHong Zhang         /* perform dense axpy */
703e900a757SHong Zhang         valtmp = ada[j];
7041d633602SHong Zhang         for (k=0; k<pnz; k++){
705e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
7061d633602SHong Zhang         }
7071d633602SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
7081d633602SHong Zhang       }
7091d633602SHong Zhang 
7101d633602SHong Zhang       /* off-diagonal portion of A */
7111d633602SHong Zhang       anz = aoi[i+1] - aoi[i];
7121d633602SHong Zhang       aoj = ao->j + aoi[i];
7131d633602SHong Zhang       aoa = ao->a + aoi[i];
7141d633602SHong Zhang       for (j=0; j<anz; j++) {
7151d633602SHong Zhang         row = aoj[j];
7161d633602SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
7171d633602SHong Zhang         pj  = pj_oth + pi_oth[row];
7181d633602SHong Zhang         pa  = pa_oth + pi_oth[row];
7191d633602SHong Zhang 
7201d633602SHong Zhang         /* perform dense axpy */
721e900a757SHong Zhang         valtmp = aoa[j];
7221d633602SHong Zhang         for (k=0; k<pnz; k++){
723e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
7241d633602SHong Zhang         }
7251d633602SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
7261d633602SHong Zhang       }
7271d633602SHong Zhang 
728f9473cf7SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
729f9473cf7SHong Zhang       /*--------------------------------------------------------------*/
7301d633602SHong Zhang       apnz = api[i+1] - api[i];
731f9473cf7SHong Zhang       /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */
732f9473cf7SHong Zhang       pnz = po->i[i+1] - po->i[i];
733e900a757SHong Zhang       poJ = po->j + po->i[i];
734a9d06231SHong Zhang       pA  = po->a + po->i[i];
7351d633602SHong Zhang       for (j=0; j<pnz; j++){
736e900a757SHong Zhang         row    = poJ[j];
737e900a757SHong Zhang         cj     = coj + coi[row];
738e900a757SHong Zhang         ca     = coa + coi[row];
739e900a757SHong Zhang         valtmp = pA[j];
7401d633602SHong Zhang         /* perform sparse axpy */
7411d633602SHong Zhang         nextap = 0;
7421d633602SHong Zhang         for (k=0; nextap<apnz; k++) {
7431d633602SHong Zhang           if (cj[k]==apJ[nextap]) { /* global column index */
744e900a757SHong Zhang             ca[k] += valtmp*apa[cj[k]]; nextap++;
7451d633602SHong Zhang           }
7461d633602SHong Zhang         }
7471d633602SHong Zhang         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
7481d633602SHong Zhang       }
749f9473cf7SHong Zhang 
750f9473cf7SHong Zhang       /* put the value into Cd (diagonal part) */
751f9473cf7SHong Zhang       pnz = pd->i[i+1] - pd->i[i];
752e900a757SHong Zhang       pdJ = pd->j + pd->i[i];
753a9d06231SHong Zhang       pA  = pd->a + pd->i[i];
754f9473cf7SHong Zhang       for (j=0; j<pnz; j++){
755e900a757SHong Zhang         row    = pdJ[j];
756e900a757SHong Zhang         cj     = bj + bi[row];
757e900a757SHong Zhang         ca     = ba + bi[row];
758e900a757SHong Zhang         valtmp = pA[j];
759f9473cf7SHong Zhang         /* perform sparse axpy */
760f9473cf7SHong Zhang         nextap = 0;
761f9473cf7SHong Zhang         for (k=0; nextap<apnz; k++) {
762f9473cf7SHong Zhang           if (cj[k]==apJ[nextap]) { /* global column index */
763e900a757SHong Zhang             ca[k] += valtmp*apa[cj[k]];
764a9d06231SHong Zhang             nextap++;
765f9473cf7SHong Zhang           }
766f9473cf7SHong Zhang         }
767f9473cf7SHong Zhang         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
768f9473cf7SHong Zhang       }
769f9473cf7SHong Zhang 
770f9473cf7SHong Zhang       /* zero the current row of A*P */
7711d633602SHong Zhang       for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0;
7721d633602SHong Zhang     }
7730496f32cSHong Zhang   } else if (sparse_axpy == 2){/* Perform two sparse axpy */
774a9d06231SHong Zhang     /*----------------------------------------------------*/
7751d633602SHong Zhang     /* malloc apa to store sparse row A[i,:]*P */
776d6ab1dc0SHong Zhang     ierr = PetscMalloc((ptap->rmax+1)*sizeof(MatScalar),&apa);CHKERRQ(ierr);
777d6ab1dc0SHong Zhang     ierr = PetscMemzero(apa,ptap->rmax*sizeof(MatScalar));CHKERRQ(ierr);
77842c57489SHong Zhang 
779a9d06231SHong Zhang     pA=pa_loc;
78042c57489SHong Zhang     for (i=0; i<am; i++) {
781f4a743e1SHong Zhang       /* form i-th sparse row of A*P */
78282412ba7SHong Zhang       apnz = api[i+1] - api[i];
78382412ba7SHong Zhang       apJ  = apj + api[i];
78442c57489SHong Zhang       /* diagonal portion of A */
78582412ba7SHong Zhang       anz = adi[i+1] - adi[i];
7861d633602SHong Zhang       adj = ad->j + adi[i];
7871d633602SHong Zhang       ada = ad->a + adi[i];
78882412ba7SHong Zhang       for (j=0; j<anz; j++) {
7891d633602SHong Zhang         row = adj[j];
79082412ba7SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
79182412ba7SHong Zhang         pj  = pj_loc + pi_loc[row];
79282412ba7SHong Zhang         pa  = pa_loc + pi_loc[row];
793e900a757SHong Zhang         valtmp = ada[j];
794f4a743e1SHong Zhang         nextp  = 0;
79582412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
79682412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
797e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
79842c57489SHong Zhang           }
79942c57489SHong Zhang         }
800dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
80142c57489SHong Zhang       }
80242c57489SHong Zhang       /* off-diagonal portion of A */
80382412ba7SHong Zhang       anz = aoi[i+1] - aoi[i];
8041d633602SHong Zhang       aoj = ao->j + aoi[i];
8051d633602SHong Zhang       aoa = ao->a + aoi[i];
80682412ba7SHong Zhang       for (j=0; j<anz; j++) {
8071d633602SHong Zhang         row = aoj[j];
80882412ba7SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
80982412ba7SHong Zhang         pj  = pj_oth + pi_oth[row];
81082412ba7SHong Zhang         pa  = pa_oth + pi_oth[row];
811e900a757SHong Zhang         valtmp = aoa[j];
812f4a743e1SHong Zhang         nextp  = 0;
81382412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
81482412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
815e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
81642c57489SHong Zhang           }
81742c57489SHong Zhang         }
818dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
81942c57489SHong Zhang       }
82042c57489SHong Zhang 
821a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
822a9d06231SHong Zhang       /*--------------------------------------------------------------*/
82382412ba7SHong Zhang       pnz = pi_loc[i+1] - pi_loc[i];
824f9473cf7SHong Zhang       pJ  = pj_loc + pi_loc[i];
82582412ba7SHong Zhang       for (j=0; j<pnz; j++) {
82642c57489SHong Zhang         nextap = 0;
827f9473cf7SHong Zhang         row    = pJ[j]; /* global index */
82882412ba7SHong Zhang         if (row < pcstart || row >=pcend) { /* put the value into Co */
829e900a757SHong Zhang           row = *poJ;
830e900a757SHong Zhang           cj  = coj + coi[row];
831e900a757SHong Zhang           ca  = coa + coi[row]; poJ++;
83282412ba7SHong Zhang         } else {                            /* put the value into Cd */
833e900a757SHong Zhang           row  = *pdJ;
834e900a757SHong Zhang           cj   = bj + bi[row];
835e900a757SHong Zhang           ca   = ba + bi[row]; pdJ++;
83642c57489SHong Zhang         }
837e900a757SHong Zhang         valtmp = pA[j];
83882412ba7SHong Zhang         for (k=0; nextap<apnz; k++) {
839e900a757SHong Zhang           if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++];
84042c57489SHong Zhang         }
841dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
842de0260b3SHong Zhang       }
8430496f32cSHong Zhang       pA += pnz;
84442c57489SHong Zhang       /* zero the current row info for A*P */
84582412ba7SHong Zhang       ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr);
84642c57489SHong Zhang     }
847a9d06231SHong Zhang   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"sparse_axpy only takes values 0, 1 and 2");
84842c57489SHong Zhang   ierr = PetscFree(apa);CHKERRQ(ierr);
84942c57489SHong Zhang 
850a9d06231SHong Zhang   /* 3) send and recv matrix values coa */
851a9d06231SHong Zhang   /*------------------------------------*/
85242c57489SHong Zhang   buf_ri = merge->buf_ri;
85342c57489SHong Zhang   buf_rj = merge->buf_rj;
85442c57489SHong Zhang   len_s  = merge->len_s;
855fc42d0c8SSatish Balay   ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
85642c57489SHong Zhang   ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
85742c57489SHong Zhang 
858a9d06231SHong Zhang   ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr);
85942c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++){
86042c57489SHong Zhang     if (!len_s[proc]) continue;
861de0260b3SHong Zhang     i = merge->owners_co[proc];
862de0260b3SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
86342c57489SHong Zhang     k++;
86442c57489SHong Zhang   }
8650c468ba9SBarry Smith   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
8660c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
86742c57489SHong Zhang 
868a9d06231SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
86942c57489SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
87082412ba7SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
87142c57489SHong Zhang 
872a9d06231SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
873a9d06231SHong Zhang   /*------------------------------------------------------*/
874687f1162SBarry Smith   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
87542c57489SHong Zhang   for (k=0; k<merge->nrecv; k++){
87642c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
87742c57489SHong Zhang     nrows       = *(buf_ri_k[k]);
87842c57489SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
87982412ba7SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
88042c57489SHong Zhang   }
88142c57489SHong Zhang 
882de0260b3SHong Zhang   for (i=0; i<cm; i++) {
88382412ba7SHong Zhang     row = owners[rank] + i; /* global row index of C_seq */
88442c57489SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
885de0260b3SHong Zhang     ba_i = ba + bi[i];
88682412ba7SHong Zhang     bnz  = bi[i+1] - bi[i];
88742c57489SHong Zhang     /* add received vals into ba */
88842c57489SHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
88942c57489SHong Zhang       /* i-th row */
89042c57489SHong Zhang       if (i == *nextrow[k]) {
89182412ba7SHong Zhang         cnz = *(nextci[k]+1) - *nextci[k];
89282412ba7SHong Zhang         cj  = buf_rj[k] + *(nextci[k]);
89382412ba7SHong Zhang         ca  = abuf_r[k] + *(nextci[k]);
89482412ba7SHong Zhang         nextcj = 0;
89582412ba7SHong Zhang         for (j=0; nextcj<cnz; j++){
89682412ba7SHong Zhang           if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */
89782412ba7SHong Zhang             ba_i[j] += ca[nextcj++];
89842c57489SHong Zhang           }
89942c57489SHong Zhang         }
90082412ba7SHong Zhang         nextrow[k]++; nextci[k]++;
901c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
90242c57489SHong Zhang       }
90342c57489SHong Zhang     }
90482412ba7SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
90542c57489SHong Zhang   }
90642c57489SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
90742c57489SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
90842c57489SHong Zhang 
909de0260b3SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
910c05d87d6SBarry Smith   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
91142c57489SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
9120572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
91342c57489SHong Zhang   PetscFunctionReturn(0);
91442c57489SHong Zhang }
915