xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision 29825010b655ce3f8e1ce5f7b74757bc50ec8a87)
182412ba7SHong Zhang 
242c57489SHong Zhang /*
342c57489SHong Zhang   Defines projective product routines where A is a MPIAIJ matrix
442c57489SHong Zhang           C = P^T * A * P
542c57489SHong Zhang */
642c57489SHong Zhang 
7c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>   /*I "petscmat.h" I*/
8c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
9c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h>
10c6db04a5SJed Brown #include <petscbt.h>
1142c57489SHong Zhang 
12fe615159SHong Zhang /*
13fe615159SHong Zhang #define DEBUG_MATPTAP
14fe615159SHong Zhang  */
15fe615159SHong Zhang 
1609573ac7SBarry Smith extern PetscErrorCode MatDestroy_MPIAIJ(Mat);
17cc6cb787SHong Zhang #undef __FUNCT__
18f8487c73SHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_PtAP"
19f8487c73SHong Zhang PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A)
20cc6cb787SHong Zhang {
21cc6cb787SHong Zhang   PetscErrorCode       ierr;
22f8487c73SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data;
23f8487c73SHong Zhang   Mat_PtAPMPI          *ptap=a->ptap;
24cc6cb787SHong Zhang 
25cc6cb787SHong Zhang   PetscFunctionBegin;
26f8487c73SHong Zhang   if (ptap){
27f8487c73SHong Zhang     ierr = PetscFree2(ptap->startsj,ptap->startsj_r);CHKERRQ(ierr);
28f8487c73SHong Zhang     ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
29a1a4e84aSHong Zhang     ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
30a1a4e84aSHong Zhang     ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
31a1a4e84aSHong Zhang     ierr = PetscFree(ptap->api);CHKERRQ(ierr);
32a1a4e84aSHong Zhang     ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
33f8487c73SHong Zhang   }
34f8487c73SHong Zhang   if (ptap->merge) {
35f8487c73SHong Zhang     Mat_Merge_SeqsToMPI  *merge=ptap->merge;
36cc6cb787SHong Zhang     ierr = PetscFree(merge->id_r);CHKERRQ(ierr);
37cc6cb787SHong Zhang     ierr = PetscFree(merge->len_s);CHKERRQ(ierr);
38cc6cb787SHong Zhang     ierr = PetscFree(merge->len_r);CHKERRQ(ierr);
39cc6cb787SHong Zhang     ierr = PetscFree(merge->bi);CHKERRQ(ierr);
40cc6cb787SHong Zhang     ierr = PetscFree(merge->bj);CHKERRQ(ierr);
41c05d87d6SBarry Smith     ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr);
42cc6cb787SHong Zhang     ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr);
43c05d87d6SBarry Smith     ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr);
44cc6cb787SHong Zhang     ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr);
4505b42c5fSBarry Smith     ierr = PetscFree(merge->coi);CHKERRQ(ierr);
4605b42c5fSBarry Smith     ierr = PetscFree(merge->coj);CHKERRQ(ierr);
4705b42c5fSBarry Smith     ierr = PetscFree(merge->owners_co);CHKERRQ(ierr);
486bf464f9SBarry Smith     ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr);
49dce485f0SBarry Smith     ierr = merge->destroy(A);CHKERRQ(ierr);
50f8487c73SHong Zhang     ierr = PetscFree(ptap->merge);CHKERRQ(ierr);
51bf0cc555SLisandro Dalcin   }
52f8487c73SHong Zhang   ierr = PetscFree(ptap);CHKERRQ(ierr);
53cc6cb787SHong Zhang   PetscFunctionReturn(0);
54cc6cb787SHong Zhang }
55cc6cb787SHong Zhang 
56cc6cb787SHong Zhang #undef __FUNCT__
57cc6cb787SHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatPtAP"
58f0c0a3a6SBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M)
59f0c0a3a6SBarry Smith {
60cc6cb787SHong Zhang   PetscErrorCode       ierr;
61cc6cb787SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
62776b82aeSLisandro Dalcin   PetscContainer       container;
63cc6cb787SHong Zhang 
64cc6cb787SHong Zhang   PetscFunctionBegin;
65cc6cb787SHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"MatMergeSeqsToMPI",(PetscObject *)&container);CHKERRQ(ierr);
66cc6cb787SHong Zhang   if (container) {
67776b82aeSLisandro Dalcin     ierr  = PetscContainerGetPointer(container,(void **)&merge);CHKERRQ(ierr);
68e7e72b3dSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
69dce485f0SBarry Smith   ierr = (*merge->duplicate)(A,op,M);CHKERRQ(ierr);
70dce485f0SBarry Smith   (*M)->ops->destroy   = merge->destroy;   /* =MatDestroy_MPIAIJ, *M doesn't duplicate A's container! */
71dce485f0SBarry Smith   (*M)->ops->duplicate = merge->duplicate; /* =MatDuplicate_ MPIAIJ */
72cc6cb787SHong Zhang   PetscFunctionReturn(0);
73cc6cb787SHong Zhang }
74cc6cb787SHong Zhang 
7542c57489SHong Zhang #undef __FUNCT__
767a7894deSKris Buschelman #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ"
777a7894deSKris Buschelman PetscErrorCode MatPtAPSymbolic_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
7842c57489SHong Zhang {
7942c57489SHong Zhang   PetscErrorCode ierr;
8042c57489SHong Zhang 
8142c57489SHong Zhang   PetscFunctionBegin;
8265e19b50SBarry 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);
837a7894deSKris Buschelman   ierr = (*P->ops->ptapsymbolic_mpiaij)(A,P,fill,C);CHKERRQ(ierr);
847a7894deSKris Buschelman   PetscFunctionReturn(0);
857a7894deSKris Buschelman }
867a7894deSKris Buschelman 
877a7894deSKris Buschelman #undef __FUNCT__
887a7894deSKris Buschelman #define __FUNCT__ "MatPtAPNumeric_MPIAIJ"
897a7894deSKris Buschelman PetscErrorCode MatPtAPNumeric_MPIAIJ(Mat A,Mat P,Mat C)
907a7894deSKris Buschelman {
917a7894deSKris Buschelman   PetscErrorCode ierr;
927a7894deSKris Buschelman 
937a7894deSKris Buschelman   PetscFunctionBegin;
9465e19b50SBarry 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);
957a7894deSKris Buschelman   ierr = (*P->ops->ptapnumeric_mpiaij)(A,P,C);CHKERRQ(ierr);
9642c57489SHong Zhang   PetscFunctionReturn(0);
9742c57489SHong Zhang }
9842c57489SHong Zhang 
9942c57489SHong Zhang #undef __FUNCT__
10042c57489SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ"
10142c57489SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
10242c57489SHong Zhang {
10342c57489SHong Zhang   PetscErrorCode       ierr;
10477584fe6SHong Zhang   Mat                  Cmpi;
105f8487c73SHong Zhang   Mat_PtAPMPI          *ptap;
106a1a86e44SBarry Smith   PetscFreeSpaceList   free_space=PETSC_NULL,current_space=PETSC_NULL;
107f8487c73SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
10842c57489SHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
10942c57489SHong Zhang   Mat_SeqAIJ           *p_loc,*p_oth;
110ded4f561SHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ;
11177584fe6SHong Zhang   PetscInt             *adi=ad->i,*aj,*aoi=ao->i,nnz;
11213f74950SBarry Smith   PetscInt             nlnk,*lnk,*owners_co,*coi,*coj,i,k,pnz,row;
113d16ca5f0SHong Zhang   PetscInt             am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n;
11442c57489SHong Zhang   PetscBT              lnkbt;
1157adad957SLisandro Dalcin   MPI_Comm             comm=((PetscObject)A)->comm;
116529bc5dcSHong Zhang   PetscMPIInt          size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
11742c57489SHong Zhang   PetscInt             **buf_rj,**buf_ri,**buf_ri_k;
11842c57489SHong Zhang   PetscInt             len,proc,*dnz,*onz,*owners;
119ded4f561SHong Zhang   PetscInt             nzi,*bi,*bj;
12042c57489SHong Zhang   PetscInt             nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
121ded4f561SHong Zhang   MPI_Request          *swaits,*rwaits;
122ded4f561SHong Zhang   MPI_Status           *sstatus,rstatus;
12342c57489SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
12477584fe6SHong Zhang   PetscInt             *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0;
125d16ca5f0SHong Zhang   PetscReal            afill=1.0,afill_tmp;
126dadf0e6bSHong 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;
12708cb4532SHong Zhang   PetscInt             rstart = P->cmap->rstart,rmax;
12808cb4532SHong Zhang   PetscScalar          *vals;
12942c57489SHong Zhang 
13042c57489SHong Zhang   PetscFunctionBegin;
13177584fe6SHong Zhang   ierr = PetscGetTime(&t00);CHKERRQ(ierr);
13283445d95SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
13383445d95SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
13483445d95SHong Zhang 
13577584fe6SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
136f8487c73SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
137f8487c73SHong Zhang   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
138f8487c73SHong Zhang   ptap->reuse    = MAT_INITIAL_MATRIX;
1396c6e00e0SHong Zhang 
14077584fe6SHong Zhang 
1416c6e00e0SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
142a1a4e84aSHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
143fe615159SHong Zhang #if defined(DEBUG_MATPTAP)
144fe615159SHong Zhang    if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_oth is done...\n",rank);
145fe615159SHong Zhang #endif
146fe615159SHong Zhang 
1476c6e00e0SHong Zhang   /* get P_loc by taking all local rows of P */
148a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
149fe615159SHong Zhang #if defined(DEBUG_MATPTAP)
150fe615159SHong Zhang   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
151fe615159SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_loc is done...\n",rank);
152fe615159SHong Zhang #endif
1536c6e00e0SHong Zhang 
154a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
155a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
1566c6e00e0SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
1576c6e00e0SHong Zhang   pi_oth = p_oth->i; pj_oth = p_oth->j;
15842c57489SHong Zhang 
15977584fe6SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
16077584fe6SHong Zhang   time_matupdate += tf-t0;
16177584fe6SHong Zhang 
162fd0ff01cSHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
163fd0ff01cSHong Zhang   /*-------------------------------------------------------------------*/
16477584fe6SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
16577584fe6SHong Zhang   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&api);CHKERRQ(ierr);
16682412ba7SHong Zhang   api[0] = 0;
16783445d95SHong Zhang 
16883445d95SHong Zhang   /* create and initialize a linked list */
16983445d95SHong Zhang   nlnk = pN+1;
170fe615159SHong Zhang #if defined(DEBUG_MATPTAP)
171fe615159SHong Zhang   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
172fe615159SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call PetscLLCreate(), pN %D, P_rmax %D %D; Annz %D\n",rank,pN,p_loc->rmax,p_oth->rmax,adi[am]+aoi[am]);
173fe615159SHong Zhang #endif
17483445d95SHong Zhang   ierr = PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
175fe615159SHong Zhang #if defined(DEBUG_MATPTAP)
176fe615159SHong Zhang   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
177fe615159SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call PetscLLCreate() is done\n",rank);
178fe615159SHong Zhang #endif
17983445d95SHong Zhang 
180d16ca5f0SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
181d16ca5f0SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
182f4a743e1SHong Zhang   current_space = free_space;
183f4a743e1SHong Zhang 
184f4a743e1SHong Zhang   for (i=0; i<am; i++) {
18582412ba7SHong Zhang     apnz = 0;
186f4a743e1SHong Zhang     /* diagonal portion of A */
187ded4f561SHong Zhang     nzi = adi[i+1] - adi[i];
18877584fe6SHong Zhang     aj  = ad->j + adi[i];
189ded4f561SHong Zhang     for (j=0; j<nzi; j++){
19077584fe6SHong Zhang       row  = aj[j];
19182412ba7SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
192ded4f561SHong Zhang       Jptr = pj_loc + pi_loc[row];
19383445d95SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
194dadf0e6bSHong Zhang       ierr = PetscLLAddSorted(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
19582412ba7SHong Zhang       apnz += nlnk;
196f4a743e1SHong Zhang     }
197f4a743e1SHong Zhang     /* off-diagonal portion of A */
198ded4f561SHong Zhang     nzi = aoi[i+1] - aoi[i];
19977584fe6SHong Zhang     aj  = ao->j + aoi[i];
200ded4f561SHong Zhang     for (j=0; j<nzi; j++){
20177584fe6SHong Zhang       row = aj[j];
20282412ba7SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
203ded4f561SHong Zhang       Jptr  = pj_oth + pi_oth[row];
204dadf0e6bSHong Zhang       ierr = PetscLLAddSorted(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
20582412ba7SHong Zhang       apnz += nlnk;
206f4a743e1SHong Zhang     }
207f4a743e1SHong Zhang 
20882412ba7SHong Zhang     api[i+1] = api[i] + apnz;
20977584fe6SHong Zhang     if (ap_rmax < apnz) ap_rmax = apnz;
210f4a743e1SHong Zhang 
21183445d95SHong Zhang     /* if free space is not available, double the total space in the list */
21282412ba7SHong Zhang     if (current_space->local_remaining<apnz) {
2132ba1acd5SHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
214f2b054eeSHong Zhang       nspacedouble++;
215f4a743e1SHong Zhang     }
216f4a743e1SHong Zhang 
217f4a743e1SHong Zhang     /* Copy data into free space, then initialize lnk */
21882412ba7SHong Zhang     ierr = PetscLLClean(pN,pN,apnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
21982412ba7SHong Zhang     current_space->array           += apnz;
22082412ba7SHong Zhang     current_space->local_used      += apnz;
22182412ba7SHong Zhang     current_space->local_remaining -= apnz;
222f4a743e1SHong Zhang   }
22382412ba7SHong Zhang   /* Allocate space for apj, initialize apj, and */
224f4a743e1SHong Zhang   /* destroy list of free space and other temporary array(s) */
22577584fe6SHong Zhang   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&apj);CHKERRQ(ierr);
22677584fe6SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
227d16ca5f0SHong Zhang   afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]);
228d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
229f4a743e1SHong Zhang 
23077584fe6SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
23177584fe6SHong Zhang   time_AP += tf-t0;
232fe615159SHong Zhang #if defined(DEBUG_MATPTAP)
233fe615159SHong Zhang   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
234fe615159SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] AP is done\n",rank);
235fe615159SHong Zhang #endif
23677584fe6SHong Zhang 
237ded4f561SHong Zhang   /* determine symbolic Co=(p->B)^T*AP - send to others */
238ded4f561SHong Zhang   /*----------------------------------------------------*/
23977584fe6SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
240de0260b3SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
24142c57489SHong Zhang 
242ded4f561SHong Zhang   /* then, compute symbolic Co = (p->B)^T*AP */
243d0f46423SBarry Smith   pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors
24483445d95SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
245de0260b3SHong Zhang   ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr);
246de0260b3SHong Zhang   coi[0] = 0;
24742c57489SHong Zhang 
248d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */
249d16ca5f0SHong Zhang   nnz           = fill*(poti[pon] + api[am]);
250a1a86e44SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);
25142c57489SHong Zhang   current_space = free_space;
25242c57489SHong Zhang 
253de0260b3SHong Zhang   for (i=0; i<pon; i++) {
254ded4f561SHong Zhang     nnz = 0;
25582412ba7SHong Zhang     pnz = poti[i+1] - poti[i];
25677584fe6SHong Zhang     ptJ = potj + poti[i];
25777584fe6SHong Zhang     for (j=0; j<pnz; j++){
25877584fe6SHong Zhang       row  = ptJ[j]; /* row of AP == col of Pot */
25982412ba7SHong Zhang       apnz = api[row+1] - api[row];
260ded4f561SHong Zhang       Jptr = apj + api[row];
26183445d95SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
262dadf0e6bSHong Zhang       ierr = PetscLLAddSorted(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
263ded4f561SHong Zhang       nnz += nlnk;
26442c57489SHong Zhang     }
26542c57489SHong Zhang 
26683445d95SHong Zhang     /* If free space is not available, double the total space in the list */
267ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
2684238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
269d16ca5f0SHong Zhang       nspacedouble++;
27042c57489SHong Zhang     }
27142c57489SHong Zhang 
27242c57489SHong Zhang     /* Copy data into free space, and zero out denserows */
273ded4f561SHong Zhang     ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
274ded4f561SHong Zhang     current_space->array           += nnz;
275ded4f561SHong Zhang     current_space->local_used      += nnz;
276ded4f561SHong Zhang     current_space->local_remaining -= nnz;
277ded4f561SHong Zhang     coi[i+1] = coi[i] + nnz;
27842c57489SHong Zhang   }
279de0260b3SHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr);
280a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
281d16ca5f0SHong Zhang   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]);
282d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
283de0260b3SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
28442c57489SHong Zhang 
285fe615159SHong Zhang #if defined(DEBUG_MATPTAP)
286fe615159SHong Zhang   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
287fe615159SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] local Co is done\n",rank);
288fe615159SHong Zhang #endif
289fe615159SHong Zhang 
290ded4f561SHong Zhang   /* send j-array (coj) of Co to other processors */
291ded4f561SHong Zhang   /*----------------------------------------------*/
29242c57489SHong Zhang   /* determine row ownership */
29383445d95SHong Zhang   ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr);
29426283091SBarry Smith   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
2957a2fc3feSBarry Smith   merge->rowmap->n = pn;
2967a2fc3feSBarry Smith   merge->rowmap->bs = 1;
29726283091SBarry Smith   ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
2987a2fc3feSBarry Smith   owners = merge->rowmap->range;
29942c57489SHong Zhang 
30042c57489SHong Zhang   /* determine the number of messages to send, their lengths */
30142c57489SHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr);
30283445d95SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
30342c57489SHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr);
30442c57489SHong Zhang   len_s = merge->len_s;
30542c57489SHong Zhang   merge->nsend = 0;
30683445d95SHong Zhang 
307de0260b3SHong Zhang   ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr);
308de0260b3SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
309de0260b3SHong Zhang 
31083445d95SHong Zhang   proc = 0;
311de0260b3SHong Zhang   for (i=0; i<pon; i++){
312de0260b3SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
313de0260b3SHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
314de0260b3SHong Zhang     len_s[proc] += coi[i+1] - coi[i];
31583445d95SHong Zhang   }
316de0260b3SHong Zhang 
317de0260b3SHong Zhang   len   = 0;  /* max length of buf_si[] */
318de0260b3SHong Zhang   owners_co[0] = 0;
31942c57489SHong Zhang   for (proc=0; proc<size; proc++){
320de0260b3SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
32183445d95SHong Zhang     if (len_si[proc]){
32242c57489SHong Zhang       merge->nsend++;
32383445d95SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
32442c57489SHong Zhang       len += len_si[proc];
32542c57489SHong Zhang     }
32642c57489SHong Zhang   }
32742c57489SHong Zhang 
328ded4f561SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
32942c57489SHong Zhang   ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
33042c57489SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
33142c57489SHong Zhang 
332ded4f561SHong Zhang   /* post the Irecv and Isend of coj */
333529bc5dcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
334529bc5dcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
335ded4f561SHong Zhang   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr);
33642c57489SHong Zhang   for (proc=0, k=0; proc<size; proc++){
33742c57489SHong Zhang     if (!len_s[proc]) continue;
338de0260b3SHong Zhang     i = owners_co[proc];
339529bc5dcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
34042c57489SHong Zhang     k++;
34142c57489SHong Zhang   }
34242c57489SHong Zhang 
343ded4f561SHong Zhang   /* receives and sends of coj are complete */
344ded4f561SHong Zhang   ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr);
34577584fe6SHong Zhang   for (i=0; i<merge->nrecv; i++){
346ded4f561SHong Zhang     ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr);
347ded4f561SHong Zhang   }
348ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
3490c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
35082412ba7SHong Zhang 
351ded4f561SHong Zhang   /* send and recv coi */
352ded4f561SHong Zhang   /*-------------------*/
353529bc5dcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
354529bc5dcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
35542c57489SHong Zhang   ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr);
35642c57489SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
35742c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++){
35842c57489SHong Zhang     if (!len_s[proc]) continue;
35942c57489SHong Zhang     /* form outgoing message for i-structure:
36042c57489SHong Zhang          buf_si[0]:                 nrows to be sent
36142c57489SHong Zhang                [1:nrows]:           row index (global)
36242c57489SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
36342c57489SHong Zhang     */
36442c57489SHong Zhang     /*-------------------------------------------*/
36542c57489SHong Zhang     nrows = len_si[proc]/2 - 1;
36642c57489SHong Zhang     buf_si_i    = buf_si + nrows+1;
36742c57489SHong Zhang     buf_si[0]   = nrows;
36842c57489SHong Zhang     buf_si_i[0] = 0;
36942c57489SHong Zhang     nrows = 0;
370de0260b3SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++){
371ded4f561SHong Zhang       nzi = coi[i+1] - coi[i];
372ded4f561SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
373de0260b3SHong Zhang       buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */
37442c57489SHong Zhang       nrows++;
37542c57489SHong Zhang     }
376529bc5dcSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
37742c57489SHong Zhang     k++;
37842c57489SHong Zhang     buf_si += len_si[proc];
37942c57489SHong Zhang   }
380ded4f561SHong Zhang   i = merge->nrecv;
381ded4f561SHong Zhang   while (i--) {
382ded4f561SHong Zhang     ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr);
383ded4f561SHong Zhang   }
384ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
3850c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
386f2b054eeSHong Zhang   /*
387ae15b995SBarry Smith   ierr = PetscInfo2(A,"nsend: %d, nrecv: %d\n",merge->nsend,merge->nrecv);CHKERRQ(ierr);
38842c57489SHong Zhang   for (i=0; i<merge->nrecv; i++){
389ae15b995SBarry 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);
39042c57489SHong Zhang   }
391f2b054eeSHong Zhang   */
39242c57489SHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
39342c57489SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
394ded4f561SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
395ded4f561SHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
39642c57489SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
39742c57489SHong Zhang 
39877584fe6SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
39977584fe6SHong Zhang   time_Co += tf-t0;
40077584fe6SHong Zhang 
401fe615159SHong Zhang #if defined(DEBUG_MATPTAP)
402fe615159SHong Zhang   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
403fe615159SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Co is done\n",rank);
404fe615159SHong Zhang #endif
405fe615159SHong Zhang 
406de0260b3SHong Zhang   /* compute the local portion of C (mpi mat) */
407de0260b3SHong Zhang   /*------------------------------------------*/
40877584fe6SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
409ded4f561SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
410ded4f561SHong Zhang 
41142c57489SHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
41242c57489SHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
41342c57489SHong Zhang   bi[0] = 0;
41442c57489SHong Zhang 
415d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
416d16ca5f0SHong Zhang   nnz           = fill*(pi_loc[pm] + api[am]);
417a1a86e44SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);
41842c57489SHong Zhang   current_space = free_space;
41942c57489SHong Zhang 
420687f1162SBarry Smith   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
42142c57489SHong Zhang   for (k=0; k<merge->nrecv; k++){
42242c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
42342c57489SHong Zhang     nrows       = *buf_ri_k[k];
42442c57489SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
42542c57489SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
42642c57489SHong Zhang   }
42742c57489SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
42808cb4532SHong Zhang   rmax = 0;
42942c57489SHong Zhang   for (i=0; i<pn; i++) {
430ded4f561SHong Zhang     /* add pdt[i,:]*AP into lnk */
431ded4f561SHong Zhang     nnz = 0;
432ded4f561SHong Zhang     pnz = pdti[i+1] - pdti[i];
43377584fe6SHong Zhang     ptJ = pdtj + pdti[i];
43477584fe6SHong Zhang     for (j=0; j<pnz; j++){
43577584fe6SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
436ded4f561SHong Zhang       apnz = api[row+1] - api[row];
437ded4f561SHong Zhang       Jptr = apj + api[row];
438ded4f561SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
439dadf0e6bSHong Zhang       ierr = PetscLLAddSorted(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
440ded4f561SHong Zhang       nnz += nlnk;
441ded4f561SHong Zhang     }
44242c57489SHong Zhang     /* add received col data into lnk */
44342c57489SHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
44442c57489SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
445ded4f561SHong Zhang         nzi = *(nextci[k]+1) - *nextci[k];
446ded4f561SHong Zhang         Jptr  = buf_rj[k] + *nextci[k];
447dadf0e6bSHong Zhang         ierr = PetscLLAddSorted(nzi,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
448ded4f561SHong Zhang         nnz += nlnk;
44942c57489SHong Zhang         nextrow[k]++; nextci[k]++;
45042c57489SHong Zhang       }
45142c57489SHong Zhang     }
45242c57489SHong Zhang 
45342c57489SHong Zhang     /* if free space is not available, make more free space */
454ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
4554238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
456d16ca5f0SHong Zhang       nspacedouble++;
45742c57489SHong Zhang     }
45842c57489SHong Zhang     /* copy data into free space, then initialize lnk */
459ded4f561SHong Zhang     ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
460ded4f561SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
461ded4f561SHong Zhang     current_space->array           += nnz;
462ded4f561SHong Zhang     current_space->local_used      += nnz;
463ded4f561SHong Zhang     current_space->local_remaining -= nnz;
464ded4f561SHong Zhang     bi[i+1] = bi[i] + nnz;
46508cb4532SHong Zhang     if (nnz > rmax) rmax = nnz;
46642c57489SHong Zhang   }
467ded4f561SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
4680572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
46942c57489SHong Zhang 
47042c57489SHong Zhang   ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
471a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
472d16ca5f0SHong Zhang   afill_tmp = (PetscReal)bi[pn]/(pi_loc[pm] + api[am]);
473d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
47442c57489SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
47542c57489SHong Zhang 
47677584fe6SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
47777584fe6SHong Zhang   time_Cd += tf-t0;
47877584fe6SHong Zhang 
47908cb4532SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part - check runex55_SA ?  */
48008cb4532SHong Zhang   /*------------------------------------------------------------------------------------------------------*/
48108cb4532SHong Zhang   ierr = PetscMalloc((rmax+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
48208cb4532SHong Zhang   ierr = PetscMemzero(vals,rmax*sizeof(PetscScalar));CHKERRQ(ierr);
48308cb4532SHong Zhang   /*
48408cb4532SHong Zhang   ierr = PetscSynchronizedPrintf(comm, "[%d] pn %D, pN %D, cstart %D, rmax %D\n", rank,pn,pN,rstart,rmax);CHKERRQ(ierr);
48508cb4532SHong Zhang   ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr);
48608cb4532SHong Zhang    */
48708cb4532SHong Zhang 
48877584fe6SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
48977584fe6SHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
49077584fe6SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
49177584fe6SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
49242c57489SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
49308cb4532SHong Zhang   ierr = MatSetBlockSize(Cmpi,1);CHKERRQ(ierr);
49408cb4532SHong Zhang   for (i=0; i<pn; i++){
49508cb4532SHong Zhang     row = i + rstart;
49608cb4532SHong Zhang     nnz = bi[i+1] - bi[i];
49708cb4532SHong Zhang     Jptr = bj + bi[i];
49808cb4532SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
49908cb4532SHong Zhang   }
50008cb4532SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
50108cb4532SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
50208cb4532SHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
50342c57489SHong Zhang 
50442c57489SHong Zhang   merge->bi            = bi;
50542c57489SHong Zhang   merge->bj            = bj;
506de0260b3SHong Zhang   merge->coi           = coi;
507de0260b3SHong Zhang   merge->coj           = coj;
50842c57489SHong Zhang   merge->buf_ri        = buf_ri;
50942c57489SHong Zhang   merge->buf_rj        = buf_rj;
510de0260b3SHong Zhang   merge->owners_co     = owners_co;
51177584fe6SHong Zhang   merge->destroy       = Cmpi->ops->destroy;
51277584fe6SHong Zhang   merge->duplicate     = Cmpi->ops->duplicate;
513cc6cb787SHong Zhang 
51477584fe6SHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
51508cb4532SHong Zhang   /* Cmpi->assembled      = PETSC_FALSE; */
51677584fe6SHong Zhang   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_PtAP;
51777584fe6SHong Zhang   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP;
51842c57489SHong Zhang 
51977584fe6SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
52077584fe6SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
521f8487c73SHong Zhang   c->ptap        = ptap;
52277584fe6SHong Zhang   ptap->api      = api;
52377584fe6SHong Zhang   ptap->apj      = apj;
524f8487c73SHong Zhang   ptap->merge    = merge;
525*29825010SHong Zhang   ptap->apnz_max = ap_rmax;
526f8487c73SHong Zhang 
52777584fe6SHong Zhang   *C = Cmpi;
528f2b054eeSHong Zhang #if defined(PETSC_USE_INFO)
529f2b054eeSHong Zhang   if (bi[pn] != 0) {
53077584fe6SHong Zhang     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
53177584fe6SHong Zhang     ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
532f2b054eeSHong Zhang   } else {
53377584fe6SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
534f2b054eeSHong Zhang   }
535f2b054eeSHong Zhang #endif
53677584fe6SHong Zhang   ierr = PetscGetTime(&tff);CHKERRQ(ierr);
53777584fe6SHong Zhang   etime += tff - t00;
53877584fe6SHong Zhang   /*
53977584fe6SHong Zhang   PetscInt prid=0;
54077584fe6SHong Zhang   if (rank == prid){
54177584fe6SHong 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);
54277584fe6SHong Zhang   }
54377584fe6SHong Zhang    */
544fe615159SHong Zhang #if defined(DEBUG_MATPTAP)
545fe615159SHong Zhang   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
546fe615159SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] exit  MatPtAPSymbolic_MPIAIJ_MPIAIJ\n",rank);
547fe615159SHong Zhang #endif
54842c57489SHong Zhang   PetscFunctionReturn(0);
54942c57489SHong Zhang }
55042c57489SHong Zhang 
55142c57489SHong Zhang #undef __FUNCT__
55242c57489SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ"
55342c57489SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
55442c57489SHong Zhang {
55542c57489SHong Zhang   PetscErrorCode       ierr;
55642c57489SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
557f8487c73SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
55842c57489SHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
559de0260b3SHong Zhang   Mat_SeqAIJ           *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
56042c57489SHong Zhang   Mat_SeqAIJ           *p_loc,*p_oth;
561f8487c73SHong Zhang   Mat_PtAPMPI          *ptap;
5621d633602SHong Zhang   PetscInt             *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp;
56382412ba7SHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj;
56482412ba7SHong Zhang   PetscInt             i,j,k,anz,pnz,apnz,nextap,row,*cj;
565e900a757SHong Zhang   MatScalar            *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp;
566d0f46423SBarry Smith   PetscInt             am=A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
5677adad957SLisandro Dalcin   MPI_Comm             comm=((PetscObject)C)->comm;
56842c57489SHong Zhang   PetscMPIInt          size,rank,taga,*len_s;
56982412ba7SHong Zhang   PetscInt             *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
57042c57489SHong Zhang   PetscInt             **buf_ri,**buf_rj;
57150f5bf86SHong Zhang   PetscInt             cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */
57242c57489SHong Zhang   MPI_Request          *s_waits,*r_waits;
57342c57489SHong Zhang   MPI_Status           *status;
57482412ba7SHong Zhang   MatScalar            **abuf_r,*ba_i,*pA,*coa,*ba;
57582412ba7SHong Zhang   PetscInt             *api,*apj,*coi,*coj;
576d0f46423SBarry Smith   PetscInt             *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend;
577a9d06231SHong Zhang   PetscInt             sparse_axpy;
578f9473cf7SHong 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;
57942c57489SHong Zhang 
58042c57489SHong Zhang   PetscFunctionBegin;
581f8487c73SHong Zhang   /* ierr = MPI_Barrier(comm);CHKERRQ(ierr); */
582f9473cf7SHong Zhang   ierr = PetscGetTime(&t00);CHKERRQ(ierr);
58342c57489SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
58442c57489SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
58542c57489SHong Zhang 
586f9473cf7SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
587f8487c73SHong Zhang   ptap  = c->ptap;
588f8487c73SHong Zhang   merge = ptap->merge;
5896c6e00e0SHong Zhang 
590a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
591f9473cf7SHong Zhang   /*--------------------------------------------------*/
592f8487c73SHong Zhang   if (ptap->reuse == MAT_INITIAL_MATRIX){
593f8487c73SHong Zhang     ptap->reuse = MAT_REUSE_MATRIX;
5946c6e00e0SHong Zhang   } else { /* update numerical values of P_oth and P_loc */
595a1a4e84aSHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
596a1a4e84aSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
5976c6e00e0SHong Zhang   }
598f8487c73SHong Zhang 
599f9473cf7SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
600f9473cf7SHong Zhang   time_matupdate += tf-t0;
60142c57489SHong Zhang 
602f9473cf7SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
603f9473cf7SHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
604f9473cf7SHong Zhang   /*--------------------------------------------------------------*/
60542c57489SHong Zhang   /* get data from symbolic products */
606a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
607a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
608a9d06231SHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a;
60942c57489SHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
61042c57489SHong Zhang 
611de0260b3SHong Zhang   coi = merge->coi; coj = merge->coj;
612de0260b3SHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr);
613de0260b3SHong Zhang   ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr);
614de0260b3SHong Zhang 
615de0260b3SHong Zhang   bi     = merge->bi; bj = merge->bj;
6167a2fc3feSBarry Smith   owners = merge->rowmap->range;
617de0260b3SHong Zhang   ierr   = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
618de0260b3SHong Zhang   ierr   = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr);
619f9473cf7SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
620f9473cf7SHong Zhang   time_malloc += tf-t0;
62142c57489SHong Zhang 
622a1a4e84aSHong Zhang   api = ptap->api; apj = ptap->apj;
6230496f32cSHong Zhang   /* flag 'sparse_axpy' determines which implementations to be used:
6240496f32cSHong Zhang        0: do dense axpy in MatPtAPNumeric() - fastest, but requires storage of a dense array apa; (default)
6250496f32cSHong Zhang        1: do one sparse axpy - uses same memory as sparse_axpy=0 and might execute less flops
6260496f32cSHong Zhang           (apnz vs. cnz in the outerproduct), slower than case '0' when cnz is not too large than apnz;
6270496f32cSHong Zhang        2: do two sparse axpy in MatPtAPNumeric() - slowest, uses a sparse array apa */
628a9d06231SHong Zhang   /* set default sparse_axpy */
629e900a757SHong Zhang   sparse_axpy = 0;
6300496f32cSHong Zhang   ierr = PetscOptionsGetInt(PETSC_NULL,"-matptap_sparseaxpy",&sparse_axpy,PETSC_NULL);CHKERRQ(ierr);
631a9d06231SHong Zhang   if (sparse_axpy == 0){  /* Do not perform sparse axpy */
632a9d06231SHong Zhang     /*--------------------------------------------------*/
633a9d06231SHong Zhang     /* malloc apa to store dense row A[i,:]*P */
634a9d06231SHong Zhang     ierr = PetscMalloc((P->cmap->N)*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
635a9d06231SHong Zhang     ierr = PetscMemzero(apa,P->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
636a9d06231SHong Zhang 
637a9d06231SHong Zhang     for (i=0; i<am; i++) {
638a9d06231SHong Zhang       ierr = PetscGetTime(&t0);CHKERRQ(ierr);
639a9d06231SHong Zhang       /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */
640a9d06231SHong Zhang       /*------------------------------------------------------------*/
641a9d06231SHong Zhang       apJ = apj + api[i];
642a9d06231SHong Zhang 
643a9d06231SHong Zhang       /* diagonal portion of A */
644a9d06231SHong Zhang       anz = adi[i+1] - adi[i];
645a9d06231SHong Zhang       adj = ad->j + adi[i];
646a9d06231SHong Zhang       ada = ad->a + adi[i];
647a9d06231SHong Zhang       for (j=0; j<anz; j++) {
648a9d06231SHong Zhang         row = adj[j];
649a9d06231SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
650a9d06231SHong Zhang         pj  = pj_loc + pi_loc[row];
651a9d06231SHong Zhang         pa  = pa_loc + pi_loc[row];
652a9d06231SHong Zhang 
653a9d06231SHong Zhang         /* perform dense axpy */
654e900a757SHong Zhang         valtmp = ada[j];
655a9d06231SHong Zhang         for (k=0; k<pnz; k++){
656e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
657a9d06231SHong Zhang         }
658a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
659a9d06231SHong Zhang       }
660a9d06231SHong Zhang 
661a9d06231SHong Zhang       /* off-diagonal portion of A */
662a9d06231SHong Zhang       anz = aoi[i+1] - aoi[i];
663a9d06231SHong Zhang       aoj = ao->j + aoi[i];
664a9d06231SHong Zhang       aoa = ao->a + aoi[i];
665a9d06231SHong Zhang       for (j=0; j<anz; j++) {
666a9d06231SHong Zhang         row = aoj[j];
667a9d06231SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
668a9d06231SHong Zhang         pj  = pj_oth + pi_oth[row];
669a9d06231SHong Zhang         pa  = pa_oth + pi_oth[row];
670a9d06231SHong Zhang 
671a9d06231SHong Zhang         /* perform dense axpy */
672e900a757SHong Zhang         valtmp = aoa[j];
673a9d06231SHong Zhang         for (k=0; k<pnz; k++){
674e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
675a9d06231SHong Zhang         }
676a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
677a9d06231SHong Zhang       }
678a9d06231SHong Zhang       ierr = PetscGetTime(&tf);CHKERRQ(ierr);
679a9d06231SHong Zhang       time_Cseq0 += tf - t0;
680a9d06231SHong Zhang 
681a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
682a9d06231SHong Zhang       /*--------------------------------------------------------------*/
683a9d06231SHong Zhang       ierr = PetscGetTime(&t0);CHKERRQ(ierr);
684a9d06231SHong Zhang       apnz = api[i+1] - api[i];
685a9d06231SHong Zhang       /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */
686a9d06231SHong Zhang       pnz = po->i[i+1] - po->i[i];
687e900a757SHong Zhang       poJ = po->j + po->i[i];
688a9d06231SHong Zhang       pA  = po->a + po->i[i];
689a9d06231SHong Zhang       for (j=0; j<pnz; j++){
690e900a757SHong Zhang         row = poJ[j];
691e900a757SHong Zhang         cnz = coi[row+1] - coi[row];
692e900a757SHong Zhang         cj  = coj + coi[row];
693e900a757SHong Zhang         ca  = coa + coi[row];
694a9d06231SHong Zhang         /* perform dense axpy */
695e900a757SHong Zhang         valtmp = pA[j];
696a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
697e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
698a9d06231SHong Zhang         }
699a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
700a9d06231SHong Zhang       }
701a9d06231SHong Zhang 
702a9d06231SHong Zhang       /* put the value into Cd (diagonal part) */
703a9d06231SHong Zhang       pnz = pd->i[i+1] - pd->i[i];
704e900a757SHong Zhang       pdJ = pd->j + pd->i[i];
705a9d06231SHong Zhang       pA  = pd->a + pd->i[i];
706a9d06231SHong Zhang       for (j=0; j<pnz; j++){
707e900a757SHong Zhang         row = pdJ[j];
708e900a757SHong Zhang         cnz = bi[row+1] - bi[row];
709e900a757SHong Zhang         cj  = bj + bi[row];
710e900a757SHong Zhang         ca  = ba + bi[row];
711a9d06231SHong Zhang         /* perform dense axpy */
712e900a757SHong Zhang         valtmp = pA[j];
713a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
714e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
715a9d06231SHong Zhang         }
716a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
717a9d06231SHong Zhang       }
718a9d06231SHong Zhang 
719a9d06231SHong Zhang       /* zero the current row of A*P */
720a9d06231SHong Zhang       for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0;
721a9d06231SHong Zhang       ierr = PetscGetTime(&tf);CHKERRQ(ierr);
722a9d06231SHong Zhang       time_Cseq1 += tf - t0;
723a9d06231SHong Zhang     }
724a9d06231SHong Zhang   } else if (sparse_axpy == 1){  /* Perform one sparse axpy */
725a9d06231SHong Zhang     /*------------------------------------------------------*/
7261d633602SHong Zhang     /* malloc apa to store dense row A[i,:]*P */
7271d633602SHong Zhang     ierr = PetscMalloc((P->cmap->N)*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
7281d633602SHong Zhang     ierr = PetscMemzero(apa,P->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
7291d633602SHong Zhang 
7301d633602SHong Zhang     for (i=0; i<am; i++) {
731f9473cf7SHong Zhang       ierr = PetscGetTime(&t0);CHKERRQ(ierr);
732f9473cf7SHong Zhang       /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */
733f9473cf7SHong Zhang       /*------------------------------------------------------------*/
7341d633602SHong Zhang       apJ = apj + api[i];
7351d633602SHong Zhang 
7361d633602SHong Zhang       /* diagonal portion of A */
7371d633602SHong Zhang       anz = adi[i+1] - adi[i];
7381d633602SHong Zhang       adj = ad->j + adi[i];
7391d633602SHong Zhang       ada = ad->a + adi[i];
7401d633602SHong Zhang       for (j=0; j<anz; j++) {
7411d633602SHong Zhang         row = adj[j];
7421d633602SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
7431d633602SHong Zhang         pj  = pj_loc + pi_loc[row];
7441d633602SHong Zhang         pa  = pa_loc + pi_loc[row];
7451d633602SHong Zhang 
7461d633602SHong Zhang         /* perform dense axpy */
747e900a757SHong Zhang         valtmp = ada[j];
7481d633602SHong Zhang         for (k=0; k<pnz; k++){
749e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
7501d633602SHong Zhang         }
7511d633602SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
7521d633602SHong Zhang       }
7531d633602SHong Zhang 
7541d633602SHong Zhang       /* off-diagonal portion of A */
7551d633602SHong Zhang       anz = aoi[i+1] - aoi[i];
7561d633602SHong Zhang       aoj = ao->j + aoi[i];
7571d633602SHong Zhang       aoa = ao->a + aoi[i];
7581d633602SHong Zhang       for (j=0; j<anz; j++) {
7591d633602SHong Zhang         row = aoj[j];
7601d633602SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
7611d633602SHong Zhang         pj  = pj_oth + pi_oth[row];
7621d633602SHong Zhang         pa  = pa_oth + pi_oth[row];
7631d633602SHong Zhang 
7641d633602SHong Zhang         /* perform dense axpy */
765e900a757SHong Zhang         valtmp = aoa[j];
7661d633602SHong Zhang         for (k=0; k<pnz; k++){
767e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
7681d633602SHong Zhang         }
7691d633602SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
7701d633602SHong Zhang       }
771f9473cf7SHong Zhang       ierr = PetscGetTime(&tf);CHKERRQ(ierr);
772f9473cf7SHong Zhang       time_Cseq0 += tf - t0;
7731d633602SHong Zhang 
774f9473cf7SHong Zhang       ierr = PetscGetTime(&t0);CHKERRQ(ierr);
775f9473cf7SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
776f9473cf7SHong Zhang       /*--------------------------------------------------------------*/
7771d633602SHong Zhang       apnz = api[i+1] - api[i];
778f9473cf7SHong Zhang       /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */
779f9473cf7SHong Zhang       pnz = po->i[i+1] - po->i[i];
780e900a757SHong Zhang       poJ = po->j + po->i[i];
781a9d06231SHong Zhang       pA  = po->a + po->i[i];
7821d633602SHong Zhang       for (j=0; j<pnz; j++){
783e900a757SHong Zhang         row    = poJ[j];
784e900a757SHong Zhang         cj     = coj + coi[row];
785e900a757SHong Zhang         ca     = coa + coi[row];
786e900a757SHong Zhang         valtmp = pA[j];
7871d633602SHong Zhang         /* perform sparse axpy */
7881d633602SHong Zhang         nextap = 0;
7891d633602SHong Zhang         for (k=0; nextap<apnz; k++) {
7901d633602SHong Zhang           if (cj[k]==apJ[nextap]) { /* global column index */
791e900a757SHong Zhang             ca[k] += valtmp*apa[cj[k]]; nextap++;
7921d633602SHong Zhang           }
7931d633602SHong Zhang         }
7941d633602SHong Zhang         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
7951d633602SHong Zhang       }
796f9473cf7SHong Zhang 
797f9473cf7SHong Zhang       /* put the value into Cd (diagonal part) */
798f9473cf7SHong Zhang       pnz = pd->i[i+1] - pd->i[i];
799e900a757SHong Zhang       pdJ = pd->j + pd->i[i];
800a9d06231SHong Zhang       pA  = pd->a + pd->i[i];
801f9473cf7SHong Zhang       for (j=0; j<pnz; j++){
802e900a757SHong Zhang         row    = pdJ[j];
803e900a757SHong Zhang         cj     = bj + bi[row];
804e900a757SHong Zhang         ca     = ba + bi[row];
805e900a757SHong Zhang         valtmp = pA[j];
806f9473cf7SHong Zhang         /* perform sparse axpy */
807f9473cf7SHong Zhang         nextap = 0;
808f9473cf7SHong Zhang         for (k=0; nextap<apnz; k++) {
809f9473cf7SHong Zhang           if (cj[k]==apJ[nextap]) { /* global column index */
810e900a757SHong Zhang             ca[k] += valtmp*apa[cj[k]];
811a9d06231SHong Zhang             nextap++;
812f9473cf7SHong Zhang           }
813f9473cf7SHong Zhang         }
814f9473cf7SHong Zhang         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
815f9473cf7SHong Zhang       }
816f9473cf7SHong Zhang 
817f9473cf7SHong Zhang       /* zero the current row of A*P */
8181d633602SHong Zhang       for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0;
819f9473cf7SHong Zhang       ierr = PetscGetTime(&tf);CHKERRQ(ierr);
820f9473cf7SHong Zhang       time_Cseq1 += tf - t0;
8211d633602SHong Zhang     }
8220496f32cSHong Zhang   } else if (sparse_axpy == 2){/* Perform two sparse axpy */
823a9d06231SHong Zhang     /*----------------------------------------------------*/
8241d633602SHong Zhang     /* malloc apa to store sparse row A[i,:]*P */
825*29825010SHong Zhang     ierr = PetscMalloc((ptap->apnz_max+1)*sizeof(MatScalar),&apa);CHKERRQ(ierr);
826*29825010SHong Zhang     ierr = PetscMemzero(apa,ptap->apnz_max*sizeof(MatScalar));CHKERRQ(ierr);
82742c57489SHong Zhang 
828a9d06231SHong Zhang     pA=pa_loc;
82942c57489SHong Zhang     for (i=0; i<am; i++) {
830a9d06231SHong Zhang       ierr = PetscGetTime(&t0);CHKERRQ(ierr);
831f4a743e1SHong Zhang       /* form i-th sparse row of A*P */
83282412ba7SHong Zhang       apnz = api[i+1] - api[i];
83382412ba7SHong Zhang       apJ  = apj + api[i];
83442c57489SHong Zhang       /* diagonal portion of A */
83582412ba7SHong Zhang       anz = adi[i+1] - adi[i];
8361d633602SHong Zhang       adj = ad->j + adi[i];
8371d633602SHong Zhang       ada = ad->a + adi[i];
83882412ba7SHong Zhang       for (j=0; j<anz; j++) {
8391d633602SHong Zhang         row = adj[j];
84082412ba7SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
84182412ba7SHong Zhang         pj  = pj_loc + pi_loc[row];
84282412ba7SHong Zhang         pa  = pa_loc + pi_loc[row];
843e900a757SHong Zhang         valtmp = ada[j];
844f4a743e1SHong Zhang         nextp  = 0;
84582412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
84682412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
847e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
84842c57489SHong Zhang           }
84942c57489SHong Zhang         }
850dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
85142c57489SHong Zhang       }
85242c57489SHong Zhang       /* off-diagonal portion of A */
85382412ba7SHong Zhang       anz = aoi[i+1] - aoi[i];
8541d633602SHong Zhang       aoj = ao->j + aoi[i];
8551d633602SHong Zhang       aoa = ao->a + aoi[i];
85682412ba7SHong Zhang       for (j=0; j<anz; j++) {
8571d633602SHong Zhang         row = aoj[j];
85882412ba7SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
85982412ba7SHong Zhang         pj  = pj_oth + pi_oth[row];
86082412ba7SHong Zhang         pa  = pa_oth + pi_oth[row];
861e900a757SHong Zhang         valtmp = aoa[j];
862f4a743e1SHong Zhang         nextp  = 0;
86382412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
86482412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
865e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
86642c57489SHong Zhang           }
86742c57489SHong Zhang         }
868dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
86942c57489SHong Zhang       }
870a9d06231SHong Zhang       ierr = PetscGetTime(&tf);CHKERRQ(ierr);
871a9d06231SHong Zhang       time_Cseq0 += tf - t0;
87242c57489SHong Zhang 
873a9d06231SHong Zhang       ierr = PetscGetTime(&t0);CHKERRQ(ierr);
874a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
875a9d06231SHong Zhang       /*--------------------------------------------------------------*/
87682412ba7SHong Zhang       pnz = pi_loc[i+1] - pi_loc[i];
877f9473cf7SHong Zhang       pJ  = pj_loc + pi_loc[i];
87882412ba7SHong Zhang       for (j=0; j<pnz; j++) {
87942c57489SHong Zhang         nextap = 0;
880f9473cf7SHong Zhang         row    = pJ[j]; /* global index */
88182412ba7SHong Zhang         if (row < pcstart || row >=pcend) { /* put the value into Co */
882e900a757SHong Zhang           row = *poJ;
883e900a757SHong Zhang           cj  = coj + coi[row];
884e900a757SHong Zhang           ca  = coa + coi[row]; poJ++;
88582412ba7SHong Zhang         } else {                            /* put the value into Cd */
886e900a757SHong Zhang           row  = *pdJ;
887e900a757SHong Zhang           cj   = bj + bi[row];
888e900a757SHong Zhang           ca   = ba + bi[row]; pdJ++;
88942c57489SHong Zhang         }
890e900a757SHong Zhang         valtmp = pA[j];
89182412ba7SHong Zhang         for (k=0; nextap<apnz; k++) {
892e900a757SHong Zhang           if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++];
89342c57489SHong Zhang         }
894dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
895de0260b3SHong Zhang       }
8960496f32cSHong Zhang       pA += pnz;
89742c57489SHong Zhang       /* zero the current row info for A*P */
89882412ba7SHong Zhang       ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr);
899a9d06231SHong Zhang       ierr = PetscGetTime(&tf);CHKERRQ(ierr);
900a9d06231SHong Zhang       time_Cseq1 += tf - t0;
90142c57489SHong Zhang     }
902a9d06231SHong Zhang   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"sparse_axpy only takes values 0, 1 and 2");
90342c57489SHong Zhang   ierr = PetscFree(apa);CHKERRQ(ierr);
90442c57489SHong Zhang 
905a9d06231SHong Zhang   /* 3) send and recv matrix values coa */
906a9d06231SHong Zhang   /*------------------------------------*/
90742c57489SHong Zhang   buf_ri = merge->buf_ri;
90842c57489SHong Zhang   buf_rj = merge->buf_rj;
90942c57489SHong Zhang   len_s  = merge->len_s;
910fc42d0c8SSatish Balay   ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
91142c57489SHong Zhang   ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
91242c57489SHong Zhang 
913a9d06231SHong Zhang   ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr);
91442c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++){
91542c57489SHong Zhang     if (!len_s[proc]) continue;
916de0260b3SHong Zhang     i = merge->owners_co[proc];
917de0260b3SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
91842c57489SHong Zhang     k++;
91942c57489SHong Zhang   }
9200c468ba9SBarry Smith   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
9210c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
92242c57489SHong Zhang 
923a9d06231SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
92442c57489SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
92582412ba7SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
92642c57489SHong Zhang 
927a9d06231SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
928a9d06231SHong Zhang   /*------------------------------------------------------*/
929f9473cf7SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
930687f1162SBarry Smith   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
93142c57489SHong Zhang 
93242c57489SHong Zhang   for (k=0; k<merge->nrecv; k++){
93342c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
93442c57489SHong Zhang     nrows       = *(buf_ri_k[k]);
93542c57489SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
93682412ba7SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
93742c57489SHong Zhang   }
93842c57489SHong Zhang 
939de0260b3SHong Zhang   for (i=0; i<cm; i++) {
94082412ba7SHong Zhang     row = owners[rank] + i; /* global row index of C_seq */
94142c57489SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
942de0260b3SHong Zhang     ba_i = ba + bi[i];
94382412ba7SHong Zhang     bnz  = bi[i+1] - bi[i];
94442c57489SHong Zhang     /* add received vals into ba */
94542c57489SHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
94642c57489SHong Zhang       /* i-th row */
94742c57489SHong Zhang       if (i == *nextrow[k]) {
94882412ba7SHong Zhang         cnz = *(nextci[k]+1) - *nextci[k];
94982412ba7SHong Zhang         cj  = buf_rj[k] + *(nextci[k]);
95082412ba7SHong Zhang         ca  = abuf_r[k] + *(nextci[k]);
95182412ba7SHong Zhang         nextcj = 0;
95282412ba7SHong Zhang         for (j=0; nextcj<cnz; j++){
95382412ba7SHong Zhang           if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */
95482412ba7SHong Zhang             ba_i[j] += ca[nextcj++];
95542c57489SHong Zhang           }
95642c57489SHong Zhang         }
95782412ba7SHong Zhang         nextrow[k]++; nextci[k]++;
95842c57489SHong Zhang       }
95942c57489SHong Zhang     }
96082412ba7SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
961dc0b31edSSatish Balay     ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
96242c57489SHong Zhang   }
96342c57489SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
96442c57489SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
965f9473cf7SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
966f9473cf7SHong Zhang   time_setvals += tf-t0;
96742c57489SHong Zhang 
968de0260b3SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
969c05d87d6SBarry Smith   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
97042c57489SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
9710572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
972f9473cf7SHong Zhang 
973f9473cf7SHong Zhang   ierr = PetscGetTime(&tff);CHKERRQ(ierr);
974f9473cf7SHong Zhang   etime += tff - t00;
975f9473cf7SHong Zhang   /*
976a9d06231SHong Zhang   PetscInt prid=0;
977a9d06231SHong Zhang   if (rank == prid){
978f9473cf7SHong 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);
979a9d06231SHong Zhang   }
980f9473cf7SHong Zhang    */
98142c57489SHong Zhang   PetscFunctionReturn(0);
98242c57489SHong Zhang }
983