xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision ae5f08671782159de19e7bb960506b6c59548140)
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);}
46a560ef98SHong Zhang 
478403a639SHong Zhang     if (merge) { /* used by alg_ptap */
48cc6cb787SHong Zhang       ierr = PetscFree(merge->id_r);CHKERRQ(ierr);
49cc6cb787SHong Zhang       ierr = PetscFree(merge->len_s);CHKERRQ(ierr);
50cc6cb787SHong Zhang       ierr = PetscFree(merge->len_r);CHKERRQ(ierr);
51cc6cb787SHong Zhang       ierr = PetscFree(merge->bi);CHKERRQ(ierr);
52cc6cb787SHong Zhang       ierr = PetscFree(merge->bj);CHKERRQ(ierr);
53c05d87d6SBarry Smith       ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr);
54cc6cb787SHong Zhang       ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr);
55c05d87d6SBarry Smith       ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr);
56cc6cb787SHong Zhang       ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr);
57445158ffSHong Zhang       ierr = PetscFree(merge->coi);CHKERRQ(ierr);
58445158ffSHong Zhang       ierr = PetscFree(merge->coj);CHKERRQ(ierr);
5905b42c5fSBarry Smith       ierr = PetscFree(merge->owners_co);CHKERRQ(ierr);
606bf464f9SBarry Smith       ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr);
61f8487c73SHong Zhang       ierr = PetscFree(ptap->merge);CHKERRQ(ierr);
62bf0cc555SLisandro Dalcin     }
63a560ef98SHong Zhang 
64a560ef98SHong Zhang     ierr = ptap->destroy(A);CHKERRQ(ierr);
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"
87150d2497SBarry Smith PETSC_INTERN 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);
966c4ed002SBarry 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);
976c4ed002SBarry 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 
99c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,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 
115*ae5f0867Sstefano_zampini #if defined(PETSC_HAVE_HYPRE)
116*ae5f0867Sstefano_zampini PETSC_INTERN PetscErrorCode MatPtAPSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat*);
117*ae5f0867Sstefano_zampini #endif
118*ae5f0867Sstefano_zampini 
11942c57489SHong Zhang #undef __FUNCT__
120aa690a28SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ"
121aa690a28SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
1220d3441aeSHong Zhang {
1230d3441aeSHong Zhang   PetscErrorCode      ierr;
1240d3441aeSHong Zhang   Mat_PtAPMPI         *ptap;
125681d504bSHong Zhang   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
1262259aa2eSHong Zhang   MPI_Comm            comm;
1272259aa2eSHong Zhang   PetscMPIInt         size,rank;
12815a3b8e2SHong Zhang   Mat                 Cmpi;
12915a3b8e2SHong Zhang   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
130d0e9a2f7SHong Zhang   PetscInt            am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
131d0e9a2f7SHong Zhang   PetscInt            *lnk,i,k,pnz,row,nsend;
13215a3b8e2SHong Zhang   PetscBT             lnkbt;
133f671be37SHong Zhang   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
13415a3b8e2SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
135d0e9a2f7SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners,nzi,nspacedouble;
13615a3b8e2SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
13715a3b8e2SHong Zhang   MPI_Request         *swaits,*rwaits;
13815a3b8e2SHong Zhang   MPI_Status          *sstatus,rstatus;
139445158ffSHong Zhang   PetscLayout         rowmap;
140445158ffSHong Zhang   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
141445158ffSHong Zhang   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
142936ad1eaSHong Zhang   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0,Crmax,*aj,*ai,*pi;
14367a12041SHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*c_loc,*c_oth;
14467a12041SHong Zhang   PetscScalar         *apv;
14578d30b94SHong Zhang   PetscTable          ta;
146*ae5f0867Sstefano_zampini #if defined(PETSC_HAVE_HYPRE)
147*ae5f0867Sstefano_zampini   const char          *algTypes[3] = {"scalable","nonscalable","hypre"};
148*ae5f0867Sstefano_zampini   PetscInt            nalg = 3;
149*ae5f0867Sstefano_zampini #else
15078d30b94SHong Zhang   const char          *algTypes[2] = {"scalable","nonscalable"};
151*ae5f0867Sstefano_zampini   PetscInt            nalg = 2;
152*ae5f0867Sstefano_zampini #endif
15378d30b94SHong Zhang   PetscInt            alg = 1; /* set default algorithm */
154aa690a28SHong Zhang #if defined(PETSC_USE_INFO)
1558cb82516SHong Zhang   PetscReal           apfill;
156aa690a28SHong Zhang #endif
157681d504bSHong Zhang #if defined(PTAP_PROFILE)
158681d504bSHong Zhang   PetscLogDouble      t0,t1,t11,t12,t2,t3,t4;
159681d504bSHong Zhang #endif
16067a12041SHong Zhang 
16167a12041SHong Zhang   PetscFunctionBegin;
16267a12041SHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
16367a12041SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
16467a12041SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
165*ae5f0867Sstefano_zampini 
166*ae5f0867Sstefano_zampini   /* create symbolic parallel matrix Cmpi */
167*ae5f0867Sstefano_zampini   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
168*ae5f0867Sstefano_zampini   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
169*ae5f0867Sstefano_zampini 
170*ae5f0867Sstefano_zampini   ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
171*ae5f0867Sstefano_zampini   ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[1],&alg,NULL);CHKERRQ(ierr);
172*ae5f0867Sstefano_zampini   ierr = PetscOptionsEnd();CHKERRQ(ierr);
173*ae5f0867Sstefano_zampini 
174*ae5f0867Sstefano_zampini   if (alg == 1) {
175*ae5f0867Sstefano_zampini     /* Do dense axpy in MatPtAPNumeric_MPIAIJ_MPIAIJ() */
176*ae5f0867Sstefano_zampini     Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ;
177*ae5f0867Sstefano_zampini   } else if (alg == 0) {
178*ae5f0867Sstefano_zampini     Cmpi->ops->ptapnumeric = 0; /* not done yet */
179*ae5f0867Sstefano_zampini     SETERRQ(comm,PETSC_ERR_ARG_SIZ,"Not done yet");
180*ae5f0867Sstefano_zampini #if defined(PETSC_HAVE_HYPRE)
181*ae5f0867Sstefano_zampini   } else if (alg == 2) {
182*ae5f0867Sstefano_zampini     /* Use boomerAMGBuildCoarseOperator */
183*ae5f0867Sstefano_zampini     ierr = MatDestroy(&Cmpi);CHKERRQ(ierr);
184*ae5f0867Sstefano_zampini     ierr = MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);CHKERRQ(ierr);
185*ae5f0867Sstefano_zampini     PetscFunctionReturn(0);
186*ae5f0867Sstefano_zampini #endif
187*ae5f0867Sstefano_zampini   } else {
188*ae5f0867Sstefano_zampini     SETERRQ1(comm,PETSC_ERR_PLIB,"Unknown algorithm number %D\n",alg);
189*ae5f0867Sstefano_zampini   }
190*ae5f0867Sstefano_zampini 
1918cb82516SHong Zhang #if defined(PTAP_PROFILE)
19280bb4639SHong Zhang   ierr = PetscTime(&t0);CHKERRQ(ierr);
1938cb82516SHong Zhang #endif
1948cb82516SHong Zhang 
19515a3b8e2SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
19615a3b8e2SHong Zhang   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
19715a3b8e2SHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
19815a3b8e2SHong Zhang 
19915a3b8e2SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
20015a3b8e2SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
20115a3b8e2SHong Zhang   /* get P_loc by taking all local rows of P */
20215a3b8e2SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
20315a3b8e2SHong Zhang 
20467a12041SHong Zhang   /* (0) compute Rd = Pd^T, Ro = Po^T  */
20567a12041SHong Zhang   /* --------------------------------- */
206de817e96SHong Zhang   ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
207de817e96SHong Zhang   ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
2088cb82516SHong Zhang #if defined(PTAP_PROFILE)
209445158ffSHong Zhang   ierr = PetscTime(&t11);CHKERRQ(ierr);
2108cb82516SHong Zhang #endif
21115a3b8e2SHong Zhang 
21267a12041SHong Zhang   /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */
21367a12041SHong Zhang   /* ----------------------------------------------------------------- */
21467a12041SHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
21567a12041SHong Zhang   p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
216445158ffSHong Zhang 
2179a6dcab7SHong Zhang   /* create and initialize a linked list */
218c373ccc6SHong Zhang   ierr = PetscTableCreate(pN,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */
2194b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta);
2204b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta);
22178d30b94SHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */
222d0e9a2f7SHong 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); */
22378d30b94SHong Zhang 
22478d30b94SHong Zhang   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
225445158ffSHong Zhang 
2268cb82516SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
227f91af8c7SBarry Smith   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr);
228445158ffSHong Zhang   current_space = free_space;
22967a12041SHong Zhang   nspacedouble  = 0;
23067a12041SHong Zhang 
23167a12041SHong Zhang   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
23267a12041SHong Zhang   api[0] = 0;
233445158ffSHong Zhang   for (i=0; i<am; i++) {
23467a12041SHong Zhang     /* diagonal portion: Ad[i,:]*P */
23567a12041SHong Zhang     ai = ad->i; pi = p_loc->i;
23667a12041SHong Zhang     nzi = ai[i+1] - ai[i];
23767a12041SHong Zhang     aj  = ad->j + ai[i];
238445158ffSHong Zhang     for (j=0; j<nzi; j++) {
239445158ffSHong Zhang       row  = aj[j];
24067a12041SHong Zhang       pnz  = pi[row+1] - pi[row];
24167a12041SHong Zhang       Jptr = p_loc->j + pi[row];
242445158ffSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
243445158ffSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
244445158ffSHong Zhang     }
24567a12041SHong Zhang     /* off-diagonal portion: Ao[i,:]*P */
24667a12041SHong Zhang     ai = ao->i; pi = p_oth->i;
24767a12041SHong Zhang     nzi = ai[i+1] - ai[i];
24867a12041SHong Zhang     aj  = ao->j + ai[i];
249445158ffSHong Zhang     for (j=0; j<nzi; j++) {
250445158ffSHong Zhang       row  = aj[j];
25167a12041SHong Zhang       pnz  = pi[row+1] - pi[row];
25267a12041SHong Zhang       Jptr = p_oth->j + pi[row];
253445158ffSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
254445158ffSHong Zhang     }
255445158ffSHong Zhang     apnz     = lnk[0];
256445158ffSHong Zhang     api[i+1] = api[i] + apnz;
257445158ffSHong Zhang     if (ap_rmax < apnz) ap_rmax = apnz;
258445158ffSHong Zhang 
259445158ffSHong Zhang     /* if free space is not available, double the total space in the list */
260445158ffSHong Zhang     if (current_space->local_remaining<apnz) {
261f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
262445158ffSHong Zhang       nspacedouble++;
263445158ffSHong Zhang     }
264445158ffSHong Zhang 
265445158ffSHong Zhang     /* Copy data into free space, then initialize lnk */
266445158ffSHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
267445158ffSHong Zhang 
268445158ffSHong Zhang     current_space->array           += apnz;
269445158ffSHong Zhang     current_space->local_used      += apnz;
270445158ffSHong Zhang     current_space->local_remaining -= apnz;
271445158ffSHong Zhang   }
272681d504bSHong Zhang   /* Allocate space for apj and apv, initialize apj, and */
273445158ffSHong Zhang   /* destroy list of free space and other temporary array(s) */
274445158ffSHong Zhang   ierr   = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr);
275445158ffSHong Zhang   ierr   = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
2769a6dcab7SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
2779a6dcab7SHong Zhang 
278aa690a28SHong Zhang   /* Create AP_loc for reuse */
279445158ffSHong Zhang   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr);
280aa690a28SHong Zhang 
281aa690a28SHong Zhang #if defined(PETSC_USE_INFO)
282aa690a28SHong Zhang   apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1);
283aa690a28SHong Zhang   ptap->AP_loc->info.mallocs           = nspacedouble;
284aa690a28SHong Zhang   ptap->AP_loc->info.fill_ratio_given  = fill;
285aa690a28SHong Zhang   ptap->AP_loc->info.fill_ratio_needed = apfill;
286aa690a28SHong Zhang 
287aa690a28SHong Zhang   if (api[am]) {
288aa690a28SHong Zhang     ierr = PetscInfo3(ptap->AP_loc,"AP_loc reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)apfill);CHKERRQ(ierr);
289aa690a28SHong Zhang     ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr);
290aa690a28SHong Zhang   } else {
291aa690a28SHong Zhang     ierr = PetscInfo(ptap->AP_loc,"AP_loc is empty \n");CHKERRQ(ierr);
292aa690a28SHong Zhang   }
293aa690a28SHong Zhang #endif
294aa690a28SHong Zhang 
2958cb82516SHong Zhang #if defined(PTAP_PROFILE)
296445158ffSHong Zhang   ierr = PetscTime(&t12);CHKERRQ(ierr);
2978cb82516SHong Zhang #endif
29815a3b8e2SHong Zhang 
299681d504bSHong Zhang   /* (2-1) compute symbolic Co = Ro*AP_loc  */
300681d504bSHong Zhang   /* ------------------------------------ */
30167a12041SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
3028cb82516SHong Zhang #if defined(PTAP_PROFILE)
30380bb4639SHong Zhang   ierr = PetscTime(&t1);CHKERRQ(ierr);
3048cb82516SHong Zhang #endif
30515a3b8e2SHong Zhang 
30667a12041SHong Zhang   /* (3) send coj of C_oth to other processors  */
30767a12041SHong Zhang   /* ------------------------------------------ */
30815a3b8e2SHong Zhang   /* determine row ownership */
309445158ffSHong Zhang   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
310445158ffSHong Zhang   rowmap->n  = pn;
311445158ffSHong Zhang   rowmap->bs = 1;
312445158ffSHong Zhang   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
313445158ffSHong Zhang   owners = rowmap->range;
31415a3b8e2SHong Zhang 
31515a3b8e2SHong Zhang   /* determine the number of messages to send, their lengths */
3168cb82516SHong Zhang   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
3178cb82516SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
31815a3b8e2SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
31915a3b8e2SHong Zhang 
32067a12041SHong Zhang   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
32167a12041SHong Zhang   coi   = c_oth->i; coj = c_oth->j;
32267a12041SHong Zhang   con   = ptap->C_oth->rmap->n;
32315a3b8e2SHong Zhang   proc  = 0;
32467a12041SHong Zhang   for (i=0; i<con; i++) {
32515a3b8e2SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
32615a3b8e2SHong Zhang     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
32715a3b8e2SHong Zhang     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
32815a3b8e2SHong Zhang   }
32915a3b8e2SHong Zhang 
33015a3b8e2SHong Zhang   len          = 0; /* max length of buf_si[], see (4) */
33115a3b8e2SHong Zhang   owners_co[0] = 0;
33267a12041SHong Zhang   nsend        = 0;
33315a3b8e2SHong Zhang   for (proc=0; proc<size; proc++) {
33415a3b8e2SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
33515a3b8e2SHong Zhang     if (len_s[proc]) {
336445158ffSHong Zhang       nsend++;
33715a3b8e2SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
33815a3b8e2SHong Zhang       len         += len_si[proc];
33915a3b8e2SHong Zhang     }
34015a3b8e2SHong Zhang   }
34115a3b8e2SHong Zhang 
34215a3b8e2SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
343445158ffSHong Zhang   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
344445158ffSHong Zhang   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
34515a3b8e2SHong Zhang 
34615a3b8e2SHong Zhang   /* post the Irecv and Isend of coj */
34715a3b8e2SHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
348445158ffSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
349445158ffSHong Zhang   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
35015a3b8e2SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
35115a3b8e2SHong Zhang     if (!len_s[proc]) continue;
35215a3b8e2SHong Zhang     i    = owners_co[proc];
35315a3b8e2SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
35415a3b8e2SHong Zhang     k++;
35515a3b8e2SHong Zhang   }
35615a3b8e2SHong Zhang 
357681d504bSHong Zhang   /* (2-2) compute symbolic C_loc = Rd*AP_loc */
358681d504bSHong Zhang   /* ---------------------------------------- */
359681d504bSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr);
360681d504bSHong Zhang   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
361681d504bSHong Zhang 
362681d504bSHong Zhang   /* receives coj are complete */
363445158ffSHong Zhang   for (i=0; i<nrecv; i++) {
364445158ffSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
36515a3b8e2SHong Zhang   }
36615a3b8e2SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
367445158ffSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
36815a3b8e2SHong Zhang 
36978d30b94SHong Zhang   /* add received column indices into ta to update Crmax */
37078d30b94SHong Zhang   for (k=0; k<nrecv; k++) {/* k-th received message */
37178d30b94SHong Zhang     Jptr = buf_rj[k];
37278d30b94SHong Zhang     for (j=0; j<len_r[k]; j++) {
37378d30b94SHong Zhang       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
37478d30b94SHong Zhang     }
37578d30b94SHong Zhang   }
37678d30b94SHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
37778d30b94SHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
3789a6dcab7SHong Zhang 
37915a3b8e2SHong Zhang   /* (4) send and recv coi */
38015a3b8e2SHong Zhang   /*-----------------------*/
38115a3b8e2SHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
382445158ffSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
38315a3b8e2SHong Zhang   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
38415a3b8e2SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
38515a3b8e2SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
38615a3b8e2SHong Zhang     if (!len_s[proc]) continue;
38715a3b8e2SHong Zhang     /* form outgoing message for i-structure:
38815a3b8e2SHong Zhang          buf_si[0]:                 nrows to be sent
38915a3b8e2SHong Zhang                [1:nrows]:           row index (global)
39015a3b8e2SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
39115a3b8e2SHong Zhang     */
39215a3b8e2SHong Zhang     /*-------------------------------------------*/
39315a3b8e2SHong Zhang     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
39415a3b8e2SHong Zhang     buf_si_i    = buf_si + nrows+1;
39515a3b8e2SHong Zhang     buf_si[0]   = nrows;
39615a3b8e2SHong Zhang     buf_si_i[0] = 0;
39715a3b8e2SHong Zhang     nrows       = 0;
39815a3b8e2SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
39915a3b8e2SHong Zhang       nzi = coi[i+1] - coi[i];
40015a3b8e2SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
40115a3b8e2SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
40215a3b8e2SHong Zhang       nrows++;
40315a3b8e2SHong Zhang     }
40415a3b8e2SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
40515a3b8e2SHong Zhang     k++;
40615a3b8e2SHong Zhang     buf_si += len_si[proc];
40715a3b8e2SHong Zhang   }
408681d504bSHong Zhang   for (i=0; i<nrecv; i++) {
409445158ffSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
41015a3b8e2SHong Zhang   }
41115a3b8e2SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
412445158ffSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
41315a3b8e2SHong Zhang 
4148cb82516SHong Zhang   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
41515a3b8e2SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
41615a3b8e2SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
41715a3b8e2SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
4188cb82516SHong Zhang #if defined(PTAP_PROFILE)
41980bb4639SHong Zhang   ierr = PetscTime(&t2);CHKERRQ(ierr);
4208cb82516SHong Zhang #endif
42167a12041SHong Zhang   /* (5) compute the local portion of Cmpi      */
42267a12041SHong Zhang   /* ------------------------------------------ */
42378d30b94SHong Zhang   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
42478d30b94SHong Zhang   ierr          = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr);
42515a3b8e2SHong Zhang   current_space = free_space;
42615a3b8e2SHong Zhang 
427445158ffSHong Zhang   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
428445158ffSHong Zhang   for (k=0; k<nrecv; k++) {
42915a3b8e2SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
43015a3b8e2SHong Zhang     nrows       = *buf_ri_k[k];
43115a3b8e2SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
43215a3b8e2SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
43315a3b8e2SHong Zhang   }
434445158ffSHong Zhang 
43578d30b94SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
43678d30b94SHong Zhang   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
43715a3b8e2SHong Zhang   for (i=0; i<pn; i++) {
43867a12041SHong Zhang     /* add C_loc into Cmpi */
43967a12041SHong Zhang     nzi  = c_loc->i[i+1] - c_loc->i[i];
44067a12041SHong Zhang     Jptr = c_loc->j + c_loc->i[i];
44167a12041SHong Zhang     ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
44215a3b8e2SHong Zhang 
44315a3b8e2SHong Zhang     /* add received col data into lnk */
444445158ffSHong Zhang     for (k=0; k<nrecv; k++) { /* k-th received message */
44515a3b8e2SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
44615a3b8e2SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
44715a3b8e2SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
44815a3b8e2SHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
44915a3b8e2SHong Zhang         nextrow[k]++; nextci[k]++;
45015a3b8e2SHong Zhang       }
45115a3b8e2SHong Zhang     }
452d0e9a2f7SHong Zhang     nzi = lnk[0];
4538cb82516SHong Zhang 
45415a3b8e2SHong Zhang     /* copy data into free space, then initialize lnk */
455d0e9a2f7SHong Zhang     ierr = PetscLLCondensedClean(pN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
456d0e9a2f7SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr);
45715a3b8e2SHong Zhang   }
45815a3b8e2SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
45915a3b8e2SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
460445158ffSHong Zhang   ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr);
4618cb82516SHong Zhang #if defined(PTAP_PROFILE)
46280bb4639SHong Zhang   ierr = PetscTime(&t3);CHKERRQ(ierr);
4638cb82516SHong Zhang #endif
46480bb4639SHong Zhang 
465*ae5f0867Sstefano_zampini   /* local sizes and preallocation */
46615a3b8e2SHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
46715a3b8e2SHong Zhang   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
46815a3b8e2SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
46915a3b8e2SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
47015a3b8e2SHong Zhang 
471445158ffSHong Zhang   /* members in merge */
472445158ffSHong Zhang   ierr = PetscFree(id_r);CHKERRQ(ierr);
473445158ffSHong Zhang   ierr = PetscFree(len_r);CHKERRQ(ierr);
474445158ffSHong Zhang   ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr);
475445158ffSHong Zhang   ierr = PetscFree(buf_ri);CHKERRQ(ierr);
476445158ffSHong Zhang   ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr);
477445158ffSHong Zhang   ierr = PetscFree(buf_rj);CHKERRQ(ierr);
478445158ffSHong Zhang   ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr);
47915a3b8e2SHong Zhang 
48015a3b8e2SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
48115a3b8e2SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
48215a3b8e2SHong Zhang   c->ptap         = ptap;
4831a47ec54SHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
4841a47ec54SHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
48515a3b8e2SHong Zhang 
48678d30b94SHong Zhang   if (alg == 1) {
4878cb82516SHong Zhang     ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
48878d30b94SHong Zhang   }
4892259aa2eSHong Zhang 
4901a47ec54SHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
4911a47ec54SHong Zhang   Cmpi->assembled        = PETSC_FALSE;
4921a47ec54SHong Zhang   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
493aa690a28SHong Zhang   Cmpi->ops->duplicate   = MatDuplicate_MPIAIJ_MatPtAP;
4941a47ec54SHong Zhang   *C                     = Cmpi;
4951a47ec54SHong Zhang 
4968cb82516SHong Zhang #if defined(PTAP_PROFILE)
49780bb4639SHong Zhang   ierr = PetscTime(&t4);CHKERRQ(ierr);
498e9c1f85fSHong Zhang   if (rank == 1) {
499445158ffSHong 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);
500e9c1f85fSHong Zhang   }
5018cb82516SHong Zhang #endif
5020d3441aeSHong Zhang   PetscFunctionReturn(0);
5030d3441aeSHong Zhang }
5040d3441aeSHong Zhang 
5050d3441aeSHong Zhang #undef __FUNCT__
506aa690a28SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ"
507aa690a28SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
5080d3441aeSHong Zhang {
5090d3441aeSHong Zhang   PetscErrorCode    ierr;
5100d3441aeSHong Zhang   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
5110d3441aeSHong Zhang   Mat_SeqAIJ        *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
5128cb82516SHong Zhang   Mat_SeqAIJ        *ap,*p_loc,*p_oth,*c_seq;
5130d3441aeSHong Zhang   Mat_PtAPMPI       *ptap = c->ptap;
5149ce11a7cSHong Zhang   Mat               AP_loc,C_loc,C_oth;
5150d3441aeSHong Zhang   PetscInt          i,rstart,rend,cm,ncols,row;
5168cb82516SHong Zhang   PetscInt          *api,*apj,am = A->rmap->n,j,col,apnz;
5178cb82516SHong Zhang   PetscScalar       *apa;
5180d3441aeSHong Zhang   const PetscInt    *cols;
5190d3441aeSHong Zhang   const PetscScalar *vals;
5208cb82516SHong Zhang #if defined(PTAP_PROFILE)
5218cb82516SHong Zhang   PetscMPIInt       rank;
5228cb82516SHong Zhang   MPI_Comm          comm;
5230d3441aeSHong Zhang   PetscLogDouble    t0,t1,t2,t3,t4,eR,eAP,eCseq,eCmpi;
5248cb82516SHong Zhang #endif
5250d3441aeSHong Zhang 
5260d3441aeSHong Zhang   PetscFunctionBegin;
5270d3441aeSHong Zhang   ierr = MatZeroEntries(C);CHKERRQ(ierr);
5280d3441aeSHong Zhang 
529e497d3c8SHong Zhang   /* 1) get R = Pd^T,Ro = Po^T */
5308cb82516SHong Zhang #if defined(PTAP_PROFILE)
5310d3441aeSHong Zhang   ierr = PetscTime(&t0);CHKERRQ(ierr);
5328cb82516SHong Zhang #endif
533748c7196SHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
5340d3441aeSHong Zhang     ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
5350d3441aeSHong Zhang     ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
536748c7196SHong Zhang   }
5378cb82516SHong Zhang #if defined(PTAP_PROFILE)
5380d3441aeSHong Zhang   ierr = PetscTime(&t1);CHKERRQ(ierr);
5390d3441aeSHong Zhang   eR = t1 - t0;
5408cb82516SHong Zhang #endif
5410d3441aeSHong Zhang 
542e497d3c8SHong Zhang   /* 2) get AP_loc */
5430d3441aeSHong Zhang   AP_loc = ptap->AP_loc;
5448cb82516SHong Zhang   ap = (Mat_SeqAIJ*)AP_loc->data;
5450d3441aeSHong Zhang 
546e497d3c8SHong Zhang   /* 2-1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
5470d3441aeSHong Zhang   /*-----------------------------------------------------*/
548748c7196SHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
549748c7196SHong Zhang     /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */
5500d3441aeSHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
5510d3441aeSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
552e497d3c8SHong Zhang   }
5530d3441aeSHong Zhang 
5548cb82516SHong Zhang   /* 2-2) compute numeric A_loc*P - dominating part */
5558cb82516SHong Zhang   /* ---------------------------------------------- */
5560d3441aeSHong Zhang   /* get data from symbolic products */
5578cb82516SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
5588cb82516SHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
5598cb82516SHong Zhang   apa   = ptap->apa;
560681d504bSHong Zhang   api   = ap->i;
561681d504bSHong Zhang   apj   = ap->j;
562e497d3c8SHong Zhang   for (i=0; i<am; i++) {
563445158ffSHong Zhang     /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */
564e497d3c8SHong Zhang     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
565e497d3c8SHong Zhang     apnz = api[i+1] - api[i];
566e497d3c8SHong Zhang     for (j=0; j<apnz; j++) {
567e497d3c8SHong Zhang       col = apj[j+api[i]];
568e497d3c8SHong Zhang       ap->a[j+ap->i[i]] = apa[col];
569e497d3c8SHong Zhang       apa[col] = 0.0;
5700d3441aeSHong Zhang     }
571aa690a28SHong Zhang     ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
5720d3441aeSHong Zhang   }
5738cb82516SHong Zhang #if defined(PTAP_PROFILE)
5740d3441aeSHong Zhang   ierr = PetscTime(&t2);CHKERRQ(ierr);
5750d3441aeSHong Zhang   eAP = t2 - t1;
5768cb82516SHong Zhang #endif
5770d3441aeSHong Zhang 
5788cb82516SHong Zhang   /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */
579445158ffSHong Zhang   ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr);
580445158ffSHong Zhang   ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr);
5819ce11a7cSHong Zhang   C_loc = ptap->C_loc;
5829ce11a7cSHong Zhang   C_oth = ptap->C_oth;
5838cb82516SHong Zhang #if defined(PTAP_PROFILE)
5840d3441aeSHong Zhang   ierr = PetscTime(&t3);CHKERRQ(ierr);
5850d3441aeSHong Zhang   eCseq = t3 - t2;
5868cb82516SHong Zhang #endif
5870d3441aeSHong Zhang 
5880d3441aeSHong Zhang   /* add C_loc and Co to to C */
5890d3441aeSHong Zhang   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
5900d3441aeSHong Zhang 
5910d3441aeSHong Zhang   /* C_loc -> C */
5920d3441aeSHong Zhang   cm    = C_loc->rmap->N;
5939ce11a7cSHong Zhang   c_seq = (Mat_SeqAIJ*)C_loc->data;
5948cb82516SHong Zhang   cols = c_seq->j;
5958cb82516SHong Zhang   vals = c_seq->a;
5960d3441aeSHong Zhang   for (i=0; i<cm; i++) {
5979ce11a7cSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
5980d3441aeSHong Zhang     row = rstart + i;
5990d3441aeSHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
6008cb82516SHong Zhang     cols += ncols; vals += ncols;
6010d3441aeSHong Zhang   }
6020d3441aeSHong Zhang 
6030d3441aeSHong Zhang   /* Co -> C, off-processor part */
6049ce11a7cSHong Zhang   cm = C_oth->rmap->N;
6059ce11a7cSHong Zhang   c_seq = (Mat_SeqAIJ*)C_oth->data;
6068cb82516SHong Zhang   cols = c_seq->j;
6078cb82516SHong Zhang   vals = c_seq->a;
6089ce11a7cSHong Zhang   for (i=0; i<cm; i++) {
6099ce11a7cSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
6100d3441aeSHong Zhang     row = p->garray[i];
6110d3441aeSHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
6128cb82516SHong Zhang     cols += ncols; vals += ncols;
6130d3441aeSHong Zhang   }
6140d3441aeSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6150d3441aeSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6168cb82516SHong Zhang #if defined(PTAP_PROFILE)
6170d3441aeSHong Zhang   ierr = PetscTime(&t4);CHKERRQ(ierr);
6180d3441aeSHong Zhang   eCmpi = t4 - t3;
6190d3441aeSHong Zhang 
6208cb82516SHong Zhang   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
6218cb82516SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
6220d3441aeSHong Zhang   if (rank==1) {
6230d3441aeSHong 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);
6240d3441aeSHong Zhang   }
6258cb82516SHong Zhang #endif
626748c7196SHong Zhang   ptap->reuse = MAT_REUSE_MATRIX;
6270d3441aeSHong Zhang   PetscFunctionReturn(0);
6280d3441aeSHong Zhang }
6290d3441aeSHong Zhang 
6300d3441aeSHong Zhang #undef __FUNCT__
6318403a639SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap"
6328403a639SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(Mat A,Mat P,PetscReal fill,Mat *C)
63342c57489SHong Zhang {
63442c57489SHong Zhang   PetscErrorCode      ierr;
63577584fe6SHong Zhang   Mat                 Cmpi;
636f8487c73SHong Zhang   Mat_PtAPMPI         *ptap;
6370298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
638f8487c73SHong Zhang   Mat_MPIAIJ          *a        =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
63942c57489SHong Zhang   Mat_SeqAIJ          *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
64042c57489SHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth;
641ded4f561SHong Zhang   PetscInt            *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ;
64277584fe6SHong Zhang   PetscInt            *adi=ad->i,*aj,*aoi=ao->i,nnz;
643a914927fSHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
644d16ca5f0SHong Zhang   PetscInt            am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n;
64542c57489SHong Zhang   PetscBT             lnkbt;
646ce94432eSBarry Smith   MPI_Comm            comm;
647a914927fSHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0;
64842c57489SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
64942c57489SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
65024ecddacSHong Zhang   PetscInt            nzi,*pti,*ptj;
65142c57489SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
652ded4f561SHong Zhang   MPI_Request         *swaits,*rwaits;
653ded4f561SHong Zhang   MPI_Status          *sstatus,rstatus;
65442c57489SHong Zhang   Mat_Merge_SeqsToMPI *merge;
65577584fe6SHong Zhang   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0;
656d16ca5f0SHong Zhang   PetscReal           afill=1.0,afill_tmp;
65742c57489SHong Zhang 
65842c57489SHong Zhang   PetscFunctionBegin;
659ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
66083445d95SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
66183445d95SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
66283445d95SHong Zhang 
663f8487c73SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
664b00a9115SJed Brown   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
665b00a9115SJed Brown   ierr        = PetscNew(&merge);CHKERRQ(ierr);
6669f4155fbSHong Zhang   ptap->merge = merge;
667f8487c73SHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
6686c6e00e0SHong Zhang 
6696c6e00e0SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
670b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
671fe615159SHong Zhang 
6726c6e00e0SHong Zhang   /* get P_loc by taking all local rows of P */
673a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
6746c6e00e0SHong Zhang 
675a1a4e84aSHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
676a1a4e84aSHong Zhang   p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
6776c6e00e0SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
6786c6e00e0SHong Zhang   pi_oth = p_oth->i; pj_oth = p_oth->j;
67942c57489SHong Zhang 
6802addb229SHong Zhang   /* (1) compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth (api,apj) */
6812addb229SHong Zhang   /*--------------------------------------------------------------------------*/
682854ce69bSBarry Smith   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
68382412ba7SHong Zhang   api[0] = 0;
68483445d95SHong Zhang 
68583445d95SHong Zhang   /* create and initialize a linked list */
686a914927fSHong Zhang   ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr);
68783445d95SHong Zhang 
688a914927fSHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) -OOM for ex56, np=8k on Intrepid! */
689f91af8c7SBarry Smith   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr);
690f4a743e1SHong Zhang   current_space = free_space;
691f4a743e1SHong Zhang 
692f4a743e1SHong Zhang   for (i=0; i<am; i++) {
693f4a743e1SHong Zhang     /* diagonal portion of A */
694ded4f561SHong Zhang     nzi = adi[i+1] - adi[i];
69577584fe6SHong Zhang     aj  = ad->j + adi[i];
696ded4f561SHong Zhang     for (j=0; j<nzi; j++) {
69777584fe6SHong Zhang       row  = aj[j];
69882412ba7SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
699ded4f561SHong Zhang       Jptr = pj_loc + pi_loc[row];
70083445d95SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
701a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
702f4a743e1SHong Zhang     }
703f4a743e1SHong Zhang     /* off-diagonal portion of A */
704ded4f561SHong Zhang     nzi = aoi[i+1] - aoi[i];
70577584fe6SHong Zhang     aj  = ao->j + aoi[i];
706ded4f561SHong Zhang     for (j=0; j<nzi; j++) {
70777584fe6SHong Zhang       row  = aj[j];
70882412ba7SHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
709ded4f561SHong Zhang       Jptr = pj_oth + pi_oth[row];
710a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
711f4a743e1SHong Zhang     }
712a914927fSHong Zhang     apnz     = lnk[0];
71382412ba7SHong Zhang     api[i+1] = api[i] + apnz;
71477584fe6SHong Zhang     if (ap_rmax < apnz) ap_rmax = apnz;
715f4a743e1SHong Zhang 
71683445d95SHong Zhang     /* if free space is not available, double the total space in the list */
71782412ba7SHong Zhang     if (current_space->local_remaining<apnz) {
718f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
719f2b054eeSHong Zhang       nspacedouble++;
720f4a743e1SHong Zhang     }
721f4a743e1SHong Zhang 
722f4a743e1SHong Zhang     /* Copy data into free space, then initialize lnk */
723a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
7242205254eSKarl Rupp 
72582412ba7SHong Zhang     current_space->array           += apnz;
72682412ba7SHong Zhang     current_space->local_used      += apnz;
72782412ba7SHong Zhang     current_space->local_remaining -= apnz;
728f4a743e1SHong Zhang   }
729a914927fSHong Zhang 
73082412ba7SHong Zhang   /* Allocate space for apj, initialize apj, and */
731f4a743e1SHong Zhang   /* destroy list of free space and other temporary array(s) */
732854ce69bSBarry Smith   ierr      = PetscMalloc1(api[am]+1,&apj);CHKERRQ(ierr);
73377584fe6SHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
734118727c9SMark F. Adams   afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1);
735d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
736f4a743e1SHong Zhang 
7372addb229SHong Zhang   /* (2) determine symbolic Co=(p->B)^T*AP - send to others (coi,coj)*/
7382addb229SHong Zhang   /*-----------------------------------------------------------------*/
739de0260b3SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
74042c57489SHong Zhang 
741ded4f561SHong Zhang   /* then, compute symbolic Co = (p->B)^T*AP */
742d0f46423SBarry Smith   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
74383445d95SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
744854ce69bSBarry Smith   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
745de0260b3SHong Zhang   coi[0] = 0;
74642c57489SHong Zhang 
747d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */
748f91af8c7SBarry Smith   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(poti[pon],api[am]));
74922d28d08SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
75042c57489SHong Zhang   current_space = free_space;
75142c57489SHong Zhang 
752de0260b3SHong Zhang   for (i=0; i<pon; i++) {
75382412ba7SHong Zhang     pnz = poti[i+1] - poti[i];
75477584fe6SHong Zhang     ptJ = potj + poti[i];
75577584fe6SHong Zhang     for (j=0; j<pnz; j++) {
75677584fe6SHong Zhang       row  = ptJ[j]; /* row of AP == col of Pot */
75782412ba7SHong Zhang       apnz = api[row+1] - api[row];
758ded4f561SHong Zhang       Jptr = apj + api[row];
75983445d95SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
760a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
76142c57489SHong Zhang     }
762a914927fSHong Zhang     nnz = lnk[0];
76342c57489SHong Zhang 
76483445d95SHong Zhang     /* If free space is not available, double the total space in the list */
765ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
766f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
767d16ca5f0SHong Zhang       nspacedouble++;
76842c57489SHong Zhang     }
76942c57489SHong Zhang 
77042c57489SHong Zhang     /* Copy data into free space, and zero out denserows */
771a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
7722205254eSKarl Rupp 
773ded4f561SHong Zhang     current_space->array           += nnz;
774ded4f561SHong Zhang     current_space->local_used      += nnz;
775ded4f561SHong Zhang     current_space->local_remaining -= nnz;
7762205254eSKarl Rupp 
777ded4f561SHong Zhang     coi[i+1] = coi[i] + nnz;
77842c57489SHong Zhang   }
7796b911d16SHong Zhang 
7806b911d16SHong Zhang   ierr      = PetscMalloc1(coi[pon],&coj);CHKERRQ(ierr);
781a1a86e44SBarry Smith   ierr      = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
782118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1);
783d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
784de0260b3SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
78542c57489SHong Zhang 
7866b911d16SHong Zhang   /* (3) send j-array (coj) of Co to other processors */
7876b911d16SHong Zhang   /*--------------------------------------------------*/
7886b911d16SHong Zhang   ierr = PetscCalloc1(size,&merge->len_s);CHKERRQ(ierr);
7896b911d16SHong Zhang   len_s        = merge->len_s;
7906b911d16SHong Zhang   merge->nsend = 0;
7916b911d16SHong Zhang 
7926b911d16SHong Zhang 
79342c57489SHong Zhang   /* determine row ownership */
79426283091SBarry Smith   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
7957a2fc3feSBarry Smith   merge->rowmap->n  = pn;
7967a2fc3feSBarry Smith   merge->rowmap->bs = 1;
7972205254eSKarl Rupp 
79826283091SBarry Smith   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
7997a2fc3feSBarry Smith   owners = merge->rowmap->range;
80042c57489SHong Zhang 
80142c57489SHong Zhang   /* determine the number of messages to send, their lengths */
802dcca6d9dSJed Brown   ierr = PetscMalloc2(size,&len_si,size,&sstatus);CHKERRQ(ierr);
80383445d95SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
804854ce69bSBarry Smith   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
805de0260b3SHong Zhang 
80683445d95SHong Zhang   proc = 0;
807de0260b3SHong Zhang   for (i=0; i<pon; i++) {
808de0260b3SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
8092addb229SHong Zhang     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
810ee6e4b5bSHong Zhang     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
81183445d95SHong Zhang   }
812de0260b3SHong Zhang 
8132addb229SHong Zhang   len          = 0; /* max length of buf_si[], see (4) */
814de0260b3SHong Zhang   owners_co[0] = 0;
81542c57489SHong Zhang   for (proc=0; proc<size; proc++) {
816de0260b3SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
817ee6e4b5bSHong Zhang     if (len_s[proc]) {
81842c57489SHong Zhang       merge->nsend++;
8192addb229SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
82042c57489SHong Zhang       len         += len_si[proc];
82142c57489SHong Zhang     }
82242c57489SHong Zhang   }
82342c57489SHong Zhang 
824ded4f561SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
8250298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
82642c57489SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
82742c57489SHong Zhang 
828ded4f561SHong Zhang   /* post the Irecv and Isend of coj */
829529bc5dcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
830529bc5dcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
831854ce69bSBarry Smith   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
83242c57489SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
83342c57489SHong Zhang     if (!len_s[proc]) continue;
834de0260b3SHong Zhang     i    = owners_co[proc];
835529bc5dcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
83642c57489SHong Zhang     k++;
83742c57489SHong Zhang   }
83842c57489SHong Zhang 
839ded4f561SHong Zhang   /* receives and sends of coj are complete */
84077584fe6SHong Zhang   for (i=0; i<merge->nrecv; i++) {
841c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
842ded4f561SHong Zhang   }
843ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
8440c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
84582412ba7SHong Zhang 
8466b911d16SHong Zhang   /* (4) send and recv coi */
8476b911d16SHong Zhang   /*-----------------------*/
848529bc5dcSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
849529bc5dcSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
850854ce69bSBarry Smith   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
85142c57489SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
85242c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
85342c57489SHong Zhang     if (!len_s[proc]) continue;
85442c57489SHong Zhang     /* form outgoing message for i-structure:
85542c57489SHong Zhang          buf_si[0]:                 nrows to be sent
85642c57489SHong Zhang                [1:nrows]:           row index (global)
85742c57489SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
85842c57489SHong Zhang     */
85942c57489SHong Zhang     /*-------------------------------------------*/
8602addb229SHong Zhang     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
86142c57489SHong Zhang     buf_si_i    = buf_si + nrows+1;
86242c57489SHong Zhang     buf_si[0]   = nrows;
86342c57489SHong Zhang     buf_si_i[0] = 0;
86442c57489SHong Zhang     nrows       = 0;
865de0260b3SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
866ded4f561SHong Zhang       nzi = coi[i+1] - coi[i];
867ded4f561SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
868de0260b3SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
86942c57489SHong Zhang       nrows++;
87042c57489SHong Zhang     }
871529bc5dcSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
87242c57489SHong Zhang     k++;
87342c57489SHong Zhang     buf_si += len_si[proc];
87442c57489SHong Zhang   }
875ded4f561SHong Zhang   i = merge->nrecv;
876ded4f561SHong Zhang   while (i--) {
877c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
878ded4f561SHong Zhang   }
879ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
8800c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
881a914927fSHong Zhang 
88224ecddacSHong Zhang   ierr = PetscFree2(len_si,sstatus);CHKERRQ(ierr);
88342c57489SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
884ded4f561SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
88542c57489SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
88642c57489SHong Zhang 
8876b911d16SHong Zhang   /* (5) compute the local portion of C (mpi mat) */
8886b911d16SHong Zhang   /*----------------------------------------------*/
889ded4f561SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
890ded4f561SHong Zhang 
89124ecddacSHong Zhang   /* allocate pti array and free space for accumulating nonzero column info */
892854ce69bSBarry Smith   ierr   = PetscMalloc1(pn+1,&pti);CHKERRQ(ierr);
89324ecddacSHong Zhang   pti[0] = 0;
89442c57489SHong Zhang 
895d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
896f91af8c7SBarry Smith   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(pi_loc[pm],api[am]));
89722d28d08SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
89842c57489SHong Zhang   current_space = free_space;
89942c57489SHong Zhang 
900dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
90142c57489SHong Zhang   for (k=0; k<merge->nrecv; k++) {
90242c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
90342c57489SHong Zhang     nrows       = *buf_ri_k[k];
90442c57489SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
90542c57489SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
90642c57489SHong Zhang   }
90742c57489SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
908936ad1eaSHong Zhang 
90942c57489SHong Zhang   for (i=0; i<pn; i++) {
910ded4f561SHong Zhang     /* add pdt[i,:]*AP into lnk */
911ded4f561SHong Zhang     pnz = pdti[i+1] - pdti[i];
91277584fe6SHong Zhang     ptJ = pdtj + pdti[i];
91377584fe6SHong Zhang     for (j=0; j<pnz; j++) {
91477584fe6SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
915ded4f561SHong Zhang       apnz = api[row+1] - api[row];
916ded4f561SHong Zhang       Jptr = apj + api[row];
917ded4f561SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
918a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
919ded4f561SHong Zhang     }
920a914927fSHong Zhang 
92142c57489SHong Zhang     /* add received col data into lnk */
92242c57489SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
92342c57489SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
924ded4f561SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
925ded4f561SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
926a914927fSHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
92742c57489SHong Zhang         nextrow[k]++; nextci[k]++;
92842c57489SHong Zhang       }
92942c57489SHong Zhang     }
930a914927fSHong Zhang     nnz = lnk[0];
93142c57489SHong Zhang 
93242c57489SHong Zhang     /* if free space is not available, make more free space */
933ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
934f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
935d16ca5f0SHong Zhang       nspacedouble++;
93642c57489SHong Zhang     }
93742c57489SHong Zhang     /* copy data into free space, then initialize lnk */
938a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
939ded4f561SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
9402205254eSKarl Rupp 
941ded4f561SHong Zhang     current_space->array           += nnz;
942ded4f561SHong Zhang     current_space->local_used      += nnz;
943ded4f561SHong Zhang     current_space->local_remaining -= nnz;
9442205254eSKarl Rupp 
94524ecddacSHong Zhang     pti[i+1] = pti[i] + nnz;
94642c57489SHong Zhang   }
947ded4f561SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
9480572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
94942c57489SHong Zhang 
950854ce69bSBarry Smith   ierr      = PetscMalloc1(pti[pn]+1,&ptj);CHKERRQ(ierr);
95124ecddacSHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,ptj);CHKERRQ(ierr);
95224ecddacSHong Zhang   afill_tmp = (PetscReal)pti[pn]/(pi_loc[pm] + api[am]+1);
953d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
95442c57489SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
95542c57489SHong Zhang 
9566b911d16SHong Zhang   /* (6) create symbolic parallel matrix Cmpi */
9576b911d16SHong Zhang   /*------------------------------------------*/
95877584fe6SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
95977584fe6SHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
96033d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
96177584fe6SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
96277584fe6SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
96342c57489SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
964a2f3521dSMark F. Adams 
965b091e509SHong Zhang   merge->bi        = pti;      /* Cseq->i */
966b091e509SHong Zhang   merge->bj        = ptj;      /* Cseq->j */
967b091e509SHong Zhang   merge->coi       = coi;      /* Co->i   */
968b091e509SHong Zhang   merge->coj       = coj;      /* Co->j   */
96942c57489SHong Zhang   merge->buf_ri    = buf_ri;
97042c57489SHong Zhang   merge->buf_rj    = buf_rj;
971de0260b3SHong Zhang   merge->owners_co = owners_co;
97277584fe6SHong Zhang   merge->destroy   = Cmpi->ops->destroy;
97342c57489SHong Zhang 
97477584fe6SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
97577584fe6SHong Zhang   c           = (Mat_MPIAIJ*)Cmpi->data;
976f8487c73SHong Zhang   c->ptap     = ptap;
97777584fe6SHong Zhang   ptap->api   = api;
97877584fe6SHong Zhang   ptap->apj   = apj;
9798403a639SHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
9808403a639SHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
9818403a639SHong Zhang 
9828403a639SHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
9838403a639SHong Zhang   Cmpi->assembled        = PETSC_FALSE;
9848403a639SHong Zhang   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
9858403a639SHong Zhang   Cmpi->ops->duplicate   = MatDuplicate_MPIAIJ_MatPtAP;
9868403a639SHong Zhang   Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap;
98777584fe6SHong Zhang   *C                     = Cmpi;
988414894bdSHong Zhang 
989414894bdSHong Zhang   /* flag 'scalable' determines which implementations to be used:
990414894bdSHong Zhang        0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa;
991414894bdSHong Zhang        1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */
992414894bdSHong Zhang   /* set default scalable */
993da0a95b2SSatish Balay   ptap->scalable = PETSC_FALSE; /* PETSC_TRUE; */
9942205254eSKarl Rupp 
995c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)Cmpi)->options,((PetscObject)Cmpi)->prefix,"-matptap_scalable",&ptap->scalable,NULL);CHKERRQ(ierr);
996414894bdSHong Zhang   if (!ptap->scalable) {  /* Do dense axpy */
9971795a4d1SJed Brown     ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
998414894bdSHong Zhang   } else {
9991795a4d1SJed Brown     ierr = PetscCalloc1(ap_rmax+1,&ptap->apa);CHKERRQ(ierr);
1000414894bdSHong Zhang   }
1001414894bdSHong Zhang 
1002f2b054eeSHong Zhang #if defined(PETSC_USE_INFO)
100324ecddacSHong Zhang   if (pti[pn] != 0) {
100457622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
100557622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1006f2b054eeSHong Zhang   } else {
100777584fe6SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1008f2b054eeSHong Zhang   }
1009f2b054eeSHong Zhang #endif
101042c57489SHong Zhang   PetscFunctionReturn(0);
101142c57489SHong Zhang }
101242c57489SHong Zhang 
101342c57489SHong Zhang #undef __FUNCT__
10148403a639SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap"
10158403a639SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap(Mat A,Mat P,Mat C)
101642c57489SHong Zhang {
101742c57489SHong Zhang   PetscErrorCode      ierr;
1018f8487c73SHong Zhang   Mat_MPIAIJ          *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
101942c57489SHong Zhang   Mat_SeqAIJ          *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
1020de0260b3SHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
102142c57489SHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth;
1022f8487c73SHong Zhang   Mat_PtAPMPI         *ptap;
10239f4155fbSHong Zhang   Mat_Merge_SeqsToMPI *merge;
10241d633602SHong Zhang   PetscInt            *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp;
102582412ba7SHong Zhang   PetscInt            *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj;
102682412ba7SHong Zhang   PetscInt            i,j,k,anz,pnz,apnz,nextap,row,*cj;
1027e900a757SHong Zhang   MatScalar           *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp;
1028d0f46423SBarry Smith   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1029ce94432eSBarry Smith   MPI_Comm            comm;
103042c57489SHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
103182412ba7SHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
103242c57489SHong Zhang   PetscInt            **buf_ri,**buf_rj;
103350f5bf86SHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
103442c57489SHong Zhang   MPI_Request         *s_waits,*r_waits;
103542c57489SHong Zhang   MPI_Status          *status;
103682412ba7SHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
103782412ba7SHong Zhang   PetscInt            *api,*apj,*coi,*coj;
1038d0f46423SBarry Smith   PetscInt            *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend;
10396ce94e8fSHong Zhang   PetscBool           scalable;
104038571addSHong Zhang #if defined(PTAP_PROFILE)
1041e497d3c8SHong Zhang   PetscLogDouble t0,t1,t2,eP,t3,t4,et2_AP=0.0,ePtAP=0.0,t2_0,t2_1,t2_2;
104238571addSHong Zhang #endif
104342c57489SHong Zhang 
104442c57489SHong Zhang   PetscFunctionBegin;
1045ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
104638571addSHong Zhang #if defined(PTAP_PROFILE)
10478563dfccSBarry Smith   ierr = PetscTime(&t0);CHKERRQ(ierr);
104838571addSHong Zhang #endif
104942c57489SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
105042c57489SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
105142c57489SHong Zhang 
1052f8487c73SHong Zhang   ptap = c->ptap;
1053ce94432eSBarry 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");
1054f8487c73SHong Zhang   merge    = ptap->merge;
1055414894bdSHong Zhang   apa      = ptap->apa;
10566ce94e8fSHong Zhang   scalable = ptap->scalable;
10576c6e00e0SHong Zhang 
1058a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
10596b911d16SHong Zhang   /*-----------------------------------------------------*/
1060f8487c73SHong Zhang   if (ptap->reuse == MAT_INITIAL_MATRIX) {
10619f4155fbSHong Zhang     /* P_oth and P_loc are obtained in MatPtASymbolic(), skip calling MatGetBrowsOfAoCols() and MatMPIAIJGetLocalMat() */
1062f8487c73SHong Zhang     ptap->reuse = MAT_REUSE_MATRIX;
10636c6e00e0SHong Zhang   } else { /* update numerical values of P_oth and P_loc */
1064b7f45c76SHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
1065a1a4e84aSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
10666c6e00e0SHong Zhang   }
106738571addSHong Zhang #if defined(PTAP_PROFILE)
10688563dfccSBarry Smith   ierr = PetscTime(&t1);CHKERRQ(ierr);
106905d62848SHong Zhang   eP = t1-t0;
107038571addSHong Zhang #endif
107105d62848SHong Zhang   /*
107205d62848SHong Zhang   printf("[%d] Ad: %d, %d; Ao: %d, %d; P_loc: %d, %d; P_oth %d, %d;\n",rank,
107305d62848SHong Zhang          a->A->rmap->N,a->A->cmap->N,a->B->rmap->N,a->B->cmap->N,
107405d62848SHong Zhang          ptap->P_loc->rmap->N,ptap->P_loc->cmap->N,
107505d62848SHong Zhang          ptap->P_oth->rmap->N,ptap->P_oth->cmap->N);
107605d62848SHong Zhang    */
1077f8487c73SHong Zhang 
1078f9473cf7SHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
1079f9473cf7SHong Zhang   /*--------------------------------------------------------------*/
108042c57489SHong Zhang   /* get data from symbolic products */
1081a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
1082a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
108389ae1891SBarry Smith   pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a;
108442c57489SHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
108542c57489SHong Zhang 
1086de0260b3SHong Zhang   coi  = merge->coi; coj = merge->coj;
10871795a4d1SJed Brown   ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
1088de0260b3SHong Zhang 
1089de0260b3SHong Zhang   bi     = merge->bi; bj = merge->bj;
10907a2fc3feSBarry Smith   owners = merge->rowmap->range;
10911795a4d1SJed Brown   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);  /* ba: Cseq->a */
109242c57489SHong Zhang 
1093a1a4e84aSHong Zhang   api = ptap->api; apj = ptap->apj;
1094d9824c15SHong Zhang 
10959516a85cSHong Zhang   if (!scalable) { /* Do dense axpy on apa (length of pN, stores A[i,:]*P) - nonscalable, but faster (could take 1/3 scalable time) */
1096b38d1a2dSJed Brown     ierr = PetscInfo(C,"Using non-scalable dense axpy\n");CHKERRQ(ierr);
10978cb82516SHong Zhang #if defined(PTAP_PROFILE)
109805d62848SHong Zhang     ierr = PetscTime(&t1);CHKERRQ(ierr);
10998cb82516SHong Zhang #endif
1100a9d06231SHong Zhang     for (i=0; i<am; i++) {
1101b091e509SHong Zhang #if defined(PTAP_PROFILE)
11028563dfccSBarry Smith       ierr = PetscTime(&t2_0);CHKERRQ(ierr);
1103b091e509SHong Zhang #endif
1104a9d06231SHong Zhang       /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */
1105a9d06231SHong Zhang       /*------------------------------------------------------------*/
1106a9d06231SHong Zhang       apJ = apj + api[i];
1107a9d06231SHong Zhang 
1108a9d06231SHong Zhang       /* diagonal portion of A */
1109a9d06231SHong Zhang       anz = adi[i+1] - adi[i];
1110a9d06231SHong Zhang       adj = ad->j + adi[i];
1111a9d06231SHong Zhang       ada = ad->a + adi[i];
1112a9d06231SHong Zhang       for (j=0; j<anz; j++) {
1113a9d06231SHong Zhang         row = adj[j];
1114a9d06231SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
1115a9d06231SHong Zhang         pj  = pj_loc + pi_loc[row];
1116a9d06231SHong Zhang         pa  = pa_loc + pi_loc[row];
1117a9d06231SHong Zhang 
1118a9d06231SHong Zhang         /* perform dense axpy */
1119e900a757SHong Zhang         valtmp = ada[j];
1120a9d06231SHong Zhang         for (k=0; k<pnz; k++) {
1121e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
1122a9d06231SHong Zhang         }
1123a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
1124a9d06231SHong Zhang       }
1125a9d06231SHong Zhang 
1126a9d06231SHong Zhang       /* off-diagonal portion of A */
1127a9d06231SHong Zhang       anz = aoi[i+1] - aoi[i];
1128a9d06231SHong Zhang       aoj = ao->j + aoi[i];
1129a9d06231SHong Zhang       aoa = ao->a + aoi[i];
1130a9d06231SHong Zhang       for (j=0; j<anz; j++) {
1131a9d06231SHong Zhang         row = aoj[j];
1132a9d06231SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
1133a9d06231SHong Zhang         pj  = pj_oth + pi_oth[row];
1134a9d06231SHong Zhang         pa  = pa_oth + pi_oth[row];
1135a9d06231SHong Zhang 
1136a9d06231SHong Zhang         /* perform dense axpy */
1137e900a757SHong Zhang         valtmp = aoa[j];
1138a9d06231SHong Zhang         for (k=0; k<pnz; k++) {
1139e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
1140a9d06231SHong Zhang         }
1141a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
1142a9d06231SHong Zhang       }
1143b091e509SHong Zhang #if defined(PTAP_PROFILE)
11448563dfccSBarry Smith       ierr    = PetscTime(&t2_1);CHKERRQ(ierr);
1145b091e509SHong Zhang       et2_AP += t2_1 - t2_0;
1146b091e509SHong Zhang #endif
1147a9d06231SHong Zhang 
1148a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
1149a9d06231SHong Zhang       /*--------------------------------------------------------------*/
1150a9d06231SHong Zhang       apnz = api[i+1] - api[i];
1151a9d06231SHong Zhang       /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */
1152a9d06231SHong Zhang       pnz = po->i[i+1] - po->i[i];
1153e900a757SHong Zhang       poJ = po->j + po->i[i];
1154a9d06231SHong Zhang       pA  = po->a + po->i[i];
1155a9d06231SHong Zhang       for (j=0; j<pnz; j++) {
1156e900a757SHong Zhang         row = poJ[j];
1157e900a757SHong Zhang         cnz = coi[row+1] - coi[row];
1158e900a757SHong Zhang         cj  = coj + coi[row];
1159e900a757SHong Zhang         ca  = coa + coi[row];
1160a9d06231SHong Zhang         /* perform dense axpy */
1161e900a757SHong Zhang         valtmp = pA[j];
1162a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
1163e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
1164a9d06231SHong Zhang         }
1165a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1166a9d06231SHong Zhang       }
1167a9d06231SHong Zhang       /* put the value into Cd (diagonal part) */
1168a9d06231SHong Zhang       pnz = pd->i[i+1] - pd->i[i];
1169e900a757SHong Zhang       pdJ = pd->j + pd->i[i];
1170a9d06231SHong Zhang       pA  = pd->a + pd->i[i];
1171a9d06231SHong Zhang       for (j=0; j<pnz; j++) {
1172e900a757SHong Zhang         row = pdJ[j];
1173e900a757SHong Zhang         cnz = bi[row+1] - bi[row];
1174e900a757SHong Zhang         cj  = bj + bi[row];
1175e900a757SHong Zhang         ca  = ba + bi[row];
1176a9d06231SHong Zhang         /* perform dense axpy */
1177e900a757SHong Zhang         valtmp = pA[j];
1178a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
1179e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
1180a9d06231SHong Zhang         }
1181a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1182a9d06231SHong Zhang       }
11838403a639SHong Zhang 
1184a9d06231SHong Zhang       /* zero the current row of A*P */
1185a9d06231SHong Zhang       for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0;
1186b091e509SHong Zhang #if defined(PTAP_PROFILE)
11878563dfccSBarry Smith       ierr      = PetscTime(&t2_2);CHKERRQ(ierr);
118805d62848SHong Zhang       ePtAP += t2_2 - t2_1;
1189b091e509SHong Zhang #endif
1190a9d06231SHong Zhang     }
119138571addSHong Zhang   } else { /* Do sparse axpy on apa (length of ap_rmax, stores A[i,:]*P) - scalable, but slower */
1192b38d1a2dSJed Brown     ierr = PetscInfo(C,"Using scalable sparse axpy\n");CHKERRQ(ierr);
119338571addSHong Zhang     /*-----------------------------------------------------------------------------------------*/
1194a9d06231SHong Zhang     pA=pa_loc;
119542c57489SHong Zhang     for (i=0; i<am; i++) {
1196b091e509SHong Zhang #if defined(PTAP_PROFILE)
11978563dfccSBarry Smith       ierr = PetscTime(&t2_0);CHKERRQ(ierr);
1198b091e509SHong Zhang #endif
1199f4a743e1SHong Zhang       /* form i-th sparse row of A*P */
120082412ba7SHong Zhang       apnz = api[i+1] - api[i];
120182412ba7SHong Zhang       apJ  = apj + api[i];
120242c57489SHong Zhang       /* diagonal portion of A */
120382412ba7SHong Zhang       anz = adi[i+1] - adi[i];
12041d633602SHong Zhang       adj = ad->j + adi[i];
12051d633602SHong Zhang       ada = ad->a + adi[i];
120682412ba7SHong Zhang       for (j=0; j<anz; j++) {
12071d633602SHong Zhang         row    = adj[j];
120882412ba7SHong Zhang         pnz    = pi_loc[row+1] - pi_loc[row];
120982412ba7SHong Zhang         pj     = pj_loc + pi_loc[row];
121082412ba7SHong Zhang         pa     = pa_loc + pi_loc[row];
1211e900a757SHong Zhang         valtmp = ada[j];
1212f4a743e1SHong Zhang         nextp  = 0;
121382412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
121482412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
1215e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
121642c57489SHong Zhang           }
121742c57489SHong Zhang         }
1218dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
121942c57489SHong Zhang       }
122042c57489SHong Zhang       /* off-diagonal portion of A */
122182412ba7SHong Zhang       anz = aoi[i+1] - aoi[i];
12221d633602SHong Zhang       aoj = ao->j + aoi[i];
12231d633602SHong Zhang       aoa = ao->a + aoi[i];
122482412ba7SHong Zhang       for (j=0; j<anz; j++) {
12251d633602SHong Zhang         row    = aoj[j];
122682412ba7SHong Zhang         pnz    = pi_oth[row+1] - pi_oth[row];
122782412ba7SHong Zhang         pj     = pj_oth + pi_oth[row];
122882412ba7SHong Zhang         pa     = pa_oth + pi_oth[row];
1229e900a757SHong Zhang         valtmp = aoa[j];
1230f4a743e1SHong Zhang         nextp  = 0;
123182412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
123282412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
1233e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
123442c57489SHong Zhang           }
123542c57489SHong Zhang         }
1236dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
123742c57489SHong Zhang       }
1238b091e509SHong Zhang #if defined(PTAP_PROFILE)
12398563dfccSBarry Smith       ierr    = PetscTime(&t2_1);CHKERRQ(ierr);
1240b091e509SHong Zhang       et2_AP += t2_1 - t2_0;
1241b091e509SHong Zhang #endif
124242c57489SHong Zhang 
1243a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
1244a9d06231SHong Zhang       /*--------------------------------------------------------------*/
124582412ba7SHong Zhang       pnz = pi_loc[i+1] - pi_loc[i];
1246f9473cf7SHong Zhang       pJ  = pj_loc + pi_loc[i];
124782412ba7SHong Zhang       for (j=0; j<pnz; j++) {
124842c57489SHong Zhang         nextap = 0;
1249f9473cf7SHong Zhang         row    = pJ[j]; /* global index */
125082412ba7SHong Zhang         if (row < pcstart || row >=pcend) { /* put the value into Co */
1251e900a757SHong Zhang           row = *poJ;
1252e900a757SHong Zhang           cj  = coj + coi[row];
1253e900a757SHong Zhang           ca  = coa + coi[row]; poJ++;
125482412ba7SHong Zhang         } else {                            /* put the value into Cd */
1255e900a757SHong Zhang           row = *pdJ;
1256e900a757SHong Zhang           cj  = bj + bi[row];
1257e900a757SHong Zhang           ca  = ba + bi[row]; pdJ++;
125842c57489SHong Zhang         }
1259e900a757SHong Zhang         valtmp = pA[j];
126082412ba7SHong Zhang         for (k=0; nextap<apnz; k++) {
1261e900a757SHong Zhang           if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++];
126242c57489SHong Zhang         }
1263dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
1264de0260b3SHong Zhang       }
12650496f32cSHong Zhang       pA += pnz;
126642c57489SHong Zhang       /* zero the current row info for A*P */
126782412ba7SHong Zhang       ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr);
1268b091e509SHong Zhang #if defined(PTAP_PROFILE)
12698563dfccSBarry Smith       ierr      = PetscTime(&t2_2);CHKERRQ(ierr);
127005d62848SHong Zhang       ePtAP += t2_2 - t2_1;
1271b091e509SHong Zhang #endif
127242c57489SHong Zhang     }
1273d9824c15SHong Zhang   }
127438571addSHong Zhang #if defined(PTAP_PROFILE)
12758563dfccSBarry Smith   ierr = PetscTime(&t2);CHKERRQ(ierr);
127638571addSHong Zhang #endif
127742c57489SHong Zhang 
1278a9d06231SHong Zhang   /* 3) send and recv matrix values coa */
1279a9d06231SHong Zhang   /*------------------------------------*/
128042c57489SHong Zhang   buf_ri = merge->buf_ri;
128142c57489SHong Zhang   buf_rj = merge->buf_rj;
128242c57489SHong Zhang   len_s  = merge->len_s;
1283fc42d0c8SSatish Balay   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
128442c57489SHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
128542c57489SHong Zhang 
1286dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
128742c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
128842c57489SHong Zhang     if (!len_s[proc]) continue;
1289de0260b3SHong Zhang     i    = merge->owners_co[proc];
1290de0260b3SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
129142c57489SHong Zhang     k++;
129242c57489SHong Zhang   }
12930c468ba9SBarry Smith   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
12940c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
129542c57489SHong Zhang 
1296a9d06231SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
129742c57489SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
129882412ba7SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
129938571addSHong Zhang #if defined(PTAP_PROFILE)
13008563dfccSBarry Smith   ierr = PetscTime(&t3);CHKERRQ(ierr);
130138571addSHong Zhang #endif
130242c57489SHong Zhang 
1303a9d06231SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1304a9d06231SHong Zhang   /*------------------------------------------------------*/
1305dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
130642c57489SHong Zhang   for (k=0; k<merge->nrecv; k++) {
130742c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
130842c57489SHong Zhang     nrows       = *(buf_ri_k[k]);
130942c57489SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
131082412ba7SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
131142c57489SHong Zhang   }
131242c57489SHong Zhang 
1313de0260b3SHong Zhang   for (i=0; i<cm; i++) {
131482412ba7SHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
131542c57489SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1316de0260b3SHong Zhang     ba_i = ba + bi[i];
131782412ba7SHong Zhang     bnz  = bi[i+1] - bi[i];
131842c57489SHong Zhang     /* add received vals into ba */
131942c57489SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
132042c57489SHong Zhang       /* i-th row */
132142c57489SHong Zhang       if (i == *nextrow[k]) {
132282412ba7SHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
132382412ba7SHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
132482412ba7SHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
132582412ba7SHong Zhang         nextcj = 0;
132682412ba7SHong Zhang         for (j=0; nextcj<cnz; j++) {
132782412ba7SHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
132882412ba7SHong Zhang             ba_i[j] += ca[nextcj++];
132942c57489SHong Zhang           }
133042c57489SHong Zhang         }
133182412ba7SHong Zhang         nextrow[k]++; nextci[k]++;
1332c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
133342c57489SHong Zhang       }
133442c57489SHong Zhang     }
133582412ba7SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
133642c57489SHong Zhang   }
133742c57489SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
133842c57489SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
133942c57489SHong Zhang 
1340de0260b3SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1341c05d87d6SBarry Smith   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
134242c57489SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
13430572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
134438571addSHong Zhang #if defined(PTAP_PROFILE)
13458563dfccSBarry Smith   ierr = PetscTime(&t4);CHKERRQ(ierr);
13469516a85cSHong Zhang   if (rank==1) {
134705d62848SHong 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);
13489516a85cSHong Zhang   }
134938571addSHong Zhang #endif
135042c57489SHong Zhang   PetscFunctionReturn(0);
135142c57489SHong Zhang }
1352