xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision 9f91fb7a2a458dd4b36f4ddc8f1de7c61b64afc6)
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);
16f8487c73SHong Zhang PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A)
17cc6cb787SHong Zhang {
18cc6cb787SHong Zhang   PetscErrorCode ierr;
19f8487c73SHong Zhang   Mat_MPIAIJ     *a=(Mat_MPIAIJ*)A->data;
20f8487c73SHong Zhang   Mat_PtAPMPI    *ptap=a->ptap;
21cc6cb787SHong Zhang 
22cc6cb787SHong Zhang   PetscFunctionBegin;
23f8487c73SHong Zhang   if (ptap) {
24c8b0795cSMark F. Adams     Mat_Merge_SeqsToMPI *merge=ptap->merge;
25b7f45c76SHong Zhang     ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
26f8487c73SHong Zhang     ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
27a1a4e84aSHong Zhang     ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
28a1a4e84aSHong Zhang     ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
29c5af039cSHong Zhang     ierr = MatDestroy(&ptap->A_loc);CHKERRQ(ierr); /* used by MatTransposeMatMult() */
3005d62848SHong Zhang     ierr = MatDestroy(&ptap->Rd);CHKERRQ(ierr);
3105d62848SHong Zhang     ierr = MatDestroy(&ptap->Ro);CHKERRQ(ierr);
328403a639SHong Zhang     if (ptap->AP_loc) { /* used by alg_rap */
33681d504bSHong Zhang       Mat_SeqAIJ *ap = (Mat_SeqAIJ*)(ptap->AP_loc)->data;
34681d504bSHong Zhang       ierr = PetscFree(ap->i);CHKERRQ(ierr);
35681d504bSHong Zhang       ierr = PetscFree2(ap->j,ap->a);CHKERRQ(ierr);
360d3441aeSHong Zhang       ierr = MatDestroy(&ptap->AP_loc);CHKERRQ(ierr);
378403a639SHong Zhang     } else { /* used by alg_ptap */
388403a639SHong Zhang       ierr = PetscFree(ptap->api);CHKERRQ(ierr);
398403a639SHong Zhang       ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
40681d504bSHong Zhang     }
412259aa2eSHong Zhang     ierr = MatDestroy(&ptap->C_loc);CHKERRQ(ierr);
422259aa2eSHong Zhang     ierr = MatDestroy(&ptap->C_oth);CHKERRQ(ierr);
43414894bdSHong Zhang     if (ptap->apa) {ierr = PetscFree(ptap->apa);CHKERRQ(ierr);}
44a560ef98SHong Zhang 
458403a639SHong Zhang     if (merge) { /* used by alg_ptap */
46cc6cb787SHong Zhang       ierr = PetscFree(merge->id_r);CHKERRQ(ierr);
47cc6cb787SHong Zhang       ierr = PetscFree(merge->len_s);CHKERRQ(ierr);
48cc6cb787SHong Zhang       ierr = PetscFree(merge->len_r);CHKERRQ(ierr);
49cc6cb787SHong Zhang       ierr = PetscFree(merge->bi);CHKERRQ(ierr);
50cc6cb787SHong Zhang       ierr = PetscFree(merge->bj);CHKERRQ(ierr);
51c05d87d6SBarry Smith       ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr);
52cc6cb787SHong Zhang       ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr);
53c05d87d6SBarry Smith       ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr);
54cc6cb787SHong Zhang       ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr);
55445158ffSHong Zhang       ierr = PetscFree(merge->coi);CHKERRQ(ierr);
56445158ffSHong Zhang       ierr = PetscFree(merge->coj);CHKERRQ(ierr);
5705b42c5fSBarry Smith       ierr = PetscFree(merge->owners_co);CHKERRQ(ierr);
586bf464f9SBarry Smith       ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr);
59f8487c73SHong Zhang       ierr = PetscFree(ptap->merge);CHKERRQ(ierr);
60bf0cc555SLisandro Dalcin     }
61a560ef98SHong Zhang 
62a560ef98SHong Zhang     ierr = ptap->destroy(A);CHKERRQ(ierr);
63f8487c73SHong Zhang     ierr = PetscFree(ptap);CHKERRQ(ierr);
64c8b0795cSMark F. Adams   }
65cc6cb787SHong Zhang   PetscFunctionReturn(0);
66cc6cb787SHong Zhang }
67cc6cb787SHong Zhang 
68aa690a28SHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M)
691a47ec54SHong Zhang {
701a47ec54SHong Zhang   PetscErrorCode ierr;
711a47ec54SHong Zhang   Mat_MPIAIJ     *a     = (Mat_MPIAIJ*)A->data;
721a47ec54SHong Zhang   Mat_PtAPMPI    *ptap  = a->ptap;
731a47ec54SHong Zhang 
741a47ec54SHong Zhang   PetscFunctionBegin;
751a47ec54SHong Zhang   ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr);
761a47ec54SHong Zhang   (*M)->ops->destroy   = ptap->destroy;
771a47ec54SHong Zhang   (*M)->ops->duplicate = ptap->duplicate;
781a47ec54SHong Zhang   PetscFunctionReturn(0);
791a47ec54SHong Zhang }
801a47ec54SHong Zhang 
81150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
8242c57489SHong Zhang {
8342c57489SHong Zhang   PetscErrorCode ierr;
848403a639SHong Zhang   PetscBool      rap=PETSC_TRUE; /* do R=P^T locally, then C=R*A*P */
8567a12041SHong Zhang   MPI_Comm       comm;
8642c57489SHong Zhang 
8742c57489SHong Zhang   PetscFunctionBegin;
8867a12041SHong Zhang   /* check if matrix local sizes are compatible */
8967a12041SHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
906c4ed002SBarry 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);
916c4ed002SBarry 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);
9267a12041SHong Zhang 
93c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-matrap",&rap,NULL);CHKERRQ(ierr);
94cf3ca8ceSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
953ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
968403a639SHong Zhang     if (rap) { /* do R=P^T locally, then C=R*A*P */
97cf3ca8ceSHong Zhang       ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);
988403a639SHong Zhang     } else {       /* do P^T*A*P */
998403a639SHong Zhang       ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(A,P,fill,C);CHKERRQ(ierr);
1000d3441aeSHong Zhang     }
1013ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
1027a7894deSKris Buschelman   }
1033ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
1048403a639SHong Zhang   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
1053ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
10642c57489SHong Zhang   PetscFunctionReturn(0);
10742c57489SHong Zhang }
10842c57489SHong Zhang 
109ae5f0867Sstefano_zampini #if defined(PETSC_HAVE_HYPRE)
110ae5f0867Sstefano_zampini PETSC_INTERN PetscErrorCode MatPtAPSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat*);
111ae5f0867Sstefano_zampini #endif
112ae5f0867Sstefano_zampini 
113e953e47cSHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,PetscReal fill,Mat *C)
1140d3441aeSHong Zhang {
1150d3441aeSHong Zhang   PetscErrorCode      ierr;
1160d3441aeSHong Zhang   Mat_PtAPMPI         *ptap;
117681d504bSHong Zhang   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
1182259aa2eSHong Zhang   MPI_Comm            comm;
1192259aa2eSHong Zhang   PetscMPIInt         size,rank;
12015a3b8e2SHong Zhang   Mat                 Cmpi;
12115a3b8e2SHong Zhang   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
122d0e9a2f7SHong Zhang   PetscInt            am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
123d0e9a2f7SHong Zhang   PetscInt            *lnk,i,k,pnz,row,nsend;
12415a3b8e2SHong Zhang   PetscBT             lnkbt;
125f671be37SHong Zhang   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
12615a3b8e2SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
127d0e9a2f7SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners,nzi,nspacedouble;
12815a3b8e2SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
12915a3b8e2SHong Zhang   MPI_Request         *swaits,*rwaits;
13015a3b8e2SHong Zhang   MPI_Status          *sstatus,rstatus;
131445158ffSHong Zhang   PetscLayout         rowmap;
132445158ffSHong Zhang   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
133445158ffSHong Zhang   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
134936ad1eaSHong Zhang   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0,Crmax,*aj,*ai,*pi;
13567a12041SHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*c_loc,*c_oth;
13667a12041SHong Zhang   PetscScalar         *apv;
13778d30b94SHong Zhang   PetscTable          ta;
138aa690a28SHong Zhang #if defined(PETSC_USE_INFO)
1398cb82516SHong Zhang   PetscReal           apfill;
140aa690a28SHong Zhang #endif
141681d504bSHong Zhang #if defined(PTAP_PROFILE)
142681d504bSHong Zhang   PetscLogDouble      t0,t1,t11,t12,t2,t3,t4;
143681d504bSHong Zhang #endif
14467a12041SHong Zhang 
14567a12041SHong Zhang   PetscFunctionBegin;
14667a12041SHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
14767a12041SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
14867a12041SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
149*9f91fb7aSHong Zhang   if (!rank) printf(" MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable...\n");
150ae5f0867Sstefano_zampini 
151ae5f0867Sstefano_zampini   /* create symbolic parallel matrix Cmpi */
152ae5f0867Sstefano_zampini   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
153ae5f0867Sstefano_zampini   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
154ae5f0867Sstefano_zampini 
155e953e47cSHong Zhang #if defined(PTAP_PROFILE)
156e953e47cSHong Zhang   ierr = PetscTime(&t0);CHKERRQ(ierr);
157e953e47cSHong Zhang #endif
158e953e47cSHong Zhang 
159e953e47cSHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
160e953e47cSHong Zhang   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
161e953e47cSHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
162e953e47cSHong Zhang 
163e953e47cSHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
164e953e47cSHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
165e953e47cSHong Zhang   /* get P_loc by taking all local rows of P */
166e953e47cSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
167e953e47cSHong Zhang 
168e953e47cSHong Zhang   /* (0) compute Rd = Pd^T, Ro = Po^T  */
169e953e47cSHong Zhang   /* --------------------------------- */
170e953e47cSHong Zhang   ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
171e953e47cSHong Zhang   ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
172e953e47cSHong Zhang #if defined(PTAP_PROFILE)
173e953e47cSHong Zhang   ierr = PetscTime(&t11);CHKERRQ(ierr);
174e953e47cSHong Zhang #endif
175e953e47cSHong Zhang 
176e953e47cSHong Zhang   /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */
177e953e47cSHong Zhang   /* ----------------------------------------------------------------- */
178e953e47cSHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
179e953e47cSHong Zhang   p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
180e953e47cSHong Zhang 
181e953e47cSHong Zhang   /* create and initialize a linked list */
182e953e47cSHong Zhang   ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */
183e953e47cSHong Zhang   MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta);
184e953e47cSHong Zhang   MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta);
185e953e47cSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */
186e953e47cSHong 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); */
187e953e47cSHong Zhang 
188e953e47cSHong Zhang   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
189e953e47cSHong Zhang 
190e953e47cSHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
191e953e47cSHong Zhang   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr);
192e953e47cSHong Zhang   current_space = free_space;
193e953e47cSHong Zhang   nspacedouble  = 0;
194e953e47cSHong Zhang 
195e953e47cSHong Zhang   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
196e953e47cSHong Zhang   api[0] = 0;
197e953e47cSHong Zhang   for (i=0; i<am; i++) {
198e953e47cSHong Zhang     /* diagonal portion: Ad[i,:]*P */
199e953e47cSHong Zhang     ai = ad->i; pi = p_loc->i;
200e953e47cSHong Zhang     nzi = ai[i+1] - ai[i];
201e953e47cSHong Zhang     aj  = ad->j + ai[i];
202e953e47cSHong Zhang     for (j=0; j<nzi; j++) {
203e953e47cSHong Zhang       row  = aj[j];
204e953e47cSHong Zhang       pnz  = pi[row+1] - pi[row];
205e953e47cSHong Zhang       Jptr = p_loc->j + pi[row];
206e953e47cSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
207e953e47cSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
208e953e47cSHong Zhang     }
209e953e47cSHong Zhang     /* off-diagonal portion: Ao[i,:]*P */
210e953e47cSHong Zhang     ai = ao->i; pi = p_oth->i;
211e953e47cSHong Zhang     nzi = ai[i+1] - ai[i];
212e953e47cSHong Zhang     aj  = ao->j + ai[i];
213e953e47cSHong Zhang     for (j=0; j<nzi; j++) {
214e953e47cSHong Zhang       row  = aj[j];
215e953e47cSHong Zhang       pnz  = pi[row+1] - pi[row];
216e953e47cSHong Zhang       Jptr = p_oth->j + pi[row];
217e953e47cSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
218e953e47cSHong Zhang     }
219e953e47cSHong Zhang     apnz     = lnk[0];
220e953e47cSHong Zhang     api[i+1] = api[i] + apnz;
221e953e47cSHong Zhang     if (ap_rmax < apnz) ap_rmax = apnz;
222e953e47cSHong Zhang 
223e953e47cSHong Zhang     /* if free space is not available, double the total space in the list */
224e953e47cSHong Zhang     if (current_space->local_remaining<apnz) {
225e953e47cSHong Zhang       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
226e953e47cSHong Zhang       nspacedouble++;
227e953e47cSHong Zhang     }
228e953e47cSHong Zhang 
229e953e47cSHong Zhang     /* Copy data into free space, then initialize lnk */
230e953e47cSHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
231e953e47cSHong Zhang 
232e953e47cSHong Zhang     current_space->array           += apnz;
233e953e47cSHong Zhang     current_space->local_used      += apnz;
234e953e47cSHong Zhang     current_space->local_remaining -= apnz;
235e953e47cSHong Zhang   }
236e953e47cSHong Zhang   /* Allocate space for apj and apv, initialize apj, and */
237e953e47cSHong Zhang   /* destroy list of free space and other temporary array(s) */
238e953e47cSHong Zhang   ierr   = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr);
239e953e47cSHong Zhang   ierr   = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
240e953e47cSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
241e953e47cSHong Zhang 
242e953e47cSHong Zhang   /* Create AP_loc for reuse */
243e953e47cSHong Zhang   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr);
244e953e47cSHong Zhang 
245e953e47cSHong Zhang #if defined(PETSC_USE_INFO)
246e953e47cSHong Zhang   apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1);
247e953e47cSHong Zhang   ptap->AP_loc->info.mallocs           = nspacedouble;
248e953e47cSHong Zhang   ptap->AP_loc->info.fill_ratio_given  = fill;
249e953e47cSHong Zhang   ptap->AP_loc->info.fill_ratio_needed = apfill;
250e953e47cSHong Zhang 
251e953e47cSHong Zhang   if (api[am]) {
252e953e47cSHong Zhang     ierr = PetscInfo3(ptap->AP_loc,"AP_loc reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)apfill);CHKERRQ(ierr);
253e953e47cSHong Zhang     ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr);
254e953e47cSHong Zhang   } else {
255e953e47cSHong Zhang     ierr = PetscInfo(ptap->AP_loc,"AP_loc is empty \n");CHKERRQ(ierr);
256e953e47cSHong Zhang   }
257e953e47cSHong Zhang #endif
258e953e47cSHong Zhang 
259e953e47cSHong Zhang #if defined(PTAP_PROFILE)
260e953e47cSHong Zhang   ierr = PetscTime(&t12);CHKERRQ(ierr);
261e953e47cSHong Zhang #endif
262e953e47cSHong Zhang 
263e953e47cSHong Zhang   /* (2-1) compute symbolic Co = Ro*AP_loc  */
264e953e47cSHong Zhang   /* ------------------------------------ */
265e953e47cSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
266e953e47cSHong Zhang #if defined(PTAP_PROFILE)
267e953e47cSHong Zhang   ierr = PetscTime(&t1);CHKERRQ(ierr);
268e953e47cSHong Zhang #endif
269e953e47cSHong Zhang 
270e953e47cSHong Zhang   /* (3) send coj of C_oth to other processors  */
271e953e47cSHong Zhang   /* ------------------------------------------ */
272e953e47cSHong Zhang   /* determine row ownership */
273e953e47cSHong Zhang   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
274e953e47cSHong Zhang   rowmap->n  = pn;
275e953e47cSHong Zhang   rowmap->bs = 1;
276e953e47cSHong Zhang   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
277e953e47cSHong Zhang   owners = rowmap->range;
278e953e47cSHong Zhang 
279e953e47cSHong Zhang   /* determine the number of messages to send, their lengths */
280e953e47cSHong Zhang   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
281e953e47cSHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
282e953e47cSHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
283e953e47cSHong Zhang 
284e953e47cSHong Zhang   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
285e953e47cSHong Zhang   coi   = c_oth->i; coj = c_oth->j;
286e953e47cSHong Zhang   con   = ptap->C_oth->rmap->n;
287e953e47cSHong Zhang   proc  = 0;
288e953e47cSHong Zhang   for (i=0; i<con; i++) {
289e953e47cSHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
290e953e47cSHong Zhang     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
291e953e47cSHong Zhang     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
292e953e47cSHong Zhang   }
293e953e47cSHong Zhang 
294e953e47cSHong Zhang   len          = 0; /* max length of buf_si[], see (4) */
295e953e47cSHong Zhang   owners_co[0] = 0;
296e953e47cSHong Zhang   nsend        = 0;
297e953e47cSHong Zhang   for (proc=0; proc<size; proc++) {
298e953e47cSHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
299e953e47cSHong Zhang     if (len_s[proc]) {
300e953e47cSHong Zhang       nsend++;
301e953e47cSHong Zhang       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
302e953e47cSHong Zhang       len         += len_si[proc];
303e953e47cSHong Zhang     }
304e953e47cSHong Zhang   }
305e953e47cSHong Zhang 
306e953e47cSHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
307e953e47cSHong Zhang   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
308e953e47cSHong Zhang   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
309e953e47cSHong Zhang 
310e953e47cSHong Zhang   /* post the Irecv and Isend of coj */
311e953e47cSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
312e953e47cSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
313e953e47cSHong Zhang   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
314e953e47cSHong Zhang   for (proc=0, k=0; proc<size; proc++) {
315e953e47cSHong Zhang     if (!len_s[proc]) continue;
316e953e47cSHong Zhang     i    = owners_co[proc];
317e953e47cSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
318e953e47cSHong Zhang     k++;
319e953e47cSHong Zhang   }
320e953e47cSHong Zhang 
321e953e47cSHong Zhang   /* (2-2) compute symbolic C_loc = Rd*AP_loc */
322e953e47cSHong Zhang   /* ---------------------------------------- */
323e953e47cSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr);
324e953e47cSHong Zhang   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
325e953e47cSHong Zhang 
326e953e47cSHong Zhang   /* receives coj are complete */
327e953e47cSHong Zhang   for (i=0; i<nrecv; i++) {
328e953e47cSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
329e953e47cSHong Zhang   }
330e953e47cSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
331e953e47cSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
332e953e47cSHong Zhang 
333e953e47cSHong Zhang   /* add received column indices into ta to update Crmax */
334e953e47cSHong Zhang   for (k=0; k<nrecv; k++) {/* k-th received message */
335e953e47cSHong Zhang     Jptr = buf_rj[k];
336e953e47cSHong Zhang     for (j=0; j<len_r[k]; j++) {
337e953e47cSHong Zhang       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
338e953e47cSHong Zhang     }
339e953e47cSHong Zhang   }
340e953e47cSHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
341e953e47cSHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
342e953e47cSHong Zhang 
343e953e47cSHong Zhang   /* (4) send and recv coi */
344e953e47cSHong Zhang   /*-----------------------*/
345e953e47cSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
346e953e47cSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
347e953e47cSHong Zhang   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
348e953e47cSHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
349e953e47cSHong Zhang   for (proc=0,k=0; proc<size; proc++) {
350e953e47cSHong Zhang     if (!len_s[proc]) continue;
351e953e47cSHong Zhang     /* form outgoing message for i-structure:
352e953e47cSHong Zhang          buf_si[0]:                 nrows to be sent
353e953e47cSHong Zhang                [1:nrows]:           row index (global)
354e953e47cSHong Zhang                [nrows+1:2*nrows+1]: i-structure index
355e953e47cSHong Zhang     */
356e953e47cSHong Zhang     /*-------------------------------------------*/
357e953e47cSHong Zhang     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
358e953e47cSHong Zhang     buf_si_i    = buf_si + nrows+1;
359e953e47cSHong Zhang     buf_si[0]   = nrows;
360e953e47cSHong Zhang     buf_si_i[0] = 0;
361e953e47cSHong Zhang     nrows       = 0;
362e953e47cSHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
363e953e47cSHong Zhang       nzi = coi[i+1] - coi[i];
364e953e47cSHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
365e953e47cSHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
366e953e47cSHong Zhang       nrows++;
367e953e47cSHong Zhang     }
368e953e47cSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
369e953e47cSHong Zhang     k++;
370e953e47cSHong Zhang     buf_si += len_si[proc];
371e953e47cSHong Zhang   }
372e953e47cSHong Zhang   for (i=0; i<nrecv; i++) {
373e953e47cSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
374e953e47cSHong Zhang   }
375e953e47cSHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
376e953e47cSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
377e953e47cSHong Zhang 
378e953e47cSHong Zhang   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
379e953e47cSHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
380e953e47cSHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
381e953e47cSHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
382e953e47cSHong Zhang #if defined(PTAP_PROFILE)
383e953e47cSHong Zhang   ierr = PetscTime(&t2);CHKERRQ(ierr);
384e953e47cSHong Zhang #endif
385e953e47cSHong Zhang   /* (5) compute the local portion of Cmpi      */
386e953e47cSHong Zhang   /* ------------------------------------------ */
387e953e47cSHong Zhang   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
388e953e47cSHong Zhang   ierr          = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr);
389e953e47cSHong Zhang   current_space = free_space;
390e953e47cSHong Zhang 
391e953e47cSHong Zhang   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
392e953e47cSHong Zhang   for (k=0; k<nrecv; k++) {
393e953e47cSHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
394e953e47cSHong Zhang     nrows       = *buf_ri_k[k];
395e953e47cSHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
396e953e47cSHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
397e953e47cSHong Zhang   }
398e953e47cSHong Zhang 
399e953e47cSHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
400e953e47cSHong Zhang   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
401e953e47cSHong Zhang   for (i=0; i<pn; i++) {
402e953e47cSHong Zhang     /* add C_loc into Cmpi */
403e953e47cSHong Zhang     nzi  = c_loc->i[i+1] - c_loc->i[i];
404e953e47cSHong Zhang     Jptr = c_loc->j + c_loc->i[i];
405e953e47cSHong Zhang     ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
406e953e47cSHong Zhang 
407e953e47cSHong Zhang     /* add received col data into lnk */
408e953e47cSHong Zhang     for (k=0; k<nrecv; k++) { /* k-th received message */
409e953e47cSHong Zhang       if (i == *nextrow[k]) { /* i-th row */
410e953e47cSHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
411e953e47cSHong Zhang         Jptr = buf_rj[k] + *nextci[k];
412e953e47cSHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
413e953e47cSHong Zhang         nextrow[k]++; nextci[k]++;
414e953e47cSHong Zhang       }
415e953e47cSHong Zhang     }
416e953e47cSHong Zhang     nzi = lnk[0];
417e953e47cSHong Zhang 
418e953e47cSHong Zhang     /* copy data into free space, then initialize lnk */
419e953e47cSHong Zhang     ierr = PetscLLCondensedClean(pN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
420e953e47cSHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr);
421e953e47cSHong Zhang   }
422e953e47cSHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
423e953e47cSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
424e953e47cSHong Zhang   ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr);
425e953e47cSHong Zhang #if defined(PTAP_PROFILE)
426e953e47cSHong Zhang   ierr = PetscTime(&t3);CHKERRQ(ierr);
427e953e47cSHong Zhang #endif
428e953e47cSHong Zhang 
429e953e47cSHong Zhang   /* local sizes and preallocation */
430e953e47cSHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
431e953e47cSHong Zhang   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
432e953e47cSHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
433e953e47cSHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
434e953e47cSHong Zhang 
435e953e47cSHong Zhang   /* members in merge */
436e953e47cSHong Zhang   ierr = PetscFree(id_r);CHKERRQ(ierr);
437e953e47cSHong Zhang   ierr = PetscFree(len_r);CHKERRQ(ierr);
438e953e47cSHong Zhang   ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr);
439e953e47cSHong Zhang   ierr = PetscFree(buf_ri);CHKERRQ(ierr);
440e953e47cSHong Zhang   ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr);
441e953e47cSHong Zhang   ierr = PetscFree(buf_rj);CHKERRQ(ierr);
442e953e47cSHong Zhang   ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr);
443e953e47cSHong Zhang 
444e953e47cSHong Zhang   /* attach the supporting struct to Cmpi for reuse */
445e953e47cSHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
446e953e47cSHong Zhang   c->ptap         = ptap;
447e953e47cSHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
448e953e47cSHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
449e953e47cSHong Zhang 
450*9f91fb7aSHong Zhang   //if (alg == 1) {
451e953e47cSHong Zhang   ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
452*9f91fb7aSHong Zhang   //}
453e953e47cSHong Zhang 
454e953e47cSHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
455e953e47cSHong Zhang   Cmpi->assembled        = PETSC_FALSE;
456e953e47cSHong Zhang   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
457e953e47cSHong Zhang   Cmpi->ops->duplicate   = MatDuplicate_MPIAIJ_MatPtAP;
458e953e47cSHong Zhang   *C                     = Cmpi;
459e953e47cSHong Zhang 
460e953e47cSHong Zhang #if defined(PTAP_PROFILE)
461e953e47cSHong Zhang   ierr = PetscTime(&t4);CHKERRQ(ierr);
462e953e47cSHong Zhang   if (rank == 1) {
463e953e47cSHong 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);
464e953e47cSHong Zhang   }
465e953e47cSHong Zhang #endif
466e953e47cSHong Zhang   PetscFunctionReturn(0);
467e953e47cSHong Zhang }
468e953e47cSHong Zhang 
469e953e47cSHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
470e953e47cSHong Zhang {
471e953e47cSHong Zhang   PetscErrorCode      ierr;
472e953e47cSHong Zhang   Mat_PtAPMPI         *ptap;
473e953e47cSHong Zhang   Mat_MPIAIJ          *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
474e953e47cSHong Zhang   MPI_Comm            comm;
475e953e47cSHong Zhang   PetscMPIInt         size,rank;
476e953e47cSHong Zhang   Mat                 Cmpi;
477e953e47cSHong Zhang   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
478e953e47cSHong Zhang   PetscInt            am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
479e953e47cSHong Zhang   PetscInt            *lnk,i,k,pnz,row,nsend;
480e953e47cSHong Zhang   PetscBT             lnkbt;
481e953e47cSHong Zhang   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
482e953e47cSHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
483e953e47cSHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners,nzi,nspacedouble;
484e953e47cSHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
485e953e47cSHong Zhang   MPI_Request         *swaits,*rwaits;
486e953e47cSHong Zhang   MPI_Status          *sstatus,rstatus;
487e953e47cSHong Zhang   PetscLayout         rowmap;
488e953e47cSHong Zhang   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
489e953e47cSHong Zhang   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
490e953e47cSHong Zhang   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0,Crmax,*aj,*ai,*pi;
491e953e47cSHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*c_loc,*c_oth;
492e953e47cSHong Zhang   PetscScalar         *apv;
493e953e47cSHong Zhang   PetscTable          ta;
494e953e47cSHong Zhang #if defined(PETSC_HAVE_HYPRE)
495e953e47cSHong Zhang   const char          *algTypes[3] = {"scalable","nonscalable","hypre"};
496e953e47cSHong Zhang   PetscInt            nalg = 3;
497e953e47cSHong Zhang #else
498e953e47cSHong Zhang   const char          *algTypes[2] = {"scalable","nonscalable"};
499e953e47cSHong Zhang   PetscInt            nalg = 2;
500e953e47cSHong Zhang #endif
501e953e47cSHong Zhang   PetscInt            alg = 1; /* set default algorithm */
502e953e47cSHong Zhang #if defined(PETSC_USE_INFO)
503e953e47cSHong Zhang   PetscReal           apfill;
504e953e47cSHong Zhang #endif
505e953e47cSHong Zhang #if defined(PTAP_PROFILE)
506e953e47cSHong Zhang   PetscLogDouble      t0,t1,t11,t12,t2,t3,t4;
507e953e47cSHong Zhang #endif
508e953e47cSHong Zhang 
509e953e47cSHong Zhang   PetscFunctionBegin;
510e953e47cSHong Zhang   /* pick an algorithm */
511ae5f0867Sstefano_zampini   ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
512ae5f0867Sstefano_zampini   ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[1],&alg,NULL);CHKERRQ(ierr);
513ae5f0867Sstefano_zampini   ierr = PetscOptionsEnd();CHKERRQ(ierr);
514ae5f0867Sstefano_zampini 
515e953e47cSHong Zhang   if (alg == 0) {
516e953e47cSHong Zhang     ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(A,P,fill,C);CHKERRQ(ierr);
517e953e47cSHong Zhang     (*C)->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ; /* _scalable; not done yet */
518e953e47cSHong Zhang     PetscFunctionReturn(0);
519e953e47cSHong Zhang 
520ae5f0867Sstefano_zampini #if defined(PETSC_HAVE_HYPRE)
521ae5f0867Sstefano_zampini   } else if (alg == 2) {
522ae5f0867Sstefano_zampini     /* Use boomerAMGBuildCoarseOperator */
523ae5f0867Sstefano_zampini     ierr = MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);CHKERRQ(ierr);
524ae5f0867Sstefano_zampini     PetscFunctionReturn(0);
525ae5f0867Sstefano_zampini #endif
526ae5f0867Sstefano_zampini   }
527ae5f0867Sstefano_zampini 
5288cb82516SHong Zhang #if defined(PTAP_PROFILE)
52980bb4639SHong Zhang   ierr = PetscTime(&t0);CHKERRQ(ierr);
5308cb82516SHong Zhang #endif
5318cb82516SHong Zhang 
532e953e47cSHong Zhang   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
533e953e47cSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
534e953e47cSHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
535e953e47cSHong Zhang 
536e953e47cSHong Zhang   /* create symbolic parallel matrix Cmpi */
537e953e47cSHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
538e953e47cSHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
539e953e47cSHong Zhang 
540e953e47cSHong Zhang   /* Do dense axpy in MatPtAPNumeric_MPIAIJ_MPIAIJ() */
541e953e47cSHong Zhang   Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ;
542e953e47cSHong Zhang 
54315a3b8e2SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
54415a3b8e2SHong Zhang   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
54515a3b8e2SHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
54615a3b8e2SHong Zhang 
54715a3b8e2SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
54815a3b8e2SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
54915a3b8e2SHong Zhang   /* get P_loc by taking all local rows of P */
55015a3b8e2SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
55115a3b8e2SHong Zhang 
55267a12041SHong Zhang   /* (0) compute Rd = Pd^T, Ro = Po^T  */
55367a12041SHong Zhang   /* --------------------------------- */
554de817e96SHong Zhang   ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr);
555de817e96SHong Zhang   ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr);
5568cb82516SHong Zhang #if defined(PTAP_PROFILE)
557445158ffSHong Zhang   ierr = PetscTime(&t11);CHKERRQ(ierr);
5588cb82516SHong Zhang #endif
55915a3b8e2SHong Zhang 
56067a12041SHong Zhang   /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */
56167a12041SHong Zhang   /* ----------------------------------------------------------------- */
56267a12041SHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
56367a12041SHong Zhang   p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
564445158ffSHong Zhang 
5659a6dcab7SHong Zhang   /* create and initialize a linked list */
56645d00d1dSHong Zhang   ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */
5674b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta);
5684b38b95cSHong Zhang   MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta);
56978d30b94SHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */
570d0e9a2f7SHong 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); */
57178d30b94SHong Zhang 
57278d30b94SHong Zhang   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
573445158ffSHong Zhang 
5748cb82516SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
575f91af8c7SBarry Smith   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr);
576445158ffSHong Zhang   current_space = free_space;
57767a12041SHong Zhang   nspacedouble  = 0;
57867a12041SHong Zhang 
57967a12041SHong Zhang   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
58067a12041SHong Zhang   api[0] = 0;
581445158ffSHong Zhang   for (i=0; i<am; i++) {
58267a12041SHong Zhang     /* diagonal portion: Ad[i,:]*P */
58367a12041SHong Zhang     ai = ad->i; pi = p_loc->i;
58467a12041SHong Zhang     nzi = ai[i+1] - ai[i];
58567a12041SHong Zhang     aj  = ad->j + ai[i];
586445158ffSHong Zhang     for (j=0; j<nzi; j++) {
587445158ffSHong Zhang       row  = aj[j];
58867a12041SHong Zhang       pnz  = pi[row+1] - pi[row];
58967a12041SHong Zhang       Jptr = p_loc->j + pi[row];
590445158ffSHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
591445158ffSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
592445158ffSHong Zhang     }
59367a12041SHong Zhang     /* off-diagonal portion: Ao[i,:]*P */
59467a12041SHong Zhang     ai = ao->i; pi = p_oth->i;
59567a12041SHong Zhang     nzi = ai[i+1] - ai[i];
59667a12041SHong Zhang     aj  = ao->j + ai[i];
597445158ffSHong Zhang     for (j=0; j<nzi; j++) {
598445158ffSHong Zhang       row  = aj[j];
59967a12041SHong Zhang       pnz  = pi[row+1] - pi[row];
60067a12041SHong Zhang       Jptr = p_oth->j + pi[row];
601445158ffSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
602445158ffSHong Zhang     }
603445158ffSHong Zhang     apnz     = lnk[0];
604445158ffSHong Zhang     api[i+1] = api[i] + apnz;
605445158ffSHong Zhang     if (ap_rmax < apnz) ap_rmax = apnz;
606445158ffSHong Zhang 
607445158ffSHong Zhang     /* if free space is not available, double the total space in the list */
608445158ffSHong Zhang     if (current_space->local_remaining<apnz) {
609f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
610445158ffSHong Zhang       nspacedouble++;
611445158ffSHong Zhang     }
612445158ffSHong Zhang 
613445158ffSHong Zhang     /* Copy data into free space, then initialize lnk */
614445158ffSHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
615445158ffSHong Zhang 
616445158ffSHong Zhang     current_space->array           += apnz;
617445158ffSHong Zhang     current_space->local_used      += apnz;
618445158ffSHong Zhang     current_space->local_remaining -= apnz;
619445158ffSHong Zhang   }
620681d504bSHong Zhang   /* Allocate space for apj and apv, initialize apj, and */
621445158ffSHong Zhang   /* destroy list of free space and other temporary array(s) */
622445158ffSHong Zhang   ierr   = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr);
623445158ffSHong Zhang   ierr   = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
6249a6dcab7SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
6259a6dcab7SHong Zhang 
626aa690a28SHong Zhang   /* Create AP_loc for reuse */
627445158ffSHong Zhang   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr);
628aa690a28SHong Zhang 
629aa690a28SHong Zhang #if defined(PETSC_USE_INFO)
630aa690a28SHong Zhang   apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1);
631aa690a28SHong Zhang   ptap->AP_loc->info.mallocs           = nspacedouble;
632aa690a28SHong Zhang   ptap->AP_loc->info.fill_ratio_given  = fill;
633aa690a28SHong Zhang   ptap->AP_loc->info.fill_ratio_needed = apfill;
634aa690a28SHong Zhang 
635aa690a28SHong Zhang   if (api[am]) {
636aa690a28SHong Zhang     ierr = PetscInfo3(ptap->AP_loc,"AP_loc reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)apfill);CHKERRQ(ierr);
637aa690a28SHong Zhang     ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr);
638aa690a28SHong Zhang   } else {
639aa690a28SHong Zhang     ierr = PetscInfo(ptap->AP_loc,"AP_loc is empty \n");CHKERRQ(ierr);
640aa690a28SHong Zhang   }
641aa690a28SHong Zhang #endif
642aa690a28SHong Zhang 
6438cb82516SHong Zhang #if defined(PTAP_PROFILE)
644445158ffSHong Zhang   ierr = PetscTime(&t12);CHKERRQ(ierr);
6458cb82516SHong Zhang #endif
64615a3b8e2SHong Zhang 
647681d504bSHong Zhang   /* (2-1) compute symbolic Co = Ro*AP_loc  */
648681d504bSHong Zhang   /* ------------------------------------ */
64967a12041SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr);
6508cb82516SHong Zhang #if defined(PTAP_PROFILE)
65180bb4639SHong Zhang   ierr = PetscTime(&t1);CHKERRQ(ierr);
6528cb82516SHong Zhang #endif
65315a3b8e2SHong Zhang 
65467a12041SHong Zhang   /* (3) send coj of C_oth to other processors  */
65567a12041SHong Zhang   /* ------------------------------------------ */
65615a3b8e2SHong Zhang   /* determine row ownership */
657445158ffSHong Zhang   ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr);
658445158ffSHong Zhang   rowmap->n  = pn;
659445158ffSHong Zhang   rowmap->bs = 1;
660445158ffSHong Zhang   ierr   = PetscLayoutSetUp(rowmap);CHKERRQ(ierr);
661445158ffSHong Zhang   owners = rowmap->range;
66215a3b8e2SHong Zhang 
66315a3b8e2SHong Zhang   /* determine the number of messages to send, their lengths */
6648cb82516SHong Zhang   ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr);
6658cb82516SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
66615a3b8e2SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
66715a3b8e2SHong Zhang 
66867a12041SHong Zhang   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
66967a12041SHong Zhang   coi   = c_oth->i; coj = c_oth->j;
67067a12041SHong Zhang   con   = ptap->C_oth->rmap->n;
67115a3b8e2SHong Zhang   proc  = 0;
67267a12041SHong Zhang   for (i=0; i<con; i++) {
67315a3b8e2SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
67415a3b8e2SHong Zhang     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
67515a3b8e2SHong Zhang     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
67615a3b8e2SHong Zhang   }
67715a3b8e2SHong Zhang 
67815a3b8e2SHong Zhang   len          = 0; /* max length of buf_si[], see (4) */
67915a3b8e2SHong Zhang   owners_co[0] = 0;
68067a12041SHong Zhang   nsend        = 0;
68115a3b8e2SHong Zhang   for (proc=0; proc<size; proc++) {
68215a3b8e2SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
68315a3b8e2SHong Zhang     if (len_s[proc]) {
684445158ffSHong Zhang       nsend++;
68515a3b8e2SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
68615a3b8e2SHong Zhang       len         += len_si[proc];
68715a3b8e2SHong Zhang     }
68815a3b8e2SHong Zhang   }
68915a3b8e2SHong Zhang 
69015a3b8e2SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
691445158ffSHong Zhang   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr);
692445158ffSHong Zhang   ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr);
69315a3b8e2SHong Zhang 
69415a3b8e2SHong Zhang   /* post the Irecv and Isend of coj */
69515a3b8e2SHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
696445158ffSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
697445158ffSHong Zhang   ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr);
69815a3b8e2SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
69915a3b8e2SHong Zhang     if (!len_s[proc]) continue;
70015a3b8e2SHong Zhang     i    = owners_co[proc];
70115a3b8e2SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
70215a3b8e2SHong Zhang     k++;
70315a3b8e2SHong Zhang   }
70415a3b8e2SHong Zhang 
705681d504bSHong Zhang   /* (2-2) compute symbolic C_loc = Rd*AP_loc */
706681d504bSHong Zhang   /* ---------------------------------------- */
707681d504bSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr);
708681d504bSHong Zhang   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
709681d504bSHong Zhang 
710681d504bSHong Zhang   /* receives coj are complete */
711445158ffSHong Zhang   for (i=0; i<nrecv; i++) {
712445158ffSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
71315a3b8e2SHong Zhang   }
71415a3b8e2SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
715445158ffSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
71615a3b8e2SHong Zhang 
71778d30b94SHong Zhang   /* add received column indices into ta to update Crmax */
71878d30b94SHong Zhang   for (k=0; k<nrecv; k++) {/* k-th received message */
71978d30b94SHong Zhang     Jptr = buf_rj[k];
72078d30b94SHong Zhang     for (j=0; j<len_r[k]; j++) {
72178d30b94SHong Zhang       ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr);
72278d30b94SHong Zhang     }
72378d30b94SHong Zhang   }
72478d30b94SHong Zhang   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
72578d30b94SHong Zhang   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
7269a6dcab7SHong Zhang 
72715a3b8e2SHong Zhang   /* (4) send and recv coi */
72815a3b8e2SHong Zhang   /*-----------------------*/
72915a3b8e2SHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
730445158ffSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
73115a3b8e2SHong Zhang   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
73215a3b8e2SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
73315a3b8e2SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
73415a3b8e2SHong Zhang     if (!len_s[proc]) continue;
73515a3b8e2SHong Zhang     /* form outgoing message for i-structure:
73615a3b8e2SHong Zhang          buf_si[0]:                 nrows to be sent
73715a3b8e2SHong Zhang                [1:nrows]:           row index (global)
73815a3b8e2SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
73915a3b8e2SHong Zhang     */
74015a3b8e2SHong Zhang     /*-------------------------------------------*/
74115a3b8e2SHong Zhang     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
74215a3b8e2SHong Zhang     buf_si_i    = buf_si + nrows+1;
74315a3b8e2SHong Zhang     buf_si[0]   = nrows;
74415a3b8e2SHong Zhang     buf_si_i[0] = 0;
74515a3b8e2SHong Zhang     nrows       = 0;
74615a3b8e2SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
74715a3b8e2SHong Zhang       nzi = coi[i+1] - coi[i];
74815a3b8e2SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
74915a3b8e2SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
75015a3b8e2SHong Zhang       nrows++;
75115a3b8e2SHong Zhang     }
75215a3b8e2SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
75315a3b8e2SHong Zhang     k++;
75415a3b8e2SHong Zhang     buf_si += len_si[proc];
75515a3b8e2SHong Zhang   }
756681d504bSHong Zhang   for (i=0; i<nrecv; i++) {
757445158ffSHong Zhang     ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
75815a3b8e2SHong Zhang   }
75915a3b8e2SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
760445158ffSHong Zhang   if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);}
76115a3b8e2SHong Zhang 
7628cb82516SHong Zhang   ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr);
76315a3b8e2SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
76415a3b8e2SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
76515a3b8e2SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
7668cb82516SHong Zhang #if defined(PTAP_PROFILE)
76780bb4639SHong Zhang   ierr = PetscTime(&t2);CHKERRQ(ierr);
7688cb82516SHong Zhang #endif
76967a12041SHong Zhang   /* (5) compute the local portion of Cmpi      */
77067a12041SHong Zhang   /* ------------------------------------------ */
77178d30b94SHong Zhang   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
77278d30b94SHong Zhang   ierr          = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr);
77315a3b8e2SHong Zhang   current_space = free_space;
77415a3b8e2SHong Zhang 
775445158ffSHong Zhang   ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr);
776445158ffSHong Zhang   for (k=0; k<nrecv; k++) {
77715a3b8e2SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
77815a3b8e2SHong Zhang     nrows       = *buf_ri_k[k];
77915a3b8e2SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
78015a3b8e2SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
78115a3b8e2SHong Zhang   }
782445158ffSHong Zhang 
78378d30b94SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
78478d30b94SHong Zhang   ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr);
78515a3b8e2SHong Zhang   for (i=0; i<pn; i++) {
78667a12041SHong Zhang     /* add C_loc into Cmpi */
78767a12041SHong Zhang     nzi  = c_loc->i[i+1] - c_loc->i[i];
78867a12041SHong Zhang     Jptr = c_loc->j + c_loc->i[i];
78967a12041SHong Zhang     ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
79015a3b8e2SHong Zhang 
79115a3b8e2SHong Zhang     /* add received col data into lnk */
792445158ffSHong Zhang     for (k=0; k<nrecv; k++) { /* k-th received message */
79315a3b8e2SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
79415a3b8e2SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
79515a3b8e2SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
79615a3b8e2SHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
79715a3b8e2SHong Zhang         nextrow[k]++; nextci[k]++;
79815a3b8e2SHong Zhang       }
79915a3b8e2SHong Zhang     }
800d0e9a2f7SHong Zhang     nzi = lnk[0];
8018cb82516SHong Zhang 
80215a3b8e2SHong Zhang     /* copy data into free space, then initialize lnk */
803d0e9a2f7SHong Zhang     ierr = PetscLLCondensedClean(pN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
804d0e9a2f7SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr);
80515a3b8e2SHong Zhang   }
80615a3b8e2SHong Zhang   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
80715a3b8e2SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
808445158ffSHong Zhang   ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr);
8098cb82516SHong Zhang #if defined(PTAP_PROFILE)
81080bb4639SHong Zhang   ierr = PetscTime(&t3);CHKERRQ(ierr);
8118cb82516SHong Zhang #endif
81280bb4639SHong Zhang 
813ae5f0867Sstefano_zampini   /* local sizes and preallocation */
81415a3b8e2SHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
81515a3b8e2SHong Zhang   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
81615a3b8e2SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
81715a3b8e2SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
81815a3b8e2SHong Zhang 
819445158ffSHong Zhang   /* members in merge */
820445158ffSHong Zhang   ierr = PetscFree(id_r);CHKERRQ(ierr);
821445158ffSHong Zhang   ierr = PetscFree(len_r);CHKERRQ(ierr);
822445158ffSHong Zhang   ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr);
823445158ffSHong Zhang   ierr = PetscFree(buf_ri);CHKERRQ(ierr);
824445158ffSHong Zhang   ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr);
825445158ffSHong Zhang   ierr = PetscFree(buf_rj);CHKERRQ(ierr);
826445158ffSHong Zhang   ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr);
82715a3b8e2SHong Zhang 
82815a3b8e2SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
82915a3b8e2SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
83015a3b8e2SHong Zhang   c->ptap         = ptap;
8311a47ec54SHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
8321a47ec54SHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
83315a3b8e2SHong Zhang 
83478d30b94SHong Zhang   if (alg == 1) {
8358cb82516SHong Zhang     ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
83678d30b94SHong Zhang   }
8372259aa2eSHong Zhang 
8381a47ec54SHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
8391a47ec54SHong Zhang   Cmpi->assembled        = PETSC_FALSE;
8401a47ec54SHong Zhang   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
841aa690a28SHong Zhang   Cmpi->ops->duplicate   = MatDuplicate_MPIAIJ_MatPtAP;
8421a47ec54SHong Zhang   *C                     = Cmpi;
8431a47ec54SHong Zhang 
8448cb82516SHong Zhang #if defined(PTAP_PROFILE)
84580bb4639SHong Zhang   ierr = PetscTime(&t4);CHKERRQ(ierr);
846e9c1f85fSHong Zhang   if (rank == 1) {
847445158ffSHong 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);
848e9c1f85fSHong Zhang   }
8498cb82516SHong Zhang #endif
8500d3441aeSHong Zhang   PetscFunctionReturn(0);
8510d3441aeSHong Zhang }
8520d3441aeSHong Zhang 
853aa690a28SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
8540d3441aeSHong Zhang {
8550d3441aeSHong Zhang   PetscErrorCode    ierr;
8560d3441aeSHong Zhang   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
8570d3441aeSHong Zhang   Mat_SeqAIJ        *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
8588cb82516SHong Zhang   Mat_SeqAIJ        *ap,*p_loc,*p_oth,*c_seq;
8590d3441aeSHong Zhang   Mat_PtAPMPI       *ptap = c->ptap;
8609ce11a7cSHong Zhang   Mat               AP_loc,C_loc,C_oth;
8610d3441aeSHong Zhang   PetscInt          i,rstart,rend,cm,ncols,row;
8628cb82516SHong Zhang   PetscInt          *api,*apj,am = A->rmap->n,j,col,apnz;
8638cb82516SHong Zhang   PetscScalar       *apa;
8640d3441aeSHong Zhang   const PetscInt    *cols;
8650d3441aeSHong Zhang   const PetscScalar *vals;
8668cb82516SHong Zhang #if defined(PTAP_PROFILE)
8678cb82516SHong Zhang   PetscMPIInt       rank;
8688cb82516SHong Zhang   MPI_Comm          comm;
8690d3441aeSHong Zhang   PetscLogDouble    t0,t1,t2,t3,t4,eR,eAP,eCseq,eCmpi;
8708cb82516SHong Zhang #endif
8710d3441aeSHong Zhang 
8720d3441aeSHong Zhang   PetscFunctionBegin;
8730d3441aeSHong Zhang   ierr = MatZeroEntries(C);CHKERRQ(ierr);
8740d3441aeSHong Zhang 
875e497d3c8SHong Zhang   /* 1) get R = Pd^T,Ro = Po^T */
8768cb82516SHong Zhang #if defined(PTAP_PROFILE)
8770d3441aeSHong Zhang   ierr = PetscTime(&t0);CHKERRQ(ierr);
8788cb82516SHong Zhang #endif
879748c7196SHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
8800d3441aeSHong Zhang     ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr);
8810d3441aeSHong Zhang     ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr);
882748c7196SHong Zhang   }
8838cb82516SHong Zhang #if defined(PTAP_PROFILE)
8840d3441aeSHong Zhang   ierr = PetscTime(&t1);CHKERRQ(ierr);
8850d3441aeSHong Zhang   eR = t1 - t0;
8868cb82516SHong Zhang #endif
8870d3441aeSHong Zhang 
888e497d3c8SHong Zhang   /* 2) get AP_loc */
8890d3441aeSHong Zhang   AP_loc = ptap->AP_loc;
8908cb82516SHong Zhang   ap = (Mat_SeqAIJ*)AP_loc->data;
8910d3441aeSHong Zhang 
892e497d3c8SHong Zhang   /* 2-1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
8930d3441aeSHong Zhang   /*-----------------------------------------------------*/
894748c7196SHong Zhang   if (ptap->reuse == MAT_REUSE_MATRIX) {
895748c7196SHong Zhang     /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */
8960d3441aeSHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
8970d3441aeSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
898e497d3c8SHong Zhang   }
8990d3441aeSHong Zhang 
9008cb82516SHong Zhang   /* 2-2) compute numeric A_loc*P - dominating part */
9018cb82516SHong Zhang   /* ---------------------------------------------- */
9020d3441aeSHong Zhang   /* get data from symbolic products */
9038cb82516SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
9048cb82516SHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
9058cb82516SHong Zhang   apa   = ptap->apa;
906681d504bSHong Zhang   api   = ap->i;
907681d504bSHong Zhang   apj   = ap->j;
908e497d3c8SHong Zhang   for (i=0; i<am; i++) {
909445158ffSHong Zhang     /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */
910e497d3c8SHong Zhang     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
911e497d3c8SHong Zhang     apnz = api[i+1] - api[i];
912e497d3c8SHong Zhang     for (j=0; j<apnz; j++) {
913e497d3c8SHong Zhang       col = apj[j+api[i]];
914e497d3c8SHong Zhang       ap->a[j+ap->i[i]] = apa[col];
915e497d3c8SHong Zhang       apa[col] = 0.0;
9160d3441aeSHong Zhang     }
917aa690a28SHong Zhang     ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
9180d3441aeSHong Zhang   }
9198cb82516SHong Zhang #if defined(PTAP_PROFILE)
9200d3441aeSHong Zhang   ierr = PetscTime(&t2);CHKERRQ(ierr);
9210d3441aeSHong Zhang   eAP = t2 - t1;
9228cb82516SHong Zhang #endif
9230d3441aeSHong Zhang 
9248cb82516SHong Zhang   /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */
925445158ffSHong Zhang   ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr);
926445158ffSHong Zhang   ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr);
9279ce11a7cSHong Zhang   C_loc = ptap->C_loc;
9289ce11a7cSHong Zhang   C_oth = ptap->C_oth;
9298cb82516SHong Zhang #if defined(PTAP_PROFILE)
9300d3441aeSHong Zhang   ierr = PetscTime(&t3);CHKERRQ(ierr);
9310d3441aeSHong Zhang   eCseq = t3 - t2;
9328cb82516SHong Zhang #endif
9330d3441aeSHong Zhang 
9340d3441aeSHong Zhang   /* add C_loc and Co to to C */
9350d3441aeSHong Zhang   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
9360d3441aeSHong Zhang 
9370d3441aeSHong Zhang   /* C_loc -> C */
9380d3441aeSHong Zhang   cm    = C_loc->rmap->N;
9399ce11a7cSHong Zhang   c_seq = (Mat_SeqAIJ*)C_loc->data;
9408cb82516SHong Zhang   cols = c_seq->j;
9418cb82516SHong Zhang   vals = c_seq->a;
9420d3441aeSHong Zhang   for (i=0; i<cm; i++) {
9439ce11a7cSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
9440d3441aeSHong Zhang     row = rstart + i;
9450d3441aeSHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
9468cb82516SHong Zhang     cols += ncols; vals += ncols;
9470d3441aeSHong Zhang   }
9480d3441aeSHong Zhang 
9490d3441aeSHong Zhang   /* Co -> C, off-processor part */
9509ce11a7cSHong Zhang   cm = C_oth->rmap->N;
9519ce11a7cSHong Zhang   c_seq = (Mat_SeqAIJ*)C_oth->data;
9528cb82516SHong Zhang   cols = c_seq->j;
9538cb82516SHong Zhang   vals = c_seq->a;
9549ce11a7cSHong Zhang   for (i=0; i<cm; i++) {
9559ce11a7cSHong Zhang     ncols = c_seq->i[i+1] - c_seq->i[i];
9560d3441aeSHong Zhang     row = p->garray[i];
9570d3441aeSHong Zhang     ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr);
9588cb82516SHong Zhang     cols += ncols; vals += ncols;
9590d3441aeSHong Zhang   }
9600d3441aeSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9610d3441aeSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9628cb82516SHong Zhang #if defined(PTAP_PROFILE)
9630d3441aeSHong Zhang   ierr = PetscTime(&t4);CHKERRQ(ierr);
9640d3441aeSHong Zhang   eCmpi = t4 - t3;
9650d3441aeSHong Zhang 
9668cb82516SHong Zhang   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
9678cb82516SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
9680d3441aeSHong Zhang   if (rank==1) {
9690d3441aeSHong 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);
9700d3441aeSHong Zhang   }
9718cb82516SHong Zhang #endif
972748c7196SHong Zhang   ptap->reuse = MAT_REUSE_MATRIX;
9730d3441aeSHong Zhang   PetscFunctionReturn(0);
9740d3441aeSHong Zhang }
9750d3441aeSHong Zhang 
9768403a639SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(Mat A,Mat P,PetscReal fill,Mat *C)
97742c57489SHong Zhang {
97842c57489SHong Zhang   PetscErrorCode      ierr;
97977584fe6SHong Zhang   Mat                 Cmpi;
980f8487c73SHong Zhang   Mat_PtAPMPI         *ptap;
9810298fd71SBarry Smith   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
982f8487c73SHong Zhang   Mat_MPIAIJ          *a        =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
98342c57489SHong Zhang   Mat_SeqAIJ          *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
98442c57489SHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth;
985ded4f561SHong Zhang   PetscInt            *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ;
98677584fe6SHong Zhang   PetscInt            *adi=ad->i,*aj,*aoi=ao->i,nnz;
987a914927fSHong Zhang   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
988d16ca5f0SHong Zhang   PetscInt            am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n;
98942c57489SHong Zhang   PetscBT             lnkbt;
990ce94432eSBarry Smith   MPI_Comm            comm;
991a914927fSHong Zhang   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0;
99242c57489SHong Zhang   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
99342c57489SHong Zhang   PetscInt            len,proc,*dnz,*onz,*owners;
99424ecddacSHong Zhang   PetscInt            nzi,*pti,*ptj;
99542c57489SHong Zhang   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
996ded4f561SHong Zhang   MPI_Request         *swaits,*rwaits;
997ded4f561SHong Zhang   MPI_Status          *sstatus,rstatus;
99842c57489SHong Zhang   Mat_Merge_SeqsToMPI *merge;
99977584fe6SHong Zhang   PetscInt            *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0;
1000d16ca5f0SHong Zhang   PetscReal           afill=1.0,afill_tmp;
100142c57489SHong Zhang 
100242c57489SHong Zhang   PetscFunctionBegin;
1003ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
100483445d95SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
100583445d95SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
100683445d95SHong Zhang 
1007f8487c73SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
1008b00a9115SJed Brown   ierr        = PetscNew(&ptap);CHKERRQ(ierr);
1009b00a9115SJed Brown   ierr        = PetscNew(&merge);CHKERRQ(ierr);
10109f4155fbSHong Zhang   ptap->merge = merge;
1011f8487c73SHong Zhang   ptap->reuse = MAT_INITIAL_MATRIX;
10126c6e00e0SHong Zhang 
10136c6e00e0SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
1014b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
1015fe615159SHong Zhang 
10166c6e00e0SHong Zhang   /* get P_loc by taking all local rows of P */
1017a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
10186c6e00e0SHong Zhang 
1019a1a4e84aSHong Zhang   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
1020a1a4e84aSHong Zhang   p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
10216c6e00e0SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
10226c6e00e0SHong Zhang   pi_oth = p_oth->i; pj_oth = p_oth->j;
102342c57489SHong Zhang 
10242addb229SHong Zhang   /* (1) compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth (api,apj) */
10252addb229SHong Zhang   /*--------------------------------------------------------------------------*/
1026854ce69bSBarry Smith   ierr   = PetscMalloc1(am+1,&api);CHKERRQ(ierr);
102782412ba7SHong Zhang   api[0] = 0;
102883445d95SHong Zhang 
102983445d95SHong Zhang   /* create and initialize a linked list */
1030a914927fSHong Zhang   ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr);
103183445d95SHong Zhang 
1032a914927fSHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) -OOM for ex56, np=8k on Intrepid! */
1033f91af8c7SBarry Smith   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr);
1034f4a743e1SHong Zhang   current_space = free_space;
1035f4a743e1SHong Zhang 
1036f4a743e1SHong Zhang   for (i=0; i<am; i++) {
1037f4a743e1SHong Zhang     /* diagonal portion of A */
1038ded4f561SHong Zhang     nzi = adi[i+1] - adi[i];
103977584fe6SHong Zhang     aj  = ad->j + adi[i];
1040ded4f561SHong Zhang     for (j=0; j<nzi; j++) {
104177584fe6SHong Zhang       row  = aj[j];
104282412ba7SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
1043ded4f561SHong Zhang       Jptr = pj_loc + pi_loc[row];
104483445d95SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
1045a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1046f4a743e1SHong Zhang     }
1047f4a743e1SHong Zhang     /* off-diagonal portion of A */
1048ded4f561SHong Zhang     nzi = aoi[i+1] - aoi[i];
104977584fe6SHong Zhang     aj  = ao->j + aoi[i];
1050ded4f561SHong Zhang     for (j=0; j<nzi; j++) {
105177584fe6SHong Zhang       row  = aj[j];
105282412ba7SHong Zhang       pnz  = pi_oth[row+1] - pi_oth[row];
1053ded4f561SHong Zhang       Jptr = pj_oth + pi_oth[row];
1054a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1055f4a743e1SHong Zhang     }
1056a914927fSHong Zhang     apnz     = lnk[0];
105782412ba7SHong Zhang     api[i+1] = api[i] + apnz;
105877584fe6SHong Zhang     if (ap_rmax < apnz) ap_rmax = apnz;
1059f4a743e1SHong Zhang 
106083445d95SHong Zhang     /* if free space is not available, double the total space in the list */
106182412ba7SHong Zhang     if (current_space->local_remaining<apnz) {
1062f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1063f2b054eeSHong Zhang       nspacedouble++;
1064f4a743e1SHong Zhang     }
1065f4a743e1SHong Zhang 
1066f4a743e1SHong Zhang     /* Copy data into free space, then initialize lnk */
1067a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
10682205254eSKarl Rupp 
106982412ba7SHong Zhang     current_space->array           += apnz;
107082412ba7SHong Zhang     current_space->local_used      += apnz;
107182412ba7SHong Zhang     current_space->local_remaining -= apnz;
1072f4a743e1SHong Zhang   }
1073a914927fSHong Zhang 
107482412ba7SHong Zhang   /* Allocate space for apj, initialize apj, and */
1075f4a743e1SHong Zhang   /* destroy list of free space and other temporary array(s) */
1076854ce69bSBarry Smith   ierr      = PetscMalloc1(api[am]+1,&apj);CHKERRQ(ierr);
107777584fe6SHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
1078118727c9SMark F. Adams   afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1);
1079d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1080f4a743e1SHong Zhang 
10812addb229SHong Zhang   /* (2) determine symbolic Co=(p->B)^T*AP - send to others (coi,coj)*/
10822addb229SHong Zhang   /*-----------------------------------------------------------------*/
1083de0260b3SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
108442c57489SHong Zhang 
1085ded4f561SHong Zhang   /* then, compute symbolic Co = (p->B)^T*AP */
1086d0f46423SBarry Smith   pon    = (p->B)->cmap->n; /* total num of rows to be sent to other processors
108783445d95SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
1088854ce69bSBarry Smith   ierr   = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr);
1089de0260b3SHong Zhang   coi[0] = 0;
109042c57489SHong Zhang 
1091d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */
1092f91af8c7SBarry Smith   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(poti[pon],api[am]));
109322d28d08SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
109442c57489SHong Zhang   current_space = free_space;
109542c57489SHong Zhang 
1096de0260b3SHong Zhang   for (i=0; i<pon; i++) {
109782412ba7SHong Zhang     pnz = poti[i+1] - poti[i];
109877584fe6SHong Zhang     ptJ = potj + poti[i];
109977584fe6SHong Zhang     for (j=0; j<pnz; j++) {
110077584fe6SHong Zhang       row  = ptJ[j]; /* row of AP == col of Pot */
110182412ba7SHong Zhang       apnz = api[row+1] - api[row];
1102ded4f561SHong Zhang       Jptr = apj + api[row];
110383445d95SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1104a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
110542c57489SHong Zhang     }
1106a914927fSHong Zhang     nnz = lnk[0];
110742c57489SHong Zhang 
110883445d95SHong Zhang     /* If free space is not available, double the total space in the list */
1109ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
1110f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1111d16ca5f0SHong Zhang       nspacedouble++;
111242c57489SHong Zhang     }
111342c57489SHong Zhang 
111442c57489SHong Zhang     /* Copy data into free space, and zero out denserows */
1115a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
11162205254eSKarl Rupp 
1117ded4f561SHong Zhang     current_space->array           += nnz;
1118ded4f561SHong Zhang     current_space->local_used      += nnz;
1119ded4f561SHong Zhang     current_space->local_remaining -= nnz;
11202205254eSKarl Rupp 
1121ded4f561SHong Zhang     coi[i+1] = coi[i] + nnz;
112242c57489SHong Zhang   }
11236b911d16SHong Zhang 
11246b911d16SHong Zhang   ierr      = PetscMalloc1(coi[pon],&coj);CHKERRQ(ierr);
1125a1a86e44SBarry Smith   ierr      = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
1126118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1);
1127d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
1128de0260b3SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
112942c57489SHong Zhang 
11306b911d16SHong Zhang   /* (3) send j-array (coj) of Co to other processors */
11316b911d16SHong Zhang   /*--------------------------------------------------*/
11326b911d16SHong Zhang   ierr = PetscCalloc1(size,&merge->len_s);CHKERRQ(ierr);
11336b911d16SHong Zhang   len_s        = merge->len_s;
11346b911d16SHong Zhang   merge->nsend = 0;
11356b911d16SHong Zhang 
11366b911d16SHong Zhang 
113742c57489SHong Zhang   /* determine row ownership */
113826283091SBarry Smith   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
11397a2fc3feSBarry Smith   merge->rowmap->n  = pn;
11407a2fc3feSBarry Smith   merge->rowmap->bs = 1;
11412205254eSKarl Rupp 
114226283091SBarry Smith   ierr   = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
11437a2fc3feSBarry Smith   owners = merge->rowmap->range;
114442c57489SHong Zhang 
114542c57489SHong Zhang   /* determine the number of messages to send, their lengths */
1146dcca6d9dSJed Brown   ierr = PetscMalloc2(size,&len_si,size,&sstatus);CHKERRQ(ierr);
114783445d95SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
1148854ce69bSBarry Smith   ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr);
1149de0260b3SHong Zhang 
115083445d95SHong Zhang   proc = 0;
1151de0260b3SHong Zhang   for (i=0; i<pon; i++) {
1152de0260b3SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
11532addb229SHong Zhang     len_si[proc]++;               /* num of rows in Co(=Pt*AP) to be sent to [proc] */
1154ee6e4b5bSHong Zhang     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
115583445d95SHong Zhang   }
1156de0260b3SHong Zhang 
11572addb229SHong Zhang   len          = 0; /* max length of buf_si[], see (4) */
1158de0260b3SHong Zhang   owners_co[0] = 0;
115942c57489SHong Zhang   for (proc=0; proc<size; proc++) {
1160de0260b3SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1161ee6e4b5bSHong Zhang     if (len_s[proc]) {
116242c57489SHong Zhang       merge->nsend++;
11632addb229SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
116442c57489SHong Zhang       len         += len_si[proc];
116542c57489SHong Zhang     }
116642c57489SHong Zhang   }
116742c57489SHong Zhang 
1168ded4f561SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
11690298fd71SBarry Smith   ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
117042c57489SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
117142c57489SHong Zhang 
1172ded4f561SHong Zhang   /* post the Irecv and Isend of coj */
1173529bc5dcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
1174529bc5dcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
1175854ce69bSBarry Smith   ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr);
117642c57489SHong Zhang   for (proc=0, k=0; proc<size; proc++) {
117742c57489SHong Zhang     if (!len_s[proc]) continue;
1178de0260b3SHong Zhang     i    = owners_co[proc];
1179529bc5dcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
118042c57489SHong Zhang     k++;
118142c57489SHong Zhang   }
118242c57489SHong Zhang 
1183ded4f561SHong Zhang   /* receives and sends of coj are complete */
118477584fe6SHong Zhang   for (i=0; i<merge->nrecv; i++) {
1185c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1186ded4f561SHong Zhang   }
1187ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
11880c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
118982412ba7SHong Zhang 
11906b911d16SHong Zhang   /* (4) send and recv coi */
11916b911d16SHong Zhang   /*-----------------------*/
1192529bc5dcSHong Zhang   ierr   = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
1193529bc5dcSHong Zhang   ierr   = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
1194854ce69bSBarry Smith   ierr   = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr);
119542c57489SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
119642c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
119742c57489SHong Zhang     if (!len_s[proc]) continue;
119842c57489SHong Zhang     /* form outgoing message for i-structure:
119942c57489SHong Zhang          buf_si[0]:                 nrows to be sent
120042c57489SHong Zhang                [1:nrows]:           row index (global)
120142c57489SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
120242c57489SHong Zhang     */
120342c57489SHong Zhang     /*-------------------------------------------*/
12042addb229SHong Zhang     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
120542c57489SHong Zhang     buf_si_i    = buf_si + nrows+1;
120642c57489SHong Zhang     buf_si[0]   = nrows;
120742c57489SHong Zhang     buf_si_i[0] = 0;
120842c57489SHong Zhang     nrows       = 0;
1209de0260b3SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1210ded4f561SHong Zhang       nzi = coi[i+1] - coi[i];
1211ded4f561SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1212de0260b3SHong Zhang       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
121342c57489SHong Zhang       nrows++;
121442c57489SHong Zhang     }
1215529bc5dcSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
121642c57489SHong Zhang     k++;
121742c57489SHong Zhang     buf_si += len_si[proc];
121842c57489SHong Zhang   }
1219ded4f561SHong Zhang   i = merge->nrecv;
1220ded4f561SHong Zhang   while (i--) {
1221c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
1222ded4f561SHong Zhang   }
1223ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
12240c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
1225a914927fSHong Zhang 
122624ecddacSHong Zhang   ierr = PetscFree2(len_si,sstatus);CHKERRQ(ierr);
122742c57489SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
1228ded4f561SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
122942c57489SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
123042c57489SHong Zhang 
12316b911d16SHong Zhang   /* (5) compute the local portion of C (mpi mat) */
12326b911d16SHong Zhang   /*----------------------------------------------*/
1233ded4f561SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
1234ded4f561SHong Zhang 
123524ecddacSHong Zhang   /* allocate pti array and free space for accumulating nonzero column info */
1236854ce69bSBarry Smith   ierr   = PetscMalloc1(pn+1,&pti);CHKERRQ(ierr);
123724ecddacSHong Zhang   pti[0] = 0;
123842c57489SHong Zhang 
1239d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1240f91af8c7SBarry Smith   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(pi_loc[pm],api[am]));
124122d28d08SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr);
124242c57489SHong Zhang   current_space = free_space;
124342c57489SHong Zhang 
1244dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
124542c57489SHong Zhang   for (k=0; k<merge->nrecv; k++) {
124642c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
124742c57489SHong Zhang     nrows       = *buf_ri_k[k];
124842c57489SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
124942c57489SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
125042c57489SHong Zhang   }
125142c57489SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
1252936ad1eaSHong Zhang 
125342c57489SHong Zhang   for (i=0; i<pn; i++) {
1254ded4f561SHong Zhang     /* add pdt[i,:]*AP into lnk */
1255ded4f561SHong Zhang     pnz = pdti[i+1] - pdti[i];
125677584fe6SHong Zhang     ptJ = pdtj + pdti[i];
125777584fe6SHong Zhang     for (j=0; j<pnz; j++) {
125877584fe6SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
1259ded4f561SHong Zhang       apnz = api[row+1] - api[row];
1260ded4f561SHong Zhang       Jptr = apj + api[row];
1261ded4f561SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
1262a914927fSHong Zhang       ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr);
1263ded4f561SHong Zhang     }
1264a914927fSHong Zhang 
126542c57489SHong Zhang     /* add received col data into lnk */
126642c57489SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
126742c57489SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
1268ded4f561SHong Zhang         nzi  = *(nextci[k]+1) - *nextci[k];
1269ded4f561SHong Zhang         Jptr = buf_rj[k] + *nextci[k];
1270a914927fSHong Zhang         ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr);
127142c57489SHong Zhang         nextrow[k]++; nextci[k]++;
127242c57489SHong Zhang       }
127342c57489SHong Zhang     }
1274a914927fSHong Zhang     nnz = lnk[0];
127542c57489SHong Zhang 
127642c57489SHong Zhang     /* if free space is not available, make more free space */
1277ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
1278f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);CHKERRQ(ierr);
1279d16ca5f0SHong Zhang       nspacedouble++;
128042c57489SHong Zhang     }
128142c57489SHong Zhang     /* copy data into free space, then initialize lnk */
1282a914927fSHong Zhang     ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
1283ded4f561SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
12842205254eSKarl Rupp 
1285ded4f561SHong Zhang     current_space->array           += nnz;
1286ded4f561SHong Zhang     current_space->local_used      += nnz;
1287ded4f561SHong Zhang     current_space->local_remaining -= nnz;
12882205254eSKarl Rupp 
128924ecddacSHong Zhang     pti[i+1] = pti[i] + nnz;
129042c57489SHong Zhang   }
1291ded4f561SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
12920572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
129342c57489SHong Zhang 
1294854ce69bSBarry Smith   ierr      = PetscMalloc1(pti[pn]+1,&ptj);CHKERRQ(ierr);
129524ecddacSHong Zhang   ierr      = PetscFreeSpaceContiguous(&free_space,ptj);CHKERRQ(ierr);
129624ecddacSHong Zhang   afill_tmp = (PetscReal)pti[pn]/(pi_loc[pm] + api[am]+1);
1297d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
129842c57489SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
129942c57489SHong Zhang 
13006b911d16SHong Zhang   /* (6) create symbolic parallel matrix Cmpi */
13016b911d16SHong Zhang   /*------------------------------------------*/
130277584fe6SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
130377584fe6SHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
130433d57670SJed Brown   ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
130577584fe6SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
130677584fe6SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
130742c57489SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1308a2f3521dSMark F. Adams 
1309b091e509SHong Zhang   merge->bi        = pti;      /* Cseq->i */
1310b091e509SHong Zhang   merge->bj        = ptj;      /* Cseq->j */
1311b091e509SHong Zhang   merge->coi       = coi;      /* Co->i   */
1312b091e509SHong Zhang   merge->coj       = coj;      /* Co->j   */
131342c57489SHong Zhang   merge->buf_ri    = buf_ri;
131442c57489SHong Zhang   merge->buf_rj    = buf_rj;
1315de0260b3SHong Zhang   merge->owners_co = owners_co;
131677584fe6SHong Zhang   merge->destroy   = Cmpi->ops->destroy;
131742c57489SHong Zhang 
131877584fe6SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
131977584fe6SHong Zhang   c           = (Mat_MPIAIJ*)Cmpi->data;
1320f8487c73SHong Zhang   c->ptap     = ptap;
132177584fe6SHong Zhang   ptap->api   = api;
132277584fe6SHong Zhang   ptap->apj   = apj;
13238403a639SHong Zhang   ptap->duplicate = Cmpi->ops->duplicate;
13248403a639SHong Zhang   ptap->destroy   = Cmpi->ops->destroy;
13258403a639SHong Zhang 
13268403a639SHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
13278403a639SHong Zhang   Cmpi->assembled        = PETSC_FALSE;
13288403a639SHong Zhang   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
13298403a639SHong Zhang   Cmpi->ops->duplicate   = MatDuplicate_MPIAIJ_MatPtAP;
13308403a639SHong Zhang   Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap;
133177584fe6SHong Zhang   *C                     = Cmpi;
1332414894bdSHong Zhang 
1333414894bdSHong Zhang   /* flag 'scalable' determines which implementations to be used:
1334414894bdSHong Zhang        0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa;
1335414894bdSHong Zhang        1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */
1336414894bdSHong Zhang   /* set default scalable */
1337da0a95b2SSatish Balay   ptap->scalable = PETSC_FALSE; /* PETSC_TRUE; */
13382205254eSKarl Rupp 
1339c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)Cmpi)->options,((PetscObject)Cmpi)->prefix,"-matptap_scalable",&ptap->scalable,NULL);CHKERRQ(ierr);
1340414894bdSHong Zhang   if (!ptap->scalable) {  /* Do dense axpy */
13411795a4d1SJed Brown     ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr);
1342414894bdSHong Zhang   } else {
13431795a4d1SJed Brown     ierr = PetscCalloc1(ap_rmax+1,&ptap->apa);CHKERRQ(ierr);
1344414894bdSHong Zhang   }
1345414894bdSHong Zhang 
1346f2b054eeSHong Zhang #if defined(PETSC_USE_INFO)
134724ecddacSHong Zhang   if (pti[pn] != 0) {
134857622a8eSBarry Smith     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
134957622a8eSBarry Smith     ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
1350f2b054eeSHong Zhang   } else {
135177584fe6SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
1352f2b054eeSHong Zhang   }
1353f2b054eeSHong Zhang #endif
135442c57489SHong Zhang   PetscFunctionReturn(0);
135542c57489SHong Zhang }
135642c57489SHong Zhang 
13578403a639SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap(Mat A,Mat P,Mat C)
135842c57489SHong Zhang {
135942c57489SHong Zhang   PetscErrorCode      ierr;
1360f8487c73SHong Zhang   Mat_MPIAIJ          *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
136142c57489SHong Zhang   Mat_SeqAIJ          *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
1362de0260b3SHong Zhang   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
136342c57489SHong Zhang   Mat_SeqAIJ          *p_loc,*p_oth;
1364f8487c73SHong Zhang   Mat_PtAPMPI         *ptap;
13659f4155fbSHong Zhang   Mat_Merge_SeqsToMPI *merge;
13661d633602SHong Zhang   PetscInt            *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp;
136782412ba7SHong Zhang   PetscInt            *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj;
136882412ba7SHong Zhang   PetscInt            i,j,k,anz,pnz,apnz,nextap,row,*cj;
1369e900a757SHong Zhang   MatScalar           *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp;
1370d0f46423SBarry Smith   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1371ce94432eSBarry Smith   MPI_Comm            comm;
137242c57489SHong Zhang   PetscMPIInt         size,rank,taga,*len_s;
137382412ba7SHong Zhang   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
137442c57489SHong Zhang   PetscInt            **buf_ri,**buf_rj;
137550f5bf86SHong Zhang   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
137642c57489SHong Zhang   MPI_Request         *s_waits,*r_waits;
137742c57489SHong Zhang   MPI_Status          *status;
137882412ba7SHong Zhang   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
137982412ba7SHong Zhang   PetscInt            *api,*apj,*coi,*coj;
1380d0f46423SBarry Smith   PetscInt            *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend;
13816ce94e8fSHong Zhang   PetscBool           scalable;
138238571addSHong Zhang #if defined(PTAP_PROFILE)
1383e497d3c8SHong Zhang   PetscLogDouble t0,t1,t2,eP,t3,t4,et2_AP=0.0,ePtAP=0.0,t2_0,t2_1,t2_2;
138438571addSHong Zhang #endif
138542c57489SHong Zhang 
138642c57489SHong Zhang   PetscFunctionBegin;
1387ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr);
138838571addSHong Zhang #if defined(PTAP_PROFILE)
13898563dfccSBarry Smith   ierr = PetscTime(&t0);CHKERRQ(ierr);
139038571addSHong Zhang #endif
139142c57489SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
139242c57489SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
139342c57489SHong Zhang 
1394f8487c73SHong Zhang   ptap = c->ptap;
1395ce94432eSBarry 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");
1396f8487c73SHong Zhang   merge    = ptap->merge;
1397414894bdSHong Zhang   apa      = ptap->apa;
13986ce94e8fSHong Zhang   scalable = ptap->scalable;
13996c6e00e0SHong Zhang 
1400a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
14016b911d16SHong Zhang   /*-----------------------------------------------------*/
1402f8487c73SHong Zhang   if (ptap->reuse == MAT_INITIAL_MATRIX) {
14039f4155fbSHong Zhang     /* P_oth and P_loc are obtained in MatPtASymbolic(), skip calling MatGetBrowsOfAoCols() and MatMPIAIJGetLocalMat() */
1404f8487c73SHong Zhang     ptap->reuse = MAT_REUSE_MATRIX;
14056c6e00e0SHong Zhang   } else { /* update numerical values of P_oth and P_loc */
1406b7f45c76SHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
1407a1a4e84aSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
14086c6e00e0SHong Zhang   }
140938571addSHong Zhang #if defined(PTAP_PROFILE)
14108563dfccSBarry Smith   ierr = PetscTime(&t1);CHKERRQ(ierr);
141105d62848SHong Zhang   eP = t1-t0;
141238571addSHong Zhang #endif
141305d62848SHong Zhang   /*
141405d62848SHong Zhang   printf("[%d] Ad: %d, %d; Ao: %d, %d; P_loc: %d, %d; P_oth %d, %d;\n",rank,
141505d62848SHong Zhang          a->A->rmap->N,a->A->cmap->N,a->B->rmap->N,a->B->cmap->N,
141605d62848SHong Zhang          ptap->P_loc->rmap->N,ptap->P_loc->cmap->N,
141705d62848SHong Zhang          ptap->P_oth->rmap->N,ptap->P_oth->cmap->N);
141805d62848SHong Zhang    */
1419f8487c73SHong Zhang 
1420f9473cf7SHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
1421f9473cf7SHong Zhang   /*--------------------------------------------------------------*/
142242c57489SHong Zhang   /* get data from symbolic products */
1423a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
1424a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
142589ae1891SBarry Smith   pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a;
142642c57489SHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
142742c57489SHong Zhang 
1428de0260b3SHong Zhang   coi  = merge->coi; coj = merge->coj;
14291795a4d1SJed Brown   ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr);
1430de0260b3SHong Zhang 
1431de0260b3SHong Zhang   bi     = merge->bi; bj = merge->bj;
14327a2fc3feSBarry Smith   owners = merge->rowmap->range;
14331795a4d1SJed Brown   ierr   = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr);  /* ba: Cseq->a */
143442c57489SHong Zhang 
1435a1a4e84aSHong Zhang   api = ptap->api; apj = ptap->apj;
1436d9824c15SHong Zhang 
14379516a85cSHong Zhang   if (!scalable) { /* Do dense axpy on apa (length of pN, stores A[i,:]*P) - nonscalable, but faster (could take 1/3 scalable time) */
1438b38d1a2dSJed Brown     ierr = PetscInfo(C,"Using non-scalable dense axpy\n");CHKERRQ(ierr);
14398cb82516SHong Zhang #if defined(PTAP_PROFILE)
144005d62848SHong Zhang     ierr = PetscTime(&t1);CHKERRQ(ierr);
14418cb82516SHong Zhang #endif
1442a9d06231SHong Zhang     for (i=0; i<am; i++) {
1443b091e509SHong Zhang #if defined(PTAP_PROFILE)
14448563dfccSBarry Smith       ierr = PetscTime(&t2_0);CHKERRQ(ierr);
1445b091e509SHong Zhang #endif
1446a9d06231SHong Zhang       /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */
1447a9d06231SHong Zhang       /*------------------------------------------------------------*/
1448a9d06231SHong Zhang       apJ = apj + api[i];
1449a9d06231SHong Zhang 
1450a9d06231SHong Zhang       /* diagonal portion of A */
1451a9d06231SHong Zhang       anz = adi[i+1] - adi[i];
1452a9d06231SHong Zhang       adj = ad->j + adi[i];
1453a9d06231SHong Zhang       ada = ad->a + adi[i];
1454a9d06231SHong Zhang       for (j=0; j<anz; j++) {
1455a9d06231SHong Zhang         row = adj[j];
1456a9d06231SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
1457a9d06231SHong Zhang         pj  = pj_loc + pi_loc[row];
1458a9d06231SHong Zhang         pa  = pa_loc + pi_loc[row];
1459a9d06231SHong Zhang 
1460a9d06231SHong Zhang         /* perform dense axpy */
1461e900a757SHong Zhang         valtmp = ada[j];
1462a9d06231SHong Zhang         for (k=0; k<pnz; k++) {
1463e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
1464a9d06231SHong Zhang         }
1465a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
1466a9d06231SHong Zhang       }
1467a9d06231SHong Zhang 
1468a9d06231SHong Zhang       /* off-diagonal portion of A */
1469a9d06231SHong Zhang       anz = aoi[i+1] - aoi[i];
1470a9d06231SHong Zhang       aoj = ao->j + aoi[i];
1471a9d06231SHong Zhang       aoa = ao->a + aoi[i];
1472a9d06231SHong Zhang       for (j=0; j<anz; j++) {
1473a9d06231SHong Zhang         row = aoj[j];
1474a9d06231SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
1475a9d06231SHong Zhang         pj  = pj_oth + pi_oth[row];
1476a9d06231SHong Zhang         pa  = pa_oth + pi_oth[row];
1477a9d06231SHong Zhang 
1478a9d06231SHong Zhang         /* perform dense axpy */
1479e900a757SHong Zhang         valtmp = aoa[j];
1480a9d06231SHong Zhang         for (k=0; k<pnz; k++) {
1481e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
1482a9d06231SHong Zhang         }
1483a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
1484a9d06231SHong Zhang       }
1485b091e509SHong Zhang #if defined(PTAP_PROFILE)
14868563dfccSBarry Smith       ierr    = PetscTime(&t2_1);CHKERRQ(ierr);
1487b091e509SHong Zhang       et2_AP += t2_1 - t2_0;
1488b091e509SHong Zhang #endif
1489a9d06231SHong Zhang 
1490a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
1491a9d06231SHong Zhang       /*--------------------------------------------------------------*/
1492a9d06231SHong Zhang       apnz = api[i+1] - api[i];
1493a9d06231SHong Zhang       /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */
1494a9d06231SHong Zhang       pnz = po->i[i+1] - po->i[i];
1495e900a757SHong Zhang       poJ = po->j + po->i[i];
1496a9d06231SHong Zhang       pA  = po->a + po->i[i];
1497a9d06231SHong Zhang       for (j=0; j<pnz; j++) {
1498e900a757SHong Zhang         row = poJ[j];
1499e900a757SHong Zhang         cnz = coi[row+1] - coi[row];
1500e900a757SHong Zhang         cj  = coj + coi[row];
1501e900a757SHong Zhang         ca  = coa + coi[row];
1502a9d06231SHong Zhang         /* perform dense axpy */
1503e900a757SHong Zhang         valtmp = pA[j];
1504a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
1505e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
1506a9d06231SHong Zhang         }
1507a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1508a9d06231SHong Zhang       }
1509a9d06231SHong Zhang       /* put the value into Cd (diagonal part) */
1510a9d06231SHong Zhang       pnz = pd->i[i+1] - pd->i[i];
1511e900a757SHong Zhang       pdJ = pd->j + pd->i[i];
1512a9d06231SHong Zhang       pA  = pd->a + pd->i[i];
1513a9d06231SHong Zhang       for (j=0; j<pnz; j++) {
1514e900a757SHong Zhang         row = pdJ[j];
1515e900a757SHong Zhang         cnz = bi[row+1] - bi[row];
1516e900a757SHong Zhang         cj  = bj + bi[row];
1517e900a757SHong Zhang         ca  = ba + bi[row];
1518a9d06231SHong Zhang         /* perform dense axpy */
1519e900a757SHong Zhang         valtmp = pA[j];
1520a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
1521e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
1522a9d06231SHong Zhang         }
1523a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
1524a9d06231SHong Zhang       }
15258403a639SHong Zhang 
1526a9d06231SHong Zhang       /* zero the current row of A*P */
1527a9d06231SHong Zhang       for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0;
1528b091e509SHong Zhang #if defined(PTAP_PROFILE)
15298563dfccSBarry Smith       ierr      = PetscTime(&t2_2);CHKERRQ(ierr);
153005d62848SHong Zhang       ePtAP += t2_2 - t2_1;
1531b091e509SHong Zhang #endif
1532a9d06231SHong Zhang     }
153338571addSHong Zhang   } else { /* Do sparse axpy on apa (length of ap_rmax, stores A[i,:]*P) - scalable, but slower */
1534b38d1a2dSJed Brown     ierr = PetscInfo(C,"Using scalable sparse axpy\n");CHKERRQ(ierr);
153538571addSHong Zhang     /*-----------------------------------------------------------------------------------------*/
1536a9d06231SHong Zhang     pA=pa_loc;
153742c57489SHong Zhang     for (i=0; i<am; i++) {
1538b091e509SHong Zhang #if defined(PTAP_PROFILE)
15398563dfccSBarry Smith       ierr = PetscTime(&t2_0);CHKERRQ(ierr);
1540b091e509SHong Zhang #endif
1541f4a743e1SHong Zhang       /* form i-th sparse row of A*P */
154282412ba7SHong Zhang       apnz = api[i+1] - api[i];
154382412ba7SHong Zhang       apJ  = apj + api[i];
154442c57489SHong Zhang       /* diagonal portion of A */
154582412ba7SHong Zhang       anz = adi[i+1] - adi[i];
15461d633602SHong Zhang       adj = ad->j + adi[i];
15471d633602SHong Zhang       ada = ad->a + adi[i];
154882412ba7SHong Zhang       for (j=0; j<anz; j++) {
15491d633602SHong Zhang         row    = adj[j];
155082412ba7SHong Zhang         pnz    = pi_loc[row+1] - pi_loc[row];
155182412ba7SHong Zhang         pj     = pj_loc + pi_loc[row];
155282412ba7SHong Zhang         pa     = pa_loc + pi_loc[row];
1553e900a757SHong Zhang         valtmp = ada[j];
1554f4a743e1SHong Zhang         nextp  = 0;
155582412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
155682412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
1557e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
155842c57489SHong Zhang           }
155942c57489SHong Zhang         }
1560dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
156142c57489SHong Zhang       }
156242c57489SHong Zhang       /* off-diagonal portion of A */
156382412ba7SHong Zhang       anz = aoi[i+1] - aoi[i];
15641d633602SHong Zhang       aoj = ao->j + aoi[i];
15651d633602SHong Zhang       aoa = ao->a + aoi[i];
156682412ba7SHong Zhang       for (j=0; j<anz; j++) {
15671d633602SHong Zhang         row    = aoj[j];
156882412ba7SHong Zhang         pnz    = pi_oth[row+1] - pi_oth[row];
156982412ba7SHong Zhang         pj     = pj_oth + pi_oth[row];
157082412ba7SHong Zhang         pa     = pa_oth + pi_oth[row];
1571e900a757SHong Zhang         valtmp = aoa[j];
1572f4a743e1SHong Zhang         nextp  = 0;
157382412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
157482412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
1575e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
157642c57489SHong Zhang           }
157742c57489SHong Zhang         }
1578dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
157942c57489SHong Zhang       }
1580b091e509SHong Zhang #if defined(PTAP_PROFILE)
15818563dfccSBarry Smith       ierr    = PetscTime(&t2_1);CHKERRQ(ierr);
1582b091e509SHong Zhang       et2_AP += t2_1 - t2_0;
1583b091e509SHong Zhang #endif
158442c57489SHong Zhang 
1585a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
1586a9d06231SHong Zhang       /*--------------------------------------------------------------*/
158782412ba7SHong Zhang       pnz = pi_loc[i+1] - pi_loc[i];
1588f9473cf7SHong Zhang       pJ  = pj_loc + pi_loc[i];
158982412ba7SHong Zhang       for (j=0; j<pnz; j++) {
159042c57489SHong Zhang         nextap = 0;
1591f9473cf7SHong Zhang         row    = pJ[j]; /* global index */
159282412ba7SHong Zhang         if (row < pcstart || row >=pcend) { /* put the value into Co */
1593e900a757SHong Zhang           row = *poJ;
1594e900a757SHong Zhang           cj  = coj + coi[row];
1595e900a757SHong Zhang           ca  = coa + coi[row]; poJ++;
159682412ba7SHong Zhang         } else {                            /* put the value into Cd */
1597e900a757SHong Zhang           row = *pdJ;
1598e900a757SHong Zhang           cj  = bj + bi[row];
1599e900a757SHong Zhang           ca  = ba + bi[row]; pdJ++;
160042c57489SHong Zhang         }
1601e900a757SHong Zhang         valtmp = pA[j];
160282412ba7SHong Zhang         for (k=0; nextap<apnz; k++) {
1603e900a757SHong Zhang           if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++];
160442c57489SHong Zhang         }
1605dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
1606de0260b3SHong Zhang       }
16070496f32cSHong Zhang       pA += pnz;
160842c57489SHong Zhang       /* zero the current row info for A*P */
160982412ba7SHong Zhang       ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr);
1610b091e509SHong Zhang #if defined(PTAP_PROFILE)
16118563dfccSBarry Smith       ierr      = PetscTime(&t2_2);CHKERRQ(ierr);
161205d62848SHong Zhang       ePtAP += t2_2 - t2_1;
1613b091e509SHong Zhang #endif
161442c57489SHong Zhang     }
1615d9824c15SHong Zhang   }
161638571addSHong Zhang #if defined(PTAP_PROFILE)
16178563dfccSBarry Smith   ierr = PetscTime(&t2);CHKERRQ(ierr);
161838571addSHong Zhang #endif
161942c57489SHong Zhang 
1620a9d06231SHong Zhang   /* 3) send and recv matrix values coa */
1621a9d06231SHong Zhang   /*------------------------------------*/
162242c57489SHong Zhang   buf_ri = merge->buf_ri;
162342c57489SHong Zhang   buf_rj = merge->buf_rj;
162442c57489SHong Zhang   len_s  = merge->len_s;
1625fc42d0c8SSatish Balay   ierr   = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
162642c57489SHong Zhang   ierr   = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
162742c57489SHong Zhang 
1628dcca6d9dSJed Brown   ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr);
162942c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++) {
163042c57489SHong Zhang     if (!len_s[proc]) continue;
1631de0260b3SHong Zhang     i    = merge->owners_co[proc];
1632de0260b3SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
163342c57489SHong Zhang     k++;
163442c57489SHong Zhang   }
16350c468ba9SBarry Smith   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
16360c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
163742c57489SHong Zhang 
1638a9d06231SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
163942c57489SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
164082412ba7SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
164138571addSHong Zhang #if defined(PTAP_PROFILE)
16428563dfccSBarry Smith   ierr = PetscTime(&t3);CHKERRQ(ierr);
164338571addSHong Zhang #endif
164442c57489SHong Zhang 
1645a9d06231SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
1646a9d06231SHong Zhang   /*------------------------------------------------------*/
1647dcca6d9dSJed Brown   ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr);
164842c57489SHong Zhang   for (k=0; k<merge->nrecv; k++) {
164942c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
165042c57489SHong Zhang     nrows       = *(buf_ri_k[k]);
165142c57489SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
165282412ba7SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
165342c57489SHong Zhang   }
165442c57489SHong Zhang 
1655de0260b3SHong Zhang   for (i=0; i<cm; i++) {
165682412ba7SHong Zhang     row  = owners[rank] + i; /* global row index of C_seq */
165742c57489SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1658de0260b3SHong Zhang     ba_i = ba + bi[i];
165982412ba7SHong Zhang     bnz  = bi[i+1] - bi[i];
166042c57489SHong Zhang     /* add received vals into ba */
166142c57489SHong Zhang     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
166242c57489SHong Zhang       /* i-th row */
166342c57489SHong Zhang       if (i == *nextrow[k]) {
166482412ba7SHong Zhang         cnz    = *(nextci[k]+1) - *nextci[k];
166582412ba7SHong Zhang         cj     = buf_rj[k] + *(nextci[k]);
166682412ba7SHong Zhang         ca     = abuf_r[k] + *(nextci[k]);
166782412ba7SHong Zhang         nextcj = 0;
166882412ba7SHong Zhang         for (j=0; nextcj<cnz; j++) {
166982412ba7SHong Zhang           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
167082412ba7SHong Zhang             ba_i[j] += ca[nextcj++];
167142c57489SHong Zhang           }
167242c57489SHong Zhang         }
167382412ba7SHong Zhang         nextrow[k]++; nextci[k]++;
1674c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
167542c57489SHong Zhang       }
167642c57489SHong Zhang     }
167782412ba7SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
167842c57489SHong Zhang   }
167942c57489SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
168042c57489SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
168142c57489SHong Zhang 
1682de0260b3SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
1683c05d87d6SBarry Smith   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
168442c57489SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
16850572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
168638571addSHong Zhang #if defined(PTAP_PROFILE)
16878563dfccSBarry Smith   ierr = PetscTime(&t4);CHKERRQ(ierr);
16889516a85cSHong Zhang   if (rank==1) {
168905d62848SHong 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);
16909516a85cSHong Zhang   }
169138571addSHong Zhang #endif
169242c57489SHong Zhang   PetscFunctionReturn(0);
169342c57489SHong Zhang }
1694