xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision 13f7495037498d59e0d9e83d0e810f28e283c435)
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 
742c57489SHong Zhang #include "src/mat/impls/aij/seq/aij.h"   /*I "petscmat.h" I*/
842c57489SHong Zhang #include "src/mat/utils/freespace.h"
942c57489SHong Zhang #include "src/mat/impls/aij/mpi/mpiaij.h"
1042c57489SHong Zhang #include "petscbt.h"
1142c57489SHong Zhang 
1242c57489SHong Zhang EXTERN PetscErrorCode MatDestroy_MPIAIJ(Mat);
1342c57489SHong Zhang #undef __FUNCT__
1442c57489SHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_PtAP"
1542c57489SHong Zhang PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A)
1642c57489SHong Zhang {
1742c57489SHong Zhang   PetscErrorCode       ierr;
18a6b2eed2SHong Zhang   Mat_MatMatMultMPI    *mult;
1942c57489SHong Zhang   PetscObjectContainer container;
2042c57489SHong Zhang 
2142c57489SHong Zhang   PetscFunctionBegin;
2242c57489SHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"MatPtAPMPI",(PetscObject *)&container);CHKERRQ(ierr);
2342c57489SHong Zhang   if (container) {
24a6b2eed2SHong Zhang     ierr  = PetscObjectContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
25a6b2eed2SHong Zhang     ierr = MatDestroy(mult->B_loc);CHKERRQ(ierr);
26a6b2eed2SHong Zhang     ierr = MatDestroy(mult->B_oth);CHKERRQ(ierr);
27a6b2eed2SHong Zhang     if (mult->abi) ierr = PetscFree(mult->abi);CHKERRQ(ierr);
28a6b2eed2SHong Zhang     if (mult->abj) ierr = PetscFree(mult->abj);CHKERRQ(ierr);
29a6b2eed2SHong Zhang     if (mult->startsj) ierr = PetscFree(mult->startsj);CHKERRQ(ierr);
30a6b2eed2SHong Zhang     if (mult->bufa) ierr = PetscFree(mult->bufa);CHKERRQ(ierr);
3142c57489SHong Zhang     ierr = PetscObjectContainerDestroy(container);CHKERRQ(ierr);
3242c57489SHong Zhang     ierr = PetscObjectCompose((PetscObject)A,"MatPtAPMPI",0);CHKERRQ(ierr);
3342c57489SHong Zhang   }
34a6b2eed2SHong Zhang   ierr = PetscFree(mult);CHKERRQ(ierr);
3542c57489SHong Zhang 
3642c57489SHong Zhang   ierr = MatDestroy_MPIAIJ(A);CHKERRQ(ierr);
3742c57489SHong Zhang   PetscFunctionReturn(0);
3842c57489SHong Zhang }
3942c57489SHong Zhang 
4042c57489SHong Zhang #undef __FUNCT__
4142c57489SHong Zhang #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ"
4242c57489SHong Zhang PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
4342c57489SHong Zhang {
4442c57489SHong Zhang   PetscErrorCode       ierr;
45a6b2eed2SHong Zhang   Mat_MatMatMultMPI    *mult;
4642c57489SHong Zhang   PetscObjectContainer container;
4742c57489SHong Zhang 
4842c57489SHong Zhang   PetscFunctionBegin;
4942c57489SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
5042c57489SHong Zhang     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
51a6b2eed2SHong Zhang     ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr);
52a6b2eed2SHong Zhang     mult->B_loc=PETSC_NULL; mult->B_oth=PETSC_NULL;
53a6b2eed2SHong Zhang     mult->abi=PETSC_NULL; mult->abj=PETSC_NULL;
54a6b2eed2SHong Zhang     mult->abnz_max = 0; /* symbolic A*P is not done yet */
5542c57489SHong Zhang 
5642c57489SHong Zhang     /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
57a6b2eed2SHong Zhang     ierr = MatGetBrowsOfAoCols(A,P,scall,&mult->startsj,&mult->bufa,&mult->B_oth);CHKERRQ(ierr);
5842c57489SHong Zhang 
5942c57489SHong Zhang     /* get P_loc by taking all local rows of P */
60a6b2eed2SHong Zhang     ierr = MatGetLocalMat(P,scall,&mult->B_loc);CHKERRQ(ierr);
6142c57489SHong Zhang 
6242c57489SHong Zhang     /* attach the supporting struct to P for reuse */
6342c57489SHong Zhang     P->ops->destroy  = MatDestroy_MPIAIJ_PtAP;
6442c57489SHong Zhang     ierr = PetscObjectContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
65a6b2eed2SHong Zhang     ierr = PetscObjectContainerSetPointer(container,mult);CHKERRQ(ierr);
6642c57489SHong Zhang     ierr = PetscObjectCompose((PetscObject)P,"MatPtAPMPI",(PetscObject)container);CHKERRQ(ierr);
6742c57489SHong Zhang 
6842c57489SHong Zhang     /* now, compute symbolic local P^T*A*P */
6942c57489SHong Zhang     ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);/* numeric product is computed as well */
7042c57489SHong Zhang     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
7142c57489SHong Zhang   } else if (scall == MAT_REUSE_MATRIX){
7242c57489SHong Zhang     ierr = PetscObjectQuery((PetscObject)P,"MatPtAPMPI",(PetscObject *)&container);CHKERRQ(ierr);
7342c57489SHong Zhang     if (container) {
74a6b2eed2SHong Zhang       ierr  = PetscObjectContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr);
7542c57489SHong Zhang     } else {
7642c57489SHong Zhang       SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix P does not posses an object container");
7742c57489SHong Zhang     }
7842c57489SHong Zhang 
7942c57489SHong Zhang     /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
80a6b2eed2SHong Zhang     ierr = MatGetBrowsOfAoCols(A,P,scall,&mult->startsj,&mult->bufa,&mult->B_oth);CHKERRQ(ierr);
8142c57489SHong Zhang 
8242c57489SHong Zhang     /* get P_loc by taking all local rows of P */
83a6b2eed2SHong Zhang     ierr = MatGetLocalMat(P,scall,&mult->B_loc);CHKERRQ(ierr);
8442c57489SHong Zhang 
8542c57489SHong Zhang   } else {
8642c57489SHong Zhang     SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",scall);
8742c57489SHong Zhang   }
8842c57489SHong Zhang   /* now, compute numeric local P^T*A*P */
8942c57489SHong Zhang   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
9042c57489SHong Zhang   ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ(A,P,*C);CHKERRQ(ierr);
9142c57489SHong Zhang   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
9242c57489SHong Zhang 
9342c57489SHong Zhang   PetscFunctionReturn(0);
9442c57489SHong Zhang }
9542c57489SHong Zhang 
9642c57489SHong Zhang #undef __FUNCT__
9742c57489SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ"
9842c57489SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
9942c57489SHong Zhang {
10042c57489SHong Zhang   PetscErrorCode       ierr;
101de0260b3SHong Zhang   Mat                  P_loc,P_oth,B_mpi;
102de0260b3SHong Zhang   Mat_MatMatMultMPI    *ap;
10342c57489SHong Zhang   PetscObjectContainer container;
10442c57489SHong Zhang   FreeSpaceList        free_space=PETSC_NULL,current_space=PETSC_NULL;
10583445d95SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data;
10642c57489SHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
10742c57489SHong Zhang   Mat_SeqAIJ           *p_loc,*p_oth;
108ded4f561SHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ;
10983445d95SHong Zhang   PetscInt             *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,nnz;
110*13f74950SBarry Smith   PetscInt             nlnk,*lnk,*owners_co,*coi,*coj,i,k,pnz,row;
111a6b2eed2SHong Zhang   PetscInt             am=A->m,pN=P->N,pn=P->n;
11242c57489SHong Zhang   PetscBT              lnkbt;
11342c57489SHong Zhang   MPI_Comm             comm=A->comm;
114ded4f561SHong Zhang   PetscMPIInt          size,rank,tag,*len_si,*len_s,*len_ri;
11542c57489SHong Zhang   PetscInt             **buf_rj,**buf_ri,**buf_ri_k;
11642c57489SHong Zhang   PetscInt             len,proc,*dnz,*onz,*owners;
117ded4f561SHong Zhang   PetscInt             nzi,*bi,*bj;
11842c57489SHong Zhang   PetscInt             nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
119ded4f561SHong Zhang   MPI_Request          *swaits,*rwaits;
120ded4f561SHong Zhang   MPI_Status           *sstatus,rstatus;
12142c57489SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
122ded4f561SHong Zhang   PetscInt             *api,*apj,*Jptr,apnz,*prmap=p->garray,pon;
123*13f74950SBarry Smith   PetscMPIInt          j;
12442c57489SHong Zhang 
12542c57489SHong Zhang   PetscFunctionBegin;
12642c57489SHong Zhang   ierr = PetscObjectQuery((PetscObject)P,"MatPtAPMPI",(PetscObject *)&container);CHKERRQ(ierr);
12742c57489SHong Zhang   if (container) {
128de0260b3SHong Zhang     ierr  = PetscObjectContainerGetPointer(container,(void **)&ap);CHKERRQ(ierr);
12942c57489SHong Zhang   } else {
13042c57489SHong Zhang     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix P does not posses an object container");
13142c57489SHong Zhang   }
13283445d95SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
13383445d95SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
13483445d95SHong Zhang 
13542c57489SHong Zhang   /* get data from P_loc and P_oth */
136a6b2eed2SHong Zhang   P_loc=ap->B_loc; P_oth=ap->B_oth;
13742c57489SHong Zhang   p_loc=(Mat_SeqAIJ*)P_loc->data; p_oth=(Mat_SeqAIJ*)P_oth->data;
13842c57489SHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pi_oth=p_oth->i; pj_oth=p_oth->j;
13942c57489SHong Zhang 
14083445d95SHong Zhang   /* first, compute symbolic AP = A_loc*P */
141de0260b3SHong Zhang   /*--------------------------------------*/
14282412ba7SHong Zhang   ierr  = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr);
14382412ba7SHong Zhang   ap->abi = api;
14482412ba7SHong Zhang   api[0] = 0;
14583445d95SHong Zhang 
14683445d95SHong Zhang   /* create and initialize a linked list */
14783445d95SHong Zhang   nlnk = pN+1;
14883445d95SHong Zhang   ierr = PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
14983445d95SHong Zhang 
15082412ba7SHong Zhang   /* Initial FreeSpace size is fill*nnz(A) */
151de0260b3SHong Zhang   ierr = GetMoreSpace((PetscInt)(fill*(adi[am]+aoi[am])),&free_space);CHKERRQ(ierr);
152f4a743e1SHong Zhang   current_space = free_space;
153f4a743e1SHong Zhang 
154f4a743e1SHong Zhang   for (i=0;i<am;i++) {
15582412ba7SHong Zhang     apnz = 0;
156f4a743e1SHong Zhang     /* diagonal portion of A */
157ded4f561SHong Zhang     nzi = adi[i+1] - adi[i];
158ded4f561SHong Zhang     for (j=0; j<nzi; j++){
15982412ba7SHong Zhang       row = *adj++;
16082412ba7SHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
161ded4f561SHong Zhang       Jptr  = pj_loc + pi_loc[row];
16283445d95SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
163ded4f561SHong Zhang       ierr = PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
16482412ba7SHong Zhang       apnz += nlnk;
165f4a743e1SHong Zhang     }
166f4a743e1SHong Zhang     /* off-diagonal portion of A */
167ded4f561SHong Zhang     nzi = aoi[i+1] - aoi[i];
168ded4f561SHong Zhang     for (j=0; j<nzi; j++){
16982412ba7SHong Zhang       row = *aoj++;
17082412ba7SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
171ded4f561SHong Zhang       Jptr  = pj_oth + pi_oth[row];
172ded4f561SHong Zhang       ierr = PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
17382412ba7SHong Zhang       apnz += nlnk;
174f4a743e1SHong Zhang     }
175f4a743e1SHong Zhang 
17682412ba7SHong Zhang     api[i+1] = api[i] + apnz;
17782412ba7SHong Zhang     if (ap->abnz_max < apnz) ap->abnz_max = apnz;
178f4a743e1SHong Zhang 
17983445d95SHong Zhang     /* if free space is not available, double the total space in the list */
18082412ba7SHong Zhang     if (current_space->local_remaining<apnz) {
181f4a743e1SHong Zhang       ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
182f4a743e1SHong Zhang     }
183f4a743e1SHong Zhang 
184f4a743e1SHong Zhang     /* Copy data into free space, then initialize lnk */
18582412ba7SHong Zhang     ierr = PetscLLClean(pN,pN,apnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
18682412ba7SHong Zhang     current_space->array           += apnz;
18782412ba7SHong Zhang     current_space->local_used      += apnz;
18882412ba7SHong Zhang     current_space->local_remaining -= apnz;
189f4a743e1SHong Zhang   }
19082412ba7SHong Zhang   /* Allocate space for apj, initialize apj, and */
191f4a743e1SHong Zhang   /* destroy list of free space and other temporary array(s) */
19282412ba7SHong Zhang   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ap->abj);CHKERRQ(ierr);
19382412ba7SHong Zhang   apj = ap->abj;
194de0260b3SHong Zhang   ierr = MakeSpaceContiguous(&free_space,ap->abj);CHKERRQ(ierr);
195f4a743e1SHong Zhang 
196ded4f561SHong Zhang   /* determine symbolic Co=(p->B)^T*AP - send to others */
197ded4f561SHong Zhang   /*----------------------------------------------------*/
198de0260b3SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
19942c57489SHong Zhang 
200ded4f561SHong Zhang   /* then, compute symbolic Co = (p->B)^T*AP */
201de0260b3SHong Zhang   pon = (p->B)->n; /* total num of rows to be sent to other processors
20283445d95SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
203de0260b3SHong Zhang   ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr);
204de0260b3SHong Zhang   coi[0] = 0;
20542c57489SHong Zhang 
20682412ba7SHong Zhang   /* set initial free space to be 3*pon*max( nnz(AP) per row) */
20782412ba7SHong Zhang   nnz           = 3*pon*ap->abnz_max + 1;
208de0260b3SHong Zhang   ierr          = GetMoreSpace(nnz,&free_space);
20942c57489SHong Zhang   current_space = free_space;
21042c57489SHong Zhang 
211de0260b3SHong Zhang   for (i=0; i<pon; i++) {
212ded4f561SHong Zhang     nnz  = 0;
21382412ba7SHong Zhang     pnz = poti[i+1] - poti[i];
21482412ba7SHong Zhang     j     = pnz;
215de0260b3SHong Zhang     ptJ   = potj + poti[i+1];
21683445d95SHong Zhang     while (j){/* assume cols are almost in increasing order, starting from its end saves computation */
21783445d95SHong Zhang       j--; ptJ--;
21882412ba7SHong Zhang       row  = *ptJ; /* row of AP == col of Pot */
21982412ba7SHong Zhang       apnz = api[row+1] - api[row];
220ded4f561SHong Zhang       Jptr   = apj + api[row];
22183445d95SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
222ded4f561SHong Zhang       ierr = PetscLLAdd(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
223ded4f561SHong Zhang       nnz += nlnk;
22442c57489SHong Zhang     }
22542c57489SHong Zhang 
22683445d95SHong Zhang     /* If free space is not available, double the total space in the list */
227ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
22842c57489SHong Zhang       ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
22942c57489SHong Zhang     }
23042c57489SHong Zhang 
23142c57489SHong Zhang     /* Copy data into free space, and zero out denserows */
232ded4f561SHong Zhang     ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
233ded4f561SHong Zhang     current_space->array           += nnz;
234ded4f561SHong Zhang     current_space->local_used      += nnz;
235ded4f561SHong Zhang     current_space->local_remaining -= nnz;
236ded4f561SHong Zhang     coi[i+1] = coi[i] + nnz;
23742c57489SHong Zhang   }
238de0260b3SHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr);
239de0260b3SHong Zhang   ierr = MakeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
240de0260b3SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
24142c57489SHong Zhang 
242ded4f561SHong Zhang   /* send j-array (coj) of Co to other processors */
243ded4f561SHong Zhang   /*----------------------------------------------*/
24442c57489SHong Zhang   /* determine row ownership */
24583445d95SHong Zhang   ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr);
24642c57489SHong Zhang   ierr = PetscMapCreate(comm,&merge->rowmap);CHKERRQ(ierr);
24742c57489SHong Zhang   ierr = PetscMapSetLocalSize(merge->rowmap,pn);CHKERRQ(ierr);
24842c57489SHong Zhang   ierr = PetscMapSetType(merge->rowmap,MAP_MPI);CHKERRQ(ierr);
24942c57489SHong Zhang   ierr = PetscMapGetGlobalRange(merge->rowmap,&owners);CHKERRQ(ierr);
25042c57489SHong Zhang 
25142c57489SHong Zhang   /* determine the number of messages to send, their lengths */
25242c57489SHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr);
25383445d95SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
25442c57489SHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr);
25542c57489SHong Zhang   len_s = merge->len_s;
25642c57489SHong Zhang   merge->nsend = 0;
25783445d95SHong Zhang 
258de0260b3SHong Zhang   ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr);
259de0260b3SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
260de0260b3SHong Zhang 
26183445d95SHong Zhang   proc = 0;
262de0260b3SHong Zhang   for (i=0; i<pon; i++){
263de0260b3SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
264de0260b3SHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
265de0260b3SHong Zhang     len_s[proc] += coi[i+1] - coi[i];
26683445d95SHong Zhang   }
267de0260b3SHong Zhang 
268de0260b3SHong Zhang   len   = 0;  /* max length of buf_si[] */
269de0260b3SHong Zhang   owners_co[0] = 0;
27042c57489SHong Zhang   for (proc=0; proc<size; proc++){
271de0260b3SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
27283445d95SHong Zhang     if (len_si[proc]){
27342c57489SHong Zhang       merge->nsend++;
27483445d95SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
27542c57489SHong Zhang       len += len_si[proc];
27642c57489SHong Zhang     }
27742c57489SHong Zhang   }
27842c57489SHong Zhang 
279ded4f561SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
28042c57489SHong Zhang   ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
28142c57489SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
28242c57489SHong Zhang 
283ded4f561SHong Zhang   /* post the Irecv and Isend of coj */
284ded4f561SHong Zhang   ierr = PetscObjectGetNewTag((PetscObject)merge->rowmap,&tag);CHKERRQ(ierr);
285ded4f561SHong Zhang   ierr = PetscPostIrecvInt(comm,tag,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
28642c57489SHong Zhang 
287ded4f561SHong Zhang   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr);
28842c57489SHong Zhang 
28942c57489SHong Zhang   for (proc=0, k=0; proc<size; proc++){
29042c57489SHong Zhang     if (!len_s[proc]) continue;
291de0260b3SHong Zhang     i = owners_co[proc];
292ded4f561SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tag,comm,swaits+k);CHKERRQ(ierr);
29342c57489SHong Zhang     k++;
29442c57489SHong Zhang   }
29542c57489SHong Zhang 
296ded4f561SHong Zhang   /* receives and sends of coj are complete */
297ded4f561SHong Zhang   ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr);
298ded4f561SHong Zhang   i = merge->nrecv;
299ded4f561SHong Zhang   while (i--) {
300ded4f561SHong Zhang     ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr);
301ded4f561SHong Zhang   }
302ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
303ded4f561SHong Zhang   if (merge->nsend){
304ded4f561SHong Zhang     ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);
30582412ba7SHong Zhang   }
30682412ba7SHong Zhang 
307ded4f561SHong Zhang   /* send and recv coi */
308ded4f561SHong Zhang   /*-------------------*/
309ded4f561SHong Zhang   ierr = PetscPostIrecvInt(comm,tag,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
31042c57489SHong Zhang 
31142c57489SHong Zhang   ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr);
31242c57489SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
31342c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++){
31442c57489SHong Zhang     if (!len_s[proc]) continue;
31542c57489SHong Zhang     /* form outgoing message for i-structure:
31642c57489SHong Zhang          buf_si[0]:                 nrows to be sent
31742c57489SHong Zhang                [1:nrows]:           row index (global)
31842c57489SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
31942c57489SHong Zhang     */
32042c57489SHong Zhang     /*-------------------------------------------*/
32142c57489SHong Zhang     nrows = len_si[proc]/2 - 1;
32242c57489SHong Zhang     buf_si_i    = buf_si + nrows+1;
32342c57489SHong Zhang     buf_si[0]   = nrows;
32442c57489SHong Zhang     buf_si_i[0] = 0;
32542c57489SHong Zhang     nrows = 0;
326de0260b3SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++){
327ded4f561SHong Zhang       nzi = coi[i+1] - coi[i];
328ded4f561SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
329de0260b3SHong Zhang       buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */
33042c57489SHong Zhang       nrows++;
33142c57489SHong Zhang     }
332ded4f561SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tag,comm,swaits+k);CHKERRQ(ierr);
33342c57489SHong Zhang     k++;
33442c57489SHong Zhang     buf_si += len_si[proc];
33542c57489SHong Zhang   }
336ded4f561SHong Zhang   i = merge->nrecv;
337ded4f561SHong Zhang   while (i--) {
338ded4f561SHong Zhang     ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr);
339ded4f561SHong Zhang   }
340ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
341ded4f561SHong Zhang   if (merge->nsend){
342ded4f561SHong Zhang     ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);
343ded4f561SHong Zhang   }
34442c57489SHong Zhang 
34542c57489SHong Zhang   ierr = PetscLogInfo((PetscObject)A,"MatMerge_SeqsToMPI: nsend: %d, nrecv: %d\n",merge->nsend,merge->nrecv);CHKERRQ(ierr);
34642c57489SHong Zhang   for (i=0; i<merge->nrecv; i++){
34742c57489SHong Zhang     ierr = PetscLogInfo((PetscObject)A,"MatMerge_SeqsToMPI:   recv len_ri=%d, len_rj=%d from [%d]\n",len_ri[i],merge->len_r[i],merge->id_r[i]);CHKERRQ(ierr);
34842c57489SHong Zhang   }
34942c57489SHong Zhang 
35042c57489SHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
35142c57489SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
352ded4f561SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
353ded4f561SHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
35442c57489SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
35542c57489SHong Zhang 
356de0260b3SHong Zhang   /* compute the local portion of C (mpi mat) */
357de0260b3SHong Zhang   /*------------------------------------------*/
358ded4f561SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
359ded4f561SHong Zhang 
36042c57489SHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
36142c57489SHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
36242c57489SHong Zhang   bi[0] = 0;
36342c57489SHong Zhang 
364ded4f561SHong Zhang   /* set initial free space to be 3*pn*max( nnz(AP) per row) */
365ded4f561SHong Zhang   nnz           = 3*pn*ap->abnz_max + 1;
366ded4f561SHong Zhang   ierr          = GetMoreSpace(nnz,&free_space);
36742c57489SHong Zhang   current_space = free_space;
36842c57489SHong Zhang 
36942c57489SHong Zhang   ierr = PetscMalloc((3*merge->nrecv+1)*sizeof(PetscInt**),&buf_ri_k);CHKERRQ(ierr);
37042c57489SHong Zhang   nextrow = buf_ri_k + merge->nrecv;
37142c57489SHong Zhang   nextci  = nextrow + merge->nrecv;
37242c57489SHong Zhang   for (k=0; k<merge->nrecv; k++){
37342c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
37442c57489SHong Zhang     nrows       = *buf_ri_k[k];
37542c57489SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
37642c57489SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
37742c57489SHong Zhang   }
37842c57489SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
37942c57489SHong Zhang   for (i=0; i<pn; i++) {
380ded4f561SHong Zhang     /* add pdt[i,:]*AP into lnk */
381ded4f561SHong Zhang     nnz = 0;
382ded4f561SHong Zhang     pnz  = pdti[i+1] - pdti[i];
383ded4f561SHong Zhang     j    = pnz;
384ded4f561SHong Zhang     ptJ  = pdtj + pdti[i+1];
385ded4f561SHong Zhang     while (j){/* assume cols are almost in increasing order, starting from its end saves computation */
386ded4f561SHong Zhang       j--; ptJ--;
387ded4f561SHong Zhang       row  = *ptJ; /* row of AP == col of Pt */
388ded4f561SHong Zhang       apnz = api[row+1] - api[row];
389ded4f561SHong Zhang       Jptr   = apj + api[row];
390ded4f561SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
391ded4f561SHong Zhang       ierr = PetscLLAdd(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
392ded4f561SHong Zhang       nnz += nlnk;
393ded4f561SHong Zhang     }
39442c57489SHong Zhang     /* add received col data into lnk */
39542c57489SHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
39642c57489SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
397ded4f561SHong Zhang         nzi = *(nextci[k]+1) - *nextci[k];
398ded4f561SHong Zhang         Jptr  = buf_rj[k] + *nextci[k];
399ded4f561SHong Zhang         ierr = PetscLLAdd(nzi,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
400ded4f561SHong Zhang         nnz += nlnk;
40142c57489SHong Zhang         nextrow[k]++; nextci[k]++;
40242c57489SHong Zhang       }
40342c57489SHong Zhang     }
40442c57489SHong Zhang 
40542c57489SHong Zhang     /* if free space is not available, make more free space */
406ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
40742c57489SHong Zhang       ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
40842c57489SHong Zhang     }
40942c57489SHong Zhang     /* copy data into free space, then initialize lnk */
410ded4f561SHong Zhang     ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
411ded4f561SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
412ded4f561SHong Zhang     current_space->array           += nnz;
413ded4f561SHong Zhang     current_space->local_used      += nnz;
414ded4f561SHong Zhang     current_space->local_remaining -= nnz;
415ded4f561SHong Zhang     bi[i+1] = bi[i] + nnz;
41642c57489SHong Zhang   }
417ded4f561SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
41842c57489SHong Zhang   ierr = PetscFree(buf_ri_k);CHKERRQ(ierr);
41942c57489SHong Zhang 
42042c57489SHong Zhang   ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
42142c57489SHong Zhang   ierr = MakeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
42242c57489SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
42342c57489SHong Zhang 
42442c57489SHong Zhang   /* create symbolic parallel matrix B_mpi */
42542c57489SHong Zhang   /*---------------------------------------*/
42642c57489SHong Zhang   ierr = MatCreate(comm,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE,&B_mpi);CHKERRQ(ierr);
42742c57489SHong Zhang   ierr = MatSetType(B_mpi,MATMPIAIJ);CHKERRQ(ierr);
42842c57489SHong Zhang   ierr = MatMPIAIJSetPreallocation(B_mpi,0,dnz,0,onz);CHKERRQ(ierr);
42942c57489SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
43042c57489SHong Zhang 
43182412ba7SHong Zhang   /* B_mpi is not ready for use - assembly will be done by MatPtAPNumeric() */
43242c57489SHong Zhang   B_mpi->assembled     = PETSC_FALSE;
43342c57489SHong Zhang   B_mpi->ops->destroy  = MatDestroy_MPIAIJ_SeqsToMPI;
43442c57489SHong Zhang   merge->bi            = bi;
43542c57489SHong Zhang   merge->bj            = bj;
436de0260b3SHong Zhang   merge->coi           = coi;
437de0260b3SHong Zhang   merge->coj           = coj;
43842c57489SHong Zhang   merge->buf_ri        = buf_ri;
43942c57489SHong Zhang   merge->buf_rj        = buf_rj;
440de0260b3SHong Zhang   merge->owners_co     = owners_co;
44142c57489SHong Zhang 
44242c57489SHong Zhang   /* attach the supporting struct to B_mpi for reuse */
44342c57489SHong Zhang   ierr = PetscObjectContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
44442c57489SHong Zhang   ierr = PetscObjectContainerSetPointer(container,merge);CHKERRQ(ierr);
44542c57489SHong Zhang   ierr = PetscObjectCompose((PetscObject)B_mpi,"MatMergeSeqsToMPI",(PetscObject)container);CHKERRQ(ierr);
44642c57489SHong Zhang   *C = B_mpi;
44742c57489SHong Zhang 
44842c57489SHong Zhang   PetscFunctionReturn(0);
44942c57489SHong Zhang }
45042c57489SHong Zhang 
45142c57489SHong Zhang #undef __FUNCT__
45242c57489SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ"
45342c57489SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
45442c57489SHong Zhang {
45542c57489SHong Zhang   PetscErrorCode       ierr;
45642c57489SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
457de0260b3SHong Zhang   Mat_MatMatMultMPI    *ap;
45842c57489SHong Zhang   PetscObjectContainer cont_merge,cont_ptap;
459de0260b3SHong Zhang   Mat_MPIAIJ     *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data;
46042c57489SHong Zhang   Mat_SeqAIJ     *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
461de0260b3SHong Zhang   Mat_SeqAIJ     *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
46242c57489SHong Zhang   Mat_SeqAIJ     *p_loc,*p_oth;
46382412ba7SHong Zhang   PetscInt       *adi=ad->i,*aoi=ao->i,*adj=ad->j,*aoj=ao->j,*apJ,nextp,flops=0;
46482412ba7SHong Zhang   PetscInt       *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj;
46582412ba7SHong Zhang   PetscInt       i,j,k,anz,pnz,apnz,nextap,row,*cj;
46682412ba7SHong Zhang   MatScalar      *ada=ad->a,*aoa=ao->a,*apa,*pa,*ca,*pa_loc,*pa_oth;
46782412ba7SHong Zhang   PetscInt       am=A->m,cm=C->m,pon=(p->B)->n;
46842c57489SHong Zhang   MPI_Comm       comm=C->comm;
46942c57489SHong Zhang   PetscMPIInt    size,rank,taga,*len_s;
47082412ba7SHong Zhang   PetscInt       *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
47142c57489SHong Zhang   PetscInt       **buf_ri,**buf_rj;
47250f5bf86SHong Zhang   PetscInt       cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */
47342c57489SHong Zhang   MPI_Request    *s_waits,*r_waits;
47442c57489SHong Zhang   MPI_Status     *status;
47582412ba7SHong Zhang   MatScalar      **abuf_r,*ba_i,*pA,*coa,*ba;
47682412ba7SHong Zhang   PetscInt       *api,*apj,*coi,*coj;
47782412ba7SHong Zhang   PetscInt       *poJ=po->j,*pdJ=pd->j,pcstart=p->cstart,pcend=p->cend;
47842c57489SHong Zhang 
47942c57489SHong Zhang   PetscFunctionBegin;
48042c57489SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
48142c57489SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
48242c57489SHong Zhang 
48342c57489SHong Zhang   ierr = PetscObjectQuery((PetscObject)C,"MatMergeSeqsToMPI",(PetscObject *)&cont_merge);CHKERRQ(ierr);
48442c57489SHong Zhang   if (cont_merge) {
48542c57489SHong Zhang     ierr  = PetscObjectContainerGetPointer(cont_merge,(void **)&merge);CHKERRQ(ierr);
48642c57489SHong Zhang   } else {
48742c57489SHong Zhang     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix C does not posses an object container");
48842c57489SHong Zhang   }
48942c57489SHong Zhang   ierr = PetscObjectQuery((PetscObject)P,"MatPtAPMPI",(PetscObject *)&cont_ptap);CHKERRQ(ierr);
49042c57489SHong Zhang   if (cont_ptap) {
491de0260b3SHong Zhang     ierr  = PetscObjectContainerGetPointer(cont_ptap,(void **)&ap);CHKERRQ(ierr);
49242c57489SHong Zhang   } else {
49342c57489SHong Zhang     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix P does not posses an object container");
49442c57489SHong Zhang   }
49542c57489SHong Zhang 
49642c57489SHong Zhang   /* get data from symbolic products */
497de0260b3SHong Zhang   p_loc = (Mat_SeqAIJ*)(ap->B_loc)->data;
498de0260b3SHong Zhang   p_oth = (Mat_SeqAIJ*)(ap->B_oth)->data;
49982412ba7SHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a,pA=pa_loc;
50042c57489SHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
50142c57489SHong Zhang 
502de0260b3SHong Zhang   coi = merge->coi; coj = merge->coj;
503de0260b3SHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr);
504de0260b3SHong Zhang   ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr);
505de0260b3SHong Zhang 
506de0260b3SHong Zhang   bi = merge->bi; bj = merge->bj;
507de0260b3SHong Zhang   ierr  = PetscMapGetGlobalRange(merge->rowmap,&owners);CHKERRQ(ierr);
508de0260b3SHong Zhang   ierr  = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
509de0260b3SHong Zhang   ierr  = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr);
51042c57489SHong Zhang 
511f4a743e1SHong Zhang   /* get data from symbolic A*P */
512de0260b3SHong Zhang   ierr = PetscMalloc((ap->abnz_max+1)*sizeof(MatScalar),&apa);CHKERRQ(ierr);
513de0260b3SHong Zhang   ierr = PetscMemzero(apa,ap->abnz_max*sizeof(MatScalar));CHKERRQ(ierr);
51442c57489SHong Zhang 
515f4a743e1SHong Zhang   /* compute numeric C_seq=P_loc^T*A_loc*P */
51682412ba7SHong Zhang   api = ap->abi; apj = ap->abj;
51742c57489SHong Zhang   for (i=0;i<am;i++) {
518f4a743e1SHong Zhang     /* form i-th sparse row of A*P */
51982412ba7SHong Zhang     apnz = api[i+1] - api[i];
52082412ba7SHong Zhang     apJ  = apj + api[i];
52142c57489SHong Zhang     /* diagonal portion of A */
52282412ba7SHong Zhang     anz  = adi[i+1] - adi[i];
52382412ba7SHong Zhang     for (j=0;j<anz;j++) {
52482412ba7SHong Zhang       row = *adj++;
52582412ba7SHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
52682412ba7SHong Zhang       pj  = pj_loc + pi_loc[row];
52782412ba7SHong Zhang       pa  = pa_loc + pi_loc[row];
528f4a743e1SHong Zhang       nextp = 0;
52982412ba7SHong Zhang       for (k=0; nextp<pnz; k++) {
53082412ba7SHong Zhang         if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
53182412ba7SHong Zhang           apa[k] += (*ada)*pa[nextp++];
53242c57489SHong Zhang         }
53342c57489SHong Zhang       }
53482412ba7SHong Zhang       flops += 2*pnz;
53542c57489SHong Zhang       ada++;
53642c57489SHong Zhang     }
53742c57489SHong Zhang     /* off-diagonal portion of A */
53882412ba7SHong Zhang     anz  = aoi[i+1] - aoi[i];
53982412ba7SHong Zhang     for (j=0; j<anz; j++) {
54082412ba7SHong Zhang       row = *aoj++;
54182412ba7SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
54282412ba7SHong Zhang       pj  = pj_oth + pi_oth[row];
54382412ba7SHong Zhang       pa  = pa_oth + pi_oth[row];
544f4a743e1SHong Zhang       nextp = 0;
54582412ba7SHong Zhang       for (k=0; nextp<pnz; k++) {
54682412ba7SHong Zhang         if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
54782412ba7SHong Zhang           apa[k] += (*aoa)*pa[nextp++];
54842c57489SHong Zhang         }
54942c57489SHong Zhang       }
55082412ba7SHong Zhang       flops += 2*pnz;
55142c57489SHong Zhang       aoa++;
55242c57489SHong Zhang     }
55342c57489SHong Zhang 
55482412ba7SHong Zhang     /* Compute P_loc[i,:]^T*AP[i,:] using outer product */
55582412ba7SHong Zhang     pnz = pi_loc[i+1] - pi_loc[i];
55682412ba7SHong Zhang     for (j=0; j<pnz; j++) {
55742c57489SHong Zhang       nextap = 0;
55882412ba7SHong Zhang       row    = *pJ++; /* global index */
55982412ba7SHong Zhang       if (row < pcstart || row >=pcend) { /* put the value into Co */
56082412ba7SHong Zhang         cj  = coj + coi[*poJ];
56182412ba7SHong Zhang         ca  = coa + coi[*poJ++];
56282412ba7SHong Zhang       } else {                            /* put the value into Cd */
56382412ba7SHong Zhang         cj   = bj + bi[*pdJ];
56482412ba7SHong Zhang         ca   = ba + bi[*pdJ++];
56542c57489SHong Zhang       }
56682412ba7SHong Zhang       for (k=0; nextap<apnz; k++) {
56782412ba7SHong Zhang         if (cj[k]==apJ[nextap]) ca[k] += (*pA)*apa[nextap++];
56842c57489SHong Zhang       }
56982412ba7SHong Zhang       flops += 2*apnz;
57082412ba7SHong Zhang       pA++;
571de0260b3SHong Zhang     }
572de0260b3SHong Zhang 
57342c57489SHong Zhang     /* zero the current row info for A*P */
57482412ba7SHong Zhang     ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr);
57542c57489SHong Zhang   }
57642c57489SHong Zhang   ierr = PetscFree(apa);CHKERRQ(ierr);
57742c57489SHong Zhang 
578de0260b3SHong Zhang   /* send and recv matrix values */
579de0260b3SHong Zhang   /*-----------------------------*/
58042c57489SHong Zhang   buf_ri = merge->buf_ri;
58142c57489SHong Zhang   buf_rj = merge->buf_rj;
58242c57489SHong Zhang   len_s  = merge->len_s;
58342c57489SHong Zhang   ierr = PetscObjectGetNewTag((PetscObject)merge->rowmap,&taga);CHKERRQ(ierr);
58442c57489SHong Zhang   ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
58542c57489SHong Zhang 
58642c57489SHong Zhang   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&s_waits);CHKERRQ(ierr);
58742c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++){
58842c57489SHong Zhang     if (!len_s[proc]) continue;
589de0260b3SHong Zhang     i = merge->owners_co[proc];
590de0260b3SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
59142c57489SHong Zhang     k++;
59242c57489SHong Zhang   }
59382412ba7SHong Zhang   ierr = PetscMalloc(size*sizeof(MPI_Status),&status);CHKERRQ(ierr);
59442c57489SHong Zhang   ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);
59542c57489SHong Zhang   ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);
59642c57489SHong Zhang   ierr = PetscFree(status);CHKERRQ(ierr);
59742c57489SHong Zhang 
59842c57489SHong Zhang   ierr = PetscFree(s_waits);CHKERRQ(ierr);
59942c57489SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
60082412ba7SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
60142c57489SHong Zhang 
60282412ba7SHong Zhang   /* insert local and received values into C */
60382412ba7SHong Zhang   /*-----------------------------------------*/
60442c57489SHong Zhang   ierr = PetscMalloc((3*merge->nrecv+1)*sizeof(PetscInt**),&buf_ri_k);CHKERRQ(ierr);
60542c57489SHong Zhang   nextrow   = buf_ri_k + merge->nrecv;
60682412ba7SHong Zhang   nextci = nextrow + merge->nrecv;
60742c57489SHong Zhang 
60842c57489SHong Zhang   for (k=0; k<merge->nrecv; k++){
60942c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
61042c57489SHong Zhang     nrows       = *(buf_ri_k[k]);
61142c57489SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
61282412ba7SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
61342c57489SHong Zhang   }
61442c57489SHong Zhang 
615de0260b3SHong Zhang   for (i=0; i<cm; i++) {
61682412ba7SHong Zhang     row = owners[rank] + i; /* global row index of C_seq */
61742c57489SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
618de0260b3SHong Zhang     ba_i = ba + bi[i];
61982412ba7SHong Zhang     bnz  = bi[i+1] - bi[i];
62042c57489SHong Zhang     /* add received vals into ba */
62142c57489SHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
62242c57489SHong Zhang       /* i-th row */
62342c57489SHong Zhang       if (i == *nextrow[k]) {
62482412ba7SHong Zhang         cnz = *(nextci[k]+1) - *nextci[k];
62582412ba7SHong Zhang         cj  = buf_rj[k] + *(nextci[k]);
62682412ba7SHong Zhang         ca  = abuf_r[k] + *(nextci[k]);
62782412ba7SHong Zhang         nextcj = 0;
62882412ba7SHong Zhang         for (j=0; nextcj<cnz; j++){
62982412ba7SHong Zhang           if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */
63082412ba7SHong Zhang             ba_i[j] += ca[nextcj++];
63142c57489SHong Zhang           }
63242c57489SHong Zhang         }
63382412ba7SHong Zhang         nextrow[k]++; nextci[k]++;
63442c57489SHong Zhang       }
63542c57489SHong Zhang     }
63682412ba7SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
63782412ba7SHong Zhang     flops += 2*cnz;
63842c57489SHong Zhang   }
63982412ba7SHong Zhang   ierr = MatSetOption(C,MAT_COLUMNS_SORTED);CHKERRQ(ierr); /* -cause delay? */
64042c57489SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
64142c57489SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
64242c57489SHong Zhang 
643de0260b3SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
64442c57489SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
64542c57489SHong Zhang   ierr = PetscFree(buf_ri_k);CHKERRQ(ierr);
64642c57489SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
64742c57489SHong Zhang 
64842c57489SHong Zhang   PetscFunctionReturn(0);
64942c57489SHong Zhang }
650