xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision 9a6dcab793401456f60aaffcc860c1d95018a977)
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);
9667a12041SHong Zhang   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) {
9767a12041SHong Zhang     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);
9867a12041SHong Zhang   }
9967a12041SHong Zhang   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) {
10067a12041SHong Zhang     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);
10167a12041SHong Zhang   }
10267a12041SHong Zhang 
1038403a639SHong Zhang   ierr = PetscOptionsGetBool(NULL,"-matrap",&rap,NULL);CHKERRQ(ierr);
104cf3ca8ceSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
1053ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
1068403a639SHong Zhang     if (rap) { /* do R=P^T locally, then C=R*A*P */
107cf3ca8ceSHong Zhang       ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);
1088403a639SHong Zhang     } else {       /* do P^T*A*P */
1098403a639SHong Zhang       ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(A,P,fill,C);CHKERRQ(ierr);
1100d3441aeSHong Zhang     }
1113ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
1127a7894deSKris Buschelman   }
1133ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
1148403a639SHong Zhang   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
1153ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
11642c57489SHong Zhang   PetscFunctionReturn(0);
11742c57489SHong Zhang }
11842c57489SHong Zhang 
1198cb82516SHong Zhang /* requires array of size pN, thus nonscalable version */
12042c57489SHong Zhang #undef __FUNCT__
121aa690a28SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ"
122aa690a28SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
1230d3441aeSHong Zhang {
1240d3441aeSHong Zhang   PetscErrorCode      ierr;
1250d3441aeSHong Zhang   Mat_PtAPMPI         *ptap;
126681d504bSHong Zhang   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
1272259aa2eSHong Zhang   MPI_Comm            comm;
1282259aa2eSHong Zhang   PetscMPIInt         size,rank;
12915a3b8e2SHong Zhang   Mat                 Cmpi;
13015a3b8e2SHong Zhang   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
13167a12041SHong Zhang   PetscInt            am=A->rmap->n,pm=P->rmap->n;
13267a12041SHong Zhang   PetscInt            *lnk,i,k,pnz,row,nnz;
133445158ffSHong Zhang   PetscInt            pN=P->cmap->N,pn=P->cmap->n;
13415a3b8e2SHong Zhang   PetscBT             lnkbt;
135f671be37SHong Zhang   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
13615a3b8e2SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
13715a3b8e2SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
13867a12041SHong Zhang   PetscInt            nzi,nspacedouble;
13915a3b8e2SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
14015a3b8e2SHong Zhang   MPI_Request         *swaits,*rwaits;
14115a3b8e2SHong Zhang   MPI_Status          *sstatus,rstatus;
142445158ffSHong Zhang   PetscLayout         rowmap;
143445158ffSHong Zhang   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
144f671be37SHong Zhang   PetscInt            nsend;
145445158ffSHong Zhang   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
146*9a6dcab7SHong Zhang   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0,Prmax,APrmax;
14767a12041SHong Zhang   PetscInt            rmax,*aj,*ai,*pi;
14867a12041SHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*c_loc,*c_oth;
14967a12041SHong Zhang   PetscScalar         *apv;
150*9a6dcab7SHong Zhang   PetscTable          ta,Cta;
151aa690a28SHong Zhang #if defined(PETSC_USE_INFO)
1528cb82516SHong Zhang   PetscReal           apfill;
153aa690a28SHong Zhang #endif
154681d504bSHong Zhang #if defined(PTAP_PROFILE)
155681d504bSHong Zhang   PetscLogDouble      t0,t1,t11,t12,t2,t3,t4;
156681d504bSHong Zhang #endif
15767a12041SHong Zhang 
15867a12041SHong Zhang   PetscFunctionBegin;
15967a12041SHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
16067a12041SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
16167a12041SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1628cb82516SHong Zhang #if defined(PTAP_PROFILE)
16380bb4639SHong Zhang   ierr = PetscTime(&t0);CHKERRQ(ierr);
1648cb82516SHong Zhang #endif
1658cb82516SHong Zhang 
16615a3b8e2SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
16715a3b8e2SHong Zhang   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
16815a3b8e2SHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
16915a3b8e2SHong Zhang 
17015a3b8e2SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
17115a3b8e2SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
17215a3b8e2SHong Zhang   /* get P_loc by taking all local rows of P */
17315a3b8e2SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
17415a3b8e2SHong Zhang 
17567a12041SHong Zhang   /* (0) compute Rd = Pd^T, Ro = Po^T  */
17667a12041SHong Zhang   /* --------------------------------- */
177de817e96SHong Zhang   ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
178de817e96SHong Zhang   ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
1798cb82516SHong Zhang #if defined(PTAP_PROFILE)
180445158ffSHong Zhang   ierr = PetscTime(&t11);CHKERRQ(ierr);
1818cb82516SHong Zhang #endif
18215a3b8e2SHong Zhang 
18367a12041SHong Zhang   /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */
18467a12041SHong Zhang   /* ----------------------------------------------------------------- */
18567a12041SHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
18667a12041SHong Zhang   p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
187445158ffSHong Zhang 
188*9a6dcab7SHong Zhang   /* create and initialize a linked list */
189*9a6dcab7SHong Zhang   Prmax = p_loc->rmax+p_oth->rmax;
190*9a6dcab7SHong Zhang #if 0
191*9a6dcab7SHong Zhang   printf("[%d] p_loc->rmax %d + p_oth->rmax %d = %d\n",rank,p_loc->rmax, p_oth->rmax,Prmax);
192*9a6dcab7SHong Zhang #endif
193*9a6dcab7SHong Zhang   ierr = PetscTableCreate(2*Prmax,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc */
194*9a6dcab7SHong Zhang   ierr = PetscTableCreate(2*Prmax,pN,&Cta);CHKERRQ(ierr);/* for compute Cmpi */
195*9a6dcab7SHong Zhang   for (row=0; row<ptap->P_loc->rmap->N; row++) {
196*9a6dcab7SHong Zhang     nzi = p_loc->i[row+1] - p_loc->i[row];
197*9a6dcab7SHong Zhang     for (j=0; j<nzi; j++) {
198*9a6dcab7SHong Zhang       Jptr = j + p_loc->j + p_loc->i[row];
199*9a6dcab7SHong Zhang       ierr = PetscTableAdd(ta,*Jptr+1,1,INSERT_VALUES);CHKERRQ(ierr); /* key must be >0 */
200*9a6dcab7SHong Zhang     }
201*9a6dcab7SHong Zhang   }
202*9a6dcab7SHong Zhang   for (row=0; row<ptap->P_oth->rmap->N; row++) {
203*9a6dcab7SHong Zhang     nzi = p_oth->i[row+1] - p_oth->i[row];
204*9a6dcab7SHong Zhang     for (j=0; j<nzi; j++) {
205*9a6dcab7SHong Zhang       Jptr = j + p_oth->j + p_oth->i[row];
206*9a6dcab7SHong Zhang       ierr = PetscTableAdd(ta,*Jptr+1,1,INSERT_VALUES);CHKERRQ(ierr); /* key must be >0 */
207*9a6dcab7SHong Zhang     }
208*9a6dcab7SHong Zhang   }
209*9a6dcab7SHong Zhang   ierr = PetscTableGetCount(ta,&Prmax);CHKERRQ(ierr); /* Prmax = nnz(sum of Prows) */
210*9a6dcab7SHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
211*9a6dcab7SHong Zhang #if 0
212*9a6dcab7SHong Zhang   //printf("[%d] Prmax %d, pN %d\n",rank,Prmax,pN);
213*9a6dcab7SHong Zhang   printf("[%d] p_loc->rmax %d + p_oth->rmax %d = %d; Prmax %d\n",rank,p_loc->rmax, p_oth->rmax,p_loc->rmax+p_oth->rmax,Prmax);
214*9a6dcab7SHong Zhang #endif
215*9a6dcab7SHong Zhang   ierr = PetscLLCondensedCreate(Prmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
216445158ffSHong Zhang 
2178cb82516SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
21867a12041SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ad->i[am]+ao->i[am]+p_loc->i[pm])),&free_space);CHKERRQ(ierr);
219445158ffSHong Zhang   current_space = free_space;
22067a12041SHong Zhang   nspacedouble  = 0;
22167a12041SHong Zhang 
22267a12041SHong Zhang   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
22367a12041SHong Zhang   api[0] = 0;
224445158ffSHong Zhang   for (i=0; i<am; i++) {
22567a12041SHong Zhang     /* diagonal portion: Ad[i,:]*P */
22667a12041SHong Zhang     ai = ad->i; pi = p_loc->i;
22767a12041SHong Zhang     nzi = ai[i+1] - ai[i];
22867a12041SHong Zhang     aj  = ad->j + ai[i];
229445158ffSHong Zhang     for (j=0; j<nzi; j++) {
230445158ffSHong Zhang       row  = aj[j];
23167a12041SHong Zhang       pnz  = pi[row+1] - pi[row];
23267a12041SHong Zhang       Jptr = p_loc->j + pi[row];
233445158ffSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
234445158ffSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
235445158ffSHong Zhang     }
23667a12041SHong Zhang     /* off-diagonal portion: Ao[i,:]*P */
23767a12041SHong Zhang     ai = ao->i; pi = p_oth->i;
23867a12041SHong Zhang     nzi = ai[i+1] - ai[i];
23967a12041SHong Zhang     aj  = ao->j + ai[i];
240445158ffSHong Zhang     for (j=0; j<nzi; j++) {
241445158ffSHong Zhang       row  = aj[j];
24267a12041SHong Zhang       pnz  = pi[row+1] - pi[row];
24367a12041SHong Zhang       Jptr = p_oth->j + pi[row];
244445158ffSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
245445158ffSHong Zhang     }
246445158ffSHong Zhang     apnz     = lnk[0];
247445158ffSHong Zhang     api[i+1] = api[i] + apnz;
248445158ffSHong Zhang     if (ap_rmax < apnz) ap_rmax = apnz;
249445158ffSHong Zhang 
250445158ffSHong Zhang     /* if free space is not available, double the total space in the list */
251445158ffSHong Zhang     if (current_space->local_remaining<apnz) {
252445158ffSHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
253445158ffSHong Zhang       nspacedouble++;
254445158ffSHong Zhang     }
255445158ffSHong Zhang 
256445158ffSHong Zhang     /* Copy data into free space, then initialize lnk */
257445158ffSHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
258445158ffSHong Zhang 
259*9a6dcab7SHong Zhang     /* add column indices into Cta */
260*9a6dcab7SHong Zhang     for (j=0; j<apnz; j++) {
261*9a6dcab7SHong Zhang       ierr = PetscTableAdd(Cta,current_space->array[j]+1,1,INSERT_VALUES);CHKERRQ(ierr); /* key must be >0 */
262*9a6dcab7SHong Zhang     }
263*9a6dcab7SHong Zhang 
264445158ffSHong Zhang     current_space->array           += apnz;
265445158ffSHong Zhang     current_space->local_used      += apnz;
266445158ffSHong Zhang     current_space->local_remaining -= apnz;
267445158ffSHong Zhang   }
268681d504bSHong Zhang   /* Allocate space for apj and apv, initialize apj, and */
269445158ffSHong Zhang   /* destroy list of free space and other temporary array(s) */
270445158ffSHong Zhang   ierr   = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr);
271445158ffSHong Zhang   ierr   = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
272*9a6dcab7SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
273*9a6dcab7SHong Zhang 
274*9a6dcab7SHong Zhang   ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr); /* rm later! */
275445158ffSHong Zhang 
276aa690a28SHong Zhang   /* Create AP_loc for reuse */
277445158ffSHong Zhang   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr);
278*9a6dcab7SHong Zhang #if 0
279*9a6dcab7SHong Zhang   Mat_SeqAIJ *ap_loc       = (Mat_SeqAIJ*)(ptap->AP_loc)->data;
280*9a6dcab7SHong Zhang   if (rank==1) printf("[%d] pN %d, AP_loc->rmax %d\n",rank,pN,ap_loc->rmax);
281*9a6dcab7SHong Zhang   if (ap_loc->rmax == 0) {
282*9a6dcab7SHong Zhang     printf("[%d] ap_loc->rmax == 0, empty AP_loc?, Prmax %d\n",rank,Prmax);
283*9a6dcab7SHong Zhang   }
284*9a6dcab7SHong Zhang #endif
285aa690a28SHong Zhang 
286aa690a28SHong Zhang #if defined(PETSC_USE_INFO)
287aa690a28SHong Zhang   apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1);
288aa690a28SHong Zhang   ptap->AP_loc->info.mallocs           = nspacedouble;
289aa690a28SHong Zhang   ptap->AP_loc->info.fill_ratio_given  = fill;
290aa690a28SHong Zhang   ptap->AP_loc->info.fill_ratio_needed = apfill;
291aa690a28SHong Zhang 
292aa690a28SHong Zhang   if (api[am]) {
293aa690a28SHong Zhang     ierr = PetscInfo3(ptap->AP_loc,"AP_loc reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)apfill);CHKERRQ(ierr);
294aa690a28SHong Zhang     ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr);
295aa690a28SHong Zhang   } else {
296aa690a28SHong Zhang     ierr = PetscInfo(ptap->AP_loc,"AP_loc is empty \n");CHKERRQ(ierr);
297aa690a28SHong Zhang   }
298aa690a28SHong Zhang #endif
299aa690a28SHong Zhang 
3008cb82516SHong Zhang #if defined(PTAP_PROFILE)
301445158ffSHong Zhang   ierr = PetscTime(&t12);CHKERRQ(ierr);
3028cb82516SHong Zhang #endif
30315a3b8e2SHong Zhang 
304681d504bSHong Zhang   /* (2-1) compute symbolic Co = Ro*AP_loc  */
305681d504bSHong Zhang   /* ------------------------------------ */
30667a12041SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
30767a12041SHong Zhang   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
3088cb82516SHong Zhang #if defined(PTAP_PROFILE)
30980bb4639SHong Zhang   ierr = PetscTime(&t1);CHKERRQ(ierr);
3108cb82516SHong Zhang #endif
31115a3b8e2SHong Zhang 
31267a12041SHong Zhang   /* (3) send coj of C_oth to other processors  */
31367a12041SHong Zhang   /* ------------------------------------------ */
31415a3b8e2SHong Zhang   /* determine row ownership */
315445158ffSHong Zhang   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
316445158ffSHong Zhang   rowmap->n  = pn;
317445158ffSHong Zhang   rowmap->bs = 1;
318445158ffSHong Zhang   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
319445158ffSHong Zhang   owners = rowmap->range;
32015a3b8e2SHong Zhang 
32115a3b8e2SHong Zhang   /* determine the number of messages to send, their lengths */
3228cb82516SHong Zhang   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
3238cb82516SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
32415a3b8e2SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
32515a3b8e2SHong Zhang 
32667a12041SHong Zhang   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
32767a12041SHong Zhang   coi   = c_oth->i; coj = c_oth->j;
32867a12041SHong Zhang   con   = ptap->C_oth->rmap->n;
32915a3b8e2SHong Zhang   proc  = 0;
33067a12041SHong Zhang   for (i=0; i<con; i++) {
33115a3b8e2SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
33215a3b8e2SHong Zhang     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
33315a3b8e2SHong Zhang     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
33415a3b8e2SHong Zhang   }
33515a3b8e2SHong Zhang 
33615a3b8e2SHong Zhang   len          = 0; /* max length of buf_si[], see (4) */
33715a3b8e2SHong Zhang   owners_co[0] = 0;
33867a12041SHong Zhang   nsend        = 0;
33915a3b8e2SHong Zhang   for (proc=0; proc<size; proc++) {
34015a3b8e2SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
34115a3b8e2SHong Zhang     if (len_s[proc]) {
342445158ffSHong Zhang       nsend++;
34315a3b8e2SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
34415a3b8e2SHong Zhang       len         += len_si[proc];
34515a3b8e2SHong Zhang     }
34615a3b8e2SHong Zhang   }
34715a3b8e2SHong Zhang 
34815a3b8e2SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
349445158ffSHong Zhang   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
350445158ffSHong Zhang   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
35115a3b8e2SHong Zhang 
35215a3b8e2SHong Zhang   /* post the Irecv and Isend of coj */
35315a3b8e2SHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
354445158ffSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
355445158ffSHong Zhang   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
35615a3b8e2SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
35715a3b8e2SHong Zhang     if (!len_s[proc]) continue;
35815a3b8e2SHong Zhang     i    = owners_co[proc];
35915a3b8e2SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
36015a3b8e2SHong Zhang     k++;
36115a3b8e2SHong Zhang   }
36215a3b8e2SHong Zhang 
363681d504bSHong Zhang   /* (2-2) compute symbolic C_loc = Rd*AP_loc */
364681d504bSHong Zhang   /* ---------------------------------------- */
365681d504bSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr);
366681d504bSHong Zhang   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
367681d504bSHong Zhang 
368681d504bSHong Zhang   /* receives coj are complete */
369445158ffSHong Zhang   for (i=0; i<nrecv; i++) {
370445158ffSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
37115a3b8e2SHong Zhang   }
37215a3b8e2SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
373445158ffSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
37415a3b8e2SHong Zhang 
375*9a6dcab7SHong Zhang   ierr = PetscTableGetCount(Cta,&APrmax);CHKERRQ(ierr); /* Prmax = nnz(sum of Prows) */
376*9a6dcab7SHong Zhang   /* printf("[%d] APrmax %d, pN %d\n",rank,APrmax,pN); */
377*9a6dcab7SHong Zhang 
37815a3b8e2SHong Zhang   /* (4) send and recv coi */
37915a3b8e2SHong Zhang   /*-----------------------*/
38015a3b8e2SHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
381445158ffSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
38215a3b8e2SHong Zhang   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
38315a3b8e2SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
38415a3b8e2SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
38515a3b8e2SHong Zhang     if (!len_s[proc]) continue;
38615a3b8e2SHong Zhang     /* form outgoing message for i-structure:
38715a3b8e2SHong Zhang          buf_si[0]:                 nrows to be sent
38815a3b8e2SHong Zhang                [1:nrows]:           row index (global)
38915a3b8e2SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
39015a3b8e2SHong Zhang     */
39115a3b8e2SHong Zhang     /*-------------------------------------------*/
39215a3b8e2SHong Zhang     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
39315a3b8e2SHong Zhang     buf_si_i    = buf_si + nrows+1;
39415a3b8e2SHong Zhang     buf_si[0]   = nrows;
39515a3b8e2SHong Zhang     buf_si_i[0] = 0;
39615a3b8e2SHong Zhang     nrows       = 0;
39715a3b8e2SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
39815a3b8e2SHong Zhang       nzi = coi[i+1] - coi[i];
39915a3b8e2SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
40015a3b8e2SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
40115a3b8e2SHong Zhang       nrows++;
40215a3b8e2SHong Zhang     }
40315a3b8e2SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
40415a3b8e2SHong Zhang     k++;
40515a3b8e2SHong Zhang     buf_si += len_si[proc];
40615a3b8e2SHong Zhang   }
407681d504bSHong Zhang   for (i=0; i<nrecv; i++) {
408445158ffSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
40915a3b8e2SHong Zhang   }
41015a3b8e2SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
411445158ffSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
41215a3b8e2SHong Zhang 
4138cb82516SHong Zhang   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
41415a3b8e2SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
41515a3b8e2SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
41615a3b8e2SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
4178cb82516SHong Zhang #if defined(PTAP_PROFILE)
41880bb4639SHong Zhang   ierr = PetscTime(&t2);CHKERRQ(ierr);
4198cb82516SHong Zhang #endif
42067a12041SHong Zhang   /* (5) compute the local portion of Cmpi      */
42167a12041SHong Zhang   /* ------------------------------------------ */
422*9a6dcab7SHong Zhang   ierr = PetscTableDestroy(&Cta);CHKERRQ(ierr);
423*9a6dcab7SHong Zhang 
424681d504bSHong Zhang   /* set initial free space to be pN */
4258cb82516SHong Zhang   ierr          = PetscFreeSpaceGet(pN,&free_space);CHKERRQ(ierr); /* non-scalable version */
42615a3b8e2SHong Zhang   current_space = free_space;
42715a3b8e2SHong Zhang 
428445158ffSHong Zhang   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
429445158ffSHong Zhang   for (k=0; k<nrecv; k++) {
43015a3b8e2SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
43115a3b8e2SHong Zhang     nrows       = *buf_ri_k[k];
43215a3b8e2SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
43315a3b8e2SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
43415a3b8e2SHong Zhang   }
43515a3b8e2SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
436445158ffSHong Zhang 
43715a3b8e2SHong Zhang   rmax = 0;
43815a3b8e2SHong Zhang   for (i=0; i<pn; i++) {
43967a12041SHong Zhang     /* add C_loc into Cmpi */
44067a12041SHong Zhang     nzi  = c_loc->i[i+1] - c_loc->i[i];
44167a12041SHong Zhang     Jptr = c_loc->j + c_loc->i[i];
44267a12041SHong Zhang     ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
44315a3b8e2SHong Zhang 
44415a3b8e2SHong Zhang     /* add received col data into lnk */
445445158ffSHong Zhang     for (k=0; k<nrecv; k++) { /* k-th received message */
44615a3b8e2SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
44715a3b8e2SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
44815a3b8e2SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
44915a3b8e2SHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
45015a3b8e2SHong Zhang         nextrow[k]++; nextci[k]++;
45115a3b8e2SHong Zhang       }
45215a3b8e2SHong Zhang     }
45315a3b8e2SHong Zhang     nnz = lnk[0];
4548cb82516SHong Zhang 
45515a3b8e2SHong Zhang     /* copy data into free space, then initialize lnk */
45615a3b8e2SHong Zhang     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
45715a3b8e2SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
45815a3b8e2SHong Zhang     if (nnz > rmax) rmax = nnz;
45915a3b8e2SHong Zhang   }
46015a3b8e2SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
46115a3b8e2SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
462445158ffSHong Zhang   ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr);
4638cb82516SHong Zhang #if defined(PTAP_PROFILE)
46480bb4639SHong Zhang   ierr = PetscTime(&t3);CHKERRQ(ierr);
4658cb82516SHong Zhang #endif
46680bb4639SHong Zhang 
46715a3b8e2SHong Zhang   /* (6) create symbolic parallel matrix Cmpi */
46815a3b8e2SHong Zhang   /*------------------------------------------*/
46915a3b8e2SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
47015a3b8e2SHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
47115a3b8e2SHong Zhang   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
47215a3b8e2SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
47315a3b8e2SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
47415a3b8e2SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
47515a3b8e2SHong Zhang 
476445158ffSHong Zhang   /* members in merge */
477445158ffSHong Zhang   ierr = PetscFree(id_r);CHKERRQ(ierr);
478445158ffSHong Zhang   ierr = PetscFree(len_r);CHKERRQ(ierr);
479445158ffSHong Zhang   ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr);
480445158ffSHong Zhang   ierr = PetscFree(buf_ri);CHKERRQ(ierr);
481445158ffSHong Zhang   ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr);
482445158ffSHong Zhang   ierr = PetscFree(buf_rj);CHKERRQ(ierr);
483445158ffSHong Zhang   ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr);
48415a3b8e2SHong Zhang 
48515a3b8e2SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
48615a3b8e2SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
48715a3b8e2SHong Zhang   c->ptap         = ptap;
48815a3b8e2SHong Zhang   ptap->rmax      = ap_rmax;
4891a47ec54SHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
4901a47ec54SHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
49115a3b8e2SHong Zhang 
4928cb82516SHong Zhang   /* Do dense axpy in MatPtAPNumeric_MPIAIJ_MPIAIJ_new() */
4938cb82516SHong Zhang   ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
4942259aa2eSHong Zhang 
4951a47ec54SHong Zhang 
4961a47ec54SHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
4971a47ec54SHong Zhang   Cmpi->assembled        = PETSC_FALSE;
4981a47ec54SHong Zhang   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
499aa690a28SHong Zhang   Cmpi->ops->duplicate   = MatDuplicate_MPIAIJ_MatPtAP;
500aa690a28SHong Zhang   Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ;
5011a47ec54SHong Zhang   *C                     = Cmpi;
5021a47ec54SHong Zhang 
5038cb82516SHong Zhang #if defined(PTAP_PROFILE)
50480bb4639SHong Zhang   ierr = PetscTime(&t4);CHKERRQ(ierr);
505e9c1f85fSHong Zhang   if (rank == 1) {
506445158ffSHong 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);
507e9c1f85fSHong Zhang   }
5088cb82516SHong Zhang #endif
5090d3441aeSHong Zhang   PetscFunctionReturn(0);
5100d3441aeSHong Zhang }
5110d3441aeSHong Zhang 
5120d3441aeSHong Zhang #undef __FUNCT__
513aa690a28SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ"
514aa690a28SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
5150d3441aeSHong Zhang {
5160d3441aeSHong Zhang   PetscErrorCode    ierr;
5170d3441aeSHong Zhang   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
5180d3441aeSHong Zhang   Mat_SeqAIJ        *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
5198cb82516SHong Zhang   Mat_SeqAIJ        *ap,*p_loc,*p_oth,*c_seq;
5200d3441aeSHong Zhang   Mat_PtAPMPI       *ptap = c->ptap;
5219ce11a7cSHong Zhang   Mat               AP_loc,C_loc,C_oth;
5220d3441aeSHong Zhang   PetscInt          i,rstart,rend,cm,ncols,row;
5238cb82516SHong Zhang   PetscInt          *api,*apj,am = A->rmap->n,j,col,apnz;
5248cb82516SHong Zhang   PetscScalar       *apa;
5250d3441aeSHong Zhang   const PetscInt    *cols;
5260d3441aeSHong Zhang   const PetscScalar *vals;
5278cb82516SHong Zhang #if defined(PTAP_PROFILE)
5288cb82516SHong Zhang   PetscMPIInt       rank;
5298cb82516SHong Zhang   MPI_Comm          comm;
5300d3441aeSHong Zhang   PetscLogDouble    t0,t1,t2,t3,t4,eR,eAP,eCseq,eCmpi;
5318cb82516SHong Zhang #endif
5320d3441aeSHong Zhang 
5330d3441aeSHong Zhang   PetscFunctionBegin;
5340d3441aeSHong Zhang   ierr = MatZeroEntries(C);CHKERRQ(ierr);
5350d3441aeSHong Zhang 
536e497d3c8SHong Zhang   /* 1) get R = Pd^T,Ro = Po^T */
5378cb82516SHong Zhang #if defined(PTAP_PROFILE)
5380d3441aeSHong Zhang   ierr = PetscTime(&t0);CHKERRQ(ierr);
5398cb82516SHong Zhang #endif
540748c7196SHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
5410d3441aeSHong Zhang     ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
5420d3441aeSHong Zhang     ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
543748c7196SHong Zhang   }
5448cb82516SHong Zhang #if defined(PTAP_PROFILE)
5450d3441aeSHong Zhang   ierr = PetscTime(&t1);CHKERRQ(ierr);
5460d3441aeSHong Zhang   eR = t1 - t0;
5478cb82516SHong Zhang #endif
5480d3441aeSHong Zhang 
549e497d3c8SHong Zhang   /* 2) get AP_loc */
5500d3441aeSHong Zhang   AP_loc = ptap->AP_loc;
5518cb82516SHong Zhang   ap = (Mat_SeqAIJ*)AP_loc->data;
5520d3441aeSHong Zhang 
553e497d3c8SHong Zhang   /* 2-1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
5540d3441aeSHong Zhang   /*-----------------------------------------------------*/
555748c7196SHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
556748c7196SHong Zhang     /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */
5570d3441aeSHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
5580d3441aeSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
559e497d3c8SHong Zhang   }
5600d3441aeSHong Zhang 
5618cb82516SHong Zhang   /* 2-2) compute numeric A_loc*P - dominating part */
5628cb82516SHong Zhang   /* ---------------------------------------------- */
5630d3441aeSHong Zhang   /* get data from symbolic products */
5648cb82516SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
5658cb82516SHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
5668cb82516SHong Zhang   apa   = ptap->apa;
567681d504bSHong Zhang   api   = ap->i;
568681d504bSHong Zhang   apj   = ap->j;
569e497d3c8SHong Zhang   for (i=0; i<am; i++) {
570445158ffSHong Zhang     /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */
571e497d3c8SHong Zhang     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
572e497d3c8SHong Zhang     apnz = api[i+1] - api[i];
573e497d3c8SHong Zhang     for (j=0; j<apnz; j++) {
574e497d3c8SHong Zhang       col = apj[j+api[i]];
575e497d3c8SHong Zhang       ap->a[j+ap->i[i]] = apa[col];
576e497d3c8SHong Zhang       apa[col] = 0.0;
5770d3441aeSHong Zhang     }
578aa690a28SHong Zhang     ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
5790d3441aeSHong Zhang   }
5808cb82516SHong Zhang #if defined(PTAP_PROFILE)
5810d3441aeSHong Zhang   ierr = PetscTime(&t2);CHKERRQ(ierr);
5820d3441aeSHong Zhang   eAP = t2 - t1;
5838cb82516SHong Zhang #endif
5840d3441aeSHong Zhang 
5858cb82516SHong Zhang   /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */
586445158ffSHong Zhang   ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr);
587445158ffSHong Zhang   ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr);
5889ce11a7cSHong Zhang   C_loc = ptap->C_loc;
5899ce11a7cSHong Zhang   C_oth = ptap->C_oth;
5908cb82516SHong Zhang #if defined(PTAP_PROFILE)
5910d3441aeSHong Zhang   ierr = PetscTime(&t3);CHKERRQ(ierr);
5920d3441aeSHong Zhang   eCseq = t3 - t2;
5938cb82516SHong Zhang #endif
5940d3441aeSHong Zhang 
5950d3441aeSHong Zhang   /* add C_loc and Co to to C */
5960d3441aeSHong Zhang   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
5970d3441aeSHong Zhang 
5980d3441aeSHong Zhang   /* C_loc -> C */
5990d3441aeSHong Zhang   cm    = C_loc->rmap->N;
6009ce11a7cSHong Zhang   c_seq = (Mat_SeqAIJ*)C_loc->data;
6018cb82516SHong Zhang   cols = c_seq->j;
6028cb82516SHong Zhang   vals = c_seq->a;
6030d3441aeSHong Zhang   for (i=0; i<cm; i++) {
6049ce11a7cSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
6050d3441aeSHong Zhang     row = rstart + i;
6060d3441aeSHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
6078cb82516SHong Zhang     cols += ncols; vals += ncols;
6080d3441aeSHong Zhang   }
6090d3441aeSHong Zhang 
6100d3441aeSHong Zhang   /* Co -> C, off-processor part */
6119ce11a7cSHong Zhang   cm = C_oth->rmap->N;
6129ce11a7cSHong Zhang   c_seq = (Mat_SeqAIJ*)C_oth->data;
6138cb82516SHong Zhang   cols = c_seq->j;
6148cb82516SHong Zhang   vals = c_seq->a;
6159ce11a7cSHong Zhang   for (i=0; i<cm; i++) {
6169ce11a7cSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
6170d3441aeSHong Zhang     row = p->garray[i];
6180d3441aeSHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
6198cb82516SHong Zhang     cols += ncols; vals += ncols;
6200d3441aeSHong Zhang   }
6210d3441aeSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6220d3441aeSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6238cb82516SHong Zhang #if defined(PTAP_PROFILE)
6240d3441aeSHong Zhang   ierr = PetscTime(&t4);CHKERRQ(ierr);
6250d3441aeSHong Zhang   eCmpi = t4 - t3;
6260d3441aeSHong Zhang 
6278cb82516SHong Zhang   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
6288cb82516SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
6290d3441aeSHong Zhang   if (rank==1) {
6300d3441aeSHong 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);
6310d3441aeSHong Zhang   }
6328cb82516SHong Zhang #endif
633748c7196SHong Zhang   ptap->reuse = MAT_REUSE_MATRIX;
6340d3441aeSHong Zhang   PetscFunctionReturn(0);
6350d3441aeSHong Zhang }
6360d3441aeSHong Zhang 
6370d3441aeSHong Zhang #undef __FUNCT__
6388403a639SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap"
6398403a639SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(Mat A,Mat P,PetscReal fill,Mat *C)
64042c57489SHong Zhang {
64142c57489SHong Zhang   PetscErrorCode      ierr;
64277584fe6SHong Zhang   Mat                 Cmpi;
643f8487c73SHong Zhang   Mat_PtAPMPI         *ptap;
6440298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
645f8487c73SHong Zhang   Mat_MPIAIJ          *a        =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
64642c57489SHong Zhang   Mat_SeqAIJ          *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
64742c57489SHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth;
648ded4f561SHong Zhang   PetscInt            *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ;
64977584fe6SHong Zhang   PetscInt            *adi=ad->i,*aj,*aoi=ao->i,nnz;
650a914927fSHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
651d16ca5f0SHong Zhang   PetscInt            am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n;
65242c57489SHong Zhang   PetscBT             lnkbt;
653ce94432eSBarry Smith   MPI_Comm            comm;
654a914927fSHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0;
65542c57489SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
65642c57489SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
65724ecddacSHong Zhang   PetscInt            nzi,*pti,*ptj;
65842c57489SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
659ded4f561SHong Zhang   MPI_Request         *swaits,*rwaits;
660ded4f561SHong Zhang   MPI_Status          *sstatus,rstatus;
66142c57489SHong Zhang   Mat_Merge_SeqsToMPI *merge;
66277584fe6SHong Zhang   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0;
663d16ca5f0SHong Zhang   PetscReal           afill=1.0,afill_tmp;
664a914927fSHong Zhang   PetscInt            rmax;
66542c57489SHong Zhang 
66642c57489SHong Zhang   PetscFunctionBegin;
667ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
66883445d95SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
66983445d95SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
67083445d95SHong Zhang 
671f8487c73SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
672b00a9115SJed Brown   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
673b00a9115SJed Brown   ierr        = PetscNew(&merge);CHKERRQ(ierr);
6749f4155fbSHong Zhang   ptap->merge = merge;
675f8487c73SHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
6766c6e00e0SHong Zhang 
6776c6e00e0SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
678b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
679fe615159SHong Zhang 
6806c6e00e0SHong Zhang   /* get P_loc by taking all local rows of P */
681a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
6826c6e00e0SHong Zhang 
683a1a4e84aSHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
684a1a4e84aSHong Zhang   p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
6856c6e00e0SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
6866c6e00e0SHong Zhang   pi_oth = p_oth->i; pj_oth = p_oth->j;
68742c57489SHong Zhang 
6882addb229SHong Zhang   /* (1) compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth (api,apj) */
6892addb229SHong Zhang   /*--------------------------------------------------------------------------*/
690854ce69bSBarry Smith   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
69182412ba7SHong Zhang   api[0] = 0;
69283445d95SHong Zhang 
69383445d95SHong Zhang   /* create and initialize a linked list */
694a914927fSHong Zhang   ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr);
69583445d95SHong Zhang 
696a914927fSHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) -OOM for ex56, np=8k on Intrepid! */
697d16ca5f0SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
698f4a743e1SHong Zhang   current_space = free_space;
699f4a743e1SHong Zhang 
700f4a743e1SHong Zhang   for (i=0; i<am; i++) {
701f4a743e1SHong Zhang     /* diagonal portion of A */
702ded4f561SHong Zhang     nzi = adi[i+1] - adi[i];
70377584fe6SHong Zhang     aj  = ad->j + adi[i];
704ded4f561SHong Zhang     for (j=0; j<nzi; j++) {
70577584fe6SHong Zhang       row  = aj[j];
70682412ba7SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
707ded4f561SHong Zhang       Jptr = pj_loc + pi_loc[row];
70883445d95SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
709a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
710f4a743e1SHong Zhang     }
711f4a743e1SHong Zhang     /* off-diagonal portion of A */
712ded4f561SHong Zhang     nzi = aoi[i+1] - aoi[i];
71377584fe6SHong Zhang     aj  = ao->j + aoi[i];
714ded4f561SHong Zhang     for (j=0; j<nzi; j++) {
71577584fe6SHong Zhang       row  = aj[j];
71682412ba7SHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
717ded4f561SHong Zhang       Jptr = pj_oth + pi_oth[row];
718a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
719f4a743e1SHong Zhang     }
720a914927fSHong Zhang     apnz     = lnk[0];
72182412ba7SHong Zhang     api[i+1] = api[i] + apnz;
72277584fe6SHong Zhang     if (ap_rmax < apnz) ap_rmax = apnz;
723f4a743e1SHong Zhang 
72483445d95SHong Zhang     /* if free space is not available, double the total space in the list */
72582412ba7SHong Zhang     if (current_space->local_remaining<apnz) {
7262ba1acd5SHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
727f2b054eeSHong Zhang       nspacedouble++;
728f4a743e1SHong Zhang     }
729f4a743e1SHong Zhang 
730f4a743e1SHong Zhang     /* Copy data into free space, then initialize lnk */
731a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
7322205254eSKarl Rupp 
73382412ba7SHong Zhang     current_space->array           += apnz;
73482412ba7SHong Zhang     current_space->local_used      += apnz;
73582412ba7SHong Zhang     current_space->local_remaining -= apnz;
736f4a743e1SHong Zhang   }
737a914927fSHong Zhang 
73882412ba7SHong Zhang   /* Allocate space for apj, initialize apj, and */
739f4a743e1SHong Zhang   /* destroy list of free space and other temporary array(s) */
740854ce69bSBarry Smith   ierr      = PetscMalloc1(api[am]+1,&apj);CHKERRQ(ierr);
74177584fe6SHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
742118727c9SMark F. Adams   afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1);
743d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
744f4a743e1SHong Zhang 
7452addb229SHong Zhang   /* (2) determine symbolic Co=(p->B)^T*AP - send to others (coi,coj)*/
7462addb229SHong Zhang   /*-----------------------------------------------------------------*/
747de0260b3SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
74842c57489SHong Zhang 
749ded4f561SHong Zhang   /* then, compute symbolic Co = (p->B)^T*AP */
750d0f46423SBarry Smith   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
75183445d95SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
752854ce69bSBarry Smith   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
753de0260b3SHong Zhang   coi[0] = 0;
75442c57489SHong Zhang 
755d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */
756d16ca5f0SHong Zhang   nnz           = fill*(poti[pon] + api[am]);
75722d28d08SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
75842c57489SHong Zhang   current_space = free_space;
75942c57489SHong Zhang 
760de0260b3SHong Zhang   for (i=0; i<pon; i++) {
76182412ba7SHong Zhang     pnz = poti[i+1] - poti[i];
76277584fe6SHong Zhang     ptJ = potj + poti[i];
76377584fe6SHong Zhang     for (j=0; j<pnz; j++) {
76477584fe6SHong Zhang       row  = ptJ[j]; /* row of AP == col of Pot */
76582412ba7SHong Zhang       apnz = api[row+1] - api[row];
766ded4f561SHong Zhang       Jptr = apj + api[row];
76783445d95SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
768a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
76942c57489SHong Zhang     }
770a914927fSHong Zhang     nnz = lnk[0];
77142c57489SHong Zhang 
77283445d95SHong Zhang     /* If free space is not available, double the total space in the list */
773ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
7744238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
775d16ca5f0SHong Zhang       nspacedouble++;
77642c57489SHong Zhang     }
77742c57489SHong Zhang 
77842c57489SHong Zhang     /* Copy data into free space, and zero out denserows */
779a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
7802205254eSKarl Rupp 
781ded4f561SHong Zhang     current_space->array           += nnz;
782ded4f561SHong Zhang     current_space->local_used      += nnz;
783ded4f561SHong Zhang     current_space->local_remaining -= nnz;
7842205254eSKarl Rupp 
785ded4f561SHong Zhang     coi[i+1] = coi[i] + nnz;
78642c57489SHong Zhang   }
7876b911d16SHong Zhang 
7886b911d16SHong Zhang   ierr      = PetscMalloc1(coi[pon],&coj);CHKERRQ(ierr);
789a1a86e44SBarry Smith   ierr      = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
790118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1);
791d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
792de0260b3SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
79342c57489SHong Zhang 
7946b911d16SHong Zhang   /* (3) send j-array (coj) of Co to other processors */
7956b911d16SHong Zhang   /*--------------------------------------------------*/
7966b911d16SHong Zhang   ierr = PetscCalloc1(size,&merge->len_s);CHKERRQ(ierr);
7976b911d16SHong Zhang   len_s        = merge->len_s;
7986b911d16SHong Zhang   merge->nsend = 0;
7996b911d16SHong Zhang 
8006b911d16SHong Zhang 
80142c57489SHong Zhang   /* determine row ownership */
80226283091SBarry Smith   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
8037a2fc3feSBarry Smith   merge->rowmap->n  = pn;
8047a2fc3feSBarry Smith   merge->rowmap->bs = 1;
8052205254eSKarl Rupp 
80626283091SBarry Smith   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
8077a2fc3feSBarry Smith   owners = merge->rowmap->range;
80842c57489SHong Zhang 
80942c57489SHong Zhang   /* determine the number of messages to send, their lengths */
810dcca6d9dSJed Brown   ierr = PetscMalloc2(size,&len_si,size,&sstatus);CHKERRQ(ierr);
81183445d95SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
812854ce69bSBarry Smith   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
813de0260b3SHong Zhang 
81483445d95SHong Zhang   proc = 0;
815de0260b3SHong Zhang   for (i=0; i<pon; i++) {
816de0260b3SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
8172addb229SHong Zhang     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
818ee6e4b5bSHong Zhang     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
81983445d95SHong Zhang   }
820de0260b3SHong Zhang 
8212addb229SHong Zhang   len          = 0; /* max length of buf_si[], see (4) */
822de0260b3SHong Zhang   owners_co[0] = 0;
82342c57489SHong Zhang   for (proc=0; proc<size; proc++) {
824de0260b3SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
825ee6e4b5bSHong Zhang     if (len_s[proc]) {
82642c57489SHong Zhang       merge->nsend++;
8272addb229SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
82842c57489SHong Zhang       len         += len_si[proc];
82942c57489SHong Zhang     }
83042c57489SHong Zhang   }
83142c57489SHong Zhang 
832ded4f561SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
8330298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
83442c57489SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
83542c57489SHong Zhang 
836ded4f561SHong Zhang   /* post the Irecv and Isend of coj */
837529bc5dcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
838529bc5dcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
839854ce69bSBarry Smith   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
84042c57489SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
84142c57489SHong Zhang     if (!len_s[proc]) continue;
842de0260b3SHong Zhang     i    = owners_co[proc];
843529bc5dcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
84442c57489SHong Zhang     k++;
84542c57489SHong Zhang   }
84642c57489SHong Zhang 
847ded4f561SHong Zhang   /* receives and sends of coj are complete */
84877584fe6SHong Zhang   for (i=0; i<merge->nrecv; i++) {
849c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
850ded4f561SHong Zhang   }
851ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
8520c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
85382412ba7SHong Zhang 
8546b911d16SHong Zhang   /* (4) send and recv coi */
8556b911d16SHong Zhang   /*-----------------------*/
856529bc5dcSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
857529bc5dcSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
858854ce69bSBarry Smith   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
85942c57489SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
86042c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
86142c57489SHong Zhang     if (!len_s[proc]) continue;
86242c57489SHong Zhang     /* form outgoing message for i-structure:
86342c57489SHong Zhang          buf_si[0]:                 nrows to be sent
86442c57489SHong Zhang                [1:nrows]:           row index (global)
86542c57489SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
86642c57489SHong Zhang     */
86742c57489SHong Zhang     /*-------------------------------------------*/
8682addb229SHong Zhang     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
86942c57489SHong Zhang     buf_si_i    = buf_si + nrows+1;
87042c57489SHong Zhang     buf_si[0]   = nrows;
87142c57489SHong Zhang     buf_si_i[0] = 0;
87242c57489SHong Zhang     nrows       = 0;
873de0260b3SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
874ded4f561SHong Zhang       nzi = coi[i+1] - coi[i];
875ded4f561SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
876de0260b3SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
87742c57489SHong Zhang       nrows++;
87842c57489SHong Zhang     }
879529bc5dcSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
88042c57489SHong Zhang     k++;
88142c57489SHong Zhang     buf_si += len_si[proc];
88242c57489SHong Zhang   }
883ded4f561SHong Zhang   i = merge->nrecv;
884ded4f561SHong Zhang   while (i--) {
885c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
886ded4f561SHong Zhang   }
887ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
8880c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
889a914927fSHong Zhang 
89024ecddacSHong Zhang   ierr = PetscFree2(len_si,sstatus);CHKERRQ(ierr);
89142c57489SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
892ded4f561SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
89342c57489SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
89442c57489SHong Zhang 
8956b911d16SHong Zhang   /* (5) compute the local portion of C (mpi mat) */
8966b911d16SHong Zhang   /*----------------------------------------------*/
897ded4f561SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
898ded4f561SHong Zhang 
89924ecddacSHong Zhang   /* allocate pti array and free space for accumulating nonzero column info */
900854ce69bSBarry Smith   ierr   = PetscMalloc1(pn+1,&pti);CHKERRQ(ierr);
90124ecddacSHong Zhang   pti[0] = 0;
90242c57489SHong Zhang 
903d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
904d16ca5f0SHong Zhang   nnz           = fill*(pi_loc[pm] + api[am]);
90522d28d08SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
90642c57489SHong Zhang   current_space = free_space;
90742c57489SHong Zhang 
908dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
90942c57489SHong Zhang   for (k=0; k<merge->nrecv; k++) {
91042c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
91142c57489SHong Zhang     nrows       = *buf_ri_k[k];
91242c57489SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
91342c57489SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
91442c57489SHong Zhang   }
91542c57489SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
91608cb4532SHong Zhang   rmax = 0;
91742c57489SHong Zhang   for (i=0; i<pn; i++) {
918ded4f561SHong Zhang     /* add pdt[i,:]*AP into lnk */
919ded4f561SHong Zhang     pnz = pdti[i+1] - pdti[i];
92077584fe6SHong Zhang     ptJ = pdtj + pdti[i];
92177584fe6SHong Zhang     for (j=0; j<pnz; j++) {
92277584fe6SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
923ded4f561SHong Zhang       apnz = api[row+1] - api[row];
924ded4f561SHong Zhang       Jptr = apj + api[row];
925ded4f561SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
926a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
927ded4f561SHong Zhang     }
928a914927fSHong Zhang 
92942c57489SHong Zhang     /* add received col data into lnk */
93042c57489SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
93142c57489SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
932ded4f561SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
933ded4f561SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
934a914927fSHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
93542c57489SHong Zhang         nextrow[k]++; nextci[k]++;
93642c57489SHong Zhang       }
93742c57489SHong Zhang     }
938a914927fSHong Zhang     nnz = lnk[0];
93942c57489SHong Zhang 
94042c57489SHong Zhang     /* if free space is not available, make more free space */
941ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
9424238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
943d16ca5f0SHong Zhang       nspacedouble++;
94442c57489SHong Zhang     }
94542c57489SHong Zhang     /* copy data into free space, then initialize lnk */
946a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
947ded4f561SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
9482205254eSKarl Rupp 
949ded4f561SHong Zhang     current_space->array           += nnz;
950ded4f561SHong Zhang     current_space->local_used      += nnz;
951ded4f561SHong Zhang     current_space->local_remaining -= nnz;
9522205254eSKarl Rupp 
95324ecddacSHong Zhang     pti[i+1] = pti[i] + nnz;
95408cb4532SHong Zhang     if (nnz > rmax) rmax = nnz;
95542c57489SHong Zhang   }
956ded4f561SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
9570572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
95842c57489SHong Zhang 
959854ce69bSBarry Smith   ierr      = PetscMalloc1(pti[pn]+1,&ptj);CHKERRQ(ierr);
96024ecddacSHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,ptj);CHKERRQ(ierr);
96124ecddacSHong Zhang   afill_tmp = (PetscReal)pti[pn]/(pi_loc[pm] + api[am]+1);
962d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
96342c57489SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
96442c57489SHong Zhang 
9656b911d16SHong Zhang   /* (6) create symbolic parallel matrix Cmpi */
9666b911d16SHong Zhang   /*------------------------------------------*/
96777584fe6SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
96877584fe6SHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
96933d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
97077584fe6SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
97177584fe6SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
97242c57489SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
973a2f3521dSMark F. Adams 
974b091e509SHong Zhang   merge->bi        = pti;      /* Cseq->i */
975b091e509SHong Zhang   merge->bj        = ptj;      /* Cseq->j */
976b091e509SHong Zhang   merge->coi       = coi;      /* Co->i   */
977b091e509SHong Zhang   merge->coj       = coj;      /* Co->j   */
97842c57489SHong Zhang   merge->buf_ri    = buf_ri;
97942c57489SHong Zhang   merge->buf_rj    = buf_rj;
980de0260b3SHong Zhang   merge->owners_co = owners_co;
98177584fe6SHong Zhang   merge->destroy   = Cmpi->ops->destroy;
98242c57489SHong Zhang 
98377584fe6SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
98477584fe6SHong Zhang   c           = (Mat_MPIAIJ*)Cmpi->data;
985f8487c73SHong Zhang   c->ptap     = ptap;
98677584fe6SHong Zhang   ptap->api   = api;
98777584fe6SHong Zhang   ptap->apj   = apj;
988d6ab1dc0SHong Zhang   ptap->rmax  = ap_rmax;
9898403a639SHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
9908403a639SHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
9918403a639SHong Zhang 
9928403a639SHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
9938403a639SHong Zhang   Cmpi->assembled        = PETSC_FALSE;
9948403a639SHong Zhang   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
9958403a639SHong Zhang   Cmpi->ops->duplicate   = MatDuplicate_MPIAIJ_MatPtAP;
9968403a639SHong Zhang   Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap;
99777584fe6SHong Zhang   *C                     = Cmpi;
998414894bdSHong Zhang 
999414894bdSHong Zhang   /* flag 'scalable' determines which implementations to be used:
1000414894bdSHong Zhang        0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa;
1001414894bdSHong Zhang        1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */
1002414894bdSHong Zhang   /* set default scalable */
1003da0a95b2SSatish Balay   ptap->scalable = PETSC_FALSE; /* PETSC_TRUE; */
10042205254eSKarl Rupp 
10050298fd71SBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)Cmpi)->prefix,"-matptap_scalable",&ptap->scalable,NULL);CHKERRQ(ierr);
1006414894bdSHong Zhang   if (!ptap->scalable) {  /* Do dense axpy */
10071795a4d1SJed Brown     ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
1008414894bdSHong Zhang   } else {
10091795a4d1SJed Brown     ierr = PetscCalloc1(ap_rmax+1,&ptap->apa);CHKERRQ(ierr);
1010414894bdSHong Zhang   }
1011414894bdSHong Zhang 
1012f2b054eeSHong Zhang #if defined(PETSC_USE_INFO)
101324ecddacSHong Zhang   if (pti[pn] != 0) {
101457622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
101557622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1016f2b054eeSHong Zhang   } else {
101777584fe6SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1018f2b054eeSHong Zhang   }
1019f2b054eeSHong Zhang #endif
102042c57489SHong Zhang   PetscFunctionReturn(0);
102142c57489SHong Zhang }
102242c57489SHong Zhang 
102342c57489SHong Zhang #undef __FUNCT__
10248403a639SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap"
10258403a639SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap(Mat A,Mat P,Mat C)
102642c57489SHong Zhang {
102742c57489SHong Zhang   PetscErrorCode      ierr;
1028f8487c73SHong Zhang   Mat_MPIAIJ          *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
102942c57489SHong Zhang   Mat_SeqAIJ          *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
1030de0260b3SHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
103142c57489SHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth;
1032f8487c73SHong Zhang   Mat_PtAPMPI         *ptap;
10339f4155fbSHong Zhang   Mat_Merge_SeqsToMPI *merge;
10341d633602SHong Zhang   PetscInt            *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp;
103582412ba7SHong Zhang   PetscInt            *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj;
103682412ba7SHong Zhang   PetscInt            i,j,k,anz,pnz,apnz,nextap,row,*cj;
1037e900a757SHong Zhang   MatScalar           *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp;
1038d0f46423SBarry Smith   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1039ce94432eSBarry Smith   MPI_Comm            comm;
104042c57489SHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
104182412ba7SHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
104242c57489SHong Zhang   PetscInt            **buf_ri,**buf_rj;
104350f5bf86SHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
104442c57489SHong Zhang   MPI_Request         *s_waits,*r_waits;
104542c57489SHong Zhang   MPI_Status          *status;
104682412ba7SHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
104782412ba7SHong Zhang   PetscInt            *api,*apj,*coi,*coj;
1048d0f46423SBarry Smith   PetscInt            *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend;
10496ce94e8fSHong Zhang   PetscBool           scalable;
105038571addSHong Zhang #if defined(PTAP_PROFILE)
1051e497d3c8SHong Zhang   PetscLogDouble t0,t1,t2,eP,t3,t4,et2_AP=0.0,ePtAP=0.0,t2_0,t2_1,t2_2;
105238571addSHong Zhang #endif
105342c57489SHong Zhang 
105442c57489SHong Zhang   PetscFunctionBegin;
1055ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
105638571addSHong Zhang #if defined(PTAP_PROFILE)
10578563dfccSBarry Smith   ierr = PetscTime(&t0);CHKERRQ(ierr);
105838571addSHong Zhang #endif
105942c57489SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
106042c57489SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
106142c57489SHong Zhang 
1062f8487c73SHong Zhang   ptap = c->ptap;
1063ce94432eSBarry 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");
1064f8487c73SHong Zhang   merge    = ptap->merge;
1065414894bdSHong Zhang   apa      = ptap->apa;
10666ce94e8fSHong Zhang   scalable = ptap->scalable;
10676c6e00e0SHong Zhang 
1068a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
10696b911d16SHong Zhang   /*-----------------------------------------------------*/
1070f8487c73SHong Zhang   if (ptap->reuse == MAT_INITIAL_MATRIX) {
10719f4155fbSHong Zhang     /* P_oth and P_loc are obtained in MatPtASymbolic(), skip calling MatGetBrowsOfAoCols() and MatMPIAIJGetLocalMat() */
1072f8487c73SHong Zhang     ptap->reuse = MAT_REUSE_MATRIX;
10736c6e00e0SHong Zhang   } else { /* update numerical values of P_oth and P_loc */
1074b7f45c76SHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
1075a1a4e84aSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
10766c6e00e0SHong Zhang   }
107738571addSHong Zhang #if defined(PTAP_PROFILE)
10788563dfccSBarry Smith   ierr = PetscTime(&t1);CHKERRQ(ierr);
107905d62848SHong Zhang   eP = t1-t0;
108038571addSHong Zhang #endif
108105d62848SHong Zhang   /*
108205d62848SHong Zhang   printf("[%d] Ad: %d, %d; Ao: %d, %d; P_loc: %d, %d; P_oth %d, %d;\n",rank,
108305d62848SHong Zhang          a->A->rmap->N,a->A->cmap->N,a->B->rmap->N,a->B->cmap->N,
108405d62848SHong Zhang          ptap->P_loc->rmap->N,ptap->P_loc->cmap->N,
108505d62848SHong Zhang          ptap->P_oth->rmap->N,ptap->P_oth->cmap->N);
108605d62848SHong Zhang    */
1087f8487c73SHong Zhang 
1088f9473cf7SHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
1089f9473cf7SHong Zhang   /*--------------------------------------------------------------*/
109042c57489SHong Zhang   /* get data from symbolic products */
1091a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
1092a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
1093a9d06231SHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a;
109442c57489SHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
109542c57489SHong Zhang 
1096de0260b3SHong Zhang   coi  = merge->coi; coj = merge->coj;
10971795a4d1SJed Brown   ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
1098de0260b3SHong Zhang 
1099de0260b3SHong Zhang   bi     = merge->bi; bj = merge->bj;
11007a2fc3feSBarry Smith   owners = merge->rowmap->range;
11011795a4d1SJed Brown   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);  /* ba: Cseq->a */
110242c57489SHong Zhang 
1103a1a4e84aSHong Zhang   api = ptap->api; apj = ptap->apj;
1104d9824c15SHong Zhang 
11059516a85cSHong Zhang   if (!scalable) { /* Do dense axpy on apa (length of pN, stores A[i,:]*P) - nonscalable, but faster (could take 1/3 scalable time) */
1106b38d1a2dSJed Brown     ierr = PetscInfo(C,"Using non-scalable dense axpy\n");CHKERRQ(ierr);
11078cb82516SHong Zhang #if defined(PTAP_PROFILE)
110805d62848SHong Zhang     ierr = PetscTime(&t1);CHKERRQ(ierr);
11098cb82516SHong Zhang #endif
1110a9d06231SHong Zhang     for (i=0; i<am; i++) {
1111b091e509SHong Zhang #if defined(PTAP_PROFILE)
11128563dfccSBarry Smith       ierr = PetscTime(&t2_0);CHKERRQ(ierr);
1113b091e509SHong Zhang #endif
1114a9d06231SHong Zhang       /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */
1115a9d06231SHong Zhang       /*------------------------------------------------------------*/
1116a9d06231SHong Zhang       apJ = apj + api[i];
1117a9d06231SHong Zhang 
1118a9d06231SHong Zhang       /* diagonal portion of A */
1119a9d06231SHong Zhang       anz = adi[i+1] - adi[i];
1120a9d06231SHong Zhang       adj = ad->j + adi[i];
1121a9d06231SHong Zhang       ada = ad->a + adi[i];
1122a9d06231SHong Zhang       for (j=0; j<anz; j++) {
1123a9d06231SHong Zhang         row = adj[j];
1124a9d06231SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
1125a9d06231SHong Zhang         pj  = pj_loc + pi_loc[row];
1126a9d06231SHong Zhang         pa  = pa_loc + pi_loc[row];
1127a9d06231SHong Zhang 
1128a9d06231SHong Zhang         /* perform dense axpy */
1129e900a757SHong Zhang         valtmp = ada[j];
1130a9d06231SHong Zhang         for (k=0; k<pnz; k++) {
1131e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
1132a9d06231SHong Zhang         }
1133a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
1134a9d06231SHong Zhang       }
1135a9d06231SHong Zhang 
1136a9d06231SHong Zhang       /* off-diagonal portion of A */
1137a9d06231SHong Zhang       anz = aoi[i+1] - aoi[i];
1138a9d06231SHong Zhang       aoj = ao->j + aoi[i];
1139a9d06231SHong Zhang       aoa = ao->a + aoi[i];
1140a9d06231SHong Zhang       for (j=0; j<anz; j++) {
1141a9d06231SHong Zhang         row = aoj[j];
1142a9d06231SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
1143a9d06231SHong Zhang         pj  = pj_oth + pi_oth[row];
1144a9d06231SHong Zhang         pa  = pa_oth + pi_oth[row];
1145a9d06231SHong Zhang 
1146a9d06231SHong Zhang         /* perform dense axpy */
1147e900a757SHong Zhang         valtmp = aoa[j];
1148a9d06231SHong Zhang         for (k=0; k<pnz; k++) {
1149e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
1150a9d06231SHong Zhang         }
1151a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
1152a9d06231SHong Zhang       }
1153b091e509SHong Zhang #if defined(PTAP_PROFILE)
11548563dfccSBarry Smith       ierr    = PetscTime(&t2_1);CHKERRQ(ierr);
1155b091e509SHong Zhang       et2_AP += t2_1 - t2_0;
1156b091e509SHong Zhang #endif
1157a9d06231SHong Zhang 
1158a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
1159a9d06231SHong Zhang       /*--------------------------------------------------------------*/
1160a9d06231SHong Zhang       apnz = api[i+1] - api[i];
1161a9d06231SHong Zhang       /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */
1162a9d06231SHong Zhang       pnz = po->i[i+1] - po->i[i];
1163e900a757SHong Zhang       poJ = po->j + po->i[i];
1164a9d06231SHong Zhang       pA  = po->a + po->i[i];
1165a9d06231SHong Zhang       for (j=0; j<pnz; j++) {
1166e900a757SHong Zhang         row = poJ[j];
1167e900a757SHong Zhang         cnz = coi[row+1] - coi[row];
1168e900a757SHong Zhang         cj  = coj + coi[row];
1169e900a757SHong Zhang         ca  = coa + coi[row];
1170a9d06231SHong Zhang         /* perform dense axpy */
1171e900a757SHong Zhang         valtmp = pA[j];
1172a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
1173e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
1174a9d06231SHong Zhang         }
1175a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1176a9d06231SHong Zhang       }
1177a9d06231SHong Zhang       /* put the value into Cd (diagonal part) */
1178a9d06231SHong Zhang       pnz = pd->i[i+1] - pd->i[i];
1179e900a757SHong Zhang       pdJ = pd->j + pd->i[i];
1180a9d06231SHong Zhang       pA  = pd->a + pd->i[i];
1181a9d06231SHong Zhang       for (j=0; j<pnz; j++) {
1182e900a757SHong Zhang         row = pdJ[j];
1183e900a757SHong Zhang         cnz = bi[row+1] - bi[row];
1184e900a757SHong Zhang         cj  = bj + bi[row];
1185e900a757SHong Zhang         ca  = ba + bi[row];
1186a9d06231SHong Zhang         /* perform dense axpy */
1187e900a757SHong Zhang         valtmp = pA[j];
1188a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
1189e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
1190a9d06231SHong Zhang         }
1191a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1192a9d06231SHong Zhang       }
11938403a639SHong Zhang 
1194a9d06231SHong Zhang       /* zero the current row of A*P */
1195a9d06231SHong Zhang       for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0;
1196b091e509SHong Zhang #if defined(PTAP_PROFILE)
11978563dfccSBarry Smith       ierr      = PetscTime(&t2_2);CHKERRQ(ierr);
119805d62848SHong Zhang       ePtAP += t2_2 - t2_1;
1199b091e509SHong Zhang #endif
1200a9d06231SHong Zhang     }
120138571addSHong Zhang   } else { /* Do sparse axpy on apa (length of ap_rmax, stores A[i,:]*P) - scalable, but slower */
1202b38d1a2dSJed Brown     ierr = PetscInfo(C,"Using scalable sparse axpy\n");CHKERRQ(ierr);
120338571addSHong Zhang     /*-----------------------------------------------------------------------------------------*/
1204a9d06231SHong Zhang     pA=pa_loc;
120542c57489SHong Zhang     for (i=0; i<am; i++) {
1206b091e509SHong Zhang #if defined(PTAP_PROFILE)
12078563dfccSBarry Smith       ierr = PetscTime(&t2_0);CHKERRQ(ierr);
1208b091e509SHong Zhang #endif
1209f4a743e1SHong Zhang       /* form i-th sparse row of A*P */
121082412ba7SHong Zhang       apnz = api[i+1] - api[i];
121182412ba7SHong Zhang       apJ  = apj + api[i];
121242c57489SHong Zhang       /* diagonal portion of A */
121382412ba7SHong Zhang       anz = adi[i+1] - adi[i];
12141d633602SHong Zhang       adj = ad->j + adi[i];
12151d633602SHong Zhang       ada = ad->a + adi[i];
121682412ba7SHong Zhang       for (j=0; j<anz; j++) {
12171d633602SHong Zhang         row    = adj[j];
121882412ba7SHong Zhang         pnz    = pi_loc[row+1] - pi_loc[row];
121982412ba7SHong Zhang         pj     = pj_loc + pi_loc[row];
122082412ba7SHong Zhang         pa     = pa_loc + pi_loc[row];
1221e900a757SHong Zhang         valtmp = ada[j];
1222f4a743e1SHong Zhang         nextp  = 0;
122382412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
122482412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
1225e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
122642c57489SHong Zhang           }
122742c57489SHong Zhang         }
1228dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
122942c57489SHong Zhang       }
123042c57489SHong Zhang       /* off-diagonal portion of A */
123182412ba7SHong Zhang       anz = aoi[i+1] - aoi[i];
12321d633602SHong Zhang       aoj = ao->j + aoi[i];
12331d633602SHong Zhang       aoa = ao->a + aoi[i];
123482412ba7SHong Zhang       for (j=0; j<anz; j++) {
12351d633602SHong Zhang         row    = aoj[j];
123682412ba7SHong Zhang         pnz    = pi_oth[row+1] - pi_oth[row];
123782412ba7SHong Zhang         pj     = pj_oth + pi_oth[row];
123882412ba7SHong Zhang         pa     = pa_oth + pi_oth[row];
1239e900a757SHong Zhang         valtmp = aoa[j];
1240f4a743e1SHong Zhang         nextp  = 0;
124182412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
124282412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
1243e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
124442c57489SHong Zhang           }
124542c57489SHong Zhang         }
1246dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
124742c57489SHong Zhang       }
1248b091e509SHong Zhang #if defined(PTAP_PROFILE)
12498563dfccSBarry Smith       ierr    = PetscTime(&t2_1);CHKERRQ(ierr);
1250b091e509SHong Zhang       et2_AP += t2_1 - t2_0;
1251b091e509SHong Zhang #endif
125242c57489SHong Zhang 
1253a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
1254a9d06231SHong Zhang       /*--------------------------------------------------------------*/
125582412ba7SHong Zhang       pnz = pi_loc[i+1] - pi_loc[i];
1256f9473cf7SHong Zhang       pJ  = pj_loc + pi_loc[i];
125782412ba7SHong Zhang       for (j=0; j<pnz; j++) {
125842c57489SHong Zhang         nextap = 0;
1259f9473cf7SHong Zhang         row    = pJ[j]; /* global index */
126082412ba7SHong Zhang         if (row < pcstart || row >=pcend) { /* put the value into Co */
1261e900a757SHong Zhang           row = *poJ;
1262e900a757SHong Zhang           cj  = coj + coi[row];
1263e900a757SHong Zhang           ca  = coa + coi[row]; poJ++;
126482412ba7SHong Zhang         } else {                            /* put the value into Cd */
1265e900a757SHong Zhang           row = *pdJ;
1266e900a757SHong Zhang           cj  = bj + bi[row];
1267e900a757SHong Zhang           ca  = ba + bi[row]; pdJ++;
126842c57489SHong Zhang         }
1269e900a757SHong Zhang         valtmp = pA[j];
127082412ba7SHong Zhang         for (k=0; nextap<apnz; k++) {
1271e900a757SHong Zhang           if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++];
127242c57489SHong Zhang         }
1273dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
1274de0260b3SHong Zhang       }
12750496f32cSHong Zhang       pA += pnz;
127642c57489SHong Zhang       /* zero the current row info for A*P */
127782412ba7SHong Zhang       ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr);
1278b091e509SHong Zhang #if defined(PTAP_PROFILE)
12798563dfccSBarry Smith       ierr      = PetscTime(&t2_2);CHKERRQ(ierr);
128005d62848SHong Zhang       ePtAP += t2_2 - t2_1;
1281b091e509SHong Zhang #endif
128242c57489SHong Zhang     }
1283d9824c15SHong Zhang   }
128438571addSHong Zhang #if defined(PTAP_PROFILE)
12858563dfccSBarry Smith   ierr = PetscTime(&t2);CHKERRQ(ierr);
128638571addSHong Zhang #endif
128742c57489SHong Zhang 
1288a9d06231SHong Zhang   /* 3) send and recv matrix values coa */
1289a9d06231SHong Zhang   /*------------------------------------*/
129042c57489SHong Zhang   buf_ri = merge->buf_ri;
129142c57489SHong Zhang   buf_rj = merge->buf_rj;
129242c57489SHong Zhang   len_s  = merge->len_s;
1293fc42d0c8SSatish Balay   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
129442c57489SHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
129542c57489SHong Zhang 
1296dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
129742c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
129842c57489SHong Zhang     if (!len_s[proc]) continue;
1299de0260b3SHong Zhang     i    = merge->owners_co[proc];
1300de0260b3SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
130142c57489SHong Zhang     k++;
130242c57489SHong Zhang   }
13030c468ba9SBarry Smith   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
13040c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
130542c57489SHong Zhang 
1306a9d06231SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
130742c57489SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
130882412ba7SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
130938571addSHong Zhang #if defined(PTAP_PROFILE)
13108563dfccSBarry Smith   ierr = PetscTime(&t3);CHKERRQ(ierr);
131138571addSHong Zhang #endif
131242c57489SHong Zhang 
1313a9d06231SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1314a9d06231SHong Zhang   /*------------------------------------------------------*/
1315dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
131642c57489SHong Zhang   for (k=0; k<merge->nrecv; k++) {
131742c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
131842c57489SHong Zhang     nrows       = *(buf_ri_k[k]);
131942c57489SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
132082412ba7SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
132142c57489SHong Zhang   }
132242c57489SHong Zhang 
1323de0260b3SHong Zhang   for (i=0; i<cm; i++) {
132482412ba7SHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
132542c57489SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1326de0260b3SHong Zhang     ba_i = ba + bi[i];
132782412ba7SHong Zhang     bnz  = bi[i+1] - bi[i];
132842c57489SHong Zhang     /* add received vals into ba */
132942c57489SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
133042c57489SHong Zhang       /* i-th row */
133142c57489SHong Zhang       if (i == *nextrow[k]) {
133282412ba7SHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
133382412ba7SHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
133482412ba7SHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
133582412ba7SHong Zhang         nextcj = 0;
133682412ba7SHong Zhang         for (j=0; nextcj<cnz; j++) {
133782412ba7SHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
133882412ba7SHong Zhang             ba_i[j] += ca[nextcj++];
133942c57489SHong Zhang           }
134042c57489SHong Zhang         }
134182412ba7SHong Zhang         nextrow[k]++; nextci[k]++;
1342c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
134342c57489SHong Zhang       }
134442c57489SHong Zhang     }
134582412ba7SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
134642c57489SHong Zhang   }
134742c57489SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
134842c57489SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
134942c57489SHong Zhang 
1350de0260b3SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1351c05d87d6SBarry Smith   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
135242c57489SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
13530572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
135438571addSHong Zhang #if defined(PTAP_PROFILE)
13558563dfccSBarry Smith   ierr = PetscTime(&t4);CHKERRQ(ierr);
13569516a85cSHong Zhang   if (rank==1) {
135705d62848SHong 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);
13589516a85cSHong Zhang   }
135938571addSHong Zhang #endif
136042c57489SHong Zhang   PetscFunctionReturn(0);
136142c57489SHong Zhang }
1362