xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision dadf0e6bf3fbf8619c99f30f0124ab4e08bd7289)
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 
1209573ac7SBarry Smith extern PetscErrorCode MatDestroy_MPIAIJ(Mat);
13cc6cb787SHong Zhang #undef __FUNCT__
14f8487c73SHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_PtAP"
15f8487c73SHong Zhang PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A)
16cc6cb787SHong Zhang {
17cc6cb787SHong Zhang   PetscErrorCode       ierr;
18f8487c73SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data;
19f8487c73SHong Zhang   Mat_PtAPMPI          *ptap=a->ptap;
20cc6cb787SHong Zhang 
21cc6cb787SHong Zhang   PetscFunctionBegin;
22f8487c73SHong Zhang   if (ptap){
23f8487c73SHong Zhang     ierr = PetscFree2(ptap->startsj,ptap->startsj_r);CHKERRQ(ierr);
24f8487c73SHong Zhang     ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
25a1a4e84aSHong Zhang     ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
26a1a4e84aSHong Zhang     ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
27a1a4e84aSHong Zhang     ierr = PetscFree(ptap->api);CHKERRQ(ierr);
28a1a4e84aSHong Zhang     ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
29f8487c73SHong Zhang   }
30f8487c73SHong Zhang   if (ptap->merge) {
31f8487c73SHong Zhang     Mat_Merge_SeqsToMPI  *merge=ptap->merge;
32cc6cb787SHong Zhang     ierr = PetscFree(merge->id_r);CHKERRQ(ierr);
33cc6cb787SHong Zhang     ierr = PetscFree(merge->len_s);CHKERRQ(ierr);
34cc6cb787SHong Zhang     ierr = PetscFree(merge->len_r);CHKERRQ(ierr);
35cc6cb787SHong Zhang     ierr = PetscFree(merge->bi);CHKERRQ(ierr);
36cc6cb787SHong Zhang     ierr = PetscFree(merge->bj);CHKERRQ(ierr);
37c05d87d6SBarry Smith     ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr);
38cc6cb787SHong Zhang     ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr);
39c05d87d6SBarry Smith     ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr);
40cc6cb787SHong Zhang     ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr);
4105b42c5fSBarry Smith     ierr = PetscFree(merge->coi);CHKERRQ(ierr);
4205b42c5fSBarry Smith     ierr = PetscFree(merge->coj);CHKERRQ(ierr);
4305b42c5fSBarry Smith     ierr = PetscFree(merge->owners_co);CHKERRQ(ierr);
446bf464f9SBarry Smith     ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr);
45dce485f0SBarry Smith     ierr = merge->destroy(A);CHKERRQ(ierr);
46f8487c73SHong Zhang     ierr = PetscFree(ptap->merge);CHKERRQ(ierr);
47bf0cc555SLisandro Dalcin   }
48f8487c73SHong Zhang   ierr = PetscFree(ptap);CHKERRQ(ierr);
49cc6cb787SHong Zhang   PetscFunctionReturn(0);
50cc6cb787SHong Zhang }
51cc6cb787SHong Zhang 
52cc6cb787SHong Zhang #undef __FUNCT__
53cc6cb787SHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatPtAP"
54f0c0a3a6SBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M)
55f0c0a3a6SBarry Smith {
56cc6cb787SHong Zhang   PetscErrorCode       ierr;
57cc6cb787SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
58776b82aeSLisandro Dalcin   PetscContainer       container;
59cc6cb787SHong Zhang 
60cc6cb787SHong Zhang   PetscFunctionBegin;
61cc6cb787SHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"MatMergeSeqsToMPI",(PetscObject *)&container);CHKERRQ(ierr);
62cc6cb787SHong Zhang   if (container) {
63776b82aeSLisandro Dalcin     ierr  = PetscContainerGetPointer(container,(void **)&merge);CHKERRQ(ierr);
64e7e72b3dSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
65dce485f0SBarry Smith   ierr = (*merge->duplicate)(A,op,M);CHKERRQ(ierr);
66dce485f0SBarry Smith   (*M)->ops->destroy   = merge->destroy;   /* =MatDestroy_MPIAIJ, *M doesn't duplicate A's container! */
67dce485f0SBarry Smith   (*M)->ops->duplicate = merge->duplicate; /* =MatDuplicate_ MPIAIJ */
68cc6cb787SHong Zhang   PetscFunctionReturn(0);
69cc6cb787SHong Zhang }
70cc6cb787SHong Zhang 
7142c57489SHong Zhang #undef __FUNCT__
727a7894deSKris Buschelman #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ"
737a7894deSKris Buschelman PetscErrorCode MatPtAPSymbolic_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
7442c57489SHong Zhang {
7542c57489SHong Zhang   PetscErrorCode ierr;
7642c57489SHong Zhang 
7742c57489SHong Zhang   PetscFunctionBegin;
7865e19b50SBarry 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);
797a7894deSKris Buschelman   ierr = (*P->ops->ptapsymbolic_mpiaij)(A,P,fill,C);CHKERRQ(ierr);
807a7894deSKris Buschelman   PetscFunctionReturn(0);
817a7894deSKris Buschelman }
827a7894deSKris Buschelman 
837a7894deSKris Buschelman #undef __FUNCT__
847a7894deSKris Buschelman #define __FUNCT__ "MatPtAPNumeric_MPIAIJ"
857a7894deSKris Buschelman PetscErrorCode MatPtAPNumeric_MPIAIJ(Mat A,Mat P,Mat C)
867a7894deSKris Buschelman {
877a7894deSKris Buschelman   PetscErrorCode ierr;
887a7894deSKris Buschelman 
897a7894deSKris Buschelman   PetscFunctionBegin;
9065e19b50SBarry 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);
917a7894deSKris Buschelman   ierr = (*P->ops->ptapnumeric_mpiaij)(A,P,C);CHKERRQ(ierr);
9242c57489SHong Zhang   PetscFunctionReturn(0);
9342c57489SHong Zhang }
9442c57489SHong Zhang 
9542c57489SHong Zhang #undef __FUNCT__
9642c57489SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ"
9742c57489SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
9842c57489SHong Zhang {
9942c57489SHong Zhang   PetscErrorCode       ierr;
10077584fe6SHong Zhang   Mat                  Cmpi;
101f8487c73SHong Zhang   Mat_PtAPMPI          *ptap;
102a1a86e44SBarry Smith   PetscFreeSpaceList   free_space=PETSC_NULL,current_space=PETSC_NULL;
103f8487c73SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
10442c57489SHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
10542c57489SHong Zhang   Mat_SeqAIJ           *p_loc,*p_oth;
106ded4f561SHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ;
10777584fe6SHong Zhang   PetscInt             *adi=ad->i,*aj,*aoi=ao->i,nnz;
10813f74950SBarry Smith   PetscInt             nlnk,*lnk,*owners_co,*coi,*coj,i,k,pnz,row;
109d16ca5f0SHong Zhang   PetscInt             am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n;
11042c57489SHong Zhang   PetscBT              lnkbt;
1117adad957SLisandro Dalcin   MPI_Comm             comm=((PetscObject)A)->comm;
112ded4f561SHong Zhang   PetscMPIInt          size,rank,tag,*len_si,*len_s,*len_ri;
11342c57489SHong Zhang   PetscInt             **buf_rj,**buf_ri,**buf_ri_k;
11442c57489SHong Zhang   PetscInt             len,proc,*dnz,*onz,*owners;
115ded4f561SHong Zhang   PetscInt             nzi,*bi,*bj;
11642c57489SHong Zhang   PetscInt             nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
117ded4f561SHong Zhang   MPI_Request          *swaits,*rwaits;
118ded4f561SHong Zhang   MPI_Status           *sstatus,rstatus;
11942c57489SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
12077584fe6SHong Zhang   PetscInt             *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0;
121d16ca5f0SHong Zhang   PetscReal            afill=1.0,afill_tmp;
122*dadf0e6bSHong Zhang   PetscLogDouble       t0,tf,etime=0.0,t00,tff,time_matupdate=0.0,time_Cd=0.0,time_AP=0.0,time_Co=0.0;
12342c57489SHong Zhang 
12442c57489SHong Zhang   PetscFunctionBegin;
12577584fe6SHong Zhang   ierr = PetscGetTime(&t00);CHKERRQ(ierr);
12683445d95SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
12783445d95SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
12883445d95SHong Zhang 
12977584fe6SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
130f8487c73SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
131f8487c73SHong Zhang   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
132f8487c73SHong Zhang   ptap->reuse    = MAT_INITIAL_MATRIX;
1336c6e00e0SHong Zhang 
13477584fe6SHong Zhang 
1356c6e00e0SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
136a1a4e84aSHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
1376c6e00e0SHong Zhang   /* get P_loc by taking all local rows of P */
138a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
1396c6e00e0SHong Zhang 
140a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
141a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
1426c6e00e0SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
1436c6e00e0SHong Zhang   pi_oth = p_oth->i; pj_oth = p_oth->j;
14442c57489SHong Zhang 
14577584fe6SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
14677584fe6SHong Zhang   time_matupdate += tf-t0;
14777584fe6SHong Zhang 
148fd0ff01cSHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
149fd0ff01cSHong Zhang   /*-------------------------------------------------------------------*/
15077584fe6SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
15177584fe6SHong Zhang   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&api);CHKERRQ(ierr);
15282412ba7SHong Zhang   api[0] = 0;
15383445d95SHong Zhang 
15483445d95SHong Zhang   /* create and initialize a linked list */
15583445d95SHong Zhang   nlnk = pN+1;
15683445d95SHong Zhang   ierr = PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
15783445d95SHong Zhang 
158d16ca5f0SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
159d16ca5f0SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
160f4a743e1SHong Zhang   current_space = free_space;
161f4a743e1SHong Zhang 
162f4a743e1SHong Zhang   for (i=0; i<am; i++) {
16382412ba7SHong Zhang     apnz = 0;
164f4a743e1SHong Zhang     /* diagonal portion of A */
165ded4f561SHong Zhang     nzi = adi[i+1] - adi[i];
16677584fe6SHong Zhang     aj  = ad->j + adi[i];
167ded4f561SHong Zhang     for (j=0; j<nzi; j++){
16877584fe6SHong Zhang       row  = aj[j];
16982412ba7SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
170ded4f561SHong Zhang       Jptr = pj_loc + pi_loc[row];
17183445d95SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
172*dadf0e6bSHong Zhang       ierr = PetscLLAddSorted(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
17382412ba7SHong Zhang       apnz += nlnk;
174f4a743e1SHong Zhang     }
175f4a743e1SHong Zhang     /* off-diagonal portion of A */
176ded4f561SHong Zhang     nzi = aoi[i+1] - aoi[i];
17777584fe6SHong Zhang     aj  = ao->j + aoi[i];
178ded4f561SHong Zhang     for (j=0; j<nzi; j++){
17977584fe6SHong Zhang       row = aj[j];
18082412ba7SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
181ded4f561SHong Zhang       Jptr  = pj_oth + pi_oth[row];
182*dadf0e6bSHong Zhang       ierr = PetscLLAddSorted(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
18382412ba7SHong Zhang       apnz += nlnk;
184f4a743e1SHong Zhang     }
185f4a743e1SHong Zhang 
18682412ba7SHong Zhang     api[i+1] = api[i] + apnz;
18777584fe6SHong Zhang     if (ap_rmax < apnz) ap_rmax = apnz;
188f4a743e1SHong Zhang 
18983445d95SHong Zhang     /* if free space is not available, double the total space in the list */
19082412ba7SHong Zhang     if (current_space->local_remaining<apnz) {
1912ba1acd5SHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
192f2b054eeSHong Zhang       nspacedouble++;
193f4a743e1SHong Zhang     }
194f4a743e1SHong Zhang 
195f4a743e1SHong Zhang     /* Copy data into free space, then initialize lnk */
19682412ba7SHong Zhang     ierr = PetscLLClean(pN,pN,apnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
19782412ba7SHong Zhang     current_space->array           += apnz;
19882412ba7SHong Zhang     current_space->local_used      += apnz;
19982412ba7SHong Zhang     current_space->local_remaining -= apnz;
200f4a743e1SHong Zhang   }
20182412ba7SHong Zhang   /* Allocate space for apj, initialize apj, and */
202f4a743e1SHong Zhang   /* destroy list of free space and other temporary array(s) */
20377584fe6SHong Zhang   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&apj);CHKERRQ(ierr);
20477584fe6SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
205d16ca5f0SHong Zhang   afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]);
206d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
207f4a743e1SHong Zhang 
20877584fe6SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
20977584fe6SHong Zhang   time_AP += tf-t0;
21077584fe6SHong Zhang 
211ded4f561SHong Zhang   /* determine symbolic Co=(p->B)^T*AP - send to others */
212ded4f561SHong Zhang   /*----------------------------------------------------*/
21377584fe6SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
214de0260b3SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
21542c57489SHong Zhang 
216ded4f561SHong Zhang   /* then, compute symbolic Co = (p->B)^T*AP */
217d0f46423SBarry Smith   pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors
21883445d95SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
219de0260b3SHong Zhang   ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr);
220de0260b3SHong Zhang   coi[0] = 0;
22142c57489SHong Zhang 
222d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */
223d16ca5f0SHong Zhang   nnz           = fill*(poti[pon] + api[am]);
224a1a86e44SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);
22542c57489SHong Zhang   current_space = free_space;
22642c57489SHong Zhang 
227de0260b3SHong Zhang   for (i=0; i<pon; i++) {
228ded4f561SHong Zhang     nnz = 0;
22982412ba7SHong Zhang     pnz = poti[i+1] - poti[i];
23077584fe6SHong Zhang     ptJ = potj + poti[i];
23177584fe6SHong Zhang     for (j=0; j<pnz; j++){
23277584fe6SHong Zhang       row  = ptJ[j]; /* row of AP == col of Pot */
23382412ba7SHong Zhang       apnz = api[row+1] - api[row];
234ded4f561SHong Zhang       Jptr = apj + api[row];
23583445d95SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
236*dadf0e6bSHong Zhang       ierr = PetscLLAddSorted(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
237ded4f561SHong Zhang       nnz += nlnk;
23842c57489SHong Zhang     }
23942c57489SHong Zhang 
24083445d95SHong Zhang     /* If free space is not available, double the total space in the list */
241ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
2424238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
243d16ca5f0SHong Zhang       nspacedouble++;
24442c57489SHong Zhang     }
24542c57489SHong Zhang 
24642c57489SHong Zhang     /* Copy data into free space, and zero out denserows */
247ded4f561SHong Zhang     ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
248ded4f561SHong Zhang     current_space->array           += nnz;
249ded4f561SHong Zhang     current_space->local_used      += nnz;
250ded4f561SHong Zhang     current_space->local_remaining -= nnz;
251ded4f561SHong Zhang     coi[i+1] = coi[i] + nnz;
25242c57489SHong Zhang   }
253de0260b3SHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr);
254a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
255d16ca5f0SHong Zhang   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]);
256d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
257de0260b3SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
25842c57489SHong Zhang 
259ded4f561SHong Zhang   /* send j-array (coj) of Co to other processors */
260ded4f561SHong Zhang   /*----------------------------------------------*/
26142c57489SHong Zhang   /* determine row ownership */
26283445d95SHong Zhang   ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr);
26326283091SBarry Smith   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
2647a2fc3feSBarry Smith   merge->rowmap->n = pn;
2657a2fc3feSBarry Smith   merge->rowmap->bs = 1;
26626283091SBarry Smith   ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
2677a2fc3feSBarry Smith   owners = merge->rowmap->range;
26842c57489SHong Zhang 
26942c57489SHong Zhang   /* determine the number of messages to send, their lengths */
27042c57489SHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr);
27183445d95SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
27242c57489SHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr);
27342c57489SHong Zhang   len_s = merge->len_s;
27442c57489SHong Zhang   merge->nsend = 0;
27583445d95SHong Zhang 
276de0260b3SHong Zhang   ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr);
277de0260b3SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
278de0260b3SHong Zhang 
27983445d95SHong Zhang   proc = 0;
280de0260b3SHong Zhang   for (i=0; i<pon; i++){
281de0260b3SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
282de0260b3SHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
283de0260b3SHong Zhang     len_s[proc] += coi[i+1] - coi[i];
28483445d95SHong Zhang   }
285de0260b3SHong Zhang 
286de0260b3SHong Zhang   len   = 0;  /* max length of buf_si[] */
287de0260b3SHong Zhang   owners_co[0] = 0;
28842c57489SHong Zhang   for (proc=0; proc<size; proc++){
289de0260b3SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
29083445d95SHong Zhang     if (len_si[proc]){
29142c57489SHong Zhang       merge->nsend++;
29283445d95SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
29342c57489SHong Zhang       len += len_si[proc];
29442c57489SHong Zhang     }
29542c57489SHong Zhang   }
29642c57489SHong Zhang 
297ded4f561SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
29842c57489SHong Zhang   ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
29942c57489SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
30042c57489SHong Zhang 
301ded4f561SHong Zhang   /* post the Irecv and Isend of coj */
302fc42d0c8SSatish Balay   ierr = PetscCommGetNewTag(comm,&tag);CHKERRQ(ierr);
303ded4f561SHong Zhang   ierr = PetscPostIrecvInt(comm,tag,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
30442c57489SHong Zhang 
305ded4f561SHong Zhang   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr);
30642c57489SHong Zhang 
30742c57489SHong Zhang   for (proc=0, k=0; proc<size; proc++){
30842c57489SHong Zhang     if (!len_s[proc]) continue;
309de0260b3SHong Zhang     i = owners_co[proc];
310ded4f561SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tag,comm,swaits+k);CHKERRQ(ierr);
31142c57489SHong Zhang     k++;
31242c57489SHong Zhang   }
31342c57489SHong Zhang 
314ded4f561SHong Zhang   /* receives and sends of coj are complete */
315ded4f561SHong Zhang   ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr);
31677584fe6SHong Zhang   for (i=0; i<merge->nrecv; i++){
317ded4f561SHong Zhang     ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr);
318ded4f561SHong Zhang   }
319ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
3200c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
32182412ba7SHong Zhang 
322ded4f561SHong Zhang   /* send and recv coi */
323ded4f561SHong Zhang   /*-------------------*/
324ded4f561SHong Zhang   ierr = PetscPostIrecvInt(comm,tag,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
32542c57489SHong Zhang 
32642c57489SHong Zhang   ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr);
32742c57489SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
32842c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++){
32942c57489SHong Zhang     if (!len_s[proc]) continue;
33042c57489SHong Zhang     /* form outgoing message for i-structure:
33142c57489SHong Zhang          buf_si[0]:                 nrows to be sent
33242c57489SHong Zhang                [1:nrows]:           row index (global)
33342c57489SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
33442c57489SHong Zhang     */
33542c57489SHong Zhang     /*-------------------------------------------*/
33642c57489SHong Zhang     nrows = len_si[proc]/2 - 1;
33742c57489SHong Zhang     buf_si_i    = buf_si + nrows+1;
33842c57489SHong Zhang     buf_si[0]   = nrows;
33942c57489SHong Zhang     buf_si_i[0] = 0;
34042c57489SHong Zhang     nrows = 0;
341de0260b3SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++){
342ded4f561SHong Zhang       nzi = coi[i+1] - coi[i];
343ded4f561SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
344de0260b3SHong Zhang       buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */
34542c57489SHong Zhang       nrows++;
34642c57489SHong Zhang     }
347ded4f561SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tag,comm,swaits+k);CHKERRQ(ierr);
34842c57489SHong Zhang     k++;
34942c57489SHong Zhang     buf_si += len_si[proc];
35042c57489SHong Zhang   }
351ded4f561SHong Zhang   i = merge->nrecv;
352ded4f561SHong Zhang   while (i--) {
353ded4f561SHong Zhang     ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr);
354ded4f561SHong Zhang   }
355ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
3560c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
357f2b054eeSHong Zhang   /*
358ae15b995SBarry Smith   ierr = PetscInfo2(A,"nsend: %d, nrecv: %d\n",merge->nsend,merge->nrecv);CHKERRQ(ierr);
35942c57489SHong Zhang   for (i=0; i<merge->nrecv; i++){
360ae15b995SBarry 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);
36142c57489SHong Zhang   }
362f2b054eeSHong Zhang   */
36342c57489SHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
36442c57489SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
365ded4f561SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
366ded4f561SHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
36742c57489SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
36842c57489SHong Zhang 
36977584fe6SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
37077584fe6SHong Zhang   time_Co += tf-t0;
37177584fe6SHong Zhang 
372de0260b3SHong Zhang   /* compute the local portion of C (mpi mat) */
373de0260b3SHong Zhang   /*------------------------------------------*/
37477584fe6SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
375ded4f561SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
376ded4f561SHong Zhang 
37742c57489SHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
37842c57489SHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
37942c57489SHong Zhang   bi[0] = 0;
38042c57489SHong Zhang 
381d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
382d16ca5f0SHong Zhang   nnz           = fill*(pi_loc[pm] + api[am]);
383a1a86e44SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);
38442c57489SHong Zhang   current_space = free_space;
38542c57489SHong Zhang 
386687f1162SBarry Smith   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
38742c57489SHong Zhang   for (k=0; k<merge->nrecv; k++){
38842c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
38942c57489SHong Zhang     nrows       = *buf_ri_k[k];
39042c57489SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
39142c57489SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
39242c57489SHong Zhang   }
39342c57489SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
39442c57489SHong Zhang   for (i=0; i<pn; i++) {
395ded4f561SHong Zhang     /* add pdt[i,:]*AP into lnk */
396ded4f561SHong Zhang     nnz = 0;
397ded4f561SHong Zhang     pnz = pdti[i+1] - pdti[i];
39877584fe6SHong Zhang     ptJ = pdtj + pdti[i];
39977584fe6SHong Zhang     for (j=0; j<pnz; j++){
40077584fe6SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
401ded4f561SHong Zhang       apnz = api[row+1] - api[row];
402ded4f561SHong Zhang       Jptr = apj + api[row];
403ded4f561SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
404*dadf0e6bSHong Zhang       ierr = PetscLLAddSorted(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
405ded4f561SHong Zhang       nnz += nlnk;
406ded4f561SHong Zhang     }
40742c57489SHong Zhang     /* add received col data into lnk */
40842c57489SHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
40942c57489SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
410ded4f561SHong Zhang         nzi = *(nextci[k]+1) - *nextci[k];
411ded4f561SHong Zhang         Jptr  = buf_rj[k] + *nextci[k];
412*dadf0e6bSHong Zhang         ierr = PetscLLAddSorted(nzi,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
413ded4f561SHong Zhang         nnz += nlnk;
41442c57489SHong Zhang         nextrow[k]++; nextci[k]++;
41542c57489SHong Zhang       }
41642c57489SHong Zhang     }
41742c57489SHong Zhang 
41842c57489SHong Zhang     /* if free space is not available, make more free space */
419ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
4204238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
421d16ca5f0SHong Zhang       nspacedouble++;
42242c57489SHong Zhang     }
42342c57489SHong Zhang     /* copy data into free space, then initialize lnk */
424ded4f561SHong Zhang     ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
425ded4f561SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
426ded4f561SHong Zhang     current_space->array           += nnz;
427ded4f561SHong Zhang     current_space->local_used      += nnz;
428ded4f561SHong Zhang     current_space->local_remaining -= nnz;
429ded4f561SHong Zhang     bi[i+1] = bi[i] + nnz;
43042c57489SHong Zhang   }
431ded4f561SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
4320572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
43342c57489SHong Zhang 
43442c57489SHong Zhang   ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
435a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
436d16ca5f0SHong Zhang   afill_tmp = (PetscReal)bi[pn]/(pi_loc[pm] + api[am]);
437d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
43842c57489SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
43942c57489SHong Zhang 
44077584fe6SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
44177584fe6SHong Zhang   time_Cd += tf-t0;
44277584fe6SHong Zhang 
44377584fe6SHong Zhang   /* create symbolic parallel matrix Cmpi */
44442c57489SHong Zhang   /*---------------------------------------*/
44577584fe6SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
44677584fe6SHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
44777584fe6SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
44877584fe6SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
44942c57489SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
45042c57489SHong Zhang 
45142c57489SHong Zhang   merge->bi            = bi;
45242c57489SHong Zhang   merge->bj            = bj;
453de0260b3SHong Zhang   merge->coi           = coi;
454de0260b3SHong Zhang   merge->coj           = coj;
45542c57489SHong Zhang   merge->buf_ri        = buf_ri;
45642c57489SHong Zhang   merge->buf_rj        = buf_rj;
457de0260b3SHong Zhang   merge->owners_co     = owners_co;
45877584fe6SHong Zhang   merge->destroy       = Cmpi->ops->destroy;
45977584fe6SHong Zhang   merge->duplicate     = Cmpi->ops->duplicate;
460cc6cb787SHong Zhang 
46177584fe6SHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
46277584fe6SHong Zhang   Cmpi->assembled      = PETSC_FALSE;
46377584fe6SHong Zhang   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_PtAP;
46477584fe6SHong Zhang   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP;
46577584fe6SHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
46642c57489SHong Zhang 
46777584fe6SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
46877584fe6SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
469f8487c73SHong Zhang   c->ptap        = ptap;
47077584fe6SHong Zhang   ptap->api      = api;
47177584fe6SHong Zhang   ptap->apj      = apj;
472f8487c73SHong Zhang   ptap->merge    = merge;
47377584fe6SHong Zhang   ptap->abnz_max = ap_rmax;
474f8487c73SHong Zhang 
47577584fe6SHong Zhang   *C = Cmpi;
476f2b054eeSHong Zhang #if defined(PETSC_USE_INFO)
477f2b054eeSHong Zhang   if (bi[pn] != 0) {
47877584fe6SHong Zhang     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
47977584fe6SHong Zhang     ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
480f2b054eeSHong Zhang   } else {
48177584fe6SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
482f2b054eeSHong Zhang   }
483f2b054eeSHong Zhang #endif
48477584fe6SHong Zhang   ierr = PetscGetTime(&tff);CHKERRQ(ierr);
48577584fe6SHong Zhang   etime += tff - t00;
48677584fe6SHong Zhang   /*
48777584fe6SHong Zhang   PetscInt prid=0;
48877584fe6SHong Zhang   if (rank == prid){
48977584fe6SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] PtAPSym time %g = matsetup %g + AP %g + Co %g + Cd %g\n",rank,etime,time_matupdate,time_AP,time_Co,time_Cd);
49077584fe6SHong Zhang   }
49177584fe6SHong Zhang    */
49242c57489SHong Zhang   PetscFunctionReturn(0);
49342c57489SHong Zhang }
49442c57489SHong Zhang 
49542c57489SHong Zhang #undef __FUNCT__
49642c57489SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ"
49742c57489SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
49842c57489SHong Zhang {
49942c57489SHong Zhang   PetscErrorCode       ierr;
50042c57489SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
501f8487c73SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
50242c57489SHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
503de0260b3SHong Zhang   Mat_SeqAIJ           *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
50442c57489SHong Zhang   Mat_SeqAIJ           *p_loc,*p_oth;
505f8487c73SHong Zhang   Mat_PtAPMPI          *ptap;
5061d633602SHong Zhang   PetscInt             *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp;
50782412ba7SHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj;
50882412ba7SHong Zhang   PetscInt             i,j,k,anz,pnz,apnz,nextap,row,*cj;
509e900a757SHong Zhang   MatScalar            *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp;
510d0f46423SBarry Smith   PetscInt             am=A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
5117adad957SLisandro Dalcin   MPI_Comm             comm=((PetscObject)C)->comm;
51242c57489SHong Zhang   PetscMPIInt          size,rank,taga,*len_s;
51382412ba7SHong Zhang   PetscInt             *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
51442c57489SHong Zhang   PetscInt             **buf_ri,**buf_rj;
51550f5bf86SHong Zhang   PetscInt             cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */
51642c57489SHong Zhang   MPI_Request          *s_waits,*r_waits;
51742c57489SHong Zhang   MPI_Status           *status;
51882412ba7SHong Zhang   MatScalar            **abuf_r,*ba_i,*pA,*coa,*ba;
51982412ba7SHong Zhang   PetscInt             *api,*apj,*coi,*coj;
520d0f46423SBarry Smith   PetscInt             *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend;
521a9d06231SHong Zhang   PetscInt             sparse_axpy;
522f9473cf7SHong Zhang   PetscLogDouble       t0,tf,etime=0.0,t00,tff,time_matupdate=0.0,time_malloc=0.0,time_Cseq0=0.0,time_Cseq1=0.0,time_setvals=0.0;
52342c57489SHong Zhang 
52442c57489SHong Zhang   PetscFunctionBegin;
525f8487c73SHong Zhang   /* ierr = MPI_Barrier(comm);CHKERRQ(ierr); */
526f9473cf7SHong Zhang   ierr = PetscGetTime(&t00);CHKERRQ(ierr);
52742c57489SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
52842c57489SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
52942c57489SHong Zhang 
530f9473cf7SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
531f8487c73SHong Zhang   ptap  = c->ptap;
532f8487c73SHong Zhang   merge = ptap->merge;
5336c6e00e0SHong Zhang 
534a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
535f9473cf7SHong Zhang   /*--------------------------------------------------*/
536f8487c73SHong Zhang   if (ptap->reuse == MAT_INITIAL_MATRIX){
537f8487c73SHong Zhang     ptap->reuse = MAT_REUSE_MATRIX;
5386c6e00e0SHong Zhang   } else { /* update numerical values of P_oth and P_loc */
539a1a4e84aSHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
540a1a4e84aSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
5416c6e00e0SHong Zhang   }
542f8487c73SHong Zhang 
543f9473cf7SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
544f9473cf7SHong Zhang   time_matupdate += tf-t0;
54542c57489SHong Zhang 
546f9473cf7SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
547f9473cf7SHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
548f9473cf7SHong Zhang   /*--------------------------------------------------------------*/
54942c57489SHong Zhang   /* get data from symbolic products */
550a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
551a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
552a9d06231SHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a;
55342c57489SHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
55442c57489SHong Zhang 
555de0260b3SHong Zhang   coi = merge->coi; coj = merge->coj;
556de0260b3SHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr);
557de0260b3SHong Zhang   ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr);
558de0260b3SHong Zhang 
559de0260b3SHong Zhang   bi     = merge->bi; bj = merge->bj;
5607a2fc3feSBarry Smith   owners = merge->rowmap->range;
561de0260b3SHong Zhang   ierr   = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
562de0260b3SHong Zhang   ierr   = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr);
563f9473cf7SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
564f9473cf7SHong Zhang   time_malloc += tf-t0;
56542c57489SHong Zhang 
566a1a4e84aSHong Zhang   api = ptap->api; apj = ptap->apj;
5670496f32cSHong Zhang   /* flag 'sparse_axpy' determines which implementations to be used:
5680496f32cSHong Zhang        0: do dense axpy in MatPtAPNumeric() - fastest, but requires storage of a dense array apa; (default)
5690496f32cSHong Zhang        1: do one sparse axpy - uses same memory as sparse_axpy=0 and might execute less flops
5700496f32cSHong Zhang           (apnz vs. cnz in the outerproduct), slower than case '0' when cnz is not too large than apnz;
5710496f32cSHong Zhang        2: do two sparse axpy in MatPtAPNumeric() - slowest, uses a sparse array apa */
572a9d06231SHong Zhang   /* set default sparse_axpy */
573e900a757SHong Zhang   sparse_axpy = 0;
5740496f32cSHong Zhang   ierr = PetscOptionsGetInt(PETSC_NULL,"-matptap_sparseaxpy",&sparse_axpy,PETSC_NULL);CHKERRQ(ierr);
575a9d06231SHong Zhang   if (sparse_axpy == 0){  /* Do not perform sparse axpy */
576a9d06231SHong Zhang     /*--------------------------------------------------*/
577a9d06231SHong Zhang     /* malloc apa to store dense row A[i,:]*P */
578a9d06231SHong Zhang     ierr = PetscMalloc((P->cmap->N)*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
579a9d06231SHong Zhang     ierr = PetscMemzero(apa,P->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
580a9d06231SHong Zhang 
581a9d06231SHong Zhang     for (i=0; i<am; i++) {
582a9d06231SHong Zhang       ierr = PetscGetTime(&t0);CHKERRQ(ierr);
583a9d06231SHong Zhang       /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */
584a9d06231SHong Zhang       /*------------------------------------------------------------*/
585a9d06231SHong Zhang       apJ = apj + api[i];
586a9d06231SHong Zhang 
587a9d06231SHong Zhang       /* diagonal portion of A */
588a9d06231SHong Zhang       anz = adi[i+1] - adi[i];
589a9d06231SHong Zhang       adj = ad->j + adi[i];
590a9d06231SHong Zhang       ada = ad->a + adi[i];
591a9d06231SHong Zhang       for (j=0; j<anz; j++) {
592a9d06231SHong Zhang         row = adj[j];
593a9d06231SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
594a9d06231SHong Zhang         pj  = pj_loc + pi_loc[row];
595a9d06231SHong Zhang         pa  = pa_loc + pi_loc[row];
596a9d06231SHong Zhang 
597a9d06231SHong Zhang         /* perform dense axpy */
598e900a757SHong Zhang         valtmp = ada[j];
599a9d06231SHong Zhang         for (k=0; k<pnz; k++){
600e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
601a9d06231SHong Zhang         }
602a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
603a9d06231SHong Zhang       }
604a9d06231SHong Zhang 
605a9d06231SHong Zhang       /* off-diagonal portion of A */
606a9d06231SHong Zhang       anz = aoi[i+1] - aoi[i];
607a9d06231SHong Zhang       aoj = ao->j + aoi[i];
608a9d06231SHong Zhang       aoa = ao->a + aoi[i];
609a9d06231SHong Zhang       for (j=0; j<anz; j++) {
610a9d06231SHong Zhang         row = aoj[j];
611a9d06231SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
612a9d06231SHong Zhang         pj  = pj_oth + pi_oth[row];
613a9d06231SHong Zhang         pa  = pa_oth + pi_oth[row];
614a9d06231SHong Zhang 
615a9d06231SHong Zhang         /* perform dense axpy */
616e900a757SHong Zhang         valtmp = aoa[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       ierr = PetscGetTime(&tf);CHKERRQ(ierr);
623a9d06231SHong Zhang       time_Cseq0 += tf - t0;
624a9d06231SHong Zhang 
625a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
626a9d06231SHong Zhang       /*--------------------------------------------------------------*/
627a9d06231SHong Zhang       ierr = PetscGetTime(&t0);CHKERRQ(ierr);
628a9d06231SHong Zhang       apnz = api[i+1] - api[i];
629a9d06231SHong Zhang       /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */
630a9d06231SHong Zhang       pnz = po->i[i+1] - po->i[i];
631e900a757SHong Zhang       poJ = po->j + po->i[i];
632a9d06231SHong Zhang       pA  = po->a + po->i[i];
633a9d06231SHong Zhang       for (j=0; j<pnz; j++){
634e900a757SHong Zhang         row = poJ[j];
635e900a757SHong Zhang         cnz = coi[row+1] - coi[row];
636e900a757SHong Zhang         cj  = coj + coi[row];
637e900a757SHong Zhang         ca  = coa + coi[row];
638a9d06231SHong Zhang         /* perform dense axpy */
639e900a757SHong Zhang         valtmp = pA[j];
640a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
641e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
642a9d06231SHong Zhang         }
643a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
644a9d06231SHong Zhang       }
645a9d06231SHong Zhang 
646a9d06231SHong Zhang       /* put the value into Cd (diagonal part) */
647a9d06231SHong Zhang       pnz = pd->i[i+1] - pd->i[i];
648e900a757SHong Zhang       pdJ = pd->j + pd->i[i];
649a9d06231SHong Zhang       pA  = pd->a + pd->i[i];
650a9d06231SHong Zhang       for (j=0; j<pnz; j++){
651e900a757SHong Zhang         row = pdJ[j];
652e900a757SHong Zhang         cnz = bi[row+1] - bi[row];
653e900a757SHong Zhang         cj  = bj + bi[row];
654e900a757SHong Zhang         ca  = ba + bi[row];
655a9d06231SHong Zhang         /* perform dense axpy */
656e900a757SHong Zhang         valtmp = pA[j];
657a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
658e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
659a9d06231SHong Zhang         }
660a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
661a9d06231SHong Zhang       }
662a9d06231SHong Zhang 
663a9d06231SHong Zhang       /* zero the current row of A*P */
664a9d06231SHong Zhang       for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0;
665a9d06231SHong Zhang       ierr = PetscGetTime(&tf);CHKERRQ(ierr);
666a9d06231SHong Zhang       time_Cseq1 += tf - t0;
667a9d06231SHong Zhang     }
668a9d06231SHong Zhang   } else if (sparse_axpy == 1){  /* Perform one sparse axpy */
669a9d06231SHong Zhang     /*------------------------------------------------------*/
6701d633602SHong Zhang     /* malloc apa to store dense row A[i,:]*P */
6711d633602SHong Zhang     ierr = PetscMalloc((P->cmap->N)*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
6721d633602SHong Zhang     ierr = PetscMemzero(apa,P->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
6731d633602SHong Zhang 
6741d633602SHong Zhang     for (i=0; i<am; i++) {
675f9473cf7SHong Zhang       ierr = PetscGetTime(&t0);CHKERRQ(ierr);
676f9473cf7SHong Zhang       /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */
677f9473cf7SHong Zhang       /*------------------------------------------------------------*/
6781d633602SHong Zhang       apJ = apj + api[i];
6791d633602SHong Zhang 
6801d633602SHong Zhang       /* diagonal portion of A */
6811d633602SHong Zhang       anz = adi[i+1] - adi[i];
6821d633602SHong Zhang       adj = ad->j + adi[i];
6831d633602SHong Zhang       ada = ad->a + adi[i];
6841d633602SHong Zhang       for (j=0; j<anz; j++) {
6851d633602SHong Zhang         row = adj[j];
6861d633602SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
6871d633602SHong Zhang         pj  = pj_loc + pi_loc[row];
6881d633602SHong Zhang         pa  = pa_loc + pi_loc[row];
6891d633602SHong Zhang 
6901d633602SHong Zhang         /* perform dense axpy */
691e900a757SHong Zhang         valtmp = ada[j];
6921d633602SHong Zhang         for (k=0; k<pnz; k++){
693e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
6941d633602SHong Zhang         }
6951d633602SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
6961d633602SHong Zhang       }
6971d633602SHong Zhang 
6981d633602SHong Zhang       /* off-diagonal portion of A */
6991d633602SHong Zhang       anz = aoi[i+1] - aoi[i];
7001d633602SHong Zhang       aoj = ao->j + aoi[i];
7011d633602SHong Zhang       aoa = ao->a + aoi[i];
7021d633602SHong Zhang       for (j=0; j<anz; j++) {
7031d633602SHong Zhang         row = aoj[j];
7041d633602SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
7051d633602SHong Zhang         pj  = pj_oth + pi_oth[row];
7061d633602SHong Zhang         pa  = pa_oth + pi_oth[row];
7071d633602SHong Zhang 
7081d633602SHong Zhang         /* perform dense axpy */
709e900a757SHong Zhang         valtmp = aoa[j];
7101d633602SHong Zhang         for (k=0; k<pnz; k++){
711e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
7121d633602SHong Zhang         }
7131d633602SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
7141d633602SHong Zhang       }
715f9473cf7SHong Zhang       ierr = PetscGetTime(&tf);CHKERRQ(ierr);
716f9473cf7SHong Zhang       time_Cseq0 += tf - t0;
7171d633602SHong Zhang 
718f9473cf7SHong Zhang       ierr = PetscGetTime(&t0);CHKERRQ(ierr);
719f9473cf7SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
720f9473cf7SHong Zhang       /*--------------------------------------------------------------*/
7211d633602SHong Zhang       apnz = api[i+1] - api[i];
722f9473cf7SHong Zhang       /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */
723f9473cf7SHong Zhang       pnz = po->i[i+1] - po->i[i];
724e900a757SHong Zhang       poJ = po->j + po->i[i];
725a9d06231SHong Zhang       pA  = po->a + po->i[i];
7261d633602SHong Zhang       for (j=0; j<pnz; j++){
727e900a757SHong Zhang         row    = poJ[j];
728e900a757SHong Zhang         cj     = coj + coi[row];
729e900a757SHong Zhang         ca     = coa + coi[row];
730e900a757SHong Zhang         valtmp = pA[j];
7311d633602SHong Zhang         /* perform sparse axpy */
7321d633602SHong Zhang         nextap = 0;
7331d633602SHong Zhang         for (k=0; nextap<apnz; k++) {
7341d633602SHong Zhang           if (cj[k]==apJ[nextap]) { /* global column index */
735e900a757SHong Zhang             ca[k] += valtmp*apa[cj[k]]; nextap++;
7361d633602SHong Zhang           }
7371d633602SHong Zhang         }
7381d633602SHong Zhang         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
7391d633602SHong Zhang       }
740f9473cf7SHong Zhang 
741f9473cf7SHong Zhang       /* put the value into Cd (diagonal part) */
742f9473cf7SHong Zhang       pnz = pd->i[i+1] - pd->i[i];
743e900a757SHong Zhang       pdJ = pd->j + pd->i[i];
744a9d06231SHong Zhang       pA  = pd->a + pd->i[i];
745f9473cf7SHong Zhang       for (j=0; j<pnz; j++){
746e900a757SHong Zhang         row    = pdJ[j];
747e900a757SHong Zhang         cj     = bj + bi[row];
748e900a757SHong Zhang         ca     = ba + bi[row];
749e900a757SHong Zhang         valtmp = pA[j];
750f9473cf7SHong Zhang         /* perform sparse axpy */
751f9473cf7SHong Zhang         nextap = 0;
752f9473cf7SHong Zhang         for (k=0; nextap<apnz; k++) {
753f9473cf7SHong Zhang           if (cj[k]==apJ[nextap]) { /* global column index */
754e900a757SHong Zhang             ca[k] += valtmp*apa[cj[k]];
755a9d06231SHong Zhang             nextap++;
756f9473cf7SHong Zhang           }
757f9473cf7SHong Zhang         }
758f9473cf7SHong Zhang         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
759f9473cf7SHong Zhang       }
760f9473cf7SHong Zhang 
761f9473cf7SHong Zhang       /* zero the current row of A*P */
7621d633602SHong Zhang       for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0;
763f9473cf7SHong Zhang       ierr = PetscGetTime(&tf);CHKERRQ(ierr);
764f9473cf7SHong Zhang       time_Cseq1 += tf - t0;
7651d633602SHong Zhang     }
7660496f32cSHong Zhang   } else if (sparse_axpy == 2){/* Perform two sparse axpy */
767a9d06231SHong Zhang     /*----------------------------------------------------*/
7681d633602SHong Zhang     /* malloc apa to store sparse row A[i,:]*P */
769f8487c73SHong Zhang     ierr = PetscMalloc((ptap->abnz_max+1)*sizeof(MatScalar),&apa);CHKERRQ(ierr);
770f8487c73SHong Zhang     ierr = PetscMemzero(apa,ptap->abnz_max*sizeof(MatScalar));CHKERRQ(ierr);
77142c57489SHong Zhang 
772a9d06231SHong Zhang     pA=pa_loc;
77342c57489SHong Zhang     for (i=0; i<am; i++) {
774a9d06231SHong Zhang       ierr = PetscGetTime(&t0);CHKERRQ(ierr);
775f4a743e1SHong Zhang       /* form i-th sparse row of A*P */
77682412ba7SHong Zhang       apnz = api[i+1] - api[i];
77782412ba7SHong Zhang       apJ  = apj + api[i];
77842c57489SHong Zhang       /* diagonal portion of A */
77982412ba7SHong Zhang       anz = adi[i+1] - adi[i];
7801d633602SHong Zhang       adj = ad->j + adi[i];
7811d633602SHong Zhang       ada = ad->a + adi[i];
78282412ba7SHong Zhang       for (j=0; j<anz; j++) {
7831d633602SHong Zhang         row = adj[j];
78482412ba7SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
78582412ba7SHong Zhang         pj  = pj_loc + pi_loc[row];
78682412ba7SHong Zhang         pa  = pa_loc + pi_loc[row];
787e900a757SHong Zhang         valtmp = ada[j];
788f4a743e1SHong Zhang         nextp  = 0;
78982412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
79082412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
791e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
79242c57489SHong Zhang           }
79342c57489SHong Zhang         }
794dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
79542c57489SHong Zhang       }
79642c57489SHong Zhang       /* off-diagonal portion of A */
79782412ba7SHong Zhang       anz = aoi[i+1] - aoi[i];
7981d633602SHong Zhang       aoj = ao->j + aoi[i];
7991d633602SHong Zhang       aoa = ao->a + aoi[i];
80082412ba7SHong Zhang       for (j=0; j<anz; j++) {
8011d633602SHong Zhang         row = aoj[j];
80282412ba7SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
80382412ba7SHong Zhang         pj  = pj_oth + pi_oth[row];
80482412ba7SHong Zhang         pa  = pa_oth + pi_oth[row];
805e900a757SHong Zhang         valtmp = aoa[j];
806f4a743e1SHong Zhang         nextp  = 0;
80782412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
80882412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
809e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
81042c57489SHong Zhang           }
81142c57489SHong Zhang         }
812dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
81342c57489SHong Zhang       }
814a9d06231SHong Zhang       ierr = PetscGetTime(&tf);CHKERRQ(ierr);
815a9d06231SHong Zhang       time_Cseq0 += tf - t0;
81642c57489SHong Zhang 
817a9d06231SHong Zhang       ierr = PetscGetTime(&t0);CHKERRQ(ierr);
818a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
819a9d06231SHong Zhang       /*--------------------------------------------------------------*/
82082412ba7SHong Zhang       pnz = pi_loc[i+1] - pi_loc[i];
821f9473cf7SHong Zhang       pJ  = pj_loc + pi_loc[i];
82282412ba7SHong Zhang       for (j=0; j<pnz; j++) {
82342c57489SHong Zhang         nextap = 0;
824f9473cf7SHong Zhang         row    = pJ[j]; /* global index */
82582412ba7SHong Zhang         if (row < pcstart || row >=pcend) { /* put the value into Co */
826e900a757SHong Zhang           row = *poJ;
827e900a757SHong Zhang           cj  = coj + coi[row];
828e900a757SHong Zhang           ca  = coa + coi[row]; poJ++;
82982412ba7SHong Zhang         } else {                            /* put the value into Cd */
830e900a757SHong Zhang           row  = *pdJ;
831e900a757SHong Zhang           cj   = bj + bi[row];
832e900a757SHong Zhang           ca   = ba + bi[row]; pdJ++;
83342c57489SHong Zhang         }
834e900a757SHong Zhang         valtmp = pA[j];
83582412ba7SHong Zhang         for (k=0; nextap<apnz; k++) {
836e900a757SHong Zhang           if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++];
83742c57489SHong Zhang         }
838dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
839de0260b3SHong Zhang       }
8400496f32cSHong Zhang       pA += pnz;
84142c57489SHong Zhang       /* zero the current row info for A*P */
84282412ba7SHong Zhang       ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr);
843a9d06231SHong Zhang       ierr = PetscGetTime(&tf);CHKERRQ(ierr);
844a9d06231SHong Zhang       time_Cseq1 += tf - t0;
84542c57489SHong Zhang     }
846a9d06231SHong Zhang   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"sparse_axpy only takes values 0, 1 and 2");
84742c57489SHong Zhang   ierr = PetscFree(apa);CHKERRQ(ierr);
84842c57489SHong Zhang 
849a9d06231SHong Zhang   /* 3) send and recv matrix values coa */
850a9d06231SHong Zhang   /*------------------------------------*/
85142c57489SHong Zhang   buf_ri = merge->buf_ri;
85242c57489SHong Zhang   buf_rj = merge->buf_rj;
85342c57489SHong Zhang   len_s  = merge->len_s;
854fc42d0c8SSatish Balay   ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
85542c57489SHong Zhang   ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
85642c57489SHong Zhang 
857a9d06231SHong Zhang   ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr);
85842c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++){
85942c57489SHong Zhang     if (!len_s[proc]) continue;
860de0260b3SHong Zhang     i = merge->owners_co[proc];
861de0260b3SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
86242c57489SHong Zhang     k++;
86342c57489SHong Zhang   }
8640c468ba9SBarry Smith   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
8650c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
86642c57489SHong Zhang 
867a9d06231SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
86842c57489SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
86982412ba7SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
87042c57489SHong Zhang 
871a9d06231SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
872a9d06231SHong Zhang   /*------------------------------------------------------*/
873f9473cf7SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
874687f1162SBarry Smith   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
87542c57489SHong Zhang 
87642c57489SHong Zhang   for (k=0; k<merge->nrecv; k++){
87742c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
87842c57489SHong Zhang     nrows       = *(buf_ri_k[k]);
87942c57489SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
88082412ba7SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
88142c57489SHong Zhang   }
88242c57489SHong Zhang 
883de0260b3SHong Zhang   for (i=0; i<cm; i++) {
88482412ba7SHong Zhang     row = owners[rank] + i; /* global row index of C_seq */
88542c57489SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
886de0260b3SHong Zhang     ba_i = ba + bi[i];
88782412ba7SHong Zhang     bnz  = bi[i+1] - bi[i];
88842c57489SHong Zhang     /* add received vals into ba */
88942c57489SHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
89042c57489SHong Zhang       /* i-th row */
89142c57489SHong Zhang       if (i == *nextrow[k]) {
89282412ba7SHong Zhang         cnz = *(nextci[k]+1) - *nextci[k];
89382412ba7SHong Zhang         cj  = buf_rj[k] + *(nextci[k]);
89482412ba7SHong Zhang         ca  = abuf_r[k] + *(nextci[k]);
89582412ba7SHong Zhang         nextcj = 0;
89682412ba7SHong Zhang         for (j=0; nextcj<cnz; j++){
89782412ba7SHong Zhang           if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */
89882412ba7SHong Zhang             ba_i[j] += ca[nextcj++];
89942c57489SHong Zhang           }
90042c57489SHong Zhang         }
90182412ba7SHong Zhang         nextrow[k]++; nextci[k]++;
90242c57489SHong Zhang       }
90342c57489SHong Zhang     }
90482412ba7SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
905dc0b31edSSatish Balay     ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
90642c57489SHong Zhang   }
90742c57489SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
90842c57489SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
909f9473cf7SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
910f9473cf7SHong Zhang   time_setvals += tf-t0;
91142c57489SHong Zhang 
912de0260b3SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
913c05d87d6SBarry Smith   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
91442c57489SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
9150572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
916f9473cf7SHong Zhang 
917f9473cf7SHong Zhang   ierr = PetscGetTime(&tff);CHKERRQ(ierr);
918f9473cf7SHong Zhang   etime += tff - t00;
919f9473cf7SHong Zhang   /*
920a9d06231SHong Zhang   PetscInt prid=0;
921a9d06231SHong Zhang   if (rank == prid){
922f9473cf7SHong Zhang    ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] PtAPNum time %g = matupdate %g + malloc %g + Cseq %g + %g + setvals %g\n",rank,etime,time_matupdate,time_malloc,time_Cseq0,time_Cseq1,time_setvals);
923a9d06231SHong Zhang   }
924f9473cf7SHong Zhang    */
92542c57489SHong Zhang   PetscFunctionReturn(0);
92642c57489SHong Zhang }
927