xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision 6c4ed00291fe44f94936dd2f04c02ab3c442e77c)
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>
118563dfccSBarry Smith #include <petsctime.h>
1242c57489SHong Zhang 
13f671be37SHong Zhang /* #define PTAP_PROFILE */
1424ecddacSHong Zhang 
1509573ac7SBarry Smith extern PetscErrorCode MatDestroy_MPIAIJ(Mat);
16cc6cb787SHong Zhang #undef __FUNCT__
17f8487c73SHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_PtAP"
18f8487c73SHong Zhang PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A)
19cc6cb787SHong Zhang {
20cc6cb787SHong Zhang   PetscErrorCode ierr;
21f8487c73SHong Zhang   Mat_MPIAIJ     *a=(Mat_MPIAIJ*)A->data;
22f8487c73SHong Zhang   Mat_PtAPMPI    *ptap=a->ptap;
23cc6cb787SHong Zhang 
24cc6cb787SHong Zhang   PetscFunctionBegin;
25f8487c73SHong Zhang   if (ptap) {
26c8b0795cSMark F. Adams     Mat_Merge_SeqsToMPI *merge=ptap->merge;
27b7f45c76SHong Zhang     ierr = PetscFree2(ptap->startsj_s,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);
31c5af039cSHong Zhang     ierr = MatDestroy(&ptap->A_loc);CHKERRQ(ierr); /* used by MatTransposeMatMult() */
3205d62848SHong Zhang     ierr = MatDestroy(&ptap->Rd);CHKERRQ(ierr);
3305d62848SHong Zhang     ierr = MatDestroy(&ptap->Ro);CHKERRQ(ierr);
348403a639SHong Zhang     if (ptap->AP_loc) { /* used by alg_rap */
35681d504bSHong Zhang       Mat_SeqAIJ *ap = (Mat_SeqAIJ*)(ptap->AP_loc)->data;
36681d504bSHong Zhang       ierr = PetscFree(ap->i);CHKERRQ(ierr);
37681d504bSHong Zhang       ierr = PetscFree2(ap->j,ap->a);CHKERRQ(ierr);
380d3441aeSHong Zhang       ierr = MatDestroy(&ptap->AP_loc);CHKERRQ(ierr);
398403a639SHong Zhang     } else { /* used by alg_ptap */
408403a639SHong Zhang       ierr = PetscFree(ptap->api);CHKERRQ(ierr);
418403a639SHong Zhang       ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
42681d504bSHong Zhang     }
432259aa2eSHong Zhang     ierr = MatDestroy(&ptap->C_loc);CHKERRQ(ierr);
442259aa2eSHong Zhang     ierr = MatDestroy(&ptap->C_oth);CHKERRQ(ierr);
45414894bdSHong Zhang     if (ptap->apa) {ierr = PetscFree(ptap->apa);CHKERRQ(ierr);}
468403a639SHong Zhang     if (merge) { /* used by alg_ptap */
47cc6cb787SHong Zhang       ierr = PetscFree(merge->id_r);CHKERRQ(ierr);
48cc6cb787SHong Zhang       ierr = PetscFree(merge->len_s);CHKERRQ(ierr);
49cc6cb787SHong Zhang       ierr = PetscFree(merge->len_r);CHKERRQ(ierr);
50cc6cb787SHong Zhang       ierr = PetscFree(merge->bi);CHKERRQ(ierr);
51cc6cb787SHong Zhang       ierr = PetscFree(merge->bj);CHKERRQ(ierr);
52c05d87d6SBarry Smith       ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr);
53cc6cb787SHong Zhang       ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr);
54c05d87d6SBarry Smith       ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr);
55cc6cb787SHong Zhang       ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr);
56445158ffSHong Zhang       ierr = PetscFree(merge->coi);CHKERRQ(ierr);
57445158ffSHong Zhang       ierr = PetscFree(merge->coj);CHKERRQ(ierr);
5805b42c5fSBarry Smith       ierr = PetscFree(merge->owners_co);CHKERRQ(ierr);
596bf464f9SBarry Smith       ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr);
60dce485f0SBarry Smith       ierr = merge->destroy(A);CHKERRQ(ierr);
61f8487c73SHong Zhang       ierr = PetscFree(ptap->merge);CHKERRQ(ierr);
62445158ffSHong Zhang     } else {
63445158ffSHong Zhang       ierr = MatDestroy_MPIAIJ(A);CHKERRQ(ierr);
64bf0cc555SLisandro Dalcin     }
65f8487c73SHong Zhang     ierr = PetscFree(ptap);CHKERRQ(ierr);
66c8b0795cSMark F. Adams   }
67cc6cb787SHong Zhang   PetscFunctionReturn(0);
68cc6cb787SHong Zhang }
69cc6cb787SHong Zhang 
70cc6cb787SHong Zhang #undef __FUNCT__
71aa690a28SHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatPtAP"
72aa690a28SHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M)
731a47ec54SHong Zhang {
741a47ec54SHong Zhang   PetscErrorCode ierr;
751a47ec54SHong Zhang   Mat_MPIAIJ     *a     = (Mat_MPIAIJ*)A->data;
761a47ec54SHong Zhang   Mat_PtAPMPI    *ptap  = a->ptap;
771a47ec54SHong Zhang 
781a47ec54SHong Zhang   PetscFunctionBegin;
791a47ec54SHong Zhang   ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr);
801a47ec54SHong Zhang   (*M)->ops->destroy   = ptap->destroy;
811a47ec54SHong Zhang   (*M)->ops->duplicate = ptap->duplicate;
821a47ec54SHong Zhang   PetscFunctionReturn(0);
831a47ec54SHong Zhang }
841a47ec54SHong Zhang 
8542c57489SHong Zhang #undef __FUNCT__
86cf3ca8ceSHong Zhang #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ"
87cf3ca8ceSHong Zhang PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
8842c57489SHong Zhang {
8942c57489SHong Zhang   PetscErrorCode ierr;
908403a639SHong Zhang   PetscBool      rap=PETSC_TRUE; /* do R=P^T locally, then C=R*A*P */
9167a12041SHong Zhang   MPI_Comm       comm;
9242c57489SHong Zhang 
9342c57489SHong Zhang   PetscFunctionBegin;
9467a12041SHong Zhang   /* check if matrix local sizes are compatible */
9567a12041SHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
96*6c4ed002SBarry Smith   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, Arow (%D, %D) != Prow (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend);
97*6c4ed002SBarry Smith   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, Acol (%D, %D) != Prow (%D,%D)",A->cmap->rstart,A->cmap->rend,P->rmap->rstart,P->rmap->rend);
9867a12041SHong Zhang 
998403a639SHong Zhang   ierr = PetscOptionsGetBool(NULL,"-matrap",&rap,NULL);CHKERRQ(ierr);
100cf3ca8ceSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
1013ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
1028403a639SHong Zhang     if (rap) { /* do R=P^T locally, then C=R*A*P */
103cf3ca8ceSHong Zhang       ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);
1048403a639SHong Zhang     } else {       /* do P^T*A*P */
1058403a639SHong Zhang       ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(A,P,fill,C);CHKERRQ(ierr);
1060d3441aeSHong Zhang     }
1073ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
1087a7894deSKris Buschelman   }
1093ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
1108403a639SHong Zhang   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
1113ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
11242c57489SHong Zhang   PetscFunctionReturn(0);
11342c57489SHong Zhang }
11442c57489SHong Zhang 
11542c57489SHong Zhang #undef __FUNCT__
116aa690a28SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ"
117aa690a28SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
1180d3441aeSHong Zhang {
1190d3441aeSHong Zhang   PetscErrorCode      ierr;
1200d3441aeSHong Zhang   Mat_PtAPMPI         *ptap;
121681d504bSHong Zhang   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
1222259aa2eSHong Zhang   MPI_Comm            comm;
1232259aa2eSHong Zhang   PetscMPIInt         size,rank;
12415a3b8e2SHong Zhang   Mat                 Cmpi;
12515a3b8e2SHong Zhang   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
126d0e9a2f7SHong Zhang   PetscInt            am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
127d0e9a2f7SHong Zhang   PetscInt            *lnk,i,k,pnz,row,nsend;
12815a3b8e2SHong Zhang   PetscBT             lnkbt;
129f671be37SHong Zhang   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
13015a3b8e2SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
131d0e9a2f7SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners,nzi,nspacedouble;
13215a3b8e2SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
13315a3b8e2SHong Zhang   MPI_Request         *swaits,*rwaits;
13415a3b8e2SHong Zhang   MPI_Status          *sstatus,rstatus;
135445158ffSHong Zhang   PetscLayout         rowmap;
136445158ffSHong Zhang   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
137445158ffSHong Zhang   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
138d0e9a2f7SHong Zhang   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0,Crmax,rmax,*aj,*ai,*pi;
13967a12041SHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*c_loc,*c_oth;
14067a12041SHong Zhang   PetscScalar         *apv;
14178d30b94SHong Zhang   PetscTable          ta;
14278d30b94SHong Zhang   const char          *algTypes[2] = {"scalable","nonscalable"};
14378d30b94SHong Zhang   PetscInt            alg=1; /* set default algorithm */
144aa690a28SHong Zhang #if defined(PETSC_USE_INFO)
1458cb82516SHong Zhang   PetscReal           apfill;
146aa690a28SHong Zhang #endif
147681d504bSHong Zhang #if defined(PTAP_PROFILE)
148681d504bSHong Zhang   PetscLogDouble      t0,t1,t11,t12,t2,t3,t4;
149681d504bSHong Zhang #endif
15067a12041SHong Zhang 
15167a12041SHong Zhang   PetscFunctionBegin;
15267a12041SHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
15367a12041SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
15467a12041SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1558cb82516SHong Zhang #if defined(PTAP_PROFILE)
15680bb4639SHong Zhang   ierr = PetscTime(&t0);CHKERRQ(ierr);
1578cb82516SHong Zhang #endif
1588cb82516SHong Zhang 
15915a3b8e2SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
16015a3b8e2SHong Zhang   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
16115a3b8e2SHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
16215a3b8e2SHong Zhang 
16315a3b8e2SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
16415a3b8e2SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
16515a3b8e2SHong Zhang   /* get P_loc by taking all local rows of P */
16615a3b8e2SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
16715a3b8e2SHong Zhang 
16867a12041SHong Zhang   /* (0) compute Rd = Pd^T, Ro = Po^T  */
16967a12041SHong Zhang   /* --------------------------------- */
170de817e96SHong Zhang   ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
171de817e96SHong Zhang   ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
1728cb82516SHong Zhang #if defined(PTAP_PROFILE)
173445158ffSHong Zhang   ierr = PetscTime(&t11);CHKERRQ(ierr);
1748cb82516SHong Zhang #endif
17515a3b8e2SHong Zhang 
17667a12041SHong Zhang   /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */
17767a12041SHong Zhang   /* ----------------------------------------------------------------- */
17867a12041SHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
17967a12041SHong Zhang   p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
180445158ffSHong Zhang 
1819a6dcab7SHong Zhang   /* create and initialize a linked list */
182d0e9a2f7SHong Zhang   Crmax = 5*(p_loc->rmax+p_oth->rmax + (PetscInt)(1.e-2*pN)); /* expected Crmax */
183d0e9a2f7SHong Zhang   if (Crmax > pN) Crmax=pN;
184d0e9a2f7SHong Zhang   ierr = PetscTableCreate(Crmax,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */
1854b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta);
1864b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta);
18778d30b94SHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */
188d0e9a2f7SHong Zhang   /* printf("[%d] est %d, Crmax %d; pN %d\n",rank,5*(p_loc->rmax+p_oth->rmax + (PetscInt)(1.e-2*pN)),Crmax,pN); */
18978d30b94SHong Zhang 
19078d30b94SHong Zhang   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
191445158ffSHong Zhang 
1928cb82516SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
19367a12041SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ad->i[am]+ao->i[am]+p_loc->i[pm])),&free_space);CHKERRQ(ierr);
194445158ffSHong Zhang   current_space = free_space;
19567a12041SHong Zhang   nspacedouble  = 0;
19667a12041SHong Zhang 
19767a12041SHong Zhang   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
19867a12041SHong Zhang   api[0] = 0;
199445158ffSHong Zhang   for (i=0; i<am; i++) {
20067a12041SHong Zhang     /* diagonal portion: Ad[i,:]*P */
20167a12041SHong Zhang     ai = ad->i; pi = p_loc->i;
20267a12041SHong Zhang     nzi = ai[i+1] - ai[i];
20367a12041SHong Zhang     aj  = ad->j + ai[i];
204445158ffSHong Zhang     for (j=0; j<nzi; j++) {
205445158ffSHong Zhang       row  = aj[j];
20667a12041SHong Zhang       pnz  = pi[row+1] - pi[row];
20767a12041SHong Zhang       Jptr = p_loc->j + pi[row];
208445158ffSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
209445158ffSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
210445158ffSHong Zhang     }
21167a12041SHong Zhang     /* off-diagonal portion: Ao[i,:]*P */
21267a12041SHong Zhang     ai = ao->i; pi = p_oth->i;
21367a12041SHong Zhang     nzi = ai[i+1] - ai[i];
21467a12041SHong Zhang     aj  = ao->j + ai[i];
215445158ffSHong Zhang     for (j=0; j<nzi; j++) {
216445158ffSHong Zhang       row  = aj[j];
21767a12041SHong Zhang       pnz  = pi[row+1] - pi[row];
21867a12041SHong Zhang       Jptr = p_oth->j + pi[row];
219445158ffSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
220445158ffSHong Zhang     }
221445158ffSHong Zhang     apnz     = lnk[0];
222445158ffSHong Zhang     api[i+1] = api[i] + apnz;
223445158ffSHong Zhang     if (ap_rmax < apnz) ap_rmax = apnz;
224445158ffSHong Zhang 
225445158ffSHong Zhang     /* if free space is not available, double the total space in the list */
226445158ffSHong Zhang     if (current_space->local_remaining<apnz) {
227445158ffSHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
228445158ffSHong Zhang       nspacedouble++;
229445158ffSHong Zhang     }
230445158ffSHong Zhang 
231445158ffSHong Zhang     /* Copy data into free space, then initialize lnk */
232445158ffSHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
233445158ffSHong Zhang 
234445158ffSHong Zhang     current_space->array           += apnz;
235445158ffSHong Zhang     current_space->local_used      += apnz;
236445158ffSHong Zhang     current_space->local_remaining -= apnz;
237445158ffSHong Zhang   }
238681d504bSHong Zhang   /* Allocate space for apj and apv, initialize apj, and */
239445158ffSHong Zhang   /* destroy list of free space and other temporary array(s) */
240445158ffSHong Zhang   ierr   = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr);
241445158ffSHong Zhang   ierr   = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
2429a6dcab7SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
2439a6dcab7SHong Zhang 
244aa690a28SHong Zhang   /* Create AP_loc for reuse */
245445158ffSHong Zhang   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr);
246aa690a28SHong Zhang 
247aa690a28SHong Zhang #if defined(PETSC_USE_INFO)
248aa690a28SHong Zhang   apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1);
249aa690a28SHong Zhang   ptap->AP_loc->info.mallocs           = nspacedouble;
250aa690a28SHong Zhang   ptap->AP_loc->info.fill_ratio_given  = fill;
251aa690a28SHong Zhang   ptap->AP_loc->info.fill_ratio_needed = apfill;
252aa690a28SHong Zhang 
253aa690a28SHong Zhang   if (api[am]) {
254aa690a28SHong Zhang     ierr = PetscInfo3(ptap->AP_loc,"AP_loc reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)apfill);CHKERRQ(ierr);
255aa690a28SHong Zhang     ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr);
256aa690a28SHong Zhang   } else {
257aa690a28SHong Zhang     ierr = PetscInfo(ptap->AP_loc,"AP_loc is empty \n");CHKERRQ(ierr);
258aa690a28SHong Zhang   }
259aa690a28SHong Zhang #endif
260aa690a28SHong Zhang 
2618cb82516SHong Zhang #if defined(PTAP_PROFILE)
262445158ffSHong Zhang   ierr = PetscTime(&t12);CHKERRQ(ierr);
2638cb82516SHong Zhang #endif
26415a3b8e2SHong Zhang 
265681d504bSHong Zhang   /* (2-1) compute symbolic Co = Ro*AP_loc  */
266681d504bSHong Zhang   /* ------------------------------------ */
26767a12041SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
26867a12041SHong Zhang   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
2698cb82516SHong Zhang #if defined(PTAP_PROFILE)
27080bb4639SHong Zhang   ierr = PetscTime(&t1);CHKERRQ(ierr);
2718cb82516SHong Zhang #endif
27215a3b8e2SHong Zhang 
27367a12041SHong Zhang   /* (3) send coj of C_oth to other processors  */
27467a12041SHong Zhang   /* ------------------------------------------ */
27515a3b8e2SHong Zhang   /* determine row ownership */
276445158ffSHong Zhang   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
277445158ffSHong Zhang   rowmap->n  = pn;
278445158ffSHong Zhang   rowmap->bs = 1;
279445158ffSHong Zhang   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
280445158ffSHong Zhang   owners = rowmap->range;
28115a3b8e2SHong Zhang 
28215a3b8e2SHong Zhang   /* determine the number of messages to send, their lengths */
2838cb82516SHong Zhang   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
2848cb82516SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
28515a3b8e2SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
28615a3b8e2SHong Zhang 
28767a12041SHong Zhang   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
28867a12041SHong Zhang   coi   = c_oth->i; coj = c_oth->j;
28967a12041SHong Zhang   con   = ptap->C_oth->rmap->n;
29015a3b8e2SHong Zhang   proc  = 0;
29167a12041SHong Zhang   for (i=0; i<con; i++) {
29215a3b8e2SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
29315a3b8e2SHong Zhang     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
29415a3b8e2SHong Zhang     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
29515a3b8e2SHong Zhang   }
29615a3b8e2SHong Zhang 
29715a3b8e2SHong Zhang   len          = 0; /* max length of buf_si[], see (4) */
29815a3b8e2SHong Zhang   owners_co[0] = 0;
29967a12041SHong Zhang   nsend        = 0;
30015a3b8e2SHong Zhang   for (proc=0; proc<size; proc++) {
30115a3b8e2SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
30215a3b8e2SHong Zhang     if (len_s[proc]) {
303445158ffSHong Zhang       nsend++;
30415a3b8e2SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
30515a3b8e2SHong Zhang       len         += len_si[proc];
30615a3b8e2SHong Zhang     }
30715a3b8e2SHong Zhang   }
30815a3b8e2SHong Zhang 
30915a3b8e2SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
310445158ffSHong Zhang   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
311445158ffSHong Zhang   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
31215a3b8e2SHong Zhang 
31315a3b8e2SHong Zhang   /* post the Irecv and Isend of coj */
31415a3b8e2SHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
315445158ffSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
316445158ffSHong Zhang   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
31715a3b8e2SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
31815a3b8e2SHong Zhang     if (!len_s[proc]) continue;
31915a3b8e2SHong Zhang     i    = owners_co[proc];
32015a3b8e2SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
32115a3b8e2SHong Zhang     k++;
32215a3b8e2SHong Zhang   }
32315a3b8e2SHong Zhang 
324681d504bSHong Zhang   /* (2-2) compute symbolic C_loc = Rd*AP_loc */
325681d504bSHong Zhang   /* ---------------------------------------- */
326681d504bSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr);
327681d504bSHong Zhang   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
328681d504bSHong Zhang 
329681d504bSHong Zhang   /* receives coj are complete */
330445158ffSHong Zhang   for (i=0; i<nrecv; i++) {
331445158ffSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
33215a3b8e2SHong Zhang   }
33315a3b8e2SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
334445158ffSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
33515a3b8e2SHong Zhang 
33678d30b94SHong Zhang   /* add received column indices into ta to update Crmax */
33778d30b94SHong Zhang   for (k=0; k<nrecv; k++) {/* k-th received message */
33878d30b94SHong Zhang     Jptr = buf_rj[k];
33978d30b94SHong Zhang     for (j=0; j<len_r[k]; j++) {
34078d30b94SHong Zhang       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
34178d30b94SHong Zhang     }
34278d30b94SHong Zhang   }
34378d30b94SHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
34478d30b94SHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
3459a6dcab7SHong Zhang 
34615a3b8e2SHong Zhang   /* (4) send and recv coi */
34715a3b8e2SHong Zhang   /*-----------------------*/
34815a3b8e2SHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
349445158ffSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
35015a3b8e2SHong Zhang   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
35115a3b8e2SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
35215a3b8e2SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
35315a3b8e2SHong Zhang     if (!len_s[proc]) continue;
35415a3b8e2SHong Zhang     /* form outgoing message for i-structure:
35515a3b8e2SHong Zhang          buf_si[0]:                 nrows to be sent
35615a3b8e2SHong Zhang                [1:nrows]:           row index (global)
35715a3b8e2SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
35815a3b8e2SHong Zhang     */
35915a3b8e2SHong Zhang     /*-------------------------------------------*/
36015a3b8e2SHong Zhang     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
36115a3b8e2SHong Zhang     buf_si_i    = buf_si + nrows+1;
36215a3b8e2SHong Zhang     buf_si[0]   = nrows;
36315a3b8e2SHong Zhang     buf_si_i[0] = 0;
36415a3b8e2SHong Zhang     nrows       = 0;
36515a3b8e2SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
36615a3b8e2SHong Zhang       nzi = coi[i+1] - coi[i];
36715a3b8e2SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
36815a3b8e2SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
36915a3b8e2SHong Zhang       nrows++;
37015a3b8e2SHong Zhang     }
37115a3b8e2SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
37215a3b8e2SHong Zhang     k++;
37315a3b8e2SHong Zhang     buf_si += len_si[proc];
37415a3b8e2SHong Zhang   }
375681d504bSHong Zhang   for (i=0; i<nrecv; i++) {
376445158ffSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
37715a3b8e2SHong Zhang   }
37815a3b8e2SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
379445158ffSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
38015a3b8e2SHong Zhang 
3818cb82516SHong Zhang   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
38215a3b8e2SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
38315a3b8e2SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
38415a3b8e2SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
3858cb82516SHong Zhang #if defined(PTAP_PROFILE)
38680bb4639SHong Zhang   ierr = PetscTime(&t2);CHKERRQ(ierr);
3878cb82516SHong Zhang #endif
38867a12041SHong Zhang   /* (5) compute the local portion of Cmpi      */
38967a12041SHong Zhang   /* ------------------------------------------ */
39078d30b94SHong Zhang   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
39178d30b94SHong Zhang   ierr          = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr);
39215a3b8e2SHong Zhang   current_space = free_space;
39315a3b8e2SHong Zhang 
394445158ffSHong Zhang   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
395445158ffSHong Zhang   for (k=0; k<nrecv; k++) {
39615a3b8e2SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
39715a3b8e2SHong Zhang     nrows       = *buf_ri_k[k];
39815a3b8e2SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
39915a3b8e2SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
40015a3b8e2SHong Zhang   }
401445158ffSHong Zhang 
40278d30b94SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
40378d30b94SHong Zhang   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
40415a3b8e2SHong Zhang   rmax = 0;
40515a3b8e2SHong Zhang   for (i=0; i<pn; i++) {
40667a12041SHong Zhang     /* add C_loc into Cmpi */
40767a12041SHong Zhang     nzi  = c_loc->i[i+1] - c_loc->i[i];
40867a12041SHong Zhang     Jptr = c_loc->j + c_loc->i[i];
40967a12041SHong Zhang     ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
41015a3b8e2SHong Zhang 
41115a3b8e2SHong Zhang     /* add received col data into lnk */
412445158ffSHong Zhang     for (k=0; k<nrecv; k++) { /* k-th received message */
41315a3b8e2SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
41415a3b8e2SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
41515a3b8e2SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
41615a3b8e2SHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
41715a3b8e2SHong Zhang         nextrow[k]++; nextci[k]++;
41815a3b8e2SHong Zhang       }
41915a3b8e2SHong Zhang     }
420d0e9a2f7SHong Zhang     nzi = lnk[0];
4218cb82516SHong Zhang 
42215a3b8e2SHong Zhang     /* copy data into free space, then initialize lnk */
423d0e9a2f7SHong Zhang     ierr = PetscLLCondensedClean(pN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
424d0e9a2f7SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr);
425d0e9a2f7SHong Zhang     if (nzi > rmax) rmax = nzi;
42615a3b8e2SHong Zhang   }
42715a3b8e2SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
42815a3b8e2SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
429445158ffSHong Zhang   ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr);
4308cb82516SHong Zhang #if defined(PTAP_PROFILE)
43180bb4639SHong Zhang   ierr = PetscTime(&t3);CHKERRQ(ierr);
4328cb82516SHong Zhang #endif
43380bb4639SHong Zhang 
43415a3b8e2SHong Zhang   /* (6) create symbolic parallel matrix Cmpi */
43515a3b8e2SHong Zhang   /*------------------------------------------*/
43615a3b8e2SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
43715a3b8e2SHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
43815a3b8e2SHong Zhang   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
43915a3b8e2SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
44015a3b8e2SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
44115a3b8e2SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
44215a3b8e2SHong Zhang 
443445158ffSHong Zhang   /* members in merge */
444445158ffSHong Zhang   ierr = PetscFree(id_r);CHKERRQ(ierr);
445445158ffSHong Zhang   ierr = PetscFree(len_r);CHKERRQ(ierr);
446445158ffSHong Zhang   ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr);
447445158ffSHong Zhang   ierr = PetscFree(buf_ri);CHKERRQ(ierr);
448445158ffSHong Zhang   ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr);
449445158ffSHong Zhang   ierr = PetscFree(buf_rj);CHKERRQ(ierr);
450445158ffSHong Zhang   ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr);
45115a3b8e2SHong Zhang 
45215a3b8e2SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
45315a3b8e2SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
45415a3b8e2SHong Zhang   c->ptap         = ptap;
45515a3b8e2SHong Zhang   ptap->rmax      = ap_rmax;
4561a47ec54SHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
4571a47ec54SHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
45815a3b8e2SHong Zhang 
45978d30b94SHong Zhang   ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
46078d30b94SHong Zhang   ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,2,algTypes[1],&alg,NULL);CHKERRQ(ierr);
46178d30b94SHong Zhang   ierr = PetscOptionsEnd();CHKERRQ(ierr);
46278d30b94SHong Zhang 
46378d30b94SHong Zhang   if (alg == 1) {
46478d30b94SHong Zhang     /* Do dense axpy in MatPtAPNumeric_MPIAIJ_MPIAIJ() */
4658cb82516SHong Zhang     ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
46678d30b94SHong Zhang     Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ;
46778d30b94SHong Zhang   } else {
468e55be12dSHong Zhang     Cmpi->ops->ptapnumeric = 0; /* not done yet */
46978d30b94SHong Zhang     SETERRQ(comm,PETSC_ERR_ARG_SIZ,"Not done yet");
47078d30b94SHong Zhang   }
4712259aa2eSHong Zhang 
4721a47ec54SHong Zhang 
4731a47ec54SHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
4741a47ec54SHong Zhang   Cmpi->assembled        = PETSC_FALSE;
4751a47ec54SHong Zhang   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
476aa690a28SHong Zhang   Cmpi->ops->duplicate   = MatDuplicate_MPIAIJ_MatPtAP;
4771a47ec54SHong Zhang   *C                     = Cmpi;
4781a47ec54SHong Zhang 
4798cb82516SHong Zhang #if defined(PTAP_PROFILE)
48080bb4639SHong Zhang   ierr = PetscTime(&t4);CHKERRQ(ierr);
481e9c1f85fSHong Zhang   if (rank == 1) {
482445158ffSHong Zhang     printf("PtAPSym: %g + %g + %g + %g + %g + %g = %g\n",t11-t0,t1-t11,t12-t11,t2-t2,t3-t2,t4-t3,t4-t0);
483e9c1f85fSHong Zhang   }
4848cb82516SHong Zhang #endif
4850d3441aeSHong Zhang   PetscFunctionReturn(0);
4860d3441aeSHong Zhang }
4870d3441aeSHong Zhang 
4880d3441aeSHong Zhang #undef __FUNCT__
489aa690a28SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ"
490aa690a28SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
4910d3441aeSHong Zhang {
4920d3441aeSHong Zhang   PetscErrorCode    ierr;
4930d3441aeSHong Zhang   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
4940d3441aeSHong Zhang   Mat_SeqAIJ        *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
4958cb82516SHong Zhang   Mat_SeqAIJ        *ap,*p_loc,*p_oth,*c_seq;
4960d3441aeSHong Zhang   Mat_PtAPMPI       *ptap = c->ptap;
4979ce11a7cSHong Zhang   Mat               AP_loc,C_loc,C_oth;
4980d3441aeSHong Zhang   PetscInt          i,rstart,rend,cm,ncols,row;
4998cb82516SHong Zhang   PetscInt          *api,*apj,am = A->rmap->n,j,col,apnz;
5008cb82516SHong Zhang   PetscScalar       *apa;
5010d3441aeSHong Zhang   const PetscInt    *cols;
5020d3441aeSHong Zhang   const PetscScalar *vals;
5038cb82516SHong Zhang #if defined(PTAP_PROFILE)
5048cb82516SHong Zhang   PetscMPIInt       rank;
5058cb82516SHong Zhang   MPI_Comm          comm;
5060d3441aeSHong Zhang   PetscLogDouble    t0,t1,t2,t3,t4,eR,eAP,eCseq,eCmpi;
5078cb82516SHong Zhang #endif
5080d3441aeSHong Zhang 
5090d3441aeSHong Zhang   PetscFunctionBegin;
5100d3441aeSHong Zhang   ierr = MatZeroEntries(C);CHKERRQ(ierr);
5110d3441aeSHong Zhang 
512e497d3c8SHong Zhang   /* 1) get R = Pd^T,Ro = Po^T */
5138cb82516SHong Zhang #if defined(PTAP_PROFILE)
5140d3441aeSHong Zhang   ierr = PetscTime(&t0);CHKERRQ(ierr);
5158cb82516SHong Zhang #endif
516748c7196SHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
5170d3441aeSHong Zhang     ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
5180d3441aeSHong Zhang     ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
519748c7196SHong Zhang   }
5208cb82516SHong Zhang #if defined(PTAP_PROFILE)
5210d3441aeSHong Zhang   ierr = PetscTime(&t1);CHKERRQ(ierr);
5220d3441aeSHong Zhang   eR = t1 - t0;
5238cb82516SHong Zhang #endif
5240d3441aeSHong Zhang 
525e497d3c8SHong Zhang   /* 2) get AP_loc */
5260d3441aeSHong Zhang   AP_loc = ptap->AP_loc;
5278cb82516SHong Zhang   ap = (Mat_SeqAIJ*)AP_loc->data;
5280d3441aeSHong Zhang 
529e497d3c8SHong Zhang   /* 2-1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
5300d3441aeSHong Zhang   /*-----------------------------------------------------*/
531748c7196SHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
532748c7196SHong Zhang     /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */
5330d3441aeSHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
5340d3441aeSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
535e497d3c8SHong Zhang   }
5360d3441aeSHong Zhang 
5378cb82516SHong Zhang   /* 2-2) compute numeric A_loc*P - dominating part */
5388cb82516SHong Zhang   /* ---------------------------------------------- */
5390d3441aeSHong Zhang   /* get data from symbolic products */
5408cb82516SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
5418cb82516SHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
5428cb82516SHong Zhang   apa   = ptap->apa;
543681d504bSHong Zhang   api   = ap->i;
544681d504bSHong Zhang   apj   = ap->j;
545e497d3c8SHong Zhang   for (i=0; i<am; i++) {
546445158ffSHong Zhang     /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */
547e497d3c8SHong Zhang     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
548e497d3c8SHong Zhang     apnz = api[i+1] - api[i];
549e497d3c8SHong Zhang     for (j=0; j<apnz; j++) {
550e497d3c8SHong Zhang       col = apj[j+api[i]];
551e497d3c8SHong Zhang       ap->a[j+ap->i[i]] = apa[col];
552e497d3c8SHong Zhang       apa[col] = 0.0;
5530d3441aeSHong Zhang     }
554aa690a28SHong Zhang     ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
5550d3441aeSHong Zhang   }
5568cb82516SHong Zhang #if defined(PTAP_PROFILE)
5570d3441aeSHong Zhang   ierr = PetscTime(&t2);CHKERRQ(ierr);
5580d3441aeSHong Zhang   eAP = t2 - t1;
5598cb82516SHong Zhang #endif
5600d3441aeSHong Zhang 
5618cb82516SHong Zhang   /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */
562445158ffSHong Zhang   ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr);
563445158ffSHong Zhang   ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr);
5649ce11a7cSHong Zhang   C_loc = ptap->C_loc;
5659ce11a7cSHong Zhang   C_oth = ptap->C_oth;
5668cb82516SHong Zhang #if defined(PTAP_PROFILE)
5670d3441aeSHong Zhang   ierr = PetscTime(&t3);CHKERRQ(ierr);
5680d3441aeSHong Zhang   eCseq = t3 - t2;
5698cb82516SHong Zhang #endif
5700d3441aeSHong Zhang 
5710d3441aeSHong Zhang   /* add C_loc and Co to to C */
5720d3441aeSHong Zhang   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
5730d3441aeSHong Zhang 
5740d3441aeSHong Zhang   /* C_loc -> C */
5750d3441aeSHong Zhang   cm    = C_loc->rmap->N;
5769ce11a7cSHong Zhang   c_seq = (Mat_SeqAIJ*)C_loc->data;
5778cb82516SHong Zhang   cols = c_seq->j;
5788cb82516SHong Zhang   vals = c_seq->a;
5790d3441aeSHong Zhang   for (i=0; i<cm; i++) {
5809ce11a7cSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
5810d3441aeSHong Zhang     row = rstart + i;
5820d3441aeSHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
5838cb82516SHong Zhang     cols += ncols; vals += ncols;
5840d3441aeSHong Zhang   }
5850d3441aeSHong Zhang 
5860d3441aeSHong Zhang   /* Co -> C, off-processor part */
5879ce11a7cSHong Zhang   cm = C_oth->rmap->N;
5889ce11a7cSHong Zhang   c_seq = (Mat_SeqAIJ*)C_oth->data;
5898cb82516SHong Zhang   cols = c_seq->j;
5908cb82516SHong Zhang   vals = c_seq->a;
5919ce11a7cSHong Zhang   for (i=0; i<cm; i++) {
5929ce11a7cSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
5930d3441aeSHong Zhang     row = p->garray[i];
5940d3441aeSHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
5958cb82516SHong Zhang     cols += ncols; vals += ncols;
5960d3441aeSHong Zhang   }
5970d3441aeSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5980d3441aeSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5998cb82516SHong Zhang #if defined(PTAP_PROFILE)
6000d3441aeSHong Zhang   ierr = PetscTime(&t4);CHKERRQ(ierr);
6010d3441aeSHong Zhang   eCmpi = t4 - t3;
6020d3441aeSHong Zhang 
6038cb82516SHong Zhang   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
6048cb82516SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
6050d3441aeSHong Zhang   if (rank==1) {
6060d3441aeSHong Zhang     ierr = PetscPrintf(MPI_COMM_SELF," R %g, AP %g, Cseq %g, Cmpi %g = %g\n", eR,eAP,eCseq,eCmpi,eR+eAP+eCseq+eCmpi);CHKERRQ(ierr);
6070d3441aeSHong Zhang   }
6088cb82516SHong Zhang #endif
609748c7196SHong Zhang   ptap->reuse = MAT_REUSE_MATRIX;
6100d3441aeSHong Zhang   PetscFunctionReturn(0);
6110d3441aeSHong Zhang }
6120d3441aeSHong Zhang 
6130d3441aeSHong Zhang #undef __FUNCT__
6148403a639SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap"
6158403a639SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(Mat A,Mat P,PetscReal fill,Mat *C)
61642c57489SHong Zhang {
61742c57489SHong Zhang   PetscErrorCode      ierr;
61877584fe6SHong Zhang   Mat                 Cmpi;
619f8487c73SHong Zhang   Mat_PtAPMPI         *ptap;
6200298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
621f8487c73SHong Zhang   Mat_MPIAIJ          *a        =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
62242c57489SHong Zhang   Mat_SeqAIJ          *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
62342c57489SHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth;
624ded4f561SHong Zhang   PetscInt            *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ;
62577584fe6SHong Zhang   PetscInt            *adi=ad->i,*aj,*aoi=ao->i,nnz;
626a914927fSHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
627d16ca5f0SHong Zhang   PetscInt            am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n;
62842c57489SHong Zhang   PetscBT             lnkbt;
629ce94432eSBarry Smith   MPI_Comm            comm;
630a914927fSHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0;
63142c57489SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
63242c57489SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
63324ecddacSHong Zhang   PetscInt            nzi,*pti,*ptj;
63442c57489SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
635ded4f561SHong Zhang   MPI_Request         *swaits,*rwaits;
636ded4f561SHong Zhang   MPI_Status          *sstatus,rstatus;
63742c57489SHong Zhang   Mat_Merge_SeqsToMPI *merge;
63877584fe6SHong Zhang   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0;
639d16ca5f0SHong Zhang   PetscReal           afill=1.0,afill_tmp;
640a914927fSHong Zhang   PetscInt            rmax;
64142c57489SHong Zhang 
64242c57489SHong Zhang   PetscFunctionBegin;
643ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
64483445d95SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
64583445d95SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
64683445d95SHong Zhang 
647f8487c73SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
648b00a9115SJed Brown   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
649b00a9115SJed Brown   ierr        = PetscNew(&merge);CHKERRQ(ierr);
6509f4155fbSHong Zhang   ptap->merge = merge;
651f8487c73SHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
6526c6e00e0SHong Zhang 
6536c6e00e0SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
654b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
655fe615159SHong Zhang 
6566c6e00e0SHong Zhang   /* get P_loc by taking all local rows of P */
657a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
6586c6e00e0SHong Zhang 
659a1a4e84aSHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
660a1a4e84aSHong Zhang   p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
6616c6e00e0SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
6626c6e00e0SHong Zhang   pi_oth = p_oth->i; pj_oth = p_oth->j;
66342c57489SHong Zhang 
6642addb229SHong Zhang   /* (1) compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth (api,apj) */
6652addb229SHong Zhang   /*--------------------------------------------------------------------------*/
666854ce69bSBarry Smith   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
66782412ba7SHong Zhang   api[0] = 0;
66883445d95SHong Zhang 
66983445d95SHong Zhang   /* create and initialize a linked list */
670a914927fSHong Zhang   ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr);
67183445d95SHong Zhang 
672a914927fSHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) -OOM for ex56, np=8k on Intrepid! */
673d16ca5f0SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
674f4a743e1SHong Zhang   current_space = free_space;
675f4a743e1SHong Zhang 
676f4a743e1SHong Zhang   for (i=0; i<am; i++) {
677f4a743e1SHong Zhang     /* diagonal portion of A */
678ded4f561SHong Zhang     nzi = adi[i+1] - adi[i];
67977584fe6SHong Zhang     aj  = ad->j + adi[i];
680ded4f561SHong Zhang     for (j=0; j<nzi; j++) {
68177584fe6SHong Zhang       row  = aj[j];
68282412ba7SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
683ded4f561SHong Zhang       Jptr = pj_loc + pi_loc[row];
68483445d95SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
685a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
686f4a743e1SHong Zhang     }
687f4a743e1SHong Zhang     /* off-diagonal portion of A */
688ded4f561SHong Zhang     nzi = aoi[i+1] - aoi[i];
68977584fe6SHong Zhang     aj  = ao->j + aoi[i];
690ded4f561SHong Zhang     for (j=0; j<nzi; j++) {
69177584fe6SHong Zhang       row  = aj[j];
69282412ba7SHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
693ded4f561SHong Zhang       Jptr = pj_oth + pi_oth[row];
694a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
695f4a743e1SHong Zhang     }
696a914927fSHong Zhang     apnz     = lnk[0];
69782412ba7SHong Zhang     api[i+1] = api[i] + apnz;
69877584fe6SHong Zhang     if (ap_rmax < apnz) ap_rmax = apnz;
699f4a743e1SHong Zhang 
70083445d95SHong Zhang     /* if free space is not available, double the total space in the list */
70182412ba7SHong Zhang     if (current_space->local_remaining<apnz) {
7022ba1acd5SHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
703f2b054eeSHong Zhang       nspacedouble++;
704f4a743e1SHong Zhang     }
705f4a743e1SHong Zhang 
706f4a743e1SHong Zhang     /* Copy data into free space, then initialize lnk */
707a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
7082205254eSKarl Rupp 
70982412ba7SHong Zhang     current_space->array           += apnz;
71082412ba7SHong Zhang     current_space->local_used      += apnz;
71182412ba7SHong Zhang     current_space->local_remaining -= apnz;
712f4a743e1SHong Zhang   }
713a914927fSHong Zhang 
71482412ba7SHong Zhang   /* Allocate space for apj, initialize apj, and */
715f4a743e1SHong Zhang   /* destroy list of free space and other temporary array(s) */
716854ce69bSBarry Smith   ierr      = PetscMalloc1(api[am]+1,&apj);CHKERRQ(ierr);
71777584fe6SHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
718118727c9SMark F. Adams   afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1);
719d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
720f4a743e1SHong Zhang 
7212addb229SHong Zhang   /* (2) determine symbolic Co=(p->B)^T*AP - send to others (coi,coj)*/
7222addb229SHong Zhang   /*-----------------------------------------------------------------*/
723de0260b3SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
72442c57489SHong Zhang 
725ded4f561SHong Zhang   /* then, compute symbolic Co = (p->B)^T*AP */
726d0f46423SBarry Smith   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
72783445d95SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
728854ce69bSBarry Smith   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
729de0260b3SHong Zhang   coi[0] = 0;
73042c57489SHong Zhang 
731d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */
732d16ca5f0SHong Zhang   nnz           = fill*(poti[pon] + api[am]);
73322d28d08SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
73442c57489SHong Zhang   current_space = free_space;
73542c57489SHong Zhang 
736de0260b3SHong Zhang   for (i=0; i<pon; i++) {
73782412ba7SHong Zhang     pnz = poti[i+1] - poti[i];
73877584fe6SHong Zhang     ptJ = potj + poti[i];
73977584fe6SHong Zhang     for (j=0; j<pnz; j++) {
74077584fe6SHong Zhang       row  = ptJ[j]; /* row of AP == col of Pot */
74182412ba7SHong Zhang       apnz = api[row+1] - api[row];
742ded4f561SHong Zhang       Jptr = apj + api[row];
74383445d95SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
744a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
74542c57489SHong Zhang     }
746a914927fSHong Zhang     nnz = lnk[0];
74742c57489SHong Zhang 
74883445d95SHong Zhang     /* If free space is not available, double the total space in the list */
749ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
7504238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
751d16ca5f0SHong Zhang       nspacedouble++;
75242c57489SHong Zhang     }
75342c57489SHong Zhang 
75442c57489SHong Zhang     /* Copy data into free space, and zero out denserows */
755a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
7562205254eSKarl Rupp 
757ded4f561SHong Zhang     current_space->array           += nnz;
758ded4f561SHong Zhang     current_space->local_used      += nnz;
759ded4f561SHong Zhang     current_space->local_remaining -= nnz;
7602205254eSKarl Rupp 
761ded4f561SHong Zhang     coi[i+1] = coi[i] + nnz;
76242c57489SHong Zhang   }
7636b911d16SHong Zhang 
7646b911d16SHong Zhang   ierr      = PetscMalloc1(coi[pon],&coj);CHKERRQ(ierr);
765a1a86e44SBarry Smith   ierr      = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
766118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1);
767d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
768de0260b3SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
76942c57489SHong Zhang 
7706b911d16SHong Zhang   /* (3) send j-array (coj) of Co to other processors */
7716b911d16SHong Zhang   /*--------------------------------------------------*/
7726b911d16SHong Zhang   ierr = PetscCalloc1(size,&merge->len_s);CHKERRQ(ierr);
7736b911d16SHong Zhang   len_s        = merge->len_s;
7746b911d16SHong Zhang   merge->nsend = 0;
7756b911d16SHong Zhang 
7766b911d16SHong Zhang 
77742c57489SHong Zhang   /* determine row ownership */
77826283091SBarry Smith   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
7797a2fc3feSBarry Smith   merge->rowmap->n  = pn;
7807a2fc3feSBarry Smith   merge->rowmap->bs = 1;
7812205254eSKarl Rupp 
78226283091SBarry Smith   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
7837a2fc3feSBarry Smith   owners = merge->rowmap->range;
78442c57489SHong Zhang 
78542c57489SHong Zhang   /* determine the number of messages to send, their lengths */
786dcca6d9dSJed Brown   ierr = PetscMalloc2(size,&len_si,size,&sstatus);CHKERRQ(ierr);
78783445d95SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
788854ce69bSBarry Smith   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
789de0260b3SHong Zhang 
79083445d95SHong Zhang   proc = 0;
791de0260b3SHong Zhang   for (i=0; i<pon; i++) {
792de0260b3SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
7932addb229SHong Zhang     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
794ee6e4b5bSHong Zhang     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
79583445d95SHong Zhang   }
796de0260b3SHong Zhang 
7972addb229SHong Zhang   len          = 0; /* max length of buf_si[], see (4) */
798de0260b3SHong Zhang   owners_co[0] = 0;
79942c57489SHong Zhang   for (proc=0; proc<size; proc++) {
800de0260b3SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
801ee6e4b5bSHong Zhang     if (len_s[proc]) {
80242c57489SHong Zhang       merge->nsend++;
8032addb229SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
80442c57489SHong Zhang       len         += len_si[proc];
80542c57489SHong Zhang     }
80642c57489SHong Zhang   }
80742c57489SHong Zhang 
808ded4f561SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
8090298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
81042c57489SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
81142c57489SHong Zhang 
812ded4f561SHong Zhang   /* post the Irecv and Isend of coj */
813529bc5dcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
814529bc5dcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
815854ce69bSBarry Smith   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
81642c57489SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
81742c57489SHong Zhang     if (!len_s[proc]) continue;
818de0260b3SHong Zhang     i    = owners_co[proc];
819529bc5dcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
82042c57489SHong Zhang     k++;
82142c57489SHong Zhang   }
82242c57489SHong Zhang 
823ded4f561SHong Zhang   /* receives and sends of coj are complete */
82477584fe6SHong Zhang   for (i=0; i<merge->nrecv; i++) {
825c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
826ded4f561SHong Zhang   }
827ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
8280c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
82982412ba7SHong Zhang 
8306b911d16SHong Zhang   /* (4) send and recv coi */
8316b911d16SHong Zhang   /*-----------------------*/
832529bc5dcSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
833529bc5dcSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
834854ce69bSBarry Smith   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
83542c57489SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
83642c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
83742c57489SHong Zhang     if (!len_s[proc]) continue;
83842c57489SHong Zhang     /* form outgoing message for i-structure:
83942c57489SHong Zhang          buf_si[0]:                 nrows to be sent
84042c57489SHong Zhang                [1:nrows]:           row index (global)
84142c57489SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
84242c57489SHong Zhang     */
84342c57489SHong Zhang     /*-------------------------------------------*/
8442addb229SHong Zhang     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
84542c57489SHong Zhang     buf_si_i    = buf_si + nrows+1;
84642c57489SHong Zhang     buf_si[0]   = nrows;
84742c57489SHong Zhang     buf_si_i[0] = 0;
84842c57489SHong Zhang     nrows       = 0;
849de0260b3SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
850ded4f561SHong Zhang       nzi = coi[i+1] - coi[i];
851ded4f561SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
852de0260b3SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
85342c57489SHong Zhang       nrows++;
85442c57489SHong Zhang     }
855529bc5dcSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
85642c57489SHong Zhang     k++;
85742c57489SHong Zhang     buf_si += len_si[proc];
85842c57489SHong Zhang   }
859ded4f561SHong Zhang   i = merge->nrecv;
860ded4f561SHong Zhang   while (i--) {
861c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
862ded4f561SHong Zhang   }
863ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
8640c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
865a914927fSHong Zhang 
86624ecddacSHong Zhang   ierr = PetscFree2(len_si,sstatus);CHKERRQ(ierr);
86742c57489SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
868ded4f561SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
86942c57489SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
87042c57489SHong Zhang 
8716b911d16SHong Zhang   /* (5) compute the local portion of C (mpi mat) */
8726b911d16SHong Zhang   /*----------------------------------------------*/
873ded4f561SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
874ded4f561SHong Zhang 
87524ecddacSHong Zhang   /* allocate pti array and free space for accumulating nonzero column info */
876854ce69bSBarry Smith   ierr   = PetscMalloc1(pn+1,&pti);CHKERRQ(ierr);
87724ecddacSHong Zhang   pti[0] = 0;
87842c57489SHong Zhang 
879d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
880d16ca5f0SHong Zhang   nnz           = fill*(pi_loc[pm] + api[am]);
88122d28d08SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
88242c57489SHong Zhang   current_space = free_space;
88342c57489SHong Zhang 
884dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
88542c57489SHong Zhang   for (k=0; k<merge->nrecv; k++) {
88642c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
88742c57489SHong Zhang     nrows       = *buf_ri_k[k];
88842c57489SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
88942c57489SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
89042c57489SHong Zhang   }
89142c57489SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
89208cb4532SHong Zhang   rmax = 0;
89342c57489SHong Zhang   for (i=0; i<pn; i++) {
894ded4f561SHong Zhang     /* add pdt[i,:]*AP into lnk */
895ded4f561SHong Zhang     pnz = pdti[i+1] - pdti[i];
89677584fe6SHong Zhang     ptJ = pdtj + pdti[i];
89777584fe6SHong Zhang     for (j=0; j<pnz; j++) {
89877584fe6SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
899ded4f561SHong Zhang       apnz = api[row+1] - api[row];
900ded4f561SHong Zhang       Jptr = apj + api[row];
901ded4f561SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
902a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
903ded4f561SHong Zhang     }
904a914927fSHong Zhang 
90542c57489SHong Zhang     /* add received col data into lnk */
90642c57489SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
90742c57489SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
908ded4f561SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
909ded4f561SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
910a914927fSHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
91142c57489SHong Zhang         nextrow[k]++; nextci[k]++;
91242c57489SHong Zhang       }
91342c57489SHong Zhang     }
914a914927fSHong Zhang     nnz = lnk[0];
91542c57489SHong Zhang 
91642c57489SHong Zhang     /* if free space is not available, make more free space */
917ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
9184238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
919d16ca5f0SHong Zhang       nspacedouble++;
92042c57489SHong Zhang     }
92142c57489SHong Zhang     /* copy data into free space, then initialize lnk */
922a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
923ded4f561SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
9242205254eSKarl Rupp 
925ded4f561SHong Zhang     current_space->array           += nnz;
926ded4f561SHong Zhang     current_space->local_used      += nnz;
927ded4f561SHong Zhang     current_space->local_remaining -= nnz;
9282205254eSKarl Rupp 
92924ecddacSHong Zhang     pti[i+1] = pti[i] + nnz;
93008cb4532SHong Zhang     if (nnz > rmax) rmax = nnz;
93142c57489SHong Zhang   }
932ded4f561SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
9330572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
93442c57489SHong Zhang 
935854ce69bSBarry Smith   ierr      = PetscMalloc1(pti[pn]+1,&ptj);CHKERRQ(ierr);
93624ecddacSHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,ptj);CHKERRQ(ierr);
93724ecddacSHong Zhang   afill_tmp = (PetscReal)pti[pn]/(pi_loc[pm] + api[am]+1);
938d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
93942c57489SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
94042c57489SHong Zhang 
9416b911d16SHong Zhang   /* (6) create symbolic parallel matrix Cmpi */
9426b911d16SHong Zhang   /*------------------------------------------*/
94377584fe6SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
94477584fe6SHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
94533d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
94677584fe6SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
94777584fe6SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
94842c57489SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
949a2f3521dSMark F. Adams 
950b091e509SHong Zhang   merge->bi        = pti;      /* Cseq->i */
951b091e509SHong Zhang   merge->bj        = ptj;      /* Cseq->j */
952b091e509SHong Zhang   merge->coi       = coi;      /* Co->i   */
953b091e509SHong Zhang   merge->coj       = coj;      /* Co->j   */
95442c57489SHong Zhang   merge->buf_ri    = buf_ri;
95542c57489SHong Zhang   merge->buf_rj    = buf_rj;
956de0260b3SHong Zhang   merge->owners_co = owners_co;
95777584fe6SHong Zhang   merge->destroy   = Cmpi->ops->destroy;
95842c57489SHong Zhang 
95977584fe6SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
96077584fe6SHong Zhang   c           = (Mat_MPIAIJ*)Cmpi->data;
961f8487c73SHong Zhang   c->ptap     = ptap;
96277584fe6SHong Zhang   ptap->api   = api;
96377584fe6SHong Zhang   ptap->apj   = apj;
964d6ab1dc0SHong Zhang   ptap->rmax  = ap_rmax;
9658403a639SHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
9668403a639SHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
9678403a639SHong Zhang 
9688403a639SHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
9698403a639SHong Zhang   Cmpi->assembled        = PETSC_FALSE;
9708403a639SHong Zhang   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
9718403a639SHong Zhang   Cmpi->ops->duplicate   = MatDuplicate_MPIAIJ_MatPtAP;
9728403a639SHong Zhang   Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap;
97377584fe6SHong Zhang   *C                     = Cmpi;
974414894bdSHong Zhang 
975414894bdSHong Zhang   /* flag 'scalable' determines which implementations to be used:
976414894bdSHong Zhang        0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa;
977414894bdSHong Zhang        1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */
978414894bdSHong Zhang   /* set default scalable */
979da0a95b2SSatish Balay   ptap->scalable = PETSC_FALSE; /* PETSC_TRUE; */
9802205254eSKarl Rupp 
9810298fd71SBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)Cmpi)->prefix,"-matptap_scalable",&ptap->scalable,NULL);CHKERRQ(ierr);
982414894bdSHong Zhang   if (!ptap->scalable) {  /* Do dense axpy */
9831795a4d1SJed Brown     ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
984414894bdSHong Zhang   } else {
9851795a4d1SJed Brown     ierr = PetscCalloc1(ap_rmax+1,&ptap->apa);CHKERRQ(ierr);
986414894bdSHong Zhang   }
987414894bdSHong Zhang 
988f2b054eeSHong Zhang #if defined(PETSC_USE_INFO)
98924ecddacSHong Zhang   if (pti[pn] != 0) {
99057622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
99157622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
992f2b054eeSHong Zhang   } else {
99377584fe6SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
994f2b054eeSHong Zhang   }
995f2b054eeSHong Zhang #endif
99642c57489SHong Zhang   PetscFunctionReturn(0);
99742c57489SHong Zhang }
99842c57489SHong Zhang 
99942c57489SHong Zhang #undef __FUNCT__
10008403a639SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap"
10018403a639SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap(Mat A,Mat P,Mat C)
100242c57489SHong Zhang {
100342c57489SHong Zhang   PetscErrorCode      ierr;
1004f8487c73SHong Zhang   Mat_MPIAIJ          *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
100542c57489SHong Zhang   Mat_SeqAIJ          *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
1006de0260b3SHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
100742c57489SHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth;
1008f8487c73SHong Zhang   Mat_PtAPMPI         *ptap;
10099f4155fbSHong Zhang   Mat_Merge_SeqsToMPI *merge;
10101d633602SHong Zhang   PetscInt            *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp;
101182412ba7SHong Zhang   PetscInt            *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj;
101282412ba7SHong Zhang   PetscInt            i,j,k,anz,pnz,apnz,nextap,row,*cj;
1013e900a757SHong Zhang   MatScalar           *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp;
1014d0f46423SBarry Smith   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1015ce94432eSBarry Smith   MPI_Comm            comm;
101642c57489SHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
101782412ba7SHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
101842c57489SHong Zhang   PetscInt            **buf_ri,**buf_rj;
101950f5bf86SHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
102042c57489SHong Zhang   MPI_Request         *s_waits,*r_waits;
102142c57489SHong Zhang   MPI_Status          *status;
102282412ba7SHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
102382412ba7SHong Zhang   PetscInt            *api,*apj,*coi,*coj;
1024d0f46423SBarry Smith   PetscInt            *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend;
10256ce94e8fSHong Zhang   PetscBool           scalable;
102638571addSHong Zhang #if defined(PTAP_PROFILE)
1027e497d3c8SHong Zhang   PetscLogDouble t0,t1,t2,eP,t3,t4,et2_AP=0.0,ePtAP=0.0,t2_0,t2_1,t2_2;
102838571addSHong Zhang #endif
102942c57489SHong Zhang 
103042c57489SHong Zhang   PetscFunctionBegin;
1031ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
103238571addSHong Zhang #if defined(PTAP_PROFILE)
10338563dfccSBarry Smith   ierr = PetscTime(&t0);CHKERRQ(ierr);
103438571addSHong Zhang #endif
103542c57489SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
103642c57489SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
103742c57489SHong Zhang 
1038f8487c73SHong Zhang   ptap = c->ptap;
1039ce94432eSBarry Smith   if (!ptap) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_INCOMP,"MatPtAP() has not been called to create matrix C yet, cannot use MAT_REUSE_MATRIX");
1040f8487c73SHong Zhang   merge    = ptap->merge;
1041414894bdSHong Zhang   apa      = ptap->apa;
10426ce94e8fSHong Zhang   scalable = ptap->scalable;
10436c6e00e0SHong Zhang 
1044a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
10456b911d16SHong Zhang   /*-----------------------------------------------------*/
1046f8487c73SHong Zhang   if (ptap->reuse == MAT_INITIAL_MATRIX) {
10479f4155fbSHong Zhang     /* P_oth and P_loc are obtained in MatPtASymbolic(), skip calling MatGetBrowsOfAoCols() and MatMPIAIJGetLocalMat() */
1048f8487c73SHong Zhang     ptap->reuse = MAT_REUSE_MATRIX;
10496c6e00e0SHong Zhang   } else { /* update numerical values of P_oth and P_loc */
1050b7f45c76SHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
1051a1a4e84aSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
10526c6e00e0SHong Zhang   }
105338571addSHong Zhang #if defined(PTAP_PROFILE)
10548563dfccSBarry Smith   ierr = PetscTime(&t1);CHKERRQ(ierr);
105505d62848SHong Zhang   eP = t1-t0;
105638571addSHong Zhang #endif
105705d62848SHong Zhang   /*
105805d62848SHong Zhang   printf("[%d] Ad: %d, %d; Ao: %d, %d; P_loc: %d, %d; P_oth %d, %d;\n",rank,
105905d62848SHong Zhang          a->A->rmap->N,a->A->cmap->N,a->B->rmap->N,a->B->cmap->N,
106005d62848SHong Zhang          ptap->P_loc->rmap->N,ptap->P_loc->cmap->N,
106105d62848SHong Zhang          ptap->P_oth->rmap->N,ptap->P_oth->cmap->N);
106205d62848SHong Zhang    */
1063f8487c73SHong Zhang 
1064f9473cf7SHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
1065f9473cf7SHong Zhang   /*--------------------------------------------------------------*/
106642c57489SHong Zhang   /* get data from symbolic products */
1067a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
1068a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
1069a9d06231SHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a;
107042c57489SHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
107142c57489SHong Zhang 
1072de0260b3SHong Zhang   coi  = merge->coi; coj = merge->coj;
10731795a4d1SJed Brown   ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
1074de0260b3SHong Zhang 
1075de0260b3SHong Zhang   bi     = merge->bi; bj = merge->bj;
10767a2fc3feSBarry Smith   owners = merge->rowmap->range;
10771795a4d1SJed Brown   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);  /* ba: Cseq->a */
107842c57489SHong Zhang 
1079a1a4e84aSHong Zhang   api = ptap->api; apj = ptap->apj;
1080d9824c15SHong Zhang 
10819516a85cSHong Zhang   if (!scalable) { /* Do dense axpy on apa (length of pN, stores A[i,:]*P) - nonscalable, but faster (could take 1/3 scalable time) */
1082b38d1a2dSJed Brown     ierr = PetscInfo(C,"Using non-scalable dense axpy\n");CHKERRQ(ierr);
10838cb82516SHong Zhang #if defined(PTAP_PROFILE)
108405d62848SHong Zhang     ierr = PetscTime(&t1);CHKERRQ(ierr);
10858cb82516SHong Zhang #endif
1086a9d06231SHong Zhang     for (i=0; i<am; i++) {
1087b091e509SHong Zhang #if defined(PTAP_PROFILE)
10888563dfccSBarry Smith       ierr = PetscTime(&t2_0);CHKERRQ(ierr);
1089b091e509SHong Zhang #endif
1090a9d06231SHong Zhang       /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */
1091a9d06231SHong Zhang       /*------------------------------------------------------------*/
1092a9d06231SHong Zhang       apJ = apj + api[i];
1093a9d06231SHong Zhang 
1094a9d06231SHong Zhang       /* diagonal portion of A */
1095a9d06231SHong Zhang       anz = adi[i+1] - adi[i];
1096a9d06231SHong Zhang       adj = ad->j + adi[i];
1097a9d06231SHong Zhang       ada = ad->a + adi[i];
1098a9d06231SHong Zhang       for (j=0; j<anz; j++) {
1099a9d06231SHong Zhang         row = adj[j];
1100a9d06231SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
1101a9d06231SHong Zhang         pj  = pj_loc + pi_loc[row];
1102a9d06231SHong Zhang         pa  = pa_loc + pi_loc[row];
1103a9d06231SHong Zhang 
1104a9d06231SHong Zhang         /* perform dense axpy */
1105e900a757SHong Zhang         valtmp = ada[j];
1106a9d06231SHong Zhang         for (k=0; k<pnz; k++) {
1107e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
1108a9d06231SHong Zhang         }
1109a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
1110a9d06231SHong Zhang       }
1111a9d06231SHong Zhang 
1112a9d06231SHong Zhang       /* off-diagonal portion of A */
1113a9d06231SHong Zhang       anz = aoi[i+1] - aoi[i];
1114a9d06231SHong Zhang       aoj = ao->j + aoi[i];
1115a9d06231SHong Zhang       aoa = ao->a + aoi[i];
1116a9d06231SHong Zhang       for (j=0; j<anz; j++) {
1117a9d06231SHong Zhang         row = aoj[j];
1118a9d06231SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
1119a9d06231SHong Zhang         pj  = pj_oth + pi_oth[row];
1120a9d06231SHong Zhang         pa  = pa_oth + pi_oth[row];
1121a9d06231SHong Zhang 
1122a9d06231SHong Zhang         /* perform dense axpy */
1123e900a757SHong Zhang         valtmp = aoa[j];
1124a9d06231SHong Zhang         for (k=0; k<pnz; k++) {
1125e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
1126a9d06231SHong Zhang         }
1127a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
1128a9d06231SHong Zhang       }
1129b091e509SHong Zhang #if defined(PTAP_PROFILE)
11308563dfccSBarry Smith       ierr    = PetscTime(&t2_1);CHKERRQ(ierr);
1131b091e509SHong Zhang       et2_AP += t2_1 - t2_0;
1132b091e509SHong Zhang #endif
1133a9d06231SHong Zhang 
1134a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
1135a9d06231SHong Zhang       /*--------------------------------------------------------------*/
1136a9d06231SHong Zhang       apnz = api[i+1] - api[i];
1137a9d06231SHong Zhang       /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */
1138a9d06231SHong Zhang       pnz = po->i[i+1] - po->i[i];
1139e900a757SHong Zhang       poJ = po->j + po->i[i];
1140a9d06231SHong Zhang       pA  = po->a + po->i[i];
1141a9d06231SHong Zhang       for (j=0; j<pnz; j++) {
1142e900a757SHong Zhang         row = poJ[j];
1143e900a757SHong Zhang         cnz = coi[row+1] - coi[row];
1144e900a757SHong Zhang         cj  = coj + coi[row];
1145e900a757SHong Zhang         ca  = coa + coi[row];
1146a9d06231SHong Zhang         /* perform dense axpy */
1147e900a757SHong Zhang         valtmp = pA[j];
1148a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
1149e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
1150a9d06231SHong Zhang         }
1151a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1152a9d06231SHong Zhang       }
1153a9d06231SHong Zhang       /* put the value into Cd (diagonal part) */
1154a9d06231SHong Zhang       pnz = pd->i[i+1] - pd->i[i];
1155e900a757SHong Zhang       pdJ = pd->j + pd->i[i];
1156a9d06231SHong Zhang       pA  = pd->a + pd->i[i];
1157a9d06231SHong Zhang       for (j=0; j<pnz; j++) {
1158e900a757SHong Zhang         row = pdJ[j];
1159e900a757SHong Zhang         cnz = bi[row+1] - bi[row];
1160e900a757SHong Zhang         cj  = bj + bi[row];
1161e900a757SHong Zhang         ca  = ba + bi[row];
1162a9d06231SHong Zhang         /* perform dense axpy */
1163e900a757SHong Zhang         valtmp = pA[j];
1164a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
1165e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
1166a9d06231SHong Zhang         }
1167a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1168a9d06231SHong Zhang       }
11698403a639SHong Zhang 
1170a9d06231SHong Zhang       /* zero the current row of A*P */
1171a9d06231SHong Zhang       for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0;
1172b091e509SHong Zhang #if defined(PTAP_PROFILE)
11738563dfccSBarry Smith       ierr      = PetscTime(&t2_2);CHKERRQ(ierr);
117405d62848SHong Zhang       ePtAP += t2_2 - t2_1;
1175b091e509SHong Zhang #endif
1176a9d06231SHong Zhang     }
117738571addSHong Zhang   } else { /* Do sparse axpy on apa (length of ap_rmax, stores A[i,:]*P) - scalable, but slower */
1178b38d1a2dSJed Brown     ierr = PetscInfo(C,"Using scalable sparse axpy\n");CHKERRQ(ierr);
117938571addSHong Zhang     /*-----------------------------------------------------------------------------------------*/
1180a9d06231SHong Zhang     pA=pa_loc;
118142c57489SHong Zhang     for (i=0; i<am; i++) {
1182b091e509SHong Zhang #if defined(PTAP_PROFILE)
11838563dfccSBarry Smith       ierr = PetscTime(&t2_0);CHKERRQ(ierr);
1184b091e509SHong Zhang #endif
1185f4a743e1SHong Zhang       /* form i-th sparse row of A*P */
118682412ba7SHong Zhang       apnz = api[i+1] - api[i];
118782412ba7SHong Zhang       apJ  = apj + api[i];
118842c57489SHong Zhang       /* diagonal portion of A */
118982412ba7SHong Zhang       anz = adi[i+1] - adi[i];
11901d633602SHong Zhang       adj = ad->j + adi[i];
11911d633602SHong Zhang       ada = ad->a + adi[i];
119282412ba7SHong Zhang       for (j=0; j<anz; j++) {
11931d633602SHong Zhang         row    = adj[j];
119482412ba7SHong Zhang         pnz    = pi_loc[row+1] - pi_loc[row];
119582412ba7SHong Zhang         pj     = pj_loc + pi_loc[row];
119682412ba7SHong Zhang         pa     = pa_loc + pi_loc[row];
1197e900a757SHong Zhang         valtmp = ada[j];
1198f4a743e1SHong Zhang         nextp  = 0;
119982412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
120082412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
1201e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
120242c57489SHong Zhang           }
120342c57489SHong Zhang         }
1204dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
120542c57489SHong Zhang       }
120642c57489SHong Zhang       /* off-diagonal portion of A */
120782412ba7SHong Zhang       anz = aoi[i+1] - aoi[i];
12081d633602SHong Zhang       aoj = ao->j + aoi[i];
12091d633602SHong Zhang       aoa = ao->a + aoi[i];
121082412ba7SHong Zhang       for (j=0; j<anz; j++) {
12111d633602SHong Zhang         row    = aoj[j];
121282412ba7SHong Zhang         pnz    = pi_oth[row+1] - pi_oth[row];
121382412ba7SHong Zhang         pj     = pj_oth + pi_oth[row];
121482412ba7SHong Zhang         pa     = pa_oth + pi_oth[row];
1215e900a757SHong Zhang         valtmp = aoa[j];
1216f4a743e1SHong Zhang         nextp  = 0;
121782412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
121882412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
1219e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
122042c57489SHong Zhang           }
122142c57489SHong Zhang         }
1222dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
122342c57489SHong Zhang       }
1224b091e509SHong Zhang #if defined(PTAP_PROFILE)
12258563dfccSBarry Smith       ierr    = PetscTime(&t2_1);CHKERRQ(ierr);
1226b091e509SHong Zhang       et2_AP += t2_1 - t2_0;
1227b091e509SHong Zhang #endif
122842c57489SHong Zhang 
1229a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
1230a9d06231SHong Zhang       /*--------------------------------------------------------------*/
123182412ba7SHong Zhang       pnz = pi_loc[i+1] - pi_loc[i];
1232f9473cf7SHong Zhang       pJ  = pj_loc + pi_loc[i];
123382412ba7SHong Zhang       for (j=0; j<pnz; j++) {
123442c57489SHong Zhang         nextap = 0;
1235f9473cf7SHong Zhang         row    = pJ[j]; /* global index */
123682412ba7SHong Zhang         if (row < pcstart || row >=pcend) { /* put the value into Co */
1237e900a757SHong Zhang           row = *poJ;
1238e900a757SHong Zhang           cj  = coj + coi[row];
1239e900a757SHong Zhang           ca  = coa + coi[row]; poJ++;
124082412ba7SHong Zhang         } else {                            /* put the value into Cd */
1241e900a757SHong Zhang           row = *pdJ;
1242e900a757SHong Zhang           cj  = bj + bi[row];
1243e900a757SHong Zhang           ca  = ba + bi[row]; pdJ++;
124442c57489SHong Zhang         }
1245e900a757SHong Zhang         valtmp = pA[j];
124682412ba7SHong Zhang         for (k=0; nextap<apnz; k++) {
1247e900a757SHong Zhang           if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++];
124842c57489SHong Zhang         }
1249dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
1250de0260b3SHong Zhang       }
12510496f32cSHong Zhang       pA += pnz;
125242c57489SHong Zhang       /* zero the current row info for A*P */
125382412ba7SHong Zhang       ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr);
1254b091e509SHong Zhang #if defined(PTAP_PROFILE)
12558563dfccSBarry Smith       ierr      = PetscTime(&t2_2);CHKERRQ(ierr);
125605d62848SHong Zhang       ePtAP += t2_2 - t2_1;
1257b091e509SHong Zhang #endif
125842c57489SHong Zhang     }
1259d9824c15SHong Zhang   }
126038571addSHong Zhang #if defined(PTAP_PROFILE)
12618563dfccSBarry Smith   ierr = PetscTime(&t2);CHKERRQ(ierr);
126238571addSHong Zhang #endif
126342c57489SHong Zhang 
1264a9d06231SHong Zhang   /* 3) send and recv matrix values coa */
1265a9d06231SHong Zhang   /*------------------------------------*/
126642c57489SHong Zhang   buf_ri = merge->buf_ri;
126742c57489SHong Zhang   buf_rj = merge->buf_rj;
126842c57489SHong Zhang   len_s  = merge->len_s;
1269fc42d0c8SSatish Balay   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
127042c57489SHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
127142c57489SHong Zhang 
1272dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
127342c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
127442c57489SHong Zhang     if (!len_s[proc]) continue;
1275de0260b3SHong Zhang     i    = merge->owners_co[proc];
1276de0260b3SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
127742c57489SHong Zhang     k++;
127842c57489SHong Zhang   }
12790c468ba9SBarry Smith   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
12800c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
128142c57489SHong Zhang 
1282a9d06231SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
128342c57489SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
128482412ba7SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
128538571addSHong Zhang #if defined(PTAP_PROFILE)
12868563dfccSBarry Smith   ierr = PetscTime(&t3);CHKERRQ(ierr);
128738571addSHong Zhang #endif
128842c57489SHong Zhang 
1289a9d06231SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1290a9d06231SHong Zhang   /*------------------------------------------------------*/
1291dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
129242c57489SHong Zhang   for (k=0; k<merge->nrecv; k++) {
129342c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
129442c57489SHong Zhang     nrows       = *(buf_ri_k[k]);
129542c57489SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
129682412ba7SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
129742c57489SHong Zhang   }
129842c57489SHong Zhang 
1299de0260b3SHong Zhang   for (i=0; i<cm; i++) {
130082412ba7SHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
130142c57489SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1302de0260b3SHong Zhang     ba_i = ba + bi[i];
130382412ba7SHong Zhang     bnz  = bi[i+1] - bi[i];
130442c57489SHong Zhang     /* add received vals into ba */
130542c57489SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
130642c57489SHong Zhang       /* i-th row */
130742c57489SHong Zhang       if (i == *nextrow[k]) {
130882412ba7SHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
130982412ba7SHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
131082412ba7SHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
131182412ba7SHong Zhang         nextcj = 0;
131282412ba7SHong Zhang         for (j=0; nextcj<cnz; j++) {
131382412ba7SHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
131482412ba7SHong Zhang             ba_i[j] += ca[nextcj++];
131542c57489SHong Zhang           }
131642c57489SHong Zhang         }
131782412ba7SHong Zhang         nextrow[k]++; nextci[k]++;
1318c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
131942c57489SHong Zhang       }
132042c57489SHong Zhang     }
132182412ba7SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
132242c57489SHong Zhang   }
132342c57489SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
132442c57489SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
132542c57489SHong Zhang 
1326de0260b3SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1327c05d87d6SBarry Smith   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
132842c57489SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
13290572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
133038571addSHong Zhang #if defined(PTAP_PROFILE)
13318563dfccSBarry Smith   ierr = PetscTime(&t4);CHKERRQ(ierr);
13329516a85cSHong Zhang   if (rank==1) {
133305d62848SHong Zhang     ierr = PetscPrintf(MPI_COMM_SELF,"  [%d] PtAPNum %g/P + %g/PtAP( %g/A*P + %g/Pt*AP ) + %g/comm + %g/Cloc = %g\n\n",rank,eP,t2-t1,et2_AP,ePtAP,t3-t2,t4-t3,t4-t0);CHKERRQ(ierr);
13349516a85cSHong Zhang   }
133538571addSHong Zhang #endif
133642c57489SHong Zhang   PetscFunctionReturn(0);
133742c57489SHong Zhang }
1338