xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision 83445d950c9699bb132b1db906a3c0e4ae82194c)
142c57489SHong Zhang /*
242c57489SHong Zhang   Defines projective product routines where A is a MPIAIJ matrix
342c57489SHong Zhang           C = P^T * A * P
442c57489SHong Zhang */
542c57489SHong Zhang 
642c57489SHong Zhang #include "src/mat/impls/aij/seq/aij.h"   /*I "petscmat.h" I*/
742c57489SHong Zhang #include "src/mat/utils/freespace.h"
842c57489SHong Zhang #include "src/mat/impls/aij/mpi/mpiaij.h"
942c57489SHong Zhang #include "petscbt.h"
1042c57489SHong Zhang 
1142c57489SHong Zhang EXTERN PetscErrorCode MatDestroy_MPIAIJ(Mat);
1242c57489SHong Zhang #undef __FUNCT__
1342c57489SHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_PtAP"
1442c57489SHong Zhang PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A)
1542c57489SHong Zhang {
1642c57489SHong Zhang   PetscErrorCode       ierr;
1742c57489SHong Zhang   Mat_MatMatMultMPI    *ptap;
1842c57489SHong Zhang   PetscObjectContainer container;
1942c57489SHong Zhang 
2042c57489SHong Zhang   PetscFunctionBegin;
2142c57489SHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"MatPtAPMPI",(PetscObject *)&container);CHKERRQ(ierr);
2242c57489SHong Zhang   if (container) {
2342c57489SHong Zhang     ierr  = PetscObjectContainerGetPointer(container,(void **)&ptap);CHKERRQ(ierr);
2442c57489SHong Zhang     ierr = MatDestroy(ptap->B_loc);CHKERRQ(ierr);
2542c57489SHong Zhang     ierr = MatDestroy(ptap->B_oth);CHKERRQ(ierr);
2642c57489SHong Zhang     ierr = ISDestroy(ptap->isrowb);CHKERRQ(ierr);
2742c57489SHong Zhang     ierr = ISDestroy(ptap->iscolb);CHKERRQ(ierr);
2842c57489SHong Zhang     if (ptap->abi) ierr = PetscFree(ptap->abi);CHKERRQ(ierr);
2942c57489SHong Zhang     if (ptap->abj) ierr = PetscFree(ptap->abj);CHKERRQ(ierr);
3042c57489SHong Zhang 
3142c57489SHong Zhang     ierr = PetscObjectContainerDestroy(container);CHKERRQ(ierr);
3242c57489SHong Zhang     ierr = PetscObjectCompose((PetscObject)A,"MatPtAPMPI",0);CHKERRQ(ierr);
3342c57489SHong Zhang   }
3442c57489SHong Zhang   ierr = PetscFree(ptap);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;
4542c57489SHong Zhang   Mat_MatMatMultMPI    *ptap;
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);
5142c57489SHong Zhang     ierr = PetscNew(Mat_MatMatMultMPI,&ptap);CHKERRQ(ierr);
5242c57489SHong Zhang     ptap->abi=PETSC_NULL; ptap->abj=PETSC_NULL;
5342c57489SHong Zhang     ptap->abnz_max = 0; /* symbolic A*P is not done yet */
5442c57489SHong Zhang 
5542c57489SHong Zhang     /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
5642c57489SHong Zhang     ierr = MatGetBrowsOfAoCols(A,P,scall,&ptap->isrowb,&ptap->iscolb,&ptap->brstart,&ptap->B_oth);CHKERRQ(ierr);
5742c57489SHong Zhang 
5842c57489SHong Zhang     /* get P_loc by taking all local rows of P */
5942c57489SHong Zhang     ierr = MatGetLocalMat(P,scall,&ptap->B_loc);CHKERRQ(ierr);
6042c57489SHong Zhang 
6142c57489SHong Zhang     /* attach the supporting struct to P for reuse */
6242c57489SHong Zhang     P->ops->destroy  = MatDestroy_MPIAIJ_PtAP;
6342c57489SHong Zhang     ierr = PetscObjectContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
6442c57489SHong Zhang     ierr = PetscObjectContainerSetPointer(container,ptap);CHKERRQ(ierr);
6542c57489SHong Zhang     ierr = PetscObjectCompose((PetscObject)P,"MatPtAPMPI",(PetscObject)container);CHKERRQ(ierr);
6642c57489SHong Zhang 
6742c57489SHong Zhang     /* now, compute symbolic local P^T*A*P */
6842c57489SHong Zhang     ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);/* numeric product is computed as well */
6942c57489SHong Zhang     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
7042c57489SHong Zhang   } else if (scall == MAT_REUSE_MATRIX){
7142c57489SHong Zhang     ierr = PetscObjectQuery((PetscObject)P,"MatPtAPMPI",(PetscObject *)&container);CHKERRQ(ierr);
7242c57489SHong Zhang     if (container) {
7342c57489SHong Zhang       ierr  = PetscObjectContainerGetPointer(container,(void **)&ptap);CHKERRQ(ierr);
7442c57489SHong Zhang     } else {
7542c57489SHong Zhang       SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix P does not posses an object container");
7642c57489SHong Zhang     }
7742c57489SHong Zhang 
7842c57489SHong Zhang     /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
7942c57489SHong Zhang     ierr = MatGetBrowsOfAoCols(A,P,scall,&ptap->isrowb,&ptap->iscolb,&ptap->brstart,&ptap->B_oth);CHKERRQ(ierr);
8042c57489SHong Zhang 
8142c57489SHong Zhang     /* get P_loc by taking all local rows of P */
8242c57489SHong Zhang     ierr = MatGetLocalMat(P,scall,&ptap->B_loc);CHKERRQ(ierr);
8342c57489SHong Zhang 
8442c57489SHong Zhang   } else {
8542c57489SHong Zhang     SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",scall);
8642c57489SHong Zhang   }
8742c57489SHong Zhang   /* now, compute numeric local P^T*A*P */
8842c57489SHong Zhang   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
8942c57489SHong Zhang   ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ(A,P,*C);CHKERRQ(ierr);
9042c57489SHong Zhang   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
9142c57489SHong Zhang 
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;
10042c57489SHong Zhang   Mat                  P_loc,P_oth;
10142c57489SHong Zhang   Mat_MatMatMultMPI    *ptap;
10242c57489SHong Zhang   PetscObjectContainer container;
10342c57489SHong Zhang   FreeSpaceList        free_space=PETSC_NULL,current_space=PETSC_NULL;
104*83445d95SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data;
10542c57489SHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
10642c57489SHong Zhang   Mat_SeqAIJ           *p_loc,*p_oth;
107*83445d95SHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pti,*ptj,*ptJ,*pjj;
108*83445d95SHong Zhang   PetscInt             *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,nnz;
109*83445d95SHong Zhang   PetscInt             pm=P->m,pn=P->n,nlnk,*lnk,*ci,*cj,*cji;
11042c57489SHong Zhang   PetscInt             aN=A->N,am=A->m,pN=P->N;
111*83445d95SHong Zhang   PetscInt             i,j,k,ptnzi,arow,prow,pnzj,cnzi;
11242c57489SHong Zhang   PetscBT              lnkbt;
113*83445d95SHong Zhang   PetscInt             prstart,prend;
11442c57489SHong Zhang   MPI_Comm             comm=A->comm;
11542c57489SHong Zhang   Mat                  B_mpi;
11642c57489SHong Zhang   PetscMPIInt          size,rank,tagi,tagj,*len_s,*len_si,*len_ri;
11742c57489SHong Zhang   PetscInt             **buf_rj,**buf_ri,**buf_ri_k;
11842c57489SHong Zhang   PetscInt             len,proc,*dnz,*onz,*owners;
11942c57489SHong Zhang   PetscInt             anzi,*bi,*bj,bnzi,nspacedouble=0;
12042c57489SHong Zhang   PetscInt             nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
12142c57489SHong Zhang   MPI_Request          *si_waits,*sj_waits,*ri_waits,*rj_waits;
12242c57489SHong Zhang   MPI_Status           *status;
12342c57489SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
124*83445d95SHong Zhang   PetscInt             *apsymi,*apsymj,*apj,apnzj,*rmap=p->garray,tnrows;
125*83445d95SHong Zhang   /*  PetscInt             tnzrows,pcstart=p->cstart,pcend=p->cend; */
12642c57489SHong Zhang 
12742c57489SHong Zhang   PetscFunctionBegin;
12842c57489SHong Zhang   ierr = PetscObjectQuery((PetscObject)P,"MatPtAPMPI",(PetscObject *)&container);CHKERRQ(ierr);
12942c57489SHong Zhang   if (container) {
13042c57489SHong Zhang     ierr  = PetscObjectContainerGetPointer(container,(void **)&ptap);CHKERRQ(ierr);
13142c57489SHong Zhang   } else {
13242c57489SHong Zhang     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix P does not posses an object container");
13342c57489SHong Zhang   }
134*83445d95SHong Zhang 
135*83445d95SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
136*83445d95SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
137*83445d95SHong Zhang 
138*83445d95SHong Zhang   /* compute symbolic C_seq = P_loc^T * A_loc * P */
139*83445d95SHong Zhang   /*----------------------------------------------*/
14042c57489SHong Zhang   /* get data from P_loc and P_oth */
14142c57489SHong Zhang   P_loc=ptap->B_loc; P_oth=ptap->B_oth; prstart=ptap->brstart;
14242c57489SHong Zhang   p_loc=(Mat_SeqAIJ*)P_loc->data; p_oth=(Mat_SeqAIJ*)P_oth->data;
14342c57489SHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pi_oth=p_oth->i; pj_oth=p_oth->j;
14442c57489SHong Zhang   prend = prstart + pm;
14542c57489SHong Zhang 
146*83445d95SHong Zhang   /* first, compute symbolic AP = A_loc*P */
147f4a743e1SHong Zhang   ierr  = PetscMalloc((am+2)*sizeof(PetscInt),&apsymi);CHKERRQ(ierr);
148f4a743e1SHong Zhang   ptap->abi = apsymi;
149f4a743e1SHong Zhang   apsymi[0] = 0;
150*83445d95SHong Zhang 
151*83445d95SHong Zhang   /* create and initialize a linked list */
152*83445d95SHong Zhang   nlnk = pN+1;
153*83445d95SHong Zhang   ierr = PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
154*83445d95SHong Zhang 
155f4a743e1SHong Zhang   /* Initial FreeSpace size is 2*nnz(A) */
156f4a743e1SHong Zhang   ierr = GetMoreSpace((PetscInt)(2*(adi[am]+aoi[am])),&free_space);CHKERRQ(ierr);
157f4a743e1SHong Zhang   current_space = free_space;
158f4a743e1SHong Zhang 
159f4a743e1SHong Zhang   for (i=0;i<am;i++) {
160f4a743e1SHong Zhang     apnzj = 0;
161f4a743e1SHong Zhang     /* diagonal portion of A */
162f4a743e1SHong Zhang     anzi = adi[i+1] - adi[i];
163f4a743e1SHong Zhang     for (j=0; j<anzi; j++){
164*83445d95SHong Zhang       prow = *adj++;
165f4a743e1SHong Zhang       pnzj = pi_loc[prow+1] - pi_loc[prow];
166f4a743e1SHong Zhang       pjj  = pj_loc + pi_loc[prow];
167*83445d95SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
168*83445d95SHong Zhang       ierr = PetscLLAdd(pnzj,pjj,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
169*83445d95SHong Zhang       apnzj += nlnk;
170f4a743e1SHong Zhang     }
171f4a743e1SHong Zhang     /* off-diagonal portion of A */
172f4a743e1SHong Zhang     anzi = aoi[i+1] - aoi[i];
173f4a743e1SHong Zhang     for (j=0; j<anzi; j++){
174*83445d95SHong Zhang       prow = *aoj++;
175f4a743e1SHong Zhang       pnzj = pi_oth[prow+1] - pi_oth[prow];
176f4a743e1SHong Zhang       pjj  = pj_oth + pi_oth[prow];
177*83445d95SHong Zhang       ierr = PetscLLAdd(pnzj,pjj,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
178*83445d95SHong Zhang       apnzj += nlnk;
179f4a743e1SHong Zhang     }
180f4a743e1SHong Zhang 
181f4a743e1SHong Zhang     apsymi[i+1] = apsymi[i] + apnzj;
182f4a743e1SHong Zhang     if (ptap->abnz_max < apnzj) ptap->abnz_max = apnzj;
183f4a743e1SHong Zhang 
184*83445d95SHong Zhang     /* if free space is not available, double the total space in the list */
185f4a743e1SHong Zhang     if (current_space->local_remaining<apnzj) {
186f4a743e1SHong Zhang       ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
187f4a743e1SHong Zhang     }
188f4a743e1SHong Zhang 
189f4a743e1SHong Zhang     /* Copy data into free space, then initialize lnk */
190*83445d95SHong Zhang     ierr = PetscLLClean(pN,pN,apnzj,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
191f4a743e1SHong Zhang     current_space->array           += apnzj;
192f4a743e1SHong Zhang     current_space->local_used      += apnzj;
193f4a743e1SHong Zhang     current_space->local_remaining -= apnzj;
194f4a743e1SHong Zhang   }
195f4a743e1SHong Zhang   /* Allocate space for apsymj, initialize apsymj, and */
196f4a743e1SHong Zhang   /* destroy list of free space and other temporary array(s) */
197f4a743e1SHong Zhang   ierr = PetscMalloc((apsymi[am]+1)*sizeof(PetscInt),&ptap->abj);CHKERRQ(ierr);
198*83445d95SHong Zhang   apsymj = ptap->abj;
199f4a743e1SHong Zhang   ierr = MakeSpaceContiguous(&free_space,ptap->abj);CHKERRQ(ierr);
200f4a743e1SHong Zhang 
20142c57489SHong Zhang   /* get ij structure of P_loc^T */
20242c57489SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(P_loc,&pti,&ptj);CHKERRQ(ierr);
20342c57489SHong Zhang 
204*83445d95SHong Zhang   /* then, compute symbolic C_seq = P_loc^T*AP */
20542c57489SHong Zhang   /* Allocate ci array, arrays for fill computation and */
20642c57489SHong Zhang   /* free space for accumulating nonzero column info */
20742c57489SHong Zhang   ierr = PetscMalloc((pN+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
20842c57489SHong Zhang   ci[0] = 0;
20942c57489SHong Zhang 
210*83445d95SHong Zhang   /* tnzrows = 0; */
211*83445d95SHong Zhang   tnrows = (p->B)->N; /* total num of rows to be sent to other processors
212*83445d95SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
21342c57489SHong Zhang 
21442c57489SHong Zhang 
21542c57489SHong Zhang   /* Set initial free space to be nnz(A) scaled by fill*P->N/P->M. */
21642c57489SHong Zhang   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
21742c57489SHong Zhang   nnz           = adi[am] + aoi[am];
21842c57489SHong Zhang   ierr          = GetMoreSpace((PetscInt)(fill*nnz*pN/aN+1),&free_space);
21942c57489SHong Zhang   current_space = free_space;
22042c57489SHong Zhang 
22142c57489SHong Zhang   /* determine symbolic info for each row of C: */
22242c57489SHong Zhang   for (i=0; i<pN; i++) {
22342c57489SHong Zhang     cnzi  = 0;
224*83445d95SHong Zhang     ptnzi = pti[i+1] - pti[i];
225*83445d95SHong Zhang     if (ptnzi){
226*83445d95SHong Zhang     j     = ptnzi;
227*83445d95SHong Zhang     ptJ   = ptj + pti[i+1];
228*83445d95SHong Zhang     while (j){/* assume cols are almost in increasing order, starting from its end saves computation */
229*83445d95SHong Zhang       j--; ptJ--;
230*83445d95SHong Zhang       arow  = *ptJ; /* row of AP == col of Pt */
231*83445d95SHong Zhang       apnzj = apsymi[arow+1] - apsymi[arow];
232*83445d95SHong Zhang       apj   = apsymj + apsymi[arow];
233*83445d95SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
234*83445d95SHong Zhang       ierr = PetscLLAdd(apnzj,apj,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
23542c57489SHong Zhang       cnzi += nlnk;
23642c57489SHong Zhang     }
23742c57489SHong Zhang 
238*83445d95SHong Zhang     /* If free space is not available, double the total space in the list */
23942c57489SHong Zhang     if (current_space->local_remaining<cnzi) {
24042c57489SHong Zhang       ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
24142c57489SHong Zhang     }
24242c57489SHong Zhang 
24342c57489SHong Zhang     /* Copy data into free space, and zero out denserows */
24442c57489SHong Zhang     ierr = PetscLLClean(pN,pN,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
24542c57489SHong Zhang     current_space->array           += cnzi;
24642c57489SHong Zhang     current_space->local_used      += cnzi;
24742c57489SHong Zhang     current_space->local_remaining -= cnzi;
248*83445d95SHong Zhang     /* tnzrows++; */
249*83445d95SHong Zhang     } /* if (ptnzi) */
25042c57489SHong Zhang     ci[i+1] = ci[i] + cnzi;
251*83445d95SHong Zhang 
25242c57489SHong Zhang   }
253*83445d95SHong Zhang   /* printf("[%d] tnzrows-pn: %d,tnrows: %d, pN: %d; pcstart/end: %d, %d\n",rank,tnzrows-pn,tnrows,pN,pcstart,pcend);*/
254*83445d95SHong Zhang 
25542c57489SHong Zhang   /* Clean up. */
25642c57489SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P_loc,&pti,&ptj);CHKERRQ(ierr);
25742c57489SHong Zhang 
25842c57489SHong Zhang   /* Allocate space for cj, initialize cj, and */
25942c57489SHong Zhang   /* destroy list of free space and other temporary array(s) */
26042c57489SHong Zhang   ierr = PetscMalloc((ci[pN]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
26142c57489SHong Zhang   ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
26242c57489SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
26342c57489SHong Zhang 
26442c57489SHong Zhang   /* add C_seq into mpi C              */
26542c57489SHong Zhang   /*-----------------------------------*/
26642c57489SHong Zhang   /* determine row ownership */
267*83445d95SHong Zhang   ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr);
26842c57489SHong Zhang   ierr = PetscMapCreate(comm,&merge->rowmap);CHKERRQ(ierr);
26942c57489SHong Zhang   ierr = PetscMapSetLocalSize(merge->rowmap,pn);CHKERRQ(ierr);
27042c57489SHong Zhang   ierr = PetscMapSetType(merge->rowmap,MAP_MPI);CHKERRQ(ierr);
27142c57489SHong Zhang   ierr = PetscMapGetGlobalRange(merge->rowmap,&owners);CHKERRQ(ierr);
27242c57489SHong Zhang 
27342c57489SHong Zhang   /* determine the number of messages to send, their lengths */
27442c57489SHong Zhang   /*---------------------------------------------------------*/
27542c57489SHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr);
276*83445d95SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
27742c57489SHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr);
27842c57489SHong Zhang   len_s  = merge->len_s;
279*83445d95SHong Zhang   len = 0;  /* max length of buf_si[] */
28042c57489SHong Zhang   merge->nsend = 0;
281*83445d95SHong Zhang 
282*83445d95SHong Zhang   proc = 0;
283*83445d95SHong Zhang   for (i=0; i<tnrows; i++){
284*83445d95SHong Zhang     while (rmap[i] >= owners[proc+1]) proc++;
285*83445d95SHong Zhang     len_si[proc]++;
286*83445d95SHong Zhang   }
28742c57489SHong Zhang   for (proc=0; proc<size; proc++){
28842c57489SHong Zhang     len_s[proc] = 0;
289*83445d95SHong Zhang     if (len_si[proc]){
29042c57489SHong Zhang       merge->nsend++;
291*83445d95SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
292*83445d95SHong Zhang       len_s[proc] = ci[owners[proc+1]] - ci[owners[proc]]; /* num of col indices to be sent to [proc] */
29342c57489SHong Zhang       len += len_si[proc];
29442c57489SHong Zhang     }
29542c57489SHong Zhang   }
29642c57489SHong Zhang 
29742c57489SHong Zhang   /* determine the number and length of messages to receive for ij-structure */
29842c57489SHong Zhang   /*-------------------------------------------------------------------------*/
29942c57489SHong Zhang   ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
30042c57489SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
30142c57489SHong Zhang 
30242c57489SHong Zhang   /* post the Irecv of j-structure */
30342c57489SHong Zhang   /*-------------------------------*/
30442c57489SHong Zhang   ierr = PetscObjectGetNewTag((PetscObject)merge->rowmap,&tagj);CHKERRQ(ierr);
30542c57489SHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rj_waits);CHKERRQ(ierr);
30642c57489SHong Zhang 
30742c57489SHong Zhang   /* post the Isend of j-structure */
30842c57489SHong Zhang   /*--------------------------------*/
30942c57489SHong Zhang   ierr = PetscMalloc((2*merge->nsend+1)*sizeof(MPI_Request),&si_waits);CHKERRQ(ierr);
31042c57489SHong Zhang   sj_waits = si_waits + merge->nsend;
31142c57489SHong Zhang 
31242c57489SHong Zhang   for (proc=0, k=0; proc<size; proc++){
31342c57489SHong Zhang     if (!len_s[proc]) continue;
31442c57489SHong Zhang     i = owners[proc];
31542c57489SHong Zhang     ierr = MPI_Isend(cj+ci[i],len_s[proc],MPIU_INT,proc,tagj,comm,sj_waits+k);CHKERRQ(ierr);
31642c57489SHong Zhang     k++;
31742c57489SHong Zhang   }
31842c57489SHong Zhang 
31942c57489SHong Zhang   /* receives and sends of j-structure are complete */
32042c57489SHong Zhang   /*------------------------------------------------*/
32142c57489SHong Zhang   ierr = PetscMalloc(size*sizeof(MPI_Status),&status);CHKERRQ(ierr);
32242c57489SHong Zhang   ierr = MPI_Waitall(merge->nrecv,rj_waits,status);CHKERRQ(ierr);
32342c57489SHong Zhang   ierr = MPI_Waitall(merge->nsend,sj_waits,status);CHKERRQ(ierr);
32442c57489SHong Zhang 
32542c57489SHong Zhang   /* send and recv i-structure */
32642c57489SHong Zhang   /*---------------------------*/
32742c57489SHong Zhang   ierr = PetscObjectGetNewTag((PetscObject)merge->rowmap,&tagi);CHKERRQ(ierr);
32842c57489SHong Zhang   ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&ri_waits);CHKERRQ(ierr);
32942c57489SHong Zhang 
33042c57489SHong Zhang   ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr);
33142c57489SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
33242c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++){
33342c57489SHong Zhang     if (!len_s[proc]) continue;
33442c57489SHong Zhang     /* form outgoing message for i-structure:
33542c57489SHong Zhang          buf_si[0]:                 nrows to be sent
33642c57489SHong Zhang                [1:nrows]:           row index (global)
33742c57489SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
33842c57489SHong Zhang     */
33942c57489SHong Zhang     /*-------------------------------------------*/
34042c57489SHong Zhang     nrows = len_si[proc]/2 - 1;
34142c57489SHong Zhang     buf_si_i    = buf_si + nrows+1;
34242c57489SHong Zhang     buf_si[0]   = nrows;
34342c57489SHong Zhang     buf_si_i[0] = 0;
34442c57489SHong Zhang     nrows = 0;
34542c57489SHong Zhang     for (i=owners[proc]; i<owners[proc+1]; i++){
34642c57489SHong Zhang       anzi = ci[i+1] - ci[i];
34742c57489SHong Zhang       if (anzi) {
34842c57489SHong Zhang         buf_si_i[nrows+1] = buf_si_i[nrows] + anzi; /* i-structure */
34942c57489SHong Zhang         buf_si[nrows+1] = i-owners[proc]; /* local row index */
35042c57489SHong Zhang         nrows++;
35142c57489SHong Zhang       }
35242c57489SHong Zhang     }
35342c57489SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,si_waits+k);CHKERRQ(ierr);
35442c57489SHong Zhang     k++;
35542c57489SHong Zhang     buf_si += len_si[proc];
35642c57489SHong Zhang   }
35742c57489SHong Zhang 
35842c57489SHong Zhang   ierr = MPI_Waitall(merge->nrecv,ri_waits,status);CHKERRQ(ierr);
35942c57489SHong Zhang   ierr = MPI_Waitall(merge->nsend,si_waits,status);CHKERRQ(ierr);
36042c57489SHong Zhang 
36142c57489SHong Zhang   ierr = PetscLogInfo((PetscObject)A,"MatMerge_SeqsToMPI: nsend: %d, nrecv: %d\n",merge->nsend,merge->nrecv);CHKERRQ(ierr);
36242c57489SHong Zhang   for (i=0; i<merge->nrecv; i++){
36342c57489SHong 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);
36442c57489SHong Zhang   }
36542c57489SHong Zhang 
36642c57489SHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
36742c57489SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
36842c57489SHong Zhang   ierr = PetscFree(rj_waits);CHKERRQ(ierr);
36942c57489SHong Zhang   ierr = PetscFree(si_waits);CHKERRQ(ierr);
37042c57489SHong Zhang   ierr = PetscFree(ri_waits);CHKERRQ(ierr);
37142c57489SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
37242c57489SHong Zhang   ierr = PetscFree(status);CHKERRQ(ierr);
37342c57489SHong Zhang 
37442c57489SHong Zhang   /* compute a local seq matrix in each processor */
37542c57489SHong Zhang   /*----------------------------------------------*/
37642c57489SHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
37742c57489SHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
37842c57489SHong Zhang   bi[0] = 0;
37942c57489SHong Zhang 
38042c57489SHong Zhang   /* create and initialize a linked list */
38142c57489SHong Zhang   nlnk = pN+1;
38242c57489SHong Zhang   ierr = PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
38342c57489SHong Zhang 
38442c57489SHong Zhang   /* initial FreeSpace size is 2*(num of local nnz(C_seq)) */
385*83445d95SHong Zhang   free_space=PETSC_NULL; current_space=PETSC_NULL;
38642c57489SHong Zhang   len = 0;
38742c57489SHong Zhang   len  = ci[owners[rank+1]] - ci[owners[rank]];
38842c57489SHong Zhang   ierr = GetMoreSpace((PetscInt)(2*len+1),&free_space);CHKERRQ(ierr);
38942c57489SHong Zhang   current_space = free_space;
39042c57489SHong Zhang 
39142c57489SHong Zhang   /* determine symbolic info for each local row */
39242c57489SHong Zhang   ierr = PetscMalloc((3*merge->nrecv+1)*sizeof(PetscInt**),&buf_ri_k);CHKERRQ(ierr);
39342c57489SHong Zhang   nextrow = buf_ri_k + merge->nrecv;
39442c57489SHong Zhang   nextci  = nextrow + merge->nrecv;
39542c57489SHong Zhang   for (k=0; k<merge->nrecv; k++){
39642c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
39742c57489SHong Zhang     nrows = *buf_ri_k[k];
39842c57489SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
39942c57489SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
40042c57489SHong Zhang   }
40142c57489SHong Zhang 
40242c57489SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
40342c57489SHong Zhang   len = 0;
40442c57489SHong Zhang   for (i=0; i<pn; i++) {
40542c57489SHong Zhang     bnzi   = 0;
40642c57489SHong Zhang     /* add local non-zero cols of this proc's C_seq into lnk */
40742c57489SHong Zhang     arow   = owners[rank] + i;
40842c57489SHong Zhang     anzi   = ci[arow+1] - ci[arow];
40942c57489SHong Zhang     cji    = cj + ci[arow];
41042c57489SHong Zhang     ierr = PetscLLAdd(anzi,cji,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
41142c57489SHong Zhang     bnzi += nlnk;
41242c57489SHong Zhang     /* add received col data into lnk */
41342c57489SHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
41442c57489SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
41542c57489SHong Zhang         anzi = *(nextci[k]+1) - *nextci[k];
41642c57489SHong Zhang         cji  = buf_rj[k] + *nextci[k];
41742c57489SHong Zhang         ierr = PetscLLAdd(anzi,cji,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
41842c57489SHong Zhang         bnzi += nlnk;
41942c57489SHong Zhang         nextrow[k]++; nextci[k]++;
42042c57489SHong Zhang       }
42142c57489SHong Zhang     }
42242c57489SHong Zhang     if (len < bnzi) len = bnzi;  /* =max(bnzi) */
42342c57489SHong Zhang 
42442c57489SHong Zhang     /* if free space is not available, make more free space */
42542c57489SHong Zhang     if (current_space->local_remaining<bnzi) {
42642c57489SHong Zhang       ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
42742c57489SHong Zhang       nspacedouble++;
42842c57489SHong Zhang     }
42942c57489SHong Zhang     /* copy data into free space, then initialize lnk */
43042c57489SHong Zhang     ierr = PetscLLClean(pN,pN,bnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
43142c57489SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],bnzi,current_space->array,dnz,onz);CHKERRQ(ierr);
43242c57489SHong Zhang 
43342c57489SHong Zhang     current_space->array           += bnzi;
43442c57489SHong Zhang     current_space->local_used      += bnzi;
43542c57489SHong Zhang     current_space->local_remaining -= bnzi;
43642c57489SHong Zhang 
43742c57489SHong Zhang     bi[i+1] = bi[i] + bnzi;
43842c57489SHong Zhang   }
43942c57489SHong Zhang 
44042c57489SHong Zhang   ierr = PetscFree(buf_ri_k);CHKERRQ(ierr);
44142c57489SHong Zhang 
44242c57489SHong Zhang   ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
44342c57489SHong Zhang   ierr = MakeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
44442c57489SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
44542c57489SHong Zhang 
44642c57489SHong Zhang   /* create symbolic parallel matrix B_mpi */
44742c57489SHong Zhang   /*---------------------------------------*/
44842c57489SHong Zhang   ierr = MatCreate(comm,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE,&B_mpi);CHKERRQ(ierr);
44942c57489SHong Zhang   ierr = MatSetType(B_mpi,MATMPIAIJ);CHKERRQ(ierr);
45042c57489SHong Zhang   ierr = MatMPIAIJSetPreallocation(B_mpi,0,dnz,0,onz);CHKERRQ(ierr);
45142c57489SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
45242c57489SHong Zhang   /* ierr = MatSetOption(B_mpi,MAT_COLUMNS_SORTED);CHKERRQ(ierr); -cause delay? */
45342c57489SHong Zhang 
45442c57489SHong Zhang   /* B_mpi is not ready for use - assembly will be done by MatMerge_SeqsToMPINumeric() */
45542c57489SHong Zhang   B_mpi->assembled     = PETSC_FALSE;
45642c57489SHong Zhang   B_mpi->ops->destroy  = MatDestroy_MPIAIJ_SeqsToMPI;
45742c57489SHong Zhang   merge->bi            = bi;
45842c57489SHong Zhang   merge->bj            = bj;
45942c57489SHong Zhang   merge->ci            = ci;
46042c57489SHong Zhang   merge->cj            = cj;
46142c57489SHong Zhang   merge->buf_ri        = buf_ri;
46242c57489SHong Zhang   merge->buf_rj        = buf_rj;
46342c57489SHong Zhang   merge->C_seq         = PETSC_NULL;
46442c57489SHong Zhang 
46542c57489SHong Zhang   /* attach the supporting struct to B_mpi for reuse */
46642c57489SHong Zhang   ierr = PetscObjectContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
46742c57489SHong Zhang   ierr = PetscObjectContainerSetPointer(container,merge);CHKERRQ(ierr);
46842c57489SHong Zhang   ierr = PetscObjectCompose((PetscObject)B_mpi,"MatMergeSeqsToMPI",(PetscObject)container);CHKERRQ(ierr);
46942c57489SHong Zhang   *C = B_mpi;
47042c57489SHong Zhang 
47142c57489SHong Zhang   PetscFunctionReturn(0);
47242c57489SHong Zhang }
47342c57489SHong Zhang 
47442c57489SHong Zhang #undef __FUNCT__
47542c57489SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ"
47642c57489SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
47742c57489SHong Zhang {
47842c57489SHong Zhang   PetscErrorCode       ierr;
47942c57489SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
48042c57489SHong Zhang   Mat_MatMatMultMPI    *ptap;
48142c57489SHong Zhang   PetscObjectContainer cont_merge,cont_ptap;
48242c57489SHong Zhang 
48342c57489SHong Zhang   PetscInt       flops=0;
48442c57489SHong Zhang   Mat_MPIAIJ     *a=(Mat_MPIAIJ*)A->data;
48542c57489SHong Zhang   Mat_SeqAIJ     *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
48642c57489SHong Zhang   Mat_SeqAIJ     *p_loc,*p_oth;
487f4a743e1SHong Zhang   PetscInt       *adi=ad->i,*aoi=ao->i,*apj,cstart=a->cstart,cend=a->cend,col;
48842c57489SHong Zhang   PetscInt       *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pjj;
48942c57489SHong Zhang   PetscInt       i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
49042c57489SHong Zhang   PetscInt       *cjj;
491f4a743e1SHong Zhang   MatScalar      *ada=ad->a,*aoa=ao->a,*apa,*paj,*cseqa,*caj;
49242c57489SHong Zhang   MatScalar      *pa_loc,*pA,*pa_oth;
49342c57489SHong Zhang   PetscInt       am=A->m,cN=C->N;
494f4a743e1SHong Zhang   PetscInt       nextp,*adj=ad->j,*aoj=ao->j;
49542c57489SHong Zhang   MPI_Comm             comm=C->comm;
49642c57489SHong Zhang   PetscMPIInt          size,rank,taga,*len_s;
49742c57489SHong Zhang   PetscInt             *owners;
49842c57489SHong Zhang   PetscInt             proc;
49942c57489SHong Zhang   PetscInt             **buf_ri,**buf_rj;
50042c57489SHong Zhang   PetscInt             cseqnzi,*bj_i,*bi,*bj,cseqrow,bnzi,nextcseqj; /* bi, bj, ba: for C (mpi mat) */
50142c57489SHong Zhang   PetscInt             nrows,**buf_ri_k,**nextrow,**nextcseqi;
50242c57489SHong Zhang   MPI_Request          *s_waits,*r_waits;
50342c57489SHong Zhang   MPI_Status           *status;
50442c57489SHong Zhang   MatScalar            **abuf_r,*ba_i;
50542c57489SHong Zhang   PetscInt             *cseqi,*cseqj;
50642c57489SHong Zhang   PetscInt             *cseqj_tmp;
50742c57489SHong Zhang   MatScalar            *cseqa_tmp;
508f4a743e1SHong Zhang   PetscInt             stages[2];
50942c57489SHong Zhang   PetscInt             *apsymi,*apsymj;
51042c57489SHong Zhang 
51142c57489SHong Zhang   PetscFunctionBegin;
512f4a743e1SHong Zhang   ierr = PetscLogStageRegister(&stages[0],"NumPtAP_local");CHKERRQ(ierr);
513f4a743e1SHong Zhang   ierr = PetscLogStageRegister(&stages[1],"NumPtAP_Comm");CHKERRQ(ierr);
51442c57489SHong Zhang 
51542c57489SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
51642c57489SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
51742c57489SHong Zhang 
51842c57489SHong Zhang   ierr = PetscObjectQuery((PetscObject)C,"MatMergeSeqsToMPI",(PetscObject *)&cont_merge);CHKERRQ(ierr);
51942c57489SHong Zhang   if (cont_merge) {
52042c57489SHong Zhang     ierr  = PetscObjectContainerGetPointer(cont_merge,(void **)&merge);CHKERRQ(ierr);
52142c57489SHong Zhang   } else {
52242c57489SHong Zhang     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix C does not posses an object container");
52342c57489SHong Zhang   }
52442c57489SHong Zhang   ierr = PetscObjectQuery((PetscObject)P,"MatPtAPMPI",(PetscObject *)&cont_ptap);CHKERRQ(ierr);
52542c57489SHong Zhang   if (cont_ptap) {
52642c57489SHong Zhang     ierr  = PetscObjectContainerGetPointer(cont_ptap,(void **)&ptap);CHKERRQ(ierr);
52742c57489SHong Zhang   } else {
52842c57489SHong Zhang     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix P does not posses an object container");
52942c57489SHong Zhang   }
53042c57489SHong Zhang 
53142c57489SHong Zhang   /* get data from symbolic products */
53242c57489SHong Zhang   p_loc=(Mat_SeqAIJ*)(ptap->B_loc)->data;
53342c57489SHong Zhang   p_oth=(Mat_SeqAIJ*)(ptap->B_oth)->data;
53442c57489SHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a; pA=pa_loc;
53542c57489SHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
53642c57489SHong Zhang 
53742c57489SHong Zhang   cseqi = merge->ci; cseqj=merge->cj;
53842c57489SHong Zhang   ierr  = PetscMalloc((cseqi[cN]+1)*sizeof(MatScalar),&cseqa);CHKERRQ(ierr);
53942c57489SHong Zhang   ierr  = PetscMemzero(cseqa,cseqi[cN]*sizeof(MatScalar));CHKERRQ(ierr);
54042c57489SHong Zhang 
541f4a743e1SHong Zhang   /* get data from symbolic A*P */
542f4a743e1SHong Zhang   ierr = PetscMalloc((ptap->abnz_max+1)*sizeof(MatScalar),&apa);CHKERRQ(ierr);
543f4a743e1SHong Zhang   ierr = PetscMemzero(apa,ptap->abnz_max*sizeof(MatScalar));CHKERRQ(ierr);
54442c57489SHong Zhang 
545f4a743e1SHong Zhang   /* compute numeric C_seq=P_loc^T*A_loc*P */
546f4a743e1SHong Zhang   ierr = PetscLogStagePush(stages[0]);CHKERRQ(ierr);
547f4a743e1SHong Zhang   apsymi = ptap->abi; apsymj = ptap->abj;
54842c57489SHong Zhang   for (i=0;i<am;i++) {
549f4a743e1SHong Zhang     /* form i-th sparse row of A*P */
550f4a743e1SHong Zhang     apnzj = apsymi[i+1] - apsymi[i];
551f4a743e1SHong Zhang     apj   = apsymj + apsymi[i];
55242c57489SHong Zhang     /* diagonal portion of A */
55342c57489SHong Zhang     anzi  = adi[i+1] - adi[i];
55442c57489SHong Zhang     for (j=0;j<anzi;j++) {
55542c57489SHong Zhang       prow = *adj;
55642c57489SHong Zhang       adj++;
55742c57489SHong Zhang       pnzj = pi_loc[prow+1] - pi_loc[prow];
55842c57489SHong Zhang       pjj  = pj_loc + pi_loc[prow];
55942c57489SHong Zhang       paj  = pa_loc + pi_loc[prow];
560f4a743e1SHong Zhang       nextp = 0;
561f4a743e1SHong Zhang       for (k=0; nextp<pnzj; k++) {
562f4a743e1SHong Zhang         if (apj[k] == pjj[nextp]) { /* col of AP == col of P */
563f4a743e1SHong Zhang           apa[k] += (*ada)*paj[nextp++];
56442c57489SHong Zhang         }
56542c57489SHong Zhang       }
56642c57489SHong Zhang       flops += 2*pnzj;
56742c57489SHong Zhang       ada++;
56842c57489SHong Zhang     }
56942c57489SHong Zhang     /* off-diagonal portion of A */
57042c57489SHong Zhang     anzi  = aoi[i+1] - aoi[i];
57142c57489SHong Zhang     for (j=0;j<anzi;j++) {
57242c57489SHong Zhang       col = a->garray[*aoj];
57342c57489SHong Zhang       if (col < cstart){
57442c57489SHong Zhang         prow = *aoj;
57542c57489SHong Zhang       } else if (col >= cend){
57642c57489SHong Zhang         prow = *aoj;
57742c57489SHong Zhang       } else {
57842c57489SHong Zhang         SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Off-diagonal portion of A has wrong column map");
57942c57489SHong Zhang       }
58042c57489SHong Zhang       aoj++;
58142c57489SHong Zhang       pnzj = pi_oth[prow+1] - pi_oth[prow];
58242c57489SHong Zhang       pjj  = pj_oth + pi_oth[prow];
58342c57489SHong Zhang       paj  = pa_oth + pi_oth[prow];
584f4a743e1SHong Zhang       nextp = 0;
585f4a743e1SHong Zhang       for (k=0; nextp<pnzj; k++) {
586f4a743e1SHong Zhang         if (apj[k] == pjj[nextp]) { /* col of AP == col of P */
587f4a743e1SHong Zhang           apa[k] += (*aoa)*paj[nextp++];
58842c57489SHong Zhang         }
58942c57489SHong Zhang       }
59042c57489SHong Zhang       flops += 2*pnzj;
59142c57489SHong Zhang       aoa++;
59242c57489SHong Zhang     }
59342c57489SHong Zhang 
59442c57489SHong Zhang     /* Compute P_loc[i,:]^T*AP[i,:] using outer product */
59542c57489SHong Zhang     pnzi = pi_loc[i+1] - pi_loc[i];
59642c57489SHong Zhang     for (j=0;j<pnzi;j++) {
59742c57489SHong Zhang       nextap = 0;
59842c57489SHong Zhang       crow   = *pJ++;
59942c57489SHong Zhang       cjj    = cseqj + cseqi[crow];
60042c57489SHong Zhang       caj    = cseqa + cseqi[crow];
60142c57489SHong Zhang       /* Perform sparse axpy operation.  Note cjj includes apj. */
60242c57489SHong Zhang       for (k=0;nextap<apnzj;k++) {
60342c57489SHong Zhang         if (cjj[k]==apj[nextap]) caj[k] += (*pA)*apa[nextap++];
60442c57489SHong Zhang       }
60542c57489SHong Zhang       flops += 2*apnzj;
60642c57489SHong Zhang       pA++;
60742c57489SHong Zhang     }
60842c57489SHong Zhang     /* zero the current row info for A*P */
60942c57489SHong Zhang     ierr = PetscMemzero(apa,apnzj*sizeof(MatScalar));CHKERRQ(ierr);
61042c57489SHong Zhang   }
61142c57489SHong Zhang 
61242c57489SHong Zhang   ierr = PetscFree(apa);CHKERRQ(ierr);
61342c57489SHong Zhang   ierr = PetscLogStagePop();
61442c57489SHong Zhang 
61542c57489SHong Zhang   bi     = merge->bi;
61642c57489SHong Zhang   bj     = merge->bj;
61742c57489SHong Zhang   buf_ri = merge->buf_ri;
61842c57489SHong Zhang   buf_rj = merge->buf_rj;
61942c57489SHong Zhang 
62042c57489SHong Zhang   ierr   = PetscMalloc(size*sizeof(MPI_Status),&status);CHKERRQ(ierr);
62142c57489SHong Zhang   ierr   = PetscMapGetGlobalRange(merge->rowmap,&owners);CHKERRQ(ierr);
62242c57489SHong Zhang 
62342c57489SHong Zhang   /* send and recv matrix values */
62442c57489SHong Zhang   /*-----------------------------*/
625f4a743e1SHong Zhang   ierr = PetscLogStagePush(stages[1]);CHKERRQ(ierr);
62642c57489SHong Zhang   len_s  = merge->len_s;
62742c57489SHong Zhang   ierr = PetscObjectGetNewTag((PetscObject)merge->rowmap,&taga);CHKERRQ(ierr);
62842c57489SHong Zhang   ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
62942c57489SHong Zhang 
63042c57489SHong Zhang   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&s_waits);CHKERRQ(ierr);
63142c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++){
63242c57489SHong Zhang     if (!len_s[proc]) continue;
63342c57489SHong Zhang     i = owners[proc];
63442c57489SHong Zhang     ierr = MPI_Isend(cseqa+cseqi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
63542c57489SHong Zhang     k++;
63642c57489SHong Zhang   }
63742c57489SHong Zhang 
63842c57489SHong Zhang   ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);
63942c57489SHong Zhang   ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);
64042c57489SHong Zhang   ierr = PetscFree(status);CHKERRQ(ierr);
64142c57489SHong Zhang 
64242c57489SHong Zhang   ierr = PetscFree(s_waits);CHKERRQ(ierr);
64342c57489SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
64442c57489SHong Zhang 
64542c57489SHong Zhang   /* insert mat values of mpimat */
64642c57489SHong Zhang   /*----------------------------*/
64742c57489SHong Zhang   ierr = PetscMalloc(cN*sizeof(MatScalar),&ba_i);CHKERRQ(ierr);
64842c57489SHong Zhang   ierr = PetscMalloc((3*merge->nrecv+1)*sizeof(PetscInt**),&buf_ri_k);CHKERRQ(ierr);
64942c57489SHong Zhang   nextrow = buf_ri_k + merge->nrecv;
65042c57489SHong Zhang   nextcseqi  = nextrow + merge->nrecv;
65142c57489SHong Zhang 
65242c57489SHong Zhang   for (k=0; k<merge->nrecv; k++){
65342c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
65442c57489SHong Zhang     nrows = *(buf_ri_k[k]);
65542c57489SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
65642c57489SHong Zhang     nextcseqi[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
65742c57489SHong Zhang   }
65842c57489SHong Zhang 
65942c57489SHong Zhang   /* set values of ba */
66042c57489SHong Zhang   for (i=0; i<C->m; i++) {
66142c57489SHong Zhang     cseqrow = owners[rank] + i; /* global row index of C_seq */
66242c57489SHong Zhang     bj_i = bj+bi[i];  /* col indices of the i-th row of C */
66342c57489SHong Zhang     bnzi = bi[i+1] - bi[i];
66442c57489SHong Zhang     ierr = PetscMemzero(ba_i,bnzi*sizeof(MatScalar));CHKERRQ(ierr);
66542c57489SHong Zhang 
66642c57489SHong Zhang     /* add local non-zero vals of this proc's C_seq into ba */
66742c57489SHong Zhang     cseqnzi = cseqi[cseqrow+1] - cseqi[cseqrow];
66842c57489SHong Zhang     cseqj_tmp = cseqj + cseqi[cseqrow];
66942c57489SHong Zhang     cseqa_tmp = cseqa + cseqi[cseqrow];
67042c57489SHong Zhang     nextcseqj = 0;
67142c57489SHong Zhang     for (j=0; nextcseqj<cseqnzi; j++){
67242c57489SHong Zhang       if (*(bj_i + j) == cseqj_tmp[nextcseqj]){ /* bcol == cseqcol */
67342c57489SHong Zhang         ba_i[j] += cseqa_tmp[nextcseqj++];
67442c57489SHong Zhang       }
67542c57489SHong Zhang     }
67642c57489SHong Zhang 
67742c57489SHong Zhang     /* add received vals into ba */
67842c57489SHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
67942c57489SHong Zhang       /* i-th row */
68042c57489SHong Zhang       if (i == *nextrow[k]) {
68142c57489SHong Zhang         cseqnzi   = *(nextcseqi[k]+1) - *nextcseqi[k];
68242c57489SHong Zhang         cseqj_tmp = buf_rj[k] + *(nextcseqi[k]);
68342c57489SHong Zhang         cseqa_tmp = abuf_r[k] + *(nextcseqi[k]);
68442c57489SHong Zhang         nextcseqj = 0;
68542c57489SHong Zhang         for (j=0; nextcseqj<cseqnzi; j++){
68642c57489SHong Zhang           if (*(bj_i + j) == cseqj_tmp[nextcseqj]){ /* bcol == cseqcol */
68742c57489SHong Zhang             ba_i[j] += cseqa_tmp[nextcseqj++];
68842c57489SHong Zhang           }
68942c57489SHong Zhang         }
69042c57489SHong Zhang         nextrow[k]++; nextcseqi[k]++;
69142c57489SHong Zhang       }
69242c57489SHong Zhang     }
69342c57489SHong Zhang     ierr = MatSetValues(C,1,&cseqrow,bnzi,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
69442c57489SHong Zhang     flops += 2*bnzi;
69542c57489SHong Zhang   }
69642c57489SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
69742c57489SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
69842c57489SHong Zhang   ierr = PetscLogStagePop();CHKERRQ(ierr);
69942c57489SHong Zhang 
70042c57489SHong Zhang   ierr = PetscFree(cseqa);CHKERRQ(ierr);
70142c57489SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
70242c57489SHong Zhang   ierr = PetscFree(ba_i);CHKERRQ(ierr);
70342c57489SHong Zhang   ierr = PetscFree(buf_ri_k);CHKERRQ(ierr);
70442c57489SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
70542c57489SHong Zhang 
70642c57489SHong Zhang   PetscFunctionReturn(0);
70742c57489SHong Zhang }
708